cuttingSurfaceBase.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) 2018-2020 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::cuttingSurfaceBase
28
29Description
30 Base for creating a MeshedSurface by performing some type of cell
31 cutting/intersection.
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 cuttingSurfaceBase.C
38
39\*---------------------------------------------------------------------------*/
40
41#ifndef cuttingSurfaceBase_H
42#define cuttingSurfaceBase_H
43
44#include "bitSet.H"
45#include "faceList.H"
46#include "MeshedSurface.H"
47#include "MeshedSurfacesFwd.H"
48
49// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50
51namespace Foam
52{
53
54// Forward Declarations
55class primitiveMesh;
56
57/*---------------------------------------------------------------------------*\
58 Class cuttingSurfaceBase Declaration
59\*---------------------------------------------------------------------------*/
62:
63 public meshedSurface
64{
65protected:
66
67 //- Typedef for convenience
68 typedef meshedSurface Mesh;
69
70
71 // Protected Data
72
73 //- List of the cells cut
75
76
77 // Protected Member Functions
78
79 //- Walk cell cuts to create faces
80 //
81 // \tparam EdgeOrientIntersect
82 // Parameter (edge&), returns bool.
83 // Orient edge for a consistent positive gradient.
84 // Checks for edge intersection (true|false).
85 //
86 // \tparam EdgeAlphaIntersect
87 // Parameter (const edge&), returns scalar.
88 // Determine alpha [0-1] for an intersecting edge.
89 // No guarantees when used with non-intersecting edges.
90 //
91 // \param cellCuts [in] The cells to walk.
92 template<class EdgeOrientIntersect, class EdgeAlphaIntersect>
93 void walkCellCuts
94 (
95 const primitiveMesh& mesh,
96 const bitSet& cellCuts,
97 const EdgeOrientIntersect& edgeOrientIntersect,
98 const EdgeAlphaIntersect& edgeAlphaIntersect,
99 const bool triangulate,
100 label nFaceCuts = 0
101 );
102
103
104 //- Cut mesh, restricted to a list of cells
105 virtual void performCut
106 (
107 const primitiveMesh& mesh,
108 const bool triangulate,
109 const labelUList& cellIdLabels
110 );
111
112 //- Cut mesh, restricted to a list of cells
113 virtual void performCut
114 (
115 const primitiveMesh& mesh,
116 const bool triangulate,
117 const bitSet& cellSelectionMask = bitSet()
118 );
119
120 //- Cut mesh, restricted to a list of cells
121 // Reclaim memory for cellSelectionMask
122 virtual void performCut
123 (
124 const primitiveMesh& mesh,
125 const bool triangulate,
126 bitSet&& cellSelectionMask
127 ) = 0;
128
129 //- Remap action on triangulation or cleanup
130 virtual void remapFaces(const labelUList& faceMap);
131
132
133 //- Check and warn if bounding boxes do not intersect
134 static void checkOverlap
135 (
136 const word callerName,
137 const boundBox& meshBounds,
138 const boundBox& userBounds
139 );
140
141
142 //- Define cell selection from bounding-box and zones.
143 //
144 // \param userBounds Optionally user-specified bounding box
145 // \param zoneNames Optionally user-specified zone names
146 // \param meshBounds [out] The effective mesh bounds after applying
147 // the user-specified zone names
148 //
149 // \return A set of nCells size with the selected cells or an empty
150 // set if no bounding-box or zones were specified.
151 static bitSet cellSelection
152 (
153 const polyMesh& mesh,
154 const boundBox& userBounds,
155 const wordRes& zoneNames,
156 boundBox& meshBounds
157 );
158
159 //- Define cell selection from bounding-box and zones.
160 //
161 // \param userBounds Optionally user-specified bounding box
162 // \param zoneNames Optionally user-specified zone names
163 // \param callerName The caller name for warnings
164 // \param warn Check and warn if the bounding box does not
165 // overlap with the mesh (or submesh)
166 //
167 // \return A set of nCells size with the selected cells or an empty
168 // set if no bounding-box or zones were specified.
169 static bitSet cellSelection
170 (
171 const polyMesh& mesh,
172 const boundBox& userBounds,
173 const wordRes& zoneNames,
174 const word callerName,
175 const bool warn
176 );
177
178
179public:
180
181 //- Debug information
182 static int debug;
183
184
185 // Constructors
186
187 //- Default construct
188 cuttingSurfaceBase() = default;
189
190
191 //- Destructors
192 virtual ~cuttingSurfaceBase() = default;
193
194
195 // Member Functions
196
197 //- The mesh cells cut
198 const labelList& meshCells() const
199 {
200 return meshCells_;
201 }
202
203 //- The mesh cells cut
205 {
206 return meshCells_;
207 }
208
209 //- Have any cells been cut?
210 bool cut() const
211 {
212 return meshCells_.size();
213 }
214
215
216 // Member Operators
217
218 //- Copy assignment
219 void operator=(const cuttingSurfaceBase& rhs);
220};
221
222// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223
224} // End namespace Foam
225
226// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227
228#ifdef NoRepository
230#endif
231
232// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233
234#endif
235
236// ************************************************************************* //
virtual label triangulate()
Triangulate in-place, returning the number of triangles added.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
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
Base for creating a MeshedSurface by performing some type of cell cutting/intersection.
void walkCellCuts(const primitiveMesh &mesh, const bitSet &cellCuts, const EdgeOrientIntersect &edgeOrientIntersect, const EdgeAlphaIntersect &edgeAlphaIntersect, const bool triangulate, label nFaceCuts=0)
Walk cell cuts to create faces.
static void checkOverlap(const word callerName, const boundBox &meshBounds, const boundBox &userBounds)
Check and warn if bounding boxes do not intersect.
virtual void performCut(const primitiveMesh &mesh, const bool triangulate, const labelUList &cellIdLabels)
Cut mesh, restricted to a list of cells.
bool cut() const
Have any cells been cut?
static bitSet cellSelection(const polyMesh &mesh, const boundBox &userBounds, const wordRes &zoneNames, boundBox &meshBounds)
Define cell selection from bounding-box and zones.
void operator=(const cuttingSurfaceBase &rhs)
Copy assignment.
labelList meshCells_
List of the cells cut.
labelList & meshCells()
The mesh cells cut.
cuttingSurfaceBase()=default
Default construct.
virtual void remapFaces(const labelUList &faceMap)
Remap action on triangulation or cleanup.
static int debug
Debug information.
meshedSurface Mesh
Typedef for convenience.
const labelList & meshCells() const
The mesh cells cut.
virtual ~cuttingSurfaceBase()=default
Destructors.
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.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)