VTK  9.3.0
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
27#ifndef vtkAppendFilter_h
28#define vtkAppendFilter_h
29
30#include "vtkFiltersCoreModule.h" // For export macro
32
33VTK_ABI_NAMESPACE_BEGIN
36
37class VTKFILTERSCORE_EXPORT vtkAppendFilter : public vtkUnstructuredGridAlgorithm
38{
39public:
42 void PrintSelf(ostream& os, vtkIndent indent) override;
43
45
49 vtkDataSet* GetInput() { return this->GetInput(0); }
51
53
58 vtkGetMacro(MergePoints, vtkTypeBool);
59 vtkSetMacro(MergePoints, vtkTypeBool);
60 vtkBooleanMacro(MergePoints, vtkTypeBool);
62
64
71 vtkSetClampMacro(Tolerance, double, 0.0, VTK_DOUBLE_MAX);
72 vtkGetMacro(Tolerance, double);
74
76
81 vtkSetMacro(ToleranceIsAbsolute, bool);
82 vtkGetMacro(ToleranceIsAbsolute, bool);
83 vtkBooleanMacro(ToleranceIsAbsolute, bool);
85
90
96
98
103 vtkSetClampMacro(OutputPointsPrecision, int, SINGLE_PRECISION, DEFAULT_PRECISION);
104 vtkGetMacro(OutputPointsPrecision, int);
106
107protected:
110
111 // Usual data generation method
114 int FillInputPortInformation(int port, vtkInformation* info) override;
115
116 // list of data sets to append together.
117 // Here as a convenience. It is a copy of the input array.
119
120 // If true we will attempt to merge points. Must also not have
121 // ghost cells defined.
123
125 double Tolerance;
126
127 // If true, tolerance is used as is. If false, tolerance is multiplied by
128 // the diagonal of the bounding box of the input.
130
131private:
132 vtkAppendFilter(const vtkAppendFilter&) = delete;
133 void operator=(const vtkAppendFilter&) = delete;
134
135 // Get all input data sets that have points, cells, or both.
136 // Caller must delete the returned vtkDataSetCollection.
137 vtkDataSetCollection* GetNonEmptyInputs(vtkInformationVector** inputVector);
138
139 void AppendArrays(int attributesType, vtkInformationVector** inputVector, vtkIdType* globalIds,
140 vtkUnstructuredGrid* output, vtkIdType totalNumberOfElements);
141};
142
143VTK_ABI_NAMESPACE_END
144#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:53
a simple class to control print indentation
Definition vtkIndent.h:29
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