triangle2D.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
26Class
27 Foam::triangle2D
28
29Description
30 2-D triangle and queries
31
32SourceFiles
33 triangle2DI.H
34 triangle2D.C
35
36\*---------------------------------------------------------------------------*/
37
38#ifndef triangle2D_H
39#define triangle2D_H
40
41#include "vector.H"
42#include "vector2D.H"
43
44// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45
46namespace Foam
47{
48
49/*---------------------------------------------------------------------------*\
50 Class triangle2D Declaration
51\*---------------------------------------------------------------------------*/
53class triangle2D
54:
55 public FixedList<vector2D, 3>
56{
57 // Private Data
58
59 //- Area
60 scalar area_;
61
62 //- Helper for intersections
63 static FixedList<vector2D, 7> work_;
64
65
66public:
68 static int debug;
69
70 //- Relative tolerance
71 static scalar relTol;
72
73 //- Absolute tolerance
74 static scalar absTol;
75
76 //- Construct from 3 2-D points
78 (
79 const vector2D& a,
80 const vector2D& b,
81 const vector2D& c,
82 const bool orient = false
83 );
84
85 //- Construct from 3 3-D points and axes
87 (
88 const vector& a3d,
89 const vector& b3d,
90 const vector& c3d,
91 const vector& axis1,
92 const vector& axis2,
93 const bool orient = false
94 );
95
96
97 // Member Functions
98
99 //- Returns:
100 //- 1 if points are ordered in anti-clockwise direction
101 //- -1 if points are ordered in clockwise direction
102 //- 0 if the triangle has collapsed to a line
103 inline label order() const;
104
105 //- Write the triangle in OBJ format
106 inline static void writeOBJ
107 (
108 Ostream& os,
109 const triangle2D& tri,
110 label offset
111 );
112
113 //- Return the number of similar points between two triangles
114 inline static label nClosePoints
115 (
116 const triangle2D& triA,
117 const triangle2D& triB
118 );
119
120 //- Return the signed area
121 inline static scalar area
122 (
123 const vector2D& a,
124 const vector2D& b,
125 const vector2D& c
126 );
127
128 //- Set the intersection between a line and segment
129 //- Return true if lines intersect
130 inline static bool lineSegmentIntersectionPoint
131 (
132 const vector2D& lp1,
133 const vector2D& lp2,
134 const vector2D& sp1,
135 const vector2D& sp2,
137 );
138
139 //- Set the intersection between two lines
140 //- Return true if lines intersect
141 inline static bool lineIntersectionPoint
142 (
143 const vector2D& a,
144 const vector2D& b,
145 const vector2D& c,
146 const vector2D& d,
148 );
149
150 //- Return true if lines ab and cd intersect
151 inline static bool lineIntersects
152 (
153 const vector2D& a,
154 const vector2D& b,
155 const vector2D& c,
156 const vector2D& d
157 );
158
159 //- Snap [this] triangle's points to those of triB if they are within
160 //- absTol
161 // Returns the number of snapped points
162 label snapClosePoints(const triangle2D& triB);
163
164 //- Return the intersection centre and area
165 void interArea
166 (
167 const triangle2D& triB,
169 scalar& area
170 ) const;
171
172 //- Return the intersection area
173 scalar interArea(const triangle2D& triB) const;
174
175 //- Return true if triB overlaps
176 bool overlaps(const triangle2D& triB) const;
177
178 //- Return the triangle area
179 inline scalar area() const noexcept;
180
181 //- Return the triangle centre
182 inline vector2D centre() const;
183
184 //- Return true if tri is within this triangle
185 inline bool contains(const triangle2D& tri) const;
186
187 //- Return true if triB is the same as this triangle
188 inline bool isSame(const triangle2D& triB) const;
189
190 //- Return true if t point p is inside this triangle
191 inline bool pointInside(const vector2D& p) const;
192};
193
194
195// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196
197} // End namespace Foam
198
199// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200
201#include "triangle2DI.H"
202
203// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204
205#endif
206
207// ************************************************************************* //
A 1D vector of objects of type <T> with a fixed length <N>.
Definition: FixedList.H:81
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Foam::intersection.
Definition: intersection.H:53
2-D triangle and queries
Definition: triangle2D.H:55
static scalar absTol
Absolute tolerance.
Definition: triangle2D.H:73
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
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 overlaps(const triangle2D &triB) const
Return true if triB overlaps.
Definition: triangle2D.C:276
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)
Namespace for OpenFOAM.
const direction noexcept
Definition: Scalar.H:223
volScalarField & b
Definition: createFields.H:27