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-------------------------------------------------------------------------------
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 "MovingPhaseModel.H"
29
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
46template<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 (
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
76template<class BasePhaseModel>
78{
79 BasePhaseModel::correct();
80}
81
82
83template<class BasePhaseModel>
86{
87 return tmp<surfaceScalarField>(phi_);
88}
89
90
91template<class BasePhaseModel>
94{
95 return phi_;
96}
97
98
99template<class BasePhaseModel>
102{
103 return tmp<surfaceScalarField>(alphaPhi_);
104}
105
106
107template<class BasePhaseModel>
110{
111 return alphaPhi_;
112}
113
114
115template<class BasePhaseModel>
118{
119 return tmp<volVectorField>(U_);
120}
121
122
123template<class BasePhaseModel>
125diffNo() const
126{
128 (
130 (
132 U_.mesh().time().timeName(),
133 U_.mesh()
134 ),
135 U_.mesh(),
137 );
138}
139
140
141// ************************************************************************* //
twoPhaseSystem & fluid
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Class which represents a moving fluid phase. Holds the velocity, fluxes and turbulence model....
virtual tmp< surfaceScalarField > diffNo() const
Diffusion number.
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.
virtual tmp< surfaceScalarField > phi() const
Constant access the volumetric flux.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
const word & name() const
The name of this phase.
Definition: phaseModel.H:114
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the finiteVolume matrix for implicit and explicit sources.
word timeName
Definition: getTimeIndex.H:3
const dimensionSet dimless
Dimensionless.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59