Go to the documentation of this file.
55 template<
class TrackingData>
63 template<
class TrackingData>
77 template<
class TrackingData>
82 const label patchFacei,
83 const point& faceCentre,
93 template<
class TrackingData>
98 const label patchFacei,
99 const point& faceCentre,
107 index_ = (
f.size() - index_) %
f.size();
113 template<
class TrackingData>
124 template<
class TrackingData>
128 const label thisCelli,
129 const label neighbourFacei,
145 if (neighbourInfo.
index() == -2)
150 else if (neighbourInfo.
index() == -1)
185 label v0 =
f[neighbourInfo.
index()];
186 label v1 =
f[(neighbourInfo.
index() + 1) %
f.size()];
198 n_ = neighbourInfo.
n();
205 template<
class TrackingData>
209 const label thisFacei,
210 const label neighbourCelli,
227 if (neighbourInfo.
index() >= 0)
237 neighbourInfo.
index()
246 n_ = neighbourInfo.
n();
253 template<
class TrackingData>
270 index_ = neighbourInfo.
index();
272 n_ = neighbourInfo.
n();
279 template<
class TrackingData>
292 inline bool Foam::directionInfo::operator==
297 return index_ == rhs.index_ && n_ == rhs.n_;
301 inline bool Foam::directionInfo::operator!=
306 return !(*
this == rhs);
static constexpr const zero Zero
Global zero (0)
void enterDomain(const polyMesh &, const polyPatch &, const label patchFacei, const point &faceCentre, TrackingData &td)
Reverse of leaveDomain.
bool equal(const directionInfo &, TrackingData &td) const
Test for equality, with TrackingData.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
static bool test(const UList< face > &faces)
bool updateFace(const polyMesh &, const label thisFacei, const label neighbourCelli, const directionInfo &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring cell.
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Mesh consisting of general polyhedral cells.
bool valid(TrackingData &td) const
Changed or contains original (invalid) value.
void transform(const polyMesh &, const tensor &, TrackingData &td)
Apply rotation matrix to any coordinates.
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.
const labelListList & faceEdges() const
Holds direction in which to split cell (in fact a local coordinate axes). Information is a label and ...
virtual const faceList & faces() const
Return raw faces.
const std::string patch
OpenFOAM patch number as a std::string.
void leaveDomain(const polyMesh &, const polyPatch &, const label patchFacei, const point &faceCentre, TrackingData &td)
Convert any absolute coordinates into relative to (patch)face.
const dimensionedScalar e
Elementary charge.
directionInfo()
Default construct, index=-1, vector::zero.
A face is a list of labels corresponding to mesh vertices.
bool sameGeometry(const polyMesh &, const directionInfo &, const scalar, TrackingData &td) const
Check for identical geometrical data (eg, cyclics checking)
bool updateCell(const polyMesh &, const label thisCelli, const label neighbourFacei, const directionInfo &neighbourInfo, const scalar tol, TrackingData &td)
Influence of neighbouring face.