Abstract profile - Gamut improvment (#6665)
* Gamut control when the chosen primaries are different from working profile * Gamut control abstract * Gamut label and history * Change to Wx Wz * Fixed crash if y primaries are set to zero * Fomated with Astylert ImProcFunctions::workingtrc and Color::primaries_to_xyz * Fixed black becomes green wit gamt abstract profile * Harmonize types in color.cc * Try to fix Multiplication result converted to larger type
This commit is contained in:
parent
21afbaf90b
commit
7afdfc1e2b
@ -1421,6 +1421,7 @@ HISTORY_MSG_ICM_AINTENT;Abstract profile intent
|
||||
HISTORY_MSG_ICM_BLUX;Primaries Blue X
|
||||
HISTORY_MSG_ICM_BLUY;Primaries Blue Y
|
||||
HISTORY_MSG_ICM_FBW;Black and White
|
||||
HISTORY_MSG_ICM_GAMUT;Gamut control
|
||||
HISTORY_MSG_ICM_GREX;Primaries Green X
|
||||
HISTORY_MSG_ICM_GREY;Primaries Green Y
|
||||
HISTORY_MSG_ICM_OUTPUT_PRIMARIES;Output - Primaries
|
||||
@ -2529,6 +2530,7 @@ TP_ICM_DCPILLUMINANT;Illuminant
|
||||
TP_ICM_DCPILLUMINANT_INTERPOLATED;Interpolated
|
||||
TP_ICM_DCPILLUMINANT_TOOLTIP;Select which embedded DCP illuminant to employ. Default is 'interpolated' which is a mix between the two based on white balance. The setting is only available if a dual-illuminant DCP with interpolation support is selected.
|
||||
TP_ICM_FBW;Black-and-White
|
||||
TP_ICM_GAMUT;Gamut control
|
||||
TP_ICM_ILLUMPRIM_TOOLTIP;Choose the illuminant closest to the shooting conditions.\nChanges can only be made when the 'Destination primaries' selection is set to 'Custom (sliders)'.
|
||||
TP_ICM_INPUTCAMERA;Camera standard
|
||||
TP_ICM_INPUTCAMERAICC;Auto-matched camera profile
|
||||
|
@ -24,6 +24,7 @@
|
||||
#include "sleef.h"
|
||||
#include "opthelper.h"
|
||||
#include "iccstore.h"
|
||||
#include <iostream>
|
||||
|
||||
using namespace std;
|
||||
|
||||
@ -1910,6 +1911,152 @@ void Color::Lch2Luv(float c, float h, float &u, float &v)
|
||||
v = c * sincosval.y;
|
||||
}
|
||||
|
||||
void Color::primaries_to_xyz(double p[6], double Wx, double Wz, double *pxyz)
|
||||
{
|
||||
//calculate Xr, Xg, Xb, Yr, Yb, Tg, Zr,Zg Zb
|
||||
double Wy = 1.0;
|
||||
double Xr = p[0] / p[1];
|
||||
double Yr = 1.0;
|
||||
double Zr = (1.0 - p[0] - p[1]) / p[1];
|
||||
double Xg = p[2] / p[3];
|
||||
double Yg = 1.0;
|
||||
double Zg = (1.0 - p[2] - p[3]) / p[3];
|
||||
double Xb = p[4] / p[5];
|
||||
double Yb = 1.0;
|
||||
double Zb = (1.0 - p[4] - p[5]) / p[5];
|
||||
|
||||
using Triple = std::array<double, 3>;
|
||||
|
||||
using Matrix = std::array<Triple, 3>;
|
||||
|
||||
Matrix input_prim;
|
||||
Matrix inv_input_prim = {};
|
||||
input_prim[0][0] = Xr;
|
||||
input_prim[0][1] = Yr;
|
||||
input_prim[0][2] = Zr;
|
||||
input_prim[1][0] = Xg;
|
||||
input_prim[1][1] = Yg;
|
||||
input_prim[1][2] = Zg;
|
||||
input_prim[2][0] = Xb;
|
||||
input_prim[2][1] = Yb;
|
||||
input_prim[2][2] = Zb;
|
||||
|
||||
//invert matrix
|
||||
if (!rtengine::invertMatrix(input_prim, inv_input_prim)) {
|
||||
std::cout << "Matrix is not invertible, skipping" << std::endl;
|
||||
}
|
||||
|
||||
//white point D50 used by LCMS
|
||||
double Wdx = 0.96420;
|
||||
double Wdy = 1.0;
|
||||
double Wdz = 0.82490;
|
||||
|
||||
double Sr = Wx * inv_input_prim [0][0] + Wy * inv_input_prim [1][0] + Wz * inv_input_prim [2][0];
|
||||
double Sg = Wx * inv_input_prim [0][1] + Wy * inv_input_prim [1][1] + Wz * inv_input_prim [2][1];
|
||||
double Sb = Wx * inv_input_prim [0][2] + Wy * inv_input_prim [1][2] + Wz * inv_input_prim [2][2];
|
||||
|
||||
//XYZ matrix for primaries and temp
|
||||
Matrix mat_xyz = {};
|
||||
mat_xyz[0][0] = Sr * Xr;
|
||||
mat_xyz[0][1] = Sr * Yr;
|
||||
mat_xyz[0][2] = Sr * Zr;
|
||||
mat_xyz[1][0] = Sg * Xg;
|
||||
mat_xyz[1][1] = Sg * Yg;
|
||||
mat_xyz[1][2] = Sg * Zg;
|
||||
mat_xyz[2][0] = Sb * Xb;
|
||||
mat_xyz[2][1] = Sb * Yb;
|
||||
mat_xyz[2][2] = Sb * Zb;
|
||||
|
||||
//chromatic adaptation Bradford
|
||||
Matrix MaBradford = {};
|
||||
MaBradford[0][0] = 0.8951;
|
||||
MaBradford[0][1] = -0.7502;
|
||||
MaBradford[0][2] = 0.0389;
|
||||
MaBradford[1][0] = 0.2664;
|
||||
MaBradford[1][1] = 1.7135;
|
||||
MaBradford[1][2] = -0.0685;
|
||||
MaBradford[2][0] = -0.1614;
|
||||
MaBradford[2][1] = 0.0367;
|
||||
MaBradford[2][2] = 1.0296;
|
||||
|
||||
Matrix Ma_oneBradford = {};
|
||||
Ma_oneBradford[0][0] = 0.9869929;
|
||||
Ma_oneBradford[0][1] = 0.4323053;
|
||||
Ma_oneBradford[0][2] = -0.0085287;
|
||||
Ma_oneBradford[1][0] = -0.1470543;
|
||||
Ma_oneBradford[1][1] = 0.5183603;
|
||||
Ma_oneBradford[1][2] = 0.0400428;
|
||||
Ma_oneBradford[2][0] = 0.1599627;
|
||||
Ma_oneBradford[2][1] = 0.0492912;
|
||||
Ma_oneBradford[2][2] = 0.9684867;
|
||||
|
||||
//R G B source
|
||||
double Rs = Wx * MaBradford[0][0] + Wy * MaBradford[1][0] + Wz * MaBradford[2][0];
|
||||
double Gs = Wx * MaBradford[0][1] + Wy * MaBradford[1][1] + Wz * MaBradford[2][1];
|
||||
double Bs = Wx * MaBradford[0][2] + Wy * MaBradford[1][2] + Wz * MaBradford[2][2];
|
||||
|
||||
// R G B destination
|
||||
double Rd = Wdx * MaBradford[0][0] + Wdy * MaBradford[1][0] + Wdz * MaBradford[2][0];
|
||||
double Gd = Wdx * MaBradford[0][1] + Wdy * MaBradford[1][1] + Wdz * MaBradford[2][1];
|
||||
double Bd = Wdx * MaBradford[0][2] + Wdy * MaBradford[1][2] + Wdz * MaBradford[2][2];
|
||||
|
||||
//cone destination
|
||||
Matrix cone_dest_sourc = {};
|
||||
cone_dest_sourc [0][0] = Rd / Rs;
|
||||
cone_dest_sourc [0][1] = 0.;
|
||||
cone_dest_sourc [0][2] = 0.;
|
||||
cone_dest_sourc [1][0] = 0.;
|
||||
cone_dest_sourc [1][1] = Gd / Gs;
|
||||
cone_dest_sourc [1][2] = 0.;
|
||||
cone_dest_sourc [2][0] = 0.;
|
||||
cone_dest_sourc [2][1] = 0.;
|
||||
cone_dest_sourc [2][2] = Bd / Bs;
|
||||
|
||||
//cone dest
|
||||
Matrix cone_ma_one = {};
|
||||
|
||||
for (int i = 0; i < 3; ++i) {
|
||||
for (int j = 0; j < 3; ++j) {
|
||||
cone_ma_one[i][j] = 0;
|
||||
|
||||
for (int k = 0; k < 3; ++k) {
|
||||
cone_ma_one[i][j] += cone_dest_sourc [i][k] * Ma_oneBradford[k][j];
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
//generate adaptation bradford matrix
|
||||
Matrix adapt_chroma = {};
|
||||
|
||||
for (int i = 0; i < 3; ++i) {
|
||||
for (int j = 0; j < 3; ++j) {
|
||||
adapt_chroma [i][j] = 0;
|
||||
|
||||
for (int k = 0; k < 3; ++k) {
|
||||
adapt_chroma[i][j] += MaBradford[i][k] * cone_ma_one[k][j];
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
Matrix mat_xyz_brad = {};
|
||||
|
||||
for (int i = 0; i < 3; ++i) {
|
||||
for (int j = 0; j < 3; ++j) {
|
||||
mat_xyz_brad[i][j] = 0;
|
||||
|
||||
for (int k = 0; k < 3; ++k) {
|
||||
mat_xyz_brad[i][j] += mat_xyz[i][k] * adapt_chroma[k][j];
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
//push result in pxyz
|
||||
for (int i = 0; i < 3; ++i) {
|
||||
for (int j = 0; j < 3; ++j) {
|
||||
pxyz[i * 3 + j] = mat_xyz_brad[i][j];
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
/*
|
||||
* Gamut mapping algorithm
|
||||
@ -1931,13 +2078,19 @@ void Color::Lch2Luv(float c, float h, float &u, float &v)
|
||||
* columns of the matrix p=xyz_rgb are RGB tristimulus primaries in XYZ
|
||||
* c is the color fixed on the boundary; and m=0 for c=0, m=1 for c=255
|
||||
*/
|
||||
|
||||
void Color::gamutmap(float &X, float Y, float &Z, const double p[3][3])
|
||||
{
|
||||
float u = 4 * X / (X + 15 * Y + 3 * Z) - u0;
|
||||
float v = 9 * Y / (X + 15 * Y + 3 * Z) - v0;
|
||||
|
||||
float epsil = 0.0001f;
|
||||
float intermXYZ = X + 15 * Y + 3 * Z;
|
||||
if(intermXYZ <= 0.f) {
|
||||
intermXYZ = epsil;
|
||||
}
|
||||
|
||||
float u = 4 * X / (intermXYZ) - u0;
|
||||
float v = 9 * Y / (intermXYZ) - v0;
|
||||
float lam[3][2];
|
||||
float lam_min = 1.0;
|
||||
float lam_min = 1.0f;
|
||||
|
||||
for (int c = 0; c < 3; c++)
|
||||
for (int m = 0; m < 2; m++) {
|
||||
@ -1955,17 +2108,24 @@ void Color::gamutmap(float &X, float Y, float &Z, const double p[3][3])
|
||||
p[0][c] * (5 * Y * p[1][c1] + m * 65535 * p[1][c1] * p[2][c2] + Y * p[2][c1] - m * 65535 * p[1][c2] * p[2][c1]) +
|
||||
m * 65535 * p[0][c2] * (p[1][c1] * p[2][c] - p[1][c] * p[2][c1])));
|
||||
|
||||
if (lam[c][m] < lam_min && lam[c][m] > 0) {
|
||||
if (lam[c][m] < lam_min && lam[c][m] > 0.f) {
|
||||
lam_min = lam[c][m];
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
u = u * lam_min + u0;
|
||||
v = v * lam_min + v0;
|
||||
u = u * (double) lam_min + u0;
|
||||
v = v * (double) lam_min + v0;
|
||||
|
||||
X = (9 * u * Y) / (4 * v);
|
||||
Z = (12 - 3 * u - 20 * v) * Y / (4 * v);
|
||||
float intermuv = 12 - 3 * u - 20 * v;
|
||||
if(intermuv < 0.f) {
|
||||
intermuv = 0.f;
|
||||
}
|
||||
Z = (intermuv) * Y / (4 * v);
|
||||
|
||||
|
||||
|
||||
}
|
||||
|
||||
void Color::skinredfloat ( float J, float h, float sres, float Sp, float dred, float protect_red, int sk, float rstprotection, float ko, float &s)
|
||||
|
@ -1847,6 +1847,13 @@ static inline void Lab2XYZ(vfloat L, vfloat a, vfloat b, vfloat &x, vfloat &y, v
|
||||
*/
|
||||
static void gamutmap(float &X, float Y, float &Z, const double p[3][3]);
|
||||
|
||||
/**
|
||||
* @brief Convert primaries in XYZ values in function of illuminant
|
||||
* @param p primaries red, gree, blue
|
||||
* @param Wx Wy white for illuminant
|
||||
* @param pxyz return matrix XYZ
|
||||
*/
|
||||
static void primaries_to_xyz (double p[6], double Wx, double Wz, double *pxyz);
|
||||
|
||||
/**
|
||||
* @brief Get HSV's hue from the Lab's hue
|
||||
|
@ -416,6 +416,8 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
{
|
||||
const TMatrix wprof = ICCStore::getInstance()->workingSpaceMatrix(params->icm.workingProfile);
|
||||
|
||||
double wprofprim[3][3];//store primaries to XYZ
|
||||
bool gamutcontrol = params->icm.gamut;
|
||||
const float toxyz[3][3] = {
|
||||
{
|
||||
static_cast<float>(wprof[0][0] / ((normalizeIn ? 65535.0 : 1.0))), //I have suppressed / Color::D50x
|
||||
@ -440,57 +442,64 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
if (settings->verbose) {
|
||||
printf("profile not accepted\n");
|
||||
}
|
||||
|
||||
return;
|
||||
}
|
||||
|
||||
if (mul == -5 && gampos == 2.4 && slpos == 12.92310) {//must be change if we change settings RT sRGB
|
||||
//only in this case we can shortcut..all process..no gamut control..because we reduce...leads to very small differences, but big speedup
|
||||
#ifdef _OPENMP
|
||||
#pragma omp parallel for schedule(dynamic, 16) if (multiThread)
|
||||
#pragma omp parallel for schedule(dynamic, 16) if (multiThread)
|
||||
#endif
|
||||
|
||||
for (int i = 0; i < ch; ++i)
|
||||
for (int j = 0; j < cw; ++j) {
|
||||
float r = src->r(i, j);
|
||||
float g = src->g(i, j);
|
||||
float b = src->b(i, j);
|
||||
r = (Color::igammatab_srgb[r]) / 65535.f;
|
||||
g = (Color::igammatab_srgb[g]) / 65535.f;
|
||||
b = (Color::igammatab_srgb[b]) / 65535.f;
|
||||
dst->r(i, j) = r;
|
||||
dst->g(i, j) = g;
|
||||
dst->b(i, j) = b;
|
||||
}
|
||||
return;
|
||||
for (int i = 0; i < ch; ++i)
|
||||
for (int j = 0; j < cw; ++j) {
|
||||
float r = src->r(i, j);
|
||||
float g = src->g(i, j);
|
||||
float b = src->b(i, j);
|
||||
r = (Color::igammatab_srgb[r]) / 65535.f;
|
||||
g = (Color::igammatab_srgb[g]) / 65535.f;
|
||||
b = (Color::igammatab_srgb[b]) / 65535.f;
|
||||
dst->r(i, j) = r;
|
||||
dst->g(i, j) = g;
|
||||
dst->b(i, j) = b;
|
||||
}
|
||||
|
||||
return;
|
||||
|
||||
}
|
||||
|
||||
if (mul == 1 ||(params->icm.wprim == ColorManagementParams::Primaries::DEFAULT && params->icm.will == ColorManagementParams::Illuminant::DEFAULT)) {//shortcut and speedup when no call primaries and illuminant - no gamut control...in this case be careful
|
||||
|
||||
if (mul == 1 || (params->icm.wprim == ColorManagementParams::Primaries::DEFAULT && params->icm.will == ColorManagementParams::Illuminant::DEFAULT)) { //shortcut and speedup when no call primaries and illuminant - no gamut control...in this case be careful
|
||||
GammaValues g_a; //gamma parameters
|
||||
double pwr = 1.0 / static_cast<double>(gampos);
|
||||
Color::calcGamma(pwr, slpos, g_a); // call to calcGamma with selected gamma and slope
|
||||
|
||||
#ifdef _OPENMP
|
||||
# pragma omp parallel for schedule(dynamic,16) if (multiThread)
|
||||
# pragma omp parallel for schedule(dynamic,16) if (multiThread)
|
||||
#endif
|
||||
|
||||
for (int y = 0; y < ch; ++y) {
|
||||
int x = 0;
|
||||
#ifdef __SSE2__
|
||||
|
||||
for (; x < cw - 3; x += 4) {
|
||||
STVFU(dst->r(y,x), F2V(65536.f) * gammalog(LVFU(src->r(y,x)), F2V(gampos), F2V(slpos), F2V(g_a[3]), F2V(g_a[4])));
|
||||
STVFU(dst->g(y,x), F2V(65536.f) * gammalog(LVFU(src->g(y,x)), F2V(gampos), F2V(slpos), F2V(g_a[3]), F2V(g_a[4])));
|
||||
STVFU(dst->b(y,x), F2V(65536.f) * gammalog(LVFU(src->b(y,x)), F2V(gampos), F2V(slpos), F2V(g_a[3]), F2V(g_a[4])));
|
||||
}
|
||||
STVFU(dst->r(y, x), F2V(65536.f) * gammalog(LVFU(src->r(y, x)), F2V(gampos), F2V(slpos), F2V(g_a[3]), F2V(g_a[4])));
|
||||
STVFU(dst->g(y, x), F2V(65536.f) * gammalog(LVFU(src->g(y, x)), F2V(gampos), F2V(slpos), F2V(g_a[3]), F2V(g_a[4])));
|
||||
STVFU(dst->b(y, x), F2V(65536.f) * gammalog(LVFU(src->b(y, x)), F2V(gampos), F2V(slpos), F2V(g_a[3]), F2V(g_a[4])));
|
||||
}
|
||||
|
||||
#endif
|
||||
|
||||
for (; x < cw; ++x) {
|
||||
dst->r(y,x) = 65536.f * gammalog(src->r(y,x), gampos, slpos, g_a[3], g_a[4]);
|
||||
dst->g(y,x) = 65536.f * gammalog(src->g(y,x), gampos, slpos, g_a[3], g_a[4]);
|
||||
dst->b(y,x) = 65536.f * gammalog(src->b(y,x), gampos, slpos, g_a[3], g_a[4]);
|
||||
dst->r(y, x) = 65536.f * gammalog(src->r(y, x), gampos, slpos, g_a[3], g_a[4]);
|
||||
dst->g(y, x) = 65536.f * gammalog(src->g(y, x), gampos, slpos, g_a[3], g_a[4]);
|
||||
dst->b(y, x) = 65536.f * gammalog(src->b(y, x), gampos, slpos, g_a[3], g_a[4]);
|
||||
}
|
||||
}
|
||||
|
||||
return;
|
||||
}
|
||||
|
||||
|
||||
|
||||
float redxx = params->icm.redx;
|
||||
float redyy = params->icm.redy;
|
||||
@ -498,6 +507,7 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
float bluyy = params->icm.bluy;
|
||||
float grexx = params->icm.grex;
|
||||
float greyy = params->icm.grey;
|
||||
float epsil = 0.0001f;
|
||||
|
||||
if (prim == 12) {//convert datas area to xy
|
||||
float redgraphx = params->icm.labgridcieALow;
|
||||
@ -519,11 +529,34 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
greyy = 0.55f * (gregraphy + 1.f) - 0.1f;
|
||||
greyy = rtengine::LIM(greyy, 0.5f, 1.f);
|
||||
}
|
||||
|
||||
//fixed crash when there is no space or too small..just a line...Possible if bx, by aligned with Gx,Gy Rx,Ry
|
||||
float ac = (greyy - redyy) / (grexx - redxx);
|
||||
//fix crash if user select 0 for redyy, bluyy, greyy
|
||||
if (redyy == 0.f) {
|
||||
redyy = epsil;
|
||||
}
|
||||
|
||||
if (bluyy == 0.f) {
|
||||
bluyy = epsil;
|
||||
}
|
||||
|
||||
if (greyy == 0.f) {
|
||||
greyy = epsil;
|
||||
}
|
||||
|
||||
//fix crash if grexx - redxx = 0
|
||||
float grered = 1.f;
|
||||
grered = grexx - redxx;
|
||||
|
||||
if (grered == 0.f) {
|
||||
grered = epsil;
|
||||
}
|
||||
|
||||
float ac = (greyy - redyy) / grered;
|
||||
float bc = greyy - ac * grexx;
|
||||
float yc = ac * bluxx + bc;
|
||||
if ((bluyy < yc + 0.0004f) && (bluyy > yc - 0.0004f)) {//under 0.0004 in some case crash because space too small
|
||||
|
||||
if ((bluyy < yc + 0.0004f) && (bluyy > yc - 0.0004f)) { //under 0.0004 in some case crash because space too small
|
||||
return;
|
||||
}
|
||||
|
||||
@ -564,7 +597,6 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
}
|
||||
|
||||
case ColorManagementParams::Primaries::ACES_P0: {
|
||||
profile = "ACESp0";
|
||||
break;
|
||||
}
|
||||
|
||||
@ -593,11 +625,13 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
break;
|
||||
}
|
||||
}
|
||||
|
||||
if (settings->verbose && prim != 0) {
|
||||
printf("prim=%i Profile Destination=%s\n", prim, profile.c_str());
|
||||
}
|
||||
|
||||
if (settings->verbose && prim != 0) {
|
||||
printf("prim=%i Profile Destination=%s\n", prim, profile.c_str());
|
||||
}
|
||||
|
||||
cmsHTRANSFORM hTransform = nullptr;
|
||||
|
||||
if (transform) {
|
||||
hTransform = transform;
|
||||
} else {
|
||||
@ -622,7 +656,9 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
|
||||
};
|
||||
double tempv4 = 5003.;
|
||||
float p[6]; //primaries
|
||||
double p[6]; //primaries
|
||||
double Wx = 1.0;
|
||||
double Wz = 1.0;
|
||||
|
||||
//primaries for 10 working profiles ==> output profiles
|
||||
if (profile == "WideGamut") {
|
||||
@ -633,6 +669,9 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
p[4] = 0.1570;
|
||||
p[5] = 0.0180;
|
||||
illum = toUnderlying(ColorManagementParams::Illuminant::D50);
|
||||
Wx = 0.964295676;
|
||||
Wz = 0.825104603;
|
||||
|
||||
} else if (profile == "Adobe RGB") {
|
||||
p[0] = 0.6400; //Adobe primaries
|
||||
p[1] = 0.3300;
|
||||
@ -642,6 +681,9 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
p[5] = 0.0600;
|
||||
tempv4 = 6504.;
|
||||
illum = toUnderlying(ColorManagementParams::Illuminant::D65);
|
||||
Wx = 0.95045471;
|
||||
Wz = 1.08905029;
|
||||
|
||||
} else if (profile == "sRGB") {
|
||||
p[0] = 0.6400; // sRGB primaries
|
||||
p[1] = 0.3300;
|
||||
@ -651,6 +693,9 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
p[5] = 0.0600;
|
||||
tempv4 = 6504.;
|
||||
illum = toUnderlying(ColorManagementParams::Illuminant::D65);
|
||||
Wx = 0.95045471;
|
||||
Wz = 1.08905029;
|
||||
|
||||
} else if (profile == "BruceRGB") {
|
||||
p[0] = 0.6400; // Bruce primaries
|
||||
p[1] = 0.3300;
|
||||
@ -660,7 +705,10 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
p[5] = 0.0600;
|
||||
tempv4 = 6504.;
|
||||
illum = toUnderlying(ColorManagementParams::Illuminant::D65);
|
||||
} else if (profile == "Beta RGB") {
|
||||
Wx = 0.95045471;
|
||||
Wz = 1.08905029;
|
||||
|
||||
} else if (profile == "Beta RGB") {
|
||||
p[0] = 0.6888; // Beta primaries
|
||||
p[1] = 0.3112;
|
||||
p[2] = 0.1986;
|
||||
@ -668,6 +716,9 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
p[4] = 0.1265;
|
||||
p[5] = 0.0352;
|
||||
illum = toUnderlying(ColorManagementParams::Illuminant::D50);
|
||||
Wx = 0.964295676;
|
||||
Wz = 0.825104603;
|
||||
|
||||
} else if (profile == "BestRGB") {
|
||||
p[0] = 0.7347; // Best primaries
|
||||
p[1] = 0.2653;
|
||||
@ -676,6 +727,9 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
p[4] = 0.1300;
|
||||
p[5] = 0.0350;
|
||||
illum = toUnderlying(ColorManagementParams::Illuminant::D50);
|
||||
Wx = 0.964295676;
|
||||
Wz = 0.825104603;
|
||||
|
||||
} else if (profile == "Rec2020") {
|
||||
p[0] = 0.7080; // Rec2020 primaries
|
||||
p[1] = 0.2920;
|
||||
@ -685,6 +739,9 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
p[5] = 0.0460;
|
||||
tempv4 = 6504.;
|
||||
illum = toUnderlying(ColorManagementParams::Illuminant::D65);
|
||||
Wx = 0.95045471;
|
||||
Wz = 1.08905029;
|
||||
|
||||
} else if (profile == "ACESp0") {
|
||||
p[0] = 0.7347; // ACES P0 primaries
|
||||
p[1] = 0.2653;
|
||||
@ -694,6 +751,9 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
p[5] = -0.0770;
|
||||
tempv4 = 6004.;
|
||||
illum = toUnderlying(ColorManagementParams::Illuminant::D60);
|
||||
Wx = 0.952646075;
|
||||
Wz = 1.008825184;
|
||||
|
||||
} else if (profile == "ACESp1") {
|
||||
p[0] = 0.713; // ACES P1 primaries
|
||||
p[1] = 0.293;
|
||||
@ -703,6 +763,9 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
p[5] = 0.044;
|
||||
tempv4 = 6004.;
|
||||
illum = toUnderlying(ColorManagementParams::Illuminant::D60);
|
||||
Wx = 0.952646075;
|
||||
Wz = 1.008825184;
|
||||
|
||||
} else if (profile == "ProPhoto") {
|
||||
p[0] = 0.7347; //ProPhoto and default primaries
|
||||
p[1] = 0.2653;
|
||||
@ -711,8 +774,11 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
p[4] = 0.0366;
|
||||
p[5] = 0.0001;
|
||||
illum = toUnderlying(ColorManagementParams::Illuminant::D50);
|
||||
Wx = 0.964295676;
|
||||
Wz = 0.825104603;
|
||||
|
||||
} else if (profile == "Custom") {
|
||||
p[0] = redxx;
|
||||
p[0] = redxx;
|
||||
p[1] = redyy;
|
||||
p[2] = grexx;
|
||||
p[3] = greyy;
|
||||
@ -743,11 +809,13 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
gammaParams[3] = 1. / slpos;
|
||||
gammaParams[5] = 0.0;
|
||||
gammaParams[6] = 0.0;
|
||||
// printf("ga0=%f ga1=%f ga2=%f ga3=%f ga4=%f\n", ga0, ga1, ga2, ga3, ga4);
|
||||
// printf("ga0=%f ga1=%f ga2=%f ga3=%f ga4=%f\n", ga0, ga1, ga2, ga3, ga4);
|
||||
|
||||
// 7 parameters for smoother curves
|
||||
cmsCIExyY xyD;
|
||||
|
||||
Glib::ustring ills = "D50";
|
||||
|
||||
switch (ColorManagementParams::Illuminant(illum)) {
|
||||
case ColorManagementParams::Illuminant::DEFAULT:
|
||||
case ColorManagementParams::Illuminant::STDA:
|
||||
@ -802,55 +870,95 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
cmsWhitePointFromTemp(&xyD, tempv4);
|
||||
|
||||
switch (ColorManagementParams::Illuminant(illum)) {
|
||||
case ColorManagementParams::Illuminant::DEFAULT:
|
||||
case ColorManagementParams::Illuminant::D55:
|
||||
case ColorManagementParams::Illuminant::DEFAULT: {
|
||||
break;
|
||||
}
|
||||
|
||||
case ColorManagementParams::Illuminant::D55: {
|
||||
Wx = 0.956565934;
|
||||
Wz = 0.920253249;
|
||||
break;
|
||||
}
|
||||
|
||||
case ColorManagementParams::Illuminant::D80: {
|
||||
Wx = 0.950095542;
|
||||
Wz = 1.284213976;
|
||||
break;
|
||||
}
|
||||
|
||||
case ColorManagementParams::Illuminant::D41: {
|
||||
Wx = 0.991488263;
|
||||
Wz = 0.631604625;
|
||||
break;
|
||||
}
|
||||
|
||||
case ColorManagementParams::Illuminant::D50: {
|
||||
xyD = {0.3457, 0.3585, 1.0}; // near LCMS values but not perfect... it's a compromise!!
|
||||
Wx = 0.964295676;
|
||||
Wz = 0.825104603;
|
||||
break;
|
||||
}
|
||||
|
||||
case ColorManagementParams::Illuminant::D60: {
|
||||
Wx = 0.952646075;
|
||||
Wz = 1.008825184;
|
||||
xyD = {0.32168, 0.33767, 1.0};
|
||||
break;
|
||||
}
|
||||
|
||||
case ColorManagementParams::Illuminant::D65: {
|
||||
Wx = 0.95045471;
|
||||
Wz = 1.08905029;
|
||||
xyD = {0.312700492, 0.329000939, 1.0};
|
||||
break;
|
||||
}
|
||||
|
||||
case ColorManagementParams::Illuminant::D120: {
|
||||
Wx = 0.979182;
|
||||
Wz = 1.623623;
|
||||
xyD = {0.269669, 0.28078, 1.0};
|
||||
break;
|
||||
}
|
||||
|
||||
case ColorManagementParams::Illuminant::STDA: {
|
||||
Wx = 1.098500393;
|
||||
Wz = 0.355848714;
|
||||
xyD = {0.447573, 0.407440, 1.0};
|
||||
ills = "stdA 2875K";
|
||||
break;
|
||||
}
|
||||
|
||||
case ColorManagementParams::Illuminant::TUNGSTEN_2000K: {
|
||||
Wx = 1.274335;
|
||||
Wz = 0.145233;
|
||||
xyD = {0.526591, 0.41331, 1.0};
|
||||
ills = "Tungsten 2000K";
|
||||
break;
|
||||
}
|
||||
|
||||
case ColorManagementParams::Illuminant::TUNGSTEN_1500K: {
|
||||
Wx = 1.489921;
|
||||
Wz = 0.053826;
|
||||
xyD = {0.585703, 0.393157, 1.0};
|
||||
ills = "Tungsten 1500K";
|
||||
break;
|
||||
}
|
||||
}
|
||||
|
||||
double wprofpri[9];
|
||||
|
||||
if (gamutcontrol) {
|
||||
//xyz in functiuon primaries and illuminant
|
||||
Color::primaries_to_xyz(p, Wx, Wz, wprofpri);
|
||||
|
||||
for (int i = 0; i < 3; ++i) {
|
||||
for (int j = 0; j < 3; ++j) {
|
||||
wprofprim[i][j] = (double) wprofpri[j * 3 + i];
|
||||
//xyz in TMatrix format
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
//D41 0.377984 0.381229
|
||||
//D55 0.332424 0.347426
|
||||
//D80 0.293755 0.309185
|
||||
@ -869,7 +977,7 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
cmsWriteTag(oprofdef, cmsSigGreenTRCTag, GammaTRC[1]);
|
||||
cmsWriteTag(oprofdef, cmsSigBlueTRCTag, GammaTRC[2]);
|
||||
|
||||
//to read XYZ values and illuminant
|
||||
//to read XYZ values and illuminant
|
||||
if (rtengine::settings->verbose) {
|
||||
cmsCIEXYZ *redT = static_cast<cmsCIEXYZ*>(cmsReadTag(oprofdef, cmsSigRedMatrixColumnTag));
|
||||
cmsCIEXYZ *greenT = static_cast<cmsCIEXYZ*>(cmsReadTag(oprofdef, cmsSigGreenMatrixColumnTag));
|
||||
@ -881,6 +989,7 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
}
|
||||
|
||||
cmsFreeToneCurve(GammaTRC[0]);
|
||||
|
||||
if (oprofdef) {
|
||||
constexpr cmsUInt32Number flags = cmsFLAGS_NOOPTIMIZE | cmsFLAGS_NOCACHE | cmsFLAGS_BLACKPOINTCOMPENSATION | cmsFLAGS_GAMUTCHECK;
|
||||
const cmsHPROFILE iprof = ICCStore::getInstance()->getXYZProfile();
|
||||
@ -889,7 +998,10 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
lcmsMutex->unlock();
|
||||
}
|
||||
}
|
||||
|
||||
if (hTransform) {
|
||||
|
||||
|
||||
#ifdef _OPENMP
|
||||
#pragma omp parallel if (multiThread)
|
||||
#endif
|
||||
@ -901,19 +1013,30 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
#pragma omp for schedule(dynamic, 16) nowait
|
||||
#endif
|
||||
|
||||
for (int i = 0; i < ch; ++i) {
|
||||
for (int i = 0; i < ch; ++i)
|
||||
{
|
||||
float *p = pBuf.data;
|
||||
|
||||
for (int j = 0; j < cw; ++j) {
|
||||
const float r = src->r(i, j);
|
||||
const float g = src->g(i, j);
|
||||
const float b = src->b(i, j);
|
||||
float X = toxyz[0][0] * r + toxyz[0][1] * g + toxyz[0][2] * b;
|
||||
float Y = toxyz[1][0] * r + toxyz[1][1] * g + toxyz[1][2] * b;
|
||||
float Z = toxyz[2][0] * r + toxyz[2][1] * g + toxyz[2][2] * b;
|
||||
|
||||
*(p++) = toxyz[0][0] * r + toxyz[0][1] * g + toxyz[0][2] * b;
|
||||
*(p++) = toxyz[1][0] * r + toxyz[1][1] * g + toxyz[1][2] * b;
|
||||
*(p++) = toxyz[2][0] * r + toxyz[2][1] * g + toxyz[2][2] * b;
|
||||
if (gamutcontrol) {
|
||||
Color::gamutmap(X, Y, Z, wprofprim);//gamut control
|
||||
}
|
||||
|
||||
*(p++) = X;
|
||||
*(p++) = Y;
|
||||
*(p++) = Z;
|
||||
}
|
||||
|
||||
p = pBuf.data;
|
||||
cmsDoTransform(hTransform, p, p, cw);
|
||||
|
||||
for (int j = 0; j < cw; ++j) {
|
||||
dst->r(i, j) = *(p++) * normalize;
|
||||
dst->g(i, j) = *(p++) * normalize;
|
||||
@ -921,10 +1044,12 @@ void ImProcFunctions::workingtrc(const Imagefloat* src, Imagefloat* dst, int cw,
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
if (!keepTransForm) {
|
||||
cmsDeleteTransform(hTransform);
|
||||
hTransform = nullptr;
|
||||
}
|
||||
|
||||
transform = hTransform;
|
||||
}
|
||||
}
|
||||
|
@ -2276,6 +2276,7 @@ ColorManagementParams::ColorManagementParams() :
|
||||
bluy(0.0001),
|
||||
preser(0.),
|
||||
fbw(false),
|
||||
gamut(false),
|
||||
labgridcieALow(0.51763),//Prophoto red = (0.7347+0.1) * 1.81818 - 1
|
||||
labgridcieBLow(-0.33582),
|
||||
labgridcieAHigh(-0.75163),//Prophoto blue
|
||||
@ -2322,6 +2323,7 @@ bool ColorManagementParams::operator ==(const ColorManagementParams& other) cons
|
||||
&& labgridcieWy == other.labgridcieWy
|
||||
&& preser == other.preser
|
||||
&& fbw == other.fbw
|
||||
&& gamut == other.gamut
|
||||
&& aRendIntent == other.aRendIntent
|
||||
&& outputProfile == other.outputProfile
|
||||
&& outputIntent == other.outputIntent
|
||||
@ -7190,6 +7192,7 @@ int ProcParams::save(const Glib::ustring& fname, const Glib::ustring& fname2, bo
|
||||
saveToKeyfile(!pedited || pedited->icm.labgridcieWy, "Color Management", "LabGridcieWy", icm.labgridcieWy, keyFile);
|
||||
saveToKeyfile(!pedited || pedited->icm.preser, "Color Management", "Preser", icm.preser, keyFile);
|
||||
saveToKeyfile(!pedited || pedited->icm.fbw, "Color Management", "Fbw", icm.fbw, keyFile);
|
||||
saveToKeyfile(!pedited || pedited->icm.gamut, "Color Management", "Gamut", icm.gamut, keyFile);
|
||||
saveToKeyfile(!pedited || pedited->icm.outputProfile, "Color Management", "OutputProfile", icm.outputProfile, keyFile);
|
||||
saveToKeyfile(
|
||||
!pedited || pedited->icm.aRendIntent,
|
||||
@ -9455,6 +9458,7 @@ int ProcParams::load(const Glib::ustring& fname, ParamsEdited* pedited)
|
||||
assignFromKeyfile(keyFile, "Color Management", "Bluy", pedited, icm.bluy, pedited->icm.bluy);
|
||||
assignFromKeyfile(keyFile, "Color Management", "Preser", pedited, icm.preser, pedited->icm.preser);
|
||||
assignFromKeyfile(keyFile, "Color Management", "Fbw", pedited, icm.fbw, pedited->icm.fbw);
|
||||
assignFromKeyfile(keyFile, "Color Management", "Gamut", pedited, icm.gamut, pedited->icm.gamut);
|
||||
assignFromKeyfile(keyFile, "Color Management", "LabGridcieALow", pedited, icm.labgridcieALow, pedited->icm.labgridcieALow);
|
||||
assignFromKeyfile(keyFile, "Color Management", "LabGridcieBLow", pedited, icm.labgridcieBLow, pedited->icm.labgridcieBLow);
|
||||
assignFromKeyfile(keyFile, "Color Management", "LabGridcieAHigh", pedited, icm.labgridcieAHigh, pedited->icm.labgridcieAHigh);
|
||||
|
@ -1936,6 +1936,7 @@ struct ColorManagementParams {
|
||||
double bluy;
|
||||
double preser;
|
||||
bool fbw;
|
||||
bool gamut;
|
||||
double labgridcieALow;
|
||||
double labgridcieBLow;
|
||||
double labgridcieAHigh;
|
||||
|
@ -62,6 +62,7 @@ ICMPanel::ICMPanel() : FoldableToolPanel(this, "icm", M("TP_ICM_LABEL")), iuncha
|
||||
EvICMpreser = m->newEvent(LUMINANCECURVE, "HISTORY_MSG_ICM_PRESER");
|
||||
EvICMLabGridciexy = m->newEvent(LUMINANCECURVE, "HISTORY_MSG_ICL_LABGRIDCIEXY");
|
||||
EvICMfbw = m->newEvent(LUMINANCECURVE, "HISTORY_MSG_ICM_FBW");
|
||||
EvICMgamut = m->newEvent(LUMINANCECURVE, "HISTORY_MSG_ICM_GAMUT");
|
||||
isBatchMode = lastToneCurve = lastApplyLookTable = lastApplyBaselineExposureOffset = lastApplyHueSatMap = false;
|
||||
|
||||
ipDialog = Gtk::manage(new MyFileChooserButton(M("TP_ICM_INPUTDLGLABEL"), Gtk::FILE_CHOOSER_ACTION_OPEN));
|
||||
@ -263,8 +264,12 @@ ICMPanel::ICMPanel() : FoldableToolPanel(this, "icm", M("TP_ICM_LABEL")), iuncha
|
||||
wprimBox->pack_start(*wprim, Gtk::PACK_EXPAND_WIDGET);
|
||||
fbw = Gtk::manage(new Gtk::CheckButton((M("TP_ICM_FBW"))));
|
||||
fbw->set_active(true);
|
||||
gamut = Gtk::manage(new Gtk::CheckButton((M("TP_ICM_GAMUT"))));
|
||||
gamut->set_active(false);
|
||||
|
||||
trcProfVBox->pack_start(*wprimBox, Gtk::PACK_EXPAND_WIDGET);
|
||||
trcProfVBox->pack_start(*fbw, Gtk::PACK_EXPAND_WIDGET);
|
||||
trcProfVBox->pack_start(*gamut, Gtk::PACK_EXPAND_WIDGET);
|
||||
|
||||
neutral = Gtk::manage (new Gtk::Button (M ("TP_ICM_NEUTRAL")));
|
||||
setExpandAlignProperties (neutral, true, false, Gtk::ALIGN_FILL, Gtk::ALIGN_START);
|
||||
@ -468,6 +473,7 @@ ICMPanel::ICMPanel() : FoldableToolPanel(this, "icm", M("TP_ICM_LABEL")), iuncha
|
||||
wprimconn = wprim->signal_changed().connect(sigc::mem_fun(*this, &ICMPanel::wprimChanged));
|
||||
|
||||
fbwconn = fbw->signal_toggled().connect(sigc::mem_fun(*this, &ICMPanel::fbwChanged));
|
||||
gamutconn = gamut->signal_toggled().connect(sigc::mem_fun(*this, &ICMPanel::gamutChanged));
|
||||
obpcconn = obpc->signal_toggled().connect(sigc::mem_fun(*this, &ICMPanel::oBPCChanged));
|
||||
tcurveconn = ckbToneCurve->signal_toggled().connect(sigc::mem_fun(*this, &ICMPanel::toneCurveChanged));
|
||||
ltableconn = ckbApplyLookTable->signal_toggled().connect(sigc::mem_fun(*this, &ICMPanel::applyLookTableChanged));
|
||||
@ -513,6 +519,7 @@ void ICMPanel::neutral_pressed ()
|
||||
wSlope->setValue(defPar.workingTRCSlope);//12.92
|
||||
preser->setValue(defPar.preser);
|
||||
fbw->set_active(defPar.fbw);
|
||||
gamut->set_active(defPar.gamut);
|
||||
wTRC->set_active(toUnderlying(ColorManagementParams::WorkingTrc::NONE));//reset to none
|
||||
will->set_active(toUnderlying(ColorManagementParams::Illuminant::DEFAULT));//reset to default - after wprim
|
||||
}
|
||||
@ -765,6 +772,7 @@ void ICMPanel::read(const ProcParams* pp, const ParamsEdited* pedited)
|
||||
|
||||
ConnectionBlocker obpcconn_(obpcconn);
|
||||
ConnectionBlocker fbwconn_(fbwconn);
|
||||
ConnectionBlocker gamutconn_(gamutconn);
|
||||
ConnectionBlocker ipc_(ipc);
|
||||
ConnectionBlocker tcurveconn_(tcurveconn);
|
||||
ConnectionBlocker ltableconn_(ltableconn);
|
||||
@ -838,6 +846,7 @@ void ICMPanel::read(const ProcParams* pp, const ParamsEdited* pedited)
|
||||
|
||||
obpc->set_active(pp->icm.outputBPC);
|
||||
fbw->set_active(pp->icm.fbw);
|
||||
gamut->set_active(pp->icm.gamut);
|
||||
ckbToneCurve->set_active(pp->icm.toneCurve);
|
||||
lastToneCurve = pp->icm.toneCurve;
|
||||
ckbApplyLookTable->set_active(pp->icm.applyLookTable);
|
||||
@ -862,6 +871,7 @@ void ICMPanel::read(const ProcParams* pp, const ParamsEdited* pedited)
|
||||
iunchanged->set_active(!pedited->icm.inputProfile);
|
||||
obpc->set_inconsistent(!pedited->icm.outputBPC);
|
||||
fbw->set_inconsistent(!pedited->icm.fbw);
|
||||
gamut->set_inconsistent(!pedited->icm.gamut);
|
||||
ckbToneCurve->set_inconsistent(!pedited->icm.toneCurve);
|
||||
ckbApplyLookTable->set_inconsistent(!pedited->icm.applyLookTable);
|
||||
ckbApplyBaselineExposureOffset->set_inconsistent(!pedited->icm.applyBaselineExposureOffset);
|
||||
@ -920,6 +930,7 @@ void ICMPanel::read(const ProcParams* pp, const ParamsEdited* pedited)
|
||||
willulab->set_sensitive(false);
|
||||
wprim->set_sensitive(false);
|
||||
fbw->set_sensitive(false);
|
||||
gamut->set_sensitive(false);
|
||||
wprimlab->set_sensitive(false);
|
||||
riaHBox->set_sensitive(false);
|
||||
redFrame->hide();
|
||||
@ -931,6 +942,7 @@ void ICMPanel::read(const ProcParams* pp, const ParamsEdited* pedited)
|
||||
willulab->set_sensitive(true);
|
||||
wprim->set_sensitive(true);
|
||||
fbw->set_sensitive(true);
|
||||
gamut->set_sensitive(true);
|
||||
wprimlab->set_sensitive(true);
|
||||
if (ColorManagementParams::Primaries(wprim->get_active_row_number()) == ColorManagementParams::Primaries::DEFAULT) {
|
||||
redFrame->hide();
|
||||
@ -973,6 +985,7 @@ void ICMPanel::read(const ProcParams* pp, const ParamsEdited* pedited)
|
||||
willulab->set_sensitive(true);
|
||||
wprim->set_sensitive(true);
|
||||
fbw->set_sensitive(true);
|
||||
gamut->set_sensitive(true);
|
||||
wprimlab->set_sensitive(true);
|
||||
wGamma->set_sensitive(false);
|
||||
wSlope->set_sensitive(false);
|
||||
@ -990,6 +1003,7 @@ void ICMPanel::read(const ProcParams* pp, const ParamsEdited* pedited)
|
||||
willulab->set_sensitive(true);
|
||||
wprim->set_sensitive(true);
|
||||
fbw->set_sensitive(true);
|
||||
gamut->set_sensitive(true);
|
||||
wprimlab->set_sensitive(true);
|
||||
wGamma->set_sensitive(false);
|
||||
wSlope->set_sensitive(false);
|
||||
@ -1007,6 +1021,7 @@ void ICMPanel::read(const ProcParams* pp, const ParamsEdited* pedited)
|
||||
willulab->set_sensitive(true);
|
||||
wprim->set_sensitive(true);
|
||||
fbw->set_sensitive(true);
|
||||
gamut->set_sensitive(true);
|
||||
wprimlab->set_sensitive(true);
|
||||
redFrame->show();
|
||||
wGamma->set_sensitive(false);
|
||||
@ -1025,6 +1040,7 @@ void ICMPanel::read(const ProcParams* pp, const ParamsEdited* pedited)
|
||||
willulab->set_sensitive(true);
|
||||
wprim->set_sensitive(true);
|
||||
fbw->set_sensitive(true);
|
||||
gamut->set_sensitive(true);
|
||||
wprimlab->set_sensitive(true);
|
||||
riaHBox->set_sensitive(true);
|
||||
if (ColorManagementParams::Primaries(wprim->get_active_row_number()) == ColorManagementParams::Primaries::DEFAULT) {
|
||||
@ -1042,6 +1058,7 @@ void ICMPanel::read(const ProcParams* pp, const ParamsEdited* pedited)
|
||||
willulab->set_sensitive(true);
|
||||
wprim->set_sensitive(true);
|
||||
fbw->set_sensitive(true);
|
||||
gamut->set_sensitive(true);
|
||||
wprimlab->set_sensitive(true);
|
||||
wGamma->set_sensitive(false);
|
||||
wSlope->set_sensitive(false);
|
||||
@ -1143,6 +1160,7 @@ void ICMPanel::write(ProcParams* pp, ParamsEdited* pedited)
|
||||
pp->icm.applyHueSatMap = ckbApplyHueSatMap->get_active();
|
||||
pp->icm.outputBPC = obpc->get_active();
|
||||
pp->icm.fbw = fbw->get_active();
|
||||
pp->icm.gamut = gamut->get_active();
|
||||
pp->icm.workingTRCGamma = wGamma->getValue();
|
||||
pp->icm.workingTRCSlope = wSlope->getValue();
|
||||
pp->icm.redx = redx->getValue();
|
||||
@ -1162,6 +1180,7 @@ void ICMPanel::write(ProcParams* pp, ParamsEdited* pedited)
|
||||
pedited->icm.aRendIntent = aRendIntent->getSelected() < 4;
|
||||
pedited->icm.outputBPC = !obpc->get_inconsistent();
|
||||
pedited->icm.fbw = !fbw->get_inconsistent();
|
||||
pedited->icm.gamut = !gamut->get_inconsistent();
|
||||
pedited->icm.dcpIlluminant = dcpIll->get_active_text() != M("GENERAL_UNCHANGED");
|
||||
pedited->icm.toneCurve = !ckbToneCurve->get_inconsistent();
|
||||
pedited->icm.applyLookTable = !ckbApplyLookTable->get_inconsistent();
|
||||
@ -1268,6 +1287,7 @@ void ICMPanel::wtrcinChanged()
|
||||
willulab->set_sensitive(false);
|
||||
wprim->set_sensitive(false);
|
||||
fbw->set_sensitive(false);
|
||||
gamut->set_sensitive(false);
|
||||
wprimlab->set_sensitive(false);
|
||||
redFrame->hide();
|
||||
riaHBox->set_sensitive(false);
|
||||
@ -1278,6 +1298,7 @@ void ICMPanel::wtrcinChanged()
|
||||
will->set_sensitive(false);
|
||||
wprim->set_sensitive(true);
|
||||
fbw->set_sensitive(true);
|
||||
gamut->set_sensitive(true);
|
||||
wprimlab->set_sensitive(true);
|
||||
willulab->set_sensitive(true);
|
||||
if (ColorManagementParams::Primaries(wprim->get_active_row_number()) == ColorManagementParams::Primaries::DEFAULT) {
|
||||
@ -1311,6 +1332,7 @@ void ICMPanel::wtrcinChanged()
|
||||
willulab->set_sensitive(true);
|
||||
wprim->set_sensitive(true);
|
||||
fbw->set_sensitive(true);
|
||||
gamut->set_sensitive(true);
|
||||
wprimlab->set_sensitive(true);
|
||||
wGamma->set_sensitive(false);
|
||||
wSlope->set_sensitive(false);
|
||||
@ -1336,6 +1358,7 @@ void ICMPanel::wtrcinChanged()
|
||||
willulab->set_sensitive(true);
|
||||
wprim->set_sensitive(true);
|
||||
fbw->set_sensitive(true);
|
||||
gamut->set_sensitive(true);
|
||||
wGamma->set_sensitive(false);
|
||||
wSlope->set_sensitive(false);
|
||||
riaHBox->set_sensitive(true);
|
||||
@ -1362,6 +1385,7 @@ void ICMPanel::wtrcinChanged()
|
||||
willulab->set_sensitive(true);
|
||||
wprim->set_sensitive(true);
|
||||
fbw->set_sensitive(true);
|
||||
gamut->set_sensitive(true);
|
||||
wprimlab->set_sensitive(true);
|
||||
wGamma->set_sensitive(false);
|
||||
wSlope->set_sensitive(false);
|
||||
@ -1389,6 +1413,7 @@ void ICMPanel::wtrcinChanged()
|
||||
willulab->set_sensitive(true);
|
||||
wprim->set_sensitive(true);
|
||||
fbw->set_sensitive(true);
|
||||
gamut->set_sensitive(true);
|
||||
wprimlab->set_sensitive(true);
|
||||
wGamma->set_sensitive(false);
|
||||
wSlope->set_sensitive(false);
|
||||
@ -1416,6 +1441,7 @@ void ICMPanel::wtrcinChanged()
|
||||
willulab->set_sensitive(true);
|
||||
wprim->set_sensitive(true);
|
||||
fbw->set_sensitive(true);
|
||||
gamut->set_sensitive(true);
|
||||
wprimlab->set_sensitive(true);
|
||||
wGamma->set_sensitive(false);
|
||||
wSlope->set_sensitive(false);
|
||||
@ -2023,6 +2049,33 @@ void ICMPanel::fbwChanged()
|
||||
}
|
||||
}
|
||||
|
||||
void ICMPanel::gamutChanged()
|
||||
{
|
||||
if (multiImage) {
|
||||
if (gamut->get_inconsistent()) {
|
||||
gamut->set_inconsistent(false);
|
||||
gamutconn.block(true);
|
||||
gamut->set_active(false);
|
||||
gamutconn.block(false);
|
||||
} else if (lastgamut) {
|
||||
gamut->set_inconsistent(true);
|
||||
}
|
||||
|
||||
lastgamut = gamut->get_active();
|
||||
}
|
||||
|
||||
if (listener) {
|
||||
if (gamut->get_inconsistent()) {
|
||||
listener->panelChanged(EvICMgamut, M("GENERAL_UNCHANGED"));
|
||||
} else if (fbw->get_active()) {
|
||||
listener->panelChanged(EvICMgamut, M("GENERAL_ENABLED"));
|
||||
} else {
|
||||
listener->panelChanged(EvICMgamut, M("GENERAL_DISABLED"));
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
void ICMPanel::setRawMeta(bool raw, const rtengine::FramesData* pMeta)
|
||||
{
|
||||
|
||||
|
@ -81,6 +81,8 @@ protected:
|
||||
bool lastfbw;
|
||||
sigc::connection fbwconn;
|
||||
bool isBatchMode;
|
||||
bool lastgamut;
|
||||
sigc::connection gamutconn;
|
||||
|
||||
private:
|
||||
rtengine::ProcEvent EvICMprimariMethod;
|
||||
@ -107,6 +109,7 @@ private:
|
||||
rtengine::ProcEvent EvICMpreser;
|
||||
rtengine::ProcEvent EvICMLabGridciexy;
|
||||
rtengine::ProcEvent EvICMfbw;
|
||||
rtengine::ProcEvent EvICMgamut;
|
||||
LabGrid *labgridcie;
|
||||
IdleRegister idle_register;
|
||||
|
||||
@ -121,6 +124,7 @@ private:
|
||||
Gtk::Box* iVBox;
|
||||
Gtk::Box* wTRCBox;
|
||||
Gtk::CheckButton* fbw;
|
||||
Gtk::CheckButton* gamut;
|
||||
|
||||
Gtk::CheckButton* obpc;
|
||||
Gtk::RadioButton* inone;
|
||||
@ -196,6 +200,7 @@ public:
|
||||
void aiChanged(int n);
|
||||
void oBPCChanged();
|
||||
void fbwChanged();
|
||||
void gamutChanged();
|
||||
void ipChanged();
|
||||
void ipSelectionChanged();
|
||||
void dcpIlluminantChanged();
|
||||
|
@ -452,6 +452,7 @@ void ParamsEdited::set(bool v)
|
||||
icm.bluy = v;
|
||||
icm.preser = v;
|
||||
icm.fbw = v;
|
||||
icm.gamut = v;
|
||||
icm.labgridcieALow = v;
|
||||
icm.labgridcieBLow = v;
|
||||
icm.labgridcieAHigh = v;
|
||||
@ -1871,6 +1872,7 @@ void ParamsEdited::initFrom(const std::vector<rtengine::procparams::ProcParams>&
|
||||
icm.labgridcieWy = icm.labgridcieWy && p.icm.labgridcieWy == other.icm.labgridcieWy;
|
||||
icm.preser = icm.preser && p.icm.preser == other.icm.preser;
|
||||
icm.fbw = icm.fbw && p.icm.fbw == other.icm.fbw;
|
||||
icm.gamut = icm.gamut && p.icm.gamut == other.icm.gamut;
|
||||
icm.aRendIntent = icm.aRendIntent && p.icm.aRendIntent == other.icm.aRendIntent;
|
||||
icm.workingTRC = icm.workingTRC && p.icm.workingTRC == other.icm.workingTRC;
|
||||
icm.will = icm.will && p.icm.will == other.icm.will;
|
||||
@ -6377,6 +6379,10 @@ void ParamsEdited::combine(rtengine::procparams::ProcParams& toEdit, const rteng
|
||||
toEdit.icm.fbw = mods.icm.fbw;
|
||||
}
|
||||
|
||||
if (icm.gamut) {
|
||||
toEdit.icm.gamut = mods.icm.gamut;
|
||||
}
|
||||
|
||||
if (icm.labgridcieALow) {
|
||||
toEdit.icm.labgridcieALow = mods.icm.labgridcieALow;
|
||||
}
|
||||
|
@ -1230,6 +1230,7 @@ struct ColorManagementParamsEdited {
|
||||
bool bluy;
|
||||
bool preser;
|
||||
bool fbw;
|
||||
bool gamut;
|
||||
bool labgridcieALow;
|
||||
bool labgridcieBLow;
|
||||
bool labgridcieAHigh;
|
||||
|
Loading…
x
Reference in New Issue
Block a user