MovingPhaseModel.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 "MovingPhaseModel.H"
29 
30 #include "multiphaseInterSystem.H"
31 
33 #include "slipFvPatchFields.H"
35 
36 #include "fvmDdt.H"
37 #include "fvmDiv.H"
38 
39 #include "fvmSup.H"
40 #include "fvcDdt.H"
41 #include "fvcDiv.H"
42 #include "surfaceInterpolate.H"
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 template<class BasePhaseModel>
48 (
50  const word& phaseName
51 )
52 :
53  BasePhaseModel(fluid, phaseName),
54  U_(fluid.mesh().lookupObject<volVectorField>("U")),
55  phi_(fluid.mesh().lookupObject<surfaceScalarField>("phi")),
56  alphaPhi_
57  (
58  IOobject
59  (
60  IOobject::groupName
61  (
62  "alphaPhi",
63  multiphaseInter::phaseModel::name()
64  ),
65  fluid.mesh().time().timeName(),
66  fluid.mesh()
67  ),
68  fluid.mesh(),
69  dimensionedScalar(dimensionSet(0, 3, -1, 0, 0), Zero)
70  )
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
76 template<class BasePhaseModel>
78 {
80 }
81 
82 
83 template<class BasePhaseModel>
86 {
87  return tmp<surfaceScalarField>(phi_);
88 }
89 
90 
91 template<class BasePhaseModel>
94 {
95  return phi_;
96 }
97 
98 
99 template<class BasePhaseModel>
102 {
103  return tmp<surfaceScalarField>(alphaPhi_);
104 }
105 
106 
107 template<class BasePhaseModel>
110 {
111  return alphaPhi_;
112 }
113 
114 
115 template<class BasePhaseModel>
118 {
119  return tmp<volVectorField>(U_);
120 }
121 
122 
123 template<class BasePhaseModel>
125 diffNo() const
126 {
128  (
129  IOobject
130  (
132  U_.mesh().time().timeName(),
133  U_.mesh()
134  ),
135  U_.mesh(),
137  );
138 }
139 
140 
141 // ************************************************************************* //
virtual tmp< volVectorField > U() const
Access const reference to U.
twoPhaseSystem & fluid
MovingPhaseModel(const multiphaseInterSystem &fluid, const word &phaseName)
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
const word & name() const
The name of this phase.
Definition: phaseModel.H:122
const dimensionSet dimless
Dimensionless.
virtual tmp< surfaceScalarField > diffNo() const
Diffusion number.
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
Calculate the matrix for the first temporal derivative.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
Calculate the divergence of the given field.
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
Calculate the matrix for the divergence of the given field and flux.
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
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Calculate the finiteVolume matrix for implicit and explicit sources.
virtual tmp< surfaceScalarField > phi() const
Constant access the volumetric flux.
virtual tmp< surfaceScalarField > alphaPhi() const
Constant access the volumetric flux of the phase.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127