cuttingPlane.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::cuttingPlane
26 
27 Description
28  Constructs plane through mesh.
29 
30  No attempt at resolving degenerate cases. Since the cut faces are
31  usually quite ugly, they will always be triangulated.
32 
33 Note
34  When the cutting plane coincides with a mesh face, the cell edge on the
35  positive side of the plane is taken.
36 
37 SourceFiles
38  cuttingPlane.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef cuttingPlane_H
43 #define cuttingPlane_H
44 
45 #include "plane.H"
46 #include "pointField.H"
47 #include "faceList.H"
48 #include "MeshedSurface.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 
55 class primitiveMesh;
56 
57 /*---------------------------------------------------------------------------*\
58  Class cuttingPlane Declaration
59 \*---------------------------------------------------------------------------*/
60 
61 class cuttingPlane
62 :
63  public plane,
64  public MeshedSurface<face>
65 {
66  //- Private typedef for convenience
68 
69 
70  // Private data
71 
72  //- List of cells cut by the plane
73  labelList cutCells_;
74 
75 
76  // Private Member Functions
77 
78  //- Determine cut cells, possibly restricted to a list of cells
79  void calcCutCells
80  (
81  const primitiveMesh&,
82  const scalarField& dotProducts,
83  const labelUList& cellIdLabels = labelUList::null()
84  );
85 
86  //- Determine intersection points (cutPoints).
87  void intersectEdges
88  (
89  const primitiveMesh&,
90  const scalarField& dotProducts,
91  List<label>& edgePoint
92  );
93 
94  //- Walk circumference of cell, starting from startEdgeI crossing
95  // only cut edges. Record cutPoint labels in faceVerts.
96  static bool walkCell
97  (
98  const primitiveMesh&,
99  const labelUList& edgePoint,
100  const label celli,
101  const label startEdgeI,
102  DynamicList<label>& faceVerts
103  );
104 
105  //- Determine cuts for all cut cells.
106  void walkCellCuts
107  (
108  const primitiveMesh& mesh,
109  const bool triangulate,
110  const labelUList& edgePoint
111  );
112 
113 
114 protected:
115 
116  // Constructors
117 
118  //- Construct plane description without cutting
119  cuttingPlane(const plane&);
120 
121 
122  // Protected Member Functions
123 
124  //- Recut mesh with existing planeDesc, restricted to a list of cells
125  void reCut
126  (
127  const primitiveMesh&,
128  const bool triangulate,
129  const labelUList& cellIdLabels = labelUList::null()
130  );
131 
132  //- Remap action on triangulation or cleanup
133  virtual void remapFaces(const labelUList& faceMap);
134 
135 
136 public:
137 
138  // Constructors
139 
140  //- Construct from plane and mesh reference,
141  // possibly restricted to a list of cells
143  (
144  const plane&,
145  const primitiveMesh&,
146  const bool triangulate,
147  const labelUList& cellIdLabels = labelUList::null()
148  );
149 
150 
151  // Member Functions
152 
153  //- Return plane used
154  const plane& planeDesc() const
155  {
156  return static_cast<const plane&>(*this);
157  }
158 
159  //- Return List of cells cut by the plane
160  const labelList& cutCells() const
161  {
162  return cutCells_;
163  }
164 
165  //- Return true or false to question: have any cells been cut?
166  bool cut() const
167  {
168  return cutCells_.size();
169  }
170 
171  //- Sample the cell field
172  template<class Type>
173  tmp<Field<Type>> sample(const Field<Type>&) const;
174 
175  template<class Type>
176  tmp<Field<Type>> sample(const tmp<Field<Type>>&) const;
177 
178 
179  // Member Operators
180 
181  void operator=(const cuttingPlane&);
182 };
183 
184 
185 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186 
187 } // End namespace Foam
188 
189 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190 
191 #ifdef NoRepository
192  #include "cuttingPlaneTemplates.C"
193 #endif
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 #endif
198 
199 // ************************************************************************* //
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
void reCut(const primitiveMesh &, const bool triangulate, const labelUList &cellIdLabels=labelUList::null())
Recut mesh with existing planeDesc, restricted to a list of cells.
Definition: cuttingPlane.C:340
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:74
A surface geometry mesh with zone information, not to be confused with the similarly named surfaceMes...
Definition: MeshedSurface.H:72
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
void operator=(const cuttingPlane &)
Definition: cuttingPlane.C:410
Geometric class that creates a 2D plane and can return the intersection point between a line and the ...
Definition: plane.H:60
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Cutting plane sampling functionality.
Constructs plane through mesh.
Definition: cuttingPlane.H:60
dynamicFvMesh & mesh
virtual void remapFaces(const labelUList &faceMap)
Remap action on triangulation or cleanup.
Definition: cuttingPlane.C:365
virtual label triangulate()
Triangulate in-place, returning the number of triangles added.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:61
static const UList< T > & null()
Return a null UList.
Definition: UListI.H:51
bool cut() const
Return true or false to question: have any cells been cut?
Definition: cuttingPlane.H:165
const labelList & cutCells() const
Return List of cells cut by the plane.
Definition: cuttingPlane.H:159
cuttingPlane(const plane &)
Construct plane description without cutting.
Definition: cuttingPlane.C:387
const plane & planeDesc() const
Return plane used.
Definition: cuttingPlane.H:153
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< Field< Type > > sample(const Field< Type > &) const
Sample the cell field.
Namespace for OpenFOAM.