StationaryPhaseModel.H
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) 2015-2018 OpenFOAM Foundation
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
26Class
27 Foam::StationaryPhaseModel
28
29Description
30 Class which represents a stationary (and therefore probably solid) phase.
31 Generates, but does not store, zero velocity and flux field and turbulent
32 qauantities. Throws an error when non-const access is requested to the
33 motion fields or when the momentum equation is requested. Usage must
34 must protect against such calls.
35
36See also
37 MovingPhaseModel
38
39SourceFiles
40 StationaryPhaseModel.C
41
42\*---------------------------------------------------------------------------*/
43
44#ifndef StationaryPhaseModel_H
45#define StationaryPhaseModel_H
46
47#include "phaseModel.H"
48
49// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50
51namespace Foam
52{
53
54/*---------------------------------------------------------------------------*\
55 Class StationaryPhaseModel Declaration
56\*---------------------------------------------------------------------------*/
57
58template<class BasePhaseModel>
60:
61 public BasePhaseModel
62{
63private:
64
65 // Private member functions
66
67 //- Create a zero geometric field
68 template<class Type, template<class> class PatchField, class GeoMesh>
70 (
71 const word& name,
72 const dimensionSet& dims,
73 const bool cache = false
74 ) const;
75
76 //- Create a zero vol field
77 template<class Type>
79 (
80 const word& name,
81 const dimensionSet& dims,
82 const bool cache = false
83 ) const;
84
85 //- Create a zero surface field
86 template<class Type>
88 (
89 const word& name,
90 const dimensionSet& dims,
91 const bool cache = false
92 ) const;
93
94
95public:
96
97 // Constructors
98
100 (
101 const phaseSystem& fluid,
102 const word& phaseName,
103 const label index
104 );
105
106
107 //- Destructor
108 virtual ~StationaryPhaseModel();
109
110
111 // Member Functions
112
113 // Momentum
114
115 //- Return whether the phase is stationary
116 virtual bool stationary() const;
117
118 //- Return the momentum equation
119 virtual tmp<fvVectorMatrix> UEqn();
120
121 //- Return the momentum equation for the face-based algorithm
122 virtual tmp<fvVectorMatrix> UfEqn();
123
124 //- Return the velocity
125 virtual tmp<volVectorField> U() const;
126
127 //- Access the velocity
128 virtual volVectorField& URef();
129
130 //- Return the volumetric flux
131 virtual tmp<surfaceScalarField> phi() const;
132
133 //- Access the volumetric flux
134 virtual surfaceScalarField& phiRef();
135
136 //- Return the volumetric flux of the phase
137 virtual tmp<surfaceScalarField> alphaPhi() const;
138
139 //- Access the volumetric flux of the phase
141
142 //- Return the mass flux of the phase
143 virtual tmp<surfaceScalarField> alphaRhoPhi() const;
144
145 //- Access the mass flux of the phase
147
148 //- Return the substantive acceleration
149 virtual tmp<volVectorField> DUDt() const;
150
151 //- Return the substantive acceleration on the faces
152 virtual tmp<surfaceScalarField> DUDtf() const;
153
154 //- Return the continuity error
155 virtual tmp<volScalarField> continuityError() const;
156
157 //- Return the continuity error due to the flow field
159
160 //- Return the continuity error due to any sources
162
163 //- Return the phase kinetic energy
164 virtual tmp<volScalarField> K() const;
165
166
167 // Compressibility (variable density)
168
169 //- Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
170 virtual tmp<volScalarField> divU() const;
171
172 //- Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
173 virtual void divU(tmp<volScalarField> divU);
174
175
176 // Turbulence
177
178 //- Return the turbulent dynamic viscosity
179 virtual tmp<volScalarField> mut() const;
180
181 //- Return the effective dynamic viscosity
182 virtual tmp<volScalarField> muEff() const;
183
184 //- Return the turbulent kinematic viscosity
185 virtual tmp<volScalarField> nut() const;
186
187 //- Return the effective kinematic viscosity
188 virtual tmp<volScalarField> nuEff() const;
189
190 //- Effective thermal turbulent diffusivity for temperature
191 // of mixture for patch [J/m/s/K]
192 using BasePhaseModel::kappaEff;
193
194 //- Return the effective thermal conductivity
195 virtual tmp<volScalarField> kappaEff() const;
196
197 //- Return the effective thermal conductivity on a patch
198 virtual tmp<scalarField> kappaEff(const label patchi) const;
199
200 //- Effective thermal turbulent diffusivity of mixture [kg/m/s]
201 using BasePhaseModel::alphaEff;
202
203 //- Return the effective thermal diffusivity
204 virtual tmp<volScalarField> alphaEff() const;
205
206 //- Return the effective thermal conductivity on a patch
207 virtual tmp<scalarField> alphaEff(const label patchi) const;
208
209 //- Return the turbulent kinetic energy
210 virtual tmp<volScalarField> k() const;
211
212 //- Return the phase-pressure'
213 // (derivative of phase-pressure w.r.t. phase-fraction)
214 virtual tmp<volScalarField> pPrime() const;
215};
216
217
218// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219
220} // End namespace Foam
221
222// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223
224#ifdef NoRepository
225 #include "StationaryPhaseModel.C"
226#endif
227
228// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229
230#endif
231
232// ************************************************************************* //
twoPhaseSystem & fluid
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:49
Class which represents a stationary (and therefore probably solid) phase. Generates,...
virtual tmp< volScalarField > continuityErrorSources() const
Return the continuity error due to any sources.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
virtual tmp< volScalarField > muEff() const
Return the effective dynamic viscosity.
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
virtual tmp< volVectorField > U() const
Return the velocity.
virtual tmp< volScalarField > mut() const
Return the turbulent dynamic viscosity.
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
virtual tmp< volScalarField > alphaEff() const
Return the effective thermal diffusivity.
virtual tmp< surfaceScalarField > alphaPhi() const
Return the volumetric flux of the phase.
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
virtual tmp< volScalarField > continuityErrorFlow() const
Return the continuity error due to the flow field.
virtual tmp< volScalarField > nuEff() const
Return the effective kinematic viscosity.
virtual ~StationaryPhaseModel()
Destructor.
virtual tmp< volScalarField > divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux.
virtual tmp< fvVectorMatrix > UfEqn()
Return the momentum equation for the face-based algorithm.
virtual tmp< volScalarField > nut() const
Return the turbulent kinematic viscosity.
virtual volVectorField & URef()
Access the velocity.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
virtual tmp< volScalarField > kappaEff() const
Return the effective thermal conductivity.
virtual bool stationary() const
Return whether the phase is stationary.
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:76
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Definition: zeroField.H:56
Namespace for OpenFOAM.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59