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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2018 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::cuttingPlane
29
30Description
31 Constructs cutting plane through a mesh.
32
33 No attempt at resolving degenerate cases.
34 Since the cut faces can be quite ugly, they will often be triangulated.
35
36SourceFiles
37 cuttingPlane.C
38
39\*---------------------------------------------------------------------------*/
40
41#ifndef cuttingPlane_H
42#define cuttingPlane_H
43
44#include "plane.H"
45#include "cuttingSurfaceBase.H"
46
47// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48
49namespace Foam
50{
51
52// Forward declarations
53class polyMesh;
54
55/*---------------------------------------------------------------------------*\
56 Class cuttingPlane Declaration
57\*---------------------------------------------------------------------------*/
59class cuttingPlane
60:
61 public plane,
63{
64 // Private Member Functions
65
66 //- Determine cut cells, possibly restricted to a list of cells
67 //
68 // \param sides [out] For each mesh point, the encoded side of the
69 // plane (0=BACK, 1=ONPLANE, 2=FRONT).
70 // \param cellCuts [in,out] On input an empty set (ie, no restriction)
71 // or subsetted cells. On output, the cells cut according to the
72 // plane-sides detection.
73 //
74 // \return number of faces cut
75 label calcCellCuts
76 (
77 const primitiveMesh& mesh,
78 PackedList<2>& sides,
80 );
81
82
83protected:
84
85 // Constructors
86
87 //- Construct from a plane description without any cutting
88 cuttingPlane(const plane& pln);
89
90
91 // Protected Member Functions
92
93 //- Cut mesh, restricted to a list of cells
95
96 //- Cut mesh, restricted to a list of cells
97 // Reclaim memory for cellIdLabels
98 virtual void performCut
99 (
100 const primitiveMesh& mesh,
101 const bool triangulate,
102 bitSet&& cellIdLabels
103 );
104
105
106 //- Check and warn if bounding boxes do not intersect,
107 //- and if the plane does not intersect the bounding boxes
108 void checkOverlap
109 (
110 const word callerName,
111 const boundBox& meshBounds,
112 const boundBox& userBounds
113 ) const;
114
115 //- Define cell selection from bounding-box and zones.
116 //
117 // \param userBounds Optionally user-specified bounding box
118 // \param zoneNames Optionally user-specified zone names
119 // \param warn Check and warn if the plane does not
120 // intersect with the bounds of the mesh (or submesh) or
121 // if the bounding box does not overlap with the mesh (or submesh)
122 //
123 // \return A set of nCells size with the selected cells or an empty
124 // set if no bounding-box or zones were specified.
126 (
127 const polyMesh& mesh,
128 const boundBox& userBounds,
129 const wordRes& zoneNames,
130 const word callerName,
131 const bool warn
132 ) const;
133
134
135public:
136
137 // Constructors
138
139 //- Construct from plane and mesh reference,
140 //- possibly restricted to a list of cells
142 (
143 const plane& pln,
144 const primitiveMesh& mesh,
145 const bool triangulate,
146 const bitSet& cellIdLabels = bitSet()
147 );
148
149 //- Construct from plane and mesh reference,
150 //- possibly restricted to a list of cells
152 (
153 const plane& pln,
154 const primitiveMesh& mesh,
155 const bool triangulate,
156 bitSet&& cellIdLabels
157 );
158
159 //- Construct from plane and mesh reference,
160 //- possibly restricted to a list of cells
162 (
163 const plane& pln,
164 const primitiveMesh& mesh,
165 const bool triangulate,
166 const labelUList& cellIdLabels
167 );
168
169
170 // Member Functions
171
172 //- The plane used.
173 const plane& planeDesc() const
174 {
175 return *this;
176 }
177
178 //- The plane used. Non-const access.
180 {
181 return *this;
182 }
183
184
185 // Member Operators
186
187 //- Copy assignment
188 void operator=(const cuttingPlane& rhs);
189};
190
191
192// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193
194} // End namespace Foam
195
196// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
197
198#endif
199
200// ************************************************************************* //
virtual label triangulate()
Triangulate in-place, returning the number of triangles added.
A dynamic list of packed unsigned integers, with the number of bits per item specified by the <Width>...
Definition: PackedList.H:129
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
Description of cuts across cells.
Definition: cellCuts.H:113
Constructs cutting plane through a mesh.
Definition: cuttingPlane.H:62
void checkOverlap(const word callerName, const boundBox &meshBounds, const boundBox &userBounds) const
bitSet cellSelection(const polyMesh &mesh, const boundBox &userBounds, const wordRes &zoneNames, const word callerName, const bool warn) const
Define cell selection from bounding-box and zones.
virtual void performCut(const primitiveMesh &mesh, const bool triangulate, bitSet &&cellIdLabels)
Cut mesh, restricted to a list of cells.
Definition: cuttingPlane.C:83
void operator=(const cuttingPlane &rhs)
Copy assignment.
Definition: cuttingPlane.C:149
const plane & planeDesc() const
The plane used.
Definition: cuttingPlane.H:172
plane & planeDesc()
The plane used. Non-const access.
Definition: cuttingPlane.H:178
Base for creating a MeshedSurface by performing some type of cell cutting/intersection.
virtual void performCut(const primitiveMesh &mesh, const bool triangulate, const labelUList &cellIdLabels)
Cut mesh, restricted to a list of cells.
Geometric class that creates a 3D plane and can return the intersection point between a line and the ...
Definition: plane.H:95
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:79
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
Namespace for OpenFOAM.