1) quantization stepsizes stored as float instead of shifted integers -> fixes a pb of precision when using very small stepsizes. 2) bug fixed when decoding until bitplane 0 -> r-value (1/2) added to the coefficient.
This commit is contained in:
parent
7f8f47566f
commit
7ee36c3a4c
|
@ -250,12 +250,13 @@ int floorlog2(int a)
|
||||||
return l;
|
return l;
|
||||||
}
|
}
|
||||||
|
|
||||||
void encode_stepsize(int stepsize, int numbps, int *expn, int *mant)
|
void encode_stepsize(float stepsize, int numbps, int *expn, int *mant)
|
||||||
{
|
{
|
||||||
int p, n;
|
int p, n, stepsizeTMP;
|
||||||
p = floorlog2(stepsize) - 13;
|
stepsizeTMP=(int) floor(stepsize * 8192.0);
|
||||||
n = 11 - floorlog2(stepsize);
|
p = floorlog2(stepsizeTMP) - 13;
|
||||||
*mant = (n < 0 ? stepsize >> -n : stepsize << n) & 0x7ff;
|
n = 11 - floorlog2(stepsizeTMP);
|
||||||
|
*mant = (n < 0 ? stepsizeTMP >> -n : stepsizeTMP << n) & 0x7ff;
|
||||||
*expn = numbps - p;
|
*expn = numbps - p;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
@ -264,7 +265,7 @@ void calc_explicit_stepsizes(j2k_tccp_t * tccp, int prec)
|
||||||
int numbands, bandno;
|
int numbands, bandno;
|
||||||
numbands = 3 * tccp->numresolutions - 2;
|
numbands = 3 * tccp->numresolutions - 2;
|
||||||
for (bandno = 0; bandno < numbands; bandno++) {
|
for (bandno = 0; bandno < numbands; bandno++) {
|
||||||
double stepsize;
|
float stepsize;
|
||||||
|
|
||||||
int resno, level, orient, gain;
|
int resno, level, orient, gain;
|
||||||
resno = bandno == 0 ? 0 : (bandno - 1) / 3 + 1;
|
resno = bandno == 0 ? 0 : (bandno - 1) / 3 + 1;
|
||||||
|
@ -280,7 +281,7 @@ void calc_explicit_stepsizes(j2k_tccp_t * tccp, int prec)
|
||||||
double norm = dwt_norms_97[orient][level];
|
double norm = dwt_norms_97[orient][level];
|
||||||
stepsize = (1 << (gain + 1)) / norm;
|
stepsize = (1 << (gain + 1)) / norm;
|
||||||
}
|
}
|
||||||
encode_stepsize((int) floor(stepsize * 8192.0), prec + gain,
|
encode_stepsize(stepsize, prec + gain,
|
||||||
&tccp->stepsizes[bandno].expn,
|
&tccp->stepsizes[bandno].expn,
|
||||||
&tccp->stepsizes[bandno].mant);
|
&tccp->stepsizes[bandno].mant);
|
||||||
}
|
}
|
||||||
|
|
|
@ -535,7 +535,7 @@ void t1_dec_clnpass(int w, int h, int bpno, int orient, int cblksty)
|
||||||
}
|
}
|
||||||
} /* VSC and BYPASS by Antonin */
|
} /* VSC and BYPASS by Antonin */
|
||||||
|
|
||||||
double t1_getwmsedec(int nmsedec, int compno, int level, int orient, int bpno, int qmfbid, double stepsize, int numcomps) //mod fixed_quality
|
double t1_getwmsedec(int nmsedec, int compno, int level, int orient, int bpno, int qmfbid, float stepsize, int numcomps) //mod fixed_quality
|
||||||
{
|
{
|
||||||
double w1, w2, wmsedec;
|
double w1, w2, wmsedec;
|
||||||
if (qmfbid == 1) {
|
if (qmfbid == 1) {
|
||||||
|
@ -550,7 +550,7 @@ double t1_getwmsedec(int nmsedec, int compno, int level, int orient, int bpno, i
|
||||||
return wmsedec;
|
return wmsedec;
|
||||||
}
|
}
|
||||||
|
|
||||||
void t1_encode_cblk(tcd_cblk_t * cblk, int orient, int compno, int level, int qmfbid, double stepsize, int cblksty, int numcomps, tcd_tile_t * tile) //mod fixed_quality
|
void t1_encode_cblk(tcd_cblk_t * cblk, int orient, int compno, int level, int qmfbid, float stepsize, int cblksty, int numcomps, tcd_tile_t * tile) //mod fixed_quality
|
||||||
{
|
{
|
||||||
int i, j;
|
int i, j;
|
||||||
int w, h;
|
int w, h;
|
||||||
|
@ -723,19 +723,17 @@ void t1_decode_cblk(tcd_cblk_t * cblk, int orient, int roishift,
|
||||||
else
|
else
|
||||||
mqc_init_dec(seg->data, seg->len);
|
mqc_init_dec(seg->data, seg->len);
|
||||||
// ddA
|
// ddA
|
||||||
|
|
||||||
if (bpno==0) cblk->lastbp=1; // Add Antonin : quantizbug1
|
|
||||||
|
|
||||||
for (passno = 0; passno < seg->numpasses; passno++) {
|
for (passno = 0; passno < seg->numpasses; passno++) {
|
||||||
switch (passtype) {
|
switch (passtype) {
|
||||||
case 0:
|
case 0:
|
||||||
t1_dec_sigpass(w, h, bpno, orient, type, cblksty);
|
t1_dec_sigpass(w, h, bpno+1, orient, type, cblksty);
|
||||||
break;
|
break;
|
||||||
case 1:
|
case 1:
|
||||||
t1_dec_refpass(w, h, bpno, type, cblksty);
|
t1_dec_refpass(w, h, bpno+1, type, cblksty);
|
||||||
break;
|
break;
|
||||||
case 2:
|
case 2:
|
||||||
t1_dec_clnpass(w, h, bpno, orient, cblksty);
|
t1_dec_clnpass(w, h, bpno+1, orient, cblksty);
|
||||||
break;
|
break;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
@ -819,7 +817,7 @@ void t1_encode_cblks(tcd_tile_t * tile, j2k_tcp_t * tcp)
|
||||||
tilec->
|
tilec->
|
||||||
x0)],
|
x0)],
|
||||||
8192 * 8192 /
|
8192 * 8192 /
|
||||||
band->stepsize) >> (13 - T1_NMSEDEC_FRACBITS);
|
((int) floor(band->stepsize * 8192))) >> (13 - T1_NMSEDEC_FRACBITS);
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
@ -896,34 +894,22 @@ void t1_decode_cblks(tcd_tile_t * tile, j2k_tcp_t * tcp)
|
||||||
if (tcp->tccps[compno].qmfbid == 1) {
|
if (tcp->tccps[compno].qmfbid == 1) {
|
||||||
for (j = 0; j < cblk->y1 - cblk->y0; j++) {
|
for (j = 0; j < cblk->y1 - cblk->y0; j++) {
|
||||||
for (i = 0; i < cblk->x1 - cblk->x0; i++) {
|
for (i = 0; i < cblk->x1 - cblk->x0; i++) {
|
||||||
tilec->data[x + i +
|
int tmp=t1_data[j][i];
|
||||||
(y + j) * (tilec->x1 -
|
if (tmp>>1==0) tilec->data[x + i + (y + j) * (tilec->x1 - tilec->x0)] = 0;
|
||||||
tilec->x0)] = t1_data[j][i];
|
else tilec->data[x + i + (y + j) * (tilec->x1 - tilec->x0)] = tmp<0?((tmp>>1) | 0x80000000)+1:(tmp>>1);
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
} else { /* if (tcp->tccps[compno].qmfbid == 0) */
|
} else { /* if (tcp->tccps[compno].qmfbid == 0) */
|
||||||
|
|
||||||
for (j = 0; j < cblk->y1 - cblk->y0; j++) {
|
for (j = 0; j < cblk->y1 - cblk->y0; j++) {
|
||||||
for (i = 0; i < cblk->x1 - cblk->x0; i++) {
|
for (i = 0; i < cblk->x1 - cblk->x0; i++) {
|
||||||
if (t1_data[j][i] == 0) {
|
float tmp=t1_data[j][i] * band->stepsize * 4096.0;
|
||||||
tilec->data[x + i +
|
int tmp2;
|
||||||
(y + j) * (tilec->x1 - tilec->x0)] = 0;
|
if (t1_data[j][i]>>1 == 0) {
|
||||||
|
tilec->data[x + i + (y + j) * (tilec->x1 - tilec->x0)] = 0;
|
||||||
} else {
|
} else {
|
||||||
|
tmp2=((int) (floor(fabs(tmp)))) + ((int) floor(fabs(tmp*2))%2);
|
||||||
// Add antonin : quantizbug1
|
tilec->data[x + i + (y + j) * (tilec->x1 - tilec->x0)] = ((tmp<0)?-tmp2:tmp2);
|
||||||
|
|
||||||
t1_data[j][i]<<=1;
|
|
||||||
|
|
||||||
//if (cblk->lastbp)
|
|
||||||
|
|
||||||
t1_data[j][i]+=t1_data[j][i]>0?1:-1;
|
|
||||||
|
|
||||||
// ddA
|
|
||||||
tilec->data[x + i +
|
|
||||||
(y + j) * (tilec->x1 -
|
|
||||||
tilec->
|
|
||||||
x0)] =
|
|
||||||
fix_mul(t1_data[j][i] << 12, band->stepsize); //Mod Antonin : quantizbug1 (before : << 13)
|
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
|
@ -87,7 +87,7 @@ void tcd_dump(tcd_image_t * img, int curtileno)
|
||||||
tcd_band_t *band = &res->bands[bandno];
|
tcd_band_t *band = &res->bands[bandno];
|
||||||
fprintf(stdout, " band {\n");
|
fprintf(stdout, " band {\n");
|
||||||
fprintf(stdout,
|
fprintf(stdout,
|
||||||
" x0=%d, y0=%d, x1=%d, y1=%d, stepsize=%d, numbps=%d\n",
|
" x0=%d, y0=%d, x1=%d, y1=%d, stepsize=%f, numbps=%d\n",
|
||||||
band->x0, band->y0, band->x1, band->y1,
|
band->x0, band->y0, band->x1, band->y1,
|
||||||
band->stepsize, band->numbps);
|
band->stepsize, band->numbps);
|
||||||
for (precno = 0; precno < res->pw * res->ph; precno++) {
|
for (precno = 0; precno < res->pw * res->ph; precno++) {
|
||||||
|
@ -270,9 +270,7 @@ void tcd_malloc_encode(j2k_image_t * img, j2k_cp_t * cp, int curtileno)
|
||||||
tccp->qmfbid ==
|
tccp->qmfbid ==
|
||||||
0 ? dwt_getgain_real(band->bandno) : dwt_getgain(band->bandno);
|
0 ? dwt_getgain_real(band->bandno) : dwt_getgain(band->bandno);
|
||||||
numbps = img->comps[compno].prec + gain;
|
numbps = img->comps[compno].prec + gain;
|
||||||
band->stepsize =
|
band->stepsize = (1.0 + ss->mant / 2048.0) * pow(2.0, numbps - ss->expn);
|
||||||
(int) floor((1.0 + ss->mant / 2048.0) *
|
|
||||||
pow(2.0, numbps - ss->expn) * 8192.0);
|
|
||||||
band->numbps = ss->expn + tccp->numgbits - 1; /* WHY -1 ? */
|
band->numbps = ss->expn + tccp->numgbits - 1; /* WHY -1 ? */
|
||||||
|
|
||||||
band->precincts =
|
band->precincts =
|
||||||
|
@ -514,9 +512,7 @@ void tcd_init_encode(j2k_image_t * img, j2k_cp_t * cp, int curtileno)
|
||||||
tccp->qmfbid ==
|
tccp->qmfbid ==
|
||||||
0 ? dwt_getgain_real(band->bandno) : dwt_getgain(band->bandno);
|
0 ? dwt_getgain_real(band->bandno) : dwt_getgain(band->bandno);
|
||||||
numbps = img->comps[compno].prec + gain;
|
numbps = img->comps[compno].prec + gain;
|
||||||
band->stepsize =
|
band->stepsize = (1.0 + ss->mant / 2048.0) * pow(2.0, numbps - ss->expn);
|
||||||
(int) floor((1.0 + ss->mant / 2048.0) *
|
|
||||||
pow(2.0, numbps - ss->expn) * 8192.0);
|
|
||||||
band->numbps = ss->expn + tccp->numgbits - 1; /* WHY -1 ? */
|
band->numbps = ss->expn + tccp->numgbits - 1; /* WHY -1 ? */
|
||||||
|
|
||||||
for (precno = 0; precno < res->pw * res->ph; precno++) {
|
for (precno = 0; precno < res->pw * res->ph; precno++) {
|
||||||
|
@ -718,13 +714,11 @@ void tcd_init(j2k_image_t * img, j2k_cp_t * cp)
|
||||||
|
|
||||||
ss = &tccp->stepsizes[resno ==
|
ss = &tccp->stepsizes[resno ==
|
||||||
0 ? 0 : 3 * (resno - 1) + bandno + 1];
|
0 ? 0 : 3 * (resno - 1) + bandno + 1];
|
||||||
gain =
|
gain =
|
||||||
tccp->qmfbid ==
|
tccp->qmfbid ==
|
||||||
0 ? dwt_getgain_real(band->bandno) : dwt_getgain(band->bandno);
|
0 ? dwt_getgain_real(band->bandno) : dwt_getgain(band->bandno);
|
||||||
numbps = img->comps[compno].prec + gain;
|
numbps = img->comps[compno].prec + gain;
|
||||||
band->stepsize =
|
band->stepsize = (1.0 + ss->mant / 2048.0) * pow(2.0, numbps - ss->expn);
|
||||||
(int) floor((1.0 + ss->mant / 2048.0) *
|
|
||||||
pow(2.0, numbps - ss->expn) * 8192.0);
|
|
||||||
band->numbps = ss->expn + tccp->numgbits - 1; /* WHY -1 ? */
|
band->numbps = ss->expn + tccp->numgbits - 1; /* WHY -1 ? */
|
||||||
|
|
||||||
band->precincts =
|
band->precincts =
|
||||||
|
@ -777,10 +771,6 @@ void tcd_init(j2k_image_t * img, j2k_cp_t * cp)
|
||||||
cblk->y0 = int_max(cblkystart, prc->y0);
|
cblk->y0 = int_max(cblkystart, prc->y0);
|
||||||
cblk->x1 = int_min(cblkxend, prc->x1);
|
cblk->x1 = int_min(cblkxend, prc->x1);
|
||||||
cblk->y1 = int_min(cblkyend, prc->y1);
|
cblk->y1 = int_min(cblkyend, prc->y1);
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
cblk->lastbp = 0; // Add Antonin : quantizbug1
|
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
@ -1571,25 +1561,16 @@ int tcd_decode_tile(unsigned char *src, int len, int tileno)
|
||||||
for (i = res->x0; i < res->x1; i++) {
|
for (i = res->x0; i < res->x1; i++) {
|
||||||
|
|
||||||
int v;
|
int v;
|
||||||
|
float tmp = (tilec->data[i - res->x0 + (j - res->y0) * tw])/8192.0;
|
||||||
double tmp =
|
int tmp2;
|
||||||
(double) tilec->data[i - res->x0 + (j - res->y0) * tw];
|
|
||||||
if (tcd_tcp->tccps[compno].qmfbid == 1) {
|
if (tcd_tcp->tccps[compno].qmfbid == 1) {
|
||||||
v = (int) tmp;
|
v = tilec->data[i - res->x0 + (j - res->y0) * tw];
|
||||||
} else {
|
} else {
|
||||||
|
tmp2=((int) (floor(fabs(tmp)))) + ((int) floor(fabs(tmp*2))%2);
|
||||||
//v = (int) tmp >> 13;
|
v = ((tmp<0)?-tmp2:tmp2);
|
||||||
|
}
|
||||||
//Mod antonin : multbug1
|
|
||||||
v =
|
|
||||||
(int) ((fabs(tmp / 8192.0) >=
|
|
||||||
floor(fabs(tmp / 8192.0)) +
|
|
||||||
0.5) ? fabs(tmp / 8192.0) + 1.0 : fabs(tmp / 8192.0));
|
|
||||||
|
|
||||||
v = (tmp < 0) ? -v : v;
|
|
||||||
|
|
||||||
//doM
|
|
||||||
}
|
|
||||||
v += adjust;
|
v += adjust;
|
||||||
|
|
||||||
tcd_img->comps[compno].data[(i - offset_x) +
|
tcd_img->comps[compno].data[(i - offset_x) +
|
||||||
|
|
|
@ -57,7 +57,6 @@ typedef struct {
|
||||||
typedef struct {
|
typedef struct {
|
||||||
int x0, y0, x1, y1; /* dimension of the code-blocks : left upper corner (x0, y0) right low corner (x1,y1) */
|
int x0, y0, x1, y1; /* dimension of the code-blocks : left upper corner (x0, y0) right low corner (x1,y1) */
|
||||||
int numbps;
|
int numbps;
|
||||||
int lastbp; /* Add antonin : quantizbug1 */
|
|
||||||
int numlenbits;
|
int numlenbits;
|
||||||
int len; /* length */
|
int len; /* length */
|
||||||
int numpasses; /* number of pass already done for the code-blocks */
|
int numpasses; /* number of pass already done for the code-blocks */
|
||||||
|
@ -84,7 +83,7 @@ typedef struct {
|
||||||
int bandno;
|
int bandno;
|
||||||
tcd_precinct_t *precincts; /* precinct information */
|
tcd_precinct_t *precincts; /* precinct information */
|
||||||
int numbps;
|
int numbps;
|
||||||
int stepsize;
|
float stepsize;
|
||||||
} tcd_band_t;
|
} tcd_band_t;
|
||||||
|
|
||||||
typedef struct {
|
typedef struct {
|
||||||
|
|
Loading…
Reference in New Issue