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-2022 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 "multiphaseInterSystem.H"
31 
32 #include "fvcDdt.H"
33 #include "fvcDiv.H"
34 #include "surfaceInterpolate.H"
35 
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 template<class BasePhaseModel>
41 (
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", multiphaseInter::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
64  (
65  "alphaPhi",
66  multiphaseInter::phaseModel::name()
67  ),
68  fluid.mesh().time().timeName(),
69  fluid.mesh()
70  ),
71  fluid.mesh(),
72  dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero)
73  )
74 {}
75 
76 
77 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
79 template<class BasePhaseModel>
81 {
83 }
84 
85 
86 template<class BasePhaseModel>
89 {
91  (
92  IOobject
93  (
95  U_.mesh().time().timeName(),
96  U_.mesh()
97  ),
98  U_.mesh(),
100  );
101 }
102 
103 
104 template<class BasePhaseModel>
107 {
108  phi_ = dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero);
109  return phi_;
110 }
111 
112 
113 template<class BasePhaseModel>
116 {
118  (
119  IOobject
120  (
122  (
123  "alphaPhi",
125  ),
126  U_.mesh().time().timeName(),
127  U_.mesh()
128  ),
129  U_.mesh(),
131  );
132 }
133 
134 
135 template<class BasePhaseModel>
138 {
139  alphaPhi_ = dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero);
140  return alphaPhi_;
141 }
142 
143 
144 template<class BasePhaseModel>
147 {
149  (
150  IOobject
151  (
153  U_.mesh().time().timeName(),
154  U_.mesh()
155  ),
156  U_.mesh(),
158  );
159 }
160 
161 
162 template<class BasePhaseModel>
164 ::diffNo() const
165 {
166  tmp<surfaceScalarField> tkapparhoCpbyDelta
167  (
168  sqr(U_.mesh().surfaceInterpolation::deltaCoeffs())
169  *fvc::interpolate(this->kappa().ref())
170  /fvc::interpolate((this->Cp()*this->rho())())
171  );
172 
173  return tkapparhoCpbyDelta;
174 }
175 
176 
177 // ************************************************************************* //
twoPhaseSystem & fluid
virtual tmp< volVectorField > U() const
Access const reference to U.
virtual tmp< surfaceScalarField > diffNo() const
Maximum diffusion number.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const word & name() const
The name of this phase.
Definition: phaseModel.H:122
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Calculate the first temporal derivative.
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
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual tmp< surfaceScalarField > phi() const
Constant access the volumetric flux. Return zero field.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
Calculate the divergence of the given field.
const volScalarField & Cp
Definition: EEqn.H:7
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2312/OpenFOAM-v2312/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
StaticPhaseModel(const multiphaseInterSystem &fluid, const word &phaseName)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
virtual tmp< surfaceScalarField > alphaPhi() const
Constant access the volumetric flux of the phase.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity