Go to the documentation of this file.
58 if (orient && area_ < 0)
120 const triangle2D& triB,
125 const triangle2D& triA = *
this;
138 static int nInter = 0;
139 OFstream
os(
"intersection-tris-"+
Foam::name(nInter++)+
".obj");
151 static FixedList<vector2D, 7> work2;
164 for (
int i0 = 0; i0 <= 2; ++i0)
168 Info<<
"\nstarting points:";
169 for (
int i = 0; i < nPoint2; ++i)
177 const label i1 = (i0 + 1) % 3;
184 for (
int j = 0; j < nPoint2; ++j)
187 const vector2D& p1 = work2[(j+1) % nPoint2];
189 if (triangle2D(c0,
c1,
p0).order() == 1)
191 if (triangle2D(c0,
c1, p1).order() == 1)
193 work_[nPoint++] = p1;
197 if (lineIntersectionPoint(
p0, p1, c0,
c1, intersection))
199 work_[nPoint++] = intersection;
205 if (triangle2D(c0,
c1, p1).order() == 1)
207 if (lineIntersectionPoint(
p0, p1, c0,
c1, intersection))
209 work_[nPoint++] = intersection;
211 work_[nPoint++] = p1;
224 for (
int i = 0; i < nPoint; ++i)
226 os <<
"v " << work_[i].x() <<
" " << work_[i].y() <<
" 0" <<
nl;
229 for (
int i = 0; i < nPoint; ++i)
231 os <<
" " << (i + 1);
236 Info<<
"Intersection points:" <<
endl;
237 for (
int i = 0; i < nPoint; ++i)
248 for (
int i = 0; i < nPoint - 1; ++i)
250 twoArea += work_[i].x()*work_[i+1].y();
251 twoArea -= work_[i].y()*work_[i+1].x();
254 twoArea += work_[nPoint-1].x()*work_[0].y();
255 twoArea -= work_[nPoint-1].y()*work_[0].x();
257 centre += work_[nPoint - 1];
258 centre /= scalar(nPoint);
270 interArea(triB, dummyCentre,
area);
281 for (
int i = 0; i < 3; ++i)
283 int i1 = (i + 1) % 3;
285 for (
int j = 0; j < 3; ++j)
287 int j1 = (j + 1) % 3;
289 if (lineIntersects(triA[i], triA[i1], triB[j], triB[
j1]))
297 (nClosePoints(triA, triB) == 3)
bool isClose(const Vector2D< Cmpt > &b, const scalar tol=1e-10) const
virtual const fileName & name() const
A class for managing temporary objects.
static constexpr const zero Zero
Ostream & endl(Ostream &os)
const dimensionedScalar c1
const dimensionedScalar b
bool overlaps(const triangle2D &triB) const
OBJstream os(runTime.globalPath()/outputName)
Vector< scalar > vector
A scalar version of the templated Vector.
void interArea(const triangle2D &triB, vector2D ¢re, scalar &area) const
Output to file stream, using an OSstream.
bool match(const UList< wordRe > &patterns, const std::string &text)
Vector2D< scalar > vector2D
A 2D vector of scalars obtained from the generic Vector2D.
bool contains(const triangle2D &tri) const
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A 1D vector of objects of type <T> with a fixed length <N>.
const dimensionedScalar e
label snapClosePoints(const triangle2D &triB)
const dimensionedScalar c
word name(const expressions::valueTypeCode typeCode)
const volScalarField & p0
triangle2D(const vector2D &a, const vector2D &b, const vector2D &c, const bool orient=false)
dimensionedScalar j1(const dimensionedScalar &ds)