turbulentMixingLengthDissipationRateInletFvPatchScalarField.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-2022 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 #include "momentumTransportModel.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
42 (
43  const fvPatch& p,
45 )
46 :
47  inletOutletFvPatchScalarField(p, iF),
48  mixingLength_(0.0),
49  kName_("k")
50 {
51  this->refValue() = 0.0;
52  this->refGrad() = 0.0;
53  this->valueFraction() = 0.0;
54 }
55 
56 
59 (
60  const fvPatch& p,
62  const dictionary& dict
63 )
64 :
65  inletOutletFvPatchScalarField(p, iF),
66  mixingLength_(dict.lookup<scalar>("mixingLength")),
67  kName_(dict.lookupOrDefault<word>("k", "k"))
68 {
69  this->phiName_ = dict.lookupOrDefault<word>("phi", "phi");
70 
71  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
72 
73  this->refValue() = 0.0;
74  this->refGrad() = 0.0;
75  this->valueFraction() = 0.0;
76 }
77 
78 
81 (
83  const fvPatch& p,
85  const fvPatchFieldMapper& mapper
86 )
87 :
88  inletOutletFvPatchScalarField(ptf, p, iF, mapper),
89  mixingLength_(ptf.mixingLength_),
90  kName_(ptf.kName_)
91 {}
92 
93 
96 (
99 )
100 :
101  inletOutletFvPatchScalarField(ptf, iF),
102  mixingLength_(ptf.mixingLength_),
103  kName_(ptf.kName_)
104 {}
105 
106 
107 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108 
110 {
111  if (updated())
112  {
113  return;
114  }
115 
116  // Lookup Cmu corresponding to the turbulence model selected
117  const momentumTransportModel& turbModel =
118  db().lookupObject<momentumTransportModel>
119  (
121  (
122  momentumTransportModel::typeName,
123  internalField().group()
124  )
125  );
126 
127  const scalar Cmu =
128  turbModel.coeffDict().lookupOrDefault<scalar>("Cmu", 0.09);
129 
130  const scalar Cmu75 = pow(Cmu, 0.75);
131 
132  const fvPatchScalarField& kp =
133  patch().lookupPatchField<volScalarField, scalar>(kName_);
134 
135  const fvsPatchScalarField& phip =
136  patch().lookupPatchField<surfaceScalarField, scalar>(this->phiName_);
137 
138  this->refValue() = Cmu75*kp*sqrt(kp)/mixingLength_;
139  this->valueFraction() = neg(phip);
140 
141  inletOutletFvPatchScalarField::updateCoeffs();
142 }
143 
144 
146 (
147  Ostream& os
148 ) const
149 {
151  writeEntry(os, "mixingLength", mixingLength_);
152  writeEntry(os, "phi", this->phiName_);
153  writeEntry(os, "k", kName_);
154  writeEntry(os, "value", *this);
155 }
156 
157 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 
161 (
164 );
165 
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 } // End namespace Foam
170 
171 // ************************************************************************* //
Foam::surfaceFields.
const char *const group
Group name for atomic constants.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
This boundary condition provides a turbulence dissipation, (epsilon) inlet condition based on a spec...
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
virtual label size() const
Return size.
Definition: fvPatch.H:157
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
virtual const dictionary & coeffDict() const =0
Const access to the coefficients dictionary.
turbulentMixingLengthDissipationRateInletFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Abstract base class for turbulence models (RAS, LES and laminar).
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:258
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864