TTK
Loading...
Searching...
No Matches
ttkIdentifierRandomizer.cpp
Go to the documentation of this file.
1#include <Shuffle.h>
3
4#include <vtkCellData.h>
5#include <vtkDataArray.h>
6#include <vtkDataSet.h>
7#include <vtkInformation.h>
8#include <vtkObjectFactory.h>
9#include <vtkPointData.h>
10
11#include <ttkMacros.h>
12#include <ttkUtils.h>
13
14#include <map>
15#include <numeric>
16#include <random>
17
19
21
22 this->SetNumberOfInputPorts(1);
23 this->SetNumberOfOutputPorts(1);
24
25 this->setDebugMsgPrefix("IdentifierRandomizer");
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
46template <typename T>
47int shuffleScalarFieldValues(const T *const inputField,
48 T *const outputField,
49 const int nValues,
50 const int seed,
51 const bool compactRange,
52 const int nThreads = 1) {
53
54 // copy input field into vector
55 std::vector<T> inputValues(inputField, inputField + nValues);
56
57 // reduce the copy
58 TTK_PSORT(nThreads, inputValues.begin(), inputValues.end());
59 const auto last = std::unique(inputValues.begin(), inputValues.end());
60 inputValues.erase(last, inputValues.end());
61
62 // copy the range of values
63 std::vector<T> shuffledValues(inputValues.size());
64 if(compactRange) {
65 std::iota(shuffledValues.begin(), shuffledValues.end(), T{});
66 } else {
67 std::copy(inputValues.begin(), inputValues.end(), shuffledValues.begin());
68 }
69
70 // shuffle them using the seed
71 std::mt19937 random_engine{};
72 random_engine.seed(seed);
73 // use the Fisher-Yates algorithm instead of std::shuffle, whose
74 // results are platform-dependent
75 ttk::shuffle(shuffledValues, random_engine);
76
77 // link original value to shuffled value correspondence
78 std::map<T, T> originalToShuffledValues{};
79 for(size_t i = 0; i < inputValues.size(); ++i) {
80 originalToShuffledValues[inputValues[i]] = shuffledValues[i];
81 }
82
83// write shuffled values inside the output scalar field
84#ifdef TTK_ENABLE_OPENMP
85#pragma omp parallel for num_threads(nThreads)
86#endif // TTK_ENABLE_OPENMP
87 for(int i = 0; i < nValues; ++i) {
88 outputField[i] = originalToShuffledValues[inputField[i]];
89 }
90
91 TTK_FORCE_USE(nThreads);
92 return 1;
93}
94
96 vtkInformationVector **inputVector,
97 vtkInformationVector *outputVector) {
98
99 ttk::Timer t;
100
101 bool isPointData = false;
102
103 vtkDataSet *input = vtkDataSet::GetData(inputVector[0]);
104 vtkDataSet *output = vtkDataSet::GetData(outputVector);
105
106 // use a pointer-base copy for the input data -- to adapt if your wrapper does
107 // not produce an output of the type of the input.
108 output->ShallowCopy(input);
109
110 // in the following, the target scalar field of the input is replaced in the
111 // variable 'output' with the result of the computation.
112 // if your wrapper produces an output of the same type of the input, you
113 // should proceed in the same way.
114 vtkDataArray *inputScalarField = this->GetInputArrayToProcess(0, inputVector);
115
116 if(!inputScalarField) {
117 printErr("Could not retrieve mandatory input array :(");
118 return 0;
119 }
120
121 if(input->GetPointData()->GetArray(inputScalarField->GetName())
122 == inputScalarField) {
123 isPointData = true;
124 }
125
126 this->printMsg("Shuffling " + std::string{isPointData ? "vertex" : "cell"}
127 + " field `" + std::string{inputScalarField->GetName()}
128 + "'...");
129
130 // allocate the memory for the output scalar field
131 vtkSmartPointer<vtkDataArray> const outputArray
132 = vtkSmartPointer<vtkDataArray>::Take(inputScalarField->NewInstance());
133 outputArray->SetName(inputScalarField->GetName());
134 outputArray->SetNumberOfComponents(1);
135 outputArray->SetNumberOfTuples(inputScalarField->GetNumberOfTuples());
136
137 switch(outputArray->GetDataType()) {
138 vtkTemplateMacro(shuffleScalarFieldValues(
139 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(inputScalarField)),
140 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(outputArray)),
141 outputArray->GetNumberOfTuples(), this->RandomSeed, this->CompactRange,
142 this->threadNumber_));
143 }
144
145 if(isPointData)
146 output->GetPointData()->AddArray(outputArray);
147 else
148 output->GetCellData()->AddArray(outputArray);
149
150 printMsg("Processed " + std::to_string(outputArray->GetNumberOfTuples())
151 + (isPointData ? " vertices." : " cells."),
152 1, t.getElapsedTime(), 1);
153
155
156 return 1;
157}
#define TTK_FORCE_USE(x)
Force the compiler to use the function/method parameter.
Definition BaseClass.h:57
#define ttkNotUsed(x)
Mark function/method parameters that are not used in the function body at all.
Definition BaseClass.h:47
#define TTK_PSORT(NTHREADS,...)
Parallel sort macro.
Definition OpenMP.h:46
static vtkInformationIntegerKey * SAME_DATA_TYPE_AS_INPUT_PORT()
TTK VTK-filter that randomly shuffles segmentation identifiers.
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
void setDebugMsgPrefix(const std::string &prefix)
Definition Debug.h:364
int printErr(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:149
double getElapsedTime()
Definition Timer.h:15
void shuffle(std::vector< T > &toShuffle, U &&rng)
Platform-independent alternative to std::shuffle implementing the Fisher-Yates shuffle algorithm.
Definition Shuffle.h:13
vtkStandardNewMacro(ttkIdentifierRandomizer)
int shuffleScalarFieldValues(const T *const inputField, T *const outputField, const int nValues, const int seed, const bool compactRange, const int nThreads=1)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)