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 | Copyright (C) 2011-2016 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_("undefined-k")
50 {
51  this->refValue() = 0.0;
52  this->refGrad() = 0.0;
53  this->valueFraction() = 0.0;
54 }
55 
58 (
60  const fvPatch& p,
62  const fvPatchFieldMapper& mapper
63 )
64 :
65  inletOutletFvPatchScalarField(ptf, p, iF, mapper),
66  mixingLength_(ptf.mixingLength_),
67  kName_(ptf.kName_)
68 {}
69 
72 (
73  const fvPatch& p,
75  const dictionary& dict
76 )
77 :
78  inletOutletFvPatchScalarField(p, iF),
79  mixingLength_(readScalar(dict.lookup("mixingLength"))),
80  kName_(dict.lookupOrDefault<word>("k", "k"))
81 {
82  this->phiName_ = dict.lookupOrDefault<word>("phi", "phi");
83 
84  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
85 
86  this->refValue() = 0.0;
87  this->refGrad() = 0.0;
88  this->valueFraction() = 0.0;
89 }
90 
93 (
95 )
96 :
97  inletOutletFvPatchScalarField(ptf),
98  mixingLength_(ptf.mixingLength_),
99  kName_(ptf.kName_)
100 {}
101 
104 (
107 )
108 :
109  inletOutletFvPatchScalarField(ptf, iF),
110  mixingLength_(ptf.mixingLength_),
111  kName_(ptf.kName_)
112 {}
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
118 {
119  if (updated())
120  {
121  return;
122  }
123 
124  // Lookup Cmu corresponding to the turbulence model selected
125  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
126  (
128  (
130  internalField().group()
131  )
132  );
133 
134  const scalar Cmu =
135  turbModel.coeffDict().lookupOrDefault<scalar>("Cmu", 0.09);
136 
137  const scalar Cmu25 = pow(Cmu, 0.25);
138 
139  const fvPatchScalarField& kp =
140  patch().lookupPatchField<volScalarField, scalar>(kName_);
141 
142  const fvsPatchScalarField& phip =
143  patch().lookupPatchField<surfaceScalarField, scalar>(this->phiName_);
144 
145  this->refValue() = sqrt(kp)/(Cmu25*mixingLength_);
146  this->valueFraction() = 1.0 - pos(phip);
147 
148  inletOutletFvPatchScalarField::updateCoeffs();
149 }
150 
151 
153 (
154  Ostream& os
155 ) const
156 {
158  os.writeKeyword("mixingLength")
159  << mixingLength_ << token::END_STATEMENT << nl;
160  os.writeKeyword("phi") << this->phiName_ << token::END_STATEMENT << nl;
161  os.writeKeyword("k") << kName_ << token::END_STATEMENT << nl;
162  writeEntry("value", os);
163 }
164 
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
169 (
172 );
173 
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 } // End namespace Foam
178 
179 // ************************************************************************* //
Foam::surfaceFields.
const char *const group
Group name for atomic constants.
This boundary condition provides a turbulence specific dissipation, (omega) inlet condition based on...
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:65
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
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
dimensionedScalar pos(const dimensionedScalar &ds)
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 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:161
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
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 operator=(const UList< Type > &)
Definition: fvPatchField.C:396
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
turbulentMixingLengthFrequencyInletFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:363
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451