46void Foam::thresholdCellFaces::calculate
48 const scalarField&
field,
49 const scalar lowerThreshold,
50 const scalar upperThreshold,
51 const bool triangulate
90 DynamicList<face>
surfFaces(0.5 * mesh_.nFaces());
91 DynamicList<label> surfCells(
surfFaces.size());
93 labelList oldToNewPoints(origPoints.size(), -1);
97 for (label facei = 0; facei < mesh_.nInternalFaces(); ++facei)
102 if (
field[own[facei]] > lowerThreshold)
104 if (
field[nei[facei]] < lowerThreshold)
109 else if (
field[nei[facei]] > lowerThreshold)
115 if (
field[own[facei]] < upperThreshold)
117 if (
field[nei[facei]] > upperThreshold)
122 else if (
field[nei[facei]] < upperThreshold)
130 const face&
f = origFaces[facei];
132 for (
const label pointi :
f)
134 if (oldToNewPoints[pointi] == -1)
136 oldToNewPoints[pointi] =
nPoints++;
151 surfFace =
f.reverseFace();
178 const polyPatch&
p =
bMesh[patchi];
192 label facei =
p.start();
199 field[own[facei]] > lowerThreshold
200 &&
field[own[facei]] < upperThreshold
203 const face&
f = origFaces[facei];
204 for (
const label pointi :
f)
206 if (oldToNewPoints[pointi] == -1)
208 oldToNewPoints[pointi] =
nPoints++;
212 label
cellId = own[facei];
232 zone.size() =
surfFaces.size() - zone.start();
248 forAll(oldToNewPoints, pointi)
250 if (oldToNewPoints[pointi] >= 0)
252 surfPoints[oldToNewPoints[pointi]] = origPoints[pointi];
262 meshCells_.transfer(surfCells);
272 const scalar lowerThreshold,
273 const scalar upperThreshold,
274 const bool triangulate
279 if (lowerThreshold > upperThreshold)
282 << lowerThreshold <<
" > " << upperThreshold <<
endl;
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
surfZoneList & storedZones()
const surfZoneList & surfZones() const
pointField & storedPoints()
virtual label triangulate()
const List< face > & surfFaces() const
List< face > & storedFaces()
label nFaces() const noexcept
void size(const label n)
Older name for setAddressableSize.
static bool & parRun() noexcept
Test if this a parallel run.
Mesh consisting of general polyhedral cells.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
virtual const faceList & faces() const
Return raw faces.
virtual const labelList & faceOwner() const
Return face owner.
virtual const labelList & faceNeighbour() const
Return face neighbour.
virtual const pointField & points() const
Return raw points.
Selects the mesh cell faces specified by a threshold value. Non-triangulated by default.
thresholdCellFaces(const polyMesh &mesh, const scalarField &field, const scalar lowerThreshold, const scalar upperThreshold, const bool triangulate=false)
Construct from mesh, field and threshold values.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define WarningInFunction
Report a warning using Foam::Warning.
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values within a list.
List< label > labelList
A List of labels.
List< surfZone > surfZoneList
List of surfZone.
List< face > faceList
List of faces.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
PrimitivePatch< List< face >, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points).
vectorField pointField
pointField is a vectorField.
#define forAll(list, i)
Loop across all elements in list.