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 
109 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
110 
111 void Foam::PDRblock::addDefaultPatches()
112 {
113  // Default boundaries with patchi == shapeFacei
114  patches_.resize(6);
115  for (label patchi=0; patchi < 6; ++patchi)
116  {
117  boundaryEntry& bentry = patches_.emplace_set(patchi);
118 
119  bentry.name_ = "patch" + Foam::name(patchi);
120  bentry.type_ = "patch";
121  bentry.size_ = 0;
122  bentry.faces_ = labelList(one{}, patchi);
123  }
124 }
125 
126 
127 void Foam::PDRblock::adjustSizes()
128 {
129  // Adjust i-j-k addressing
130  sizes().x() = grid_.x().nCells();
131  sizes().y() = grid_.y().nCells();
132  sizes().z() = grid_.z().nCells();
133 
134  if (sizes().x() <= 0 || sizes().y() <= 0 || sizes().z() <= 0)
135  {
136  // Sanity check. Silently disallow bad sizing
137  ijkMesh::clear();
138 
139  grid_.x().clear();
140  grid_.y().clear();
141  grid_.z().clear();
142 
143  bounds_.reset();
144  edgeLimits_.min() = 0;
145  edgeLimits_.max() = 0;
146  return;
147  }
148 
149  // Adjust boundBox
150  bounds_ = bounds(grid_.x(), grid_.y(), grid_.z());
151 
152  // Min/max edge lengths
153  edgeLimits_.clear();
154 
155  edgeLimits_.add(grid_.x().edgeLimits());
156  edgeLimits_.add(grid_.y().edgeLimits());
157  edgeLimits_.add(grid_.z().edgeLimits());
158 }
159 
160 
161 void Foam::PDRblock::readGridControl
162 (
163  const direction cmpt,
164  const dictionary& dict,
165  const scalar scaleFactor,
166  expansionType expandType
167 )
168 {
169  gridControl& ctrl = control_[cmpt];
170 
171  // Begin/end nodes for each segment
172  scalarList& knots = static_cast<scalarList&>(ctrl);
173 
174  // The number of division per segment
175  labelList& divisions = ctrl.divisions_;
176 
177  // The expansion ratio per segment
178  scalarList& expansion = ctrl.expansion_;
179 
180  expansion.clear(); // expansion is optional
181 
182  Log << "Reading grid control for "
183  << vector::componentNames[cmpt] << " direction" << nl;
184 
185  // Points
186  // ~~~~~~
187 
188  dict.readEntry("points", knots);
189 
190  const label nSegments = (knots.size()-1);
191 
192  if (nSegments < 1)
193  {
195  << "Must define at least two control points:"
196  << " in " << vector::componentNames[cmpt]
197  << " direction" << nl
198  << exit(FatalIOError);
199  }
200 
201  // Apply scaling
202  if (scaleFactor > 0)
203  {
204  for (scalar& pt : knots)
205  {
206  pt *= scaleFactor;
207  }
208  }
209 
210  // Do points increase monotonically?
211  checkMonotonic(cmpt, knots);
212 
213  Log << " points : " << flatOutput(knots) << nl;
214 
215 
216  // Divisions
217  // ~~~~~~~~~
218 
219  dict.readEntry("nCells", divisions);
220 
221  label nTotalDivs = 0;
222  for (const label ndiv : divisions)
223  {
224  nTotalDivs += ndiv;
225 
226  if (ndiv < 1)
227  {
229  << "Negative or zero divisions:"
230  << " in " << vector::componentNames[cmpt]
231  << " direction" << nl
232  << exit(FatalIOError);
233  }
234  }
235 
236  if (divisions.size() != nSegments)
237  {
239  << "Mismatch in number of divisions and number of points:"
240  << " in " << vector::componentNames[cmpt]
241  << " direction" << nl
242  << exit(FatalIOError);
243  }
244  else if (!nTotalDivs)
245  {
247  << "No grid divisions at all:"
248  << " in " << vector::componentNames[cmpt]
249  << " direction" << nl
250  << exit(FatalIOError);
251  }
252 
253  Log << " nCells : " << flatOutput(divisions) << nl;
254 
255 
256  // Expansion ratios
257  // ~~~~~~~~~~~~~~~~
258 
259  dict.readIfPresent("ratios", expansion);
260 
261  // Correct expansion ratios - negative is the same as inverse.
262  for (scalar& expRatio : expansion)
263  {
264  if (expRatio < 0)
265  {
266  expRatio = 1.0/(-expRatio);
267  }
268  }
269 
270  if (expansion.size() && divisions.size() != expansion.size())
271  {
273  << "Mismatch in number of expansion ratios and divisions:"
274  << " in " << vector::componentNames[cmpt]
275  << " direction" << nl
276  << exit(FatalIOError);
277  }
278 
279  if (expansion.empty())
280  {
281  expansion.resize(nSegments, scalar(1));
282 
283  if (expandType != expansionType::EXPAND_UNIFORM)
284  {
285  expandType = expansionType::EXPAND_UNIFORM;
286 
287  Log << "Warning: no 'ratios', use uniform spacing" << nl;
288  }
289  }
290  else
291  {
292  switch (expandType)
293  {
294  case expansionType::EXPAND_UNIFORM:
295  {
296  expansion = scalar(1);
297  break;
298  }
299 
300  case expansionType::EXPAND_RATIO:
301  {
302  Log << " ratios : " << flatOutput(expansion) << nl;
303  break;
304  }
305 
306  case expansionType::EXPAND_RELATIVE:
307  {
308  Log << " relative : " << flatOutput(expansion) << nl;
309 
310  auto divIter = divisions.cbegin();
311 
312  for (scalar& expRatio : expansion)
313  {
314  expRatio = relativeToGeometricRatio(expRatio, *divIter);
315  ++divIter;
316  }
317 
318  Log << " ratios : " << flatOutput(expansion) << nl;
319  break;
320  }
321  }
322  }
323 
324 
325  // Define the grid points
326 
327  scalarList& gridPoint = grid_[cmpt];
328  gridPoint.resize(nTotalDivs+1);
329 
330  label start = 0;
331 
332  for (label segmenti=0; segmenti < nSegments; ++segmenti)
333  {
334  const label nDiv = divisions[segmenti];
335  const scalar expRatio = expansion[segmenti];
336 
337  SubList<scalar> subPoint(gridPoint, nDiv+1, start);
338  start += nDiv;
339 
340  subPoint[0] = knots[segmenti];
341  subPoint[nDiv] = knots[segmenti+1];
342 
343  const scalar dist = (subPoint.last() - subPoint.first());
344 
345  if
346  (
347  expandType == expansionType::EXPAND_UNIFORM
348  || equal(expRatio, 1)
349  )
350  {
351  for (label i=1; i < nDiv; ++i)
352  {
353  subPoint[i] = (subPoint[0] + (dist * i)/nDiv);
354  }
355  }
356  else
357  {
358  // Geometric expansion factor from the expansion ratio
359  const scalar expFact = calcGexp(expRatio, nDiv);
360 
361  for (label i=1; i < nDiv; ++i)
362  {
363  subPoint[i] =
364  (
365  subPoint[0]
366  + dist * (1.0 - pow(expFact, i))/(1.0 - pow(expFact, nDiv))
367  );
368  }
369  }
370  }
371 }
372 
373 
374 void Foam::PDRblock::readBoundary(const dictionary& dict)
375 {
376  Log << "Reading boundary entries" << nl;
377 
378  PtrList<entry> patchEntries;
379 
380  if (dict.found("boundary"))
381  {
382  PtrList<entry> newEntries(dict.lookup("boundary"));
383  patchEntries.transfer(newEntries);
384  }
385 
386  patches_.clear();
387  patches_.resize(patchEntries.size());
388 
389  // Hex cell has 6 sides
390  const labelRange validFaceId(0, 6);
391 
392  // Track which sides have been handled
393  labelHashSet handled;
394 
395  forAll(patchEntries, patchi)
396  {
397  if (!patchEntries.set(patchi))
398  {
399  continue;
400  }
401 
402  const entry& e = patchEntries[patchi];
403 
404  if (!e.isDict())
405  {
407  << "Entry " << e
408  << " in boundary section is not a valid dictionary."
409  << exit(FatalIOError);
410  }
411 
412  const dictionary& patchDict = e.dict();
413 
414  const word& patchName = e.keyword();
415 
416  const word patchType(patchDict.get<word>("type"));
417 
418  labelList faceLabels(patchDict.get<labelList>("faces"));
419 
420  labelHashSet patchFaces(faceLabels);
421 
422  if (faceLabels.size() != patchFaces.size())
423  {
425  << "Patch: " << patchName
426  << ": Ignoring duplicate face ids" << nl;
427  }
428 
429  label nErased = patchFaces.filterKeys(validFaceId);
430  if (nErased)
431  {
433  << "Patch: " << patchName << ": Ignoring " << nErased
434  << " faces with invalid identifiers" << nl;
435  }
436 
437  nErased = patchFaces.erase(handled);
438  if (nErased)
439  {
441  << "Patch: " << patchName << ": Ignoring " << nErased
442  << " faces, which were already handled" << nl;
443  }
444 
445  // Mark these faces as having been handled
446  handled += patchFaces;
447 
448 
449  // Save information for later access during mesh creation.
450 
451  boundaryEntry& bentry = patches_.emplace_set(patchi);
452 
453  bentry.name_ = patchName;
454  bentry.type_ = patchType;
455  bentry.size_ = 0;
456  bentry.faces_ = patchFaces.sortedToc();
457 
458  // Count patch faces
459  for (const label shapeFacei : bentry.faces_)
460  {
461  bentry.size_ += nBoundaryFaces(shapeFacei);
462  }
463  }
464 
465 
466  // Deal with unhandled faces
467 
468  labelHashSet missed(identity(6));
469  missed.erase(handled);
470 
471  if (missed.size())
472  {
473  boundaryEntry& bentry = patches_.emplace_back();
474 
475  bentry.name_ = "defaultFaces";
476  bentry.type_ = emptyPolyPatch::typeName;
477  bentry.size_ = 0;
478  bentry.faces_ = missed.sortedToc();
479 
480  // Count patch faces
481  for (const label shapeFacei : bentry.faces_)
482  {
483  bentry.size_ += nBoundaryFaces(shapeFacei);
484  }
485 
486 
487  // Use name/type from "defaultPatch" entry if available
488  // - be generous with error handling
489 
490  const dictionary* defaultEntry = dict.findDict("defaultPatch");
491  if (defaultEntry)
492  {
493  defaultEntry->readIfPresent("name", bentry.name_);
494  defaultEntry->readIfPresent("type", bentry.type_);
495  }
496  }
497 
498  if (verbose_)
499  {
500  label patchi = 0;
501  for (const boundaryEntry& bentry : patches_)
502  {
503  Info<< " patch: " << patchi
504  << " (size: " << bentry.size_
505  << " type: " << bentry.type_
506  << ") name: " << bentry.name_
507  << " faces: " << flatOutput(bentry.faces_) << nl;
508 
509  ++patchi;
510  }
511  }
512 }
513 
515 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
516 
518 :
519  PDRblock(dictionary::null, false)
520 {}
522 
524 (
525  const UList<scalar>& xgrid,
526  const UList<scalar>& ygrid,
527  const UList<scalar>& zgrid
528 )
529 :
530  PDRblock(dictionary::null, false)
531 {
532  addDefaultPatches();
533  reset(xgrid, ygrid, zgrid);
534 }
535 
536 
537 Foam::PDRblock::PDRblock(const boundBox& box, const labelVector& nCells)
538 :
539  PDRblock(dictionary::null, false)
540 {
541  addDefaultPatches();
542  reset(box, nCells);
543 }
544 
545 
546 Foam::PDRblock::PDRblock(const dictionary& dict, bool verboseOutput)
547 :
548  ijkMesh(),
549  meshDict_(dict),
550  grid_(),
551  outer_(),
552  bounds_(),
553  patches_(),
554  edgeLimits_(0,0),
555  verbose_(verboseOutput)
556 {
557  if (!dict.isNullDict())
558  {
559  read(dict);
560  }
561 }
562 
564 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
565 
567 {
568  const scalar scaleFactor(dict.getOrDefault<scalar>("scale", -1));
569 
570  expansionType expandType
571  (
572  expansionNames_.getOrDefault
573  (
574  "expansion",
575  dict,
576  expansionType::EXPAND_RATIO
577  )
578  );
579 
580  readGridControl(0, dict.subDict("x"), scaleFactor, expandType);
581  readGridControl(1, dict.subDict("y"), scaleFactor, expandType);
582  readGridControl(2, dict.subDict("z"), scaleFactor, expandType);
583 
584  adjustSizes();
585 
586  readBoundary(dict);
587 
588  // Outer treatment: (none | extend | box | sphere)
589  outer_.clear();
590 
591  const dictionary* outerDictPtr = dict.findDict("outer");
592  if (outerDictPtr)
593  {
594  outer_.read(*outerDictPtr);
595  }
596  outer_.report(Info);
597 
598  return true;
599 }
601 
603 (
604  const UList<scalar>& xgrid,
605  const UList<scalar>& ygrid,
606  const UList<scalar>& zgrid
607 )
608 {
609  static_cast<scalarList&>(grid_.x()) = xgrid;
610  static_cast<scalarList&>(grid_.y()) = ygrid;
611  static_cast<scalarList&>(grid_.z()) = zgrid;
612 
613  #ifdef FULLDEBUG
614  for (direction cmpt=0; cmpt < vector::nComponents; ++cmpt)
615  {
616  checkMonotonic(cmpt, grid_[cmpt]);
617  }
618  #endif
619 
620  adjustSizes();
621 
622  // Adjust boundaries
623  for (boundaryEntry& bentry : patches_)
624  {
625  bentry.size_ = 0;
626 
627  // Count patch faces
628  for (const label shapeFacei : bentry.faces_)
629  {
630  bentry.size_ += nBoundaryFaces(shapeFacei);
631  }
632  }
633 }
634 
635 
636 void Foam::PDRblock::reset(const boundBox& box, const labelVector& nCells)
637 {
638  grid_.x().reset(box.min().x(), box.max().x(), nCells.x());
639  grid_.y().reset(box.min().y(), box.max().y(), nCells.y());
640  grid_.z().reset(box.min().z(), box.max().z(), nCells.z());
641 
642  adjustSizes();
643 
644  // Adjust boundaries
645  for (boundaryEntry& bentry : patches_)
646  {
647  bentry.size_ = 0;
648 
649  // Count patch faces
650  for (const label shapeFacei : bentry.faces_)
651  {
652  bentry.size_ += nBoundaryFaces(shapeFacei);
653  }
654  }
655 }
656 
657 
658 bool Foam::PDRblock::findCell(const point& pt, labelVector& pos) const
659 {
660  // Out-of-bounds is handled explicitly, for efficiency and consistency,
661  // but principally to ensure that findLower() returns a valid
662  // result when the point is to the right of the bounds.
663 
664  // Since findLower returns the lower index, it corresponds to the
665  // cell in which the point is found
666 
667  if (!bounds_.contains(pt))
668  {
669  return false;
670  }
671 
672  for (direction cmpt=0; cmpt < labelVector::nComponents; ++cmpt)
673  {
674  // Binary search
675  pos[cmpt] = findLower(grid_[cmpt], pt[cmpt]);
676  }
677 
678  return true;
679 }
680 
681 
682 bool Foam::PDRblock::gridIndex
683 (
684  const point& pt,
685  labelVector& pos,
686  const scalar relTol
687 ) const
688 {
689  const scalar tol = relTol * edgeLimits_.min();
690 
691  for (direction cmpt=0; cmpt < labelVector::nComponents; ++cmpt)
692  {
693  // Linear search
694  pos[cmpt] = grid_[cmpt].findIndex(pt[cmpt], tol);
695 
696  if (pos[cmpt] < 0) return false;
697  }
698 
699  return true;
700 }
701 
702 
703 Foam::labelVector Foam::PDRblock::findCell(const point& pt) const
704 {
706 
707  if (findCell(pt, pos))
708  {
709  return pos;
710  }
711 
712  return labelVector(-1,-1,-1);
713 }
715 
716 Foam::labelVector Foam::PDRblock::gridIndex
717 (
718  const point& pt,
719  const scalar relTol
720 ) const
721 {
723 
724  if (gridIndex(pt, pos, relTol))
725  {
726  return pos;
727  }
728 
729  return labelVector(-1,-1,-1);
730 }
731 
732 
734 {
735  return grading(control_);
736 }
737 
738 
740 {
741  switch (cmpt)
742  {
743  case vector::X :
744  case vector::Y :
745  case vector::Z :
746  {
747  return control_[cmpt].grading();
748  break;
749  }
750 
751  default :
753  << "Not gridControl for direction " << label(cmpt) << endl
754  << exit(FatalError);
755  break;
756  }
757 
758  return gradingDescriptors();
759 }
760 
761 
762 // ************************************************************************* //
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:153
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:608
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:563
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:514
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:130
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
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:730
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:637
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:600
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 it is a 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 ...