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-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 "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::adjustSizes()
118 {
119  // Adjust i-j-k addressing
120  sizes().x() = grid_.x().nCells();
121  sizes().y() = grid_.y().nCells();
122  sizes().z() = grid_.z().nCells();
123 
124  if (sizes().x() <= 0 || sizes().y() <= 0 || sizes().z() <= 0)
125  {
126  // Sanity check. Silently disallow bad sizing
127  ijkMesh::clear();
128 
129  grid_.x().clear();
130  grid_.y().clear();
131  grid_.z().clear();
132 
133  bounds_ = boundBox::invertedBox;
134  edgeLimits_.min() = 0;
135  edgeLimits_.max() = 0;
136  return;
137  }
138 
139  // Adjust boundBox
140  bounds_ = bounds(grid_.x(), grid_.y(), grid_.z());
141 
142  // Min/max edge lengths
143  edgeLimits_.clear();
144 
145  edgeLimits_.add(grid_.x().edgeLimits());
146  edgeLimits_.add(grid_.y().edgeLimits());
147  edgeLimits_.add(grid_.z().edgeLimits());
148 }
149 
150 
151 void Foam::PDRblock::readGridControl
152 (
153  const direction cmpt,
154  const dictionary& dict,
155  const scalar scaleFactor,
156  expansionType expandType
157 )
158 {
159  gridControl& ctrl = control_[cmpt];
160 
161  // Begin/end nodes for each segment
162  scalarList& knots = static_cast<scalarList&>(ctrl);
163 
164  // The number of division per segment
165  labelList& divisions = ctrl.divisions_;
166 
167  // The expansion ratio per segment
168  scalarList& expansion = ctrl.expansion_;
169 
170  expansion.clear(); // expansion is optional
171 
172  Log << "Reading grid control for "
173  << vector::componentNames[cmpt] << " direction" << nl;
174 
175  // Points
176  // ~~~~~~
177 
178  dict.readEntry("points", knots);
179 
180  const label nSegments = (knots.size()-1);
181 
182  if (nSegments < 1)
183  {
185  << "Must define at least two control points:"
186  << " in " << vector::componentNames[cmpt]
187  << " direction" << nl
188  << exit(FatalIOError);
189  }
190 
191  // Apply scaling
192  if (scaleFactor > 0)
193  {
194  for (scalar& pt : knots)
195  {
196  pt *= scaleFactor;
197  }
198  }
199 
200  // Do points increase monotonically?
201  checkMonotonic(cmpt, knots);
202 
203  Log << " points : " << flatOutput(knots) << nl;
204 
205 
206  // Divisions
207  // ~~~~~~~~~
208 
209  dict.readEntry("nCells", divisions);
210 
211  label nTotalDivs = 0;
212  for (const label ndiv : divisions)
213  {
214  nTotalDivs += ndiv;
215 
216  if (ndiv < 1)
217  {
219  << "Negative or zero divisions:"
220  << " in " << vector::componentNames[cmpt]
221  << " direction" << nl
222  << exit(FatalIOError);
223  }
224  }
225 
226  if (divisions.size() != nSegments)
227  {
229  << "Mismatch in number of divisions and number of points:"
230  << " in " << vector::componentNames[cmpt]
231  << " direction" << nl
232  << exit(FatalIOError);
233  }
234  else if (!nTotalDivs)
235  {
237  << "No grid divisions at all:"
238  << " in " << vector::componentNames[cmpt]
239  << " direction" << nl
240  << exit(FatalIOError);
241  }
242 
243  Log << " nCells : " << flatOutput(divisions) << nl;
244 
245 
246  // Expansion ratios
247  // ~~~~~~~~~~~~~~~~
248 
249  dict.readIfPresent("ratios", expansion);
250 
251  // Correct expansion ratios - negative is the same as inverse.
252  for (scalar& expRatio : expansion)
253  {
254  if (expRatio < 0)
255  {
256  expRatio = 1.0/(-expRatio);
257  }
258  }
259 
260  if (expansion.size() && divisions.size() != expansion.size())
261  {
263  << "Mismatch in number of expansion ratios and divisions:"
264  << " in " << vector::componentNames[cmpt]
265  << " direction" << nl
266  << exit(FatalIOError);
267  }
268 
269  if (expansion.empty())
270  {
271  expansion.resize(nSegments, scalar(1));
272 
273  if (expandType != expansionType::EXPAND_UNIFORM)
274  {
275  expandType = expansionType::EXPAND_UNIFORM;
276 
277  Log << "Warning: no 'ratios', use uniform spacing" << nl;
278  }
279  }
280  else
281  {
282  switch (expandType)
283  {
284  case expansionType::EXPAND_UNIFORM:
285  {
286  expansion = scalar(1);
287  break;
288  }
289 
290  case expansionType::EXPAND_RATIO:
291  {
292  Log << " ratios : " << flatOutput(expansion) << nl;
293  break;
294  }
295 
296  case expansionType::EXPAND_RELATIVE:
297  {
298  Log << " relative : " << flatOutput(expansion) << nl;
299 
300  auto divIter = divisions.cbegin();
301 
302  for (scalar& expRatio : expansion)
303  {
304  expRatio = relativeToGeometricRatio(expRatio, *divIter);
305  ++divIter;
306  }
307 
308  Log << " ratios : " << flatOutput(expansion) << nl;
309  break;
310  }
311  }
312  }
313 
314 
315  // Define the grid points
316 
317  scalarList& gridPoint = grid_[cmpt];
318  gridPoint.resize(nTotalDivs+1);
319 
320  label start = 0;
321 
322  for (label segmenti=0; segmenti < nSegments; ++segmenti)
323  {
324  const label nDiv = divisions[segmenti];
325  const scalar expRatio = expansion[segmenti];
326 
327  SubList<scalar> subPoint(gridPoint, nDiv+1, start);
328  start += nDiv;
329 
330  subPoint[0] = knots[segmenti];
331  subPoint[nDiv] = knots[segmenti+1];
332 
333  const scalar dist = (subPoint.last() - subPoint.first());
334 
335  if
336  (
337  expandType == expansionType::EXPAND_UNIFORM
338  || equal(expRatio, 1)
339  )
340  {
341  for (label i=1; i < nDiv; ++i)
342  {
343  subPoint[i] = (subPoint[0] + (dist * i)/nDiv);
344  }
345  }
346  else
347  {
348  // Geometric expansion factor from the expansion ratio
349  const scalar expFact = calcGexp(expRatio, nDiv);
350 
351  for (label i=1; i < nDiv; ++i)
352  {
353  subPoint[i] =
354  (
355  subPoint[0]
356  + dist * (1.0 - pow(expFact, i))/(1.0 - pow(expFact, nDiv))
357  );
358  }
359  }
360  }
361 }
362 
363 
364 void Foam::PDRblock::readBoundary(const dictionary& dict)
365 {
366  Log << "Reading boundary entries" << nl;
367 
368  PtrList<entry> patchEntries;
369 
370  if (dict.found("boundary"))
371  {
372  PtrList<entry> newEntries(dict.lookup("boundary"));
373  patchEntries.transfer(newEntries);
374  }
375 
376  patches_.clear();
377  patches_.resize(patchEntries.size());
378 
379  // Hex cell has 6 sides
380  const labelRange validFaceId(0, 6);
381 
382  // Track which sides have been handled
383  labelHashSet handled;
384 
385  forAll(patchEntries, patchi)
386  {
387  if (!patchEntries.set(patchi))
388  {
389  continue;
390  }
391 
392  const entry& e = patchEntries[patchi];
393 
394  if (!e.isDict())
395  {
397  << "Entry " << e
398  << " in boundary section is not a valid dictionary."
399  << exit(FatalIOError);
400  }
401 
402  const dictionary& patchDict = e.dict();
403 
404  const word& patchName = e.keyword();
405 
406  const word patchType(patchDict.get<word>("type"));
407 
408  labelList faceLabels(patchDict.get<labelList>("faces"));
409 
410  labelHashSet patchFaces(faceLabels);
411 
412  if (faceLabels.size() != patchFaces.size())
413  {
415  << "Patch: " << patchName
416  << ": Ignoring duplicate face ids" << nl;
417  }
418 
419  label nErased = patchFaces.filterKeys(validFaceId);
420  if (nErased)
421  {
423  << "Patch: " << patchName << ": Ignoring " << nErased
424  << " faces with invalid identifiers" << nl;
425  }
426 
427  nErased = patchFaces.erase(handled);
428  if (nErased)
429  {
431  << "Patch: " << patchName << ": Ignoring " << nErased
432  << " faces, which were already handled" << nl;
433  }
434 
435  // Mark these faces as having been handled
436  handled += patchFaces;
437 
438 
439  // Save information for later access during mesh creation.
440 
441  patches_.set(patchi, new boundaryEntry());
442 
443  boundaryEntry& bentry = patches_[patchi];
444 
445  bentry.name_ = patchName;
446  bentry.type_ = patchType;
447  bentry.size_ = 0;
448  bentry.faces_ = patchFaces.sortedToc();
449 
450  // Count patch faces
451  for (const label shapeFacei : bentry.faces_)
452  {
453  bentry.size_ += nBoundaryFaces(shapeFacei);
454  }
455  }
456 
457 
458  // Deal with unhandled faces
459 
460  labelHashSet missed(identity(6));
461  missed.erase(handled);
462 
463  if (missed.size())
464  {
465  patches_.append(new boundaryEntry());
466 
467  boundaryEntry& bentry = patches_.last();
468 
469  bentry.name_ = "defaultFaces";
470  bentry.type_ = emptyPolyPatch::typeName;
471  bentry.size_ = 0;
472  bentry.faces_ = missed.sortedToc();
473 
474  // Count patch faces
475  for (const label shapeFacei : bentry.faces_)
476  {
477  bentry.size_ += nBoundaryFaces(shapeFacei);
478  }
479 
480 
481  // Use name/type from "defaultPatch" entry if available
482  // - be generous with error handling
483 
484  const dictionary* defaultEntry = dict.findDict("defaultPatch");
485  if (defaultEntry)
486  {
487  defaultEntry->readIfPresent("name", bentry.name_);
488  defaultEntry->readIfPresent("type", bentry.type_);
489  }
490  }
491 
492  if (verbose_)
493  {
494  label patchi = 0;
495  for (const boundaryEntry& bentry : patches_)
496  {
497  Info<< " patch: " << patchi
498  << " (size: " << bentry.size_
499  << " type: " << bentry.type_
500  << ") name: " << bentry.name_
501  << " faces: " << flatOutput(bentry.faces_) << nl;
502 
503  ++patchi;
504  }
505  }
506 }
507 
508 
509 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
510 
512 :
513  PDRblock(dictionary::null, false)
514 {}
515 
516 
518 (
519  const UList<scalar>& xgrid,
520  const UList<scalar>& ygrid,
521  const UList<scalar>& zgrid
522 )
523 :
524  PDRblock(dictionary::null, false)
525 {
526  // Default boundaries with patchi == shapeFacei
527  patches_.resize(6);
528  for (label patchi=0; patchi < 6; ++patchi)
529  {
530  patches_.set(patchi, new boundaryEntry());
531 
532  boundaryEntry& bentry = patches_[patchi];
533 
534  bentry.name_ = "patch" + Foam::name(patchi);
535  bentry.type_ = "patch";
536  bentry.size_ = 0;
537  bentry.faces_ = labelList(one{}, patchi);
538  }
539 
540  reset(xgrid, ygrid, zgrid);
541 }
542 
543 
544 Foam::PDRblock::PDRblock(const dictionary& dict, bool verboseOutput)
545 :
546  ijkMesh(),
547  meshDict_(dict),
548  grid_(),
549  outer_(),
550  bounds_(),
551  patches_(),
552  edgeLimits_(0,0),
553  verbose_(verboseOutput)
554 {
555  if (!dict.isNullDict())
556  {
557  read(dict);
558  }
559 }
560 
561 
562 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
563 
565 {
566  const scalar scaleFactor(dict.getOrDefault<scalar>("scale", -1));
567 
568  expansionType expandType
569  (
570  expansionNames_.getOrDefault
571  (
572  "expansion",
573  dict,
574  expansionType::EXPAND_RATIO
575  )
576  );
577 
578  readGridControl(0, dict.subDict("x"), scaleFactor, expandType);
579  readGridControl(1, dict.subDict("y"), scaleFactor, expandType);
580  readGridControl(2, dict.subDict("z"), scaleFactor, expandType);
581 
582  adjustSizes();
583 
584  readBoundary(dict);
585 
586  // Outer treatment: (none | extend | box | sphere)
587  outer_.clear();
588 
589  const dictionary* outerDictPtr = dict.findDict("outer");
590  if (outerDictPtr)
591  {
592  outer_.read(*outerDictPtr);
593  }
594  outer_.report(Info);
595 
596  return true;
597 }
598 
599 
601 (
602  const UList<scalar>& xgrid,
603  const UList<scalar>& ygrid,
604  const UList<scalar>& zgrid
605 )
606 {
607  static_cast<scalarList&>(grid_.x()) = xgrid;
608  static_cast<scalarList&>(grid_.y()) = ygrid;
609  static_cast<scalarList&>(grid_.z()) = zgrid;
610 
611  #ifdef FULLDEBUG
612  for (direction cmpt=0; cmpt < vector::nComponents; ++cmpt)
613  {
614  checkMonotonic(cmpt, grid_[cmpt]);
615  }
616  #endif
617 
618  adjustSizes();
619 
620  // Adjust boundaries
621  for (boundaryEntry& bentry : patches_)
622  {
623  bentry.size_ = 0;
624 
625  // Count patch faces
626  for (const label shapeFacei : bentry.faces_)
627  {
628  bentry.size_ += nBoundaryFaces(shapeFacei);
629  }
630  }
631 }
632 
633 
634 bool Foam::PDRblock::findCell(const point& pt, labelVector& pos) const
635 {
636  // Out-of-bounds is handled explicitly, for efficiency and consistency,
637  // but principally to ensure that findLower() returns a valid
638  // result when the point is to the right of the bounds.
639 
640  // Since findLower returns the lower index, it corresponds to the
641  // cell in which the point is found
642 
643  if (!bounds_.contains(pt))
644  {
645  return false;
646  }
647 
648  for (direction cmpt=0; cmpt < labelVector::nComponents; ++cmpt)
649  {
650  // Binary search
651  pos[cmpt] = findLower(grid_[cmpt], pt[cmpt]);
652  }
653 
654  return true;
655 }
656 
657 
658 bool Foam::PDRblock::gridIndex
659 (
660  const point& pt,
661  labelVector& pos,
662  const scalar relTol
663 ) const
664 {
665  const scalar tol = relTol * edgeLimits_.min();
666 
667  for (direction cmpt=0; cmpt < labelVector::nComponents; ++cmpt)
668  {
669  // Linear search
670  pos[cmpt] = grid_[cmpt].findIndex(pt[cmpt], tol);
671 
672  if (pos[cmpt] < 0) return false;
673  }
674 
675  return true;
676 }
677 
678 
679 Foam::labelVector Foam::PDRblock::findCell(const point& pt) const
680 {
682 
683  if (findCell(pt, pos))
684  {
685  return pos;
686  }
687 
688  return labelVector(-1,-1,-1);
689 }
690 
691 
692 Foam::labelVector Foam::PDRblock::gridIndex
693 (
694  const point& pt,
695  const scalar relTol
696 ) const
697 {
699 
700  if (gridIndex(pt, pos, relTol))
701  {
702  return pos;
703  }
704 
705  return labelVector(-1,-1,-1);
706 }
707 
708 
710 {
711  return grading(control_);
712 }
713 
714 
716 {
717  switch (cmpt)
718  {
719  case vector::X :
720  case vector::Y :
721  case vector::Z :
722  {
723  return control_[cmpt].grading();
724  break;
725  }
726 
727  default :
729  << "Not gridControl for direction " << label(cmpt) << endl
730  << exit(FatalError);
731  break;
732  }
733 
734  return gradingDescriptors();
735 }
736 
737 
738 // ************************************************************************* //
Foam::dictionary::findDict
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionaryI.H:127
Foam::PDRblock::expansionNames_
const static Enum< expansionType > expansionNames_
Named enumerations for the expansion type.
Definition: PDRblock.H:170
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::scalarList
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Log
#define Log
Definition: PDRblock.C:35
Foam::Vector< scalar >::Z
Definition: Vector.H:81
Foam::Vector< scalar >::Y
Definition: Vector.H:81
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::componentNames
static const char *const componentNames[]
Definition: VectorSpace.H:114
Foam::List::resize
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
Foam::PDRblock::grading
Vector< gradingDescriptors > grading() const
Equivalent edge grading descriptors in (x,y,z) directions.
Definition: PDRblock.C:709
Foam::PDRblock::null
static const PDRblock & null()
Return a PDRblock reference to a nullObject.
Definition: PDRblock.C:109
Foam::boundBox::invertedBox
static const boundBox invertedBox
A large inverted boundBox: min/max == +/- ROOTVGREAT.
Definition: boundBox.H:86
Foam::fac::ndiv
tmp< GeometricField< Type, faPatchField, areaMesh > > ndiv(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facNDiv.C:50
Foam::one
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:61
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::gradingDescriptors
List of gradingDescriptor for the sections of a block with additional IO functionality.
Definition: gradingDescriptors.H:58
Foam::calcGexp
scalar calcGexp(const scalar expRatio, const label nDiv)
Calculate the geometric expansion factor from the expansion ratio.
Definition: lineDivide.C:37
Foam::ijkAddressing::clear
void clear()
Reset to (0,0,0) sizing.
Definition: ijkAddressingI.H:97
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
gradingDescriptors.H
Foam::dictionary::null
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:392
Foam::labelVector
Vector< label > labelVector
Vector of labels.
Definition: labelVector.H:51
Foam::PDRblock::reset
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:601
Foam::ijkMesh
A simple i-j-k (row-major order) to linear addressing for a rectilinear mesh. Since the underlying me...
Definition: ijkMesh.H:57
Foam::findLower
label findLower(const ListType &input, const T &val, const label start, const ComparePredicate &comp)
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::Vector< scalar >::X
Definition: Vector.H:81
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:302
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
PDRblock.H
Foam::dictionary::lookup
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:386
dict
dictionary dict
Definition: searchingEngine.H:14
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:123
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::flatOutput
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:216
Foam::PDRblock
A single block x-y-z rectilinear mesh addressable as i,j,k with simplified creation....
Definition: PDRblock.H:153
emptyPolyPatch.H
Foam::dictionary::read
bool read(Istream &is)
Read dictionary from Istream. Discards the header.
Definition: dictionaryIO.C:141
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
reset
meshPtr reset(new Foam::fvMesh(Foam::IOobject(regionName, runTime.timeName(), runTime, Foam::IOobject::MUST_READ), false))
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::PDRblock::expansionType
expansionType
The expansion type.
Definition: PDRblock.H:162
Foam::VectorSpace::max
static const Form max
Definition: VectorSpace.H:117
Foam::Vector< scalar >
Foam::List< scalar >
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::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::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
Foam::direction
uint8_t direction
Definition: direction.H:52
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
ListOps.H
Various functions to operate on Lists.
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::equal
bool equal(const T &s1, const T &s2)
Compare two values for equality.
Definition: doubleFloat.H:46
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::PDRblock::read
bool read(const dictionary &dict)
Read dictionary.
Definition: PDRblock.C:564
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:405
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::nComponents
static constexpr direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:101
Foam::PDRblock::PDRblock
PDRblock()
Default construct, zero-size, inverted bounds etc.
Definition: PDRblock.C:511
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:177