conformalVoronoiMesh.C
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) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2020-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 "conformalVoronoiMesh.H"
30 #include "initialPointsMethod.H"
31 #include "relaxationModel.H"
32 #include "faceAreaWeightModel.H"
33 #include "meshSearch.H"
34 #include "vectorTools.H"
35 #include "IOmanip.H"
36 #include "indexedCellChecks.H"
37 #include "controlMeshRefinement.H"
38 #include "smoothAlignmentSolver.H"
39 #include "OBJstream.H"
40 #include "indexedVertexOps.H"
41 #include "DelaunayMeshTools.H"
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47  defineTypeNameAndDebug(conformalVoronoiMesh, 0);
48 }
49 
50 const Foam::Enum
51 <
53 >
55 ({
56  { dualMeshPointType::internal, "internal" },
57  { dualMeshPointType::surface, "surface" },
58  { dualMeshPointType::featureEdge, "featureEdge" },
59  { dualMeshPointType::featurePoint, "featurePoint" },
60  { dualMeshPointType::constrained, "constrained" },
61 });
62 
63 
64 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
65 
66 void Foam::conformalVoronoiMesh::cellSizeMeshOverlapsBackground() const
67 {
68  const cellShapeControlMesh& cellSizeMesh =
69  cellShapeControl_.shapeControlMesh();
70 
71  DynamicList<Foam::point> pts(number_of_vertices());
72 
73  for
74  (
75  Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
76  vit != finite_vertices_end();
77  ++vit
78  )
79  {
80  if (vit->internalOrBoundaryPoint() && !vit->referred())
81  {
82  pts.append(topoint(vit->point()));
83  }
84  }
85 
86  boundBox bb(pts);
87 
88  boundBox cellSizeMeshBb = cellSizeMesh.bounds();
89 
90  bool fullyContained = cellSizeMeshBb.contains(bb);
91 
92  if (!fullyContained)
93  {
94  Pout<< "Triangulation not fully contained in cell size mesh." << endl
95  << "Cell Size Mesh Bounds = " << cellSizeMeshBb << endl
96  << "foamyHexMesh Bounds = " << bb << endl;
97  }
98 
99  Info<< "Triangulation is "
100  << (returnReduceAnd(fullyContained) ? "fully" : "not fully")
101  << " contained in the cell size mesh"
102  << endl;
103 }
104 
105 
106 void Foam::conformalVoronoiMesh::insertInternalPoints
107 (
108  List<Point>& points,
109  bool distribute
110 )
111 {
112  const label nPoints = returnReduce(points.size(), sumOp<label>());
113 
114  Info<< " " << nPoints << " points to insert..." << endl;
115 
116  if (Pstream::parRun() && distribute)
117  {
118  List<Foam::point> transferPoints(points.size());
119 
120  forAll(points, pI)
121  {
122  transferPoints[pI] = topoint(points[pI]);
123  }
124 
125  // Send the points that are not on this processor to the appropriate
126  // place
128  (
129  decomposition_().distributePoints(transferPoints)
130  );
131 
132  transferPoints.clear();
133 
134  map().distribute(points);
135  }
136 
137  label preReinsertionSize(number_of_vertices());
138 
139  insert(points.begin(), points.end());
140 
141  const label nInserted = returnReduce
142  (
143  label(number_of_vertices()) - preReinsertionSize,
144  sumOp<label>()
145  );
146 
147  Info<< " " << nInserted << " points inserted"
148  << ", failed to insert " << nPoints - nInserted
149  << " ("
150  << 100.0*(nPoints - nInserted)/(nInserted + SMALL)
151  << " %)"<< endl;
152 
153  for
154  (
155  Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
156  vit != finite_vertices_end();
157  ++vit
158  )
159  {
161  {
162  vit->index() = getNewVertexIndex();
163  vit->type() = Vb::vtInternal;
164  }
165  }
166 }
167 
168 
169 Foam::Map<Foam::label> Foam::conformalVoronoiMesh::insertPointPairs
170 (
171  List<Vb>& vertices,
172  bool distribute,
173  bool reIndex
174 )
175 {
176  if (Pstream::parRun() && distribute)
177  {
178  autoPtr<mapDistribute> mapDist =
179  decomposition_().distributePoints(vertices);
180 
181  // Re-index the point pairs if one or both have been distributed.
182  // If both, remove
183 
184  // If added a point, then need to know its point pair
185  // If one moved, then update procIndex locally
186 
187  forAll(vertices, vI)
188  {
189  vertices[vI].procIndex() = Pstream::myProcNo();
190  }
191  }
192 
193  label preReinsertionSize(number_of_vertices());
194 
195  Map<label> oldToNewIndices =
196  this->DelaunayMesh<Delaunay>::insertPoints(vertices, reIndex);
197 
198  const label nReinserted = returnReduce
199  (
200  label(number_of_vertices()) - preReinsertionSize,
201  sumOp<label>()
202  );
203 
204  Info<< " Reinserted " << nReinserted << " vertices out of "
205  << returnReduce(vertices.size(), sumOp<label>())
206  << endl;
207 
208  return oldToNewIndices;
209 }
210 
211 
212 void Foam::conformalVoronoiMesh::insertSurfacePointPairs
213 (
214  const pointIndexHitAndFeatureList& surfaceHits,
215  const fileName fName,
216  DynamicList<Vb>& pts
217 )
218 {
219  forAll(surfaceHits, i)
220  {
221  vectorField norm(1);
222 
223  const pointIndexHit surfaceHit = surfaceHits[i].first();
224  const label featureIndex = surfaceHits[i].second();
225 
226  allGeometry_[featureIndex].getNormal
227  (
228  List<pointIndexHit>(1, surfaceHit),
229  norm
230  );
231 
232  const vector& normal = norm[0];
233 
234  const Foam::point& surfacePt = surfaceHit.hitPoint();
235 
237  geometryToConformTo_.meshableSide(featureIndex, surfaceHit);
238 
239  if (meshableSide == extendedFeatureEdgeMesh::BOTH)
240  {
241  createBafflePointPair
242  (
243  pointPairDistance(surfacePt),
244  surfacePt,
245  normal,
246  true,
247  pts
248  );
249  }
250  else if (meshableSide == extendedFeatureEdgeMesh::INSIDE)
251  {
252  createPointPair
253  (
254  pointPairDistance(surfacePt),
255  surfacePt,
256  normal,
257  true,
258  pts
259  );
260  }
261  else if (meshableSide == extendedFeatureEdgeMesh::OUTSIDE)
262  {
263  createPointPair
264  (
265  pointPairDistance(surfacePt),
266  surfacePt,
267  -normal,
268  true,
269  pts
270  );
271  }
272  else
273  {
275  << meshableSide << ", bad"
276  << endl;
277  }
278  }
279 
280  if (foamyHexMeshControls().objOutput() && !fName.empty())
281  {
282  DelaunayMeshTools::writeOBJ(time().path()/fName, pts);
283  }
284 }
285 
286 
287 void Foam::conformalVoronoiMesh::insertEdgePointGroups
288 (
289  const pointIndexHitAndFeatureList& edgeHits,
290  const fileName fName,
291  DynamicList<Vb>& pts
292 )
293 {
294  forAll(edgeHits, i)
295  {
296  if (edgeHits[i].first().hit())
297  {
298  const extendedFeatureEdgeMesh& feMesh
299  (
300  geometryToConformTo_.features()[edgeHits[i].second()]
301  );
302 
303 // const bool isBaffle =
304 // geometryToConformTo_.isBaffleFeature(edgeHits[i].second());
305 
306  createEdgePointGroup
307  (
308  feMesh,
309  edgeHits[i].first(),
310  pts
311  );
312  }
313  }
314 
315  if (foamyHexMeshControls().objOutput() && !fName.empty())
316  {
317  DelaunayMeshTools::writeOBJ(time().path()/fName, pts);
318  }
319 }
320 
321 
322 bool Foam::conformalVoronoiMesh::nearFeaturePt(const Foam::point& pt) const
323 {
324  scalar exclusionRangeSqr = featurePointExclusionDistanceSqr(pt);
325 
326  pointIndexHit info;
327  label featureHit;
328 
329  geometryToConformTo_.findFeaturePointNearest
330  (
331  pt,
332  exclusionRangeSqr,
333  info,
334  featureHit
335  );
336 
337  return info.hit();
338 }
339 
340 
341 bool Foam::conformalVoronoiMesh::surfacePtNearFeatureEdge
342 (
343  const Foam::point& pt
344 ) const
345 {
346  scalar exclusionRangeSqr = surfacePtExclusionDistanceSqr(pt);
347 
348  pointIndexHit info;
349  label featureHit;
350 
351  geometryToConformTo_.findEdgeNearest
352  (
353  pt,
354  exclusionRangeSqr,
355  info,
356  featureHit
357  );
358 
359  return info.hit();
360 }
361 
362 
363 void Foam::conformalVoronoiMesh::insertInitialPoints()
364 {
365  Info<< nl << "Inserting initial points" << endl;
366 
367  timeCheck("Before initial points call");
368 
369  List<Point> initPts = initialPointsMethod_->initialPoints();
370 
371  timeCheck("After initial points call");
372 
373  // Assume that the initial points method made the correct decision for
374  // which processor each point should be on, so give distribute = false
375  insertInternalPoints(initPts);
376 
377  if (initialPointsMethod_->fixInitialPoints())
378  {
379  for
380  (
381  Finite_vertices_iterator vit = finite_vertices_begin();
382  vit != finite_vertices_end();
383  ++vit
384  )
385  {
386  vit->fixed() = true;
387  }
388  }
389 
390  if (foamyHexMeshControls().objOutput())
391  {
393  (
394  time().path()/"initialPoints.obj",
395  *this,
397  );
398  }
399 }
400 
401 
402 void Foam::conformalVoronoiMesh::distribute()
403 {
404  if (!Pstream::parRun())
405  {
406  return;
407  }
408 
409  DynamicList<Foam::point> points(number_of_vertices());
410  DynamicList<Foam::indexedVertexEnum::vertexType> types
411  (
412  number_of_vertices()
413  );
414  DynamicList<scalar> sizes(number_of_vertices());
415  DynamicList<tensor> alignments(number_of_vertices());
416 
417  for
418  (
419  Finite_vertices_iterator vit = finite_vertices_begin();
420  vit != finite_vertices_end();
421  ++vit
422  )
423  {
424  if (vit->real())
425  {
426  points.append(topoint(vit->point()));
427  types.append(vit->type());
428  sizes.append(vit->targetCellSize());
429  alignments.append(vit->alignment());
430  }
431  }
432 
433  autoPtr<mapDistribute> mapDist =
435 
436  mapDist().distribute(types);
437  mapDist().distribute(sizes);
438  mapDist().distribute(alignments);
439 
440  // Reset the entire tessellation
442 
443  Info<< nl << " Inserting distributed tessellation" << endl;
444 
445  // Internal points have to be inserted first
446 
447  DynamicList<Vb> verticesToInsert(points.size());
448 
449  forAll(points, pI)
450  {
451  verticesToInsert.append
452  (
453  Vb
454  (
455  toPoint(points[pI]),
456  -1,
457  types[pI],
459  )
460  );
461 
462  verticesToInsert.last().targetCellSize() = sizes[pI];
463  verticesToInsert.last().alignment() = alignments[pI];
464  }
465 
466  this->rangeInsertWithInfo
467  (
468  verticesToInsert.begin(),
469  verticesToInsert.end(),
470  true
471  );
472 
473  Info<< " Total number of vertices after redistribution "
474  << returnReduce
475  (
476  label(number_of_vertices()), sumOp<label>()
477  )
478  << endl;
479 }
480 
481 
482 void Foam::conformalVoronoiMesh::buildCellSizeAndAlignmentMesh()
483 {
484  controlMeshRefinement meshRefinement
485  (
486  cellShapeControl_
487  );
488 
489  smoothAlignmentSolver meshAlignmentSmoother
490  (
491  cellShapeControl_.shapeControlMesh()
492  );
493 
494  meshRefinement.initialMeshPopulation(decomposition_);
495 
496  cellShapeControlMesh& cellSizeMesh = cellShapeControl_.shapeControlMesh();
497 
498  if (Pstream::parRun())
499  {
500  if (!distributeBackground(cellSizeMesh))
501  {
502  // Synchronise the cell size mesh if it has not been distributed
503  cellSizeMesh.distribute(decomposition_());
504  }
505  }
506 
507  const dictionary& motionControlDict
508  = foamyHexMeshControls().foamyHexMeshDict().subDict("motionControl");
509 
510  const label nMaxIter =
511  motionControlDict.get<label>("maxRefinementIterations");
512 
513  Info<< "Maximum number of refinement iterations : " << nMaxIter << endl;
514 
515  for (label i = 0; i < nMaxIter; ++i)
516  {
517  label nAdded = meshRefinement.refineMesh(decomposition_);
518  //cellShapeControl_.refineMesh(decomposition_);
519  reduce(nAdded, sumOp<label>());
520 
521  if (Pstream::parRun())
522  {
523  cellSizeMesh.distribute(decomposition_());
524  }
525 
526  Info<< " Iteration " << i
527  << " Added = " << nAdded << " points"
528  << endl;
529 
530  if (nAdded == 0)
531  {
532  break;
533  }
534  }
535 
536  if (Pstream::parRun())
537  {
538  // Need to distribute the cell size mesh to cover the background mesh
539  if (!distributeBackground(cellSizeMesh))
540  {
541  cellSizeMesh.distribute(decomposition_());
542  }
543  }
544 
545  meshAlignmentSmoother.smoothAlignments
546  (
547  motionControlDict.get<label>("maxSmoothingIterations")
548  );
549 
550  Info<< "Background cell size and alignment mesh:" << endl;
551  cellSizeMesh.printInfo(Info);
552 
553  Info<< "Triangulation is "
554  << (cellSizeMesh.is_valid() ? "valid" : "not valid!" )
555  << endl;
556 
557  if (foamyHexMeshControls().writeCellShapeControlMesh())
558  {
559  //cellSizeMesh.writeTriangulation();
560  cellSizeMesh.write();
561  }
562 
563  if (foamyHexMeshControls().printVertexInfo())
564  {
565  cellSizeMesh.printVertexInfo(Info);
566  }
567 
568 // Info<< "Estimated number of cells in final mesh = "
569 // << returnReduce
570 // (
571 // cellSizeMesh.estimateCellCount(decomposition_),
572 // sumOp<label>()
573 // )
574 // << endl;
575 }
576 
577 
578 void Foam::conformalVoronoiMesh::setVertexSizeAndAlignment()
579 {
580  Info<< nl << "Calculating target cell alignment and size" << endl;
581 
582  for
583  (
584  Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
585  vit != finite_vertices_end();
586  vit++
587  )
588  {
589  if (vit->internalOrBoundaryPoint())
590  {
591  pointFromPoint pt = topoint(vit->point());
592 
593  cellShapeControls().cellSizeAndAlignment
594  (
595  pt,
596  vit->targetCellSize(),
597  vit->alignment()
598  );
599  }
600  }
601 }
602 
603 
604 Foam::face Foam::conformalVoronoiMesh::buildDualFace
605 (
606  const Delaunay::Finite_edges_iterator& eit
607 ) const
608 {
609  Cell_circulator ccStart = incident_cells(*eit);
610  Cell_circulator cc1 = ccStart;
611  Cell_circulator cc2 = cc1;
612 
613  // Advance the second circulator so that it always stays on the next
614  // cell around the edge;
615  cc2++;
616 
617  DynamicList<label> verticesOnFace;
618 
619  label nUniqueVertices = 0;
620 
621  do
622  {
623  if
624  (
625  cc1->hasFarPoint() || cc2->hasFarPoint()
626  || is_infinite(cc1) || is_infinite(cc2)
627  )
628  {
629  Cell_handle c = eit->first;
630  Vertex_handle vA = c->vertex(eit->second);
631  Vertex_handle vB = c->vertex(eit->third);
632 
633 // DelaunayMeshTools::drawDelaunayCell(Pout, cc1);
634 // DelaunayMeshTools::drawDelaunayCell(Pout, cc2);
635 
637  << "Dual face uses circumcenter defined by a "
638  << "Delaunay tetrahedron with no internal "
639  << "or boundary points. Defining Delaunay edge ends: "
640  << vA->info() << " "
641  << vB->info() << nl
642  <<endl;//<< exit(FatalError);
643  }
644  else
645  {
646  label cc1I = cc1->cellIndex();
647  label cc2I = cc2->cellIndex();
648 
649  if (cc1I != cc2I)
650  {
651  if (!verticesOnFace.found(cc1I))
652  {
653  nUniqueVertices++;
654  }
655 
656  verticesOnFace.append(cc1I);
657  }
658  }
659 
660  cc1++;
661  cc2++;
662 
663  } while (cc1 != ccStart);
664 
665  verticesOnFace.shrink();
666 
667  if (verticesOnFace.size() >= 3 && nUniqueVertices < 3)
668  {
669  // There are not enough unique vertices on this face to
670  // justify its size, it may have a form like:
671 
672  // Vertices:
673  // A B
674  // A B
675 
676  // Face:
677  // ABAB
678 
679  // Setting the size to be below 3, so that it will not be
680  // created
681 
682  verticesOnFace.setSize(nUniqueVertices);
683  }
684 
685  return face(verticesOnFace);
686 }
687 
688 
689 Foam::label Foam::conformalVoronoiMesh::maxFilterCount
690 (
691  const Delaunay::Finite_edges_iterator& eit
692 ) const
693 {
694  Cell_circulator ccStart = incident_cells(*eit);
695  Cell_circulator cc = ccStart;
696 
697  label maxFC = 0;
698 
699  do
700  {
701  if (cc->hasFarPoint())
702  {
703  Cell_handle c = eit->first;
704  Vertex_handle vA = c->vertex(eit->second);
705  Vertex_handle vB = c->vertex(eit->third);
706 
708  << "Dual face uses circumcenter defined by a "
709  << "Delaunay tetrahedron with no internal "
710  << "or boundary points. Defining Delaunay edge ends: "
711  << topoint(vA->point()) << " "
712  << topoint(vB->point()) << nl
713  << exit(FatalError);
714  }
715 
716  if (cc->filterCount() > maxFC)
717  {
718  maxFC = cc->filterCount();
719  }
720 
721  cc++;
722 
723  } while (cc != ccStart);
724 
725  return maxFC;
726 }
727 
728 
729 bool Foam::conformalVoronoiMesh::ownerAndNeighbour
730 (
731  Vertex_handle vA,
732  Vertex_handle vB,
733  label& owner,
734  label& neighbour
735 ) const
736 {
737  bool reverse = false;
738 
739  owner = -1;
740  neighbour = -1;
741 
742  label dualCellIndexA = vA->index();
743 
744  if (!vA->internalOrBoundaryPoint() || vA->referred())
745  {
746  if (!vA->constrained())
747  {
748  dualCellIndexA = -1;
749  }
750  }
751 
752  label dualCellIndexB = vB->index();
753 
754  if (!vB->internalOrBoundaryPoint() || vB->referred())
755  {
756  if (!vB->constrained())
757  {
758  dualCellIndexB = -1;
759  }
760  }
761 
762  if (dualCellIndexA == -1 && dualCellIndexB == -1)
763  {
765  << "Attempting to create a face joining "
766  << "two unindexed dual cells "
767  << exit(FatalError);
768  }
769  else if (dualCellIndexA == -1 || dualCellIndexB == -1)
770  {
771  // boundary face, find which is the owner
772 
773  if (dualCellIndexA == -1)
774  {
775  owner = dualCellIndexB;
776 
777  reverse = true;
778  }
779  else
780  {
781  owner = dualCellIndexA;
782  }
783  }
784  else
785  {
786  // internal face, find the lower cell to be the owner
787 
788  if (dualCellIndexB > dualCellIndexA)
789  {
790  owner = dualCellIndexA;
791  neighbour = dualCellIndexB;
792  }
793  else
794  {
795  owner = dualCellIndexB;
796  neighbour = dualCellIndexA;
797 
798  // reverse face order to correctly orientate normal
799  reverse = true;
800  }
801  }
802 
803  return reverse;
804 }
805 
806 
807 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
808 
809 Foam::conformalVoronoiMesh::conformalVoronoiMesh
810 (
811  const Time& runTime,
812  const dictionary& foamyHexMeshDict,
813  const fileName& decompDictFile
814 )
815 :
816  DistributedDelaunayMesh<Delaunay>(runTime),
817  runTime_(runTime),
818  rndGen_(64293*Pstream::myProcNo()),
819  foamyHexMeshControls_(foamyHexMeshDict),
820  allGeometry_
821  (
822  IOobject
823  (
824  "cvSearchableSurfaces",
825  runTime_.constant(),
826  "triSurface",
827  runTime_,
828  IOobject::MUST_READ,
829  IOobject::NO_WRITE
830  ),
831  foamyHexMeshDict.subDict("geometry"),
832  foamyHexMeshDict.getOrDefault("singleRegionName", true)
833  ),
834  geometryToConformTo_
835  (
836  runTime_,
837  rndGen_,
838  allGeometry_,
839  foamyHexMeshDict.subDict("surfaceConformation")
840  ),
841  decomposition_
842  (
843  Pstream::parRun()
844  ? new backgroundMeshDecomposition
845  (
846  runTime_,
847  rndGen_,
848  geometryToConformTo_,
849  foamyHexMeshControls().foamyHexMeshDict().subDict
850  (
851  "backgroundMeshDecomposition"
852  ),
853  decompDictFile
854  )
855  : nullptr
856  ),
857  cellShapeControl_
858  (
859  runTime_,
860  foamyHexMeshControls_,
861  allGeometry_,
862  geometryToConformTo_
863  ),
864  limitBounds_(),
865  ptPairs_(*this),
866  ftPtConformer_(*this),
867  edgeLocationTreePtr_(),
868  surfacePtLocationTreePtr_(),
869  surfaceConformationVertices_(),
870  initialPointsMethod_
871  (
872  initialPointsMethod::New
873  (
874  foamyHexMeshDict.subDict("initialPoints"),
875  runTime_,
876  rndGen_,
877  geometryToConformTo_,
878  cellShapeControl_,
879  decomposition_
880  )
881  ),
882  relaxationModel_
883  (
884  relaxationModel::New
885  (
886  foamyHexMeshDict.subDict("motionControl"),
887  runTime_
888  )
889  ),
890  faceAreaWeightModel_
891  (
892  faceAreaWeightModel::New
893  (
894  foamyHexMeshDict.subDict("motionControl")
895  )
896  )
897 {}
898 
899 
900 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
901 
903 {}
904 
905 
906 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
907 
909 {
910  if (foamyHexMeshControls().objOutput())
911  {
912  geometryToConformTo_.writeFeatureObj("foamyHexMesh");
913  }
914 
915  buildCellSizeAndAlignmentMesh();
916 
917  insertInitialPoints();
918 
919  insertFeaturePoints(true);
920 
921  setVertexSizeAndAlignment();
922 
923  cellSizeMeshOverlapsBackground();
924 
925  // Improve the guess that the backgroundMeshDecomposition makes with the
926  // initial positions. Use before building the surface conformation to
927  // better balance the surface conformation load.
928  distributeBackground(*this);
929 
930  buildSurfaceConformation();
931 
932  // The introduction of the surface conformation may have distorted the
933  // balance of vertices, distribute if necessary.
934  distributeBackground(*this);
935 
936  if (Pstream::parRun())
937  {
938  sync(decomposition_().procBounds());
939  }
940 
941  // Do not store the surface conformation until after it has been
942  // (potentially) redistributed.
943  storeSurfaceConformation();
944 
945  // Report any Delaunay vertices that do not think that they are in the
946  // domain the processor they are on.
947  // reportProcessorOccupancy();
948 
949  cellSizeMeshOverlapsBackground();
950 
951  if (foamyHexMeshControls().printVertexInfo())
952  {
953  printVertexInfo(Info);
954  }
955 
956  if (foamyHexMeshControls().objOutput())
957  {
959  (
960  time().path()/"internalPoints_" + time().timeName() + ".obj",
961  *this,
964  );
965  }
966 }
967 
968 
970 {
971  if (Pstream::parRun())
972  {
973  decomposition_.reset
974  (
975  new backgroundMeshDecomposition
976  (
977  runTime_,
978  rndGen_,
979  geometryToConformTo_,
980  foamyHexMeshControls().foamyHexMeshDict().subDict
981  (
982  "backgroundMeshDecomposition"
983  )
984  )
985  );
986  }
987 
988  insertInitialPoints();
989 
990  insertFeaturePoints();
991 
992  // Improve the guess that the backgroundMeshDecomposition makes with the
993  // initial positions. Use before building the surface conformation to
994  // better balance the surface conformation load.
995  distributeBackground(*this);
996 
997  buildSurfaceConformation();
998 
999  // The introduction of the surface conformation may have distorted the
1000  // balance of vertices, distribute if necessary.
1001  distributeBackground(*this);
1002 
1003  if (Pstream::parRun())
1004  {
1005  sync(decomposition_().procBounds());
1006  }
1007 
1008  cellSizeMeshOverlapsBackground();
1009 
1010  if (foamyHexMeshControls().printVertexInfo())
1011  {
1012  printVertexInfo(Info);
1013  }
1014 }
1015 
1016 
1018 {
1019  timeCheck("Start of move");
1020 
1021  scalar relaxation = relaxationModel_->relaxation();
1022 
1023  Info<< nl << "Relaxation = " << relaxation << endl;
1024 
1025  pointField dualVertices(number_of_finite_cells());
1026 
1027  this->resetCellCount();
1028 
1029  // Find the dual point of each tetrahedron and assign it an index.
1030  for
1031  (
1032  Delaunay::Finite_cells_iterator cit = finite_cells_begin();
1033  cit != finite_cells_end();
1034  ++cit
1035  )
1036  {
1037  cit->cellIndex() = Cb::ctUnassigned;
1038 
1039  if (cit->anyInternalOrBoundaryDualVertex())
1040  {
1041  cit->cellIndex() = getNewCellIndex();
1042 
1043  dualVertices[cit->cellIndex()] = cit->dual();
1044  }
1045 
1046  if (cit->hasFarPoint())
1047  {
1048  cit->cellIndex() = Cb::ctFar;
1049  }
1050  }
1051 
1052  dualVertices.setSize(cellCount());
1053 
1054  setVertexSizeAndAlignment();
1055 
1056  timeCheck("Determined sizes and alignments");
1057 
1058  Info<< nl << "Determining vertex displacements" << endl;
1059 
1060  vectorField cartesianDirections(3);
1061 
1062  cartesianDirections[0] = vector(1, 0, 0);
1063  cartesianDirections[1] = vector(0, 1, 0);
1064  cartesianDirections[2] = vector(0, 0, 1);
1065 
1066  vectorField displacementAccumulator
1067  (
1068  number_of_vertices(),
1069  Zero
1070  );
1071 
1072  bitSet pointToBeRetained(number_of_vertices(), true);
1073 
1074  DynamicList<Point> pointsToInsert(number_of_vertices());
1075 
1076  for
1077  (
1078  Delaunay::Finite_edges_iterator eit = finite_edges_begin();
1079  eit != finite_edges_end();
1080  ++eit
1081  )
1082  {
1083  Cell_handle c = eit->first;
1084  Vertex_handle vA = c->vertex(eit->second);
1085  Vertex_handle vB = c->vertex(eit->third);
1086 
1087  if
1088  (
1089  (
1090  vA->internalPoint() && !vA->referred()
1091  && vB->internalOrBoundaryPoint()
1092  )
1093  || (
1094  vB->internalPoint() && !vB->referred()
1095  && vA->internalOrBoundaryPoint()
1096  )
1097  )
1098  {
1099  pointFromPoint dVA = topoint(vA->point());
1100  pointFromPoint dVB = topoint(vB->point());
1101 
1102  Field<vector> alignmentDirsA
1103  (
1104  vA->alignment().T() & cartesianDirections
1105  );
1106  Field<vector> alignmentDirsB
1107  (
1108  vB->alignment().T() & cartesianDirections
1109  );
1110 
1111  Field<vector> alignmentDirs(alignmentDirsA);
1112 
1113  forAll(alignmentDirsA, aA)
1114  {
1115  const vector& a = alignmentDirsA[aA];
1116 
1117  scalar maxDotProduct = 0.0;
1118 
1119  forAll(alignmentDirsB, aB)
1120  {
1121  const vector& b = alignmentDirsB[aB];
1122 
1123  const scalar dotProduct = a & b;
1124 
1125  if (mag(dotProduct) > maxDotProduct)
1126  {
1127  maxDotProduct = mag(dotProduct);
1128 
1129  alignmentDirs[aA] = a + sign(dotProduct)*b;
1130 
1131  alignmentDirs[aA].normalise();
1132  }
1133  }
1134  }
1135 
1136  vector rAB = dVA - dVB;
1137 
1138  scalar rABMag = mag(rAB);
1139 
1140  if (rABMag < SMALL)
1141  {
1142  // Removal of close points
1143 
1144  if
1145  (
1146  vA->internalPoint() && !vA->referred() && !vA->fixed()
1147  && vB->internalPoint() && !vB->referred() && !vB->fixed()
1148  )
1149  {
1150  // Only insert a point at the midpoint of
1151  // the short edge if neither attached
1152  // point has already been identified to be
1153  // removed.
1154 
1155  if
1156  (
1157  pointToBeRetained.test(vA->index())
1158  && pointToBeRetained.test(vB->index())
1159  )
1160  {
1161  const Foam::point pt(0.5*(dVA + dVB));
1162 
1163  if (internalPointIsInside(pt))
1164  {
1165  pointsToInsert.append(toPoint(pt));
1166  }
1167  }
1168  }
1169 
1170  if (vA->internalPoint() && !vA->referred() && !vA->fixed())
1171  {
1172  pointToBeRetained.unset(vA->index());
1173  }
1174 
1175  if (vB->internalPoint() && !vB->referred() && !vB->fixed())
1176  {
1177  pointToBeRetained.unset(vB->index());
1178  }
1179 
1180  // Do not consider this Delaunay edge any further
1181 
1182  continue;
1183  }
1184 
1185  forAll(alignmentDirs, aD)
1186  {
1187  vector& alignmentDir = alignmentDirs[aD];
1188 
1189  scalar dotProd = rAB & alignmentDir;
1190 
1191  if (dotProd < 0)
1192  {
1193  // swap the direction of the alignment so that has the
1194  // same sense as rAB
1195  alignmentDir *= -1;
1196  dotProd *= -1;
1197  }
1198 
1199  const scalar alignmentDotProd = dotProd/rABMag;
1200 
1201  if
1202  (
1203  alignmentDotProd
1204  > foamyHexMeshControls().cosAlignmentAcceptanceAngle()
1205  )
1206  {
1207  scalar targetCellSize =
1209 
1210  scalar targetFaceArea = sqr(targetCellSize);
1211 
1212  const vector originalAlignmentDir = alignmentDir;
1213 
1214  // Update cell size and face area
1215  cellShapeControls().aspectRatio().updateCellSizeAndFaceArea
1216  (
1217  alignmentDir,
1218  targetFaceArea,
1219  targetCellSize
1220  );
1221 
1222  // Vector to move end points around middle of vector
1223  // to align edge (i.e. dual face normal) with alignment
1224  // directions.
1225  vector delta = alignmentDir - 0.5*rAB;
1226 
1227  face dualFace = buildDualFace(eit);
1228 
1229 // Pout<< dualFace << endl;
1230 // Pout<< " " << vA->info() << endl;
1231 // Pout<< " " << vB->info() << endl;
1232 
1233  const scalar faceArea = dualFace.mag(dualVertices);
1234 
1235  // Update delta vector
1236  cellShapeControls().aspectRatio().updateDeltaVector
1237  (
1238  originalAlignmentDir,
1239  targetCellSize,
1240  rABMag,
1241  delta
1242  );
1243 
1244 // if (targetFaceArea == 0)
1245 // {
1246 // Pout<< vA->info() << vB->info();
1247 //
1248 // Cell_handle ch = locate(vA->point());
1249 // if (is_infinite(ch))
1250 // {
1251 // Pout<< "vA " << vA->targetCellSize() << endl;
1252 // }
1253 //
1254 // ch = locate(vB->point());
1255 // if (is_infinite(ch))
1256 // {
1257 // Pout<< "vB " << vB->targetCellSize() << endl;
1258 // }
1259 // }
1260 
1261  delta *= faceAreaWeightModel_->faceAreaWeight
1262  (
1263  faceArea/targetFaceArea
1264  );
1265 
1266  if
1267  (
1268  (
1269  (vA->internalPoint() && vB->internalPoint())
1270  && (!vA->referred() || !vB->referred())
1271 // ||
1272 // (
1273 // vA->referredInternalPoint()
1274 // && vB->referredInternalPoint()
1275 // )
1276  )
1277  && rABMag
1278  > foamyHexMeshControls().insertionDistCoeff()
1279  *targetCellSize
1280  && faceArea
1281  > foamyHexMeshControls().faceAreaRatioCoeff()
1282  *targetFaceArea
1283  && alignmentDotProd
1284  > foamyHexMeshControls().cosInsertionAcceptanceAngle()
1285  )
1286  {
1287  // Point insertion
1288  if
1289  (
1290  !geometryToConformTo_.findSurfaceAnyIntersection
1291  (
1292  dVA,
1293  dVB
1294  )
1295  )
1296  {
1297  const Foam::point newPt(0.5*(dVA + dVB));
1298 
1299  // Prevent insertions spanning surfaces
1300  if (internalPointIsInside(newPt))
1301  {
1302  if (Pstream::parRun())
1303  {
1304  if
1305  (
1306  decomposition().positionOnThisProcessor
1307  (
1308  newPt
1309  )
1310  )
1311  {
1312  pointsToInsert.append(toPoint(newPt));
1313  }
1314  }
1315  else
1316  {
1317  pointsToInsert.append(toPoint(newPt));
1318  }
1319  }
1320  }
1321  }
1322  else if
1323  (
1324  (
1325  (vA->internalPoint() && !vA->referred())
1326  || (vB->internalPoint() && !vB->referred())
1327  )
1328  && rABMag
1329  < foamyHexMeshControls().removalDistCoeff()
1330  *targetCellSize
1331  )
1332  {
1333  // Point removal
1334  if
1335  (
1336  (
1337  vA->internalPoint()
1338  && !vA->referred()
1339  && !vA->fixed()
1340  )
1341  &&
1342  (
1343  vB->internalPoint()
1344  && !vB->referred()
1345  && !vB->fixed()
1346  )
1347  )
1348  {
1349  // Only insert a point at the midpoint of
1350  // the short edge if neither attached
1351  // point has already been identified to be
1352  // removed.
1353  if
1354  (
1355  pointToBeRetained.test(vA->index())
1356  && pointToBeRetained.test(vB->index())
1357  )
1358  {
1359  const Foam::point pt(0.5*(dVA + dVB));
1360 
1361  if (internalPointIsInside(pt))
1362  {
1363  pointsToInsert.append(toPoint(pt));
1364  }
1365  }
1366  }
1367 
1368  if
1369  (
1370  vA->internalPoint()
1371  && !vA->referred()
1372  && !vA->fixed()
1373  )
1374  {
1375  pointToBeRetained.unset(vA->index());
1376  }
1377 
1378  if
1379  (
1380  vB->internalPoint()
1381  && !vB->referred()
1382  && !vB->fixed()
1383  )
1384  {
1385  pointToBeRetained.unset(vB->index());
1386  }
1387  }
1388  else
1389  {
1390  if
1391  (
1392  vA->internalPoint()
1393  && !vA->referred()
1394  && !vA->fixed()
1395  )
1396  {
1397  if (vB->fixed())
1398  {
1399  displacementAccumulator[vA->index()] += 2*delta;
1400  }
1401  else
1402  {
1403  displacementAccumulator[vA->index()] += delta;
1404  }
1405  }
1406 
1407  if
1408  (
1409  vB->internalPoint()
1410  && !vB->referred()
1411  && !vB->fixed()
1412  )
1413  {
1414  if (vA->fixed())
1415  {
1416  displacementAccumulator[vB->index()] -= 2*delta;
1417  }
1418  else
1419  {
1420  displacementAccumulator[vB->index()] -= delta;
1421  }
1422  }
1423  }
1424  }
1425  }
1426  }
1427  }
1428 
1429  Info<< "Limit displacements" << endl;
1430 
1431  // Limit displacements that pierce, or get too close to the surface
1432  for
1433  (
1434  Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
1435  vit != finite_vertices_end();
1436  ++vit
1437  )
1438  {
1439  if (vit->internalPoint() && !vit->referred() && !vit->fixed())
1440  {
1441  if (pointToBeRetained.test(vit->index()))
1442  {
1443  limitDisplacement
1444  (
1445  vit,
1446  displacementAccumulator[vit->index()]
1447  );
1448  }
1449  }
1450  }
1451 
1452  vector totalDisp = gSum(displacementAccumulator);
1453  scalar totalDist = gSum(mag(displacementAccumulator));
1454 
1455  displacementAccumulator *= relaxation;
1456 
1457  Info<< "Sum displacements" << endl;
1458 
1459  label nPointsToRetain = 0;
1460  label nPointsToRemove = 0;
1461 
1462  for
1463  (
1464  Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
1465  vit != finite_vertices_end();
1466  ++vit
1467  )
1468  {
1469  if (vit->internalPoint() && !vit->referred() && !vit->fixed())
1470  {
1471  if (pointToBeRetained.test(vit->index()))
1472  {
1473  // Convert vit->point() to FOAM vector (double) to do addition,
1474  // avoids memory increase because a record of the constructions
1475  // would be kept otherwise.
1476  // See [email protected]:
1477  // "Memory issue with openSUSE 11.3, exact kernel, adding
1478  // points/vectors"
1479  // 14/1/2011.
1480  // Only necessary if using an exact constructions kernel
1481  // (extended precision)
1482  Foam::point pt
1483  (
1484  topoint(vit->point())
1485  + displacementAccumulator[vit->index()]
1486  );
1487 
1488  if (internalPointIsInside(pt))
1489  {
1490  pointsToInsert.append(toPoint(pt));
1491  nPointsToRemove++;
1492  }
1493 
1494  nPointsToRetain++;
1495  }
1496  }
1497  }
1498 
1499  pointsToInsert.shrink();
1500 
1501  Info<< returnReduce
1502  (
1503  nPointsToRetain - nPointsToRemove,
1504  sumOp<label>()
1505  )
1506  << " internal points are outside the domain. "
1507  << "They will not be inserted." << endl;
1508 
1509  // Save displacements to file.
1510  if (foamyHexMeshControls().objOutput() && time().writeTime())
1511  {
1512  Info<< "Writing point displacement vectors to file." << endl;
1513  OFstream str
1514  (
1515  time().path()/"displacements_" + runTime_.timeName() + ".obj"
1516  );
1517 
1518  for
1519  (
1520  Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
1521  vit != finite_vertices_end();
1522  ++vit
1523  )
1524  {
1525  if (vit->internalPoint() && !vit->referred())
1526  {
1527  if (pointToBeRetained.test(vit->index()))
1528  {
1529  meshTools::writeOBJ(str, topoint(vit->point()));
1530 
1531  str << "vn "
1532  << displacementAccumulator[vit->index()][0] << " "
1533  << displacementAccumulator[vit->index()][1] << " "
1534  << displacementAccumulator[vit->index()][2] << " "
1535  << endl;
1536  }
1537  }
1538  }
1539  }
1540 
1541  // Remove the entire tessellation
1543 
1544  timeCheck("Displacement calculated");
1545 
1546  Info<< nl<< "Inserting displaced tessellation" << endl;
1547 
1548  insertFeaturePoints(true);
1549 
1550  insertInternalPoints(pointsToInsert, true);
1551 
1552  // Fix points that have not been significantly displaced
1553 // for
1554 // (
1555 // Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
1556 // vit != finite_vertices_end();
1557 // ++vit
1558 // )
1559 // {
1560 // if (vit->internalPoint())
1561 // {
1562 // if
1563 // (
1564 // mag(displacementAccumulator[vit->index()])
1565 // < 0.1*targetCellSize(topoint(vit->point()))
1566 // )
1567 // {
1568 // vit->setVertexFixed();
1569 // }
1570 // }
1571 // }
1572 
1573  timeCheck("Internal points inserted");
1574 
1575  {
1576  // Check that no index is shared between any of the local points
1577  labelHashSet usedIndices;
1578  for
1579  (
1580  Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
1581  vit != finite_vertices_end();
1582  ++vit
1583  )
1584  {
1585  if (!vit->referred() && !usedIndices.insert(vit->index()))
1586  {
1588  << "Index already used! Could not insert: " << nl
1589  << vit->info()
1590  << abort(FatalError);
1591  }
1592  }
1593  }
1594 
1595  conformToSurface();
1596 
1597  if (foamyHexMeshControls().objOutput())
1598  {
1600  (
1601  time().path()/"internalPoints_" + time().timeName() + ".obj",
1602  *this,
1604  );
1605 
1606  if (reconformToSurface())
1607  {
1609  (
1610  time().path()/"boundaryPoints_" + time().timeName() + ".obj",
1611  *this
1612  );
1613 
1615  (
1616  time().path()/"internalBoundaryPoints_" + time().timeName()
1617  + ".obj",
1618  *this,
1621  );
1622 
1624  (
1625  time().path()/"externalBoundaryPoints_" + time().timeName()
1626  + ".obj",
1627  *this,
1630  );
1631 
1632  OBJstream multipleIntersections
1633  (
1634  "multipleIntersections_"
1635  + time().timeName()
1636  + ".obj"
1637  );
1638 
1639  for
1640  (
1641  Delaunay::Finite_edges_iterator eit = finite_edges_begin();
1642  eit != finite_edges_end();
1643  ++eit
1644  )
1645  {
1646  Cell_handle c = eit->first;
1647  Vertex_handle vA = c->vertex(eit->second);
1648  Vertex_handle vB = c->vertex(eit->third);
1649 
1650  Foam::point ptA(topoint(vA->point()));
1651  Foam::point ptB(topoint(vB->point()));
1652 
1653  List<pointIndexHit> surfHits;
1654  labelList hitSurfaces;
1655 
1656  geometryToConformTo().findSurfaceAllIntersections
1657  (
1658  ptA,
1659  ptB,
1660  surfHits,
1661  hitSurfaces
1662  );
1663 
1664  if
1665  (
1666  surfHits.size() >= 2
1667  || (
1668  surfHits.size() == 0
1669  && (
1670  (vA->externalBoundaryPoint()
1671  && vB->internalBoundaryPoint())
1672  || (vB->externalBoundaryPoint()
1673  && vA->internalBoundaryPoint())
1674  )
1675  )
1676  )
1677  {
1678  multipleIntersections.writeLine(ptA, ptB);
1679  }
1680  }
1681  }
1682  }
1683 
1684  timeCheck("After conformToSurface");
1685 
1686  if (foamyHexMeshControls().printVertexInfo())
1687  {
1688  printVertexInfo(Info);
1689  }
1690 
1691  if (time().writeTime())
1692  {
1693  writeMesh(time().timeName());
1694  }
1695 
1696  Info<< nl
1697  << "Total displacement = " << totalDisp << nl
1698  << "Total distance = " << totalDist << nl
1699  << endl;
1700 }
1701 
1702 
1703 void Foam::conformalVoronoiMesh::checkCoPlanarCells() const
1704 {
1705  typedef CGAL::Exact_predicates_exact_constructions_kernel Kexact;
1706  typedef CGAL::Point_3<Kexact> PointExact;
1707 
1708  if (!is_valid())
1709  {
1710  Pout<< "Triangulation is invalid!" << endl;
1711  }
1712 
1713  OFstream str("badCells.obj");
1714 
1715  label badCells = 0;
1716 
1717  for
1718  (
1719  Delaunay::Finite_cells_iterator cit = finite_cells_begin();
1720  cit != finite_cells_end();
1721  ++cit
1722  )
1723  {
1724  const scalar quality = foamyHexMeshChecks::coplanarTet(cit, 1e-16);
1725 
1726  if (quality == 0)
1727  {
1728  Pout<< "COPLANAR: " << cit->info() << nl
1729  << " quality = " << quality << nl
1730  << " dual = " << topoint(cit->dual()) << endl;
1731 
1732  DelaunayMeshTools::drawDelaunayCell(str, cit, badCells++);
1733 
1734  FixedList<PointExact, 4> cellVerticesExact(PointExact(0,0,0));
1735  forAll(cellVerticesExact, vI)
1736  {
1737  cellVerticesExact[vI] = PointExact
1738  (
1739  cit->vertex(vI)->point().x(),
1740  cit->vertex(vI)->point().y(),
1741  cit->vertex(vI)->point().z()
1742  );
1743  }
1744 
1745  PointExact synchronisedDual = CGAL::circumcenter<Kexact>
1746  (
1747  cellVerticesExact[0],
1748  cellVerticesExact[1],
1749  cellVerticesExact[2],
1750  cellVerticesExact[3]
1751  );
1752 
1753  Foam::point exactPt
1754  (
1755  CGAL::to_double(synchronisedDual.x()),
1756  CGAL::to_double(synchronisedDual.y()),
1757  CGAL::to_double(synchronisedDual.z())
1758  );
1759 
1760  Info<< "inexact = " << cit->dual() << nl
1761  << "exact = " << exactPt << endl;
1762  }
1763  }
1764 
1765  Pout<< "There are " << badCells << " bad cells out of "
1766  << number_of_finite_cells() << endl;
1767 
1768 
1769  label nNonGabriel = 0;
1770  for
1771  (
1772  Delaunay::Finite_facets_iterator fit = finite_facets_begin();
1773  fit != finite_facets_end();
1774  ++fit
1775  )
1776  {
1777  if (!is_Gabriel(*fit))
1778  {
1779  nNonGabriel++;//Pout<< "Non-gabriel face" << endl;
1780  }
1781  }
1782 
1783  Pout<< "There are " << nNonGabriel << " non-Gabriel faces out of "
1784  << number_of_finite_facets() << endl;
1785 }
1786 
1787 
1788 // ************************************************************************* //
scalar delta
dimensionedScalar sign(const dimensionedScalar &ds)
CGAL::Delaunay_triangulation_3< K, Tds, CompactLocator > Delaunay
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
~conformalVoronoiMesh()
Destructor.
pointFromPoint topoint(const Point &P)
void writeOBJ(const fileName &fName, const List< Foam::point > &points)
Write list of points to file.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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:600
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:526
dimensionedSymmTensor sqr(const dimensionedVector &dv)
srcOptions insert("case", fileName(rootDirSource/caseDirSource))
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:529
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
Definition: pointIndexHit.H:58
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1586
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
bool distribute(const boundBox &bb)
InfoProxy< IOstream > info() const noexcept
Return info proxy, used to print IOstream information to a stream.
Definition: IOstream.H:532
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:196
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition: UPstream.H:1611
void drawDelaunayCell(Ostream &os, const CellHandle &c, label offset=0)
Draws a tet cell to an output stream. The offset is supplied as the tet.
void reduce(T &value, [[maybe_unused]] BinaryOp bop, [[maybe_unused]] const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce)
An indexed form of CGAL::Triangulation_vertex_base_3<K> used to keep track of the Delaunay vertices i...
Definition: indexedVertex.H:59
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:286
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
void move()
Move the vertices according to the controller, re-conforming to the.
pointField vertices(const blockVertexList &bvl)
word timeName
Definition: getTimeIndex.H:3
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
bool uninitialised(const VertexType &v)
Type gSum(const FieldField< Field, Type > &f)
const pointField & points
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
label nPoints
bool returnReduceAnd(const bool value, const int communicator=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
sideVolumeType
Normals point to the outside.
void reset()
Clear the entire triangulation.
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
PointFrompoint toPoint(const Foam::point &p)
void reverse(UList< T > &list, const label n)
Reverse the first n elements of the list.
Definition: UListI.H:514
Vector< scalar > vector
Definition: vector.H:57
errorManip< error > abort(error &err)
Definition: errorManip.H:139
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:385
const wordList surface
Standard surface field types (scalar, vector, tensor, etc)
Istream and Ostream manipulators taking arguments.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
defineTypeNameAndDebug(combustionModel, 0)
static const Enum< dualMeshPointType > dualMeshPointTypeNames_
void writeBoundaryPoints(const fileName &fName, const Triangulation &t)
Write the boundary Delaunay points to an OBJ file.
const wordList internal
Standard dimensioned field types (scalar, vector, tensor, etc)
Foam::scalar averageCellSize(const VertexType &vA, const VertexType &vB)
Return the target cell size from that stored on a pair of Delaunay vertices,.
#define WarningInFunction
Report a warning using Foam::Warning.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
const dimensionedScalar c
Speed of light in a vacuum.
messageStream Info
Information stream (stdout output on master, null elsewhere)
cellShapeControlMesh & shapeControlMesh()
Field< vector > vectorField
Specialisation of Field<T> for vector.
scalar coplanarTet(Cell &c, const scalar tol=1e-12)
Map< label > insertPoints(const List< Vb > &vertices, const bool reIndex)
Insert the list of vertices (calls rangeInsertWithInfo)
iterator end() noexcept
Return an iterator to end traversing the UList.
Definition: UListI.H:429
List< label > labelList
A List of labels.
Definition: List.H:61
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Namespace for OpenFOAM.
A HashTable to objects of type <T> with a label key.
const pointField & pts
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127