Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::PDRblock
28
29Description
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
61SourceFiles
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
81namespace Foam
82{
83
84// Forward Declarations
85class IOobject;
86class blockMesh;
87class polyMesh;
89
90/*---------------------------------------------------------------------------*\
91 Class PDRblock Declaration
92\*---------------------------------------------------------------------------*/
93
94class PDRblock
95:
96 public ijkMesh
97{
98public:
99
100 // Data Types
101
102 //- The expansion type
103 enum expansionType : uint8_t
104 {
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 good() const noexcept;
130
131 //- The number of cells in this direction.
132 inline label nCells() const noexcept;
133
134 //- The number of points in this direction.
135 inline label nPoints() const noexcept;
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 //- Return value clamped to min/max limits.
179 // If the range is invalid, always return the value.
180 inline const scalar& clamp(const scalar& val) const;
182 //- Return value clamped to min/max limits.
183 const scalar& clip(const scalar& val) const { return clamp(val); }
184 };
185
186
187 //- The begin/end nodes for each segment,
188 //- with divisions and expansion for each segment
189 // Not normally used outside of PDRblock
190 struct gridControl
191 :
192 public scalarList
193 {
194 //- The number of division per segment
196
197 //- The expansion ratio per segment
199
200 //- Total number of cells in this direction
201 label nCells() const;
202
203 //- Return edge grading descriptors for the locations
204 // \see Foam::gradingDescriptor
206
207 //- Resize lists
208 void resize(label len);
209
210 //- Add point/divisions/expand to end of list (push_back)
211 void append(const scalar p, label nDiv, scalar expRatio=1);
212
213 //- Add point/divisions/expand to front of list (push_front)
214 void prepend(const scalar p, label nDiv, scalar expRatio=1);
215
216 //- Write as dictionary contents for specified vector direction
217 void writeDict(Ostream& os, const direction cmpt) const;
218 };
219
220
221private:
222
223 // Private Classes
224
225 //- Extracted patch settings
226 struct boundaryEntry
227 {
228 //- The patch name
229 word name_;
230
231 //- The patch type
232 word type_;
233
234 //- The patch size
235 label size_;
236
237 //- The associated block face ids [0,5]
238 labelList faces_;
239 };
240
241 //- The begin/end nodes for each segment,
242 //- with divisions and expansion for each segment
243 struct outerControl
244 {
245 //- The control type
246 enum controlType : uint8_t
247 {
248 OUTER_NONE = 0,
249 OUTER_EXTEND,
250 OUTER_BOX,
251 OUTER_SPHERE
252 };
253
254 //- Named enumerations for the control type
255 const static Enum<controlType> controlNames_;
256
257 //- The control type
258 controlType type_;
260 //- The expansion type
261 expansionType expandType_;
262
263 //- True if on the ground
264 bool onGround_;
265
266 //- Relative size(s) for the outer region
267 vector2D relSize_;
268
269 //- Number of cells in outer region
270 // Generally only single component is used
271 labelVector2D nCells_;
272
273 //- Expansion ratio(s) for the outer region
274 vector2D expansion_;
275
276
277 // Constructors
278
279 //- Default construct. NONE
280 outerControl();
281
282
283 // Member Functions
284
285 //- Reset to default (NONE) values
286 void clear();
287
288 //- Is enabled (not NONE)
289 bool active() const;
290
291 //- Project on to sphere (is SPHERE)
292 bool isSphere() const;
293
294 //- Is the outer region on the ground?
295 bool onGround() const;
297 //- Define that the outer region is on the ground or not
298 // \return the old value
299 bool onGround(const bool on);
300
301 //- Read content from dictionary
302 void read(const dictionary& dict);
304 //- Report information about outer region
305 void report(Ostream& os) const;
306 };
307
309 // Private Data
310
311 //- Reference to mesh dictionary
312 const dictionary& meshDict_;
313
314 //- The grid controls in (i,j,k / x,y,z) directions.
315 Vector<gridControl> control_;
316
317 //- The grid points in all (i,j,k / x,y,z) directions,
318 //- after applying the internal subdivisions.
319 Vector<location> grid_;
320
321 //- Control for the outer-region (if any)
322 outerControl outer_;
323
324 //- The mesh bounding box
325 boundBox bounds_;
326
327 //- The boundary patch information
328 PtrList<boundaryEntry> patches_;
329
330 //- The min/max edge lengths
331 scalarMinMax edgeLimits_;
332
333 //- Verbosity
334 bool verbose_;
335
336
337 // Private Member Functions
338
339 //- Check that points increase monotonically
340 static bool checkMonotonic
341 (
342 const direction cmpt,
343 const UList<scalar>& pts
344 );
345
346 //- Add zero-sized patches with default naming/types
347 void addDefaultPatches();
348
349 //- Adjust sizing for updated grid points
350 void adjustSizes();
351
352 //- Read and define grid points in given direction
353 void readGridControl
354 (
355 const direction cmpt,
356 const dictionary& dict,
357 const scalar scaleFactor = -1,
359 );
360
361 //- Read "boundary" information
362 void readBoundary(const dictionary& dict);
363
364 //- Populate point field for the block
365 void createPoints(pointField& pts) const;
366
367 //- Add internal faces to lists.
368 // Lists must be properly sized!
369 // \return the number of faces added
370 label addInternalFaces
371 (
372 faceList::iterator& faceIter,
373 labelList::iterator& ownIter,
374 labelList::iterator& neiIter
375 ) const;
376
377 //- Add boundary faces for the shape face to lists
378 // Lists must be properly sized!
379 // \return the number of faces added
380 label addBoundaryFaces
381 (
382 const direction shapeFacei,
383 faceList::iterator& faceIter,
384 labelList::iterator& ownIter
385 ) const;
386
387 //- Obtain i,j,k index for cell enclosing this location
388 // \return false for out-of-bounds
389 bool findCell(const point& pt, labelVector& pos) const;
390
391 //- Obtain i,j,k grid index for point location
392 // \return false for out-of-bounds and off-grid
393 bool gridIndex
394 (
395 const point& pt,
397 const scalar tol
398 ) const;
399
400 //- The bounding box of the grid points
401 static boundBox bounds
402 (
403 const scalarList& x,
404 const scalarList& y,
405 const scalarList& z
406 );
407
408 //- Equivalent edge grading descriptors in (x,y,z) directions.
410 (
411 const Vector<gridControl>& ctrl
412 );
413
414 //- Mesh sizes based on the controls
415 static labelVector sizes
416 (
417 const Vector<gridControl>& ctrl
418 );
419
420
421 // Mesh Generation
422
423 //- Create a blockMesh
424 autoPtr<blockMesh> createBlockMesh(const IOobject& io) const;
425
426 //- Create polyMesh via blockMesh
427 autoPtr<polyMesh> meshBlockMesh(const IOobject& io) const;
428
429
430public:
431
432 // Static Member Functions
433
434 //- Return a null PDRblock (reference to a nullObject).
435 static const PDRblock& null() noexcept
436 {
438 }
439
440
441 // Constructors
442
443 //- Default construct, zero-size, inverted bounds etc
444 PDRblock();
445
446 //- Construct from components
448 (
449 const UList<scalar>& xgrid,
450 const UList<scalar>& ygrid,
451 const UList<scalar>& zgrid
452 );
453
454 //- Construct from cube with specified griding
455 PDRblock(const boundBox& box, const labelVector& nCells);
456
457 //- Construct from dictionary
458 explicit PDRblock(const dictionary& dict, bool verboseOutput=false);
459
460
461 // Member Functions
462
463 //- Read dictionary
464 bool read(const dictionary& dict);
465
466 //- Reset grid locations and mesh i-j-k sizing
467 void reset
468 (
469 const UList<scalar>& xgrid,
470 const UList<scalar>& ygrid,
471 const UList<scalar>& zgrid
472 );
473
474 //- Reset cube and mesh i-j-k sizing
475 void reset(const boundBox& box, const labelVector& nCells);
476
477
478 // Access
479
480 //- The grid point locations in the i,j,k (x,y,z) directions.
481 const Vector<location>& grid() const noexcept { return grid_; }
482
483 //- Equivalent edge grading descriptors in (x,y,z) directions.
485
486 //- Equivalent edge grading descriptors in specified (x,y,z) direction.
487 gradingDescriptors grading(const direction cmpt) const;
488
489
490 // Mesh Information
491
492 //- Mesh sizing as per ijkMesh
493 using ijkMesh::sizes;
494
495 //- The mesh bounding box
496 const boundBox& bounds() const noexcept { return bounds_; }
497
498 //- The min/max edge length
499 const scalarMinMax& edgeLimits() const noexcept { return edgeLimits_; }
500
501 //- Cell size in x-direction at i position.
502 inline scalar dx(const label i) const;
503
504 //- Cell size in x-direction at i position.
505 inline scalar dx(const labelVector& ijk) const;
506
507 //- Cell size in y-direction at j position.
508 inline scalar dy(const label j) const;
509
510 //- Cell size in y-direction at j position.
511 inline scalar dy(const labelVector& ijk) const;
512
513 //- Cell size in z-direction at k position.
514 inline scalar dz(const label k) const;
515
516 //- Cell size in z-direction at k position.
517 inline scalar dz(const labelVector& ijk) const;
518
519 //- Cell dimensions at i,j,k position.
520 inline vector span(const label i, const label j, const label k) const;
521
522 //- Cell dimensions at i,j,k position.
523 inline vector span(const labelVector& ijk) const;
524
525 //- Grid point at i,j,k position.
526 inline point grid(const label i, const label j, const label k) const;
527
528 //- Grid point at i,j,k position.
529 inline point grid(const labelVector& ijk) const;
530
531 //- Cell centre at i,j,k position.
532 inline point C(const label i, const label j, const label k) const;
533
534 //- Cell centre at i,j,k position.
535 inline point C(const labelVector& ijk) const;
536
537 //- Cell volume at i,j,k position.
538 inline scalar V(const label i, const label j, const label k) const;
539
540 //- Cell volume at i,j,k position.
541 inline scalar V(const labelVector& ijk) const;
542
543 //- Characteristic cell size at i,j,k position.
544 // This is the cubic root of the volume
545 inline scalar width(const label i, const label j, const label k) const;
546
547 //- Characteristic cell size at i,j,k position.
548 // This is the cubic root of the volume
549 inline scalar width(const labelVector& ijk) const;
550
551
552 // Searching
553
554 //- Return i,j,k index for cell enclosing this location
555 // The value (-1,-1,-1) is returned for out-of-bounds (not found).
556 labelVector findCell(const point& pt) const;
557
558 //- Obtain i,j,k grid index for point location within specified
559 // relative tolerance of the min edge length
560 // The value (-1,-1,-1) is returned for out-of-bounds (not found).
561 // and off-grid
562 labelVector gridIndex(const point& pt, const scalar relTol=0.01) const;
563
564
565 // Mesh Generation
566
567 //- Output content for an equivalent blockMeshDict
568 // Optionally generate header/footer content
569 Ostream& blockMeshDict(Ostream& os, const bool withHeader=false) const;
570
571 //- Content for an equivalent blockMeshDict
573
574 //- Write an equivalent blockMeshDict
575 void writeBlockMeshDict(const IOobject& io) const;
576
577 //- Create polyMesh for grid definition and patch information
578 autoPtr<polyMesh> mesh(const IOobject& io) const;
579
580 //- Create polyMesh for inner-mesh only,
581 //- ignore any outer block definitions
582 autoPtr<polyMesh> innerMesh(const IOobject& io) const;
583};
584
585
586// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
587
588} // End namespace Foam
589
590// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
591
592#include "PDRblockI.H"
593
594// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
595
596#endif
597
598// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
scalar y
label k
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Grid locations in an axis direction.
Definition PDRblock.H:184
scalar length() const
The difference between min/max values, zero for an empty list.
Definition PDRblockI.H:66
label findIndex(const scalar p, const scalar tol) const
Find the grid index, within the given tolerance.
void checkIndex(const label i) const
Check that element index is within valid range.
Definition PDRblockI.H:72
const scalar & clamp(const scalar &val) const
Return value clamped to min/max limits.
Definition PDRblockI.H:120
bool good() const noexcept
The location list is valid if it contains 2 or more points.
Definition PDRblockI.H:24
label findCell(const scalar p) const
Find the cell index enclosing this location.
scalar centre() const
Mid-point location, zero for an empty list.
Definition PDRblockI.H:60
scalarMinMax edgeLimits() const
Return min/max edge lengths.
label nPoints() const noexcept
The number of points in this direction.
Definition PDRblockI.H:36
void reset(const scalar low, const scalar upp, const label nCells)
Reset with min/max range and number of divisions.
const scalar & min() const
The first() value is considered the min value.
Definition PDRblockI.H:48
label nCells() const noexcept
The number of cells in this direction.
Definition PDRblockI.H:30
const scalar & clip(const scalar &val) const
Return value clamped to min/max limits.
Definition PDRblock.H:286
gradingDescriptors grading() const
Return edge grading descriptors for the locations.
const scalar & max() const
The last() value is considered the max value.
Definition PDRblockI.H:54
scalar width(const label i) const
Cell size at element position.
Definition PDRblockI.H:84
bool contains(const scalar p) const
True if the location is within the range.
Definition PDRblockI.H:42
A single block x-y-z rectilinear mesh addressable as i,j,k with simplified creation....
Definition PDRblock.H:152
scalar V(const label i, const label j, const label k) const
Cell volume at i,j,k position.
Definition PDRblockI.H:240
const scalarMinMax & edgeLimits() const noexcept
The min/max edge length.
Definition PDRblock.H:746
static const Enum< expansionType > expansionNames_
Named enumerations for the expansion type.
Definition PDRblock.H:170
scalar dx(const label i) const
Cell size in x-direction at i position.
Definition PDRblockI.H:140
point C(const label i, const label j, const label k) const
Cell centre at i,j,k position.
Definition PDRblockI.H:217
bool read(const dictionary &dict)
Read dictionary.
Definition PDRblock.C:563
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:600
vector span(const label i, const label j, const label k) const
Cell dimensions at i,j,k position.
Definition PDRblockI.H:177
PDRblock()
Default construct, zero-size, inverted bounds etc.
Definition PDRblock.C:514
autoPtr< polyMesh > innerMesh(const IOobject &io) const
Create polyMesh for inner-mesh only, ignore any outer block definitions.
const Vector< location > & grid() const noexcept
The grid point locations in the i,j,k (x,y,z) directions.
Definition PDRblock.H:718
scalar dy(const label j) const
Cell size in y-direction at j position.
Definition PDRblockI.H:152
const boundBox & bounds() const noexcept
The mesh bounding box.
Definition PDRblock.H:741
static const PDRblock & null() noexcept
Return a null PDRblock (reference to a nullObject).
Definition PDRblock.H:656
scalar dz(const label k) const
Cell size in z-direction at k position.
Definition PDRblockI.H:164
Vector< gradingDescriptors > grading() const
Equivalent edge grading descriptors in (x,y,z) directions.
Definition PDRblock.C:730
expansionType
The expansion type.
Definition PDRblock.H:161
@ EXPAND_RELATIVE
Relative expansion ratio.
Definition PDRblock.H:164
@ EXPAND_UNIFORM
Uniform expansion (ie, no expansion).
Definition PDRblock.H:162
@ EXPAND_RATIO
End/start ratio.
Definition PDRblock.H:163
void writeBlockMeshDict(const IOobject &io) const
Write an equivalent blockMeshDict.
scalar width(const label i, const label j, const label k) const
Characteristic cell size at i,j,k position.
Definition PDRblockI.H:257
dictionary blockMeshDict() const
Content for an equivalent blockMeshDict.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
UList(const UList< scalar > &) noexcept=default
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition Vector.H:61
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A multi-block mesh generator.
Definition blockMesh.H:164
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
List of gradingDescriptor for the sections of a block with additional IO functionality.
void checkIndex(const label i, const label j, const label k, const bool allowExtra=false) const
Check indices are within ni,nj,nk range.
const labelVector & sizes() const noexcept
The (i,j,k) addressing dimensions.
void clear()
Reset to (0,0,0) sizing.
A simple i-j-k (row-major order) to linear addressing for a rectilinear mesh. Since the underlying me...
Definition ijkMesh.H:56
ijkMesh()
Construct zero-sized.
Definition ijkMeshI.H:23
label nCells() const
The number of mesh cells (nx*ny*nz) in the i-j-k mesh.
Definition ijkMeshI.H:61
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
const auto & io
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
label findIndex(const ListType &input, typename ListType::const_reference val, const label start=0)
Deprecated(2017-10) search for first occurrence of the given element.
Definition ListOps.H:517
const T & NullObjectRef() noexcept
Const reference (of type T) to the nullObject.
Definition nullObject.H:228
List< label > labelList
A List of labels.
Definition List.H:62
Vector< label > labelVector
Vector of labels.
Definition labelVector.H:47
Vector2D< scalar > vector2D
A 2D vector of scalars obtained from the generic Vector2D.
Definition vector2D.H:56
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition MinMax.H:97
Vector2D< label > labelVector2D
A 2D vector of labels obtained from the generic Vector2D.
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
const direction noexcept
Definition scalarImpl.H:265
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
dictionary dict
const pointField & pts
The begin/end nodes for each segment, with divisions and expansion for each segment.
Definition PDRblock.H:299
gradingDescriptors grading() const
Return edge grading descriptors for the locations.
scalarList expansion_
The expansion ratio per segment.
Definition PDRblock.H:308
labelList divisions_
The number of division per segment.
Definition PDRblock.H:303
void prepend(const scalar p, label nDiv, scalar expRatio=1)
Add point/divisions/expand to front of list (push_front).
void append(const scalar p, label nDiv, scalar expRatio=1)
Add point/divisions/expand to end of list (push_back).
void resize(label len)
Resize lists.
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.