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