VTK  9.2.5
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 
63 #ifndef vtkExtractCells_h
64 #define vtkExtractCells_h
65 
66 #include "vtkFiltersExtractionModule.h" // For export macro
68 
69 class vtkIdList;
70 class vtkExtractCellsSTLCloak;
71 
72 class VTKFILTERSEXTRACTION_EXPORT vtkExtractCells : public vtkUnstructuredGridAlgorithm
73 {
74 public:
76 
80  void PrintSelf(ostream& os, vtkIndent indent) override;
81  static vtkExtractCells* New();
83 
90 
96 
102 
104 
107  void SetCellIds(const vtkIdType* ptr, vtkIdType numValues);
108  void AddCellIds(const vtkIdType* ptr, vtkIdType numValues);
110 
112 
118  vtkSetMacro(ExtractAllCells, bool);
119  vtkGetMacro(ExtractAllCells, bool);
120  vtkBooleanMacro(ExtractAllCells, bool);
122 
124 
129  vtkSetMacro(AssumeSortedAndUniqueIds, bool);
130  vtkGetMacro(AssumeSortedAndUniqueIds, bool);
131  vtkBooleanMacro(AssumeSortedAndUniqueIds, bool);
133 protected:
135  ~vtkExtractCells() override;
136 
139  bool Copy(vtkDataSet* input, vtkUnstructuredGrid* output);
140 
141  vtkExtractCellsSTLCloak* CellList = nullptr;
142  vtkIdType SubSetUGridCellArraySize = 0;
143  vtkIdType SubSetUGridFacesArraySize = 0;
144  bool ExtractAllCells = false;
145  bool AssumeSortedAndUniqueIds = false;
146 
147 private:
148  vtkExtractCells(const vtkExtractCells&) = delete;
149  void operator=(const vtkExtractCells&) = delete;
150 };
151 
152 #endif
abstract class to specify dataset behavior
Definition: vtkDataSet.h:172
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.
static vtkExtractCells * New()
Standard methods for construction, type info, and printing.
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.
list of point or cell ids
Definition: vtkIdList.h:143
a simple class to control print indentation
Definition: vtkIndent.h:119
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
@ info
Definition: vtkX3D.h:382
@ port
Definition: vtkX3D.h:453
int vtkIdType
Definition: vtkType.h:332