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-2017 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  << "Turbulence intensity should be specified as a fraction 0-1 "
82  "of the mean velocity\n"
83  " value given is " << intensity_ << nl
84  << " on patch " << this->patch().name()
85  << " of field " << this->internalField().name()
86  << " in file " << this->internalField().objectPath()
87  << exit(FatalError);
88  }
89 
90  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
91 
92  this->refValue() = 0.0;
93  this->refGrad() = 0.0;
94  this->valueFraction() = 0.0;
95 }
96 
99 (
101 )
102 :
103  inletOutletFvPatchScalarField(ptf),
104  intensity_(ptf.intensity_),
105  UName_(ptf.UName_)
106 {}
107 
108 
111 (
114 )
115 :
116  inletOutletFvPatchScalarField(ptf, iF),
117  intensity_(ptf.intensity_),
118  UName_(ptf.UName_)
119 {}
120 
121 
122 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
123 
126 {
127  if (updated())
128  {
129  return;
130  }
131 
132  const fvPatchVectorField& Up =
133  patch().lookupPatchField<volVectorField, vector>(UName_);
134 
135  const fvsPatchScalarField& phip =
136  patch().lookupPatchField<surfaceScalarField, scalar>(this->phiName_);
137 
138  this->refValue() = 1.5*sqr(intensity_)*magSqr(Up);
139  this->valueFraction() = 1.0 - pos0(phip);
140 
141  inletOutletFvPatchScalarField::updateCoeffs();
142 }
143 
144 
146 (
147  Ostream& os
148 ) const
149 {
151  os.writeKeyword("intensity") << intensity_ << token::END_STATEMENT << nl;
152  writeEntryIfDifferent<word>(os, "U", "U", UName_);
153  writeEntryIfDifferent<word>(os, "phi", "phi", this->phiName_);
154  writeEntry("value", os);
155 }
156 
157 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 
160 namespace Foam
161 {
163  (
166  );
167 }
168 
169 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
Foam::surfaceFields.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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
This boundary condition provides a turbulent kinetic energy condition, based on user-supplied turbule...
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Macros for easy insertion into run-time selection tables.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual label size() const
Return size.
Definition: fvPatch.H:161
dimensioned< scalar > magSqr(const dimensioned< Type > &)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
volScalarField scalarField(fieldObject, mesh)
dimensionedScalar pos0(const dimensionedScalar &ds)
static const char nl
Definition: Ostream.H:262
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
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.
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