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. More...
Classes | |
struct | gridControl |
class | location |
Public Types | |
enum | expansionType : uint8_t { EXPAND_UNIFORM, EXPAND_RATIO, EXPAND_RELATIVE } |
Public Member Functions | |
PDRblock () | |
PDRblock (const UList< scalar > &xgrid, const UList< scalar > &ygrid, const UList< scalar > &zgrid) | |
PDRblock (const dictionary &dict, bool verboseOutput=false) | |
bool | read (const dictionary &dict) |
void | reset (const UList< scalar > &xgrid, const UList< scalar > &ygrid, const UList< scalar > &zgrid) |
const Vector< location > & | grid () const |
Vector< gradingDescriptors > | grading () const |
gradingDescriptors | grading (const direction cmpt) const |
const boundBox & | bounds () const |
const scalarMinMax & | edgeLimits () const |
scalar | dx (const label i) const |
scalar | dx (const labelVector &ijk) const |
scalar | dy (const label j) const |
scalar | dy (const labelVector &ijk) const |
scalar | dz (const label k) const |
scalar | dz (const labelVector &ijk) const |
vector | span (const label i, const label j, const label k) const |
vector | span (const labelVector &ijk) const |
point | grid (const label i, const label j, const label k) const |
point | grid (const labelVector &ijk) const |
point | C (const label i, const label j, const label k) const |
point | C (const labelVector &ijk) const |
scalar | V (const label i, const label j, const label k) const |
scalar | V (const labelVector &ijk) const |
scalar | width (const label i, const label j, const label k) const |
scalar | width (const labelVector &ijk) const |
labelVector | findCell (const point &pt) const |
labelVector | gridIndex (const point &pt, const scalar relTol=0.01) const |
Ostream & | blockMeshDict (Ostream &os, const bool withHeader=false) const |
dictionary | blockMeshDict () const |
void | writeBlockMeshDict (const IOobject &io) const |
autoPtr< polyMesh > | mesh (const IOobject &io) const |
autoPtr< polyMesh > | innerMesh (const IOobject &io) const |
![]() | |
ijkMesh () | |
ijkMesh (const labelVector &ijk) | |
ijkMesh (const label nx, const label ny, const label nz) | |
label | nPoints () const |
label | nCells () const |
label | nFaces () const |
label | nInternalFaces () const |
label | nBoundaryFaces () const |
label | nBoundaryFaces (const direction shapeFacei) const |
label | cellLabel (const label i, const label j, const label k) const |
label | cellLabel (const labelVector &ijk) const |
label | pointLabel (const label i, const label j, const label k) const |
label | pointLabel (const labelVector &ijk) const |
hexCell | vertLabels (const label i, const label j, const label k) const |
hexCell | vertLabels (const labelVector &ijk) const |
![]() | |
ijkAddressing () | |
ijkAddressing (const labelVector &ijk) | |
ijkAddressing (const label ni, const label nj, const label nk) | |
bool | empty () const |
const labelVector & | sizes () const |
labelVector & | sizes () |
label | size () const |
const label & | size (const vector::components cmpt) const |
void | clear () |
void | reset (const label ni, const label nj, const label nk) |
void | reset (const labelVector &newSizes) |
label | index (const label i, const label j, const label k) const |
label | index (const labelVector &ijk) const |
labelVector | index (const label idx) const |
void | checkIndex (const label i, const label j, const label k, const bool allowExtra=false) const |
void | checkIndex (const labelVector &ijk, const bool allowExtra=false) const |
void | checkSizes () const |
void | checkSizes (const labelVector &other) const |
void | checkSizes (const label nTotal) const |
Static Public Member Functions | |
static const PDRblock & | null () |
Static Public Attributes | |
const static Enum< expansionType > | expansionNames_ |
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.
Dictionary controls
Property | Description | Required | Default |
---|---|---|---|
x | X-direction grid specification | yes | |
y | Y-direction grid specification | yes | |
z | Z-direction grid specification | yes | |
scale | Point scaling | no | 1.0 |
expansion | Type of expansion (ratio/relative) | no | ratio |
boundary | Boundary patches | yes | |
defaultPatch | Default patch specification | no |
Grid coordinate controls
Property | Description | Required | Default |
---|---|---|---|
points | Locations defining the mesh segment | yes | |
nCells | Divisions per mesh segment | yes | |
ratios | Expansion values per segment | no |
A negative expansion value is trapped and treated as its reciprocal. by default, the expansion is as per blockMesh and represents the ratio of end-size / start-size for the section. Alternatively, the relative size can be given.
Definition at line 149 of file PDRblock.H.
enum expansionType : uint8_t |
Enumerator | |
---|---|
EXPAND_UNIFORM | Uniform expansion (ie, no expansion) |
EXPAND_RATIO | End/start ratio. |
EXPAND_RELATIVE | Relative expansion ratio. |
Definition at line 158 of file PDRblock.H.
PDRblock | ( | ) |
Definition at line 508 of file PDRblock.C.
PDRblock | ( | const UList< scalar > & | xgrid, |
const UList< scalar > & | ygrid, | ||
const UList< scalar > & | zgrid | ||
) |
Definition at line 515 of file PDRblock.C.
References Foam::name(), and reset().
|
explicit |
Definition at line 541 of file PDRblock.C.
References dict, and PDRblock::read().
|
static |
Definition at line 106 of file PDRblock.C.
bool read | ( | const dictionary & | dict | ) |
Definition at line 561 of file PDRblock.C.
References dict, dictionary::findDict(), dictionary::getOrDefault(), Foam::Info, dictionary::read(), and dictionary::subDict().
Referenced by PDRblock::PDRblock().
void reset | ( | const UList< scalar > & | xgrid, |
const UList< scalar > & | ygrid, | ||
const UList< scalar > & | zgrid | ||
) |
Definition at line 598 of file PDRblock.C.
References VectorSpace< Vector< scalar >, scalar, 3 >::nComponents.
|
inline |
Definition at line 141 of file PDRblockI.H.
Foam::Vector< Foam::gradingDescriptors > grading | ( | ) | const |
Definition at line 706 of file PDRblock.C.
Foam::gradingDescriptors grading | ( | const direction | cmpt | ) | const |
Definition at line 712 of file PDRblock.C.
References Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, Vector< scalar >::X, Vector< scalar >::Y, and Vector< scalar >::Z.
|
inline |
Definition at line 153 of file PDRblockI.H.
|
inline |
Definition at line 147 of file PDRblockI.H.
|
inline |
Definition at line 159 of file PDRblockI.H.
Referenced by PDRblock::span(), and PDRblock::V().
|
inline |
Definition at line 165 of file PDRblockI.H.
References Vector< Cmpt >::x().
|
inline |
Definition at line 171 of file PDRblockI.H.
Referenced by PDRblock::span(), and PDRblock::V().
|
inline |
Definition at line 177 of file PDRblockI.H.
References Vector< Cmpt >::y().
|
inline |
Definition at line 183 of file PDRblockI.H.
References k.
Referenced by PDRblock::span(), and PDRblock::V().
|
inline |
Definition at line 189 of file PDRblockI.H.
References Vector< Cmpt >::z().
|
inline |
Definition at line 196 of file PDRblockI.H.
References PDRblock::dx(), PDRblock::dy(), PDRblock::dz(), and k.
|
inline |
Definition at line 206 of file PDRblockI.H.
References PDRblock::dx(), PDRblock::dy(), and PDRblock::dz().
|
inline |
Definition at line 213 of file PDRblockI.H.
References k.
|
inline |
Definition at line 223 of file PDRblockI.H.
References Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().
|
inline |
Definition at line 236 of file PDRblockI.H.
References k.
|
inline |
Definition at line 246 of file PDRblockI.H.
References Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().
|
inline |
Definition at line 259 of file PDRblockI.H.
References PDRblock::dx(), PDRblock::dy(), PDRblock::dz(), and k.
Referenced by PDRblock::width().
|
inline |
Definition at line 269 of file PDRblockI.H.
References PDRblock::dx(), PDRblock::dy(), PDRblock::dz(), Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().
|
inline |
Definition at line 276 of file PDRblockI.H.
References Foam::cbrt(), k, and PDRblock::V().
|
inline |
Definition at line 286 of file PDRblockI.H.
References Foam::cbrt(), and PDRblock::V().
Foam::labelVector findCell | ( | const point & | pt | ) | const |
Definition at line 676 of file PDRblock.C.
References Foam::pos().
Foam::labelVector gridIndex | ( | const point & | pt, |
const scalar | relTol = 0.01 |
||
) | const |
Definition at line 690 of file PDRblock.C.
References Foam::pos().
Foam::Ostream & blockMeshDict | ( | Ostream & | os, |
const bool | withHeader = false |
||
) | const |
Definition at line 231 of file PDRblockBlockMesh.C.
References DynamicList::append(), Foam::begIndentList(), boundBox::centre(), Foam::constant::electromagnetic::e, Foam::endIndentList(), hex, Foam::indent(), HashSet::insert(), boundBox::max(), Foam::max(), boundBox::min(), Foam::min(), Foam::New(), Foam::nl, os(), Foam::outputIndent(), boundBox::points(), Foam::projGeomName, ref(), Foam::relativeToGeometricRatio(), Foam::serializeFace(), Foam::serializeHex(), Foam::serializeProjectEdge(), Foam::serializeProjectFace(), Foam::serializeProjectPoints(), boundBox::span(), IOobject::writeHeader(), Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().
Foam::dictionary blockMeshDict | ( | ) | const |
Definition at line 770 of file PDRblockBlockMesh.C.
References os(), and List::transfer().
void writeBlockMeshDict | ( | const IOobject & | io | ) | const |
Definition at line 782 of file PDRblockBlockMesh.C.
References IOobject::db(), IOstream::defaultPrecision(), Foam::endl(), Foam::Info, IOobject::local(), Foam::max(), OFstream::name(), IOobject::name(), Foam::nl, IOobject::NO_READ, IOobject::NO_WRITE, IOobject::objectPath(), os(), OSstream::precision(), TimePaths::relativePath(), TimePaths::system(), objectRegistry::time(), IOobject::writeEndDivider(), and IOobject::writeHeader().
Foam::autoPtr< Foam::polyMesh > mesh | ( | const IOobject & | io | ) | const |
Definition at line 374 of file PDRblockCreate.C.
References Foam::Info, and Foam::nl.
Foam::autoPtr< Foam::polyMesh > innerMesh | ( | const IOobject & | io | ) | const |
Definition at line 297 of file PDRblockCreate.C.
References polyMesh::addPatches(), IOobject::AUTO_WRITE, polyMesh::boundaryMesh(), meshPtr, autoPtr::New(), polyPatch::New(), nPoints, patches, PtrList::set(), and IOobject::writeOpt().
|
static |
Definition at line 166 of file PDRblock.H.
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.