VTK  9.3.1
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 
27 #ifndef vtkExtractCells_h
28 #define vtkExtractCells_h
29 
30 #include "vtkFiltersCoreModule.h" // For export macro
31 #include "vtkSmartPointer.h" // For vtkSmartPointer
33 
34 VTK_ABI_NAMESPACE_BEGIN
35 class vtkIdList;
36 class vtkExtractCellsIdList;
37 
38 class VTKFILTERSCORE_EXPORT vtkExtractCells : public vtkUnstructuredGridAlgorithm
39 {
40 public:
42 
46  void PrintSelf(ostream& os, vtkIndent indent) override;
47  static vtkExtractCells* New();
49 
55  void SetCellList(vtkIdList* l);
56 
61  void AddCellList(vtkIdList* l);
62 
67  void AddCellRange(vtkIdType from, vtkIdType to);
68 
70 
73  void SetCellIds(const vtkIdType* ptr, vtkIdType numValues);
74  void AddCellIds(const vtkIdType* ptr, vtkIdType numValues);
76 
78 
84  vtkSetMacro(ExtractAllCells, bool);
85  vtkGetMacro(ExtractAllCells, bool);
86  vtkBooleanMacro(ExtractAllCells, bool);
88 
90 
95  vtkSetMacro(AssumeSortedAndUniqueIds, bool);
96  vtkGetMacro(AssumeSortedAndUniqueIds, bool);
97  vtkBooleanMacro(AssumeSortedAndUniqueIds, bool);
99 
101 
106  vtkSetMacro(PassThroughCellIds, bool);
107  vtkGetMacro(PassThroughCellIds, bool);
108  vtkBooleanMacro(PassThroughCellIds, bool);
110 
112 
117  vtkSetMacro(OutputPointsPrecision, int);
118  vtkGetMacro(OutputPointsPrecision, int);
120 
122 
130  vtkSetClampMacro(BatchSize, unsigned int, 1, VTK_INT_MAX);
131  vtkGetMacro(BatchSize, unsigned int);
133 protected:
134  vtkExtractCells();
135  ~vtkExtractCells() override;
136 
138  int FillInputPortInformation(int port, vtkInformation* info) override;
139  bool Copy(vtkDataSet* input, vtkUnstructuredGrid* output);
140 
142  bool ExtractAllCells = false;
143  bool AssumeSortedAndUniqueIds = false;
144  bool PassThroughCellIds = true;
145  int OutputPointsPrecision = DEFAULT_PRECISION;
146  unsigned int BatchSize = 1000;
147 
148 private:
149  vtkExtractCells(const vtkExtractCells&) = delete;
150  void operator=(const vtkExtractCells&) = delete;
151 };
152 
153 VTK_ABI_NAMESPACE_END
154 #endif
virtual int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
This is called by the superclass.
Store vtkAlgorithm input/output information.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:52
#define VTK_INT_MAX
Definition: vtkType.h:144
static vtkUnstructuredGridAlgorithm * New()
vtkSmartPointer< vtkExtractCellsIdList > CellList
int vtkIdType
Definition: vtkType.h:315
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
a simple class to control print indentation
Definition: vtkIndent.h:28
list of point or cell ids
Definition: vtkIdList.h:22
dataset represents arbitrary combinations of all possible cell types
Superclass for algorithms that produce only unstructured grid as output.
Store zero or more vtkInformation instances.
subset a vtkDataSet to create a vtkUnstructuredGrid