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) 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 "phaseSystem.H"
31 #include "diameterModel.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(phaseModel, 0);
38  defineRunTimeSelectionTable(phaseModel, phaseSystem);
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 
81 {
83  return autoPtr<phaseModel>(nullptr);
84 }
85 
86 
87 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
88 
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
96 {
97  return name_;
98 }
99 
100 
102 {
103  return name_;
104 }
105 
106 
107 Foam::label Foam::phaseModel::index() const
108 {
109  return index_;
110 }
111 
112 
114 {
115  return fluid_;
116 }
117 
118 
120 {
121  return residualAlpha_;
122 }
123 
124 
125 Foam::scalar Foam::phaseModel::alphaMax() const
126 {
127  return alphaMax_;
128 }
129 
130 
132 {
133  return diameterModel_().d();
134 }
135 
136 
138 {
139  return diameterModel_;
140 }
141 
142 
144 {
145  diameterModel_->correct();
146 }
147 
148 
150 {}
151 
152 
154 {}
155 
156 
158 {}
159 
160 
162 {}
163 
164 
166 {
167  return diameterModel_->read(fluid_.subDict(name_));
168 }
169 
170 
172 {
173  surfaceScalarField::Boundary& alphaPhiBf = alphaPhi.boundaryFieldRef();
174  const volScalarField::Boundary& alphaBf = boundaryField();
175  const surfaceScalarField::Boundary& phiBf = phi()().boundaryField();
176 
177  forAll(alphaPhiBf, patchi)
178  {
179  fvsPatchScalarField& alphaPhip = alphaPhiBf[patchi];
180 
181  if (!alphaPhip.coupled())
182  {
183  alphaPhip = phiBf[patchi]*alphaBf[patchi];
184  }
185  }
186 }
187 
188 
189 // ************************************************************************* //
Foam::phaseModel::correctTurbulence
virtual void correctTurbulence()
Correct the turbulence.
Definition: phaseModel.C:157
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:107
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:131
Foam::phaseModel::alphaMax
scalar alphaMax() const
Return the maximum phase-fraction (e.g. packing limit)
Definition: phaseModel.C:125
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::phaseModel::correctEnergyTransport
virtual void correctEnergyTransport()
Correct the energy transport.
Definition: phaseModel.C:161
Foam::phaseModel::keyword
const word & keyword() const
Return the name of the phase for use as the keyword in PtrDictionary.
Definition: phaseModel.C:101
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:143
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::phaseModel::~phaseModel
virtual ~phaseModel()
Destructor.
Definition: phaseModel.C:89
Foam::phaseModel::correctThermo
virtual void correctThermo()
Correct the thermodynamics.
Definition: phaseModel.C:153
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:436
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:137
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:95
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:171
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::phaseModel::read
virtual bool read()
Read phase properties dictionary.
Definition: phaseModel.C:165
Foam::phaseModel::fluid
const phaseSystem & fluid() const
Return the system to which this phase belongs.
Definition: phaseModel.C:113
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:80
Foam::phaseModel::residualAlpha
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
Definition: phaseModel.C:119
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:149
Foam::phaseModel::phaseModel
phaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
Definition: phaseModel.C:45