Takes a set of patches to start MeshWave from. More...
Public Member Functions | |
patchDataWave (const polyMesh &mesh, const labelHashSet &patchIDs, const UPtrList< Field< Type > > &initialPatchValuePtrs, bool correctWalls=true) | |
Construct from mesh, information on patches to initialize and flag. More... | |
virtual | ~patchDataWave () |
Destructor. More... | |
virtual void | correct () |
Correct for mesh geom/topo changes. More... | |
const scalarField & | distance () const |
scalarField & | distance () |
Non const access so we can 'transfer' contents for efficiency. More... | |
const FieldField< Field, scalar > & | patchDistance () const |
FieldField< Field, scalar > & | patchDistance () |
const Field< Type > & | cellData () const |
Field< Type > & | cellData () |
const FieldField< Field, Type > & | patchData () const |
FieldField< Field, Type > & | patchData () |
label | nUnset () const |
![]() | |
ClassName ("cellDistFuncs") | |
cellDistFuncs (const polyMesh &mesh) | |
Construct from mesh. More... | |
const polyMesh & | mesh () const |
Access mesh. More... | |
labelHashSet | getPatchIDs (const wordReList &patchNames) const |
Return the set of patch IDs corresponding to the given names. More... | |
template<class Type > | |
labelHashSet | getPatchIDs () const |
Get patchIDs of/derived off certain type (e.g. 'processorPolyPatch') More... | |
scalar | smallestDist (const point &p, const polyPatch &patch, const label nWallFaces, const labelList &wallFaces, label &meshFaceI) const |
Calculate smallest true distance (and face index) More... | |
label | getPointNeighbours (const primitivePatch &, const label patchFaceI, labelList &) const |
Get faces sharing point with face on patch. More... | |
label | maxPatchSize (const labelHashSet &patchIDs) const |
Size of largest patch (out of supplied subset of patches) More... | |
label | sumPatchSize (const labelHashSet &patchIDs) const |
Sum of patch sizes (out of supplied subset of patches). More... | |
void | correctBoundaryFaceCells (const labelHashSet &patchIDs, scalarField &wallDistCorrected, Map< label > &nearestFace) const |
Correct all cells connected to boundary (via face). Sets values in. More... | |
void | correctBoundaryPointCells (const labelHashSet &patchIDs, scalarField &wallDistCorrected, Map< label > &nearestFace) const |
Correct all cells connected to wall (via point). Sets values in. More... | |
template<class Type > | |
Foam::labelHashSet | getPatchIDs () const |
Private Types | |
typedef TransferType::dataType | Type |
Private Member Functions | |
void | setChangedFaces (const labelHashSet &patchIDs, labelList &, List< TransferType > &) const |
Set initial set of changed faces. More... | |
label | getValues (const MeshWave< TransferType > &) |
Copy MeshWave values into *this. More... | |
Private Attributes | |
labelHashSet | patchIDs_ |
Current patch subset (stored as patchIDs) More... | |
const UPtrList< Field< Type > > & | initialPatchValuePtrs_ |
Reference to initial extra data at patch faces. More... | |
bool | correctWalls_ |
Do accurate distance calculation for near-wall cells. More... | |
label | nUnset_ |
Number of cells/faces unset after MeshWave has finished. More... | |
scalarField | distance_ |
Distance at cell centres. More... | |
FieldField< Field, scalar > | patchDistance_ |
Distance at patch faces. More... | |
Field< Type > | cellData_ |
Extra data at cell centres. More... | |
FieldField< Field, Type > | patchData_ |
Extra data at patch faces. More... | |
Takes a set of patches to start MeshWave from.
Holds after construction distance at cells and distance at patches (like patchWave), but also additional transported data. It is used, for example, in the y+ calculation.
Definition at line 63 of file patchDataWave.H.
|
private |
Definition at line 70 of file patchDataWave.H.
patchDataWave | ( | const polyMesh & | mesh, |
const labelHashSet & | patchIDs, | ||
const UPtrList< Field< Type > > & | initialPatchValuePtrs, | ||
bool | correctWalls = true |
||
) |
Construct from mesh, information on patches to initialize and flag.
whether or not to correct wall. Calculate for all cells. correctWalls : correct wall (face&point) cells for correct distance, searching neighbours.
Definition at line 172 of file patchDataWave.C.
References correct().
|
virtual |
Destructor.
Definition at line 196 of file patchDataWave.C.
|
private |
Set initial set of changed faces.
Definition at line 33 of file patchDataWave.C.
References polyPatch::faceCentres(), forAll, HashTable::found(), mesh, and polyPatch::start().
|
private |
Copy MeshWave values into *this.
Definition at line 75 of file patchDataWave.C.
References MeshWave< Type, TrackingData >::allCellInfo(), MeshWave< Type, TrackingData >::allFaceInfo(), MeshWave< Type, TrackingData >::data(), forAll, Foam::mag(), mesh, scalarField(), Foam::sqrt(), and polyPatch::start().
|
virtual |
Correct for mesh geom/topo changes.
Definition at line 204 of file patchDataWave.C.
References MeshWave< Type, TrackingData >::allFaceInfo(), forAll, and mesh.
|
inline |
Definition at line 145 of file patchDataWave.H.
References patchDataWave< TransferType >::distance_.
Referenced by wallDistData< TransferType >::correct(), and meshWave::correct().
|
inline |
Non const access so we can 'transfer' contents for efficiency.
Definition at line 151 of file patchDataWave.H.
References patchDataWave< TransferType >::distance_.
|
inline |
Definition at line 156 of file patchDataWave.H.
References patchDataWave< TransferType >::patchDistance_.
Referenced by wallDistData< TransferType >::correct(), and meshWave::correct().
|
inline |
Definition at line 161 of file patchDataWave.H.
References patchDataWave< TransferType >::patchDistance_.
Definition at line 166 of file patchDataWave.H.
References patchDataWave< TransferType >::cellData_.
Referenced by wallDistData< TransferType >::correct(), and meshWave::correct().
Definition at line 171 of file patchDataWave.H.
References patchDataWave< TransferType >::cellData_.
|
inline |
Definition at line 176 of file patchDataWave.H.
References patchDataWave< TransferType >::patchData_.
Referenced by wallDistData< TransferType >::correct(), and meshWave::correct().
|
inline |
Definition at line 181 of file patchDataWave.H.
References patchDataWave< TransferType >::patchData_.
|
inline |
Definition at line 186 of file patchDataWave.H.
References patchDataWave< TransferType >::nUnset_.
Referenced by wallDistData< TransferType >::correct(), and meshWave::correct().
|
private |
Current patch subset (stored as patchIDs)
Definition at line 76 of file patchDataWave.H.
Reference to initial extra data at patch faces.
Definition at line 79 of file patchDataWave.H.
|
private |
Do accurate distance calculation for near-wall cells.
Definition at line 82 of file patchDataWave.H.
|
private |
Number of cells/faces unset after MeshWave has finished.
Definition at line 89 of file patchDataWave.H.
Referenced by patchDataWave< TransferType >::nUnset().
|
private |
Distance at cell centres.
Definition at line 92 of file patchDataWave.H.
Referenced by patchDataWave< TransferType >::distance().
|
private |
Distance at patch faces.
Definition at line 95 of file patchDataWave.H.
Referenced by patchDataWave< TransferType >::patchDistance().
Extra data at cell centres.
Definition at line 98 of file patchDataWave.H.
Referenced by patchDataWave< TransferType >::cellData().
|
private |
Extra data at patch faces.
Definition at line 101 of file patchDataWave.H.
Referenced by patchDataWave< TransferType >::patchData().
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.