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 -------------------------------------------------------------------------------
10 License
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 
51 inline Foam::label Foam::triangle2D::nClosePoints
52 (
53  const triangle2D& triA,
54  const triangle2D& triB
55 )
56 {
57  label nClose = 0;
58 
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 
78 inline 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,
162  vector2D& intersection
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 
193 inline Foam::scalar Foam::triangle2D::area() const noexcept
194 {
195  return area_;
196 }
197 
198 
199 inline Foam::label Foam::triangle2D::order() const
200 {
201  return mag(area_) < SMALL ? 0 : sign(area_);
202 }
203 
204 
205 inline 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 
214 inline 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 
225 inline 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 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Foam::triangle2D::order
label order() const
Definition: triangle2DI.H:192
Foam::triangle2D::lineIntersects
static bool lineIntersects(const vector2D &a, const vector2D &b, const vector2D &c, const vector2D &d)
Definition: triangle2DI.H:163
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::triangle2D::area
scalar area() const noexcept
Definition: triangle2DI.H:186
Foam::Vector2D::isClose
bool isClose(const Vector2D< Cmpt > &b, const scalar tol=1e-10) const
Definition: Vector2DI.H:128
s
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))
Definition: gmvOutputSpray.H:25
Foam::triangle2D::isSame
bool isSame(const triangle2D &triB) const
Definition: triangle2DI.H:207
Foam::triangle2D::centre
vector2D centre() const
Definition: triangle2DI.H:85
Foam::triangle2D::nClosePoints
static label nClosePoints(const triangle2D &triA, const triangle2D &triB)
Definition: triangle2DI.H:45
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::triangle2D::writeOBJ
static void writeOBJ(Ostream &os, const triangle2D &tri, label offset)
Definition: triangle2DI.H:26
Foam::Vector2D< scalar >
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:159
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
OFstream.H
Foam::Vector2D::y
const Cmpt & y() const
Definition: Vector2DI.H:66
Foam::constant::physicoChemical::b
const dimensionedScalar b
Definition: createFields.H:27
Foam::Vector2D::x
const Cmpt & x() const
Definition: Vector2DI.H:59
Foam::triangle2D::lineIntersectionPoint
static bool lineIntersectionPoint(const vector2D &a, const vector2D &b, const vector2D &c, const vector2D &d, vector2D &intersection)
Definition: triangle2DI.H:150
Foam::triangle2D
2-D triangle and queries
Definition: triangle2D.H:48
Foam::triangle2D::lineSegmentIntersectionPoint
static bool lineSegmentIntersectionPoint(const vector2D &lp1, const vector2D &lp2, const vector2D &sp1, const vector2D &sp2, vector2D &intersection)
Definition: triangle2DI.H:94
os
OBJstream os(runTime.globalPath()/outputName)
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:49
Foam::stringOps::match
bool match(const UList< wordRe > &patterns, const std::string &text)
Definition: stringOps.H:72
DebugInfo
#define DebugInfo
Definition: messageStream.H:436
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::triangle2D::pointInside
bool pointInside(const vector2D &p) const
Definition: triangle2DI.H:218
Foam::vector2D
Vector2D< scalar > vector2D
A 2D vector of scalars obtained from the generic Vector2D.
Definition: vector2D.H:47
Foam::triangle2D::contains
bool contains(const triangle2D &tri) const
Definition: triangle2DI.H:198
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::FixedList< bool, 3 >
Foam::intersection
Foam::intersection.
Definition: intersection.H:48
Foam::constant::universal::c
const dimensionedScalar c
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:52
Foam::Vector2D::perp
scalar perp(const Vector2D< Cmpt > &b) const
Definition: Vector2DI.H:120