All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
indexedCellI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 template<class Gt, class Cb>
28 (
29  const Foam::globalIndex& globalDelaunayVertexIndices
30 ) const
31 {
32  Foam::tetCell tVGI;
33 
34  for (int i = 0; i < 4; i++)
35  {
36  Vertex_handle v = this->vertex(i);
37 
38  // Finding the global index of each Delaunay vertex
39  tVGI[i] = globalDelaunayVertexIndices.toGlobal
40  (
42  v->index()
43  );
44  }
45 
46  return tVGI;
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 template<class Gt, class Cb>
54 :
55  Cb(),
56  index_(ctUnassigned),
57  filterCount_(0)
58 {}
59 
60 
61 template<class Gt, class Cb>
63 (
65 )
66 :
67  Cb(v0, v1, v2, v3),
68  index_(ctUnassigned),
69  filterCount_(0)
70 {}
71 
72 
73 template<class Gt, class Cb>
75 (
76  Vertex_handle v0,
77  Vertex_handle v1,
78  Vertex_handle v2,
79  Vertex_handle v3,
80  Cell_handle n0,
81  Cell_handle n1,
82  Cell_handle n2,
83  Cell_handle n3
84 )
85 :
86  Cb(v0, v1, v2, v3, n0, n1, n2, n3),
87  index_(ctUnassigned),
88  filterCount_(0)
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 
94 template<class Gt, class Cb>
96 {
97  return index_;
98 }
99 
100 
101 template<class Gt, class Cb>
103 {
104  return index_;
105 }
106 
107 
108 #ifdef CGAL_INEXACT
109 
110  template<class Gt, class Cb>
112  {
113  return reinterpret_cast<const Foam::point&>(this->circumcenter());
114  }
115 
116 #else
117 
118  template<class Gt, class Cb>
120  {
121  const typename Gt::Point_3& P = this->circumcenter();
122 
123  return Foam::point
124  (
125  CGAL::to_double(P.x()),
126  CGAL::to_double(P.y()),
127  CGAL::to_double(P.z())
128  );
129  }
130 
131 #endif
132 
133 
134 template<class Gt, class Cb>
136 {
137  return index_ == ctUnassigned;
138 }
139 
140 
141 template<class Gt, class Cb>
143 {
144  return filterCount_;
145 }
146 
147 
148 template<class Gt, class Cb>
150 {
151  return filterCount_;
152 }
153 
154 
155 template<class Gt, class Cb>
157 {
158  return
159  (
160  (
161  this->vertex(0)->real()
162  || this->vertex(1)->real()
163  || this->vertex(2)->real()
164  || this->vertex(3)->real()
165  )
166  &&
167  !(
168  this->vertex(0)->farPoint()
169  || this->vertex(1)->farPoint()
170  || this->vertex(2)->farPoint()
171  || this->vertex(3)->farPoint()
172  )
173  );
174 }
175 
176 
177 template<class Gt, class Cb>
179 {
180  return
181  (
182  this->vertex(0)->farPoint()
183  || this->vertex(1)->farPoint()
184  || this->vertex(2)->farPoint()
185  || this->vertex(3)->farPoint()
186  );
187 }
188 
189 
190 template<class Gt, class Cb>
192 {
193  return
194  (
195  this->vertex(0)->referred()
196  || this->vertex(1)->referred()
197  || this->vertex(2)->referred()
198  || this->vertex(3)->referred()
199  );
200 }
201 
202 
203 template<class Gt, class Cb>
205 {
206  return
207  (
208  this->vertex(0)->featurePoint()
209  || this->vertex(1)->featurePoint()
210  || this->vertex(2)->featurePoint()
211  || this->vertex(3)->featurePoint()
212  );
213 }
214 
215 
216 template<class Gt, class Cb>
218 {
219  return
220  (
221  this->vertex(0)->seedPoint()
222  || this->vertex(1)->seedPoint()
223  || this->vertex(2)->seedPoint()
224  || this->vertex(3)->seedPoint()
225  );
226 }
227 
228 
229 template<class Gt, class Cb>
231 {
232  return
233  (
234  this->vertex(0)->internalPoint()
235  || this->vertex(1)->internalPoint()
236  || this->vertex(2)->internalPoint()
237  || this->vertex(3)->internalPoint()
238  );
239 }
240 
241 
242 template<class Gt, class Cb>
244 {
245  return
246  (
247  this->vertex(0)->boundaryPoint()
248  || this->vertex(1)->boundaryPoint()
249  || this->vertex(2)->boundaryPoint()
250  || this->vertex(3)->boundaryPoint()
251  );
252 }
253 
254 
255 template<class Gt, class Cb>
257 {
258  return
259  (
260  this->vertex(0)->constrained()
261  || this->vertex(1)->constrained()
262  || this->vertex(2)->constrained()
263  || this->vertex(3)->constrained()
264  );
265 }
266 
267 
268 template<class Gt, class Cb>
270 {
271  return
272  (
273  !this->hasFarPoint()
274  &&
275  (
276  this->vertex(0)->referred()
277  || this->vertex(1)->referred()
278  || this->vertex(2)->referred()
279  || this->vertex(3)->referred()
280  )
281  &&
282  (
283  this->vertex(0)->real()
284  || this->vertex(1)->real()
285  || this->vertex(2)->real()
286  || this->vertex(3)->real()
287  )
288  );
289 }
290 
291 
292 template<class Gt, class Cb>
294 {
295  Foam::label lowestProc = -1;
296 
297  for (int i = 0; i < 4; ++i)
298  {
299  if (this->vertex(i)->referred())
300  {
301  lowestProc = min(lowestProc, this->vertex(i)->procIndex());
302  }
303  }
304 
305  return lowestProc;
306 }
307 
308 
309 template<class Gt, class Cb>
311 (
312  const Foam::globalIndex& globalDelaunayVertexIndices
313 ) const
314 {
315  // tetVertexGlobalIndices
316  Foam::tetCell tVGI
317  = unsortedVertexGlobalIndices(globalDelaunayVertexIndices);
318 
319  // bubble sort
320  for (int i = 0; i < tVGI.size(); i++)
321  {
322  for (int j = tVGI.size() - 1 ; j > i; j--)
323  {
324  if (tVGI[j - 1] > tVGI[j])
325  {
326  Foam::Swap(tVGI[j - 1], tVGI[j]);
327  }
328  }
329  }
330 
331  return tVGI;
332 }
333 
334 
335 template<class Gt, class Cb>
338 (
339  const Foam::globalIndex& globalDelaunayVertexIndices
340 ) const
341 {
342  // tetVertexGlobalIndices
343  Foam::tetCell tVGI
344  = unsortedVertexGlobalIndices(globalDelaunayVertexIndices);
345 
347 
348  // bubble sort
349  for (int i = 0; i < tVGI.size(); i++)
350  {
351  for (int j = tVGI.size() - 1 ; j > i; j--)
352  {
353  if (tVGI[j - 1] > tVGI[j])
354  {
355  Foam::Swap(tVGI[j - 1], tVGI[j]);
356  Foam::Swap(vertexMap[j - 1], vertexMap[j]);
357  }
358  }
359  }
360 
361  for (int i = 0; i < 4; i++)
362  {
363  tVGI[i] = vertexMap[i];
364  }
365 
366  return tVGI;
367 }
368 
369 
370 template<class Gt, class Cb>
372 {
373  return
374  (
375  this->vertex(0)->internalOrBoundaryPoint()
376  || this->vertex(1)->internalOrBoundaryPoint()
377  || this->vertex(2)->internalOrBoundaryPoint()
378  || this->vertex(3)->internalOrBoundaryPoint()
379  );
380 }
381 
382 
383 template<class Gt, class Cb>
385 {
386  return
387  (
388  this->vertex(0)->internalOrBoundaryPoint()
389  || this->vertex(0)->externalBoundaryPoint()
390  || this->vertex(1)->internalOrBoundaryPoint()
391  || this->vertex(1)->externalBoundaryPoint()
392  || this->vertex(2)->internalOrBoundaryPoint()
393  || this->vertex(2)->externalBoundaryPoint()
394  || this->vertex(3)->internalOrBoundaryPoint()
395  || this->vertex(3)->externalBoundaryPoint()
396  );
397 }
398 
399 
400 template<class Gt, class Cb>
402 {
403 // return
404 // (
405 // this->vertex(0)->boundaryPoint()
406 // && this->vertex(1)->boundaryPoint()
407 // && this->vertex(2)->boundaryPoint()
408 // && this->vertex(3)->boundaryPoint()
409 // );
410 
411  return
412  (
413  (
414  this->vertex(0)->internalBoundaryPoint()
415  || this->vertex(1)->internalBoundaryPoint()
416  || this->vertex(2)->internalBoundaryPoint()
417  || this->vertex(3)->internalBoundaryPoint()
418  )
419  && (
420  this->vertex(0)->externalBoundaryPoint()
421  || this->vertex(1)->externalBoundaryPoint()
422  || this->vertex(2)->externalBoundaryPoint()
423  || this->vertex(3)->externalBoundaryPoint()
424  )
425  );
426 
427 // Foam::label nBoundaryPoints = 0;
428 //
429 // for (Foam::label i = 0; i < 4; ++i)
430 // {
431 // Vertex_handle v = this->vertex(i);
432 //
433 // if (v->boundaryPoint())
434 // {
435 // nBoundaryPoints++;
436 // }
437 // }
438 //
439 // return (nBoundaryPoints > 1);
440 }
441 
442 
443 template<class Gt, class Cb>
445 {
446  return
447  (
448  (
449  this->vertex(0)->internalBaffleSurfacePoint()
450  || this->vertex(1)->internalBaffleSurfacePoint()
451  || this->vertex(2)->internalBaffleSurfacePoint()
452  || this->vertex(3)->internalBaffleSurfacePoint()
453  )
454  && (
455  this->vertex(0)->externalBaffleSurfacePoint()
456  || this->vertex(1)->externalBaffleSurfacePoint()
457  || this->vertex(2)->externalBaffleSurfacePoint()
458  || this->vertex(3)->externalBaffleSurfacePoint()
459  )
460  );
461 }
462 
463 
464 template<class Gt, class Cb>
466 {
467  return
468  (
469  (
470  this->vertex(0)->internalBaffleEdgePoint()
471  || this->vertex(1)->internalBaffleEdgePoint()
472  || this->vertex(2)->internalBaffleEdgePoint()
473  || this->vertex(3)->internalBaffleEdgePoint()
474  )
475  && (
476  this->vertex(0)->externalBaffleEdgePoint()
477  || this->vertex(1)->externalBaffleEdgePoint()
478  || this->vertex(2)->externalBaffleEdgePoint()
479  || this->vertex(3)->externalBaffleEdgePoint()
480  )
481  );
482 }
483 
484 
485 template<class Gt, class Cb>
487 {
488  return
489  (
490  this->vertex(0)->featureEdgePoint()
491  && this->vertex(1)->featureEdgePoint()
492  && this->vertex(2)->featureEdgePoint()
493  && this->vertex(3)->featureEdgePoint()
494  );
495 // (
496 // this->vertex(0)->featureEdgePoint()
497 // || this->vertex(1)->featureEdgePoint()
498 // || this->vertex(2)->featureEdgePoint()
499 // || this->vertex(3)->featureEdgePoint()
500 // )
501 // && (
502 // (
503 // this->vertex(0)->featureEdgePoint()
504 // || this->vertex(0)->featurePoint()
505 // )
506 // && (
507 // this->vertex(1)->featureEdgePoint()
508 // || this->vertex(1)->featurePoint()
509 // )
510 // && (
511 // this->vertex(2)->featureEdgePoint()
512 // || this->vertex(2)->featurePoint()
513 // )
514 // && (
515 // this->vertex(3)->featureEdgePoint()
516 // || this->vertex(3)->featurePoint()
517 // )
518 // )
519 // );
520 }
521 
522 
523 template<class Gt, class Cb>
525 {
526  return
527  (
528  this->vertex(0)->featurePoint()
529  && this->vertex(1)->featurePoint()
530  && this->vertex(2)->featurePoint()
531  && this->vertex(3)->featurePoint()
532  );
533 }
534 
535 
536 template<class Gt, class Cb>
538 {
539  return
540  (
541  this->vertex(0)->nearProcBoundary()
542  || this->vertex(1)->nearProcBoundary()
543  || this->vertex(2)->nearProcBoundary()
544  || this->vertex(3)->nearProcBoundary()
545  );
546 }
547 
548 
549 template<class Gt, class Cb>
551 {
552  Foam::label nMasters = 0;
553  Foam::label nSlaves = 0;
554 
555  Vertex_handle vM[2];
556  Vertex_handle vS[2];
557 
558  for (Foam::label i = 0; i < 4; ++i)
559  {
560  Vertex_handle v = this->vertex(i);
561 
562  if (v->internalBoundaryPoint())
563  {
564  vM[nMasters] = v;
565  nMasters++;
566  }
567 
568  if (v->externalBoundaryPoint())
569  {
570  vS[nSlaves] = v;
571  nSlaves++;
572  }
573  }
574 
575  Foam::label nPairs = 0;
576 
577  if (nMasters == 2 && nSlaves == 2)
578  {
581 
582  if
583  (
584  vM[0]->type() == vS[0]->index()
585  && vM[0]->index() == vS[0]->type()
586  )
587  {
588  vp0 = reinterpret_cast<const Foam::point&>(vM[0]->point())
589  - reinterpret_cast<const Foam::point&>(vS[0]->point());
590  nPairs++;
591  }
592  else if
593  (
594  vM[0]->type() == vS[1]->index()
595  && vM[0]->index() == vS[1]->type()
596  )
597  {
598  vp0 = reinterpret_cast<const Foam::point&>(vM[0]->point())
599  - reinterpret_cast<const Foam::point&>(vS[1]->point());
600  nPairs++;
601  }
602 
603  if
604  (
605  vM[1]->type() == vS[0]->index()
606  && vM[1]->index() == vS[0]->type()
607  )
608  {
609  vp1 = reinterpret_cast<const Foam::point&>(vM[1]->point())
610  - reinterpret_cast<const Foam::point&>(vS[0]->point());
611  nPairs++;
612  }
613  else if
614  (
615  vM[1]->type() == vS[1]->index()
616  && vM[1]->index() == vS[1]->type()
617  )
618  {
619  vp1 = reinterpret_cast<const Foam::point&>(vM[1]->point())
620  - reinterpret_cast<const Foam::point&>(vS[1]->point());
621  nPairs++;
622  }
623 
624  if (nPairs == 2)
625  {
626  if (Foam::vectorTools::areParallel(vp0, vp1))
627  {
628  Foam::Pout<< "PARALLEL" << Foam::endl;
629 
630  return true;
631  }
632  }
633  }
634 
635  return false;
636 }
637 
638 
639 template<class Gt, class Cb>
641 {
642  int featureVertex = -1;
643  for (int i = 0; i < 4; ++i)
644  {
645  if (this->vertex(i)->constrained())
646  {
647  featureVertex = i;
648  }
649  }
650 
651  // Pick cell with a face attached to an infinite cell
652 
653  if (featureVertex != -1)
654  {
655  Vertex_handle v1 =
656  this->vertex(Tds::vertex_triple_index(featureVertex, 0));
657  Vertex_handle v2 =
658  this->vertex(Tds::vertex_triple_index(featureVertex, 1));
659  Vertex_handle v3 =
660  this->vertex(Tds::vertex_triple_index(featureVertex, 2));
661 
662  if (v1->internalBoundaryPoint())
663  {
664  if
665  (
666  v2->externalBoundaryPoint()
667  && v3->externalBoundaryPoint()
668  )
669  {
670  return true;
671  }
672  }
673  else if (v2->internalBoundaryPoint())
674  {
675  if
676  (
677  v1->externalBoundaryPoint()
678  && v3->externalBoundaryPoint()
679  )
680  {
681  return true;
682  }
683  }
684  else if (v3->internalBoundaryPoint())
685  {
686  if
687  (
688  v1->externalBoundaryPoint()
689  && v2->externalBoundaryPoint()
690  )
691  {
692  return true;
693  }
694  }
695  }
696 
697  return false;
698 }
699 
700 
701 template<class Gt, class Cb>
703 {
704  int featureVertex = -1;
705  for (int i = 0; i < 4; ++i)
706  {
707  if (this->vertex(i)->constrained())
708  {
709  featureVertex = i;
710  }
711  }
712 
713  // Pick cell with a face attached to an infinite cell
714 
715  if (featureVertex != -1)
716  {
717  Vertex_handle v1 =
718  this->vertex(Tds::vertex_triple_index(featureVertex, 0));
719  Vertex_handle v2 =
720  this->vertex(Tds::vertex_triple_index(featureVertex, 1));
721  Vertex_handle v3 =
722  this->vertex(Tds::vertex_triple_index(featureVertex, 2));
723 
724  if (v1->externalBoundaryPoint())
725  {
726  if
727  (
728  v2->internalBoundaryPoint()
729  && v3->internalBoundaryPoint()
730  )
731  {
732  return true;
733  }
734  }
735  else if (v2->externalBoundaryPoint())
736  {
737  if
738  (
739  v1->internalBoundaryPoint()
740  && v3->internalBoundaryPoint()
741  )
742  {
743  return true;
744  }
745  }
746  else if (v3->externalBoundaryPoint())
747  {
748  if
749  (
750  v1->internalBoundaryPoint()
751  && v2->internalBoundaryPoint()
752  )
753  {
754  return true;
755  }
756  }
757  }
758 
759  return false;
760 }
761 
762 
763 // * * * * * * * * * * * * * * * Friend Functions * * * * * * * * * * * * * //
bool baffleSurfaceDualVertex() const
Definition: indexedCellI.H:444
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::label vertexLowestProc() const
Definition: indexedCellI.H:293
CGAL::indexedCell< K > Cb
bool parallelDualVertex() const
Does the Dual vertex form part of a processor patch.
Definition: indexedCellI.H:269
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
Foam::label & cellIndex()
Definition: indexedCellI.H:95
Cb::Cell_handle Cell_handle
Definition: indexedCell.H:119
bool hasFarPoint() const
Does the Delaunay cell have a far point.
Definition: indexedCellI.H:178
bool featurePointInternalCell() const
Definition: indexedCellI.H:702
bool featureEdgeDualVertex() const
A dual vertex on a feature edge will result from this Delaunay cell.
Definition: indexedCellI.H:486
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const Foam::point dual()
Definition: indexedCellI.H:119
An indexed form of CGAL::Triangulation_cell_base_3<K> used to keep track of the Delaunay cells (tets)...
Definition: indexedCell.H:56
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
bool baffleEdgeDualVertex() const
Definition: indexedCellI.H:465
bool featurePointExternalCell() const
Definition: indexedCellI.H:640
Cb::Vertex_handle Vertex_handle
Definition: indexedCell.H:118
bool unassigned() const
Definition: indexedCellI.H:135
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
bool anyInternalOrBoundaryDualVertex() const
Is the Delaunay cell real or referred (or mixed), i.e. all vertices.
Definition: indexedCellI.H:384
void Swap(T &a, T &b)
Definition: Swap.H:43
bool hasInternalPoint() const
Definition: indexedCellI.H:230
bool hasConstrainedPoint() const
Definition: indexedCellI.H:256
static const zero Zero
Definition: zero.H:97
A tetrahedral cell primitive.
Definition: tetCell.H:58
bool hasBoundaryPoint() const
Definition: indexedCellI.H:243
bool hasFeaturePoint() const
Does the Delaunay cell have a feature point.
Definition: indexedCellI.H:204
label size() const
Return the number of elements in the FixedList.
Definition: FixedListI.H:419
bool nearProcBoundary() const
Definition: indexedCellI.H:537
bool potentialCoplanarCell() const
Definition: indexedCellI.H:550
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
bool real() const
Is the Delaunay cell real, i.e. any real vertex.
Definition: indexedCellI.H:156
vector point
Point is a vector.
Definition: point.H:41
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
bool areParallel(const Vector< T > &a, const Vector< T > &b, const T &tolerance=small)
Test if a and b are parallel: a^b = 0.
Definition: vectorTools.H:54
Foam::FixedList< Foam::label, 4 > globallyOrderedCellVertices(const Foam::globalIndex &globalDelaunayVertexIndices) const
Using the globalIndex object, return a list of four vertices with.
Definition: indexedCellI.H:338
bool featurePointDualVertex() const
A dual vertex on a feature point will result from this Delaunay cell.
Definition: indexedCellI.H:524
bool hasSeedPoint() const
Does the Delaunay cell have a seed point.
Definition: indexedCellI.H:217
bool hasReferredPoint() const
Does the Delaunay cell have a referred point.
Definition: indexedCellI.H:191
bool internalOrBoundaryDualVertex() const
Is the Delaunay cell part of the final dual mesh, i.e. any vertex.
Definition: indexedCellI.H:371
Foam::tetCell vertexGlobalIndices(const Foam::globalIndex &globalDelaunayVertexIndices) const
Using the globalIndex object, return a list of four (sorted) global.
Definition: indexedCellI.H:311
bool boundaryDualVertex() const
A dual vertex on the boundary will result from a Delaunay cell with.
Definition: indexedCellI.H:401