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-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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"
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33namespace Foam
34{
35namespace multiphaseInter
36{
39}
40}
41
42// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43
45(
47 const word& phaseName
48)
49:
51 (
53 (
54 IOobject::groupName("alpha", phaseName),
55 fluid.mesh().time().timeName(),
56 fluid.mesh(),
57 IOobject::READ_IF_PRESENT,
58 IOobject::AUTO_WRITE
59 ),
60 fluid.mesh(),
62 ),
63 fluid_(fluid),
64 name_(phaseName)
65{}
66
67
68// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
69
72(
74 const word& phaseName
75)
76{
77 const dictionary& dict = fluid.subDict(phaseName);
78
79 const word modelType(dict.get<word>("type"));
80
81 Info<< "Selecting phaseModel for "
82 << phaseName << ": " << modelType << endl;
83
84 auto* ctorPtr = multiphaseInterSystemConstructorTable(modelType);
85
86 if (!ctorPtr)
87 {
89 (
90 dict,
91 "phaseModel",
92 modelType,
93 *multiphaseInterSystemConstructorTablePtr_
94 ) << exit(FatalIOError);
95 }
96
97 return ctorPtr(fluid, phaseName);
98}
99
100
101// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102
105{
106 return fluid_;
107}
108
109
111{
112 thermo().correct();
113}
114
115
117{
118 // do nothing
119}
120
121
123{
124 return thermo().rho();
125}
126
127
130{
131 return thermo().rho(patchI);
132}
133
134
136{
137 return thermo().hc();
138}
139
140
142{
143 return thermo().Cp();
144}
145
146
148(
149 const scalarField& p,
150 const scalarField& T,
151 const label patchI
152) const
153{
154 return (thermo().Cp(p, T, patchI));
155}
156
157
159{
160 return thermo().Cv();
161}
162
163
165(
166 const scalarField& p,
167 const scalarField& T,
168 const label patchI
169) const
170{
171 return thermo().Cv(p, T, patchI);
172}
173
174
176{
177 return thermo().gamma();
178}
179
180
182(
183 const scalarField& p,
184 const scalarField& T,
185 const label patchI
186) const
187{
188 return thermo().gamma(p, T, patchI);
189}
190
191
193{
194 return thermo().Cpv();
195}
196
197
199(
200 const scalarField& p,
201 const scalarField& T,
202 const label patchI
203) const
204{
205 return thermo().Cpv(p, T, patchI);
206}
207
208
210{
211 return thermo().CpByCpv();
212}
213
214
216(
217 const scalarField& p,
218 const scalarField& T,
219 const label patchI
220) const
221{
222 return thermo().CpByCpv(p, T, patchI);
223}
224
225
227{
228 return thermo().alpha();
229}
230
231
234{
235 return thermo().alpha(patchI);
236}
237
238
240{
241 return thermo().kappa();
242}
243
244
247{
248 return thermo().kappa(patchI);
249}
250
251
254{
255 return thermo().alphahe();
256}
257
258
261{
262 return thermo().alphahe(patchI);
263}
264
265
267(
268 const volScalarField& kappat
269) const
270{
271 tmp<volScalarField> kappaEff(kappa() + kappat);
272 kappaEff.ref().rename("kappaEff" + name_);
273 return kappaEff;
274}
275
276
278(
279 const scalarField& kappat,
280 const label patchI
281) const
282{
283 return (kappa(patchI) + kappat);
284}
285
286
288(
289 const volScalarField& alphat
290) const
291{
292 return (thermo().alpha() + alphat);
293}
294
295
297(
298 const scalarField& alphat,
299 const label patchI
300) const
301{
302 return (thermo().alpha(patchI) + alphat);
303}
304
305
307{
308 return thermo().mu();
309}
310
311
314{
315 return thermo().mu(patchi);
316}
317
318
320{
321 return thermo().nu();
322}
323
324
327{
328 return thermo().nu(patchi);
329}
330
331
333{
334 return true;
335}
336
337
338// ************************************************************************* //
twoPhaseSystem & fluid
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
tmp< volScalarField > alphaEff() const
Effective thermal turbulent diffusivity of mixture [kg/m/s].
tmp< volScalarField > Cp() const
Return phase Cp.
Definition: phaseModel.C:141
virtual void correct()
Correct phase thermo.
Definition: phaseModel.C:110
tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
Definition: phaseModel.C:253
virtual tmp< volScalarField > mu() const
Return the mixture dymanic viscosity.
Definition: phaseModel.C:306
virtual tmp< volScalarField > nu() const
Return the mixture kinematic viscosity.
Definition: phaseModel.C:319
const multiphaseInterSystem & fluid() const
Return the system to which this phase belongs.
Definition: phaseModel.C:104
tmp< volScalarField > CpByCpv() const
Heat capacity ratio for phase [].
Definition: phaseModel.C:209
virtual void correctTurbulence()
Correct the turbulence.
Definition: phaseModel.C:116
tmp< volScalarField > hc() const
Chemical enthalpy for phase [J/kg].
Definition: phaseModel.C:135
tmp< volScalarField > gamma() const
Gamma = Cp/Cv of phase[].
Definition: phaseModel.C:175
tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of phase [J/m/s/K].
Definition: phaseModel.C:239
tmp< volScalarField > rho() const
Return the phase density.
Definition: phaseModel.C:122
const volScalarField & alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
Definition: phaseModel.C:226
tmp< volScalarField > Cv() const
Return Cv of the phase.
Definition: phaseModel.C:158
virtual bool read()
Read phase properties dictionary.
Definition: phaseModel.C:332
tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume for phase [J/kg/K].
Definition: phaseModel.C:192
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
const volScalarField & T
const volScalarField & Cp
Definition: EEqn.H:7
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
word timeName
Definition: getTimeIndex.H:3
kappaEff
Definition: TEqn.H:10
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
IOerror FatalIOError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
volScalarField & alpha
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict