34 template<
class TrackingData>
35 inline bool Foam::externalPointEdgePoint::update
38 const externalPointEdgePoint&
w2,
43 const scalar dist2 = pt.distSqr(
w2.origin());
49 origin_ =
w2.origin();
54 const scalar
diff = distSqr_ - dist2;
62 if ((
diff < SMALL) || ((distSqr_ > SMALL) && (
diff/distSqr_ < tol)))
71 origin_ =
w2.origin();
78 template<
class TrackingData>
79 inline bool Foam::externalPointEdgePoint::update
81 const externalPointEdgePoint&
w2,
89 distSqr_ =
w2.distSqr();
90 origin_ =
w2.origin();
95 const scalar
diff = distSqr_ -
w2.distSqr();
103 if ((
diff < SMALL) || ((distSqr_ > SMALL) && (
diff/distSqr_ < tol)))
111 distSqr_ =
w2.distSqr();
112 origin_ =
w2.origin();
141 template<
class TrackingData>
149 template<
class TrackingData>
165 if ((distSqr() > SMALL) && ((
diff/distSqr()) < tol))
177 template<
class TrackingData>
180 const polyPatch&
patch,
181 const label patchPointi,
190 template<
class TrackingData>
201 template<
class TrackingData>
205 const label patchPointi,
215 template<
class TrackingData>
226 return update(td.points_[pointi], edgeInfo, tol, td);
230 template<
class TrackingData>
240 return update(td.points_[pointi], newPointInfo, tol, td);
244 template<
class TrackingData>
252 return update(newPointInfo, tol, td);
256 template<
class TrackingData>
268 return update(
e.centre(td.points_), pointInfo, tol, td);
272 template<
class TrackingData>
285 inline bool Foam::externalPointEdgePoint::operator==
290 return origin_ == rhs.origin_ && distSqr_ == rhs.distSqr_;
294 inline bool Foam::externalPointEdgePoint::operator!=
299 return !(*
this == rhs);
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
bool updatePoint(const polyMesh &mesh, const label pointi, const label edgeI, const externalPointEdgePoint &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on point.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
bool valid(TrackingData &td) const
Changed or contains original (invalid) value.
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
const dimensionedScalar e
Elementary charge.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
void transform(const tensor &rotTensor, TrackingData &td)
Apply rotation matrix to origin.
externalPointEdgePoint()
Default construct.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
bool sameGeometry(const externalPointEdgePoint &, const scalar tol, TrackingData &td) const
Check for identical geometrical data (eg, cyclics checking)
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
vector point
Point is a vector.
const std::string patch
OpenFOAM patch number as a std::string.
Mesh consisting of general polyhedral cells.
bool equal(const externalPointEdgePoint &, TrackingData &td) const
Test for equality, with TrackingData.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
A patch is a list of labels that address the faces in the global face list.
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Tensor of scalars, i.e. Tensor<scalar>.
void leaveDomain(const polyPatch &patch, const label patchPointi, const point &pos, TrackingData &td)
Convert origin to relative vector to leaving point (= point coordinate)
bool updateEdge(const polyMesh &mesh, const label edgeI, const label pointi, const externalPointEdgePoint &pointInfo, const scalar tol, TrackingData &td)
Influence of point on edge.
void enterDomain(const polyPatch &patch, const label patchPointi, const point &pos, TrackingData &td)
Convert relative origin to absolute by adding entering point.