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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2017-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "fvPatchFieldMapper.H"
32 #include "surfaceFields.H"
33 #include "volFields.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
42 )
43 :
44  inletOutletFvPatchScalarField(p, iF),
45  intensity_(0.0),
46  UName_("U")
47 {
48  this->refValue() = 0.0;
49  this->refGrad() = 0.0;
50  this->valueFraction() = 0.0;
51 }
52 
55 (
57  const fvPatch& p,
59  const fvPatchFieldMapper& mapper
60 )
61 :
62  inletOutletFvPatchScalarField(ptf, p, iF, mapper),
63  intensity_(ptf.intensity_),
64  UName_(ptf.UName_)
65 {}
66 
69 (
70  const fvPatch& p,
72  const dictionary& dict
73 )
74 :
75  inletOutletFvPatchScalarField(p, iF),
76  intensity_(dict.get<scalar>("intensity")),
77  UName_(dict.getOrDefault<word>("U", "U"))
78 {
79  this->patchType() = dict.getOrDefault<word>("patchType", word::null);
80  this->phiName_ = dict.getOrDefault<word>("phi", "phi");
81 
82  if (intensity_ < 0 || intensity_ > 1)
83  {
85  << "Turbulence intensity should be specified as a fraction 0-1 "
86  "of the mean velocity\n"
87  " value given is " << intensity_ << nl
88  << " on patch " << this->patch().name()
89  << " of field " << this->internalField().name()
90  << " in file " << this->internalField().objectPath()
91  << exit(FatalError);
92  }
93 
94  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
95 
96  this->refValue() = 0.0;
97  this->refGrad() = 0.0;
98  this->valueFraction() = 0.0;
99 }
100 
103 (
105 )
106 :
107  inletOutletFvPatchScalarField(ptf),
108  intensity_(ptf.intensity_),
109  UName_(ptf.UName_)
110 {}
111 
112 
115 (
118 )
119 :
120  inletOutletFvPatchScalarField(ptf, iF),
121  intensity_(ptf.intensity_),
122  UName_(ptf.UName_)
123 {}
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
130 {
131  if (updated())
132  {
133  return;
134  }
135 
136  const fvPatchVectorField& Up =
137  patch().lookupPatchField<volVectorField, vector>(UName_);
138 
139  const fvsPatchScalarField& phip =
140  patch().lookupPatchField<surfaceScalarField, scalar>(this->phiName_);
141 
142  this->refValue() = 1.5*sqr(intensity_)*magSqr(Up);
143  this->valueFraction() = 1.0 - pos0(phip);
144 
145  inletOutletFvPatchScalarField::updateCoeffs();
146 }
147 
148 
150 (
151  Ostream& os
152 ) const
153 {
155  os.writeEntry("intensity", intensity_);
156  os.writeEntryIfDifferent<word>("U", "U", UName_);
157  os.writeEntryIfDifferent<word>("phi", "phi", this->phiName_);
158  writeEntry("value", os);
159 }
160 
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 namespace Foam
165 {
167  (
170  );
171 }
172 
173 // ************************************************************************* //
Foam::fvPatchField< vector >
volFields.H
Foam::fvPatchField< scalar >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
turbulentIntensityKineticEnergyInletFvPatchScalarField.H
Foam::Ostream::writeEntryIfDifferent
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:248
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::turbulentIntensityKineticEnergyInletFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: turbulentIntensityKineticEnergyInletFvPatchScalarField.C:129
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::turbulentIntensityKineticEnergyInletFvPatchScalarField::turbulentIntensityKineticEnergyInletFvPatchScalarField
turbulentIntensityKineticEnergyInletFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: turbulentIntensityKineticEnergyInletFvPatchScalarField.C:39
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::turbulentIntensityKineticEnergyInletFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: turbulentIntensityKineticEnergyInletFvPatchScalarField.C:150
Foam::turbulentIntensityKineticEnergyInletFvPatchScalarField
This boundary condition provides a turbulent kinetic energy condition, based on user-supplied turbule...
Definition: turbulentIntensityKineticEnergyInletFvPatchScalarField.H:122
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54