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 changedFacesInfo[changedFacei] =
100 data.boundaryField()[
patchi][patchFacei] ...,
101 mesh.boundaryMesh()[
patchi][patchFacei],
103 mesh.Cf().boundaryField()[
patchi][patchFacei],
109 List<FvWallInfoType> internalFaceInfo(mesh.nInternalFaces());
110 List<List<FvWallInfoType>> patchFaceInfo
112 FvFaceCellWave<FvWallInfoType, TrackingData>::template
113 sizesListList<List<List<FvWallInfoType>>>
115 FvFaceCellWave<FvWallInfoType, TrackingData>::template
116 listListSizes(mesh.boundary()),
120 List<FvWallInfoType> cellInfo(mesh.nCells());
127 FvFaceCellWave<FvWallInfoType, TrackingData>
wave
135 wave.setFaceInfo(changedPatchAndFaces, changedFacesInfo);
139 wave.iterate(mesh.globalData().nTotalCells() + 1);
146 wave.iterate(nCorrections - 1);
153 forAll(internalInfo, internali)
155 const bool valid = internalInfo[internali].valid(td);
161 distance.primitiveFieldRef()[internali] =
162 internalInfo[internali].dist(td);
164 (void)std::initializer_list<nil>
166 data.primitiveFieldRef()[internali] =
167 internalInfo[internali].data(td),
176 const bool valid = patchFaceInfo[
patchi][patchFacei].valid(td);
182 distance.boundaryFieldRef()[
patchi][patchFacei] =
183 patchFaceInfo[
patchi][patchFacei].dist(td) + small;
185 (void)std::initializer_list<nil>
187 data.boundaryFieldRef()[
patchi][patchFacei] =
188 patchFaceInfo[
patchi][patchFacei].data(td),
199 template<
template<
class>
class PatchField,
class GeoMesh>
204 const scalar minFaceFraction,
205 GeometricField<scalar, PatchField, GeoMesh>& distance
209 wave<FvWallInfo<wallPoint>>
215 FvFaceCellWave<FvWallInfo<wallPoint>>::defaultTrackingData_
220 template<
template<
class>
class PatchField,
class GeoMesh>
225 const scalar minFaceFraction,
226 const label nCorrections,
227 GeometricField<scalar, PatchField, GeoMesh>& distance
230 wave<FvWallInfo<wallFace>>
236 FvFaceCellWave<FvWallInfo<wallFace>>::defaultTrackingData_
241 template<
template<
class>
class PatchField,
class GeoMesh>
246 const scalar minFaceFraction,
247 const label nCorrections,
248 GeometricField<scalar, PatchField, GeoMesh>& distance
251 const List<labelPair> changedPatchAndFaces =
255 wave<FvWallInfo<wallPoint>>
258 changedPatchAndFaces,
261 FvFaceCellWave<FvWallInfo<wallPoint>>::defaultTrackingData_
264 wave<FvWallInfo<wallFace>>
267 changedPatchAndFaces,
270 FvFaceCellWave<FvWallInfo<wallFace>>::defaultTrackingData_
279 template<
class>
class WallLocation,
281 template<
class>
class PatchField,
289 const scalar minFaceFraction,
290 GeometricField<scalar, PatchField, GeoMesh>& distance,
291 GeometricField<DataType, PatchField, GeoMesh>& data,
296 wave<FvWallInfo<WallLocation<wallPoint>>, TrackingData>
310 template<
class>
class WallLocation,
312 template<
class>
class PatchField,
320 const scalar minFaceFraction,
321 const label nCorrections,
322 GeometricField<scalar, PatchField, GeoMesh>& distance,
323 GeometricField<DataType, PatchField, GeoMesh>& data,
327 wave<FvWallInfo<WallLocation<wallFace>>, TrackingData>
341 template<
class>
class WallLocation,
343 template<
class>
class PatchField,
351 const scalar minFaceFraction,
352 const label nCorrections,
353 GeometricField<scalar, PatchField, GeoMesh>& distance,
354 GeometricField<DataType, PatchField, GeoMesh>& data,
358 const List<labelPair> changedPatchAndFaces =
362 wave<FvWallInfo<WallLocation<wallPoint>>, TrackingData>
365 changedPatchAndFaces,
372 wave<FvWallInfo<WallLocation<wallFace>>, TrackingData>
375 changedPatchAndFaces,
393 template<
class WallLocation>
403 template<
class>
class PatchField,
411 const scalar minFaceFraction,
412 GeometricField<scalar, PatchField, GeoMesh>& distance,
413 GeometricField<DataType, PatchField, GeoMesh>& data,
418 calculate<WallLocationDataType<DataType>::template
type>
433 template<
class>
class PatchField,
441 const scalar minFaceFraction,
442 const label nCorrections,
443 GeometricField<scalar, PatchField, GeoMesh>& distance,
444 GeometricField<DataType, PatchField, GeoMesh>& data,
448 correct<WallLocationDataType<DataType>::template
type>
463 template<
class>
class PatchField,
471 const scalar minFaceFraction,
472 const label nCorrections,
473 GeometricField<scalar, PatchField, GeoMesh>& distance,
474 GeometricField<DataType, PatchField, GeoMesh>& data,
479 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.