treeDataCell.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  Copyright (C) 2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::treeDataCell
29 
30 Description
31  Encapsulation of data needed to search in/for cells. Used to find the
32  cell containing a point (e.g. cell-cell mapping).
33 
34 SourceFiles
35  treeDataCell.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef Foam_treeDataCell_H
40 #define Foam_treeDataCell_H
41 
42 #include "polyMesh.H"
43 #include "treeBoundBoxList.H"
44 #include "volumeType.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 // Forward Declarations
52 template<class Type> class indexedOctree;
53 
54 /*---------------------------------------------------------------------------*\
55  Class treeDataCell Declaration
56 \*---------------------------------------------------------------------------*/
57 
58 class treeDataCell
59 {
60  // Private Data
61 
62  //- Reference to the mesh
63  const polyMesh& mesh_;
64 
65  //- Subset of cells to work on
66  const labelList cellLabels_;
67 
68  //- Use subset of cells (cellLabels)
69  const bool useSubset_;
70 
71  //- Whether to precalculate and store cell bounding box
72  const bool cacheBb_;
73 
74  //- How to decide if point is inside cell
75  const polyMesh::cellDecomposition decompMode_;
76 
77  //- Cell bounding boxes (valid only if cacheBb_)
78  treeBoundBoxList bbs_;
79 
80 
81  // Private Member Functions
82 
83  //- Get cell bounding box at specified index
84  inline treeBoundBox getBounds(const label index) const;
85 
86  //- Initialise all member data
87  void update();
88 
89 public:
90 
91 
92  class findNearestOp
93  {
94  const indexedOctree<treeDataCell>& tree_;
95 
96  public:
97 
99 
100  void operator()
101  (
102  const labelUList& indices,
103  const point& sample,
104 
105  scalar& nearestDistSqr,
106  label& minIndex,
107  point& nearestPoint
108  ) const;
109 
110  void operator()
111  (
112  const labelUList& indices,
113  const linePointRef& ln,
114 
115  treeBoundBox& tightest,
116  label& minIndex,
117  point& linePoint,
118  point& nearestPoint
119  ) const;
120  };
121 
122 
123  class findIntersectOp
124  {
125  const indexedOctree<treeDataCell>& tree_;
126 
127  public:
128 
130 
131  bool operator()
132  (
133  const label index,
134  const point& start,
135  const point& end,
136  point& intersectionPoint
137  ) const;
138  };
139 
140 
141  // Declare name of the class
142  ClassNameNoDebug("treeDataCell");
143 
144 
145  // Constructors (cachable versions)
146 
147  //- Construct from mesh, using all cells in mesh.
149  (
150  const bool cacheBb,
151  const polyMesh& mesh,
153  );
154 
155  //- Construct from mesh, copying subset of cells.
157  (
158  const bool cacheBb,
159  const polyMesh& mesh,
160  const labelUList& cellLabels,
162  );
163 
164  //- Construct from mesh, moving subset of cells
166  (
167  const bool cacheBb,
168  const polyMesh& mesh,
171  );
172 
173 
174  // Constructors (non-caching)
175 
176  //- Construct from mesh, using all cells in mesh.
178  (
179  const polyMesh& mesh,
181  )
182  :
183  treeDataCell(false, mesh, decompMode)
184  {}
185 
186  //- Construct from mesh, copying subset of cells.
188  (
189  const polyMesh& mesh,
190  const labelUList& cellLabels,
192  )
193  :
195  {}
196 
197  //- Construct from mesh, moving subset of cells
199  (
200  const polyMesh& mesh,
203  )
204  :
205  treeDataCell(false, mesh, std::move(cellLabels), decompMode)
206  {}
207 
208 
209  // Static Functions
210 
211  //- Calculate and return bounding boxes for all mesh cells
212  static treeBoundBoxList boxes
213  (
214  const primitiveMesh& mesh
215  );
216 
217  //- Calculate and return bounding boxes for specified mesh cells
218  static treeBoundBoxList boxes
219  (
220  const primitiveMesh& mesh,
221  const labelUList& cellIds
222  );
223 
224  //- Return bounding box of specified mesh cells
225  static treeBoundBox bounds
226  (
227  const primitiveMesh& mesh,
228  const labelUList& cellIds
229  );
230 
231 
232  // Member Functions
233 
234  //- Object dimension == 3 (volume element)
235  int nDim() const noexcept { return 3; }
236 
237  //- Return bounding box for the specified cell indices
238  treeBoundBox bounds(const labelUList& indices) const;
239 
240 
241  // Access
242 
243  //- Reference to the supporting mesh
244  const polyMesh& mesh() const noexcept { return mesh_; }
245 
246  //- The cell decomposition mode used
248  {
249  return decompMode_;
250  }
251 
252  //- The subset of cell ids to use
253  const labelList& cellLabels() const noexcept { return cellLabels_; }
254 
255  //- Use a subset of cells
256  bool useSubset() const noexcept { return useSubset_; }
257 
258  //- Is the effective cell selection empty?
259  bool empty() const noexcept
260  {
261  return useSubset_ ? cellLabels_.empty() : !mesh_.nCells();
262  }
263 
264  //- The size of the cell selection
265  label size() const noexcept
266  {
267  return useSubset_ ? cellLabels_.size() : mesh_.nCells();
268  }
269 
270  //- Map from shape index to original (non-subset) cell label
271  label objectIndex(const label index) const
272  {
273  return useSubset_ && index >= 0 ? cellLabels_[index] : index;
274  }
275 
276  //TBD //- Cell at specified shape index
280  //- Representative point (cell centre) at shape index
281  const point& centre(const label index) const
282  {
283  return mesh_.cellCentres()[objectIndex(index)];
284  }
285 
286  //- Representative point cloud for contained shapes.
287  //- One point per shape, corresponding to the cell centres.
288  tmp<pointField> centres() const;
289 
290 
291  // Search
293  //- Get type (inside,outside,mixed,unknown) of point w.r.t. surface.
294  // Only makes sense for closed surfaces.
296  (
298  const point&
299  ) const
300  {
303  }
304 
305  //- Does (bb of) shape at index overlap searchBox
306  bool overlaps
307  (
308  const label index,
309  const treeBoundBox& searchBox
310  ) const;
311 
312  //- Does shape at index contain sample
313  bool contains
314  (
315  const label index,
316  const point& sample
317  ) const;
319  //- Calculates nearest (to sample) point in shape.
320  // Returns actual point and distance (squared)
321  void findNearest
322  (
323  const labelUList& indices,
324  const point& sample,
325 
326  scalar& nearestDistSqr,
327  label& nearestIndex,
328  point& nearestPoint
329  ) const;
330 };
331 
332 
333 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
334 
335 } // End namespace Foam
336 
337 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
338 
339 
340 #endif
341 
342 // ************************************************************************* //
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
label size() const noexcept
The size of the cell selection.
Definition: treeDataCell.H:310
A line primitive.
Definition: line.H:52
tmp< pointField > centres() const
Representative point cloud for contained shapes. One point per shape, corresponding to the cell centr...
Definition: treeDataCell.C:246
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:75
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
const polyMesh & mesh() const noexcept
Reference to the supporting mesh.
Definition: treeDataCell.H:279
bool overlaps(const label index, const treeBoundBox &searchBox) const
Does (bb of) shape at index overlap searchBox.
Definition: treeDataCell.C:258
An enumeration wrapper for classification of a location as being inside/outside of a volume...
Definition: volumeType.H:55
static treeBoundBoxList boxes(const primitiveMesh &mesh)
Calculate and return bounding boxes for all mesh cells.
Definition: treeDataCell.C:85
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:103
Unknown state.
Definition: volumeType.H:64
ClassNameNoDebug("treeDataCell")
findIntersectOp(const indexedOctree< treeDataCell > &tree)
Definition: treeDataCell.C:294
const labelList & cellLabels() const noexcept
The subset of cell ids to use.
Definition: treeDataCell.H:292
Tree tree(triangles.begin(), triangles.end())
volumeType getVolumeType(const indexedOctree< treeDataCell > &, const point &) const
Get type (inside,outside,mixed,unknown) of point w.r.t. surface.
Definition: treeDataCell.H:350
static treeBoundBox bounds(const primitiveMesh &mesh, const labelUList &cellIds)
Return bounding box of specified mesh cells.
Definition: treeDataCell.C:105
void findNearest(const labelUList &indices, const point &sample, scalar &nearestDistSqr, label &nearestIndex, point &nearestPoint) const
Calculates nearest (to sample) point in shape.
Definition: treeDataCell.C:303
const vectorField & cellCentres() const
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: POSIX.C:1237
const point & centre(const label index) const
Representative point (cell centre) at shape index.
Definition: treeDataCell.H:330
const direction noexcept
Definition: Scalar.H:258
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
mesh update()
Encapsulation of data needed to search in/for cells. Used to find the cell containing a point (e...
Definition: treeDataCell.H:53
treeDataCell(const bool cacheBb, const polyMesh &mesh, polyMesh::cellDecomposition decompMode)
Construct from mesh, using all cells in mesh.
Definition: treeDataCell.C:159
bool empty() const noexcept
Is the effective cell selection empty?
Definition: treeDataCell.H:302
polyMesh::cellDecomposition decompMode() const noexcept
The cell decomposition mode used.
Definition: treeDataCell.H:284
findNearestOp(const indexedOctree< treeDataCell > &tree)
Definition: treeDataCell.C:285
Non-pointer based hierarchical recursive searching.
label objectIndex(const label index) const
Map from shape index to original (non-subset) cell label.
Definition: treeDataCell.H:318
label nCells() const noexcept
Number of mesh cells.
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:90
int nDim() const noexcept
Object dimension == 3 (volume element)
Definition: treeDataCell.H:266
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
bool contains(const label index, const point &sample) const
Does shape at index contain sample.
Definition: treeDataCell.C:273
Namespace for OpenFOAM.
bool useSubset() const noexcept
Use a subset of cells.
Definition: treeDataCell.H:297