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  g_(fluid.mesh().lookupObject<uniformDimensionedVectorField>("g")),
76  thermophysicalTransport_
77  (
79  <
80  phaseCompressible::momentumTransportModel,
81  transportThermoModel
82  >::New(this->momentumTransport_, this->thermo_)
83  )
84 {}
85 
86 
87 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
88 
89 template<class BasePhaseModel>
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
96 template<class BasePhaseModel>
98 {
99  BasePhaseModel::correctThermo();
100 
101  this->thermo_->correct();
102 }
103 
104 
105 template<class BasePhaseModel>
107 {
108  return false;
109 }
110 
111 
112 template<class BasePhaseModel>
115 {
116  BasePhaseModel::predictThermophysicalTransport();
117  thermophysicalTransport_->predict();
118 }
119 
120 
121 template<class BasePhaseModel>
124 {
125  BasePhaseModel::correctThermophysicalTransport();
126  thermophysicalTransport_->correct();
127 }
128 
129 
130 template<class BasePhaseModel>
133 {
134  return thermophysicalTransport_->kappaEff(patchi);
135 }
136 
137 
138 template<class BasePhaseModel>
141 {
142  return thermophysicalTransport_->divq(he);
143 }
144 
145 
146 template<class BasePhaseModel>
149 {
150  return thermophysicalTransport_->divj(Yi);
151 }
152 
153 
154 template<class BasePhaseModel>
157 {
158  const volScalarField& alpha = *this;
159  const volScalarField& rho = this->rho();
160 
161  const tmp<volVectorField> tU(this->U());
162  const volVectorField& U(tU());
163 
164  const tmp<surfaceScalarField> talphaRhoPhi(this->alphaRhoPhi());
165  const surfaceScalarField& alphaRhoPhi(talphaRhoPhi());
166 
167  const tmp<volScalarField> tcontErr(this->continuityError());
168  const volScalarField& contErr(tcontErr());
169 
170  tmp<volScalarField> tK(this->K());
171  const volScalarField& K(tK());
172 
173  volScalarField& he = this->thermo_->he();
174 
175  tmp<fvScalarMatrix> tEEqn
176  (
177  fvm::ddt(alpha, rho, he)
178  + fvm::div(alphaRhoPhi, he)
179  - fvm::Sp(contErr, he)
180 
181  + fvc::ddt(alpha, rho, K) + fvc::div(alphaRhoPhi, K)
182  - contErr*K
183 
184  + this->divq(he)
185  ==
186  alpha*rho*(U&g_)
187  + alpha*this->Qdot()
188  );
189 
190  // Add the appropriate pressure-work term
191  if (he.name() == this->thermo_->phasePropertyName("e"))
192  {
193  tEEqn.ref() += filterPressureWork
194  (
195  fvc::div
196  (
197  fvc::absolute(alphaRhoPhi, alpha, rho, U),
198  this->fluidThermo().p()/rho
199  )
200  + (fvc::ddt(alpha) - contErr/rho)*this->fluidThermo().p()
201  );
202  }
203  else if (this->thermo_->dpdt())
204  {
205  tEEqn.ref() -= filterPressureWork(alpha*this->fluid().dpdt());
206  }
207 
208  return tEEqn;
209 }
210 
211 
212 // ************************************************************************* //
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 void predictThermophysicalTransport()
Predict the energy transport e.g. alphat.
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction.
virtual void correctThermophysicalTransport()
Correct the energy transport e.g. alphat.
virtual bool isothermal() const
Return whether the phase is isothermal.
virtual tmp< fvScalarMatrix > heEqn()
Return the enthalpy equation.
virtual tmp< scalarField > kappaEff(const label patchi) const
Return the effective thermal conductivity on a patch.
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Generic GeometricField class.
const word & name() const
Return name.
Definition: IOobject.H:310
Templated base class for multiphase thermophysical transport models.
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
Abstract base class for turbulence models (RAS, LES and laminar).
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.
label patchi
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))
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
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 > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
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:64
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
thermo he()
scalar Qdot
Definition: solveChemistry.H:2
volScalarField & p