triangle2DI.H
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 "OFstream.H"
29
30// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
31
33(
34 Ostream& os,
35 const triangle2D& tri,
36 label offset
37)
38{
39 os << "v " << tri[0].x() << " " << tri[0].y() << " 0" << nl
40 << "v " << tri[1].x() << " " << tri[1].y() << " 0" << nl
41 << "v " << tri[2].x() << " " << tri[2].y() << " 0" << nl
42 << "l"
43 << " " << 1 + offset
44 << " " << 2 + offset
45 << " " << 3 + offset
46 << " " << 1 + offset
47 << endl;
48}
49
50
52(
53 const triangle2D& triA,
54 const triangle2D& triB
55)
56{
57 label nClose = 0;
58
59 FixedList<bool, 3> match(true);
60
61 for (auto& a : triA)
62 {
63 forAll(triB, tb)
64 {
65 if (match[tb] && a.isClose(triB[tb]))
66 {
67 match[tb] = false;
68 ++nClose;
69 break;
70 }
71 }
72 }
73
74 return nClose;
75}
76
77
78inline Foam::scalar Foam::triangle2D::area
79(
80 const vector2D& a,
81 const vector2D& b,
82 const vector2D& c
83)
84{
85 const vector2D e1(a - c);
86 const vector2D e2(b - c);
87
88 return 0.5*e1.perp(e2);
89}
90
91
93{
94 const triangle2D& tri = *this;
95
96 return (1/3.0)*(tri[0] + tri[1] + tri[2]);
97}
98
99
101(
102 const vector2D& lp1,
103 const vector2D& lp2,
104 const vector2D& sp1,
105 const vector2D& sp2,
107)
108{
109 const vector2D r(lp2 - lp1);
110 const vector2D s(sp2 - sp1);
111
112 const scalar rcs = r.perp(s);
113
114 if (mag(rcs) > ROOTVSMALL)
115 {
116 const scalar u = (sp1 - lp1).perp(r)/rcs;
117
118 if (u >= -relTol && u <= 1+relTol)
119 {
120 intersection = sp1 + u*s;
121 return true;
122 }
123 }
124
125/*
126 // Collapsed to line
127 if (mag(triangle2D(lp1, lp2, sp1).area()) < absTol)
128 {
129 intersection = sp1;
130 return true;
131 }
132
133 if (mag(triangle2D(lp1, lp2, sp2).area()) < absTol)
134 {
135 intersection = sp2;
136 return true;
137 }
138*/
139
140 if (debug)
141 {
142 OFstream os("bad-intersection.obj");
143 os << "v " << lp1.x() << " " << lp1.y() << " 0" << nl
144 << "v " << lp2.x() << " " << lp2.y() << " 0" << nl
145 << "v " << sp1.x() << " " << sp1.y() << " 0" << nl
146 << "v " << sp2.x() << " " << sp2.y() << " 0" << nl
147 << "l 1 2"
148 << "l 3 4"
149 << endl;
150 }
151
152 return false;
153}
154
155
157(
158 const vector2D& a,
159 const vector2D& b,
160 const vector2D& c,
161 const vector2D& d,
163)
164{
165 return lineSegmentIntersectionPoint(c, d, a, b, intersection);
166}
167
168
170(
171 const vector2D& a,
172 const vector2D& b,
173 const vector2D& c,
174 const vector2D& d
175)
176{
177 if
178 (
179 (triangle2D(a, c, d).order() == triangle2D(b, c, d).order())
180 || (triangle2D(a, b, c).order() == triangle2D(a, b, d).order())
181 )
182 {
183 return false;
184 }
185
186
187 DebugInfo<< "line " << a << b << " intersects " << c << d << endl;
188
189 return true;
190}
191
192
193inline Foam::scalar Foam::triangle2D::area() const noexcept
194{
195 return area_;
196}
197
198
199inline Foam::label Foam::triangle2D::order() const
200{
201 return mag(area_) < SMALL ? 0 : sign(area_);
202}
203
204
205inline bool Foam::triangle2D::contains(const triangle2D& tri) const
206{
207 return
208 pointInside(tri[0])
209 && pointInside(tri[1])
210 && pointInside(tri[2]);
211}
212
213
214inline bool Foam::triangle2D::isSame(const triangle2D& triB) const
215{
216 const triangle2D& triA = *this;
217
218 return
219 triA[0].isClose(triB[0])
220 && triA[1].isClose(triB[1])
221 && triA[2].isClose(triB[2]);
222}
223
224
225inline bool Foam::triangle2D::pointInside(const vector2D& p) const
226{
227 const triangle2D& triA = *this;
228
229 for (int i = 0; i < 3; ++i)
230 {
231 if ((triA[(i + 1) % 3] - triA[i]).perp(p - triA[i]) < 0)
232 {
233 return false;
234 }
235 }
236
237 return true;
238}
239
240
241// ************************************************************************* //
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
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
const Cmpt & y() const
Access to the vector y component.
Definition: Vector2DI.H:73
scalar perp(const Vector2D< Cmpt > &b) const
Perp dot product (dot product with perpendicular vector)
Definition: Vector2DI.H:136
const Cmpt & x() const
Access to the vector x component.
Definition: Vector2DI.H:66
Foam::intersection.
Definition: intersection.H:53
2-D triangle and queries
Definition: triangle2D.H:55
static void writeOBJ(Ostream &os, const triangle2D &tri, label offset)
Write the triangle in OBJ format.
Definition: triangle2DI.H:33
bool isSame(const triangle2D &triB) const
Return true if triB is the same as this triangle.
Definition: triangle2DI.H:214
bool pointInside(const vector2D &p) const
Return true if t point p is inside this triangle.
Definition: triangle2DI.H:225
static bool lineIntersects(const vector2D &a, const vector2D &b, const vector2D &c, const vector2D &d)
Return true if lines ab and cd intersect.
Definition: triangle2DI.H:170
static label nClosePoints(const triangle2D &triA, const triangle2D &triB)
Return the number of similar points between two triangles.
Definition: triangle2DI.H:52
label order() const
Definition: triangle2DI.H:199
scalar area() const noexcept
Return the triangle area.
Definition: triangle2DI.H:193
static bool lineIntersectionPoint(const vector2D &a, const vector2D &b, const vector2D &c, const vector2D &d, vector2D &intersection)
Definition: triangle2DI.H:157
static bool lineSegmentIntersectionPoint(const vector2D &lp1, const vector2D &lp2, const vector2D &sp1, const vector2D &sp2, vector2D &intersection)
Definition: triangle2DI.H:101
vector2D centre() const
Return the triangle centre.
Definition: triangle2DI.H:92
bool contains(const triangle2D &tri) const
Return true if tri is within this triangle.
Definition: triangle2DI.H:205
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define DebugInfo
Report an information message using Foam::Info.
dimensionedScalar sign(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const direction noexcept
Definition: Scalar.H:223
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & b
Definition: createFields.H:27
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333