turbulentMixingLengthFrequencyInletFvPatchScalarField.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  mixingLength_(0.0),
52  kName_("undefined-k")
53 {
54  this->refValue() = Zero;
55  this->refGrad() = Zero;
56  this->valueFraction() = 0.0;
57 }
58 
61 (
63  const fvPatch& p,
65  const fvPatchFieldMapper& mapper
66 )
67 :
68  inletOutletFvPatchScalarField(ptf, p, iF, mapper),
69  mixingLength_(ptf.mixingLength_),
70  kName_(ptf.kName_)
71 {}
72 
75 (
76  const fvPatch& p,
78  const dictionary& dict
79 )
80 :
81  inletOutletFvPatchScalarField(p, iF),
82  mixingLength_(dict.get<scalar>("mixingLength")),
83  kName_(dict.getOrDefault<word>("k", "k"))
84 {
85  this->phiName_ = dict.getOrDefault<word>("phi", "phi");
86 
87  this->readValueEntry(dict, IOobjectOption::MUST_READ);
88 
89  this->refValue() = Zero;
90  this->refGrad() = Zero;
91  this->valueFraction() = 0.0;
92 }
93 
96 (
98 )
99 :
100  inletOutletFvPatchScalarField(ptf),
101  mixingLength_(ptf.mixingLength_),
102  kName_(ptf.kName_)
103 {}
104 
107 (
110 )
111 :
112  inletOutletFvPatchScalarField(ptf, iF),
113  mixingLength_(ptf.mixingLength_),
114  kName_(ptf.kName_)
115 {}
116 
117 
118 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119 
121 {
122  if (updated())
123  {
124  return;
125  }
126 
127  // Lookup Cmu corresponding to the turbulence model selected
128  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
129  (
131  (
133  internalField().group()
134  )
135  );
136 
137  const scalar Cmu =
138  turbModel.coeffDict().getOrDefault<scalar>("Cmu", 0.09);
139 
140  const scalar Cmu25 = pow(Cmu, 0.25);
141 
142  const auto& kp =
143  patch().lookupPatchField<volScalarField>(kName_);
144 
145  const auto& phip =
146  patch().lookupPatchField<surfaceScalarField>(this->phiName_);
147 
148  this->refValue() = sqrt(kp)/(Cmu25*mixingLength_);
149  this->valueFraction() = neg(phip);
150 
151  inletOutletFvPatchScalarField::updateCoeffs();
152 }
153 
154 
156 (
157  Ostream& os
158 ) const
159 {
161  os.writeEntry("mixingLength", mixingLength_);
162  os.writeEntry("phi", this->phiName_);
163  os.writeEntry("k", kName_);
165 }
166 
167 
168 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 
171 (
173  turbulentMixingLengthFrequencyInletFvPatchScalarField
174 );
175 
176 
177 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 
179 } // End namespace Foam
180 
181 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Foam::surfaceFields.
dictionary dict
This boundary condition provides a turbulence specific dissipation, (omega) inlet condition based on...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
dimensionedScalar sqrt(const dimensionedScalar &ds)
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:375
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:81
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.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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
OBJstream os(runTime.globalPath()/outputName)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
turbulentMixingLengthFrequencyInletFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127