TTK
Loading...
Searching...
No Matches
ttkArrayPreconditioning.cpp
Go to the documentation of this file.
2#include <Triangulation.h>
4#include <ttkMacros.h>
5#include <ttkUtils.h>
6
7#include <regex>
8#include <vtkCommand.h>
9#include <vtkDataArray.h>
10#include <vtkDataSet.h>
11#include <vtkIdTypeArray.h>
12#include <vtkInformation.h>
13#include <vtkIntArray.h>
14#include <vtkPointData.h>
15
17
19 this->SetNumberOfInputPorts(1);
20 this->SetNumberOfOutputPorts(1);
21
22 // ensure that modifying the selection re-triggers the filter
23 // (c.f. vtkPassSelectedArrays.cxx)
24 this->ArraySelection->AddObserver(
25 vtkCommand::ModifiedEvent, this, &ttkArrayPreconditioning::Modified);
26}
27
29 vtkInformation *info) {
30 if(port == 0) {
31 info->Set(vtkAlgorithm::INPUT_REQUIRED_DATA_TYPE(), "vtkDataSet");
32 return 1;
33 }
34 return 0;
35}
36
38 vtkInformation *info) {
39 if(port == 0) {
41 return 1;
42 }
43 return 0;
44}
45
47 vtkInformationVector **inputVector,
48 vtkInformationVector *outputVector) {
49
50 auto input = vtkDataSet::GetData(inputVector[0]);
51 auto output = vtkDataSet::GetData(outputVector);
52 ttk::Timer tm{};
53
54 int const keepGoing = checkEmptyMPIInput<vtkDataSet>(input);
55 if(keepGoing < 2) {
56 return keepGoing;
57 }
58#ifdef TTK_ENABLE_MPI
59 const auto triangulation{this->GetTriangulation(input)};
60#endif
61 output->ShallowCopy(input);
62
63 auto pointData = input->GetPointData();
64 size_t const nVertices = input->GetNumberOfPoints();
65
66 std::vector<vtkDataArray *> scalarArrays{};
67
68 if(SelectFieldsWithRegexp) {
69 // select all input point data arrays whose name is matching the regexp
70 const auto n = pointData->GetNumberOfArrays();
71 for(int i = 0; i < n; ++i) {
72 auto array = pointData->GetArray(i);
73 if(array != nullptr && array->GetName() != nullptr
74 && std::regex_match(array->GetName(), std::regex(RegexpString))) {
75 scalarArrays.emplace_back(array);
76 }
77 }
78 } else {
79 // get all selected input point data arrays
80 for(int i = 0; i < pointData->GetNumberOfArrays(); ++i) {
81 auto array = pointData->GetArray(i);
82 if(array != nullptr && array->GetName() != nullptr
83 && ArraySelection->ArrayIsEnabled(array->GetName())) {
84 scalarArrays.emplace_back(array);
85 }
86 }
87 }
88
89#ifdef TTK_ENABLE_MPI
90 if(GlobalOrder) {
91 if(ttk::isRunningWithMPI()) {
92
93 this->preconditionTriangulation(triangulation);
94#ifdef TTK_ENABLE_MPI_TIME
95 ttk::Timer t_mpi;
96 ttk::startMPITimer(t_mpi, ttk::MPIrank_, ttk::MPIsize_);
97#endif
98 // add the order array for every scalar array, except the ghostcells, the
99 // rankarray and the global ids
100 for(auto scalarArray : scalarArrays) {
101 int status = 0;
102 std::string arrayName = std::string(scalarArray->GetName());
103 if(arrayName != "GlobalPointIds" && arrayName != "vtkGhostType"
104 && arrayName != "RankArray") {
105 this->printMsg("Arrayname: " + arrayName);
106 vtkNew<ttkSimplexIdTypeArray> orderArray{};
107 orderArray->SetName(
109 orderArray->SetNumberOfComponents(1);
110 orderArray->SetNumberOfTuples(nVertices);
111
112 this->printMsg(std::to_string(scalarArray->GetDataType()));
113 ttkTypeMacroAT(scalarArray->GetDataType(), triangulation->getType(),
115 static_cast<const T1 *>(triangulation->getData()),
117 ttkUtils::GetPointer<T0>(scalarArray), nVertices)));
118
119 // On error cancel filter execution
120 if(status != 1)
121 return 0;
122 output->GetPointData()->AddArray(orderArray);
123 triangulation->setIsOrderArrayGlobal(
124 ttkUtils::GetVoidPointer(scalarArray), true);
125 }
126 }
127#ifdef TTK_ENABLE_MPI_TIME
128 double elapsedTime
129 = ttk::endMPITimer(t_mpi, ttk::MPIrank_, ttk::MPIsize_);
130 if(ttk::MPIrank_ == 0) {
131 printMsg("Array preconditioning performed using "
132 + std::to_string(ttk::MPIsize_)
133 + " MPI processes lasted: " + std::to_string(elapsedTime));
134 }
135#endif
136 this->printMsg("Preconditioned selected scalar arrays", 1.0,
137 tm.getElapsedTime(), this->threadNumber_);
138 return 1;
139 } else {
140 this->printMsg("Necessary arrays are present, TTK is built with MPI "
141 "support, but not run with mpirun. Running sequentially.");
142 }
143 }
144#endif
145#ifdef TTK_ENABLE_MPI_TIME
146 ttk::Timer t_mpi;
147 ttk::startMPITimer(t_mpi, ttk::MPIrank_, ttk::MPIsize_);
148#endif
149 for(auto scalarArray : scalarArrays) {
150 vtkNew<ttkSimplexIdTypeArray> orderArray{};
151 orderArray->SetName(
153 orderArray->SetNumberOfComponents(1);
154 orderArray->SetNumberOfTuples(nVertices);
155
156 switch(scalarArray->GetDataType()) {
157 vtkTemplateMacro(ttk::sortVertices(
158 nVertices, static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(scalarArray)),
159 static_cast<int *>(nullptr),
160 static_cast<ttk::SimplexId *>(ttkUtils::GetVoidPointer(orderArray)),
161 this->threadNumber_));
162 }
163#ifdef TTK_ENABLE_MPI
164 triangulation->setIsOrderArrayGlobal(
165 ttkUtils::GetVoidPointer(scalarArray), false);
166#endif // TTK_ENABLE_MPI
167 output->GetPointData()->AddArray(orderArray);
168 this->printMsg("Generated order array for scalar array `"
169 + std::string{scalarArray->GetName()} + "'");
170 }
171#ifdef TTK_ENABLE_MPI_TIME
172 double elapsedTime = ttk::endMPITimer(t_mpi, ttk::MPIrank_, ttk::MPIsize_);
173 if(ttk::MPIrank_ == 0) {
174 printMsg("Array preconditioning performed using "
175 + std::to_string(ttk::MPIsize_)
176 + " MPI processes lasted: " + std::to_string(elapsedTime));
177 }
178#endif
179 this->printMsg("Preconditioned selected scalar arrays", 1.0,
180 tm.getElapsedTime(), this->threadNumber_);
181
182 return 1;
183}
#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()
int checkEmptyMPIInput(inputType *input)
This method tests whether the input is a nullptr. If the computation is being done on multiple proces...
ttk::Triangulation * GetTriangulation(vtkDataSet *dataSet)
static std::string GetOrderArrayName(vtkDataArray *const array)
TTK VTK-filter to generate order fields.
int FillInputPortInformation(int port, vtkInformation *info) override
int FillOutputPortInformation(int port, vtkInformation *info) override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
static DT * GetPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.h:59
int preconditionTriangulation(AbstractTriangulation *triangulation)
int processScalarArray(const triangulationType *triangulation, ttk::SimplexId *orderArray, const DT *scalarArray, const size_t nVerts) const
COMMON_EXPORTS int MPIsize_
Definition BaseClass.cpp:10
int SimplexId
Identifier type for simplices of any dimension.
Definition DataTypes.h:22
void sortVertices(const size_t nVerts, const scalarType *const scalars, const idType *const offsets, SimplexId *const order, const int nThreads)
Sort vertices according to scalars disambiguated by offsets.
COMMON_EXPORTS int MPIrank_
Definition BaseClass.cpp:9
vtkStandardNewMacro(ttkArrayPreconditioning)
#define ttkTypeMacroAT(group0, group1, call)
Definition ttkMacros.h:272
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/| (_) |"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)