36 #ifndef fvPatchDistWave_H 37 #define fvPatchDistWave_H 58 const scalar minFaceFraction
71 template<
class>
class PatchField,
79 const label nCorrections,
86 template<
template<
class>
class PatchField,
class GeoMesh>
91 const scalar minFaceFraction,
96 template<
template<
class>
class PatchField,
class GeoMesh>
101 const scalar minFaceFraction,
102 const label nCorrections,
107 template<
template<
class>
class PatchField,
class GeoMesh>
112 const scalar minFaceFraction,
113 const label nCorrections,
120 template<
class>
class WallInfoData,
121 template<
class>
class PatchField,
123 class TrackingData =
int 129 const scalar minFaceFraction,
132 <
typename WallInfoData<wallPoint>::dataType, PatchField, GeoMesh>&
141 template<
class>
class WallInfoData,
142 template<
class>
class PatchField,
144 class TrackingData =
int 150 const scalar minFaceFraction,
151 const label nCorrections,
154 <
typename WallInfoData<wallPoint>::dataType, PatchField, GeoMesh>&
163 template<
class>
class WallInfoData,
164 template<
class>
class PatchField,
166 class TrackingData =
int 172 const scalar minFaceFraction,
173 const label nCorrections,
176 <
typename WallInfoData<wallPoint>::dataType, PatchField, GeoMesh>&
Wave propagation of information through grid. Every iteration information goes through one layer of c...
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void correct(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, const label nCorrections, GeometricField< scalar, PatchField, GeoMesh > &distance)
Correct distance data from patches.
Holds information (coordinate and normal) regarding nearest wall point.
scalar distance(const vector &p1, const vector &p2)
Generic GeometricField class.
label calculate(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, GeometricField< scalar, PatchField, GeoMesh > &distance)
Calculate distance data from patches.
label calculateAndCorrect(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, const label nCorrections, GeometricField< scalar, PatchField, GeoMesh > &distance)
Calculate and correct distance data from patches.
Database for solution and other reduced data.
Takes a set of patches to start FvFaceCellWave from and computed the distance at patches and possibly...
Mesh data needed to do the Finite Volume discretisation.
List< labelPair > getChangedPatchAndFaces(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction)
Get initial set of changed faces.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
label wave(const fvMesh &mesh, const List< labelPair > &changedPatchAndFaces, const label nCorrections, GeometricField< scalar, PatchField, GeoMesh > &distance, TrackingData &td, GeometricField< DataType, PatchField, GeoMesh > &... data)
Wave distance (and maybe additional) data from faces. If nCorrections is.