Go to the documentation of this file.
88 class gradingDescriptors;
126 inline bool valid()
const;
129 inline label
nCells()
const;
135 inline bool contains(
const scalar
p)
const;
138 inline const scalar&
min()
const;
141 inline const scalar&
max()
const;
144 inline scalar
centre()
const;
147 inline scalar
length()
const;
153 inline scalar
width(
const label i)
const;
157 inline scalar
C(
const label i)
const;
173 label
findIndex(
const scalar
p,
const scalar tol)
const;
178 inline const scalar&
clip(
const scalar& val)
const;
206 void append(
const scalar
p, label nDiv, scalar expRatio=1);
209 void prepend(
const scalar
p, label nDiv, scalar expRatio=1);
241 enum controlType : uint8_t
287 bool isSphere()
const;
290 bool onGround()
const;
294 bool onGround(
const bool on);
335 static bool checkMonotonic
349 const scalar scaleFactor = -1,
362 label addInternalFaces
364 faceList::iterator& faceIter,
365 labelList::iterator& ownIter,
366 labelList::iterator& neiIter
372 label addBoundaryFaces
375 faceList::iterator& faceIter,
376 labelList::iterator& ownIter
485 inline scalar
dx(
const label i)
const;
491 inline scalar
dy(
const label j)
const;
497 inline scalar
dz(
const label
k)
const;
503 inline vector span(
const label i,
const label j,
const label
k)
const;
509 inline point grid(
const label i,
const label j,
const label
k)
const;
515 inline point C(
const label i,
const label j,
const label
k)
const;
521 inline scalar
V(
const label i,
const label j,
const label
k)
const;
528 inline scalar
width(
const label i,
const label j,
const label
k)
const;
const static Enum< expansionType > expansionNames_
Named enumerations for the expansion type.
scalar width(const label i, const label j, const label k) const
Characteristic cell size at i,j,k position.
vectorField pointField
pointField is a vectorField.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
scalar dz(const label k) const
Cell size in z-direction at k position.
List< scalar > scalarList
A List of scalars.
autoPtr< polyMesh > innerMesh(const IOobject &io) const
Uniform expansion (ie, no expansion)
bool valid() const
The location list is valid if it contains 2 or more points.
A class for handling words, derived from Foam::string.
gradingDescriptors grading() const
Return edge grading descriptors for the locations.
scalar dy(const label j) const
Cell size in y-direction at j position.
Vector< gradingDescriptors > grading() const
Equivalent edge grading descriptors in (x,y,z) directions.
bool contains(const scalar p) const
True if the location is within the range.
void checkIndex(const label i) const
Check that element index is within valid range.
const labelVector & sizes() const
The (i,j,k) addressing dimensions.
scalarMinMax edgeLimits() const
Return min/max edge lengths.
const scalar & min() const
The first() value is considered the min value.
void prepend(const scalar p, label nDiv, scalar expRatio=1)
Add point/divisions/expand to front of list (push_front)
List of gradingDescriptor for the sections of a block with additional IO functionality.
void clear()
Reset to (0,0,0) sizing.
const boundBox & bounds() const
The mesh bounding box.
const scalar & clip(const scalar &val) const
void resize(label len)
Resize lists.
Vector< label > labelVector
Vector of labels.
void append(const scalar p, label nDiv, scalar expRatio=1)
Add point/divisions/expand to end of list (push_back)
void reset(const UList< scalar > &xgrid, const UList< scalar > &ygrid, const UList< scalar > &zgrid)
Reset grid locations and mesh i-j-k sizing.
label findCell(const scalar p) const
Find the cell index enclosing this location.
autoPtr< polyMesh > mesh(const IOobject &io) const
Create polyMesh for grid definition and patch information.
scalar dx(const label i) const
Cell size in x-direction at i position.
A simple i-j-k (row-major order) to linear addressing for a rectilinear mesh. Since the underlying me...
scalar length() const
The difference between min/max values, zero for an empty list.
const scalarMinMax & edgeLimits() const
The min/max edge length.
point C(const label i, const label j, const label k) const
Cell centre at i,j,k position.
scalar centre() const
Mid-point location, zero for an empty list.
MinMax< scalar > scalarMinMax
A scalar min/max range.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
scalar V(const label i, const label j, const label k) const
Cell volume at i,j,k position.
OBJstream os(runTime.globalPath()/outputName)
const Vector< location > & grid() const
The grid point locations in the i,j,k (x,y,z) directions.
label nPoints() const
The number of points in this direction.
Vector< scalar > vector
A scalar version of the templated Vector.
label nCells() const
Total number of cells in this direction.
A single block x-y-z rectilinear mesh addressable as i,j,k with simplified creation....
scalar width(const label i) const
Cell size at element position.
vector span(const label i, const label j, const label k) const
Cell dimensions at i,j,k position.
gradingDescriptors grading() const
Return edge grading descriptors for the locations.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
expansionType
The expansion type.
Vector2D< scalar > vector2D
A 2D vector of scalars obtained from the generic Vector2D.
Vector2D< label > labelVector2D
A 2D vector of labels obtained from the generic Vector2D.
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
label findIndex(const scalar p, const scalar tol) const
Find the grid index, within the given tolerance.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
label k
Boltzmann constant.
void writeDict(Ostream &os, const direction cmpt) const
Write as dictionary contents for specified vector direction.
A bounding box defined in terms of min/max extrema points.
void writeBlockMeshDict(const IOobject &io) const
Write an equivalent blockMeshDict.
scalarList expansion_
The expansion ratio per segment.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
vector point
Point is a vector.
label nCells() const
The number of cells in this direction.
Relative expansion ratio.
dictionary blockMeshDict() const
Content for an equivalent blockMeshDict.
bool read(const dictionary &dict)
Read dictionary.
const scalar & max() const
The last() value is considered the max value.
labelList divisions_
The number of division per segment.
PDRblock()
Default construct, zero-size, inverted bounds etc.
Grid locations in an axis direction.
dimensionedScalar pos(const dimensionedScalar &ds)
scalar C(const label i) const
Cell centre at element position.