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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "phaseModel.H"
29 #include "phaseSystem.H"
30 #include "diameterModel.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(phaseModel, 0);
37  defineRunTimeSelectionTable(phaseModel, phaseSystem);
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const phaseSystem& fluid,
46  const word& phaseName,
47  const label index
48 )
49 :
51  (
52  IOobject
53  (
54  IOobject::groupName("alpha", phaseName),
55  fluid.mesh().time().timeName(),
56  fluid.mesh(),
57  IOobject::READ_IF_PRESENT,
58  IOobject::AUTO_WRITE
59  ),
60  fluid.mesh(),
61  dimensionedScalar("zero", dimless, 0)
62  ),
63 
64  fluid_(fluid),
65  name_(phaseName),
66  index_(index),
67  residualAlpha_
68  (
69  "residualAlpha",
70  dimless,
71  fluid.subDict(phaseName)
72  ),
73  alphaMax_(fluid.subDict(phaseName).lookupOrDefault<scalar>("alphaMax", 1))
74 {
75  diameterModel_ = diameterModel::New(fluid.subDict(phaseName), *this);
76 }
77 
78 
80 {
82  return autoPtr<phaseModel>(nullptr);
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
87 
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 
95 {
96  return name_;
97 }
98 
99 
101 {
102  return name_;
103 }
104 
105 
107 {
108  return index_;
109 }
110 
111 
113 {
114  return fluid_;
115 }
116 
117 
119 {
120  return residualAlpha_;
121 }
122 
123 
124 Foam::scalar Foam::phaseModel::alphaMax() const
125 {
126  return alphaMax_;
127 }
128 
129 
131 {
132  return diameterModel_().d();
133 }
134 
135 
137 {
138  return diameterModel_;
139 }
140 
141 
143 {
144  diameterModel_->correct();
145 }
146 
147 
149 {}
150 
151 
153 {}
154 
155 
157 {}
158 
159 
161 {}
162 
163 
165 {
166  return diameterModel_->read(fluid_.subDict(name_));
167 }
168 
169 
171 {
172  surfaceScalarField::Boundary& alphaPhiBf = alphaPhi.boundaryFieldRef();
173  const volScalarField::Boundary& alphaBf = boundaryField();
174  const surfaceScalarField::Boundary& phiBf = phi()().boundaryField();
175 
176  forAll(alphaPhiBf, patchi)
177  {
178  fvsPatchScalarField& alphaPhip = alphaPhiBf[patchi];
179 
180  if (!alphaPhip.coupled())
181  {
182  alphaPhip = phiBf[patchi]*alphaBf[patchi];
183  }
184  }
185 }
186 
187 
188 // ************************************************************************* //
Foam::phaseModel::correctTurbulence
virtual void correctTurbulence()
Correct the turbulence.
Definition: phaseModel.C:156
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::phaseModel::index
label index() const
Return the index of the phase.
Definition: phaseModel.C:106
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::phaseModel::d
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
Definition: phaseModel.C:130
Foam::phaseModel::alphaMax
scalar alphaMax() const
Return the maximum phase-fraction (e.g. packing limit)
Definition: phaseModel.C:124
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::phaseModel::correctEnergyTransport
virtual void correctEnergyTransport()
Correct the energy transport.
Definition: phaseModel.C:160
Foam::phaseModel::keyword
const word & keyword() const
Return the name of the phase for use as the keyword in PtrDictionary.
Definition: phaseModel.C:100
phaseModel.H
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
Foam::phaseModel::correct
virtual void correct()
Correct the phase properties.
Definition: phaseModel.C:142
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::phaseModel::~phaseModel
virtual ~phaseModel()
Destructor.
Definition: phaseModel.C:88
Foam::phaseModel::correctThermo
virtual void correctThermo()
Correct the thermodynamics.
Definition: phaseModel.C:152
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:419
Foam::label
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:62
alphaPhi
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::phaseModel::dPtr
const autoPtr< diameterModel > & dPtr() const
Return const-reference to diameterModel of the phase.
Definition: phaseModel.C:136
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::phaseModel::name
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:94
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::phaseModel::correctInflowOutflow
void correctInflowOutflow(surfaceScalarField &alphaPhi) const
Ensure that the flux at inflow/outflow BCs is preserved.
Definition: phaseModel.C:170
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:752
Foam::phaseModel::read
virtual bool read()
Read phase properties dictionary.
Definition: phaseModel.C:164
Foam::phaseModel::fluid
const phaseSystem & fluid() const
Return the system to which this phase belongs.
Definition: phaseModel.C:112
Foam::phaseSystem
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:69
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::phaseModel::clone
autoPtr< phaseModel > clone() const
Return clone.
Definition: phaseModel.C:79
Foam::phaseModel::residualAlpha
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
Definition: phaseModel.C:118
Foam::fvsPatchField::coupled
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvsPatchField.H:310
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::phaseModel::correctKinematics
virtual void correctKinematics()
Correct the kinematics.
Definition: phaseModel.C:148
Foam::phaseModel::phaseModel
phaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
Definition: phaseModel.C:44