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-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 "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 :
50  (
51  referencePhase
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(),
64  )
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 
79  fluid_(fluid),
80  name_(phaseName),
81  index_(index),
82  residualAlpha_
83  (
84  "residualAlpha",
85  dimless,
86  fluid.subDict(phaseName).lookup("residualAlpha")
87  ),
88  alphaMax_(fluid.subDict(phaseName).lookupOrDefault("alphaMax", 1.0))
89 {
90  diameterModel_ = diameterModel::New(fluid.subDict(phaseName), *this);
91 }
92 
93 
95 {
97  return autoPtr<phaseModel>(nullptr);
98 }
99 
100 
101 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
102 
104 {}
105 
106 
107 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108 
110 {
111  return name_;
112 }
113 
114 
116 {
117  return name_;
118 }
119 
120 
122 {
123  return index_;
124 }
125 
126 
128 {
129  return fluid_;
130 }
131 
132 
134 {
135  return residualAlpha_;
136 }
137 
138 
139 Foam::scalar Foam::phaseModel::alphaMax() const
140 {
141  return alphaMax_;
142 }
143 
144 
146 {
147  return diameterModel_().d();
148 }
149 
150 
152 {
153  return diameterModel_;
154 }
155 
156 
158 {
159  diameterModel_->correct();
160 }
161 
162 
164 {}
165 
166 
168 {}
169 
170 
172 {}
173 
174 
176 {}
177 
178 
180 {}
181 
182 
184 {}
185 
186 
188 {}
189 
190 
192 {}
193 
194 
196 {}
197 
198 
200 {}
201 
202 
204 {
205  return diameterModel_->read(fluid_.subDict(name_));
206 }
207 
208 
210 {
211  surfaceScalarField::Boundary& alphaPhiBf = alphaPhi.boundaryFieldRef();
212  const volScalarField::Boundary& alphaBf = boundaryField();
213  const surfaceScalarField::Boundary& phiBf = phiRef().boundaryField();
214 
215  forAll(alphaPhiBf, patchi)
216  {
217  fvsPatchScalarField& alphaPhip = alphaPhiBf[patchi];
218 
219  if (!alphaPhip.coupled())
220  {
221  alphaPhip = phiBf[patchi]*alphaBf[patchi];
222  }
223  }
224 }
225 
226 
227 // ************************************************************************* //
#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
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:998
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:82
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:199
virtual void correctKinematics()
Correct the kinematics.
Definition: phaseModel.C:167
const phaseSystem & fluid() const
Return the system to which this phase belongs.
Definition: phaseModel.C:127
virtual void correctReactions()
Correct the reactions.
Definition: phaseModel.C:175
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
Definition: phaseModel.C:133
virtual void correctSpecies()
Correct the species concentrations.
Definition: phaseModel.C:179
virtual void correctThermo()
Correct the thermodynamics.
Definition: phaseModel.C:171
virtual void correct()
Correct the phase properties.
Definition: phaseModel.C:157
virtual ~phaseModel()
Destructor.
Definition: phaseModel.C:103
virtual void predictThermophysicalTransport()
Predict the energy transport.
Definition: phaseModel.C:187
virtual void correctMomentumTransport()
Correct the momentumTransport.
Definition: phaseModel.C:191
scalar alphaMax() const
Return the maximum phase-fraction (e.g. packing limit)
Definition: phaseModel.C:139
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
Definition: phaseModel.C:145
virtual void correctThermophysicalTransport()
Correct the energy transport.
Definition: phaseModel.C:195
label index() const
Return the index of the phase.
Definition: phaseModel.C:121
const word & keyword() const
Return the name of the phase for use as the keyword in PtrDictionary.
Definition: phaseModel.C:115
phaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
Definition: phaseModel.C:42
const autoPtr< diameterModel > & dPtr() const
Return const-reference to diameterModel of the phase.
Definition: phaseModel.C:151
autoPtr< phaseModel > clone() const
Return clone.
Definition: phaseModel.C:94
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
Definition: phaseModel.C:209
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:109
virtual void predictMomentumTransport()
Predict the momentumTransport.
Definition: phaseModel.C:183
virtual void correctContinuityError(const volScalarField &source)
Correct the continuity error.
Definition: phaseModel.C:163
virtual bool read()
Read phase properties dictionary.
Definition: phaseModel.C:203
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:353
label patchi
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
const dimensionSet dimless
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47