StaticPhaseModel.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) 2017 OpenCFD Ltd.
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 "StaticPhaseModel.H"
29 
30 #include "phaseSystem.H"
31 
32 #include "fvcDdt.H"
33 #include "fvcDiv.H"
34 #include "surfaceInterpolate.H"
35 
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 template<class BasePhaseModel>
41 (
42  const phaseSystem& fluid,
43  const word& phaseName
44 )
45 :
46  BasePhaseModel(fluid, phaseName),
47  U_(fluid.mesh().lookupObject<volVectorField>("U")),
48  phi_
49  (
50  IOobject
51  (
52  IOobject::groupName("phi", phaseModel::name()),
53  fluid.mesh().time().timeName(),
54  fluid.mesh()
55  ),
56  fluid.mesh(),
57  dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero)
58  ),
59  alphaPhi_
60  (
61  IOobject
62  (
63  IOobject::groupName("alphaPhi", phaseModel::name()),
64  fluid.mesh().time().timeName(),
65  fluid.mesh()
66  ),
67  fluid.mesh(),
68  dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero)
69  )
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
75 template<class BasePhaseModel>
77 {
79 }
80 
81 
82 template<class BasePhaseModel>
85 {
87  (
88  IOobject
89  (
90  IOobject::groupName("phi", phaseModel::name()),
91  U_.mesh().time().timeName(),
92  U_.mesh()
93  ),
94  U_.mesh(),
95  dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero)
96  );
97 }
98 
99 
100 template<class BasePhaseModel>
103 {
104  phi_ = dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero);
105  return phi_;
106 }
107 
108 
109 template<class BasePhaseModel>
112 {
114  (
116  (
117  IOobject
118  (
119  IOobject::groupName("alphaPhi", phaseModel::name()),
120  U_.mesh().time().timeName(),
121  U_.mesh()
122  ),
123  U_.mesh(),
124  dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero)
125  )
126  );
127 }
128 
129 
130 template<class BasePhaseModel>
133 {
134  alphaPhi_ = dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero);
135  return alphaPhi_;
136 }
137 
138 
139 template<class BasePhaseModel>
142 {
143  return tmp<volVectorField>
144  (
145  new volVectorField
146  (
147  IOobject
148  (
149  IOobject::groupName("U", phaseModel::name()),
150  U_.mesh().time().timeName(),
151  U_.mesh()
152  ),
153  U_.mesh(),
155  )
156  );
157 }
158 
159 
160 template<class BasePhaseModel>
162 ::diffNo() const
163 {
164  tmp<surfaceScalarField> tkapparhoCpbyDelta
165  (
166  sqr(U_.mesh().surfaceInterpolation::deltaCoeffs())
167  *fvc::interpolate(this->kappa().ref())
168  /fvc::interpolate((this->Cp()*this->rho())())
169  );
170 
171  return tkapparhoCpbyDelta;
172 }
173 
174 
175 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimVelocity
const dimensionSet dimVelocity
fvcDiv.H
Calculate the divergence of the given field.
Foam::dimensionSet
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:108
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
Foam::constant::electromagnetic::kappa
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Foam::StaticPhaseModel::phi
virtual tmp< surfaceScalarField > phi() const
Constant access the volumetric flux. Return zero field.
Definition: StaticPhaseModel.C:84
rho
rho
Definition: readInitialConditions.H:88
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::StaticPhaseModel::correct
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
Definition: StaticPhaseModel.C:76
correct
fvOptions correct(rho)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::StaticPhaseModel::diffNo
virtual tmp< surfaceScalarField > diffNo() const
Maximum diffusion number.
Definition: StaticPhaseModel.C:162
Foam::StaticPhaseModel::U
virtual tmp< volVectorField > U() const
Access const reference to U.
Definition: StaticPhaseModel.C:141
StaticPhaseModel.H
Foam::StaticPhaseModel::StaticPhaseModel
StaticPhaseModel(const phaseSystem &fluid, const word &phaseName)
Definition: StaticPhaseModel.C:41
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::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Cp
const volScalarField & Cp
Definition: EEqn.H:7
Foam::StaticPhaseModel::alphaPhi
virtual tmp< surfaceScalarField > alphaPhi() const
Constant access the volumetric flux of the phase.
Definition: StaticPhaseModel.C:111
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
fvcDdt.H
Calculate the first temporal derivative.
Foam::phaseSystem
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:66
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.