phaseModel.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) 2017-2021 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "phaseModel.H"
29 #include "phaseSystem.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(phaseModel, 0);
36  defineRunTimeSelectionTable(phaseModel, phaseSystem);
37 }
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42 (
43  const phaseSystem& fluid,
44  const word& phaseName
45 )
46 :
48  (
49  IOobject
50  (
51  IOobject::groupName("alpha", phaseName),
52  fluid.mesh().time().timeName(),
53  fluid.mesh(),
54  IOobject::READ_IF_PRESENT,
55  IOobject::AUTO_WRITE
56  ),
57  fluid.mesh(),
59  ),
60  fluid_(fluid),
61  name_(phaseName)
62 {}
63 
64 
65 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
66 
69 (
70  const phaseSystem& fluid,
71  const word& phaseName
72 )
73 {
74  const dictionary& dict = fluid.subDict(phaseName);
75 
76  const word modelType(dict.get<word>("type"));
77 
78  Info<< "Selecting phaseModel for "
79  << phaseName << ": " << modelType << endl;
80 
81  auto* ctorPtr = phaseSystemConstructorTable(modelType);
82 
83  if (!ctorPtr)
84  {
86  (
87  dict,
88  "phaseModel",
89  modelType,
90  *phaseSystemConstructorTablePtr_
91  ) << exit(FatalIOError);
92  }
93 
94  return ctorPtr(fluid, phaseName);
95 }
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
101 {
102  return fluid_;
103 }
104 
105 
107 {
108  thermo().correct();
109 }
110 
111 
113 {
114  // do nothing
115 }
116 
117 
119 {
120  return thermo().rho();
121 }
122 
123 
125 {
126  return thermo().rho(patchI);
127 }
128 
129 
131 {
132  return thermo().hc();
133 }
134 
135 
137 {
138  return thermo().Cp();
139 }
140 
141 
143 (
144  const scalarField& p,
145  const scalarField& T,
146  const label patchI
147 ) const
148 {
149  return (thermo().Cp(p, T, patchI));
150 }
151 
152 
154 {
155  return thermo().Cv();
156 }
157 
158 
160 (
161  const scalarField& p,
162  const scalarField& T,
163  const label patchI
164 ) const
165 {
166  return thermo().Cv(p, T, patchI);
167 }
168 
169 
171 {
172  return thermo().gamma();
173 }
174 
175 
177 (
178  const scalarField& p,
179  const scalarField& T,
180  const label patchI
181 ) const
182 {
183  return thermo().gamma(p, T, patchI);
184 }
185 
186 
188 {
189  return thermo().Cpv();
190 }
191 
192 
194 (
195  const scalarField& p,
196  const scalarField& T,
197  const label patchI
198 ) const
199 {
200  return thermo().Cpv(p, T, patchI);
201 }
202 
203 
205 {
206  return thermo().CpByCpv();
207 }
208 
209 
211 (
212  const scalarField& p,
213  const scalarField& T,
214  const label patchI
215 ) const
216 {
217  return thermo().CpByCpv(p, T, patchI);
218 }
219 
220 
222 {
223  return thermo().alpha();
224 }
225 
226 
227 const Foam::scalarField& Foam::phaseModel::alpha(const label patchI) const
228 {
229  return thermo().alpha(patchI);
230 }
231 
232 
234 {
235  return thermo().kappa();
236 }
237 
238 
240 {
241  return thermo().kappa(patchI);
242 }
243 
244 
246 {
247  return thermo().alphahe();
248 }
249 
250 
252 {
253  return thermo().alphahe(patchI);
254 }
255 
256 
258 (
259  const volScalarField& kappat
260 ) const
261 {
262  tmp<volScalarField> kappaEff(kappa() + kappat);
263  kappaEff.ref().rename("kappaEff" + name_);
264  return kappaEff;
265 }
266 
267 
269 (
270  const scalarField& kappat,
271  const label patchI
272 ) const
273 {
274  return (kappa(patchI) + kappat);
275 }
276 
277 
279 (
280  const volScalarField& alphat
281 ) const
282 {
283  return (thermo().alpha() + alphat);
284 }
285 
286 
288 (
289  const scalarField& alphat,
290  const label patchI
291 ) const
292 {
293  return (thermo().alpha(patchI) + alphat);
294 }
295 
296 
298 {
299  return thermo().mu();
300 }
301 
302 
304 {
305  return thermo().mu(patchi);
306 }
307 
308 
310 {
311  return thermo().nu();
312 }
313 
314 
316 {
317  return thermo().nu(patchi);
318 }
319 
320 
322 {
323  return true;
324 }
325 
326 
327 // ************************************************************************* //
Foam::phaseModel::mu
virtual tmp< volScalarField > mu() const
Return the mixture dymanic viscosity.
Definition: phaseModel.C:297
phaseModel.H
Foam::phaseModel::correctTurbulence
virtual void correctTurbulence()
Correct the turbulence.
Definition: phaseModel.C:112
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::phaseModel::Cv
tmp< volScalarField > Cv() const
Return Cv of the phase.
Definition: phaseModel.C:153
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::phaseModel::kappa
const dimensionedScalar & kappa() const
Definition: phaseModel.H:157
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::phaseModel::kappaEff
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent diffusivity for temperature.
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::basicThermo::Cp
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure [J/kg/K].
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
thermo
psiReactionThermo & thermo
Definition: createFields.H:28
kappaEff
kappaEff
Definition: TEqn.H:10
Foam::phaseModel::correct
void correct()
Correct the phase properties.
Definition: phaseModel.C:215
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
Foam::basicThermo::gamma
virtual tmp< volScalarField > gamma() const =0
Gamma = Cp/Cv [].
Foam::constant::electromagnetic::kappa
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Foam::basicThermo::hc
virtual tmp< volScalarField > hc() const =0
Chemical enthalpy [J/kg].
Foam::psiThermo::rho
virtual tmp< volScalarField > rho() const
Density [kg/m^3] - uses current value of pressure.
Definition: psiThermo.C:146
Foam::phaseModel::Cp
const dimensionedScalar & Cp() const
Definition: phaseModel.H:162
Foam::phaseModel::phaseModel
phaseModel(const word &phaseName, const dictionary &phaseDict, const fvMesh &mesh)
Definition: phaseModel.C:39
Foam::phaseModel::CpByCpv
tmp< volScalarField > CpByCpv() const
Heat capacity ratio for phase [].
Definition: phaseModel.C:204
Foam::basicThermo::correct
virtual void correct()=0
Update properties.
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::basicThermo::alphahe
virtual tmp< volScalarField > alphahe() const =0
Thermal diffusivity for energy of mixture [kg/m/s].
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::phaseModel::Cpv
tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume for phase [J/kg/K].
Definition: phaseModel.C:187
Foam::phaseModel::nu
const dimensionedScalar & nu() const
Definition: phaseModel.H:152
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::basicThermo::kappa
virtual tmp< volScalarField > kappa() const =0
Thermal diffusivity for temperature of mixture [J/m/s/K].
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::basicThermo::CpByCpv
virtual tmp< volScalarField > CpByCpv() const =0
Heat capacity ratio [].
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::phaseModel::gamma
tmp< volScalarField > gamma() const
Gamma = Cp/Cv of phase[].
Definition: phaseModel.C:170
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::psiThermo::mu
virtual tmp< volScalarField > mu() const
Dynamic viscosity of mixture [kg/m/s].
Definition: psiThermo.C:178
Foam::fluidThermo::nu
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
Definition: fluidThermo.C:103
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::basicThermo::alpha
virtual const volScalarField & alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
Definition: basicThermo.C:626
Foam::phaseModel::alpha
const volScalarField & alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
Definition: phaseModel.C:221
Foam::phaseModel::read
virtual bool read()
Read phase properties dictionary.
Definition: phaseModel.C:321
Foam::basicThermo::Cpv
virtual tmp< volScalarField > Cpv() const =0
Heat capacity at constant pressure/volume [J/kg/K].
Foam::phaseModel::hc
tmp< volScalarField > hc() const
Chemical enthalpy for phase [J/kg].
Definition: phaseModel.C:130
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::phaseModel::fluid
const phaseSystem & fluid() const
Return the system to which this phase belongs.
Definition: phaseModel.C:100
Foam::phaseModel::alphaEff
virtual tmp< volScalarField > alphaEff() const =0
Effective thermal turbulent diffusivity of mixture [kg/m/s].
Foam::phaseModel::alphahe
tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
Definition: phaseModel.C:245
Foam::phaseSystem
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:66
Foam::phaseModel::rho
const dimensionedScalar & rho() const
Definition: phaseModel.H:167
Foam::basicThermo::Cv
virtual tmp< volScalarField > Cv() const =0
Heat capacity at constant volume [J/kg/K].
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::phaseModel::New
static autoPtr< phaseModel > New(const phaseSystem &fluid, const word &phaseName)
Definition: phaseModel.C:69
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189