VTK  9.6.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
26
27#ifndef vtkExtractCells_h
28#define vtkExtractCells_h
29
30#include "vtkFiltersCoreModule.h" // For export macro
31#include "vtkSmartPointer.h" // For vtkSmartPointer
33#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
34
35VTK_ABI_NAMESPACE_BEGIN
36class vtkIdList;
37class vtkExtractCellsIdList;
38
40{
41public:
43
47 void PrintSelf(ostream& os, vtkIndent indent) override;
50
57
63
69
71
74 void SetCellIds(const vtkIdType* ptr, vtkIdType numValues);
75 void AddCellIds(const vtkIdType* ptr, vtkIdType numValues);
77
79
85 vtkSetMacro(ExtractAllCells, bool);
86 vtkGetMacro(ExtractAllCells, bool);
87 vtkBooleanMacro(ExtractAllCells, bool);
89
91
96 vtkSetMacro(AssumeSortedAndUniqueIds, bool);
97 vtkGetMacro(AssumeSortedAndUniqueIds, bool);
98 vtkBooleanMacro(AssumeSortedAndUniqueIds, bool);
100
102
107 vtkSetMacro(PassThroughCellIds, bool);
108 vtkGetMacro(PassThroughCellIds, bool);
109 vtkBooleanMacro(PassThroughCellIds, bool);
111
113
118 vtkSetMacro(OutputPointsPrecision, int);
119 vtkGetMacro(OutputPointsPrecision, int);
121
123
131 vtkSetClampMacro(BatchSize, unsigned int, 1, VTK_INT_MAX);
132 vtkGetMacro(BatchSize, unsigned int);
134protected:
137
139 int FillInputPortInformation(int port, vtkInformation* info) override;
140 bool Copy(vtkDataSet* input, vtkUnstructuredGrid* output);
141
143 bool ExtractAllCells = false;
147 unsigned int BatchSize = 1000;
148
149private:
150 vtkExtractCells(const vtkExtractCells&) = delete;
151 void operator=(const vtkExtractCells&) = delete;
152};
153
154VTK_ABI_NAMESPACE_END
155#endif
abstract class to specify dataset behavior
Definition vtkDataSet.h:57
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.
unsigned int BatchSize
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
static vtkExtractCells * New()
Standard methods for construction, type info, and printing.
list of point or cell ids
Definition vtkIdList.h:24
a simple class to control print indentation
Definition vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Hold a reference to a vtkObjectBase instance.
dataset represents arbitrary combinations of all possible cell types
int vtkIdType
Definition vtkType.h:368
#define VTK_INT_MAX
Definition vtkType.h:197
#define VTK_MARSHALAUTO