AnisothermalPhaseModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "AnisothermalPhaseModel.H"
27 #include "phaseSystem.H"
28 #include "fvcMeshPhi.H"
29 #include "fvcDdt.H"
30 #include "fvmDiv.H"
31 #include "fvmSup.H"
32 
33 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
34 
35 template<class BasePhaseModel>
38 (
39  const tmp<volScalarField>& pressureWork
40 ) const
41 {
42  const volScalarField& alpha = *this;
43 
44  scalar pressureWorkAlphaLimit =
45  this->thermo_->properties()
46  .lookupOrDefault("pressureWorkAlphaLimit", 0.0);
47 
48  if (pressureWorkAlphaLimit > 0)
49  {
50  return
51  (
52  max(alpha - pressureWorkAlphaLimit, scalar(0))
53  /max(alpha - pressureWorkAlphaLimit, pressureWorkAlphaLimit)
54  )*pressureWork;
55  }
56  else
57  {
58  return pressureWork;
59  }
60 }
61 
62 
63 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64 
65 template<class BasePhaseModel>
67 (
68  const phaseSystem& fluid,
69  const word& phaseName,
70  const bool referencePhase,
71  const label index
72 )
73 :
74  BasePhaseModel(fluid, phaseName, referencePhase, index)
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
79 
80 template<class BasePhaseModel>
82 {}
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
87 template<class BasePhaseModel>
89 {
90  BasePhaseModel::correctThermo();
91 
92  this->thermo_->correct();
93 }
94 
95 
96 template<class BasePhaseModel>
98 {
99  return false;
100 }
101 
102 
103 template<class BasePhaseModel>
106 {
107  const volScalarField& alpha = *this;
108  const volScalarField& rho = this->rho();
109 
110  const tmp<volVectorField> tU(this->U());
111  const volVectorField& U(tU());
112 
113  const tmp<surfaceScalarField> talphaRhoPhi(this->alphaRhoPhi());
114  const surfaceScalarField& alphaRhoPhi(talphaRhoPhi());
115 
116  const tmp<volScalarField> tcontErr(this->continuityError());
117  const volScalarField& contErr(tcontErr());
118 
119  tmp<volScalarField> tK(this->K());
120  const volScalarField& K(tK());
121 
122  volScalarField& he = this->thermo_->he();
123 
124  tmp<fvScalarMatrix> tEEqn
125  (
126  fvm::ddt(alpha, rho, he)
127  + fvm::div(alphaRhoPhi, he)
128  - fvm::Sp(contErr, he)
129 
130  + fvc::ddt(alpha, rho, K) + fvc::div(alphaRhoPhi, K)
131  - contErr*K
132 
133  + this->divq(he)
134  ==
135  alpha*this->Qdot()
136  );
137 
138  // Add the appropriate pressure-work term
139  if (he.name() == this->thermo_->phasePropertyName("e"))
140  {
141  tEEqn.ref() += filterPressureWork
142  (
143  fvc::div
144  (
145  fvc::absolute(alphaRhoPhi, alpha, rho, U),
146  this->thermo().p()/rho
147  )
148  + (fvc::ddt(alpha) - contErr/rho)*this->thermo().p()
149  );
150  }
151  else if (this->thermo_->dpdt())
152  {
153  tEEqn.ref() -= filterPressureWork(alpha*this->fluid().dpdt());
154  }
155 
156  return tEEqn;
157 }
158 
159 
160 // ************************************************************************* //
Class which represents a phase for which the temperature (strictly energy) varies....
virtual ~AnisothermalPhaseModel()
Destructor.
AnisothermalPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
virtual void correctThermo()
Correct the thermodynamics.
virtual bool isothermal() const
Return whether the phase is isothermal.
virtual tmp< fvScalarMatrix > heEqn()
Return the enthalpy equation.
Generic GeometricField class.
const word & name() const
Return name.
Definition: IOobject.H:310
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
Calculate the first temporal derivative.
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for implicit and explicit sources.
K
Definition: pEqn.H:75
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
thermo he()
scalar Qdot
Definition: solveChemistry.H:2
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31