36 template<
class FvWallInfoType,
class TrackingData>
43 return wave.cellInfo();
46 template<
class FvWallInfoType,
class TrackingData>
53 return wave.internalFaceInfo();
66 template<
class>
class PatchField,
73 const List<labelPair>& changedPatchAndFaces,
74 const label nCorrections,
75 GeometricField<scalar, PatchField, GeoMesh>& distance,
77 GeometricField<DataType, PatchField, GeoMesh>& ... data
86 if (!
calculate && nCorrections == 0)
return 0;
89 List<FvWallInfoType> changedFacesInfo(changedPatchAndFaces.size());
90 forAll(changedPatchAndFaces, changedFacei)
93 changedPatchAndFaces[changedFacei].first();
94 const label patchFacei =
95 changedPatchAndFaces[changedFacei].second();
97 const label polyFacei =
98 mesh.polyFacesBf()[
patchi][patchFacei];
100 changedFacesInfo[changedFacei] =
103 data.boundaryField()[
patchi][patchFacei] ...,
104 mesh.faces()[polyFacei],
106 mesh.Cf().boundaryField()[
patchi][patchFacei],
112 List<FvWallInfoType> internalFaceInfo(mesh.nInternalFaces());
113 List<List<FvWallInfoType>> patchFaceInfo
115 FvFaceCellWave<FvWallInfoType, TrackingData>::template
116 sizesListList<List<List<FvWallInfoType>>>
118 FvFaceCellWave<FvWallInfoType, TrackingData>::template
119 listListSizes(mesh.boundary()),
123 List<FvWallInfoType> cellInfo(mesh.nCells());
130 FvFaceCellWave<FvWallInfoType, TrackingData>
wave
138 wave.setFaceInfo(changedPatchAndFaces, changedFacesInfo);
142 wave.iterate(mesh.globalData().nTotalCells() + 1);
149 wave.iterate(nCorrections - 1);
156 forAll(internalInfo, internali)
158 const bool valid = internalInfo[internali].valid(td);
164 distance.primitiveFieldRef()[internali] =
165 internalInfo[internali].dist(td);
167 (void)std::initializer_list<nil>
169 data.primitiveFieldRef()[internali] =
170 internalInfo[internali].data(td),
179 const bool valid = patchFaceInfo[
patchi][patchFacei].valid(td);
185 distance.boundaryFieldRef()[
patchi][patchFacei] =
186 patchFaceInfo[
patchi][patchFacei].dist(td) + small;
188 (void)std::initializer_list<nil>
190 data.boundaryFieldRef()[
patchi][patchFacei] =
191 patchFaceInfo[
patchi][patchFacei].data(td),
202 template<
template<
class>
class PatchField,
class GeoMesh>
207 const scalar minFaceFraction,
208 GeometricField<scalar, PatchField, GeoMesh>& distance
212 wave<FvWallInfo<wallPoint>>
218 FvFaceCellWave<FvWallInfo<wallPoint>>::defaultTrackingData_
223 template<
template<
class>
class PatchField,
class GeoMesh>
228 const scalar minFaceFraction,
229 const label nCorrections,
230 GeometricField<scalar, PatchField, GeoMesh>& distance
233 wave<FvWallInfo<wallFace>>
239 FvFaceCellWave<FvWallInfo<wallFace>>::defaultTrackingData_
244 template<
template<
class>
class PatchField,
class GeoMesh>
249 const scalar minFaceFraction,
250 const label nCorrections,
251 GeometricField<scalar, PatchField, GeoMesh>& distance
254 const List<labelPair> changedPatchAndFaces =
258 wave<FvWallInfo<wallPoint>>
261 changedPatchAndFaces,
264 FvFaceCellWave<FvWallInfo<wallPoint>>::defaultTrackingData_
267 wave<FvWallInfo<wallFace>>
270 changedPatchAndFaces,
273 FvFaceCellWave<FvWallInfo<wallFace>>::defaultTrackingData_
282 template<
class>
class WallLocation,
284 template<
class>
class PatchField,
292 const scalar minFaceFraction,
293 GeometricField<scalar, PatchField, GeoMesh>& distance,
294 GeometricField<DataType, PatchField, GeoMesh>& data,
299 wave<FvWallInfo<WallLocation<wallPoint>>, TrackingData>
313 template<
class>
class WallLocation,
315 template<
class>
class PatchField,
323 const scalar minFaceFraction,
324 const label nCorrections,
325 GeometricField<scalar, PatchField, GeoMesh>& distance,
326 GeometricField<DataType, PatchField, GeoMesh>& data,
330 wave<FvWallInfo<WallLocation<wallFace>>, TrackingData>
344 template<
class>
class WallLocation,
346 template<
class>
class PatchField,
354 const scalar minFaceFraction,
355 const label nCorrections,
356 GeometricField<scalar, PatchField, GeoMesh>& distance,
357 GeometricField<DataType, PatchField, GeoMesh>& data,
361 const List<labelPair> changedPatchAndFaces =
365 wave<FvWallInfo<WallLocation<wallPoint>>, TrackingData>
368 changedPatchAndFaces,
375 wave<FvWallInfo<WallLocation<wallFace>>, TrackingData>
378 changedPatchAndFaces,
396 template<
class WallLocation>
406 template<
class>
class PatchField,
414 const scalar minFaceFraction,
415 GeometricField<scalar, PatchField, GeoMesh>& distance,
416 GeometricField<DataType, PatchField, GeoMesh>& data,
421 calculate<WallLocationDataType<DataType>::template
type>
436 template<
class>
class PatchField,
444 const scalar minFaceFraction,
445 const label nCorrections,
446 GeometricField<scalar, PatchField, GeoMesh>& distance,
447 GeometricField<DataType, PatchField, GeoMesh>& data,
451 correct<WallLocationDataType<DataType>::template
type>
466 template<
class>
class PatchField,
474 const scalar minFaceFraction,
475 const label nCorrections,
476 GeometricField<scalar, PatchField, GeoMesh>& distance,
477 GeometricField<DataType, PatchField, GeoMesh>& data,
482 calculateAndCorrect<WallLocationDataType<DataType>::template
type>
#define forAll(list, i)
Loop across all elements in list.
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Generic GeometricField class.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Holds information regarding nearest wall point. Used in wall distance calculation.
Takes a set of patches to start FvFaceCellWave from and computed the distance at patches and possibly...
bool valid(const PtrList< ModelType > &l)
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.
const List< FvWallInfoType > & getInternalInfo(const volScalarField &distance, FvFaceCellWave< FvWallInfoType, TrackingData > &wave)
void correct(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, const label nCorrections, GeometricField< scalar, PatchField, GeoMesh > &distance)
Correct distance data from patches.
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.
List< labelPair > getChangedPatchAndFaces(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction)
Get initial set of changed faces.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.