Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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
31
32#include "fvcDdt.H"
33#include "fvcDiv.H"
34#include "surfaceInterpolate.H"
35
36
37// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38
39template<class BasePhaseModel>
41(
43 const word& phaseName
44)
45:
46 BasePhaseModel(fluid, phaseName),
47 U_(fluid.mesh().lookupObject<volVectorField>("U")),
48 phi_
49 (
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 (
62 (
63 IOobject::groupName
64 (
65 "alphaPhi",
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
79template<class BasePhaseModel>
82 BasePhaseModel::correct();
83}
84
85
86template<class BasePhaseModel>
89{
91 (
94 U_.mesh(),
96 );
97}
98
99
100template<class BasePhaseModel>
103{
104 phi_ = Zero;
105 return phi_;
106}
107
108
109template<class BasePhaseModel>
112{
114 (
116 (
117 "alphaPhi",
119 ),
121 U_.mesh(),
123 );
124}
125
126
127template<class BasePhaseModel>
130{
131 alphaPhi_ = Zero;
132 return alphaPhi_;
133}
134
135
136template<class BasePhaseModel>
139{
141 (
144 U_.mesh(),
146 );
147}
148
149
150template<class BasePhaseModel>
152::diffNo() const
153{
154 tmp<surfaceScalarField> tkapparhoCpbyDelta
155 (
156 sqr(U_.mesh().surfaceInterpolation::deltaCoeffs())
157 *fvc::interpolate(this->kappa().ref())
158 /fvc::interpolate((this->Cp()*this->rho())())
159 );
160
161 return tkapparhoCpbyDelta;
162}
163
164
165// ************************************************************************* //
twoPhaseSystem & fluid
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvsPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Class which represents a static fluid phase.
virtual void correct()
Correct the phase properties other than the thermo and turbulence.
virtual tmp< volVectorField > U() const
Access const reference to U.
virtual tmp< surfaceScalarField > alphaPhi() const
Constant access the volumetric flux of the phase.
StaticPhaseModel(const multiphaseInterSystem &fluid, const word &phaseName)
virtual tmp< surfaceScalarField > phi() const
Constant access the volumetric flux. Return zero field.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
const word & name() const
The name of this phase.
Definition phaseModel.H:122
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
dynamicFvMesh & mesh
Calculate the first temporal derivative.
Calculate the divergence of the given field.
word timeName
Definition getTimeIndex.H:3
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.
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimVelocity
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.