CR2RE Pipeline Reference Manual 1.6.10
hdrl_resample.c
1/*
2 * This file is part of the HDRLDEMO pipeline
3 * Copyright (C) 2020 European Southern Observatory
4 *
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
18 */
19
20#ifdef HAVE_CONFIG_H
21#include <config.h>
22#endif
23
24/*-----------------------------------------------------------------------------
25 Includes
26 -----------------------------------------------------------------------------*/
27#include <cpl.h>
28#include <string.h>
29#include <math.h>
30
31#include "hdrl_resample.h"
32#include "hdrl_utils.h"
33
34#ifndef _OPENMP
35#define omp_get_max_threads() 1
36#define omp_get_thread_num() 0
37#else
38#include <omp.h>
39#endif
40
41#include <sys/time.h>
42
43/*-----------------------------------------------------------------------------
44 Defines
45 -----------------------------------------------------------------------------*/
46
47/* Not very elegant but it removes compiler warnings */
48#if !defined(_POSIX_C_SOURCE)
49#define _POSIX_C_SOURCE 200112L
50#elif (_POSIX_C_SOURCE - 0) < 200112L
51#undef _POSIX_C_SOURCE
52#define _POSIX_C_SOURCE 200112L
53#endif
54
55
56/* maximum keyword length for FITS headers, including '\0' */
57#define KEYWORD_LENGTH 81
58
59/* Default field margin (in percent), if the user does not specify any. 5
60 * percent is also used in the software package swarp */
61#define FIELDMARGIN 5.
62
63/* use bits 0-52 for the value (the pixel table row), this allows to convert
64 * pixel tables with up to 9e15 pixels into a pixel grid */
65#define PT_IDX_MASK 0x1FFFFFFFFFFFFFll
66
67/* use bits 53-62 to store the thread ID, this allows parallelization with up to
68 * 1024 cores */
69#define XMAP_BITMASK 0x3FFll /* 1023 */
70#define XMAP_LSHIFT 53ll
71
72typedef struct {
73 double crpix1, crpix2; /* the reference point */
74 double crval1, crval2; /* the coordinate values at the ref. point */
75 double cd11, cd12, cd21, cd22; /* the four elements of the CDi_j matrix */
76 double cddet; /* the determinant of the CDi_j matrix */
77} hdrl_resample_smallwcs;
78
79typedef struct {
80 unsigned int npix; /* number of pixels in this grid point */
81 cpl_size *pix; /* the row number(s) in the pixel table */
82} hdrl_resample_pixels_ext;
83
84typedef struct {
85 /* The pixel grid array, elements can be:
86 0: empty
87 positive: row_number in the pixel table
88 negative: -(i_ext+1) in the extension array,
89 bits 53-62 contain the map index */
90 cpl_size *pix;
91 cpl_size nx; /* horizontal spatial size */
92 cpl_size ny; /* vertical spatial size */
93 cpl_size nz; /* size in dispersion direction */
94 unsigned short nmaps; /* number of extension maps */
95 cpl_size *nxmap; /* number of filled pixels in the extension maps */
96 cpl_size *nxalloc; /* number of allocated pixels in the extension maps */
97 hdrl_resample_pixels_ext **xmaps; /* the extension maps */
98} hdrl_resample_pixgrid;
99
100/*----------------------------------------------------------------------------*/
153/*----------------------------------------------------------------------------*/
154
157/*----------------------------------------------------------------------------*/
176/*----------------------------------------------------------------------------*/
177
178/*----------------------------------------------------------------------------*/
197/*----------------------------------------------------------------------------*/
198
201/*-----------------------------------------------------------------------------
202 RESAMPLE Parameter Definition
203 -----------------------------------------------------------------------------*/
204
205typedef struct {
206 HDRL_PARAMETER_HEAD;
207 hdrl_resample_outgrid method;
208 double delta_ra ; /* step size in right ascension [deg] */
209 double delta_dec ; /* step size in declination [deg] */
210 double delta_lambda ; /* step size in wavelength direction [m] */
211 const cpl_wcs* wcs ; /* World Coordinate System */
212 cpl_boolean recalc_limits;
213 double ra_min; /* Minimal Right ascension [deg] */
214 double ra_max; /* Maximal Right ascension [deg] */
215 double dec_min; /* Minimal Declination [deg] */
216 double dec_max; /* Maximal Declination [deg] */
217 double lambda_min; /* Minimal wavelength [m] */
218 double lambda_max; /* Maximal wavelength [m] */
219 double fieldmargin; /* Field margin to add [percent] */
220} hdrl_resample_outgrid_parameter;
221
222/* Parameter type */
223static hdrl_parameter_typeobj hdrl_resample_outgrid_parameter_type = {
224 HDRL_PARAMETER_RESAMPLE_OUTGRID, /* type */
225 (hdrl_alloc *)&cpl_malloc, /* fp_alloc */
226 (hdrl_free *)&cpl_free, /* fp_free */
227 NULL, /* fp_destroy */
228 sizeof(hdrl_resample_outgrid_parameter), /* obj_size */
229};
230
231typedef struct {
232 HDRL_PARAMETER_HEAD;
233 hdrl_resample_method method;
234 int loop_distance;
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;
248
249/* Parameter type */
250static hdrl_parameter_typeobj hdrl_resample_method_parameter_type = {
251 HDRL_PARAMETER_RESAMPLE_METHOD, /* type */
252 (hdrl_alloc *)&cpl_malloc, /* fp_alloc */
253 (hdrl_free *)&cpl_free, /* fp_free */
254 NULL, /* fp_destroy */
255 sizeof(hdrl_resample_method_parameter), /* obj_size */
256};
257
259/*-----------------------------------------------------------------------------
260 Functions
261 -----------------------------------------------------------------------------*/
262
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);
268
269static cpl_error_code
270hdrl_resample_inputtable_verify(const cpl_table *ResTable);
271
272/*----------------------------------------------------------------------------*/
278/*----------------------------------------------------------------------------*/
279
280static cpl_error_code
281hdrl_resample_wcs_print(const cpl_wcs *wcs)
282{
283 cpl_ensure_code(wcs, CPL_ERROR_NULL_INPUT);
284
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);
289
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);
293
294 cpl_msg_debug(cpl_func, "NAXIS: %d", naxis);
295
296 int testerr = 0;
297
298 cpl_msg_indent_more();
299 /* Check NAXIS */
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));
303 }
304 cpl_msg_indent_less();
305
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);
314
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);
321
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);
326 if (ctype) {
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));
329 }
330
331 if (cunit) {
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));
334 }
335 cpl_msg_indent_less();
336
337 /* Is it a 3D cube or a 2D image */
338 cpl_size cdncol = cpl_matrix_get_ncol(cd);
339 if (cdncol == 3) {
340
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);
348
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);
356
357 cpl_msg_debug(cpl_func, "CRPIX3: %g", crpix3);
358 cpl_msg_debug(cpl_func, "CRVAL3: %g", crval3);
359
360 if (ctype) {
361 cpl_msg_debug(cpl_func, "CTYPE3: %s", cpl_array_get_string(ctype, 2));
362 }
363
364 if (cunit) {
365 cpl_msg_debug(cpl_func, "CUNIT3: %s", cpl_array_get_string(cunit, 2));
366 }
367
368 cpl_msg_indent_less();
369 }
370
371 return cpl_error_get_code();
372}
373
374/*----------------------------------------------------------------------------*/
380/*----------------------------------------------------------------------------*/
381static cpl_error_code
382hdrl_resample_outgrid_parameter_print(hdrl_resample_outgrid_parameter* p)
383{
384
385 cpl_ensure_code(p, CPL_ERROR_NULL_INPUT);
386 /* Content of the outgrid parameter structure */
387
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);
400
401 /* World Coordinate System */
402 hdrl_resample_wcs_print(p->wcs);
403 cpl_msg_indent_less();
404 return cpl_error_get_code();
405
406}
407/*----------------------------------------------------------------------------*/
413/*----------------------------------------------------------------------------*/
414static cpl_error_code
415hdrl_resample_method_parameter_print(hdrl_resample_method_parameter* p)
416{
417 if (hdrl_resample_parameter_method_verify((hdrl_parameter*)p) != CPL_ERROR_NONE) {
418 return cpl_error_get_code();
419 }
420
421 /* Content of the outgrid parameter structure */
422
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);
456 }
457
458 cpl_msg_indent_less();
459 return cpl_error_get_code();
460}
461
462/*----------------------------------------------------------------------------*/
475/*----------------------------------------------------------------------------*/
476cpl_error_code
477hdrl_wcs_xy_to_radec(const cpl_wcs *wcs, double x, double y, double *ra,
478 double *dec) {
479 cpl_ensure_code(wcs && ra && dec, CPL_ERROR_NULL_INPUT);
480 double * xy = NULL;
481 double * radec = NULL;
482 cpl_matrix * from = NULL;
483 cpl_matrix * to = NULL;
484 cpl_array * status = NULL;
485
486 /* Load up the information */
487 int naxis = cpl_wcs_get_image_naxis(wcs);
488 from = cpl_matrix_new(1, naxis);
489 xy = cpl_matrix_get_data(from);
490 xy[0] = x;
491 xy[1] = y;
492
493 /* Call the conversion routine */
494
495 cpl_wcs_convert(wcs, from, &to, &status, CPL_WCS_PHYS2WORLD);
496
497 /* Pass it back now */
498
499 radec = cpl_matrix_get_data(to);
500 *ra = radec[0];
501 *dec = radec[1];
502
503 /* Tidy and exit */
504
505 cpl_matrix_delete(from);
506 cpl_matrix_delete(to);
507 cpl_array_delete(status);
508 return cpl_error_get_code();
509}
510
511
512/*----------------------------------------------------------------------------*/
521/*----------------------------------------------------------------------------*/
522static double
523hdrl_resample_pfits_get_crpix(const cpl_propertylist *aHeaders, unsigned int aAxis)
524{
525 cpl_errorstate prestate = cpl_errorstate_get();
526 cpl_ensure(aHeaders, CPL_ERROR_NULL_INPUT, 0);
527
528 char keyword[KEYWORD_LENGTH];
529 snprintf(keyword, KEYWORD_LENGTH, "CRPIX%u", aAxis);
530 const double value = cpl_propertylist_get_double(aHeaders, keyword);
531 /* default to 0.0 as per FITS Standard v3.0 */
532 cpl_ensure(cpl_errorstate_is_equal(prestate), cpl_error_get_code(), 0.0);
533 return value;
534}
535
536/*----------------------------------------------------------------------------*/
545/*----------------------------------------------------------------------------*/
546static double
547hdrl_resample_pfits_get_crval(const cpl_propertylist *aHeaders, unsigned int aAxis)
548{
549 cpl_errorstate prestate = cpl_errorstate_get();
550 cpl_ensure(aHeaders, CPL_ERROR_NULL_INPUT, 0);
551
552 char keyword[KEYWORD_LENGTH];
553 snprintf(keyword, KEYWORD_LENGTH, "CRVAL%u", aAxis);
554 const double value = cpl_propertylist_get_double(aHeaders, keyword);
555 /* default to 0.0 as per FITS Standard v3.0 */
556 cpl_ensure(cpl_errorstate_is_equal(prestate), cpl_error_get_code(), 0.0);
557 return value;
558}
559
560/*----------------------------------------------------------------------------*/
570/*----------------------------------------------------------------------------*/
571static double
572hdrl_resample_pfits_get_cd(const cpl_propertylist *aHeaders, unsigned int aAxisI,
573 unsigned int aAxisJ)
574{
575 cpl_errorstate prestate = cpl_errorstate_get();
576 cpl_ensure(aHeaders, CPL_ERROR_NULL_INPUT, 0);
577
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);
581 /* default to 0.0 as per FITS Standard v3.0 */
582 cpl_ensure(cpl_errorstate_is_equal(prestate), cpl_error_get_code(), 0.0);
583 return value;
584}
585
586/*----------------------------------------------------------------------------*/
603/*----------------------------------------------------------------------------*/
604static hdrl_resample_smallwcs *
605hdrl_resample_smallwcs_new(cpl_propertylist *aHeader)
606{
607
608 cpl_ensure(aHeader, CPL_ERROR_NULL_INPUT, NULL);
609 hdrl_resample_smallwcs *wcs = cpl_calloc(1, sizeof(hdrl_resample_smallwcs));
610
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)) {
617 /* all these headers default to 0.0 following the FITS *
618 * Standard, so we can ignore any errors set here */
619 cpl_errorstate_set(prestate);
620 }
621
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.) {
629 /* FITS Standard says to handle the CD matrix like the PC *
630 * matrix in this case, with 1 for the diagonal elements */
631 wcs->cd11 = wcs->cd22 = wcs->cddet = 1.;
632 cpl_errorstate_set(prestate); /* not a real error */
633 }
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);
637 }
638
639 return wcs;
640} /* hdrl_resample_smallwcs_new() */
641
642/*---------------------------------------------------------------------------*/
648/*---------------------------------------------------------------------------*/
649static void
650hdrl_resample_pixgrid_delete(hdrl_resample_pixgrid *aGrid)
651{
652 if (!aGrid) {
653 return;
654 }
655 cpl_free(aGrid->pix);
656 aGrid->pix = NULL;
657 /* clean up extension maps */
658 unsigned short ix;
659 for (ix = 0; ix < aGrid->nmaps; ix++) {
660 cpl_size iext; /* index in this extension map */
661 for (iext = 0; iext < aGrid->nxalloc[ix]; iext++) {
662 cpl_free(aGrid->xmaps[ix][iext].pix);
663 } /* for iext (all allocated pixels in this extension map) */
664
665 cpl_free(aGrid->xmaps[ix]);
666 } /* for ix (all extension maps) */
667 cpl_free(aGrid->xmaps);
668 aGrid->xmaps = NULL;
669 cpl_free(aGrid->nxmap);
670 aGrid->nxmap = NULL;
671 cpl_free(aGrid->nxalloc);
672 aGrid->nxalloc = NULL;
673 cpl_free(aGrid);
674} /* hdrl_resample_pixgrid_delete() */
675
676
677/*---------------------------------------------------------------------------*/
733/*---------------------------------------------------------------------------*/
734hdrl_resample_result *
735hdrl_resample_compute(const cpl_table *ResTable,
736 hdrl_parameter *method,
737 hdrl_parameter *outputgrid,
738 const cpl_wcs* wcs)
739{
740 cpl_ensure(ResTable && method, CPL_ERROR_NULL_INPUT, NULL);
741 cpl_ensure(wcs != NULL, CPL_ERROR_NULL_INPUT, NULL);
742
743 if (hdrl_resample_inputtable_verify(ResTable)) {
744 return NULL;
745 }
746 if (hdrl_resample_parameter_method_verify((hdrl_parameter *)method)) {
747 return NULL;
748 }
749 if (hdrl_resample_parameter_outgrid_verify((hdrl_parameter *)outputgrid)) {
750 return NULL;
751 }
752
753 if (hdrl_resample_inputtable_verify(ResTable)) {
754 return NULL;
755 }
756
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;
761
762 /*Assign the wcs*/
763 aParams_outputgrid->wcs = wcs;
764
765 /* Recalculate the limits if the user did not specify any */
766 cpl_boolean recalc_limits = aParams_outputgrid->recalc_limits;
767
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);
775
776 /* We have the rare case that the image spans over ra = 0.*/
777 if(ramax - ramin > 180){
778 const double *ra = cpl_table_get_data_double_const(ResTable,
779 HDRL_RESAMPLE_TABLE_RA);
780
781 /* set to extreme values for a start */
782 ramin = 0.;
783 ramax = 360.;
784 cpl_size nrow = cpl_table_get_nrow(ResTable);
785
786 for (cpl_size i = 0; i < nrow; i++) {
787 if (ra[i] > ramin && ra[i] <= 180.) ramin = ra[i]; /* get the maximum */
788 if (ra[i] < ramax && ra[i] > 180.) ramax = ra[i]; /* get the minimum */
789 }
790 }
791
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;
798 }
799
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);
803
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);
807
808 /* create cube and cast to generic pointer to save code duplication */
809 hdrl_resample_result *cube = NULL;
810 hdrl_resample_pixgrid *grid = NULL;
811
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,
817 &grid);
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);
822
823 /* Cleanup wcs for 2D / 3D case */
824
825 if(hdrl_imagelist_get_size(cube->himlist) == 1) { /* 2D case */
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;
834 }
835
836
837 hdrl_resample_pixgrid_delete(grid);
838 return cube;
839} /* hdrl_resample_compute() */
840
841/*---------------------------------------------------------------------------*/
847/*---------------------------------------------------------------------------*/
848void
849hdrl_resample_result_delete(hdrl_resample_result *aCube)
850{
851 /* if the cube does not exists at all, we don't need to do anything */
852 if (!aCube) {
853 return;
854 }
855
856 /* checks for the existence of the sub-images *
857 * are done in the CPL functions */
858
859 hdrl_imagelist_delete(aCube->himlist);
860 aCube->himlist = NULL;
861 /* delete the FITS header, too */
862 cpl_propertylist_delete(aCube->header);
863 aCube->header = NULL;
864
865 cpl_free(aCube);
866} /* hdrl_resample_result_delete() */
867
868/*----------------------------------------------------------------------------*/
895/*----------------------------------------------------------------------------*/
896static cpl_error_code
897hdrl_resample_wcs_projplane_from_celestial(hdrl_resample_outgrid_parameter *
898 aParams_outputgrid, double aRA,
899 double aDEC, double *aX, double *aY)
900{
901 cpl_ensure_code(aParams_outputgrid && aX && aY, CPL_ERROR_NULL_INPUT);
902 /* make sure that the header represents TAN */
903 //const char * type1 = aParams_outputgrid->limits.ctype1;
904 //const char * type2 = aParams_outputgrid->limits.ctype2;
905 /*
906 cpl_ensure_code(type1 && type2 && !strncmp(type1, "RA---TAN", 9) &&
907 !strncmp(type2, "DEC--TAN", 9),
908 CPL_ERROR_UNSUPPORTED_MODE);
909 */
910
911 int err = 0;
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);
915
916 /* spherical coordinate shift/translation */
917 double a = aRA / CPL_MATH_DEG_RAD, /* RA in radians */
918 d = aDEC / CPL_MATH_DEG_RAD, /* DEC in radians */
919 /* alpha_p and delta_p in Paper II (in radians) */
920 //ap = aParams_outputgrid->limits.crval1 / CPL_MATH_DEG_RAD,
921 //dp = aParams_outputgrid->limits.crval2 / 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);
929 /* spherical deprojection */
930 *aX = R_theta * sin(phi),
931 *aY = -R_theta * cos(phi);
932
933 return CPL_ERROR_NONE;
934} /* hdrl_resample_wcs_projplane_from_celestial() */
935
936/*---------------------------------------------------------------------------*/
954/*---------------------------------------------------------------------------*/
955static inline void
956hdrl_resample_wcs_pixel_from_celestial_fast(hdrl_resample_smallwcs *aWCS,
957 double aRA, double aDEC, double *aX, double *aY)
958{
959 if(!aWCS) return;
960 /* spherical coordinate shift/translation */
961 /*Calabretta & Greisen 2002 A&A 395, 1077 (Paper II) */
962 /* aRA is alpha in Paper II, aDEC is delta, aWCS->crval1 is alpha_p, *
963 * aWCS->crval2 is delta_p, all of them in units of radians, eq (5), arg=atan2 */
964
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);
972 /* spherical deprojection */
973 double x = R_theta * sin(phi),
974 y = -R_theta * cos(phi);
975 /* inverse linear transformation */
976 *aX = (aWCS->cd22 * x - aWCS->cd12 * y) / aWCS->cddet + aWCS->crpix1;
977 *aY = (aWCS->cd11 * y - aWCS->cd21 * x) / aWCS->cddet + aWCS->crpix2;
978} /* hdrl_resample_wcs_pixel_from_celestial_fast() */
979
980
981/*---------------------------------------------------------------------------*/
994/*---------------------------------------------------------------------------*/
995static cpl_error_code
996hdrl_resample_compute_size(hdrl_resample_outgrid_parameter *aParams_outputgrid,
997 int *aX, int *aY, int *aZ)
998{
999 cpl_ensure_code(aParams_outputgrid && aX && aY && aZ, CPL_ERROR_NULL_INPUT);
1000 const char func[] = "hdrl_resample_compute_size"; /* pretend to be in that fct... */
1001 double x1, y1, x2, y2;
1002
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,
1008 &x1, &y1);
1009 hdrl_resample_wcs_projplane_from_celestial(aParams_outputgrid, ramax, decmax,
1010 &x2, &y2);
1011 *aX = lround(fabs(x2 - x1) / aParams_outputgrid->delta_ra) + 1;
1012 *aY = lround(fabs(y2 - y1) / aParams_outputgrid->delta_dec) + 1;
1013
1014 double lmin = aParams_outputgrid->lambda_min;
1015 double lmax = aParams_outputgrid->lambda_max;
1016
1017 *aZ = (int)ceil((lmax - lmin) / aParams_outputgrid->delta_lambda) + 1;
1018
1019 cpl_msg_debug(func, "Output cube size %d x %d x %d (fit to data)",
1020 *aX, *aY, *aZ);
1021 return CPL_ERROR_NONE;
1022} /* hdrl_resample_compute_size() */
1023
1024/*---------------------------------------------------------------------------*/
1040/*---------------------------------------------------------------------------*/
1041static void
1042hdrl_resample_pixgrid_add(hdrl_resample_pixgrid *aGrid, cpl_size aIndex,
1043 cpl_size aRow, unsigned short aXIdx)
1044{
1045
1046 if (aIndex < 0 || aGrid == NULL) {
1047 return;
1048 }
1049
1050 if (aGrid->pix[aIndex] == 0 && aRow > 0) {
1051 /* First pixel is stored directly. */
1052 aGrid->pix[aIndex] = aRow;
1053 } else if (aGrid->pix[aIndex] == 0 && aRow == 0) {
1054 /* Special case: we cannot put "0" into the main map. */
1055 cpl_size iext = aGrid->nxmap[aXIdx]++;
1056 if (aGrid->nxmap[aXIdx] > aGrid->nxalloc[aXIdx]) {
1057 /* double the number of allocated entries */
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));
1062 }
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) {
1068 /* When a second pixel is added, put both to the extension map. */
1069 cpl_size iext = aGrid->nxmap[aXIdx]++;
1070 if (aGrid->nxmap[aXIdx] > aGrid->nxalloc[aXIdx]) {
1071 /* double the number of allocated entries */
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));
1076 }
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));
1082 } else {
1083 /* Append additional pixels to the extension map. */
1084 cpl_size iext = (-aGrid->pix[aIndex] - 1) & PT_IDX_MASK;
1085 /* index of the new entry in this grid point */
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;
1091 }
1092} /* hdrl_resample_pixgrid_add() */
1093
1094/*---------------------------------------------------------------------------*/
1103/*---------------------------------------------------------------------------*/
1104static inline cpl_size
1105hdrl_resample_pixgrid_get_count(hdrl_resample_pixgrid *aGrid, cpl_size aIndex)
1106{
1107 if (aIndex < 0 || aGrid == NULL) {
1108 return 0;
1109 }
1110 /* get entry in pixel grid */
1111 cpl_size p = aGrid->pix[aIndex];
1112 if (p == 0) { /* points nowhere --> no pixels */
1113 return 0;
1114 }
1115 if (p > 0) { /* points to pixel table --> 1 pixel */
1116 return 1;
1117 }
1118 /* p is negative, so points to an extension map, get its index */
1119 unsigned short ix = (-p >> XMAP_LSHIFT) & XMAP_BITMASK;
1120 p = (-p - 1) & PT_IDX_MASK;
1121 /* the npix component now gives the number of pixels in this grid index */
1122 return aGrid->xmaps[ix][p].npix;
1123} /* hdrl_resample_pixgrid_get_count() */
1124
1125/*---------------------------------------------------------------------------*/
1139/*---------------------------------------------------------------------------*/
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)
1143{
1144 if (aGrid == NULL) {
1145 return -1;
1146 }
1147 if (!aAllowOutside &&
1148 (aX < 0 || aX >= aGrid->nx || aY < 0 || aY >= aGrid->ny ||
1149 aZ < 0 || aZ >= aGrid->nz)) {
1150 return -1;
1151 }
1152 if (aX < 0) {
1153 aX = 0;
1154 }
1155 if (aX >= aGrid->nx) {
1156 aX = aGrid->nx - 1;
1157 }
1158 if (aY < 0) {
1159 aY = 0;
1160 }
1161 if (aY >= aGrid->ny) {
1162 aY = aGrid->ny - 1;
1163 }
1164 if (aZ < 0) {
1165 aZ = 0;
1166 }
1167 if (aZ >= aGrid->nz) {
1168 aZ = aGrid->nz - 1;
1169 }
1170 return aX + aGrid->nx * (aY + aGrid->ny * aZ);
1171} /* hdrl_resample_pixgrid_get_index() */
1172
1173/*---------------------------------------------------------------------------*/
1183/*---------------------------------------------------------------------------*/
1184static hdrl_resample_pixgrid *
1185hdrl_resample_pixgrid_new(cpl_size aSizeX, cpl_size aSizeY, cpl_size aSizeZ,
1186 unsigned short aNMaps)
1187{
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));
1195 /* extension maps for possibly multiple threads */
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));
1200 return pixels;
1201} /* hdrl_resample_pixgrid_new() */
1202
1203/*---------------------------------------------------------------------------*/
1236/*---------------------------------------------------------------------------*/
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)
1240{
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);
1244 if (nrow == 0) {
1245 cpl_msg_error(cpl_func, "Invalid pixel table (no entries?)");
1246 cpl_error_set(cpl_func, CPL_ERROR_NULL_INPUT);
1247 return NULL;
1248 }
1249 cpl_ensure(aXSize > 0 && aYSize > 0 && aZSize > 0, CPL_ERROR_ILLEGAL_INPUT,
1250 NULL);
1251
1252
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);
1256
1257 hdrl_resample_smallwcs *wcs = hdrl_resample_smallwcs_new(aHeader);
1258
1259 /* get all (relevant) table columns for easy pointer access */
1260
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);
1269 cpl_free(wcs);
1270 return NULL;
1271 }
1272
1273 wcs->crval1 /= CPL_MATH_DEG_RAD; /* convert to radians before calling... */
1274 wcs->crval2 /= CPL_MATH_DEG_RAD; /* ...hdrl_resample_wcs_pixel_from_celestial_fast() */
1275
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;
1281
1282 /* check for the selected pixels in the pixel table, only those *
1283 * are used to construct the pixel grid; since constructing the *
1284 * array of selected pixels costs significant amounts of time, *
1285 * do that only when not all pixels are selected! */
1286 cpl_array *asel = NULL;
1287 const cpl_size *sel = NULL;
1288 cpl_size nsel = cpl_table_count_selected(ResTable);
1289 if (nsel < nrow) {
1290 asel = cpl_table_where_selected(ResTable);
1291 sel = cpl_array_get_data_cplsize_const(asel);
1292 }
1293
1294 /* can use at most XMAP_BITMASK threads so that the bitmask does not *
1295 * overflow, but ensure that we are not using more cores than available... */
1296 int nth = omp_get_max_threads() > XMAP_BITMASK ? XMAP_BITMASK
1297 : omp_get_max_threads();
1298 /* prepare the ranges for the different threads, store them in arrays */
1299 cpl_array *az1 = cpl_array_new(nth, CPL_TYPE_INT),
1300 *az2 = cpl_array_new(nth, CPL_TYPE_INT);
1301 if (aZSize < nth) {
1302 /* pre-fill arrays with values that cause the threads to do nothing */
1303 cpl_array_fill_window_int(az1, aZSize, nth, -1);
1304 cpl_array_fill_window_int(az2, aZSize, nth, -2);
1305 }
1306 /* now fill the (first) ones with real ranges */
1307 double base = nth > aZSize ? 1. : (double)aZSize / nth;
1308 cpl_size ith;
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));
1312 } /* for */
1313 /* make sure that we don't lose pixels at the edges of the wavelength *
1314 * range, put them into the extreme threads by making their ranges larger; *
1315 * set the relevant array entries to something close to the largest value, *
1316 * that we can still add as an integer (to compute the z-range) */
1317 cpl_array_set_int(az1, 0, -INT_MAX / 2 + 1);
1318 cpl_array_set_int(az2, ith - 1, INT_MAX / 2 - 1);
1319
1320 /* create the pixel grid with extension maps for threads */
1321 hdrl_resample_pixgrid *grid = hdrl_resample_pixgrid_new(aXSize, aYSize, aZSize, nth);
1322
1323 /* parallel region to fill the pixel grid */
1324
1325 struct timeval tv1, tv2;
1326
1327 cpl_msg_debug(cpl_func,"Starting parallel loop in hdrl_resample_pixgrid_create");
1328 gettimeofday(&tv1, NULL);
1329
1330#pragma omp parallel num_threads(nth) /* default(none) as req. by Ralf */ \
1331 shared(aXSize, aYSize, aZSize, az1, az2, cd33, crpix3, crval3, grid, \
1332 lbda, nsel, sel, showprogress, \
1333 timeinit, timeprogress, wcs, xpos, ypos)
1334
1335 {
1336 /* split the work up into threads, for non-overlapping wavelength ranges */
1337 unsigned short ithread = omp_get_thread_num(); /* index of this thread */
1338 int z1 = cpl_array_get_int(az1, ithread, NULL),
1339 z2 = cpl_array_get_int(az2, ithread, NULL),
1340 zrange = z2 - z1 + 1;
1341
1342 /* check if we actually need to enter the (parallel) loop, i.e. *
1343 * if the current thread is handling any wavelength planes */
1344
1345 /* now the actual parallel loop */
1346 cpl_size isel;
1347 for (isel = 0 ; zrange > 0 && isel < nsel; isel++) {
1348#pragma omp master /* only output progress from the master thread */
1349 if (showprogress && !((isel+1) % 1000000ll)) { /* output before every millionth entry */
1350 double timenow = cpl_test_get_walltime();
1351 if (timenow - timeprogress > 30.) { /* and more than half a minute passed */
1352 timeprogress = timenow;
1353 double percent = 100. * (isel + 1.) / nsel,
1354 elapsed = timeprogress - timeinit,
1355 remaining = (100. - percent) * elapsed / percent;
1356 /* overwritable only exists for INFO mode, but we check *
1357 * above that we want this only for DEBUG mode output... */
1358 cpl_msg_info_overwritable(cpl_func, "pixel grid creation is %.1f%% "
1359 "complete, %.1fs elapsed, ~%.1fs remaining",
1360 percent, elapsed, remaining);
1361 } /* if: 1/2 min passed */
1362 } /* if: want debug output */
1363
1364 /* either use the index from the array of selected rows *
1365 * or the row numberdirectly (for a fully selected table) */
1366 cpl_size n = sel ? sel[isel] : isel;
1367 int z = -1;
1368
1369 z = lround((lbda[n] - crval3) / cd33 + crpix3) - 1;
1370
1371 if (z < z1 || z > z2) {
1372 continue; /* skip this entry, one of the other threads handles it */
1373 }
1374
1375 /* determine the pixel coordinates in the grid (indices, starting at 0) */
1376 double xpx = 0., ypx = 0.;
1377
1378 hdrl_resample_wcs_pixel_from_celestial_fast(wcs, (xpos[n]) / CPL_MATH_DEG_RAD,
1379 (ypos[n]) / CPL_MATH_DEG_RAD,
1380 &xpx, &ypx);
1381
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);
1385
1386
1387 /* write the pixel values to the correct place in the grid */
1388 hdrl_resample_pixgrid_add(grid, idx, n, ithread);
1389 } /* for isel (all selected pixel table rows) */
1390
1391 /* Clean up the possibly too many allocations; this is not strictly *
1392 * needed but nice to only consume as much memory as we need. */
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];
1397 } /* omp parallel */
1398
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));
1403
1404 cpl_array_delete(asel);
1405 cpl_free(wcs);
1406 cpl_array_delete(az1);
1407 cpl_array_delete(az2);
1408
1409
1410 cpl_size idx, npix = 0;
1411 for (idx = 0; idx < aXSize * aYSize * aZSize; idx++) {
1412 npix += hdrl_resample_pixgrid_get_count(grid, idx);
1413 }
1414 cpl_size nxmap = 0;
1415 for (idx = 0; idx < (cpl_size)grid->nmaps; idx++) {
1416 nxmap += grid->nxmap[idx];
1417 }
1418 if (npix != nsel) {
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);
1424 cpl_free(msg);
1425 }
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,
1433 cpufini - cpuinit);
1434
1435 return grid;
1436} /* hdrl_resample_pixgrid_create() */
1437
1438/*----------------------------------------------------------------------------*/
1458/*----------------------------------------------------------------------------*/
1459static cpl_error_code
1460hdrl_resample_wcs_get_scales(hdrl_resample_outgrid_parameter *aParams_outputgrid,
1461 double *aXScale, double *aYScale)
1462{
1463 cpl_ensure_code(aParams_outputgrid && aXScale && aYScale, CPL_ERROR_NULL_INPUT);
1464
1465 cpl_errorstate prestate = cpl_errorstate_get();
1466
1467 const cpl_matrix *cd = cpl_wcs_get_cd(aParams_outputgrid->wcs);
1468 //double cd11 = aParams_outputgrid->limits.cd11;
1469 //double cd22 = aParams_outputgrid->limits.cd22;
1470 //double cd12 = aParams_outputgrid->limits.cd12;
1471 //double cd21 = aParams_outputgrid->limits.cd21;
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);
1476
1477
1478 double det = cd11 * cd22 - cd12 * cd21;
1479 cpl_ensure_code(cpl_errorstate_is_equal(prestate), cpl_error_get_code());
1480
1481 if (det < 0.) {
1482 cd12 *= -1;
1483 cd11 *= -1;
1484 }
1485 if (cd12 == 0. && cd21 == 0.) { /* matrix without rotation */
1486 *aXScale = cd11;
1487 *aYScale = cd22;
1488 return CPL_ERROR_NONE;
1489 }
1490 *aXScale = sqrt(cd11*cd11 + cd12*cd12); /* only the absolute value */
1491 *aYScale = sqrt(cd22*cd22 + cd21*cd21);
1492 return CPL_ERROR_NONE;
1493} /* hdrl_resample_wcs_get_scales() */
1494/*---------------------------------------------------------------------------*/
1503/*---------------------------------------------------------------------------*/
1504static inline const cpl_size *
1505hdrl_resample_pixgrid_get_rows(hdrl_resample_pixgrid *aGrid, cpl_size aIndex)
1506{
1507
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);
1511 /* get entry in pixel grid */
1512 cpl_size p = aGrid->pix[aIndex];
1513 if (p == 0) { /* points nowhere --> no array */
1514 return NULL;
1515 }
1516 if (p > 0) { /* points to pixel table */
1517 return aGrid->pix + aIndex;
1518 }
1519 /* p is negative, so points to an extension map, get its array */
1520 unsigned short ix = (-p >> XMAP_LSHIFT) & XMAP_BITMASK;
1521 p = (-p - 1) & PT_IDX_MASK;
1522 /* the pix component provides the array of pixel table rows in this index */
1523 return aGrid->xmaps[ix][p].pix;
1524} /* hdrl_resample_pixgrid_get_rows() */
1525
1526
1527/*---------------------------------------------------------------------------*/
1544/*---------------------------------------------------------------------------*/
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)
1549{
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);
1558
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);
1563
1564 cpl_wcs *wcscpl = cpl_wcs_new_from_propertylist(aCube->header);
1565
1566 if (hdrl_resample_inputtable_verify(ResTable)) {
1567 return CPL_ERROR_ILLEGAL_INPUT;
1568 }
1569
1570
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);
1577
1578 /* If our data was astrometrically calibrated, we need to scale the *
1579 * data units to the pixel size in all three dimensions so that the *
1580 * radius computation works again. *
1581 * Otherwise dx~5.6e-5deg won't contribute to the weighting at all. */
1582
1583 double xnorm = 1., ynorm = 1., znorm = 1. ;
1584
1585
1586 hdrl_resample_wcs_get_scales(aParams_outputgrid, &xnorm, &ynorm);
1587 //TODO: we should check that xnorm, ynorm, znorm are not zero
1588 xnorm = 1. / xnorm;
1589 ynorm = 1. / ynorm;
1590 //znorm = 1. / aParams_outputgrid->limits.cd33;
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);
1594 }
1595/*
1596 cpl_msg_debug(cpl_func, "Norm factors from wcs - xnorm: %g ynorm: %g znorm: %g",
1597 xnorm, ynorm, znorm);
1598 need to use the real coordinate offset for celestial spherical
1599 // We are now working with the full astrometric solution
1600 //ptxoff = hdrl_resample_pfits_get_crval(ResTable->header, 1);
1601 //ptyoff = hdrl_resample_pfits_get_crval(ResTable->header, 2);
1602*/
1603
1604 struct timeval tv1, tv2;
1605 cpl_msg_debug(cpl_func,"Starting parallel loop in hdrl_resample_cube_nearest");
1606 gettimeofday(&tv1, NULL);
1607
1608#pragma omp parallel for default(none) /* as req. by Ralf */ \
1609 shared(aCube, aGrid, cd33, crpix3, crval3, data, dq, lbda, \
1610 stat, stdout, wcscpl, \
1611 xnorm, xpos, ynorm, ypos, znorm) collapse(2)
1612
1613 for (cpl_size l = 0; l < aGrid->nz; l++) {
1614 for (cpl_size i = 0; i < aGrid->nx; i++) {
1615 double *pdata = cpl_image_get_data_double(hdrl_image_get_image(hdrl_imagelist_get(aCube->himlist, l)));
1616 double *pstat = cpl_image_get_data_double(hdrl_image_get_error(hdrl_imagelist_get(aCube->himlist, l)));
1617 cpl_binary *pdq = cpl_mask_get_data(hdrl_image_get_mask(hdrl_imagelist_get(aCube->himlist, l)));
1618 /* wavelength of center of current grid cell (l is index starting at 0) */
1619 double lambda = (l + 1. - crpix3) * cd33 + crval3;
1620
1621 cpl_size j;
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);
1626
1627 /* x and y position of center of current grid cell (i, j start at 0) */
1628 double x = 0., y = 0.;
1629
1630 // We are now working with the full astrometric solution
1631 // hdrl_resample_wcs_celestial_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
1632 hdrl_wcs_xy_to_radec(wcscpl, i + 1, j + 1, &x, &y);
1633
1634 if (n_rows == 1) {
1635 if ((cpl_binary)dq[rows[0]] == CPL_BINARY_0) {
1636 /* if there is only one pixel in the cell, just use it */
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]];
1640 } else {
1641 pdq[i + j * aGrid->nx] = CPL_BINARY_1;
1642 }
1643 } else if (n_rows >= 2) {
1644 /* loop through all available values and take the closest one */
1645 cpl_size n, nbest = -1;
1646 double dbest = FLT_MAX; /* some unlikely large value to start with*/
1647 for (n = 0; n < n_rows; n++) {
1648 /* do not use bad pixels */
1649 if((cpl_binary)dq[rows[n]] != CPL_BINARY_0) {
1650 continue;
1651 }
1652 /* the differences for the cell center and the current pixel */
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);
1657
1658 /* Not strictly necessary for NN, but still scale the RA *
1659 * distance properly, see hdrl_resample_cube_weighted(). */
1660 dx *= cos(y * CPL_MATH_RAD_DEG);
1661
1662 if (dthis < dbest) {
1663 nbest = n;
1664 dbest = dthis;
1665 }
1666 }
1667 if (nbest >= 0) { /* We found a good nearest neighbor */
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]];
1671 }
1672 } else {
1673 /* npix == 0: do nothing, pixel stays zero */
1674 pdq[i + j * aGrid->nx] = CPL_BINARY_1;
1675 }
1676 } /* for j (y direction) */
1677 } /* for i (x direction) */
1678 } /* for l (wavelength planes) */
1679
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));
1684
1685 /* Make sure that the bpm of the image and the error are in sinc as we are
1686 * working with pointers */
1687 cpl_size size = hdrl_imagelist_get_size(aCube->himlist);
1688 for(cpl_size i = 0; i < size; i++){
1689 /* sync image and error bpm ignoring what is in error before */
1690 cpl_image_reject_from_mask(hdrl_image_get_error(hdrl_imagelist_get(aCube->himlist, i)),
1691 hdrl_image_get_mask(hdrl_imagelist_get(aCube->himlist, i) ));
1692 }
1693
1694 cpl_wcs_delete(wcscpl);
1695
1696 return CPL_ERROR_NONE;
1697} /* hdrl_resample_cube_nearest() */
1698
1699/*---------------------------------------------------------------------------*/
1712/*---------------------------------------------------------------------------*/
1713static inline double
1714hdrl_resample_weight_function_renka(double r, double r_c)
1715{
1716 if (r == 0) {
1717 return FLT_MAX;
1718 } else if (r >= r_c) {
1719 return DBL_MIN;
1720 } else {
1721 double p = (r_c - r) / (r_c * r);
1722 return p*p;
1723 }
1724} /* hdrl_resample_weight_function_renka() */
1725
1726/*---------------------------------------------------------------------------*/
1745/*---------------------------------------------------------------------------*/
1746static inline double
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)
1750{
1751 /* compute the three terms in the numerator: if the offset *
1752 * plus the output halfsize is less than the input halfsize, *
1753 * then that side is fully contained in the input pixel */
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;
1757 /* any negative value means that the input pixel is completely *
1758 * outside the target voxel, so it doesn't contribute */
1759 if (x <= 0 || y <= 0 || z <= 0) {
1760 return 0.;
1761 }
1762 /* any value > input size means this dimension of the input pixel *
1763 * is completely inside the target voxel, so that's the limit! */
1764 return (x > xin ? xin : x) * (y > yin ? yin : y) * (z > zin ? zin : z)
1765 / (xin * yin * zin);
1766} /* hdrl_resample_weight_function_drizzle() */
1767
1768
1769/*---------------------------------------------------------------------------*/
1780/*---------------------------------------------------------------------------*/
1781static inline double
1782hdrl_resample_weight_function_linear(double r)
1783{
1784 return r == 0 ? FLT_MAX : 1. / r;
1785} /* hdrl_resample_weight_function_linear() */
1786
1787/*---------------------------------------------------------------------------*/
1798/*---------------------------------------------------------------------------*/
1799static inline double
1800hdrl_resample_weight_function_quadratic(double r2)
1801{
1802 return r2 == 0 ? FLT_MAX : 1. / r2;
1803} /* hdrl_resample_weight_function_quadratic() */
1804
1805/*---------------------------------------------------------------------------*/
1812/*---------------------------------------------------------------------------*/
1813static inline double
1814hdrl_resample_weight_function_sinc(double r)
1815{
1816 return fabs(r) < DBL_EPSILON ? 1. : sin(CPL_MATH_PI * r) / (CPL_MATH_PI * r);
1817} /* hdrl_resample_weight_function_sinc() */
1818
1819/*---------------------------------------------------------------------------*/
1830/*---------------------------------------------------------------------------*/
1831static inline double
1832hdrl_resample_weight_function_lanczos(double dx, double dy, double dz,
1833 unsigned int ld, unsigned int lks)
1834{
1835 /* Adding 0.5 as for a loop distance of 0 the weight should only drop to 0 if
1836 * the distance is larger then half the pixel */
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);
1841} /* hdrl_resample_weight_function_lanczos() */
1842
1843/*---------------------------------------------------------------------------*/
1865/*---------------------------------------------------------------------------*/
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)
1872{
1873
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);
1882
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);
1886
1887 hdrl_resample_smallwcs *wcs = hdrl_resample_smallwcs_new(aCube->header);
1888 cpl_wcs *wcscpl = cpl_wcs_new_from_propertylist(aCube->header);
1889
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);
1896
1897 /* If our data was astrometrically calibrated, we need to scale the *
1898 * data units to the pixel size in all three dimensions so that the *
1899 * radius computation works again. *
1900 * Otherwise dx~5.6e-5deg won't contribute to the weighting at all. */
1901
1902 double xnorm = 1., ynorm = 1., znorm = 1. ;
1903
1904 hdrl_resample_wcs_get_scales(aParams_outputgrid, &xnorm, &ynorm);
1905 xnorm = 1. / xnorm;
1906 ynorm = 1. / ynorm;
1907 //znorm = 1. / aParams_outputgrid->limits.cd33;
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);
1911 }
1912/*
1913 cpl_msg_debug(cpl_func, "Norm factors from wcs - xnorm: %g ynorm: %g znorm: %g",
1914 xnorm, ynorm, znorm);
1915 need to use the real coordinate offset for celestial spherical
1916 // We are now working with the full astrometric solution
1917 //ptxoff = hdrl_resample_pfits_get_crval(ResTable->header, 1);
1918 //ptyoff = hdrl_resample_pfits_get_crval(ResTable->header, 2);
1919*/
1920
1921
1922 /* scale the input critical radius by the voxel radius */
1923 double renka_rc = aParams_method->renka_critical_radius /*beware of rotation! */
1924 * sqrt((wcs->cd11*xnorm)*(wcs->cd11*xnorm)
1925 + (wcs->cd22*ynorm)*(wcs->cd22*ynorm)
1926 + (cd33*znorm)*(cd33*znorm));
1927
1928 /* loop distance (to take into account surrounding pixels) verification */
1929 int ld = aParams_method->loop_distance;
1930 if (ld < 0) {
1931 ld = 0;
1932 cpl_msg_debug(cpl_func, "Overriding loop distance ld=%d", ld);
1933 }
1934
1935 /* Lankzos kernel size (lks) verification */
1936 int lks = aParams_method->lanczos_kernel_size;
1937 if (lks <= 0) {
1938 lks = 1;
1939 cpl_msg_debug(cpl_func, "Overriding lanczos kernel size lks=%d", lks);
1940 }
1941
1942 /* Should 1/variance be used as an additional weight */
1943 cpl_boolean wght = aParams_method->use_errorweights;
1944
1945 /* pixel sizes in all three directions, scaled by pixfrac, and *
1946 * output pixel sizes (absolute values), as needed for drizzle */
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);
1951
1952 struct timeval tv1, tv2;
1953 cpl_msg_debug(cpl_func,"Starting parallel loop in hdrl_resample_cube_weighted");
1954 gettimeofday(&tv1, NULL);
1955
1956#pragma omp parallel for default(none) /* as req. by Ralf */ \
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)
1961
1962 for (cpl_size l = 0; l < aGrid->nz; l++) {
1963 for (cpl_size i = 0; i < aGrid->nx; i++) {
1964 double *pdata = cpl_image_get_data_double(hdrl_image_get_image(hdrl_imagelist_get(aCube->himlist, l)));
1965 double *pstat = cpl_image_get_data_double(hdrl_image_get_error(hdrl_imagelist_get(aCube->himlist, l)));
1966 cpl_binary *pdq = cpl_mask_get_data(hdrl_image_get_mask(hdrl_imagelist_get(aCube->himlist, l)));
1967
1968 /* wavelength of center of current grid cell (l is index starting at 0) */
1969 double lambda = (l + 1. - crpix3) * cd33 + crval3;
1970 double zout2 = zout; /* correct the output pixel size for log-lambda */
1971
1972 for (cpl_size j = 0; j < aGrid->ny; j++) {
1973 /* x and y position of center of current grid cell (i, j start at 0) */
1974 double x, y;
1975
1976 // We are now working with the full astrometric solution
1977 // hdrl_resample_wcs_celestial_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
1978 hdrl_wcs_xy_to_radec(wcscpl, i + 1, j + 1, &x, &y);
1979 //cpl_msg_debug(cpl_func,"case1 i %d j %d x %10.7f y %10.7f ",i+1,j+1,x,y);
1980
1981 double sumdata = 0, sumstat = 0, sumweight = 0;
1982 double flux=0;
1983 cpl_size npoints = 0;
1984
1985 /* loop through surrounding cells and their contained pixels */
1986 cpl_size i2;
1987 for (i2 = i - ld; i2 <= i + ld; i2++) {
1988 cpl_size j2;
1989 for (j2 = j - ld; j2 <= j + ld; j2++) {
1990 cpl_size l2;
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);
1995
1996 const cpl_size *rows2 = hdrl_resample_pixgrid_get_rows(aGrid, idx2);
1997 cpl_size n;
1998 for (n = 0; n < n_rows2; n++) {
1999 if (dq[rows2[n]]) { /* exclude all bad pixels */
2000
2001 continue;
2002 }
2003
2004 double dx = fabs(x - (xpos[rows2[n]])),
2005 dy = fabs(y - (ypos[rows2[n]])),
2006 dlambda = fabs(lambda - lbda[rows2[n]]),
2007 r2 = 0;
2008
2009 /* Since the distances of RA in degrees get larger the *
2010 * closer we get to the celestial pole, we have to *
2011 * compensate for that by multiplying the distance in *
2012 * RA by cos(delta), to make it comparable to the *
2013 * distances in pixels for the different kernels below. */
2014
2015 // We are now working with the full astrometric solution
2016 dx *= cos(y * CPL_MATH_RAD_DEG);
2017
2018 if (aParams_method->method != HDRL_RESAMPLE_METHOD_DRIZZLE) {
2019 dx *= xnorm;
2020 dy *= ynorm;
2021 dlambda *= znorm;
2022 r2 = dx * dx + dy * dy + dlambda * dlambda;
2023 //cpl_msg_debug(cpl_func,"n %lld dx %10.7f dy %10.7f dlambda %10.7f",n,dx,dy,dlambda);
2024 }
2025 double weight = 0.;
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,
2030 xout, yout, zout2,
2031 dx, dy, dlambda);
2032 } else if (aParams_method->method == HDRL_RESAMPLE_METHOD_LINEAR) {
2033
2034 weight = hdrl_resample_weight_function_linear(sqrt(r2));
2035 //cpl_msg_debug(cpl_func,"r2=%10.7f weight=%10.7f",r2,weight);
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,
2040 dlambda, ld, lks);
2041 }
2042
2043 if (wght && stat[rows2[n]] > 0.) { /* User wants to weight by 1/variance */
2044 /* apply it on top of the weight computed here */
2045 weight /= (stat[rows2[n]] * stat[rows2[n]]);
2046 }
2047
2048 sumweight += weight;
2049 sumdata += data[rows2[n]] * weight;
2050 flux += data[rows2[n]];
2051 sumstat += stat[rows2[n]] * stat[rows2[n]] * weight * weight;
2052 npoints++;
2053
2054 } /* for n (all pixels in grid cell) */
2055 } /* for l2 (lambda direction) */
2056 } /* for j2 (y direction) */
2057 } /* for i2 (x direction) */
2058
2059
2060 /* if no points were found, we cannot divide by the summed weight *
2061 * and don't need to set the output pixel value (it's 0 already), *
2062 * only set the relevant Euro3D bad pixel flag
2063 * In some cases only sumweight * sumweight is really zero so this
2064 * check was additionally added for the error propagation part*/
2065
2066 if (!npoints || !isnormal(sumweight) || !isnormal(sumweight * sumweight)) {
2067 pdq[i + j * aGrid->nx] = CPL_BINARY_1;
2068 continue;
2069 }
2070
2071 /* divide results by weight of summed pixels */
2072 sumdata /= sumweight;
2073 sumstat /= sumweight * sumweight;
2074
2075 pdata[i + j * aGrid->nx] = sumdata;
2076 /* Going back from variance to errors */
2077 pstat[i + j * aGrid->nx] = sqrt(sumstat);
2078 pdq[i + j * aGrid->nx] = CPL_BINARY_0; /* now we can mark it as good */
2079
2080 } /* for j (y direction) */
2081 } /* for i (x direction) */
2082 } /* for l (wavelength planes) */
2083
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));
2088
2089 /* Make sure that the bpm of the image and the error are in sinc as we are
2090 * working with pointers */
2091 cpl_size size = hdrl_imagelist_get_size(aCube->himlist);
2092 for(cpl_size i = 0; i < size; i++){
2093 /* sync image and error bpm ignoring what is in error before */
2094 cpl_image_reject_from_mask(hdrl_image_get_error(hdrl_imagelist_get(aCube->himlist, i)),
2095 hdrl_image_get_mask(hdrl_imagelist_get(aCube->himlist, i) ));
2096 }
2097
2098 cpl_free(wcs);
2099 cpl_wcs_delete(wcscpl);
2100
2101 return CPL_ERROR_NONE;
2102} /* hdrl_resample_cube_weighted() */
2103/*---------------------------------------------------------------------------*/
2122/*---------------------------------------------------------------------------*/
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)
2126{
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);
2130
2131 cpl_ensure_code(cube && aParams_outputgrid, CPL_ERROR_NULL_INPUT);
2132 /* TODO: cube->header may be already allocated. In this case why a new and
2133 * not an over-write of the content?
2134 */
2135 cube->header = cpl_propertylist_new ();
2136 hdrl_wcs_to_propertylist (aParams_outputgrid->wcs, cube->header,
2137 FALSE);
2138
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");
2143
2144
2145
2146 /* set NAXIS for later handling of the WCS */
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);
2151 /* if pixel table was astrometrically calibrated, use its WCS headers *
2152 * Axis 1: x or RA, axis 2: y or DEC, axis 3: lambda
2153 *
2154 * */
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.);
2161
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;
2166
2167 /* Following swarp we put CRPIX and CRVAL to the center of the field */
2168 cpl_propertylist_update_double(cube->header, "CRPIX1", (double)((xsize + 1) / 2.));
2169 cpl_propertylist_update_double(cube->header, "CRPIX2", (double)((ysize + 1) / 2.));
2170
2171 if(ramax - ramin < 180) {
2172 /* To be checked: Both values are in 0 - 180 or 180 - 360 */
2173 cpl_propertylist_update_double(cube->header, "CRVAL1", (ramin + ramax) / 2.);
2174 } else {
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);
2179 } else {
2180 cpl_propertylist_update_double(cube->header, "CRVAL1", ramax + (diff1 + diff2) / 2.);
2181 }
2182 }
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);
2189 /* fill in empty cross-terms of the CDi_j matrix */
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.);
2194
2195 return cpl_error_get_code();
2196}
2197
2198/*---------------------------------------------------------------------------*/
2252/*---------------------------------------------------------------------------*/
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)
2258{
2259 cpl_ensure(ResTable && aParams_method && aParams_outputgrid && aGrid,
2260 CPL_ERROR_NULL_INPUT, NULL);
2261
2262 /* compute or set the size of the output grid depending on *
2263 * the inputs and the data available in the pixel table */
2264
2265
2266 /* compute output sizes; wavelength is different in that it is *
2267 * more useful to contain partly empty areas within the field *
2268 * for the extreme ends, so use ceil() */
2269 int xsize = 0, ysize = 0, zsize = 0;
2270
2271 hdrl_resample_compute_size(aParams_outputgrid, &xsize, &ysize, &zsize);
2272
2273 /* Following swarp for x and y : Add a margin in field size */
2274
2275 xsize *= (100. + aParams_outputgrid->fieldmargin)/100.;
2276 ysize *= (100. + aParams_outputgrid->fieldmargin)/100.;
2277
2278 cpl_ensure(xsize > 0 && ysize > 0 && zsize > 0, CPL_ERROR_ILLEGAL_OUTPUT,
2279 NULL);
2280
2281 double time = cpl_test_get_walltime();
2282
2283 /* create the structure for the output datacube */
2284 hdrl_resample_result *cube = cpl_calloc(1, sizeof(hdrl_resample_result));
2285
2286 hdrl_resampling_set_outputgrid (xsize, ysize, zsize, cube, aParams_outputgrid);
2287
2288 if (aParams_method->method < HDRL_RESAMPLE_METHOD_NONE) {
2289 /* fill the cube for the data */
2290 cube->himlist = hdrl_imagelist_new();
2291
2292 for (cpl_size i = 0; i < zsize; i++) {
2293 hdrl_image * image = hdrl_image_new(xsize, ysize);
2294
2295 /* Set all pixels a priori to bad - do not use pointers to keep the
2296 * bpm of the data and error image in sync*/
2297 for (cpl_size j = 1; j <= xsize; j++) {
2298 for (cpl_size k = 1; k <= ysize; k++) {
2299 hdrl_image_reject(image, j, k);
2300 }
2301 }
2302
2303 hdrl_imagelist_set(cube->himlist, image, i);
2304 } /* for i (all wavelength planes) */
2305 } /* if method not HDRL_RESAMPLE_METHOD_NONE */
2306
2307 /* convert the pixel table into a pixel grid */
2308 hdrl_resample_pixgrid *grid = hdrl_resample_pixgrid_create(ResTable, cube->header,
2309 xsize, ysize, zsize);
2310 if (!grid) {
2312 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND, "Could not create"
2313 " pixel grid!");
2314 if (aGrid) {
2315 *aGrid = NULL;
2316 }
2317 return NULL;
2318 } /* if not grid */
2319
2320 double timeinit = cpl_test_get_walltime(),
2321 cpuinit = cpl_test_get_cputime();
2322
2323 /* do the resampling */
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);
2329 break;
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);
2337 break;
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 "
2342 "distance ld=%d)",
2343 aParams_method->method == HDRL_RESAMPLE_METHOD_LINEAR
2344 ? "linear"
2345 : (aParams_method->method == HDRL_RESAMPLE_METHOD_QUADRATIC
2346 ? "quadratic"
2347 : "lanczos"),
2348 aParams_method->loop_distance);
2349 rc = hdrl_resample_cube_weighted(cube, ResTable, grid, aParams_method,
2350 aParams_outputgrid);
2351 break;
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,
2360 aParams_method,
2361 aParams_outputgrid);
2362 break;
2363 case HDRL_RESAMPLE_METHOD_NONE:
2364 /* cpl_msg_debug(cpl_func, "Method %d (no resampling)", aParams_method->method); */
2365 break;
2366 default:
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;
2371 }
2372
2373 double timefini = cpl_test_get_walltime(),
2374 cpufini = cpl_test_get_cputime();
2375
2376 /* now that we have resampled we can either remove the pixel grid or save it */
2377 if (aGrid) {
2378 *aGrid = grid;
2379 } else {
2380 hdrl_resample_pixgrid_delete(grid);
2381 }
2382
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());
2390 return NULL;
2391 }
2392
2393 return cube;
2394} /* hdrl_resample_cube() */
2395
2396/*---------------------------------------------------------------------------*/
2406/*---------------------------------------------------------------------------*/
2407cpl_error_code
2408hdrl_wcs_to_propertylist(const cpl_wcs * wcs, cpl_propertylist * header,
2409 cpl_boolean only2d)
2410{
2411 cpl_ensure_code(wcs && header, CPL_ERROR_NULL_INPUT);
2412 int err = 0;
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);
2417
2418 const cpl_matrix *cd = cpl_wcs_get_cd(wcs);
2419
2420 const cpl_array *dims = cpl_wcs_get_image_dims(wcs);
2421 int naxis = cpl_wcs_get_image_naxis(wcs);
2422
2423
2424 /* Check NAXIS */
2425 for (cpl_size i = 0; i < naxis; i++) {
2426 if (i == 0) {
2427 cpl_propertylist_update_int(header, "NAXIS", naxis);
2428 }
2429 char * buf = cpl_sprintf("NAXIS%lld", i + 1);
2430 cpl_propertylist_update_int(header, buf,
2431 cpl_array_get_int(dims, i, &err));
2432 cpl_free(buf);
2433 }
2434
2435/* Make sure to have the right NAXIS keywords if 2D is forced */
2436 if (only2d == TRUE) {
2437 cpl_propertylist_update_int(header, "NAXIS", 2);
2438
2439 if(cpl_propertylist_has(header, "NAXIS3")){
2440 cpl_propertylist_erase(header, "NAXIS3");
2441 }
2442 }
2443
2444 /* for 2D images */
2445 if (crval) {
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));
2450 }
2451
2452 if (crpix) {
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));
2457 }
2458
2459 if (ctype) {
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));
2464 }
2465
2466 if (cunit) {
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));
2471 }
2472
2473 if (cd) {
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);
2482 }
2483
2484 /* for 3D cubes */
2485 if (only2d == FALSE && cpl_array_get_size(crval) > 2) {
2486
2487 if (crval) {
2488 cpl_propertylist_update_double(header, "CRVAL3",
2489 cpl_array_get_double(crval, 2, &err));
2490 }
2491
2492 if (crpix) {
2493 cpl_propertylist_update_double(header, "CRPIX3",
2494 cpl_array_get_double(crpix, 2, &err));
2495 }
2496
2497 if (ctype) {
2498 cpl_propertylist_update_string(header, "CTYPE3",
2499 cpl_array_get_string(ctype, 2));
2500 }
2501
2502 if (cunit) {
2503 cpl_propertylist_update_string(header, "CUNIT3",
2504 cpl_array_get_string(cunit, 2));
2505 }
2506
2507 if (cd) {
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);
2518 }
2519 } /* if 3D */
2520 return CPL_ERROR_NONE;
2521}
2522
2523/*----------------------------------------------------------------------------*/
2530/*----------------------------------------------------------------------------*/
2531
2532static cpl_error_code
2533hdrl_resample_create_table(cpl_table** tab, const cpl_size size)
2534{
2535
2536 cpl_ensure_code(tab, CPL_ERROR_NULL_INPUT);
2537 cpl_ensure_code(size > 0 , CPL_ERROR_ILLEGAL_INPUT);
2538
2539 *tab = cpl_table_new(size);
2540 /*
2541 * MUSE_PIXTABLE_XPOS "xpos" *< x-position column of a MUSE pixel table
2542 * MUSE_PIXTABLE_YPOS "ypos" *< y-position column of a MUSE pixel table
2543 * MUSE_PIXTABLE_LAMBDA "lambda" *< wavelength column of a MUSE pixel table
2544 * MUSE_PIXTABLE_DATA "data" *< data column of a MUSE pixel table
2545 * MUSE_PIXTABLE_DQ "dq" *< data quality column of a MUSE pixel table
2546 * MUSE_PIXTABLE_STAT "stat" *< error column of a MUSE pixel table
2547 */
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);
2560
2561 /* init column values */
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.);
2568
2569 return cpl_error_get_code();
2570}
2571
2572/*----------------------------------------------------------------------------*/
2580/*----------------------------------------------------------------------------*/
2581cpl_table*
2582hdrl_resample_image_to_table(const hdrl_image * hima, const cpl_wcs* wcs)
2583{
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;
2588 cpl_table* tab;
2589 ilist = hdrl_imagelist_new();
2590 hdrl_imagelist_set(ilist, (hdrl_image *)hima, 0);
2591
2592 tab = hdrl_resample_imagelist_to_table(ilist, wcs);
2593
2594 /* cleanup memory */
2595 hdrl_imagelist_unset(ilist, 0);
2596 hdrl_imagelist_delete(ilist);
2597
2598
2599 return tab;
2600}
2601
2602/*----------------------------------------------------------------------------*/
2610/*----------------------------------------------------------------------------*/
2611cpl_table*
2612hdrl_resample_imagelist_to_table(const hdrl_imagelist* himlist,
2613 const cpl_wcs* wcs)
2614{
2615
2616 cpl_ensure(himlist, CPL_ERROR_NULL_INPUT, NULL);
2617 cpl_ensure(wcs, CPL_ERROR_NULL_INPUT, NULL);
2618
2619 cpl_msg_debug(cpl_func, "Converting Dataset to table");
2620
2621 cpl_size naxis1 = hdrl_imagelist_get_size_x(himlist);
2622 cpl_size naxis2 = hdrl_imagelist_get_size_y(himlist);
2623 cpl_size naxis3 = hdrl_imagelist_get_size(himlist);
2624
2625 cpl_msg_debug(cpl_func,"Dataset dimentions (x, y, l): (%lld, %lld, %lld)",
2626 naxis1, naxis2, naxis3);
2627
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);
2631
2632 int testerr = 0;
2633 double crpix3 = 0.;
2634 double crval3 = 0.;
2635 double cdelt3 = 0.;
2636
2637 if (naxis3 > 1) { /* We have a cube */
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); /* CD3_3 */
2641/*
2642 cpl_msg_debug(cpl_func,"crpix3: %g crval3: %g cdelt3: %g", crpix3,
2643 crval3, cdelt3);
2644*/
2645 }
2646
2647 cpl_size tab_size = naxis1 * naxis2 * naxis3;
2648 cpl_table* tab = NULL;
2649 /* Prefill the full table */
2650 hdrl_resample_create_table(&tab, tab_size);
2651
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);
2658
2659 struct timeval tv1, tv2;
2660 cpl_msg_debug(cpl_func,"Starting parallel loop in hdrl_imagelist_to_table");
2661 gettimeofday(&tv1, NULL);
2662
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;
2668
2669 /* Fill the data */
2670 const hdrl_image* hima = hdrl_imagelist_get_const(himlist, k);
2671 const cpl_image* imadata = hdrl_image_get_image_const(hima);
2672 const cpl_image* imaerrs = hdrl_image_get_error_const(hima);
2673 const cpl_mask* imamask = hdrl_image_get_mask_const(hima);
2674 const double* pimadata = cpl_image_get_data_double_const(imadata);
2675
2676 if (imaerrs) { /* Fill the errors */
2677 pimaerr = cpl_image_get_data_double_const(imaerrs);
2678 }
2679 if (imamask) { /* Fill the bpm */
2680 pimabpm = cpl_mask_get_data_const(imamask);
2681 }
2682
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],
2688 &ptabypos[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];
2693 /* Insert only good pixels */
2694 if (!isfinite(pimadata[j_naxis1 + i]) ||
2695 ptabbpm[raw] != CPL_BINARY_0 ){
2696 ptabbpm[raw] = CPL_BINARY_1;
2697
2698 }
2699 }
2700 }
2701 }
2702
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));
2707
2708 return tab;
2709}
2710
2711/*----------------------------------------------------------------------------*/
2732/*----------------------------------------------------------------------------*/
2733hdrl_parameter*
2735 const double delta_dec)
2736{
2737
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;
2743
2744 p->recalc_limits = CPL_TRUE;
2745
2746 /* This function asks to recalculate the limits in the hdrl_resample_compute
2747 * function - therefore we put dummy values for the moment. */
2748
2749 p->dec_min = 0.1;
2750 p->dec_max = 0.2;
2751 p->ra_min = 0.1;
2752 p->ra_max = 0.2;
2753
2754 /* in case of 2D sets some defaults dummy values for 3rd dimension */
2755 p->lambda_min = 0;
2756 p->lambda_max = 0;
2757 p->delta_lambda = 1;
2758
2759 /* Default field margin in percent taken from swarp. */
2760 p->fieldmargin = FIELDMARGIN;
2761
2762 p->wcs = NULL;
2763
2764 if (hdrl_resample_parameter_outgrid_verify((hdrl_parameter*)p) != CPL_ERROR_NONE) {
2765 cpl_free(p);
2766 return NULL;
2767 }
2768 return (hdrl_parameter *)p;
2769
2770}
2771
2772/*----------------------------------------------------------------------------*/
2794/*----------------------------------------------------------------------------*/
2795hdrl_parameter*
2797 const double delta_dec, const double delta_lambda)
2798{
2799
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;
2806
2807 p->recalc_limits = CPL_TRUE;
2808 /* This function asks to recalculate the limits in the hdrl_resample_compute
2809 * function - therefore we put dummy values for the moment. */
2810
2811 p->dec_min = 0.1;
2812 p->dec_max = 0.2;
2813 p->ra_min = 0.1;
2814 p->ra_max = 0.2;
2815 p->lambda_min = 0.;
2816 p->lambda_max = 0.;
2817
2818 /* Default field margin in percent taken from swarp. */
2819 p->fieldmargin = FIELDMARGIN;
2820
2821 p->wcs = NULL;
2822
2823 if (hdrl_resample_parameter_outgrid_verify((hdrl_parameter*)p) != CPL_ERROR_NONE) {
2824 cpl_free(p);
2825 return NULL;
2826 }
2827 return (hdrl_parameter *)p;
2828
2829}
2830
2831/*----------------------------------------------------------------------------*/
2854/*----------------------------------------------------------------------------*/
2855hdrl_parameter*
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)
2863{
2864
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;
2870
2871
2872 p->recalc_limits = CPL_FALSE; /*This function takes the limits from the user*/
2873 p->dec_min = dec_min;
2874 p->dec_max = dec_max;
2875 p->ra_min = ra_min;
2876 p->ra_max = ra_max;
2877
2878 /* in case of 2D sets some defaults dummy values for 3rd dimension */
2879 p->lambda_min = 0.;
2880 p->lambda_max = 0.;
2881 p->delta_lambda = 1.;
2882
2883 p->fieldmargin = fieldmargin;
2884
2885 p->wcs = NULL;
2886
2887 if (hdrl_resample_parameter_outgrid_verify((hdrl_parameter*)p) != CPL_ERROR_NONE) {
2888 cpl_free(p);
2889 return NULL;
2890 }
2891 return (hdrl_parameter *)p;
2892
2893}
2894
2895/*----------------------------------------------------------------------------*/
2923/*----------------------------------------------------------------------------*/
2924hdrl_parameter*
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)
2935{
2936
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;
2943
2944 p->recalc_limits = CPL_FALSE; /*This function takes the limits from the user*/
2945
2946 p->dec_min = dec_min;
2947 p->dec_max = dec_max;
2948 p->ra_min = ra_min;
2949 p->ra_max = ra_max;
2950
2951 p->lambda_min = lambda_min;
2952 p->lambda_max = lambda_max;
2953
2954 p->fieldmargin = fieldmargin;
2955
2956 p->wcs = NULL;
2957
2958 if (hdrl_resample_parameter_outgrid_verify((hdrl_parameter*)p) != CPL_ERROR_NONE) {
2959 cpl_free(p);
2960 return NULL;
2961 }
2962 return (hdrl_parameter *)p;
2963
2964}
2965
2966/*----------------------------------------------------------------------------*/
2988/*----------------------------------------------------------------------------*/
2989hdrl_parameter*
2991 cpl_boolean use_errorweights,
2992 const double critical_radius)
2993{
2994
2995 hdrl_resample_method_parameter * p = (hdrl_resample_method_parameter *)
2996 hdrl_parameter_new(&hdrl_resample_method_parameter_type);
2997
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;
3002
3003 /* fill rest with dummy input */
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;
3008
3009 if (hdrl_resample_parameter_method_verify((hdrl_parameter*)p) != CPL_ERROR_NONE) {
3010 cpl_free(p);
3011 return NULL;
3012 }
3013 return (hdrl_parameter*) p;
3014}
3015
3016/*----------------------------------------------------------------------------*/
3035/*----------------------------------------------------------------------------*/
3036hdrl_parameter*
3038 cpl_boolean use_errorweights)
3039{
3040
3041 hdrl_resample_method_parameter * p = (hdrl_resample_method_parameter *)
3042 hdrl_parameter_new(&hdrl_resample_method_parameter_type);
3043
3044
3045 p->method = HDRL_RESAMPLE_METHOD_LINEAR;
3046 p->loop_distance = loop_distance;
3047 p->use_errorweights = use_errorweights;
3048
3049 /* fill rest with dummy input */
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;
3055
3056 if (hdrl_resample_parameter_method_verify((hdrl_parameter*)p) != CPL_ERROR_NONE) {
3057 cpl_free(p);
3058 return NULL;
3059 }
3060 return (hdrl_parameter*) p;
3061}
3062
3063/*----------------------------------------------------------------------------*/
3083/*----------------------------------------------------------------------------*/
3084hdrl_parameter*
3086 cpl_boolean use_errorweights)
3087{
3088
3089 hdrl_resample_method_parameter * p = (hdrl_resample_method_parameter *)
3090 hdrl_parameter_new(&hdrl_resample_method_parameter_type);
3091
3092 p->method = HDRL_RESAMPLE_METHOD_QUADRATIC;
3093 p->loop_distance = loop_distance;
3094 p->use_errorweights = use_errorweights;
3095
3096 /* fill rest with dummy input */
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;
3102
3103 if (hdrl_resample_parameter_method_verify((hdrl_parameter*)p) != CPL_ERROR_NONE) {
3104 cpl_free(p);
3105 return NULL;
3106 }
3107 return (hdrl_parameter*) p;
3108}
3109
3110/*----------------------------------------------------------------------------*/
3125/*----------------------------------------------------------------------------*/
3126hdrl_parameter*
3128{
3129
3130 hdrl_resample_method_parameter * p = (hdrl_resample_method_parameter *)
3131 hdrl_parameter_new(&hdrl_resample_method_parameter_type);
3132
3133
3134 p->method = HDRL_RESAMPLE_METHOD_NEAREST;
3135 p->loop_distance = 0;
3136 p->use_errorweights = CPL_FALSE;
3137
3138 /* fill rest with dummy input */
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;
3144
3145 if (hdrl_resample_parameter_method_verify((hdrl_parameter*)p) != CPL_ERROR_NONE) {
3146 cpl_free(p);
3147 return NULL;
3148 }
3149 return (hdrl_parameter*) p;
3150}
3151
3152/*----------------------------------------------------------------------------*/
3173/*----------------------------------------------------------------------------*/
3174hdrl_parameter*
3176 cpl_boolean use_errorweights,
3177 const int kernel_size)
3178{
3179
3180 hdrl_resample_method_parameter * p = (hdrl_resample_method_parameter *)
3181 hdrl_parameter_new(&hdrl_resample_method_parameter_type);
3182
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;
3187 /* fill rest with dummy input */
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;
3192
3193 if (hdrl_resample_parameter_method_verify((hdrl_parameter*)p) != CPL_ERROR_NONE) {
3194 cpl_free(p);
3195 return NULL;
3196 }
3197 return (hdrl_parameter*) p;
3198}
3199
3200/*----------------------------------------------------------------------------*/
3228/*----------------------------------------------------------------------------*/
3229hdrl_parameter*
3231 cpl_boolean use_errorweights,
3232 const double pix_frac_x,
3233 const double pix_frac_y,
3234 const double pix_frac_lambda)
3235{
3236
3237 hdrl_resample_method_parameter * p = (hdrl_resample_method_parameter *)
3238 hdrl_parameter_new(&hdrl_resample_method_parameter_type);
3239
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;
3246
3247 /* fill rest with dummy input */
3248 p->renka_critical_radius = 0.1;
3249 p->lanczos_kernel_size = 2;
3250
3251 if (hdrl_resample_parameter_method_verify((hdrl_parameter*)p) != CPL_ERROR_NONE) {
3252 cpl_free(p);
3253 return NULL;
3254 }
3255 return (hdrl_parameter*) p;
3256
3257}
3258/*----------------------------------------------------------------------------*/
3267/*----------------------------------------------------------------------------*/
3268cpl_error_code
3270
3271 const hdrl_resample_outgrid_parameter * param_loc =
3272 (const hdrl_resample_outgrid_parameter *)hp ;
3273
3274 cpl_error_ensure(hp != NULL, CPL_ERROR_NULL_INPUT,
3275 return CPL_ERROR_NULL_INPUT, "NULL Input Parameters");
3276
3277 cpl_error_ensure(hdrl_resample_parameter_outgrid_check(hp),
3278 CPL_ERROR_ILLEGAL_INPUT, return CPL_ERROR_ILLEGAL_INPUT,
3279 "Here we expect a resample outgrid parameter") ;
3280
3281 cpl_error_ensure(
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");
3286
3287/*
3288 * The wcs is filled later on by the compute function so it can not be checked
3289 * at this stage
3290 *
3291 cpl_error_ensure(param_loc->wcs != NULL, CPL_ERROR_NULL_INPUT,
3292 return CPL_ERROR_NULL_INPUT, "WCS must be defined");
3293*/
3294
3295
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");
3299
3300 cpl_error_ensure(param_loc->delta_dec > 0, CPL_ERROR_ILLEGAL_INPUT,
3301 return CPL_ERROR_ILLEGAL_INPUT, "declination stepsize "
3302 "must be > 0");
3303
3304 cpl_error_ensure(param_loc->delta_lambda > 0, CPL_ERROR_ILLEGAL_INPUT,
3305 return CPL_ERROR_ILLEGAL_INPUT, "wavelength stepsize "
3306 "must be > 0");
3307
3308
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");
3312
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");
3316
3317
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");
3321
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");
3325
3326 cpl_error_ensure(param_loc->fieldmargin >= 0., CPL_ERROR_ILLEGAL_INPUT,
3327 return CPL_ERROR_ILLEGAL_INPUT, "The field margin must"
3328 " be >= 0.");
3329
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");
3338
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");
3343
3344 return CPL_ERROR_NONE ;
3345}
3346
3347/*----------------------------------------------------------------------------*/
3356/*----------------------------------------------------------------------------*/
3357cpl_error_code
3358hdrl_resample_parameter_method_verify(const hdrl_parameter * hp){
3359
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");
3364
3365 cpl_error_ensure(hdrl_resample_parameter_method_check(hp),
3366 CPL_ERROR_ILLEGAL_INPUT, return CPL_ERROR_ILLEGAL_INPUT,
3367 "Here we expect a resample method parameter") ;
3368
3369 /* checks on parameter methods */
3370 cpl_error_ensure(
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");
3379
3380 /* checks on common parameter elements */
3381 cpl_error_ensure(param_loc->loop_distance >= 0, CPL_ERROR_ILLEGAL_INPUT,
3382 return CPL_ERROR_ILLEGAL_INPUT, "The loop distance must "
3383 "be >=0");
3384
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");
3389
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:
3395 break;
3396 case HDRL_RESAMPLE_METHOD_RENKA:
3397 /* checks on peculiar parameter elements */
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");
3401 break ;
3402
3403 case HDRL_RESAMPLE_METHOD_DRIZZLE:
3404 /* checks on peculiar parameter elements */
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 "
3408 "must be > 0");
3409
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 "
3413 "must be > 0");
3414
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 "
3418 "must be > 0");
3419 break;
3420
3421 case HDRL_RESAMPLE_METHOD_LANCZOS:
3422 /* checks on peculiar parameter elements */
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");
3426 break;
3427 }
3428
3429 return CPL_ERROR_NONE ;
3430}
3431/*----------------------------------------------------------------------------*/
3440/*----------------------------------------------------------------------------*/
3441int
3442hdrl_resample_parameter_outgrid_check(const hdrl_parameter * self){
3443 /* Check if the method is of proper type */
3444 cpl_ensure_code(self != NULL, CPL_ERROR_NULL_INPUT);
3445
3446 return hdrl_parameter_check_type(self, &hdrl_resample_outgrid_parameter_type);
3447}
3448/*----------------------------------------------------------------------------*/
3457/*----------------------------------------------------------------------------*/
3458int
3459hdrl_resample_parameter_method_check(const hdrl_parameter * self){
3460 /* Check if the method is of proper type */
3461 cpl_ensure_code(self != NULL, CPL_ERROR_NULL_INPUT);
3462 return hdrl_parameter_check_type(self, &hdrl_resample_method_parameter_type);
3463}
3464
3465
3466/*----------------------------------------------------------------------------*/
3473/*----------------------------------------------------------------------------*/
3474static cpl_error_code
3475hdrl_resample_inputtable_verify(const cpl_table *ResTable){
3476
3477 cpl_error_ensure(ResTable != NULL, CPL_ERROR_NULL_INPUT,
3478 return CPL_ERROR_NULL_INPUT, "No Table as input");
3479
3480 /* Check the existence of all columns */
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");
3505
3506 /* Check the format of all columns */
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");
3531
3532 return cpl_error_get_code();
3533}
3535/* EOF */
3536
cpl_mask * hdrl_image_get_mask(hdrl_image *himg)
get cpl bad pixel mask from image
Definition: hdrl_image.c:157
cpl_image * hdrl_image_get_error(hdrl_image *himg)
get error as cpl image
Definition: hdrl_image.c:131
const cpl_mask * hdrl_image_get_mask_const(const hdrl_image *himg)
get cpl bad pixel mask from image
Definition: hdrl_image.c:175
const cpl_image * hdrl_image_get_error_const(const hdrl_image *himg)
get error as cpl image
Definition: hdrl_image.c:144
cpl_image * hdrl_image_get_image(hdrl_image *himg)
get data as cpl image
Definition: hdrl_image.c:105
cpl_error_code hdrl_image_reject(hdrl_image *self, cpl_size xpos, cpl_size ypos)
mark pixel as bad
Definition: hdrl_image.c:427
const cpl_image * hdrl_image_get_image_const(const hdrl_image *himg)
get data as cpl image
Definition: hdrl_image.c:118
hdrl_image * hdrl_image_new(cpl_size nx, cpl_size ny)
create new zero filled hdrl image
Definition: hdrl_image.c:311
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().