VTK  9.1.0
vtkExtractCells.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkExtractCells.h
5
6 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7 All rights reserved.
8 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9
10 This software is distributed WITHOUT ANY WARRANTY; without even
11 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12 PURPOSE. See the above copyright notice for more information.
13
14=========================================================================*/
15/*----------------------------------------------------------------------------
16 Copyright (c) Sandia Corporation
17 See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details.
18----------------------------------------------------------------------------*/
19
60#ifndef vtkExtractCells_h
61#define vtkExtractCells_h
62
63#include "vtkFiltersExtractionModule.h" // For export macro
65
66class vtkIdList;
67class vtkExtractCellsSTLCloak;
68
69class VTKFILTERSEXTRACTION_EXPORT vtkExtractCells : public vtkUnstructuredGridAlgorithm
70{
71public:
73
77 void PrintSelf(ostream& os, vtkIndent indent) override;
80
87
93
99
101
104 void SetCellIds(const vtkIdType* ptr, vtkIdType numValues);
105 void AddCellIds(const vtkIdType* ptr, vtkIdType numValues);
107
109
115 vtkSetMacro(ExtractAllCells, bool);
116 vtkGetMacro(ExtractAllCells, bool);
117 vtkBooleanMacro(ExtractAllCells, bool);
119
121
126 vtkSetMacro(AssumeSortedAndUniqueIds, bool);
127 vtkGetMacro(AssumeSortedAndUniqueIds, bool);
128 vtkBooleanMacro(AssumeSortedAndUniqueIds, bool);
130protected:
133
135 int FillInputPortInformation(int port, vtkInformation* info) override;
136 bool Copy(vtkDataSet* input, vtkUnstructuredGrid* output);
137
138 vtkExtractCellsSTLCloak* CellList = nullptr;
139 vtkIdType SubSetUGridCellArraySize = 0;
140 vtkIdType SubSetUGridFacesArraySize = 0;
141 bool ExtractAllCells = false;
142 bool AssumeSortedAndUniqueIds = false;
143
144private:
145 vtkExtractCells(const vtkExtractCells&) = delete;
146 void operator=(const vtkExtractCells&) = delete;
147};
148
149#endif
abstract class to specify dataset behavior
Definition: vtkDataSet.h:166
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.
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:140
a simple class to control print indentation
Definition: vtkIndent.h:113
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
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