00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021 #define CPL_TYPE_ADD(a) CONCAT2X(a, CPL_TYPE)
00022 #define CPL_TYPE_ADD_CONST(a) CONCAT3X(a, CPL_TYPE, const)
00023
00024
00025
00043
00044 static cpl_error_code
00045 CPL_TYPE_ADD(hawki_geom_img_offset_saa)(cpl_image * final,
00046 cpl_image * contrib,
00047 int rmin, int rmax,
00048 double start_x,
00049 double start_y,
00050 const cpl_imagelist * ilist,
00051 const cpl_bivector * offs,
00052 const cpl_vector * xyprofile,
00053 int tabsperpix)
00054
00055 {
00056
00057 const int nima = cpl_imagelist_get_size(ilist);
00058 const cpl_image * img1 = cpl_imagelist_get_const(ilist, 0);
00059 const int sizex = cpl_image_get_size_x(img1);
00060 const int sizey = cpl_image_get_size_y(img1);
00061 const int nx = cpl_image_get_size_x(final);
00062 const int ny = cpl_image_get_size_y(final);
00063 const int rtot = rmin + rmax;
00064 CPL_TYPE * po = CPL_TYPE_ADD(cpl_image_get_data)(final);
00065 const double * offs_x = cpl_bivector_get_x_data_const(offs);
00066 const double * offs_y = cpl_bivector_get_y_data_const(offs);
00067 const double * interp_kernel = cpl_vector_get_data_const(xyprofile);
00068 double * acc = (double*)cpl_malloc(nx * nima * sizeof(double));
00069 int * pcontrib = cpl_image_get_data_int(contrib);
00070 const CPL_TYPE ** ppi = cpl_malloc(nima * sizeof(CPL_TYPE*));
00071 const cpl_binary ** ppibpm = cpl_malloc(nima * sizeof(cpl_binary*));
00072 int * ncontrib = (int*)cpl_calloc(sizeof(int), nx);
00073 int * offset_ij = (int*)cpl_malloc(2 * nima * sizeof(int));
00074 double * rsc = (double*)cpl_malloc(nima * 8 * sizeof(double));
00075 int j, p;
00076 int bpmflag = 0;
00077
00078 for (p=0; p < nima; p++) {
00079
00080 const double offset_x = start_x - offs_x[p];
00081 const double offset_y = start_y - offs_y[p];
00082 const double suboff_x = offset_x - floor(offset_x);
00083 const double suboff_y = offset_y - floor(offset_y);
00084
00085
00086 const int tabx = (int)(0.5 + suboff_x * (double)(tabsperpix));
00087 const int taby = (int)(0.5 + suboff_y * (double)(tabsperpix));
00088
00089 const double sumrs
00090 = ( interp_kernel[ tabsperpix + tabx]
00091 + interp_kernel[ tabx]
00092 + interp_kernel[ tabsperpix - tabx]
00093 + interp_kernel[2 * tabsperpix - tabx])
00094 * ( interp_kernel[ tabsperpix + taby]
00095 + interp_kernel[ taby]
00096 + interp_kernel[ tabsperpix - taby]
00097 + interp_kernel[2 * tabsperpix - taby]);
00098
00099
00100 offset_ij[0 + p * 2] = (int)floor(offset_x);
00101 offset_ij[1 + p * 2] = (int)floor(offset_y);
00102
00103 ppi[p] = CPL_TYPE_ADD_CONST(cpl_image_get_data)
00104 (cpl_imagelist_get_const(ilist, p));
00105 ppi[p] += offset_ij[0 + p * 2];
00106
00107 if(cpl_image_count_rejected(cpl_imagelist_get_const(ilist, p)) > 0)
00108 bpmflag = 1;
00109
00110
00111
00112
00113
00114 rsc[0 + p * 8] = interp_kernel[ tabsperpix + tabx];
00115 rsc[1 + p * 8] = interp_kernel[ tabx];
00116 rsc[2 + p * 8] = interp_kernel[ tabsperpix - tabx];
00117 rsc[3 + p * 8] = interp_kernel[2 * tabsperpix - tabx];
00118 rsc[4 + p * 8] = interp_kernel[ tabsperpix + taby] / sumrs;
00119 rsc[5 + p * 8] = interp_kernel[ taby] / sumrs;
00120 rsc[6 + p * 8] = interp_kernel[ tabsperpix - taby] / sumrs;
00121 rsc[7 + p * 8] = interp_kernel[2 * tabsperpix - taby] / sumrs;
00122
00123 }
00124
00125
00126 if(bpmflag)
00127 {
00128 for (p=0; p < nima; p++)
00129 {
00130 const cpl_mask * mask;
00131 mask = cpl_image_get_bpm_const(cpl_imagelist_get_const(ilist, p));
00132
00133 if(mask == NULL)
00134 {
00135 ppibpm[p] = NULL;
00136 }
00137 else
00138 {
00139 ppibpm[p] = cpl_mask_get_data_const
00140 (cpl_image_get_bpm_const(cpl_imagelist_get_const(ilist, p)));
00141 ppibpm[p] += offset_ij[0 + p * 2];
00142 }
00143 }
00144 }
00145
00146
00147 for (j=0; j < ny; j++) {
00148 int i;
00149
00150 for (p=0; p < nima; p++) {
00151
00152 const int py = j + offset_ij[1 + p * 2];
00153
00154
00155 if (py > 1 && py < (sizey-2)) {
00156
00157 const CPL_TYPE * pi = ppi[p] + py * sizex;
00158 const cpl_binary * pibpm;
00159
00160 if(bpmflag && ppibpm != NULL)
00161 pibpm = ppibpm[p] + py * sizex;
00162
00163
00164 const int i0 = CX_MAX(0, 2 - offset_ij[0 + p * 2]);
00165
00166 const int i1 = CX_MIN(sizex - 2 - offset_ij[0 + p * 2], nx);
00167
00168 for (i = i0; i < i1; i++) {
00169
00170
00171
00172
00173
00174
00175 if(bpmflag && pibpm != NULL)
00176 {
00177 if((pibpm[i - 1 - sizex] && rsc[4 + 8*p]*rsc[0 + 8*p]) ||
00178 (pibpm[i - sizex] && rsc[4 + 8*p]*rsc[1 + 8*p]) ||
00179 (pibpm[i + 1 - sizex] && rsc[4 + 8*p]*rsc[2 + 8*p]) ||
00180 (pibpm[i + 2 - sizex] && rsc[4 + 8*p]*rsc[3 + 8*p]) ||
00181 (pibpm[i - 1] && rsc[5 + 8*p]*rsc[0 + 8*p]) ||
00182 (pibpm[i ] && rsc[5 + 8*p]*rsc[1 + 8*p]) ||
00183 (pibpm[i + 1] && rsc[5 + 8*p]*rsc[2 + 8*p]) ||
00184 (pibpm[i + 2] && rsc[5 + 8*p]*rsc[3 + 8*p]) ||
00185 (pibpm[i - 1 + sizex] && rsc[6 + 8*p]*rsc[0 + 8*p]) ||
00186 (pibpm[i + sizex] && rsc[6 + 8*p]*rsc[1 + 8*p]) ||
00187 (pibpm[i + 1 + sizex] && rsc[6 + 8*p]*rsc[2 + 8*p]) ||
00188 (pibpm[i + 2 + sizex] && rsc[6 + 8*p]*rsc[3 + 8*p]) ||
00189 (pibpm[i - 1 + sizex * 2] && rsc[7 + 8*p]*rsc[0 + 8*p]) ||
00190 (pibpm[i + sizex * 2] && rsc[7 + 8*p]*rsc[1 + 8*p]) ||
00191 (pibpm[i + 1 + sizex * 2] && rsc[7 + 8*p]*rsc[2 + 8*p]) ||
00192 (pibpm[i + 2 + sizex * 2] && rsc[7 + 8*p]*rsc[3 + 8*p]))
00193 {
00194 continue;
00195 }
00196 }
00197
00198
00199
00200
00201
00202 acc[ncontrib[i] + i * nima] =
00203 rsc[4 + 8*p] * ( rsc[0 + 8*p] * pi[i - 1 - sizex] +
00204 rsc[1 + 8*p] * pi[i - sizex] +
00205 rsc[2 + 8*p] * pi[i + 1 - sizex] +
00206 rsc[3 + 8*p] * pi[i + 2 - sizex] ) +
00207 rsc[5 + 8*p] * ( rsc[0 + 8*p] * pi[i - 1] +
00208 rsc[1 + 8*p] * pi[i ] +
00209 rsc[2 + 8*p] * pi[i + 1] +
00210 rsc[3 + 8*p] * pi[i + 2] ) +
00211 rsc[6 + 8*p] * ( rsc[0 + 8*p] * pi[i - 1 + sizex] +
00212 rsc[1 + 8*p] * pi[i + sizex] +
00213 rsc[2 + 8*p] * pi[i + 1 + sizex] +
00214 rsc[3 + 8*p] * pi[i + 2 + sizex] ) +
00215 rsc[7 + 8*p] * ( rsc[0 + 8*p] * pi[i - 1 + sizex * 2] +
00216 rsc[1 + 8*p] * pi[i + sizex * 2] +
00217 rsc[2 + 8*p] * pi[i + 1 + sizex * 2] +
00218 rsc[3 + 8*p] * pi[i + 2 + sizex * 2] );
00219
00220 ncontrib[i]++;
00221 }
00222 }
00223 }
00224 for (i=0; i < nx; i++) {
00225 if (ncontrib[i] > rtot) {
00226 double finpix = 0.0;
00227 double * acci = acc + i * nima;
00228 hawki_geom_img_get_min_max_double(acci, ncontrib[i],
00229 rmin, rmax);
00230 for (p=rmin; p<(ncontrib[i]-rmax); p++)
00231 finpix += acci[p];
00232 finpix /= (double)(ncontrib[i]-rtot);
00233
00234 po[i+j*nx] = (CPL_TYPE)finpix;
00235 pcontrib[i+j*nx] = ncontrib[i] - rtot;
00236 } else {
00237 cpl_image_reject(final, i+1, j+1);
00238 }
00239 ncontrib[i] = 0;
00240 }
00241 }
00242
00243 cpl_free(offset_ij);
00244 cpl_free(rsc);
00245 cpl_free(acc);
00246 cpl_free(ncontrib);
00247 cpl_free((void*)ppi);
00248
00249 return CPL_ERROR_NONE;
00250 }
00251
00252
00253
00269
00270 static cpl_error_code
00271 CPL_TYPE_ADD(hawki_geom_img_offset_saa_all)(cpl_image * final,
00272 cpl_image * contrib,
00273 double start_x,
00274 double start_y,
00275 const cpl_imagelist * ilist,
00276 const cpl_bivector * offs,
00277 const cpl_vector * xyprofile,
00278 int tabsperpix)
00279
00280 {
00281
00282 const int nima = cpl_imagelist_get_size(ilist);
00283 const cpl_image * img1 = cpl_imagelist_get_const(ilist, 0);
00284 const int sizex = cpl_image_get_size_x(img1);
00285 const int sizey = cpl_image_get_size_y(img1);
00286 const int nx = cpl_image_get_size_x(final);
00287 const int ny = cpl_image_get_size_y(final);
00288
00289 CPL_TYPE * po = CPL_TYPE_ADD(cpl_image_get_data)(final);
00290 const double * offs_x = cpl_bivector_get_x_data_const(offs);
00291 const double * offs_y = cpl_bivector_get_y_data_const(offs);
00292 const double * interp_kernel = cpl_vector_get_data_const(xyprofile);
00293
00294 int * pcontrib = cpl_image_get_data_int(contrib);
00295
00296 unsigned char * pcb = cpl_calloc(1, nx * ny);
00297 const unsigned char * pbjr = pcb;
00298
00299
00300
00301
00302 double rsc[8];
00303 int ninter = 0;
00304 int nadd = 0;
00305 int i, j, p;
00306
00307
00308 for (p = 0; p < nima; p++) {
00309 const cpl_image * imgp = cpl_imagelist_get_const(ilist, p);
00310 const CPL_TYPE * pij = CPL_TYPE_ADD_CONST(cpl_image_get_data)(imgp);
00311 const double offset_x = start_x - offs_x[p];
00312 const double offset_y = start_y - offs_y[p];
00313 const double suboff_x = offset_x - floor(offset_x);
00314 const double suboff_y = offset_y - floor(offset_y);
00315 const int offset_i = (int)floor(offset_x);
00316 const int offset_j = (int)floor(offset_y);
00317
00318 const int tabx = (int)(0.5 + suboff_x * (double)(tabsperpix));
00319 const int taby = (int)(0.5 + suboff_y * (double)(tabsperpix));
00320
00321 CPL_TYPE * poj = po;
00322 int * pcj = pcontrib;
00323 unsigned char * pbj = pcb;
00324
00325 if (tabx == 0 && taby == 0) {
00326
00327
00328
00329
00330 const int px0 = CX_MAX(0, offset_i);
00331 const int py0 = CX_MAX(0, offset_j);
00332
00333 const int px1 = CX_MIN(sizex, nx + offset_i);
00334 const int py1 = CX_MIN(sizey, ny + offset_j);
00335
00336 int px, py;
00337
00338
00339
00340 poj += (py0 - offset_j) * nx - offset_i;
00341 pcj += (py0 - offset_j) * nx - offset_i;
00342 pbj += (py0 - offset_j) * nx - offset_i;
00343 pij += py0 * sizex;
00344
00345 nadd++;
00346 for (py = py0; py < py1; py++, poj += nx, pcj += nx, pbj += nx,
00347 pij += sizex) {
00348
00349 for (px = px0; px < px1; px++) {
00350
00351 poj[px] += pij[px];
00352
00353 if (++pbj[px] == 0) pcj[px] += 256;
00354 }
00355 }
00356 } else {
00357
00358 const double sumrs
00359 = ( interp_kernel[ tabsperpix + tabx]
00360 + interp_kernel[ tabx]
00361 + interp_kernel[ tabsperpix - tabx]
00362 + interp_kernel[2 * tabsperpix - tabx])
00363 * ( interp_kernel[ tabsperpix + taby]
00364 + interp_kernel[ taby]
00365 + interp_kernel[ tabsperpix - taby]
00366 + interp_kernel[2 * tabsperpix - taby]);
00367
00368
00369
00370
00371
00372
00373 const int px0 = CX_MAX(2, offset_i);
00374 const int py0 = CX_MAX(2, offset_j);
00375
00376 const int px1 = CX_MIN(sizex - 2, nx + offset_i);
00377 const int py1 = CX_MIN(sizey - 2, ny + offset_j);
00378
00379 int px, py;
00380
00381 ninter++;
00382
00383
00384
00385
00386
00387
00388
00389 rsc[0] = interp_kernel[ tabsperpix + tabx];
00390 rsc[1] = interp_kernel[ tabx];
00391 rsc[2] = interp_kernel[ tabsperpix - tabx];
00392 rsc[3] = interp_kernel[2 * tabsperpix - tabx];
00393 rsc[4] = interp_kernel[ tabsperpix + taby] / sumrs;
00394 rsc[5] = interp_kernel[ taby] / sumrs;
00395 rsc[6] = interp_kernel[ tabsperpix - taby] / sumrs;
00396 rsc[7] = interp_kernel[2 * tabsperpix - taby] / sumrs;
00397
00398
00399
00400
00401 poj += (py0 - offset_j) * nx - offset_i;
00402 pcj += (py0 - offset_j) * nx - offset_i;
00403 pbj += (py0 - offset_j) * nx - offset_i;
00404 pij += py0 * sizex;
00405
00406 for (py = py0; py < py1; py++, poj += nx, pcj += nx, pbj += nx, pij += sizex) {
00407
00408 for (px = px0; px < px1; px++) {
00409
00410
00411 poj[px] +=
00412 rsc[4] * ( rsc[0] * pij[px - 1 - sizex] +
00413 rsc[1] * pij[px - sizex] +
00414 rsc[2] * pij[px + 1 - sizex] +
00415 rsc[3] * pij[px + 2 - sizex] ) +
00416 rsc[5] * ( rsc[0] * pij[px - 1] +
00417 rsc[1] * pij[px ] +
00418 rsc[2] * pij[px + 1] +
00419 rsc[3] * pij[px + 2] ) +
00420 rsc[6] * ( rsc[0] * pij[px - 1 + sizex] +
00421 rsc[1] * pij[px + sizex] +
00422 rsc[2] * pij[px + 1 + sizex] +
00423 rsc[3] * pij[px + 2 + sizex] ) +
00424 rsc[7] * ( rsc[0] * pij[px - 1 + sizex * 2] +
00425 rsc[1] * pij[px + sizex * 2] +
00426 rsc[2] * pij[px + 1 + sizex * 2] +
00427 rsc[3] * pij[px + 2 + sizex * 2] );
00428
00429 if (++pbj[px] == 0) pcj[px] += 256;
00430 }
00431 }
00432 }
00433 }
00434
00435 for (j = 0; j < ny; j++, po += nx, pcontrib += nx, pbjr += nx) {
00436 for (i = 0; i < nx; i++) {
00437 pcontrib[i] += (int)pbjr[i];
00438 if (pcontrib[i] > 0) {
00439 po[i] /= (CPL_TYPE)pcontrib[i];
00440 } else {
00441 cpl_image_reject(final, i+1, j+1);
00442 }
00443 }
00444 }
00445
00446 cpl_free(pcb);
00447
00448 return CPL_ERROR_NONE;
00449 }
00450
00451
00452 #undef CPL_TYPE_ADD
00453 #undef CPL_TYPE_ADD_CONST