PDRblockBlockMesh.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) 2020-2021 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "PDRblock.H"
29 #include "ListOps.H"
30 #include "gradingDescriptors.H"
31 #include "objectRegistry.H"
32 #include "Time.H"
33 #include "IOdictionary.H"
34 #include "Fstream.H"
35 #include "OTstream.H"
36 #include "edgeHashes.H"
37 
38 #include "cellModel.H"
39 #include "blockMesh.H"
40 #include "polyMesh.H"
41 #include "searchableSphere.H"
42 
43 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 // Name for the projection geometry
49 static const word projKeyword("project");
50 
51 // Name for the projection geometry
52 static const word projGeomName("sphere");
53 
54 
55 //- Calculate geometric ratio from relative ratio
56 inline scalar relativeToGeometricRatio
57 (
58  const scalar expRatio,
59  const label nDiv
60 )
61 {
62  return nDiv > 1 ? pow(expRatio, (nDiv - 1)) : 1.0;
63 }
64 
65 
66 // Output space-separated flat list. No size prefix.
67 template<class T>
68 static Ostream& outputFlatList(Ostream& os, const UList<T>& list)
69 {
71  label i = 0;
72  for (const label val : list)
73  {
74  if (i++) os << token::SPACE;
75  os << val;
76  }
78 
79  return os;
80 }
81 
82 
83 // Begin indent list
84 static inline Ostream& begIndentList(Ostream& os)
85 {
87  return os;
88 }
89 
90 // End indent list
91 static inline Ostream& endIndentList(Ostream& os)
92 {
94  return os;
95 }
96 
97 
98 // Output list contents (newline separated) indented.
99 template<class T>
100 static Ostream& outputIndent(Ostream& os, const UList<T>& list)
101 {
102  for (const T& val : list)
103  {
104  os << indent << val << nl;
105  }
106  return os;
107 }
108 
109 
110 //
111 static Ostream& serializeHex
112 (
113  Ostream& os,
114  const labelUList& hexVerts,
115  const labelVector& hexCount,
116  const Vector<gradingDescriptors> hexGrade,
117  const word& zoneName = word::null
118 )
119 {
121  outputFlatList(os, hexVerts);
122 
123  if (!zoneName.empty())
124  {
125  os << token::SPACE << zoneName;
126  }
127 
128  os << token::SPACE << hexCount << nl
129  << indent << word("edgeGrading") << nl;
130 
131  begIndentList(os);
132 
133  // Grading (x/y/z)
134  for (const gradingDescriptors& gds : hexGrade)
135  {
136  begIndentList(os);
137  outputIndent(os, gds);
138  endIndentList(os) << nl;
139  }
140 
141  endIndentList(os) << nl;
142  return os;
143 }
144 
145 
146 // Generate list with entries:
147 //
148 // project (x y z) (geometry)
149 //
150 static Ostream& serializeProjectPoints
151 (
152  Ostream& os,
153  const UList<point>& list
154 )
155 {
156  for (const point& p : list)
157  {
159  << p
160  << token::SPACE
162  }
163 
164  return os;
165 }
166 
167 
168 // Generate entry:
169 //
170 // project (beg end) (geometry)
171 //
172 static Ostream& serializeProjectEdge
173 (
174  Ostream& os,
175  const edge& e
176 )
177 {
179 
180  if (e.sorted())
181  {
182  os << e.first() << token::SPACE << e.second();
183  }
184  else
185  {
186  os << e.second() << token::SPACE << e.first();
187  }
188 
189  os << token::SPACE
191 
192  return os;
193 }
194 
195 
196 // Generate entry:
197 //
198 // (0 1 2 ..)
199 //
200 static Ostream& serializeFace
201 (
202  Ostream& os,
203  const face& list
204 )
205 {
206  os << indent;
207  outputFlatList(os, list);
208  os << nl;
209  return os;
210 }
211 
212 
213 // Generate entry:
214 //
215 // project (0 1 2 ..) geometry
216 //
217 static Ostream& serializeProjectFace
218 (
219  Ostream& os,
220  const face& list
221 )
222 {
224  outputFlatList(os, list);
225  os << token::SPACE << projGeomName << nl;
226 
227  return os;
228 }
229 
230 } // namespace Foam
231 
232 
233 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
234 
236 (
237  Ostream& os,
238  const bool withHeader
239 ) const
240 {
241  if (withHeader)
242  {
243  // Use dummy time for fake objectRegistry
244  autoPtr<Time> dummyTimePtr(Time::New());
245 
246  IOdictionary iodict
247  (
248  IOobject
249  (
250  "blockMeshDict",
251  dummyTimePtr->system(),
252  *dummyTimePtr,
253  IOobject::NO_READ,
254  IOobject::NO_WRITE,
255  false // no register
256  )
257  );
258 
259  iodict.writeHeader(os);
260  }
261 
262  const cellModel& hex = cellModel::ref(cellModel::HEX);
263 
264  // The mesh topology will normally be an O-grid with a central (inner)
265  // block and 6 outer blocks.
266  // The inner block is described by (0 1 2 3 4 5 6 7).
267  // The additional points for the outer region: (8 9 19 11 12 13 14 15)
268 
269  // The outer blocks will be addressed according to their
270  // placement w.r.t. the inner block, and defined such that they retain
271  // the same global orientation as the inner block.
272  // For example, face 0 of all blocks will be on the logical x-min side
273  // of the mesh.
274 
275  // List of hex vertices, ordered with outer blocks first to allow
276  // direct addressing by their logical position.
278  hexVerts
279  ({
280  {8, 0, 3, 11, 12, 4, 7, 15}, // x-min block (face 0)
281  {1, 9, 10, 2, 5, 13, 14, 6}, // x-max block (face 1)
282  {8, 9, 1, 0, 12, 13, 5, 4}, // y-min block (face 2)
283  {3, 2, 10, 11, 7, 6, 14, 15}, // y-max block (face 3)
284  {8, 9, 10, 11, 0, 1, 2, 3}, // z-min block (face 4)
285  {4, 5, 6, 7, 12, 13, 14, 15}, // z-max block (face 5)
286  {0, 1, 2, 3, 4, 5, 6, 7}, // Inner box description
287  {8, 9, 10, 11, 12, 13, 14, 15}, // Outer box description
288  });
289 
290  // The face or the logical block index (for the O-grid)
291  enum faceIndex
292  {
293  X_Min = 0,
294  X_Max = 1,
295  Y_Min = 2,
296  Y_Max = 3,
297  Z_Min = 4,
298  Z_Max = 5,
299  Inner_Block = 6, // Inner block description
300  Outer_Block = 7, // Outer bounding box description
301  };
302 
303  // Lists of block/face for outside and ground faces
304  DynamicList<labelPair> outerFaces(8);
305  DynamicList<labelPair> groundFaces(8);
306 
307  // We handle a few fixed topology configurations
308 
309  enum blockTopologyType
310  {
311  INNER_ONLY = 0, // No outer region
312  EXTENDED = 1, // Outer created by extending inner region
313  CLIP_BOTTOM = 5, // Outer O-grid on 5 sides
314  FULL_OUTER = 6 // Outer O-grid on 6 sides
315  };
316 
317 
318 
319  // Expansion ratios need conversion from relative to geometric
320  const bool useRelToGeom =
321  (expansionType::EXPAND_RATIO == outer_.expandType_);
322 
323 
324  // Physical dimensions
325 
326  Vector<gridControl> ctrl(control_);
327  boundBox innerCorners(bounds(ctrl.x(), ctrl.y(), ctrl.z()));
328 
329  boundBox outerCorners;
330 
331 
332  point radialCentre(innerCorners.centre());
333  vector radialSizes(0.5*innerCorners.span());
334 
335  blockTopologyType outerTopology = INNER_ONLY;
336 
337 
338  if (outer_.active())
339  {
340  outerTopology = FULL_OUTER;
341 
342  // Convert from relative size
343  radialSizes.x() *= outer_.relSize_.x();
344  radialSizes.y() *= outer_.relSize_.y();
345  radialSizes.z() *= min(outer_.relSize_.x(), outer_.relSize_.y());
346 
347  if (outer_.onGround())
348  {
349  outerTopology = CLIP_BOTTOM;
350  radialCentre.z() = innerCorners.min().z();
351  radialSizes.z() *= 2;
352  }
353 
354  // Corners for a box
355  outerCorners.min() = radialCentre - radialSizes;
356  outerCorners.max() = radialCentre + radialSizes;
357 
358  if (outer_.onGround())
359  {
360  outerCorners.min().z() = innerCorners.min().z();
361  }
362 
363 
364  if (outer_.isSphere())
365  {
366  // For spheroid projection, don't trust that blockMesh does it
367  // properly. Give some reasonable estimates of the corners.
368 
369  // Use dummy Time for objectRegistry
370  autoPtr<Time> dummyTimePtr(Time::New());
371 
372  const IOobject io
373  (
374  "sphere",
375  *dummyTimePtr,
376  IOobject::NO_READ,
377  IOobject::NO_WRITE,
378  false // do not register
379  );
380 
381  searchableSphere sphere(io, radialCentre, radialSizes);
382 
383  pointField queries(2);
384  queries[0] = outerCorners.min();
385  queries[1] = outerCorners.max();
386 
387  List<pointIndexHit> hits;
388  sphere.findNearest
389  (
390  queries,
391  scalarField(2, GREAT),
392  hits
393  );
394 
395  outerCorners.min() = hits[0].hitPoint();
396  outerCorners.max() = hits[1].hitPoint();
397  }
398  else if (outerControl::OUTER_EXTEND == outer_.type_)
399  {
400  outerTopology = EXTENDED;
401 
402  // Extend the inner block
403  label outerCount;
404  scalar expRatio;
405 
406  outerCount = outer_.nCells_.x();
407  expRatio = outer_.expansion_.x();
408  if (useRelToGeom)
409  {
410  expRatio = relativeToGeometricRatio(expRatio, outerCount);
411  }
412 
413  ctrl.x().prepend(outerCorners.min().x(), outerCount, -expRatio);
414  ctrl.x().append(outerCorners.max().x(), outerCount, expRatio);
415 
416 
417  outerCount = outer_.nCells_.y();
418  expRatio = outer_.expansion_.y();
419  if (useRelToGeom)
420  {
421  expRatio = relativeToGeometricRatio(expRatio, outerCount);
422  }
423 
424  ctrl.y().prepend(outerCorners.min().y(), outerCount, -expRatio);
425  ctrl.y().append(outerCorners.max().y(), outerCount, expRatio);
426 
427  outerCount = max(outer_.nCells_.x(), outer_.nCells_.y());
428  expRatio = min(outer_.expansion_.x(), outer_.expansion_.y());
429  if (useRelToGeom)
430  {
431  expRatio = relativeToGeometricRatio(expRatio, outerCount);
432  }
433 
434  if (!outer_.onGround())
435  {
436  ctrl.z().prepend(outerCorners.min().z(), outerCount, -expRatio);
437  }
438  ctrl.z().append(outerCorners.max().z(), outerCount, expRatio);
439 
440  // Update corners
441  innerCorners = bounds(ctrl.x(), ctrl.y(), ctrl.z());
442  outerCorners = innerCorners;
443  }
444  }
445 
446 
447  const Vector<gradingDescriptors> innerGrading(grading(ctrl));
448  const labelVector innerCount(sizes(ctrl));
449 
450  labelVector hexCount;
452 
453 
454  const label radialCount = outer_.nCells_.x();
455  scalar expRatio = outer_.expansion_.x();
456 
457  if (useRelToGeom)
458  {
459  expRatio = relativeToGeometricRatio(expRatio, radialCount);
460  }
461 
462  const gradingDescriptors radialInward
463  (
464  gradingDescriptor{-expRatio}
465  );
466 
467  const gradingDescriptors radialOutward
468  (
469  gradingDescriptor{expRatio}
470  );
471 
472 
473  if (EXTENDED == outerTopology)
474  {
475  // The inner block is extended to become the outer faces
476  outerFaces.append
477  ({
478  labelPair(Inner_Block, X_Min),
479  labelPair(Inner_Block, X_Max),
480  labelPair(Inner_Block, Y_Min),
481  labelPair(Inner_Block, Y_Max),
482  labelPair(Inner_Block, Z_Max)
483  });
484 
485  // The ground faces vs outside faces
486  if (outer_.onGround())
487  {
488  groundFaces.append
489  (
490  labelPair(Inner_Block, Z_Min)
491  );
492  }
493  else
494  {
495  outerFaces.append
496  (
497  labelPair(Inner_Block, Z_Min)
498  );
499  }
500  }
501  else if (CLIP_BOTTOM == outerTopology || FULL_OUTER == outerTopology)
502  {
503  // The outside faces
504  outerFaces.append
505  ({
506  labelPair(X_Min, X_Min),
507  labelPair(X_Max, X_Max),
508  labelPair(Y_Min, Y_Min),
509  labelPair(Y_Max, Y_Max),
510  labelPair(Z_Max, Z_Max)
511  });
512 
513  // The ground faces
514  if (CLIP_BOTTOM == outerTopology)
515  {
516  groundFaces.append
517  ({
518  labelPair(X_Min, Z_Min),
519  labelPair(X_Max, Z_Min),
520  labelPair(Y_Min, Z_Min),
521  labelPair(Y_Max, Z_Min),
522  // Note: {Z_Min, Z_Min} will not exist
523  labelPair(Inner_Block, Z_Min)
524  });
525  }
526  else
527  {
528  outerFaces.append
529  (
530  labelPair(Z_Min, Z_Min)
531  );
532  }
533  }
534 
535 
536  if (outer_.isSphere())
537  {
538  os.beginBlock("geometry");
539  {
540  os.beginBlock(projGeomName);
541  {
542  os.writeEntry("type", "sphere");
543  os.writeEntry("origin", radialCentre);
544  os.writeEntry("radius", radialSizes);
545  }
546  os.endBlock();
547  }
548  os.endBlock();
549  }
550 
551  // vertices
552  {
553  os << nl << word("vertices") << nl;
554  begIndentList(os);
555 
556  pointField corners(innerCorners.points());
557 
558  // inner
559  outputIndent(os, corners);
560 
561  // outer
562  if (CLIP_BOTTOM == outerTopology || FULL_OUTER == outerTopology)
563  {
564  corners = outerCorners.points();
565 
566  if (outer_.isSphere())
567  {
568  serializeProjectPoints(os, corners);
569  }
570  else
571  {
572  outputIndent(os, corners);
573  }
574  }
575 
576  endIndentList(os) << token::END_STATEMENT << nl;
577  }
578 
579 
580  // blocks
581  {
582  word innerZoneName = "inner";
583  if (INNER_ONLY == outerTopology || EXTENDED == outerTopology)
584  {
585  innerZoneName.clear();
586  }
587 
588  os << nl << word("blocks") << nl;
589  begIndentList(os);
590 
591  // Inner block
592  hexCount = innerCount;
593 
595  (
596  os,
597  hexVerts[Inner_Block],
598  hexCount,
599  innerGrading,
600  innerZoneName
601  );
602 
603  // outer
604  if (CLIP_BOTTOM == outerTopology || FULL_OUTER == outerTopology)
605  {
606  // Radial direction = X
607  hexCount = innerCount;
608  hexGrade = innerGrading;
609 
610  hexCount.x() = radialCount;
611 
612  // Face 0: x-min
613  {
614  hexGrade.x() = radialInward;
615  serializeHex(os, hexVerts[X_Min], hexCount, hexGrade);
616  }
617  // Face 1: x-max
618  {
619  hexGrade.x() = radialOutward;
620  serializeHex(os, hexVerts[X_Max], hexCount, hexGrade);
621  }
622 
623 
624  // Radial direction = Y
625  hexCount = innerCount;
626  hexGrade = innerGrading;
627 
628  hexCount.y() = radialCount;
629 
630  // Face 2: y-min
631  {
632  hexGrade.y() = radialInward;
633  serializeHex(os, hexVerts[Y_Min], hexCount, hexGrade);
634  }
635  // Face 3: y-max
636  {
637  hexGrade.y() = radialOutward;
638  serializeHex(os, hexVerts[Y_Max], hexCount, hexGrade);
639  }
640 
641 
642  // Radial direction = Z
643  hexCount = innerCount;
644  hexGrade = innerGrading;
645 
646  hexCount.z() = radialCount;
647 
648  // Face 4: z-min
649  if (!outer_.onGround())
650  {
651  hexGrade.z() = radialInward;
652  serializeHex(os, hexVerts[Z_Min], hexCount, hexGrade);
653  }
654  // Face 5: z-max
655  {
656  hexGrade.z() = radialOutward;
657  serializeHex(os, hexVerts[Z_Max], hexCount, hexGrade);
658  }
659  }
660 
661  endIndentList(os) << token::END_STATEMENT << nl;
662  }
663 
664 
665  // edges
666  {
667  os << nl << word("edges") << nl;
668  begIndentList(os);
669 
670  if (outer_.isSphere() && outerFaces.size())
671  {
672  // Edges for the outer face of the block
673  edgeHashSet projEdges(32);
674 
675  for (const labelPair& pr : outerFaces)
676  {
677  projEdges.insert
678  (
679  hex.face(pr.second(), hexVerts[pr.first()]).edges()
680  );
681  }
682 
683  for (const edge& e : projEdges.sortedToc())
684  {
686  }
687  }
688 
689  endIndentList(os) << token::END_STATEMENT << nl;
690  }
691 
692 
693  // faces
694  {
695  os << nl << word("faces") << nl;
696  begIndentList(os);
697 
698  if (outer_.isSphere() && outerFaces.size())
699  {
700  for (const labelPair& pr : outerFaces)
701  {
703  (
704  os,
705  hex.face(pr.second(), hexVerts[pr.first()])
706  );
707  }
708  }
709 
710  endIndentList(os) << token::END_STATEMENT << nl;
711  }
712 
713 
714  // boundary
715  {
716  os << nl << word("boundary") << nl;
717  begIndentList(os);
718 
719  // outer
720  {
721  os.beginBlock("outer");
722  os.writeEntry("type", word("patch"));
723 
724  os << indent << word("faces") << nl;
725  begIndentList(os);
726 
727  for (const labelPair& pr : outerFaces)
728  {
730  (
731  os,
732  hex.face(pr.second(), hexVerts[pr.first()])
733  );
734  }
735 
736  endIndentList(os) << token::END_STATEMENT << nl;
737 
738  os.endBlock();
739  }
740 
741  if (outer_.onGround())
742  {
743  os.beginBlock("ground");
744  os.writeEntry("type", word("wall"));
745 
746  os << indent << word("faces") << nl;
747  begIndentList(os);
748 
749  for (const labelPair& pr : groundFaces)
750  {
752  (
753  os,
754  hex.face(pr.second(), hexVerts[pr.first()])
755  );
756  }
757 
758  endIndentList(os) << token::END_STATEMENT << nl;
759  os.endBlock();
760  }
761 
762  endIndentList(os) << token::END_STATEMENT << nl;
763  }
764 
765 
766  if (withHeader)
767  {
768  IOobject::writeEndDivider(os);
769  }
770 
771  return os;
772 }
773 
774 
776 {
777  OTstream os;
778  blockMeshDict(os);
779 
780  ITstream is;
781  is.transfer(os.tokens());
782 
783  return dictionary(is);
784 }
785 
786 
788 {
789  IOdictionary iodict
790  (
791  IOobject
792  (
793  io.name(),
794  io.db().time().system(),
795  io.local(),
796  io.db(),
799  false // no register
800  )
801  );
802 
803  OFstream os(iodict.objectPath());
804 
805  Info<< nl
806  << "Generate blockMeshDict: "
807  << iodict.db().time().relativePath(os.name()) << endl;
808 
809  // Set precision for points to 10
811 
812  iodict.writeHeader(os);
813 
814  // Just like writeData, but without copying beforehand
815  this->blockMeshDict(os);
816 
817  iodict.writeEndDivider(os);
818 }
819 
820 
822 Foam::PDRblock::createBlockMesh(const IOobject& io) const
823 {
824  IOdictionary iodict
825  (
826  IOobject
827  (
828  "blockMeshDict.PDRblockMesh",
829  io.db().time().system(),
830  io.local(),
831  io.db(),
834  false // no register
835  ),
836  blockMeshDict()
837  );
838 
839  return autoPtr<blockMesh>::New(iodict);
840 }
841 
842 
844 Foam::PDRblock::meshBlockMesh(const IOobject& io) const
845 {
846  const bool oldVerbose = blockMesh::verboseOutput;
847  blockMesh::verboseOutput = false;
848 
849  autoPtr<polyMesh> meshPtr(createBlockMesh(io)->mesh(io));
850 
851  blockMesh::verboseOutput = oldVerbose;
852 
853  // This is a bit ugly.
854  // For extend, we still wish to have an 'inner' cellZone,
855  // but we meshed the entirety.
856 
857  if
858  (
859  outerControl::OUTER_EXTEND == outer_.type_
860  && meshPtr->cellZones().empty()
861  )
862  {
863  const boundBox innerBox
864  (
865  bounds(control_.x(), control_.y(), control_.z())
866  );
867 
868  const label nZoneCellsMax =
869  (
870  control_.x().nCells()
871  * control_.y().nCells()
872  * control_.z().nCells()
873  );
874 
875 
876  polyMesh& pmesh = *meshPtr;
877 
878  List<cellZone*> cz(1);
879  cz[0] = new cellZone
880  (
881  "inner",
882  labelList(nZoneCellsMax),
883  0, // zonei
884  pmesh.cellZones()
885  );
886 
887  cellZone& innerZone = *(cz[0]);
888 
889  const vectorField& cc = pmesh.cellCentres();
890 
891  label nZoneCells = 0;
892 
893  for
894  (
895  label celli = 0;
896  celli < cc.size() && nZoneCells < nZoneCellsMax;
897  ++celli
898  )
899  {
900  if (innerBox.contains(cc[celli]))
901  {
902  innerZone[nZoneCells] = celli;
903  ++nZoneCells;
904  }
905  }
906 
907  innerZone.resize(nZoneCells);
908 
909  pmesh.pointZones().clear();
910  pmesh.faceZones().clear();
911  pmesh.cellZones().clear();
912  pmesh.addZones(List<pointZone*>(), List<faceZone*>(), cz);
913  }
914 
915  return meshPtr;
916 }
917 
918 
919 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
hex
const cellModel & hex
Definition: createBlockMesh.H:1
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
meshPtr
Foam::autoPtr< Foam::fvMesh > meshPtr(nullptr)
Foam::projKeyword
static const word projKeyword("project")
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::OFstream::name
virtual const fileName & name() const
Read/write access to the name of the stream.
Definition: OSstream.H:107
searchableSphere.H
Foam::cellModel::modelNames
static const Enum< modelType > modelNames
Names of commonly used cellModels corresponding to modelType.
Definition: cellModel.H:92
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:55
Foam::cellModel::HEX
hex
Definition: cellModel.H:81
Foam::IOobject::writeEndDivider
static Ostream & writeEndDivider(Ostream &os)
Write the standard end file divider.
Definition: IOobjectWriteHeader.C:142
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
objectRegistry.H
Foam::serializeProjectEdge
static Ostream & serializeProjectEdge(Ostream &os, const edge &e)
Definition: PDRblockBlockMesh.C:173
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
ref
rDeltaT ref()
Foam::gradingDescriptors
List of gradingDescriptor for the sections of a block with additional IO functionality.
Definition: gradingDescriptors.H:58
blockMesh.H
polyMesh.H
Foam::HashSet< edge, Hash< edge > >
Foam::incrIndent
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:346
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
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::projGeomName
static const word projGeomName("sphere")
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
gradingDescriptors.H
Foam::labelPair
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:54
Foam::boundBox::span
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:127
sphere
Specialization of rigidBody to construct a sphere given the mass and radius.
Foam::endIndentList
static Ostream & endIndentList(Ostream &os)
Definition: PDRblockBlockMesh.C:91
Foam::TimePaths::relativePath
fileName relativePath(const fileName &input, const bool caseTag=false) const
Definition: TimePathsI.H:87
Foam::Field< vector >
Foam::ITstream
An input stream of tokens.
Definition: ITstream.H:52
cellModel.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::boundBox::centre
point centre() const
The centre (midpoint) of the bounding box.
Definition: boundBoxI.H:115
Foam::IOobject::local
const fileName & local() const noexcept
Definition: IOobjectI.H:208
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:492
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:511
Foam::begIndentList
static Ostream & begIndentList(Ostream &os)
Definition: PDRblockBlockMesh.C:84
Foam::List::transfer
void transfer(List< T > &list)
Definition: List.C:456
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
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
PDRblock.H
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::serializeHex
static Ostream & serializeHex(Ostream &os, const labelUList &hexVerts, const labelVector &hexCount, const Vector< gradingDescriptors > hexGrade, const word &zoneName=word::null)
Definition: PDRblockBlockMesh.C:112
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::relativeToGeometricRatio
scalar relativeToGeometricRatio(const scalar expRatio, const label nDiv)
Calculate geometric ratio from relative ratio.
Definition: PDRblock.C:65
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::decrIndent
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:353
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:339
Foam::serializeProjectFace
static Ostream & serializeProjectFace(Ostream &os, const face &list)
Definition: PDRblockBlockMesh.C:218
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::OSstream::precision
virtual int precision() const
Get precision of output field.
Definition: OSstream.C:326
Foam::IOobject::name
const word & name() const noexcept
Return name.
Definition: IOobjectI.H:65
Foam::IOstream::defaultPrecision
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:342
IOdictionary.H
Time.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
Foam::TimePaths::system
const word & system() const
Return system name.
Definition: TimePathsI.H:102
edgeHashes.H
Foam::gradingDescriptor
Handles the specification for grading within a section of a block.
Definition: gradingDescriptor.H:72
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::Pair< label >
Foam::OTstream
A simple output token stream that can be used to build token lists. Always UNCOMPRESSED.
Definition: OTstream.H:56
Fstream.H
Foam::Vector< label >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::FixedList
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: HashTable.H:104
Foam::token::SPACE
Space [isspace].
Definition: token.H:125
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:103
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
Foam::serializeProjectPoints
static Ostream & serializeProjectPoints(Ostream &os, const UList< point > &list)
Definition: PDRblockBlockMesh.C:151
Foam::word::null
static const word null
An empty word.
Definition: word.H:80
Foam::outputFlatList
static Ostream & outputFlatList(Ostream &os, const UList< T > &list)
Definition: PDRblockBlockMesh.C:68
Foam::boundBox
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
Foam::token::END_LIST
End list [isseparator].
Definition: token.H:156
Foam::PDRblock::writeBlockMeshDict
void writeBlockMeshDict(const IOobject &io) const
Write an equivalent blockMeshDict.
Definition: PDRblockBlockMesh.C:787
Foam::outputIndent
static Ostream & outputIndent(Ostream &os, const UList< T > &list)
Definition: PDRblockBlockMesh.C:100
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::IOobject::writeHeader
bool writeHeader(Ostream &os) const
Write header with current type()
Definition: IOobjectWriteHeader.C:277
ListOps.H
Various functions to operate on Lists.
Foam::cellModel
Maps a geometry to a set of cell primitives.
Definition: cellModel.H:72
Foam::serializeFace
static Ostream & serializeFace(Ostream &os, const face &list)
Definition: PDRblockBlockMesh.C:201
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::IOobject::objectPath
fileName objectPath() const
The complete path + object name.
Definition: IOobjectI.H:214
Foam::boundBox::points
tmp< pointField > points() const
Corner points in an order corresponding to a 'hex' cell.
Definition: boundBox.C:118
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::token::BEGIN_LIST
Begin list [isseparator].
Definition: token.H:155
OTstream.H
Foam::PDRblock::blockMeshDict
dictionary blockMeshDict() const
Content for an equivalent blockMeshDict.
Definition: PDRblockBlockMesh.C:775
Foam::objectRegistry::time
const Time & time() const noexcept
Return time registry.
Definition: objectRegistry.H:178
Foam::blockMesh::verboseOutput
static bool verboseOutput
The default verbosity (true)
Definition: blockMesh.H:339
Foam::IOobject::db
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition: IOobject.C:487
Foam::searchableSphere
Searching on general spheroid.
Definition: searchableSphere.H:92