Loading...
Searching...
No Matches
turbulentIntensityKineticEnergyInletFvPatchScalarField.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-2017 OpenFOAM Foundation
9 Copyright (C) 2017-2020 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 "fvPatchFieldMapper.H"
32#include "surfaceFields.H"
33#include "volFields.H"
34
35// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36
39(
40 const fvPatch& p,
42)
43:
44 inletOutletFvPatchScalarField(p, iF),
45 intensity_(0.0),
46 UName_("U")
47{
48 this->refValue() = Zero;
49 this->refGrad() = Zero;
50 this->valueFraction() = 0.0;
51}
52
55(
57 const fvPatch& p,
59 const fvPatchFieldMapper& mapper
61:
62 inletOutletFvPatchScalarField(ptf, p, iF, mapper),
63 intensity_(ptf.intensity_),
64 UName_(ptf.UName_)
65{}
66
69(
70 const fvPatch& p,
72 const dictionary& dict
73)
74:
75 inletOutletFvPatchScalarField(p, iF),
76 intensity_(dict.get<scalar>("intensity")),
77 UName_(dict.getOrDefault<word>("U", "U"))
78{
80 this->phiName_ = dict.getOrDefault<word>("phi", "phi");
81
82 if (intensity_ < 0 || intensity_ > 1)
83 {
85 << "Turbulence intensity should be specified as a fraction 0-1 "
86 "of the mean velocity\n"
87 " value given is " << intensity_ << nl
88 << " on patch " << this->patch().name()
89 << " of field " << this->internalField().name()
90 << " in file " << this->internalField().objectPath()
91 << exit(FatalError);
92 }
93
94 this->readValueEntry(dict, IOobjectOption::MUST_READ);
95
96 this->refValue() = Zero;
97 this->refGrad() = Zero;
98 this->valueFraction() = 0.0;
99}
100
103(
105)
107 inletOutletFvPatchScalarField(ptf),
108 intensity_(ptf.intensity_),
109 UName_(ptf.UName_)
110{}
111
112
115(
118)
119:
120 inletOutletFvPatchScalarField(ptf, iF),
121 intensity_(ptf.intensity_),
122 UName_(ptf.UName_)
123{}
124
125
126// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127
130{
131 if (updated())
132 {
133 return;
134 }
135
136 const auto& Up =
137 patch().lookupPatchField<volVectorField>(UName_);
138
139 const auto& phip =
140 patch().lookupPatchField<surfaceScalarField>(this->phiName_);
141
142 this->refValue() = 1.5*sqr(intensity_)*magSqr(Up);
143 this->valueFraction() = neg(phip);
144
145 inletOutletFvPatchScalarField::updateCoeffs();
146}
147
148
150(
151 Ostream& os
152) const
153{
155 os.writeEntry("intensity", intensity_);
156 os.writeEntryIfDifferent<word>("U", "U", UName_);
157 os.writeEntryIfDifferent<word>("phi", "phi", this->phiName_);
160
161
162// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163
164namespace Foam
165{
167 (
169 turbulentIntensityKineticEnergyInletFvPatchScalarField
170 );
171}
172
173// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
@ MUST_READ
Reading required.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
virtual void readDict(const dictionary &dict)
Read dictionary entries.
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
This boundary condition provides a turbulent kinetic energy condition, based on user-supplied turbule...
turbulentIntensityKineticEnergyInletFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensionedScalar neg(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
Foam::surfaceFields.