42void Foam::faceCollapser::insert
44 const labelList& elems,
45 const label excludeElem,
51 if (elems[i] != excludeElem)
60Foam::label Foam::faceCollapser::findEdge
70 label edgeI = edgeLabels[i];
72 const edge&
e = edges[edgeI];
76 (
e[0] == v0 &&
e[1] == v1)
77 || (
e[0] == v1 &&
e[1] == v0)
85 <<
" and " << v1 <<
" in edge labels " << edgeLabels
95void Foam::faceCollapser::filterFace
102 const face&
f = mesh_.faces()[facei];
103 const labelList& fEdges = mesh_.faceEdges()[facei];
120 label edgeI =
findEdge(mesh_.edges(), fEdges, v0, v1);
123 splitEdges.find(edgeI);
125 if (edgeFnd != splitEdges.end())
133 if (v0 == mesh_.edges()[edgeI].start())
137 newFace.append(extraVerts[i]);
144 newFace.append(extraVerts[i]);
149 face newF(newFace.shrink());
160 if (mesh_.isInternalFace(facei))
162 nei = mesh_.faceNeighbour()[facei];
166 patchi = mesh_.boundaryMesh().whichPatch(facei);
170 label zoneID = mesh_.faceZones().whichZone(facei);
172 bool zoneFlip =
false;
176 const faceZone& fZone = mesh_.faceZones()[zoneID];
178 zoneFlip = fZone.flipMap()[fZone.whichFace(facei)];
187 mesh_.faceOwner()[facei],
203Foam::faceCollapser::faceCollapser(
const polyMesh&
mesh)
220 const edgeList& edges = mesh_.edges();
221 const faceList& faces = mesh_.faces();
242 const face&
f = faces[facei];
244 const label fpA = fpStart[i];
245 const label fpB = fpEnd[i];
250 Pout<<
"Face:" <<
f <<
" collapsed to fp:" << fpA <<
' ' << fpB
251 <<
" with points:" << pA <<
' ' << pB
263 dist[fpB] =
magSqr(pB - pA);
270 label fp =
f.fcIndex(fpMin1);
279 if (w <= dist[fpMin1])
282 w = dist[fpMin1] + 1
e-6*(dist[fpB] - dist[fpA]);
286 pA +
Foam::sqrt(w / (dist[fpB] - dist[fpA]))*(pB - pA)
289 Pout<<
"Adapting position of vertex " <<
f[fp] <<
" on face "
290 <<
f <<
" from " << near.
point() <<
" to " << newPoint
315 fp =
f.fcIndex(fpMin1);
332 if (w <= dist[fpMin1])
335 w = dist[fpMin1] + 1
e-6*(dist[fpB] - dist[fpA]);
339 pA +
Foam::sqrt(w / (dist[fpB] - dist[fpA]))*(pB - pA)
342 Pout<<
"Adapting position of vertex " <<
f[fp] <<
" on face "
343 <<
f <<
" from " << near.point() <<
" to " << newPoint
346 near.setPoint(newPoint);
374 if (dist.indices()[dist.size()-1] != fpB)
376 OFstream str(
"conflictingFace.obj");
380 <<
"Trying to collapse face:" << facei <<
" vertices:" <<
f
381 <<
" to edges between vertices " <<
f[fpA] <<
" and "
382 <<
f[fpB] <<
" but " <<
f[fpB] <<
" does not seem to be the"
383 <<
" vertex furthest away from " <<
f[fpA] <<
endl
384 <<
"Dumped conflicting face to obj file conflictingFace.obj"
390 Pout<<
"Face:" <<
f <<
" fpA:" << fpA <<
" fpB:" << fpB <<
nl;
395 label fp = dist.indices()[i];
397 Pout<<
" fp:" << fp <<
" distance:" << dist[i] <<
nl;
402 const labelList& fEdges = mesh_.faceEdges()[facei];
415 label sorted0 = sortedFp[fp];
416 label sorted1 = sortedFp[fp1];
424 if (sorted0 < sorted1)
428 for (label j = sorted0+1; j < sorted1; j++)
430 edgePoints.append(
f[dist.indices()[j]]);
437 for (label j = sorted1+1; j < sorted0; j++)
439 edgePoints.append(
f[dist.indices()[j]]);
443 if (edgePoints.size())
447 label edgeI =
findEdge(edges, fEdges,
f[fp],
f[fp1]);
449 const edge&
e = edges[edgeI];
451 if (fpToFp1 == (
f[fp] ==
e.start()))
453 splitEdges.insert(edgeI, edgePoints);
458 splitEdges.insert(edgeI, edgePoints);
462 insert(edgeFaces[edgeI], facei, affectedFaces);
469 Pout<<
"Split edge:" << iter.key()
470 <<
" verts:" << mesh_.edges()[iter.key()]
473 const labelList& edgePoints = iter.val();
477 Pout<<
" " << edgePoints[i] <<
nl;
493 affectedFaces.erase(facei);
501 for (
const label facei : affectedFaces)
503 filterFace(splitEdges, facei, meshMod);
labelList faceLabels(nFaceLabels)
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
void append(const T &val)
Copy append an element to the end of this list.
DynamicList< T, SizeMin > & shrink()
Calls shrink_to_fit() and returns a reference to the DynamicList.
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
bool erase(const iterator &iter)
Erase an entry specified by given iterator.
void append(const T &val)
Append an element at the end of the list.
A HashTable to objects of type <T> with a label key.
typename parent_type::const_iterator const_iterator
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
void setPoint(const point_type &p)
Set the point.
const point_type & point() const noexcept
Return the point, no checks.
A list that is sorted upon construction or when explicitly requested with the sort() method.
const labelList & indices() const noexcept
Return the list of sorted indices. Updated every sort.
void sort()
Forward (stable) sort the list (if changed after construction).
label rcIndex(const label i) const noexcept
The reverse circular index. The previous index in the list which returns to the last at the beginning...
void size(const label n)
Older name for setAddressableSize.
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
scalar distSqr(const Vector< Cmpt > &v2) const
The L2-norm distance squared from another vector. The magSqr() of the difference.
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
void setRefinement(const labelList &faceLabels, const labelList &fpA, const labelList &fpB, polyTopoChange &) const
Collapse faces along endpoints. Play commands into.
A subset of mesh faces organised as a primitive patch.
A face is a list of labels corresponding to mesh vertices.
PointHit< Point > nearestDist(const Point &p) const
Return nearest distance to line from a given point.
Mesh consisting of general polyhedral cells.
Class describing modification of a face.
Class describing modification of a point.
Class containing data for face removal.
Direct mesh changes based on v1.3 polyTopoChange syntax.
label setAction(const topoAction &action)
For compatibility with polyTopoChange: set topological action.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
void set(List< bool > &bools, const labelUList &locations)
Set the listed locations (assign 'true').
List< edge > edgeList
List of edge.
List< labelList > labelListList
List of labelList.
PointHit< point > pointHit
A PointHit with a 3D point.
List< label > labelList
A List of labels.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
void reverse(UList< T > &list, const label n)
Reverse the first n elements of the list.
List< face > faceList
List of faces.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
line< point, const point & > linePointRef
A line using referred points.
errorManip< error > abort(error &err)
vector point
Point is a vector.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a).
nonInt insert("surfaceSum(((S|magSf)*S)")
#define forAll(list, i)
Loop across all elements in list.
#define forAllReverse(list, i)
Reverse loop across all elements in list.
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.