Public Member Functions | Private Types | Private Member Functions | Private Attributes
patchDataWave< TransferType > Class Template Reference

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

Inheritance diagram for patchDataWave< TransferType >:
Inheritance graph
[legend]
Collaboration diagram for patchDataWave< TransferType >:
Collaboration graph
[legend]

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...
 
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...
 
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< TypecellData_
 Extra data at cell centres. More...
 
FieldField< Field, TypepatchData_
 Extra data at patch faces. More...
 

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.

Member Typedef Documentation

◆ Type

typedef TransferType::dataType Type
private

Definition at line 70 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().

Here is the call graph for this function:

◆ ~patchDataWave()

~patchDataWave
virtual

Destructor.

Definition at line 196 of file patchDataWave.C.

Member Function Documentation

◆ setChangedFaces()

void setChangedFaces ( const labelHashSet patchIDs,
labelList changedFaces,
List< TransferType > &  faceDist 
) const
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().

Here is the call graph for this function:

◆ getValues()

Foam::label getValues ( const MeshWave< TransferType > &  waveInfo)
private

◆ correct()

void correct
virtual

Correct for mesh geom/topo changes.

Definition at line 204 of file patchDataWave.C.

References MeshWave< Type, TrackingData >::allFaceInfo(), forAll, and mesh.

Here is the call graph for this function:

◆ distance() [1/2]

const scalarField& distance ( ) const
inline

Definition at line 145 of file patchDataWave.H.

References patchDataWave< TransferType >::distance_.

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 151 of file patchDataWave.H.

References patchDataWave< TransferType >::distance_.

◆ patchDistance() [1/2]

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

Definition at line 156 of file patchDataWave.H.

References patchDataWave< TransferType >::patchDistance_.

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 161 of file patchDataWave.H.

References patchDataWave< TransferType >::patchDistance_.

◆ cellData() [1/2]

const Field<Type>& cellData ( ) const
inline

Definition at line 166 of file patchDataWave.H.

References patchDataWave< TransferType >::cellData_.

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 171 of file patchDataWave.H.

References patchDataWave< TransferType >::cellData_.

◆ patchData() [1/2]

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

Definition at line 176 of file patchDataWave.H.

References patchDataWave< TransferType >::patchData_.

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 181 of file patchDataWave.H.

References patchDataWave< TransferType >::patchData_.

◆ nUnset()

label nUnset ( ) const
inline

Definition at line 186 of file patchDataWave.H.

References patchDataWave< TransferType >::nUnset_.

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

Here is the caller graph for this function:

Field Documentation

◆ patchIDs_

labelHashSet patchIDs_
private

Current patch subset (stored as patchIDs)

Definition at line 76 of file patchDataWave.H.

◆ initialPatchValuePtrs_

const UPtrList<Field<Type> >& initialPatchValuePtrs_
private

Reference to initial extra data at patch faces.

Definition at line 79 of file patchDataWave.H.

◆ correctWalls_

bool correctWalls_
private

Do accurate distance calculation for near-wall cells.

Definition at line 82 of file patchDataWave.H.

◆ nUnset_

label nUnset_
private

Number of cells/faces unset after MeshWave has finished.

Definition at line 89 of file patchDataWave.H.

Referenced by patchDataWave< TransferType >::nUnset().

◆ distance_

scalarField distance_
private

Distance at cell centres.

Definition at line 92 of file patchDataWave.H.

Referenced by patchDataWave< TransferType >::distance().

◆ patchDistance_

FieldField<Field, scalar> patchDistance_
private

Distance at patch faces.

Definition at line 95 of file patchDataWave.H.

Referenced by patchDataWave< TransferType >::patchDistance().

◆ cellData_

Field<Type> cellData_
private

Extra data at cell centres.

Definition at line 98 of file patchDataWave.H.

Referenced by patchDataWave< TransferType >::cellData().

◆ patchData_

FieldField<Field, Type> patchData_
private

Extra data at patch faces.

Definition at line 101 of file patchDataWave.H.

Referenced by patchDataWave< TransferType >::patchData().


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