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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-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::backgroundMeshDecomposition
28 
29 Description
30  Store a background polyMesh to use for the decomposition of space and
31  queries for parallel conformalVoronoiMesh.
32 
33  The requirements are:
34 
35  - To have a decomposition of space which can quickly interrogate an
36  arbitrary location from any processor to reliably and unambiguously
37  determine which processor owns the space that the point is in, i.e. as
38  the vertices move, or need inserted as part of the surface conformation,
39  send them to the correct proc.
40 
41  - To be able to be dynamically built, refined and redistributed to other
42  procs the partitioning as the meshing progresses to balance the load.
43 
44  - To be able to query whether a sphere (the circumsphere of a Delaunay tet)
45  overlaps any part of the space defined by the structure, and whether a
46  ray (Voronoi edge) penetrates any part of the space defined by the
47  structure, this is what determines if points get referred to a processor.
48 
49 SourceFiles
50  backgroundMeshDecompositionI.H
51  backgroundMeshDecomposition.C
52 
53 \*---------------------------------------------------------------------------*/
54 
55 #ifndef backgroundMeshDecomposition_H
56 #define backgroundMeshDecomposition_H
57 
58 #include "fvMesh.H"
59 #include "hexRef8.H"
60 #include "cellSet.H"
61 #include "meshTools.H"
62 #include "polyTopoChange.H"
63 #include "mapPolyMesh.H"
64 #include "decompositionMethod.H"
65 #include "fvMeshDistribute.H"
66 #include "removeCells.H"
67 #include "mapDistributePolyMesh.H"
68 #include "globalIndex.H"
69 #include "treeBoundBox.H"
70 #include "primitivePatch.H"
71 #include "labelList.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, specifying 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  //- merge distance required by fvMeshDistribute
131  scalar mergeDist_;
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  //- No copy construct
182  (
184  ) = delete;
185 
186  //- No copy assignment
187  void operator=(const backgroundMeshDecomposition&) = delete;
188 
189 
190 public:
191 
192  //- Runtime type information
193  ClassName("backgroundMeshDecomposition");
194 
195 
196  // Constructors
197 
198  //- Construct from components in foamyHexMesh operation
200  (
201  const Time& runTime,
202  Random& rndGen,
203  const conformationSurfaces& geometryToConformTo,
204  const dictionary& coeffsDict,
205  const fileName& decompDictFile = ""
206  );
207 
208 
209  //- Destructor
210  ~backgroundMeshDecomposition() = default;
211 
212 
213  // Member Functions
214 
215  //- Build a mapDistribute for the supplied destination processor data
216  static autoPtr<mapDistribute> buildMap(const List<label>& toProc);
217 
218  //- Redistribute the background mesh based on a supplied weight field,
219  // returning a map to use to redistribute vertices.
221  (
222  volScalarField& cellWeights
223  );
224 
225  //- Distribute supplied the points to the appropriate processor
226  template<class PointType>
228 
229  //- Is the given position inside the domain of this decomposition
230  bool positionOnThisProcessor(const point& pt) const;
231 
232  //- Are the given positions inside the domain of this decomposition
234 
235  //- Does the given box overlap the faces of the boundary of this
236  // processor
237  bool overlapsThisProcessor(const treeBoundBox& box) const;
238 
239  //- Does the given sphere overlap the faces of the boundary of this
240  // processor
242  (
243  const point& centre,
244  const scalar radiusSqr
245  ) const;
246 
247  //- Find nearest intersection of line between start and end, (exposing
248  // underlying indexedOctree)
250  (
251  const point& start,
252  const point& end
253  ) const;
254 
255  //- Find any intersection of line between start and end, (exposing
256  // underlying indexedOctree)
258  (
259  const point& start,
260  const point& end
261  ) const;
262 
263  //- What processor is the given position on?
264  template<class PointType>
266 
267  //- What is the nearest processor to the given position?
269 
270  //- Which processors are intersected by the line segment, returns all
271  // processors whose boundary patch is intersected by the sphere. By
272  // default this does not return the processor that the query is
273  // launched from, it is assumed that the point is on that processor.
274  // The index data member of the pointIndexHit is replaced with the
275  // processor index.
277  (
278  const List<point>& starts,
279  const List<point>& ends,
280  bool includeOwnProcessor = false
281  ) const;
282 
284  (
285  const point& centre,
286  const scalar& radiusSqr
287  ) const;
288 
290  (
291  const point& centre,
292  const scalar radiusSqr
293  ) const;
294 
295 // //- Which processors overlap the given sphere, returns all processors
296 // // whose boundary patch is touched by the sphere or whom the sphere
297 // // is inside. By default this does not return the processor that the
298 // // query is launched from, it is assumed that the point is on that
299 // // processor.
300 // labelListList overlapsProcessors
301 // (
302 // const List<point>& centres,
303 // const List<scalar>& radiusSqrs,
304 // const Delaunay& T,
305 // bool includeOwnProcessor
306 // ) const;
307 
308  // Access
309 
310  //- Return access to the underlying mesh
311  inline const fvMesh& mesh() const;
312 
313  //- Return access to the underlying tree
314  inline const indexedOctree<treeDataBPatch>& tree() const;
315 
316  //- Return the boundBox of this processor
317  inline const treeBoundBox& procBounds() const;
318 
319  //- Return the cell level of the underlying mesh
320  inline const labelList& cellLevel() const;
321 
322  //- Return the point level of the underlying mesh
323  inline const labelList& pointLevel() const;
324 
325 };
326 
327 
328 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
329 
330 } // End namespace Foam
331 
332 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
333 
335 
336 #ifdef NoRepository
338 #endif
339 
340 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
341 
342 #endif
343 
344 // ************************************************************************* //
bool positionOnThisProcessor(const point &pt) const
Is the given position inside the domain of this decomposition.
PrimitivePatch<::Foam::List< face >, const pointField > bPatch
A class for handling file names.
Definition: fileName.H:72
const fvMesh & mesh() const
Return access to the underlying mesh.
~backgroundMeshDecomposition()=default
Destructor.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const indexedOctree< treeDataBPatch > & tree() const
Return access to the underlying tree.
ClassName("backgroundMeshDecomposition")
Runtime type information.
Random rndGen
Definition: createFields.H:23
engineTime & runTime
This class describes the interaction of an object (often a face) and a point. It carries the info of ...
Definition: pointIndexHit.H:44
labelList overlapProcessors(const point &centre, const scalar radiusSqr) const
An enumeration wrapper for classification of a location as being inside/outside of a volume...
Definition: volumeType.H:55
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
static autoPtr< mapDistribute > buildMap(const List< label > &toProc)
Build a mapDistribute for the supplied destination processor data.
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:62
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.
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 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:52
Random number generator.
Definition: Random.H:55
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
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.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
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
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.
const treeBoundBox & procBounds() const
Return the boundBox of this processor.
Namespace for OpenFOAM.
bool overlapsThisProcessor(const treeBoundBox &box) const
Does the given box overlap the faces of the boundary of this.
const pointField & pts