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) 2015-2019 OpenFOAM Foundation
9  Copyright (C) 2019-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "phaseModel.H"
30 #include "phaseSystem.H"
31 #include "diameterModel.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(phaseModel, 0);
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const phaseSystem& fluid,
47  const word& phaseName,
48  const label index
49 )
50 :
52  (
53  IOobject
54  (
55  IOobject::groupName("alpha", phaseName),
56  fluid.mesh().time().timeName(),
57  fluid.mesh(),
58  IOobject::READ_IF_PRESENT,
59  IOobject::AUTO_WRITE
60  ),
61  fluid.mesh(),
62  dimensionedScalar("zero", dimless, 0)
63  ),
64 
65  fluid_(fluid),
66  name_(phaseName),
67  index_(index),
68  residualAlpha_
69  (
70  "residualAlpha",
71  dimless,
72  fluid.subDict(phaseName)
73  ),
74  alphaMax_(fluid.subDict(phaseName).getOrDefault<scalar>("alphaMax", 1))
75 {
76  diameterModel_ = diameterModel::New(fluid.subDict(phaseName), *this);
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
81 
84 {
86  return nullptr;
87 }
88 
89 
92 (
93  const phaseSystem& fluid,
94  const word& phaseName,
95  const label index
96 )
97 {
98  const dictionary& dict = fluid.subDict(phaseName);
99 
100  const word modelType(dict.get<word>("type"));
101 
102  Info<< "Selecting phaseModel for "
103  << phaseName << ": " << modelType << endl;
104 
105  auto* ctorPtr = phaseSystemConstructorTable(modelType);
106 
107  if (!ctorPtr)
108  {
110  (
111  dict,
112  "phaseModel",
113  modelType,
114  *phaseSystemConstructorTablePtr_
115  ) << abort(FatalIOError);
116  }
117 
118  return ctorPtr(fluid, phaseName, index);
119 }
120 
121 
122 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
123 
125 {}
126 
127 
128 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129 
130 const Foam::word& Foam::phaseModel::name() const
131 {
132  return name_;
133 }
134 
137 {
138  return name_;
139 }
140 
142 Foam::label Foam::phaseModel::index() const
143 {
144  return index_;
145 }
146 
149 {
150  return fluid_;
151 }
152 
155 {
156  return residualAlpha_;
157 }
158 
159 
160 Foam::scalar Foam::phaseModel::alphaMax() const
161 {
162  return alphaMax_;
163 }
164 
167 {
168  return diameterModel_().d();
169 }
170 
171 
173 {
174  return diameterModel_;
175 }
176 
179 {
180  diameterModel_->correct();
181 }
182 
183 
185 {}
186 
187 
189 {}
190 
191 
193 {}
194 
195 
197 {}
198 
199 
201 {
202  return diameterModel_->read(fluid_.subDict(name_));
203 }
204 
205 
207 {
209  const volScalarField::Boundary& alphaBf = boundaryField();
210  const surfaceScalarField::Boundary& phiBf = phi()().boundaryField();
211 
212  forAll(alphaPhiBf, patchi)
213  {
214  fvsPatchScalarField& alphaPhip = alphaPhiBf[patchi];
215 
216  if (!alphaPhip.coupled())
217  {
218  alphaPhip = phiBf[patchi]*alphaBf[patchi];
219  }
220  }
221 }
222 
223 
224 // ************************************************************************* //
twoPhaseSystem & fluid
fvsPatchField< scalar > fvsPatchScalarField
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
Definition: phaseModel.C:233
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
autoPtr< phaseModel > clone() const
Return clone.
Definition: phaseModel.C:201
const phaseSystem & fluid() const
Return the system to which this phase belongs.
Definition: phaseModel.C:141
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const word & keyword() const
Definition: phaseModel.H:171
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
Type of boundary fields.
const dimensionSet dimless
Dimensionless.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
word timeName
Definition: getTimeIndex.H:3
dynamicFvMesh & mesh
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:69
scalar alphaMax() const
Return the maximum phase-fraction (e.g. packing limit)
Definition: phaseModel.C:153
A class for handling words, derived from Foam::string.
Definition: word.H:63
static autoPtr< diameterModel > New(const dictionary &diameterProperties, const phaseModel &phase)
Definition: diameterModel.C:50
static autoPtr< phaseModel > New(const phaseSystem &fluid, const word &phaseName, const label index)
Definition: phaseModel.C:85
errorManip< error > abort(error &err)
Definition: errorManip.H:139
void correct()
Correct the phase properties.
Definition: phaseModel.C:209
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
virtual void correctThermo()
Correct the thermodynamics.
Definition: phaseModel.C:181
virtual void correctKinematics()
Correct the kinematics.
Definition: phaseModel.C:177
label index() const
Return the index of the phase.
Definition: phaseModel.C:135
const autoPtr< diameterModel > & dPtr() const
Return const-reference to diameterModel of the phase.
Definition: phaseModel.C:165
phaseModel(const word &phaseName, const dictionary &phaseDict, const fvMesh &mesh)
Definition: phaseModel.C:33
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
Definition: phaseModel.C:147
const word & name() const
Definition: phaseModel.H:166
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual void correctEnergyTransport()
Correct the energy transport.
Definition: phaseModel.C:189
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
virtual bool read()
Read phase properties dictionary.
Definition: phaseModel.C:193
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
tmp< volScalarField > d() const
Definition: phaseModel.C:251
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:635
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
virtual void correctTurbulence()
Correct the turbulence.
Definition: phaseModel.C:185
Namespace for OpenFOAM.
virtual ~phaseModel()
Destructor.
Definition: phaseModel.C:195
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...