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 -------------------------------------------------------------------------------
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::treeDataCell
28 
29 Description
30  Encapsulation of data needed to search in/for cells. Used to find the
31  cell containing a point (e.g. cell-cell mapping).
32 
33 SourceFiles
34  treeDataCell.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef treeDataCell_H
39 #define treeDataCell_H
40 
41 #include "polyMesh.H"
42 #include "treeBoundBoxList.H"
43 #include "volumeType.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 // Forward declaration of classes
51 template<class Type> class indexedOctree;
52 
53 /*---------------------------------------------------------------------------*\
54  Class treeDataCell Declaration
55 \*---------------------------------------------------------------------------*/
56 
57 class treeDataCell
58 {
59  // Private data
60 
61  const polyMesh& mesh_;
62 
63  //- Subset of cells to work on
64  const labelList cellLabels_;
65 
66  //- Whether to precalculate and store cell bounding box
67  const bool cacheBb_;
68 
69  //- How to decide if point is inside cell
70  const polyMesh::cellDecomposition decompMode_;
71 
72  //- Cell bounding boxes (valid only if cacheBb_)
73  treeBoundBoxList bbs_;
74 
75 
76  // Private Member Functions
77 
78  //- Calculate cell bounding box
79  treeBoundBox calcCellBb(const label celli) const;
80 
81  //- Initialise all member data
82  void update();
83 
84 public:
85 
86 
87  class findNearestOp
88  {
89  const indexedOctree<treeDataCell>& tree_;
90 
91  public:
92 
94 
95  void operator()
96  (
97  const labelUList& indices,
98  const point& sample,
99 
100  scalar& nearestDistSqr,
101  label& minIndex,
102  point& nearestPoint
103  ) const;
104 
105  void operator()
106  (
107  const labelUList& indices,
108  const linePointRef& ln,
109 
110  treeBoundBox& tightest,
111  label& minIndex,
112  point& linePoint,
113  point& nearestPoint
114  ) const;
115  };
116 
117 
118  class findIntersectOp
119  {
120  const indexedOctree<treeDataCell>& tree_;
121 
122  public:
123 
125 
126  bool operator()
127  (
128  const label index,
129  const point& start,
130  const point& end,
131  point& intersectionPoint
132  ) const;
133  };
134 
135 
136  // Declare name of the class and its debug switch
137  ClassName("treeDataCell");
138 
139 
140  // Constructors
141 
142  //- Construct from mesh, copying subset of cells.
144  (
145  const bool cacheBb,
146  const polyMesh& mesh,
147  const labelUList& cellLabels,
149  );
150 
151  //- Construct from mesh, moving subset of cells
153  (
154  const bool cacheBb,
155  const polyMesh& mesh,
158  );
159 
160  //- Construct from mesh. Uses all cells in mesh.
162  (
163  const bool cacheBb,
164  const polyMesh& mesh,
166  );
167 
168 
169  // Member Functions
170 
171  // Access
172 
173  inline const labelList& cellLabels() const
174  {
175  return cellLabels_;
176  }
177 
178  inline const polyMesh& mesh() const
179  {
180  return mesh_;
181  }
182 
184  {
185  return decompMode_;
186  }
187 
188  inline label size() const
189  {
190  return cellLabels_.size();
191  }
192 
193  //- Get representative point cloud for all shapes inside
194  // (one point per shape)
195  pointField shapePoints() const;
196 
197 
198  // Search
199 
200  //- Get type (inside,outside,mixed,unknown) of point w.r.t. surface.
201  // Only makes sense for closed surfaces.
203  (
205  const point&
206  ) const
207  {
209  return volumeType::UNKNOWN;
210  }
211 
212  //- Does (bb of) shape at index overlap bb
213  bool overlaps
214  (
215  const label index,
216  const treeBoundBox& sampleBb
217  ) const;
218 
219  //- Does shape at index contain sample
220  bool contains
221  (
222  const label index,
223  const point& sample
224  ) const;
225 };
226 
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 } // End namespace Foam
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 
235 #endif
236 
237 // ************************************************************************* //
Foam::polyMesh::cellDecomposition
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:100
Foam::treeDataCell::shapePoints
pointField shapePoints() const
Get representative point cloud for all shapes inside.
Definition: treeDataCell.C:159
Foam::treeBoundBox
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:86
Foam::treeDataCell::findNearestOp::findNearestOp
findNearestOp(const indexedOctree< treeDataCell > &tree)
Definition: treeDataCell.C:140
Foam::treeDataCell::mesh
const polyMesh & mesh() const
Definition: treeDataCell.H:177
polyMesh.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::treeDataCell::getVolumeType
volumeType getVolumeType(const indexedOctree< treeDataCell > &, const point &) const
Get type (inside,outside,mixed,unknown) of point w.r.t. surface.
Definition: treeDataCell.H:202
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
volumeType.H
Foam::treeDataCell
Encapsulation of data needed to search in/for cells. Used to find the cell containing a point (e....
Definition: treeDataCell.H:56
Foam::volumeType
An enumeration wrapper for classification of a location as being inside/outside of a volume.
Definition: volumeType.H:60
Foam::Field< vector >
Foam::volumeType::UNKNOWN
Unknown state.
Definition: volumeType.H:67
Foam::treeDataCell::cellLabels
const labelList & cellLabels() const
Definition: treeDataCell.H:172
Foam::indexedOctree
Non-pointer based hierarchical recursive searching.
Definition: treeDataEdge.H:50
Foam::treeDataCell::size
label size() const
Definition: treeDataCell.H:187
Foam::treeDataCell::decompMode
polyMesh::cellDecomposition decompMode() const
Definition: treeDataCell.H:182
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::treeDataCell::overlaps
bool overlaps(const label index, const treeBoundBox &sampleBb) const
Does (bb of) shape at index overlap bb.
Definition: treeDataCell.C:173
Foam::treeDataCell::ClassName
ClassName("treeDataCell")
treeBoundBoxList.H
Foam::treeDataCell::findNearestOp
Definition: treeDataCell.H:86
Foam::treeDataCell::treeDataCell
treeDataCell(const bool cacheBb, const polyMesh &mesh, const labelUList &cellLabels, const polyMesh::cellDecomposition decompMode)
Construct from mesh, copying subset of cells.
Definition: treeDataCell.C:90
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::UList< label >
Foam::treeDataCell::contains
bool contains(const label index, const point &sample) const
Does shape at index contain sample.
Definition: treeDataCell.C:188
Foam::line
A line primitive.
Definition: line.H:53
Foam::ln
bool ln(const fileName &src, const fileName &dst)
Create a softlink. dst should not exist. Returns true if successful.
Definition: MSwindows.C:925
Foam::treeDataCell::findIntersectOp::findIntersectOp
findIntersectOp(const indexedOctree< treeDataCell > &tree)
Definition: treeDataCell.C:149
Foam::treeDataCell::findIntersectOp
Definition: treeDataCell.H:117
sample
Minimal example by using system/controlDict.functions: