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-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
30 #include "volFields.H"
32 #include "fvPatchFieldMapper.H"
33 #include "momentOfInertia.H"
34 #include "Fstream.H"
35 #include "globalIndex.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 Foam::label Foam::turbulentDFSEMInletFvPatchVectorField::seedIterMax_ = 1000;
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
43 void Foam::turbulentDFSEMInletFvPatchVectorField::writeEddyOBJ() const
44 {
45  {
46  // Output the bounding box
47  OFstream os(db().time().path()/"eddyBox.obj");
48 
49  const polyPatch& pp = this->patch().patch();
50  const labelList& boundaryPoints = pp.boundaryPoints();
51  const pointField& localPoints = pp.localPoints();
52 
53  vector offset = patchNormal_*maxSigmaX_;
54  forAll(boundaryPoints, i)
55  {
56  point p = localPoints[boundaryPoints[i]];
57  p += offset;
58  os << "v " << p.x() << " " << p.y() << " " << p.z() << nl;
59  }
60 
61  forAll(boundaryPoints, i)
62  {
63  point p = localPoints[boundaryPoints[i]];
64  p -= offset;
65  os << "v " << p.x() << " " << p.y() << " " << p.z() << nl;
66  }
67 
68  // Draw lines between points
69  // Note: need to order to avoid crossing patch
70  //const label nPoint = boundaryPoints.size();
71  //
72  //forAll(boundaryPoints, i)
73  //{
74  // label i1 = i;
75  // label i2 = (i + 1) % nPoint;
76  // os << "l " << i1 << " " << i2 << nl;
77  //}
78  //
79  //forAll(boundaryPoints, i)
80  //{
81  // label i1 = i + nPoint;
82  // label i2 = ((i + 1) % nPoint) + nPoint;
83  // os << "l " << i1 << " " << i2 << nl;
84  //}
85  }
86 
87  {
88  const Time& time = db().time();
89  OFstream os
90  (
91  time.path()/"eddies_" + Foam::name(time.timeIndex()) + ".obj"
92  );
93 
94  label pointOffset = 0;
95  forAll(eddies_, eddyI)
96  {
97  const eddy& e = eddies_[eddyI];
98  pointOffset += e.writeSurfaceOBJ(pointOffset, patchNormal_, os);
99  }
100  }
101 }
102 
103 
104 void Foam::turbulentDFSEMInletFvPatchVectorField::writeLumleyCoeffs() const
105 {
106  // Output list of xi vs eta
107 
108  // Before interpolation/raw data
109  if (interpolateR_)
110  {
111  fileName valsFile
112  (
113  fileHandler().filePath
114  (
115  fileName
116  (
117  db().time().path()
118  /db().time().caseConstant()
119  /"boundaryData"
120  /this->patch().name()
121  /"0"
122  /"R"
123  )
124  )
125  );
126 
127  autoPtr<ISstream> isPtr
128  (
129  fileHandler().NewIFstream
130  (
131  valsFile
132  )
133  );
134 
135  Field<symmTensor> Rexp(isPtr());
136 
137  OFstream os(db().time().path()/"lumley_input.out");
138 
139  os << "# xi" << token::TAB << "eta" << endl;
140 
141  forAll(Rexp, faceI)
142  {
143  // Normalised anisotropy tensor
144  symmTensor devR = dev(Rexp[faceI]/(tr(Rexp[faceI])));
145 
146  // Second tensor invariant
147  scalar ii = min(0, invariantII(devR));
148 
149  // Third tensor invariant
150  scalar iii = invariantIII(devR);
151 
152  // xi, eta
153  // See Pope - characterization of Reynolds-stress anisotropy
154  scalar xi = cbrt(0.5*iii);
155  scalar eta = sqrt(-ii/3.0);
156  os << xi << token::TAB << eta << token::TAB
157  << ii << token::TAB << iii << endl;
158  }
159  }
160 
161  // After interpolation
162  {
163  OFstream os(db().time().path()/"lumley_interpolated.out");
164 
165  os << "# xi" << token::TAB << "eta" << endl;
166 
167  forAll(R_, faceI)
168  {
169  // Normalised anisotropy tensor
170  symmTensor devR = dev(R_[faceI]/(tr(R_[faceI])));
171 
172  // Second tensor invariant
173  scalar ii = min(0, invariantII(devR));
174 
175  // Third tensor invariant
176  scalar iii = invariantIII(devR);
177 
178  // xi, eta
179  // See Pope - characterization of Reynolds-stress anisotropy
180  scalar xi = cbrt(0.5*iii);
181  scalar eta = sqrt(-ii/3.0);
182  os << xi << token::TAB << eta << token::TAB
183  << ii << token::TAB << iii << endl;
184  }
185  }
186 }
187 
188 
190 Foam::turbulentDFSEMInletFvPatchVectorField::patchMapper() const
191 {
192  // Initialise interpolation (2D planar interpolation by triangulation)
193  if (mapperPtr_.empty())
194  {
195  // Reread values and interpolate
196  fileName samplePointsFile
197  (
198  this->db().time().path()
199  /this->db().time().caseConstant()
200  /"boundaryData"
201  /this->patch().name()
202  /"points"
203  );
204 
205  pointField samplePoints((IFstream(samplePointsFile)()));
206 
207  if (debug)
208  {
210  << " Read " << samplePoints.size() << " sample points from "
211  << samplePointsFile << endl;
212  }
213 
214 
215  // tbd: run-time selection
216  bool nearestOnly =
217  (
218  !mapMethod_.empty()
219  && mapMethod_ != "planarInterpolation"
220  );
221 
222  // Allocate the interpolator
223  mapperPtr_.reset
224  (
225  new pointToPointPlanarInterpolation
226  (
227  samplePoints,
228  this->patch().patch().faceCentres(),
229  perturb_,
230  nearestOnly
231  )
232  );
233  }
234 
235  return *mapperPtr_;
236 }
237 
238 
239 void Foam::turbulentDFSEMInletFvPatchVectorField::initialisePatch()
240 {
241  const vectorField nf(patch().nf());
242 
243  // Patch normal points into domain
244  patchNormal_ = -gAverage(nf);
245 
246  // Check that patch is planar
247  scalar error = max(magSqr(patchNormal_ + nf));
248 
249  if (error > SMALL)
250  {
252  << "Patch " << patch().name() << " is not planar"
253  << endl;
254  }
255 
256  patchNormal_ /= mag(patchNormal_) + ROOTVSMALL;
257 
258 
259  // Decompose the patch faces into triangles to enable point search
260 
261  const polyPatch& patch = this->patch().patch();
262  const pointField& points = patch.points();
263 
264  // Triangulate the patch faces and create addressing
265  DynamicList<label> triToFace(2*patch.size());
266  DynamicList<scalar> triMagSf(2*patch.size());
267  DynamicList<face> triFace(2*patch.size());
268  DynamicList<face> tris(5);
269 
270  // Set zero value at the start of the tri area list
271  triMagSf.append(0.0);
272 
273  forAll(patch, faceI)
274  {
275  const face& f = patch[faceI];
276 
277  tris.clear();
278  f.triangles(points, tris);
279 
280  forAll(tris, i)
281  {
282  triToFace.append(faceI);
283  triFace.append(tris[i]);
284  triMagSf.append(tris[i].mag(points));
285  }
286  }
287 
288  forAll(sumTriMagSf_, i)
289  {
290  sumTriMagSf_[i] = 0.0;
291  }
292 
293  sumTriMagSf_[Pstream::myProcNo() + 1] = sum(triMagSf);
294 
295  Pstream::listCombineGather(sumTriMagSf_, maxEqOp<scalar>());
296  Pstream::listCombineScatter(sumTriMagSf_);
297 
298  for (label i = 1; i < triMagSf.size(); i++)
299  {
300  triMagSf[i] += triMagSf[i-1];
301  }
302 
303  // Transfer to persistent storage
304  triFace_.transfer(triFace);
305  triToFace_.transfer(triToFace);
306  triCumulativeMagSf_.transfer(triMagSf);
307 
308  // Convert sumTriMagSf_ into cumulative sum of areas per proc
309  for (label i = 1; i < sumTriMagSf_.size(); i++)
310  {
311  sumTriMagSf_[i] += sumTriMagSf_[i-1];
312  }
313 
314  // Global patch area (over all processors)
315  patchArea_ = sumTriMagSf_.last();
316 
317  // Local patch bounds (this processor)
318  patchBounds_ = boundBox(patch.localPoints(), false);
319  patchBounds_.inflate(0.1);
320 
321  // Determine if all eddies spawned from a single processor
322  singleProc_ = patch.size() == returnReduce(patch.size(), sumOp<label>());
323  reduce(singleProc_, orOp<bool>());
324 }
325 
326 
327 void Foam::turbulentDFSEMInletFvPatchVectorField::initialiseEddyBox()
328 {
329  const scalarField& magSf = patch().magSf();
330 
331  //const scalarField cellDx(Foam::sqrt(magSf));
332  const scalarField cellDx(max(Foam::sqrt(magSf), 2/patch().deltaCoeffs()));
333 
334  // Inialise eddy box extents
335  forAll(*this, faceI)
336  {
337  scalar& s = sigmax_[faceI];
338 
339  // Length scale in x direction (based on eq. 14)
340  s = mag(L_[faceI]);
341  s = min(s, kappa_*delta_);
342 
343  // Allow eddies to be smaller than the mesh scale as suggested by
344  // the reference?
345  // s = min(s, nCellPerEddy_*cellDx[faceI]);
346  s = max(s, nCellPerEddy_*cellDx[faceI]);
347  }
348 
349  // Maximum extent across all processors
350  maxSigmaX_ = gMax(sigmax_);
351 
352  // Eddy box volume
353  v0_ = 2*gSum(magSf)*maxSigmaX_;
354 
355  if (debug)
356  {
357  Info<< "Patch: " << patch().patch().name() << " eddy box:" << nl
358  << " volume : " << v0_ << nl
359  << " maxSigmaX : " << maxSigmaX_ << nl
360  << endl;
361  }
362 }
363 
364 
365 Foam::pointIndexHit Foam::turbulentDFSEMInletFvPatchVectorField::setNewPosition
366 (
367  const bool global
368 )
369 {
370  // Initialise to miss
371  pointIndexHit pos(false, vector::max, -1);
372 
373  const polyPatch& patch = this->patch().patch();
374  const pointField& points = patch.points();
375 
376  if (global)
377  {
378  scalar areaFraction = rndGen_.globalPosition<scalar>(0, patchArea_);
379 
380  // Determine which processor to use
381  label procI = 0;
382  forAllReverse(sumTriMagSf_, i)
383  {
384  if (areaFraction >= sumTriMagSf_[i])
385  {
386  procI = i;
387  break;
388  }
389  }
390 
391  if (Pstream::myProcNo() == procI)
392  {
393  // Find corresponding decomposed face triangle
394  label triI = 0;
395  scalar offset = sumTriMagSf_[procI];
396  forAllReverse(triCumulativeMagSf_, i)
397  {
398  if (areaFraction > triCumulativeMagSf_[i] + offset)
399  {
400  triI = i;
401  break;
402  }
403  }
404 
405  // Find random point in triangle
406  const face& tf = triFace_[triI];
407  const triPointRef tri(points[tf[0]], points[tf[1]], points[tf[2]]);
408 
409  pos.setHit();
410  pos.setIndex(triToFace_[triI]);
411  pos.rawPoint() = tri.randomPoint(rndGen_);
412  }
413  }
414  else
415  {
416  // Find corresponding decomposed face triangle on local processor
417  label triI = 0;
418  scalar maxAreaLimit = triCumulativeMagSf_.last();
419  scalar areaFraction = rndGen_.position<scalar>(0, maxAreaLimit);
420 
421  forAllReverse(triCumulativeMagSf_, i)
422  {
423  if (areaFraction > triCumulativeMagSf_[i])
424  {
425  triI = i;
426  break;
427  }
428  }
429 
430  // Find random point in triangle
431  const face& tf = triFace_[triI];
432  const triPointRef tri(points[tf[0]], points[tf[1]], points[tf[2]]);
433 
434  pos.setHit();
435  pos.setIndex(triToFace_[triI]);
436  pos.rawPoint() = tri.randomPoint(rndGen_);
437  }
438 
439  return pos;
440 }
441 
442 
443 void Foam::turbulentDFSEMInletFvPatchVectorField::initialiseEddies()
444 {
445  DynamicList<eddy> eddies(size());
446 
447  // Initialise eddy properties
448  scalar sumVolEddy = 0;
449  scalar sumVolEddyAllProc = 0;
450 
451  while (sumVolEddyAllProc/v0_ < d_)
452  {
453  bool search = true;
454  label iter = 0;
455 
456  while (search && iter++ < seedIterMax_)
457  {
458  // Get new parallel consistent position
459  pointIndexHit pos(setNewPosition(true));
460  label faceI = pos.index();
461 
462  // Note: only 1 processor will pick up this face
463  if (faceI != -1)
464  {
465  eddy e
466  (
467  faceI,
468  pos.hitPoint(),
469  rndGen_.position<scalar>(-maxSigmaX_, maxSigmaX_),
470  sigmax_[faceI],
471  R_[faceI],
472  rndGen_
473  );
474 
475  // If eddy valid, patchFaceI is non-zero
476  if (e.patchFaceI() != -1)
477  {
478  eddies.append(e);
479  sumVolEddy += e.volume();
480  search = false;
481  }
482  }
483  // else eddy on remote processor
484 
485  reduce(search, andOp<bool>());
486  }
487 
488 
489  sumVolEddyAllProc = returnReduce(sumVolEddy, sumOp<scalar>());
490  }
491  eddies_.transfer(eddies);
492 
493  nEddy_ = eddies_.size();
494 
495  if (debug)
496  {
497  Pout<< "Patch:" << patch().patch().name();
498 
499  if (Pstream::parRun())
500  {
501  Pout<< " processor:" << Pstream::myProcNo();
502  }
503 
504  Pout<< " seeded:" << nEddy_ << " eddies" << endl;
505  }
506 
507  reduce(nEddy_, sumOp<label>());
508 
509  if (nEddy_ > 0)
510  {
511  Info<< "Turbulent DFSEM patch: " << patch().name()
512  << " seeded " << nEddy_ << " eddies with total volume "
513  << sumVolEddyAllProc
514  << endl;
515  }
516  else
517  {
519  << "Patch: " << patch().patch().name()
520  << " on field " << internalField().name()
521  << ": No eddies seeded - please check your set-up" << endl;
522  }
523 }
524 
525 
526 void Foam::turbulentDFSEMInletFvPatchVectorField::convectEddies
527 (
528  const scalar deltaT
529 )
530 {
531  // Note: all operations applied to local processor only
532 
533  label nRecycled = 0;
534 
535  forAll(eddies_, eddyI)
536  {
537  eddy& e = eddies_[eddyI];
538  e.move(deltaT*(UMean_ & patchNormal_));
539 
540  const scalar position0 = e.x();
541 
542  // Check to see if eddy has exited downstream box plane
543  if (position0 > maxSigmaX_)
544  {
545  bool search = true;
546  label iter = 0;
547 
548  while (search && iter++ < seedIterMax_)
549  {
550  // Spawn new eddy with new random properties (intensity etc)
551  pointIndexHit pos(setNewPosition(false));
552  label faceI = pos.index();
553 
554  e = eddy
555  (
556  faceI,
557  pos.hitPoint(),
558  position0 - floor(position0/maxSigmaX_)*maxSigmaX_,
559  sigmax_[faceI],
560  R_[faceI],
561  rndGen_
562  );
563 
564  if (e.patchFaceI() != -1)
565  {
566  search = false;
567  }
568  }
569 
570  nRecycled++;
571  }
572  }
573 
574  reduce(nRecycled, sumOp<label>());
575 
576  if (debug && nRecycled > 0)
577  {
578  Info<< "Patch: " << patch().patch().name() << " recycled "
579  << nRecycled << " eddies" << endl;
580  }
581 }
582 
583 
584 Foam::vector Foam::turbulentDFSEMInletFvPatchVectorField::uDashEddy
585 (
586  const List<eddy>& eddies,
587  const point& patchFaceCf
588 ) const
589 {
590  vector uDash(Zero);
591 
592  forAll(eddies, k)
593  {
594  const eddy& e = eddies[k];
595  uDash += e.uDash(patchFaceCf, patchNormal_);
596  }
597 
598  return uDash;
599 }
600 
601 
602 void Foam::turbulentDFSEMInletFvPatchVectorField::calcOverlappingProcEddies
603 (
604  List<List<eddy>>& overlappingEddies
605 ) const
606 {
607  int oldTag = UPstream::msgType();
608  UPstream::msgType() = oldTag + 1;
609 
610  List<boundBox> patchBBs(Pstream::nProcs());
611  patchBBs[Pstream::myProcNo()] = patchBounds_;
612  Pstream::gatherList(patchBBs);
613  Pstream::scatterList(patchBBs);
614 
615  // Per processor indices into all segments to send
616  List<DynamicList<label>> dynSendMap(Pstream::nProcs());
617 
618  forAll(eddies_, i)
619  {
620  // Collect overlapping eddies
621  const eddy& e = eddies_[i];
622 
623  // Eddy bounds
624  point x = e.position(patchNormal_);
625  boundBox ebb = e.bounds();
626  ebb.min() += x;
627  ebb.max() += x;
628 
629  forAll(patchBBs, procI)
630  {
631  // Not including intersection with local patch
632  if (procI != Pstream::myProcNo())
633  {
634  if (ebb.overlaps(patchBBs[procI]))
635  {
636  dynSendMap[procI].append(i);
637  }
638  }
639  }
640  }
641 
642  labelListList sendMap(Pstream::nProcs());
643  forAll(sendMap, procI)
644  {
645  sendMap[procI].transfer(dynSendMap[procI]);
646  }
647 
648  // Send the number of eddies for local processors to receive
649  labelListList sendSizes(Pstream::nProcs());
650  sendSizes[Pstream::myProcNo()].setSize(Pstream::nProcs());
651  forAll(sendMap, procI)
652  {
653  sendSizes[Pstream::myProcNo()][procI] = sendMap[procI].size();
654  }
655  Pstream::gatherList(sendSizes);
656  Pstream::scatterList(sendSizes);
657 
658  // Determine order of receiving
659  labelListList constructMap(Pstream::nProcs());
660 
661  // Local segment first
662  constructMap[Pstream::myProcNo()] = identity
663  (
664  sendMap[Pstream::myProcNo()].size()
665  );
666 
667  label segmentI = constructMap[Pstream::myProcNo()].size();
668  forAll(constructMap, procI)
669  {
670  if (procI != Pstream::myProcNo())
671  {
672  // What I need to receive is what other processor is sending to me
673  label nRecv = sendSizes[procI][Pstream::myProcNo()];
674  constructMap[procI].setSize(nRecv);
675 
676  for (label i = 0; i < nRecv; i++)
677  {
678  constructMap[procI][i] = segmentI++;
679  }
680  }
681  }
682 
683  mapDistribute map(segmentI, std::move(sendMap), std::move(constructMap));
684 
685  PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
686 
687  for (label domain = 0; domain < Pstream::nProcs(); domain++)
688  {
689  const labelList& sendElems = map.subMap()[domain];
690 
691  if (domain != Pstream::myProcNo() && sendElems.size())
692  {
693  List<eddy> subEddies(UIndirectList<eddy>(eddies_, sendElems));
694 
695  UOPstream toDomain(domain, pBufs);
696 
697  toDomain<< subEddies;
698  }
699  }
700 
701  // Start receiving
702  pBufs.finishedSends();
703 
704  // Consume
705  for (label domain = 0; domain < Pstream::nProcs(); domain++)
706  {
707  const labelList& recvElems = map.constructMap()[domain];
708 
709  if (domain != Pstream::myProcNo() && recvElems.size())
710  {
711  UIPstream str(domain, pBufs);
712  {
713  str >> overlappingEddies[domain];
714  }
715  }
716  }
717 
718  // Restore tag
719  UPstream::msgType() = oldTag;
720 }
721 
722 
723 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
724 
727 (
728  const fvPatch& p,
730 )
731 :
733  delta_(Zero),
734  d_(Zero),
735  kappa_(Zero),
736 
737  perturb_(1e-5),
738  mapMethod_("nearestCell"),
739  mapperPtr_(nullptr),
740  interpolateR_(false),
741  interpolateL_(false),
742  interpolateU_(false),
743  R_(),
744  L_(),
745  U_(),
746  UMean_(Zero),
747 
748  patchArea_(-1),
749  triFace_(),
750  triToFace_(),
751  triCumulativeMagSf_(),
752  sumTriMagSf_(Pstream::nProcs() + 1, Zero),
753 
754  eddies_(Zero),
755  nCellPerEddy_(5),
756  patchNormal_(Zero),
757  v0_(Zero),
758  rndGen_(Pstream::myProcNo()),
759  sigmax_(size(), Zero),
760  maxSigmaX_(Zero),
761  nEddy_(Zero),
762  curTimeIndex_(-1),
763  patchBounds_(boundBox::invertedBox),
764  singleProc_(false),
765  writeEddies_(false)
766 {}
767 
768 
771 (
773  const fvPatch& p,
775  const fvPatchFieldMapper& mapper
776 )
777 :
778  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
779  delta_(ptf.delta_),
780  d_(ptf.d_),
781  kappa_(ptf.kappa_),
782 
783  perturb_(ptf.perturb_),
784  mapMethod_(ptf.mapMethod_),
785  mapperPtr_(nullptr),
786  interpolateR_(ptf.interpolateR_),
787  interpolateL_(ptf.interpolateL_),
788  interpolateU_(ptf.interpolateU_),
789  R_(ptf.R_, mapper),
790  L_(ptf.L_, mapper),
791  U_(ptf.U_, mapper),
792  UMean_(ptf.UMean_),
793 
794  patchArea_(ptf.patchArea_),
795  triFace_(ptf.triFace_),
796  triToFace_(ptf.triToFace_),
797  triCumulativeMagSf_(ptf.triCumulativeMagSf_),
798  sumTriMagSf_(ptf.sumTriMagSf_),
799 
800  eddies_(ptf.eddies_),
801  nCellPerEddy_(ptf.nCellPerEddy_),
802  patchNormal_(ptf.patchNormal_),
803  v0_(ptf.v0_),
804  rndGen_(ptf.rndGen_),
805  sigmax_(ptf.sigmax_, mapper),
806  maxSigmaX_(ptf.maxSigmaX_),
807  nEddy_(Zero),
808  curTimeIndex_(-1),
809  patchBounds_(ptf.patchBounds_),
810  singleProc_(ptf.singleProc_),
811  writeEddies_(ptf.writeEddies_)
812 {}
813 
814 
817 (
818  const fvPatch& p,
820  const dictionary& dict
821 )
822 :
824  delta_(dict.get<scalar>("delta")),
825  d_(dict.getOrDefault<scalar>("d", 1.0)),
826  kappa_(dict.getOrDefault<scalar>("kappa", 0.41)),
827 
828  perturb_(dict.getOrDefault<scalar>("perturb", 1e-5)),
829  mapMethod_(dict.getOrDefault<word>("mapMethod", "nearestCell")),
830  mapperPtr_(nullptr),
831  interpolateR_(dict.getOrDefault<bool>("interpolateR", false)),
832  interpolateL_(dict.getOrDefault<bool>("interpolateL", false)),
833  interpolateU_(dict.getOrDefault<bool>("interpolateU", false)),
834  R_(interpolateOrRead<symmTensor>("R", dict, interpolateR_)),
835  L_(interpolateOrRead<scalar>("L", dict, interpolateL_)),
836  U_(interpolateOrRead<vector>("U", dict, interpolateU_)),
837  UMean_(Zero),
838 
839  patchArea_(-1),
840  triFace_(),
841  triToFace_(),
842  triCumulativeMagSf_(),
843  sumTriMagSf_(Pstream::nProcs() + 1, Zero),
844 
845  eddies_(),
846  nCellPerEddy_(dict.getOrDefault<label>("nCellPerEddy", 5)),
847  patchNormal_(Zero),
848  v0_(Zero),
849  rndGen_(),
850  sigmax_(size(), Zero),
851  maxSigmaX_(Zero),
852  nEddy_(Zero),
853  curTimeIndex_(-1),
854  patchBounds_(boundBox::invertedBox),
855  singleProc_(false),
856  writeEddies_(dict.getOrDefault<bool>("writeEddies", false))
857 {
858  eddy::debug = debug;
859 
860  checkStresses(R_);
861 
862  // Set UMean as patch area average value
863  UMean_ = gSum(U_*patch().magSf())/(gSum(patch().magSf()) + ROOTVSMALL);
864 }
865 
866 
869 (
871 )
872 :
874  delta_(ptf.delta_),
875  d_(ptf.d_),
876  kappa_(ptf.kappa_),
877 
878  perturb_(ptf.perturb_),
879  mapMethod_(ptf.mapMethod_),
880  mapperPtr_(nullptr),
881  interpolateR_(ptf.interpolateR_),
882  interpolateL_(ptf.interpolateL_),
883  interpolateU_(ptf.interpolateU_),
884  R_(ptf.R_),
885  L_(ptf.L_),
886  U_(ptf.U_),
887  UMean_(ptf.UMean_),
888 
889  patchArea_(ptf.patchArea_),
890  triFace_(ptf.triFace_),
891  triToFace_(ptf.triToFace_),
892  triCumulativeMagSf_(ptf.triCumulativeMagSf_),
893  sumTriMagSf_(ptf.sumTriMagSf_),
894 
895  eddies_(ptf.eddies_),
896  nCellPerEddy_(ptf.nCellPerEddy_),
897  patchNormal_(ptf.patchNormal_),
898  v0_(ptf.v0_),
899  rndGen_(ptf.rndGen_),
900  sigmax_(ptf.sigmax_),
901  maxSigmaX_(ptf.maxSigmaX_),
902  nEddy_(Zero),
903  curTimeIndex_(-1),
904  patchBounds_(ptf.patchBounds_),
905  singleProc_(ptf.singleProc_),
906  writeEddies_(ptf.writeEddies_)
907 {}
908 
909 
912 (
915 )
916 :
918  delta_(ptf.delta_),
919  d_(ptf.d_),
920  kappa_(ptf.kappa_),
921 
922  perturb_(ptf.perturb_),
923  mapMethod_(ptf.mapMethod_),
924  mapperPtr_(nullptr),
925  interpolateR_(ptf.interpolateR_),
926  interpolateL_(ptf.interpolateL_),
927  interpolateU_(ptf.interpolateU_),
928  R_(ptf.R_),
929  L_(ptf.L_),
930  U_(ptf.U_),
931  UMean_(ptf.UMean_),
932 
933  patchArea_(ptf.patchArea_),
934  triFace_(ptf.triFace_),
935  triToFace_(ptf.triToFace_),
936  triCumulativeMagSf_(ptf.triCumulativeMagSf_),
937  sumTriMagSf_(ptf.sumTriMagSf_),
938 
939  eddies_(ptf.eddies_),
940  nCellPerEddy_(ptf.nCellPerEddy_),
941  patchNormal_(ptf.patchNormal_),
942  v0_(ptf.v0_),
943  rndGen_(ptf.rndGen_),
944  sigmax_(ptf.sigmax_),
945  maxSigmaX_(ptf.maxSigmaX_),
946  nEddy_(Zero),
947  curTimeIndex_(-1),
948  patchBounds_(ptf.patchBounds_),
949  singleProc_(ptf.singleProc_),
950  writeEddies_(ptf.writeEddies_)
951 {}
952 
953 
954 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
955 
957 (
958  const symmTensorField& Rf
959 )
960 {
961  // Perform checks of the stress tensor based on Cholesky decomposition
962  // constraints
963 
964  forAll(Rf, facei)
965  {
966  const symmTensor& R = Rf[facei];
967 
968  if (R.xx() <= 0)
969  {
971  << "Reynolds stress " << R << " at face " << facei
972  << " does not obey the constraint: R_xx > 0"
973  << exit(FatalError);
974  }
975 
976  scalar a_xx = sqrt(R.xx());
977 
978  scalar a_xy = R.xy()/a_xx;
979 
980  scalar a_yy_2 = R.yy() - sqr(a_xy);
981 
982  if (a_yy_2 < 0)
983  {
985  << "Reynolds stress " << R << " at face " << facei
986  << " leads to an invalid Cholesky decomposition due to the "
987  << "constraint R_yy - sqr(a_xy) >= 0"
988  << exit(FatalError);
989  }
990 
991  scalar a_yy = Foam::sqrt(a_yy_2);
992 
993  scalar a_xz = R.xz()/a_xx;
994 
995  scalar a_yz = (R.yz() - a_xy*a_xz)*a_yy;
996 
997  scalar a_zz_2 = R.zz() - sqr(a_xz) - sqr(a_yz);
998 
999  if (a_zz_2 < 0)
1000  {
1002  << "Reynolds stress " << R << " at face " << facei
1003  << " leads to an invalid Cholesky decomposition due to the "
1004  << "constraint R_zz - sqr(a_xz) - sqr(a_yz) >= 0"
1005  << exit(FatalError);
1006  }
1007 
1008  scalar a_zz = Foam::sqrt(a_zz_2);
1009 
1010  if (debug)
1011  {
1012  Pout<< "R: " << R
1013  << " a_xx:" << a_xx << " a_xy:" << a_xy << " a_xz:" << a_xy
1014  << " a_yy:" << a_yy << " a_yz:" << a_yz
1015  << " a_zz:" << a_zz
1016  << endl;
1017  }
1018  }
1019 
1020  return true;
1021 }
1022 
1023 
1026  const fvPatchFieldMapper& m
1027 )
1028 {
1030 
1031  // Clear interpolator
1032  mapperPtr_.clear();
1033  R_.autoMap(m);
1034  L_.autoMap(m);
1035  U_.autoMap(m);
1036 
1037  sigmax_.autoMap(m);
1038 }
1039 
1040 
1043  const fvPatchVectorField& ptf,
1044  const labelList& addr
1045 )
1046 {
1048 
1049  const turbulentDFSEMInletFvPatchVectorField& dfsemptf =
1050  refCast<const turbulentDFSEMInletFvPatchVectorField>(ptf);
1051 
1052  R_.rmap(dfsemptf.R_, addr);
1053  L_.rmap(dfsemptf.L_, addr);
1054  U_.rmap(dfsemptf.U_, addr);
1055 
1056  // Clear interpolator
1057  mapperPtr_.clear();
1058 
1059  sigmax_.rmap(dfsemptf.sigmax_, addr);
1060 }
1061 
1062 
1064 {
1065  if (updated())
1066  {
1067  return;
1068  }
1069 
1070  if (curTimeIndex_ == -1)
1071  {
1072  initialisePatch();
1073 
1074  initialiseEddyBox();
1075 
1076  initialiseEddies();
1077  }
1078 
1079 
1080  if (curTimeIndex_ != db().time().timeIndex())
1081  {
1082  if (debug)
1083  {
1084  label n = eddies_.size();
1085  Info<< "Number of eddies: " << returnReduce(n, sumOp<label>())
1086  << endl;
1087  }
1088 
1089  const scalar deltaT = db().time().deltaTValue();
1090 
1091  // Move eddies using mean velocity
1092  convectEddies(deltaT);
1093 
1094  // Set velocity
1095  vectorField& U = *this;
1096  //U = UMean_;
1097  U = U_;
1098 
1099  const pointField& Cf = patch().Cf();
1100 
1101  // Apply second part of normalisation coefficient
1102  // Note: factor of 2 required to match reference stresses?
1103  const scalar FACTOR = 2;
1104  const scalar c = FACTOR*Foam::sqrt(10*v0_)/Foam::sqrt(scalar(nEddy_));
1105 
1106  // In parallel, need to collect all eddies that will interact with
1107  // local faces
1108 
1109  if (singleProc_ || !Pstream::parRun())
1110  {
1111  forAll(U, faceI)
1112  {
1113  U[faceI] += c*uDashEddy(eddies_, Cf[faceI]);
1114  }
1115  }
1116  else
1117  {
1118  // Process local eddy contributions
1119  forAll(U, faceI)
1120  {
1121  U[faceI] += c*uDashEddy(eddies_, Cf[faceI]);
1122  }
1123 
1124  // Add contributions from overlapping eddies
1125  List<List<eddy>> overlappingEddies(Pstream::nProcs());
1126  calcOverlappingProcEddies(overlappingEddies);
1127 
1128  forAll(overlappingEddies, procI)
1129  {
1130  const List<eddy>& eddies = overlappingEddies[procI];
1131 
1132  if (eddies.size())
1133  {
1134  //Pout<< "Applying " << eddies.size()
1135  // << " eddies from processor " << procI << endl;
1136 
1137  forAll(U, faceI)
1138  {
1139  U[faceI] += c*uDashEddy(eddies, Cf[faceI]);
1140  }
1141  }
1142  }
1143  }
1144 
1145  // Re-scale to ensure correct flow rate
1146  scalar fCorr =
1147  gSum((UMean_ & patchNormal_)*patch().magSf())
1148  /gSum(U & -patch().Sf());
1149 
1150  U *= fCorr;
1151 
1152  if (debug)
1153  {
1154  Info<< "Patch:" << patch().patch().name()
1155  << " min/max(U):" << gMin(U) << ", " << gMax(U) << endl;
1156  }
1157 
1158  curTimeIndex_ = db().time().timeIndex();
1159 
1160  if (writeEddies_)
1161  {
1162  writeEddyOBJ();
1163  }
1164 
1165  if (debug && db().time().writeTime())
1166  {
1167  writeLumleyCoeffs();
1168  }
1169  }
1170 
1171  fixedValueFvPatchVectorField::updateCoeffs();
1172 }
1173 
1174 
1176 {
1178  writeEntry("value", os);
1179  os.writeEntry("delta", delta_);
1180  os.writeEntryIfDifferent<scalar>("d", 1.0, d_);
1181  os.writeEntryIfDifferent<scalar>("kappa", 0.41, kappa_);
1182  os.writeEntryIfDifferent<scalar>("perturb", 1e-5, perturb_);
1183  os.writeEntryIfDifferent<label>("nCellPerEddy", 5, nCellPerEddy_);
1184  os.writeEntryIfDifferent("writeEddies", false, writeEddies_);
1185 
1186  if (!interpolateR_)
1187  {
1188  R_.writeEntry("R", os);
1189  }
1190 
1191  if (!interpolateL_)
1192  {
1193  L_.writeEntry("L", os);
1194  }
1195 
1196  if (!interpolateU_)
1197  {
1198  U_.writeEntry("U", os);
1199  }
1200 
1201  if (!mapMethod_.empty())
1202  {
1204  (
1205  "mapMethod",
1206  "nearestCell",
1207  mapMethod_
1208  );
1209  }
1210 }
1211 
1212 
1213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1214 
1215 namespace Foam
1216 {
1218  (
1221  );
1222 }
1223 
1224 
1225 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField< vector >
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
volFields.H
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
Foam::Ostream::writeEntryIfDifferent
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:231
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::turbulentDFSEMInletFvPatchVectorField::checkStresses
static bool checkStresses(const symmTensorField &Rf)
Helper function to check that Reynold stresses are valid.
Definition: turbulentDFSEMInletFvPatchVectorField.C:957
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:316
Foam::SymmTensor< scalar >
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::UPstream::commsTypes::nonBlocking
Foam::eddy::debug
static int debug
Definition: eddy.H:143
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
s
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))
Definition: gmvOutputSpray.H:25
Foam::invariantIII
Cmpt invariantIII(const SymmTensor< Cmpt > &st)
Return the 3rd invariant of a symmetric tensor.
Definition: SymmTensorI.H:431
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
momentOfInertia.H
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::Pstream::scatterList
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
Definition: gatherScatterList.C:215
globalIndex.H
Foam::boundBox::invertedBox
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition: boundBox.H:86
Foam::UPstream::nProcs
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:426
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:414
Foam::pointToPointPlanarInterpolation
Interpolates between two sets of unstructured points using 2D Delaunay triangulation....
Definition: pointToPointPlanarInterpolation.H:53
Foam::turbulentDFSEMInletFvPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: turbulentDFSEMInletFvPatchVectorField.C:1175
Foam::fileHandler
const fileOperation & fileHandler()
Get current file handler.
Definition: fileOperation.C:1181
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:182
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
fvPatchFieldMapper.H
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
Foam::turbulentDFSEMInletFvPatchVectorField::turbulentDFSEMInletFvPatchVectorField
turbulentDFSEMInletFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: turbulentDFSEMInletFvPatchVectorField.C:727
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::sumOp
Definition: ops.H:213
Foam::fixedValueFvPatchField
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Definition: fixedValueFvPatchField.H:80
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
R
#define R(A, B, C, D, E, F, K, M)
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:55
turbulentDFSEMInletFvPatchVectorField.H
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Field< symmTensor >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Foam::symmTensor
SymmTensor< scalar > symmTensor
SymmTensor of scalars.
Definition: symmTensor.H:50
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::fvPatchField::rmap
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:303
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::invariantII
Cmpt invariantII(const SymmTensor< Cmpt > &st)
Return the 2nd invariant of a symmetric tensor.
Definition: SymmTensorI.H:419
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:444
Foam::turbulentDFSEMInletFvPatchVectorField::rmap
virtual void rmap(const fvPatchVectorField &ptf, const labelList &addr)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: turbulentDFSEMInletFvPatchVectorField.C:1042
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:290
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
U
U
Definition: pEqn.H:72
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:491
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:52
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:45
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Fstream.H
Input/output from file streams.
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::max
static const Vector< scalar > max
Definition: VectorSpace.H:117
f
labelList f(nPoints)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:62
Foam::List< label >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::token::TAB
Tab [isspace].
Definition: token.H:113
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::turbulentDFSEMInletFvPatchVectorField
Velocity boundary condition including synthesised eddies for use with LES and DES turbulent flows.
Definition: turbulentDFSEMInletFvPatchVectorField.H:179
points
const pointField & points
Definition: gmvOutputHeader.H:1
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
path
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
Foam::search
fileName search(const word &file, const fileName &directory)
Recursively search the given directory for the file.
Definition: fileName.C:576
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:219
Foam::turbulentDFSEMInletFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: turbulentDFSEMInletFvPatchVectorField.C:1063
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
forAllReverse
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: stdFoam.H:303
Foam::turbulentDFSEMInletFvPatchVectorField::autoMap
virtual void autoMap(const fvPatchFieldMapper &m)
Map (and resize as needed) from self given a mapping object.
Definition: turbulentDFSEMInletFvPatchVectorField.C:1025
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:51
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:432
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:122
triFace
face triFace(3)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::fvPatchField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: fvPatchField.C:240
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::triPointRef
triangle< point, const point & > triPointRef
Definition: triangle.H:78
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:177