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-2019 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 )
36 :
38  ranGen_(label(0)),
39  fluctuationScale_(Zero),
40  referenceField_(p.size()),
41  alpha_(0.1),
42  curTimeIndex_(-1)
43 {}
44 
45 
46 template<class Type>
48 (
49  const fvPatch& p,
51  const dictionary& dict
52 )
53 :
54  fixedValueFvPatchField<Type>(p, iF, dict, false),
55  ranGen_(label(0)),
56  fluctuationScale_(pTraits<Type>(dict.lookup("fluctuationScale"))),
57  referenceField_("referenceField", dict, p.size()),
58  alpha_(dict.lookupOrDefault<scalar>("alpha", 0.1)),
59  curTimeIndex_(-1)
60 {
61  if (dict.found("value"))
62  {
64  (
65  Field<Type>("value", dict, p.size())
66  );
67  }
68  else
69  {
71  }
72 }
73 
74 
75 template<class Type>
77 (
79  const fvPatch& p,
81  const fvPatchFieldMapper& mapper
82 )
83 :
84  fixedValueFvPatchField<Type>(ptf, p, iF, mapper),
85  ranGen_(label(0)),
86  fluctuationScale_(ptf.fluctuationScale_),
87  referenceField_(mapper(ptf.referenceField_)),
88  alpha_(ptf.alpha_),
89  curTimeIndex_(-1)
90 {}
91 
92 
93 template<class Type>
95 (
97 )
98 :
100  ranGen_(ptf.ranGen_),
101  fluctuationScale_(ptf.fluctuationScale_),
102  referenceField_(ptf.referenceField_),
103  alpha_(ptf.alpha_),
104  curTimeIndex_(-1)
105 {}
106 
107 
108 template<class Type>
110 (
113 )
114 :
116  ranGen_(ptf.ranGen_),
117  fluctuationScale_(ptf.fluctuationScale_),
118  referenceField_(ptf.referenceField_),
119  alpha_(ptf.alpha_),
120  curTimeIndex_(-1)
121 {}
122 
123 
124 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125 
126 template<class Type>
128 (
129  const fvPatchFieldMapper& m
130 )
131 {
133  m(referenceField_, referenceField_);
134 }
135 
136 
137 template<class Type>
139 (
140  const fvPatchField<Type>& ptf,
141  const labelList& addr
142 )
143 {
145 
146  const turbulentInletFvPatchField<Type>& tiptf =
147  refCast<const turbulentInletFvPatchField<Type>>(ptf);
148 
149  referenceField_.rmap(tiptf.referenceField_, addr);
150 }
151 
152 
153 template<class Type>
155 {
156  if (this->updated())
157  {
158  return;
159  }
160 
161  if (curTimeIndex_ != this->db().time().timeIndex())
162  {
163  Field<Type>& patchField = *this;
164 
165  Field<Type> randomField(this->size());
166 
167  forAll(patchField, facei)
168  {
169  randomField[facei] = ranGen_.sample01<Type>();
170  }
171 
172  // Correction-factor to compensate for the loss of RMS fluctuation
173  // due to the temporal correlation introduced by the alpha parameter.
174  scalar rmsCorr = sqrt(12*(2*alpha_ - sqr(alpha_)))/alpha_;
175 
176  patchField =
177  (1 - alpha_)*patchField
178  + alpha_*
179  (
180  referenceField_
181  + rmsCorr*cmptMultiply
182  (
183  randomField - 0.5*pTraits<Type>::one,
184  fluctuationScale_
185  )*mag(referenceField_)
186  );
187 
188  curTimeIndex_ = this->db().time().timeIndex();
189  }
190 
192 }
193 
194 
195 template<class Type>
197 {
199  writeEntry(os, "fluctuationScale", fluctuationScale_);
200  writeEntry(os, "referenceField", referenceField_);
201  writeEntry(os, "alpha", alpha_);
202  writeEntry(os, "value", *this);
203 }
204 
205 
206 // ************************************************************************* //
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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.
Pre-declare SubField and related Field type.
Definition: Field.H:56
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
turbulentInletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
This boundary condition generates a fluctuating inlet condition by adding a random component to a ref...
Foam::fvPatchFieldMapper.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
static const zero Zero
Definition: zero.H:97
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
label timeIndex
Definition: getTimeIndex.H:4
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField & p