PDRblock.H
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 Class
27  Foam::PDRblock
28 
29 Description
30  A single block x-y-z rectilinear mesh addressable as i,j,k with
31  simplified creation. Some of the input is similar to blockMeshDict,
32  but since this specialization is for a single-block that is aligned
33  with the x-y-z directions, it provides a different means of specifying
34  the mesh.
35 
36  Dictionary controls
37  \table
38  Property | Description | Required | Default
39  x | X-direction grid specification | yes |
40  y | Y-direction grid specification | yes |
41  z | Z-direction grid specification | yes |
42  scale | Point scaling | no | 1.0
43  expansion | Type of expansion (ratio/relative) | no | ratio
44  boundary | Boundary patches | yes |
45  defaultPatch | Default patch specification | no |
46  \endtable
47 
48  Grid coordinate controls
49  \table
50  Property| Description | Required | Default
51  points | Locations defining the mesh segment | yes |
52  nCells | Divisions per mesh segment | yes |
53  ratios | Expansion values per segment | no |
54  \endtable
55 
56  A negative expansion value is trapped and treated as its reciprocal.
57  by default, the expansion is as per blockMesh and represents the ratio
58  of end-size / start-size for the section.
59  Alternatively, the relative size can be given.
60 
61 SourceFiles
62  PDRblockI.H
63  PDRblock.C
64  PDRblockCreate.C
65 
66 \*---------------------------------------------------------------------------*/
67 
68 #ifndef Foam_PDRblock_H
69 #define Foam_PDRblock_H
70 
71 #include "ijkMesh.H"
72 #include "boundBox.H"
73 #include "pointField.H"
74 #include "faceList.H"
75 #include "Enum.H"
76 #include "vector2D.H"
77 #include "labelVector2D.H"
78 
79 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
80 
81 namespace Foam
82 {
83 
84 // Forward Declarations
85 class IOobject;
86 class blockMesh;
87 class polyMesh;
88 class gradingDescriptors;
89 
90 /*---------------------------------------------------------------------------*\
91  Class PDRblock Declaration
92 \*---------------------------------------------------------------------------*/
93 
94 class PDRblock
95 :
96  public ijkMesh
97 {
98 public:
99 
100  // Data Types
101 
102  //- The expansion type
103  enum expansionType : uint8_t
104  {
106  EXPAND_RATIO,
108  };
109 
110  //- Named enumerations for the expansion type
111  const static Enum<expansionType> expansionNames_;
112 
113 
114  // Public Classes
115 
116  //- Grid locations in an axis direction.
117  // The number of points is one larger than the number of elements
118  // it represents
119  class location
120  :
121  public scalarList
122  {
123  public:
124 
125  //- Reset with min/max range and number of divisions
126  void reset(const scalar low, const scalar upp, const label nCells);
127 
128  //- The location list is valid if it contains 2 or more points
129  inline bool valid() const;
130 
131  //- The number of cells in this direction.
132  inline label nCells() const;
133 
134  //- The number of points in this direction.
135  inline label nPoints() const;
136 
137  //- True if the location is within the range
138  inline bool contains(const scalar p) const;
139 
140  //- The first() value is considered the min value.
141  inline const scalar& min() const;
142 
143  //- The last() value is considered the max value.
144  inline const scalar& max() const;
145 
146  //- Mid-point location, zero for an empty list.
147  inline scalar centre() const;
148 
149  //- The difference between min/max values, zero for an empty list.
150  inline scalar length() const;
151 
152  //- Check that element index is within valid range.
153  inline void checkIndex(const label i) const;
154 
155  //- Cell size at element position.
156  inline scalar width(const label i) const;
157 
158  //- Cell centre at element position.
159  // Treats -1 and nCells positions like a halo cell.
160  inline scalar C(const label i) const;
161 
162  //- Return min/max edge lengths
165  //- Return edge grading descriptors for the locations
166  // \see Foam::gradingDescriptor
167  gradingDescriptors grading() const;
168 
169  //- Find the cell index enclosing this location
170  // \return -1 for out-of-bounds
171  label findCell(const scalar p) const;
172 
173  //- Find the grid index, within the given tolerance
174  // Return -1 for out-of-bounds and -2 for a point that is
175  // within bounds, but not aligned with a grid point.
176  label findIndex(const scalar p, const scalar tol) const;
177 
178  //- If out of range, return the respective min/max limits,
179  //- otherwise return the value itself.
180  // If the range is invalid, always return the value.
181  inline const scalar& clip(const scalar& val) const;
182  };
183 
184 
185  //- The begin/end nodes for each segment,
186  //- with divisions and expansion for each segment
187  // Not normally used outside of PDRblock
188  struct gridControl
189  :
190  public scalarList
191  {
192  //- The number of division per segment
194 
195  //- The expansion ratio per segment
197 
198  //- Total number of cells in this direction
199  label nCells() const;
200 
201  //- Return edge grading descriptors for the locations
202  // \see Foam::gradingDescriptor
203  gradingDescriptors grading() const;
204 
205  //- Resize lists
206  void resize(label len);
207 
208  //- Add point/divisions/expand to end of list (push_back)
209  void append(const scalar p, label nDiv, scalar expRatio=1);
210 
211  //- Add point/divisions/expand to front of list (push_front)
212  void prepend(const scalar p, label nDiv, scalar expRatio=1);
213 
214  //- Write as dictionary contents for specified vector direction
215  void writeDict(Ostream& os, const direction cmpt) const;
216  };
217 
218 
219 private:
220 
221  // Private Classes
222 
223  //- Extracted patch settings
224  struct boundaryEntry
225  {
226  //- The patch name
227  word name_;
228 
229  //- The patch type
230  word type_;
231 
232  //- The patch size
233  label size_;
234 
235  //- The associated block face ids [0,5]
236  labelList faces_;
237  };
238 
239  //- The begin/end nodes for each segment,
240  //- with divisions and expansion for each segment
241  struct outerControl
242  {
243  //- The control type
244  enum controlType : uint8_t
245  {
246  OUTER_NONE = 0,
247  OUTER_EXTEND,
248  OUTER_BOX,
249  OUTER_SPHERE
250  };
251 
252  //- Named enumerations for the control type
253  const static Enum<controlType> controlNames_;
254 
255  //- The control type
256  controlType type_;
257 
258  //- The expansion type
259  expansionType expandType_;
260 
261  //- True if on the ground
262  bool onGround_;
263 
264  //- Relative size(s) for the outer region
265  vector2D relSize_;
266 
267  //- Number of cells in outer region
268  // Generally only single component is used
269  labelVector2D nCells_;
270 
271  //- Expansion ratio(s) for the outer region
272  vector2D expansion_;
273 
274 
275  // Constructors
276 
277  //- Default construct. NONE
278  outerControl();
279 
280 
281  // Member Functions
282 
283  //- Reset to default (NONE) values
284  void clear();
285 
286  //- Is enabled (not NONE)
287  bool active() const;
288 
289  //- Project on to sphere (is SPHERE)
290  bool isSphere() const;
291 
292  //- Is the outer region on the ground?
293  bool onGround() const;
294 
295  //- Define that the outer region is on the ground or not
296  // \return the old value
297  bool onGround(const bool on);
298 
299  //- Read content from dictionary
300  void read(const dictionary& dict);
301 
302  //- Report information about outer region
303  void report(Ostream& os) const;
304  };
305 
306 
307  // Private Data
308 
309  //- Reference to mesh dictionary
310  const dictionary& meshDict_;
311 
312  //- The grid controls in (i,j,k / x,y,z) directions.
313  Vector<gridControl> control_;
314 
315  //- The grid points in all (i,j,k / x,y,z) directions,
316  //- after applying the internal subdivisions.
317  Vector<location> grid_;
318 
319  //- Control for the outer-region (if any)
320  outerControl outer_;
321 
322  //- The mesh bounding box
323  boundBox bounds_;
324 
325  //- The boundary patch information
326  PtrList<boundaryEntry> patches_;
327 
328  //- The min/max edge lengths
329  scalarMinMax edgeLimits_;
330 
331  //- Verbosity
332  bool verbose_;
333 
334 
335  // Private Member Functions
336 
337  //- Check that points increase monotonically
338  static bool checkMonotonic
339  (
340  const direction cmpt,
341  const UList<scalar>& pts
342  );
343 
344  //- Add zero-sized patches with default naming/types
345  void addDefaultPatches();
346 
347  //- Adjust sizing for updated grid points
348  void adjustSizes();
349 
350  //- Read and define grid points in given direction
351  void readGridControl
352  (
353  const direction cmpt,
354  const dictionary& dict,
355  const scalar scaleFactor = -1,
356  expansionType expandType = expansionType::EXPAND_RATIO
357  );
358 
359  //- Read "boundary" information
360  void readBoundary(const dictionary& dict);
361 
362  //- Populate point field for the block
363  void createPoints(pointField& pts) const;
364 
365  //- Add internal faces to lists.
366  // Lists must be properly sized!
367  // \return the number of faces added
368  label addInternalFaces
369  (
370  faceList::iterator& faceIter,
371  labelList::iterator& ownIter,
372  labelList::iterator& neiIter
373  ) const;
374 
375  //- Add boundary faces for the shape face to lists
376  // Lists must be properly sized!
377  // \return the number of faces added
378  label addBoundaryFaces
379  (
380  const direction shapeFacei,
381  faceList::iterator& faceIter,
382  labelList::iterator& ownIter
383  ) const;
384 
385  //- Obtain i,j,k index for cell enclosing this location
386  // \return false for out-of-bounds
387  bool findCell(const point& pt, labelVector& pos) const;
388 
389  //- Obtain i,j,k grid index for point location
390  // \return false for out-of-bounds and off-grid
391  bool gridIndex
392  (
393  const point& pt,
394  labelVector& pos,
395  const scalar tol
396  ) const;
397 
398  //- The bounding box of the grid points
399  static boundBox bounds
400  (
401  const scalarList& x,
402  const scalarList& y,
403  const scalarList& z
404  );
405 
406  //- Equivalent edge grading descriptors in (x,y,z) directions.
408  (
409  const Vector<gridControl>& ctrl
410  );
411 
412  //- Mesh sizes based on the controls
413  static labelVector sizes
414  (
415  const Vector<gridControl>& ctrl
416  );
417 
418 
419  // Mesh Generation
420 
421  //- Create a blockMesh
422  autoPtr<blockMesh> createBlockMesh(const IOobject& io) const;
423 
424  //- Create polyMesh via blockMesh
425  autoPtr<polyMesh> meshBlockMesh(const IOobject& io) const;
426 
427 
428 public:
429 
430  // Static Member Functions
431 
432  //- Return a PDRblock reference to a nullObject
433  static const PDRblock& null();
434 
435 
436  // Constructors
437 
438  //- Default construct, zero-size, inverted bounds etc
439  PDRblock();
440 
441  //- Construct from components
442  PDRblock
443  (
444  const UList<scalar>& xgrid,
445  const UList<scalar>& ygrid,
446  const UList<scalar>& zgrid
447  );
448 
449  //- Construct from cube with specified griding
450  PDRblock(const boundBox& box, const labelVector& nCells);
451 
452  //- Construct from dictionary
453  explicit PDRblock(const dictionary& dict, bool verboseOutput=false);
454 
455 
456  // Member Functions
457 
458  //- Read dictionary
459  bool read(const dictionary& dict);
460 
461  //- Reset grid locations and mesh i-j-k sizing
462  void reset
463  (
464  const UList<scalar>& xgrid,
465  const UList<scalar>& ygrid,
466  const UList<scalar>& zgrid
467  );
468 
469  //- Reset cube and mesh i-j-k sizing
470  void reset(const boundBox& box, const labelVector& nCells);
471 
472 
473  // Access
474 
475  //- The grid point locations in the i,j,k (x,y,z) directions.
476  const Vector<location>& grid() const noexcept { return grid_; }
477 
478  //- Equivalent edge grading descriptors in (x,y,z) directions.
480 
481  //- Equivalent edge grading descriptors in specified (x,y,z) direction.
482  gradingDescriptors grading(const direction cmpt) const;
483 
484 
485  // Mesh Information
486 
487  //- Mesh sizing as per ijkMesh
488  using ijkMesh::sizes;
489 
490  //- The mesh bounding box
491  const boundBox& bounds() const noexcept { return bounds_; }
492 
493  //- The min/max edge length
494  const scalarMinMax& edgeLimits() const noexcept { return edgeLimits_; }
495 
496  //- Cell size in x-direction at i position.
497  inline scalar dx(const label i) const;
498 
499  //- Cell size in x-direction at i position.
500  inline scalar dx(const labelVector& ijk) const;
501 
502  //- Cell size in y-direction at j position.
503  inline scalar dy(const label j) const;
504 
505  //- Cell size in y-direction at j position.
506  inline scalar dy(const labelVector& ijk) const;
507 
508  //- Cell size in z-direction at k position.
509  inline scalar dz(const label k) const;
510 
511  //- Cell size in z-direction at k position.
512  inline scalar dz(const labelVector& ijk) const;
513 
514  //- Cell dimensions at i,j,k position.
515  inline vector span(const label i, const label j, const label k) const;
516 
517  //- Cell dimensions at i,j,k position.
518  inline vector span(const labelVector& ijk) const;
519 
520  //- Grid point at i,j,k position.
521  inline point grid(const label i, const label j, const label k) const;
522 
523  //- Grid point at i,j,k position.
524  inline point grid(const labelVector& ijk) const;
525 
526  //- Cell centre at i,j,k position.
527  inline point C(const label i, const label j, const label k) const;
528 
529  //- Cell centre at i,j,k position.
530  inline point C(const labelVector& ijk) const;
531 
532  //- Cell volume at i,j,k position.
533  inline scalar V(const label i, const label j, const label k) const;
534 
535  //- Cell volume at i,j,k position.
536  inline scalar V(const labelVector& ijk) const;
537 
538  //- Characteristic cell size at i,j,k position.
539  // This is the cubic root of the volume
540  inline scalar width(const label i, const label j, const label k) const;
541 
542  //- Characteristic cell size at i,j,k position.
543  // This is the cubic root of the volume
544  inline scalar width(const labelVector& ijk) const;
545 
546 
547  // Searching
548 
549  //- Return i,j,k index for cell enclosing this location
550  // The value (-1,-1,-1) is returned for out-of-bounds (not found).
551  labelVector findCell(const point& pt) const;
552 
553  //- Obtain i,j,k grid index for point location within specified
554  // relative tolerance of the min edge length
555  // The value (-1,-1,-1) is returned for out-of-bounds (not found).
556  // and off-grid
557  labelVector gridIndex(const point& pt, const scalar relTol=0.01) const;
558 
559 
560  // Mesh Generation
561 
562  //- Output content for an equivalent blockMeshDict
563  // Optionally generate header/footer content
564  Ostream& blockMeshDict(Ostream& os, const bool withHeader=false) const;
565 
566  //- Content for an equivalent blockMeshDict
567  dictionary blockMeshDict() const;
568 
569  //- Write an equivalent blockMeshDict
570  void writeBlockMeshDict(const IOobject& io) const;
571 
572  //- Create polyMesh for grid definition and patch information
573  autoPtr<polyMesh> mesh(const IOobject& io) const;
574 
575  //- Create polyMesh for inner-mesh only,
576  //- ignore any outer block definitions
577  autoPtr<polyMesh> innerMesh(const IOobject& io) const;
578 };
579 
580 
581 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
582 
583 } // End namespace Foam
584 
585 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
586 
587 #include "PDRblockI.H"
588 
589 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
590 
591 #endif
592 
593 // ************************************************************************* //
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
const scalarMinMax & edgeLimits() const noexcept
The min/max edge length.
Definition: PDRblock.H:739
scalar centre() const
Mid-point location, zero for an empty list.
Definition: PDRblockI.H:60
const scalar & max() const
The last() value is considered the max value.
Definition: PDRblockI.H:54
dictionary dict
point C(const label i, const label j, const label k) const
Cell centre at i,j,k position.
Definition: PDRblockI.H:217
uint8_t direction
Definition: direction.H:48
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition: MinMax.H:97
label nPoints() const
The number of points in this direction.
Definition: PDRblockI.H:36
const scalar & clip(const scalar &val) const
If out of range, return the respective min/max limits, otherwise return the value itself...
Definition: PDRblockI.H:120
label findIndex(const scalar p, const scalar tol) const
Find the grid index, within the given tolerance.
void resize(label len)
Resize lists.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
scalar V(const label i, const label j, const label k) const
Cell volume at i,j,k position.
Definition: PDRblockI.H:240
void prepend(const scalar p, label nDiv, scalar expRatio=1)
Add point/divisions/expand to front of list (push_front)
Uniform expansion (ie, no expansion)
Definition: PDRblock.H:162
void checkIndex(const label i) const
Check that element index is within valid range.
Definition: PDRblockI.H:72
label nCells() const
The number of mesh cells (nx*ny*nz) in the i-j-k mesh.
Definition: ijkMeshI.H:61
const scalar & min() const
The first() value is considered the min value.
Definition: PDRblockI.H:48
Vector2D< scalar > vector2D
A 2D vector of scalars obtained from the generic Vector2D.
Definition: vector2D.H:56
vector span(const label i, const label j, const label k) const
Cell dimensions at i,j,k position.
Definition: PDRblockI.H:177
label * iterator
Random access iterator for traversing a UList.
Definition: UList.H:172
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
label k
Boltzmann constant.
expansionType
The expansion type.
Definition: PDRblock.H:160
dictionary blockMeshDict() const
Content for an equivalent blockMeshDict.
void writeBlockMeshDict(const IOobject &io) const
Write an equivalent blockMeshDict.
bool read(const dictionary &dict)
Read dictionary.
Definition: PDRblock.C:569
dimensionedScalar pos(const dimensionedScalar &ds)
scalar y
Vector2D< label > labelVector2D
A 2D vector of labels obtained from the generic Vector2D.
Definition: labelVector2D.H:46
PDRblock()
Default construct, zero-size, inverted bounds etc.
Definition: PDRblock.C:520
bool contains(const scalar p) const
True if the location is within the range.
Definition: PDRblockI.H:42
The begin/end nodes for each segment, with divisions and expansion for each segment.
Definition: PDRblock.H:292
scalar dz(const label k) const
Cell size in z-direction at k position.
Definition: PDRblockI.H:164
A class for handling words, derived from Foam::string.
Definition: word.H:63
static const Enum< expansionType > expansionNames_
Named enumerations for the expansion type.
Definition: PDRblock.H:170
const labelVector & sizes() const
The (i,j,k) addressing dimensions.
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
Vector< scalar > vector
Definition: vector.H:57
autoPtr< polyMesh > innerMesh(const IOobject &io) const
Create polyMesh for inner-mesh only, ignore any outer block definitions.
scalar width(const label i, const label j, const label k) const
Characteristic cell size at i,j,k position.
Definition: PDRblockI.H:257
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:99
scalarList expansion_
The expansion ratio per segment.
Definition: PDRblock.H:304
Templated 3D Vector derived from VectorSpace adding construction from 3 components, element access using x(), y() and z() member functions and the inner-product (dot-product) and cross-product operators.
Definition: Vector.H:58
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
const direction noexcept
Definition: Scalar.H:258
static const PDRblock & null()
Return a PDRblock reference to a nullObject.
Definition: PDRblock.C:106
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Vector< gradingDescriptors > grading() const
Equivalent edge grading descriptors in (x,y,z) directions.
Definition: PDRblock.C:736
scalar length() const
The difference between min/max values, zero for an empty list.
Definition: PDRblockI.H:66
void append(const scalar p, label nDiv, scalar expRatio=1)
Add point/divisions/expand to end of list (push_back)
OBJstream os(runTime.globalPath()/outputName)
label nCells() const
The number of cells in this direction.
Definition: PDRblockI.H:30
scalar dy(const label j) const
Cell size in y-direction at j position.
Definition: PDRblockI.H:152
scalarMinMax edgeLimits() const
Return min/max edge lengths.
label findCell(const scalar p) const
Find the cell index enclosing this location.
scalar width(const label i) const
Cell size at element position.
Definition: PDRblockI.H:84
void clear()
Reset to (0,0,0) sizing.
vector point
Point is a vector.
Definition: point.H:37
void reset(const UList< scalar > &xgrid, const UList< scalar > &ygrid, const UList< scalar > &zgrid)
Reset grid locations and mesh i-j-k sizing.
Definition: PDRblock.C:606
label nCells() const
Total number of cells in this direction.
void writeDict(Ostream &os, const direction cmpt) const
Write as dictionary contents for specified vector direction.
labelList divisions_
The number of division per segment.
Definition: PDRblock.H:299
End/start ratio.
Definition: PDRblock.H:163
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
bool valid() const
The location list is valid if it contains 2 or more points.
Definition: PDRblockI.H:24
Vector< label > labelVector
Vector of labels.
Definition: labelVector.H:47
gradingDescriptors grading() const
Return edge grading descriptors for the locations.
volScalarField & p
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
autoPtr< polyMesh > mesh(const IOobject &io) const
Create polyMesh for grid definition and patch information.
scalar C(const label i) const
Cell centre at element position.
Definition: PDRblockI.H:94
List of gradingDescriptor for the sections of a block with additional IO functionality.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:171
gradingDescriptors grading() const
Return edge grading descriptors for the locations.
const boundBox & bounds() const noexcept
The mesh bounding box.
Definition: PDRblock.H:734
scalar dx(const label i) const
Cell size in x-direction at i position.
Definition: PDRblockI.H:140
void reset(const scalar low, const scalar upp, const label nCells)
Reset with min/max range and number of divisions.
Namespace for OpenFOAM.
const Vector< location > & grid() const noexcept
The grid point locations in the i,j,k (x,y,z) directions.
Definition: PDRblock.H:711
Relative expansion ratio.
Definition: PDRblock.H:164
const pointField & pts