turbulentInletFvPatchField.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 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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 template<class Type>
37 (
38  const fvPatch& p,
40 )
41 :
43  ranGen_(label(0)),
44  fluctuationScale_(pTraits<Type>::zero),
45  referenceField_(p.size()),
46  alpha_(0.1),
47  curTimeIndex_(-1)
48 {}
49 
50 
51 template<class Type>
53 (
55  const fvPatch& p,
57  const fvPatchFieldMapper& mapper
58 )
59 :
60  fixedValueFvPatchField<Type>(ptf, p, iF, mapper),
61  ranGen_(label(0)),
62  fluctuationScale_(ptf.fluctuationScale_),
63  referenceField_(ptf.referenceField_, mapper),
64  alpha_(ptf.alpha_),
65  curTimeIndex_(-1)
66 {}
67 
68 
69 template<class Type>
71 (
72  const fvPatch& p,
74  const dictionary& dict
75 )
76 :
78  ranGen_(label(0)),
79  fluctuationScale_(pTraits<Type>(dict.lookup("fluctuationScale"))),
80  referenceField_("referenceField", dict, p.size()),
81  alpha_(dict.lookupOrDefault<scalar>("alpha", 0.1)),
82  curTimeIndex_(-1)
83 {
84  if (dict.found("value"))
85  {
87  (
88  Field<Type>("value", dict, p.size())
89  );
90  }
91  else
92  {
94  }
95 }
96 
97 
98 template<class Type>
100 (
102 )
103 :
105  ranGen_(ptf.ranGen_),
106  fluctuationScale_(ptf.fluctuationScale_),
107  referenceField_(ptf.referenceField_),
108  alpha_(ptf.alpha_),
109  curTimeIndex_(-1)
110 {}
111 
112 
113 template<class Type>
115 (
118 )
119 :
121  ranGen_(ptf.ranGen_),
122  fluctuationScale_(ptf.fluctuationScale_),
123  referenceField_(ptf.referenceField_),
124  alpha_(ptf.alpha_),
125  curTimeIndex_(-1)
126 {}
127 
128 
129 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
130 
131 template<class Type>
133 (
134  const fvPatchFieldMapper& m
135 )
136 {
138  referenceField_.autoMap(m);
139 }
140 
141 
142 template<class Type>
144 (
145  const fvPatchField<Type>& ptf,
146  const labelList& addr
147 )
148 {
150 
151  const turbulentInletFvPatchField<Type>& tiptf =
152  refCast<const turbulentInletFvPatchField<Type> >(ptf);
153 
154  referenceField_.rmap(tiptf.referenceField_, addr);
155 }
156 
157 
158 template<class Type>
160 {
161  if (this->updated())
162  {
163  return;
164  }
165 
166  if (curTimeIndex_ != this->db().time().timeIndex())
167  {
168  Field<Type>& patchField = *this;
169 
170  Field<Type> randomField(this->size());
171 
172  forAll(patchField, facei)
173  {
174  ranGen_.randomise(randomField[facei]);
175  }
176 
177  // Correction-factor to compensate for the loss of RMS fluctuation
178  // due to the temporal correlation introduced by the alpha parameter.
179  scalar rmsCorr = sqrt(12*(2*alpha_ - sqr(alpha_)))/alpha_;
180 
181  patchField =
182  (1 - alpha_)*patchField
183  + alpha_*
184  (
185  referenceField_
186  + rmsCorr*cmptMultiply
187  (
188  randomField - 0.5*pTraits<Type>::one,
189  fluctuationScale_
190  )*mag(referenceField_)
191  );
192 
193  curTimeIndex_ = this->db().time().timeIndex();
194  }
195 
197 }
198 
199 
200 template<class Type>
202 {
204  os.writeKeyword("fluctuationScale")
205  << fluctuationScale_ << token::END_STATEMENT << nl;
206  referenceField_.writeEntry("referenceField", os);
207  os.writeKeyword("alpha") << alpha_ << token::END_STATEMENT << nl;
208  this->writeEntry("value", os);
209 }
210 
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 } // End namespace Foam
215 
216 // ************************************************************************* //
dimensionedScalar sqrt(const dimensionedScalar &ds)
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
turbulentInletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:301
label timeIndex
Definition: getTimeIndex.H:4
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
void size(const label)
Override size to be inconsistent with allocated storage.
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.
dictionary dict
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
static const char nl
Definition: Ostream.H:260
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:291
volScalarField & p
Definition: createFields.H:51
#define forAll(list, i)
Definition: UList.H:421
virtual void write(Ostream &) const
Write.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Pre-declare SubField and related Field type.
Definition: Field.H:57
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: fvPatchField.C:230
Traits class for primitives.
Definition: pTraits.H:50
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
virtual void operator==(const fvPatchField< Type > &)
Definition: fvPatchField.C:559
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
This boundary condition generates a fluctuating inlet condition by adding a random component to a ref...