turbulentIntensityKineticEnergyInletFvPatchScalarField.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) 2011-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 
28 #include "fvPatchFieldMapper.H"
29 #include "surfaceFields.H"
30 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39  const dictionary& dict
40 )
41 :
42  inletOutletFvPatchScalarField(p, iF),
43  intensity_(dict.lookup<scalar>("intensity")),
44  UName_(dict.lookupOrDefault<word>("U", "U"))
45 {
46  this->phiName_ = dict.lookupOrDefault<word>("phi", "phi");
47 
48  if (intensity_ < 0 || intensity_ > 1)
49  {
51  << "Turbulence intensity should be specified as a fraction 0-1 "
52  "of the mean velocity\n"
53  " value given is " << intensity_ << nl
54  << " on patch " << this->patch().name()
55  << " of field " << this->internalField().name()
56  << " in file " << this->internalField().objectPath()
57  << exit(FatalError);
58  }
59 
61 
62  this->refValue() = 0.0;
63  this->refGrad() = 0.0;
64  this->valueFraction() = 0.0;
65 }
66 
67 
70 (
72  const fvPatch& p,
74  const fvPatchFieldMapper& mapper
75 )
76 :
77  inletOutletFvPatchScalarField(ptf, p, iF, mapper),
78  intensity_(ptf.intensity_),
79  UName_(ptf.UName_)
80 {}
81 
82 
85 (
88 )
89 :
90  inletOutletFvPatchScalarField(ptf, iF),
91  intensity_(ptf.intensity_),
92  UName_(ptf.UName_)
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
100 {
101  if (updated())
102  {
103  return;
104  }
105 
106  const fvPatchVectorField& Up =
107  patch().lookupPatchField<volVectorField, vector>(UName_);
108 
109  const fvsPatchScalarField& phip =
110  patch().lookupPatchField<surfaceScalarField, scalar>(this->phiName_);
111 
112  this->refValue() = 1.5*sqr(intensity_)*magSqr(Up);
113  this->valueFraction() = neg(phip);
114 
115  inletOutletFvPatchScalarField::updateCoeffs();
116 }
117 
118 
120 (
121  Ostream& os
122 ) const
123 {
125  writeEntry(os, "intensity", intensity_);
126  writeEntryIfDifferent<word>(os, "U", "U", UName_);
127  writeEntryIfDifferent<word>(os, "phi", "phi", this->phiName_);
128  writeEntry(os, "value", *this);
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133 
134 namespace Foam
135 {
137  (
140  );
141 }
142 
143 // ************************************************************************* //
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...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:82
This boundary condition provides a turbulent kinetic energy condition, based on user-supplied turbule...
turbulentIntensityKineticEnergyInletFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
VolField< vector > volVectorField
Definition: volFieldsFwd.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
SurfaceField< scalar > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar neg(const dimensionedScalar &ds)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
static const char nl
Definition: Ostream.H:260
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dictionary dict
volScalarField & p
Foam::surfaceFields.