AnisothermalPhaseModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2018 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "AnisothermalPhaseModel.H"
30 #include "phaseSystem.H"
31 
32 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
33 
34 template<class BasePhaseModel>
37 (
38  const tmp<volScalarField>& pressureWork
39 ) const
40 {
41  const volScalarField& alpha = *this;
42 
43  const scalar pressureWorkAlphaLimit =
44  this->thermo_->getOrDefault("pressureWorkAlphaLimit", scalar(0));
45 
46  if (pressureWorkAlphaLimit > 0)
47  {
48  return
49  (
50  max(alpha - pressureWorkAlphaLimit, scalar(0))
51  /max(alpha - pressureWorkAlphaLimit, pressureWorkAlphaLimit)
52  )*pressureWork;
53  }
54 
55  return pressureWork;
56 }
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
61 template<class BasePhaseModel>
63 (
64  const phaseSystem& fluid,
65  const word& phaseName,
66  const label index
67 )
68 :
69  BasePhaseModel(fluid, phaseName, index)
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
74 
75 template<class BasePhaseModel>
77 {}
78 
79 
80 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
81 
82 template<class BasePhaseModel>
84 {
86 
87  this->thermo_->correct();
88 }
89 
90 
91 template<class BasePhaseModel>
93 {
94  return false;
95 }
96 
97 
98 template<class BasePhaseModel>
101 {
102  const volScalarField& alpha = *this;
103 
104  const volVectorField U(this->U());
105  const surfaceScalarField alphaPhi(this->alphaPhi());
106  const surfaceScalarField alphaRhoPhi(this->alphaRhoPhi());
107 
108  const volScalarField contErr(this->continuityError());
109  const volScalarField K(this->K());
110 
111  volScalarField& he = this->thermo_->he();
112 
113  tmp<fvScalarMatrix> tEEqn
114  (
115  fvm::ddt(alpha, this->rho(), he)
116  + fvm::div(alphaRhoPhi, he)
117  - fvm::Sp(contErr, he)
118 
119  + fvc::ddt(alpha, this->rho(), K) + fvc::div(alphaRhoPhi, K)
120  - contErr*K
121 
123  (
125  *fvc::interpolate(this->alphaEff()),
126  he
127  )
128  ==
129  alpha*this->Qdot()
130  );
131 
132  // Add the appropriate pressure-work term
133  if (he.name() == this->thermo_->phasePropertyName("e"))
134  {
135  tEEqn.ref() += filterPressureWork
136  (
137  fvc::div(fvc::absolute(alphaPhi, alpha, U), this->thermo().p())
138  + this->thermo().p()*fvc::ddt(alpha)
139  );
140  }
141  else if (this->thermo_->dpdt())
142  {
143  tEEqn.ref() -= filterPressureWork(alpha*this->fluid().dpdt());
144  }
145 
146  return tEEqn;
147 }
148 
149 
150 // ************************************************************************* //
contErr
volScalarField::Internal contErr((fvc::ddt(rho)+fvc::div(rhoPhi) -(fvOptions(alpha1, mixture.thermo1().rho())&rho1) -(fvOptions(alpha2, mixture.thermo2().rho())&rho2))())
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::AnisothermalPhaseModel::AnisothermalPhaseModel
AnisothermalPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
Definition: AnisothermalPhaseModel.C:63
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
thermo
psiReactionThermo & thermo
Definition: createFields.H:28
Sp
zeroField Sp
Definition: alphaSuSp.H:2
Foam::AnisothermalPhaseModel::correctThermo
virtual void correctThermo()
Correct the thermodynamics.
Definition: AnisothermalPhaseModel.C:83
Foam::AnisothermalPhaseModel::~AnisothermalPhaseModel
virtual ~AnisothermalPhaseModel()
Destructor.
Definition: AnisothermalPhaseModel.C:76
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
rho
rho
Definition: readInitialConditions.H:88
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:58
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
AnisothermalPhaseModel.H
alphaPhi
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::AnisothermalPhaseModel::heEqn
virtual tmp< fvScalarMatrix > heEqn()
Return the enthalpy equation.
Definition: AnisothermalPhaseModel.C:100
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Qdot
scalar Qdot
Definition: solveChemistry.H:2
Foam::AnisothermalPhaseModel
Class which represents a phase for which the temperature (strictly energy) varies....
Definition: AnisothermalPhaseModel.H:52
correctThermo
mixture correctThermo()
Foam::dimensioned::getOrDefault
static dimensioned< Type > getOrDefault(const word &name, const dictionary &dict, const dimensionSet &dims=dimless, const Type &deflt=Type(Zero))
Construct dimensioned from dictionary, with default value.
Definition: dimensionedType.C:348
U
U
Definition: pEqn.H:72
he
volScalarField & he
Definition: YEEqn.H:52
Foam::fac::ddt
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:47
alphaEff
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
Foam::AnisothermalPhaseModel::isothermal
virtual bool isothermal() const
Return whether the phase is isothermal.
Definition: AnisothermalPhaseModel.C:92
Foam::phaseSystem
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:63
dpdt
volScalarField & dpdt
Definition: setRegionFluidFields.H:32
Foam::fac::laplacian
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:47
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190