PDRblock.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2019-2023 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 "emptyPolyPatch.H"
31 #include "gradingDescriptors.H"
32 
33 // Output when verbosity is enabled
34 #undef Log
35 #define Log if (verbose_) Info
36 
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 const Foam::Enum
41 <
43 >
45 ({
46  { expansionType::EXPAND_UNIFORM, "uniform" },
47  { expansionType::EXPAND_RATIO, "ratio" },
48  { expansionType::EXPAND_RELATIVE, "relative"}
49 });
50 
51 
52 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 
57 //- Calculate geometric expansion factor from expansion ratio
58 inline scalar calcGexp(const scalar expRatio, const label nDiv)
59 {
60  return nDiv > 1 ? pow(expRatio, 1.0/(nDiv - 1)) : 0.0;
61 }
62 
63 //- Calculate geometric ratio from relative ratio
64 inline scalar relativeToGeometricRatio
65 (
66  const scalar expRatio,
67  const label nDiv
68 )
69 {
70  return nDiv > 1 ? pow(expRatio, (nDiv - 1)) : 1.0;
71 }
72 
73 } // End namespace Foam
74 
75 
76 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
77 
78 bool Foam::PDRblock::checkMonotonic
79 (
80  const direction cmpt,
81  const UList<scalar>& pts
82 )
83 {
84  const label len = pts.size();
85 
86  if (!len)
87  {
88  return false;
89  }
90 
91  const scalar& minVal = pts[0];
92 
93  for (label i=1; i < len; ++i)
94  {
95  if (pts[i] <= minVal)
96  {
98  << "Points in " << vector::componentNames[cmpt]
99  << " direction do not increase monotonically" << nl
100  << flatOutput(pts) << nl << nl
101  << exit(FatalError);
102  }
103  }
104 
105  return true;
106 }
107 
108 
110 {
111  return NullObjectRef<PDRblock>();
112 }
113 
114 
115 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
116 
117 void Foam::PDRblock::addDefaultPatches()
118 {
119  // Default boundaries with patchi == shapeFacei
120  patches_.resize(6);
121  for (label patchi=0; patchi < 6; ++patchi)
122  {
123  boundaryEntry& bentry = patches_.emplace_set(patchi);
124 
125  bentry.name_ = "patch" + Foam::name(patchi);
126  bentry.type_ = "patch";
127  bentry.size_ = 0;
128  bentry.faces_ = labelList(one{}, patchi);
129  }
130 }
131 
132 
133 void Foam::PDRblock::adjustSizes()
134 {
135  // Adjust i-j-k addressing
136  sizes().x() = grid_.x().nCells();
137  sizes().y() = grid_.y().nCells();
138  sizes().z() = grid_.z().nCells();
139 
140  if (sizes().x() <= 0 || sizes().y() <= 0 || sizes().z() <= 0)
141  {
142  // Sanity check. Silently disallow bad sizing
143  ijkMesh::clear();
144 
145  grid_.x().clear();
146  grid_.y().clear();
147  grid_.z().clear();
148 
149  bounds_.reset();
150  edgeLimits_.min() = 0;
151  edgeLimits_.max() = 0;
152  return;
153  }
154 
155  // Adjust boundBox
156  bounds_ = bounds(grid_.x(), grid_.y(), grid_.z());
157 
158  // Min/max edge lengths
159  edgeLimits_.clear();
160 
161  edgeLimits_.add(grid_.x().edgeLimits());
162  edgeLimits_.add(grid_.y().edgeLimits());
163  edgeLimits_.add(grid_.z().edgeLimits());
164 }
165 
166 
167 void Foam::PDRblock::readGridControl
168 (
169  const direction cmpt,
170  const dictionary& dict,
171  const scalar scaleFactor,
172  expansionType expandType
173 )
174 {
175  gridControl& ctrl = control_[cmpt];
176 
177  // Begin/end nodes for each segment
178  scalarList& knots = static_cast<scalarList&>(ctrl);
179 
180  // The number of division per segment
181  labelList& divisions = ctrl.divisions_;
182 
183  // The expansion ratio per segment
184  scalarList& expansion = ctrl.expansion_;
185 
186  expansion.clear(); // expansion is optional
187 
188  Log << "Reading grid control for "
189  << vector::componentNames[cmpt] << " direction" << nl;
190 
191  // Points
192  // ~~~~~~
193 
194  dict.readEntry("points", knots);
195 
196  const label nSegments = (knots.size()-1);
197 
198  if (nSegments < 1)
199  {
201  << "Must define at least two control points:"
202  << " in " << vector::componentNames[cmpt]
203  << " direction" << nl
204  << exit(FatalIOError);
205  }
206 
207  // Apply scaling
208  if (scaleFactor > 0)
209  {
210  for (scalar& pt : knots)
211  {
212  pt *= scaleFactor;
213  }
214  }
215 
216  // Do points increase monotonically?
217  checkMonotonic(cmpt, knots);
218 
219  Log << " points : " << flatOutput(knots) << nl;
220 
221 
222  // Divisions
223  // ~~~~~~~~~
224 
225  dict.readEntry("nCells", divisions);
226 
227  label nTotalDivs = 0;
228  for (const label ndiv : divisions)
229  {
230  nTotalDivs += ndiv;
231 
232  if (ndiv < 1)
233  {
235  << "Negative or zero divisions:"
236  << " in " << vector::componentNames[cmpt]
237  << " direction" << nl
238  << exit(FatalIOError);
239  }
240  }
241 
242  if (divisions.size() != nSegments)
243  {
245  << "Mismatch in number of divisions and number of points:"
246  << " in " << vector::componentNames[cmpt]
247  << " direction" << nl
248  << exit(FatalIOError);
249  }
250  else if (!nTotalDivs)
251  {
253  << "No grid divisions at all:"
254  << " in " << vector::componentNames[cmpt]
255  << " direction" << nl
256  << exit(FatalIOError);
257  }
258 
259  Log << " nCells : " << flatOutput(divisions) << nl;
260 
261 
262  // Expansion ratios
263  // ~~~~~~~~~~~~~~~~
264 
265  dict.readIfPresent("ratios", expansion);
266 
267  // Correct expansion ratios - negative is the same as inverse.
268  for (scalar& expRatio : expansion)
269  {
270  if (expRatio < 0)
271  {
272  expRatio = 1.0/(-expRatio);
273  }
274  }
275 
276  if (expansion.size() && divisions.size() != expansion.size())
277  {
279  << "Mismatch in number of expansion ratios and divisions:"
280  << " in " << vector::componentNames[cmpt]
281  << " direction" << nl
282  << exit(FatalIOError);
283  }
284 
285  if (expansion.empty())
286  {
287  expansion.resize(nSegments, scalar(1));
288 
289  if (expandType != expansionType::EXPAND_UNIFORM)
290  {
291  expandType = expansionType::EXPAND_UNIFORM;
292 
293  Log << "Warning: no 'ratios', use uniform spacing" << nl;
294  }
295  }
296  else
297  {
298  switch (expandType)
299  {
300  case expansionType::EXPAND_UNIFORM:
301  {
302  expansion = scalar(1);
303  break;
304  }
305 
306  case expansionType::EXPAND_RATIO:
307  {
308  Log << " ratios : " << flatOutput(expansion) << nl;
309  break;
310  }
311 
312  case expansionType::EXPAND_RELATIVE:
313  {
314  Log << " relative : " << flatOutput(expansion) << nl;
315 
316  auto divIter = divisions.cbegin();
317 
318  for (scalar& expRatio : expansion)
319  {
320  expRatio = relativeToGeometricRatio(expRatio, *divIter);
321  ++divIter;
322  }
323 
324  Log << " ratios : " << flatOutput(expansion) << nl;
325  break;
326  }
327  }
328  }
329 
330 
331  // Define the grid points
332 
333  scalarList& gridPoint = grid_[cmpt];
334  gridPoint.resize(nTotalDivs+1);
335 
336  label start = 0;
337 
338  for (label segmenti=0; segmenti < nSegments; ++segmenti)
339  {
340  const label nDiv = divisions[segmenti];
341  const scalar expRatio = expansion[segmenti];
342 
343  SubList<scalar> subPoint(gridPoint, nDiv+1, start);
344  start += nDiv;
345 
346  subPoint[0] = knots[segmenti];
347  subPoint[nDiv] = knots[segmenti+1];
348 
349  const scalar dist = (subPoint.last() - subPoint.first());
350 
351  if
352  (
353  expandType == expansionType::EXPAND_UNIFORM
354  || equal(expRatio, 1)
355  )
356  {
357  for (label i=1; i < nDiv; ++i)
358  {
359  subPoint[i] = (subPoint[0] + (dist * i)/nDiv);
360  }
361  }
362  else
363  {
364  // Geometric expansion factor from the expansion ratio
365  const scalar expFact = calcGexp(expRatio, nDiv);
366 
367  for (label i=1; i < nDiv; ++i)
368  {
369  subPoint[i] =
370  (
371  subPoint[0]
372  + dist * (1.0 - pow(expFact, i))/(1.0 - pow(expFact, nDiv))
373  );
374  }
375  }
376  }
377 }
378 
379 
380 void Foam::PDRblock::readBoundary(const dictionary& dict)
381 {
382  Log << "Reading boundary entries" << nl;
383 
384  PtrList<entry> patchEntries;
385 
386  if (dict.found("boundary"))
387  {
388  PtrList<entry> newEntries(dict.lookup("boundary"));
389  patchEntries.transfer(newEntries);
390  }
391 
392  patches_.clear();
393  patches_.resize(patchEntries.size());
394 
395  // Hex cell has 6 sides
396  const labelRange validFaceId(0, 6);
397 
398  // Track which sides have been handled
399  labelHashSet handled;
400 
401  forAll(patchEntries, patchi)
402  {
403  if (!patchEntries.set(patchi))
404  {
405  continue;
406  }
407 
408  const entry& e = patchEntries[patchi];
409 
410  if (!e.isDict())
411  {
413  << "Entry " << e
414  << " in boundary section is not a valid dictionary."
415  << exit(FatalIOError);
416  }
417 
418  const dictionary& patchDict = e.dict();
419 
420  const word& patchName = e.keyword();
421 
422  const word patchType(patchDict.get<word>("type"));
423 
424  labelList faceLabels(patchDict.get<labelList>("faces"));
425 
426  labelHashSet patchFaces(faceLabels);
427 
428  if (faceLabels.size() != patchFaces.size())
429  {
431  << "Patch: " << patchName
432  << ": Ignoring duplicate face ids" << nl;
433  }
434 
435  label nErased = patchFaces.filterKeys(validFaceId);
436  if (nErased)
437  {
439  << "Patch: " << patchName << ": Ignoring " << nErased
440  << " faces with invalid identifiers" << nl;
441  }
442 
443  nErased = patchFaces.erase(handled);
444  if (nErased)
445  {
447  << "Patch: " << patchName << ": Ignoring " << nErased
448  << " faces, which were already handled" << nl;
449  }
450 
451  // Mark these faces as having been handled
452  handled += patchFaces;
453 
454 
455  // Save information for later access during mesh creation.
456 
457  boundaryEntry& bentry = patches_.emplace_set(patchi);
458 
459  bentry.name_ = patchName;
460  bentry.type_ = patchType;
461  bentry.size_ = 0;
462  bentry.faces_ = patchFaces.sortedToc();
463 
464  // Count patch faces
465  for (const label shapeFacei : bentry.faces_)
466  {
467  bentry.size_ += nBoundaryFaces(shapeFacei);
468  }
469  }
470 
471 
472  // Deal with unhandled faces
473 
474  labelHashSet missed(identity(6));
475  missed.erase(handled);
476 
477  if (missed.size())
478  {
479  boundaryEntry& bentry = patches_.emplace_back();
480 
481  bentry.name_ = "defaultFaces";
482  bentry.type_ = emptyPolyPatch::typeName;
483  bentry.size_ = 0;
484  bentry.faces_ = missed.sortedToc();
485 
486  // Count patch faces
487  for (const label shapeFacei : bentry.faces_)
488  {
489  bentry.size_ += nBoundaryFaces(shapeFacei);
490  }
491 
492 
493  // Use name/type from "defaultPatch" entry if available
494  // - be generous with error handling
495 
496  const dictionary* defaultEntry = dict.findDict("defaultPatch");
497  if (defaultEntry)
498  {
499  defaultEntry->readIfPresent("name", bentry.name_);
500  defaultEntry->readIfPresent("type", bentry.type_);
501  }
502  }
503 
504  if (verbose_)
505  {
506  label patchi = 0;
507  for (const boundaryEntry& bentry : patches_)
508  {
509  Info<< " patch: " << patchi
510  << " (size: " << bentry.size_
511  << " type: " << bentry.type_
512  << ") name: " << bentry.name_
513  << " faces: " << flatOutput(bentry.faces_) << nl;
514 
515  ++patchi;
516  }
517  }
518 }
519 
521 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
522 
524 :
525  PDRblock(dictionary::null, false)
526 {}
528 
530 (
531  const UList<scalar>& xgrid,
532  const UList<scalar>& ygrid,
533  const UList<scalar>& zgrid
534 )
535 :
536  PDRblock(dictionary::null, false)
537 {
538  addDefaultPatches();
539  reset(xgrid, ygrid, zgrid);
540 }
541 
542 
543 Foam::PDRblock::PDRblock(const boundBox& box, const labelVector& nCells)
544 :
545  PDRblock(dictionary::null, false)
546 {
547  addDefaultPatches();
548  reset(box, nCells);
549 }
550 
551 
552 Foam::PDRblock::PDRblock(const dictionary& dict, bool verboseOutput)
553 :
554  ijkMesh(),
555  meshDict_(dict),
556  grid_(),
557  outer_(),
558  bounds_(),
559  patches_(),
560  edgeLimits_(0,0),
561  verbose_(verboseOutput)
562 {
563  if (!dict.isNullDict())
564  {
565  read(dict);
566  }
567 }
568 
570 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
571 
573 {
574  const scalar scaleFactor(dict.getOrDefault<scalar>("scale", -1));
575 
576  expansionType expandType
577  (
578  expansionNames_.getOrDefault
579  (
580  "expansion",
581  dict,
582  expansionType::EXPAND_RATIO
583  )
584  );
585 
586  readGridControl(0, dict.subDict("x"), scaleFactor, expandType);
587  readGridControl(1, dict.subDict("y"), scaleFactor, expandType);
588  readGridControl(2, dict.subDict("z"), scaleFactor, expandType);
589 
590  adjustSizes();
591 
592  readBoundary(dict);
593 
594  // Outer treatment: (none | extend | box | sphere)
595  outer_.clear();
596 
597  const dictionary* outerDictPtr = dict.findDict("outer");
598  if (outerDictPtr)
599  {
600  outer_.read(*outerDictPtr);
601  }
602  outer_.report(Info);
603 
604  return true;
605 }
607 
609 (
610  const UList<scalar>& xgrid,
611  const UList<scalar>& ygrid,
612  const UList<scalar>& zgrid
613 )
614 {
615  static_cast<scalarList&>(grid_.x()) = xgrid;
616  static_cast<scalarList&>(grid_.y()) = ygrid;
617  static_cast<scalarList&>(grid_.z()) = zgrid;
618 
619  #ifdef FULLDEBUG
620  for (direction cmpt=0; cmpt < vector::nComponents; ++cmpt)
621  {
622  checkMonotonic(cmpt, grid_[cmpt]);
623  }
624  #endif
625 
626  adjustSizes();
627 
628  // Adjust boundaries
629  for (boundaryEntry& bentry : patches_)
630  {
631  bentry.size_ = 0;
632 
633  // Count patch faces
634  for (const label shapeFacei : bentry.faces_)
635  {
636  bentry.size_ += nBoundaryFaces(shapeFacei);
637  }
638  }
639 }
640 
641 
642 void Foam::PDRblock::reset(const boundBox& box, const labelVector& nCells)
643 {
644  grid_.x().reset(box.min().x(), box.max().x(), nCells.x());
645  grid_.y().reset(box.min().y(), box.max().y(), nCells.y());
646  grid_.z().reset(box.min().z(), box.max().z(), nCells.z());
647 
648  adjustSizes();
649 
650  // Adjust boundaries
651  for (boundaryEntry& bentry : patches_)
652  {
653  bentry.size_ = 0;
654 
655  // Count patch faces
656  for (const label shapeFacei : bentry.faces_)
657  {
658  bentry.size_ += nBoundaryFaces(shapeFacei);
659  }
660  }
661 }
662 
663 
664 bool Foam::PDRblock::findCell(const point& pt, labelVector& pos) const
665 {
666  // Out-of-bounds is handled explicitly, for efficiency and consistency,
667  // but principally to ensure that findLower() returns a valid
668  // result when the point is to the right of the bounds.
669 
670  // Since findLower returns the lower index, it corresponds to the
671  // cell in which the point is found
672 
673  if (!bounds_.contains(pt))
674  {
675  return false;
676  }
677 
678  for (direction cmpt=0; cmpt < labelVector::nComponents; ++cmpt)
679  {
680  // Binary search
681  pos[cmpt] = findLower(grid_[cmpt], pt[cmpt]);
682  }
683 
684  return true;
685 }
686 
687 
688 bool Foam::PDRblock::gridIndex
689 (
690  const point& pt,
691  labelVector& pos,
692  const scalar relTol
693 ) const
694 {
695  const scalar tol = relTol * edgeLimits_.min();
696 
697  for (direction cmpt=0; cmpt < labelVector::nComponents; ++cmpt)
698  {
699  // Linear search
700  pos[cmpt] = grid_[cmpt].findIndex(pt[cmpt], tol);
701 
702  if (pos[cmpt] < 0) return false;
703  }
704 
705  return true;
706 }
707 
708 
709 Foam::labelVector Foam::PDRblock::findCell(const point& pt) const
710 {
712 
713  if (findCell(pt, pos))
714  {
715  return pos;
716  }
717 
718  return labelVector(-1,-1,-1);
719 }
721 
722 Foam::labelVector Foam::PDRblock::gridIndex
723 (
724  const point& pt,
725  const scalar relTol
726 ) const
727 {
729 
730  if (gridIndex(pt, pos, relTol))
731  {
732  return pos;
733  }
734 
735  return labelVector(-1,-1,-1);
736 }
737 
738 
740 {
741  return grading(control_);
742 }
743 
744 
746 {
747  switch (cmpt)
748  {
749  case vector::X :
750  case vector::Y :
751  case vector::Z :
752  {
753  return control_[cmpt].grading();
754  break;
755  }
756 
757  default :
759  << "Not gridControl for direction " << label(cmpt) << endl
760  << exit(FatalError);
761  break;
762  }
763 
764  return gradingDescriptors();
765 }
766 
767 
768 // ************************************************************************* //
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return an entry data stream. FatalIOError if not found, or not a stream. ...
Definition: dictionary.C:367
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
scalar relativeToGeometricRatio(const scalar expRatio, const label nDiv)
Calculate geometric ratio from relative ratio.
Definition: PDRblock.C:62
label findLower(const ListType &input, const T &val, const label start, const ComparePredicate &comp)
Binary search to find the index of the last element in a sorted list that is less than value...
uint8_t direction
Definition: direction.H:46
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
bool read(Istream &is)
Read dictionary from Istream (discards the header). Reads entries until EOF or when the first token i...
Definition: dictionaryIO.C:134
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
bool equal(const T &a, const T &b)
Compare two values for equality.
Definition: label.H:164
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
label nCells() const
The number of mesh cells (nx*ny*nz) in the i-j-k mesh.
Definition: ijkMeshI.H:61
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
const Cmpt & y() const noexcept
Access to the vector y component.
Definition: Vector.H:140
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
const point & min() const noexcept
Minimum describing the bounding box.
Definition: boundBoxI.H:162
expansionType
The expansion type.
Definition: PDRblock.H:160
labelList faceLabels(nFaceLabels)
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
Various functions to operate on Lists.
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:104
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
bool read(const dictionary &dict)
Read dictionary.
Definition: PDRblock.C:569
dimensionedScalar pos(const dimensionedScalar &ds)
scalar y
const point & max() const noexcept
Maximum describing the bounding box.
Definition: boundBoxI.H:168
PDRblock()
Default construct, zero-size, inverted bounds etc.
Definition: PDRblock.C:520
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:137
static const Enum< expansionType > expansionNames_
Named enumerations for the expansion type.
Definition: PDRblock.H:170
A single block x-y-z rectilinear mesh addressable as i,j,k with simplified creation. Some of the input is similar to blockMeshDict, but since this specialization is for a single-block that is aligned with the x-y-z directions, it provides a different means of specifying the mesh.
Definition: PDRblock.H:149
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:105
static constexpr direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:109
static const PDRblock & null()
Return a PDRblock reference to a nullObject.
Definition: PDRblock.C:106
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Vector< gradingDescriptors > grading() const
Equivalent edge grading descriptors in (x,y,z) directions.
Definition: PDRblock.C:736
const Cmpt & x() const noexcept
Access to the vector x component.
Definition: Vector.H:135
scalar calcGexp(const scalar expRatio, const label nDiv)
Calculate the geometric expansion factor from the expansion ratio.
Definition: lineDivide.C:32
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
A simple i-j-k (row-major order) to linear addressing for a rectilinear mesh. Since the underlying me...
Definition: ijkMesh.H:53
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void clear()
Reset to (0,0,0) sizing.
const Cmpt & z() const noexcept
Access to the vector z component.
Definition: Vector.H:145
vector point
Point is a vector.
Definition: point.H:37
#define WarningInFunction
Report a warning using Foam::Warning.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
PtrList< volScalarField > & Y
#define Log
Definition: PDRblock.C:28
void reset(const UList< scalar > &xgrid, const UList< scalar > &ygrid, const UList< scalar > &zgrid)
Reset grid locations and mesh i-j-k sizing.
Definition: PDRblock.C:606
messageStream Info
Information stream (stdout output on master, null elsewhere)
Vector< label > labelVector
Vector of labels.
Definition: labelVector.H:47
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
List< label > labelList
A List of labels.
Definition: List.H:62
tmp< GeometricField< Type, faPatchField, areaMesh > > ndiv(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facNDiv.C:43
List of gradingDescriptor for the sections of a block with additional IO functionality.
Namespace for OpenFOAM.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225
const pointField & pts
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and a sub-dictionary) otherwise return nullptr...
Definition: dictionaryI.H:124
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...