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 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 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  //- Calculate geometric expansion factor from expansion ratio
37  inline scalar calcGexp(const scalar expRatio, const label nDiv)
38  {
39  return nDiv > 1 ? pow(expRatio, 1.0/(nDiv - 1)) : 0.0;
40  }
41 }
42 
43 
44 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
45 
46 bool Foam::PDRblock::checkMonotonic
47 (
48  const direction cmpt,
49  const UList<scalar>& pts
50 )
51 {
52  const label len = pts.size();
53 
54  if (!len)
55  {
56  return false;
57  }
58 
59  const scalar& minVal = pts[0];
60 
61  for (label i=1; i < len; ++i)
62  {
63  if (pts[i] <= minVal)
64  {
66  << "Points in " << vector::componentNames[cmpt]
67  << " direction do not increase monotonically" << nl
68  << flatOutput(pts) << nl << nl
69  << exit(FatalError);
70  }
71  }
72 
73  return true;
74 }
75 
76 
78 {
79  return NullObjectRef<PDRblock>();
80 }
81 
82 
83 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
84 
85 void Foam::PDRblock::adjustSizes()
86 {
87  // Adjust i-j-k addressing
88  sizes().x() = grid_.x().nCells();
89  sizes().y() = grid_.y().nCells();
90  sizes().z() = grid_.z().nCells();
91 
92  if (sizes().x() <= 0 || sizes().y() <= 0 || sizes().z() <= 0)
93  {
94  // Sanity check. Silently disallow bad sizing
96 
97  grid_.x().clear();
98  grid_.y().clear();
99  grid_.z().clear();
100 
101  bounds_ = boundBox::invertedBox;
102  minEdgeLen_ = Zero;
103  return;
104  }
105 
106  // Adjust boundBox
107  bounds_.min().x() = grid_.x().first();
108  bounds_.min().y() = grid_.y().first();
109  bounds_.min().z() = grid_.z().first();
110 
111  bounds_.max().x() = grid_.x().last();
112  bounds_.max().y() = grid_.y().last();
113  bounds_.max().z() = grid_.z().last();
114 
115  // Min edge length
116  minEdgeLen_ = GREAT;
117 
118  for (direction cmpt=0; cmpt < vector::nComponents; ++cmpt)
119  {
120  const label nEdge = grid_[cmpt].nCells();
121 
122  for (label edgei=0; edgei < nEdge; ++edgei)
123  {
124  const scalar len = grid_[cmpt].width(edgei);
125  minEdgeLen_ = min(minEdgeLen_, len);
126  }
127  }
128 }
129 
130 
131 void Foam::PDRblock::readGridControl
132 (
133  const direction cmpt,
134  const dictionary& dict,
135  const scalar scaleFactor
136 )
137 {
138  //- The begin/end nodes for each segment
139  scalarList knots;
140 
141  // The number of division per segment
142  labelList divisions;
143 
144  // The expansion ratio per segment
145  scalarList expansion;
146 
147  if (verbose_)
148  {
149  Info<< "Reading grid control for "
150  << vector::componentNames[cmpt] << " direction" << nl;
151  }
152 
153  // Points
154  // ~~~~~~
155 
156  dict.readEntry("points", knots);
157 
158  const label nSegments = (knots.size()-1);
159 
160  if (nSegments < 1)
161  {
163  << "Must define at least two control points:"
164  << " in " << vector::componentNames[cmpt]
165  << " direction" << nl
166  << exit(FatalIOError);
167  }
168 
169  // Apply scaling
170  if (scaleFactor > 0)
171  {
172  for (scalar& pt : knots)
173  {
174  pt *= scaleFactor;
175  }
176  }
177 
178  // Do points increase monotonically?
179  checkMonotonic(cmpt, knots);
180 
181  if (verbose_)
182  {
183  Info<< " points : " << flatOutput(knots) << nl;
184  }
185 
186 
187  // Divisions
188  // ~~~~~~~~~
189 
190  dict.readEntry("nCells", divisions);
191 
192  label nTotalDivs = 0;
193  for (const label ndiv : divisions)
194  {
195  nTotalDivs += ndiv;
196 
197  if (ndiv < 1)
198  {
200  << "Negative or zero divisions:"
201  << " in " << vector::componentNames[cmpt]
202  << " direction" << nl
203  << exit(FatalIOError);
204  }
205  }
206 
207  if (divisions.size() != nSegments)
208  {
210  << "Mismatch in number of divisions and number of points:"
211  << " in " << vector::componentNames[cmpt]
212  << " direction" << nl
213  << exit(FatalIOError);
214  }
215  else if (!nTotalDivs)
216  {
218  << "No grid divisions at all:"
219  << " in " << vector::componentNames[cmpt]
220  << " direction" << nl
221  << exit(FatalIOError);
222  }
223 
224 
225  if (verbose_)
226  {
227  Info<< " nCells : " << flatOutput(divisions) << nl;
228  }
229 
230 
231  // Expansion ratios
232  // ~~~~~~~~~~~~~~~~
233 
234  dict.readIfPresent("ratios", expansion);
235 
236  // Correct expansion ratios - negative is the same as inverse.
237  for (scalar& expRatio : expansion)
238  {
239  if (expRatio < 0)
240  {
241  expRatio = 1.0/(-expRatio);
242  }
243  }
244 
245  if (expansion.size() && expansion.size() != nSegments)
246  {
248  << "Mismatch in number of expansion ratios and divisions:"
249  << " in " << vector::componentNames[cmpt]
250  << " direction" << nl
251  << exit(FatalIOError);
252  }
253 
254  if (expansion.empty())
255  {
256  expansion.resize(nSegments, scalar(1));
257 
258  if (verbose_)
259  {
260  Info<< "Warning: 'ratios' not specified, using constant spacing"
261  << nl;
262  }
263  }
264  else
265  {
266  if (verbose_)
267  {
268  Info<< " ratios : " << flatOutput(expansion) << nl;
269  }
270  }
271 
272 
273 
274  // Define the grid points
275 
276  scalarList& gridPoint = grid_[cmpt];
277  gridPoint.resize(nTotalDivs+1);
278 
279  label start = 0;
280 
281  for (label segmenti=0; segmenti < nSegments; ++segmenti)
282  {
283  const label nDiv = divisions[segmenti];
284  const scalar expRatio = expansion[segmenti];
285 
286  SubList<scalar> subPoint(gridPoint, nDiv+1, start);
287  start += nDiv;
288 
289  subPoint[0] = knots[segmenti];
290  subPoint[nDiv] = knots[segmenti+1];
291 
292  const scalar dist = (subPoint.last() - subPoint.first());
293 
294  if (equal(expRatio, 1))
295  {
296  for (label i=1; i < nDiv; ++i)
297  {
298  subPoint[i] = (subPoint[0] + (dist * i)/nDiv);
299  }
300  }
301  else
302  {
303  // Geometric expansion factor from the expansion ratio
304  const scalar expFact = calcGexp(expRatio, nDiv);
305 
306  for (label i=1; i < nDiv; ++i)
307  {
308  subPoint[i] =
309  (
310  subPoint[0]
311  + dist * (1.0 - pow(expFact, i))/(1.0 - pow(expFact, nDiv))
312  );
313  }
314  }
315  }
316 }
317 
318 
319 void Foam::PDRblock::readBoundary(const dictionary& dict)
320 {
321  if (verbose_)
322  {
323  Info<< "Reading boundary entries" << nl;
324  }
325 
326  PtrList<entry> patchEntries;
327 
328  if (dict.found("boundary"))
329  {
330  PtrList<entry> newEntries(dict.lookup("boundary"));
331  patchEntries.transfer(newEntries);
332  }
333 
334  patches_.clear();
335  patches_.resize(patchEntries.size());
336 
337  // Hex cell has 6 sides
338  const labelRange validFaceId(0, 6);
339 
340  // Track which sides have been handled
341  labelHashSet handled;
342 
343  forAll(patchEntries, patchi)
344  {
345  if (!patchEntries.set(patchi))
346  {
347  continue;
348  }
349 
350  const entry& e = patchEntries[patchi];
351 
352  if (!e.isDict())
353  {
355  << "Entry " << e
356  << " in boundary section is not a valid dictionary."
357  << exit(FatalIOError);
358  }
359 
360  const dictionary& patchDict = e.dict();
361 
362  const word& patchName = e.keyword();
363 
364  const word patchType(patchDict.get<word>("type"));
365 
366  labelList faceLabels(patchDict.get<labelList>("faces"));
367 
368  labelHashSet patchFaces(faceLabels);
369 
370  if (faceLabels.size() != patchFaces.size())
371  {
373  << "Patch: " << patchName
374  << ": Ignoring duplicate face ids" << nl;
375  }
376 
377  label nErased = patchFaces.filterKeys(validFaceId);
378  if (nErased)
379  {
381  << "Patch: " << patchName << ": Ignoring " << nErased
382  << " faces with invalid identifiers" << nl;
383  }
384 
385  nErased = patchFaces.erase(handled);
386  if (nErased)
387  {
389  << "Patch: " << patchName << ": Ignoring " << nErased
390  << " faces, which were already handled" << nl;
391  }
392 
393  // Mark these faces as having been handled
394  handled += patchFaces;
395 
396 
397  // Save information for later access during mesh creation.
398 
399  patches_.set(patchi, new boundaryEntry());
400 
401  boundaryEntry& bentry = patches_[patchi];
402 
403  bentry.name_ = patchName;
404  bentry.type_ = patchType;
405  bentry.size_ = 0;
406  bentry.faces_ = patchFaces.sortedToc();
407 
408  // Count patch faces
409  for (const label shapeFacei : bentry.faces_)
410  {
411  bentry.size_ += nBoundaryFaces(shapeFacei);
412  }
413  }
414 
415 
416  // Deal with unhandled faces
417 
418  labelHashSet missed(identity(6));
419  missed.erase(handled);
420 
421  if (missed.size())
422  {
423  patches_.append(new boundaryEntry());
424 
425  boundaryEntry& bentry = patches_.last();
426 
427  bentry.name_ = "defaultFaces";
428  bentry.type_ = emptyPolyPatch::typeName;
429  bentry.size_ = 0;
430  bentry.faces_ = missed.sortedToc();
431 
432  // Count patch faces
433  for (const label shapeFacei : bentry.faces_)
434  {
435  bentry.size_ += nBoundaryFaces(shapeFacei);
436  }
437 
438 
439  // Use name/type from "defaultPatch" entry if available
440  // - be generous with error handling
441 
442  const dictionary* defaultEntry = dict.findDict("defaultPatch");
443  if (defaultEntry)
444  {
445  defaultEntry->readIfPresent("name", bentry.name_);
446  defaultEntry->readIfPresent("type", bentry.type_);
447  }
448  }
449 
450  if (verbose_)
451  {
452  label patchi = 0;
453  for (const boundaryEntry& bentry : patches_)
454  {
455  Info<< " patch: " << patchi
456  << " (size: " << bentry.size_
457  << " type: " << bentry.type_
458  << ") name: " << bentry.name_
459  << " faces: " << flatOutput(bentry.faces_) << nl;
460 
461  ++patchi;
462  }
463  }
464 }
465 
466 
467 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
468 
470 (
471  const UList<scalar>& xgrid,
472  const UList<scalar>& ygrid,
473  const UList<scalar>& zgrid
474 )
475 :
476  PDRblock()
477 {
478  // Default boundaries with patchi == shapeFacei
479  patches_.resize(6);
480  for (label patchi=0; patchi < 6; ++patchi)
481  {
482  patches_.set(patchi, new boundaryEntry());
483 
484  boundaryEntry& bentry = patches_[patchi];
485 
486  bentry.name_ = "patch" + Foam::name(patchi);
487  bentry.type_ = "patch";
488  bentry.size_ = 0;
489  bentry.faces_ = labelList(one(), patchi);
490  }
491 
492  reset(xgrid, ygrid, zgrid);
493 }
494 
495 
497 :
498  PDRblock()
499 {
500  verbose_ = verbose;
501 
502  read(dict);
503 }
504 
505 
506 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
507 
509 {
510  // Mark no scaling with invalid value
511  const scalar scaleFactor = dict.lookupOrDefault<scalar>("scale", -1);
512 
513  readGridControl(0, dict.subDict("x"), scaleFactor);
514  readGridControl(1, dict.subDict("y"), scaleFactor);
515  readGridControl(2, dict.subDict("z"), scaleFactor);
516 
517  adjustSizes();
518 
519  readBoundary(dict);
520 
521  return true;
522 }
523 
524 
526 (
527  const UList<scalar>& xgrid,
528  const UList<scalar>& ygrid,
529  const UList<scalar>& zgrid
530 )
531 {
532  static_cast<scalarList&>(grid_.x()) = xgrid;
533  static_cast<scalarList&>(grid_.y()) = ygrid;
534  static_cast<scalarList&>(grid_.z()) = zgrid;
535 
536  #ifdef FULLDEBUG
537  for (direction cmpt=0; cmpt < vector::nComponents; ++cmpt)
538  {
539  checkMonotonic(cmpt, grid_[cmpt]);
540  }
541  #endif
542 
543  adjustSizes();
544 
545  // Adjust boundaries
546  for (boundaryEntry& bentry : patches_)
547  {
548  bentry.size_ = 0;
549 
550  // Count patch faces
551  for (const label shapeFacei : bentry.faces_)
552  {
553  bentry.size_ += nBoundaryFaces(shapeFacei);
554  }
555  }
556 }
557 
558 
559 bool Foam::PDRblock::findCell(const point& pt, labelVector& pos) const
560 {
561  // Out-of-bounds is handled explicitly, for efficiency and consistency,
562  // but principally to ensure that findLower() returns a valid
563  // result when the point is to the right of the bounds.
564 
565  // Since findLower returns the lower index, it corresponds to the
566  // cell in which the point is found
567 
568  if (!bounds_.contains(pt))
569  {
570  return false;
571  }
572 
573  for (direction cmpt=0; cmpt < labelVector::nComponents; ++cmpt)
574  {
575  // Binary search
576  pos[cmpt] = findLower(grid_[cmpt], pt[cmpt]);
577  }
578 
579  return true;
580 }
581 
582 
583 bool Foam::PDRblock::gridIndex
584 (
585  const point& pt,
586  labelVector& pos,
587  const scalar relTol
588 ) const
589 {
590  const scalar tol = relTol * minEdgeLen_;
591 
592  for (direction cmpt=0; cmpt < labelVector::nComponents; ++cmpt)
593  {
594  // Linear search
595  pos[cmpt] = grid_[cmpt].findIndex(pt[cmpt], tol);
596 
597  if (pos[cmpt] < 0) return false;
598  }
599 
600  return true;
601 }
602 
603 
604 Foam::labelVector Foam::PDRblock::findCell(const point& pt) const
605 {
607 
608  if (findCell(pt, pos))
609  {
610  return pos;
611  }
612 
613  return labelVector(-1,-1,-1);
614 }
615 
616 
617 Foam::labelVector Foam::PDRblock::gridIndex
618 (
619  const point& pt,
620  const scalar relTol
621 ) const
622 {
624 
625  if (gridIndex(pt, pos, relTol))
626  {
627  return pos;
628  }
629 
630  return labelVector(-1,-1,-1);
631 }
632 
633 
634 // ************************************************************************* //
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:503
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::scalarList
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
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.
Definition: zero.H:128
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:359
Foam::PDRblock::null
static const PDRblock & null()
Return a PDRblock reference to a nullObject.
Definition: PDRblock.C:77
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), which can be used to avoid manipulating objects that are...
Definition: one.H:60
Foam::FatalIOError
IOerror FatalIOError
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.H:1241
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
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
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:526
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
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:523
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:419
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:141
emptyPolyPatch.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::nl
constexpr char nl
Definition: Ostream.H:372
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::start
label ListType::const_reference const label start
Definition: ListOps.H:408
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:375
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:415
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:294
Foam::PDRblock::read
bool read(const dictionary &dict)
Read dictionary.
Definition: PDRblock.C:508
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