OGS
VtkColorLookupTable Class Reference

Detailed Description

Calculates and stores a colour lookup table.

Based on a start colour and an end colour, RGB-values are interpolated and stored in vector of GeoLib::Color. If no colours are set, default values are used for start (blue) and end (red). The number of entries of the colour table can be set in the constructor, the default value is 256. If additional colours are inserted into the table using setColor() the interpolation will be calculated iteratively between set colour values. Interpolation can be linear (default) or exponential. Based on the set range of values, colour values can be retrieved using getColor().

Definition at line 34 of file VtkColorLookupTable.h.

#include <VtkColorLookupTable.h>

Inheritance diagram for VtkColorLookupTable:
[legend]
Collaboration diagram for VtkColorLookupTable:
[legend]

Public Member Functions

 vtkTypeMacro (VtkColorLookupTable, vtkLookupTable)
 
void Build () override
 Builds the colour table based on the previously set parameters. This method should only be called after all options have been set. More...
 
void setColor (double pos, DataHolderLib::Color const &color)
 
void getColor (vtkIdType indx, unsigned char rgba[4]) const
 
DataHolderLib::LUTType getInterpolationType () const
 Returns the type of interpolation used. More...
 
void setInterpolationType (DataHolderLib::LUTType type)
 Sets the type of interpolation. More...
 
void setLookupTable (DataHolderLib::ColorLookupTable const &lut)
 Imports settings of OGS lookup table class. More...
 
void writeToFile (const std::string &filename)
 Exports a color table to a file. More...
 
void SetTableValueRGBA (vtkIdType idx, unsigned char rgba[4])
 Set a value within the LUT. More...
 
void GetTableValue (vtkIdType idx, unsigned char rgba[4])
 Get a value from the LUT. More...
 

Static Public Member Functions

static VtkColorLookupTableNew ()
 Create new objects with New() because of VTKs object reference counting. More...
 

Static Public Attributes

static const int DEFAULTMINVALUE = -9999
 
static const int DEFAULTMAXVALUE = 9999
 

Protected Member Functions

 VtkColorLookupTable ()
 Constructor. More...
 
 ~VtkColorLookupTable () override
 Destructor. More...
 

Private Member Functions

unsigned char linInterpolation (unsigned char a, unsigned char b, double p) const
 Interpolates values linearly. More...
 
unsigned char expInterpolation (unsigned char a, unsigned char b, double gamma, double p) const
 Interpolates values exponentially. gamma should roughly be in [0,4), for gamma=1 interpolation is linear. More...
 

Private Attributes

std::map< double, unsigned char * > _dict
 
DataHolderLib::LUTType _type {DataHolderLib::LUTType::LINEAR}
 

Constructor & Destructor Documentation

◆ VtkColorLookupTable()

VtkColorLookupTable::VtkColorLookupTable ( )
protecteddefault

Constructor.

◆ ~VtkColorLookupTable()

VtkColorLookupTable::~VtkColorLookupTable ( )
overrideprotected

Destructor.

Definition at line 29 of file VtkColorLookupTable.cpp.

30 {
31  for (auto& it : _dict)
32  {
33  delete it.second;
34  }
35 }
std::map< double, unsigned char * > _dict

References _dict.

Member Function Documentation

◆ Build()

void VtkColorLookupTable::Build ( )
override

Builds the colour table based on the previously set parameters. This method should only be called after all options have been set.

Definition at line 53 of file VtkColorLookupTable.cpp.

