isoSurfaceBase.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) 2019-2022 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::isoSurfaceBase
28
29Description
30 Low-level components common to various iso-surface algorithms.
31
32Note
33 The interpolation samplers currently require a volField for the cell
34 values. This is largely a restriction imposed by the point algorithm
35 and may be revised in the future.
36
37SourceFiles
38 isoSurfaceBase.C
39 isoSurfaceBaseNew.C
40
41\*---------------------------------------------------------------------------*/
42
43#ifndef Foam_isoSurfaceBase_H
44#define Foam_isoSurfaceBase_H
45
46#include "isoSurfaceParams.H"
47#include "bitSet.H"
48#include "scalarField.H"
49#include "volumeType.H"
50#include "volFieldsFwd.H"
51#include "MeshedSurface.H"
52#include "MeshedSurfacesFwd.H"
53
54// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55
56namespace Foam
57{
58
59// Forward Declarations
60class polyMesh;
61class tetCell;
62
63/*---------------------------------------------------------------------------*\
64 Class isoSurfaceBase Declaration
65\*---------------------------------------------------------------------------*/
68:
69 public meshedSurface,
70 public isoSurfaceParams
71{
72public:
73
74 // Data Types
75
76 //- The type of cell/face cuts
77 enum cutType : uint8_t
78 {
79 NOTCUT = 0,
80 CUT = 0x1,
81 TETCUT = 0x2,
82 SPHERE = 0x4,
83 ANYCUT = 0xF,
85 UNVISITED = 0x10,
86 BLOCKED = 0x20,
87 SPECIAL = 0xF0
88 };
89
90
91protected:
92
93 // Protected typedefs for convenience
94 typedef meshedSurface Mesh;
95
96 // Typedef for code transition
97 typedef cutType cellCutType;
98
99
100 // Protected Data
101
102 //- Reference to mesh
103 const polyMesh& mesh_;
104
105 //- Cell values
106 const scalarField& cVals_;
107
108 //- Point values
109 const scalarField& pVals_;
110
111 //- Iso value
112 const scalar iso_;
113
114
115 // Controls, restrictions
116
117 //- Optional boundary faces to ignore.
118 // Eg, Used to exclude cyclicACMI (since duplicate faces)
120
121
122 // Sampling information
123
124 //- For every face, the original cell in mesh
126
127
128 // Protected Member Functions
129
130 //- Count the number of cuts matching the mask type
131 // Checks as bitmask or as zero.
132 static label countCutType
133 (
134 const UList<cutType>& cuts,
135 const uint8_t maskValue
136 );
137
138 //- Dummy templated interpolate method
139 template<class Type>
141 (
144 ) const
145 {
146 return nullptr;
147 }
148
149 //- No copy construct
150 isoSurfaceBase(const isoSurfaceBase&) = delete;
151
152 //- No copy assignment
153 void operator=(const isoSurfaceBase&) = delete;
154
155
156public:
157
158 // Typedefs for code transition
161
162
163 // Constructors
164
165 //- Construct with mesh, cell/point values and iso-value
167 (
168 const polyMesh& mesh,
169 const scalarField& cellValues,
171 const scalar iso,
172 const isoSurfaceParams& params = isoSurfaceParams()
173 );
174
175
176 // Selector
177
178 //- Create for specified algorithm type
179 // Currently uses hard-code lookups based in isoSurfaceParams
181 (
182 const isoSurfaceParams& params,
185 const scalar iso,
186 const bitSet& ignoreCells = bitSet()
187 );
188
189
190 // Member Functions
191
192 // Access, Edit
193
194 //- The mesh for which the iso-surface is associated
195 const polyMesh& mesh() const noexcept
196 {
197 return mesh_;
198 }
199
200 //- The mesh cell values used for creating the iso-surface
201 const scalarField& cellValues() const noexcept
202 {
203 return cVals_;
204 }
205
206 //- The mesh point values used for creating the iso-surface
207 const scalarField& pointValues() const noexcept
208 {
209 return pVals_;
210 }
211
212 //- The iso-value associated with the surface
213 scalar isoValue() const noexcept
214 {
215 return iso_;
216 }
217
218 //- For each face, the original cell in mesh
219 const labelList& meshCells() const noexcept
220 {
221 return meshCells_;
222 }
223
224 //- For each face, the original cell in mesh
226 {
227 return meshCells_;
228 }
229
230
231 // Helpers
232
233 //- Restore non-BLOCKED state to an UNVISITED state
234 static void resetCuts(UList<cutType>& cuts);
235
236 //- Mark ignoreCells as BLOCKED
237 label blockCells
238 (
239 UList<cutType>& cuts,
240 const bitSet& ignoreCells
241 ) const;
242
243 //- Mark cells inside/outside a (valid) bound box as BLOCKED
244 // The volType is INSIDE or OUTSIDE only
245 label blockCells
246 (
247 UList<cutType>& cuts,
248 const boundBox& bb,
249 const volumeType::type volType
250 ) const;
251
252
253 // Cutting
254
255 //- Set ignoreBoundaryFaces to ignore cyclics (cyclicACMI)
256 void ignoreCyclics();
257
258 //- Populate a list of candidate cell cuts using getCellCutType()
259 label calcCellCuts(List<cutType>& cuts) const;
260
261 //- Determine face cut for an individual face
262 cutType getFaceCutType(const label facei) const;
263
264 //- Cell cut for an individual cell, with special handling
265 //- for TETCUT and SPHERE cuts
266 cutType getCellCutType(const label celli) const;
267
268
269 // Sampling
270
271#undef declareIsoSurfaceInterpolateMethod
272#define declareIsoSurfaceInterpolateMethod(Type) \
273 \
274 virtual tmp<Field<Type>> \
275 interpolate \
276 ( \
277 const VolumeField<Type>& cellValues, \
278 const Field<Type>& pointValues \
279 ) const;
286};
287
288
289// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
290
291} // End namespace Foam
292
293// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
294
295#endif
296
297// ************************************************************************* //
Generic GeometricField class.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:94
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
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
Low-level components common to various iso-surface algorithms.
static void resetCuts(UList< cutType > &cuts)
Restore non-BLOCKED state to an UNVISITED state.
bitSet ignoreBoundaryFaces_
Optional boundary faces to ignore.
isoSurfaceBase(const isoSurfaceBase &)=delete
No copy construct.
static autoPtr< isoSurfaceBase > New(const isoSurfaceParams &params, const volScalarField &cellValues, const scalarField &pointValues, const scalar iso, const bitSet &ignoreCells=bitSet())
Create for specified algorithm type.
static label countCutType(const UList< cutType > &cuts, const uint8_t maskValue)
Count the number of cuts matching the mask type.
cutType
The type of cell/face cuts.
@ BLOCKED
Blocked (never cut)
@ ANYCUT
Any cut type (bitmask)
@ SPECIAL
Bitmask for specials.
@ SPHERE
All edges to cell centre cut.
@ TETCUT
Cell cut is a tet.
const scalarField & cellValues() const noexcept
The mesh cell values used for creating the iso-surface.
const scalarField & pVals_
Point values.
scalar isoValue() const noexcept
The iso-value associated with the surface.
label calcCellCuts(List< cutType > &cuts) const
Populate a list of candidate cell cuts using getCellCutType()
labelList meshCells_
For every face, the original cell in mesh.
const scalar iso_
Iso value.
void operator=(const isoSurfaceBase &)=delete
No copy assignment.
cutType getFaceCutType(const label facei) const
Determine face cut for an individual face.
void ignoreCyclics()
Set ignoreBoundaryFaces to ignore cyclics (cyclicACMI)
labelList & meshCells() noexcept
For each face, the original cell in mesh.
const polyMesh & mesh() const noexcept
The mesh for which the iso-surface is associated.
meshedSurface Mesh
cutType getCellCutType(const label celli) const
const labelList & meshCells() const noexcept
For each face, the original cell in mesh.
const scalarField & pointValues() const noexcept
The mesh point values used for creating the iso-surface.
const polyMesh & mesh_
Reference to mesh.
label blockCells(UList< cutType > &cuts, const bitSet &ignoreCells) const
Mark ignoreCells as BLOCKED.
const scalarField & cVals_
Cell values.
tmp< Field< Type > > interpolateTemplate(const VolumeField< Type > &cellValues, const Field< Type > &pointValues) const
Dummy templated interpolate method.
Preferences for controlling iso-surface algorithms.
algorithmType
The algorithm types.
filterType
The filtering (regularization) to apply.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A class for managing temporary objects.
Definition: tmp.H:65
type
Volume classification types.
Definition: volumeType.H:66
#define declareIsoSurfaceInterpolateMethod(Type)
Namespace for OpenFOAM.
const direction noexcept
Definition: Scalar.H:223