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-2020 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  (
139  IOobject::groupName("alphaPhi", this->name()),
140  this->mesh(),
142  );
143 }
144 
145 
146 template<class BasePhaseModel>
149 {
151  << "Cannot access the volumetric flux of a stationary phase"
152  << exit(FatalError);
153 
154  return const_cast<surfaceScalarField&>(surfaceScalarField::null());
155 }
156 
157 
158 template<class BasePhaseModel>
161 {
163  (
164  IOobject::groupName("alphaRhoPhi", this->name()),
165  this->mesh(),
167  );
168 }
169 
170 
171 template<class BasePhaseModel>
174 {
176  << "Cannot access the mass flux of a stationary phase"
177  << exit(FatalError);
178 
179  return const_cast<surfaceScalarField&>(surfaceScalarField::null());
180 }
181 
182 
183 template<class BasePhaseModel>
186 {
187  return volVectorField::New
188  (
189  IOobject::groupName("DUDt", this->name()),
190  this->mesh(),
192  );
193 }
194 
195 
196 template<class BasePhaseModel>
199 {
201  (
202  IOobject::groupName("DUDtf", this->name()),
203  this->mesh(),
205  );
206 }
207 
208 
209 template<class BasePhaseModel>
212 {
213  return volScalarField::New
214  (
215  IOobject::groupName("continuityError", this->name()),
216  this->mesh(),
218  );
219 }
220 
221 
222 template<class BasePhaseModel>
225 {
226  return volScalarField::New
227  (
228  IOobject::groupName("K", this->name()),
229  this->mesh(),
231  );
232 }
233 
234 
235 template<class BasePhaseModel>
238 {
239  return tmp<volScalarField>();
240 }
241 
242 
243 template<class BasePhaseModel>
245 (
246  tmp<volScalarField> divU
247 )
248 {
250  << "Cannot set the dilatation rate of a stationary phase"
251  << exit(FatalError);
252 }
253 
254 
255 template<class BasePhaseModel>
258 {
259  return this->thermo().kappa(patchi);
260 }
261 
262 
263 template<class BasePhaseModel>
266 {
267  return volScalarField::New
268  (
269  IOobject::groupName("k", this->name()),
270  this->mesh(),
272  );
273 }
274 
275 
276 template<class BasePhaseModel>
279 {
280  return volScalarField::New
281  (
282  IOobject::groupName("pPrime", this->name()),
283  this->mesh(),
285  );
286 }
287 
288 
289 template<class BasePhaseModel>
292 {
293  const volScalarField& alpha = *this;
294 
295  return -fvm::laplacian
296  (
297  fvc::interpolate(alpha)*fvc::interpolate(this->thermo().alpha()),
298  he
299  );
300 }
301 
302 
303 // ************************************************************************* //
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.
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)
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< vector >> &)
Return a temporary field constructed from name,.
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
rhoReactionThermo & thermo
Definition: createFields.H:28
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 dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
virtual tmp< scalarField > kappaEff(const label patchi) const
Return the effective thermal conductivity on a patch.
dynamicFvMesh & mesh
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.
phaseSystem & fluid
Definition: createFields.H:11
static const zero Zero
Definition: zero.H:97
const dimensionSet dimPressure
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
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.
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;.
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.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual volVectorField & URef()
Access the velocity.
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
zeroField divU
Definition: alphaSuSp.H:3
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 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.
const dimensionSet dimVelocity