39#ifndef pointEdgeStructuredWalk_H
40#define pointEdgeStructuredWalk_H
85 template<
class TrackingData>
104 const point& previousPoint,
107 const label
index = -1
134 inline bool inZone()
const;
140 template<
class TrackingData>
141 inline bool valid(TrackingData&
td)
const;
144 template<
class TrackingData>
154 template<
class TrackingData>
158 const label patchPointi,
164 template<
class TrackingData>
168 const label patchPointi,
174 template<
class TrackingData>
182 template<
class TrackingData>
195 template<
class TrackingData>
207 template<
class TrackingData>
216 template<
class TrackingData>
228 template<
class TrackingData>
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Determines length of string of edges walked to point.
bool operator==(const pointEdgeStructuredWalk &) const
Test for equality.
bool updatePoint(const polyMesh &mesh, const label pointi, const label edgeI, const pointEdgeStructuredWalk &edgeInfo, const scalar tol, TrackingData &td)
Influence of edge on point.
bool sameGeometry(const pointEdgeStructuredWalk &, const scalar tol, TrackingData &td) const
Check for identical geometrical data (eg, cyclics checking).
label index() const
Index (if any) associated with data.
pointEdgeStructuredWalk()
Default construct.
bool equal(const pointEdgeStructuredWalk &, TrackingData &) const
Test for equality, with TrackingData.
bool operator!=(const pointEdgeStructuredWalk &) const
Test for inequality.
const vector & data() const
Tracking data.
friend Istream & operator>>(Istream &, pointEdgeStructuredWalk &)
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.
scalar dist() const
The distance information.
friend Ostream & operator<<(Ostream &, const pointEdgeStructuredWalk &)
bool inZone() const
True if starting point is valid (ie, not point::max).
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.
bool updateEdge(const polyMesh &mesh, const label edgeI, const label pointi, const pointEdgeStructuredWalk &pointInfo, const scalar tol, TrackingData &td)
Influence of point on edge.
Mesh consisting of general polyhedral cells.
A patch is a list of labels that address the faces in the global face list.
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
dimensionedScalar pos(const dimensionedScalar &ds)
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
Istream & operator>>(Istream &, directionInfo &)
vector point
Point is a vector.
A template class to specify that a data type can be considered as being contiguous in memory.