VTK  9.5.20250910
vtkExtractCells.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-FileCopyrightText: Copyright (c) Sandia Corporation
3// SPDX-License-Identifier: BSD-3-Clause
4
56#ifndef vtkExtractCells_h
57#define vtkExtractCells_h
58
59#include "vtkFiltersCoreModule.h" // For export macro
60#include "vtkSmartPointer.h" // For vtkSmartPointer
62#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
63
64VTK_ABI_NAMESPACE_BEGIN
65class vtkIdList;
66class vtkExtractCellsIdList;
67
69{
70public:
72
76 void PrintSelf(ostream& os, vtkIndent indent) override;
79
86
92
98
100
103 void SetCellIds(const vtkIdType* ptr, vtkIdType numValues);
104 void AddCellIds(const vtkIdType* ptr, vtkIdType numValues);
106
108
114 vtkSetMacro(ExtractAllCells, bool);
115 vtkGetMacro(ExtractAllCells, bool);
116 vtkBooleanMacro(ExtractAllCells, bool);
118
120
125 vtkSetMacro(AssumeSortedAndUniqueIds, bool);
126 vtkGetMacro(AssumeSortedAndUniqueIds, bool);
127 vtkBooleanMacro(AssumeSortedAndUniqueIds, bool);
129
131
136 vtkSetMacro(PassThroughCellIds, bool);
137 vtkGetMacro(PassThroughCellIds, bool);
138 vtkBooleanMacro(PassThroughCellIds, bool);
140
142
147 vtkSetMacro(OutputPointsPrecision, int);
148 vtkGetMacro(OutputPointsPrecision, int);
150
152
160 vtkSetClampMacro(BatchSize, unsigned int, 1, VTK_INT_MAX);
161 vtkGetMacro(BatchSize, unsigned int);
163protected:
166
168 int FillInputPortInformation(int port, vtkInformation* info) override;
169 bool Copy(vtkDataSet* input, vtkUnstructuredGrid* output);
170
172 bool ExtractAllCells = false;
173 bool AssumeSortedAndUniqueIds = false;
174 bool PassThroughCellIds = true;
175 int OutputPointsPrecision = DEFAULT_PRECISION;
176 unsigned int BatchSize = 1000;
177
178private:
179 vtkExtractCells(const vtkExtractCells&) = delete;
180 void operator=(const vtkExtractCells&) = delete;
181};
182
183VTK_ABI_NAMESPACE_END
184#endif
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
subset a vtkDataSet to create a vtkUnstructuredGrid
void SetCellIds(const vtkIdType *ptr, vtkIdType numValues)
Another way to provide ids using a pointer to vtkIdType array.
void SetCellList(vtkIdList *l)
Set the list of cell IDs that the output vtkUnstructuredGrid will be composed of.
~vtkExtractCells() override
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for construction, type info, and printing.
vtkSmartPointer< vtkExtractCellsIdList > CellList
void AddCellList(vtkIdList *l)
Add the supplied list of cell IDs to those that will be included in the output vtkUnstructuredGrid.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
void AddCellIds(const vtkIdType *ptr, vtkIdType numValues)
Another way to provide ids using a pointer to vtkIdType array.
void AddCellRange(vtkIdType from, vtkIdType to)
Add this range of cell IDs to those that will be included in the output vtkUnstructuredGrid.
bool Copy(vtkDataSet *input, vtkUnstructuredGrid *output)
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
static vtkExtractCells * New()
Standard methods for construction, type info, and printing.
list of point or cell ids
Definition vtkIdList.h:133
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Hold a reference to a vtkObjectBase instance.
Superclass for algorithms that produce only unstructured grid as output.
dataset represents arbitrary combinations of all possible cell types
int vtkIdType
Definition vtkType.h:332
#define VTK_INT_MAX
Definition vtkType.h:161
#define VTK_MARSHALAUTO