Loading...
Searching...
No Matches
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-2024 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 const polyPatch& patch = this->patch().patch();
136 const pointField& points = patch.points();
137
138 // Triangulate the patch faces and create addressing
139 {
140 label nTris = 0;
141 for (const face& f : patch)
142 {
143 nTris += f.nTriangles();
144 }
145
146 DynamicList<labelledTri> dynTriFace(nTris);
147 DynamicList<face> tris(8); // work array
148
149 forAll(patch, facei)
150 {
151 const face& f = patch[facei];
152
153 tris.clear();
154 f.triangles(points, tris);
155
156 for (const auto& t : tris)
157 {
158 dynTriFace.emplace_back(t[0], t[1], t[2], facei);
159 }
160 }
161
162 // Transfer to persistent storage
163 triFace_.transfer(dynTriFace);
164 }
165
166
167 const label myProci = UPstream::myProcNo();
168 const label numProc = UPstream::nProcs();
169
170 // Calculate the cumulative triangle weights
171 {
172 triCumulativeMagSf_.resize_nocopy(triFace_.size()+1);
173
174 auto iter = triCumulativeMagSf_.begin();
175
176 // Set zero value at the start of the tri area/weight list
177 scalar patchArea = 0;
178 *iter++ = patchArea;
179
180 // Calculate cumulative and total area
181 for (const auto& t : triFace_)
182 {
183 patchArea += t.mag(points);
184 *iter++ = patchArea;
185 }
186
187 sumTriMagSf_.resize_nocopy(numProc+1);
188 sumTriMagSf_[0] = 0;
189
190 {
191 scalarList::subList slice(sumTriMagSf_, numProc, 1);
192 slice[myProci] = patchArea;
194 }
195
196 // Convert to cumulative sum of areas per proc
197 for (label i = 1; i < sumTriMagSf_.size(); ++i)
198 {
199 sumTriMagSf_[i] += sumTriMagSf_[i-1];
200 }
201 }
202
203 // Global patch area (over all processors)
204 patchArea_ = sumTriMagSf_.back();
205
206 // Local patch bounds (this processor)
207 patchBounds_ = boundBox(patch.localPoints(), false);
208 patchBounds_.inflate(0.1);
209
210 // Determine if all eddies spawned from a single processor
211 singleProc_ = returnReduceOr
212 (
213 patch.size() == returnReduce(patch.size(), sumOp<label>())
214 );
215}
216
217
218void Foam::turbulentDFSEMInletFvPatchVectorField::initialiseEddyBox()
219{
220 const scalarField& magSf = patch().magSf();
221
222 const scalarField L(L_->value(db().time().timeOutputValue())/Lref_);
223
224 // (PCF:Eq. 14)
225 const scalarField cellDx(max(Foam::sqrt(magSf), 2/patch().deltaCoeffs()));
226
227 // Inialise eddy box extents
228 forAll(*this, faceI)
229 {
230 scalar& s = sigmax_[faceI];
231
232 // Average length scale (SST:Eq. 24)
233 // Personal communication regarding (PCR:Eq. 14)
234 // - the min operator in Eq. 14 is a typo, and should be a max operator
235 s = min(mag(L[faceI]), kappa_*delta_);
236 s = max(s, nCellPerEddy_*cellDx[faceI]);
237 }
238
239 // Maximum extent across all processors
240 maxSigmaX_ = gMax(sigmax_);
241
242 // Eddy box volume
243 v0_ = 2*gSum(magSf)*maxSigmaX_;
244
245 if (debug)
246 {
247 Info<< "Patch: " << patch().patch().name() << " eddy box:" << nl
248 << " volume : " << v0_ << nl
249 << " maxSigmaX : " << maxSigmaX_ << nl
250 << endl;
251 }
252}
253
254
255Foam::pointIndexHit Foam::turbulentDFSEMInletFvPatchVectorField::setNewPosition
256(
257 const bool global
258)
259{
260 const polyPatch& patch = this->patch().patch();
261 const pointField& points = patch.points();
262
263 label triI = -1;
264
265 if (global)
266 {
267 const scalar areaFraction =
268 rndGen_.globalPosition<scalar>(0, patchArea_);
269
270 // Determine which processor to use
271 label proci = 0;
272 forAllReverse(sumTriMagSf_, i)
273 {
274 if (areaFraction >= sumTriMagSf_[i])
275 {
276 proci = i;
277 break;
278 }
279 }
280
281 if (UPstream::myProcNo() == proci)
282 {
283 // Find corresponding decomposed face triangle
284 triI = 0;
285 const scalar offset = sumTriMagSf_[proci];
286 forAllReverse(triCumulativeMagSf_, i)
287 {
288 if (areaFraction > triCumulativeMagSf_[i] + offset)
289 {
290 triI = i;
291 break;
292 }
293 }
294 }
295 }
296 else
297 {
298 // Find corresponding decomposed face triangle on local processor
299 triI = 0;
300 const scalar maxAreaLimit = triCumulativeMagSf_.back();
301 const scalar areaFraction = rndGen_.position<scalar>(0, maxAreaLimit);
302
303 forAllReverse(triCumulativeMagSf_, i)
304 {
305 if (areaFraction > triCumulativeMagSf_[i])
306 {
307 triI = i;
308 break;
309 }
310 }
311 }
312
313
314 if (triI >= 0)
315 {
316 return pointIndexHit
317 (
318 true,
319 // Find random point in triangle
320 triFace_[triI].tri(points).randomPoint(rndGen_),
321 triFace_[triI].index()
322 );
323 }
324
325 // No hit
326 return pointIndexHit(false, vector::max, -1);
327}
328
329
330void Foam::turbulentDFSEMInletFvPatchVectorField::initialiseEddies()
331{
332 const scalar t = db().time().timeOutputValue();
333 const symmTensorField R(R_->value(t)/sqr(Uref_));
334
335 DynamicList<eddy> eddies(size());
336
337 // Initialise eddy properties
338 scalar sumVolEddy = 0;
339 scalar sumVolEddyAllProc = 0;
340
341 while (sumVolEddyAllProc/v0_ < d_)
342 {
343 bool search = true;
344 label iter = 0;
345
346 while (search && iter++ < seedIterMax_)
347 {
348 // Get new parallel consistent position
349 pointIndexHit pos(setNewPosition(true));
350 const label patchFaceI = pos.index();
351
352 // Note: only 1 processor will pick up this face
353 if (patchFaceI != -1)
354 {
355 eddy e
356 (
357 patchFaceI,
358 pos.hitPoint(),
359 rndGen_.position<scalar>(-maxSigmaX_, maxSigmaX_),
360 sigmax_[patchFaceI],
361 R[patchFaceI],
362 rndGen_
363 );
364
365 // If eddy valid, patchFaceI is non-zero
366 if (e.patchFaceI() != -1)
367 {
368 eddies.append(e);
369 sumVolEddy += e.volume();
370 search = false;
371 }
372 }
373 // else eddy on remote processor
374
376 }
377
378
379 sumVolEddyAllProc = returnReduce(sumVolEddy, sumOp<scalar>());
380 }
381 eddies_.transfer(eddies);
382
383 nEddy_ = eddies_.size();
384
385 if (debug)
386 {
387 Pout<< "Patch:" << patch().patch().name();
388
389 if (Pstream::parRun())
390 {
391 Pout<< " processor:" << Pstream::myProcNo();
392 }
393
394 Pout<< " seeded:" << nEddy_ << " eddies" << endl;
395 }
396
397 reduce(nEddy_, sumOp<label>());
398
399 if (nEddy_ > 0)
400 {
401 Info<< "Turbulent DFSEM patch: " << patch().name()
402 << " seeded " << nEddy_ << " eddies with total volume "
403 << sumVolEddyAllProc
404 << endl;
405 }
406 else
407 {
409 << "Patch: " << patch().patch().name()
410 << " on field " << internalField().name()
411 << ": No eddies seeded - please check your set-up"
412 << endl;
413 }
414}
415
416
417void Foam::turbulentDFSEMInletFvPatchVectorField::convectEddies
418(
419 const vector& UBulk,
420 const scalar deltaT
421)
422{
423 const scalar t = db().time().timeOutputValue();
424 const symmTensorField R(R_->value(t)/sqr(Uref_));
425
426 // Note: all operations applied to local processor only
427
428 label nRecycled = 0;
429
430 forAll(eddies_, eddyI)
431 {
432 eddy& e = eddies_[eddyI];
433 e.move(deltaT*(UBulk & patchNormal_));
434
435 const scalar position0 = e.x();
436
437 // Check to see if eddy has exited downstream box plane
438 if (position0 > maxSigmaX_)
439 {
440 bool search = true;
441 label iter = 0;
442
443 while (search && iter++ < seedIterMax_)
444 {
445 // Spawn new eddy with new random properties (intensity etc)
446 pointIndexHit pos(setNewPosition(false));
447 const label patchFaceI = pos.index();
448
449 e = eddy
450 (
451 patchFaceI,
452 pos.hitPoint(),
453 position0 - floor(position0/maxSigmaX_)*maxSigmaX_,
454 sigmax_[patchFaceI],
455 R[patchFaceI],
456 rndGen_
457 );
458
459 if (e.patchFaceI() != -1)
460 {
461 search = false;
462 }
463 }
464
465 nRecycled++;
466 }
467 }
468
469 if (debug)
470 {
471 reduce(nRecycled, sumOp<label>());
472
473 if (nRecycled)
474 {
475 Info<< "Patch: " << patch().patch().name()
476 << " recycled " << nRecycled << " eddies"
477 << endl;
478 }
479 }
480}
481
482
483Foam::vector Foam::turbulentDFSEMInletFvPatchVectorField::uPrimeEddy
484(
485 const List<eddy>& eddies,
486 const point& patchFaceCf
487) const
488{
489 vector uPrime(Zero);
490
491 forAll(eddies, k)
492 {
493 const eddy& e = eddies[k];
494 uPrime += e.uPrime(patchFaceCf, patchNormal_);
495 }
496
497 return uPrime;
498}
499
500
501void Foam::turbulentDFSEMInletFvPatchVectorField::calcOverlappingProcEddies
502(
503 List<List<eddy>>& overlappingEddies
504) const
505{
506 const int oldTag = UPstream::incrMsgType();
507
509 patchBBs[Pstream::myProcNo()] = patchBounds_;
510 Pstream::allGatherList(patchBBs);
511
512 // Per processor indices into all segments to send
514
515 // Collect overlapping eddies
516 forAll(eddies_, i)
517 {
518 const eddy& e = eddies_[i];
519
520 // Eddy bounds
521 const point x(e.position(patchNormal_));
522 boundBox ebb(e.bounds());
523 ebb.min() += x;
524 ebb.max() += x;
525
526 forAll(patchBBs, proci)
527 {
528 // Not including intersection with local patch
529 if (proci != Pstream::myProcNo())
530 {
531 if (ebb.overlaps(patchBBs[proci]))
532 {
533 sendMap[proci].push_back(i);
534 }
535 }
536 }
537 }
538
539
540 PstreamBuffers pBufs;
541
542 forAll(sendMap, proci)
543 {
544 if (proci != Pstream::myProcNo() && !sendMap[proci].empty())
545 {
546 UOPstream os(proci, pBufs);
547
548 os << UIndirectList<eddy>(eddies_, sendMap[proci]);
549 }
550 }
551
552 pBufs.finishedSends();
553
554 for (const int proci : pBufs.allProcs())
555 {
556 if (proci != Pstream::myProcNo() && pBufs.recvDataCount(proci))
557 {
558 UIPstream is(proci, pBufs);
559
560 is >> overlappingEddies[proci];
561 }
562 }
564 UPstream::msgType(oldTag); // Restore tag
565}
566
567
568// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
569
572(
573 const fvPatch& p,
575)
576:
578 U_(nullptr),
579 R_(nullptr),
580 L_(nullptr),
581 delta_(1.0),
582 d_(1.0),
583 kappa_(0.41),
584 Uref_(1.0),
585 Lref_(1.0),
586 scale_(1.0),
587 m_(0.5),
588 nCellPerEddy_(5),
589
590 patchArea_(-1),
591 triFace_(),
592 triCumulativeMagSf_(),
593 sumTriMagSf_(Pstream::nProcs() + 1, Zero),
594 patchNormal_(Zero),
595 patchBounds_(),
596
597 eddies_(Zero),
598 v0_(Zero),
599 rndGen_(Pstream::myProcNo()),
600 sigmax_(size(), Zero),
601 maxSigmaX_(Zero),
602 nEddy_(0),
603 curTimeIndex_(-1),
604 singleProc_(false),
605 writeEddies_(false)
606{}
607
608
611(
613 const fvPatch& p,
615 const fvPatchFieldMapper& mapper
616)
617:
618 fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
619 U_(ptf.U_.clone(patch().patch())),
620 R_(ptf.R_.clone(patch().patch())),
621 L_(ptf.L_.clone(patch().patch())),
622 delta_(ptf.delta_),
623 d_(ptf.d_),
624 kappa_(ptf.kappa_),
625 Uref_(ptf.Uref_),
626 Lref_(ptf.Lref_),
627 scale_(ptf.scale_),
628 m_(ptf.m_),
629 nCellPerEddy_(ptf.nCellPerEddy_),
630
631 patchArea_(ptf.patchArea_),
632 triFace_(ptf.triFace_),
633 triCumulativeMagSf_(ptf.triCumulativeMagSf_),
634 sumTriMagSf_(ptf.sumTriMagSf_),
635 patchNormal_(ptf.patchNormal_),
636 patchBounds_(ptf.patchBounds_),
637
638 eddies_(ptf.eddies_),
639 v0_(ptf.v0_),
640 rndGen_(ptf.rndGen_),
641 sigmax_(ptf.sigmax_, mapper),
642 maxSigmaX_(ptf.maxSigmaX_),
643 nEddy_(ptf.nEddy_),
644 curTimeIndex_(-1),
645 singleProc_(ptf.singleProc_),
646 writeEddies_(ptf.writeEddies_)
647{}
648
649
652(
653 const fvPatch& p,
655 const dictionary& dict
656)
657:
659 U_(PatchFunction1<vector>::New(patch().patch(), "U", dict)),
660 R_(PatchFunction1<symmTensor>::New(patch().patch(), "R", dict)),
661 L_(PatchFunction1<scalar>::New(patch().patch(), "L", dict)),
662 delta_(dict.getCheck<scalar>("delta", scalarMinMax::ge(0))),
663 d_(dict.getCheckOrDefault<scalar>("d", 1, scalarMinMax::ge(SMALL))),
664 kappa_(dict.getCheckOrDefault<scalar>("kappa", 0.41, scalarMinMax::ge(0))),
665 Uref_(dict.getCheckOrDefault<scalar>("Uref", 1, scalarMinMax::ge(SMALL))),
666 Lref_(dict.getCheckOrDefault<scalar>("Lref", 1, scalarMinMax::ge(SMALL))),
667 scale_(dict.getCheckOrDefault<scalar>("scale", 1, scalarMinMax::ge(0))),
668 m_(dict.getCheckOrDefault<scalar>("m", 0.5, scalarMinMax::ge(0))),
669 nCellPerEddy_(dict.getOrDefault<label>("nCellPerEddy", 5)),
670
671 patchArea_(-1),
672 triFace_(),
673 triCumulativeMagSf_(),
674 sumTriMagSf_(Pstream::nProcs() + 1, Zero),
675 patchNormal_(Zero),
676 patchBounds_(),
677
678 eddies_(),
679 v0_(Zero),
680 rndGen_(),
681 sigmax_(size(), Zero),
682 maxSigmaX_(Zero),
683 nEddy_(0),
684 curTimeIndex_(-1),
685 singleProc_(false),
686 writeEddies_(dict.getOrDefault("writeEddies", false))
687{
689
690 const scalar t = db().time().timeOutputValue();
691 const symmTensorField R(R_->value(t)/sqr(Uref_));
692
694}
695
696
699(
702)
703:
705 U_(ptf.U_.clone(patch().patch())),
706 R_(ptf.R_.clone(patch().patch())),
707 L_(ptf.L_.clone(patch().patch())),
708 delta_(ptf.delta_),
709 d_(ptf.d_),
710 kappa_(ptf.kappa_),
711 Uref_(ptf.Uref_),
712 Lref_(ptf.Lref_),
713 scale_(ptf.scale_),
714 m_(ptf.m_),
715 nCellPerEddy_(ptf.nCellPerEddy_),
716
717 patchArea_(ptf.patchArea_),
718 triFace_(ptf.triFace_),
719 triCumulativeMagSf_(ptf.triCumulativeMagSf_),
720 sumTriMagSf_(ptf.sumTriMagSf_),
721 patchNormal_(ptf.patchNormal_),
722 patchBounds_(ptf.patchBounds_),
723
724 eddies_(ptf.eddies_),
725 v0_(ptf.v0_),
726 rndGen_(ptf.rndGen_),
727 sigmax_(ptf.sigmax_),
728 maxSigmaX_(ptf.maxSigmaX_),
729 nEddy_(ptf.nEddy_),
730 curTimeIndex_(-1),
731 singleProc_(ptf.singleProc_),
732 writeEddies_(ptf.writeEddies_)
733{}
734
735
738(
740)
742 turbulentDFSEMInletFvPatchVectorField(ptf, ptf.internalField())
743{}
744
745
746// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
747
749(
750 const symmTensorField& R
751)
752{
753 constexpr label maxDiffs = 5;
754 label nDiffs = 0;
755
756 // (S:Eq. 4a-4c)
757 forAll(R, i)
758 {
759 bool diff = false;
760
761 if (maxDiffs < nDiffs)
762 {
763 Info<< "More than " << maxDiffs << " times"
764 << " Reynolds-stress realizability checks failed."
765 << " Skipping further comparisons." << endl;
766 return;
767 }
768
769 const symmTensor& r = R[i];
770
771 if (r.xx() < 0)
772 {
774 << "Reynolds stress " << r << " at index " << i
775 << " does not obey the constraint: Rxx >= 0"
776 << endl;
777 diff = true;
778 }
779
780 if ((r.xx()*r.yy() - sqr(r.xy())) < 0)
781 {
783 << "Reynolds stress " << r << " at index " << i
784 << " does not obey the constraint: Rxx*Ryy - sqr(Rxy) >= 0"
785 << endl;
786 diff = true;
787 }
788
789 if (det(r) < 0)
790 {
792 << "Reynolds stress " << r << " at index " << i
793 << " does not obey the constraint: det(R) >= 0"
794 << endl;
795 diff = true;
796 }
797
798 if (diff)
800 ++nDiffs;
801 }
802 }
803}
804
805
807(
808 const scalarField& R
809)
810{
811 if (min(R) <= 0)
812 {
814 << "Reynolds stresses contain at least one "
815 << "nonpositive element. min(R) = " << min(R)
816 << exit(FatalError);
817 }
818}
819
820
822(
823 const fvPatchFieldMapper& m
824)
825{
827
828 if (U_)
829 {
830 U_->autoMap(m);
831 }
832 if (R_)
833 {
834 R_->autoMap(m);
835 }
836 if (L_)
837 {
838 L_->autoMap(m);
839 }
840
841 sigmax_.autoMap(m);
842}
843
844
846(
847 const fvPatchVectorField& ptf,
848 const labelList& addr
849)
850{
852
853 const auto& dfsemptf =
855
856 if (U_)
857 {
858 U_->rmap(dfsemptf.U_(), addr);
859 }
860 if (R_)
861 {
862 R_->rmap(dfsemptf.R_(), addr);
863 }
864 if (L_)
865 {
866 L_->rmap(dfsemptf.L_(), addr);
867 }
868
869 sigmax_.rmap(dfsemptf.sigmax_, addr);
870}
871
872
874{
875 if (updated())
876 {
877 return;
878 }
879
880 if (curTimeIndex_ == -1)
881 {
882 initialisePatch();
883
884 initialiseEddyBox();
885
886 initialiseEddies();
887 }
888
889
890 if (curTimeIndex_ != db().time().timeIndex())
891 {
893 U_->value(db().time().timeOutputValue())/Uref_;
894
895 // (PCR:p. 522)
896 const vector UBulk
897 (
898 gWeightedAverage(patch().magSf(), UMean())
899 );
900
901 // Move eddies using bulk velocity
902 const scalar deltaT = db().time().deltaTValue();
903 convectEddies(UBulk, deltaT);
904
905 // Set mean velocity
906 vectorField& U = *this;
907 U = UMean;
908
909 // Apply second part of normalisation coefficient
910 const scalar c =
911 scale_*Foam::pow(10*v0_, m_)/Foam::sqrt(scalar(nEddy_));
912
913 // In parallel, need to collect all eddies that will interact with
914 // local faces
915
916 const pointField& Cf = patch().Cf();
917
918 if (singleProc_ || !Pstream::parRun())
919 {
920 forAll(U, faceI)
921 {
922 U[faceI] += c*uPrimeEddy(eddies_, Cf[faceI]);
923 }
924 }
925 else
926 {
927 // Process local eddy contributions
928 forAll(U, faceI)
929 {
930 U[faceI] += c*uPrimeEddy(eddies_, Cf[faceI]);
931 }
932
933 // Add contributions from overlapping eddies
934 List<List<eddy>> overlappingEddies(Pstream::nProcs());
935 calcOverlappingProcEddies(overlappingEddies);
936
937 for (const List<eddy>& eddies : overlappingEddies)
938 {
939 if (eddies.size())
940 {
941 forAll(U, faceI)
942 {
943 U[faceI] += c*uPrimeEddy(eddies, Cf[faceI]);
944 }
945 }
946 }
947 }
948
949 // Re-scale to ensure correct flow rate
950 const scalar fCorr =
951 gSum((UBulk & patchNormal_)*patch().magSf())
952 /gSum(U & -patch().Sf());
953
954 U *= fCorr;
955
956 curTimeIndex_ = db().time().timeIndex();
957
958 if (writeEddies_)
959 {
960 writeEddyOBJ();
961 }
962
963 if (debug)
964 {
965 auto limits = gMinMax(*this);
966
967 Info<< "Magnitude of bulk velocity: " << UBulk << endl;
968
969 Info<< "Number of eddies: "
970 << returnReduce(eddies_.size(), sumOp<label>())
971 << endl;
972
973 Info<< "Patch:" << patch().patch().name()
974 << " min/max(U):" << limits.min() << ", " << limits.max()
975 << endl;
976
977 if (db().time().writeTime())
978 {
979 writeLumleyCoeffs();
980 }
982 }
983
984 fixedValueFvPatchVectorField::updateCoeffs();
985}
986
987
989{
991 os.writeEntry("delta", delta_);
992 os.writeEntryIfDifferent<scalar>("d", 1.0, d_);
993 os.writeEntryIfDifferent<scalar>("kappa", 0.41, kappa_);
994 os.writeEntryIfDifferent<scalar>("Uref", 1.0, Uref_);
995 os.writeEntryIfDifferent<scalar>("Lref", 1.0, Lref_);
996 os.writeEntryIfDifferent<scalar>("scale", 1.0, scale_);
997 os.writeEntryIfDifferent<scalar>("m", 0.5, m_);
998 os.writeEntryIfDifferent<label>("nCellPerEddy", 5, nCellPerEddy_);
999 os.writeEntryIfDifferent("writeEddies", false, writeEddies_);
1000 if (U_)
1001 {
1002 U_->writeData(os);
1003 }
1004 if (R_)
1005 {
1006 R_->writeData(os);
1007 }
1008 if (L_)
1009 {
1010 L_->writeData(os);
1011 }
1014
1015
1016// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1017
1018namespace Foam
1019{
1021 (
1023 turbulentDFSEMInletFvPatchVectorField
1024 );
1025}
1026
1027
1028// ************************************************************************* //
label k
#define R(A, B, C, D, E, F, K, M)
Macros for easy insertion into run-time selection tables.
volVectorField UMean(UMeanHeader, mesh)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
SubList< scalar > subList
Definition List.H:129
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition Ostream.H:331
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:346
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void finishedSends(const bool wait=true)
Mark the send phase as being finished.
Inter-processor communications stream.
Definition Pstream.H:59
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual communication.
const Cmpt & yy() const noexcept
Definition SymmTensor.H:154
const Cmpt & xx() const noexcept
Definition SymmTensor.H:150
const Cmpt & xy() const noexcept
Definition SymmTensor.H:151
Input inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Definition UIPstream.H:313
Output inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Definition UOPstream.H:408
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static void reduceAnd(bool &value, const int communicator=worldComm)
Logical (and) reduction (MPI_AllReduce).
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition UPstream.H:1948
static const Form max
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Class to describe eddies for the turbulentDFSEMInletFvPatchVectorField boundary condition.
Definition eddy.H:64
static int debug
Flag to activate debug statements.
Definition eddy.H:172
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition error.H:74
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
A FieldMapper for finite-volume patch fields.
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 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
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
A class for managing temporary objects.
Definition tmp.H:75
@ TAB
Tab [isspace].
Definition token.H:142
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.
turbulentDFSEMInletFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual tmp< fvPatchField< vector > > clone() const
Return a clone.
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.
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
U
Definition pEqn.H:72
volScalarField & p
auto limits
Definition setRDeltaT.H:186
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
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 dimensionedScalar c
Speed of light in a vacuum.
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar det(const dimensionedSphericalTensor &dt)
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
dimensionedScalar pos(const dimensionedScalar &ds)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
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)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
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)
Cmpt invariantII(const SymmTensor2D< Cmpt > &st)
Return the 2nd invariant of a SymmTensor2D.
messageStream Info
Information stream (stdout output on master, null elsewhere).
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
Cmpt invariantIII(const SymmTensor< Cmpt > &st)
Return the 3rd invariant of a SymmTensor.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition triad.C:373
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
Type gMax(const FieldField< Field, Type > &f)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
fvPatchField< vector > fvPatchVectorField
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition symmTensor.H:55
fileName search(const word &file, const fileName &directory)
Recursively search the given directory for the file.
Definition fileName.C:642
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
label timeIndex
labelList f(nPoints)
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition stdFoam.H:315
const vector L(dict.get< vector >("L"))