turbulentInletFvPatchField.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-2024 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 
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class Type>
32 (
33  const fvPatch& p,
35  const dictionary& dict
36 )
37 :
38  fixedValueFvPatchField<Type>(p, iF, dict, false),
39  ranGen_(label(0)),
40  fluctuationScale_(dict.lookup<Type>("fluctuationScale", unitFraction)),
41  referenceField_("referenceField", iF.dimensions(), dict, p.size()),
42  alpha_(dict.lookupOrDefault<scalar>("alpha", unitFraction, 0.1)),
43  curTimeIndex_(-1)
44 {
45  if (dict.found("value"))
46  {
48  (
49  Field<Type>("value", iF.dimensions(), dict, p.size())
50  );
51  }
52  else
53  {
55  }
56 }
57 
58 
59 template<class Type>
61 (
63  const fvPatch& p,
65  const fieldMapper& mapper
66 )
67 :
68  fixedValueFvPatchField<Type>(ptf, p, iF, mapper),
69  ranGen_(label(0)),
70  fluctuationScale_(ptf.fluctuationScale_),
71  referenceField_(mapper(ptf.referenceField_)),
72  alpha_(ptf.alpha_),
73  curTimeIndex_(-1)
74 {}
75 
76 
77 template<class Type>
79 (
82 )
83 :
84  fixedValueFvPatchField<Type>(ptf, iF),
85  ranGen_(ptf.ranGen_),
86  fluctuationScale_(ptf.fluctuationScale_),
87  referenceField_(ptf.referenceField_),
88  alpha_(ptf.alpha_),
89  curTimeIndex_(-1)
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
95 template<class Type>
97 (
98  const fvPatchField<Type>& ptf,
99  const fieldMapper& mapper
100 )
101 {
103 
104  const turbulentInletFvPatchField<Type>& tiptf =
105  refCast<const turbulentInletFvPatchField<Type>>(ptf);
106 
107  mapper(referenceField_, tiptf.referenceField_);
108 }
109 
110 
111 template<class Type>
113 (
114  const fvPatchField<Type>& ptf
115 )
116 {
118 
119  const turbulentInletFvPatchField<Type>& tiptf =
120  refCast<const turbulentInletFvPatchField<Type>>(ptf);
121 
122  referenceField_.reset(tiptf.referenceField_);
123 }
124 
125 
126 template<class Type>
128 {
129  if (this->updated())
130  {
131  return;
132  }
133 
134  if (curTimeIndex_ != this->db().time().timeIndex())
135  {
136  Field<Type>& patchField = *this;
137 
138  Field<Type> randomField(this->size());
139 
140  forAll(patchField, facei)
141  {
142  randomField[facei] = ranGen_.sample01<Type>();
143  }
144 
145  // Correction-factor to compensate for the loss of RMS fluctuation
146  // due to the temporal correlation introduced by the alpha parameter.
147  scalar rmsCorr = sqrt(12*(2*alpha_ - sqr(alpha_)))/alpha_;
148 
149  patchField =
150  (1 - alpha_)*patchField
151  + alpha_*
152  (
153  referenceField_
154  + rmsCorr*cmptMultiply
155  (
156  randomField - 0.5*pTraits<Type>::one,
157  fluctuationScale_
158  )*mag(referenceField_)
159  );
160 
161  curTimeIndex_ = this->db().time().timeIndex();
162  }
163 
165 }
166 
167 
168 template<class Type>
170 {
172  writeEntry(os, "fluctuationScale", fluctuationScale_);
173  writeEntry(os, "referenceField", referenceField_);
174  writeEntry(os, "alpha", alpha_);
175  writeEntry(os, "value", *this);
176 }
177 
178 
179 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Pre-declare SubField and related Field type.
Definition: Field.H:83
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
virtual void operator==(const fvPatchField< Type > &)
Definition: fvPatchField.C:415
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:202
virtual void reset(const fvPatchField< Type > &)
Reset the fvPatchField to the given fvPatchField.
Definition: fvPatchField.C:195
virtual void map(const fvPatchField< Type > &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:185
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Traits class for primitives.
Definition: pTraits.H:53
This boundary condition generates a fluctuating inlet condition by adding a random component to a ref...
virtual void write(Ostream &) const
Write.
turbulentInletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void reset(const fvPatchField< Type > &)
Reset the fvPatchField to the given fvPatchField.
virtual void map(const fvPatchField< Type > &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
const unitConversion unitFraction
label timeIndex
Definition: getTimeIndex.H:4
dictionary dict
volScalarField & p