All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2019 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  //- Merge distance required by fvMeshDistribute
134  scalar mergeDist_;
135 
136  //- Scale of a cell span vs cell size used to decide to refine a cell
137  scalar spanScale_;
138 
139  //- Smallest minimum cell size allowed, i.e. to avoid high initial
140  // refinement of areas of small size
141  scalar minCellSizeLimit_;
142 
143  //- Minimum normal level of refinement
144  label minLevels_;
145 
146  //- How fine should the initial sample of the volume a box be to
147  // investigate the local cell size
148  label volRes_;
149 
150  //- Allowed factor above the average cell weight before a background
151  // cell needs to be split
152  scalar maxCellWeightCoeff_;
153 
154 
155  // Private Member Functions
156 
157  void initialRefinement();
158 
159  //- Print details of the decomposed mesh
160  void printMeshData(const polyMesh& mesh) const;
161 
162  //- Estimate the number of vertices that will be in this cell, returns
163  // true if the cell is to be split because of the density ratio inside
164  // it
165  bool refineCell
166  (
167  label celli,
168  volumeType volType,
169  scalar& weightEstimate
170  ) const;
171 
172  //- Select cells for refinement at the surface of the geometry to be
173  // meshed
174  labelList selectRefinementCells
175  (
176  List<volumeType>& volumeStatus,
177  volScalarField& cellWeights
178  ) const;
179 
180  //- Build the surface patch and search tree
181  void buildPatchAndTree();
182 
183 
184 public:
185 
186  //- Runtime type information
187  ClassName("backgroundMeshDecomposition");
188 
189 
190  // Constructors
191 
192  //- Construct from components in foamyHexMesh operation
194  (
195  const Time& runTime,
196  Random& rndGen,
197  const conformationSurfaces& geometryToConformTo,
198  const dictionary& coeffsDict
199  );
200 
201  //- Disallow default bitwise copy construction
203  (
205  ) = delete;
206 
207 
208  //- Destructor
210 
211 
212  // Member Functions
213 
214  //- Build a mapDistribute for the supplied destination processor data
215  static autoPtr<mapDistribute> buildMap(const List<label>& toProc);
216 
217  //- Redistribute the background mesh based on a supplied weight field,
218  // returning a map to use to redistribute vertices.
220  (
221  volScalarField& cellWeights
222  );
223 
224  //- Distribute supplied the points to the appropriate processor
225  template<class PointType>
227 
228  //- Is the given position inside the domain of this decomposition
229  bool positionOnThisProcessor(const point& pt) const;
230 
231  //- Are the given positions inside the domain of this decomposition
232  boolList positionOnThisProcessor(const List<point>& pts) const;
233 
234  //- Does the given box overlap the faces of the boundary of this
235  // processor
236  bool overlapsThisProcessor(const treeBoundBox& box) const;
237 
238  //- Does the given sphere overlap the faces of the boundary of this
239  // processor
241  (
242  const point& centre,
243  const scalar radiusSqr
244  ) const;
245 
246  //- Find nearest intersection of line between start and end, (exposing
247  // underlying indexedOctree)
249  (
250  const point& start,
251  const point& end
252  ) const;
253 
254  //- Find any intersection of line between start and end, (exposing
255  // underlying indexedOctree)
257  (
258  const point& start,
259  const point& end
260  ) const;
261 
262  //- What processor is the given position on?
263  template<class PointType>
264  labelList processorPosition(const List<PointType>& pts) const;
265 
266  //- What is the nearest processor to the given position?
268 
269  //- Which processors are intersected by the line segment, returns all
270  // processors whose boundary patch is intersected by the sphere. By
271  // default this does not return the processor that the query is
272  // launched from, it is assumed that the point is on that processor.
273  // The index data member of the pointIndexHit is replaced with the
274  // processor index.
276  (
277  const List<point>& starts,
278  const List<point>& ends,
279  bool includeOwnProcessor = false
280  ) const;
281 
283  (
284  const point& centre,
285  const scalar& radiusSqr
286  ) const;
287 
289  (
290  const point& centre,
291  const scalar radiusSqr
292  ) const;
293 
294 
295  // Access
296 
297  //- Return access to the underlying mesh
298  inline const fvMesh& mesh() const;
299 
300  //- Return access to the underlying tree
301  inline const indexedOctree<treeDataBPatch>& tree() const;
302 
303  //- Return the boundBox of this processor
304  inline const treeBoundBox& procBounds() const;
305 
306  //- Return the cell level of the underlying mesh
307  inline const labelList& cellLevel() const;
308 
309  //- Return the point level of the underlying mesh
310  inline const labelList& pointLevel() const;
311 
312  //- Return the current decomposition method
313  inline const decompositionMethod& decomposer() const;
314 
315 
316  // Member Operators
317 
318  //- Disallow default bitwise assignment
319  void operator=(const backgroundMeshDecomposition&) = delete;
320 };
321 
322 
323 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
324 
325 } // End namespace Foam
326 
327 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
328 
330 
331 #ifdef NoRepository
333 #endif
334 
335 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
336 
337 #endif
338 
339 // ************************************************************************* //
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:158
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:52
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.