backgroundMeshDecomposition.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 | Copyright (C) 2011-2016 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  //- Random number generator
107  Random& rndGen_;
108 
109  //- Mesh stored on for this processor, specifiying the domain that it
110  // is responsible for.
111  fvMesh mesh_;
112 
113  //- Refinement object
114  hexRef8 meshCutter_;
115 
116  //- Patch containing an independent representation of the surface of the
117  // mesh of this processor
118  autoPtr<bPatch> boundaryFacesPtr_;
119 
120  //- Search tree for the boundaryFaces_ patch
122 
123  //- The bounds of all background meshes on all processors
124  treeBoundBoxList allBackgroundMeshBounds_;
125 
126  //- The overall bounds of all of the background meshes, used to test if
127  // a point that is not found on any processor is in the domain at all
128  treeBoundBox globalBackgroundBounds_;
129 
130  //- Decomposition dictionary
131  IOdictionary decomposeDict_;
132 
133  //- Decomposition method
134  autoPtr<decompositionMethod> decomposerPtr_;
135 
136  //- Merge distance required by fvMeshDistribute
137  scalar mergeDist_;
138 
139  //- Scale of a cell span vs cell size used to decide to refine a cell
140  scalar spanScale_;
141 
142  //- Smallest minimum cell size allowed, i.e. to avoid high initial
143  // refinement of areas of small size
144  scalar minCellSizeLimit_;
145 
146  //- Minimum normal level of refinement
147  label minLevels_;
148 
149  //- How fine should the initial sample of the volume a box be to
150  // investigate the local cell size
151  label volRes_;
152 
153  //- Allowed factor above the average cell weight before a background
154  // cell needs to be split
155  scalar maxCellWeightCoeff_;
156 
157 
158  // Private Member Functions
159 
160  void initialRefinement();
161 
162  //- Print details of the decomposed mesh
163  void printMeshData(const polyMesh& mesh) const;
164 
165  //- Estimate the number of vertices that will be in this cell, returns
166  // true if the cell is to be split because of the density ratio inside
167  // it
168  bool refineCell
169  (
170  label celli,
171  volumeType volType,
172  scalar& weightEstimate
173  ) const;
174 
175  //- Select cells for refinement at the surface of the geometry to be
176  // meshed
177  labelList selectRefinementCells
178  (
179  List<volumeType>& volumeStatus,
180  volScalarField& cellWeights
181  ) const;
182 
183  //- Build the surface patch and search tree
184  void buildPatchAndTree();
185 
186  //- Disallow default bitwise copy construct
188 
189  //- Disallow default bitwise assignment
190  void operator=(const backgroundMeshDecomposition&);
191 
192 
193 public:
194 
195  //- Runtime type information
196  ClassName("backgroundMeshDecomposition");
197 
198 
199  // Constructors
200 
201  //- Construct from components in foamyHexMesh operation
203  (
204  const Time& runTime,
205  Random& rndGen,
206  const conformationSurfaces& geometryToConformTo,
207  const dictionary& coeffsDict
208  );
209 
210 
211  //- Destructor
213 
214 
215  // Member Functions
216 
217  //- Build a mapDistribute for the supplied destination processor data
218  static autoPtr<mapDistribute> buildMap(const List<label>& toProc);
219 
220  //- Redistribute the background mesh based on a supplied weight field,
221  // returning a map to use to redistribute vertices.
223  (
224  volScalarField& cellWeights
225  );
226 
227  //- Distribute supplied the points to the appropriate processor
228  template<class PointType>
230 
231  //- Is the given position inside the domain of this decomposition
232  bool positionOnThisProcessor(const point& pt) const;
233 
234  //- Are the given positions inside the domain of this decomposition
235  boolList positionOnThisProcessor(const List<point>& pts) const;
236 
237  //- Does the given box overlap the faces of the boundary of this
238  // processor
239  bool overlapsThisProcessor(const treeBoundBox& box) const;
240 
241  //- Does the given sphere overlap the faces of the boundary of this
242  // processor
244  (
245  const point& centre,
246  const scalar radiusSqr
247  ) const;
248 
249  //- Find nearest intersection of line between start and end, (exposing
250  // underlying indexedOctree)
252  (
253  const point& start,
254  const point& end
255  ) const;
256 
257  //- Find any intersection of line between start and end, (exposing
258  // underlying indexedOctree)
260  (
261  const point& start,
262  const point& end
263  ) const;
264 
265  //- What processor is the given position on?
266  template<class PointType>
267  labelList processorPosition(const List<PointType>& pts) const;
268 
269  //- What is the nearest processor to the given position?
271 
272  //- Which processors are intersected by the line segment, returns all
273  // processors whose boundary patch is intersected by the sphere. By
274  // default this does not return the processor that the query is
275  // launched from, it is assumed that the point is on that processor.
276  // The index data member of the pointIndexHit is replaced with the
277  // processor index.
279  (
280  const List<point>& starts,
281  const List<point>& ends,
282  bool includeOwnProcessor = false
283  ) const;
284 
286  (
287  const point& centre,
288  const scalar& radiusSqr
289  ) const;
290 
292  (
293  const point& centre,
294  const scalar radiusSqr
295  ) const;
296 
297 // //- Which processors overlap the given sphere, returns all processors
298 // // whose boundary patch is touched by the sphere or whom the sphere
299 // // is inside. By default this does not return the processor that the
300 // // query is launched from, it is assumed that the point is on that
301 // // processor.
302 // labelListList overlapsProcessors
303 // (
304 // const List<point>& centres,
305 // const List<scalar>& radiusSqrs,
306 // const Delaunay& T,
307 // bool includeOwnProcessor
308 // ) const;
309 
310  // Access
311 
312  //- Return access to the underlying mesh
313  inline const fvMesh& mesh() const;
314 
315  //- Return access to the underlying tree
316  inline const indexedOctree<treeDataBPatch>& tree() const;
317 
318  //- Return the boundBox of this processor
319  inline const treeBoundBox& procBounds() const;
320 
321  //- Return the cell level of the underlying mesh
322  inline const labelList& cellLevel() const;
323 
324  //- Return the point level of the underlying mesh
325  inline const labelList& pointLevel() const;
326 
327  //- Return the current decomposition method
328  inline const decompositionMethod& decomposer() const;
329 };
330 
331 
332 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
333 
334 } // End namespace Foam
335 
336 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
337 
339 
340 #ifdef NoRepository
342 #endif
343 
344 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
345 
346 #endif
347 
348 // ************************************************************************* //
cachedRandom rndGen(label(0),-1)
labelList processorNearestPosition(const List< point > &pts) const
What is the nearest processor to the given position?
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 list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
ClassName("backgroundMeshDecomposition")
Runtime type information.
labelList overlapProcessors(const point &centre, const scalar radiusSqr) const
bool positionOnThisProcessor(const point &pt) const
Is the given position inside the domain of this decomposition.
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
const labelList & cellLevel() const
Return the cell level of the underlying mesh.
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
const labelList & pointLevel() const
Return the point level of the underlying mesh.
treeDataPrimitivePatch< bPatch > treeDataBPatch
const fvMesh & mesh() const
Return access to the underlying mesh.
A list of faces which address into the list of points.
const pointField & points
Refinement of (split) hexes using polyTopoChange.
Definition: hexRef8.H:64
labelList processorPosition(const List< PointType > &pts) const
What processor is the given position on?
Encapsulation of data needed to search on PrimitivePatches.
Abstract base class for decomposition.
Container with cells to refine. Refinement given as single direction.
Definition: refineCell.H:56
PrimitivePatch< face, List, const pointField > bPatch
Simple random number generator.
Definition: Random.H:49
pointIndexHit findLine(const point &start, const point &end) const
Find nearest intersection of line between start and end, (exposing.
bool overlapsThisProcessor(const treeBoundBox &box) const
Does the given box overlap the faces of the boundary of this.
const decompositionMethod & decomposer() const
Return the current decomposition method.
bool overlapsOtherProcessors(const point &centre, const scalar &radiusSqr) const
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
const indexedOctree< treeDataBPatch > & tree() const
Return access to the underlying tree.
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:53
autoPtr< mapDistribute > distributePoints(List< PointType > &points) const
Distribute supplied the points to the appropriate processor.
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.
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.
const treeBoundBox & procBounds() const
Return the boundBox of this processor.
~backgroundMeshDecomposition()
Destructor.
Namespace for OpenFOAM.
pointIndexHit findLineAny(const point &start, const point &end) const
Find any intersection of line between start and end, (exposing.