edgeI.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) 2017-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 "IOstreams.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 inline int Foam::edge::compare(const edge& a, const edge& b)
34 {
35  return labelPair::compare(a, b);
36 }
37 
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
41 inline Foam::edge::edge()
42 :
43  labelPair(-1, -1)
44 {}
45 
46 
47 inline Foam::edge::edge(const label from, const label to)
48 :
49  labelPair(from, to)
50 {}
51 
52 
53 inline Foam::edge::edge(const labelPair& pair)
54 :
55  labelPair(pair.first(), pair.second())
56 {}
57 
58 
59 inline Foam::edge::edge(const FixedList<label, 2>& list)
60 :
61  labelPair(list.get<0>(), list.get<1>())
62 {}
63 
64 
65 inline Foam::edge::edge(const label from, const label to, const bool doSort)
66 :
67  labelPair(from, to, doSort)
68 {}
69 
70 
71 inline Foam::edge::edge(const FixedList<label, 2>& list, const bool doSort)
72 :
73  labelPair(list, doSort)
74 {}
75 
76 
77 inline Foam::edge::edge
78 (
79  const labelUList& list,
80  const FixedList<label, 2>& indices
81 )
82 :
83  labelPair(list[indices.get<0>()], list[indices.get<1>()])
84 {}
85 
86 
87 inline Foam::edge::edge(Istream& is)
88 :
89  labelPair(is)
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
95 inline Foam::label Foam::edge::minVertex() const noexcept
96 {
97  return (first() < second() ? first() : second());
98 }
99 
101 inline Foam::label Foam::edge::maxVertex() const noexcept
102 {
103  return (second() < first() ? first() : second());
104 }
105 
107 inline bool Foam::edge::good() const noexcept
108 {
109  return (first() != second() && first() >= 0 && second() >= 0);
110 }
111 
112 
113 inline bool Foam::edge::contains(const label vertex) const noexcept
114 {
115  // -1: always false
116  return
117  (
118  vertex >= 0
119  && (vertex == first() || vertex == second())
120  );
121 }
122 
123 
124 inline Foam::label Foam::edge::which(const label vertex) const
125 {
126  // -1: always false
127  if (vertex >= 0)
128  {
129  if (vertex == first())
130  {
131  return 0;
132  }
133  if (vertex == second())
134  {
135  return 1;
136  }
137  }
138 
139  return -1;
140 }
141 
143 inline bool Foam::edge::connected(const edge& other) const
144 {
145  return (other.contains(first()) || other.contains(second()));
146 }
147 
148 
149 inline Foam::label Foam::edge::commonVertex(const edge& other) const
150 {
151  if (other.contains(first()))
152  {
153  return first();
154  }
155  if (other.contains(second()))
156  {
157  return second();
158  }
159 
160  // No shared vertex.
161  return -1;
162 }
163 
164 
165 inline Foam::label Foam::edge::otherVertex(const label vertex) const
166 {
167  if (vertex == first())
168  {
169  return second();
170  }
171  if (vertex == second())
172  {
173  return first();
174  }
175 
176  // The given vertex is not on the edge in the first place.
177  return -1;
178 }
179 
180 
181 inline Foam::label Foam::edge::collapse()
182 {
183  // Cannot resize FixedList, so mark duplicates with '-1'
184  // (the lower vertex is retained)
185  // catch any '-1' (eg, if called multiple times)
186 
187  label n = 2;
188  if (first() == second() || second() < 0)
189  {
190  second() = -1;
191  --n;
192  }
193  if (first() < 0)
194  {
195  --n;
196  }
197 
198  return n;
199 }
200 
202 inline Foam::edge Foam::edge::reverseEdge() const
203 {
204  return Foam::edge(second(), first());
205 }
206 
207 
208 inline void Foam::edge::clear()
209 {
210  first() = -1;
211  second() = -1;
212 }
213 
214 
215 inline Foam::label Foam::edge::count() const
216 {
217  label n = 2;
218  if (first() == second() || second() < 0)
219  {
220  --n;
221  }
222  if (first() < 0)
223  {
224  --n;
225  }
226 
227  return n;
228 }
229 
231 inline bool Foam::edge::empty() const noexcept
232 {
233  return (first() < 0 && second() < 0);
234 }
235 
236 
237 inline bool Foam::edge::insert(const label vertex)
238 {
239  if (vertex < 0)
240  {
241  // Cannot insert invalid vertex labels (use direct assignment for that)
242  return false;
243  }
244 
245  if (first() < 0)
246  {
247  // Store at first, if not duplicate of second
248  if (vertex != second())
249  {
250  first() = vertex;
251  return true;
252  }
253  }
254  else if (second() < 0)
255  {
256  // Store at second, if not duplicate of first
257  if (vertex != first())
258  {
259  second() = vertex;
260  return true;
261  }
262  }
264  return false;
265 }
266 
267 
268 template<class InputIterator>
269 inline Foam::label Foam::edge::insert
270 (
271  InputIterator begIter,
272  InputIterator endIter
273 )
274 {
275  // Available slots.
276  // Don't use count() since it has special treatment for duplicates
277  const int maxChange = ((first() < 0 ? 1 : 0) + (second() < 0 ? 1 : 0));
278 
279  int changed = 0;
280  for (; changed < maxChange && begIter != endIter; ++begIter)
281  {
282  if (insert(*begIter))
283  {
284  ++changed;
285  }
286  }
287 
288  return changed;
289 }
290 
291 
292 inline Foam::label Foam::edge::insert(std::initializer_list<label> list)
293 {
294  return insert(list.begin(), list.end());
295 }
296 
297 
298 template<unsigned N>
299 inline Foam::label Foam::edge::insert(const FixedList<label, N>& list)
300 {
301  return insert(list.begin(), list.end());
302 }
303 
305 inline Foam::label Foam::edge::insert(const labelUList& list)
306 {
307  return insert(list.begin(), list.end());
308 }
309 
310 
311 inline Foam::label Foam::edge::erase(const label vertex)
312 {
313  if (vertex < 0)
314  {
315  // Can never remove invalid point labels!
316  return 0;
317  }
318 
319  label n = 0;
320  if (vertex == first())
321  {
322  first() = -1;
323  ++n;
324  }
325 
326  // Automatically handle duplicates, which should not have been there anyhow
327  if (vertex == second())
328  {
329  second() = -1;
330  ++n;
331  }
333  return n;
334 }
335 
336 
337 template<class InputIterator>
338 inline Foam::label Foam::edge::erase
339 (
340  InputIterator begIter,
341  InputIterator endIter
342 )
343 {
344  // Occupied slots.
345  // Don't use count() since it has special treatment for duplicates
346  const int maxChange = ((first() >= 0 ? 1 : 0) + (second() >= 0 ? 1 : 0));
347 
348  int changed = 0;
349  for (; changed < maxChange && begIter != endIter; ++begIter)
350  {
351  changed += erase(*begIter);
352  }
353 
354  return changed;
355 }
356 
357 
358 inline Foam::label Foam::edge::erase(std::initializer_list<label> list)
359 {
360  return erase(list.begin(), list.end());
361 }
362 
363 
364 template<unsigned N>
365 inline Foam::label Foam::edge::erase(const FixedList<label, N>& list)
366 {
367  return erase(list.begin(), list.end());
368 }
369 
370 
371 inline Foam::label Foam::edge::erase(const labelUList& list)
372 {
373  return erase(list.begin(), list.end());
374 }
375 
376 
377 // Geometric
378 
379 inline Foam::point Foam::edge::centre(const UList<point>& pts) const
380 {
381  #ifdef FULLDEBUG
382  if (first() < 0 || second() < 0)
383  {
385  << "negative point index on edge " << *this
386  << abort(FatalError);
387  }
388  #endif
389 
390  return 0.5*(pts[first()] + pts[second()]);
391 }
392 
393 
394 inline Foam::vector Foam::edge::vec(const UList<point>& pts) const
395 {
396  #ifdef FULLDEBUG
397  if (first() < 0 || second() < 0)
398  {
400  << "negative point index on edge " << *this
401  << abort(FatalError);
402  }
403  #endif
404 
405  return pts[second()] - pts[first()];
406 }
407 
408 
409 inline Foam::vector Foam::edge::unitVec(const UList<point>& pts) const
410 {
411  #ifdef FULLDEBUG
412  if (first() < 0 || second() < 0)
413  {
415  << "negative point index on edge " << *this
416  << abort(FatalError);
417  }
418  #endif
419 
420  const vector v = (pts[second()] - pts[first()]);
421  const scalar s(Foam::mag(v));
422 
423  return s < ROOTVSMALL ? Zero : v/s;
424 }
425 
427 inline Foam::scalar Foam::edge::mag(const UList<point>& pts) const
428 {
429  return pts[first()].dist(pts[second()]);
430 }
431 
432 
433 inline Foam::scalar Foam::edge::magSqr(const UList<point>& pts) const
434 {
435  return pts[first()].distSqr(pts[second()]);
436 }
437 
438 
440 Foam::edge::box(const UList<point>& pts) const
441 {
442  #ifdef FULLDEBUG
443  if (first() < 0 || second() < 0)
444  {
446  << "negative point index on edge " << *this
447  << abort(FatalError);
448  }
449  #endif
450 
451  return linePointRef::box(pts[first()], pts[second()]);
452 }
453 
454 
455 inline Foam::linePointRef Foam::edge::line(const UList<point>& pts) const
456 {
457  #ifdef FULLDEBUG
458  if (first() < 0 || second() < 0)
459  {
461  << "negative point index on edge " << *this
462  << abort(FatalError);
463  }
464  #endif
466  return linePointRef(pts[first()], pts[second()]);
467 }
468 
469 
470 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
471 
472 inline Foam::label& Foam::edge::operator[](const label i)
473 {
474  #ifdef FULLDEBUG
475  if (i < 0 || i > 1)
476  {
478  << "Index " << i << " out of range [0,1]" << abort(FatalError);
479  }
480  #endif
481  return (i ? second() : first());
482 }
483 
484 
485 inline const Foam::label& Foam::edge::operator[](const label i) const
486 {
487  #ifdef FULLDEBUG
488  if (i < 0 || i > 1)
489  {
491  << "Index " << i << " out of range [0,1]" << abort(FatalError);
492  }
493  #endif
494  return (i ? second() : first());
495 }
496 
497 
498 // * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * * //
500 inline bool Foam::operator==(const edge& a, const edge& b)
501 {
502  return edge::compare(a,b) != 0;
503 }
504 
505 
506 inline bool Foam::operator!=(const edge& a, const edge& b)
507 {
508  return edge::compare(a,b) == 0;
509 }
510 
511 
512 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
label count() const
Return the number of unique, valid (non -1) vertex labels.
Definition: edgeI.H:208
A line primitive.
Definition: line.H:52
scalar mag(const UList< point > &pts) const
The length (L2-norm) of the edge vector.
Definition: edgeI.H:420
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
iterator begin() noexcept
Return an iterator to begin traversing the FixedList.
Definition: FixedListI.H:530
srcOptions erase("case")
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
edge reverseEdge() const
Return reverse edge as copy.
Definition: edgeI.H:195
linePointRef line(const UList< point > &pts) const
Return edge line.
Definition: edgeI.H:448
srcOptions insert("case", fileName(rootDirSource/caseDirSource))
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
iterator end() noexcept
Return an iterator to end traversing the FixedList.
Definition: FixedListI.H:554
label b() const noexcept
The second vertex.
Definition: edge.H:132
bool connected(const edge &other) const
True if the edge has at least one vertex in common with other.
Definition: edgeI.H:136
Pair< Point > box() const
The enclosing (bounding) box for the line.
Definition: lineI.H:149
bool empty() const noexcept
Return true if edge has no valid vertex labels.
Definition: edgeI.H:224
vector unitVec(const UList< point > &pts) const
Return the unit vector (from first to second).
Definition: edgeI.H:402
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
label erase(const label vertex)
Remove an existing vertex from the edge and set its location to &#39;-1&#39;. A negative vertex label never r...
Definition: edgeI.H:304
void clear()
&#39;Clears&#39; edge by setting both ends to invalid vertex labels.
Definition: edgeI.H:201
label minVertex() const noexcept
Return the smallest vertex label used by the edge.
Definition: edgeI.H:88
line< point, const point & > linePointRef
A line using referred points.
Definition: line.H:66
label commonVertex(const edge &other) const
Return vertex common with other edge or -1 on failure.
Definition: edgeI.H:142
static int compare(const Pair< label > &a, const Pair< label > &b)
Compare Pairs.
Definition: PairI.H:24
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
label which(const label vertex) const
Return local index (0,1) of vertex label in edge -1 on failure.
Definition: edgeI.H:117
label collapse()
&#39;Collapse&#39; edge by marking duplicate vertex labels as &#39;-1&#39;, the lower vertex is retained.
Definition: edgeI.H:174
vector vec(const UList< point > &pts) const
Return the vector (from first to second).
Definition: edgeI.H:387
errorManip< error > abort(error &err)
Definition: errorManip.H:139
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:321
Pair< point > box(const UList< point > &pts) const
The enclosing (bounding) box for the edge.
Definition: edgeI.H:433
const direction noexcept
Definition: Scalar.H:258
bool good() const noexcept
True if the vertices are unique and non-negative.
Definition: edgeI.H:100
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
label a() const noexcept
The first vertex.
Definition: edge.H:127
label otherVertex(const label vertex) const
Given one vertex label, return the other one.
Definition: edgeI.H:158
bool insert(const label vertex)
Fill any open slot with the vertex label (if not previously contained in the edge).
Definition: edgeI.H:230
label maxVertex() const noexcept
Return the largest vertex label used by the edge.
Definition: edgeI.H:94
edge()
Default construct, with invalid vertex labels (-1)
Definition: edgeI.H:34
scalar magSqr(const UList< point > &pts) const
The length (L2-norm) squared of the edge vector.
Definition: edgeI.H:426
bool contains(const label vertex) const noexcept
Return true if the vertex label is contained in the edge.
Definition: edgeI.H:106
point centre(const UList< point > &pts) const
Return centre point (centroid) of the edge.
Definition: edgeI.H:372
label n
bool operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:297
iterator end() noexcept
Return an iterator to end traversing the UList.
Definition: UListI.H:365
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
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))
label & operator[](const label i)
Return edge element. Index should be limited to 0/1.
Definition: edgeI.H:465
static int compare(const edge &a, const edge &b)
Compare edges.
Definition: edgeI.H:26
const pointField & pts
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133