ThermalPhaseModel.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-2025 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 "ThermalPhaseModel.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  ThermophysicalTransportPhaseModel<BasePhaseModel>
75  (
76  fluid,
77  phaseName,
78  referencePhase,
79  index
80  ),
81  g_(fluid.mesh().lookupObject<uniformDimensionedVectorField>("g"))
82 {}
83 
84 
85 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
86 
87 template<class BasePhaseModel>
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 
94 template<class BasePhaseModel>
96 {
97  BasePhaseModel::correctThermo();
98 
99  this->thermo_->correct();
100 }
101 
102 
103 template<class BasePhaseModel>
105 {
106  return false;
107 }
108 
109 
110 template<class BasePhaseModel>
113 {
114  const volScalarField& alpha = *this;
115  const volScalarField& rho = this->rho();
116 
117  const tmp<volVectorField> tU(this->U());
118  const volVectorField& U(tU());
119 
120  const tmp<surfaceScalarField> talphaRhoPhi(this->alphaRhoPhi());
121  const surfaceScalarField& alphaRhoPhi(talphaRhoPhi());
122 
123  const tmp<volScalarField> tcontErr(this->continuityError());
124  const volScalarField& contErr(tcontErr());
125 
126  tmp<volScalarField> tK(this->K());
127  const volScalarField& K(tK());
128 
129  volScalarField& he = this->thermo_->he();
130 
131  tmp<fvScalarMatrix> tEEqn
132  (
133  fvm::ddt(alpha, rho, he)
134  + fvm::div(alphaRhoPhi, he)
135  - fvm::Sp(contErr, he)
136 
137  + fvc::ddt(alpha, rho, K) + fvc::div(alphaRhoPhi, K)
138  - contErr*K
139 
140  + this->divq(he)
141  ==
142  alpha*rho*(U&g_)
143  + alpha*this->Qdot()
144  );
145 
146  // Add the appropriate pressure-work term
147  if (he.name() == this->thermo_->phasePropertyName("e"))
148  {
149  tEEqn.ref() += filterPressureWork
150  (
151  fvc::div
152  (
153  fvc::absolute(alphaRhoPhi, alpha, rho, U),
154  this->fluidThermo().p()/rho
155  )
156  + (fvc::ddt(alpha) - contErr/rho)*this->fluidThermo().p()
157  );
158  }
159  else if (this->thermo_->dpdt())
160  {
161  tEEqn.ref() -= filterPressureWork(alpha*this->fluid().dpdt());
162  }
163 
164  return tEEqn;
165 }
166 
167 
168 // ************************************************************************* //
Generic GeometricField class.
const word & name() const
Return name.
Definition: IOobject.H:307
Class which represents a phase for which the temperature (strictly energy) varies....
virtual ~ThermalPhaseModel()
Destructor.
ThermalPhaseModel(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.
Class which makes thermophysical transport modelling available to derived classes.
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:56
Class to represent a system of phases.
Definition: phaseSystem.H:74
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:197
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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 > > 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:62
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
thermo he()
scalar Qdot
Definition: solveChemistry.H:2
volScalarField & p