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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29
30// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31
32template<class Type>
34(
35 const fvPatch& p,
37)
38:
39 fixedValueFvPatchField<Type>(p, iF),
40 ranGen_(label(0)),
41 fluctuationScale_(Zero),
42 referenceField_(p.size()),
43 alpha_(0.1),
44 curTimeIndex_(-1)
45{}
46
47
48template<class Type>
50(
51 const fvPatch& p,
53 const dictionary& dict
54)
55:
56 fixedValueFvPatchField<Type>(p, iF, dict, false),
57 ranGen_(label(0)),
58 fluctuationScale_(dict.get<Type>("fluctuationScale")),
59 referenceField_("referenceField", dict, p.size()),
60 alpha_(dict.getOrDefault<scalar>("alpha", 0.1)),
61 curTimeIndex_(-1)
62{
63 if (dict.found("value"))
64 {
66 (
67 Field<Type>("value", dict, p.size())
68 );
69 }
70 else
71 {
73 }
74}
75
76
77template<class Type>
79(
81 const fvPatch& p,
83 const fvPatchFieldMapper& mapper
84)
85:
86 fixedValueFvPatchField<Type>(ptf, p, iF, mapper),
87 ranGen_(label(0)),
88 fluctuationScale_(ptf.fluctuationScale_),
89 referenceField_(ptf.referenceField_, mapper),
90 alpha_(ptf.alpha_),
91 curTimeIndex_(-1)
92{}
93
94
95template<class Type>
97(
99)
100:
101 fixedValueFvPatchField<Type>(ptf),
102 ranGen_(ptf.ranGen_),
103 fluctuationScale_(ptf.fluctuationScale_),
104 referenceField_(ptf.referenceField_),
105 alpha_(ptf.alpha_),
106 curTimeIndex_(-1)
107{}
108
109
110template<class Type>
112(
115)
116:
117 fixedValueFvPatchField<Type>(ptf, iF),
118 ranGen_(ptf.ranGen_),
119 fluctuationScale_(ptf.fluctuationScale_),
120 referenceField_(ptf.referenceField_),
121 alpha_(ptf.alpha_),
122 curTimeIndex_(-1)
123{}
124
125
126// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127
128template<class Type>
130(
131 const fvPatchFieldMapper& m
132)
133{
135 referenceField_.autoMap(m);
136}
137
138
139template<class Type>
141(
142 const fvPatchField<Type>& ptf,
143 const labelList& addr
144)
145{
147
149 refCast<const turbulentInletFvPatchField<Type>>(ptf);
150
151 referenceField_.rmap(tiptf.referenceField_, addr);
152}
153
154
155template<class Type>
157{
158 if (this->updated())
159 {
160 return;
161 }
162
163 if (curTimeIndex_ != this->db().time().timeIndex())
164 {
165 Field<Type>& patchField = *this;
166
167 Field<Type> randomField(this->size());
168
169 forAll(patchField, facei)
170 {
171 ranGen_.randomise01<Type>(randomField[facei]);
172 }
173
174 // Correction-factor to compensate for the loss of RMS fluctuation
175 // due to the temporal correlation introduced by the alpha parameter.
176 scalar rmsCorr = sqrt(12*(2*alpha_ - sqr(alpha_)))/alpha_;
177
178 patchField =
179 (1 - alpha_)*patchField
180 + alpha_*
181 (
182 referenceField_
183 + rmsCorr*cmptMultiply
184 (
185 randomField - 0.5*pTraits<Type>::one,
186 fluctuationScale_
187 )*mag(referenceField_)
188 );
189
190 curTimeIndex_ = this->db().time().timeIndex();
191 }
192
194}
195
196
197template<class Type>
199{
201 os.writeEntry("fluctuationScale", fluctuationScale_);
202 referenceField_.writeEntry("referenceField", os);
203 os.writeEntry("alpha", alpha_);
204 this->writeEntry("value", os);
205}
206
207
208// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type.
Definition: Field.H:82
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
void rmap(const atmBoundaryLayer &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:62
friend bool operator==(const refineCell &rc1, const refineCell &rc2)
Definition: refineCell.H:97
This boundary condition produces spatiotemporal-variant field by summing a set of pseudo-random numbe...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
label timeIndex
Definition: getTimeIndex.H:30
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333