conformalVoronoiMesh.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::conformalVoronoiMesh
26 
27 Description
28 
29 SourceFiles
30  conformalVoronoiMeshI.H
31  conformalVoronoiMesh.C
32  conformalVoronoiMeshZones.C
33  conformalVoronoiMeshIO.C
34  conformalVoronoiMeshConformToSurface.C
35  conformalVoronoiMeshFeaturePoints.C
36  conformalVoronoiMeshCalcDualMesh.C
37  conformalVoronoiMeshTemplates.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef conformalVoronoiMesh_H
42 #define conformalVoronoiMesh_H
43 
45 #include "searchableSurfaces.H"
46 #include "conformationSurfaces.H"
47 #include "cellShapeControl.H"
48 #include "cvControls.H"
49 #include "DynamicList.H"
50 #include "PackedBoolList.H"
51 #include "Time.H"
52 #include "polyMesh.H"
53 #include "plane.H"
54 #include "SortableList.H"
55 #include "meshTools.H"
56 #include "dynamicIndexedOctree.H"
57 #include "dynamicTreeDataPoint.H"
58 #include "indexedOctree.H"
59 #include "treeDataPoint.H"
60 #include "unitConversion.H"
61 #include "transform.H"
62 #include "volFields.H"
63 #include "fvMesh.H"
64 #include "labelPair.H"
65 #include "HashSet.H"
66 #include "memInfo.H"
67 #include "point.H"
68 #include "cellSet.H"
69 #include "wallPolyPatch.H"
70 #include "processorPolyPatch.H"
72 #include "globalIndex.H"
73 #include "pointFeatureEdgesTypes.H"
74 #include "pointConversion.H"
75 #include "Pair.H"
77 #include "featurePointConformer.H"
78 #include "pointPairs.H"
79 
80 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
81 
82 namespace Foam
83 {
84 
85 // Forward declaration of classes
86 class initialPointsMethod;
87 class relaxationModel;
88 class faceAreaWeightModel;
89 class backgroundMeshDecomposition;
90 class OBJstream;
91 
92 /*---------------------------------------------------------------------------*\
93  Class conformalVoronoiMesh Declaration
94 \*---------------------------------------------------------------------------*/
95 
97 :
98  public DistributedDelaunayMesh<Delaunay>
99 {
100 public:
102  typedef Delaunay::Vertex_handle Vertex_handle;
103  typedef Delaunay::Cell_handle Cell_handle;
104  typedef Delaunay::Edge Edge;
105  typedef Delaunay::Facet Facet;
106  typedef Delaunay::Point Point;
113 
114  // Static data
116  enum dualMeshPointType
117  {
118  internal = 0,
119  surface = 1,
122  constrained = 4
123  };
126 
127 
128 private:
129 
130  // Static data
131 
132  static const scalar searchConeAngle;
133 
134  static const scalar searchAngleOppositeSurface;
135 
136 
137  // Private data
138 
139  //- The time registry of the application
140  const Time& runTime_;
141 
142  //- Random number generator
143  mutable Random rndGen_;
144 
145  //- Controls for the conformal Voronoi meshing process
146  cvControls foamyHexMeshControls_;
147 
148  //- All geometry of the meshing process, including surfaces to be
149  // conformed to and those to be used for refinement
150  searchableSurfaces allGeometry_;
151 
152  //- The surfaces to conform to
153  conformationSurfaces geometryToConformTo_;
154 
155  //- Background mesh decomposition, only available in parallel.
157 
158  //- The cell shape control object
159  cellShapeControl cellShapeControl_;
160 
161  //- Limiting bound box before infinity begins
162  treeBoundBox limitBounds_;
163 
164  //-
165  mutable pointPairs<Delaunay> ptPairs_;
166 
167  //-
168  featurePointConformer ftPtConformer_;
169 
170  //- Search tree for edge point locations
172  edgeLocationTreePtr_;
173 
174  mutable DynamicList<Foam::point> existingEdgeLocations_;
175 
176  //- Search tree for surface point locations
178  surfacePtLocationTreePtr_;
179 
180  mutable DynamicList<Foam::point> existingSurfacePtLocations_;
181 
182  //- Store the surface and feature edge conformation locations to be
183  // reinserted
184  List<Vb> surfaceConformationVertices_;
185 
186  //- Method for inserting initial points. Runtime selectable.
187  autoPtr<initialPointsMethod> initialPointsMethod_;
188 
189  //- Relaxation coefficient model. Runtime selectable.
190  autoPtr<relaxationModel> relaxationModel_;
191 
192  //- Face area weight function. Runtime selectable.
193  autoPtr<faceAreaWeightModel> faceAreaWeightModel_;
194 
195 
196  // Private Member Functions
197 
198  inline scalar defaultCellSize() const;
199 
200  //- Return the local target cell size at the given location. Takes
201  // boolean argument to allow speed-up of queries if the point is going
202  // to be on a surface.
203  inline scalar targetCellSize(const Foam::point& pt) const;
204 
205  //- Return the target cell size from that stored on a pair of
206  // Delaunay vertices, including the possibility that one of
207  // them is not an internalOrBoundaryPoint, and so will not
208  // have valid data.
209  inline scalar averageAnyCellSize
210  (
211  const Vertex_handle& vA,
212  const Vertex_handle& vB
213  ) const;
214 
215  //- The average target cell size of a Delaunay facet, i.e., of
216  // a dual edge
217  inline scalar averageAnyCellSize
218  (
219  const Delaunay::Finite_facets_iterator& fit
220  ) const;
221 
222  //- Insert Delaunay vertices using the CGAL range insertion method,
223  // optionally check processor occupancy and distribute to other
224  // processors
225  void insertInternalPoints
226  (
228  const bool distribute = false
229  );
230 
231  Map<label> insertPointPairs
232  (
234  bool distribute,
235  bool reIndex
236  );
237 
238  //- Create a point-pair at a ppDist distance either side of
239  // surface point surfPt, in the direction n
240  inline void createPointPair
241  (
242  const scalar ppDist,
243  const Foam::point& surfPt,
244  const vector& n,
245  const bool ptPair,
246  DynamicList<Vb>& pts
247  ) const;
248 
249  inline Foam::point perturbPoint(const Foam::point& pt) const;
250 
251  //- Create a point-pair at a ppDist distance either side of
252  // surface point surfPt, in the direction n
253  inline void createBafflePointPair
254  (
255  const scalar ppDist,
256  const Foam::point& surfPt,
257  const vector& n,
258  const bool ptPair,
259  DynamicList<Vb>& pts
260  ) const;
261 
262  //- Check internal point is completely inside the meshable region
263  inline bool internalPointIsInside(const Foam::point& pt) const;
264 
265  //- Insert pairs of points on the surface with the given normals, at the
266  // specified spacing
267  void insertSurfacePointPairs
268  (
269  const pointIndexHitAndFeatureList& surfaceHits,
270  const fileName fName,
271  DynamicList<Vb>& pts
272  );
273 
274  //- Insert groups of points to conform to an edge given a list of
275  // pointIndexHits specifying the location and edge index of the point
276  // to be conformed to on the corresponding entry in featureHit
277  void insertEdgePointGroups
278  (
279  const pointIndexHitAndFeatureList& edgeHits,
280  const fileName fName,
281  DynamicList<Vb>& pts
282  );
283 
284  void createEdgePointGroupByCirculating
285  (
286  const extendedFeatureEdgeMesh& feMesh,
287  const pointIndexHit& edHit,
288  DynamicList<Vb>& pts
289  ) const;
290 
291  bool meshableRegion
292  (
293  const plane::side side,
295  ) const;
296 
297  bool regionIsInside
298  (
300  const vector& normalA,
302  const vector& normalB,
303  const vector& masterPtVec
304  ) const;
305 
306  //- Create points to conform to an external edge
307  void createExternalEdgePointGroup
308  (
309  const extendedFeatureEdgeMesh& feMesh,
310  const pointIndexHit& edHit,
311  DynamicList<Vb>& pts
312  ) const;
313 
314  //- Create points to conform to an internal edge
315  void createInternalEdgePointGroup
316  (
317  const extendedFeatureEdgeMesh& feMesh,
318  const pointIndexHit& edHit,
319  DynamicList<Vb>& pts
320  ) const;
321 
322  //- Create points to conform to a flat edge
323  void createFlatEdgePointGroup
324  (
325  const extendedFeatureEdgeMesh& feMesh,
326  const pointIndexHit& edHit,
327  DynamicList<Vb>& pts
328  ) const;
329 
330  //- Create points to conform to an open edge
331  void createOpenEdgePointGroup
332  (
333  const extendedFeatureEdgeMesh& feMesh,
334  const pointIndexHit& edHit,
335  DynamicList<Vb>& pts
336  ) const;
337 
338  //- Create points to conform to multiply connected edge
339  void createMultipleEdgePointGroup
340  (
341  const extendedFeatureEdgeMesh& feMesh,
342  const pointIndexHit& edHit,
343  DynamicList<Vb>& pts
344  ) const;
345 
346  //- Determine and insert point groups at the feature points
347  void insertFeaturePoints(bool distribute = false);
348 
349  //- Check if a location is in exclusion range around a feature point
350  bool nearFeaturePt(const Foam::point& pt) const;
351 
352  //- Check if a surface point is in exclusion range around a feature edge
353  bool surfacePtNearFeatureEdge(const Foam::point& pt) const;
354 
355  //- Insert the initial points into the triangulation, based on the
356  // initialPointsMethod
357  void insertInitialPoints();
358 
359  //- In parallel redistribute the backgroundMeshDecomposition and
360  // vertices to balance the number of vertices on each processor.
361  // Returns true if the background mesh changes as this removes all
362  // referred vertices, so the parallel interface may need rebuilt.
363  template<class Triangulation>
364  bool distributeBackground(const Triangulation& mesh);
365 
366  // Test for full containment
367  void cellSizeMeshOverlapsBackground() const;
368 
369  //-
370  void distribute();
371 
372  void buildCellSizeAndAlignmentMesh();
373 
374  //- Set the size and alignment data for each vertex
375  void setVertexSizeAndAlignment();
376 
377  //- Builds a dual face by circulating around the supplied edge.
378  face buildDualFace
379  (
380  const Delaunay::Finite_edges_iterator& eit
381  ) const;
382 
383  boolList dualFaceBoundaryPoints
384  (
385  const Delaunay::Finite_edges_iterator& eit
386  ) const;
387 
388  //- Finds the maximum filterCount of the dual vertices
389  // (Delaunay cells) that form the dual face produced by the
390  // supplied edge
391  label maxFilterCount
392  (
393  const Delaunay::Finite_edges_iterator& eit
394  ) const;
395 
396  //- Determines the owner and neighbour labels for dual cells
397  // corresponding to the dual face formed by the supplied
398  // Delaunay vertices. If the dual face is a boundary face
399  // then neighbour = -1. Returns true if the dual face
400  // created by vA -> vB needs to be reversed to be correctly
401  // orientated.
402  bool ownerAndNeighbour
403  (
404  Vertex_handle vA,
405  Vertex_handle vB,
406  label& owner,
407  label& neighbour
408  ) const;
409 
410  //- Insert the necessary point pairs to conform to the surface, either
411  // from stored results, or trigger a re-conformation
412  void conformToSurface();
413 
414  //- Decision making function for when to rebuild the surface
415  // conformation
416  bool reconformToSurface() const;
417 
418  //- Determines geometrically whether a vertex is close to a surface
419  // This is an optimisation
420  label findVerticesNearBoundaries();
421 
422  //- Create and insert the necessary point pairs to conform to the
423  // surface, then store the result
424  void buildSurfaceConformation();
425 
426  label synchroniseEdgeTrees
427  (
428  const DynamicList<label>& edgeToTreeShape,
429  pointIndexHitAndFeatureList& featureEdgeHits
430  );
431 
432  label synchroniseSurfaceTrees
433  (
434  const DynamicList<label>& surfaceToTreeShape,
435  pointIndexHitAndFeatureList& surfaceHits
436  );
437 
438  bool surfaceLocationConformsToInside
439  (
440  const pointIndexHitAndFeature& info
441  ) const;
442 
443  //- Check to see if dual cell specified by given vertex iterator
444  // intersects the boundary and hence reqires a point-pair
445  bool dualCellSurfaceAnyIntersection
446  (
447  const Delaunay::Finite_vertices_iterator& vit
448  ) const;
449 
450  //- Return all intersections
451  bool dualCellSurfaceAllIntersections
452  (
453  const Delaunay::Finite_vertices_iterator& vit,
454  pointIndexHitAndFeatureDynList& info
455  ) const;
456 
457  //- Return false if the line is entirely outside the current processor
458  // domain, true is either point is inside, or the processor domain
459  // bounadry is intersected (i.e. the points are box outside but the
460  // line cuts. The points will be moved onto the box where they
461  // intersect.
462  bool clipLineToProc
463  (
464  const Foam::point& pt,
465  Foam::point& a,
466  Foam::point& b
467  ) const;
468 
469  //- Find the "worst" protrusion of a dual cell through the surface,
470  // subject to the maxSurfaceProtrusion tolerance
471  void dualCellLargestSurfaceProtrusion
472  (
473  const Delaunay::Finite_vertices_iterator& vit,
474  pointIndexHit& surfHit,
475  label& hitSurface
476  ) const;
477 
478  void dualCellLargestSurfaceIncursion
479  (
480  const Delaunay::Finite_vertices_iterator& vit,
481  pointIndexHit& surfHit,
482  label& hitSurface
483  ) const;
484 
485  //- Write out vertex-processor occupancy information for debugging
486  void reportProcessorOccupancy();
487 
488  //- Write out debugging information about the surface conformation
489  // quality
490 // void reportSurfaceConformationQuality();
491 
492  //- Limit the displacement of a point so that it doesn't penetrate the
493  // surface to be meshed or come too close to it
494  void limitDisplacement
495  (
496  const Delaunay::Finite_vertices_iterator& vit,
497  vector& displacement,
498  label callCount = 0
499  ) const;
500 
501  //- Find angle between the normals of two close surface points.
502  scalar angleBetweenSurfacePoints(Foam::point pA, Foam::point pB) const;
503 
504  //- Check if a surface point is near another.
505  bool nearSurfacePoint
506  (
507  pointIndexHitAndFeature& pHit
508  ) const;
509 
510  //- Append a point to the surface point tree and the existing list
511  bool appendToSurfacePtTree
512  (
513  const Foam::point& pt
514  ) const;
515 
516  //- Append a point to the edge location tree and the existing list
517  bool appendToEdgeLocationTree
518  (
519  const Foam::point& pt
520  ) const;
521 
522  //- Return a list of the nearest feature edge locations
523  List<pointIndexHit> nearestFeatureEdgeLocations
524  (
525  const Foam::point& pt
526  ) const;
527 
528  //- Check if a point is near any feature edge points.
529  bool pointIsNearFeatureEdgeLocation(const Foam::point& pt) const;
530 
531  bool pointIsNearFeatureEdgeLocation
532  (
533  const Foam::point& pt,
534  pointIndexHit& info
535  ) const;
536 
537  //- Check if a point is near any surface conformation points.
538  bool pointIsNearSurfaceLocation(const Foam::point& pt) const;
539 
540  bool pointIsNearSurfaceLocation
541  (
542  const Foam::point& pt,
543  pointIndexHit& info
544  ) const;
545 
546  //- Check if a location is in the exclusion range of an existing feature
547  //- Edge conformation location
548  bool nearFeatureEdgeLocation
549  (
550  const pointIndexHit& pHit,
551  pointIndexHit& nearestEdgeHit
552  ) const;
553 
554  //- Build or rebuild the edge location tree
555  void buildEdgeLocationTree
556  (
557  const DynamicList<Foam::point>& existingEdgeLocations
558  ) const;
559 
560  //- Build or rebuild the surface point location tree
561  void buildSurfacePtLocationTree
562  (
563  const DynamicList<Foam::point>& existingSurfacePtLocations
564  ) const;
565 
566  //- Process the surface conformation locations to decide which surface
567  // and edge conformation locations to add
568  void addSurfaceAndEdgeHits
569  (
570  const Foam::point& vit,
571  const pointIndexHitAndFeatureDynList& surfaceIntersections,
572  scalar surfacePtReplaceDistCoeffSqr,
573  scalar edgeSearchDistCoeffSqr,
574  pointIndexHitAndFeatureDynList& surfaceHits,
575  pointIndexHitAndFeatureDynList& featureEdgeHits,
576  DynamicList<label>& surfaceToTreeShape,
577  DynamicList<label>& edgeToTreeShape,
578  Map<scalar>& surfacePtToEdgePtDist,
579  bool firstPass
580  ) const;
581 
582  //- Store the surface conformation with the indices offset to be
583  // relative to zero
584  void storeSurfaceConformation();
585 
586  //- Reinsert the surface conformation re-offsetting indices to be
587  // relative to new number of internal vertices
588  void reinsertSurfaceConformation();
589 
590  void checkCells();
591 
592  void checkDuals();
593 
594  void checkVertices();
595 
596  void checkCoPlanarCells() const;
597 
598  //- Dual calculation
599  void calcDualMesh
600  (
602  labelList& boundaryPts,
603  faceList& faces,
604  labelList& owner,
605  labelList& neighbour,
608  pointField& cellCentres,
609  labelList& cellToDelaunayVertex,
610  labelListList& patchToDelaunayVertex,
611  PackedBoolList& boundaryFacesToRemove
612  );
613 
614  void calcNeighbourCellCentres
615  (
616  const polyMesh& mesh,
617  const pointField& cellCentres,
618  pointField& neiCc
619  ) const;
620 
621  void selectSeparatedCoupledFaces
622  (
623  const polyMesh& mesh,
624  boolList& selected
625  ) const;
626 
627  //- From meshRefinementBaffles.C. Use insidePoint for a surface to
628  // determine the cell zone.
629  void findCellZoneInsideWalk
630  (
631  const polyMesh& mesh,
632  const labelList& locationSurfaces,
633  const labelList& faceToSurface,
634  labelList& cellToSurface
635  ) const;
636 
637  //- Calculate the cell zones from cellCentres using all closed surfaces
638  labelList calcCellZones(const pointField& cellCentres) const;
639 
640  //- Calculate the face zones
641  void calcFaceZones
642  (
643  const polyMesh& mesh,
644  const pointField& cellCentres,
645  const labelList& cellToSurface,
646  labelList& faceToSurface,
647  boolList& flipMap
648  ) const;
649 
650  //- Add zones to the polyMesh
651  void addZones(polyMesh& mesh, const pointField& cellCentres) const;
652 
653  //- Tet mesh calculation
654  void calcTetMesh
655  (
657  labelList& pointToDelaunayVertex,
658  faceList& faces,
659  labelList& owner,
660  labelList& neighbour,
663  );
664 
665  //- Determines if the dual face constructed by the Delaunay
666  // edge is a boundary face
667  inline bool isBoundaryDualFace
668  (
669  const Delaunay::Finite_edges_iterator& eit
670  ) const;
671 
672  //- Which processors are attached to the dual edge represented by this
673  // Delaunay facet
674  inline List<label> processorsAttached
675  (
676  const Delaunay::Finite_facets_iterator& fit
677  ) const;
678 
679  //- Determines if the edge constructed from the face is on
680  // a processor patch
681  inline bool isParallelDualEdge
682  (
683  const Delaunay::Finite_facets_iterator& fit
684  ) const;
685 
686  //- Determines if the dual face constructed by the Delaunay
687  // edge is a processor boundary face
688  inline bool isProcBoundaryEdge
689  (
690  const Delaunay::Finite_edges_iterator& eit
691  ) const;
692 
693  //- Merge vertices that are identical
694  void mergeIdenticalDualVertices
695  (
696  const pointField& pts,
697  labelList& boundaryPts
698  );
699 
700  label mergeIdenticalDualVertices
701  (
702  const pointField& pts,
703  Map<label>& dualPtIndexMap
704  ) const;
705 
706  //- Identify the face labels of the deferred collapse faces
707  void deferredCollapseFaceSet
708  (
709  labelList& owner,
710  labelList& neighbour,
711  const HashSet<labelPair, labelPair::Hash<>>& deferredCollapseFaces
712  ) const;
713 
714  //- Check whether the cell sizes are fine enough. Creates a polyMesh.
715  void checkCellSizing();
716 
717  //- Find all cells with a patch face that is not near the surface. The
718  // allowed offset is the fraction of the target cell size.
719  labelHashSet findOffsetPatchFaces
720  (
721  const polyMesh& mesh,
722  const scalar allowedOffset
723  ) const;
724 
725  //- Create a polyMesh and check its quality, reports which
726  // elements damage the mesh quality, allowing backtracking.
727  labelHashSet checkPolyMeshQuality(const pointField& pts) const;
728 
729  label classifyBoundaryPoint(Cell_handle cit) const;
730 
731  //- Index all of the Delaunay cells and calculate their
732  //- Dual points
733  void indexDualVertices
734  (
735  pointField& pts,
736  labelList& boundaryPts
737  );
738 
739  //- Re-index all of the Delaunay cells
740  void reindexDualVertices
741  (
742  const Map<label>& dualPtIndexMap,
743  labelList& boundaryPts
744  );
745 
746  label createPatchInfo
747  (
750  ) const;
751 
752  vector calcSharedPatchNormal(Cell_handle c1, Cell_handle c2) const;
753 
754  bool boundaryDualFace(Cell_handle c1, Cell_handle c2) const;
755 
756  //- Create all of the internal and boundary faces
757  void createFacesOwnerNeighbourAndPatches
758  (
759  const pointField& pts,
760  faceList& faces,
761  labelList& owner,
762  labelList& neighbour,
765  labelListList& patchPointPairSlaves,
766  PackedBoolList& boundaryFacesToRemove,
767  bool includeEmptyPatches = false
768  ) const;
769 
770  //- Sort the faces, owner and neighbour lists into
771  // upper-triangular order. For internal faces only, use
772  // before adding patch faces
773  void sortFaces
774  (
775  faceList& faces,
776  labelList& owner,
777  labelList& neighbour
778  ) const;
779 
780  //- Sort the processor patches so that the faces are in the same order
781  // on both processors
782  void sortProcPatches
783  (
784  List<DynamicList<face>>& patchFaces,
785  List<DynamicList<label>>& patchOwners,
786  List<DynamicList<label>>& patchPointPairSlaves,
787  labelPairPairDynListList& patchSortingIndices
788  ) const;
789 
790  //- Add the faces and owner information for the patches
791  void addPatches
792  (
793  const label nInternalFaces,
794  faceList& faces,
795  labelList& owner,
797  PackedBoolList& boundaryFacesToRemove,
798  const List<DynamicList<face>>& patchFaces,
799  const List<DynamicList<label>>& patchOwners,
800  const List<DynamicList<bool>>& indirectPatchFace
801  ) const;
802 
803  //- Remove points that are no longer used by any faces
804  void removeUnusedPoints
805  (
806  faceList& faces,
807  pointField& pts,
808  labelList& boundaryPts
809  ) const;
810 
811  //- Remove dual cells that are not used by any faces. Return compaction
812  // map.
813  labelList removeUnusedCells
814  (
815  labelList& owner,
816  labelList& neighbour
817  ) const;
818 
819  //- Create an empty fvMesh
820  autoPtr<fvMesh> createDummyMesh
821  (
822  const IOobject& io,
823  const wordList& patchNames,
825  ) const;
826 
827  //- Create a polyMesh from points.
828  autoPtr<polyMesh> createPolyMeshFromPoints(const pointField& pts) const;
829 
830  void checkProcessorPatchesMatch
831  (
833  ) const;
834 
835  void reorderPoints
836  (
838  labelList& boundaryPts,
839  faceList& faces,
840  const label nInternalFaces
841  ) const;
842 
843  //- Rotate the faces on processor patches if necessary
844  void reorderProcessorPatches
845  (
846  const word& meshName,
847  const fileName& instance,
848  const pointField& points,
849  faceList& faces,
850  const wordList& patchNames,
852  ) const;
853 
854  void writePointPairs(const fileName& fName) const;
855 
856  //- Disallow default bitwise copy construct
858 
859  //- Disallow default bitwise assignment
860  void operator=(const conformalVoronoiMesh&);
861 
862 
863 public:
864 
865  //- Runtime type information
866  ClassName("conformalVoronoiMesh");
867 
868 
869  // Constructors
870 
871  //- Construct from Time and foamyHexMeshDict
873  (
874  const Time& runTime,
875  const dictionary& foamyHexMeshDict
876  );
877 
878 
879  //- Destructor
881 
882 
883  // Member Functions
884 
885  void initialiseForMotion();
886 
888 
889  //- Move the vertices according to the controller, re-conforming to the
890  // surface as required
891  void move();
892 
893 // //- Which other processors does each sphere overlap
894 // labelListList overlapsProc
895 // (
896 // const List<Foam::point>& centres,
897 // const List<scalar>& radiusSqrs
898 // ) const;
899 
900 
901  // Access
902 
903  //- Return the Time object
904  inline const Time& time() const;
905 
906  //- Return the random number generator
907  inline Random& rndGen() const;
908 
909  //- Return the allGeometry object
910  inline const searchableSurfaces& allGeometry() const;
911 
912  //- Return the conformationSurfaces object
913  inline const conformationSurfaces& geometryToConformTo() const;
914 
915  //- Return the backgroundMeshDecomposition
916  inline const backgroundMeshDecomposition& decomposition() const;
917 
918  //- Return the cellShapeControl object
919  inline const cellShapeControl& cellShapeControls() const;
920 
921  //- Return the foamyHexMeshControls object
922  inline const cvControls& foamyHexMeshControls() const;
923 
924 
925  // Query
926 
927  //- Return the local point pair separation at the given location
928  inline scalar pointPairDistance(const Foam::point& pt) const;
929 
930  //- Return the local mixed feature point placement distance
931  inline scalar mixedFeaturePointDistance
932  (
933  const Foam::point& pt
934  ) const;
935 
936  //- Return the square of the local feature point exclusion distance
938  (
939  const Foam::point& pt
940  ) const;
941 
942  //- Return the square of the local feature edge exclusion distance
943  inline scalar featureEdgeExclusionDistanceSqr
944  (
945  const Foam::point& pt
946  ) const;
947 
948  //- Return the square of the local surface point exclusion distance
949  inline scalar surfacePtExclusionDistanceSqr
950  (
951  const Foam::point& pt
952  ) const;
953 
954  //- Return the square of the local surface search distance
955  inline scalar surfaceSearchDistanceSqr(const Foam::point& pt) const;
956 
957  //- Return the local maximum surface protrusion distance
958  inline scalar maxSurfaceProtrusion(const Foam::point& pt) const;
959 
960  //- Call the appropriate function to conform to an edge
962  (
963  const extendedFeatureEdgeMesh& feMesh,
964  const pointIndexHit& edHit,
965  DynamicList<Vb>& pts
966  ) const;
967 
968 
969  // Write
970 
971  //- Write the elapsedCpuTime and memory usage, with an optional
972  // description
973  static void timeCheck
974  (
975  const Time& runTime,
976  const string& description = string::null,
977  const bool check = true
978  );
979 
980  void timeCheck
981  (
982  const string& description = string::null
983  ) const;
984 
985  //- Prepare data and call writeMesh for polyMesh and
986  // tetDualMesh
987  void writeMesh(const fileName& instance);
988 
989  //- Write mesh to disk
990  void writeMesh
991  (
992  const word& meshName,
993  const fileName& instance,
995  labelList& boundaryPts,
996  faceList& faces,
997  labelList& owner,
998  labelList& neighbour,
999  const wordList& patchNames,
1001  const pointField& cellCentres,
1002  PackedBoolList& boundaryFacesToRemove
1003  ) const;
1004 
1005  //- Calculate and write a field of the target cell size,
1006  // target cell volume, actual cell volume and equivalent
1007  // actual cell size (cbrt(actual cell volume)).
1008  void writeCellSizes(const fvMesh& mesh) const;
1009 
1010  void writeCellAlignments(const fvMesh& mesh) const;
1011 
1012  //- Calculate and write the cell centres.
1013  void writeCellCentres(const fvMesh& mesh) const;
1014 
1015  //- Find the cellSet of the boundary cells which have points that
1016  // protrude out of the surface beyond a tolerance.
1018 };
1019 
1020 
1021 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1022 
1023 } // End namespace Foam
1024 
1025 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1026 
1027 #include "conformalVoronoiMeshI.H"
1028 
1029 #ifdef NoRepository
1031 #endif
1032 
1033 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1034 
1035 #endif
1036 
1037 // ************************************************************************* //
labelHashSet findRemainingProtrusionSet(const polyMesh &mesh) const
Find the cellSet of the boundary cells which have points that.
Delaunay::Cell_handle Cell_handle
A HashTable with keys but without contents.
Definition: HashSet.H:59
~conformalVoronoiMesh()
Destructor.
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
A class for handling file names.
Definition: fileName.H:69
Delaunay::Vertex_handle Vertex_handle
const searchableSurfaces & allGeometry() const
Return the allGeometry object.
scalar surfaceSearchDistanceSqr(const Foam::point &pt) const
Return the square of the local surface search distance.
scalar pointPairDistance(const Foam::point &pt) const
Return the local point pair separation at the given location.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
A 2-tuple for storing two objects of different types.
Definition: HashTable.H:66
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
Unit conversion functions.
const cellShapeControl & cellShapeControls() const
Return the cellShapeControl object.
const dimensionedScalar c2
Second radiation constant: default SI units: [m.K].
scalar mixedFeaturePointDistance(const Foam::point &pt) const
Return the local mixed feature point placement distance.
engineTime & runTime
static void timeCheck(const Time &runTime, const string &description=string::null, const bool check=true)
Write the elapsedCpuTime and memory usage, with an optional.
ClassName("conformalVoronoiMesh")
Runtime type information.
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
void writeCellCentres(const fvMesh &mesh) const
Calculate and write the cell centres.
List< pointIndexHitAndFeature > pointIndexHitAndFeatureList
const conformationSurfaces & geometryToConformTo() const
Return the conformationSurfaces object.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
void writeMesh(const fileName &instance)
Prepare data and call writeMesh for polyMesh and.
Random & rndGen() const
Return the random number generator.
void createEdgePointGroup(const extendedFeatureEdgeMesh &feMesh, const pointIndexHit &edHit, DynamicList< Vb > &pts) const
Call the appropriate function to conform to an edge.
void move()
Move the vertices according to the controller, re-conforming to the.
pointField vertices(const blockVertexList &bvl)
void writeCellSizes(const fvMesh &mesh) const
Calculate and write a field of the target cell size,.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
dynamicFvMesh & mesh
static char meshName[]
Definition: globalFoam.H:7
Tuple2< pointIndexHit, label > pointIndexHitAndFeature
const pointField & points
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
A class for handling words, derived from string.
Definition: word.H:59
3D tensor transformation operations.
wordList patchNames(nPatches)
Container for searchableSurfaces.
sideVolumeType
Normals point to the outside.
scalar surfacePtExclusionDistanceSqr(const Foam::point &pt) const
Return the square of the local surface point exclusion distance.
Random number generator.
Definition: Random.H:57
Controls for the conformalVoronoiMesh mesh generator.
Definition: cvControls.H:53
static const string null
An empty string.
Definition: string.H:86
DynamicList< pointIndexHitAndFeature > pointIndexHitAndFeatureDynList
const backgroundMeshDecomposition & decomposition() const
Return the backgroundMeshDecomposition.
The Delaunay vertices required to conform to a feature point can be determined upon initialisation be...
scalar maxSurfaceProtrusion(const Foam::point &pt) const
Return the local maximum surface protrusion distance.
A bit-packed bool list.
side
Side of the plane.
Definition: plane.H:65
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
void writeCellAlignments(const fvMesh &mesh) const
static const NamedEnum< dualMeshPointType, 5 > dualMeshPointTypeNames_
label n
Store a background polyMesh to use for the decomposition of space and queries for parallel conformalV...
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const cvControls & foamyHexMeshControls() const
Return the foamyHexMeshControls object.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const Time & time() const
Return the Time object.
scalar featureEdgeExclusionDistanceSqr(const Foam::point &pt) const
Return the square of the local feature edge exclusion distance.
CGAL data structures used for 3D Delaunay meshing.
PtrList< dictionary > patchDicts
Definition: readKivaGrid.H:537
scalar featurePointExclusionDistanceSqr(const Foam::point &pt) const
Return the square of the local feature point exclusion distance.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
List< DynamicList< Pair< labelPair > > > labelPairPairDynListList
Namespace for OpenFOAM.