37 #ifndef intersectionPatchToPatch_H
38 #define intersectionPatchToPatch_H
49 namespace patchToPatches
95 template<
class Po
intField>
114 const scalar magArea =
mag(
area);
115 const scalar magPArea =
mag(
p.area);
120 magArea == 0 ?
p.centre
122 : (magArea*
centre + magPArea*
p.centre)/(magArea + magPArea);
146 return os <<
p.area <<
p.centre;
152 return is >>
p.area >>
p.centre;
183 static_cast<const part&
>(a) ==
static_cast<const part&
>(
b)
196 return os << static_cast<const part&>(
c) <<
c.nbr;
202 return is >>
static_cast<part&
>(
c) >>
c.nbr;
315 virtual bool intersectFaces
321 const label srcFacei,
326 virtual void initialise
351 virtual label finalise
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
A 2-tuple for storing two objects of different types.
const Type2 & second() const
Return second.
const Type1 & first() const
Return first.
A face is a list of labels corresponding to mesh vertices.
static Tuple2< vector, point > areaAndCentreStabilised(const PointField &)
Return vector area and centre point given face points. Stabilised.
Traits class for primitives.
Class to generate coupling geometry between two primitive patches.
tmp< Field< Type > > tgtToSrc(const Field< Type > &tgtFld) const
Interpolate a target patch field to the source with no left.
bool reverse() const
Flag to indicate that the two patches are co-directional and.
Class to generate patchToPatch coupling geometry. A full geometric intersection is done between a fac...
TypeName("intersection")
Runtime type information.
const List< DynamicList< couple > > & tgtCouples() const
For each target face, the target and source areas for each.
const List< part > & srcErrorParts() const
For each source face, the area associated with mismatch.
const List< DynamicList< couple > > & srcCouples() const
For each source face, the source and target areas for each.
const List< scalar > & srcCoverage() const
Return the proportion of the source faces that are coupled.
~intersection()
Destructor.
static treeBoundBox srcBoxStatic(const face &srcFace, const pointField &srcPoints, const vectorField &srcPointNormals)
Get the bound box for a source face.
static int debugSrcFacei
Extra debugging for intersections between specific faces. Named.
const List< part > & srcEdgeParts() const
For each source edge, the non-coupled geometry associated.
intersection(const bool reverse)
Construct from components.
const List< scalar > & tgtCoverage() const
Return the proportion of the target faces that are coupled.
Triangulation of three-dimensional polygons.
A class for managing temporary objects without reference counting.
Standard boundBox + extra functionality for use in octree.
A triangular face using a FixedList of labels corresponding to mesh vertices.
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
const dimensionedScalar c
Speed of light in a vacuum.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
dimensioned< scalar > mag(const dimensioned< Type > &)
Structure to store the geometry associated with the coupling.
part nbr
The neighbour part.
friend bool operator==(const couple &a, const couple &b)
Equality comparison.
friend bool operator!=(const couple &a, const couple &b)
Inequality comparison.
friend Istream & operator>>(Istream &is, couple &c)
Input stream operator.
friend Ostream & operator<<(Ostream &os, const couple &c)
Output stream operator.
couple()
Default construct.
Structure to store the geometry associated with part of a patch.
void operator-=(const part &p)
Subtract another part from this one.
vector area
The area of this part.
friend Istream & operator>>(Istream &is, part &p)
Input stream operator.
friend Ostream & operator<<(Ostream &os, const part &p)
Output stream operator.
friend bool operator!=(const part &a, const part &b)
Inequality comparison.
point centre
The centre of this part.
part operator-() const
Negate this part.
friend bool operator==(const part &a, const part &b)
Equality comparison.
void operator+=(const part &p)
Add another part to this one.
Class to encapsulate the topology of a point within a triangle intersection.