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) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2020 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 "twoPhaseSystem.H"
31 #include "diameterModel.H"
32 #include "fvMatrix.H"
34 #include "dragModel.H"
35 #include "heatTransferModel.H"
38 #include "slipFvPatchFields.H"
40 #include "fvcFlux.H"
41 #include "surfaceInterpolate.H"
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const twoPhaseSystem& fluid,
50  const word& phaseName
51 )
52 :
54  (
55  IOobject
56  (
57  IOobject::groupName("alpha", phaseName),
58  fluid.mesh().time().timeName(),
59  fluid.mesh(),
60  IOobject::READ_IF_PRESENT,
61  IOobject::AUTO_WRITE
62  ),
63  fluid.mesh(),
65  ),
66  fluid_(fluid),
67  name_(phaseName),
68  phaseDict_
69  (
70  phaseProperties.subDict(name_)
71  ),
72  residualAlpha_
73  (
74  "residualAlpha",
75  dimless,
76  fluid.subDict(phaseName)
77  ),
78  alphaMax_(phaseDict_.getOrDefault<scalar>("alphaMax", 1)),
79  thermo_(rhoThermo::New(fluid.mesh(), name_)),
80  U_
81  (
82  IOobject
83  (
84  IOobject::groupName("U", name_),
85  fluid.mesh().time().timeName(),
86  fluid.mesh(),
87  IOobject::MUST_READ,
88  IOobject::AUTO_WRITE
89  ),
90  fluid.mesh()
91  ),
92  alphaPhi_
93  (
94  IOobject
95  (
96  IOobject::groupName("alphaPhi", name_),
97  fluid.mesh().time().timeName(),
98  fluid.mesh()
99  ),
100  fluid.mesh(),
101  dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero)
102  ),
103  alphaRhoPhi_
104  (
105  IOobject
106  (
107  IOobject::groupName("alphaRhoPhi", name_),
108  fluid.mesh().time().timeName(),
109  fluid.mesh()
110  ),
111  fluid.mesh(),
112  dimensionedScalar(dimensionSet(1, 0, -1, 0, 0), Zero)
113  )
114 {
115  alphaPhi_.setOriented();
116  alphaRhoPhi_.setOriented();
117 
118  thermo_->validate("phaseModel " + name_, "h", "e");
119 
120  const word phiName = IOobject::groupName("phi", name_);
121 
122  IOobject phiHeader
123  (
124  phiName,
125  fluid_.mesh().time().timeName(),
126  fluid_.mesh(),
128  );
129 
130  if (phiHeader.typeHeaderOk<surfaceScalarField>(true))
131  {
132  Info<< "Reading face flux field " << phiName << endl;
133 
134  phiPtr_.reset
135  (
137  (
138  IOobject
139  (
140  phiName,
141  fluid_.mesh().time().timeName(),
142  fluid_.mesh(),
145  ),
146  fluid_.mesh()
147  )
148  );
149  }
150  else
151  {
152  Info<< "Calculating face flux field " << phiName << endl;
153 
154  wordList phiTypes
155  (
156  U_.boundaryField().size(),
158  );
159 
160  forAll(U_.boundaryField(), i)
161  {
162  if
163  (
164  isA<fixedValueFvPatchVectorField>(U_.boundaryField()[i])
165  || isA<slipFvPatchVectorField>(U_.boundaryField()[i])
166  || isA<partialSlipFvPatchVectorField>(U_.boundaryField()[i])
167  )
168  {
169  phiTypes[i] = fixedValueFvsPatchScalarField::typeName;
170  }
171  }
172 
173  phiPtr_.reset
174  (
176  (
177  IOobject
178  (
179  phiName,
180  fluid_.mesh().time().timeName(),
181  fluid_.mesh(),
184  ),
185  fvc::flux(U_),
186  phiTypes
187  )
188  );
189  }
190 
191  dPtr_ = diameterModel::New
192  (
193  phaseDict_,
194  *this
195  );
196 
197  turbulence_ =
199  (
200  *this,
201  thermo_->rho(),
202  U_,
203  alphaRhoPhi_,
204  phi(),
205  *this
206  );
207 }
208 
209 
210 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
213 {}
214 
215 
216 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
217 
219 {
220  return fluid_.otherPhase(*this);
221 }
222 
223 
225 {
226  return dPtr_().d();
227 }
228 
229 
232 {
233  return *turbulence_;
234 }
235 
236 
239 {
240  return *turbulence_;
241 }
242 
243 
245 {
246  return dPtr_->correct();
247 }
248 
249 
250 bool Foam::phaseModel::read(const dictionary& phaseProperties)
251 {
252  phaseDict_ = phaseProperties.subDict(name_);
253  return dPtr_->read(phaseDict_);
254 }
255 
256 
258 {
260  const volScalarField::Boundary& alphaBf = boundaryField();
261  const surfaceScalarField::Boundary& phiBf = phi().boundaryField();
262 
263  forAll(alphaPhiBf, patchi)
264  {
265  fvsPatchScalarField& alphaPhip = alphaPhiBf[patchi];
266 
267  if (!alphaPhip.coupled())
268  {
269  alphaPhip = phiBf[patchi]*alphaBf[patchi];
270  }
271  }
272 }
273 
274 
275 // ************************************************************************* //
twoPhaseSystem & fluid
fvsPatchField< scalar > fvsPatchScalarField
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
Definition: phaseModel.C:233
const phaseModel & otherPhase() const
Return the other phase in this two-phase system.
Definition: phaseModel.C:211
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Helper class to manage multi-specie phase properties.
const PhaseCompressibleTurbulenceModel< phaseModel > & turbulence() const
Return the turbulence model.
Definition: phaseModel.C:231
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
Type of boundary fields.
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
const surfaceScalarField & phi() const
Definition: phaseModel.H:218
Class which solves the volume fraction equations for two phases.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
word timeName
Definition: getTimeIndex.H:3
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
Calculate the face-flux of the given field.
static autoPtr< PhaseCompressibleTurbulenceModel > New(const alphaField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transportModel, const word &propertiesName=turbulenceModel::propertiesName)
Return a reference to the selected turbulence model.
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
static autoPtr< diameterModel > New(const dictionary &diameterProperties, const phaseModel &phase)
Definition: diameterModel.C:50
void setOriented(bool on=true) noexcept
Set the oriented flag: on/off.
void correct()
Correct the phase properties.
Definition: phaseModel.C:209
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
phaseModel(const word &phaseName, const dictionary &phaseDict, const fvMesh &mesh)
Definition: phaseModel.C:33
Nothing to be read.
Automatically write from objectRegistry::writeObject()
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Basic thermodynamic properties based on density.
Definition: rhoThermo.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
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
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
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:23
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
virtual ~phaseModel()
Destructor.
Definition: phaseModel.C:195
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127