cuttingPlane.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 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
67 
69 
70 
71  // Private Data
72 
73  //- List of cells cut by the plane
74  labelList cutCells_;
75 
76 
77  // Private Member Functions
78 
79  //- Determine cut cells, possibly restricted to a list of cells
80  void calcCutCells
81  (
82  const primitiveMesh&,
83  const scalarField& dotProducts,
84  const labelUList& cellIdLabels = labelUList::null()
85  );
86 
87  //- Determine intersection points (cutPoints).
88  void intersectEdges
89  (
90  const primitiveMesh&,
91  const scalarField& dotProducts,
92  List<label>& edgePoint
93  );
94 
95  //- Walk circumference of cell, starting from startEdgeI crossing
96  // only cut edges. Record cutPoint labels in faceVerts.
97  static bool walkCell
98  (
99  const primitiveMesh&,
100  const labelUList& edgePoint,
101  const label celli,
102  const label startEdgeI,
103  DynamicList<label>& faceVerts
104  );
105 
106  //- Determine cuts for all cut cells.
107  void walkCellCuts
108  (
109  const primitiveMesh& mesh,
110  const bool triangulate,
111  const labelUList& edgePoint
112  );
113 
114 
115 protected:
116 
117  // Constructors
118 
119  //- Construct plane description without cutting
120  cuttingPlane(const plane&);
121 
122 
123  // Protected Member Functions
124 
125  //- Recut mesh with existing planeDesc, restricted to a list of cells
126  void reCut
127  (
128  const primitiveMesh&,
129  const bool triangulate,
130  const labelUList& cellIdLabels = labelUList::null()
131  );
132 
133  //- Remap action on triangulation or cleanup
134  virtual void remapFaces(const labelUList& faceMap);
135 
136 
137 public:
138 
139  // Constructors
140 
141  //- Construct from plane and mesh reference,
142  // possibly restricted to a list of cells
144  (
145  const plane&,
146  const primitiveMesh&,
147  const bool triangulate,
148  const labelUList& cellIdLabels = labelUList::null()
149  );
150 
151 
152  // Member Functions
153 
154  //- Return plane used
155  const plane& planeDesc() const
156  {
157  return static_cast<const plane&>(*this);
158  }
159 
160  //- Return List of cells cut by the plane
161  const labelList& cutCells() const
162  {
163  return cutCells_;
164  }
165 
166  //- Return true or false to question: have any cells been cut?
167  bool cut() const
168  {
169  return cutCells_.size();
170  }
171 
172  //- Sample the cell field
173  template<class Type>
174  tmp<Field<Type>> sample(const Field<Type>&) const;
175 
176  template<class Type>
177  tmp<Field<Type>> sample(const tmp<Field<Type>>&) const;
178 };
179 
180 
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 
183 } // End namespace Foam
184 
185 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186 
187 #ifdef NoRepository
188  #include "cuttingPlaneTemplates.C"
189 #endif
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 #endif
194 
195 // ************************************************************************* //
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:364
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:164
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:389
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:60
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:166
const labelList & cutCells() const
Return List of cells cut by the plane.
Definition: cuttingPlane.H:160
cuttingPlane(const plane &)
Construct plane description without cutting.
Definition: cuttingPlane.C:410
const plane & planeDesc() const
Return plane used.
Definition: cuttingPlane.H:154
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.