patchDataWave< TransferType > Class Template Reference

Takes a set of patches to start MeshWave from. More...

Inheritance diagram for patchDataWave< TransferType >:
Collaboration diagram for patchDataWave< TransferType >:

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 scalarFielddistance () const
 
scalarFielddistance ()
 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
 
- Public Member Functions inherited from cellDistFuncs
 ClassName ("cellDistFuncs")
 
 cellDistFuncs (const polyMesh &mesh)
 Construct from mesh. More...
 
 cellDistFuncs (const cellDistFuncs &)=delete
 Disallow default bitwise copy construction. More...
 
const polyMeshmesh () 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...
 
void operator= (const cellDistFuncs &)=delete
 Disallow default bitwise assignment. More...
 
template<class Type >
Foam::labelHashSet getPatchIDs () const
 

Detailed Description

template<class TransferType>
class Foam::patchDataWave< TransferType >

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.

See also
The patchWave class.
Source files

Definition at line 63 of file patchDataWave.H.

Constructor & Destructor Documentation

◆ patchDataWave()

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.

◆ ~patchDataWave()

~patchDataWave ( )
virtual

Destructor.

Definition at line 196 of file patchDataWave.C.

Member Function Documentation

◆ correct()

void correct ( )
virtual

Correct for mesh geom/topo changes.

Definition at line 204 of file patchDataWave.C.

References MeshWave< Type, TrackingData >::allFaceInfo(), UList< T >::data(), forAll, mesh, and HashTable< T, Key, Hash >::toc().

Here is the call graph for this function:

◆ distance() [1/2]

const scalarField& distance ( ) const
inline

Definition at line 142 of file patchDataWave.H.

Referenced by wallDistData< TransferType >::correct(), and meshWave::correct().

Here is the caller graph for this function:

◆ distance() [2/2]

scalarField& distance ( )
inline

Non const access so we can 'transfer' contents for efficiency.

Definition at line 148 of file patchDataWave.H.

◆ patchDistance() [1/2]

const FieldField<Field, scalar>& patchDistance ( ) const
inline

Definition at line 153 of file patchDataWave.H.

Referenced by wallDistData< TransferType >::correct(), and meshWave::correct().

Here is the caller graph for this function:

◆ patchDistance() [2/2]

FieldField<Field, scalar>& patchDistance ( )
inline

Definition at line 158 of file patchDataWave.H.

◆ cellData() [1/2]

const Field<Type>& cellData ( ) const
inline

Definition at line 163 of file patchDataWave.H.

Referenced by wallDistData< TransferType >::correct(), and meshWave::correct().

Here is the caller graph for this function:

◆ cellData() [2/2]

Field<Type>& cellData ( )
inline

Definition at line 168 of file patchDataWave.H.

◆ patchData() [1/2]

const FieldField<Field, Type>& patchData ( ) const
inline

Definition at line 173 of file patchDataWave.H.

Referenced by wallDistData< TransferType >::correct(), and meshWave::correct().

Here is the caller graph for this function:

◆ patchData() [2/2]

FieldField<Field, Type>& patchData ( )
inline

Definition at line 178 of file patchDataWave.H.

◆ nUnset()

label nUnset ( ) const
inline

Definition at line 183 of file patchDataWave.H.

Referenced by wallDistData< TransferType >::correct(), and meshWave::correct().

Here is the caller graph for this function:

The documentation for this class was generated from the following files: