backgroundMeshDecomposition.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) 2011-2021 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::backgroundMeshDecomposition
26 
27 Description
28  Store a background polyMesh to use for the decomposition of space and
29  queries for parallel conformalVoronoiMesh.
30 
31  The requirements are:
32 
33  - To have a decomposition of space which can quickly interrogate an
34  arbitrary location from any processor to reliably and unambiguously
35  determine which processor owns the space that the point is in, i.e. as
36  the vertices move, or need inserted as part of the surface conformation,
37  send them to the correct proc.
38 
39  - To be able to be dynamically built, refined and redistributed to other
40  procs the partitioning as the meshing progresses to balance the load.
41 
42  - To be able to query whether a sphere (the circumsphere of a Delaunay tet)
43  overlaps any part of the space defined by the structure, and whether a
44  ray (Voronoi edge) penetrates any part of the space defined by the
45  structure, this is what determines if points get referred to a processor.
46 
47 SourceFiles
48  backgroundMeshDecompositionI.H
49  backgroundMeshDecomposition.C
50 
51 \*---------------------------------------------------------------------------*/
52 
53 #ifndef backgroundMeshDecomposition_H
54 #define backgroundMeshDecomposition_H
55 
56 #include "fvMesh.H"
57 #include "hexRef8.H"
58 #include "cellSet.H"
59 #include "meshTools.H"
60 #include "polyTopoChange.H"
61 #include "mapPolyMesh.H"
62 #include "decompositionMethod.H"
63 #include "fvMeshDistribute.H"
64 #include "removeCells.H"
65 #include "mapDistributePolyMesh.H"
66 #include "globalIndex.H"
67 #include "treeBoundBox.H"
68 #include "primitivePatch.H"
69 #include "face.H"
70 #include "labelList.H"
71 #include "pointField.H"
72 #include "indexedOctree.H"
73 #include "treeDataPrimitivePatch.H"
74 #include "volumeType.H"
76 
77 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78 
79 namespace Foam
80 {
81 
84 
85 class Time;
86 class Random;
88 
89 /*---------------------------------------------------------------------------*\
90  Class backgroundMeshDecomposition Declaration
91 \*---------------------------------------------------------------------------*/
92 
94 {
95  // Private Data
96 
97  //- Method details dictionary
98  // dictionary coeffsDict_;
99 
100  //- Reference to runtime
101  const Time& runTime_;
102 
103  //- Reference to surface
104  const conformationSurfaces& geometryToConformTo_;
105 
106  //- Mesh stored on for this processor, specifying the domain that it
107  // is responsible for.
108  fvMesh mesh_;
109 
110  //- Refinement object
111  hexRef8 meshCutter_;
112 
113  //- Patch containing an independent representation of the surface of the
114  // mesh of this processor
115  autoPtr<bPatch> boundaryFacesPtr_;
116 
117  //- Search tree for the boundaryFaces_ patch
119 
120  //- The bounds of all background meshes on all processors
121  treeBoundBoxList allBackgroundMeshBounds_;
122 
123  //- The overall bounds of all of the background meshes, used to test if
124  // a point that is not found on any processor is in the domain at all
125  treeBoundBox globalBackgroundBounds_;
126 
127  //- Decomposition dictionary
128  IOdictionary decomposeDict_;
129 
130  //- Decomposition method
131  autoPtr<decompositionMethod> decomposerPtr_;
132 
133  //- Scale of a cell span vs cell size used to decide to refine a cell
134  scalar spanScale_;
135 
136  //- Smallest minimum cell size allowed, i.e. to avoid high initial
137  // refinement of areas of small size
138  scalar minCellSizeLimit_;
139 
140  //- Minimum normal level of refinement
141  label minLevels_;
142 
143  //- How fine should the initial sample of the volume a box be to
144  // investigate the local cell size
145  label volRes_;
146 
147  //- Allowed factor above the average cell weight before a background
148  // cell needs to be split
149  scalar maxCellWeightCoeff_;
150 
151 
152  // Private Member Functions
153 
154  void initialRefinement();
155 
156  //- Print details of the decomposed mesh
157  void printMeshData(const polyMesh& mesh) const;
158 
159  //- Estimate the number of vertices that will be in this cell, returns
160  // true if the cell is to be split because of the density ratio inside
161  // it
162  bool refineCell
163  (
164  label celli,
165  volumeType volType,
166  scalar& weightEstimate
167  ) const;
168 
169  //- Select cells for refinement at the surface of the geometry to be
170  // meshed
171  labelList selectRefinementCells
172  (
173  List<volumeType>& volumeStatus,
174  volScalarField& cellWeights
175  ) const;
176 
177  //- Build the surface patch and search tree
178  void buildPatchAndTree();
179 
180 
181 public:
182 
183  //- Runtime type information
184  ClassName("backgroundMeshDecomposition");
185 
186 
187  // Constructors
188 
189  //- Construct from components in foamyHexMesh operation
191  (
192  const Time& runTime,
193  Random& rndGen,
194  const conformationSurfaces& geometryToConformTo,
195  const dictionary& coeffsDict
196  );
197 
198  //- Disallow default bitwise copy construction
200  (
202  ) = delete;
203 
204 
205  //- Destructor
207 
208 
209  // Member Functions
210 
211  //- Build a mapDistribute for the supplied destination processor data
212  static autoPtr<mapDistribute> buildMap(const List<label>& toProc);
213 
214  //- Redistribute the background mesh based on a supplied weight field,
215  // returning a map to use to redistribute vertices.
217  (
218  volScalarField& cellWeights
219  );
220 
221  //- Distribute supplied the points to the appropriate processor
222  template<class PointType>
224 
225  //- Is the given position inside the domain of this decomposition
226  bool positionOnThisProcessor(const point& pt) const;
227 
228  //- Are the given positions inside the domain of this decomposition
229  boolList positionOnThisProcessor(const List<point>& pts) const;
230 
231  //- Does the given box overlap the faces of the boundary of this
232  // processor
233  bool overlapsThisProcessor(const treeBoundBox& box) const;
234 
235  //- Does the given sphere overlap the faces of the boundary of this
236  // processor
238  (
239  const point& centre,
240  const scalar radiusSqr
241  ) const;
242 
243  //- Find nearest intersection of line between start and end, (exposing
244  // underlying indexedOctree)
246  (
247  const point& start,
248  const point& end
249  ) const;
250 
251  //- Find any intersection of line between start and end, (exposing
252  // underlying indexedOctree)
254  (
255  const point& start,
256  const point& end
257  ) const;
258 
259  //- What processor is the given position on?
260  template<class PointType>
261  labelList processorPosition(const List<PointType>& pts) const;
262 
263  //- What is the nearest processor to the given position?
265 
266  //- Which processors are intersected by the line segment, returns all
267  // processors whose boundary patch is intersected by the sphere. By
268  // default this does not return the processor that the query is
269  // launched from, it is assumed that the point is on that processor.
270  // The index data member of the pointIndexHit is replaced with the
271  // processor index.
273  (
274  const List<point>& starts,
275  const List<point>& ends,
276  bool includeOwnProcessor = false
277  ) const;
278 
280  (
281  const point& centre,
282  const scalar& radiusSqr
283  ) const;
284 
286  (
287  const point& centre,
288  const scalar radiusSqr
289  ) const;
290 
291 
292  // Access
293 
294  //- Return access to the underlying mesh
295  inline const fvMesh& mesh() const;
296 
297  //- Return access to the underlying tree
298  inline const indexedOctree<treeDataBPatch>& tree() const;
299 
300  //- Return the boundBox of this processor
301  inline const treeBoundBox& procBounds() const;
302 
303  //- Return the cell level of the underlying mesh
304  inline const labelList& cellLevel() const;
305 
306  //- Return the point level of the underlying mesh
307  inline const labelList& pointLevel() const;
308 
309  //- Return the current decomposition method
310  inline const decompositionMethod& decomposer() const;
311 
312 
313  // Member Operators
314 
315  //- Disallow default bitwise assignment
316  void operator=(const backgroundMeshDecomposition&) = delete;
317 };
318 
319 
320 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
321 
322 } // End namespace Foam
323 
324 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
325 
327 
328 #ifdef NoRepository
330 #endif
331 
332 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
333 
334 #endif
335 
336 // ************************************************************************* //
bool positionOnThisProcessor(const point &pt) const
Is the given position inside the domain of this decomposition.
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
const fvMesh & mesh() const
Return access to the underlying mesh.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
const indexedOctree< treeDataBPatch > & tree() const
Return access to the underlying tree.
void operator=(const backgroundMeshDecomposition &)=delete
Disallow default bitwise assignment.
ClassName("backgroundMeshDecomposition")
Runtime type information.
engineTime & runTime
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
labelList overlapProcessors(const point &centre, const scalar radiusSqr) const
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
static autoPtr< mapDistribute > buildMap(const List< label > &toProc)
Build a mapDistribute for the supplied destination processor data.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Random rndGen(label(0))
treeDataPrimitivePatch< bPatch > treeDataBPatch
A list of faces which address into the list of points.
labelList processorNearestPosition(const List< point > &pts) const
What is the nearest processor to the given position?
const pointField & points
bool overlapsOtherProcessors(const point &centre, const scalar &radiusSqr) const
Refinement of (split) hexes using polyTopoChange.
Definition: hexRef8.H:64
pointIndexHit findLine(const point &start, const point &end) const
Find nearest intersection of line between start and end, (exposing.
Encapsulation of data needed to search on PrimitivePatches.
const decompositionMethod & decomposer() const
Return the current decomposition method.
List< List< pointIndexHit > > intersectsProcessors(const List< point > &starts, const List< point > &ends, bool includeOwnProcessor=false) const
Which processors are intersected by the line segment, returns all.
Abstract base class for decomposition.
const labelList & cellLevel() const
Return the cell level of the underlying mesh.
autoPtr< mapDistribute > distributePoints(List< PointType > &points) const
Distribute supplied the points to the appropriate processor.
Container with cells to refine. Refinement given as single direction.
Definition: refineCell.H:56
Random number generator.
Definition: Random.H:57
pointIndexHit findLineAny(const point &start, const point &end) const
Find any intersection of line between start and end, (exposing.
labelList processorPosition(const List< PointType > &pts) const
What processor is the given position on?
const labelList & pointLevel() const
Return the point level of the underlying mesh.
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:47
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
PrimitivePatch< faceList, const pointField > bPatch
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
autoPtr< mapDistributePolyMesh > distribute(volScalarField &cellWeights)
Redistribute the background mesh based on a supplied weight field,.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
CGAL data structures used for 3D Delaunay meshing.
backgroundMeshDecomposition(const Time &runTime, Random &rndGen, const conformationSurfaces &geometryToConformTo, const dictionary &coeffsDict)
Construct from components in foamyHexMesh operation.
const treeBoundBox & procBounds() const
Return the boundBox of this processor.
~backgroundMeshDecomposition()
Destructor.
Namespace for OpenFOAM.
bool overlapsThisProcessor(const treeBoundBox &box) const
Does the given box overlap the faces of the boundary of this.