turbulentIntensityKineticEnergyInletFvPatchScalarField.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-2012 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 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39 )
40 :
41  inletOutletFvPatchScalarField(p, iF),
42  intensity_(0.0),
43  UName_("U")
44 {
45  this->refValue() = 0.0;
46  this->refGrad() = 0.0;
47  this->valueFraction() = 0.0;
48 }
49 
52 (
54  const fvPatch& p,
56  const fvPatchFieldMapper& mapper
57 )
58 :
59  inletOutletFvPatchScalarField(ptf, p, iF, mapper),
60  intensity_(ptf.intensity_),
61  UName_(ptf.UName_)
62 {}
63 
66 (
67  const fvPatch& p,
69  const dictionary& dict
70 )
71 :
72  inletOutletFvPatchScalarField(p, iF),
73  intensity_(readScalar(dict.lookup("intensity"))),
74  UName_(dict.lookupOrDefault<word>("U", "U"))
75 {
76  this->phiName_ = dict.lookupOrDefault<word>("phi", "phi");
77 
78  if (intensity_ < 0 || intensity_ > 1)
79  {
81  (
82  "turbulentIntensityKineticEnergyInletFvPatchScalarField::"
83  "turbulentIntensityKineticEnergyInletFvPatchScalarField"
84  "("
85  "const fvPatch&, "
86  "const DimensionedField<scalar, volMesh>&, "
87  "const dictionary&"
88  ")"
89  ) << "Turbulence intensity should be specified as a fraction 0-1 "
90  "of the mean velocity\n"
91  " value given is " << intensity_ << nl
92  << " on patch " << this->patch().name()
93  << " of field " << this->dimensionedInternalField().name()
94  << " in file " << this->dimensionedInternalField().objectPath()
95  << exit(FatalError);
96  }
97 
98  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
99 
100  this->refValue() = 0.0;
101  this->refGrad() = 0.0;
102  this->valueFraction() = 0.0;
103 }
104 
107 (
109 )
110 :
111  inletOutletFvPatchScalarField(ptf),
112  intensity_(ptf.intensity_),
113  UName_(ptf.UName_)
114 {}
115 
116 
119 (
122 )
123 :
124  inletOutletFvPatchScalarField(ptf, iF),
125  intensity_(ptf.intensity_),
126  UName_(ptf.UName_)
127 {}
128 
129 
130 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
131 
134 {
135  if (updated())
136  {
137  return;
138  }
139 
140  const fvPatchVectorField& Up =
141  patch().lookupPatchField<volVectorField, vector>(UName_);
142 
143  const fvsPatchScalarField& phip =
144  patch().lookupPatchField<surfaceScalarField, scalar>(this->phiName_);
145 
146  this->refValue() = 1.5*sqr(intensity_)*magSqr(Up);
147  this->valueFraction() = 1.0 - pos(phip);
148 
149  inletOutletFvPatchScalarField::updateCoeffs();
150 }
151 
152 
154 (
155  Ostream& os
156 ) const
157 {
159  os.writeKeyword("intensity") << intensity_ << token::END_STATEMENT << nl;
160  writeEntryIfDifferent<word>(os, "U", "U", UName_);
161  writeEntryIfDifferent<word>(os, "phi", "phi", this->phiName_);
162  writeEntry("value", os);
163 }
164 
165 
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
167 
168 namespace Foam
169 {
171  (
174  );
175 }
176 
177 // ************************************************************************* //
Foam::surfaceFields.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
#define readScalar
Definition: doubleScalar.C:38
dimensioned< scalar > magSqr(const dimensioned< Type > &)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from string.
Definition: word.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::fvPatchFieldMapper.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:354
Namespace for OpenFOAM.
static const char nl
Definition: Ostream.H:260
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Macros for easy insertion into run-time selection tables.
This boundary condition provides a turbulent kinetic energy condition, based on user-supplied turbule...
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
rDeltaT dimensionedInternalField()
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
error FatalError
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
dimensionedScalar pos(const dimensionedScalar &ds)
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField scalarField(fieldObject, mesh)
virtual label size() const
Return size.
Definition: fvPatch.H:161
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
turbulentIntensityKineticEnergyInletFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.