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-2023 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 #include "fvcLaplacian.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class BasePhaseModel>
33 (
34  const phaseSystem& fluid,
35  const word& phaseName,
36  const bool referencePhase,
37  const label index
38 )
39 :
40  BasePhaseModel(fluid, phaseName, referencePhase, index)
41 {}
42 
43 
44 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
45 
46 template<class BasePhaseModel>
48 {}
49 
50 
51 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
52 
53 template<class BasePhaseModel>
55 {
56  return true;
57 }
58 
59 
60 template<class BasePhaseModel>
63 {
65  << "Cannot construct a momentum equation for a stationary phase"
66  << exit(FatalError);
67 
68  return tmp<fvVectorMatrix>();
69 }
70 
71 
72 template<class BasePhaseModel>
75 {
77  << "Cannot construct a momentum equation for a stationary phase"
78  << exit(FatalError);
79 
80  return tmp<fvVectorMatrix>();
81 }
82 
83 
84 template<class BasePhaseModel>
87 {
88  return volVectorField::New
89  (
90  IOobject::groupName("U", this->name()),
91  this->mesh(),
93  );
94 }
95 
96 
97 template<class BasePhaseModel>
100 {
102  << "Cannot access the velocity of a stationary phase"
103  << exit(FatalError);
104 
105  return const_cast<volVectorField&>(volVectorField::null());
106 }
107 
108 
109 template<class BasePhaseModel>
112 {
114  << "Cannot access the velocity of a stationary phase"
115  << exit(FatalError);
116 
117  return volVectorField::null();
118 }
119 
120 
121 template<class BasePhaseModel>
124 {
126  (
127  IOobject::groupName("phi", this->name()),
128  this->mesh(),
130  );
131 }
132 
133 
134 template<class BasePhaseModel>
137 {
139  << "Cannot access the flux of a stationary phase"
140  << exit(FatalError);
141 
142  return const_cast<surfaceScalarField&>(surfaceScalarField::null());
143 }
144 
145 
146 template<class BasePhaseModel>
149 {
151  << "Cannot access the flux of a stationary phase"
152  << exit(FatalError);
153 
154  return surfaceScalarField::null();
155 }
156 
157 
158 template<class BasePhaseModel>
161 {
163  << "Cannot access the face velocity of a stationary phase"
164  << abort(FatalError);
165 
167  return Uf_;
168 }
169 
170 
171 template<class BasePhaseModel>
174 {
176  << "Cannot access the face velocity of a stationary phase"
177  << exit(FatalError);
178 
179  return const_cast<surfaceVectorField&>(surfaceVectorField::null());
180 }
181 
182 
183 template<class BasePhaseModel>
186 {
188  << "Cannot access the face velocity of a stationary phase"
189  << exit(FatalError);
190 
191  return surfaceVectorField::null();
192 }
193 
194 
195 template<class BasePhaseModel>
198 {
200  (
201  IOobject::groupName("alphaPhi", this->name()),
202  this->mesh(),
204  );
205 }
206 
207 
208 template<class BasePhaseModel>
211 {
213  << "Cannot access the volumetric flux of a stationary phase"
214  << exit(FatalError);
215 
216  return const_cast<surfaceScalarField&>(surfaceScalarField::null());
217 }
218 
219 
220 template<class BasePhaseModel>
223 {
225  << "Cannot access the volumetric flux of a stationary phase"
226  << exit(FatalError);
227 
228  return surfaceScalarField::null();
229 }
230 
231 
232 template<class BasePhaseModel>
235 {
237  (
238  IOobject::groupName("alphaRhoPhi", this->name()),
239  this->mesh(),
241  );
242 }
243 
244 
245 template<class BasePhaseModel>
248 {
250  << "Cannot access the mass flux of a stationary phase"
251  << exit(FatalError);
252 
253  return const_cast<surfaceScalarField&>(surfaceScalarField::null());
254 }
255 
256 
257 template<class BasePhaseModel>
260 {
262  << "Cannot access the mass flux of a stationary phase"
263  << exit(FatalError);
264 
265  return surfaceScalarField::null();
266 }
267 
268 
269 template<class BasePhaseModel>
272 {
273  return volVectorField::New
274  (
275  IOobject::groupName("DUDt", this->name()),
276  this->mesh(),
278  );
279 }
280 
281 
282 template<class BasePhaseModel>
285 {
287  (
288  IOobject::groupName("DUDtf", this->name()),
289  this->mesh(),
291  );
292 }
293 
294 
295 template<class BasePhaseModel>
298 {
299  return volScalarField::New
300  (
301  IOobject::groupName("continuityError", this->name()),
302  this->mesh(),
304  );
305 }
306 
307 
308 template<class BasePhaseModel>
311 {
312  return volScalarField::New
313  (
314  IOobject::groupName("K", this->name()),
315  this->mesh(),
317  );
318 }
319 
320 
321 template<class BasePhaseModel>
324 {
325  static autoPtr<volScalarField> divU_;
326  return divU_;
327 }
328 
329 
330 template<class BasePhaseModel>
332 (
334 )
335 {
337  << "Cannot set the dilatation rate of a stationary phase"
338  << exit(FatalError);
339 }
340 
341 
342 template<class BasePhaseModel>
345 {
346  return this->thermo().kappa().boundaryField()[patchi];
347 }
348 
349 
350 template<class BasePhaseModel>
353 {
354  return volScalarField::New
355  (
356  IOobject::groupName("k", this->name()),
357  this->mesh(),
359  );
360 }
361 
362 
363 template<class BasePhaseModel>
366 {
367  return volScalarField::New
368  (
369  IOobject::groupName("pPrime", this->name()),
370  this->mesh(),
372  );
373 }
374 
375 
376 template<class BasePhaseModel>
379 {
380  const volScalarField& alpha = *this;
381 
382  const surfaceScalarField alphaKappa
383  (
384  alpha.name() + '*' + this->thermo().kappa().name(),
385  fvc::interpolate(alpha)*fvc::interpolate(this->thermo().kappa())
386  );
387 
388  // Return heat flux source as an implicit energy correction
389  // to the temperature gradient flux
390  return
391  -fvc::laplacian(alphaKappa, this->thermo().T())
393  (
394  alphaKappa/fvc::interpolate(this->thermo().Cpv()),
395  he
396  );
397 }
398 
399 
400 // ************************************************************************* //
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
const word & name() const
Return name.
Definition: IOobject.H:310
static word groupName(Name name, const word &group)
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
virtual tmp< surfaceScalarField > DUDtf() const
Return the substantive acceleration on the faces.
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
virtual tmp< volVectorField > U() const
Return the velocity.
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
virtual const autoPtr< volScalarField > & divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
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 const autoPtr< surfaceVectorField > & Uf() const
Return the face velocity.
virtual surfaceVectorField & UfRef()
Access the face velocity.
virtual ~StationaryPhaseModel()
Destructor.
StationaryPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
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 volVectorField & URef()
Access the velocity.
virtual tmp< volVectorField > DUDt() const
Return the substantive acceleration.
virtual bool stationary() const
Return whether the phase is stationary.
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
virtual tmp< scalarField > kappaEff(const label patchi) const
Return the effective thermal conductivity on a patch.
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
virtual const volScalarField & kappa() const =0
Thermal conductivity of mixture [W/m/K].
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Calculate the laplacian of the given field.
label patchi
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< fvMatrix< Type > > laplacianCorrection(const VolField< scalar > &gamma, const VolField< Type > &vf)
Definition: fvmLaplacian.C:340
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static const zero Zero
Definition: zero.H:97
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
const dimensionSet dimPressure
dimensionedSymmTensor sqr(const dimensionedVector &dv)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimTime
const dimensionSet dimDensity
const dimensionSet dimVolume
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const dimensionSet dimMass
const dimensionSet dimVelocity
error FatalError
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
thermo he()
fluidMulticomponentThermo & thermo
Definition: createFields.H:31