33 namespace fvPatchDistWave
36 template<
class WallInfo,
class TrackingData>
46 template<
class WallInfo,
class TrackingData>
66 template<
class>
class PatchField,
74 const label nCorrections,
86 if (!calculate && nCorrections == 0)
return 0;
90 forAll(changedPatchAndFaces, changedFacei)
93 changedPatchAndFaces[changedFacei].
first();
94 const label patchFacei =
95 changedPatchAndFaces[changedFacei].second();
97 changedFacesInfo[changedFacei] =
135 wave.
setFaceInfo(changedPatchAndFaces, changedFacesInfo);
146 wave.
iterate(nCorrections - 1);
153 forAll(internalInfo, internali)
155 const bool valid = internalInfo[internali].valid(td);
157 if (calculate || valid)
162 internalInfo[internali].dist(td);
164 (void)std::initializer_list<nil>
167 internalInfo[internali].
data(td),
176 const bool valid = patchFaceInfo[
patchi][patchFacei].valid(td);
178 if (calculate || valid)
183 patchFaceInfo[
patchi][patchFacei].dist(td) + small;
185 (void)std::initializer_list<nil>
188 patchFaceInfo[
patchi][patchFacei].data(td),
199 template<
template<
class>
class PatchField,
class GeoMesh>
204 const scalar minFaceFraction,
209 wave<FvWallInfo<wallPoint>>
220 template<
template<
class>
class PatchField,
class GeoMesh>
225 const scalar minFaceFraction,
226 const label nCorrections,
230 wave<FvWallInfo<wallFace>>
241 template<
template<
class>
class PatchField,
class GeoMesh>
246 const scalar minFaceFraction,
247 const label nCorrections,
255 wave<FvWallInfo<wallPoint>>
258 changedPatchAndFaces,
264 wave<FvWallInfo<wallFace>>
267 changedPatchAndFaces,
279 template<
class>
class WallInfoData,
280 template<
class>
class PatchField,
288 const scalar minFaceFraction,
291 <
typename WallInfoData<wallPoint>::dataType, PatchField, GeoMesh>&
297 wave<WallInfoData<wallPoint>, TrackingData>
311 template<
class>
class WallInfoData,
312 template<
class>
class PatchField,
320 const scalar minFaceFraction,
321 const label nCorrections,
324 <
typename WallInfoData<wallPoint>::dataType, PatchField, GeoMesh>&
329 wave<WallInfoData<wallFace>, TrackingData>
343 template<
class>
class WallInfoData,
344 template<
class>
class PatchField,
352 const scalar minFaceFraction,
353 const label nCorrections,
356 <
typename WallInfoData<wallPoint>::dataType, PatchField, GeoMesh>&
365 wave<WallInfoData<wallPoint>, TrackingData>
368 changedPatchAndFaces,
375 wave<WallInfoData<wallFace>, TrackingData>
378 changedPatchAndFaces,
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
#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...
List< Type > & cellInfo()
Access cellInfo.
const surfaceVectorField & Cf() const
Return face centres.
label nInternalFaces() const
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...
void size(const label)
Override size to be inconsistent with allocated storage.
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 nTotalCells() const
Return total number of cells in decomposed mesh.
scalar distance(const vector &p1, const vector &p2)
Generic GeometricField class.
Holds information regarding type of cell. Used in inside/outside determination in cellClassification...
T & first()
Return the first element of the list.
const List< WallInfo > & getInternalInfo(const volScalarField &distance, FvFaceCellWave< WallInfo, TrackingData > &wave)
virtual const pointField & points() const
Return raw points.
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.
const globalMeshData & globalData() const
Return parallel info.
List< Type > & internalFaceInfo()
Access internalFaceInfo.
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
void setFaceInfo(const List< labelPair > &changedPatchAndFaces, const List< Type > &changedFacesInfo)
Set initial changed faces.
T * data()
Return a pointer to the first data element,.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
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.
virtual label faceToCell()
Propagate from face to cell. Returns total number of cells.
const volVectorField & C() const
Return cell centres.
virtual label iterate(const label maxIter)
Iterate until no changes or maxIter reached. Returns actual.
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.
A zero-sized class without any storage. Used, for example, in HashSet.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.