54 {
55  double range[2];
56  this->GetTableRange(range);
57  const double interval = range[1] - range[0];
58  this->SetNumberOfTableValues(static_cast<vtkIdType>(ceil(interval)) + 1);
59  // const vtkIdType nColours = this->GetNumberOfTableValues();
60  if (!_dict.empty())
61  {
62  // make sure that color map starts with the first color in the
63  // dictionary
64  unsigned char startcolor[4] = {0, 0, 0, 0};
65  std::pair<std::size_t, unsigned char*> lastValue(0, startcolor);
66 
67  for (auto& it : _dict)
68  {
69  double val = (it.first < range[0])
70  ? range[0]
71  : ((it.first > range[1]) ? range[1] : it.first);
72  auto const nextIndex =
73  static_cast<std::size_t>(std::floor(val - range[0]));
74 
75  this->SetTableValueRGBA(nextIndex, it.second);
76 
77  if (nextIndex - lastValue.first > 0)
78  {
79  for (std::size_t i = lastValue.first + 1; i < nextIndex; i++)
80  {
81  unsigned char int_rgba[4];
82  double pos =
83  (i - lastValue.first) /
84  (static_cast<double>(nextIndex - lastValue.first));
85 
87  {
88  for (std::size_t j = 0; j < 4; j++)
89  {
90  int_rgba[j] = linInterpolation(
91  (lastValue.second)[j], (it.second)[j], pos);
92  }
93  }
95  {
96  for (std::size_t j = 0; j < 4; j++)
97  {
98  int_rgba[j] =
99  expInterpolation((lastValue.second)[j],
100  (it.second)[j], 0.2, pos);
101  }
102  }
103  else
104  { // no interpolation
105  for (std::size_t j = 0; j < 4; j++)
106  {
107  int_rgba[j] = (lastValue.second)[j];
108  }
109  }
110 
111  this->SetTableValueRGBA(i, int_rgba);
112  }
113  }
114 
115  lastValue.first = nextIndex;
116  lastValue.second = it.second;
117  }
118  }
119  else
120  {
121  vtkLookupTable::Build();
122  }
123 }
DataHolderLib::LUTType _type
unsigned char expInterpolation(unsigned char a, unsigned char b, double gamma, double p) const
Interpolates values exponentially. gamma should roughly be in [0,4), for gamma=1 interpolation is lin...
void SetTableValueRGBA(vtkIdType idx, unsigned char rgba[4])
Set a value within the LUT.
unsigned char linInterpolation(unsigned char a, unsigned char b, double p) const
Interpolates values linearly.

References _dict, _type, expInterpolation(), DataHolderLib::EXPONENTIAL, DataHolderLib::LINEAR, linInterpolation(), and SetTableValueRGBA().

Referenced by VtkCompositeElementSelectionFilter::GetLookupTable(), VtkCompositeColorByHeightFilter::init(), and setLookupTable().

◆ expInterpolation()

unsigned char VtkColorLookupTable::expInterpolation ( unsigned char  a,
unsigned char  b,
double  gamma,
double  p 
) const
private

Interpolates values exponentially. gamma should roughly be in [0,4), for gamma=1 interpolation is linear.

Definition at line 44 of file VtkColorLookupTable.cpp.

48 {
49  assert(gamma > 0 && gamma < 4);
50  return static_cast<unsigned char>((b - a) * pow(p, gamma) + a);
51 }
static const double p

References MathLib::p.

Referenced by Build().

◆ getColor()

void VtkColorLookupTable::getColor ( vtkIdType  indx,
unsigned char  rgba[4] 
) const

Definition at line 193 of file VtkColorLookupTable.cpp.

194 {
195  indx = ((indx < this->TableRange[0])
196  ? static_cast<vtkIdType>(this->TableRange[0])
197  : (indx >= this->TableRange[1]
198  ? static_cast<vtkIdType>(this->TableRange[1]) - 1
199  : indx));
200  indx = static_cast<std::size_t>(std::floor(
201  (indx - this->TableRange[0]) *
202  (this->NumberOfColors / (this->TableRange[1] - this->TableRange[0]))));
203 
204  unsigned char* _rgba;
205  _rgba = this->Table->GetPointer(indx * 4);
206  for (std::size_t i = 0; i < 4; i++)
207  {
208  rgba[i] = _rgba[i];
209  }
210 }

◆ getInterpolationType()

DataHolderLib::LUTType VtkColorLookupTable::getInterpolationType ( ) const
inline

Returns the type of interpolation used.

Definition at line 63 of file VtkColorLookupTable.h.

63 { return _type; }

References _type.

Referenced by VtkColorByHeightFilter::PrintSelf().

◆ GetTableValue()

void VtkColorLookupTable::GetTableValue ( vtkIdType  idx,
unsigned char  rgba[4] 
)

Get a value from the LUT.

Definition at line 171 of file VtkColorLookupTable.cpp.

172 {
173  double value[4];
174  vtkLookupTable::GetTableValue(idx, value);
175 
176  for (unsigned i = 0; i < 4; ++i)
177  {
178  rgba[i] = static_cast<unsigned char>(value[i] * 255.0);
179  }
180 }

Referenced by writeToFile().

◆ linInterpolation()

unsigned char VtkColorLookupTable::linInterpolation ( unsigned char  a,
unsigned char  b,
double  p 
) const
private

Interpolates values linearly.

Definition at line 37 of file VtkColorLookupTable.cpp.

40 {
41  return static_cast<unsigned char>(a * (1 - p) + b * p);
42 }

References MathLib::p.

Referenced by Build().

◆ New()

static VtkColorLookupTable* VtkColorLookupTable::New ( )
static

Create new objects with New() because of VTKs object reference counting.

Referenced by VtkCompositeElementSelectionFilter::GetLookupTable(), and VtkAlgorithmProperties::SetLookUpTable().

◆ setColor()

void VtkColorLookupTable::setColor ( double  pos,
DataHolderLib::Color const &  color 
)

