VTK  9.3.20240424
vtkAppendFilter.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
96#ifndef vtkAppendFilter_h
97#define vtkAppendFilter_h
98
99#include "vtkFiltersCoreModule.h" // For export macro
101
102VTK_ABI_NAMESPACE_BEGIN
105
106class VTKFILTERSCORE_EXPORT vtkAppendFilter : public vtkUnstructuredGridAlgorithm
107{
108public:
111 void PrintSelf(ostream& os, vtkIndent indent) override;
112
114
118 vtkDataSet* GetInput() { return this->GetInput(0); }
120
122
127 vtkGetMacro(MergePoints, vtkTypeBool);
128 vtkSetMacro(MergePoints, vtkTypeBool);
129 vtkBooleanMacro(MergePoints, vtkTypeBool);
131
133
140 vtkSetClampMacro(Tolerance, double, 0.0, VTK_DOUBLE_MAX);
141 vtkGetMacro(Tolerance, double);
143
145
150 vtkSetMacro(ToleranceIsAbsolute, bool);
151 vtkGetMacro(ToleranceIsAbsolute, bool);
152 vtkBooleanMacro(ToleranceIsAbsolute, bool);
154
159
165
167
172 vtkSetClampMacro(OutputPointsPrecision, int, SINGLE_PRECISION, DEFAULT_PRECISION);
173 vtkGetMacro(OutputPointsPrecision, int);
175
176protected:
179
180 // Usual data generation method
183 int FillInputPortInformation(int port, vtkInformation* info) override;
184
185 // list of data sets to append together.
186 // Here as a convenience. It is a copy of the input array.
188
189 // If true we will attempt to merge points. Must also not have
190 // ghost cells defined.
192
194 double Tolerance;
195
196 // If true, tolerance is used as is. If false, tolerance is multiplied by
197 // the diagonal of the bounding box of the input.
199
200private:
201 vtkAppendFilter(const vtkAppendFilter&) = delete;
202 void operator=(const vtkAppendFilter&) = delete;
203
204 // Get all input data sets that have points, cells, or both.
205 // Caller must delete the returned vtkDataSetCollection.
206 vtkDataSetCollection* GetNonEmptyInputs(vtkInformationVector** inputVector);
207
208 void AppendArrays(int attributesType, vtkInformationVector** inputVector, vtkIdType* globalIds,
209 vtkUnstructuredGrid* output, vtkIdType totalNumberOfElements);
210};
211
212VTK_ABI_NAMESPACE_END
213#endif
appends one or more datasets together into a single unstructured grid
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
vtkDataSetCollection * InputList
void RemoveInputData(vtkDataSet *in)
Remove a dataset from the list of data to append.
vtkDataSet * GetInput()
Get any input of this filter.
vtkTypeBool MergePoints
int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
~vtkAppendFilter() override
static vtkAppendFilter * New()
vtkDataSetCollection * GetInputList()
Returns a copy of the input array.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkDataSet * GetInput(int idx)
Get any input of this filter.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
represent and manipulate attribute data in a dataset
maintain an unordered list of dataset objects
abstract class to specify dataset behavior
Definition vtkDataSet.h:166
a simple class to control print indentation
Definition vtkIndent.h:108
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 vtkTypeBool
Definition vtkABI.h:64
int vtkIdType
Definition vtkType.h:315
#define VTK_DOUBLE_MAX
Definition vtkType.h:154