Loading...
Searching...
No Matches
lineI.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) 2011 OpenFOAM Foundation
9 Copyright (C) 2018-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "zero.H"
30#include "IOstreams.H"
31
32// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35:
36 Pair<point>(pts.a(), pts.b())
37{}
38
39
41(
42 const UList<point>& points,
43 const FixedList<label, 2>& indices
45:
46 Pair<point>(points[indices.get<0>()], points[indices.get<1>()])
47{}
48
49
50template<class Point, class PointRef>
52(
53 const Point& from,
54 const Point& to
55)
57 a_(from),
58 b_(to)
59{}
60
61
62template<class Point, class PointRef>
64(
65 const UList<Point>& points,
66 const FixedList<label, 2>& indices
67)
69 a_(points[indices.template get<0>()]),
70 b_(points[indices.template get<1>()])
71{}
72
73
74template<class Point, class PointRef>
77 is >> *this;
78}
79
80
81// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82
84{
85 return linePointRef(a(), b());
86}
87
88
89template<class Point, class PointRef>
91{
92 return 0.5*(a_ + b_);
93}
94
95
97{
98 return 0.5*(a() + b());
99}
100
101
102template<class Point, class PointRef>
103inline Foam::scalar Foam::line<Point, PointRef>::mag() const
104{
105 return ::Foam::mag(b() - a());
106}
107
108
109inline Foam::scalar Foam::linePoints::mag() const
110{
111 return ::Foam::mag(b() - a());
112}
113
114
115template<class Point, class PointRef>
116inline Foam::scalar Foam::line<Point, PointRef>::magSqr() const
117{
118 return ::Foam::magSqr(b() - a());
119}
120
121
122inline Foam::scalar Foam::linePoints::magSqr() const
123{
124 return ::Foam::magSqr(b() - a());
125}
126
127
128template<class Point, class PointRef>
130{
131 return (b() - a());
132}
133
134
136{
137 return (b() - a());
138}
139
140
141template<class Point, class PointRef>
143{
144 const Point v = (b_ - a_);
145
146 #ifdef __clang__
147 volatile // Use volatile to avoid aggressive branch optimization
148 #endif
149 const scalar s(::Foam::mag(v));
150
151 return s < ROOTVSMALL ? Zero : v/s;
152}
153
154
156{
157 return normalised(b() - a());
158}
159
160
161template<class Point, class PointRef>
163(
164 const Point& p0,
165 const Point& p1
167{
168 return Pair<Point>(min(p0, p1), max(p0, p1));
169}
170
171
172template<class Point, class PointRef>
174{
175 return line<Point, PointRef>::box(a_, b_);
176}
177
178
180{
181 return linePointRef::box(a(), b());
182}
183
184
185template<class Point, class PointRef>
187(
188 const Point& p
189) const
190{
191 Point v = vec();
192
193 Point w(p - a_);
194
195 const scalar c1 = v & w;
196
197 if (c1 <= 0)
198 {
199 return PointHit<Point>(false, a_, Foam::mag(p - a_), true);
200 }
201
202 const scalar c2 = v & v;
203
204 if (c2 <= c1)
205 {
206 return PointHit<Point>(false, b_, Foam::mag(p - b_), true);
207 }
208
209 const scalar b = c1/c2;
210
211 Point pb(a_ + b*v);
212
213 return PointHit<Point>(true, pb, Foam::mag(p - pb), false);
214}
215
216
217template<class Point, class PointRef>
219(
221 Point& thisPt,
222 Point& edgePt
223) const
224{
225 // From Mathworld Line-Line distance/(Gellert et al. 1989, p. 538).
226 Point a(end() - start());
227 Point b(edge.end() - edge.start());
228 Point c(edge.start() - start());
229
230 Point crossab = a ^ b;
231 const scalar magCrossSqr = Foam::magSqr(crossab);
232
233 if (magCrossSqr > VSMALL)
234 {
235 scalar s = ((c ^ b) & crossab)/magCrossSqr;
236 scalar t = ((c ^ a) & crossab)/magCrossSqr;
237
238 // Check for end points outside of range 0..1
239 if (s >= 0 && s <= 1 && t >= 0 && t <= 1)
240 {
241 // Both inside range 0..1
242 thisPt = start() + a*s;
243 edgePt = edge.start() + b*t;
244 }
245 else
246 {
247 // Do brute force. Distance of everything to everything.
248 // Can quite possibly be improved!
249
250 // From edge endpoints to *this
251 PointHit<Point> this0(nearestDist(edge.start()));
252 PointHit<Point> this1(nearestDist(edge.end()));
253 scalar thisDist = min(this0.distance(), this1.distance());
254
255 // From *this to edge
256 PointHit<Point> edge0(edge.nearestDist(start()));
257 PointHit<Point> edge1(edge.nearestDist(end()));
258 scalar edgeDist = min(edge0.distance(), edge1.distance());
259
260 if (thisDist < edgeDist)
261 {
262 if (this0.distance() < this1.distance())
263 {
264 thisPt = this0.point();
265 edgePt = edge.start();
266 }
267 else
268 {
269 thisPt = this1.point();
270 edgePt = edge.end();
271 }
272 }
273 else
274 {
275 if (edge0.distance() < edge1.distance())
276 {
277 thisPt = start();
278 edgePt = edge0.point();
279 }
280 else
281 {
282 thisPt = end();
283 edgePt = edge1.point();
284 }
285 }
286 }
287 }
288 else
289 {
290 // Parallel lines. Find overlap of both lines by projecting onto
291 // direction vector (now equal for both lines).
292
293 scalar edge0 = edge.start() & a;
294 scalar edge1 = edge.end() & a;
295 bool edgeOrder = edge0 < edge1;
296
297 scalar minEdge = (edgeOrder ? edge0 : edge1);
298 scalar maxEdge = (edgeOrder ? edge1 : edge0);
299 const Point& minEdgePt = (edgeOrder ? edge.start() : edge.end());
300 const Point& maxEdgePt = (edgeOrder ? edge.end() : edge.start());
301
302 scalar this0 = start() & a;
303 scalar this1 = end() & a;
304 bool thisOrder = this0 < this1;
305
306 scalar minThis = min(this0, this1);
307 scalar maxThis = max(this1, this0);
308 const Point& minThisPt = (thisOrder ? start() : end());
309 const Point& maxThisPt = (thisOrder ? end() : start());
310
311 if (maxEdge < minThis)
312 {
313 // edge completely below *this
314 edgePt = maxEdgePt;
315 thisPt = minThisPt;
316 }
317 else if (maxEdge < maxThis)
318 {
319 // maxEdge inside interval of *this
320 edgePt = maxEdgePt;
321 thisPt = nearestDist(edgePt).point();
322 }
323 else
324 {
325 // maxEdge outside. Check if minEdge inside.
326 if (minEdge < minThis)
327 {
328 // Edge completely envelops this. Take any this point and
329 // determine nearest on edge.
330 thisPt = minThisPt;
331 edgePt = edge.nearestDist(thisPt).point();
332 }
333 else if (minEdge < maxThis)
334 {
335 // minEdge inside this interval.
336 edgePt = minEdgePt;
337 thisPt = nearestDist(edgePt).point();
338 }
339 else
340 {
341 // minEdge outside this interval
342 edgePt = minEdgePt;
343 thisPt = maxThisPt;
344 }
345 }
346 }
347
348 return Foam::mag(thisPt - edgePt);
349}
350
351
352// * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * * //
353
354template<class Point, class PointRef>
355inline Foam::Istream& Foam::operator>>
356(
357 Istream& is,
359)
360{
361 is.readBegin("line");
362 is >> l.a_ >> l.b_;
363 is.readEnd("line");
364
365 is.check(FUNCTION_NAME);
366 return is;
367}
368
369
370template<class Point, class PointRef>
371inline Foam::Ostream& Foam::operator<<
372(
373 Ostream& os,
374 const line<Point, PointRef>& l
375)
376{
378 << l.a_ << token::SPACE
379 << l.b_
381 return os;
382}
383
384
385// ************************************************************************* //
CGAL::Point_3< K > Point
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
T & get() noexcept
Element access using compile-time indexing.
Definition FixedListI.H:152
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
An ordered pair of two objects of type <T> with first() and second() elements.
Definition Pair.H:66
Describes the interaction of a object and a (templated) point. It carries the info of a successful hi...
Definition pointHit.H:60
scalar distance() const noexcept
Return distance to hit.
Definition pointHit.H:169
const point_type & point() const noexcept
Return the point, no checks.
Definition pointHit.H:161
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
label & end() noexcept
The end (second/last) vertex label.
Definition edge.H:165
label start() const noexcept
The start (first) vertex label.
Definition edge.H:155
vector vec() const
Return start-to-end vector.
Definition lineI.H:128
linePointRef ln() const
Return as line reference.
Definition lineI.H:76
const point & b() const noexcept
The second vertex.
Definition line.H:122
Pair< point > box() const
The enclosing (bounding) box for the line.
Definition lineI.H:172
scalar mag() const
The magnitude (length) of the line.
Definition lineI.H:102
const point & a() const noexcept
The first vertex.
Definition line.H:117
scalar magSqr() const
The magnitude squared (length squared) of the line.
Definition lineI.H:115
vector unitVec() const
Return the unit vector (start-to-end).
Definition lineI.H:148
linePoints()=default
Default construct.
point centre() const
Return centre (centroid).
Definition lineI.H:89
A line primitive.
Definition line.H:180
PointHit< Point > nearestDist(const Point &p) const
Return nearest distance to line from a given point.
Definition lineI.H:180
PointRef end() const noexcept
The end (second) point.
Definition line.H:252
line(const Point &from, const Point &to)
Construct from two points.
Definition lineI.H:45
PointRef a() const noexcept
The first point.
Definition line.H:227
Point centre() const
Return centre (centroid).
Definition lineI.H:83
scalar mag() const
The magnitude (length) of the line.
Definition lineI.H:96
Point unitVec() const
Return the unit vector (start-to-end).
Definition lineI.H:135
Pair< Point > box() const
The enclosing (bounding) box for the line.
Definition lineI.H:166
scalar magSqr() const
The magnitude squared (length squared) of the line.
Definition lineI.H:109
PointRef b() const noexcept
The second point.
Definition line.H:232
PointRef start() const noexcept
The start (first) point.
Definition line.H:247
Point vec() const
Return start-to-end vector.
Definition lineI.H:122
@ BEGIN_LIST
Begin list [isseparator].
Definition token.H:174
@ END_LIST
End list [isseparator].
Definition token.H:175
@ SPACE
Space [isspace].
Definition token.H:144
volScalarField & p
const volScalarField & p0
Definition EEqn.H:36
OBJstream os(runTime.globalPath()/outputName)
const pointField & points
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 FUNCTION_NAME
const char * end
Definition SVGTools.H:223
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
line< point, const point & > linePointRef
A line using referred points.
Definition line.H:66
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Vector< scalar > vector
Definition vector.H:57
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
const pointField & pts
volScalarField & b