Definition at line 182 of file VtkColorLookupTable.cpp.

184 {
185  auto* dict_rgba = new unsigned char[4];
186  for (std::size_t i = 0; i < 4; i++)
187  {
188  dict_rgba[i] = color[i];
189  }
190  _dict.insert(std::pair<double, unsigned char*>(pos, dict_rgba));
191 }

References _dict.

Referenced by VtkCompositeElementSelectionFilter::GetLookupTable(), VtkCompositeColorByHeightFilter::init(), and setLookupTable().

◆ setInterpolationType()

void VtkColorLookupTable::setInterpolationType ( DataHolderLib::LUTType  type)
inline

Sets the type of interpolation.

Definition at line 66 of file VtkColorLookupTable.h.

66 { _type = type; }

References _type.

Referenced by VtkCompositeColorByHeightFilter::init(), and setLookupTable().

◆ setLookupTable()

void VtkColorLookupTable::setLookupTable ( DataHolderLib::ColorLookupTable const &  lut)

Imports settings of OGS lookup table class.

Definition at line 125 of file VtkColorLookupTable.cpp.

127 {
128  std::size_t const n_colors(lut.size());
129  for (std::size_t i = 0; i < n_colors; ++i)
130  {
131  setColor(std::get<0>(lut[i]), std::get<1>(lut[i]));
132  }
133  setInterpolationType(lut.getInterpolationType());
134  auto const range(lut.getTableRange());
135  SetTableRange(range.first, range.second);
136  Build();
137 }
void setColor(double pos, DataHolderLib::Color const &color)
void Build() override
Builds the colour table based on the previously set parameters. This method should only be called aft...
void setInterpolationType(DataHolderLib::LUTType type)
Sets the type of interpolation.

References Build(), DataHolderLib::ColorLookupTable::getInterpolationType(), DataHolderLib::ColorLookupTable::getTableRange(), setColor(), setInterpolationType(), and DataHolderLib::ColorLookupTable::size().

Referenced by VtkAlgorithmProperties::SetLookUpTable().

◆ SetTableValueRGBA()

void VtkColorLookupTable::SetTableValueRGBA ( vtkIdType  idx,
unsigned char  rgba[4] 
)

Set a value within the LUT.

Definition at line 159 of file VtkColorLookupTable.cpp.

161 {
162  double value[4];
163 
164  for (unsigned i = 0; i < 4; ++i)
165  {
166  value[i] = rgba[i] / 255.0;
167  }
168  vtkLookupTable::SetTableValue(idx, value);
169 }

Referenced by Build().

◆ vtkTypeMacro()

VtkColorLookupTable::vtkTypeMacro ( VtkColorLookupTable  ,
vtkLookupTable   
)

◆ writeToFile()

void VtkColorLookupTable::writeToFile ( const std::string &  filename)

Exports a color table to a file.

Definition at line 139 of file VtkColorLookupTable.cpp.

140 {
141  std::stringstream strout;
142  strout << "Writing color table to " << filename << " ... ";
143  std::ofstream out(filename.c_str(), std::ios::out);
144 
145  std::size_t nColors = this->GetNumberOfTableValues();
146  for (std::size_t i = 0; i < nColors; i++)
147  {
148  unsigned char rgba[4];
149  this->GetTableValue(i, rgba);
150  out << i << "\t" << rgba[0] << "\t" << rgba[1] << "\t" << rgba[2]
151  << "\n";
152  }
153 
154  strout << " done." << std::endl;
155  INFO("{:s}", strout.str());
156  out.close();
157 }
void INFO(char const *fmt, Args const &... args)
Definition: Logging.h:32
void GetTableValue(vtkIdType idx, unsigned char rgba[4])
Get a value from the LUT.

References GetTableValue(), and INFO().

Member Data Documentation

◆ _dict

std::map<double, unsigned char*> VtkColorLookupTable::_dict
private

Definition at line 94 of file VtkColorLookupTable.h.

Referenced by ~VtkColorLookupTable(), Build(), and setColor().

◆ _type

DataHolderLib::LUTType VtkColorLookupTable::_type {DataHolderLib::LUTType::LINEAR}
private

Definition at line 95 of file VtkColorLookupTable.h.

Referenced by Build(), getInterpolationType(), and setInterpolationType().

◆ DEFAULTMAXVALUE

const int VtkColorLookupTable::DEFAULTMAXVALUE = 9999
static

Definition at line 38 of file VtkColorLookupTable.h.

◆ DEFAULTMINVALUE

const int VtkColorLookupTable::DEFAULTMINVALUE = -9999
static

Definition at line 37 of file VtkColorLookupTable.h.


The documentation for this class was generated from the following files: