Loading...
Searching...
No Matches
epsilonWallFunctionFvPatchScalarField.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) 2011-2019 OpenFOAM Foundation
9 Copyright (C) 2017-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
31#include "turbulenceModel.H"
32#include "fvMatrix.H"
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
38
39// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
40
42{
43 if (master_ != -1)
44 {
45 return;
46 }
47
48 const auto& epsilon =
49 static_cast<const volScalarField&>(this->internalField());
50
51 const volScalarField::Boundary& bf = epsilon.boundaryField();
52
53 label master = -1;
54 forAll(bf, patchi)
55 {
57 {
59
60 if (master == -1)
61 {
62 master = patchi;
63 }
65 epf.master() = master;
66 }
67 }
68}
69
70
72{
73 const auto& epsilon =
74 static_cast<const volScalarField&>(this->internalField());
75
76 const volScalarField::Boundary& bf = epsilon.boundaryField();
77
78 const fvMesh& mesh = epsilon.mesh();
79
80 if (initialised_ && !mesh.changing())
81 {
82 return;
83 }
84
85 volScalarField weights
86 (
88 (
89 "weights",
90 mesh.time().timeName(),
91 mesh,
95 ),
96 mesh,
98 );
99
100 DynamicList<label> epsilonPatches(bf.size());
101 forAll(bf, patchi)
102 {
104 {
105 epsilonPatches.append(patchi);
106
107 const labelUList& faceCells = bf[patchi].patch().faceCells();
108 for (const auto& faceCell : faceCells)
109 {
110 ++weights[faceCell];
111 }
112 }
113 }
114
115 cornerWeights_.setSize(bf.size());
116
117 for (const auto& patchi : epsilonPatches)
118 {
119 const fvPatchScalarField& wf = weights.boundaryField()[patchi];
120 cornerWeights_[patchi] = 1.0/wf.patchInternalField();
121 }
122
123 G_.setSize(internalField().size(), Zero);
125
126 initialised_ = true;
127}
128
129
132(
133 const label patchi
134)
135{
136 const auto& epsilon =
137 static_cast<const volScalarField&>(this->internalField());
138
139 const volScalarField::Boundary& bf = epsilon.boundaryField();
140
141 const auto& epf =
143
144 return const_cast<epsilonWallFunctionFvPatchScalarField&>(epf);
145}
146
147
149(
151 scalarField& G0,
152 scalarField& epsilon0
153)
154{
155 // Accumulate all of the G and epsilon contributions
156 forAll(cornerWeights_, patchi)
157 {
158 if (!cornerWeights_[patchi].empty())
159 {
160 epsilonWallFunctionFvPatchScalarField& epf = epsilonPatch(patchi);
161
162 const List<scalar>& w = cornerWeights_[patchi];
163
164 epf.calculate(turbulence, w, epf.patch(), G0, epsilon0);
165 }
166 }
167
168 // Apply zero-gradient condition for epsilon
169 forAll(cornerWeights_, patchi)
170 {
171 if (!cornerWeights_[patchi].empty())
172 {
173 epsilonWallFunctionFvPatchScalarField& epf = epsilonPatch(patchi);
175 epf == scalarField(epsilon0, epf.patch().faceCells());
176 }
177 }
178}
179
180
182(
183 const turbulenceModel& turbModel,
184 const List<scalar>& cornerWeights,
185 const fvPatch& patch,
186 scalarField& G0,
187 scalarField& epsilon0
188)
189{
190 const label patchi = patch.index();
191
192 const scalar Cmu25 = pow025(wallCoeffs_.Cmu());
193 const scalar Cmu75 = pow(wallCoeffs_.Cmu(), 0.75);
194 const scalar kappa = wallCoeffs_.kappa();
195 const scalar yPlusLam = wallCoeffs_.yPlusLam();
196
197 const scalarField& y = turbModel.y()[patchi];
198
199 const labelUList& faceCells = patch.faceCells();
200
201 const tmp<scalarField> tnuw = turbModel.nu(patchi);
202 const scalarField& nuw = tnuw();
203
204 const tmp<volScalarField> tk = turbModel.k();
205 const volScalarField& k = tk();
206
207 // Calculate y-plus
208 const auto yPlus = [&](const label facei) -> scalar
209 {
210 return
211 (
212 Cmu25*y[facei]*sqrt(k[faceCells[facei]])/nuw[facei]
213 );
214 };
215
216 // Contribution from the viscous sublayer
217 const auto epsilonVis = [&](const label facei) -> scalar
218 {
219 return
220 (
221 2.0*k[faceCells[facei]]*nuw[facei]
222 / sqr(y[facei])
223 );
224 };
225
226 // Contribution from the inertial sublayer
227 const auto epsilonLog = [&](const label facei) -> scalar
228 {
229 return
230 (
231 Cmu75*pow(k[faceCells[facei]], 1.5)
232 / (kappa*y[facei])
233 );
234 };
235
236 switch (blender_)
237 {
238 case blenderType::STEPWISE:
239 {
240 forAll(faceCells, facei)
241 {
242 if (lowReCorrection_ && yPlus(facei) < yPlusLam)
243 {
244 epsilon0[faceCells[facei]] +=
245 cornerWeights[facei]
246 * epsilonVis(facei);
247 }
248 else
249 {
250 epsilon0[faceCells[facei]] +=
251 cornerWeights[facei]
252 * epsilonLog(facei);
253 }
254 }
255 break;
256 }
257
258 case blenderType::BINOMIAL:
259 {
260 forAll(faceCells, facei)
261 {
262 // (ME:Eqs. 15-16)
263 epsilon0[faceCells[facei]] +=
264 cornerWeights[facei]
265 * pow
266 (
267 pow(epsilonVis(facei), n_) + pow(epsilonLog(facei), n_),
268 scalar(1)/n_
269 );
270 }
271 break;
272 }
273
274 case blenderType::MAX:
275 {
276 forAll(faceCells, facei)
277 {
278 // (PH:Eq. 27)
279 epsilon0[faceCells[facei]] +=
280 cornerWeights[facei]
281 * max(epsilonVis(facei), epsilonLog(facei));
282 }
283 break;
284 }
285
286 case blenderType::EXPONENTIAL:
287 {
288 forAll(faceCells, facei)
289 {
290 // (PH:p. 193)
291 const scalar yPlusFace = yPlus(facei);
292 const scalar Gamma =
293 0.001*pow4(yPlusFace)/(scalar(1) + yPlusFace);
294 const scalar invGamma = scalar(1)/(Gamma + ROOTVSMALL);
295
296 epsilon0[faceCells[facei]] +=
297 cornerWeights[facei]
298 * (
299 epsilonVis(facei)*exp(-Gamma)
300 + epsilonLog(facei)*exp(-invGamma)
301 );
302 }
303 break;
304 }
305
306 case blenderType::TANH:
307 {
308 forAll(faceCells, facei)
309 {
310 // (KAS:Eqs. 33-34)
311 const scalar epsilonVisFace = epsilonVis(facei);
312 const scalar epsilonLogFace = epsilonLog(facei);
313 const scalar b1 = epsilonVisFace + epsilonLogFace;
314 const scalar b2 =
315 pow
316 (
317 pow(epsilonVisFace, 1.2) + pow(epsilonLogFace, 1.2),
318 1.0/1.2
319 );
320 const scalar phiTanh = tanh(pow4(0.1*yPlus(facei)));
321
322 epsilon0[faceCells[facei]] +=
323 cornerWeights[facei]
324 * (phiTanh*b1 + (1 - phiTanh)*b2);
325 }
326 break;
327 }
328 }
329
330 const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
331 const scalarField magGradUw(mag(Uw.snGrad()));
332
333 const tmp<scalarField> tnutw = turbModel.nut(patchi);
334 const scalarField& nutw = tnutw();
335
336 forAll(faceCells, facei)
337 {
338 if (!lowReCorrection_ || (yPlus(facei) > yPlusLam))
339 {
340 G0[faceCells[facei]] +=
341 cornerWeights[facei]
342 *(nutw[facei] + nuw[facei])
343 *magGradUw[facei]
344 *Cmu25*sqrt(k[faceCells[facei]])
345 /(kappa*y[facei]);
346 }
347 }
348}
349
350
352(
353 Ostream& os
354) const
355{
357 os.writeEntryIfDifferent<bool>("lowReCorrection", false, lowReCorrection_);
358 wallCoeffs_.writeEntries(os);
359}
360
361
362// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
363
366(
367 const fvPatch& p,
369)
370:
371 fixedValueFvPatchField<scalar>(p, iF),
373 lowReCorrection_(false),
374 initialised_(false),
375 master_(-1),
377 G_(),
378 epsilon_(),
380{}
381
382
385(
387 const fvPatch& p,
389 const fvPatchFieldMapper& mapper
390)
391:
392 fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
394 lowReCorrection_(ptf.lowReCorrection_),
395 initialised_(false),
396 master_(-1),
398 G_(),
399 epsilon_(),
401{}
402
403
406(
407 const fvPatch& p,
409 const dictionary& dict
410)
411:
412 fixedValueFvPatchField<scalar>(p, iF, dict),
413 wallFunctionBlenders(dict, blenderType::STEPWISE, scalar(2)),
414 lowReCorrection_(dict.getOrDefault("lowReCorrection", false)),
415 initialised_(false),
416 master_(-1),
417 wallCoeffs_(dict),
418 G_(),
419 epsilon_(),
421{
422 // Apply zero-gradient condition on start-up
423 this->extrapolateInternal();
424}
425
426
429(
431)
432:
433 fixedValueFvPatchField<scalar>(ewfpsf),
434 wallFunctionBlenders(ewfpsf),
435 lowReCorrection_(ewfpsf.lowReCorrection_),
436 initialised_(false),
437 master_(-1),
439 G_(),
440 epsilon_(),
442{}
443
444
447(
450)
451:
452 fixedValueFvPatchField<scalar>(ewfpsf, iF),
453 wallFunctionBlenders(ewfpsf),
454 lowReCorrection_(ewfpsf.lowReCorrection_),
455 initialised_(false),
456 master_(-1),
457 wallCoeffs_(ewfpsf.wallCoeffs_),
458 G_(),
461{}
462
463
464// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
465
467(
468 bool init
469)
470{
471 if (patch().index() == master_)
472 {
473 if (init)
474 {
475 G_ = 0.0;
476 }
477
478 return G_;
479 }
480
481 return epsilonPatch(master_).G();
482}
483
484
486(
487 bool init
488)
489{
490 if (patch().index() == master_)
491 {
492 if (init)
493 {
494 epsilon_ = 0.0;
495 }
496
497 return epsilon_;
498 }
499
500 return epsilonPatch(master_).epsilon(init);
501}
502
503
505{
506 if (updated())
507 {
508 return;
509 }
510
511 const auto& turbModel = db().lookupObject<turbulenceModel>
512 (
514 (
516 internalField().group()
517 )
518 );
519
520 setMaster();
521
522 if (patch().index() == master_)
523 {
524 createAveragingWeights();
525 calculateTurbulenceFields(turbModel, G(true), epsilon(true));
526 }
527
528 const scalarField& G0 = this->G();
529 const scalarField& epsilon0 = this->epsilon();
530
531 typedef DimensionedField<scalar, volMesh> FieldType;
532
533 FieldType& G = db().lookupObjectRef<FieldType>(turbModel.GName());
534
535 FieldType& epsilon = const_cast<FieldType&>(internalField());
536
537 forAll(*this, facei)
538 {
539 const label celli = patch().faceCells()[facei];
540
541 G[celli] = G0[celli];
542 epsilon[celli] = epsilon0[celli];
543 }
544
546}
547
548
550(
551 const scalarField& weights
552)
553{
554 if (updated())
555 {
556 return;
557 }
558
559 const auto& turbModel = db().lookupObject<turbulenceModel>
560 (
562 (
564 internalField().group()
565 )
566 );
567
568 setMaster();
569
570 if (patch().index() == master_)
571 {
572 createAveragingWeights();
573 calculateTurbulenceFields(turbModel, G(true), epsilon(true));
574 }
575
576 const scalarField& G0 = this->G();
577 const scalarField& epsilon0 = this->epsilon();
578
579 typedef DimensionedField<scalar, volMesh> FieldType;
580
581 FieldType& G = db().lookupObjectRef<FieldType>(turbModel.GName());
582
583 FieldType& epsilon = const_cast<FieldType&>(internalField());
584
585 scalarField& epsilonf = *this;
586
587 // Only set the values if the weights are > tolerance
588 forAll(weights, facei)
589 {
590 const scalar w = weights[facei];
591
592 if (w > tolerance_)
593 {
594 const label celli = patch().faceCells()[facei];
595
596 G[celli] = (1.0 - w)*G[celli] + w*G0[celli];
597 epsilon[celli] = (1.0 - w)*epsilon[celli] + w*epsilon0[celli];
598 epsilonf[facei] = epsilon[celli];
600 }
601
603}
604
605
607(
608 fvMatrix<scalar>& matrix
609)
610{
611 if (manipulatedMatrix())
612 {
613 return;
614 }
615
616 matrix.setValues(patch().faceCells(), patchInternalField());
617
619}
620
621
623(
624 fvMatrix<scalar>& matrix,
625 const Field<scalar>& weights
626)
627{
628 if (manipulatedMatrix())
629 {
630 return;
631 }
632
633 DynamicList<label> constraintCells(weights.size());
634 DynamicList<scalar> constraintValues(weights.size());
635 const labelUList& faceCells = patch().faceCells();
636
637 const DimensionedField<scalar, volMesh>& fld = internalField();
638
639 forAll(weights, facei)
640 {
641 // Only set the values if the weights are > tolerance
642 if (weights[facei] > tolerance_)
643 {
644 const label celli = faceCells[facei];
645
646 constraintCells.append(celli);
647 constraintValues.append(fld[celli]);
648 }
649 }
650
651 if (debug)
652 {
653 Pout<< "Patch: " << patch().name()
654 << ": number of constrained cells = " << constraintCells.size()
655 << " out of " << patch().size()
656 << endl;
657 }
659 matrix.setValues(constraintCells, constraintValues);
660
662}
663
664
666(
667 Ostream& os
668) const
669{
671 writeLocalEntries(os);
674
675
676// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
677
678namespace Foam
679{
681 (
683 epsilonWallFunctionFvPatchScalarField
684 );
685}
686
687
688// ************************************************************************* //
scalar y
label k
Macros for easy insertion into run-time selection tables.
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
bool empty() const noexcept
Definition UList.H:701
void size(const label n)
Definition UList.H:118
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
This boundary condition provides wall functions for the turbulent kinetic energy dissipation rate (i....
List< List< scalar > > cornerWeights_
List of averaging corner weights.
epsilonWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
scalarField epsilon_
Local copy of turbulence epsilon field.
void writeLocalEntries(Ostream &) const
Write local wall function variables.
virtual void manipulateMatrix(fvMatrix< scalar > &matrix)
Manipulate matrix.
static scalar tolerance_
Tolerance used in weighted calculations.
const bool lowReCorrection_
Apply low-Re correction term (default = no).
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void createAveragingWeights()
Create the averaging weights for cells which are bounded by multiple wall function faces.
wallFunctionCoefficients wallCoeffs_
Wall-function coefficients.
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &epsilon)
Calculate the epsilon and G.
virtual label & master()
Return non-const access to the master patch ID.
scalarField & G(bool init=false)
Return non-const access to the master's G field.
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
virtual epsilonWallFunctionFvPatchScalarField & epsilonPatch(const label patchi)
Helper function to return non-const access to an epsilon patch.
virtual void calculateTurbulenceFields(const turbulenceModel &turbulence, scalarField &G0, scalarField &epsilon0)
Main driver to calculate the turbulence fields.
scalarField & epsilon(bool init=false)
Return non-const access to the master's epsilon field.
virtual void setMaster()
Set the master patch - master is responsible for updating all wall function patches.
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
fixedValueFvPatchField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
void setValues(const labelUList &cellLabels, const Type &value)
Set solution in given cells to the specified value and eliminate the corresponding equations from the...
Definition fvMatrix.C:971
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
const objectRegistry & db() const
The associated objectRegistry.
const fvPatch & patch() const noexcept
Return the patch.
bool manipulatedMatrix() const noexcept
True if the matrix has already been manipulated.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
const DimensionedField< scalar, volMesh > & internalField() const noexcept
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
virtual const labelUList & faceCells() const
Return faceCells.
Definition fvPatch.C:107
bool changing() const noexcept
Is mesh changing (topology changing and/or moving).
Definition polyMesh.H:768
A class for managing temporary objects.
Definition tmp.H:75
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
const nearWallDist & y() const
Return the near wall distances.
The class wallFunctionBlenders is a base class that hosts common entries for various derived wall-fun...
void writeEntries(Ostream &) const
Write wall-function blending data as dictionary entries.
wallFunctionBlenders()
Default construct with default coefficients.
blenderType
Options for the blending treatment of viscous and inertial sublayers.
@ STEPWISE
"Stepwise switch (discontinuous)"
enum blenderType blender_
Blending treatment.
scalar n_
Binomial blending exponent being used when blenderType is blenderType::BINOMIAL.
volScalarField & p
dynamicFvMesh & mesh
scalar yPlus
scalar epsilon
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
constexpr const char *const group
Group name for atomic constants.
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
const dimensionedScalar epsilon0
Electric constant: default SI units: [F/m].
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar G
Newtonian constant of gravitation.
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar tanh(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar pow4(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
fvPatchField< vector > fvPatchVectorField
fvPatchField< scalar > fvPatchScalarField
dimensionedScalar pow025(const dimensionedScalar &ds)
static bool initialised_(false)
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299