phaseModel.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-2024 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 "phaseModel.H"
27 #include "phaseSystem.H"
28 #include "diameterModel.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
36 }
37 
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42 (
43  const phaseSystem& fluid,
44  const word& phaseName,
45  const bool referencePhase,
46  const label index
47 )
48 :
49  // volScalarField
50  // (
51  // referencePhase
52  // ? volScalarField
53  // (
54  // IOobject
55  // (
56  // IOobject::groupName("alpha", phaseName),
57  // fluid.mesh().time().name(),
58  // fluid.mesh(),
59  // IOobject::NO_READ,
60  // IOobject::AUTO_WRITE
61  // ),
62  // fluid.mesh(),
63  // dimensionedScalar(dimless, 0)
64  // )
65  // : volScalarField
66  // (
67  // IOobject
68  // (
69  // IOobject::groupName("alpha", phaseName),
70  // fluid.mesh().time().name(),
71  // fluid.mesh(),
72  // IOobject::MUST_READ,
73  // IOobject::AUTO_WRITE
74  // ),
75  // fluid.mesh()
76  // )
77  // ),
78  // Replaced by the following because the Clang compiler does not use
79  // std::move to transfer the result of the ternary operator
81  (
82  IOobject
83  (
84  IOobject::groupName("alpha", phaseName),
85  fluid.mesh().time().name(),
86  fluid.mesh(),
87  IOobject::NO_READ,
88  IOobject::AUTO_WRITE
89  ),
90  referencePhase
92  (
93  IOobject::groupName("alpha", phaseName),
94  fluid.mesh(),
96  )
98  (
99  new volScalarField
100  (
101  IOobject
102  (
103  IOobject::groupName("alpha", phaseName),
104  fluid.mesh().time().name(),
105  fluid.mesh(),
106  IOobject::MUST_READ,
107  IOobject::AUTO_WRITE,
108  false
109  ),
110  fluid.mesh()
111  )
112  )
113  ),
114 
115  fluid_(fluid),
116  name_(phaseName),
117  index_(index),
118  residualAlpha_
119  (
120  "residualAlpha",
121  dimless,
122  fluid.subDict(phaseName).lookup("residualAlpha")
123  ),
124  alphaMax_(fluid.subDict(phaseName).lookupOrDefault("alphaMax", 1.0))
125 {
126  diameterModel_ = diameterModel::New(fluid.subDict(phaseName), *this);
127 }
128 
129 
131 {
133  return autoPtr<phaseModel>(nullptr);
134 }
135 
136 
137 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
138 
140 {}
141 
142 
143 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
144 
146 {
147  return name_;
148 }
149 
150 
152 {
153  return name_;
154 }
155 
156 
158 {
159  return index_;
160 }
161 
162 
164 {
165  return fluid_;
166 }
167 
168 
170 {
171  return residualAlpha_;
172 }
173 
174 
175 Foam::scalar Foam::phaseModel::alphaMax() const
176 {
177  return alphaMax_;
178 }
179 
180 
182 {
183  return diameterModel_().d();
184 }
185 
186 
188 {
189  return diameterModel_();
190 }
191 
192 
194 {
195  diameterModel_->correct();
196 }
197 
198 
200 {}
201 
202 
204 {}
205 
206 
208 {}
209 
210 
212 {}
213 
214 
216 {}
217 
218 
220 {}
221 
222 
224 {}
225 
226 
228 {}
229 
230 
232 {}
233 
234 
236 {}
237 
238 
240 {
241  return diameterModel_->read(fluid_.subDict(name_));
242 }
243 
244 
246 {
247  surfaceScalarField::Boundary& alphaPhiBf = alphaPhi.boundaryFieldRef();
248  const volScalarField::Boundary& alphaBf = boundaryField();
249  const surfaceScalarField::Boundary& phiBf = phiRef().boundaryField();
250 
251  forAll(alphaPhiBf, patchi)
252  {
253  fvsPatchScalarField& alphaPhip = alphaPhiBf[patchi];
254 
255  if (!alphaPhip.coupled())
256  {
257  alphaPhip = phiBf[patchi]*alphaBf[patchi];
258  }
259  }
260 }
261 
262 
263 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricBoundaryField class.
Generic GeometricField class.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Abstract base-class for dispersed-phase particle diameter models.
Definition: diameterModel.H:52
static autoPtr< diameterModel > New(const dictionary &diameterProperties, const phaseModel &phase)
Select from dictionary and phase.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:83
virtual bool coupled() const
Return true if this patch field is coupled.
virtual void correctUf()
Correct the face velocity for moving meshes.
Definition: phaseModel.C:235
virtual void correctKinematics()
Correct the kinematics.
Definition: phaseModel.C:203
const phaseSystem & fluid() const
Return the system to which this phase belongs.
Definition: phaseModel.C:163
virtual void correctReactions()
Correct the reactions.
Definition: phaseModel.C:211
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
Definition: phaseModel.C:169
virtual void correctSpecies()
Correct the species concentrations.
Definition: phaseModel.C:215
virtual void correctThermo()
Correct the thermodynamics.
Definition: phaseModel.C:207
virtual void correct()
Correct the phase properties.
Definition: phaseModel.C:193
virtual ~phaseModel()
Destructor.
Definition: phaseModel.C:139
virtual void predictThermophysicalTransport()
Predict the energy transport.
Definition: phaseModel.C:223
virtual void correctMomentumTransport()
Correct the momentumTransport.
Definition: phaseModel.C:227
scalar alphaMax() const
Return the maximum phase-fraction (e.g. packing limit)
Definition: phaseModel.C:175
const diameterModel & diameter() const
Return a reference to the diameterModel of the phase.
Definition: phaseModel.C:187
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
Definition: phaseModel.C:181
virtual void correctThermophysicalTransport()
Correct the energy transport.
Definition: phaseModel.C:231
label index() const
Return the index of the phase.
Definition: phaseModel.C:157
const word & keyword() const
Return the name of the phase for use as the keyword in PtrDictionary.
Definition: phaseModel.C:151
phaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
Definition: phaseModel.C:42
autoPtr< phaseModel > clone() const
Return clone.
Definition: phaseModel.C:130
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
Definition: phaseModel.C:245
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
virtual void predictMomentumTransport()
Predict the momentumTransport.
Definition: phaseModel.C:219
virtual void correctContinuityError(const volScalarField &source)
Correct the continuity error.
Definition: phaseModel.C:199
virtual bool read()
Read phase properties dictionary.
Definition: phaseModel.C:239
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
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
label patchi
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Namespace for OpenFOAM.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
defineTypeNameAndDebug(combustionModel, 0)