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-2020 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 
32 // Output when verbosity is enabled
33 #undef Log
34 #define Log if (verbose_) Info
35 
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 const Foam::Enum
40 <
42 >
43 Foam::PDRblock::expansionNames_
44 ({
45  { expansionType::EXPAND_UNIFORM, "uniform" },
46  { expansionType::EXPAND_RATIO, "ratio" },
47  { expansionType::EXPAND_RELATIVE, "relative"}
48 });
49 
50 
51 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55  //- Calculate geometric expansion factor from expansion ratio
56  inline scalar calcGexp(const scalar expRatio, const label nDiv)
57  {
58  return nDiv > 1 ? pow(expRatio, 1.0/(nDiv - 1)) : 0.0;
59  }
60 
61  //- Calculate geometric ratio from relative ratio
62  inline scalar relativeToGeometricRatio
63  (
64  const scalar expRatio,
65  const label nDiv
66  )
67  {
68  return nDiv > 1 ? pow(expRatio, (nDiv - 1)) : 1.0;
69  }
70 
71 } // End namespace Foam
72 
73 
74 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
75 
76 bool Foam::PDRblock::checkMonotonic
77 (
78  const direction cmpt,
79  const UList<scalar>& pts
80 )
81 {
82  const label len = pts.size();
83 
84  if (!len)
85  {
86  return false;
87  }
88 
89  const scalar& minVal = pts[0];
90 
91  for (label i=1; i < len; ++i)
92  {
93  if (pts[i] <= minVal)
94  {
96  << "Points in " << vector::componentNames[cmpt]
97  << " direction do not increase monotonically" << nl
98  << flatOutput(pts) << nl << nl
99  << exit(FatalError);
100  }
101  }
102 
103  return true;
104 }
105 
106 
108 {
109  return NullObjectRef<PDRblock>();
110 }
111 
112 
113 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
114 
115 void Foam::PDRblock::adjustSizes()
116 {
117  // Adjust i-j-k addressing
118  sizes().x() = grid_.x().nCells();
119  sizes().y() = grid_.y().nCells();
120  sizes().z() = grid_.z().nCells();
121 
122  if (sizes().x() <= 0 || sizes().y() <= 0 || sizes().z() <= 0)
123  {
124  // Sanity check. Silently disallow bad sizing
125  ijkMesh::clear();
126 
127  grid_.x().clear();
128  grid_.y().clear();
129  grid_.z().clear();
130 
131  bounds_ = boundBox::invertedBox;
132  minEdgeLen_ = Zero;
133  return;
134  }
135 
136  // Adjust boundBox
137  bounds_.min().x() = grid_.x().first();
138  bounds_.min().y() = grid_.y().first();
139  bounds_.min().z() = grid_.z().first();
140 
141  bounds_.max().x() = grid_.x().last();
142  bounds_.max().y() = grid_.y().last();
143  bounds_.max().z() = grid_.z().last();
144 
145  // Min edge length
146  minEdgeLen_ = GREAT;
147 
148  for (direction cmpt=0; cmpt < vector::nComponents; ++cmpt)
149  {
150  const label nEdge = grid_[cmpt].nCells();
151 
152  for (label edgei=0; edgei < nEdge; ++edgei)
153  {
154  const scalar len = grid_[cmpt].width(edgei);
155  minEdgeLen_ = min(minEdgeLen_, len);
156  }
157  }
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  //- The begin/end nodes for each segment
170  scalarList knots;
171 
172  // The number of division per segment
173  labelList divisions;
174 
175  // The expansion ratio per segment
176  scalarList expansion;
177 
178  Log << "Reading grid control for "
179  << vector::componentNames[cmpt] << " direction" << nl;
180 
181  // Points
182  // ~~~~~~
183 
184  dict.readEntry("points", knots);
185 
186  const label nSegments = (knots.size()-1);
187 
188  if (nSegments < 1)
189  {
191  << "Must define at least two control points:"
192  << " in " << vector::componentNames[cmpt]
193  << " direction" << nl
194  << exit(FatalIOError);
195  }
196 
197  // Apply scaling
198  if (scaleFactor > 0)
199  {
200  for (scalar& pt : knots)
201  {
202  pt *= scaleFactor;
203  }
204  }
205 
206  // Do points increase monotonically?
207  checkMonotonic(cmpt, knots);
208 
209  Log << " points : " << flatOutput(knots) << nl;
210 
211 
212  // Divisions
213  // ~~~~~~~~~
214 
215  dict.readEntry("nCells", divisions);
216 
217  label nTotalDivs = 0;
218  for (const label ndiv : divisions)
219  {
220  nTotalDivs += ndiv;
221 
222  if (ndiv < 1)
223  {
225  << "Negative or zero divisions:"
226  << " in " << vector::componentNames[cmpt]
227  << " direction" << nl
228  << exit(FatalIOError);
229  }
230  }
231 
232  if (divisions.size() != nSegments)
233  {
235  << "Mismatch in number of divisions and number of points:"
236  << " in " << vector::componentNames[cmpt]
237  << " direction" << nl
238  << exit(FatalIOError);
239  }
240  else if (!nTotalDivs)
241  {
243  << "No grid divisions at all:"
244  << " in " << vector::componentNames[cmpt]
245  << " direction" << nl
246  << exit(FatalIOError);
247  }
248 
249  Log << " nCells : " << flatOutput(divisions) << nl;
250 
251 
252  // Expansion ratios
253  // ~~~~~~~~~~~~~~~~
254 
255  dict.readIfPresent("ratios", expansion);
256 
257  // Correct expansion ratios - negative is the same as inverse.
258  for (scalar& expRatio : expansion)
259  {
260  if (expRatio < 0)
261  {
262  expRatio = 1.0/(-expRatio);
263  }
264  }
265 
266  if (expansion.size() && divisions.size() != expansion.size())
267  {
269  << "Mismatch in number of expansion ratios and divisions:"
270  << " in " << vector::componentNames[cmpt]
271  << " direction" << nl
272  << exit(FatalIOError);
273  }
274 
275  if (expansion.empty())
276  {
277  expansion.resize(nSegments, scalar(1));
278 
279  if (expandType != expansionType::EXPAND_UNIFORM)
280  {
281  expandType = expansionType::EXPAND_UNIFORM;
282 
283  Log << "Warning: no 'ratios', use uniform spacing" << nl;
284  }
285  }
286  else
287  {
288  switch (expandType)
289  {
290  case expansionType::EXPAND_UNIFORM:
291  {
292  expansion = scalar(1);
293  break;
294  }
295 
296  case expansionType::EXPAND_RATIO:
297  {
298  Log << " ratios : " << flatOutput(expansion) << nl;
299  break;
300  }
301 
302  case expansionType::EXPAND_RELATIVE:
303  {
304  Log << " relative : " << flatOutput(expansion) << nl;
305 
306  auto divIter = divisions.cbegin();
307 
308  for (scalar& expRatio : expansion)
309  {
310  expRatio = relativeToGeometricRatio(expRatio, *divIter);
311  ++divIter;
312  }
313 
314  Log << " ratios : " << flatOutput(expansion) << nl;
315  break;
316  }
317  }
318  }
319 
320 
321  // Define the grid points
322 
323  scalarList& gridPoint = grid_[cmpt];
324  gridPoint.resize(nTotalDivs+1);
325 
326  label start = 0;
327 
328  for (label segmenti=0; segmenti < nSegments; ++segmenti)
329  {
330  const label nDiv = divisions[segmenti];
331  const scalar expRatio = expansion[segmenti];
332 
333  SubList<scalar> subPoint(gridPoint, nDiv+1, start);
334  start += nDiv;
335 
336  subPoint[0] = knots[segmenti];
337  subPoint[nDiv] = knots[segmenti+1];
338 
339  const scalar dist = (subPoint.last() - subPoint.first());
340 
341  if
342  (
343  expandType == expansionType::EXPAND_UNIFORM
344  || equal(expRatio, 1)
345  )
346  {
347  for (label i=1; i < nDiv; ++i)
348  {
349  subPoint[i] = (subPoint[0] + (dist * i)/nDiv);
350  }
351  }
352  else
353  {
354  // Geometric expansion factor from the expansion ratio
355  const scalar expFact = calcGexp(expRatio, nDiv);
356 
357  for (label i=1; i < nDiv; ++i)
358  {
359  subPoint[i] =
360  (
361  subPoint[0]
362  + dist * (1.0 - pow(expFact, i))/(1.0 - pow(expFact, nDiv))
363  );
364  }
365  }
366  }
367 }
368 
369 
370 void Foam::PDRblock::readBoundary(const dictionary& dict)
371 {
372  Log << "Reading boundary entries" << nl;
373 
374  PtrList<entry> patchEntries;
375 
376  if (dict.found("boundary"))
377  {
378  PtrList<entry> newEntries(dict.lookup("boundary"));
379  patchEntries.transfer(newEntries);
380  }
381 
382  patches_.clear();
383  patches_.resize(patchEntries.size());
384 
385  // Hex cell has 6 sides
386  const labelRange validFaceId(0, 6);
387 
388  // Track which sides have been handled
389  labelHashSet handled;
390 
391  forAll(patchEntries, patchi)
392  {
393  if (!patchEntries.set(patchi))
394  {
395  continue;
396  }
397 
398  const entry& e = patchEntries[patchi];
399 
400  if (!e.isDict())
401  {
403  << "Entry " << e
404  << " in boundary section is not a valid dictionary."
405  << exit(FatalIOError);
406  }
407 
408  const dictionary& patchDict = e.dict();
409 
410  const word& patchName = e.keyword();
411 
412  const word patchType(patchDict.get<word>("type"));
413 
414  labelList faceLabels(patchDict.get<labelList>("faces"));
415 
416  labelHashSet patchFaces(faceLabels);
417 
418  if (faceLabels.size() != patchFaces.size())
419  {
421  << "Patch: " << patchName
422  << ": Ignoring duplicate face ids" << nl;
423  }
424 
425  label nErased = patchFaces.filterKeys(validFaceId);
426  if (nErased)
427  {
429  << "Patch: " << patchName << ": Ignoring " << nErased
430  << " faces with invalid identifiers" << nl;
431  }
432 
433  nErased = patchFaces.erase(handled);
434  if (nErased)
435  {
437  << "Patch: " << patchName << ": Ignoring " << nErased
438  << " faces, which were already handled" << nl;
439  }
440 
441  // Mark these faces as having been handled
442  handled += patchFaces;
443 
444 
445  // Save information for later access during mesh creation.
446 
447  patches_.set(patchi, new boundaryEntry());
448 
449  boundaryEntry& bentry = patches_[patchi];
450 
451  bentry.name_ = patchName;
452  bentry.type_ = patchType;
453  bentry.size_ = 0;
454  bentry.faces_ = patchFaces.sortedToc();
455 
456  // Count patch faces
457  for (const label shapeFacei : bentry.faces_)
458  {
459  bentry.size_ += nBoundaryFaces(shapeFacei);
460  }
461  }
462 
463 
464  // Deal with unhandled faces
465 
466  labelHashSet missed(identity(6));
467  missed.erase(handled);
468 
469  if (missed.size())
470  {
471  patches_.append(new boundaryEntry());
472 
473  boundaryEntry& bentry = patches_.last();
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 
514 
515 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
516 
518 (
519  const UList<scalar>& xgrid,
520  const UList<scalar>& ygrid,
521  const UList<scalar>& zgrid
522 )
523 :
524  PDRblock()
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 
545 :
546  PDRblock()
547 {
548  verbose_ = verbose;
549 
550  read(dict);
551 }
552 
553 
554 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
555 
557 {
558  const scalar scaleFactor(dict.getOrDefault<scalar>("scale", -1));
559 
560  expansionType expandType
561  (
562  expansionNames_.getOrDefault
563  (
564  "expansion",
565  dict,
566  expansionType::EXPAND_RATIO
567  )
568  );
569 
570  readGridControl(0, dict.subDict("x"), scaleFactor, expandType);
571  readGridControl(1, dict.subDict("y"), scaleFactor, expandType);
572  readGridControl(2, dict.subDict("z"), scaleFactor, expandType);
573 
574  adjustSizes();
575 
576  readBoundary(dict);
577 
578  return true;
579 }
580 
581 
583 (
584  const UList<scalar>& xgrid,
585  const UList<scalar>& ygrid,
586  const UList<scalar>& zgrid
587 )
588 {
589  static_cast<scalarList&>(grid_.x()) = xgrid;
590  static_cast<scalarList&>(grid_.y()) = ygrid;
591  static_cast<scalarList&>(grid_.z()) = zgrid;
592 
593  #ifdef FULLDEBUG
594  for (direction cmpt=0; cmpt < vector::nComponents; ++cmpt)
595  {
596  checkMonotonic(cmpt, grid_[cmpt]);
597  }
598  #endif
599 
600  adjustSizes();
601 
602  // Adjust boundaries
603  for (boundaryEntry& bentry : patches_)
604  {
605  bentry.size_ = 0;
606 
607  // Count patch faces
608  for (const label shapeFacei : bentry.faces_)
609  {
610  bentry.size_ += nBoundaryFaces(shapeFacei);
611  }
612  }
613 }
614 
615 
616 bool Foam::PDRblock::findCell(const point& pt, labelVector& pos) const
617 {
618  // Out-of-bounds is handled explicitly, for efficiency and consistency,
619  // but principally to ensure that findLower() returns a valid
620  // result when the point is to the right of the bounds.
621 
622  // Since findLower returns the lower index, it corresponds to the
623  // cell in which the point is found
624 
625  if (!bounds_.contains(pt))
626  {
627  return false;
628  }
629 
630  for (direction cmpt=0; cmpt < labelVector::nComponents; ++cmpt)
631  {
632  // Binary search
633  pos[cmpt] = findLower(grid_[cmpt], pt[cmpt]);
634  }
635 
636  return true;
637 }
638 
639 
640 bool Foam::PDRblock::gridIndex
641 (
642  const point& pt,
643  labelVector& pos,
644  const scalar relTol
645 ) const
646 {
647  const scalar tol = relTol * minEdgeLen_;
648 
649  for (direction cmpt=0; cmpt < labelVector::nComponents; ++cmpt)
650  {
651  // Linear search
652  pos[cmpt] = grid_[cmpt].findIndex(pt[cmpt], tol);
653 
654  if (pos[cmpt] < 0) return false;
655  }
656 
657  return true;
658 }
659 
660 
661 Foam::labelVector Foam::PDRblock::findCell(const point& pt) const
662 {
664 
665  if (findCell(pt, pos))
666  {
667  return pos;
668  }
669 
670  return labelVector(-1,-1,-1);
671 }
672 
673 
674 Foam::labelVector Foam::PDRblock::gridIndex
675 (
676  const point& pt,
677  const scalar relTol
678 ) const
679 {
681 
682  if (gridIndex(pt, pos, relTol))
683  {
684  return pos;
685  }
686 
687  return labelVector(-1,-1,-1);
688 }
689 
690 
691 // ************************************************************************* //
Foam::dictionary::findDict
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionary.C:508
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
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:34
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::componentNames
static const char *const componentNames[]
Definition: VectorSpace.H:114
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
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: dictionary.C:364
Foam::PDRblock::null
static const PDRblock & null()
Return a PDRblock reference to a nullObject.
Definition: PDRblock.C:107
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::calcGexp
scalar calcGexp(const scalar expRatio, const label nDiv)
Calculate the geometric expansion factor from the expansion ratio.
Definition: lineDivide.C:36
Foam::ijkAddressing::clear
void clear()
Reset to (0,0,0) sizing.
Definition: ijkAddressingI.H:97
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::PDRblock::expansionType
expansionType
The expansion type.
Definition: PDRblock.H:158
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
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:583
Foam::findLower
label findLower(const ListType &input, const T &val, const label start, const ComparePredicate &comp)
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:528
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:314
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
PDRblock.H
Foam::List::resize
void resize(const label newSize)
Adjust allocated size of list.
Definition: ListI.H:139
Foam::dictionary::lookup
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:424
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:121
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::PDRblock
A single block x-y-z rectilinear mesh addressable as i,j,k with simplified creation....
Definition: PDRblock.H:149
emptyPolyPatch.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::relativeToGeometricRatio
scalar relativeToGeometricRatio(const scalar expRatio, const label nDiv)
Calculate geometric ratio from relative ratio.
Definition: PDRblock.C:63
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:372
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::flatOutput
FlatOutput< Container > flatOutput(const Container &obj, label len=0)
Global flatOutput function.
Definition: FlatOutput.H:85
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::direction
uint8_t direction
Definition: direction.H:47
x
x
Definition: LISASMDCalcMethod2.H:52
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:392
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 with label keys and label hasher.
Definition: HashSet.H:410
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:122
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:298
Foam::PDRblock::read
bool read(const dictionary &dict)
Read dictionary.
Definition: PDRblock.C:556
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:417
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()
Construct zero-size.
Definition: PDRblockI.H:30
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:177