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