31#include "hdrl_resample.h"
32#include "hdrl_utils.h"
35#define omp_get_max_threads() 1
36#define omp_get_thread_num() 0
48#if !defined(_POSIX_C_SOURCE)
49#define _POSIX_C_SOURCE 200112L
50#elif (_POSIX_C_SOURCE - 0) < 200112L
52#define _POSIX_C_SOURCE 200112L
57#define KEYWORD_LENGTH 81
65#define PT_IDX_MASK 0x1FFFFFFFFFFFFFll
69#define XMAP_BITMASK 0x3FFll
70#define XMAP_LSHIFT 53ll
73 double crpix1, crpix2;
74 double crval1, crval2;
75 double cd11, cd12, cd21, cd22;
77} hdrl_resample_smallwcs;
82} hdrl_resample_pixels_ext;
97 hdrl_resample_pixels_ext **xmaps;
98} hdrl_resample_pixgrid;
207 hdrl_resample_outgrid method;
210 double delta_lambda ;
212 cpl_boolean recalc_limits;
220} hdrl_resample_outgrid_parameter;
223static hdrl_parameter_typeobj hdrl_resample_outgrid_parameter_type = {
224 HDRL_PARAMETER_RESAMPLE_OUTGRID,
225 (hdrl_alloc *)&cpl_malloc,
226 (hdrl_free *)&cpl_free,
228 sizeof(hdrl_resample_outgrid_parameter),
233 hdrl_resample_method method;
236 cpl_boolean use_errorweights;
240 double drizzle_pix_frac_x;
241 double drizzle_pix_frac_y;
242 double drizzle_pix_frac_lambda;
244 double renka_critical_radius;
246 int lanczos_kernel_size;
247} hdrl_resample_method_parameter;
250static hdrl_parameter_typeobj hdrl_resample_method_parameter_type = {
251 HDRL_PARAMETER_RESAMPLE_METHOD,
252 (hdrl_alloc *)&cpl_malloc,
253 (hdrl_free *)&cpl_free,
255 sizeof(hdrl_resample_method_parameter),
263static hdrl_resample_result *
264hdrl_resample_cube(
const cpl_table *ResTable,
265 hdrl_resample_method_parameter *aParams_method,
266 hdrl_resample_outgrid_parameter *aParams_outputgrid,
267 hdrl_resample_pixgrid **aGrid);
270hdrl_resample_inputtable_verify(
const cpl_table *ResTable);
281hdrl_resample_wcs_print(
const cpl_wcs *wcs)
283 cpl_ensure_code(wcs, CPL_ERROR_NULL_INPUT);
285 const cpl_array *crval = cpl_wcs_get_crval(wcs);
286 const cpl_array *crpix = cpl_wcs_get_crpix(wcs);
287 const cpl_array *ctype = cpl_wcs_get_ctype(wcs);
288 const cpl_array *cunit = cpl_wcs_get_cunit(wcs);
290 const cpl_matrix *cd = cpl_wcs_get_cd(wcs);
291 const cpl_array *dims = cpl_wcs_get_image_dims(wcs);
292 int naxis = cpl_wcs_get_image_naxis(wcs);
294 cpl_msg_debug(cpl_func,
"NAXIS: %d", naxis);
298 cpl_msg_indent_more();
300 for (cpl_size i = 0; i < naxis; i++) {
301 cpl_msg_debug(cpl_func,
"NAXIS%lld: %d", i + 1,
302 cpl_array_get_int(dims, i, &testerr));
304 cpl_msg_indent_less();
306 double cd11 = cpl_matrix_get(cd, 0, 0);
307 double cd12 = cpl_matrix_get(cd, 0, 1);
308 double cd21 = cpl_matrix_get(cd, 1, 0);
309 double cd22 = cpl_matrix_get(cd, 1, 1);
310 double crpix1 = cpl_array_get_double(crpix, 0, &testerr);
311 double crpix2 = cpl_array_get_double(crpix, 1, &testerr);
312 double crval1 = cpl_array_get_double(crval, 0, &testerr);
313 double crval2 = cpl_array_get_double(crval, 1, &testerr);
315 cpl_msg_debug(cpl_func,
"1st and 2nd dimension");
316 cpl_msg_indent_more();
317 cpl_msg_debug(cpl_func,
"CD1_1: %g", cd11);
318 cpl_msg_debug(cpl_func,
"CD1_2: %g", cd12);
319 cpl_msg_debug(cpl_func,
"CD2_1: %g", cd21);
320 cpl_msg_debug(cpl_func,
"CD2_2: %g", cd22);
322 cpl_msg_debug(cpl_func,
"CRPIX1: %g", crpix1);
323 cpl_msg_debug(cpl_func,
"CRPIX2: %g", crpix2);
324 cpl_msg_debug(cpl_func,
"CRVAL1: %g", crval1);
325 cpl_msg_debug(cpl_func,
"CRVAL2: %g", crval2);
327 cpl_msg_debug(cpl_func,
"CTYPE1: %s", cpl_array_get_string(ctype, 0));
328 cpl_msg_debug(cpl_func,
"CTYPE2: %s", cpl_array_get_string(ctype, 1));
332 cpl_msg_debug(cpl_func,
"CUNIT1: %s", cpl_array_get_string(cunit, 0));
333 cpl_msg_debug(cpl_func,
"CUNIT2: %s", cpl_array_get_string(cunit, 1));
335 cpl_msg_indent_less();
338 cpl_size cdncol = cpl_matrix_get_ncol(cd);
341 double cd13 = cpl_matrix_get(cd, 0, 2);
342 double cd23 = cpl_matrix_get(cd, 1, 2);
343 double cd31 = cpl_matrix_get(cd, 2, 0);
344 double cd32 = cpl_matrix_get(cd, 2, 1);
345 double cd33 = cpl_matrix_get(cd, 2, 2);
346 double crval3 = cpl_array_get_double(crval, 2, &testerr);
347 double crpix3 = cpl_array_get_double(crpix, 2, &testerr);
349 cpl_msg_debug(cpl_func,
"3rd dimension");
350 cpl_msg_indent_more();
351 cpl_msg_debug(cpl_func,
"CD1_3: %g", cd13);
352 cpl_msg_debug(cpl_func,
"CD2_3: %g", cd23);
353 cpl_msg_debug(cpl_func,
"CD3_1: %g", cd31);
354 cpl_msg_debug(cpl_func,
"CD3_2: %g", cd32);
355 cpl_msg_debug(cpl_func,
"CD3_3: %g", cd33);
357 cpl_msg_debug(cpl_func,
"CRPIX3: %g", crpix3);
358 cpl_msg_debug(cpl_func,
"CRVAL3: %g", crval3);
361 cpl_msg_debug(cpl_func,
"CTYPE3: %s", cpl_array_get_string(ctype, 2));
365 cpl_msg_debug(cpl_func,
"CUNIT3: %s", cpl_array_get_string(cunit, 2));
368 cpl_msg_indent_less();
371 return cpl_error_get_code();
382hdrl_resample_outgrid_parameter_print(hdrl_resample_outgrid_parameter* p)
385 cpl_ensure_code(p, CPL_ERROR_NULL_INPUT);
388 cpl_msg_indent_more();
389 cpl_msg_debug(cpl_func,
"delta_ra: %g", p->delta_ra);
390 cpl_msg_debug(cpl_func,
"delta_dec: %g", p->delta_dec);
391 cpl_msg_debug(cpl_func,
"delta_lambda: %g", p->delta_lambda);
392 cpl_msg_debug(cpl_func,
"ra_min: %g", p->ra_min);
393 cpl_msg_debug(cpl_func,
"ra_max: %g", p->ra_max);
394 cpl_msg_debug(cpl_func,
"dec_min: %g", p->dec_min);
395 cpl_msg_debug(cpl_func,
"dec_max: %g", p->dec_max);
396 cpl_msg_debug(cpl_func,
"lambda_min: %g", p->lambda_min);
397 cpl_msg_debug(cpl_func,
"lambda_max: %g", p->lambda_max);
398 cpl_msg_debug(cpl_func,
"fieldmargin: %g", p->fieldmargin);
399 cpl_msg_debug(cpl_func,
"recalc_limits: %d", p->recalc_limits);
402 hdrl_resample_wcs_print(p->wcs);
403 cpl_msg_indent_less();
404 return cpl_error_get_code();
415hdrl_resample_method_parameter_print(hdrl_resample_method_parameter* p)
418 return cpl_error_get_code();
423 cpl_msg_indent_more();
424 if(p->method == HDRL_RESAMPLE_METHOD_NEAREST) {
425 cpl_msg_debug(cpl_func,
"method: %s",
"NEAREST" );
426 }
else if (p->method == HDRL_RESAMPLE_METHOD_RENKA) {
427 cpl_msg_debug(cpl_func,
"method: %s",
"RENKA" );
428 cpl_msg_debug(cpl_func,
"loop_distance: %d", p->loop_distance);
429 cpl_msg_debug(cpl_func,
"use_errorweights: %s",
430 p->use_errorweights == CPL_TRUE ?
"TRUE" :
"FALSE");
431 cpl_msg_debug(cpl_func,
"renka_critical_radius: %g", p->renka_critical_radius);
432 }
else if (p->method == HDRL_RESAMPLE_METHOD_LINEAR) {
433 cpl_msg_debug(cpl_func,
"method: %s",
"LINEAR" );
434 cpl_msg_debug(cpl_func,
"loop_distance: %d", p->loop_distance);
435 cpl_msg_debug(cpl_func,
"use_errorweights: %s",
436 p->use_errorweights == CPL_TRUE ?
"TRUE" :
"FALSE");
437 }
else if (p->method == HDRL_RESAMPLE_METHOD_QUADRATIC) {
438 cpl_msg_debug(cpl_func,
"method: %s",
"QUADRATIC)" );
439 cpl_msg_debug(cpl_func,
"loop_distance: %d", p->loop_distance);
440 cpl_msg_debug(cpl_func,
"use_errorweights: %s",
441 p->use_errorweights == CPL_TRUE ?
"TRUE" :
"FALSE");
442 }
else if (p->method == HDRL_RESAMPLE_METHOD_DRIZZLE) {
443 cpl_msg_debug(cpl_func,
"method: %s",
"DRIZZLE" );
444 cpl_msg_debug(cpl_func,
"loop_distance: %d", p->loop_distance);
445 cpl_msg_debug(cpl_func,
"use_errorweights: %s",
446 p->use_errorweights == CPL_TRUE ?
"TRUE" :
"FALSE");
447 cpl_msg_debug(cpl_func,
"drizzle_pix_frac_x: %g", p->drizzle_pix_frac_x);
448 cpl_msg_debug(cpl_func,
"drizzle_pix_frac_y: %g", p->drizzle_pix_frac_y);
449 cpl_msg_debug(cpl_func,
"drizzle_pix_frac_lambda: %g", p->drizzle_pix_frac_lambda);
450 }
else if (p->method == HDRL_RESAMPLE_METHOD_LANCZOS) {
451 cpl_msg_debug(cpl_func,
"method: %s",
"LANCZOS" );
452 cpl_msg_debug(cpl_func,
"loop_distance: %d", p->loop_distance);
453 cpl_msg_debug(cpl_func,
"use_errorweights: %s",
454 p->use_errorweights == CPL_TRUE ?
"TRUE" :
"FALSE");
455 cpl_msg_debug(cpl_func,
"lanczos_kernel_size: %d", p->lanczos_kernel_size);
458 cpl_msg_indent_less();
459 return cpl_error_get_code();
477hdrl_wcs_xy_to_radec(
const cpl_wcs *wcs,
double x,
double y,
double *ra,
479 cpl_ensure_code(wcs && ra && dec, CPL_ERROR_NULL_INPUT);
481 double * radec = NULL;
482 cpl_matrix * from = NULL;
483 cpl_matrix * to = NULL;
484 cpl_array * status = NULL;
487 int naxis = cpl_wcs_get_image_naxis(wcs);
488 from = cpl_matrix_new(1, naxis);
489 xy = cpl_matrix_get_data(from);
495 cpl_wcs_convert(wcs, from, &to, &status, CPL_WCS_PHYS2WORLD);
499 radec = cpl_matrix_get_data(to);
505 cpl_matrix_delete(from);
506 cpl_matrix_delete(to);
507 cpl_array_delete(status);
508 return cpl_error_get_code();
523hdrl_resample_pfits_get_crpix(
const cpl_propertylist *aHeaders,
unsigned int aAxis)
525 cpl_errorstate prestate = cpl_errorstate_get();
526 cpl_ensure(aHeaders, CPL_ERROR_NULL_INPUT, 0);
528 char keyword[KEYWORD_LENGTH];
529 snprintf(keyword, KEYWORD_LENGTH,
"CRPIX%u", aAxis);
530 const double value = cpl_propertylist_get_double(aHeaders, keyword);
532 cpl_ensure(cpl_errorstate_is_equal(prestate), cpl_error_get_code(), 0.0);
547hdrl_resample_pfits_get_crval(
const cpl_propertylist *aHeaders,
unsigned int aAxis)
549 cpl_errorstate prestate = cpl_errorstate_get();
550 cpl_ensure(aHeaders, CPL_ERROR_NULL_INPUT, 0);
552 char keyword[KEYWORD_LENGTH];
553 snprintf(keyword, KEYWORD_LENGTH,
"CRVAL%u", aAxis);
554 const double value = cpl_propertylist_get_double(aHeaders, keyword);
556 cpl_ensure(cpl_errorstate_is_equal(prestate), cpl_error_get_code(), 0.0);
572hdrl_resample_pfits_get_cd(
const cpl_propertylist *aHeaders,
unsigned int aAxisI,
575 cpl_errorstate prestate = cpl_errorstate_get();
576 cpl_ensure(aHeaders, CPL_ERROR_NULL_INPUT, 0);
578 char keyword[KEYWORD_LENGTH];
579 snprintf(keyword, KEYWORD_LENGTH,
"CD%u_%u", aAxisI, aAxisJ);
580 const double value = cpl_propertylist_get_double(aHeaders, keyword);
582 cpl_ensure(cpl_errorstate_is_equal(prestate), cpl_error_get_code(), 0.0);
604static hdrl_resample_smallwcs *
605hdrl_resample_smallwcs_new(cpl_propertylist *aHeader)
608 cpl_ensure(aHeader, CPL_ERROR_NULL_INPUT, NULL);
609 hdrl_resample_smallwcs *wcs = cpl_calloc(1,
sizeof(hdrl_resample_smallwcs));
611 cpl_errorstate prestate = cpl_errorstate_get();
612 wcs->crpix1 = hdrl_resample_pfits_get_crpix(aHeader, 1);
613 wcs->crpix2 = hdrl_resample_pfits_get_crpix(aHeader, 2);
614 wcs->crval1 = hdrl_resample_pfits_get_crval(aHeader, 1);
615 wcs->crval2 = hdrl_resample_pfits_get_crval(aHeader, 2);
616 if (!cpl_errorstate_is_equal(prestate)) {
619 cpl_errorstate_set(prestate);
622 prestate = cpl_errorstate_get();
623 wcs->cd11 = hdrl_resample_pfits_get_cd(aHeader, 1, 1);
624 wcs->cd22 = hdrl_resample_pfits_get_cd(aHeader, 2, 2);
625 wcs->cd12 = hdrl_resample_pfits_get_cd(aHeader, 1, 2);
626 wcs->cd21 = hdrl_resample_pfits_get_cd(aHeader, 2, 1);
627 if (!cpl_errorstate_is_equal(prestate) &&
628 wcs->cd11 == 0. && wcs->cd12 == 0. && wcs->cd21 == 0. && wcs->cd22 == 0.) {
631 wcs->cd11 = wcs->cd22 = wcs->cddet = 1.;
632 cpl_errorstate_set(prestate);
634 wcs->cddet = wcs->cd11 * wcs->cd22 - wcs->cd12 * wcs->cd21;
635 if (wcs->cddet == 0.) {
636 cpl_error_set(cpl_func, CPL_ERROR_SINGULAR_MATRIX);
650hdrl_resample_pixgrid_delete(hdrl_resample_pixgrid *aGrid)
655 cpl_free(aGrid->pix);
659 for (ix = 0; ix < aGrid->nmaps; ix++) {
661 for (iext = 0; iext < aGrid->nxalloc[ix]; iext++) {
662 cpl_free(aGrid->xmaps[ix][iext].pix);
665 cpl_free(aGrid->xmaps[ix]);
667 cpl_free(aGrid->xmaps);
669 cpl_free(aGrid->nxmap);
671 cpl_free(aGrid->nxalloc);
672 aGrid->nxalloc = NULL;
734hdrl_resample_result *
736 hdrl_parameter *method,
737 hdrl_parameter *outputgrid,
740 cpl_ensure(ResTable && method, CPL_ERROR_NULL_INPUT, NULL);
741 cpl_ensure(wcs != NULL, CPL_ERROR_NULL_INPUT, NULL);
743 if (hdrl_resample_inputtable_verify(ResTable)) {
753 if (hdrl_resample_inputtable_verify(ResTable)) {
757 hdrl_resample_outgrid_parameter *aParams_outputgrid =
758 (hdrl_resample_outgrid_parameter *) outputgrid;
759 hdrl_resample_method_parameter *aParams_method =
760 (hdrl_resample_method_parameter *) method;
763 aParams_outputgrid->wcs = wcs;
766 cpl_boolean recalc_limits = aParams_outputgrid->recalc_limits;
768 if (recalc_limits == CPL_TRUE) {
769 double ramin = cpl_table_get_column_min(ResTable, HDRL_RESAMPLE_TABLE_RA);
770 double ramax = cpl_table_get_column_max(ResTable, HDRL_RESAMPLE_TABLE_RA);
771 double decmin = cpl_table_get_column_min(ResTable, HDRL_RESAMPLE_TABLE_DEC);
772 double decmax = cpl_table_get_column_max(ResTable, HDRL_RESAMPLE_TABLE_DEC);
773 double lmin = cpl_table_get_column_min(ResTable, HDRL_RESAMPLE_TABLE_LAMBDA);
774 double lmax = cpl_table_get_column_max(ResTable, HDRL_RESAMPLE_TABLE_LAMBDA);
777 if(ramax - ramin > 180){
778 const double *ra = cpl_table_get_data_double_const(ResTable,
779 HDRL_RESAMPLE_TABLE_RA);
784 cpl_size nrow = cpl_table_get_nrow(ResTable);
786 for (cpl_size i = 0; i < nrow; i++) {
787 if (ra[i] > ramin && ra[i] <= 180.) ramin = ra[i];
788 if (ra[i] < ramax && ra[i] > 180.) ramax = ra[i];
792 aParams_outputgrid->ra_min = ramin;
793 aParams_outputgrid->ra_max = ramax;
794 aParams_outputgrid->dec_min = decmin;
795 aParams_outputgrid->dec_max = decmax;
796 aParams_outputgrid->lambda_min = lmin;
797 aParams_outputgrid->lambda_max = lmax;
800 cpl_msg_debug(cpl_func,
"Content of the outgrid parameter structure "
801 "hdrl_resample_outgrid_parameter when resampling starts:");
802 hdrl_resample_outgrid_parameter_print(aParams_outputgrid);
804 cpl_msg_debug(cpl_func,
"Content of the method parameter structure "
805 "hdrl_resample_method_parameter when resampling starts:");
806 hdrl_resample_method_parameter_print(aParams_method);
809 hdrl_resample_result *cube = NULL;
810 hdrl_resample_pixgrid *grid = NULL;
812 cpl_msg_debug(cpl_func,
"Resampling starts ...");
813 cpl_msg_indent_more();
814 cube = hdrl_resample_cube(ResTable,
815 (hdrl_resample_method_parameter *)aParams_method,
816 (hdrl_resample_outgrid_parameter *)aParams_outputgrid,
818 cpl_msg_indent_less();
819 cpl_ensure(cube, CPL_ERROR_NULL_INPUT, NULL);
820 cpl_ensure(cube->header, CPL_ERROR_NULL_INPUT, NULL);
821 cpl_ensure(cube->himlist, CPL_ERROR_NULL_INPUT, NULL);
826 cpl_propertylist *header = cpl_propertylist_new();
827 cpl_wcs* wcs_local = cpl_wcs_new_from_propertylist(cube->header);
828 hdrl_wcs_to_propertylist(wcs_local, header, CPL_TRUE);
829 cpl_wcs_delete(wcs_local);
830 cpl_propertylist_delete(cube->header);
831 cpl_propertylist_set_comment(header,
"CTYPE1",
"Gnomonic projection");
832 cpl_propertylist_set_comment(header,
"CTYPE2",
"Gnomonic projection");
833 cube->header = header;
837 hdrl_resample_pixgrid_delete(grid);
860 aCube->himlist = NULL;
862 cpl_propertylist_delete(aCube->header);
863 aCube->header = NULL;
897hdrl_resample_wcs_projplane_from_celestial(hdrl_resample_outgrid_parameter *
898 aParams_outputgrid,
double aRA,
899 double aDEC,
double *aX,
double *aY)
901 cpl_ensure_code(aParams_outputgrid && aX && aY, CPL_ERROR_NULL_INPUT);
912 const cpl_array *crval = cpl_wcs_get_crval(aParams_outputgrid->wcs);
913 double crval1 = cpl_array_get_double(crval, 0, &err);
914 double crval2 = cpl_array_get_double(crval, 1, &err);
917 double a = aRA / CPL_MATH_DEG_RAD,
918 d = aDEC / CPL_MATH_DEG_RAD,
922 ap = crval1 / CPL_MATH_DEG_RAD,
923 dp = crval2 / CPL_MATH_DEG_RAD,
924 phi = atan2(-cos(d) * sin(a - ap),
925 sin(d) * cos(dp) - cos(d) * sin(dp) * cos(a-ap))
926 + 180 / CPL_MATH_DEG_RAD,
927 theta = asin(sin(d) * sin(dp) + cos(d) * cos(dp) * cos(a-ap)),
928 R_theta = CPL_MATH_DEG_RAD / tan(theta);
930 *aX = R_theta * sin(phi),
931 *aY = -R_theta * cos(phi);
933 return CPL_ERROR_NONE;
956hdrl_resample_wcs_pixel_from_celestial_fast(hdrl_resample_smallwcs *aWCS,
957 double aRA,
double aDEC,
double *aX,
double *aY)
965 double phi = atan2(-cos(aDEC) * sin(aRA - aWCS->crval1),
966 sin(aDEC) * cos(aWCS->crval2)
967 - cos(aDEC) * sin(aWCS->crval2) * cos(aRA - aWCS->crval1))
968 + 180 / CPL_MATH_DEG_RAD,
969 theta = asin(sin(aDEC) * sin(aWCS->crval2)
970 + cos(aDEC) * cos(aWCS->crval2) * cos(aRA - aWCS->crval1)),
971 R_theta = CPL_MATH_DEG_RAD / tan(theta);
973 double x = R_theta * sin(phi),
974 y = -R_theta * cos(phi);
976 *aX = (aWCS->cd22 * x - aWCS->cd12 * y) / aWCS->cddet + aWCS->crpix1;
977 *aY = (aWCS->cd11 * y - aWCS->cd21 * x) / aWCS->cddet + aWCS->crpix2;
996hdrl_resample_compute_size(hdrl_resample_outgrid_parameter *aParams_outputgrid,
997 int *aX,
int *aY,
int *aZ)
999 cpl_ensure_code(aParams_outputgrid && aX && aY && aZ, CPL_ERROR_NULL_INPUT);
1000 const char func[] =
"hdrl_resample_compute_size";
1001 double x1, y1, x2, y2;
1003 double ramin = aParams_outputgrid->ra_min;
1004 double ramax = aParams_outputgrid->ra_max;
1005 double decmin = aParams_outputgrid->dec_min;
1006 double decmax = aParams_outputgrid->dec_max;
1007 hdrl_resample_wcs_projplane_from_celestial(aParams_outputgrid, ramin, decmin,
1009 hdrl_resample_wcs_projplane_from_celestial(aParams_outputgrid, ramax, decmax,
1011 *aX = lround(fabs(x2 - x1) / aParams_outputgrid->delta_ra) + 1;
1012 *aY = lround(fabs(y2 - y1) / aParams_outputgrid->delta_dec) + 1;
1014 double lmin = aParams_outputgrid->lambda_min;
1015 double lmax = aParams_outputgrid->lambda_max;
1017 *aZ = (int)ceil((lmax - lmin) / aParams_outputgrid->delta_lambda) + 1;
1019 cpl_msg_debug(func,
"Output cube size %d x %d x %d (fit to data)",
1021 return CPL_ERROR_NONE;
1042hdrl_resample_pixgrid_add(hdrl_resample_pixgrid *aGrid, cpl_size aIndex,
1043 cpl_size aRow,
unsigned short aXIdx)
1046 if (aIndex < 0 || aGrid == NULL) {
1050 if (aGrid->pix[aIndex] == 0 && aRow > 0) {
1052 aGrid->pix[aIndex] = aRow;
1053 }
else if (aGrid->pix[aIndex] == 0 && aRow == 0) {
1055 cpl_size iext = aGrid->nxmap[aXIdx]++;
1056 if (aGrid->nxmap[aXIdx] > aGrid->nxalloc[aXIdx]) {
1058 aGrid->nxalloc[aXIdx] = 2 * aGrid->nxmap[aXIdx];
1059 aGrid->xmaps[aXIdx] = cpl_realloc(aGrid->xmaps[aXIdx],
1060 aGrid->nxalloc[aXIdx]
1061 *
sizeof(hdrl_resample_pixels_ext));
1063 aGrid->xmaps[aXIdx][iext].npix = 1;
1064 aGrid->xmaps[aXIdx][iext].pix = cpl_malloc(
sizeof(cpl_size));
1065 aGrid->xmaps[aXIdx][iext].pix[0] = aRow;
1066 aGrid->pix[aIndex] = -(iext + 1 + ((cpl_size)aXIdx << XMAP_LSHIFT));
1067 }
else if (aGrid->pix[aIndex] > 0) {
1069 cpl_size iext = aGrid->nxmap[aXIdx]++;
1070 if (aGrid->nxmap[aXIdx] > aGrid->nxalloc[aXIdx]) {
1072 aGrid->nxalloc[aXIdx] = 2 * aGrid->nxmap[aXIdx];
1073 aGrid->xmaps[aXIdx] = cpl_realloc(aGrid->xmaps[aXIdx],
1074 aGrid->nxalloc[aXIdx]
1075 *
sizeof(hdrl_resample_pixels_ext));
1077 aGrid->xmaps[aXIdx][iext].npix = 2;
1078 aGrid->xmaps[aXIdx][iext].pix = cpl_malloc(2 *
sizeof(cpl_size));
1079 aGrid->xmaps[aXIdx][iext].pix[0] = aGrid->pix[aIndex];
1080 aGrid->xmaps[aXIdx][iext].pix[1] = aRow;
1081 aGrid->pix[aIndex] = -(iext + 1 + ((cpl_size)aXIdx << XMAP_LSHIFT));
1084 cpl_size iext = (-aGrid->pix[aIndex] - 1) & PT_IDX_MASK;
1086 unsigned int ipix = aGrid->xmaps[aXIdx][iext].npix;
1087 aGrid->xmaps[aXIdx][iext].npix++;
1088 aGrid->xmaps[aXIdx][iext].pix = cpl_realloc(aGrid->xmaps[aXIdx][iext].pix,
1089 (ipix + 1) *
sizeof(cpl_size));
1090 aGrid->xmaps[aXIdx][iext].pix[ipix] = aRow;
1104static inline cpl_size
1105hdrl_resample_pixgrid_get_count(hdrl_resample_pixgrid *aGrid, cpl_size aIndex)
1107 if (aIndex < 0 || aGrid == NULL) {
1111 cpl_size p = aGrid->pix[aIndex];
1119 unsigned short ix = (-p >> XMAP_LSHIFT) & XMAP_BITMASK;
1120 p = (-p - 1) & PT_IDX_MASK;
1122 return aGrid->xmaps[ix][p].npix;
1140static inline cpl_size
1141hdrl_resample_pixgrid_get_index(hdrl_resample_pixgrid *aGrid, cpl_size aX,
1142 cpl_size aY, cpl_size aZ, cpl_boolean aAllowOutside)
1144 if (aGrid == NULL) {
1147 if (!aAllowOutside &&
1148 (aX < 0 || aX >= aGrid->nx || aY < 0 || aY >= aGrid->ny ||
1149 aZ < 0 || aZ >= aGrid->nz)) {
1155 if (aX >= aGrid->nx) {
1161 if (aY >= aGrid->ny) {
1167 if (aZ >= aGrid->nz) {
1170 return aX + aGrid->nx * (aY + aGrid->ny * aZ);
1184static hdrl_resample_pixgrid *
1185hdrl_resample_pixgrid_new(cpl_size aSizeX, cpl_size aSizeY, cpl_size aSizeZ,
1186 unsigned short aNMaps)
1188 cpl_ensure(aSizeX > 0 && aSizeY > 0 && aSizeZ > 0 && aNMaps > 0,
1189 CPL_ERROR_ILLEGAL_INPUT, NULL);
1190 hdrl_resample_pixgrid *pixels = cpl_calloc(1,
sizeof(hdrl_resample_pixgrid));
1191 pixels->nx = aSizeX;
1192 pixels->ny = aSizeY;
1193 pixels->nz = aSizeZ;
1194 pixels->pix = cpl_calloc(aSizeX * aSizeY * aSizeZ,
sizeof(cpl_size));
1196 pixels->nmaps = aNMaps;
1197 pixels->nxalloc = cpl_calloc(aNMaps,
sizeof(cpl_size));
1198 pixels->xmaps = cpl_calloc(aNMaps,
sizeof(hdrl_resample_pixels_ext *));
1199 pixels->nxmap = cpl_calloc(aNMaps,
sizeof(cpl_size));
1237static hdrl_resample_pixgrid *
1238hdrl_resample_pixgrid_create(
const cpl_table *ResTable, cpl_propertylist *aHeader,
1239 cpl_size aXSize, cpl_size aYSize, cpl_size aZSize)
1241 cpl_ensure(ResTable, CPL_ERROR_NULL_INPUT, NULL);
1242 cpl_ensure(aHeader, CPL_ERROR_NULL_INPUT, NULL);
1243 cpl_size nrow = cpl_table_get_nrow(ResTable);
1245 cpl_msg_error(cpl_func,
"Invalid pixel table (no entries?)");
1246 cpl_error_set(cpl_func, CPL_ERROR_NULL_INPUT);
1249 cpl_ensure(aXSize > 0 && aYSize > 0 && aZSize > 0, CPL_ERROR_ILLEGAL_INPUT,
1253 double crval3 = hdrl_resample_pfits_get_crval(aHeader, 3),
1254 crpix3 = hdrl_resample_pfits_get_crpix(aHeader, 3),
1255 cd33 = hdrl_resample_pfits_get_cd(aHeader, 3, 3);
1257 hdrl_resample_smallwcs *wcs = hdrl_resample_smallwcs_new(aHeader);
1261 const double *xpos = cpl_table_get_data_double_const(ResTable, HDRL_RESAMPLE_TABLE_RA),
1262 *ypos = cpl_table_get_data_double_const(ResTable, HDRL_RESAMPLE_TABLE_DEC),
1263 *lbda = cpl_table_get_data_double_const(ResTable, HDRL_RESAMPLE_TABLE_LAMBDA);
1264 if (!xpos || !ypos || !lbda) {
1265 cpl_msg_error(cpl_func,
"Missing pixel table column (%p %p %p): %s",
1266 (
void *)xpos, (
void *)ypos, (
void *)lbda,
1267 cpl_error_get_message());
1268 cpl_error_set(cpl_func, CPL_ERROR_DATA_NOT_FOUND);
1273 wcs->crval1 /= CPL_MATH_DEG_RAD;
1274 wcs->crval2 /= CPL_MATH_DEG_RAD;
1276 double timeinit = cpl_test_get_walltime(),
1277 timeprogress = timeinit,
1278 cpuinit = cpl_test_get_cputime();
1279 cpl_boolean showprogress = cpl_msg_get_level() == CPL_MSG_DEBUG
1280 || cpl_msg_get_log_level() == CPL_MSG_DEBUG;
1286 cpl_array *asel = NULL;
1287 const cpl_size *sel = NULL;
1288 cpl_size nsel = cpl_table_count_selected(ResTable);
1290 asel = cpl_table_where_selected(ResTable);
1291 sel = cpl_array_get_data_cplsize_const(asel);
1296 int nth = omp_get_max_threads() > XMAP_BITMASK ? XMAP_BITMASK
1297 : omp_get_max_threads();
1299 cpl_array *az1 = cpl_array_new(nth, CPL_TYPE_INT),
1300 *az2 = cpl_array_new(nth, CPL_TYPE_INT);
1303 cpl_array_fill_window_int(az1, aZSize, nth, -1);
1304 cpl_array_fill_window_int(az2, aZSize, nth, -2);
1307 double base = nth > aZSize ? 1. : (double)aZSize / nth;
1309 for (ith = 0; ith < nth && ith < aZSize; ith++) {
1310 cpl_array_set_int(az1, ith, lround(base * ith));
1311 cpl_array_set_int(az2, ith, lround(base * (ith + 1) - 1));
1317 cpl_array_set_int(az1, 0, -INT_MAX / 2 + 1);
1318 cpl_array_set_int(az2, ith - 1, INT_MAX / 2 - 1);
1321 hdrl_resample_pixgrid *grid = hdrl_resample_pixgrid_new(aXSize, aYSize, aZSize, nth);
1325 struct timeval tv1, tv2;
1327 cpl_msg_debug(cpl_func,
"Starting parallel loop in hdrl_resample_pixgrid_create");
1328 gettimeofday(&tv1, NULL);
1330#pragma omp parallel num_threads(nth) \
1331 shared(aXSize, aYSize, aZSize, az1, az2, cd33, crpix3, crval3, grid, \
1332 lbda, nsel, sel, showprogress, \
1333 timeinit, timeprogress, wcs, xpos, ypos)
1337 unsigned short ithread = omp_get_thread_num();
1338 int z1 = cpl_array_get_int(az1, ithread, NULL),
1339 z2 = cpl_array_get_int(az2, ithread, NULL),
1340 zrange = z2 - z1 + 1;
1347 for (isel = 0 ; zrange > 0 && isel < nsel; isel++) {
1349 if (showprogress && !((isel+1) % 1000000ll)) {
1350 double timenow = cpl_test_get_walltime();
1351 if (timenow - timeprogress > 30.) {
1352 timeprogress = timenow;
1353 double percent = 100. * (isel + 1.) / nsel,
1354 elapsed = timeprogress - timeinit,
1355 remaining = (100. - percent) * elapsed / percent;
1358 cpl_msg_info_overwritable(cpl_func,
"pixel grid creation is %.1f%% "
1359 "complete, %.1fs elapsed, ~%.1fs remaining",
1360 percent, elapsed, remaining);
1366 cpl_size n = sel ? sel[isel] : isel;
1369 z = lround((lbda[n] - crval3) / cd33 + crpix3) - 1;
1371 if (z < z1 || z > z2) {
1376 double xpx = 0., ypx = 0.;
1378 hdrl_resample_wcs_pixel_from_celestial_fast(wcs, (xpos[n]) / CPL_MATH_DEG_RAD,
1379 (ypos[n]) / CPL_MATH_DEG_RAD,
1382 int x = lround(xpx) - 1,
1383 y = lround(ypx) - 1;
1384 cpl_size idx = hdrl_resample_pixgrid_get_index(grid, x, y, z, CPL_TRUE);
1388 hdrl_resample_pixgrid_add(grid, idx, n, ithread);
1393 grid->xmaps[ithread] = cpl_realloc(grid->xmaps[ithread],
1394 grid->nxmap[ithread]
1395 *
sizeof(hdrl_resample_pixels_ext));
1396 grid->nxalloc[ithread] = grid->nxmap[ithread];
1399 gettimeofday(&tv2, NULL);
1400 cpl_msg_debug(cpl_func,
"Wall time for hdrl_resample_pixgrid_create was %f seconds\n",
1401 (
double) (tv2.tv_usec - tv1.tv_usec) / 1000000 +
1402 (
double) (tv2.tv_sec - tv1.tv_sec));
1404 cpl_array_delete(asel);
1406 cpl_array_delete(az1);
1407 cpl_array_delete(az2);
1410 cpl_size idx, npix = 0;
1411 for (idx = 0; idx < aXSize * aYSize * aZSize; idx++) {
1412 npix += hdrl_resample_pixgrid_get_count(grid, idx);
1415 for (idx = 0; idx < (cpl_size)grid->nmaps; idx++) {
1416 nxmap += grid->nxmap[idx];
1419 char *msg = cpl_sprintf(
"Pixels got lost while creating the cube (input "
1420 "pixel table: %"CPL_SIZE_FORMAT
", output pixel grid"
1421 ": %"CPL_SIZE_FORMAT
")", nsel, npix);
1422 cpl_msg_error(cpl_func,
"%s:", msg);
1423 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_OUTPUT,
"%s!", msg);
1426 double timefini = cpl_test_get_walltime(),
1427 cpufini = cpl_test_get_cputime();
1428 cpl_msg_debug(cpl_func,
"pixel grid: %dx%dx%d, %"CPL_SIZE_FORMAT
" pixels "
1429 "total, %"CPL_SIZE_FORMAT
" (%.1f%%) in %hu extension maps; took"
1430 " %gs (wall-clock) and %gs (CPU) to create", (
int)grid->nx,
1431 (
int)grid->ny, (
int)grid->nz, npix, nxmap,
1432 (
double)nxmap / npix * 100., grid->nmaps, timefini - timeinit,
1459static cpl_error_code
1460hdrl_resample_wcs_get_scales(hdrl_resample_outgrid_parameter *aParams_outputgrid,
1461 double *aXScale,
double *aYScale)
1463 cpl_ensure_code(aParams_outputgrid && aXScale && aYScale, CPL_ERROR_NULL_INPUT);
1465 cpl_errorstate prestate = cpl_errorstate_get();
1467 const cpl_matrix *cd = cpl_wcs_get_cd(aParams_outputgrid->wcs);
1472 double cd11 = cpl_matrix_get(cd, 0, 0);
1473 double cd12 = cpl_matrix_get(cd, 0, 1);
1474 double cd21 = cpl_matrix_get(cd, 1, 0);
1475 double cd22 = cpl_matrix_get(cd, 1, 1);
1478 double det = cd11 * cd22 - cd12 * cd21;
1479 cpl_ensure_code(cpl_errorstate_is_equal(prestate), cpl_error_get_code());
1485 if (cd12 == 0. && cd21 == 0.) {
1488 return CPL_ERROR_NONE;
1490 *aXScale = sqrt(cd11*cd11 + cd12*cd12);
1491 *aYScale = sqrt(cd22*cd22 + cd21*cd21);
1492 return CPL_ERROR_NONE;
1504static inline const cpl_size *
1505hdrl_resample_pixgrid_get_rows(hdrl_resample_pixgrid *aGrid, cpl_size aIndex)
1508 cpl_ensure(aGrid , CPL_ERROR_NULL_INPUT, 0);
1509 cpl_ensure(aIndex >= 0 , CPL_ERROR_ILLEGAL_INPUT, 0);
1510 cpl_ensure(aIndex < aGrid->nx * aGrid->ny * aGrid->nz, CPL_ERROR_ILLEGAL_INPUT, 0);
1512 cpl_size p = aGrid->pix[aIndex];
1517 return aGrid->pix + aIndex;
1520 unsigned short ix = (-p >> XMAP_LSHIFT) & XMAP_BITMASK;
1521 p = (-p - 1) & PT_IDX_MASK;
1523 return aGrid->xmaps[ix][p].pix;
1545static cpl_error_code
1546hdrl_resample_cube_nearest(hdrl_resample_result *aCube,
const cpl_table *ResTable,
1547 hdrl_resample_pixgrid *aGrid,
1548 hdrl_resample_outgrid_parameter *aParams_outputgrid)
1550 cpl_ensure_code(aCube && ResTable && aGrid && aParams_outputgrid,
1551 CPL_ERROR_NULL_INPUT);
1552 cpl_ensure_code(cpl_propertylist_has(aCube->header,
"CRVAL3") == 1,
1553 CPL_ERROR_ILLEGAL_INPUT);
1554 cpl_ensure_code(cpl_propertylist_has(aCube->header,
"CRPIX3") == 1,
1555 CPL_ERROR_ILLEGAL_INPUT);
1556 cpl_ensure_code(cpl_propertylist_has(aCube->header,
"CD3_3") == 1,
1557 CPL_ERROR_ILLEGAL_INPUT);
1559 double crval3 = hdrl_resample_pfits_get_crval(aCube->header, 3);
1560 double crpix3 = hdrl_resample_pfits_get_crpix(aCube->header, 3);
1561 double cd33 = hdrl_resample_pfits_get_cd(aCube->header, 3, 3);
1562 cpl_ensure_code(cd33 != 0, CPL_ERROR_ILLEGAL_INPUT);
1564 cpl_wcs *wcscpl = cpl_wcs_new_from_propertylist(aCube->header);
1566 if (hdrl_resample_inputtable_verify(ResTable)) {
1567 return CPL_ERROR_ILLEGAL_INPUT;
1571 const double *xpos = cpl_table_get_data_double_const(ResTable, HDRL_RESAMPLE_TABLE_RA),
1572 *ypos = cpl_table_get_data_double_const(ResTable, HDRL_RESAMPLE_TABLE_DEC),
1573 *lbda = cpl_table_get_data_double_const(ResTable, HDRL_RESAMPLE_TABLE_LAMBDA),
1574 *data = cpl_table_get_data_double_const(ResTable, HDRL_RESAMPLE_TABLE_DATA),
1575 *stat = cpl_table_get_data_double_const(ResTable, HDRL_RESAMPLE_TABLE_ERRORS);
1576 const int *dq = cpl_table_get_data_int_const(ResTable, HDRL_RESAMPLE_TABLE_BPM);
1583 double xnorm = 1., ynorm = 1., znorm = 1. ;
1586 hdrl_resample_wcs_get_scales(aParams_outputgrid, &xnorm, &ynorm);
1591 const cpl_matrix *cd = cpl_wcs_get_cd(aParams_outputgrid->wcs);
1592 if (cpl_matrix_get_ncol(cd) == 3) {
1593 znorm = 1. / cpl_matrix_get(cd, 2, 2);
1604 struct timeval tv1, tv2;
1605 cpl_msg_debug(cpl_func,
"Starting parallel loop in hdrl_resample_cube_nearest");
1606 gettimeofday(&tv1, NULL);
1608#pragma omp parallel for default(none) \
1609 shared(aCube, aGrid, cd33, crpix3, crval3, data, dq, lbda, \
1610 stat, stdout, wcscpl, \
1611 xnorm, xpos, ynorm, ypos, znorm) collapse(2)
1613 for (cpl_size l = 0; l < aGrid->nz; l++) {
1614 for (cpl_size i = 0; i < aGrid->nx; i++) {
1619 double lambda = (l + 1. - crpix3) * cd33 + crval3;
1622 for (j = 0; j < aGrid->ny; j++) {
1623 cpl_size idx = hdrl_resample_pixgrid_get_index(aGrid, i, j, l, CPL_FALSE),
1624 n_rows = hdrl_resample_pixgrid_get_count(aGrid, idx);
1625 const cpl_size *rows = hdrl_resample_pixgrid_get_rows(aGrid, idx);
1628 double x = 0., y = 0.;
1632 hdrl_wcs_xy_to_radec(wcscpl, i + 1, j + 1, &x, &y);
1635 if ((cpl_binary)dq[rows[0]] == CPL_BINARY_0) {
1637 pdata[i + j * aGrid->nx] = data[rows[0]];
1638 pstat[i + j * aGrid->nx] = stat[rows[0]];
1639 pdq[i + j * aGrid->nx] = (cpl_binary)dq[rows[0]];
1641 pdq[i + j * aGrid->nx] = CPL_BINARY_1;
1643 }
else if (n_rows >= 2) {
1645 cpl_size n, nbest = -1;
1646 double dbest = FLT_MAX;
1647 for (n = 0; n < n_rows; n++) {
1649 if((cpl_binary)dq[rows[n]] != CPL_BINARY_0) {
1653 double dx = fabs(x - xpos[rows[n]]) * xnorm,
1654 dy = fabs(y - ypos[rows[n]]) * ynorm,
1655 dlambda = fabs(lambda - lbda[rows[n]]) * znorm,
1656 dthis = sqrt(dx*dx + dy*dy + dlambda*dlambda);
1660 dx *= cos(y * CPL_MATH_RAD_DEG);
1662 if (dthis < dbest) {
1668 pdata[i + j * aGrid->nx] = data[rows[nbest]];
1669 pstat[i + j * aGrid->nx] = stat[rows[nbest]];
1670 pdq[i + j * aGrid->nx] = (cpl_binary)dq[rows[nbest]];
1674 pdq[i + j * aGrid->nx] = CPL_BINARY_1;
1680 gettimeofday(&tv2, NULL);
1681 cpl_msg_debug(cpl_func,
"Wall time for hdrl_resample_cube_nearest was %f seconds\n",
1682 (
double) (tv2.tv_usec - tv1.tv_usec) / 1000000 +
1683 (
double) (tv2.tv_sec - tv1.tv_sec));
1688 for(cpl_size i = 0; i < size; i++){
1694 cpl_wcs_delete(wcscpl);
1696 return CPL_ERROR_NONE;
1714hdrl_resample_weight_function_renka(
double r,
double r_c)
1718 }
else if (r >= r_c) {
1721 double p = (r_c - r) / (r_c * r);
1747hdrl_resample_weight_function_drizzle(
double xin,
double yin,
double zin,
1748 double xout,
double yout,
double zout,
1749 double dx,
double dy,
double dz)
1754 double x = (dx + xout / 2.) <= xin / 2. ? xout : (xin + xout) / 2. - dx,
1755 y = (dy + yout / 2.) <= yin / 2. ? yout : (yin + yout) / 2. - dy,
1756 z = (dz + zout / 2.) <= zin / 2. ? zout : (zin + zout) / 2. - dz;
1759 if (x <= 0 || y <= 0 || z <= 0) {
1764 return (x > xin ? xin : x) * (y > yin ? yin : y) * (z > zin ? zin : z)
1765 / (xin * yin * zin);
1782hdrl_resample_weight_function_linear(
double r)
1784 return r == 0 ? FLT_MAX : 1. / r;
1800hdrl_resample_weight_function_quadratic(
double r2)
1802 return r2 == 0 ? FLT_MAX : 1. / r2;
1814hdrl_resample_weight_function_sinc(
double r)
1816 return fabs(r) < DBL_EPSILON ? 1. : sin(CPL_MATH_PI * r) / (CPL_MATH_PI * r);
1832hdrl_resample_weight_function_lanczos(
double dx,
double dy,
double dz,
1833 unsigned int ld,
unsigned int lks)
1837 return (fabs(dx) >= (ld + 0.5) || fabs(dy) >= (ld + 0.5) || fabs(dz) > (ld + 0.5)) ? 0.
1838 : hdrl_resample_weight_function_sinc(dx) * hdrl_resample_weight_function_sinc(dx / lks)
1839 * hdrl_resample_weight_function_sinc(dy) * hdrl_resample_weight_function_sinc(dy / lks)
1840 * hdrl_resample_weight_function_sinc(dz) * hdrl_resample_weight_function_sinc(dz / lks);
1866static cpl_error_code
1867hdrl_resample_cube_weighted(hdrl_resample_result *aCube,
1868 const cpl_table *ResTable,
1869 hdrl_resample_pixgrid *aGrid,
1870 hdrl_resample_method_parameter *aParams_method,
1871 hdrl_resample_outgrid_parameter *aParams_outputgrid)
1874 cpl_ensure_code(aCube && ResTable && aGrid && aParams_method &&
1875 aParams_outputgrid, CPL_ERROR_NULL_INPUT);
1876 cpl_ensure_code(cpl_propertylist_has(aCube->header,
"CRVAL3") == 1,
1877 CPL_ERROR_ILLEGAL_INPUT);
1878 cpl_ensure_code(cpl_propertylist_has(aCube->header,
"CRPIX3") == 1,
1879 CPL_ERROR_ILLEGAL_INPUT);
1880 cpl_ensure_code(cpl_propertylist_has(aCube->header,
"CD3_3") == 1,
1881 CPL_ERROR_ILLEGAL_INPUT);
1883 double crval3 = hdrl_resample_pfits_get_crval(aCube->header, 3),
1884 crpix3 = hdrl_resample_pfits_get_crpix(aCube->header, 3),
1885 cd33 = hdrl_resample_pfits_get_cd(aCube->header, 3, 3);
1887 hdrl_resample_smallwcs *wcs = hdrl_resample_smallwcs_new(aCube->header);
1888 cpl_wcs *wcscpl = cpl_wcs_new_from_propertylist(aCube->header);
1890 const double *xpos = cpl_table_get_data_double_const(ResTable, HDRL_RESAMPLE_TABLE_RA),
1891 *ypos = cpl_table_get_data_double_const(ResTable, HDRL_RESAMPLE_TABLE_DEC),
1892 *lbda = cpl_table_get_data_double_const(ResTable, HDRL_RESAMPLE_TABLE_LAMBDA),
1893 *data = cpl_table_get_data_double_const(ResTable, HDRL_RESAMPLE_TABLE_DATA),
1894 *stat = cpl_table_get_data_double_const(ResTable, HDRL_RESAMPLE_TABLE_ERRORS);
1895 const int *dq = cpl_table_get_data_int_const(ResTable, HDRL_RESAMPLE_TABLE_BPM);
1902 double xnorm = 1., ynorm = 1., znorm = 1. ;
1904 hdrl_resample_wcs_get_scales(aParams_outputgrid, &xnorm, &ynorm);
1908 const cpl_matrix *cd = cpl_wcs_get_cd(aParams_outputgrid->wcs);
1909 if (cpl_matrix_get_ncol(cd) == 3) {
1910 znorm = 1. / cpl_matrix_get(cd, 2, 2);
1923 double renka_rc = aParams_method->renka_critical_radius
1924 * sqrt((wcs->cd11*xnorm)*(wcs->cd11*xnorm)
1925 + (wcs->cd22*ynorm)*(wcs->cd22*ynorm)
1926 + (cd33*znorm)*(cd33*znorm));
1929 int ld = aParams_method->loop_distance;
1932 cpl_msg_debug(cpl_func,
"Overriding loop distance ld=%d", ld);
1936 int lks = aParams_method->lanczos_kernel_size;
1939 cpl_msg_debug(cpl_func,
"Overriding lanczos kernel size lks=%d", lks);
1943 cpl_boolean wght = aParams_method->use_errorweights;
1947 double xsz = aParams_method->drizzle_pix_frac_x / xnorm,
1948 ysz = aParams_method->drizzle_pix_frac_y / ynorm,
1949 zsz = aParams_method->drizzle_pix_frac_lambda / znorm,
1950 xout = fabs(wcs->cd11), yout = fabs(wcs->cd22), zout = fabs(cd33);
1952 struct timeval tv1, tv2;
1953 cpl_msg_debug(cpl_func,
"Starting parallel loop in hdrl_resample_cube_weighted");
1954 gettimeofday(&tv1, NULL);
1956#pragma omp parallel for default(none) \
1957 shared(aCube, aParams_method, aParams_outputgrid, aGrid, ResTable, cd33, crpix3, crval3, data, \
1958 dq, lbda, ld, lks, renka_rc, stat, \
1959 stdout, wcs, wcscpl, wght, xnorm, xout, xpos, xsz, ynorm, yout,\
1960 ypos, ysz, znorm, zout, zsz) collapse(2)
1962 for (cpl_size l = 0; l < aGrid->nz; l++) {
1963 for (cpl_size i = 0; i < aGrid->nx; i++) {
1969 double lambda = (l + 1. - crpix3) * cd33 + crval3;
1970 double zout2 = zout;
1972 for (cpl_size j = 0; j < aGrid->ny; j++) {
1978 hdrl_wcs_xy_to_radec(wcscpl, i + 1, j + 1, &x, &y);
1981 double sumdata = 0, sumstat = 0, sumweight = 0;
1983 cpl_size npoints = 0;
1987 for (i2 = i - ld; i2 <= i + ld; i2++) {
1989 for (j2 = j - ld; j2 <= j + ld; j2++) {
1991 for (l2 = l - ld; l2 <= l + ld; l2++) {
1992 cpl_size idx2 = hdrl_resample_pixgrid_get_index(aGrid, i2, j2, l2, CPL_FALSE);
1993 if (idx2 < 0)
continue;
1994 cpl_size n_rows2 = hdrl_resample_pixgrid_get_count(aGrid, idx2);
1996 const cpl_size *rows2 = hdrl_resample_pixgrid_get_rows(aGrid, idx2);
1998 for (n = 0; n < n_rows2; n++) {
2004 double dx = fabs(x - (xpos[rows2[n]])),
2005 dy = fabs(y - (ypos[rows2[n]])),
2006 dlambda = fabs(lambda - lbda[rows2[n]]),
2016 dx *= cos(y * CPL_MATH_RAD_DEG);
2018 if (aParams_method->method != HDRL_RESAMPLE_METHOD_DRIZZLE) {
2022 r2 = dx * dx + dy * dy + dlambda * dlambda;
2026 if (aParams_method->method == HDRL_RESAMPLE_METHOD_RENKA) {
2027 weight = hdrl_resample_weight_function_renka(sqrt(r2), renka_rc);
2028 }
else if (aParams_method->method == HDRL_RESAMPLE_METHOD_DRIZZLE) {
2029 weight = hdrl_resample_weight_function_drizzle(xsz, ysz, zsz,
2032 }
else if (aParams_method->method == HDRL_RESAMPLE_METHOD_LINEAR) {
2034 weight = hdrl_resample_weight_function_linear(sqrt(r2));
2036 }
else if (aParams_method->method == HDRL_RESAMPLE_METHOD_QUADRATIC) {
2037 weight = hdrl_resample_weight_function_quadratic(r2);
2038 }
else if (aParams_method->method == HDRL_RESAMPLE_METHOD_LANCZOS) {
2039 weight = hdrl_resample_weight_function_lanczos(dx, dy,
2043 if (wght && stat[rows2[n]] > 0.) {
2045 weight /= (stat[rows2[n]] * stat[rows2[n]]);
2048 sumweight += weight;
2049 sumdata += data[rows2[n]] * weight;
2050 flux += data[rows2[n]];
2051 sumstat += stat[rows2[n]] * stat[rows2[n]] * weight * weight;
2066 if (!npoints || !isnormal(sumweight) || !isnormal(sumweight * sumweight)) {
2067 pdq[i + j * aGrid->nx] = CPL_BINARY_1;
2072 sumdata /= sumweight;
2073 sumstat /= sumweight * sumweight;
2075 pdata[i + j * aGrid->nx] = sumdata;
2077 pstat[i + j * aGrid->nx] = sqrt(sumstat);
2078 pdq[i + j * aGrid->nx] = CPL_BINARY_0;
2084 gettimeofday(&tv2, NULL);
2085 cpl_msg_debug(cpl_func,
"Wall time for hdrl_resample_cube_weighted was %f seconds\n",
2086 (
double) (tv2.tv_usec - tv1.tv_usec) / 1000000 +
2087 (
double) (tv2.tv_sec - tv1.tv_sec));
2092 for(cpl_size i = 0; i < size; i++){
2099 cpl_wcs_delete(wcscpl);
2101 return CPL_ERROR_NONE;
2123static cpl_error_code hdrl_resampling_set_outputgrid (
2124 int xsize,
int ysize,
int zsize, hdrl_resample_result *cube,
2125 hdrl_resample_outgrid_parameter *aParams_outputgrid)
2127 cpl_ensure_code(xsize > 0, CPL_ERROR_ILLEGAL_INPUT);
2128 cpl_ensure_code(ysize > 0,CPL_ERROR_ILLEGAL_INPUT);
2129 cpl_ensure_code(zsize >= 0, CPL_ERROR_ILLEGAL_INPUT);
2131 cpl_ensure_code(cube && aParams_outputgrid, CPL_ERROR_NULL_INPUT);
2135 cube->header = cpl_propertylist_new ();
2136 hdrl_wcs_to_propertylist (aParams_outputgrid->wcs, cube->header,
2139 cpl_propertylist_update_string (cube->header,
"CTYPE1",
"RA---TAN");
2140 cpl_propertylist_update_string (cube->header,
"CTYPE2",
"DEC--TAN");
2141 cpl_propertylist_set_comment(cube->header,
"CTYPE1",
"Gnomonic projection");
2142 cpl_propertylist_set_comment(cube->header,
"CTYPE2",
"Gnomonic projection");
2147 cpl_propertylist_update_int (cube->header,
"NAXIS", 3);
2148 cpl_propertylist_update_int (cube->header,
"NAXIS1", xsize);
2149 cpl_propertylist_update_int (cube->header,
"NAXIS2", ysize);
2150 cpl_propertylist_update_int (cube->header,
"NAXIS3", zsize);
2155 cpl_propertylist_update_double (cube->header,
"CD1_1",
2156 -aParams_outputgrid->delta_ra);
2157 cpl_propertylist_update_double (cube->header,
"CD2_2",
2158 aParams_outputgrid->delta_dec);
2159 cpl_propertylist_update_double (cube->header,
"CD1_2", 0.);
2160 cpl_propertylist_update_double (cube->header,
"CD2_1", 0.);
2162 double ramin = aParams_outputgrid->ra_min;
2163 double ramax = aParams_outputgrid->ra_max;
2164 double decmin = aParams_outputgrid->dec_min;
2165 double decmax = aParams_outputgrid->dec_max;
2168 cpl_propertylist_update_double(cube->header,
"CRPIX1", (
double)((xsize + 1) / 2.));
2169 cpl_propertylist_update_double(cube->header,
"CRPIX2", (
double)((ysize + 1) / 2.));
2171 if(ramax - ramin < 180) {
2173 cpl_propertylist_update_double(cube->header,
"CRVAL1", (ramin + ramax) / 2.);
2175 double diff1 = 360. - ramax;
2176 double diff2 = ramin - 0.;
2177 if (diff1 < diff2) {
2178 cpl_propertylist_update_double(cube->header,
"CRVAL1", ramin - (diff1 + diff2) / 2);
2180 cpl_propertylist_update_double(cube->header,
"CRVAL1", ramax + (diff1 + diff2) / 2.);
2183 cpl_propertylist_update_double(cube->header,
"CRVAL2", (decmin + decmax) / 2.);
2184 cpl_propertylist_update_double (cube->header,
"CD3_3",
2185 aParams_outputgrid->delta_lambda);
2186 cpl_propertylist_update_double (cube->header,
"CRPIX3", 1.);
2187 cpl_propertylist_update_double (cube->header,
"CRVAL3",
2188 aParams_outputgrid->lambda_min);
2190 cpl_propertylist_update_double (cube->header,
"CD1_3", 0.);
2191 cpl_propertylist_update_double (cube->header,
"CD2_3", 0.);
2192 cpl_propertylist_update_double (cube->header,
"CD3_1", 0.);
2193 cpl_propertylist_update_double (cube->header,
"CD3_2", 0.);
2195 return cpl_error_get_code();
2253static hdrl_resample_result *
2254hdrl_resample_cube(
const cpl_table *ResTable,
2255 hdrl_resample_method_parameter *aParams_method,
2256 hdrl_resample_outgrid_parameter *aParams_outputgrid,
2257 hdrl_resample_pixgrid **aGrid)
2259 cpl_ensure(ResTable && aParams_method && aParams_outputgrid && aGrid,
2260 CPL_ERROR_NULL_INPUT, NULL);
2269 int xsize = 0, ysize = 0, zsize = 0;
2271 hdrl_resample_compute_size(aParams_outputgrid, &xsize, &ysize, &zsize);
2275 xsize *= (100. + aParams_outputgrid->fieldmargin)/100.;
2276 ysize *= (100. + aParams_outputgrid->fieldmargin)/100.;
2278 cpl_ensure(xsize > 0 && ysize > 0 && zsize > 0, CPL_ERROR_ILLEGAL_OUTPUT,
2281 double time = cpl_test_get_walltime();
2284 hdrl_resample_result *cube = cpl_calloc(1,
sizeof(hdrl_resample_result));
2286 hdrl_resampling_set_outputgrid (xsize, ysize, zsize, cube, aParams_outputgrid);
2288 if (aParams_method->method < HDRL_RESAMPLE_METHOD_NONE) {
2292 for (cpl_size i = 0; i < zsize; i++) {
2297 for (cpl_size j = 1; j <= xsize; j++) {
2298 for (cpl_size k = 1; k <= ysize; k++) {
2308 hdrl_resample_pixgrid *grid = hdrl_resample_pixgrid_create(ResTable, cube->header,
2309 xsize, ysize, zsize);
2312 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
"Could not create"
2320 double timeinit = cpl_test_get_walltime(),
2321 cpuinit = cpl_test_get_cputime();
2324 cpl_error_code rc = CPL_ERROR_NONE;
2325 switch (aParams_method->method) {
2326 case HDRL_RESAMPLE_METHOD_NEAREST:
2327 cpl_msg_debug(cpl_func,
"Starting resampling, using method \"nearest\"");
2328 rc = hdrl_resample_cube_nearest(cube, ResTable, grid, aParams_outputgrid);
2330 case HDRL_RESAMPLE_METHOD_RENKA:
2331 cpl_msg_debug(cpl_func,
"Starting resampling, using method \"renka\" "
2332 "(critical radius rc=%f, loop distance ld=%d)",
2333 aParams_method->renka_critical_radius,
2334 aParams_method->loop_distance);
2335 rc = hdrl_resample_cube_weighted(cube, ResTable, grid, aParams_method,
2336 aParams_outputgrid);
2338 case HDRL_RESAMPLE_METHOD_LINEAR:
2339 case HDRL_RESAMPLE_METHOD_QUADRATIC:
2340 case HDRL_RESAMPLE_METHOD_LANCZOS:
2341 cpl_msg_debug(cpl_func,
"Starting resampling, using method \"%s\" (loop "
2343 aParams_method->method == HDRL_RESAMPLE_METHOD_LINEAR
2345 : (aParams_method->method == HDRL_RESAMPLE_METHOD_QUADRATIC
2348 aParams_method->loop_distance);
2349 rc = hdrl_resample_cube_weighted(cube, ResTable, grid, aParams_method,
2350 aParams_outputgrid);
2352 case HDRL_RESAMPLE_METHOD_DRIZZLE:
2353 cpl_msg_debug(cpl_func,
"Starting resampling, using method \"drizzle\" "
2354 "(pixfrac f=%.3f,%.3f,%.3f, loop distance ld=%d)",
2355 aParams_method->drizzle_pix_frac_x,
2356 aParams_method->drizzle_pix_frac_y,
2357 aParams_method->drizzle_pix_frac_lambda,
2358 aParams_method->loop_distance);
2359 rc = hdrl_resample_cube_weighted(cube, ResTable, grid,
2361 aParams_outputgrid);
2363 case HDRL_RESAMPLE_METHOD_NONE:
2367 cpl_msg_error(cpl_func,
"Don't know this resampling method: %d",
2368 aParams_method->method);
2369 cpl_error_set(cpl_func, CPL_ERROR_UNSUPPORTED_MODE);
2370 rc = CPL_ERROR_UNSUPPORTED_MODE;
2373 double timefini = cpl_test_get_walltime(),
2374 cpufini = cpl_test_get_cputime();
2380 hdrl_resample_pixgrid_delete(grid);
2383 cpl_msg_debug(cpl_func,
"resampling took %.3fs (wall-clock) and %.3fs "
2384 "(%.3fs CPU, %d CPUs) for hdrl_resample_cube*() alone",
2385 timefini - time, timefini - timeinit, cpufini - cpuinit,
2386 omp_get_max_threads());
2387 if (rc != CPL_ERROR_NONE) {
2388 cpl_msg_error(cpl_func,
"resampling failed: %s", cpl_error_get_message());
2408hdrl_wcs_to_propertylist(
const cpl_wcs * wcs, cpl_propertylist * header,
2411 cpl_ensure_code(wcs && header, CPL_ERROR_NULL_INPUT);
2413 const cpl_array *crval = cpl_wcs_get_crval(wcs);
2414 const cpl_array *crpix = cpl_wcs_get_crpix(wcs);
2415 const cpl_array *ctype = cpl_wcs_get_ctype(wcs);
2416 const cpl_array *cunit = cpl_wcs_get_cunit(wcs);
2418 const cpl_matrix *cd = cpl_wcs_get_cd(wcs);
2420 const cpl_array *dims = cpl_wcs_get_image_dims(wcs);
2421 int naxis = cpl_wcs_get_image_naxis(wcs);
2425 for (cpl_size i = 0; i < naxis; i++) {
2427 cpl_propertylist_update_int(header,
"NAXIS", naxis);
2429 char * buf = cpl_sprintf(
"NAXIS%lld", i + 1);
2430 cpl_propertylist_update_int(header, buf,
2431 cpl_array_get_int(dims, i, &err));
2436 if (only2d == TRUE) {
2437 cpl_propertylist_update_int(header,
"NAXIS", 2);
2439 if(cpl_propertylist_has(header,
"NAXIS3")){
2440 cpl_propertylist_erase(header,
"NAXIS3");
2446 cpl_propertylist_update_double(header,
"CRVAL1",
2447 cpl_array_get_double(crval, 0, &err));
2448 cpl_propertylist_update_double(header,
"CRVAL2",
2449 cpl_array_get_double(crval, 1, &err));
2453 cpl_propertylist_update_double(header,
"CRPIX1",
2454 cpl_array_get_double(crpix, 0, &err));
2455 cpl_propertylist_update_double(header,
"CRPIX2",
2456 cpl_array_get_double(crpix, 1, &err));
2460 cpl_propertylist_update_string(header,
"CTYPE1",
2461 cpl_array_get_string(ctype, 0));
2462 cpl_propertylist_update_string(header,
"CTYPE2",
2463 cpl_array_get_string(ctype, 1));
2467 cpl_propertylist_update_string(header,
"CUNIT1",
2468 cpl_array_get_string(cunit, 0));
2469 cpl_propertylist_update_string(header,
"CUNIT2",
2470 cpl_array_get_string(cunit, 1));
2474 double cd11 = cpl_matrix_get(cd, 0, 0);
2475 double cd12 = cpl_matrix_get(cd, 0, 1);
2476 double cd21 = cpl_matrix_get(cd, 1, 0);
2477 double cd22 = cpl_matrix_get(cd, 1, 1);
2478 cpl_propertylist_update_double(header,
"CD1_1", cd11);
2479 cpl_propertylist_update_double(header,
"CD1_2", cd12);
2480 cpl_propertylist_update_double(header,
"CD2_1", cd21);
2481 cpl_propertylist_update_double(header,
"CD2_2", cd22);
2485 if (only2d == FALSE && cpl_array_get_size(crval) > 2) {
2488 cpl_propertylist_update_double(header,
"CRVAL3",
2489 cpl_array_get_double(crval, 2, &err));
2493 cpl_propertylist_update_double(header,
"CRPIX3",
2494 cpl_array_get_double(crpix, 2, &err));
2498 cpl_propertylist_update_string(header,
"CTYPE3",
2499 cpl_array_get_string(ctype, 2));
2503 cpl_propertylist_update_string(header,
"CUNIT3",
2504 cpl_array_get_string(cunit, 2));
2508 double cd13 = cpl_matrix_get(cd, 0, 2);
2509 double cd23 = cpl_matrix_get(cd, 1, 2);
2510 double cd31 = cpl_matrix_get(cd, 2, 0);
2511 double cd32 = cpl_matrix_get(cd, 2, 1);
2512 double cd33 = cpl_matrix_get(cd, 2, 2);
2513 cpl_propertylist_update_double(header,
"CD1_3", cd13);
2514 cpl_propertylist_update_double(header,
"CD2_3", cd23);
2515 cpl_propertylist_update_double(header,
"CD3_1", cd31);
2516 cpl_propertylist_update_double(header,
"CD3_2", cd32);
2517 cpl_propertylist_update_double(header,
"CD3_3", cd33);
2520 return CPL_ERROR_NONE;
2532static cpl_error_code
2533hdrl_resample_create_table(cpl_table** tab,
const cpl_size size)
2536 cpl_ensure_code(tab, CPL_ERROR_NULL_INPUT);
2537 cpl_ensure_code(size > 0 , CPL_ERROR_ILLEGAL_INPUT);
2539 *tab = cpl_table_new(size);
2548 cpl_table_new_column(*tab,
2549 HDRL_RESAMPLE_TABLE_RA, HDRL_RESAMPLE_TABLE_RA_TYPE);
2550 cpl_table_new_column(*tab,
2551 HDRL_RESAMPLE_TABLE_DEC, HDRL_RESAMPLE_TABLE_DEC_TYPE);
2552 cpl_table_new_column(*tab,
2553 HDRL_RESAMPLE_TABLE_LAMBDA, HDRL_RESAMPLE_TABLE_LAMBDA_TYPE);
2554 cpl_table_new_column(*tab,
2555 HDRL_RESAMPLE_TABLE_DATA, HDRL_RESAMPLE_TABLE_DATA_TYPE);
2556 cpl_table_new_column(*tab,
2557 HDRL_RESAMPLE_TABLE_BPM, HDRL_RESAMPLE_TABLE_BPM_TYPE);
2558 cpl_table_new_column(*tab,
2559 HDRL_RESAMPLE_TABLE_ERRORS, HDRL_RESAMPLE_TABLE_ERRORS_TYPE);
2562 cpl_table_fill_column_window_double(*tab, HDRL_RESAMPLE_TABLE_RA, 0, size, 0.);
2563 cpl_table_fill_column_window_double(*tab, HDRL_RESAMPLE_TABLE_DEC, 0, size, 0.);
2564 cpl_table_fill_column_window_double(*tab, HDRL_RESAMPLE_TABLE_LAMBDA,0, size, 0.);
2565 cpl_table_fill_column_window_double(*tab, HDRL_RESAMPLE_TABLE_DATA, 0, size, 0.);
2566 cpl_table_fill_column_window_int(*tab, HDRL_RESAMPLE_TABLE_BPM, 0, size, 0);
2567 cpl_table_fill_column_window_double(*tab, HDRL_RESAMPLE_TABLE_ERRORS,0, size, 0.);
2569 return cpl_error_get_code();
2584 cpl_ensure(hima, CPL_ERROR_NULL_INPUT, NULL);
2585 cpl_ensure(wcs, CPL_ERROR_NULL_INPUT, NULL);
2586 cpl_msg_debug(cpl_func,
"Converting Data to table");
2587 hdrl_imagelist* ilist = NULL;
2616 cpl_ensure(himlist, CPL_ERROR_NULL_INPUT, NULL);
2617 cpl_ensure(wcs, CPL_ERROR_NULL_INPUT, NULL);
2619 cpl_msg_debug(cpl_func,
"Converting Dataset to table");
2625 cpl_msg_debug(cpl_func,
"Dataset dimentions (x, y, l): (%lld, %lld, %lld)",
2626 naxis1, naxis2, naxis3);
2628 const cpl_array *crval = cpl_wcs_get_crval(wcs);
2629 const cpl_array *crpix = cpl_wcs_get_crpix(wcs);
2630 const cpl_matrix *cd = cpl_wcs_get_cd(wcs);
2638 crpix3 = cpl_array_get_double(crpix, 2, &testerr);
2639 crval3 = cpl_array_get_double(crval, 2, &testerr);
2640 cdelt3 = cpl_matrix_get(cd, 2, 2);
2647 cpl_size tab_size = naxis1 * naxis2 * naxis3;
2648 cpl_table* tab = NULL;
2650 hdrl_resample_create_table(&tab, tab_size);
2652 double* ptabxpos = cpl_table_get_data_double(tab, HDRL_RESAMPLE_TABLE_RA);
2653 double* ptabypos = cpl_table_get_data_double(tab, HDRL_RESAMPLE_TABLE_DEC);
2654 double* ptablambda = cpl_table_get_data_double(tab, HDRL_RESAMPLE_TABLE_LAMBDA);
2655 double* ptabdata = cpl_table_get_data_double(tab, HDRL_RESAMPLE_TABLE_DATA);
2656 int* ptabbpm = cpl_table_get_data_int(tab, HDRL_RESAMPLE_TABLE_BPM);
2657 double* ptaberr = cpl_table_get_data_double(tab, HDRL_RESAMPLE_TABLE_ERRORS);
2659 struct timeval tv1, tv2;
2660 cpl_msg_debug(cpl_func,
"Starting parallel loop in hdrl_imagelist_to_table");
2661 gettimeofday(&tv1, NULL);
2663 HDRL_OMP(omp parallel
for collapse(2))
2664 for(cpl_size k = 0; k < naxis3; k++) {
2665 for(cpl_size j = 0; j < naxis2; j++) {
2666 const double* pimaerr = NULL;
2667 const cpl_binary * pimabpm = NULL;
2674 const double* pimadata = cpl_image_get_data_double_const(imadata);
2677 pimaerr = cpl_image_get_data_double_const(imaerrs);
2680 pimabpm = cpl_mask_get_data_const(imamask);
2683 cpl_size k_naxis1_naxis2 = naxis1 * naxis2 * k;
2684 cpl_size j_naxis1 = j * naxis1;
2685 for(cpl_size i = 0; i < naxis1; i++) {
2686 cpl_size raw = k_naxis1_naxis2 + j_naxis1 + i;
2687 hdrl_wcs_xy_to_radec(wcs, i+1 , j+1 , &ptabxpos[raw],
2689 ptabdata[raw] = pimadata[j_naxis1 + i];
2690 if (naxis3 > 1) ptablambda[raw] = crval3 + cdelt3 * (k - crpix3 + 1.);
2691 if (imaerrs) ptaberr[raw] = pimaerr[j_naxis1 + i];
2692 if (imamask) ptabbpm[raw] = pimabpm[j_naxis1 + i];
2694 if (!isfinite(pimadata[j_naxis1 + i]) ||
2695 ptabbpm[raw] != CPL_BINARY_0 ){
2696 ptabbpm[raw] = CPL_BINARY_1;
2703 gettimeofday(&tv2, NULL);
2704 cpl_msg_debug(cpl_func,
"Wall time for hdrl_imagelist_to_table was %f seconds\n",
2705 (
double) (tv2.tv_usec - tv1.tv_usec) / 1000000 +
2706 (
double) (tv2.tv_sec - tv1.tv_sec));
2735 const double delta_dec)
2738 hdrl_resample_outgrid_parameter * p = (hdrl_resample_outgrid_parameter *)
2739 hdrl_parameter_new(&hdrl_resample_outgrid_parameter_type);
2740 p->method = HDRL_RESAMPLE_OUTGRID_2D;
2741 p->delta_ra = delta_ra;
2742 p->delta_dec = delta_dec;
2744 p->recalc_limits = CPL_TRUE;
2757 p->delta_lambda = 1;
2760 p->fieldmargin = FIELDMARGIN;
2768 return (hdrl_parameter *)p;
2797 const double delta_dec,
const double delta_lambda)
2800 hdrl_resample_outgrid_parameter * p = (hdrl_resample_outgrid_parameter *)
2801 hdrl_parameter_new(&hdrl_resample_outgrid_parameter_type);
2802 p->method = HDRL_RESAMPLE_OUTGRID_3D;
2803 p->delta_ra = delta_ra;
2804 p->delta_dec = delta_dec;
2805 p->delta_lambda = delta_lambda;
2807 p->recalc_limits = CPL_TRUE;
2819 p->fieldmargin = FIELDMARGIN;
2827 return (hdrl_parameter *)p;
2857 const double delta_dec,
2858 const double ra_min,
2859 const double ra_max,
2860 const double dec_min,
2861 const double dec_max,
2862 const double fieldmargin)
2865 hdrl_resample_outgrid_parameter * p = (hdrl_resample_outgrid_parameter *)
2866 hdrl_parameter_new(&hdrl_resample_outgrid_parameter_type);
2867 p->method = HDRL_RESAMPLE_OUTGRID_2D;
2868 p->delta_ra = delta_ra;
2869 p->delta_dec = delta_dec;
2872 p->recalc_limits = CPL_FALSE;
2873 p->dec_min = dec_min;
2874 p->dec_max = dec_max;
2881 p->delta_lambda = 1.;
2883 p->fieldmargin = fieldmargin;
2891 return (hdrl_parameter *)p;
2926 const double delta_dec,
2927 const double delta_lambda,
2928 const double ra_min,
2929 const double ra_max,
2930 const double dec_min,
2931 const double dec_max,
2932 const double lambda_min,
2933 const double lambda_max,
2934 const double fieldmargin)
2937 hdrl_resample_outgrid_parameter * p = (hdrl_resample_outgrid_parameter *)
2938 hdrl_parameter_new(&hdrl_resample_outgrid_parameter_type);
2939 p->method = HDRL_RESAMPLE_OUTGRID_3D;
2940 p->delta_ra = delta_ra;
2941 p->delta_dec = delta_dec;
2942 p->delta_lambda = delta_lambda;
2944 p->recalc_limits = CPL_FALSE;
2946 p->dec_min = dec_min;
2947 p->dec_max = dec_max;
2951 p->lambda_min = lambda_min;
2952 p->lambda_max = lambda_max;
2954 p->fieldmargin = fieldmargin;
2962 return (hdrl_parameter *)p;
2991 cpl_boolean use_errorweights,
2992 const double critical_radius)
2995 hdrl_resample_method_parameter * p = (hdrl_resample_method_parameter *)
2996 hdrl_parameter_new(&hdrl_resample_method_parameter_type);
2998 p->method = HDRL_RESAMPLE_METHOD_RENKA;
2999 p->loop_distance = loop_distance;
3000 p->use_errorweights = use_errorweights;
3001 p->renka_critical_radius = critical_radius;
3004 p->drizzle_pix_frac_x = 0.1;
3005 p->drizzle_pix_frac_y = 0.1;
3006 p->drizzle_pix_frac_lambda = 0.1;
3007 p->lanczos_kernel_size = 2;
3013 return (hdrl_parameter*) p;
3038 cpl_boolean use_errorweights)
3041 hdrl_resample_method_parameter * p = (hdrl_resample_method_parameter *)
3042 hdrl_parameter_new(&hdrl_resample_method_parameter_type);
3045 p->method = HDRL_RESAMPLE_METHOD_LINEAR;
3046 p->loop_distance = loop_distance;
3047 p->use_errorweights = use_errorweights;
3050 p->renka_critical_radius = 0.1;
3051 p->drizzle_pix_frac_x = 0.1;
3052 p->drizzle_pix_frac_y = 0.1;
3053 p->drizzle_pix_frac_lambda = 0.1;
3054 p->lanczos_kernel_size = 2;
3060 return (hdrl_parameter*) p;
3086 cpl_boolean use_errorweights)
3089 hdrl_resample_method_parameter * p = (hdrl_resample_method_parameter *)
3090 hdrl_parameter_new(&hdrl_resample_method_parameter_type);
3092 p->method = HDRL_RESAMPLE_METHOD_QUADRATIC;
3093 p->loop_distance = loop_distance;
3094 p->use_errorweights = use_errorweights;
3097 p->renka_critical_radius = 0.1;
3098 p->drizzle_pix_frac_x = 0.1;
3099 p->drizzle_pix_frac_y = 0.1;
3100 p->drizzle_pix_frac_lambda = 0.1;
3101 p->lanczos_kernel_size = 2;
3107 return (hdrl_parameter*) p;
3130 hdrl_resample_method_parameter * p = (hdrl_resample_method_parameter *)
3131 hdrl_parameter_new(&hdrl_resample_method_parameter_type);
3134 p->method = HDRL_RESAMPLE_METHOD_NEAREST;
3135 p->loop_distance = 0;
3136 p->use_errorweights = CPL_FALSE;
3139 p->renka_critical_radius = 0.1;
3140 p->drizzle_pix_frac_x = 0.1;
3141 p->drizzle_pix_frac_y = 0.1;
3142 p->drizzle_pix_frac_lambda = 0.1;
3143 p->lanczos_kernel_size = 2;
3149 return (hdrl_parameter*) p;
3176 cpl_boolean use_errorweights,
3177 const int kernel_size)
3180 hdrl_resample_method_parameter * p = (hdrl_resample_method_parameter *)
3181 hdrl_parameter_new(&hdrl_resample_method_parameter_type);
3183 p->method = HDRL_RESAMPLE_METHOD_LANCZOS;
3184 p->loop_distance = loop_distance;
3185 p->use_errorweights = use_errorweights;
3186 p->lanczos_kernel_size = kernel_size;
3188 p->renka_critical_radius = 0.1;
3189 p->drizzle_pix_frac_x = 0.1;
3190 p->drizzle_pix_frac_y = 0.1;
3191 p->drizzle_pix_frac_lambda = 0.1;
3197 return (hdrl_parameter*) p;
3231 cpl_boolean use_errorweights,
3232 const double pix_frac_x,
3233 const double pix_frac_y,
3234 const double pix_frac_lambda)
3237 hdrl_resample_method_parameter * p = (hdrl_resample_method_parameter *)
3238 hdrl_parameter_new(&hdrl_resample_method_parameter_type);
3240 p->method = HDRL_RESAMPLE_METHOD_DRIZZLE;
3241 p->loop_distance = loop_distance;
3242 p->use_errorweights = use_errorweights;
3243 p->drizzle_pix_frac_x = pix_frac_x;
3244 p->drizzle_pix_frac_y = pix_frac_y;
3245 p->drizzle_pix_frac_lambda = pix_frac_lambda;
3248 p->renka_critical_radius = 0.1;
3249 p->lanczos_kernel_size = 2;
3255 return (hdrl_parameter*) p;
3271 const hdrl_resample_outgrid_parameter * param_loc =
3272 (
const hdrl_resample_outgrid_parameter *)hp ;
3274 cpl_error_ensure(hp != NULL, CPL_ERROR_NULL_INPUT,
3275 return CPL_ERROR_NULL_INPUT,
"NULL Input Parameters");
3278 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
3279 "Here we expect a resample outgrid parameter") ;
3282 param_loc->recalc_limits == CPL_TRUE ||
3283 param_loc->recalc_limits == CPL_FALSE,
3284 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
3285 "Unsupported resample recalc_limits value");
3296 cpl_error_ensure(param_loc->delta_ra > 0, CPL_ERROR_ILLEGAL_INPUT,
3297 return CPL_ERROR_ILLEGAL_INPUT,
"right ascension "
3298 "stepsize must be > 0");
3300 cpl_error_ensure(param_loc->delta_dec > 0, CPL_ERROR_ILLEGAL_INPUT,
3301 return CPL_ERROR_ILLEGAL_INPUT,
"declination stepsize "
3304 cpl_error_ensure(param_loc->delta_lambda > 0, CPL_ERROR_ILLEGAL_INPUT,
3305 return CPL_ERROR_ILLEGAL_INPUT,
"wavelength stepsize "
3309 cpl_error_ensure(param_loc->ra_min >= 0,
3310 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
3311 "Minimum right ascension must be >= 0");
3313 cpl_error_ensure(param_loc->ra_max >= 0,
3314 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
3315 "Maximum right ascension must be >= 0");
3318 cpl_error_ensure(param_loc->lambda_min >= 0,
3319 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
3320 "Minimum wavelength must be >= 0");
3322 cpl_error_ensure(param_loc->lambda_max >= 0,
3323 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
3324 "Maximum wavelength must be >= 0");
3326 cpl_error_ensure(param_loc->fieldmargin >= 0., CPL_ERROR_ILLEGAL_INPUT,
3327 return CPL_ERROR_ILLEGAL_INPUT,
"The field margin must"
3330 cpl_error_ensure(param_loc->ra_max >= param_loc->ra_min ,
3331 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
3332 "The maximum right ascension must be >= "
3333 "the minimum right ascension");
3334 cpl_error_ensure(param_loc->dec_max >= param_loc->dec_min ,
3335 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
3336 "The maximum declination must be >= "
3337 "the minimum declination");
3339 cpl_error_ensure(param_loc->lambda_max >= param_loc->lambda_min ,
3340 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
3341 "The maximum wavelength must be >= "
3342 "the minimum wavelength");
3344 return CPL_ERROR_NONE ;
3360 const hdrl_resample_method_parameter * param_loc =
3361 (
const hdrl_resample_method_parameter *)hp ;
3362 cpl_error_ensure(hp != NULL, CPL_ERROR_NULL_INPUT,
3363 return CPL_ERROR_NULL_INPUT,
"NULL Input Parameters");
3366 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
3367 "Here we expect a resample method parameter") ;
3371 param_loc->method == HDRL_RESAMPLE_METHOD_NEAREST ||
3372 param_loc->method == HDRL_RESAMPLE_METHOD_LINEAR ||
3373 param_loc->method == HDRL_RESAMPLE_METHOD_QUADRATIC ||
3374 param_loc->method == HDRL_RESAMPLE_METHOD_LANCZOS ||
3375 param_loc->method == HDRL_RESAMPLE_METHOD_DRIZZLE ||
3376 param_loc->method == HDRL_RESAMPLE_METHOD_RENKA,
3377 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
3378 "Unsupported resample method");
3381 cpl_error_ensure(param_loc->loop_distance >= 0, CPL_ERROR_ILLEGAL_INPUT,
3382 return CPL_ERROR_ILLEGAL_INPUT,
"The loop distance must "
3385 cpl_error_ensure(param_loc->use_errorweights == CPL_TRUE ||
3386 param_loc->use_errorweights == CPL_FALSE,
3387 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
3388 "Unsupported resample use_errorweights value");
3390 switch (param_loc->method) {
3391 case HDRL_RESAMPLE_METHOD_NEAREST:
3392 case HDRL_RESAMPLE_METHOD_LINEAR:
3393 case HDRL_RESAMPLE_METHOD_QUADRATIC:
3394 case HDRL_RESAMPLE_METHOD_NONE:
3396 case HDRL_RESAMPLE_METHOD_RENKA:
3398 cpl_error_ensure(param_loc->renka_critical_radius > 0.,
3399 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
3400 "Critical radius of the Renka method must be > 0");
3403 case HDRL_RESAMPLE_METHOD_DRIZZLE:
3405 cpl_error_ensure(param_loc->drizzle_pix_frac_x > 0.,
3406 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
3407 "Drizzle down-scaling factor in x direction "
3410 cpl_error_ensure(param_loc->drizzle_pix_frac_y > 0.,
3411 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
3412 "Drizzle down-scaling factor in y direction "
3415 cpl_error_ensure(param_loc->drizzle_pix_frac_lambda > 0.,
3416 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
3417 "Drizzle down-scaling factor in z/lambda direction "
3421 case HDRL_RESAMPLE_METHOD_LANCZOS:
3423 cpl_error_ensure(param_loc->lanczos_kernel_size > 0,
3424 CPL_ERROR_ILLEGAL_INPUT,
return CPL_ERROR_ILLEGAL_INPUT,
3425 "The kernel size of the Lanczos method must be > 0");
3429 return CPL_ERROR_NONE ;
3444 cpl_ensure_code(self != NULL, CPL_ERROR_NULL_INPUT);
3446 return hdrl_parameter_check_type(self, &hdrl_resample_outgrid_parameter_type);
3461 cpl_ensure_code(self != NULL, CPL_ERROR_NULL_INPUT);
3462 return hdrl_parameter_check_type(self, &hdrl_resample_method_parameter_type);
3474static cpl_error_code
3475hdrl_resample_inputtable_verify(
const cpl_table *ResTable){
3477 cpl_error_ensure(ResTable != NULL, CPL_ERROR_NULL_INPUT,
3478 return CPL_ERROR_NULL_INPUT,
"No Table as input");
3481 cpl_error_ensure(cpl_table_has_column(ResTable, HDRL_RESAMPLE_TABLE_DATA )
3482 == 1, CPL_ERROR_INCOMPATIBLE_INPUT,
3483 return CPL_ERROR_INCOMPATIBLE_INPUT,
3484 "Missing data table column");
3485 cpl_error_ensure(cpl_table_has_column(ResTable, HDRL_RESAMPLE_TABLE_BPM )
3486 == 1, CPL_ERROR_INCOMPATIBLE_INPUT,
3487 return CPL_ERROR_INCOMPATIBLE_INPUT,
3488 "Missing bpm table column");
3489 cpl_error_ensure(cpl_table_has_column(ResTable, HDRL_RESAMPLE_TABLE_ERRORS)
3490 == 1, CPL_ERROR_INCOMPATIBLE_INPUT,
3491 return CPL_ERROR_INCOMPATIBLE_INPUT,
3492 "Missing error table column");
3493 cpl_error_ensure(cpl_table_has_column(ResTable, HDRL_RESAMPLE_TABLE_RA )
3494 == 1, CPL_ERROR_INCOMPATIBLE_INPUT,
3495 return CPL_ERROR_INCOMPATIBLE_INPUT,
3496 "Missing right ascension table column");
3497 cpl_error_ensure(cpl_table_has_column(ResTable, HDRL_RESAMPLE_TABLE_DEC )
3498 == 1, CPL_ERROR_INCOMPATIBLE_INPUT,
3499 return CPL_ERROR_INCOMPATIBLE_INPUT,
3500 "Missing declination table column");
3501 cpl_error_ensure(cpl_table_has_column(ResTable, HDRL_RESAMPLE_TABLE_LAMBDA)
3502 == 1, CPL_ERROR_INCOMPATIBLE_INPUT,
3503 return CPL_ERROR_INCOMPATIBLE_INPUT,
3504 "Missing wavelength table column");
3507 cpl_error_ensure(cpl_table_get_column_type(ResTable, HDRL_RESAMPLE_TABLE_DATA )
3508 == HDRL_RESAMPLE_TABLE_DATA_TYPE, CPL_ERROR_INCOMPATIBLE_INPUT,
3509 return CPL_ERROR_INCOMPATIBLE_INPUT,
3510 "Data table column has wrong format");
3511 cpl_error_ensure(cpl_table_get_column_type(ResTable, HDRL_RESAMPLE_TABLE_BPM )
3512 == HDRL_RESAMPLE_TABLE_BPM_TYPE, CPL_ERROR_INCOMPATIBLE_INPUT,
3513 return CPL_ERROR_INCOMPATIBLE_INPUT,
3514 "Bpm table column has wrong format");
3515 cpl_error_ensure(cpl_table_get_column_type(ResTable, HDRL_RESAMPLE_TABLE_ERRORS)
3516 == HDRL_RESAMPLE_TABLE_ERRORS_TYPE, CPL_ERROR_INCOMPATIBLE_INPUT,
3517 return CPL_ERROR_INCOMPATIBLE_INPUT,
3518 "Error table column has wrong format");
3519 cpl_error_ensure(cpl_table_get_column_type(ResTable, HDRL_RESAMPLE_TABLE_RA )
3520 == HDRL_RESAMPLE_TABLE_RA_TYPE, CPL_ERROR_INCOMPATIBLE_INPUT,
3521 return CPL_ERROR_INCOMPATIBLE_INPUT,
3522 "Right ascension table column has wrong format");
3523 cpl_error_ensure(cpl_table_get_column_type(ResTable, HDRL_RESAMPLE_TABLE_DEC )
3524 == HDRL_RESAMPLE_TABLE_DEC_TYPE, CPL_ERROR_INCOMPATIBLE_INPUT,
3525 return CPL_ERROR_INCOMPATIBLE_INPUT,
3526 "Declination table column has wrong format");
3527 cpl_error_ensure(cpl_table_get_column_type(ResTable, HDRL_RESAMPLE_TABLE_LAMBDA)
3528 == HDRL_RESAMPLE_TABLE_LAMBDA_TYPE, CPL_ERROR_INCOMPATIBLE_INPUT,
3529 return CPL_ERROR_INCOMPATIBLE_INPUT,
3530 "Wavelength table column has wrong format");
3532 return cpl_error_get_code();
cpl_mask * hdrl_image_get_mask(hdrl_image *himg)
get cpl bad pixel mask from image
cpl_image * hdrl_image_get_error(hdrl_image *himg)
get error as cpl image
const cpl_mask * hdrl_image_get_mask_const(const hdrl_image *himg)
get cpl bad pixel mask from image
const cpl_image * hdrl_image_get_error_const(const hdrl_image *himg)
get error as cpl image
cpl_image * hdrl_image_get_image(hdrl_image *himg)
get data as cpl image
cpl_error_code hdrl_image_reject(hdrl_image *self, cpl_size xpos, cpl_size ypos)
mark pixel as bad
const cpl_image * hdrl_image_get_image_const(const hdrl_image *himg)
get data as cpl image
hdrl_image * hdrl_image_new(cpl_size nx, cpl_size ny)
create new zero filled hdrl image
hdrl_image * hdrl_imagelist_unset(hdrl_imagelist *himlist, cpl_size pos)
Remove an image from an imagelist.
cpl_error_code hdrl_imagelist_set(hdrl_imagelist *himlist, hdrl_image *himg, cpl_size pos)
Insert an image into an imagelist.
cpl_size hdrl_imagelist_get_size_y(const hdrl_imagelist *himlist)
Get number of rows of images in the imagelist.
void hdrl_imagelist_delete(hdrl_imagelist *himlist)
Free all memory used by a hdrl_imagelist object including the images.
const hdrl_image * hdrl_imagelist_get_const(const hdrl_imagelist *himlist, cpl_size inum)
Get an image from a list of images.
cpl_size hdrl_imagelist_get_size(const hdrl_imagelist *himlist)
Get the number of images in the imagelist.
hdrl_imagelist * hdrl_imagelist_new(void)
Create an empty imagelist.
cpl_size hdrl_imagelist_get_size_x(const hdrl_imagelist *himlist)
Get number of colums of images in the imagelist.
hdrl_image * hdrl_imagelist_get(const hdrl_imagelist *himlist, cpl_size inum)
Get an image from a list of images.
hdrl_parameter * hdrl_resample_parameter_create_nearest(void)
Creates a resample nearest neighbor hdrl parameter object. The algorithm does not use any weighting f...
cpl_error_code hdrl_resample_parameter_method_verify(const hdrl_parameter *hp)
verify parameters have proper values
int hdrl_resample_parameter_method_check(const hdrl_parameter *self)
check method is of proper type
hdrl_parameter * hdrl_resample_parameter_create_renka(const int loop_distance, cpl_boolean use_errorweights, const double critical_radius)
Creates a resample renka hdrl parameter object. The algorithm uses a modified Shepard-like distance w...
hdrl_parameter * hdrl_resample_parameter_create_lanczos(const int loop_distance, cpl_boolean use_errorweights, const int kernel_size)
Creates a resample Lanczos hdrl parameter object. The algorithm uses a restricted SINC distance weigh...
hdrl_parameter * hdrl_resample_parameter_create_drizzle(const int loop_distance, cpl_boolean use_errorweights, const double pix_frac_x, const double pix_frac_y, const double pix_frac_lambda)
Creates a resample drizzle hdrl parameter object. The algorithm uses a drizzle-like distance weightin...
hdrl_parameter * hdrl_resample_parameter_create_quadratic(const int loop_distance, cpl_boolean use_errorweights)
Creates a resample quadratic hdrl parameter object. The algorithm uses a quadratic inverse distance w...
hdrl_parameter * hdrl_resample_parameter_create_linear(const int loop_distance, cpl_boolean use_errorweights)
Creates a resample linear hdrl parameter object. The algorithm uses a linear inverse distance weighti...
hdrl_parameter * hdrl_resample_parameter_create_outgrid3D(const double delta_ra, const double delta_dec, const double delta_lambda)
Creates a resample_outgrid hdrl parameter object for a 3 dimensional interpolation,...
hdrl_parameter * hdrl_resample_parameter_create_outgrid3D_userdef(const double delta_ra, const double delta_dec, const double delta_lambda, const double ra_min, const double ra_max, const double dec_min, const double dec_max, const double lambda_min, const double lambda_max, const double fieldmargin)
Creates a resample_outgrid hdrl parameter object for a 3 dimensional interpolation,...
int hdrl_resample_parameter_outgrid_check(const hdrl_parameter *self)
check method is of proper type
hdrl_parameter * hdrl_resample_parameter_create_outgrid2D(const double delta_ra, const double delta_dec)
Creates a resample_outgrid hdrl parameter object for a 2 dimensional interpolation,...
hdrl_parameter * hdrl_resample_parameter_create_outgrid2D_userdef(const double delta_ra, const double delta_dec, const double ra_min, const double ra_max, const double dec_min, const double dec_max, const double fieldmargin)
Creates a resample_outgrid hdrl parameter object for a 2 dimensional interpolation,...
cpl_error_code hdrl_resample_parameter_outgrid_verify(const hdrl_parameter *hp)
verify parameters have proper values
hdrl_resample_result * hdrl_resample_compute(const cpl_table *ResTable, hdrl_parameter *method, hdrl_parameter *outputgrid, const cpl_wcs *wcs)
High level resampling function.
void hdrl_resample_result_delete(hdrl_resample_result *aCube)
Deallocates the memory associated to a hdrl_resample_result object.
cpl_table * hdrl_resample_image_to_table(const hdrl_image *hima, const cpl_wcs *wcs)
Convert a hdrl image into a cpl table that can be given as input to hdrl_resample_compute()
cpl_table * hdrl_resample_imagelist_to_table(const hdrl_imagelist *himlist, const cpl_wcs *wcs)
Convert a hdrl imagelist into a cpl table that can be given as input to hdrl_resample_compute().