colour.recovery.spectral_primary_decomposition_Mallett2019

colour.recovery.spectral_primary_decomposition_Mallett2019(colourspace, cmfs=XYZ_ColourMatchingFunctions(name='CIE 1931 2 Degree Standard Observer', ...), illuminant=SpectralDistribution(name='D65', ...), metric=<function norm>, metric_args=(), optimisation_kwargs=None)[source]

Performs the spectral primary decomposition as described in Mallett and Yuksel (2019) for given RGB colourspace.

Parameters
  • colourspace (RGB_Colourspace) – RGB colourspace.

  • cmfs (XYZ_ColourMatchingFunctions, optional) – Standard observer colour matching functions.

  • illuminant (SpectralDistribution, optional) – Illuminant spectral distribution.

  • metric (unicode, optional) –

    Function to be minimised, i.e. the objective function.

    metric(basis, *metric_args) -> float

    where basis is three reflectances concatenated together, each with a shape matching shape.

  • metric_args (tuple, optional) – Additional arguments passed to metric.

  • optimisation_kwargs (dict_like, optional) – Parameters for scipy.optimize.minimize() definition.

Returns

Basis functions for given RGB colourspace.

Return type

MultiSpectralDistributions

References

[]

Notes

  • In-addition to the BT.709 primaries used by the sRGB colourspace, [] tried BT.2020, P3 D65, Adobe RGB 1998, NTSC (1987), Pal/Secam, ProPhoto RGB, and Adobe Wide Gamut RGB primaries, every one of which encompasses a larger (albeit not-always-enveloping) set of CIE L*a*b* colours than BT.709. Of these, only Pal/Secam produces a feasible basis, which is relatively unsurprising since it is very similar to BT.709, whereas the others are significantly larger.

Examples

>>> from colour.colorimetry import SpectralShape
>>> from colour.models import RGB_COLOURSPACE_PAL_SECAM
>>> from colour.utilities import numpy_print_options
>>> cmfs = (
...     MSDS_CMFS_STANDARD_OBSERVER['CIE 1931 2 Degree Standard Observer'].
...     copy().align(SpectralShape(360, 780, 10))
... )
>>> illuminant = SDS_ILLUMINANTS['D65'].copy().align(cmfs.shape)
>>> msds = spectral_primary_decomposition_Mallett2019(
...     RGB_COLOURSPACE_PAL_SECAM, cmfs, illuminant, optimisation_kwargs={
...         'options': {'ftol': 1e-5}
...     }
... )
>>> with numpy_print_options(suppress=True):
...     # Doctests skip for Python 2.x compatibility.
...     print(msds)  
[[ 360.            0.3324728...    0.3332663...    0.3342608...]
 [ 370.            0.3323307...    0.3327746...    0.3348946...]
 [ 380.            0.3341115...    0.3323995...    0.3334889...]
 [ 390.            0.3337570...    0.3298092...    0.3364336...]
 [ 400.            0.3209352...    0.3218213...    0.3572433...]
 [ 410.            0.2881025...    0.2837628...    0.4281346...]
 [ 420.            0.1836749...    0.1838893...    0.6324357...]
 [ 430.            0.0187212...    0.0529655...    0.9283132...]
 [ 440.            0.       ...    0.       ...    1.       ...]
 [ 450.            0.       ...    0.       ...    1.       ...]
 [ 460.            0.       ...    0.       ...    1.       ...]
 [ 470.            0.       ...    0.0509556...    0.9490443...]
 [ 480.            0.       ...    0.2933996...    0.7066003...]
 [ 490.            0.       ...    0.5001396...    0.4998603...]
 [ 500.            0.       ...    0.6734805...    0.3265195...]
 [ 510.            0.       ...    0.8555213...    0.1444786...]
 [ 520.            0.       ...    0.9999985...    0.0000014...]
 [ 530.            0.       ...    1.       ...    0.       ...]
 [ 540.            0.       ...    1.       ...    0.       ...]
 [ 550.            0.       ...    1.       ...    0.       ...]
 [ 560.            0.       ...    0.9924229...    0.       ...]
 [ 570.            0.       ...    0.9913344...    0.0083032...]
 [ 580.            0.0370289...    0.9145168...    0.0484542...]
 [ 590.            0.7100075...    0.2898477...    0.0001446...]
 [ 600.            1.       ...    0.       ...    0.       ...]
 [ 610.            1.       ...    0.       ...    0.       ...]
 [ 620.            1.       ...    0.       ...    0.       ...]
 [ 630.            1.       ...    0.       ...    0.       ...]
 [ 640.            0.9711347...    0.0137659...    0.0150993...]
 [ 650.            0.7996619...    0.1119379...    0.0884001...]
 [ 660.            0.6064640...    0.202815 ...    0.1907209...]
 [ 670.            0.4662959...    0.2675037...    0.2662005...]
 [ 680.            0.4010958...    0.2998989...    0.2990052...]
 [ 690.            0.3617485...    0.3208921...    0.3173592...]
 [ 700.            0.3496691...    0.3247855...    0.3255453...]
 [ 710.            0.3433979...    0.3273540...    0.329248 ...]
 [ 720.            0.3358860...    0.3345583...    0.3295556...]
 [ 730.            0.3349498...    0.3314232...    0.3336269...]
 [ 740.            0.3359954...    0.3340147...    0.3299897...]
 [ 750.            0.3310392...    0.3327595...    0.3362012...]
 [ 760.            0.3346883...    0.3314158...    0.3338957...]
 [ 770.            0.3332167...    0.333371 ...    0.3334122...]
 [ 780.            0.3319670...    0.3325476...    0.3354852...]]