Loading...
Searching...
No Matches
colourTable.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2019 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "colourTable.H"
29#include "colourTools.H"
30#include "ListOps.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34const Foam::Enum
35<
37>
39({
40 { interpolationType::RGB, "rgb" },
42 { interpolationType::DIVERGING, "diverging" },
43});
44
45
46// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47
49(
50 const List<Tuple2<scalar, vector>>& values,
51 const interpolationType interp
53:
54 table_(values),
55 interp_(interp)
56{}
57
58
60(
62 const interpolationType interp
64:
65 table_(std::move(values)),
66 interp_(interp)
67{}
68
69
71(
72 const dictionary& dict,
73 const interpolationType interp
74)
75:
76 table_(),
77 interp_(interp)
78{
79 dict.readEntry("table", table_);
80 interpolationTypeNames.readIfPresent("interpolate", dict, interp_);
81}
82
83
84// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
85
89}
90
91
92// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93
94Foam::vector Foam::colourTable::value(const scalar x) const
95{
96 if (x <= 0)
97 {
98 return table_.first().second();
99 }
100
101 if (x >= 1)
102 {
103 return table_.last().second();
104 }
105
106
107 label idx = findLower
108 (
109 table_, x, 0,
110 [](const pair_type& pr, const scalar& val)
111 {
112 // Test first element
113 return (pr.first() <= val);
114 }
115 );
116
117 if (idx == -1)
118 {
119 // Use first element only
120 return table_.first().second();
121 }
122 else if (idx == table_.size()-1)
123 {
124 // Use last element only
125 return table_.last().second();
126 }
127
128 const scalar t0 = table_[idx].first();
129 const scalar t1 = table_[idx+1].first();
130
131 const scalar s = (x - t0)/(t1 - t0);
132
133 const vector& rgb0 = table_[idx].second();
134 const vector& rgb1 = table_[idx+1].second();
135
136 if (interp_ == DIVERGING)
137 {
138 return colourTools::interpolateDiverging(s, rgb0, rgb1);
139 }
140 else if (interp_ == HSV)
141 {
142 return colourTools::interpolateHSV(s, rgb0, rgb1);
144
145 return colourTools::interpolateRGB(s, rgb0, rgb1);
146}
147
148
150Foam::colourTable::table(const label nColours) const
151{
152 List<Tuple2<scalar, vector>> lut(nColours);
153
154 for (label i=0; i < nColours; ++i)
155 {
156 const scalar x = scalar(i)/scalar(nColours-1);
157
158 lut[i] = pair_type(x, value(x));
159 }
160
161 return lut;
162}
163
164
166{
167 os.beginBlock();
168 os.writeEntry("interpolate", interpolationTypeNames[interp_]);
169 os.writeEntry("table", table_);
170 os.endBlock();
171
172 return os;
173}
174
175
176Foam::Ostream& Foam::operator<<(Ostream& os, const colourTable& tbl)
177{
178 tbl.writeDict(os);
179
180 return os;
181}
182
183
184// ************************************************************************* //
Various functions to operate on Lists.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A 2-tuple for storing two objects of dissimilar types. The container is similar in purpose to std::pa...
Definition Tuple2.H:51
const T1 & first() const noexcept
Access the first element.
Definition Tuple2.H:132
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
Base class for generating a colour table from node points.
Definition colourTable.H:76
vector value(const scalar x) const
Return the colour at x. The input is clipped to 0-1 range.
Definition colourTable.C:87
Ostream & writeDict(Ostream &os) const
Write as dictionary format.
static autoPtr< colourTable > New(Istream &is)
Read as dictionary content.
Definition colourTable.C:79
Tuple2< scalar, vector > pair_type
The data lookup type.
colourTable(const List< Tuple2< scalar, vector > > &values, const interpolationType interp=interpolationType::RGB)
Copy construct from table values.
Definition colourTable.C:42
static const Enum< interpolationType > interpolationTypeNames
Enumeration names for interpolationType.
List< Tuple2< scalar, vector > > table(const label nColours) const
Return a discrete lookup table of colours.
interpolationType
Internal interpolation type.
Definition colourTable.H:83
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
OBJstream os(runTime.globalPath()/outputName)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
vector interpolateRGB(scalar s, const vector &rgb1, const vector &rgb2)
Interpolate RGB values in RGB colourspace.
void interpolateDiverging(scalar s, const vector &rgb1, const vector &rgb2, vector &result)
Interpolate RGB values with diverging color map.
void interpolateHSV(scalar s, const vector &rgb1, const vector &rgb2, vector &result)
Interpolate RGB values in HSV colourspace.
label findLower(const ListType &input, const T &val, const label start, const ComparePredicate &comp)
Binary search to find the index of the last element in a sorted list that is less than value.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
Vector< scalar > vector
Definition vector.H:57
dictionary dict