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-2015 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
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
115 
116  //- Patch containing an independent representation of the surface of the
117  // mesh of this processor
119 
120  //- Search tree for the boundaryFaces_ patch
122 
123  //- The bounds of all background meshes on all processors
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
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
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
172  (
173  List<volumeType>& volumeStatus,
174  volScalarField& cellWeights
175  ) const;
176 
177  //- Build the surface patch and search tree
178  void buildPatchAndTree();
179 
180  //- Disallow default bitwise copy construct
182 
183  //- Disallow default bitwise assignment
185 
186 
187 public:
188 
189  //- Runtime type information
190  ClassName("backgroundMeshDecomposition");
191 
192 
193  // Constructors
194 
195  //- Construct from components in foamyHexMesh operation
197  (
198  const Time& runTime,
199  Random& rndGen,
200  const conformationSurfaces& geometryToConformTo,
201  const dictionary& coeffsDict,
202  const fileName& decompDictFile = ""
203  );
204 
205 
206  //- Destructor
208 
209 
210  // Member Functions
211 
212  //- Build a mapDistribute for the supplied destination processor data
213  static autoPtr<mapDistribute> buildMap(const List<label>& toProc);
214 
215  //- Redistribute the background mesh based on a supplied weight field,
216  // returning a map to use to redistribute vertices.
218  (
219  volScalarField& cellWeights
220  );
221 
222  //- Distribute supplied the points to the appropriate processor
223  template<typename PointType>
225 
226  //- Is the given position inside the domain of this decomposition
227  bool positionOnThisProcessor(const point& pt) const;
228 
229  //- Are the given positions inside the domain of this decomposition
230  boolList positionOnThisProcessor(const List<point>& pts) const;
231 
232  //- Does the given box overlap the faces of the boundary of this
233  // processor
234  bool overlapsThisProcessor(const treeBoundBox& box) const;
235 
236  //- Does the given sphere overlap the faces of the boundary of this
237  // processor
239  (
240  const point& centre,
241  const scalar radiusSqr
242  ) const;
243 
244  //- Find nearest intersection of line between start and end, (exposing
245  // underlying indexedOctree)
247  (
248  const point& start,
249  const point& end
250  ) const;
251 
252  //- Find any intersection of line between start and end, (exposing
253  // underlying indexedOctree)
255  (
256  const point& start,
257  const point& end
258  ) const;
259 
260  //- What processor is the given position on?
261  template<typename PointType>
262  labelList processorPosition(const List<PointType>& pts) const;
263 
264  //- What is the nearest processor to the given position?
266 
267  //- Which processors are intersected by the line segment, returns all
268  // processors whose boundary patch is intersected by the sphere. By
269  // default this does not return the processor that the query is
270  // launched from, it is assumed that the point is on that processor.
271  // The index data member of the pointIndexHit is replaced with the
272  // processor index.
274  (
275  const List<point>& starts,
276  const List<point>& ends,
277  bool includeOwnProcessor = false
278  ) const;
279 
281  (
282  const point& centre,
283  const scalar& radiusSqr
284  ) const;
285 
287  (
288  const point& centre,
289  const scalar radiusSqr
290  ) const;
291 
292 // //- Which processors overlap the given sphere, returns all processors
293 // // whose boundary patch is touched by the sphere or whom the sphere
294 // // is inside. By default this does not return the processor that the
295 // // query is launched from, it is assumed that the point is on that
296 // // processor.
297 // labelListList overlapsProcessors
298 // (
299 // const List<point>& centres,
300 // const List<scalar>& radiusSqrs,
301 // const Delaunay& T,
302 // bool includeOwnProcessor
303 // ) const;
304 
305  // Access
306 
307  //- Return access to the underlying mesh
308  inline const fvMesh& mesh() const;
309 
310  //- Return access to the underlying tree
311  inline const indexedOctree<treeDataBPatch>& tree() const;
312 
313  //- Return the boundBox of this processor
314  inline const treeBoundBox& procBounds() const;
315 
316  //- Return the cell level of the underlying mesh
317  inline const labelList& cellLevel() const;
318 
319  //- Return the point level of the underlying mesh
320  inline const labelList& pointLevel() const;
321 
322 };
323 
324 
325 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
326 
327 } // End namespace Foam
328 
329 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
330 
332 
333 #ifdef NoRepository
335 #endif
336 
337 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
338 
339 #endif
340 
341 // ************************************************************************* //
Foam::Random
Simple random number generator.
Definition: Random.H:49
Foam::backgroundMeshDecomposition::meshCutter_
hexRef8 meshCutter_
Refinement object.
Definition: backgroundMeshDecomposition.H:113
meshTools.H
Foam::backgroundMeshDecomposition::selectRefinementCells
labelList selectRefinementCells(List< volumeType > &volumeStatus, volScalarField &cellWeights) const
Select cells for refinement at the surface of the geometry to be.
Foam::backgroundMeshDecomposition::rndGen_
Random & rndGen_
Random number generator.
Definition: backgroundMeshDecomposition.H:106
backgroundMeshDecompositionTemplates.C
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::backgroundMeshDecomposition::~backgroundMeshDecomposition
~backgroundMeshDecomposition()
Destructor.
Foam::backgroundMeshDecomposition::operator=
void operator=(const backgroundMeshDecomposition &)
Disallow default bitwise assignment.
Foam::backgroundMeshDecomposition::minLevels_
label minLevels_
Minimum normal level of refinement.
Definition: backgroundMeshDecomposition.H:140
Foam::backgroundMeshDecomposition::bFTreePtr_
autoPtr< indexedOctree< treeDataBPatch > > bFTreePtr_
Search tree for the boundaryFaces_ patch.
Definition: backgroundMeshDecomposition.H:120
mapPolyMesh.H
globalIndex.H
Foam::backgroundMeshDecomposition::volRes_
label volRes_
How fine should the initial sample of the volume a box be to.
Definition: backgroundMeshDecomposition.H:144
Foam::treeBoundBox
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:75
Foam::backgroundMeshDecomposition::intersectsProcessors
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.
polyTopoChange.H
indexedOctree.H
face.H
Foam::backgroundMeshDecomposition::runTime_
const Time & runTime_
Method details dictionary.
Definition: backgroundMeshDecomposition.H:100
Foam::backgroundMeshDecomposition::boundaryFacesPtr_
autoPtr< bPatch > boundaryFacesPtr_
Patch containing an independent representation of the surface of the.
Definition: backgroundMeshDecomposition.H:117
Foam::backgroundMeshDecomposition::backgroundMeshDecomposition
backgroundMeshDecomposition(const backgroundMeshDecomposition &)
Disallow default bitwise copy construct.
Foam::backgroundMeshDecomposition::buildPatchAndTree
void buildPatchAndTree()
Build the surface patch and search tree.
decompositionMethod.H
removeCells.H
Foam::backgroundMeshDecomposition::spanScale_
scalar spanScale_
Scale of a cell span vs cell size used to decide to refine a cell.
Definition: backgroundMeshDecomposition.H:133
Foam::backgroundMeshDecomposition::mesh_
fvMesh mesh_
Mesh stored on for this processor, specifiying the domain that it.
Definition: backgroundMeshDecomposition.H:110
Foam::conformationSurfaces
Definition: conformationSurfaces.H:53
Foam::backgroundMeshDecomposition::geometryToConformTo_
const conformationSurfaces & geometryToConformTo_
Reference to surface.
Definition: backgroundMeshDecomposition.H:103
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
mapDistributePolyMesh.H
Foam::backgroundMeshDecomposition::findLineAny
pointIndexHit findLineAny(const point &start, const point &end) const
Find any intersection of line between start and end, (exposing.
volumeType.H
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
Foam::label
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
Foam::volumeType
Definition: volumeType.H:54
labelList.H
Foam::backgroundMeshDecomposition::mesh
const fvMesh & mesh() const
Return access to the underlying mesh.
Definition: backgroundMeshDecompositionI.H:28
treeBoundBox.H
Foam::backgroundMeshDecomposition::mergeDist_
scalar mergeDist_
merge distance required by fvMeshDistribute
Definition: backgroundMeshDecomposition.H:130
hexRef8.H
Foam::backgroundMeshDecomposition::refineCell
bool refineCell(label cellI, volumeType volType, scalar &weightEstimate) const
Estimate the number of vertices that will be in this cell, returns.
Foam::indexedOctree
Non-pointer based hierarchical recursive searching.
Definition: treeDataTriSurface.H:47
treeDataPrimitivePatch.H
Foam::backgroundMeshDecomposition::overlapsOtherProcessors
bool overlapsOtherProcessors(const point &centre, const scalar &radiusSqr) const
Foam::backgroundMeshDecomposition::ClassName
ClassName("backgroundMeshDecomposition")
Runtime type information.
Foam::backgroundMeshDecomposition::processorNearestPosition
labelList processorNearestPosition(const List< point > &pts) const
What is the nearest processor to the given position?
Foam::backgroundMeshDecomposition::printMeshData
void printMeshData(const polyMesh &mesh) const
Print details of the decomposed mesh.
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::backgroundMeshDecomposition::buildMap
static autoPtr< mapDistribute > buildMap(const List< label > &toProc)
Build a mapDistribute for the supplied destination processor data.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::backgroundMeshDecomposition::globalBackgroundBounds_
treeBoundBox globalBackgroundBounds_
The overall bounds of all of the background meshes, used to test if.
Definition: backgroundMeshDecomposition.H:127
backgroundMeshDecompositionI.H
Foam::backgroundMeshDecomposition::processorPosition
labelList processorPosition(const List< PointType > &pts) const
What processor is the given position on?
Foam::backgroundMeshDecomposition::overlapProcessors
labelList overlapProcessors(const point &centre, const scalar radiusSqr) const
pointField.H
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
fvMeshDistribute.H
Foam::hexRef8
Refinement of (split) hexes using polyTopoChange.
Definition: hexRef4.H:64
Foam::backgroundMeshDecomposition::findLine
pointIndexHit findLine(const point &start, const point &end) const
Find nearest intersection of line between start and end, (exposing.
Foam::treeDataBPatch
treeDataPrimitivePatch< bPatch > treeDataBPatch
Definition: backgroundMeshDecomposition.H:82
Foam::bPatch
PrimitivePatch< face, List, const pointField > bPatch
Definition: pairPatchAgglomeration.H:45
Foam::Vector< scalar >
Foam::backgroundMeshDecomposition::procBounds
const treeBoundBox & procBounds() const
Return the boundBox of this processor.
Definition: backgroundMeshDecompositionI.H:42
Foam::treeDataPrimitivePatch
Encapsulation of data needed to search on PrimitivePatches.
Definition: treeDataPrimitivePatch.H:61
Foam::List< treeBoundBox >
Foam::backgroundMeshDecomposition::tree
const indexedOctree< treeDataBPatch > & tree() const
Return access to the underlying tree.
Definition: backgroundMeshDecompositionI.H:35
Foam::backgroundMeshDecomposition
Store a background polyMesh to use for the decomposition of space and queries for parallel conformalV...
Definition: backgroundMeshDecomposition.H:92
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::backgroundMeshDecomposition::pointLevel
const labelList & pointLevel() const
Return the point level of the underlying mesh.
Definition: backgroundMeshDecompositionI.H:54
Foam::backgroundMeshDecomposition::distribute
autoPtr< mapDistributePolyMesh > distribute(volScalarField &cellWeights)
Redistribute the background mesh based on a supplied weight field,.
Foam::backgroundMeshDecomposition::minCellSizeLimit_
scalar minCellSizeLimit_
Smallest minimum cell size allowed, i.e. to avoid high initial.
Definition: backgroundMeshDecomposition.H:137
Foam::backgroundMeshDecomposition::distributePoints
autoPtr< mapDistribute > distributePoints(List< PointType > &points) const
Distribute supplied the points to the appropriate processor.
cellSet.H
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::backgroundMeshDecomposition::initialRefinement
void initialRefinement()
Foam::backgroundMeshDecomposition::cellLevel
const labelList & cellLevel() const
Return the cell level of the underlying mesh.
Definition: backgroundMeshDecompositionI.H:48
CGALTriangulation3Ddefs.H
CGAL data structures used for 3D Delaunay meshing.
Foam::backgroundMeshDecomposition::overlapsThisProcessor
bool overlapsThisProcessor(const treeBoundBox &box) const
Does the given box overlap the faces of the boundary of this.
Foam::backgroundMeshDecomposition::allBackgroundMeshBounds_
treeBoundBoxList allBackgroundMeshBounds_
The bounds of all background meshes on all processors.
Definition: backgroundMeshDecomposition.H:123
rndGen
cachedRandom rndGen(label(0), -1)
Foam::backgroundMeshDecomposition::maxCellWeightCoeff_
scalar maxCellWeightCoeff_
Allowed factor above the average cell weight before a background.
Definition: backgroundMeshDecomposition.H:148
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88
Foam::backgroundMeshDecomposition::positionOnThisProcessor
bool positionOnThisProcessor(const point &pt) const
Is the given position inside the domain of this decomposition.