Loading...
Searching...
No Matches
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, IOobjectOption::NO_READ),
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 (!this->readValueEntry(dict))
64 {
65 fvPatchField<Type>::operator=(referenceField_);
66 }
67}
68
69
70template<class Type>
72(
73 const turbulentInletFvPatchField<Type>& ptf,
74 const fvPatch& p,
75 const DimensionedField<Type, volMesh>& iF,
76 const fvPatchFieldMapper& mapper
77)
78:
79 fixedValueFvPatchField<Type>(ptf, p, iF, mapper),
80 ranGen_(label(0)),
81 fluctuationScale_(ptf.fluctuationScale_),
82 referenceField_(ptf.referenceField_, mapper),
83 alpha_(ptf.alpha_),
84 curTimeIndex_(-1)
85{}
86
87
88template<class Type>
90(
92)
93:
94 fixedValueFvPatchField<Type>(ptf),
95 ranGen_(ptf.ranGen_),
96 fluctuationScale_(ptf.fluctuationScale_),
97 referenceField_(ptf.referenceField_),
98 alpha_(ptf.alpha_),
99 curTimeIndex_(-1)
100{}
101
102
103template<class Type>
105(
108)
109:
110 fixedValueFvPatchField<Type>(ptf, iF),
111 ranGen_(ptf.ranGen_),
112 fluctuationScale_(ptf.fluctuationScale_),
113 referenceField_(ptf.referenceField_),
114 alpha_(ptf.alpha_),
115 curTimeIndex_(-1)
116{}
117
118
119// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120
121template<class Type>
123(
124 const fvPatchFieldMapper& m
125)
128 referenceField_.autoMap(m);
129}
130
131
132template<class Type>
134(
135 const fvPatchField<Type>& ptf,
136 const labelList& addr
137)
138{
140
143
144 referenceField_.rmap(tiptf.referenceField_, addr);
145}
146
147
148template<class Type>
150{
151 if (this->updated())
152 {
153 return;
154 }
155
156 if (curTimeIndex_ != this->db().time().timeIndex())
157 {
158 Field<Type>& patchField = *this;
159
160 Field<Type> randomField(this->size());
161
162 forAll(patchField, facei)
163 {
164 ranGen_.randomise01<Type>(randomField[facei]);
165 }
166
167 // Correction-factor to compensate for the loss of RMS fluctuation
168 // due to the temporal correlation introduced by the alpha parameter.
169 scalar rmsCorr = sqrt(12*(2*alpha_ - sqr(alpha_)))/alpha_;
170
171 patchField =
172 (1 - alpha_)*patchField
173 + alpha_*
174 (
175 referenceField_
176 + rmsCorr*cmptMultiply
177 (
178 randomField - 0.5*pTraits<Type>::one,
179 fluctuationScale_
180 )*mag(referenceField_)
181 );
182
183 curTimeIndex_ = this->db().time().timeIndex();
185
187}
188
189
190template<class Type>
192{
194 os.writeEntry("fluctuationScale", fluctuationScale_);
195 referenceField_.writeEntry("referenceField", os);
196 os.writeEntry("alpha", alpha_);
198}
199
200
201// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
bool get(const label i) const
Definition UList.H:868
void size(const label n)
Definition UList.H:118
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
fixedValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
const objectRegistry & db() const
The associated objectRegistry.
bool updated() const noexcept
True if the boundary condition has already been updated.
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...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual void operator=(const UList< Type > &)
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.
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
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 write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
turbulentInletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
label timeIndex
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299