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-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 * * * * * * * * * * * * * * //
33 
35 :
36  Pair<point>(pts.a(), pts.b())
37 {}
38 
39 
41 (
42  const UList<point>& points,
43  const FixedList<label, 2>& indices
44 )
45 :
46  Pair<point>(points[indices.get<0>()], points[indices.get<1>()])
47 {}
48 
49 
50 template<class Point, class PointRef>
52 (
53  const Point& from,
54  const Point& to
55 )
56 :
57  a_(from),
58  b_(to)
59 {}
60 
61 
62 template<class Point, class PointRef>
64 (
65  const UList<Point>& points,
66  const FixedList<label, 2>& indices
67 )
68 :
69  a_(points[indices.template get<0>()]),
70  b_(points[indices.template get<1>()])
71 {}
72 
73 
74 template<class Point, class PointRef>
76 {
77  is >> *this;
78 }
79 
80 
81 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82 
84 {
85  return linePointRef(a(), b());
86 }
87 
88 
89 template<class Point, class PointRef>
91 {
92  return 0.5*(a_ + b_);
93 }
94 
95 
97 {
98  return 0.5*(a() + b());
99 }
100 
101 
102 template<class Point, class PointRef>
103 inline Foam::scalar Foam::line<Point, PointRef>::mag() const
104 {
105  return ::Foam::mag(b() - a());
106 }
107 
108 
109 inline Foam::scalar Foam::linePoints::mag() const
110 {
111  return ::Foam::mag(b() - a());
112 }
113 
114 
115 template<class Point, class PointRef>
117 {
118  return (b() - a());
119 }
120 
121 
123 {
124  return (b() - a());
125 }
126 
127 
128 template<class Point, class PointRef>
130 {
131  const Point v = (b_ - a_);
132  const scalar s(::Foam::mag(v));
133 
134  return s < ROOTVSMALL ? Zero : v/s;
135 }
136 
137 
139 {
140  return normalised(b() - a());
141 }
142 
143 
144 template<class Point, class PointRef>
146 (
147  const Point& p0,
148  const Point& p1
149 )
150 {
151  return Pair<Point>(min(p0, p1), max(p0, p1));
152 }
153 
154 
155 template<class Point, class PointRef>
157 {
158  return line<Point, PointRef>::box(a_, b_);
159 }
160 
161 
163 {
164  return linePointRef::box(a(), b());
165 }
166 
167 
168 template<class Point, class PointRef>
170 (
171  const Point& p
172 ) const
173 {
174  Point v = vec();
175 
176  Point w(p - a_);
177 
178  const scalar c1 = v & w;
179 
180  if (c1 <= 0)
181  {
182  return PointHit<Point>(false, a_, Foam::mag(p - a_), true);
183  }
184 
185  const scalar c2 = v & v;
186 
187  if (c2 <= c1)
188  {
189  return PointHit<Point>(false, b_, Foam::mag(p - b_), true);
190  }
191 
192  const scalar b = c1/c2;
193 
194  Point pb(a_ + b*v);
196  return PointHit<Point>(true, pb, Foam::mag(p - pb), false);
197 }
198 
199 
200 template<class Point, class PointRef>
202 (
204  Point& thisPt,
205  Point& edgePt
206 ) const
207 {
208  // From Mathworld Line-Line distance/(Gellert et al. 1989, p. 538).
209  Point a(end() - start());
210  Point b(edge.end() - edge.start());
211  Point c(edge.start() - start());
212 
213  Point crossab = a ^ b;
214  const scalar magCrossSqr = Foam::magSqr(crossab);
215 
216  if (magCrossSqr > VSMALL)
217  {
218  scalar s = ((c ^ b) & crossab)/magCrossSqr;
219  scalar t = ((c ^ a) & crossab)/magCrossSqr;
220 
221  // Check for end points outside of range 0..1
222  if (s >= 0 && s <= 1 && t >= 0 && t <= 1)
223  {
224  // Both inside range 0..1
225  thisPt = start() + a*s;
226  edgePt = edge.start() + b*t;
227  }
228  else
229  {
230  // Do brute force. Distance of everything to everything.
231  // Can quite possibly be improved!
232 
233  // From edge endpoints to *this
234  PointHit<Point> this0(nearestDist(edge.start()));
235  PointHit<Point> this1(nearestDist(edge.end()));
236  scalar thisDist = min(this0.distance(), this1.distance());
237 
238  // From *this to edge
239  PointHit<Point> edge0(edge.nearestDist(start()));
240  PointHit<Point> edge1(edge.nearestDist(end()));
241  scalar edgeDist = min(edge0.distance(), edge1.distance());
242 
243  if (thisDist < edgeDist)
244  {
245  if (this0.distance() < this1.distance())
246  {
247  thisPt = this0.point();
248  edgePt = edge.start();
249  }
250  else
251  {
252  thisPt = this1.point();
253  edgePt = edge.end();
254  }
255  }
256  else
257  {
258  if (edge0.distance() < edge1.distance())
259  {
260  thisPt = start();
261  edgePt = edge0.point();
262  }
263  else
264  {
265  thisPt = end();
266  edgePt = edge1.point();
267  }
268  }
269  }
270  }
271  else
272  {
273  // Parallel lines. Find overlap of both lines by projecting onto
274  // direction vector (now equal for both lines).
275 
276  scalar edge0 = edge.start() & a;
277  scalar edge1 = edge.end() & a;
278  bool edgeOrder = edge0 < edge1;
279 
280  scalar minEdge = (edgeOrder ? edge0 : edge1);
281  scalar maxEdge = (edgeOrder ? edge1 : edge0);
282  const Point& minEdgePt = (edgeOrder ? edge.start() : edge.end());
283  const Point& maxEdgePt = (edgeOrder ? edge.end() : edge.start());
284 
285  scalar this0 = start() & a;
286  scalar this1 = end() & a;
287  bool thisOrder = this0 < this1;
288 
289  scalar minThis = min(this0, this1);
290  scalar maxThis = max(this1, this0);
291  const Point& minThisPt = (thisOrder ? start() : end());
292  const Point& maxThisPt = (thisOrder ? end() : start());
293 
294  if (maxEdge < minThis)
295  {
296  // edge completely below *this
297  edgePt = maxEdgePt;
298  thisPt = minThisPt;
299  }
300  else if (maxEdge < maxThis)
301  {
302  // maxEdge inside interval of *this
303  edgePt = maxEdgePt;
304  thisPt = nearestDist(edgePt).point();
305  }
306  else
307  {
308  // maxEdge outside. Check if minEdge inside.
309  if (minEdge < minThis)
310  {
311  // Edge completely envelops this. Take any this point and
312  // determine nearest on edge.
313  thisPt = minThisPt;
314  edgePt = edge.nearestDist(thisPt).point();
315  }
316  else if (minEdge < maxThis)
317  {
318  // minEdge inside this interval.
319  edgePt = minEdgePt;
320  thisPt = nearestDist(edgePt).point();
321  }
322  else
323  {
324  // minEdge outside this interval
325  edgePt = minEdgePt;
326  thisPt = maxThisPt;
327  }
328  }
329  }
330 
331  return Foam::mag(thisPt - edgePt);
332 }
333 
334 
335 // * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * * //
336 
337 template<class Point, class PointRef>
338 inline Foam::Istream& Foam::operator>>
339 (
340  Istream& is,
341  line<Point, PointRef>& l
342 )
343 {
344  is.readBegin("line");
345  is >> l.a_ >> l.b_;
346  is.readEnd("line");
347 
348  is.check(FUNCTION_NAME);
349  return is;
350 }
351 
352 
353 template<class Point, class PointRef>
354 inline Foam::Ostream& Foam::operator<<
355 (
356  Ostream& os,
357  const line<Point, PointRef>& l
358 )
359 {
361  << l.a_ << token::SPACE
362  << l.b_
363  << token::END_LIST;
364  return os;
365 }
366 
367 
368 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Point unitVec() const
Return the unit vector (start-to-end)
Definition: lineI.H:122
line(const Point &from, const Point &to)
Construct from two points.
Definition: lineI.H:45
A line primitive.
Definition: line.H:52
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
linePointRef ln() const
Return as line reference.
Definition: lineI.H:76
scalar mag() const
The magnitude (length) of the line.
Definition: lineI.H:96
bool readBegin(const char *funcName)
Begin read of data chunk, starts with &#39;(&#39;.
Definition: Istream.C:104
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Pair< point > box() const
The enclosing (bounding) box for the line.
Definition: lineI.H:155
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
scalar mag() const
The magnitude (length) of the line.
Definition: lineI.H:102
Begin list [isseparator].
Definition: token.H:158
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:673
PointRef end() const noexcept
The end (second) point.
Definition: line.H:247
Pair< Point > box() const
The enclosing (bounding) box for the line.
Definition: lineI.H:149
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
vector vec() const
Return start-to-end vector.
Definition: lineI.H:115
label start() const noexcept
The start (first) vertex label.
Definition: edge.H:147
vector unitVec() const
Return the unit vector (start-to-end)
Definition: lineI.H:131
linePoints()=default
Default construct.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: instant.H:46
line< point, const point & > linePointRef
A line using referred points.
Definition: line.H:66
const pointField & points
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:59
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Space [isspace].
Definition: token.H:128
End list [isseparator].
Definition: token.H:159
Describes the interaction of a object and a (templated) point. It carries the info of a successful hi...
Definition: pointHit.H:43
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:99
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:194
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
point centre() const
Return centre (centroid)
Definition: lineI.H:89
CGAL::Point_3< K > Point
Point centre() const
Return centre (centroid)
Definition: lineI.H:83
Point vec() const
Return start-to-end vector.
Definition: lineI.H:109
const dimensionedScalar c
Speed of light in a vacuum.
scalar distance() const noexcept
Return distance to hit.
Definition: pointHit.H:169
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
const point_type & point() const noexcept
Return the point, no checks.
Definition: pointHit.H:161
label end() const noexcept
The end (last/second) vertex label.
Definition: edge.H:152
volScalarField & p
PointRef start() const noexcept
The start (first) point.
Definition: line.H:242
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))
PointHit< Point > nearestDist(const Point &p) const
Return nearest distance to line from a given point.
Definition: lineI.H:163
const volScalarField & p0
Definition: EEqn.H:36
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
const pointField & pts
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133