StationaryPhaseModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "StationaryPhaseModel.H"
27 
28 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
29 
30 template<class BasePhaseModel>
31 template<class Type, template<class> class PatchField, class GeoMesh>
34 (
35  const word& name,
36  const dimensionSet& dims,
37  const bool cache
38 ) const
39 {
40  return tmp<GeometricField<Type, PatchField, GeoMesh>>
41  (
42  new GeometricField<Type, PatchField, GeoMesh>
43  (
44  IOobject
45  (
46  IOobject::groupName(name, this->name()),
47  this->mesh().time().timeName(),
48  this->mesh()
49  ),
50  this->mesh(),
51  dimensioned<Type>("zero", dims, pTraits<Type>::zero)
52  )
53  );
54 }
55 
56 
57 template<class BasePhaseModel>
58 template<class Type>
61 (
62  const word& name,
63  const dimensionSet& dims,
64  const bool cache
65 ) const
66 {
67  return zeroField<Type, fvPatchField, volMesh>(name, dims, cache);
68 }
69 
70 
71 template<class BasePhaseModel>
72 template<class Type>
75 (
76  const word& name,
77  const dimensionSet& dims,
78  const bool cache
79 ) const
80 {
81  return zeroField<Type, fvsPatchField, surfaceMesh>(name, dims, cache);
82 }
83 
84 
85 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
86 
87 template<class BasePhaseModel>
89 (
90  const phaseSystem& fluid,
91  const word& phaseName,
92  const label index
93 )
94 :
95  BasePhaseModel(fluid, phaseName, index)
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
100 
101 template<class BasePhaseModel>
103 {}
104 
105 
106 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107 
108 template<class BasePhaseModel>
110 {
111  return true;
112 }
113 
114 
115 template<class BasePhaseModel>
118 {
120  << "Cannot construct a momentum equation for a stationary phase"
121  << exit(FatalError);
122 
123  return tmp<fvVectorMatrix>();
124 }
125 
126 
127 template<class BasePhaseModel>
130 {
132  << "Cannot construct a momentum equation for a stationary phase"
133  << exit(FatalError);
134 
135  return tmp<fvVectorMatrix>();
136 }
137 
138 
139 template<class BasePhaseModel>
142 {
143  return zeroVolField<vector>("U", dimVelocity, true);
144 }
145 
146 
147 template<class BasePhaseModel>
150 {
152  << "Cannot access the velocity of a stationary phase"
153  << exit(FatalError);
154 
155  return const_cast<volVectorField&>(volVectorField::null());
156 }
157 
158 
159 template<class BasePhaseModel>
162 {
163  return zeroSurfaceField<scalar>("phi", dimVolume/dimTime);
164 }
165 
166 
167 template<class BasePhaseModel>
170 {
172  << "Cannot access the flux of a stationary phase"
173  << exit(FatalError);
174 
175  return const_cast<surfaceScalarField&>(surfaceScalarField::null());
176 }
177 
178 
179 template<class BasePhaseModel>
182 {
183  return zeroSurfaceField<scalar>("alphaPhi", dimVolume/dimTime);
184 }
185 
186 
187 template<class BasePhaseModel>
190 {
192  << "Cannot access the volumetric flux of a stationary phase"
193  << exit(FatalError);
194 
195  return const_cast<surfaceScalarField&>(surfaceScalarField::null());
196 }
197 
198 
199 template<class BasePhaseModel>
202 {
203  return zeroSurfaceField<scalar>("alphaRhoPhi", dimMass/dimTime);
204 }
205 
206 
207 template<class BasePhaseModel>
210 {
212  << "Cannot access the mass flux of a stationary phase"
213  << exit(FatalError);
214 
215  return const_cast<surfaceScalarField&>(surfaceScalarField::null());
216 }
217 
218 
219 template<class BasePhaseModel>
222 {
223  return zeroVolField<vector>("DUDt", dimVelocity/dimTime);
224 }
225 
226 
227 template<class BasePhaseModel>
230 {
231  return zeroSurfaceField<scalar>("DUDtf", dimVelocity*dimArea/dimTime);
232 }
233 
234 
235 template<class BasePhaseModel>
238 {
239  return zeroVolField<scalar>("continuityError", dimDensity/dimTime);
240 }
241 
242 
243 template<class BasePhaseModel>
246 {
247  return zeroVolField<scalar>("continuityErrorFlow", dimDensity/dimTime);
248 }
249 
250 
251 template<class BasePhaseModel>
254 {
255  return zeroVolField<scalar>("continuityErrorSources", dimDensity/dimTime);
256 }
257 
258 
259 template<class BasePhaseModel>
262 {
263  return zeroVolField<scalar>("K", sqr(dimVelocity));
264 }
265 
266 
267 template<class BasePhaseModel>
270 {
271  return tmp<volScalarField>();
272 }
273 
274 
275 template<class BasePhaseModel>
277 (
278  tmp<volScalarField> divU
279 )
280 {
282  << "Cannot set the dilatation rate of a stationary phase"
283  << exit(FatalError);
284 }
285 
286 
287 template<class BasePhaseModel>
290 {
291  return zeroVolField<scalar>("continuityError", dimDynamicViscosity);
292 }
293 
294 
295 template<class BasePhaseModel>
298 {
299  return this->thermo().mu();
300 }
301 
302 
303 template<class BasePhaseModel>
306 {
307  return zeroVolField<scalar>("continuityError", dimViscosity);
308 }
309 
310 
311 template<class BasePhaseModel>
314 {
315  return this->thermo().nu();
316 }
317 
318 
319 template<class BasePhaseModel>
322 {
323  return this->thermo().kappa();
324 }
325 
326 
327 template<class BasePhaseModel>
330 {
331  return this->thermo().kappa(patchi);
332 }
333 
334 
335 template<class BasePhaseModel>
338 {
339  return this->thermo().alpha();
340 }
341 
342 
343 template<class BasePhaseModel>
346 {
347  return this->thermo().alpha(patchi);
348 }
349 
350 
351 template<class BasePhaseModel>
354 {
355  return zeroVolField<scalar>("k", sqr(dimVelocity));
356 }
357 
358 
359 template<class BasePhaseModel>
362 {
363  return zeroVolField<scalar>("pPrime", dimPressure);
364 }
365 
366 
367 // ************************************************************************* //
virtual tmp< volScalarField > mut() const
Return the turbulent dynamic viscosity.
virtual tmp< surfaceScalarField > alphaPhi() const
Return the volumetric flux of the phase.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
multiphaseSystem & fluid
Definition: createFields.H:11
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual tmp< fvVectorMatrix > UfEqn()
Return the momentum equation for the face-based algorithm.
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimViscosity
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
virtual tmp< volScalarField > continuityErrorSources() const
Return the continuity error due to any sources.
rhoReactionThermo & thermo
Definition: createFields.H:28
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
const dimensionSet dimDynamicViscosity
virtual tmp< volScalarField > divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual tmp< volScalarField > nuEff() const
Return the effective kinematic viscosity.
virtual tmp< volScalarField > kappaEff() const
Return the effective thermal conductivity.
virtual tmp< volVectorField > U() const
Return the velocity.
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Class which represents a stationary (and therefore probably solid) phase. Generates, but does not store, zero velocity and flux field and turbulent qauantities. Throws an error when non-const access is requested to the motion fields or when the momentum equation is requested. Usage must must protect against such calls.
dynamicFvMesh & mesh
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
word timeName
Definition: getTimeIndex.H:3
const dimensionSet dimPressure
virtual tmp< volScalarField > continuityErrorFlow() const
Return the continuity error due to the flow field.
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimDensity
static const GeometricField< vector, fvPatchField, volMesh > & null()
Return a null geometric field.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual volVectorField & URef()
Access the velocity.
virtual tmp< volScalarField > muEff() const
Return the effective dynamic viscosity.
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
zeroField divU
Definition: alphaSuSp.H:3
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual tmp< volScalarField > alphaEff() const
Return the effective thermal diffusivity.
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
virtual ~StationaryPhaseModel()
Destructor.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
virtual tmp< volScalarField > nut() const
Return the turbulent kinematic viscosity.
virtual bool stationary() const
Return whether the phase is stationary.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
StationaryPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
const dimensionSet dimVelocity