TTK
Loading...
Searching...
No Matches
ttkHelloWorld.cpp
Go to the documentation of this file.
1#include <ttkHelloWorld.h>
2
3#include <vtkInformation.h>
4
5#include <vtkDataArray.h>
6#include <vtkDataSet.h>
7#include <vtkObjectFactory.h>
8#include <vtkPointData.h>
9#include <vtkSmartPointer.h>
10
11#include <ttkMacros.h>
12#include <ttkUtils.h>
13
14// A VTK macro that enables the instantiation of this class via ::New()
15// You do not have to modify this
17
31 this->SetNumberOfInputPorts(1);
32 this->SetNumberOfOutputPorts(1);
33}
34
42int ttkHelloWorld::FillInputPortInformation(int port, vtkInformation *info) {
43 if(port == 0) {
44 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
45 return 1;
46 }
47 return 0;
48}
49
65int ttkHelloWorld::FillOutputPortInformation(int port, vtkInformation *info) {
66 if(port == 0) {
68 return 1;
69 }
70 return 0;
71}
72
86int ttkHelloWorld::RequestData(vtkInformation *ttkNotUsed(request),
87 vtkInformationVector **inputVector,
88 vtkInformationVector *outputVector) {
89
90 // Get input object from input vector
91 // Note: has to be a vtkDataSet as required by FillInputPortInformation
92 vtkDataSet *inputDataSet = vtkDataSet::GetData(inputVector[0]);
93 if(!inputDataSet)
94 return 0;
95
96 // Get input array that will be processed
97 //
98 // Note: VTK provides abstract functionality to handle array selections, but
99 // this essential functionality is unfortunately not well documented.
100 // Before you read further, please keep in mind the the TTK developer
101 // team is not responsible for the existing VTK Api ;-)
102 //
103 // In a nutshell, prior to the RequestData execution one has to call
104 //
105 // SetInputArrayToProcess (
106 // int idx,
107 // int port,
108 // int connection,
109 // int fieldAssociation,
110 // const char *name
111 // )
112 //
113 // The parameter 'idx' is often misunderstood: lets say the filter
114 // requires n arrays, then idx enumerates them from 0 to n-1.
115 //
116 // The 'port' is the input port index at which the object is connected
117 // from which we want to get the array.
118 //
119 // The 'connection' is the connection index at that port (we have to
120 // specify this because VTK allows multiple connections at the same
121 // input port).
122 //
123 // The 'fieldAssociation' integer specifies if the array should be taken
124 // from 0: point data, 1: cell data, or 2: field data.
125 //
126 // The final parameter is the 'name' of the array.
127 //
128 // Example: SetInputArrayToProcess(3,1,0,1,"EdgeLength") will store that
129 // for the 3rd array the filter needs the cell data array named
130 // "EdgeLength" that it will retrieve from the vtkDataObject
131 // at input port 1 (first connection). During the RequestData
132 // method one can then actually retrieve the 3rd array it
133 // requires for its computation by calling
134 // GetInputArrayToProcess(3, inputVector)
135 //
136 // If this filter is run within ParaView, then the UI will automatically
137 // call SetInputArrayToProcess (see HelloWorld.xml file).
138 //
139 // During the RequestData execution one can then retrieve an actual
140 // array with the method "GetInputArrayToProcess".
141 vtkDataArray *inputArray = this->GetInputArrayToProcess(0, inputVector);
142 if(!inputArray) {
143 this->printErr("Unable to retrieve input array.");
144 return 0;
145 }
146
147 // To make sure that the selected array can be processed by this filter,
148 // one should also check that the array association and format is correct.
149 if(this->GetInputArrayAssociation(0, inputVector) != 0) {
150 this->printErr("Input array needs to be a point data array.");
151 return 0;
152 }
153 if(inputArray->GetNumberOfComponents() != 1) {
154 this->printErr("Input array needs to be a scalar array.");
155 return 0;
156 }
157
158 // If all checks pass then log which array is going to be processed.
159 this->printMsg("Starting computation...");
160 this->printMsg(" Scalar Array: " + std::string(inputArray->GetName()));
161
162 // Create an output array that has the same data type as the input array
163 // Note: vtkSmartPointers are well documented
164 // (https://vtk.org/Wiki/VTK/Tutorials/SmartPointers)
165 vtkSmartPointer<vtkDataArray> const outputArray
166 = vtkSmartPointer<vtkDataArray>::Take(inputArray->NewInstance());
167 outputArray->SetName(this->OutputArrayName.data()); // set array name
168 outputArray->SetNumberOfComponents(1); // only one component per tuple
169 outputArray->SetNumberOfTuples(inputArray->GetNumberOfTuples());
170
171 // Get ttk::triangulation of the input vtkDataSet (will create one if one does
172 // not exist already).
173 ttk::Triangulation *triangulation
174 = ttkAlgorithm::GetTriangulation(inputDataSet);
175 if(!triangulation)
176 return 0;
177
178 // Precondition the triangulation (e.g., enable fetching of vertex neighbors)
179 this->preconditionTriangulation(triangulation); // implemented in base class
180
181 // Templatize over the different input array data types and call the base code
182 int status = 0; // this integer checks if the base code returns an error
183 ttkVtkTemplateMacro(inputArray->GetDataType(), triangulation->getType(),
184 (status = this->computeAverages<VTK_TT, TTK_TT>(
185 (VTK_TT *)ttkUtils::GetVoidPointer(outputArray),
186 (VTK_TT *)ttkUtils::GetVoidPointer(inputArray),
187 (TTK_TT *)triangulation->getData())));
188
189 // On error cancel filter execution
190 if(status != 1)
191 return 0;
192
193 // Get output vtkDataSet (which was already instantiated based on the
194 // information provided by FillOutputPortInformation)
195 vtkDataSet *outputDataSet = vtkDataSet::GetData(outputVector, 0);
196
197 // make a SHALLOW copy of the input
198 outputDataSet->ShallowCopy(inputDataSet);
199
200 // add to the output point data the computed output array
201 outputDataSet->GetPointData()->AddArray(outputArray);
202
203 // return success
204 return 1;
205}
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
static vtkInformationIntegerKey * SAME_DATA_TYPE_AS_INPUT_PORT()
ttk::Triangulation * GetTriangulation(vtkDataSet *dataSet)
TTK VTK-filter that wraps the ttk::HelloWorld module.
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
int FillInputPortInformation(int port, vtkInformation *info) override
int FillOutputPortInformation(int port, vtkInformation *info) override
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
int preconditionTriangulation(ttk::AbstractTriangulation *triangulation) const
Definition HelloWorld.h:43
Triangulation is a class that provides time and memory efficient traversal methods on triangulations ...
AbstractTriangulation * getData()
Triangulation::Type getType() const
vtkStandardNewMacro(ttkHelloWorld)
#define ttkVtkTemplateMacro(dataType, triangulationType, call)
Definition ttkMacros.h:69
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)