68 #ifndef Foam_PDRblock_H 69 #define Foam_PDRblock_H 88 class gradingDescriptors;
126 void reset(
const scalar low,
const scalar upp,
const label
nCells);
138 inline
bool contains(const scalar
p) const;
141 inline const scalar&
min() const;
144 inline const scalar&
max() const;
147 inline scalar
centre() const;
150 inline scalar
length() const;
156 inline scalar
width(const label i) const;
160 inline scalar
C(const label i) const;
176 label
findIndex(const scalar
p, const scalar tol) const;
181 inline const scalar&
clip(const scalar& val) const;
209 void append(
const scalar
p, label nDiv, scalar expRatio=1);
212 void prepend(
const scalar
p, label nDiv, scalar expRatio=1);
244 enum controlType : uint8_t
253 const static Enum<controlType> controlNames_;
290 bool isSphere()
const;
293 bool onGround()
const;
297 bool onGround(
const bool on);
338 static bool checkMonotonic
345 void addDefaultPatches();
355 const scalar scaleFactor = -1,
368 label addInternalFaces
378 label addBoundaryFaces
497 inline scalar
dx(
const label i)
const;
503 inline scalar
dy(
const label j)
const;
509 inline scalar
dz(
const label
k)
const;
515 inline vector span(
const label i,
const label j,
const label
k)
const;
521 inline point grid(
const label i,
const label j,
const label
k)
const;
527 inline point C(
const label i,
const label j,
const label
k)
const;
533 inline scalar
V(
const label i,
const label j,
const label
k)
const;
540 inline scalar
width(
const label i,
const label j,
const label
k)
const;
564 Ostream&
blockMeshDict(Ostream&
os,
const bool withHeader=
false)
const;
573 autoPtr<polyMesh>
mesh(
const IOobject&
io)
const;
577 autoPtr<polyMesh>
innerMesh(
const IOobject&
io)
const;
List< scalar > scalarList
List of scalar.
scalar centre() const
Mid-point location, zero for an empty list.
const scalar & max() const
The last() value is considered the max value.
Graphite solid properties.
MinMax< scalar > scalarMinMax
A scalar min/max range.
const scalar & clip(const scalar &val) const
If out of range, return the respective min/max limits, otherwise return the value itself...
label findIndex(const scalar p, const scalar tol) const
Find the grid index, within the given tolerance.
A face is a list of labels corresponding to mesh vertices.
void resize(const label len)
Adjust allocated size of list.
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.
void append(const T &val)
Append an element at the end of the list.
Uniform expansion (ie, no expansion)
void checkIndex(const label i) const
Check that element index is within valid range.
const scalar & min() const
The first() value is considered the min value.
Vector2D< scalar > vector2D
A 2D vector of scalars obtained from the generic Vector2D.
vector span(const label i, const label j, const label k) const
Cell dimensions at i,j,k position.
label * iterator
Random access iterator for traversing a UList.
A bounding box defined in terms of min/max extrema points.
label k
Boltzmann constant.
static const List< T > & null()
Return a null List.
label nCells() const noexcept
The number of cells in this direction.
expansionType
The expansion type.
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.
dimensionedScalar pos(const dimensionedScalar &ds)
Vector2D< label > labelVector2D
A 2D vector of labels obtained from the generic Vector2D.
PDRblock()
Default construct, zero-size, inverted bounds etc.
bool contains(const scalar p) const
True if the location is within the range.
bool good() const noexcept
The location list is valid if it contains 2 or more points.
void clear()
Clear the list, i.e. set size to zero.
The begin/end nodes for each segment, with divisions and expansion for each segment.
scalar dz(const label k) const
Cell size in z-direction at k position.
A class for handling words, derived from Foam::string.
static const Enum< expansionType > expansionNames_
Named enumerations for the expansion type.
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.
autoPtr< polyMesh > innerMesh(const IOobject &io) const
Create polyMesh for inner-mesh only, ignore any outer block definitions.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
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.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
scalar length() const
The difference between min/max values, zero for an empty list.
OBJstream os(runTime.globalPath()/outputName)
scalar dy(const label j) const
Cell size in y-direction at j position.
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.
vector point
Point is a vector.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Vector< label > labelVector
Vector of labels.
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.
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...
gradingDescriptors grading() const
Return edge grading descriptors for the locations.
const boundBox & bounds() const noexcept
The mesh bounding box.
scalar dx(const label i) const
Cell size in x-direction at i position.
void reset(const scalar low, const scalar upp, const label nCells)
Reset with min/max range and number of divisions.
label nPoints() const noexcept
The number of points in this direction.
const Vector< location > & grid() const noexcept
The grid point locations in the i,j,k (x,y,z) directions.
Relative expansion ratio.