TTK
Loading...
Searching...
No Matches
ttkScalarFieldCriticalPoints.cpp
Go to the documentation of this file.
2
3#include <vtkInformation.h>
4
5#include <vtkDoubleArray.h>
6#include <vtkFloatArray.h>
7#include <vtkIdTypeArray.h>
8#include <vtkIntArray.h>
9#include <vtkNew.h>
10#include <vtkPointData.h>
11#include <vtkPolyData.h>
12#include <vtkSignedCharArray.h>
13
14#include <ttkMacros.h>
15#include <ttkUtils.h>
16
17using namespace std;
18using namespace ttk;
19
21
23
24 this->SetNumberOfInputPorts(1);
25 this->SetNumberOfOutputPorts(1);
26}
27
29
31 int port, vtkInformation *info) {
32 if(port == 0)
33 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
34 else
35 return 0;
36
37 return 1;
38}
39
41 int port, vtkInformation *info) {
42 if(port == 0)
43 info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkPolyData");
44 else
45 return 0;
46
47 return 1;
48}
49
51 vtkInformation *ttkNotUsed(request),
52 vtkInformationVector **inputVector,
53 vtkInformationVector *outputVector) {
54
55 vtkDataSet *input = vtkDataSet::GetData(inputVector[0]);
56 vtkPolyData *output = vtkPolyData::GetData(outputVector, 0);
57
59
60 int const keepGoing = checkEmptyMPIInput<Triangulation>(triangulation);
61 if(keepGoing < 2) {
62 return keepGoing;
63 }
64
65 if(VertexBoundary)
66 triangulation->preconditionBoundaryVertices();
67
68 // in the following, the target scalar field of the input is replaced in the
69 // variable 'output' with the result of the computation.
70 // if your wrapper produces an output of the same type of the input, you
71 // should proceed in the same way.
72
73 vtkDataArray *inputScalarField = this->GetInputArrayToProcess(0, inputVector);
74 if(!inputScalarField)
75 return 0;
76
77 vtkDataArray *offsetField = this->GetOrderArray(
78 input, 0, triangulation, false, 1, ForceInputOffsetScalarField);
79
80 // setting up the base layer
81 this->preconditionTriangulation(triangulation);
82 this->setOutput(&criticalPoints_);
83
84#ifdef TTK_ENABLE_MPI_TIME
85 ttk::Timer t_mpi;
86 ttk::startMPITimer(t_mpi, ttk::MPIrank_, ttk::MPIsize_);
87#endif
88
89 printMsg("Starting computation...");
90 printMsg({{" Scalar Array", inputScalarField->GetName()},
91 {" Offset Array", offsetField ? offsetField->GetName() : "None"}});
92
93 int status = 0;
95 triangulation->getType(),
96 (status = this->execute(
97 static_cast<SimplexId *>(ttkUtils::GetVoidPointer(offsetField)),
98 (TTK_TT *)triangulation->getData())));
99
100 if(status < 0)
101 return 0;
102
103#ifdef TTK_ENABLE_MPI_TIME
104 double elapsedTime = ttk::endMPITimer(t_mpi, ttk::MPIrank_, ttk::MPIsize_);
105 if(ttk::MPIrank_ == 0) {
106 printMsg("Computation performed using " + std::to_string(ttk::MPIsize_)
107 + " MPI processes lasted :" + std::to_string(elapsedTime));
108 }
109#endif
110
111 // allocate the output
112 vtkNew<vtkSignedCharArray> vertexTypes{};
113 vertexTypes->SetNumberOfComponents(1);
114 vertexTypes->SetNumberOfTuples(criticalPoints_.size());
115 vertexTypes->SetName("CriticalType");
116
117 vtkNew<vtkPoints> pointSet{};
118 pointSet->SetNumberOfPoints(criticalPoints_.size());
119
120#ifdef TTK_ENABLE_OPENMP
121#pragma omp parallel for num_threads(this->threadNumber_)
122#endif // TTK_ENABLE_OPENMP
123 for(size_t i = 0; i < criticalPoints_.size(); i++) {
124 std::array<double, 3> p{};
125 input->GetPoint(criticalPoints_[i].first, p.data());
126 pointSet->SetPoint(i, p.data());
127 vertexTypes->SetTuple1(i, (float)criticalPoints_[i].second);
128 }
129
130 ttkUtils::CellVertexFromPoints(output, pointSet);
131 output->GetPointData()->AddArray(vertexTypes);
132
133 if(VertexBoundary) {
134 vtkNew<vtkSignedCharArray> vertexBoundary{};
135 vertexBoundary->SetNumberOfComponents(1);
136 vertexBoundary->SetNumberOfTuples(criticalPoints_.size());
137 vertexBoundary->SetName("IsOnBoundary");
138
139#ifdef TTK_ENABLE_OPENMP
140#pragma omp parallel for num_threads(threadNumber_)
141#endif
142 for(size_t i = 0; i < criticalPoints_.size(); i++) {
143 vertexBoundary->SetTuple1(
144 i, (signed char)triangulation->isVertexOnBoundary(
145 criticalPoints_[i].first));
146 }
147
148 output->GetPointData()->AddArray(vertexBoundary);
149 } else {
150 output->GetPointData()->RemoveArray("IsOnBoundary");
151 }
152
153 if(VertexIds) {
154 vtkNew<ttkSimplexIdTypeArray> vertexIds{};
155 vertexIds->SetNumberOfComponents(1);
156 vertexIds->SetNumberOfTuples(criticalPoints_.size());
157 vertexIds->SetName(ttk::VertexScalarFieldName);
158 for(size_t i = 0; i < criticalPoints_.size(); i++) {
159#ifdef TTK_ENABLE_MPI
160 if(hasInitializedMPI()) {
161 vertexIds->SetTuple1(
162 i, triangulation->getVertexGlobalId(criticalPoints_[i].first));
163 } else
164#endif // TTK_ENABLE_MPI
165 {
166 vertexIds->SetTuple1(i, criticalPoints_[i].first);
167 }
168 }
169
170 output->GetPointData()->AddArray(vertexIds);
171 } else {
172 output->GetPointData()->RemoveArray(ttk::VertexScalarFieldName);
173 }
174
175 if(VertexScalars) {
176 for(SimplexId i = 0; i < input->GetPointData()->GetNumberOfArrays(); i++) {
177
178 vtkDataArray *scalarField = input->GetPointData()->GetArray(i);
179 vtkSmartPointer<vtkDataArray> const scalarArray{
180 scalarField->NewInstance()};
181
182 scalarArray->SetNumberOfComponents(scalarField->GetNumberOfComponents());
183 scalarArray->SetNumberOfTuples(criticalPoints_.size());
184 scalarArray->SetName(scalarField->GetName());
185 std::vector<double> value(scalarField->GetNumberOfComponents());
186 for(size_t j = 0; j < criticalPoints_.size(); j++) {
187 scalarField->GetTuple(criticalPoints_[j].first, value.data());
188 scalarArray->SetTuple(j, value.data());
189 }
190 output->GetPointData()->AddArray(scalarArray);
191 }
192 } else {
193 for(SimplexId i = 0; i < input->GetPointData()->GetNumberOfArrays(); i++) {
194 output->GetPointData()->RemoveArray(
195 input->GetPointData()->GetArray(i)->GetName());
196 }
197 }
198
199 return 1;
200}
#define ttkTemplateMacro(triangulationType, call)
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
ttk::Triangulation * GetTriangulation(vtkDataSet *dataSet)
vtkDataArray * GetOrderArray(vtkDataSet *const inputData, const int scalarArrayIdx, ttk::Triangulation *triangulation, const bool getGlobalOrder=false, const int orderArrayIdx=0, const bool enforceOrderArrayIdx=false)
TTK VTK-filter for the computation of critical points in PL scalar fields defined on PL manifolds.
int FillOutputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
~ttkScalarFieldCriticalPoints() override
int FillInputPortInformation(int port, vtkInformation *info) override
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
static int CellVertexFromPoints(vtkDataSet *const dataSet, vtkPoints *const points)
Definition ttkUtils.cpp:328
void preconditionTriangulation(AbstractTriangulation *triangulation)
void setOutput(std::vector< std::pair< SimplexId, char > > *criticalPoints)
Triangulation is a class that provides time and memory efficient traversal methods on triangulations ...
AbstractTriangulation * getData()
bool isVertexOnBoundary(const SimplexId &vertexId) const override
Triangulation::Type getType() const
int preconditionBoundaryVertices() override
The Topology ToolKit.
COMMON_EXPORTS int MPIsize_
Definition BaseClass.cpp:10
COMMON_EXPORTS int MPIrank_
Definition BaseClass.cpp:9
const char VertexScalarFieldName[]
default name for vertex scalar field
Definition DataTypes.h:35
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
vtkStandardNewMacro(ttkScalarFieldCriticalPoints)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)