turbulentDigitalFilterInletFvPatchField.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) 2019-2022 OpenCFD Ltd.
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
30#include "faceAreaWeightAMI.H"
32
33// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34
35template<class Type>
38{
39 const vectorField nf(this->patch().nf());
40
41 // Patch normal points into domain
42 vector patchNormal(-gAverage(nf));
43
44 // Check that patch is planar
45 const scalar error = max(magSqr(patchNormal + nf));
46
47 if (error > SMALL)
48 {
50 << "Patch " << this->patch().name() << " is not planar"
51 << endl;
52 }
53
54 return patchNormal.normalise();
55}
56
57
58template<class Type>
60{
61 L_.initialise();
62
63 AMIPtr_->calculate(this->patch().patch(), L_.patch());
64
65 patchNormal_ = calcPatchNormal();
66}
67
68
69template<class Type>
71(
72 Field<Type>& fld
73)
74{
75 Field<Type> sourceFld;
76
77 if (Pstream::master())
78 {
79 sourceFld = L_.convolve();
80 L_.shift();
81 L_.refill();
82 }
83
84 // Map two-point correlations (integral scales)
85 plusEqOp<Type> cop;
86 AMIPtr_->interpolateToSource
87 (
88 sourceFld,
89 multiplyWeightedOp<Type, plusEqOp<Type>>(cop),
90 fld,
91 UList<Type>::null()
92 );
93
94 // Map forward-stepwise method correlations if requested
95 if (L_.fsm())
96 {
97 L_.correlate(fld);
98 }
99}
100
101
102template<class Type>
104(
105 scalarField& fld
106) const
107{
108 const scalar t = this->db().time().timeOutputValue();
109 scalarField R(Rptr_->value(t));
110
111 // Lund-Wu-Squires transformation for scalar fields
112 R = Foam::sqrt(R);
113
114 // Map transformed Reynolds stresses field onto patch for scalar fields
115 fld *= R;
116}
117
118
119template<class Type>
121(
122 vectorField& fld
123) const
124{
125 const scalar t = this->db().time().timeOutputValue();
126 symmTensorField R(Rptr_->value(t));
127
128 // Lund-Wu-Squires transformation for vector fields
129 for (symmTensor& r : R)
130 {
131 // (KSJ:Eq. 5)
132 r.xx() = Foam::sqrt(r.xx());
133 r.xy() /= r.xx();
134 r.xz() /= r.xx();
135 r.yy() = Foam::sqrt(r.yy() - sqr(r.xy()));
136 r.yz() = (r.yz() - r.xy()*r.xz())/r.yy();
137 r.zz() = Foam::sqrt(r.zz() - sqr(r.xz()) - sqr(r.yz()));
138 }
139
140 // Map transformed Reynolds stresses field onto patch for vector fields
141 forAll(fld, i)
142 {
143 vector& u = fld[i];
144 const symmTensor& r = R[i];
145
146 // (KSJ:p. 658, item-e)
147 u.z() = u.x()*r.xz() + u.y()*r.yz() + u.z()*r.zz();
148 u.y() = u.x()*r.xy() + u.y()*r.yy();
149 u.x() = u.x()*r.xx();
150 }
151}
152
153
154template<class Type>
156(
157 scalarField& fld
158) const
159{
160 const scalar t = this->db().time().timeOutputValue();
161
162 fld += meanPtr_->value(t);
163}
164
165
166template<class Type>
168(
169 vectorField& fld
170) const
171{
172 const scalar t = this->db().time().timeOutputValue();
173 tmp<vectorField> tmean = meanPtr_->value(t);
174 const vectorField& mean = tmean.cref();
175
176 // Calculate flow-rate correction factor for vector fields (KCX:Eq. 8)
177 const vector bulk
178 (
179 gSum(mean*this->patch().magSf())
180 /(gSum(this->patch().magSf()) + ROOTVSMALL)
181 );
182
183 const scalar correct
184 (
185 gSum((bulk & patchNormal_)*this->patch().magSf())
186 /gSum(mean & -this->patch().Sf())
187 );
188
189 // Map mean field onto patch for vector fields
190 fld += mean;
191
192 // Correct flow rate
193 fld *= correct;
194}
195
196
197// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
198
199template<class Type>
202(
203 const fvPatch& p,
205)
206:
207 fixedValueFvPatchField<Type>(p, iF),
208 AMIPtr_(new faceAreaWeightAMI(true, false)),
209 meanPtr_(nullptr),
210 Rptr_(nullptr),
211 curTimeIndex_(-1),
212 patchNormal_(Zero),
213 L_(p)
214{}
215
216
217template<class Type>
220(
222 const fvPatch& p,
224 const fvPatchFieldMapper& mapper
225)
226:
227 fixedValueFvPatchField<Type>(ptf, p, iF, mapper),
228 AMIPtr_(ptf.AMIPtr_.clone()),
229 meanPtr_(ptf.meanPtr_.clone(this->patch().patch())),
230 Rptr_(ptf.Rptr_.clone(this->patch().patch())),
231 curTimeIndex_(ptf.curTimeIndex_),
232 patchNormal_(ptf.patchNormal_),
233 L_(p, ptf.L_)
234{}
235
236
237template<class Type>
240(
241 const fvPatch& p,
243 const dictionary& dict
244)
245:
246 fixedValueFvPatchField<Type>(p, iF, dict),
247 AMIPtr_
248 (
250 (
251 dict.getOrDefault("AMIMethod", faceAreaWeightAMI::typeName),
252 dict,
253 true // flipNormals
254 )
255 ),
256 meanPtr_
257 (
258 PatchFunction1<Type>::New
259 (
260 this->patch().patch(),
261 "mean",
262 dict
263 )
264 ),
265 Rptr_
266 (
267 PatchFunction1<TypeR>::New
268 (
269 this->patch().patch(),
270 "R",
271 dict
272 )
273 ),
274 curTimeIndex_(-1),
275 patchNormal_(Zero),
276 L_(p, dict)
277{
279
280 // Check if varying or fixed time-step computation
281 if (!L_.fsm() && this->db().time().isAdjustTimeStep())
282 {
284 << "Varying time-step computations are not "
285 << "supported by the digital filter method."
286 << endl;
287 }
288
289 const scalar t = this->db().time().timeOutputValue();
290 const Field<TypeR> R(Rptr_->value(t));
291
293}
294
295
296template<class Type>
299(
301)
302:
303 fixedValueFvPatchField<Type>(ptf),
304 AMIPtr_(ptf.AMIPtr_.clone()),
305 meanPtr_(ptf.meanPtr_.clone(this->patch().patch())),
306 Rptr_(ptf.Rptr_.clone(this->patch().patch())),
307 curTimeIndex_(ptf.curTimeIndex_),
308 patchNormal_(ptf.patchNormal_),
309 L_(ptf.L_)
310{}
311
312
313template<class Type>
316(
319)
320:
321 fixedValueFvPatchField<Type>(ptf, iF),
322 AMIPtr_(ptf.AMIPtr_.clone()),
323 meanPtr_(ptf.meanPtr_.clone(this->patch().patch())),
324 Rptr_(ptf.Rptr_.clone(this->patch().patch())),
325 curTimeIndex_(ptf.curTimeIndex_),
326 patchNormal_(ptf.patchNormal_),
327 L_(ptf.L_)
328{}
329
330
331// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
332
333template<class Type>
335(
336 const fvPatchFieldMapper& m
337)
338{
340
341 if (meanPtr_)
342 {
343 meanPtr_->autoMap(m);
344 }
345 if (Rptr_)
346 {
347 Rptr_->autoMap(m);
348 }
349}
350
351
352template<class Type>
354(
355 const fvPatchField<Type>& ptf,
356 const labelList& addr
357)
358{
360
361 const auto& dfmptf =
362 refCast<const turbulentDigitalFilterInletFvPatchField<Type>>(ptf);
363
364 if (meanPtr_)
365 {
366 meanPtr_->rmap(dfmptf.meanPtr_(), addr);
367 }
368 if (Rptr_)
369 {
370 Rptr_->rmap(dfmptf.Rptr_(), addr);
371 }
372}
373
374
375template<class Type>
377{
378 if (this->updated())
379 {
380 return;
381 }
382
383 if (curTimeIndex_ == -1)
384 {
385 initialisePatch();
386 }
387
388 if (curTimeIndex_ != this->db().time().timeIndex())
389 {
390 Field<Type>& fld = *this;
391 fld = Zero;
392
393 mapL(fld);
394
395 mapR(fld);
396
397 mapMean(fld);
398
399 curTimeIndex_ = this->db().time().timeIndex();
400 }
401
403}
404
405
406template<class Type>
408(
409 Ostream& os
410) const
411{
413
414 if (meanPtr_)
415 {
416 meanPtr_->writeData(os);
417 }
418 if (Rptr_)
419 {
420 Rptr_->writeData(os);
421 }
422 if (AMIPtr_)
423 {
424 AMIPtr_->write(os);
425 }
426 L_.write(os);
427
428 this->writeEntry("value", os);
429}
430
431
432// ************************************************************************* //
#define R(A, B, C, D, E, F, K, M)
Y[inertIndex] max(0.0)
Macros for easy insertion into run-time selection tables.
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
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
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
scalar timeOutputValue() const
Return current time value.
Definition: TimeStateI.H:31
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
Face area weighted Arbitrary Mesh Interface (AMI) method.
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
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:210
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
const Time & time() const noexcept
Return time registry.
Class to describe the integral-scale container being used in the turbulentDigitalFilterInletFvPatchFi...
static void checkStresses(const symmTensorField &R)
Check if input Reynolds stresses are valid.
Digital-filter based boundary condition for vector- and scalar-based quantities (e....
virtual void rmap(const fvPatchField< Type > &ptf, const labelList &addr)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void autoMap(const fvPatchFieldMapper &m)
Map (and resize as needed) from self given a mapping object.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
volScalarField & p
thermo correct()
OBJstream os(runTime.globalPath()/outputName)
#define WarningInFunction
Report a warning using Foam::Warning.
const std::string patch
OpenFOAM patch number as a std::string.
Type gSum(const FieldField< Field, Type > &f)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition: symmTensor.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
Type gAverage(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
label timeIndex
Definition: getTimeIndex.H:30
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333