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