turbulentDFSEMInletFvPatchVectorField.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) 2015 OpenFOAM Foundation
9 Copyright (C) 2016-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
31#include "momentOfInertia.H"
32#include "OFstream.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36Foam::label Foam::turbulentDFSEMInletFvPatchVectorField::seedIterMax_ = 1000;
37
38// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39
40void Foam::turbulentDFSEMInletFvPatchVectorField::writeEddyOBJ() const
41{
42 {
43 // Output the bounding box
44 OFstream os(db().time().path()/"eddyBox.obj");
45
46 const polyPatch& pp = this->patch().patch();
47 const labelList& boundaryPoints = pp.boundaryPoints();
48 const pointField& localPoints = pp.localPoints();
49
50 const vector offset(patchNormal_*maxSigmaX_);
51 forAll(boundaryPoints, i)
52 {
53 point p = localPoints[boundaryPoints[i]];
54 p += offset;
55 os << "v " << p.x() << " " << p.y() << " " << p.z() << nl;
56 }
57
58 forAll(boundaryPoints, i)
59 {
60 point p = localPoints[boundaryPoints[i]];
61 p -= offset;
62 os << "v " << p.x() << " " << p.y() << " " << p.z() << nl;
63 }
64 }
65
66 {
67 const Time& time = db().time();
68 OFstream os
69 (
70 time.path()/"eddies_" + Foam::name(time.timeIndex()) + ".obj"
71 );
72
73 label pointOffset = 0;
74 forAll(eddies_, eddyI)
75 {
76 const eddy& e = eddies_[eddyI];
77 pointOffset += e.writeSurfaceOBJ(pointOffset, patchNormal_, os);
78 }
79 }
80}
81
82
83void Foam::turbulentDFSEMInletFvPatchVectorField::writeLumleyCoeffs() const
84{
85 // Output list of xi vs eta
86
87 OFstream os(db().time().path()/"lumley_interpolated.out");
88
89 os << "# xi" << token::TAB << "eta" << endl;
90
91 const scalar t = db().time().timeOutputValue();
92 const symmTensorField R(R_->value(t)/sqr(Uref_));
93
94 forAll(R, faceI)
95 {
96 // Normalised anisotropy tensor
97 const symmTensor devR(dev(R[faceI]/(tr(R[faceI]))));
98
99 // Second tensor invariant
100 const scalar ii = min(0, invariantII(devR));
101
102 // Third tensor invariant
103 const scalar iii = invariantIII(devR);
104
105 // xi, eta
106 // See Pope - characterization of Reynolds-stress anisotropy
107 const scalar xi = cbrt(0.5*iii);
108 const scalar eta = sqrt(-ii/3.0);
109 os << xi << token::TAB << eta << token::TAB
110 << ii << token::TAB << iii << endl;
111 }
112}
113
114
115void Foam::turbulentDFSEMInletFvPatchVectorField::initialisePatch()
116{
117 const vectorField nf(patch().nf());
118
119 // Patch normal points into domain
120 patchNormal_ = -gAverage(nf);
121
122 // Check that patch is planar
123 const scalar error = max(magSqr(patchNormal_ + nf));
124
125 if (error > SMALL)
126 {
128 << "Patch " << patch().name() << " is not planar"
129 << endl;
130 }
131
132 patchNormal_ /= mag(patchNormal_) + ROOTVSMALL;
133
134
135 // Decompose the patch faces into triangles to enable point search
136
137 const polyPatch& patch = this->patch().patch();
138 const pointField& points = patch.points();
139
140 // Triangulate the patch faces and create addressing
141 DynamicList<label> triToFace(2*patch.size());
142 DynamicList<scalar> triMagSf(2*patch.size());
143 DynamicList<face> triFace(2*patch.size());
144 DynamicList<face> tris(5);
145
146 // Set zero value at the start of the tri area list
147 triMagSf.append(0.0);
148
149 forAll(patch, faceI)
150 {
151 const face& f = patch[faceI];
152
153 tris.clear();
154 f.triangles(points, tris);
155
156 forAll(tris, i)
157 {
158 triToFace.append(faceI);
159 triFace.append(tris[i]);
160 triMagSf.append(tris[i].mag(points));
161 }
162 }
163
164 sumTriMagSf_ = Zero;
165 sumTriMagSf_[Pstream::myProcNo() + 1] = sum(triMagSf);
166
167 Pstream::listCombineAllGather(sumTriMagSf_, maxEqOp<scalar>());
168
169 for (label i = 1; i < triMagSf.size(); ++i)
170 {
171 triMagSf[i] += triMagSf[i-1];
172 }
173
174 // Transfer to persistent storage
175 triFace_.transfer(triFace);
176 triToFace_.transfer(triToFace);
177 triCumulativeMagSf_.transfer(triMagSf);
178
179 // Convert sumTriMagSf_ into cumulative sum of areas per proc
180 for (label i = 1; i < sumTriMagSf_.size(); ++i)
181 {
182 sumTriMagSf_[i] += sumTriMagSf_[i-1];
183 }
184
185 // Global patch area (over all processors)
186 patchArea_ = sumTriMagSf_.last();
187
188 // Local patch bounds (this processor)
189 patchBounds_ = boundBox(patch.localPoints(), false);
190 patchBounds_.inflate(0.1);
191
192 // Determine if all eddies spawned from a single processor
193 singleProc_ = patch.size() == returnReduce(patch.size(), sumOp<label>());
194 reduce(singleProc_, orOp<bool>());
195}
196
197
198void Foam::turbulentDFSEMInletFvPatchVectorField::initialiseEddyBox()
199{
200 const scalarField& magSf = patch().magSf();
201
202 const scalarField L(L_->value(db().time().timeOutputValue())/Lref_);
203
204 // (PCF:Eq. 14)
205 const scalarField cellDx(max(Foam::sqrt(magSf), 2/patch().deltaCoeffs()));
206
207 // Inialise eddy box extents
208 forAll(*this, faceI)
209 {
210 scalar& s = sigmax_[faceI];
211
212 // Average length scale (SST:Eq. 24)
213 // Personal communication regarding (PCR:Eq. 14)
214 // - the min operator in Eq. 14 is a typo, and should be a max operator
215 s = min(mag(L[faceI]), kappa_*delta_);
216 s = max(s, nCellPerEddy_*cellDx[faceI]);
217 }
218
219 // Maximum extent across all processors
220 maxSigmaX_ = gMax(sigmax_);
221
222 // Eddy box volume
223 v0_ = 2*gSum(magSf)*maxSigmaX_;
224
225 if (debug)
226 {
227 Info<< "Patch: " << patch().patch().name() << " eddy box:" << nl
228 << " volume : " << v0_ << nl
229 << " maxSigmaX : " << maxSigmaX_ << nl
230 << endl;
231 }
232}
233
234
235Foam::pointIndexHit Foam::turbulentDFSEMInletFvPatchVectorField::setNewPosition
236(
237 const bool global
238)
239{
240 // Initialise to miss
241 pointIndexHit pos(false, vector::max, -1);
242
243 const polyPatch& patch = this->patch().patch();
244 const pointField& points = patch.points();
245
246 if (global)
247 {
248 const scalar areaFraction =
249 rndGen_.globalPosition<scalar>(0, patchArea_);
250
251 // Determine which processor to use
252 label procI = 0;
253 forAllReverse(sumTriMagSf_, i)
254 {
255 if (areaFraction >= sumTriMagSf_[i])
256 {
257 procI = i;
258 break;
259 }
260 }
261
262 if (Pstream::myProcNo() == procI)
263 {
264 // Find corresponding decomposed face triangle
265 label triI = 0;
266 const scalar offset = sumTriMagSf_[procI];
267 forAllReverse(triCumulativeMagSf_, i)
268 {
269 if (areaFraction > triCumulativeMagSf_[i] + offset)
270 {
271 triI = i;
272 break;
273 }
274 }
275
276 // Find random point in triangle
277 const face& tf = triFace_[triI];
278 const triPointRef tri(points[tf[0]], points[tf[1]], points[tf[2]]);
279
280 pos.setHit();
281 pos.setIndex(triToFace_[triI]);
282 pos.rawPoint() = tri.randomPoint(rndGen_);
283 }
284 }
285 else
286 {
287 // Find corresponding decomposed face triangle on local processor
288 label triI = 0;
289 const scalar maxAreaLimit = triCumulativeMagSf_.last();
290 const scalar areaFraction = rndGen_.position<scalar>(0, maxAreaLimit);
291
292 forAllReverse(triCumulativeMagSf_, i)
293 {
294 if (areaFraction > triCumulativeMagSf_[i])
295 {
296 triI = i;
297 break;
298 }
299 }
300
301 // Find random point in triangle
302 const face& tf = triFace_[triI];
303 const triPointRef tri(points[tf[0]], points[tf[1]], points[tf[2]]);
304
305 pos.setHit();
306 pos.setIndex(triToFace_[triI]);
307 pos.rawPoint() = tri.randomPoint(rndGen_);
308 }
309
310 return pos;
311}
312
313
314void Foam::turbulentDFSEMInletFvPatchVectorField::initialiseEddies()
315{
316 const scalar t = db().time().timeOutputValue();
317 const symmTensorField R(R_->value(t)/sqr(Uref_));
318
319 DynamicList<eddy> eddies(size());
320
321 // Initialise eddy properties
322 scalar sumVolEddy = 0;
323 scalar sumVolEddyAllProc = 0;
324
325 while (sumVolEddyAllProc/v0_ < d_)
326 {
327 bool search = true;
328 label iter = 0;
329
330 while (search && iter++ < seedIterMax_)
331 {
332 // Get new parallel consistent position
333 pointIndexHit pos(setNewPosition(true));
334 const label patchFaceI = pos.index();
335
336 // Note: only 1 processor will pick up this face
337 if (patchFaceI != -1)
338 {
339 eddy e
340 (
341 patchFaceI,
342 pos.hitPoint(),
343 rndGen_.position<scalar>(-maxSigmaX_, maxSigmaX_),
344 sigmax_[patchFaceI],
345 R[patchFaceI],
346 rndGen_
347 );
348
349 // If eddy valid, patchFaceI is non-zero
350 if (e.patchFaceI() != -1)
351 {
352 eddies.append(e);
353 sumVolEddy += e.volume();
354 search = false;
355 }
356 }
357 // else eddy on remote processor
358
359 reduce(search, andOp<bool>());
360 }
361
362
363 sumVolEddyAllProc = returnReduce(sumVolEddy, sumOp<scalar>());
364 }
365 eddies_.transfer(eddies);
366
367 nEddy_ = eddies_.size();
368
369 if (debug)
370 {
371 Pout<< "Patch:" << patch().patch().name();
372
373 if (Pstream::parRun())
374 {
375 Pout<< " processor:" << Pstream::myProcNo();
376 }
377
378 Pout<< " seeded:" << nEddy_ << " eddies" << endl;
379 }
380
381 reduce(nEddy_, sumOp<label>());
382
383 if (nEddy_ > 0)
384 {
385 Info<< "Turbulent DFSEM patch: " << patch().name()
386 << " seeded " << nEddy_ << " eddies with total volume "
387 << sumVolEddyAllProc
388 << endl;
389 }
390 else
391 {
393 << "Patch: " << patch().patch().name()
394 << " on field " << internalField().name()
395 << ": No eddies seeded - please check your set-up"
396 << endl;
397 }
398}
399
400
401void Foam::turbulentDFSEMInletFvPatchVectorField::convectEddies
402(
403 const vector& UBulk,
404 const scalar deltaT
405)
406{
407 const scalar t = db().time().timeOutputValue();
408 const symmTensorField R(R_->value(t)/sqr(Uref_));
409
410 // Note: all operations applied to local processor only
411
412 label nRecycled = 0;
413
414 forAll(eddies_, eddyI)
415 {
416 eddy& e = eddies_[eddyI];
417 e.move(deltaT*(UBulk & patchNormal_));
418
419 const scalar position0 = e.x();
420
421 // Check to see if eddy has exited downstream box plane
422 if (position0 > maxSigmaX_)
423 {
424 bool search = true;
425 label iter = 0;
426
427 while (search && iter++ < seedIterMax_)
428 {
429 // Spawn new eddy with new random properties (intensity etc)
430 pointIndexHit pos(setNewPosition(false));
431 const label patchFaceI = pos.index();
432
433 e = eddy
434 (
435 patchFaceI,
436 pos.hitPoint(),
437 position0 - floor(position0/maxSigmaX_)*maxSigmaX_,
438 sigmax_[patchFaceI],
439 R[patchFaceI],
440 rndGen_
441 );
442
443 if (e.patchFaceI() != -1)
444 {
445 search = false;
446 }
447 }
448
449 nRecycled++;
450 }
451 }
452
453 reduce(nRecycled, sumOp<label>());
454
455 if (debug && nRecycled > 0)
456 {
457 Info<< "Patch: " << patch().patch().name()
458 << " recycled " << nRecycled << " eddies"
459 << endl;
460 }
461}
462
463
464Foam::vector Foam::turbulentDFSEMInletFvPatchVectorField::uPrimeEddy
465(
466 const List<eddy>& eddies,
467 const point& patchFaceCf
468) const
469{
470 vector uPrime(Zero);
471
472 forAll(eddies, k)
473 {
474 const eddy& e = eddies[k];
475 uPrime += e.uPrime(patchFaceCf, patchNormal_);
476 }
477
478 return uPrime;
479}
480
481
482void Foam::turbulentDFSEMInletFvPatchVectorField::calcOverlappingProcEddies
483(
484 List<List<eddy>>& overlappingEddies
485) const
486{
487 int oldTag = UPstream::msgType();
488 UPstream::msgType() = oldTag + 1;
489
490 List<boundBox> patchBBs(Pstream::nProcs());
491 patchBBs[Pstream::myProcNo()] = patchBounds_;
492 Pstream::allGatherList(patchBBs);
493
494 // Per processor indices into all segments to send
495 List<DynamicList<label>> dynSendMap(Pstream::nProcs());
496
497 forAll(eddies_, i)
498 {
499 // Collect overlapping eddies
500 const eddy& e = eddies_[i];
501
502 // Eddy bounds
503 const point x(e.position(patchNormal_));
504 boundBox ebb(e.bounds());
505 ebb.min() += x;
506 ebb.max() += x;
507
508 forAll(patchBBs, procI)
509 {
510 // Not including intersection with local patch
511 if (procI != Pstream::myProcNo())
512 {
513 if (ebb.overlaps(patchBBs[procI]))
514 {
515 dynSendMap[procI].append(i);
516 }
517 }
518 }
519 }
520
522 forAll(sendMap, procI)
523 {
524 sendMap[procI].transfer(dynSendMap[procI]);
525 }
526
527 // Send the number of eddies for local processors to receive
528 labelListList sendSizes(Pstream::nProcs());
529 sendSizes[Pstream::myProcNo()].setSize(Pstream::nProcs());
530 forAll(sendMap, procI)
531 {
532 sendSizes[Pstream::myProcNo()][procI] = sendMap[procI].size();
533 }
534 Pstream::allGatherList(sendSizes);
535
536 // Determine order of receiving
537 labelListList constructMap(Pstream::nProcs());
538
539 // Local segment first
540 constructMap[Pstream::myProcNo()] = identity
541 (
542 sendMap[Pstream::myProcNo()].size()
543 );
544
545 label segmentI = constructMap[Pstream::myProcNo()].size();
546 forAll(constructMap, procI)
547 {
548 if (procI != Pstream::myProcNo())
549 {
550 // What I need to receive is what other processor is sending to me
551 const label nRecv = sendSizes[procI][Pstream::myProcNo()];
552 constructMap[procI].setSize(nRecv);
553
554 for (label i = 0; i < nRecv; ++i)
555 {
556 constructMap[procI][i] = segmentI++;
557 }
558 }
559 }
560
561 mapDistribute map(segmentI, std::move(sendMap), std::move(constructMap));
562
563 PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
564
565 for (const int domain : Pstream::allProcs())
566 {
567 const labelList& sendElems = map.subMap()[domain];
568
569 if (domain != Pstream::myProcNo() && sendElems.size())
570 {
571 List<eddy> subEddies(UIndirectList<eddy>(eddies_, sendElems));
572
573 UOPstream toDomain(domain, pBufs);
574
575 toDomain<< subEddies;
576 }
577 }
578
579 // Start receiving
580 pBufs.finishedSends();
581
582 // Consume
583 for (const int domain : Pstream::allProcs())
584 {
585 const labelList& recvElems = map.constructMap()[domain];
586
587 if (domain != Pstream::myProcNo() && recvElems.size())
588 {
589 UIPstream str(domain, pBufs);
590 {
591 str >> overlappingEddies[domain];
592 }
593 }
594 }
595
596 // Restore tag
597 UPstream::msgType() = oldTag;
598}
599
600
601// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
602
603Foam::turbulentDFSEMInletFvPatchVectorField::
604turbulentDFSEMInletFvPatchVectorField
605(
606 const fvPatch& p,
608)
609:
611 U_(nullptr),
612 R_(nullptr),
613 L_(nullptr),
614 delta_(1.0),
615 d_(1.0),
616 kappa_(0.41),
617 Uref_(1.0),
618 Lref_(1.0),
619 scale_(1.0),
620 m_(0.5),
621 nCellPerEddy_(5),
622
623 patchArea_(-1),
624 triFace_(),
625 triToFace_(),
626 triCumulativeMagSf_(),
627 sumTriMagSf_(Pstream::nProcs() + 1, Zero),
628 patchNormal_(Zero),
629 patchBounds_(boundBox::invertedBox),
630
631 eddies_(Zero),
632 v0_(Zero),
633 rndGen_(Pstream::myProcNo()),
634 sigmax_(size(), Zero),
635 maxSigmaX_(Zero),
636 nEddy_(0),
637 curTimeIndex_(-1),
638 singleProc_(false),
639 writeEddies_(false)
640{}
641
642
643Foam::turbulentDFSEMInletFvPatchVectorField::
644turbulentDFSEMInletFvPatchVectorField
645(
647 const fvPatch& p,
649 const fvPatchFieldMapper& mapper
650)
651:
652 fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
653 U_(ptf.U_.clone(patch().patch())),
654 R_(ptf.R_.clone(patch().patch())),
655 L_(ptf.L_.clone(patch().patch())),
656 delta_(ptf.delta_),
657 d_(ptf.d_),
658 kappa_(ptf.kappa_),
659 Uref_(ptf.Uref_),
660 Lref_(ptf.Lref_),
661 scale_(ptf.scale_),
662 m_(ptf.m_),
663 nCellPerEddy_(ptf.nCellPerEddy_),
664
665 patchArea_(ptf.patchArea_),
666 triFace_(ptf.triFace_),
667 triToFace_(ptf.triToFace_),
668 triCumulativeMagSf_(ptf.triCumulativeMagSf_),
669 sumTriMagSf_(ptf.sumTriMagSf_),
670 patchNormal_(ptf.patchNormal_),
671 patchBounds_(ptf.patchBounds_),
672
673 eddies_(ptf.eddies_),
674 v0_(ptf.v0_),
675 rndGen_(ptf.rndGen_),
676 sigmax_(ptf.sigmax_, mapper),
677 maxSigmaX_(ptf.maxSigmaX_),
678 nEddy_(ptf.nEddy_),
679 curTimeIndex_(-1),
680 singleProc_(ptf.singleProc_),
681 writeEddies_(ptf.writeEddies_)
682{}
683
684
685Foam::turbulentDFSEMInletFvPatchVectorField::
686turbulentDFSEMInletFvPatchVectorField
687(
688 const fvPatch& p,
690 const dictionary& dict
691)
692:
694 U_(PatchFunction1<vector>::New(patch().patch(), "U", dict)),
695 R_(PatchFunction1<symmTensor>::New(patch().patch(), "R", dict)),
696 L_(PatchFunction1<scalar>::New(patch().patch(), "L", dict)),
697 delta_(dict.getCheck<scalar>("delta", scalarMinMax::ge(0))),
698 d_(dict.getCheckOrDefault<scalar>("d", 1, scalarMinMax::ge(SMALL))),
699 kappa_(dict.getCheckOrDefault<scalar>("kappa", 0.41, scalarMinMax::ge(0))),
700 Uref_(dict.getCheckOrDefault<scalar>("Uref", 1, scalarMinMax::ge(SMALL))),
701 Lref_(dict.getCheckOrDefault<scalar>("Lref", 1, scalarMinMax::ge(SMALL))),
702 scale_(dict.getCheckOrDefault<scalar>("scale", 1, scalarMinMax::ge(0))),
703 m_(dict.getCheckOrDefault<scalar>("m", 0.5, scalarMinMax::ge(0))),
704 nCellPerEddy_(dict.getOrDefault<label>("nCellPerEddy", 5)),
705
706 patchArea_(-1),
707 triFace_(),
708 triToFace_(),
709 triCumulativeMagSf_(),
710 sumTriMagSf_(Pstream::nProcs() + 1, Zero),
711 patchNormal_(Zero),
712 patchBounds_(boundBox::invertedBox),
713
714 eddies_(),
715 v0_(Zero),
716 rndGen_(),
717 sigmax_(size(), Zero),
718 maxSigmaX_(Zero),
719 nEddy_(0),
720 curTimeIndex_(-1),
721 singleProc_(false),
722 writeEddies_(dict.getOrDefault("writeEddies", false))
723{
724 eddy::debug = debug;
725
726 const scalar t = db().time().timeOutputValue();
727 const symmTensorField R(R_->value(t)/sqr(Uref_));
728
730}
731
732
733Foam::turbulentDFSEMInletFvPatchVectorField::
734turbulentDFSEMInletFvPatchVectorField
735(
737)
738:
740 U_(ptf.U_.clone(patch().patch())),
741 R_(ptf.R_.clone(patch().patch())),
742 L_(ptf.L_.clone(patch().patch())),
743 delta_(ptf.delta_),
744 d_(ptf.d_),
745 kappa_(ptf.kappa_),
746 Uref_(ptf.Uref_),
747 Lref_(ptf.Lref_),
748 scale_(ptf.scale_),
749 m_(ptf.m_),
750 nCellPerEddy_(ptf.nCellPerEddy_),
751
752 patchArea_(ptf.patchArea_),
753 triFace_(ptf.triFace_),
754 triToFace_(ptf.triToFace_),
755 triCumulativeMagSf_(ptf.triCumulativeMagSf_),
756 sumTriMagSf_(ptf.sumTriMagSf_),
757 patchNormal_(ptf.patchNormal_),
758 patchBounds_(ptf.patchBounds_),
759
760 eddies_(ptf.eddies_),
761 v0_(ptf.v0_),
762 rndGen_(ptf.rndGen_),
763 sigmax_(ptf.sigmax_),
764 maxSigmaX_(ptf.maxSigmaX_),
765 nEddy_(ptf.nEddy_),
766 curTimeIndex_(-1),
767 singleProc_(ptf.singleProc_),
768 writeEddies_(ptf.writeEddies_)
769{}
770
771
772Foam::turbulentDFSEMInletFvPatchVectorField::
773turbulentDFSEMInletFvPatchVectorField
774(
777)
778:
780 U_(ptf.U_.clone(patch().patch())),
781 R_(ptf.R_.clone(patch().patch())),
782 L_(ptf.L_.clone(patch().patch())),
783 delta_(ptf.delta_),
784 d_(ptf.d_),
785 kappa_(ptf.kappa_),
786 Uref_(ptf.Uref_),
787 Lref_(ptf.Lref_),
788 scale_(ptf.scale_),
789 m_(ptf.m_),
790 nCellPerEddy_(ptf.nCellPerEddy_),
791
792 patchArea_(ptf.patchArea_),
793 triFace_(ptf.triFace_),
794 triToFace_(ptf.triToFace_),
795 triCumulativeMagSf_(ptf.triCumulativeMagSf_),
796 sumTriMagSf_(ptf.sumTriMagSf_),
797 patchNormal_(ptf.patchNormal_),
798 patchBounds_(ptf.patchBounds_),
799
800 eddies_(ptf.eddies_),
801 v0_(ptf.v0_),
802 rndGen_(ptf.rndGen_),
803 sigmax_(ptf.sigmax_),
804 maxSigmaX_(ptf.maxSigmaX_),
805 nEddy_(ptf.nEddy_),
806 curTimeIndex_(-1),
807 singleProc_(ptf.singleProc_),
808 writeEddies_(ptf.writeEddies_)
809{}
810
811
812// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
813
815(
816 const symmTensorField& R
817)
818{
819 constexpr label maxDiffs = 5;
820 label nDiffs = 0;
821
822 // (S:Eq. 4a-4c)
823 forAll(R, i)
824 {
825 bool diff = false;
826
827 if (maxDiffs < nDiffs)
828 {
829 Info<< "More than " << maxDiffs << " times"
830 << " Reynolds-stress realizability checks failed."
831 << " Skipping further comparisons." << endl;
832 return;
833 }
834
835 const symmTensor& r = R[i];
836
837 if (r.xx() < 0)
838 {
840 << "Reynolds stress " << r << " at index " << i
841 << " does not obey the constraint: Rxx >= 0"
842 << endl;
843 diff = true;
844 }
845
846 if ((r.xx()*r.yy() - sqr(r.xy())) < 0)
847 {
849 << "Reynolds stress " << r << " at index " << i
850 << " does not obey the constraint: Rxx*Ryy - sqr(Rxy) >= 0"
851 << endl;
852 diff = true;
853 }
854
855 if (det(r) < 0)
856 {
858 << "Reynolds stress " << r << " at index " << i
859 << " does not obey the constraint: det(R) >= 0"
860 << endl;
861 diff = true;
862 }
863
864 if (diff)
865 {
866 ++nDiffs;
867 }
868 }
869}
870
871
873(
874 const scalarField& R
875)
876{
877 if (min(R) <= 0)
878 {
880 << "Reynolds stresses contain at least one "
881 << "nonpositive element. min(R) = " << min(R)
882 << exit(FatalError);
883 }
884}
885
886
888(
889 const fvPatchFieldMapper& m
890)
891{
893
894 if (U_)
895 {
896 U_->autoMap(m);
897 }
898 if (R_)
899 {
900 R_->autoMap(m);
901 }
902 if (L_)
903 {
904 L_->autoMap(m);
905 }
906
907 sigmax_.autoMap(m);
908}
909
910
912(
913 const fvPatchVectorField& ptf,
914 const labelList& addr
915)
916{
918
919 const auto& dfsemptf =
920 refCast<const turbulentDFSEMInletFvPatchVectorField>(ptf);
921
922 if (U_)
923 {
924 U_->rmap(dfsemptf.U_(), addr);
925 }
926 if (R_)
927 {
928 R_->rmap(dfsemptf.R_(), addr);
929 }
930 if (L_)
931 {
932 L_->rmap(dfsemptf.L_(), addr);
933 }
934
935 sigmax_.rmap(dfsemptf.sigmax_, addr);
936}
937
938
940{
941 if (updated())
942 {
943 return;
944 }
945
946 if (curTimeIndex_ == -1)
947 {
948 initialisePatch();
949
950 initialiseEddyBox();
951
952 initialiseEddies();
953 }
954
955
956 if (curTimeIndex_ != db().time().timeIndex())
957 {
959 U_->value(db().time().timeOutputValue())/Uref_;
960
961 // (PCR:p. 522)
962 const vector UBulk
963 (
964 gSum(UMean()*patch().magSf())
965 /(gSum(patch().magSf()) + ROOTVSMALL)
966 );
967
968 // Move eddies using bulk velocity
969 const scalar deltaT = db().time().deltaTValue();
970 convectEddies(UBulk, deltaT);
971
972 // Set mean velocity
973 vectorField& U = *this;
974 U = UMean;
975
976 // Apply second part of normalisation coefficient
977 const scalar c =
978 scale_*Foam::pow(10*v0_, m_)/Foam::sqrt(scalar(nEddy_));
979
980 // In parallel, need to collect all eddies that will interact with
981 // local faces
982
983 const pointField& Cf = patch().Cf();
984
985 if (singleProc_ || !Pstream::parRun())
986 {
987 forAll(U, faceI)
988 {
989 U[faceI] += c*uPrimeEddy(eddies_, Cf[faceI]);
990 }
991 }
992 else
993 {
994 // Process local eddy contributions
995 forAll(U, faceI)
996 {
997 U[faceI] += c*uPrimeEddy(eddies_, Cf[faceI]);
998 }
999
1000 // Add contributions from overlapping eddies
1001 List<List<eddy>> overlappingEddies(Pstream::nProcs());
1002 calcOverlappingProcEddies(overlappingEddies);
1003
1004 forAll(overlappingEddies, procI)
1005 {
1006 const List<eddy>& eddies = overlappingEddies[procI];
1007
1008 if (eddies.size())
1009 {
1010 forAll(U, faceI)
1011 {
1012 U[faceI] += c*uPrimeEddy(eddies, Cf[faceI]);
1013 }
1014 }
1015 }
1016 }
1017
1018 // Re-scale to ensure correct flow rate
1019 const scalar fCorr =
1020 gSum((UBulk & patchNormal_)*patch().magSf())
1021 /gSum(U & -patch().Sf());
1022
1023 U *= fCorr;
1024
1025 curTimeIndex_ = db().time().timeIndex();
1026
1027 if (writeEddies_)
1028 {
1029 writeEddyOBJ();
1030 }
1031
1032 if (debug)
1033 {
1034 Info<< "Magnitude of bulk velocity: " << UBulk << endl;
1035
1036 label n = eddies_.size();
1037 Info<< "Number of eddies: " << returnReduce(n, sumOp<label>())
1038 << endl;
1039
1040 Info<< "Patch:" << patch().patch().name()
1041 << " min/max(U):" << gMin(U) << ", " << gMax(U)
1042 << endl;
1043
1044 if (db().time().writeTime())
1045 {
1046 writeLumleyCoeffs();
1047 }
1048 }
1049 }
1050
1051 fixedValueFvPatchVectorField::updateCoeffs();
1052}
1053
1054
1056{
1058 os.writeEntry("delta", delta_);
1059 os.writeEntryIfDifferent<scalar>("d", 1.0, d_);
1060 os.writeEntryIfDifferent<scalar>("kappa", 0.41, kappa_);
1061 os.writeEntryIfDifferent<scalar>("Uref", 1.0, Uref_);
1062 os.writeEntryIfDifferent<scalar>("Lref", 1.0, Lref_);
1063 os.writeEntryIfDifferent<scalar>("scale", 1.0, scale_);
1064 os.writeEntryIfDifferent<scalar>("m", 0.5, m_);
1065 os.writeEntryIfDifferent<label>("nCellPerEddy", 5, nCellPerEddy_);
1066 os.writeEntryIfDifferent("writeEddies", false, writeEddies_);
1067 if (U_)
1068 {
1069 U_->writeData(os);
1070 }
1071 if (R_)
1072 {
1073 R_->writeData(os);
1074 }
1075 if (L_)
1076 {
1077 L_->writeData(os);
1078 }
1079 writeEntry("value", os);
1080}
1081
1082
1083// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1084
1085namespace Foam
1086{
1088 (
1091 );
1092}
1093
1094
1095// ************************************************************************* //
label k
#define R(A, B, C, D, E, F, K, M)
label n
Macros for easy insertion into run-time selection tables.
volVectorField UMean(UMeanHeader, mesh)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
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
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:251
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:66
UPstream::rangeType allProcs() const noexcept
Range of ranks indices associated with PstreamBuffers.
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
Inter-processor communications stream.
Definition: Pstream.H:63
static void allGatherList(const List< commsStruct > &comms, List< T > &values, const int tag, const label comm)
static void listCombineAllGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
After completion all processors have the same data.
const Cmpt & xx() const
Definition: SymmTensorI.H:97
const Cmpt & xy() const
Definition: SymmTensorI.H:103
const Cmpt & yy() const
Definition: SymmTensorI.H:121
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
@ nonBlocking
"nonBlocking"
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:556
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
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 bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
static int debug
Flag to activate debug statements.
Definition: eddy.H:147
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.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
static const complex max
complex (VGREAT,VGREAT)
Definition: complex.H:293
int myProcNo() const noexcept
Return processor number.
A class for managing temporary objects.
Definition: tmp.H:65
@ TAB
Tab [isspace].
Definition: token.H:123
The turbulentDFSEMInlet is a synthesised-eddy based velocity inlet boundary condition to generate syn...
static void checkStresses(const symmTensorField &R)
Check if input Reynolds stresses are valid.
virtual void rmap(const fvPatchVectorField &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...
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
U
Definition: pEqn.H:72
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
const pointField & points
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define WarningInFunction
Report a warning using Foam::Warning.
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar det(const dimensionedSphericalTensor &dt)
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Type gSum(const FieldField< Field, Type > &f)
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i)
Definition: labelList.C:38
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
List< label > labelList
A List of labels.
Definition: List.H:66
dimensionedSymmTensor sqr(const dimensionedVector &dv)
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
Cmpt invariantII(const SymmTensor< Cmpt > &st)
Return the 2nd invariant of a SymmTensor.
Definition: SymmTensorI.H:468
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition: symmTensor.H:59
messageStream Info
Information stream (stdout output on master, null elsewhere)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
vector point
Point is a vector.
Definition: point.H:43
Cmpt invariantIII(const SymmTensor< Cmpt > &st)
Return the 3rd invariant of a SymmTensor.
Definition: SymmTensorI.H:480
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangle.H:71
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:378
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
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.
error FatalError
dimensionedScalar cbrt(const dimensionedScalar &ds)
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Type gMax(const FieldField< Field, Type > &f)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
fileName search(const word &file, const fileName &directory)
Recursively search the given directory for the file.
Definition: fileName.C:624
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
label timeIndex
Definition: getTimeIndex.H:30
face triFace(3)
labelList f(nPoints)
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:346
const vector L(dict.get< vector >("L"))