TTK
Loading...
Searching...
No Matches
ttkCinemaDarkroomColorMapping.cpp
Go to the documentation of this file.
1#include "BaseClass.h"
3
4#include <vtkImageData.h>
5#include <vtkInformation.h>
6#include <vtkObjectFactory.h>
7#include <vtkPointData.h>
8#include <vtkUnsignedCharArray.h>
9
10#include <ttkMacros.h>
11#include <ttkUtils.h>
12
13#include <cmath>
14
16
18 this->setDebugMsgPrefix("CinemaDarkroomColorMapping");
19
20 this->SetNumberOfInputPorts(1);
21 this->SetNumberOfOutputPorts(1);
22}
23
25
26template <typename DT>
27int mapScalarsToColor(unsigned char *color,
28 const std::vector<double> &colorMap,
29 const double *nanColor,
30 const DT *array,
31 const double range[2],
32 const size_t nPixels,
33 const int threadNumber) {
34 const size_t nKeys = colorMap.size() / 4;
35 const double valueDelta = range[1] - range[0];
36
37#ifdef TTK_ENABLE_OPENMP
38#pragma omp parallel for num_threads(threadNumber)
39#endif
40 for(size_t i = 0; i < nPixels; i++) {
41
42 const double value = (double)array[i];
43 if(std::isnan(value)) {
44 size_t const idx = i * 3;
45 color[idx + 0] = 255.0 * nanColor[0];
46 color[idx + 1] = 255.0 * nanColor[1];
47 color[idx + 2] = 255.0 * nanColor[2];
48 continue;
49 }
50
51 const double normalizedValue
52 = std::max(0.0, std::min(1.0, (value - range[0]) / valueDelta));
53
54 size_t ki = 0;
55 for(size_t k = 1; k < nKeys; k++) {
56 if(normalizedValue <= colorMap[k * 4]) {
57 ki = k - 1;
58 break;
59 }
60 }
61
62 double const lambda = (normalizedValue - colorMap[ki * 4])
63 / (colorMap[(ki + 1) * 4] - colorMap[ki * 4]);
64 double const lambdaInv = 1 - lambda;
65
66 size_t const idx = i * 3;
67 size_t const idx2 = ki * 4;
68 color[idx + 0]
69 = 255.0 * (lambdaInv * colorMap[idx2 + 1] + lambda * colorMap[idx2 + 5]);
70 color[idx + 1]
71 = 255.0 * (lambdaInv * colorMap[idx2 + 2] + lambda * colorMap[idx2 + 6]);
72 color[idx + 2]
73 = 255.0 * (lambdaInv * colorMap[idx2 + 3] + lambda * colorMap[idx2 + 7]);
74 }
75
76 TTK_FORCE_USE(threadNumber);
77 return 1;
78}
79
81 vtkInformation *ttkNotUsed(request),
82 vtkInformationVector **inputVector,
83 vtkInformationVector *outputVector) {
84
85 ttk::Timer timer;
86 this->printMsg("Applying Color Map", 0, 0, this->threadNumber_,
88
89 auto input = vtkImageData::GetData(inputVector[0]);
90 auto output = vtkImageData::GetData(outputVector);
91 output->ShallowCopy(input);
92
93 auto scalarArray = this->GetInputArrayToProcess(0, output);
94 if(!scalarArray || this->GetInputArrayAssociation(0, output) != 0
95 || scalarArray->GetNumberOfComponents() != 1) {
96 this->printErr("Unable to retrieve point scalar array.");
97 return 0;
98 }
99
100 size_t const nPixels = scalarArray->GetNumberOfTuples();
101
102 double range[2];
103 scalarArray->GetRange(range);
104
106 colorArray->SetName("Diffuse");
107 colorArray->SetNumberOfComponents(3);
108 colorArray->SetNumberOfTuples(nPixels);
109 output->GetPointData()->AddArray(colorArray);
110 auto colorArrayData
111 = static_cast<unsigned char *>(ttkUtils::GetVoidPointer(colorArray));
112
113 std::vector<double> manualColorMap;
114 const std::vector<double> *colorMap{nullptr};
115
116 if(this->ColorMap == -1) { // Solid
117 manualColorMap.resize(8);
118 manualColorMap[0] = 0;
119 manualColorMap[1] = this->SingleColor[0];
120 manualColorMap[2] = this->SingleColor[1];
121 manualColorMap[3] = this->SingleColor[2];
122 manualColorMap[4] = 1;
123 manualColorMap[5] = this->SingleColor[0];
124 manualColorMap[6] = this->SingleColor[1];
125 manualColorMap[7] = this->SingleColor[2];
126 colorMap = &manualColorMap;
127 } else if(this->ColorMap == -2) {
128 int const status = ttkUtils::stringListToDoubleVector(
129 this->ManualColorMap, manualColorMap);
130 if(!status || manualColorMap.size() < 8 || manualColorMap.size() % 4 != 0) {
131 this->printErr("Invalid manual color map input.");
132 return 0;
133 }
134 colorMap = &manualColorMap;
135 } else {
136 if(this->ColorMap < 0 || this->ColorMap >= (int)this->ColorMaps.size()) {
137 this->printErr("Invalid color map index: "
138 + std::to_string(this->ColorMap));
139 return 0;
140 }
141 colorMap = &this->ColorMaps[this->ColorMap];
142 }
143
144 switch(scalarArray->GetDataType()) {
145 vtkTemplateMacro(mapScalarsToColor<VTK_TT>(
146 colorArrayData, *colorMap, this->NANColor,
147 static_cast<VTK_TT *>(ttkUtils::GetVoidPointer(scalarArray)), range,
148 nPixels, this->threadNumber_));
149 }
150
151 this->printMsg(
152 "Applying Color Map", 1, timer.getElapsedTime(), this->threadNumber_);
153
154 return 1;
155}
#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
Performs color mapping of a scalar field.
~ttkCinemaDarkroomColorMapping() override
int RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
static void * GetVoidPointer(vtkDataArray *array, vtkIdType start=0)
Definition ttkUtils.cpp:226
static int stringListToDoubleVector(const std::string &iString, std::vector< double > &v)
Definition ttkUtils.cpp:133
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
int mapScalarsToColor(unsigned char *color, const std::vector< double > &colorMap, const double *nanColor, const DT *array, const double range[2], const size_t nPixels, const int threadNumber)
vtkStandardNewMacro(ttkCinemaDarkroomColorMapping)
printMsg(debug::output::BOLD+" | | | | | . \\ | | (__| | / __/| |_| / __/|__ _|"+debug::output::ENDCOLOR, debug::Priority::PERFORMANCE, debug::LineMode::NEW, stream)