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-2018 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 "turbulenceModel.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 (
61  const fvPatch& p,
63  const fvPatchFieldMapper& mapper
64 )
65 :
66  inletOutletFvPatchScalarField(ptf, p, iF, mapper),
67  mixingLength_(ptf.mixingLength_),
68  kName_(ptf.kName_)
69 {}
70 
71 
74 (
75  const fvPatch& p,
77  const dictionary& dict
78 )
79 :
80  inletOutletFvPatchScalarField(p, iF),
81  mixingLength_(readScalar(dict.lookup("mixingLength"))),
82  kName_(dict.lookupOrDefault<word>("k", "k"))
83 {
84  this->phiName_ = dict.lookupOrDefault<word>("phi", "phi");
85 
86  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
87 
88  this->refValue() = 0.0;
89  this->refGrad() = 0.0;
90  this->valueFraction() = 0.0;
91 }
92 
93 
96 (
98 )
99 :
100  inletOutletFvPatchScalarField(ptf),
101  mixingLength_(ptf.mixingLength_),
102  kName_(ptf.kName_)
103 {}
104 
105 
108 (
111 )
112 :
113  inletOutletFvPatchScalarField(ptf, iF),
114  mixingLength_(ptf.mixingLength_),
115  kName_(ptf.kName_)
116 {}
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
122 {
123  if (updated())
124  {
125  return;
126  }
127 
128  // Lookup Cmu corresponding to the turbulence model selected
129  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
130  (
132  (
134  internalField().group()
135  )
136  );
137 
138  const scalar Cmu =
139  turbModel.coeffDict().lookupOrDefault<scalar>("Cmu", 0.09);
140 
141  const scalar Cmu75 = pow(Cmu, 0.75);
142 
143  const fvPatchScalarField& kp =
144  patch().lookupPatchField<volScalarField, scalar>(kName_);
145 
146  const fvsPatchScalarField& phip =
147  patch().lookupPatchField<surfaceScalarField, scalar>(this->phiName_);
148 
149  this->refValue() = Cmu75*kp*sqrt(kp)/mixingLength_;
150  this->valueFraction() = 1.0 - pos0(phip);
151 
152  inletOutletFvPatchScalarField::updateCoeffs();
153 }
154 
155 
157 (
158  Ostream& os
159 ) const
160 {
162  os.writeKeyword("mixingLength")
163  << mixingLength_ << token::END_STATEMENT << nl;
164  os.writeKeyword("phi") << this->phiName_ << token::END_STATEMENT << nl;
165  os.writeKeyword("k") << kName_ << token::END_STATEMENT << nl;
166  writeEntry("value", os);
167 }
168 
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
173 (
176 );
177 
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 } // End namespace Foam
182 
183 // ************************************************************************* //
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:137
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
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:362
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
This boundary condition provides a turbulence dissipation, (epsilon) inlet condition based on a spec...
static const word propertiesName
Default name of the turbulence properties dictionary.
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.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
virtual label size() const
Return size.
Definition: fvPatch.H:161
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
dimensionedScalar pos0(const dimensionedScalar &ds)
static const char nl
Definition: Ostream.H:265
turbulentMixingLengthDissipationRateInletFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
virtual const dictionary & coeffDict() const =0
Const access to the coefficients dictionary.
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:395
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:576