triangle2D.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2020-2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "triangle2D.H"
29#include "OFstream.H"
30
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
34
35Foam::scalar Foam::triangle2D::relTol = 1e-8;
36
37Foam::scalar Foam::triangle2D::absTol = 1e-10;
38
39Foam::FixedList<Foam::vector2D, 7> Foam::triangle2D::work_
40(
41 vector2D::zero
42);
43
44
45// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46
48(
49 const vector2D& a,
50 const vector2D& b,
51 const vector2D& c,
52 const bool orient
53)
54:
55 FixedList<vector2D, 3>({a, b, c}),
56 area_(area(a, b, c))
57{
58 if (orient && area_ < 0)
59 {
60 // Inverted triangle
61 triangle2D& tri = *this;
62 vector2D tmp = tri[0];
63 tri[0] = tri[2];
64 tri[2] = tmp;
65
66 area_ = mag(area_);
67 }
68}
69
70
72(
73 const vector& a3d,
74 const vector& b3d,
75 const vector& c3d,
76 const vector& axis1,
77 const vector& axis2,
78 const bool orient
79)
80:
82 (
83 vector2D(axis1 & a3d, axis2 & a3d),
84 vector2D(axis1 & b3d, axis2 & b3d),
85 vector2D(axis1 & c3d, axis2 & c3d),
86 orient
87 )
88{}
89
90
91// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92
94{
95 label nSnapped = 0;
96 triangle2D& triA = *this;
97
98 FixedList<bool, 3> match(true);
99
100 for (auto& a : triA)
101 {
102 forAll(triB, tb)
103 {
104 if (match[tb] && a.isClose(triB[tb]))
105 {
106 a = triB[tb];
107 match[tb] = false;
108 ++nSnapped;
109 break;
110 }
111 }
112 }
113
114 return nSnapped;
115}
116
117
119(
120 const triangle2D& triB,
121 vector2D& centre,
122 scalar& area
123) const
124{
125 const triangle2D& triA = *this;
126
127 // Potential short-cut if the triangles are the same-ish. Happens rarely
128 // for moving mesh cases.
129 // if (nClosePoints(triA, triB) == 3)
130 // {
131 // centre = triA.centre();
132 // area = triA.area();
133 // return;
134 // }
135
136 if (debug)
137 {
138 static int nInter = 0;
139 OFstream os("intersection-tris-"+Foam::name(nInter++)+".obj");
140 writeOBJ(os, triA, 0);
141 writeOBJ(os, triB, 3);
142 Info<< "written " << os.name() << endl;
143 }
144
145
146 // Use work_ to store intersections
147
148 // Current number of intersections
149 int nPoint = 0;
150
151 static FixedList<vector2D, 7> work2;
152
153 // Clipped triangle starts as triA
154 work2[0] = triA[0];
155 work2[1] = triA[1];
156 work2[2] = triA[2];
157
158 int nPoint2 = 3;
159
160
162
163 // Cut triA using triB's edges as clipping planes
164 for (int i0 = 0; i0 <= 2; ++i0)
165 {
166 if (debug)
167 {
168 Info<< "\nstarting points:";
169 for (int i = 0; i < nPoint2; ++i)
170 {
171 Info<< work2[i];
172 }
173 Info<< endl;
174 }
175
176 // Clipping plane points
177 const label i1 = (i0 + 1) % 3;
178 const vector2D& c0 = triB[i0];
179 const vector2D& c1 = triB[i1];
180
181 nPoint = 0;
182
183 // Test against all intersection poly points
184 for (int j = 0; j < nPoint2; ++j)
185 {
186 const vector2D& p0 = work2[j];
187 const vector2D& p1 = work2[(j+1) % nPoint2];
188
189 if (triangle2D(c0, c1, p0).order() == 1)
190 {
191 if (triangle2D(c0, c1, p1).order() == 1)
192 {
193 work_[nPoint++] = p1;
194 }
195 else
196 {
197 if (lineIntersectionPoint(p0, p1, c0, c1, intersection))
198 {
199 work_[nPoint++] = intersection;
200 }
201 }
202 }
203 else
204 {
205 if (triangle2D(c0, c1, p1).order() == 1)
206 {
207 if (lineIntersectionPoint(p0, p1, c0, c1, intersection))
208 {
209 work_[nPoint++] = intersection;
210 }
211 work_[nPoint++] = p1;
212 }
213 }
214 }
215
216 work2 = work_;
217 nPoint2 = nPoint;
218 }
219
220 if (debug)
221 {
222 static int n = 0;
223 OFstream os("intersection-poly-"+Foam::name(n++)+".obj");
224 for (int i = 0; i < nPoint; ++i)
225 {
226 os << "v " << work_[i].x() << " " << work_[i].y() << " 0" << nl;
227 }
228 os << "l";
229 for (int i = 0; i < nPoint; ++i)
230 {
231 os << " " << (i + 1);
232 }
233 os << " 1" << endl;
234 Info<< "written " << os.name() << endl;
235
236 Info<< "Intersection points:" << endl;
237 for (int i = 0; i < nPoint; ++i)
238 {
239 Info<< " " << work_[i] << endl;
240 }
241 }
242
243 // Calculate the area of the clipped triangle
244 scalar twoArea = 0;
245 centre = vector2D::zero;
246 if (nPoint)
247 {
248 for (int i = 0; i < nPoint - 1; ++i)
249 {
250 twoArea += work_[i].x()*work_[i+1].y();
251 twoArea -= work_[i].y()*work_[i+1].x();
252 centre += work_[i];
253 }
254 twoArea += work_[nPoint-1].x()*work_[0].y();
255 twoArea -= work_[nPoint-1].y()*work_[0].x();
256
257 centre += work_[nPoint - 1];
258 centre /= scalar(nPoint);
259 }
260
261 area = 0.5*twoArea;
262}
263
264
265Foam::scalar Foam::triangle2D::interArea(const triangle2D& triB) const
266{
267 vector2D dummyCentre(Zero);
268 scalar area;
269
270 interArea(triB, dummyCentre, area);
271
272 return area;
273}
274
275
277{
278 const triangle2D& triA = *this;
279
280 // True if any of triB's edges intersect a triA edge
281 for (int i = 0; i < 3; ++i)
282 {
283 int i1 = (i + 1) % 3;
284
285 for (int j = 0; j < 3; ++j)
286 {
287 int j1 = (j + 1) % 3;
288
289 if (lineIntersects(triA[i], triA[i1], triB[j], triB[j1]))
290 {
291 return true;
292 }
293 }
294 }
295
296 return
297 (nClosePoints(triA, triB) == 3) // same tri
298 || triA.contains(triB) // triA contains triB
299 || triB.contains(triA); // triB contains triA
300}
301
302
303// ************************************************************************* //
label n
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: FixedList.H:81
Output to file stream, using an OSstream.
Definition: OFstream.H:57
virtual const fileName & name() const
Read/write access to the name of the stream.
Definition: OSstream.H:107
int overlaps
Flag to control which overlap calculations are performed.
Definition: PDRparams.H:97
const Cmpt & x() const
Access to the vector x component.
Definition: Vector2DI.H:66
Foam::intersection.
Definition: intersection.H:53
A class for managing temporary objects.
Definition: tmp.H:65
2-D triangle and queries
Definition: triangle2D.H:55
static scalar absTol
Absolute tolerance.
Definition: triangle2D.H:73
label order() const
Definition: triangle2DI.H:199
static int debug
Definition: triangle2D.H:67
label snapClosePoints(const triangle2D &triB)
Definition: triangle2D.C:93
void interArea(const triangle2D &triB, vector2D &centre, scalar &area) const
Return the intersection centre and area.
Definition: triangle2D.C:119
static scalar relTol
Relative tolerance.
Definition: triangle2D.H:70
bool contains(const triangle2D &tri) const
Return true if tri is within this triangle.
Definition: triangle2DI.H:205
const volScalarField & p0
Definition: EEqn.H:36
OBJstream os(runTime.globalPath()/outputName)
dimensionedScalar j1(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & b
Definition: createFields.H:27
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
A non-counting (dummy) refCount.
Definition: refCount.H:59