Loading...
Searching...
No Matches
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-2025 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>
37Foam::turbulentDigitalFilterInletFvPatchField<Type>::calcPatchNormal() const
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>
59void Foam::turbulentDigitalFilterInletFvPatchField<Type>::initialisePatch()
60{
61 L_.initialise();
62
63 AMIPtr_->calculate(this->patch().patch(), L_.patch());
64
65 patchNormal_ = calcPatchNormal();
66}
67
68
69template<class Type>
70void Foam::turbulentDigitalFilterInletFvPatchField<Type>::mapL
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>
103void Foam::turbulentDigitalFilterInletFvPatchField<Type>::mapR
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>
120void Foam::turbulentDigitalFilterInletFvPatchField<Type>::mapR
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>
155void Foam::turbulentDigitalFilterInletFvPatchField<Type>::mapMean
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>
167void Foam::turbulentDigitalFilterInletFvPatchField<Type>::mapMean
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 = gWeightedAverage(this->patch().magSf(), mean);
178
179 const scalar correct
180 (
181 gSum((bulk & patchNormal_)*this->patch().magSf())
182 /gSum(mean & -this->patch().Sf())
183 );
184
185 // Map mean field onto patch for vector fields
186 fld += mean;
187
188 // Correct flow rate
190}
191
192
193// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
194
195template<class Type>
198(
199 const fvPatch& p,
201)
202:
203 fixedValueFvPatchField<Type>(p, iF),
204 AMIPtr_(new faceAreaWeightAMI(true, false)),
205 meanPtr_(nullptr),
206 Rptr_(nullptr),
207 curTimeIndex_(-1),
208 patchNormal_(Zero),
209 L_(p)
210{}
211
212
213template<class Type>
216(
218 const fvPatch& p,
220 const fvPatchFieldMapper& mapper
221)
222:
223 fixedValueFvPatchField<Type>(ptf, p, iF, mapper),
224 AMIPtr_(ptf.AMIPtr_.clone()),
225 meanPtr_(ptf.meanPtr_.clone(this->patch().patch())),
226 Rptr_(ptf.Rptr_.clone(this->patch().patch())),
227 curTimeIndex_(ptf.curTimeIndex_),
228 patchNormal_(ptf.patchNormal_),
229 L_(p, ptf.L_)
230{}
231
232
233template<class Type>
236(
237 const fvPatch& p,
239 const dictionary& dict
240)
241:
242 fixedValueFvPatchField<Type>(p, iF, dict),
243 AMIPtr_
244 (
246 (
247 dict.getOrDefault("AMIMethod", faceAreaWeightAMI::typeName),
248 dict,
249 true // flipNormals
250 )
251 ),
252 meanPtr_
253 (
254 PatchFunction1<Type>::New
255 (
256 this->patch().patch(),
257 "mean",
258 dict
259 )
260 ),
261 Rptr_
262 (
263 PatchFunction1<TypeR>::New
264 (
265 this->patch().patch(),
266 "R",
267 dict
268 )
269 ),
270 curTimeIndex_(-1),
271 patchNormal_(Zero),
272 L_(p, dict)
273{
275
276 // Check if varying or fixed time-step computation
277 if (!L_.fsm() && this->db().time().isAdjustTimeStep())
278 {
279 WarningInFunction
280 << "Varying time-step computations are not "
281 << "supported by the digital filter method."
282 << endl;
283 }
284
285 const scalar t = this->db().time().timeOutputValue();
300 AMIPtr_(ptf.AMIPtr_.clone()),
301 meanPtr_(ptf.meanPtr_.clone(this->patch().patch())),
302 Rptr_(ptf.Rptr_.clone(this->patch().patch())),
303 curTimeIndex_(ptf.curTimeIndex_),
304 patchNormal_(ptf.patchNormal_),
305 L_(ptf.L_)
306{}
307
308
309template<class Type>
312(
315)
316:
317 fixedValueFvPatchField<Type>(ptf, iF),
318 AMIPtr_(ptf.AMIPtr_.clone()),
319 meanPtr_(ptf.meanPtr_.clone(this->patch().patch())),
320 Rptr_(ptf.Rptr_.clone(this->patch().patch())),
321 curTimeIndex_(ptf.curTimeIndex_),
322 patchNormal_(ptf.patchNormal_),
323 L_(ptf.L_)
324{}
325
326
327// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
328
329template<class Type>
331(
332 const fvPatchFieldMapper& m
333)
334{
336
337 if (meanPtr_)
338 {
339 meanPtr_->autoMap(m);
340 }
341 if (Rptr_)
343 Rptr_->autoMap(m);
344 }
345}
346
347
348template<class Type>
350(
351 const fvPatchField<Type>& ptf,
352 const labelList& addr
353)
354{
356
357 const auto& dfmptf =
359
360 if (meanPtr_)
361 {
362 meanPtr_->rmap(dfmptf.meanPtr_(), addr);
363 }
364 if (Rptr_)
366 Rptr_->rmap(dfmptf.Rptr_(), addr);
367 }
368}
369
370
371template<class Type>
373{
374 if (this->updated())
375 {
376 return;
377 }
378
379 if (curTimeIndex_ == -1)
380 {
381 initialisePatch();
382 }
383
384 if (curTimeIndex_ != this->db().time().timeIndex())
385 {
386 Field<Type>& fld = *this;
387 fld = Zero;
388
389 mapL(fld);
390
391 mapR(fld);
392
393 mapMean(fld);
394
395 curTimeIndex_ = this->db().time().timeIndex();
397
399}
400
401
402template<class Type>
404(
405 Ostream& os
406) const
407{
409
410 if (meanPtr_)
411 {
412 meanPtr_->writeData(os);
413 }
414 if (Rptr_)
415 {
416 Rptr_->writeData(os);
417 }
418 if (AMIPtr_)
419 {
420 AMIPtr_->write(os);
421 }
422 L_.write(os);
423
425}
426
427
428// ************************************************************************* //
#define R(A, B, C, D, E, F, K, M)
Macros for easy insertion into run-time selection tables.
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-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 that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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...
fixedValueFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
const objectRegistry & db() const
The associated objectRegistry.
const fvPatch & patch() const noexcept
Return the patch.
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.
static tmp< fvPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
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.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
static int debug
Flag to activate debug statements.
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 write(Ostream &) const
Write boundary condition settings.
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.
turbulentDigitalFilterInletFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual tmp< fvPatchField< Type > > clone() const
Return a clone.
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.
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Type gSum(const FieldField< Field, Type > &f)
Type gWeightedAverage(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted average of a field, using the mag() of the weights.
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Vector< scalar > vector
Definition vector.H:57
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition symmTensor.H:55
label timeIndex
dictionary dict
fvOptions correct(rho)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299