-
-
Notifications
You must be signed in to change notification settings - Fork 252
/
jakob2019.py
781 lines (646 loc) · 26.8 KB
/
jakob2019.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
# -*- coding: utf-8 -*-
"""
Jakob and Hanika (2019) - Reflectance Recovery
==============================================
Defines objects for reflectance recovery, i.e. spectral upsampling, using
*Jakob and Hanika (2019)* method:
- :func:`colour.recovery.sd_Jakob2019`
- :func:`colour.recovery.find_coefficients_Jakob2019`
- :func:`colour.recovery.XYZ_to_sd_Jakob2019`
- :class:`colour.recovery.Jakob2019Interpolator`
References
----------
- :cite:`Jakob2019` : Jakob, W., & Hanika, J. (2019). A Low‐Dimensional
Function Space for Efficient Spectral Upsampling. Computer Graphics Forum,
38(2), 147–155. doi:10.1111/cgf.13626
"""
from __future__ import division, print_function, unicode_literals
import numpy as np
import struct
from scipy.optimize import minimize
from scipy.interpolate import RegularGridInterpolator
from colour import ILLUMINANT_SDS
from colour.algebra import spow, smoothstep_function
from colour.colorimetry import (
STANDARD_OBSERVER_CMFS, SpectralDistribution, SpectralShape,
intermediate_lightness_function_CIE1976, sd_ones, sd_to_XYZ)
from colour.difference import JND_CIE1976
from colour.models import XYZ_to_xy, XYZ_to_Lab, RGB_to_XYZ
from colour.utilities import (as_float_array, domain_range_scale, full,
index_along_last_axis, to_domain_1,
runtime_warning, zeros)
__author__ = 'Colour Developers'
__copyright__ = 'Copyright (C) 2013-2020 - Colour Developers'
__license__ = 'New BSD License - https://opensource.org/licenses/BSD-3-Clause'
__maintainer__ = 'Colour Developers'
__email__ = 'colour-developers@colour-science.org'
__status__ = 'Production'
__all__ = [
'JAKOB2019_SPECTRAL_SHAPE', 'StopMinimizationEarly', 'sd_Jakob2019',
'error_function', 'dimensionalise_coefficients', 'lightness_scale',
'find_coefficients_Jakob2019', 'XYZ_to_sd_Jakob2019',
'Jakob2019Interpolator'
]
JAKOB2019_SPECTRAL_SHAPE = SpectralShape(360, 780, 5)
"""
Spectral shape for *Jakob and Hanika (2019)* method.
JAKOB2019_SPECTRAL_SHAPE : SpectralShape
"""
class StopMinimizationEarly(Exception):
"""
The exception used to stop :func:`scipy.optimize.minimize` once the
value of the minimized function is small enough. *SciPy* doesn't currently
offer a better way of doing it.
Attributes
----------
coefficients : ndarray, (3,)
Coefficients (function arguments) when this exception was raised.
error : float
Error (function value) when this exception was raised.
"""
def __init__(self, coefficients, error):
self.coefficients = coefficients
self.error = error
def sd_Jakob2019(coefficients, shape=JAKOB2019_SPECTRAL_SHAPE):
"""
Returns a spectral distribution following the spectral model given by
*Jakob and Hanika (2019)*.
Parameters
----------
coefficients : array_like
Dimensionless coefficients for *Jakob and Hanika (2019)* reflectance
spectral model.
shape : SpectralShape, optional
Shape used by the spectral distribution.
Returns
-------
SpectralDistribution
*Jakob and Hanika (2019)* spectral distribution.
Examples
--------
>>> from colour.utilities import numpy_print_options
>>> with numpy_print_options(suppress=True):
... sd_Jakob2019([-9e-05, 8.5e-02, -20], SpectralShape(400, 700, 20))
... # doctest: +ELLIPSIS
SpectralDistribution([[ 400. , 0.3143046...],
[ 420. , 0.4133320...],
[ 440. , 0.4880034...],
[ 460. , 0.5279562...],
[ 480. , 0.5319346...],
[ 500. , 0.5 ...],
[ 520. , 0.4326202...],
[ 540. , 0.3373544...],
[ 560. , 0.2353056...],
[ 580. , 0.1507665...],
[ 600. , 0.0931332...],
[ 620. , 0.0577434...],
[ 640. , 0.0367011...],
[ 660. , 0.0240879...],
[ 680. , 0.0163316...],
[ 700. , 0.0114118...]],
interpolator=SpragueInterpolator,
interpolator_kwargs={},
extrapolator=Extrapolator,
extrapolator_kwargs={...})
"""
c_0, c_1, c_2 = as_float_array(coefficients)
wl = shape.range()
U = c_0 * wl ** 2 + c_1 * wl + c_2
R = 1 / 2 + U / (2 * np.sqrt(1 + U ** 2))
name = 'Jakob (2019) - {0} (COEFF)'.format(coefficients)
return SpectralDistribution(R, wl, name=name)
def error_function(coefficients,
target,
cmfs,
illuminant,
max_error=None,
additional_data=False):
"""
Computes :math:`\\Delta E_{76}` between the target colour and the colour
defined by given spectral model, along with its gradient.
Parameters
----------
coefficients : array_like
Dimensionless coefficients for *Jakob and Hanika (2019)* reflectance
spectral model.
target : array_like, (3,)
*CIE L\\*a\\*b\\** colourspace array of the target colour.
cmfs : XYZ_ColourMatchingFunctions
Standard observer colour matching functions.
illuminant : SpectralDistribution
Illuminant spectral distribution.
max_error : float, optional
Raise ``StopMinimizationEarly`` if the error is smaller than this.
The default is *None* and the function doesn't raise anything.
additional_data : bool, optional
If *True*, some intermediate calculations are returned, for use in
correctness tests: R, XYZ and Lab.
Returns
-------
error : float
The computed :math:`\\Delta E_{76}` error.
derror : ndarray, (3,)
The gradient of error, i.e. the first derivatives of error with respect
to the input coefficients.
R : ndarray
Computed spectral reflectance.
XYZ : ndarray, (3,)
*CIE XYZ* tristimulus values corresponding to ``R``.
Lab : ndarray, (3,)
*CIE L\\*a\\*b\\** colourspace array corresponding to ``XYZ``.
Raises
------
StopMinimizationEarly
Raised when the error is below ``max_error``.
"""
c_0, c_1, c_2 = as_float_array(coefficients)
wv = np.linspace(0, 1, len(cmfs.shape))
U = c_0 * wv ** 2 + c_1 * wv + c_2
t1 = np.sqrt(1 + U ** 2)
R = 1 / 2 + U / (2 * t1)
t2 = 1 / (2 * t1) - U ** 2 / (2 * t1 ** 3)
dR = np.array([wv ** 2 * t2, wv * t2, t2])
E = illuminant.values * R
dE = illuminant.values * dR
dw = cmfs.wavelengths[1] - cmfs.wavelengths[0]
k = 1 / (np.sum(cmfs.values[:, 1] * illuminant.values) * dw)
XYZ = k * np.dot(E, cmfs.values) * dw
dXYZ = np.transpose(k * np.dot(dE, cmfs.values) * dw)
XYZ_n = sd_to_XYZ(illuminant, cmfs)
XYZ_n /= XYZ_n[1]
XYZ_XYZ_n = XYZ / XYZ_n
XYZ_f = intermediate_lightness_function_CIE1976(XYZ, XYZ_n)
dXYZ_f = np.where(
XYZ_XYZ_n[..., np.newaxis] > (24 / 116) ** 3,
1 / (3 * spow(XYZ_n[..., np.newaxis], 1 / 3) * spow(
XYZ[..., np.newaxis], 2 / 3)) * dXYZ,
(841 / 108) * dXYZ / XYZ_n[..., np.newaxis],
)
def intermediate_XYZ_to_Lab(XYZ_i, offset=16):
"""
Returns the final intermediate value for the *CIE Lab* to *CIE XYZ*
conversion.
"""
return np.array([
116 * XYZ_i[1] - offset, 500 * (XYZ_i[0] - XYZ_i[1]),
200 * (XYZ_i[1] - XYZ_i[2])
])
Lab_i = intermediate_XYZ_to_Lab(XYZ_f)
dLab_i = intermediate_XYZ_to_Lab(dXYZ_f, 0)
error = np.sqrt(np.sum((Lab_i - target) ** 2))
if max_error is not None and error <= max_error:
raise StopMinimizationEarly(coefficients, error)
derror = np.sum(
dLab_i * (Lab_i[..., np.newaxis] - target[..., np.newaxis]),
axis=0) / error
if additional_data:
return error, derror, R, XYZ, Lab_i
else:
return error, derror
def dimensionalise_coefficients(coefficients, shape):
"""
Rescales the dimensionless coefficients to given spectral shape.
A dimensionless form of the reflectance spectral model is used in the
optimisation process. Instead of the usual spectral shape, specified in
nanometers, it is normalised to the [0, 1] range. A side effect is that
computed coefficients work only with the normalised range and need to be
rescaled to regain units and be compatible with standard shapes.
Parameters
----------
coefficients : array_like, (3,)
Dimensionless coefficients.
shape : SpectralShape
Spectral distribution shape used in calculations.
Returns
-------
ndarray, (3,)
Dimensionful coefficients, with units of
:math:`\\frac{1}{\\mathrm{nm}^2}`, :math:`\\frac{1}{\\mathrm{nm}}`
and 1, respectively.
"""
cp_0, cp_1, cp_2 = coefficients
span = shape.end - shape.start
c_0 = cp_0 / span ** 2
c_1 = cp_1 / span - 2 * cp_0 * shape.start / span ** 2
c_2 = (
cp_0 * shape.start ** 2 / span ** 2 - cp_1 * shape.start / span + cp_2)
return np.array([c_0, c_1, c_2])
def lightness_scale(steps):
"""
Creates a non-linear lightness scale, as described in *Jakob and Hanika
(2019)*. The spacing between very dark and very bright (and saturated)
colours is made smaller, because in those regions coefficients tend to
change rapidly and a finer resolution is needed.
Parameters
----------
steps : int
Samples/steps count along the non-linear lightness scale.
Returns
-------
ndarray
Non-linear lightness scale.
Examples
--------
>>> lightness_scale(5) # doctest: +ELLIPSIS
array([ 0. , 0.0656127..., 0.5 , 0.9343872..., \
1. ])
"""
linear = np.linspace(0, 1, steps)
return smoothstep_function(smoothstep_function(linear))
def find_coefficients_Jakob2019(
XYZ,
cmfs=STANDARD_OBSERVER_CMFS['CIE 1931 2 Degree Standard Observer']
.copy().align(JAKOB2019_SPECTRAL_SHAPE),
illuminant=ILLUMINANT_SDS['D65'].copy().align(
JAKOB2019_SPECTRAL_SHAPE),
coefficients_0=zeros(3),
max_error=JND_CIE1976,
dimensionalise=True):
"""
Computes the coefficients for *Jakob and Hanika (2019)* reflectance
spectral model.
Parameters
----------
XYZ : array_like, (3,)
*CIE XYZ* tristimulus values to find the coefficients for.
cmfs : XYZ_ColourMatchingFunctions
Standard observer colour matching functions.
illuminant : SpectralDistribution
Illuminant spectral distribution.
coefficients_0 : array_like, (3,), optional
Starting coefficients for the solver.
max_error : float, optional
Maximal acceptable error. Set higher to save computational time.
If *None*, the solver will keep going until it is very close to the
minimum. The default is ``ACCEPTABLE_DELTA_E``.
dimensionalise : bool, optional
If *True*, returned coefficients are dimensionful and will not work
correctly if fed back as ``coefficients_0``. The default is *True*.
Returns
-------
coefficients : ndarray, (3,)
Computed coefficients that best fit the given colour.
error : float
:math:`\\Delta E_{76}` between the target colour and the colour
corresponding to the computed coefficients.
"""
shape = cmfs.shape
if illuminant.shape != shape:
runtime_warning(
'Aligning "{0}" illuminant shape to "{1}" colour matching '
'functions shape.'.format(illuminant.name, cmfs.name))
illuminant = illuminant.copy().align(cmfs.shape)
def optimize(target_o, coefficients_0_o):
"""
Minimises the error function using *L-BFGS-B* method.
"""
try:
result = minimize(
error_function,
coefficients_0_o, (target_o, cmfs, illuminant, max_error),
method='L-BFGS-B',
jac=True)
return result.x, result.fun
except StopMinimizationEarly as error:
return error.coefficients, error.error
xy_n = XYZ_to_xy(sd_to_XYZ(illuminant))
XYZ_good = full(3, 0.5)
coefficients_good = zeros(3)
divisions = 3
while divisions < 10:
XYZ_r = XYZ_good
coefficient_r = coefficients_good
keep_divisions = False
coefficients_0 = coefficient_r
for i in range(1, divisions):
XYZ_i = (XYZ - XYZ_r) * i / (divisions - 1) + XYZ_r
Lab_i = XYZ_to_Lab(XYZ_i)
coefficients_0, error = optimize(Lab_i, coefficients_0)
if error > max_error:
break
else:
XYZ_good = XYZ_i
coefficients_good = coefficients_0
keep_divisions = True
else:
break
if not keep_divisions:
divisions += 2
target = XYZ_to_Lab(XYZ, xy_n)
coefficients, error = optimize(target, coefficients_0)
if dimensionalise:
coefficients = dimensionalise_coefficients(coefficients, shape)
return coefficients, error
def XYZ_to_sd_Jakob2019(
XYZ,
cmfs=STANDARD_OBSERVER_CMFS['CIE 1931 2 Degree Standard Observer']
.copy().align(JAKOB2019_SPECTRAL_SHAPE),
illuminant=sd_ones(JAKOB2019_SPECTRAL_SHAPE),
optimisation_kwargs=None,
additional_data=False):
"""
Recovers the spectral distribution of given RGB colourspace array
using *Jakob and Hanika (2019)* method.
Parameters
----------
XYZ : array_like, (3,)
*CIE XYZ* tristimulus values to recover the spectral distribution from.
cmfs : XYZ_ColourMatchingFunctions
Standard observer colour matching functions.
illuminant : SpectralDistribution
Illuminant spectral distribution.
optimisation_kwargs : dict_like, optional
Parameters for :func:`colour.recovery.find_coefficients_Jakob2019`
definition.
additional_data : bool, optional
If *True*, ``error`` will be returned alongside ``sd``.
Returns
-------
sd : SpectralDistribution
Recovered spectral distribution.
error : float
:math:`\\Delta E_{76}` between the target colour and the colour
corresponding to the computed coefficients.
Examples
--------
>>> from colour.colorimetry import ILLUMINANTS, sd_to_XYZ_integration
>>> from colour.models import XYZ_to_sRGB
>>> from colour.utilities import numpy_print_options
>>> XYZ = np.array([0.20654008, 0.12197225, 0.05136952])
>>> cmfs = (
... STANDARD_OBSERVER_CMFS['CIE 1931 2 Degree Standard Observer'].
... copy().align(SpectralShape(360, 780, 10))
... )
>>> sd = XYZ_to_sd_Jakob2019(XYZ, cmfs)
>>> with numpy_print_options(suppress=True):
... # Doctests skip for Python 2.x compatibility.
... sd # doctest: +SKIP
SpectralDistribution([[ 360. , 0.3717653...],
[ 370. , 0.2551674...],
[ 380. , 0.1796344...],
[ 390. , 0.1324685...],
[ 400. , 0.1025098...],
[ 410. , 0.0828360...],
[ 420. , 0.0694814...],
[ 430. , 0.0601729...],
[ 440. , 0.0535745...],
[ 450. , 0.0488773...],
[ 460. , 0.0455801...],
[ 470. , 0.0433696...],
[ 480. , 0.0420540...],
[ 490. , 0.0415260...],
[ 500. , 0.0417442...],
[ 510. , 0.0427256...],
[ 520. , 0.0445487...],
[ 530. , 0.0473671...],
[ 540. , 0.0514381...],
[ 550. , 0.0571745...],
[ 560. , 0.0652386...],
[ 570. , 0.0767126...],
[ 580. , 0.0934152...],
[ 590. , 0.1184910...],
[ 600. , 0.1574567...],
[ 610. , 0.2196538...],
[ 620. , 0.3181180...],
[ 630. , 0.4598423...],
[ 640. , 0.6224910...],
[ 650. , 0.7601476...],
[ 660. , 0.8516183...],
[ 670. , 0.9061111...],
[ 680. , 0.9381461...],
[ 690. , 0.9575256...],
[ 700. , 0.9697334...],
[ 710. , 0.9777401...],
[ 720. , 0.9831864...],
[ 730. , 0.9870109...],
[ 740. , 0.9897712...],
[ 750. , 0.9918114...],
[ 760. , 0.9933506...],
[ 770. , 0.9945329...],
[ 780. , 0.9954555...]],
interpolator=SpragueInterpolator,
interpolator_kwargs={},
extrapolator=Extrapolator,
extrapolator_kwargs={...})
>>> sd_to_XYZ_integration(sd) / 100 # doctest: +ELLIPSIS
array([ 0.2065209..., 0.1220029..., 0.0513715...])
"""
XYZ = to_domain_1(XYZ)
if optimisation_kwargs is None:
optimisation_kwargs = {}
with domain_range_scale('ignore'):
coefficients, error = find_coefficients_Jakob2019(
XYZ, cmfs, illuminant, **optimisation_kwargs)
sd = sd_Jakob2019(coefficients, cmfs.shape)
sd.name = 'Jakob (2019) - {0}'.format(XYZ)
if additional_data:
return sd, error
else:
return sd
class Jakob2019Interpolator:
"""
Class for working with pre-computed lookup tables for the
*Jakob and Hanika (2019)* spectral upsampling method. It allows significant
time savings by performing the expensive numerical optimization ahead of
time and storing the results in a file.
The file format is compatible with the code and *.coeff* files in
supplemental material published alongside the article.
Attributes
----------
scale
coefficients
resolution
Methods
-------
RGB_to_coefficients
RGB_to_sd
generate
read
write
"""
def __init__(self):
self.scale = None
self.coefficients = None
self.resolution = None
def _create_cube(self):
"""
Creates a :class:`scipy.interpolate.RegularGridInterpolator` class
instance for read or generated coefficients.
"""
samples = np.linspace(0, 1, self.resolution)
axes = ([0, 1, 2], self.scale, samples, samples)
self.cube = RegularGridInterpolator(
axes, self.coefficients, bounds_error=False)
def RGB_to_coefficients(self, RGB):
"""
Look up a given *RGB* colourspace array and return corresponding
coefficients. Interpolation is used for colours not on the table grid.
Parameters
----------
RGB : ndarray, (3,)
*RGB* colourspace array.
Returns
-------
coefficients : ndarray, (3,)
Corresponding coefficients that can be passed to
:func:`colour.recovery.jakob2019.sd_Jakob2019` to obtain a
spectral distribution.
"""
RGB = as_float_array(RGB)
vmax = np.max(RGB, axis=-1)
chroma = RGB / (np.expand_dims(vmax, -1) + 1e-10)
imax = np.argmax(RGB, axis=-1)
v1 = index_along_last_axis(RGB, imax)
v2 = index_along_last_axis(chroma, (imax + 2) % 3)
v3 = index_along_last_axis(chroma, (imax + 1) % 3)
coords = np.stack([imax, v1, v2, v3], axis=-1)
return self.cube(coords).squeeze()
def RGB_to_sd(self, RGB, shape=JAKOB2019_SPECTRAL_SHAPE):
"""
Looks up a given *RGB* colourspace array and return the corresponding
spectral distribution.
Parameters
----------
RGB : ndarray, (3,)
*RGB* colourspace array.
shape : SpectralShape, optional
Shape used by the spectral distribution.
Returns
-------
sd : SpectralDistribution
Spectral distribution corresponding with the RGB* colourspace
array.
"""
sd = sd_Jakob2019(self.RGB_to_coefficients(RGB), shape)
sd.name = 'Jakob (2019) - {0} (RGB)'.format(RGB)
return sd
def generate(self,
colourspace,
cmfs,
illuminant,
resolution,
verbose=False,
print_callable=print):
"""
Creates a lookup table for a given *RGB* colourspace with given
resolution.
Parameters
----------
colourspace: RGB_Colourspace
The *RGB* colourspace to create a lookup table for.
cmfs : XYZ_ColourMatchingFunctions
Standard observer colour matching functions.
illuminant : SpectralDistribution
Illuminant spectral distribution.
resolution : int
The resolution of the lookup table. Higher values will decrease
errors but at the cost of a much longer run time. The published
*.coeff* files have a resolution of 64.
verbose : bool, optional
If *True*, the default, information about the progress is printed
to the standard output.
print_callable : callable, optional
Callable used to for verbose.
"""
illuminant_xy = XYZ_to_xy(sd_to_XYZ(illuminant))
# It could be interesting to have different resolutions for lightness
# and chromaticity, but the current file format doesn't allow it.
lightness_steps = resolution
chroma_steps = resolution
self.scale = lightness_scale(lightness_steps)
self.coefficients = np.empty((3, chroma_steps, chroma_steps,
lightness_steps, 3))
cube_indexes = np.ndindex(3, chroma_steps, chroma_steps)
# First, create a list of all the fully bright colours with the order
# matching cube_indexes.
samples = np.linspace(0, 1, chroma_steps)
yx = np.meshgrid(*[[1], samples, samples], indexing='ij')
yx = np.transpose(yx).reshape(-1, 3)
chromas = np.concatenate(
[yx, np.roll(yx, 1, axis=1),
np.roll(yx, 2, axis=1)])
# TODO: Replace this with a proper progress bar.
if verbose:
print_callable(
'{0:>6} {1:>6} {2:>6} {3:>13} {4:>13} {5:>13} {6}'.format(
'R', 'G', 'B', 'c0', 'c1', 'c2', 'Delta E'))
# TODO: Send the list to a multiprocessing pool; this takes a while.
for (i, j, k), chroma in zip(cube_indexes, chromas):
if verbose:
print_callable(
'i={0}, j={1}, k={2}, R={3:.6f}, G={4:.6f}, B={5:.6f}'.
format(i, j, k, chroma[0], chroma[1], chroma[2]))
def optimize(L, coefficients_0):
"""
Solves for a specific lightness and stores the result in the
appropriate cell.
"""
RGB = self.scale[L] * chroma
XYZ = RGB_to_XYZ(RGB, colourspace.whitepoint, illuminant_xy,
colourspace.RGB_to_XYZ_matrix)
coefficients, error = find_coefficients_Jakob2019(
XYZ,
cmfs,
illuminant,
coefficients_0,
dimensionalise=False)
if verbose:
print_callable(
'{0:.4f} {1:.4f} {2:.4f} '
'{3:>13.6f} {4:>13.6f} {5:>13.6f} {6:.6f}'.format(
RGB[0], RGB[1], RGB[2], coefficients[0],
coefficients[1], coefficients[2], error))
self.coefficients[i, L, j, k, :] = dimensionalise_coefficients(
coefficients, cmfs.shape)
return coefficients
# Starts from somewhere in the middle, similarly to how feedback
# works in "colour.recovery.find_coefficients_Jakob2019"
# definition.
middle_L = lightness_steps // 3
middle_coefficients = optimize(middle_L, (0, 0, 0))
# Goes down the lightness scale.
coefficients_0 = middle_coefficients
for L in reversed(range(0, middle_L)):
coefficients_0 = optimize(L, coefficients_0)
# Goes up the lightness scale.
coefficients_0 = middle_coefficients
for L in range(middle_L + 1, lightness_steps):
coefficients_0 = optimize(L, coefficients_0)
self.resolution = lightness_steps
self._create_cube()
def read(self, path):
"""
Loads a lookup table from a *.coeff* file.
Parameters
----------
path : unicode
Path to the file.
"""
with open(path, 'rb') as coeff_file:
if coeff_file.read(4).decode('ISO-8859-1') != 'SPEC':
raise ValueError(
'Bad magic number, this is likely not the right file type!'
)
self.resolution = struct.unpack('i', coeff_file.read(4))[0]
self.scale = np.fromfile(
coeff_file, count=self.resolution, dtype=np.float32)
self.coefficients = np.fromfile(
coeff_file,
count=3 * (self.resolution ** 3) * 3,
dtype=np.float32)
self.coefficients = self.coefficients.reshape(
3, self.resolution, self.resolution, self.resolution, 3)
self._create_cube()
def write(self, path):
"""
Writes the lookup table to a *.coeff* file.
Parameters
----------
path : unicode
Path to the file.
"""
with open(path, 'wb') as coeff_file:
coeff_file.write(b'SPEC')
coeff_file.write(struct.pack('i', self.coefficients.shape[1]))
np.float32(self.scale).tofile(coeff_file)
np.float32(self.coefficients).tofile(coeff_file)