turbulentMixingLengthDissipationRateInletFvPatchScalarField.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) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 #include "turbulenceModel.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
45 (
46  const fvPatch& p,
48 )
49 :
50  inletOutletFvPatchScalarField(p, iF),
51  kName_("k"),
52  mixingLength_(0.0),
53  Cmu_(0.0)
54 {
55  this->refValue() = Zero;
56  this->refGrad() = Zero;
57  this->valueFraction() = 0.0;
58 }
59 
60 
63 (
65  const fvPatch& p,
67  const fvPatchFieldMapper& mapper
68 )
69 :
70  inletOutletFvPatchScalarField(ptf, p, iF, mapper),
71  kName_(ptf.kName_),
72  mixingLength_(ptf.mixingLength_),
73  Cmu_(ptf.Cmu_)
74 {}
75 
76 
79 (
80  const fvPatch& p,
82  const dictionary& dict
83 )
84 :
85  inletOutletFvPatchScalarField(p, iF),
86  kName_(dict.getOrDefault<word>("k", "k")),
87  mixingLength_
88  (
89  dict.getCheck<scalar>("mixingLength", scalarMinMax::ge(SMALL))
90  ),
91  Cmu_(dict.getCheckOrDefault<scalar>("Cmu", 0.09, scalarMinMax::ge(SMALL)))
92 {
93  this->phiName_ = dict.getOrDefault<word>("phi", "phi");
94 
95  this->readValueEntry(dict, IOobjectOption::MUST_READ);
96 
97  this->refValue() = Zero;
98  this->refGrad() = Zero;
99  this->valueFraction() = 0.0;
100 }
101 
102 
105 (
107 )
108 :
109  inletOutletFvPatchScalarField(ptf),
110  kName_(ptf.kName_),
111  mixingLength_(ptf.mixingLength_),
112  Cmu_(ptf.Cmu_)
113 {}
114 
115 
118 (
121 )
122 :
123  inletOutletFvPatchScalarField(ptf, iF),
124  kName_(ptf.kName_),
125  mixingLength_(ptf.mixingLength_),
126  Cmu_(ptf.Cmu_)
127 {}
128 
129 
130 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
131 
133 {
134  if (updated())
135  {
136  return;
137  }
138 
139  // Lookup Cmu corresponding to the turbulence model selected
140  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
141  (
143  (
145  internalField().group()
146  )
147  );
148 
149  Cmu_ = turbModel.coeffDict().getOrDefault<scalar>("Cmu", Cmu_);
150 
151  const scalar Cmu75 = pow(Cmu_, 0.75);
152 
153  const auto& kp =
154  patch().lookupPatchField<volScalarField>(kName_);
155 
156  const auto& phip =
157  patch().lookupPatchField<surfaceScalarField>(this->phiName_);
158 
159  this->refValue() = (Cmu75/mixingLength_)*pow(kp, 1.5);
160  this->valueFraction() = neg(phip);
161 
162  inletOutletFvPatchScalarField::updateCoeffs();
163 }
164 
165 
167 (
168  Ostream& os
169 ) const
170 {
172  os.writeEntry("mixingLength", mixingLength_);
173  os.writeEntry("phi", this->phiName_);
174  os.writeEntry("k", kName_);
176 }
177 
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
182 (
184  turbulentMixingLengthDissipationRateInletFvPatchScalarField
185 );
186 
187 
188 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
189 
190 } // End namespace Foam
191 
192 // ************************************************************************* //
Foam::surfaceFields.
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:403
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
constexpr const char *const group
Group name for atomic constants.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
This boundary condition provides an inlet condition for turbulent kinetic energy dissipation rate...
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
fvPatchField< scalar > fvPatchScalarField
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
A FieldMapper for finite-volume patch fields.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
turbulentMixingLengthDissipationRateInletFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
OBJstream os(runTime.globalPath()/outputName)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const std::string patch
OpenFOAM patch number as a std::string.
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127