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(), "vtkDataObject");
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
95template <typename T>
97 vtkMultiBlockDataSet *input,
98 vtkMultiBlockDataSet *output,
99 const int nThreads) {
100
101 int n_blocks = input->GetNumberOfBlocks();
102 std::vector<T> inputValues;
103 for(int i = 0; i < n_blocks; i++) {
104 vtkDataSet *block = vtkDataSet::SafeDownCast(input->GetBlock(i));
105 vtkDataArray *inputScalarField = this->GetInputArrayToProcess(0, block);
106 int nValues = inputScalarField->GetNumberOfTuples();
107 const T *const inputScalarFieldPtr
108 = static_cast<T *>(ttkUtils::GetVoidPointer(inputScalarField));
109 inputValues.insert(
110 inputValues.end(), inputScalarFieldPtr, inputScalarFieldPtr + nValues);
111 }
112
113 TTK_PSORT(nThreads, inputValues.begin(), inputValues.end());
114 const auto last = std::unique(inputValues.begin(), inputValues.end());
115 inputValues.erase(last, inputValues.end());
116 std::vector<T> shuffledValues(inputValues.size());
117 if(CompactRange) {
118 std::iota(shuffledValues.begin(), shuffledValues.end(), T{});
119 } else {
120 std::copy(inputValues.begin(), inputValues.end(), shuffledValues.begin());
121 }
122 // shuffle them using the seed
123 std::mt19937 random_engine{};
124 random_engine.seed(RandomSeed);
125 // use the Fisher-Yates algorithm instead of std::shuffle, whose
126 // results are platform-dependent
127 ttk::shuffle(shuffledValues, random_engine);
128
129 // link original value to shuffled value correspondence
130 std::map<T, T> originalToShuffledValues{};
131 for(size_t i = 0; i < inputValues.size(); ++i) {
132 originalToShuffledValues[inputValues[i]] = shuffledValues[i];
133 }
134
135// write shuffled values inside the output scalar field
136#ifdef TTK_ENABLE_OPENMP
137#pragma omp parallel for num_threads(nThreads)
138#endif // TTK_ENABLE_OPENMP
139 for(int i = 0; i < n_blocks; ++i) {
140 vtkDataSet *block = vtkDataSet::SafeDownCast(input->GetBlock(i));
141 vtkDataArray *inputScalarField = this->GetInputArrayToProcess(0, block);
142 int nValues = inputScalarField->GetNumberOfTuples();
143
144 vtkSmartPointer<vtkDataArray> const outputArray
145 = vtkSmartPointer<vtkDataArray>::Take(inputScalarField->NewInstance());
146 outputArray->SetName(inputScalarField->GetName());
147 outputArray->SetNumberOfComponents(1);
148 outputArray->SetNumberOfTuples(inputScalarField->GetNumberOfTuples());
149 T *const inputArrayPtr
150 = static_cast<T *>(ttkUtils::GetVoidPointer(inputScalarField));
151 T *const outputArrayPtr
152 = static_cast<T *>(ttkUtils::GetVoidPointer(outputArray));
153 for(int j = 0; j < nValues; ++j) {
154 T newValue = originalToShuffledValues[inputArrayPtr[j]];
155 outputArrayPtr[j] = newValue;
156 }
157
158 vtkDataSet *outputBlock = vtkDataSet::SafeDownCast(output->GetBlock(i));
159 if(block->GetPointData()->GetArray(inputScalarField->GetName())) {
160 outputBlock->GetPointData()->AddArray(outputArray);
161 } else {
162 outputBlock->GetCellData()->AddArray(outputArray);
163 }
164 }
165 TTK_FORCE_USE(nThreads);
166 return 1;
167}
168
170 vtkInformationVector **inputVector,
171 vtkInformationVector *outputVector) {
172
173 ttk::Timer t;
174
175 vtkDataSet *input
176 = vtkDataSet::SafeDownCast(vtkDataSet::GetData(inputVector[0]));
177 if(input) {
178 vtkDataSet *output = vtkDataSet::GetData(outputVector);
179 output->ShallowCopy(input);
180 // use a pointer-base copy for the input data -- to adapt if your wrapper
181 // does not produce an output of the type of the input.
182
183 // in the following, the target scalar field of the input is replaced in the
184 // variable 'output' with the result of the computation.
185 // if your wrapper produces an output of the same type of the input, you
186 // should proceed in the same way.
187 vtkDataArray *inputScalarField
188 = this->GetInputArrayToProcess(0, inputVector);
189
190 if(!inputScalarField) {
191 printErr("Could not retrieve mandatory input array :(");
192 return 0;
193 }
194
195 bool isPointData = false;
196 if(input->GetPointData()->GetArray(inputScalarField->GetName())
197 == inputScalarField) {
198 isPointData = true;
199 }
200
201 this->printMsg("Shuffling " + std::string{isPointData ? "vertex" : "cell"}
202 + " field `" + std::string{inputScalarField->GetName()}
203 + "'...");
204
205 // allocate the memory for the output scalar field
206 vtkSmartPointer<vtkDataArray> const outputArray
207 = vtkSmartPointer<vtkDataArray>::Take(inputScalarField->NewInstance());
208 outputArray->SetName(inputScalarField->GetName());
209 outputArray->SetNumberOfComponents(1);
210 outputArray->SetNumberOfTuples(inputScalarField->GetNumberOfTuples());
211
212 switch(outputArray->GetDataType()) {
213 vtkTemplateMacro(shuffleScalarFieldValues(
214 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(inputScalarField)),
215 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(outputArray)),
216 outputArray->GetNumberOfTuples(), this->RandomSeed, this->CompactRange,
217 this->threadNumber_));
218 }
219
220 if(isPointData)
221 output->GetPointData()->AddArray(outputArray);
222 else
223 output->GetCellData()->AddArray(outputArray);
224
225 printMsg("Processed " + std::to_string(outputArray->GetNumberOfTuples())
226 + (isPointData ? " vertices." : " cells."),
227 1, t.getElapsedTime(), 1);
228
230
231 } else {
232 vtkMultiBlockDataSet *input_mb = vtkMultiBlockDataSet::SafeDownCast(
233 vtkMultiBlockDataSet::GetData(inputVector[0]));
234 vtkMultiBlockDataSet *output_mb
235 = vtkMultiBlockDataSet::GetData(outputVector);
236 output_mb->ShallowCopy(input_mb);
237 if(!input_mb) {
238 printMsg("Invalid input.");
239 }
240 int n_blocks = input_mb->GetNumberOfBlocks();
241 int currentType = 0;
242 // check if the multiblock input and the input scalar field are valid and
243 // retrieve the data type
244 for(int i = 0; i < n_blocks; i++) {
245 vtkDataSet *block = vtkDataSet::SafeDownCast(input_mb->GetBlock(i));
246 if(!block)
247 printMsg("Block " + std::to_string(i) + " invalid.");
248 vtkDataArray *inputScalarField = this->GetInputArrayToProcess(0, block);
249
250 if(!inputScalarField) {
251 printWrn(
252 "Block " + std::to_string(i)
253 + " does not have the required input scalar field as data array.");
254 continue;
255 }
256
257 if(i == 0) {
258 currentType = inputScalarField->GetDataType();
259 continue;
260 } else {
261 if(currentType != inputScalarField->GetDataType()) {
262 printErr("All block's input scalar field must have the same type.");
263 return 0;
264 }
265 }
266 }
267 switch(currentType) {
269 input_mb, output_mb, this->threadNumber_));
270 }
271 }
272
273 return 1;
274}
#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 shuffleScalarFieldValuesMultiBlock(vtkMultiBlockDataSet *input, vtkMultiBlockDataSet *output, const int nThreads=1)
int FillOutputPortInformation(int port, vtkInformation *info) override
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
int printWrn(const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
Definition Debug.h:159
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)