34template<
class TrackingData>
35inline bool Foam::pointEdgeCollapse::update
37 const pointEdgeCollapse&
w2,
54 if (
w2.collapseIndex_ == -1 || collapseIndex_ == -1)
60 if (
w2.collapsePriority_ < collapsePriority_)
64 else if (
w2.collapsePriority_ > collapsePriority_)
71 if (
w2.collapseIndex_ < collapseIndex_)
76 else if (
w2.collapseIndex_ == collapseIndex_)
78 bool identicalPoint = samePoint(
w2.collapsePoint_);
80 bool nearer = (
magSqr(
w2.collapsePoint_) <
magSqr(collapsePoint_));
101inline bool Foam::pointEdgeCollapse::samePoint(
const point& pt)
const
103 bool isLegal1 = (
cmptMin(collapsePoint_) < 0.5*GREAT);
104 bool isLegal2 = (
cmptMin(pt) < 0.5*GREAT);
106 if (isLegal1 && isLegal2)
108 return mag(collapsePoint_ - pt) < 1
e-9;
112 return isLegal1 == isLegal2;
121 collapsePoint_(GREAT, GREAT, GREAT),
123 collapsePriority_(-2)
129 const point& collapsePoint,
130 const label collapseIndex,
131 const label collapsePriority
134 collapsePoint_(collapsePoint),
135 collapseIndex_(collapseIndex),
142template<
class TrackingData>
145 return collapseIndex_ != -2;
149template<
class TrackingData>
153 const label patchPointi,
158 collapsePoint_ -= coord;
162template<
class TrackingData>
175template<
class TrackingData>
179 const label patchPointi,
185 collapsePoint_ += coord;
190template<
class TrackingData>
206template<
class TrackingData>
216 return update(newPointInfo, tol,
td);
221template<
class TrackingData>
229 return update(newPointInfo, tol,
td);
234template<
class TrackingData>
249template<
class TrackingData>
262inline bool Foam::pointEdgeCollapse::operator==
268 collapseIndex_ ==
rhs.collapseIndex_
269 && collapsePriority_ ==
rhs.collapsePriority_
270 && samePoint(
rhs.collapsePoint_);
274inline bool Foam::pointEdgeCollapse::operator!=
279 return !(*
this ==
rhs);
Determines length of string of edges walked to point.
bool updateEdge(const polyMesh &mesh, const label edgeI, const label pointi, const pointEdgeCollapse &pointInfo, const scalar tol, TrackingData &td)
Influence of point on edge.
pointEdgeCollapse()
Default construct.
bool equal(const pointEdgeCollapse &, TrackingData &) const
Test for equality, with TrackingData.
void leaveDomain(const polyPatch &patch, const label patchPointi, const point &pos, TrackingData &td)
Convert origin to relative vector to leaving point.
void transform(const tensor &rotTensor, TrackingData &td)
Apply rotation matrix to origin.
void enterDomain(const polyPatch &patch, const label patchPointi, const point &pos, TrackingData &td)
Convert relative origin to absolute by adding entering point.
bool valid(TrackingData &td) const
Changed or contains original (invalid) value.
label collapsePriority() const
bool updatePoint(const polyMesh &mesh, const label pointi, const label edgeI, const pointEdgeCollapse &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on point.
const point & collapsePoint() const
label collapseIndex() const
Mesh consisting of general polyhedral cells.
A patch is a list of labels that address the faces in the global face list.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
vector point
Point is a vector.
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)