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-2021 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class BasePhaseModel>
32 (
33  const phaseSystem& fluid,
34  const word& phaseName,
35  const bool referencePhase,
36  const label index
37 )
38 :
39  BasePhaseModel(fluid, phaseName, referencePhase, index)
40 {}
41 
42 
43 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
44 
45 template<class BasePhaseModel>
47 {}
48 
49 
50 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
51 
52 template<class BasePhaseModel>
54 {
55  return true;
56 }
57 
58 
59 template<class BasePhaseModel>
62 {
64  << "Cannot construct a momentum equation for a stationary phase"
65  << exit(FatalError);
66 
67  return tmp<fvVectorMatrix>();
68 }
69 
70 
71 template<class BasePhaseModel>
74 {
76  << "Cannot construct a momentum equation for a stationary phase"
77  << exit(FatalError);
78 
79  return tmp<fvVectorMatrix>();
80 }
81 
82 
83 template<class BasePhaseModel>
86 {
87  return volVectorField::New
88  (
89  IOobject::groupName("U", this->name()),
90  this->mesh(),
92  );
93 }
94 
95 
96 template<class BasePhaseModel>
99 {
101  << "Cannot access the velocity of a stationary phase"
102  << exit(FatalError);
103 
104  return const_cast<volVectorField&>(volVectorField::null());
105 }
106 
107 
108 template<class BasePhaseModel>
111 {
113  (
114  IOobject::groupName("phi", this->name()),
115  this->mesh(),
117  );
118 }
119 
120 
121 template<class BasePhaseModel>
124 {
126  << "Cannot access the flux of a stationary phase"
127  << exit(FatalError);
128 
129  return const_cast<surfaceScalarField&>(surfaceScalarField::null());
130 }
131 
132 
133 template<class BasePhaseModel>
136 {
138  << "Cannot access the face velocity of a stationary phase"
139  << abort(FatalError);
140 
141  return tmp<Foam::surfaceVectorField>();
142 }
143 
144 
145 template<class BasePhaseModel>
148 {
150  << "Cannot access the face velocity of a stationary phase"
151  << exit(FatalError);
152 
153  return const_cast<surfaceVectorField&>(surfaceVectorField::null());
154 }
155 
156 
157 template<class BasePhaseModel>
160 {
162  (
163  IOobject::groupName("alphaPhi", this->name()),
164  this->mesh(),
166  );
167 }
168 
169 
170 template<class BasePhaseModel>
173 {
175  << "Cannot access the volumetric flux of a stationary phase"
176  << exit(FatalError);
177 
178  return const_cast<surfaceScalarField&>(surfaceScalarField::null());
179 }
180 
181 
182 template<class BasePhaseModel>
185 {
187  (
188  IOobject::groupName("alphaRhoPhi", this->name()),
189  this->mesh(),
191  );
192 }
193 
194 
195 template<class BasePhaseModel>
198 {
200  << "Cannot access the mass flux of a stationary phase"
201  << exit(FatalError);
202 
203  return const_cast<surfaceScalarField&>(surfaceScalarField::null());
204 }
205 
206 
207 template<class BasePhaseModel>
210 {
211  return volVectorField::New
212  (
213  IOobject::groupName("DUDt", this->name()),
214  this->mesh(),
216  );
217 }
218 
219 
220 template<class BasePhaseModel>
223 {
225  (
226  IOobject::groupName("DUDtf", this->name()),
227  this->mesh(),
229  );
230 }
231 
232 
233 template<class BasePhaseModel>
236 {
237  return volScalarField::New
238  (
239  IOobject::groupName("continuityError", this->name()),
240  this->mesh(),
242  );
243 }
244 
245 
246 template<class BasePhaseModel>
249 {
250  return volScalarField::New
251  (
252  IOobject::groupName("K", this->name()),
253  this->mesh(),
255  );
256 }
257 
258 
259 template<class BasePhaseModel>
262 {
263  return tmp<volScalarField>();
264 }
265 
266 
267 template<class BasePhaseModel>
269 (
270  tmp<volScalarField> divU
271 )
272 {
274  << "Cannot set the dilatation rate of a stationary phase"
275  << exit(FatalError);
276 }
277 
278 
279 template<class BasePhaseModel>
282 {
283  return this->thermo().kappa(patchi);
284 }
285 
286 
287 template<class BasePhaseModel>
290 {
291  return volScalarField::New
292  (
293  IOobject::groupName("k", this->name()),
294  this->mesh(),
296  );
297 }
298 
299 
300 template<class BasePhaseModel>
303 {
304  return volScalarField::New
305  (
306  IOobject::groupName("pPrime", this->name()),
307  this->mesh(),
309  );
310 }
311 
312 
313 template<class BasePhaseModel>
316 {
317  const volScalarField& alpha = *this;
318 
319  return -fvm::laplacian
320  (
321  fvc::interpolate(alpha)*fvc::interpolate(this->thermo().alpha()),
322  he
323  );
324 }
325 
326 
327 // ************************************************************************* //
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
fluidReactionThermo & thermo
Definition: createFields.H:28
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
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:323
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 dimPressure
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< vector >> &)
Return a temporary field constructed from name,.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
virtual tmp< volScalarField > divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual tmp< volVectorField > U() const
Return the velocity.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const dimensionSet dimTime
virtual tmp< scalarField > kappaEff(const label patchi) const
Return the effective thermal conductivity on a patch.
dynamicFvMesh & mesh
virtual tmp< surfaceVectorField > Uf() const
Return the face velocity.
static word groupName(Name name, const word &group)
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
const dimensionSet dimDensity
phaseSystem & fluid
Definition: createFields.H:11
static const zero Zero
Definition: zero.H:97
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimVelocity
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
tmp< volScalarField > divU
Definition: alphaSuSp.H:9
const dimensionSet dimMass
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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.
static const GeometricField< vector, fvPatchField, volMesh > & null()
Return a null geometric field.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
StationaryPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
virtual volVectorField & URef()
Access the velocity.
const dimensionSet dimVolume
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
virtual surfaceVectorField & UfRef()
Access the face velocity.
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
virtual ~StationaryPhaseModel()
Destructor.
virtual bool stationary() const
Return whether the phase is stationary.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.