StationaryPhaseModel.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) 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
26\*---------------------------------------------------------------------------*/
27
29
30// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
31
32template<class BasePhaseModel>
33template<class Type, template<class> class PatchField, class GeoMesh>
36(
37 const word& name,
38 const dimensionSet& dims,
39 const bool cache
40) const
41{
43 (
45 (
47 (
48 IOobject::groupName(name, this->name()),
49 this->mesh().time().timeName(),
50 this->mesh()
51 ),
52 this->mesh(),
53 dimensioned<Type>(dims, Zero)
54 )
55 );
56}
57
58
59template<class BasePhaseModel>
60template<class Type>
63(
64 const word& name,
65 const dimensionSet& dims,
66 const bool cache
67) const
68{
70}
71
72
73template<class BasePhaseModel>
74template<class Type>
77(
78 const word& name,
79 const dimensionSet& dims,
80 const bool cache
81) const
82{
84}
85
86
87// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
88
89template<class BasePhaseModel>
91(
92 const phaseSystem& fluid,
93 const word& phaseName,
94 const label index
95)
96:
97 BasePhaseModel(fluid, phaseName, index)
98{}
99
100
101// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
102
103template<class BasePhaseModel>
105{}
106
107
108// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
109
110template<class BasePhaseModel>
112{
113 return true;
114}
115
116
117template<class BasePhaseModel>
120{
122 << "Cannot construct a momentum equation for a stationary phase"
123 << exit(FatalError);
124
125 return tmp<fvVectorMatrix>();
126}
127
128
129template<class BasePhaseModel>
132{
134 << "Cannot construct a momentum equation for a stationary phase"
135 << exit(FatalError);
136
137 return tmp<fvVectorMatrix>();
138}
139
140
141template<class BasePhaseModel>
144{
145 return zeroVolField<vector>("U", dimVelocity, true);
146}
147
148
149template<class BasePhaseModel>
152{
154 << "Cannot access the velocity of a stationary phase"
155 << exit(FatalError);
156
157 return const_cast<volVectorField&>(volVectorField::null());
158}
159
160
161template<class BasePhaseModel>
164{
165 return zeroSurfaceField<scalar>("phi", dimVolume/dimTime);
166}
167
168
169template<class BasePhaseModel>
172{
174 << "Cannot access the flux of a stationary phase"
175 << exit(FatalError);
176
177 return const_cast<surfaceScalarField&>(surfaceScalarField::null());
178}
179
180
181template<class BasePhaseModel>
184{
185 return zeroSurfaceField<scalar>("alphaPhi", dimVolume/dimTime);
186}
187
188
189template<class BasePhaseModel>
192{
194 << "Cannot access the volumetric flux of a stationary phase"
195 << exit(FatalError);
196
197 return const_cast<surfaceScalarField&>(surfaceScalarField::null());
198}
199
200
201template<class BasePhaseModel>
204{
205 return zeroSurfaceField<scalar>("alphaRhoPhi", dimMass/dimTime);
206}
207
208
209template<class BasePhaseModel>
212{
214 << "Cannot access the mass flux of a stationary phase"
215 << exit(FatalError);
216
217 return const_cast<surfaceScalarField&>(surfaceScalarField::null());
218}
219
220
221template<class BasePhaseModel>
224{
225 return zeroVolField<vector>("DUDt", dimVelocity/dimTime);
226}
227
228
229template<class BasePhaseModel>
232{
233 return zeroSurfaceField<scalar>("DUDtf", dimVelocity*dimArea/dimTime);
234}
235
236
237template<class BasePhaseModel>
240{
241 return zeroVolField<scalar>("continuityError", dimDensity/dimTime);
242}
243
244
245template<class BasePhaseModel>
248{
249 return zeroVolField<scalar>("continuityErrorFlow", dimDensity/dimTime);
250}
251
252
253template<class BasePhaseModel>
256{
257 return zeroVolField<scalar>("continuityErrorSources", dimDensity/dimTime);
258}
259
260
261template<class BasePhaseModel>
264{
265 return zeroVolField<scalar>("K", sqr(dimVelocity));
266}
267
268
269template<class BasePhaseModel>
272{
273 return tmp<volScalarField>();
274}
275
276
277template<class BasePhaseModel>
279(
281)
282{
284 << "Cannot set the dilatation rate of a stationary phase"
285 << exit(FatalError);
286}
287
288
289template<class BasePhaseModel>
292{
293 return zeroVolField<scalar>("continuityError", dimDynamicViscosity);
294}
295
296
297template<class BasePhaseModel>
300{
301 return this->thermo().mu();
302}
303
304
305template<class BasePhaseModel>
308{
309 return zeroVolField<scalar>("continuityError", dimViscosity);
310}
311
312
313template<class BasePhaseModel>
316{
317 return this->thermo().nu();
318}
319
320
321template<class BasePhaseModel>
324{
325 return this->thermo().kappa();
326}
327
328
329template<class BasePhaseModel>
332{
333 return this->thermo().kappa(patchi);
334}
335
336
337template<class BasePhaseModel>
340{
341 return this->thermo().alpha();
342}
343
344
345template<class BasePhaseModel>
348{
349 return this->thermo().alpha(patchi);
350}
351
352
353template<class BasePhaseModel>
356{
357 return zeroVolField<scalar>("k", sqr(dimVelocity));
358}
359
360
361template<class BasePhaseModel>
364{
365 return zeroVolField<scalar>("pPrime", dimPressure);
366}
367
368
369// ************************************************************************* //
twoPhaseSystem & fluid
Generic GeometricField class.
static const GeometricField< vector, fvPatchField, volMesh > & null()
Return a null geometric field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
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
Generic dimensioned Type class.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:76
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
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
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
zeroField divU
Definition: alphaSuSp.H:3
word timeName
Definition: getTimeIndex.H:3
const dimensionSet dimPressure
const dimensionSet dimViscosity
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
const dimensionSet dimVelocity
const dimensionSet dimDynamicViscosity
error FatalError
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51