36void Foam::patchWave::setChangedFaces
39 DynamicList<label>& changedFaces,
40 DynamicList<wallPoint>& faceDist
49 const polyPatch&
patch =
mesh.boundaryMesh()[patchi];
51 const auto areaFraction(
patch.areaFraction());
53 const auto faceCentres(
patch.faceCentres());
55 forAll(faceCentres, patchFacei)
60 || (areaFraction()[patchFacei] > 0.5)
63 changedFaces.append(
patch.start() + patchFacei);
64 faceDist.append(wallPoint(faceCentres[patchFacei], 0.0));
70 for (
const label facei : sourceIDs_)
72 changedFaces.append(facei);
78 mesh.faceCentres()[facei],
98 scalar dist =
cellInfo[celli].distSqr();
100 if (
cellInfo[celli].valid(waveInfo.data()))
106 distance_[celli] = dist;
113 forAll(patchDistance_, patchi)
120 patchDistance_.set(patchi, patchDistPtr);
124 forAll(patchField, patchFacei)
126 label meshFacei =
patch.start() + patchFacei;
128 scalar dist = faceInfo[meshFacei].distSqr();
130 if (faceInfo[meshFacei].valid(waveInfo.data()))
134 patchField[patchFacei] =
Foam::sqrt(dist) + SMALL;
138 patchField[patchFacei] = dist;
152 const polyMesh&
mesh,
154 const bool correctWalls,
160 correctWalls_(correctWalls),
162 distance_(
mesh.nCells()),
163 patchDistance_(
mesh.boundaryMesh().size()),
164 sourceIDs_(sourceIDs)
182 label nPatch = sumPatchSize(patchIDs_) + sourceIDs_.size();
188 setChangedFaces(patchIDs_, changedFaces, faceDist);
196 mesh().globalData().nTotalCells()+1
200 nUnset_ = getValues(waveInfo);
212 patchIDs_.sortedToc(),
221 correctBoundaryFaceCells
228 correctBoundaryPointCells
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void setSize(label n)
Alias for resize().
A HashTable to objects of type <T> with a label key.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
label sumPatchSize(const labelHashSet &patchIDs) const
Sum of patch sizes (out of supplied subset of patches).
const polyMesh & mesh() const
Access mesh.
void correctBoundaryFaceCells(const labelHashSet &patchIDs, scalarField &wallDistCorrected, Map< label > &nearestFace) const
Correct all cells connected to boundary (via face). Sets values in.
void correctBoundaryCells(const labelList &patchIDs, const bool doPointCells, scalarField &wallDistCorrected, Map< label > &nearestFace) const
Correct all cells connected to any of the patches in patchIDs. Sets.
void correctBoundaryPointCells(const labelHashSet &patchIDs, scalarField &wallDistCorrected, Map< label > &nearestFace) const
Correct all cells connected to wall (via point). Sets values in.
static bool useCombinedWallPatch
Use combined-wall-patches wall distance v.s. v2406 per-patch distance. Default is true.
Holds information regarding type of cell. Used in inside/outside determination in cellClassification.
virtual void correct()
Correct for mesh geom/topo changes.
virtual ~patchWave()
Destructor.
patchWave(const polyMesh &mesh, const labelHashSet &patchIDs, bool correctWalls=true, const labelList &sourceFaceIDs=labelList())
Construct from mesh and patches to initialize to 0 and flag.
Mesh consisting of general polyhedral cells.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
A patch is a list of labels that address the faces in the global face list.
const std::string patch
OpenFOAM patch number as a std::string.
List< label > labelList
A List of labels.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar sqrt(const dimensionedScalar &ds)
#define forAll(list, i)
Loop across all elements in list.