36 template<
class FvWallInfoType,
class TrackingData>
43 return wave.cellInfo();
46 template<
class FvWallInfoType,
class TrackingData>
53 return wave.internalFaceInfo();
72 const List<labelPair>& changedPatchAndFaces,
73 const label nCorrections,
74 GeometricField<scalar, GeoMesh>& distance,
76 GeometricField<DataType, GeoMesh>& ... data
85 if (!
calculate && nCorrections == 0)
return 0;
88 List<FvWallInfoType> changedFacesInfo(changedPatchAndFaces.size());
89 forAll(changedPatchAndFaces, changedFacei)
92 changedPatchAndFaces[changedFacei].first();
93 const label patchFacei =
94 changedPatchAndFaces[changedFacei].second();
96 const label polyFacei =
99 changedFacesInfo[changedFacei] =
102 data.boundaryField()[
patchi][patchFacei] ...,
112 List<List<FvWallInfoType>> patchFaceInfo
114 FvFaceCellWave<FvWallInfoType, TrackingData>::template
115 sizesListList<List<List<FvWallInfoType>>>
117 FvFaceCellWave<FvWallInfoType, TrackingData>::template
129 FvFaceCellWave<FvWallInfoType, TrackingData>
wave
137 wave.setFaceInfo(changedPatchAndFaces, changedFacesInfo);
148 wave.iterate(nCorrections - 1);
155 forAll(internalInfo, internali)
157 const bool valid = internalInfo[internali].valid(td);
163 distance.primitiveFieldRef()[internali] =
164 internalInfo[internali].dist(td);
166 (void)std::initializer_list<nil>
168 data.primitiveFieldRef()[internali] =
169 internalInfo[internali].data(td),
178 const bool valid = patchFaceInfo[
patchi][patchFacei].valid(td);
184 distance.boundaryFieldRef()[
patchi][patchFacei] =
185 patchFaceInfo[
patchi][patchFacei].dist(td) + small;
187 (void)std::initializer_list<nil>
189 data.boundaryFieldRef()[
patchi][patchFacei] =
190 patchFaceInfo[
patchi][patchFacei].data(td),
201 template<
class GeoMesh>
206 const scalar minFaceFraction,
207 GeometricField<scalar, GeoMesh>& distance
211 wave<FvWallInfo<wallPoint>>
217 FvFaceCellWave<FvWallInfo<wallPoint>>::defaultTrackingData_
222 template<
class GeoMesh>
227 const scalar minFaceFraction,
228 const label nCorrections,
229 GeometricField<scalar, GeoMesh>& distance
232 wave<FvWallInfo<wallFace>>
238 FvFaceCellWave<FvWallInfo<wallFace>>::defaultTrackingData_
243 template<
class GeoMesh>
248 const scalar minFaceFraction,
249 const label nCorrections,
250 GeometricField<scalar, GeoMesh>& distance
253 const List<labelPair> changedPatchAndFaces =
257 wave<FvWallInfo<wallPoint>>
260 changedPatchAndFaces,
263 FvFaceCellWave<FvWallInfo<wallPoint>>::defaultTrackingData_
266 wave<FvWallInfo<wallFace>>
269 changedPatchAndFaces,
272 FvFaceCellWave<FvWallInfo<wallFace>>::defaultTrackingData_
281 template<
class>
class WallLocation,
290 const scalar minFaceFraction,
291 GeometricField<scalar, GeoMesh>& distance,
292 GeometricField<DataType, GeoMesh>& data,
297 wave<FvWallInfo<WallLocation<wallPoint>>, TrackingData>
311 template<
class>
class WallLocation,
320 const scalar minFaceFraction,
321 const label nCorrections,
322 GeometricField<scalar, GeoMesh>& distance,
323 GeometricField<DataType, GeoMesh>& data,
327 wave<FvWallInfo<WallLocation<wallFace>>, TrackingData>
341 template<
class>
class WallLocation,
350 const scalar minFaceFraction,
351 const label nCorrections,
352 GeometricField<scalar, GeoMesh>& distance,
353 GeometricField<DataType, GeoMesh>& data,
357 const List<labelPair> changedPatchAndFaces =
361 wave<FvWallInfo<WallLocation<wallPoint>>, TrackingData>
364 changedPatchAndFaces,
371 wave<FvWallInfo<WallLocation<wallFace>>, TrackingData>
374 changedPatchAndFaces,
392 template<
class WallLocation>
399 template<
class DataType,
class GeoMesh,
class TrackingData>
404 const scalar minFaceFraction,
405 GeometricField<scalar, GeoMesh>& distance,
406 GeometricField<DataType, GeoMesh>& data,
411 calculate<WallLocationDataType<DataType>::template
type>
423 template<
class DataType,
class GeoMesh,
class TrackingData>
428 const scalar minFaceFraction,
429 const label nCorrections,
430 GeometricField<scalar, GeoMesh>& distance,
431 GeometricField<DataType, GeoMesh>& data,
435 correct<WallLocationDataType<DataType>::template
type>
447 template<
class DataType,
class GeoMesh,
class TrackingData>
452 const scalar minFaceFraction,
453 const label nCorrections,
454 GeometricField<scalar, GeoMesh>& distance,
455 GeometricField<DataType, GeoMesh>& data,
460 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.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
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.
const volVectorField & C() const
Return cell centres.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
const surfaceVectorField & Cf() const
Return face centres.
const GeometricBoundaryField< label, surfaceMesh > & polyFacesBf() const
Return face-poly-face addressing.
Takes a set of patches to start FvFaceCellWave from and computed the distance at patches and possibly...
label nTotalCells() const
Return total number of cells in decomposed mesh.
virtual const faceList & faces() const
Return raw faces.
const globalMeshData & globalData() const
Return parallel info.
virtual const pointField & points() const
Return raw points.
label nInternalFaces() const
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
bool valid(const PtrList< ModelType > &l)
label wave(const fvMesh &mesh, const List< labelPair > &changedPatchAndFaces, const label nCorrections, GeometricField< scalar, GeoMesh > &distance, TrackingData &td, GeometricField< DataType, 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, GeoMesh > &distance)
Correct distance data from patches.
label calculate(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, GeometricField< scalar, GeoMesh > &distance)
Calculate distance data from patches.
List< labelPair > getChangedPatchAndFaces(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction)
Get initial set of changed faces.
label calculateAndCorrect(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, const label nCorrections, GeometricField< scalar, GeoMesh > &distance)
Calculate and correct distance data from patches.
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.