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-2020 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 {
34  defineTypeNameAndDebug(phaseModel, 0);
35  defineRunTimeSelectionTable(phaseModel, phaseSystem);
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().timeName(),
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().timeName(),
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 
109 const Foam::word& Foam::phaseModel::name() const
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 
175 {}
176 
178 {}
179 
181 {}
182 
183 
185 {}
186 
187 
189 {
190  return diameterModel_->read(fluid_.subDict(name_));
191 }
192 
193 
195 {
196  surfaceScalarField::Boundary& alphaPhiBf = alphaPhi.boundaryFieldRef();
197  const volScalarField::Boundary& alphaBf = boundaryField();
198  const surfaceScalarField::Boundary& phiBf = phi()().boundaryField();
199 
200  forAll(alphaPhiBf, patchi)
201  {
202  fvsPatchScalarField& alphaPhip = alphaPhiBf[patchi];
203 
204  if (!alphaPhip.coupled())
205  {
206  alphaPhip = phiBf[patchi]*alphaBf[patchi];
207  }
208  }
209 }
210 
211 
212 // ************************************************************************* //
fvsPatchField< scalar > fvsPatchScalarField
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
label index() const
Return the index of the phase.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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 surfaceScalarField & phi() const
Return the volumetric flux.
virtual void correctKinematics()
Correct the kinematics.
autoPtr< phaseModel > clone() const
Return clone.
const word & name() const
Definition: phaseModel.H:109
static autoPtr< diameterModel > New(const dictionary &diameterProperties, const phaseModel &phase)
Generic dimensioned Type class.
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
phaseModel(const word &phaseName, const volScalarField &p, const volScalarField &T)
Construct from components.
virtual void correctSpecies()
Correct the species concentrations.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
virtual ~phaseModel()
Destructor.
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
virtual bool read()
Read phase properties dictionary.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:67
A class for handling words, derived from string.
Definition: word.H:59
const word & keyword() const
Definition: phaseModel.H:114
word timeName
Definition: getTimeIndex.H:3
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
scalar alphaMax() const
Return the maximum phase-fraction (e.g. packing limit)
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const phaseSystem & fluid() const
Return the system to which this phase belongs.
const autoPtr< diameterModel > & dPtr() const
Return const-reference to diameterModel of the phase.
virtual void correctThermo()
Correct the thermodynamics.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual void correctEnergyTransport()
Correct the energy transport.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
virtual void correctTurbulence()
Correct the turbulence.
virtual void correctReactions()
Correct the reactions.
Namespace for OpenFOAM.
virtual void correctContinuityError(const volScalarField &source)
Correct the continuity error.