ERIS Pipeline Reference Manual 1.8.15
eris_ifu_resample.c
1/* $Id: eris_resample.c,v 1.9 2013-03-26 17:00:44 jtaylor Exp $
2 *
3 * This file is part of the ERIS Pipeline
4 * Copyright (C) 2002,2003 European Southern Observatory
5 *
6 * This program is free software; you can redistribute it and/or modify
7 * it under the terms of the GNU General Public License as published by
8 * the Free Software Foundation; either version 2 of the License, or
9 * (at your option) any later version.
10 *
11 * This program is distributed in the hope that it will be useful,
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 * GNU General Public License for more details.
15 *
16 * You should have received a copy of the GNU General Public License
17 * along with this program; if not, write to the Free Software
18 * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
19 */
20
21
22#ifdef HAVE_CONFIG_H
23#include <config.h>
24#endif
25
26/*-----------------------------------------------------------------------------
27 Includes
28 -----------------------------------------------------------------------------*/
29
30
31
32#include <string.h>
33#include <hdrl.h>
34#include <eris_ifu_resample.h>
35#include <eris_ifu_dfs.h>
36#include <eris_ifu_utils.h>
37#include <eris_utils.h>
38#include <eris_pfits.h>
39#include <eris_ifu_sdp.h>
40/*----------------------------------------------------------------------------*/
64/*----------------------------------------------------------------------------*/
65#define MAX_COMBINE_CUBE_PLANE_DIAGONAL 1000
66
67
83static cpl_error_code
84eris_ifu_resample_update_table(const cpl_table* append, cpl_table** storage) {
85
86 cpl_ensure_code(append, CPL_ERROR_NULL_INPUT);
87
88 if(*storage == NULL) {
89 *storage = cpl_table_new(0);
90 cpl_table_copy_structure(*storage, append);
91 cpl_table_insert(*storage, append, 0);
92 } else {
93 cpl_size nrow = cpl_table_get_nrow(*storage);
94 cpl_table_insert(*storage, append, nrow);
95 }
96 eris_check_error_code("eris_ifu_resample_update_table");
97 return cpl_error_get_code();
98
99}
100
101
122static cpl_error_code
123eris_ifu_resample_tablesave (const char* recipe_name, const cpl_parameter *par,
124 const cpl_parameterlist *parlist, cpl_frameset *frameset,
125 cpl_table *tab_final, const char* pipefile_prefix)
126{
127 /* Save the final table if required */
128 char* param_name = cpl_sprintf("eris.%s.save-table", recipe_name);
129 par = cpl_parameterlist_find_const (parlist, param_name);
130 cpl_free(param_name);
131
132 const cpl_boolean save_table = cpl_parameter_get_bool (par);
133
134 if (save_table)
135 {
136 cpl_msg_info (cpl_func, "Saving the table before resampling");
137 cpl_propertylist *applist = cpl_propertylist_new ();
138 cpl_propertylist_update_string (applist, CPL_DFS_PRO_CATG,
139 "RESAMPLE_TABLE");
140 char* fname = cpl_sprintf("%s_resample_table.fits", pipefile_prefix);
141 cpl_dfs_save_table (frameset, NULL, parlist, frameset, NULL, tab_final,
142 NULL, recipe_name, applist, NULL, PACKAGE "/" PACKAGE_VERSION,
143 fname);
144 cpl_free(fname);
145 cpl_propertylist_delete (applist);
146 }
147 eris_check_error_code("eris_ifu_resample_tablesave");
148 return cpl_error_get_code();
149}
150
151
152
153/*----------------------------------------------------------------------------*/
169/*----------------------------------------------------------------------------*/
170
171/* set propertylist (int) value */
172static cpl_error_code
173eris_ifu_resample_parameters_get_int(const cpl_parameterlist* parlist,
174 const char* pname, int *pvalue)
175{
176 const cpl_parameter* p = cpl_parameterlist_find_const(parlist, pname);
177 cpl_boolean p_has_changed = eris_param_has_changed(p);
178
179 if ( cpl_parameter_get_default_flag(p) && p_has_changed != CPL_FALSE) {
180 *pvalue = cpl_parameter_get_int(p);
181 } else {
182 *pvalue = cpl_parameter_get_default_int(p);
183 }
184 eris_check_error_code("eris_ifu_resample_parameters_get_int");
185 return cpl_error_get_code();
186}
187
188/*----------------------------------------------------------------------------*/
205/*----------------------------------------------------------------------------*/
206
207
208static cpl_error_code
209set_double_from_parameters(const cpl_parameterlist* parlist,
210 const char* pname,
211 const double table_value,
212 double *value) {
213
214 double temp_value;
215 eris_parameters_get_double (parlist, pname, &temp_value);
216 if (temp_value != -1)
217 {
218 *value = temp_value;
219 }
220 else
221 {
222 *value = table_value;
223 }
224 eris_check_error_code("set_double_from_parameters");
225 return cpl_error_get_code();
226}
227
228
229/*----------------------------------------------------------------------------*/
246/*----------------------------------------------------------------------------*/
247
248
249static cpl_error_code
250set_int_from_parameters(const cpl_parameterlist* parlist,
251 const char* pname,
252 const int table_value,
253 int *value) {
254
255 int temp_value;
256 eris_ifu_resample_parameters_get_int (parlist, pname, &temp_value);
257 if (temp_value != -1)
258 {
259 *value = temp_value;
260 }
261 else
262 {
263 *value = table_value;
264 }
265 eris_check_error_code("set_int_from_parameters");
266 return cpl_error_get_code();
267}
268
269
270
271/*----------------------------------------------------------------------------*/
294/*----------------------------------------------------------------------------*/
295
296static hdrl_parameter *
297eris_ifu_resample_set_method (const cpl_parameterlist *parlist,
298 const char* context)
299{
300 hdrl_parameter *aParams_method = NULL;
301 const cpl_parameter *p = NULL;
302 cpl_boolean p_has_changed = 0;
303 int loop_distance = 0;
304 double critical_radius_renka = 0.;
305 int kernel_size_lanczos = 0;
306 double pix_frac_drizzle_x = 0.;
307 double pix_frac_drizzle_y= 0.;
308 double pix_frac_drizzle_l= 0.;
309 cpl_boolean use_errorweights = FALSE;
310
311 const cpl_parameter *par = NULL;
312 char* param_name = NULL;
313 /* Apply additional weights based on the error if required */
314
315 param_name = cpl_sprintf("%s.method.use-errorweights",context);
316 par = cpl_parameterlist_find_const (parlist, param_name);
317 use_errorweights = cpl_parameter_get_bool (par);
318 cpl_free(param_name);
319
321 param_name = cpl_sprintf("%s.method.loop-distance",context);
322 p = cpl_parameterlist_find_const (parlist, param_name);
323
324 p_has_changed = eris_param_has_changed (p);
325 if (cpl_parameter_get_default_flag (p) && p_has_changed != CPL_FALSE)
326 {
327 loop_distance = cpl_parameter_get_int (p);
328 }
329 else
330 {
331 loop_distance = LOOP_DISTANCE;
332 }
333 cpl_free(param_name);
334
336 param_name = cpl_sprintf("%s.method.renka.critical-radius",context);
337 set_double_from_parameters(parlist, param_name,
338 RENKA_CRITICAL_RADIUS, &critical_radius_renka);
339 cpl_free(param_name);
340
342 param_name = cpl_sprintf("%s.method.lanczos.kernel-size",context);
343 set_int_from_parameters(parlist, param_name,
344 LANCZOS_KERNEL_SIZE, &kernel_size_lanczos);
345 cpl_free(param_name);
346
352 param_name = cpl_sprintf("%s.method.drizzle.downscale-x",context);
353 set_double_from_parameters(parlist, param_name, DRIZZLE_DOWN_SCALING_FACTOR_X,
354 &pix_frac_drizzle_x);
355 cpl_free(param_name);
356
357 param_name = cpl_sprintf("%s.method.drizzle.downscale-y",context);
358 set_double_from_parameters(parlist, param_name, DRIZZLE_DOWN_SCALING_FACTOR_Y,
359 &pix_frac_drizzle_y);
360 cpl_free(param_name);
361
362 param_name = cpl_sprintf("%s.method.drizzle.downscale-z",context);
363 set_double_from_parameters(parlist, param_name, DRIZZLE_DOWN_SCALING_FACTOR_Z,
364 &pix_frac_drizzle_l);
365 cpl_free(param_name);
366
367 /* Create the right re-sampling parameter */
368 param_name = cpl_sprintf("%s.method",context);
369 par = cpl_parameterlist_find_const (parlist, param_name);
370 const char *method = cpl_parameter_get_string (par);
371
372 if (strcmp (method, "NEAREST") == 0)
373 {
375 }
376 else if (strcmp (method, "RENKA") == 0)
377 {
378 aParams_method = hdrl_resample_parameter_create_renka(loop_distance,
379 use_errorweights,
380 critical_radius_renka);
381 }
382 else if (strcmp (method, "LINEAR") == 0)
383 {
384 aParams_method = hdrl_resample_parameter_create_linear(loop_distance,
385 use_errorweights);
386 }
387 else if (strcmp (method, "QUADRATIC") == 0)
388 {
389 aParams_method = hdrl_resample_parameter_create_quadratic(loop_distance,
390 use_errorweights);
391 }
392 else if (strcmp (method, "DRIZZLE") == 0)
393 {
394 aParams_method = hdrl_resample_parameter_create_drizzle(loop_distance,
395 use_errorweights,
396 pix_frac_drizzle_x,
397 pix_frac_drizzle_y,
398 pix_frac_drizzle_l);
399 }
400 else if (strcmp (method, "LANCZOS") == 0)
401 {
402 aParams_method = hdrl_resample_parameter_create_lanczos(loop_distance,
403 use_errorweights,
404 kernel_size_lanczos);
405 }
406 else
407 {
408 aParams_method = hdrl_resample_parameter_create_lanczos(loop_distance,
409 use_errorweights,
410 kernel_size_lanczos);
411 cpl_msg_warning (cpl_func,
412 "%s is an unsupported method! Default to LANCZOS",
413 method);
414 }
415 cpl_free(param_name);
416 //eris_ifu_resample_print_method_params(aParams_method);
417 eris_check_error_code("eris_ifu_resample_set_method");
418 return aParams_method;
419}
420
421/*----------------------------------------------------------------------------*/
441/*----------------------------------------------------------------------------*/
442
443static cpl_error_code
444eris_ifu_resample_get_cdelt123(cpl_wcs *wcs, double *cdelt1, double *cdelt2,
445 double *cdelt3)
446{
447 cpl_ensure_code(wcs, CPL_ERROR_NULL_INPUT);
448
449
450 double cd1_1 = 0;
451 double cd1_2 = 0;
452 double cd2_1 = 0;
453 double cd2_2 = 0;
454 double cd3_3 = 0;
455
456 double dx=0;
457 double dy=0;
458
459 const cpl_matrix *cd = cpl_wcs_get_cd(wcs);
460 cd1_1 = cpl_matrix_get(cd, 0, 0);
461 cd1_2 = cpl_matrix_get(cd, 1, 0);
462 cd2_1 = cpl_matrix_get(cd, 0, 1);
463 cd2_2 = cpl_matrix_get(cd, 1, 1);
464
465 dx = sqrt(pow(cd1_1,2) + pow(cd2_1,2));
466 dy = sqrt(pow(cd1_2,2)+ pow(cd2_2,2));
467
468// double det = cd1_1 * cd2_2 - cd1_2 * cd2_1;
469// if (det < 0)
470// dx *= -1.0;
471
472// *cdelt1 = fabs(cpl_matrix_get(cd, 0, 0)); /* CD1_1 */
473// *cdelt2 = fabs(cpl_matrix_get(cd, 1, 1)); /* CD2_2 */
474 *cdelt1 = dx;
475 *cdelt2 = dy;
476
477 /* Changes for Cube and image */
478 cpl_size cdncol = cpl_matrix_get_ncol(cd);
479 if (cdncol == 3) {
480 cd3_3= fabs(cpl_matrix_get(cd, 2, 2)); /* CD3_3 */
481 }
482 else {
483 cd3_3 = 1.; /* CD3_3 */
484 }
485
486 *cdelt3 = cd3_3;
487
488 /* Check keys (for debug) */
489 cpl_msg_debug(cpl_func, "cdelt1: %g, cdelt2: %g, cdelt3: %g",
490 *cdelt1, *cdelt2, *cdelt3);
491 eris_check_error_code("eris_ifu_resample_get_cdelt123");
492 return cpl_error_get_code();
493}
494
495
496/*----------------------------------------------------------------------------*/
518/*----------------------------------------------------------------------------*/
519
520static hdrl_parameter *
521eris_ifu_resample_set_outputgrid (const char* recipe_name,
522 const cpl_parameterlist *parlist, cpl_table *muse_table, cpl_wcs *wcs)
523{
524 /* Should be done in a dedicated _new function */
525 cpl_error_ensure(parlist != NULL, CPL_ERROR_NULL_INPUT,
526 return NULL, "NULL Input Parameters");
527 cpl_error_ensure(muse_table != NULL, CPL_ERROR_NULL_INPUT,
528 return NULL, "NULL Input table");
529 cpl_error_ensure(wcs != NULL, CPL_ERROR_NULL_INPUT,
530 return NULL, "NULL Input wcs");
531
532 hdrl_parameter *aParams_outputgrid = NULL;
533 char* param_name = NULL;
534 char* context = NULL;
535 double ra_min = 0.;
536 double ra_max = 0.;
537 double dec_min = 0.;
538 double dec_max = 0.;
539 double lambda_min = 0.;
540 double lambda_max = 0.;
541
542 double dx = 0.;
543 double dy = 0.;
544 double dlambda = 0.;
545
546 context = cpl_sprintf("eris.%s", recipe_name);
547
548 /* init relevant parameters */
549 double ra_min_tmp =
550 cpl_table_get_column_min (muse_table, HDRL_RESAMPLE_TABLE_RA);
551 double ra_max_tmp =
552 cpl_table_get_column_max (muse_table, HDRL_RESAMPLE_TABLE_RA);
553 double dec_min_tmp =
554 cpl_table_get_column_min (muse_table, HDRL_RESAMPLE_TABLE_DEC);
555 double dec_max_tmp =
556 cpl_table_get_column_max (muse_table, HDRL_RESAMPLE_TABLE_DEC);
557 double lambda_min_tmp =
558 cpl_table_get_column_min (muse_table, HDRL_RESAMPLE_TABLE_LAMBDA);
559 double lambda_max_tmp =
560 cpl_table_get_column_max (muse_table, HDRL_RESAMPLE_TABLE_LAMBDA);
561
562 /* We have the rare case that the image spans over ra = 0.*/
563 if(ra_max_tmp - ra_min_tmp > 180){
564 const double *ra = cpl_table_get_data_double_const(muse_table,
565 HDRL_RESAMPLE_TABLE_RA);
566 /* Should we also take the bpm into account ?
567 const int *bpm = cpl_table_get_data_int_const(ResTable,
568 HDRL_RESAMPLE_TABLE_BPM);
569 */
570
571 /* set to extreme values for a start */
572 ra_min_tmp = 0.;
573 ra_max_tmp = 360.;
574 cpl_size nrow = cpl_table_get_nrow(muse_table);
575
576 for (cpl_size i = 0; i < nrow; i++) {
577 if (ra[i] > ra_min_tmp && ra[i] <= 180.) ra_min_tmp = ra[i]; /* get the maximum */
578 if (ra[i] < ra_max_tmp && ra[i] > 180.) ra_max_tmp = ra[i]; /* get the minimum */
579 }
580 }
581
582 /* Check output (for debug) */
583 cpl_msg_debug (cpl_func, "min x %10.7f", ra_min_tmp);
584 cpl_msg_debug (cpl_func, "max x %10.7f", ra_max_tmp);
585 cpl_msg_debug (cpl_func, "min y %10.7f", dec_min_tmp);
586 cpl_msg_debug (cpl_func, "max y %10.7f", dec_max_tmp);
587 cpl_msg_debug (cpl_func, "min lambda %10.7f", lambda_min_tmp);
588 cpl_msg_debug (cpl_func, "max lambda %10.7f", lambda_max_tmp);
589 /*
590 param_name = cpl_sprintf("%s.outgrid.ra-min", context);
591 set_double_from_parameters(parlist, param_name, ra_min_tmp, &ra_min);
592 cpl_free(param_name);
593 */
594 ra_min = ra_min_tmp;
595 /*
596 param_name = cpl_sprintf("%s.outgrid.ra-max", context);
597 set_double_from_parameters(parlist, param_name, ra_max_tmp, &ra_max);
598 cpl_free(param_name);
599 */
600 ra_max = ra_max_tmp;
601
602 /*
603 param_name = cpl_sprintf("%s.outgrid.dec-min", context);
604 set_double_from_parameters(parlist, param_name, dec_min_tmp, &dec_min);
605 cpl_free(param_name);
606 */
607 dec_min = dec_min_tmp;
608
609 /*
610 param_name = cpl_sprintf("%s.outgrid.dec-max", context);
611 set_double_from_parameters(parlist, param_name, dec_max_tmp, &dec_max);
612 cpl_free(param_name);
613 */
614 dec_max = dec_max_tmp;
615
616 /*
617 param_name = cpl_sprintf("%s.outgrid.lambda-min", context);
618 set_double_from_parameters(parlist, param_name, lambda_min_tmp, &lambda_min);
619 cpl_free(param_name);
620
621 param_name = cpl_sprintf("%s.outgrid.lambda-max", context);
622 set_double_from_parameters(parlist, param_name, lambda_max_tmp, &lambda_max);
623 cpl_free(param_name);
624 */
625 lambda_max = lambda_max_tmp;
626
627
628 /* Reset all delta values */
629 double cdelt1 = 0., cdelt2 = 0., cdelt3 = 0.;
630 eris_ifu_resample_get_cdelt123 (wcs, &cdelt1, &cdelt2, &cdelt3);
631 /*
632 param_name = cpl_sprintf("%s.outgrid.delta-ra", context);
633 set_double_from_parameters(parlist, param_name, cdelt1, &dx);
634 cpl_free(param_name);
635 */
636 dx = cdelt1;
637
638 /*
639 param_name = cpl_sprintf("%s.outgrid.delta-dec", context);
640 set_double_from_parameters(parlist, param_name, cdelt2, &dy);
641 cpl_free(param_name);
642 */
643 dy = cdelt2;
644
645 /*
646 param_name = cpl_sprintf("%s.outgrid.delta-lambda", context);
647 set_double_from_parameters(parlist, param_name, cdelt3, &dlambda);
648 cpl_free(param_name);
649 */
650 dlambda = cdelt3;
651
652
653 /* Assign the field margin */
654 const cpl_parameter *par = NULL;
655 param_name = cpl_sprintf("%s.fieldmargin", context);
656 par = cpl_parameterlist_find_const (parlist, param_name);
657 cpl_free(param_name);
658
659 double fieldmargin = 0.;
660 fieldmargin = cpl_parameter_get_double(par);
661
662 /* create final outgrid parameter structure */
663 int naxis = cpl_wcs_get_image_naxis(wcs);
664
665 if(naxis == 2) {
666 aParams_outputgrid =
668 ra_min, ra_max,
669 dec_min, dec_max,
670 fieldmargin);
671 } else {
672 //cpl_msg_info(cpl_func,"dx: %g dy: %g dlambda: %g", dx, dy, dlambda);
673 //cpl_msg_info(cpl_func,"ra_min: %g ra_max: %g, dec_min: %g dec_max: %g, lambda_min: %g, lambda_max: %g fieldmarging: %g",
674 // ra_min, ra_max, dec_min, dec_max, lambda_min, lambda_max, fieldmargin);
675 aParams_outputgrid =
677 ra_min, ra_max,
678 dec_min, dec_max,
679 lambda_min, lambda_max,
680 fieldmargin);
681 }
682
683 double dist_x = (ra_max-ra_min)/dx;
684 double dist_y = (dec_max-dec_min)/dy;
685
686 cpl_msg_debug(cpl_func,"size1: %g size2: %g",dist_x, dist_y);
687 double dist = sqrt(dist_x * dist_x + dist_y * dist_y);
688 cpl_msg_info(cpl_func,"Output cube diagonal size: %g", dist);
689 if ( dist > MAX_COMBINE_CUBE_PLANE_DIAGONAL) {
690 cpl_msg_info(cpl_func,"Internal threshold on max cube plane diagonal size: %d",MAX_COMBINE_CUBE_PLANE_DIAGONAL);
691 cpl_msg_warning(cpl_func,"Resampled cube plane diagonal greater than threshold. Skip cube combination.");
692 cpl_error_set(cpl_func,CPL_ERROR_INCOMPATIBLE_INPUT);
693 cpl_free(context);
694 eris_check_error_code("eris_ifu_resample_set_outputgrid");
695 return aParams_outputgrid;
696 }
697 cpl_free(context);
698 eris_check_error_code("eris_ifu_resample_set_outputgrid");
699
700 return aParams_outputgrid;
701}
702
703/*----------------------------------------------------------------------------*/
723/*----------------------------------------------------------------------------*/
724
725static hdrl_parameter *
726eris_ifu_resample_set_outputgrid2D (const char* recipe_name,
727 const cpl_parameterlist *parlist, cpl_table *muse_table, cpl_wcs *wcs)
728{
729 /* Should be done in a dedicated _new function */
730 cpl_error_ensure(parlist != NULL, CPL_ERROR_NULL_INPUT,
731 return NULL, "NULL Input Parameters");
732 cpl_error_ensure(muse_table != NULL, CPL_ERROR_NULL_INPUT,
733 return NULL, "NULL Input table");
734 cpl_error_ensure(wcs != NULL, CPL_ERROR_NULL_INPUT,
735 return NULL, "NULL Input wcs");
736
737 hdrl_parameter *aParams_outputgrid = NULL;
738 char* param_name = NULL;
739 char* context = NULL;
740 double ra_min = 0.;
741 double ra_max = 0.;
742 double dec_min = 0.;
743 double dec_max = 0.;
744
745 double dx = 0.;
746 double dy = 0.;
747
748 context = cpl_sprintf("eris.%s", recipe_name);
749 //cpl_table_save (muse_table, NULL, NULL, "muse_table.fits", CPL_IO_CREATE);
750 /* init relevant parameters */
751 double ra_min_tmp =
752 cpl_table_get_column_min (muse_table, HDRL_RESAMPLE_TABLE_RA);
753 double ra_max_tmp =
754 cpl_table_get_column_max (muse_table, HDRL_RESAMPLE_TABLE_RA);
755 double dec_min_tmp =
756 cpl_table_get_column_min (muse_table, HDRL_RESAMPLE_TABLE_DEC);
757 double dec_max_tmp =
758 cpl_table_get_column_max (muse_table, HDRL_RESAMPLE_TABLE_DEC);
759
760 //cpl_msg_info(cpl_func,"ra_max_tmp: %g ra_min_tmp: %g",ra_max_tmp, ra_min_tmp);
761 //cpl_msg_info(cpl_func,"dec_max_tmp: %g dec_min_tmp: %g",dec_max_tmp, dec_min_tmp);
762 /* We have the rare case that the image spans over ra = 0.*/
763 if(ra_max_tmp - ra_min_tmp > 180){
764 const double *ra = cpl_table_get_data_double_const(muse_table,
765 HDRL_RESAMPLE_TABLE_RA);
766 /* Should we also take the bpm into account ?
767 const int *bpm = cpl_table_get_data_int_const(ResTable,
768 HDRL_RESAMPLE_TABLE_BPM);
769 */
770
771 /* set to extreme values for a start */
772 ra_min_tmp = 0.;
773 ra_max_tmp = 360.;
774 cpl_size nrow = cpl_table_get_nrow(muse_table);
775
776 for (cpl_size i = 0; i < nrow; i++) {
777 if (ra[i] > ra_min_tmp && ra[i] <= 180.) ra_min_tmp = ra[i]; /* get the maximum */
778 if (ra[i] < ra_max_tmp && ra[i] > 180.) ra_max_tmp = ra[i]; /* get the minimum */
779 }
780 }
781
782 /* Check output (for debug) */
783 cpl_msg_debug (cpl_func, "min x %10.7f", ra_min_tmp);
784 cpl_msg_debug (cpl_func, "max x %10.7f", ra_max_tmp);
785 cpl_msg_debug (cpl_func, "min y %10.7f", dec_min_tmp);
786 cpl_msg_debug (cpl_func, "max y %10.7f", dec_max_tmp);
787
788 /*
789 param_name = cpl_sprintf("%s.outgrid.ra-min", context);
790 set_double_from_parameters(parlist, param_name, ra_min_tmp, &ra_min);
791 cpl_free(param_name);
792 */
793 ra_min = ra_min_tmp;
794 /*
795 param_name = cpl_sprintf("%s.outgrid.ra-max", context);
796 set_double_from_parameters(parlist, param_name, ra_max_tmp, &ra_max);
797 cpl_free(param_name);
798 */
799 ra_max = ra_max_tmp;
800 /*
801 param_name = cpl_sprintf("%s.outgrid.dec-min", context);
802 set_double_from_parameters(parlist, param_name, dec_min_tmp, &dec_min);
803 cpl_free(param_name);
804 */
805 dec_min = dec_min_tmp;
806
807 /*
808 param_name = cpl_sprintf("%s.outgrid.dec-max", context);
809 set_double_from_parameters(parlist, param_name, dec_max_tmp, &dec_max);
810 cpl_free(param_name);
811 */
812 dec_max = dec_max_tmp;
813
814
815
816 /* Reset all delta values */
817 double cdelt1 = 0., cdelt2 = 0., cdelt3 = 0.;
818 eris_ifu_resample_get_cdelt123 (wcs, &cdelt1, &cdelt2, &cdelt3);
819 /*
820 param_name = cpl_sprintf("%s.outgrid.delta-ra", context);
821 set_double_from_parameters(parlist, param_name, cdelt1, &dx);
822 cpl_free(param_name);
823 */
824 dx = cdelt1;
825 /*
826 param_name = cpl_sprintf("%s.outgrid.delta-dec", context);
827 set_double_from_parameters(parlist, param_name, cdelt2, &dy);
828 cpl_free(param_name);
829 */
830 dy = cdelt2;
831
832 /* Assign the field margin */
833 const cpl_parameter *par = NULL;
834 param_name = cpl_sprintf("%s.fieldmargin", context);
835 par = cpl_parameterlist_find_const (parlist, param_name);
836 cpl_free(param_name);
837
838 double fieldmargin = 0.;
839 fieldmargin = cpl_parameter_get_double(par);
840
841 /* create final outgrid parameter structure */
842 //int naxis = cpl_wcs_get_image_naxis(wcs);
843
844 aParams_outputgrid =
846 ra_min, ra_max,
847 dec_min, dec_max,
848 fieldmargin);
849 //cpl_msg_info(cpl_func,"ra_max: %g ra_min: %g dx: %g",ra_max, ra_min, dx);
850 //cpl_msg_info(cpl_func,"dec_max: %g dec_min: %g dy: %g",dec_max, dec_min, dy);
851 double dist_x = (ra_max-ra_min)/dx;
852 double dist_y = (dec_max-dec_min)/dy;
853
854 cpl_msg_debug(cpl_func,"size1: %g size2: %g",dist_x, dist_y);
855 double dist = sqrt(dist_x * dist_x + dist_y * dist_y);
856 cpl_msg_debug(cpl_func,"Output cube diagonal size: %g", dist);
857 if ( dist > MAX_COMBINE_CUBE_PLANE_DIAGONAL) {
858 cpl_msg_info(cpl_func,"Internal threshold on max cube plane diagonal size: %d",MAX_COMBINE_CUBE_PLANE_DIAGONAL);
859 cpl_msg_warning(cpl_func,"Resampled cube plane diagonal %g greater than threshold. Skip cube combination.", dist);
860 cpl_error_set(cpl_func,CPL_ERROR_INCOMPATIBLE_INPUT);
861 eris_check_error_code("eris_ifu_resample_set_outputgrid2D");
862 return aParams_outputgrid;
863 }
864 cpl_free(context);
865 eris_check_error_code("eris_ifu_resample_set_outputgrid2D");
866
867 return aParams_outputgrid;
868}
869
870
871/*----------------------------------------------------------------------------*/
899cpl_error_code
900eris_ifu_resample_save_cube(hdrl_resample_result *aCube,
901 const char *procatg,
902 const char *recipe,
903 const char *filename,
904 const cpl_parameterlist *parlist,
905 cpl_frameset *frameset,
906 cpl_boolean gen_phase3)
907{
908
909 cpl_ensure_code(aCube && aCube->header, CPL_ERROR_NULL_INPUT);
910
911// cpl_error_code rc;
912 //cpl_propertylist *header = NULL;
913 cpl_propertylist* dataHeader;
914 cpl_propertylist* errsHeader;
915 cpl_propertylist* qualHeader;
916
917
918 deqErrorType errorType = rmse;
919 deqQualityType qualityType = flag16bit;
920 char* pipe_id = cpl_sprintf("%s%s%s", PACKAGE, "/", PACKAGE_VERSION);
921
922 /* Create FITS headers for the extensions */
923 if(cpl_propertylist_has(aCube->header, "CDELT3")) {
924 cpl_propertylist_erase(aCube->header, "CDELT3");
925 }
926
927 if(!cpl_propertylist_has(aCube->header, CUNIT3)) {
928 cpl_propertylist_append_string(aCube->header,"CUNIT3","um");
929 }
930 if(!cpl_propertylist_has(aCube->header, "BUNIT")) {
931 cpl_propertylist_append_string(aCube->header,"BUNIT", UNIT_ADU);
932 }
933
934
935 dataHeader = cpl_propertylist_duplicate(aCube->header);
936 errsHeader = cpl_propertylist_duplicate(aCube->header);
937 qualHeader = cpl_propertylist_duplicate(aCube->header);
938 cpl_propertylist_update_string(qualHeader,"BUNIT", "");
939
940 eris_ifu_heades_add_hduclass_common(dataHeader, "IMAGE");
942
943 eris_ifu_heades_add_hduclass_common(errsHeader, "IMAGE");
944 eris_ifu_heades_add_hduclass_errs(errsHeader, errorType);
945
946 eris_ifu_heades_add_hduclass_common(qualHeader, "IMAGE");
947 eris_ifu_heades_add_hduclass_qual(qualHeader, qualityType);
948
949 if(cpl_propertylist_has(dataHeader,"CTYPE1")) {
950 const char* ctype1 = cpl_propertylist_get_string(dataHeader,"CTYPE1");
951 if (!cpl_propertylist_has(errsHeader,"CTYPE1"))
952 cpl_propertylist_append_string(errsHeader, "CTYPE1", ctype1);
953 if (!cpl_propertylist_has(qualHeader,"CTYPE1"))
954 cpl_propertylist_append_string(qualHeader, "CTYPE1", ctype1);
955 }
956 if(cpl_propertylist_has(dataHeader,"CTYPE2")) {
957 const char* ctype2 = cpl_propertylist_get_string(dataHeader,"CTYPE2");
958 if (!cpl_propertylist_has(errsHeader,"CTYPE2"))
959 cpl_propertylist_append_string(errsHeader, "CTYPE2",ctype2);
960 if (!cpl_propertylist_has(errsHeader,"CTYPE2"))
961 cpl_propertylist_append_string(qualHeader, "CTYPE2", ctype2);
962 }
963
964 cpl_propertylist* applist = cpl_propertylist_duplicate(aCube->header);
965 if(gen_phase3) {
966 /*
967 if(cpl_propertylist_has(aCube->header,"MJD-OBS")) {
968 double mjd_obs = cpl_propertylist_get_double(aCube->header, "MJD-OBS");
969 cpl_propertylist_append_double(header, "MJD-OBS",mjd_obs);
970 }
971 */
972 eris_ifu_sdp_properties * prop = eris_ifu_sdp_properties_collect(aCube,
973 frameset, parlist, recipe);
974
975 eris_ifu_sdp_properties_update(applist, prop);
976 //cpl_propertylist_dump(applist,stdout);
978 }
979
980 /* Create an empty primary HDU with only the header information*/
981 eris_ifu_plist_erase_wcs(applist);
982 eris_pfits_clean_header(dataHeader, CPL_TRUE);
983 eris_pfits_clean_header(errsHeader, CPL_TRUE);
984 eris_pfits_clean_header(qualHeader, CPL_TRUE);
985 cpl_propertylist_append_string(applist, "PRODCATG", "SCIENCE.CUBE.IFS");
986
987 cpl_propertylist_update_string(applist, CPL_DFS_PRO_CATG, procatg);
988 if(cpl_propertylist_has(applist,"BUNIT")) {
989 cpl_propertylist_erase(applist,"BUNIT");
990 }
991 if(cpl_propertylist_has(applist,"COMMENT")) {
992 cpl_msg_warning(cpl_func,"remove comments");
993 cpl_propertylist_erase(applist,"COMMENT");
994 }
995 if(cpl_propertylist_has(dataHeader,"COMMENT")) {
996 cpl_msg_warning(cpl_func,"remove comments");
997 cpl_propertylist_erase(dataHeader,"COMMENT");
998 }
999 cpl_dfs_save_propertylist(frameset, NULL, parlist, frameset,
1000 NULL, recipe, applist, "RADECSYS", pipe_id, filename);
1001 cpl_free(pipe_id);
1002 cpl_imagelist* data = cpl_imagelist_new();
1003 cpl_imagelist* errs = cpl_imagelist_new();
1004 cpl_imagelist* qual = cpl_imagelist_new();
1005
1006 eris_ifu_split3_hdrl_imagelist(aCube->himlist, data, errs, qual);
1007
1008 /* Write the extension units */
1009 cpl_type data_type = cpl_image_get_type(cpl_imagelist_get(data,0));
1010
1011 cpl_imagelist_save(data, filename, data_type, dataHeader, CPL_IO_EXTEND);
1012 cpl_imagelist_save(errs, filename, data_type, errsHeader, CPL_IO_EXTEND);
1013 cpl_imagelist_save(qual, filename, CPL_TYPE_INT, qualHeader, CPL_IO_EXTEND);
1014 cpl_size ix_max = cpl_imagelist_get_size(data);
1015 cpl_image * q;
1016 for (cpl_size ix=ix_max-1; ix >= 0; ix--) {
1017 cpl_imagelist_unset(data, ix);
1018 cpl_imagelist_unset(errs, ix);
1019 q=cpl_imagelist_unset(qual, ix);
1020 cpl_image_delete(q);
1021 }
1022
1023 cpl_imagelist_delete(data);
1024 cpl_imagelist_delete(errs);
1025 cpl_imagelist_delete(qual);
1026
1027 cpl_propertylist_delete(applist);
1028 /* changed up to here */
1029
1030 /* this creates a problem */
1031// eris_ifu_plist_erase_wcs(header);
1032
1033
1034// if (rc != CPL_ERROR_NONE) {
1035// return cpl_error_get_code();
1036// }
1037
1038
1039// cpl_size planes = hdrl_imagelist_get_size(aCube->himlist);
1040
1041
1042 //cpl_propertylist_delete(header);
1043 eris_ifu_free_propertylist(&dataHeader);
1044 eris_ifu_free_propertylist(&errsHeader);
1045 eris_ifu_free_propertylist(&qualHeader);
1046 eris_check_error_code("eris_ifu_resample_save_cube");
1047 return cpl_error_get_code();
1048}
1049
1050
1067static cpl_wcs *
1068eris_ifu_resample_get_wcs_from_frameset(cpl_frameset* frameset,
1069 const char* procatg) {
1070
1071 cpl_ensure(frameset, CPL_ERROR_NULL_INPUT, NULL);
1072
1073 /* Read wcs from firts raw image */
1074 cpl_frameset* in_set = NULL;
1075
1076 if ((in_set = eris_ifu_extract_frameset(frameset, procatg)) == NULL){
1077 cpl_error_set_message(cpl_func, CPL_ERROR_FILE_NOT_FOUND, "Missing RAW "
1078 "files");
1079 return NULL;
1080 }
1081
1082 cpl_frame* frame = cpl_frameset_get_position(in_set, 0);
1083 const char* fname = cpl_frame_get_filename(frame);
1084
1085
1086 cpl_wcs *wcs = NULL;
1087 cpl_propertylist* head = NULL;
1088
1089 cpl_errorstate prestate = cpl_errorstate_get();
1090 head = cpl_propertylist_load(fname, 0);
1091 wcs = cpl_wcs_new_from_propertylist(head);
1092 if (wcs == NULL) {
1093 /* Not possible to read wcs - trying from extension */
1094 cpl_errorstate_set(prestate);
1095 cpl_propertylist_delete(head);
1096 } else {
1097 cpl_propertylist_delete(head);
1098 cpl_frameset_delete(in_set);
1099 return wcs;
1100 }
1101
1102 prestate = cpl_errorstate_get();
1103 head = cpl_propertylist_load(fname, 1);
1104 wcs = cpl_wcs_new_from_propertylist(head);
1105 cpl_propertylist_delete(head);
1106 cpl_frameset_delete(in_set);
1107 eris_check_error_code("eris_ifu_resample_get_wcs_from_frameset");
1108 return wcs ;
1109
1110}
1111
1112
1140static cpl_table *
1141eris_ifu_resample_frameset_to_table(const cpl_frame *data_frm,
1142 const cpl_size data_ext_id,
1143 const cpl_frame *errs_frm,
1144 const cpl_size errs_ext_id,
1145 const cpl_frame *qual_frm,
1146 const cpl_size qual_ext_id,
1147 const cpl_boolean is_variance,
1148 const cpl_boolean subtract_bkg,
1149 const int edgetrim)
1150{
1151
1152 cpl_ensure(data_frm, CPL_ERROR_NULL_INPUT, NULL);
1153 cpl_msg_severity level = cpl_msg_get_level();
1154 cpl_msg_set_level(CPL_MSG_INFO);
1155 const char *name = cpl_frame_get_filename(data_frm);
1156 cpl_imagelist *dlist = NULL;
1157 cpl_imagelist *elist = NULL;
1158 cpl_imagelist *qlist = NULL;
1159
1160 cpl_table *tab_ext = NULL;
1161
1162 cpl_errorstate prestate = cpl_errorstate_get();
1163
1164 dlist = cpl_imagelist_load(name, CPL_TYPE_DOUBLE, data_ext_id);
1165
1166 if (dlist == NULL) {
1167 /*It was not an image nor a imagelist - reset cpl error and return NULL*/
1168 if (!cpl_errorstate_is_equal(prestate)) {
1169 cpl_errorstate_set(prestate);
1170 }
1171 return NULL;
1172 }
1173
1174 cpl_propertylist *xheader_data = cpl_propertylist_load(name,
1175 data_ext_id);
1176 cpl_wcs *wcs = cpl_wcs_new_from_propertylist(xheader_data);
1177
1178 if (errs_frm != NULL) {
1179 name = cpl_frame_get_filename(errs_frm);
1180 elist = cpl_imagelist_load(name, CPL_TYPE_DOUBLE, errs_ext_id);
1181 if(is_variance) {
1182 /* transform variance into error */
1183 cpl_imagelist_power(elist, 0.5);
1184 }
1185 }
1186 if (qual_frm != NULL) {
1187 name = cpl_frame_get_filename(qual_frm);
1188 qlist = cpl_imagelist_load(name, CPL_TYPE_INT, qual_ext_id);
1189 }
1190
1191 cpl_size size = cpl_imagelist_get_size(dlist);
1192 /* ingest pixel quality in data */
1193
1194 if (qual_frm != NULL && size > 0){
1195 for(cpl_size k = 0; k < size; k++) {
1196 cpl_image* data = cpl_imagelist_get(dlist, k);
1197 cpl_image* qual = cpl_imagelist_get(qlist, k);
1198
1199 /*we use INT_MAX instead of 1.1 as some pipeline
1200 * may use pixel codes as qualifier */
1201 cpl_mask* mask = cpl_mask_threshold_image_create(qual, 0, INT_MAX);
1202
1203 cpl_image_reject_from_mask(data, mask);
1204 cpl_mask_delete(mask);
1205 cpl_imagelist_set(dlist, data, k);
1206 }
1207 }
1208
1209 /* In case the error is passed as variance we take the square root. This could
1210 * cause on some data to add a bad pixel mask - thus we have to synchronize it
1211 * */
1212 if(elist != NULL && is_variance){
1213 for(cpl_size k = 0; k < size; k++) {
1214 cpl_image* data = cpl_imagelist_get(dlist, k);
1215 cpl_mask* data_mask = cpl_image_get_bpm(data);
1216
1217 cpl_image* error = cpl_imagelist_get(elist, k);
1218 cpl_mask* error_mask = cpl_image_get_bpm(error);
1219
1220 cpl_mask* merged = cpl_mask_duplicate(data_mask);
1221
1222 /*Add original bad pixels to previous iteration*/
1223 cpl_mask_or(merged, error_mask);
1224 cpl_image_reject_from_mask(data, merged);
1225 cpl_image_reject_from_mask(error, merged);
1226 cpl_mask_delete(merged);
1227
1228 }
1229 }
1230
1231 hdrl_imagelist* hlist = NULL;
1232 if (size > 0) {
1233 hlist = hdrl_imagelist_create(dlist, elist);
1234 }
1235 int edge_trim = edgetrim;
1236 if (edge_trim > 0) {
1237 cpl_msg_info(cpl_func, "Trim input image edges of %d pixels", edge_trim);
1238 /* trim each product plane edge of edgetrim pixels */
1239 cpl_size sx = hdrl_imagelist_get_size_x(hlist);
1240 cpl_size sy = hdrl_imagelist_get_size_y(hlist);
1241 cpl_size sz = hdrl_imagelist_get_size(hlist);
1242 if (edge_trim >= 0.5* sx) {
1243 edge_trim = 0;
1244 cpl_msg_warning(cpl_func, "edge-trim must be smaller than half image "
1245 "size. Reset to 0");
1246 }
1247 if (edge_trim >= 0.5* sy) {
1248 edge_trim = 0;
1249 cpl_msg_warning(cpl_func, "edge-trim must be smaller than half image "
1250 "size. Reset to 0");
1251 }
1252 for(cpl_size k = 0; k < sz; k++) {
1253
1254 hdrl_image* hima = hdrl_imagelist_get(hlist, k);
1255 /* Note FITS convention to loop over pixels */
1256 /* trim lower image border along X direction */
1257 for(cpl_size j = 1; j <= edge_trim; j++) {
1258 for(cpl_size i = 1; i <= sx; i++) {
1259 hdrl_image_reject(hima, i, j);
1260 }
1261 }
1262 /* trim upper image border along X direction */
1263 for(cpl_size j = sy - edge_trim + 1; j <= sy; j++) {
1264 for(cpl_size i = 1; i <= sx; i++) {
1265 hdrl_image_reject(hima, i, j);
1266 }
1267 }
1268 /* trim left image border along Y direction */
1269 for(cpl_size j = 1; j <= sy; j++) {
1270 for(cpl_size i = 1; i <= edge_trim; i++) {
1271 hdrl_image_reject(hima, i, j);
1272 }
1273 }
1274 /* trim right image border along Y direction */
1275 for(cpl_size j = 1; j <= sy; j++) {
1276 for(cpl_size i = sx - edge_trim + 1; i <= sx; i++) {
1277 hdrl_image_reject(hima, i, j);
1278 }
1279 }
1280 }
1281 }
1282
1283 /* single list are not anymore needed as data are copied to the hdrl_imagelist
1284 * therefore we free them */
1285 cpl_imagelist_delete(dlist);
1286 if (qual_frm != NULL)
1287 cpl_imagelist_delete(qlist);
1288 if (errs_frm != NULL)
1289 cpl_imagelist_delete(elist);
1290
1291 /* Do the calculation */
1292 if (size == 1){
1293 /* Subtract the background if we have an image and not a cube */
1294 cpl_msg_info(cpl_func, "Reading the image ...");
1295 hdrl_image* hima = hdrl_imagelist_get(hlist, 0);
1296
1297 /* Subtract the background on request */
1298 if(subtract_bkg == CPL_TRUE) {
1299 cpl_msg_info(cpl_func, "Subtracting the median as requested ...");
1301 }
1302
1304 tab_ext = hdrl_resample_image_to_table(hima, wcs);
1305 } else {
1306 cpl_msg_info(cpl_func, "Converting imagelist to table with hdrl function");
1307 tab_ext = hdrl_resample_imagelist_to_table(hlist, wcs);
1308 }
1309
1310 /*Cleanup the memory */
1311 if (hlist != NULL)
1312 hdrl_imagelist_delete(hlist);
1313 cpl_wcs_delete(wcs);
1314 cpl_propertylist_delete(xheader_data);
1315 eris_check_error_code("eris_ifu_resample_frameset_to_table");
1316 cpl_msg_set_level(level);
1317 return tab_ext;
1318
1319}
1320
1321
1322/*----------------------------------------------------------------------------*/
1336/*----------------------------------------------------------------------------*/
1337
1338static cpl_error_code
1339eris_ifu_resample_wcs_print(cpl_wcs *wcs)
1340{
1341 cpl_ensure_code(wcs, CPL_ERROR_NULL_INPUT);
1342
1343 const cpl_array *crval = cpl_wcs_get_crval(wcs);
1344 const cpl_array *crpix = cpl_wcs_get_crpix(wcs);
1345 const cpl_array *ctype = cpl_wcs_get_ctype(wcs);
1346 const cpl_array *cunit = cpl_wcs_get_cunit(wcs);
1347
1348 const cpl_matrix *cd = cpl_wcs_get_cd(wcs);
1349 const cpl_array *dims = cpl_wcs_get_image_dims(wcs);
1350 cpl_size naxis = cpl_wcs_get_image_naxis(wcs);
1351
1352 cpl_msg_info(cpl_func, "NAXIS: %lld", naxis);
1353 int testerr = 0;
1354
1355 cpl_msg_indent_more();
1356 /* Check NAXIS */
1357 for (cpl_size i = 0; i < naxis; i++) {
1358 cpl_msg_info(cpl_func, "NAXIS%lld: %d", i + 1,
1359 cpl_array_get_int(dims, i, &testerr));
1360 }
1361 cpl_msg_indent_less();
1362
1363 double cd11 = cpl_matrix_get(cd, 0, 0);
1364 double cd12 = cpl_matrix_get(cd, 0, 1);
1365 double cd21 = cpl_matrix_get(cd, 1, 0);
1366 double cd22 = cpl_matrix_get(cd, 1, 1);
1367 double crpix1 = cpl_array_get_double(crpix, 0, &testerr);
1368 double crpix2 = cpl_array_get_double(crpix, 1, &testerr);
1369 double crval1 = cpl_array_get_double(crval, 0, &testerr);
1370 double crval2 = cpl_array_get_double(crval, 1, &testerr);
1371
1372 cpl_msg_info(cpl_func, "1st and 2nd dimension");
1373 cpl_msg_indent_more();
1374 cpl_msg_info(cpl_func, "CD1_1: %g", cd11);
1375 cpl_msg_info(cpl_func, "CD1_2: %g", cd12);
1376 cpl_msg_info(cpl_func, "CD2_1: %g", cd21);
1377 cpl_msg_info(cpl_func, "CD2_2: %g", cd22);
1378
1379 cpl_msg_info(cpl_func, "CRPIX1: %g", crpix1);
1380 cpl_msg_info(cpl_func, "CRPIX2: %g", crpix2);
1381 cpl_msg_info(cpl_func, "CRVAL1: %f", crval1);
1382 cpl_msg_info(cpl_func, "CRVAL2: %f", crval2);
1383 if (ctype) {
1384 cpl_msg_info(cpl_func, "CTYPE1: %s", cpl_array_get_string(ctype, 0));
1385 cpl_msg_info(cpl_func, "CTYPE2: %s", cpl_array_get_string(ctype, 1));
1386 }
1387
1388 if (cunit) {
1389 cpl_msg_info(cpl_func, "CUNIT1: %s", cpl_array_get_string(cunit, 0));
1390 cpl_msg_info(cpl_func, "CUNIT2: %s", cpl_array_get_string(cunit, 1));
1391 }
1392 cpl_msg_indent_less();
1393
1394 /* Is it a 3D cube or a 2D image */
1395 cpl_size cdncol = cpl_matrix_get_ncol(cd);
1396 if (cdncol == 3) {
1397
1398 double cd13 = cpl_matrix_get(cd, 0, 2);
1399 double cd23 = cpl_matrix_get(cd, 1, 2);
1400 double cd31 = cpl_matrix_get(cd, 2, 0);
1401 double cd32 = cpl_matrix_get(cd, 2, 1);
1402 double cd33 = cpl_matrix_get(cd, 2, 2);
1403 double crval3 = cpl_array_get_double(crval, 2, &testerr);
1404 double crpix3 = cpl_array_get_double(crpix, 2, &testerr);
1405
1406 cpl_msg_info(cpl_func, "3rd dimension");
1407 cpl_msg_indent_more();
1408 cpl_msg_info(cpl_func, "CD1_3: %g", cd13);
1409 cpl_msg_info(cpl_func, "CD2_3: %g", cd23);
1410 cpl_msg_info(cpl_func, "CD3_1: %g", cd31);
1411 cpl_msg_info(cpl_func, "CD3_2: %g", cd32);
1412 cpl_msg_info(cpl_func, "CD3_3: %g", cd33);
1413
1414 cpl_msg_info(cpl_func, "CRPIX3: %g", crpix3);
1415 cpl_msg_info(cpl_func, "CRVAL3: %g", crval3);
1416
1417 if (ctype) {
1418 cpl_msg_info(cpl_func, "CTYPE3: %s", cpl_array_get_string(ctype, 2));
1419 }
1420
1421 if (cunit) {
1422 cpl_msg_info(cpl_func, "CUNIT3: %s", cpl_array_get_string(cunit, 2));
1423 }
1424
1425 cpl_msg_indent_less();
1426 }
1427 eris_check_error_code("eris_ifu_resample_wcs_print");
1428 return cpl_error_get_code();
1429}
1430
1431
1455static cpl_table*
1456eris_ifu_resample_get_table_from_frameset(const cpl_frame* data_frm,
1457 const cpl_frame* errs_frm,
1458 const cpl_frame* qual_frm,
1459 const cpl_size data_ext_id,
1460 const cpl_size errs_ext_id,
1461 const cpl_size qual_ext_id,
1462 const cpl_boolean is_variance,
1463 const cpl_boolean subtract_bkg,
1464 const int edge_trim) {
1465
1466 cpl_ensure(data_frm, CPL_ERROR_NULL_INPUT, NULL);
1467 /* Assumption - all have the same number of extensions */
1468 cpl_size next = cpl_frame_get_nextensions(data_frm);
1469 cpl_msg_info(cpl_func, "Analysing and processing file %s",
1470 cpl_frame_get_filename(data_frm));
1471
1472 /* Assumption
1473 * if data_ext_id == -1 the loop can be done from 0 to next for data, error,
1474 * and quality; only the filename changes
1475 * if data_ext_id != -1 no loop is needed and is difficult to do in the loop
1476 * */
1477
1478 cpl_msg_indent_more();
1479 cpl_table* table = NULL;
1480 if(data_ext_id != -1) {
1481 cpl_msg_info(cpl_func, "Extension: %02lld", data_ext_id);
1482 /* No Loop, only use the extension given by the user */
1483 table = eris_ifu_resample_frameset_to_table(data_frm, data_ext_id,
1484 errs_frm, errs_ext_id,
1485 qual_frm, qual_ext_id,
1486 is_variance, subtract_bkg, edge_trim);
1487 } else {
1488 /* Loop */
1489 for (cpl_size i = 0; i < next + 1; i++ ) {
1490 cpl_table * table_local = NULL;
1491 cpl_msg_info(cpl_func,"Extension: %02lld", i);
1492 table_local = eris_ifu_resample_frameset_to_table(data_frm, i,
1493 errs_frm, i,
1494 qual_frm, i,
1495 is_variance, subtract_bkg, edge_trim);
1496
1497 /* If table_local == NULL go to the next extension */
1498 if (table_local == NULL) {
1499 cpl_msg_info(cpl_func, "No siutable data found - continuing");
1500 continue;
1501 }
1502
1503 /* Now table_local is always != NULL */
1504 if (table == NULL){
1505 table = table_local;
1506 continue;
1507 } else {
1508 cpl_size nrow = cpl_table_get_nrow(table);
1509 cpl_table_insert(table, table_local, nrow);
1510 cpl_table_delete(table_local);
1511 table_local = NULL;
1512 }
1513 }
1514 }
1515 cpl_msg_indent_less();
1516 eris_check_error_code("eris_ifu_resample_get_table_from_frameset");
1517
1518 return table;
1519}
1520
1523// * @internal
1524// * @brief Find the aperture(s) with the greatest flux
1525// * @param self The aperture object
1526// * @param ind The aperture-indices in order of decreasing flux
1527// * @param nfind Number of indices to find
1528// * @return CPL_ERROR_NONE or the relevant _cpl_error_code_ on error
1529// *
1530// * nfind must be at least 1 and at most the size of the aperture object.
1531// *
1532// * The ind array must be able to hold (at least) nfind integers.
1533// * On success the first nfind elements of ind point to indices of the
1534// * aperture object.
1535// *
1536// * To find the single ind of the aperture with the maximum flux use simply:
1537// * int ind;
1538// * eris_apertures_find_max_flux(self, &ind, 1);
1539// *
1540// */
1542//static cpl_error_code eris_apertures_find_max_flux(const cpl_apertures * self,
1543// int * ind, int nfind)
1544//{
1545// const int nsize = cpl_apertures_get_size(self);
1546// int ifind;
1547
1548
1549// cpl_ensure_code(nsize > 0, cpl_error_get_code());
1550// cpl_ensure_code(ind, CPL_ERROR_NULL_INPUT);
1551// cpl_ensure_code(nfind > 0, CPL_ERROR_ILLEGAL_INPUT);
1552// cpl_ensure_code(nfind <= nsize, CPL_ERROR_ILLEGAL_INPUT);
1553
1554// for (ifind=0; ifind < nfind; ifind++) {
1555// double maxflux = -1;
1556// int maxind = -1;
1557// int i;
1558// for (i=1; i <= nsize; i++) {
1559// int k;
1560
1561// /* The flux has to be the highest among those not already found */
1562// for (k=0; k < ifind; k++) if (ind[k] == i) break;
1563
1564// if (k == ifind) {
1565// /* i has not been inserted into ind */
1566// const double flux = cpl_apertures_get_flux(self, i);
1567
1568// if (maxind < 0 || flux > maxflux) {
1569// maxind = i;
1570// maxflux = flux;
1571// }
1572// }
1573// }
1574// ind[ifind] = maxind;
1575// }
1576// eris_check_error_code("eris_apertures_find_max_flux");
1577// return CPL_ERROR_NONE;
1578
1579//}
1580
1582// * @brief filter input mask by kernel of given sizes and mode
1583// *
1584// * @param input_mask input mask
1585// * @param kernel_nx kernel X size
1586// * @param kernel_ny kernel X size
1587// * @param filter filter mode
1588// * @return mask corresponding to what filtered
1589// * @doc
1590// *
1591// * */
1592//static cpl_mask *
1593//eris_bpm_filter(
1594// const cpl_mask * input_mask,
1595// cpl_size kernel_nx,
1596// cpl_size kernel_ny,
1597// cpl_filter_mode filter)
1598//{
1599// cpl_mask * kernel = NULL;
1600// cpl_mask * filtered_mask = NULL;
1601// cpl_mask * expanded_mask = NULL;
1602// cpl_mask * expanded_filtered_mask = NULL;
1603
1604// /* Check Entries */
1605// cpl_ensure(input_mask != NULL, CPL_ERROR_NULL_INPUT, NULL);
1606// cpl_ensure(kernel_nx >= 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
1607// cpl_ensure(kernel_ny >= 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
1608// cpl_ensure(filter == CPL_FILTER_EROSION || filter == CPL_FILTER_DILATION ||
1609// filter == CPL_FILTER_OPENING || filter == CPL_FILTER_CLOSING,
1610// CPL_ERROR_ILLEGAL_INPUT, NULL);
1611
1612// /* Only odd-sized masks allowed */
1613// cpl_ensure((kernel_nx&1) == 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
1614// cpl_ensure((kernel_ny&1) == 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
1615
1616// kernel = cpl_mask_new(kernel_nx, kernel_ny);
1617// cpl_mask_not(kernel); /* All values set to unity*/
1618
1619// /* Enlarge the original mask with the kernel size and assume that outside
1620// * all pixels are good */
1621// expanded_mask = cpl_mask_new(
1622// cpl_mask_get_size_x(input_mask) + 2 * kernel_nx,
1623// cpl_mask_get_size_y(input_mask) + 2 * kernel_ny);
1624
1625// cpl_mask_copy(expanded_mask, input_mask, kernel_nx + 1, kernel_ny +1 );
1626
1627// expanded_filtered_mask = cpl_mask_new(cpl_mask_get_size_x(expanded_mask),
1628// cpl_mask_get_size_y(expanded_mask));
1629
1630// if(cpl_mask_filter(expanded_filtered_mask, expanded_mask, kernel, filter,
1631// CPL_BORDER_ZERO) != CPL_ERROR_NONE) {
1632
1633// cpl_mask_delete(kernel);
1634// cpl_mask_delete(expanded_filtered_mask);
1635// cpl_mask_delete(expanded_mask);
1636// return NULL;
1637// }
1638
1639// /* Extract the original mask from the expanded mask */
1640// filtered_mask = cpl_mask_extract(expanded_filtered_mask,
1641// kernel_nx+1, kernel_ny + 1,
1642// cpl_mask_get_size_x(input_mask) + kernel_nx,
1643// cpl_mask_get_size_y(input_mask) + kernel_ny);
1644
1645
1646// /* Free memory */
1647// cpl_mask_delete(kernel);
1648// cpl_mask_delete(expanded_filtered_mask);
1649// cpl_mask_delete(expanded_mask);
1650// eris_check_error_code("eris_bpm_filter");
1651// return filtered_mask;
1652//}
1653
1655// * @brief generate mask corresponding to post-filtered image according to method
1656// *
1657// * @param std_med_ima input image
1658// * @param pfm post filter mask method
1659// * @param pfx post filer x size
1660// * @param pfy post filer x size
1661// * @param kappa value used in kappa-sigma clipping
1662// * @return cpl_error_code
1663// * @doc
1664// *
1665// * */
1666//static cpl_mask*
1667//eris_object_mask(cpl_image* std_med_ima, const char* pfm, const int pfx,
1668// const int pfy, const double kappa)
1669//{
1670// cpl_mask* object_mask;
1671// double mad;
1672// double median;
1673// cpl_mask* nan_mask;
1674// /* flag NAN in image (due to SPIFFIER BP conventions) before computing mad */
1675// cpl_image_reject_value(std_med_ima,CPL_VALUE_NAN);
1676// /* extract mask of NANs to be used later */
1677// nan_mask = cpl_image_get_bpm(std_med_ima);
1678// //cpl_mask_save(nan_mask,"nan_mask.fits",NULL,CPL_IO_CREATE);
1679// median = cpl_image_get_mad(std_med_ima, &mad);
1680// double low_cut = median - 10 * kappa * mad * CPL_MATH_STD_MAD;
1681// double hi_cut = median + kappa * mad * CPL_MATH_STD_MAD;
1682// /* takes as low cut a very small value to be sure to get all negative pixels
1683// * as good */
1684// low_cut = -1.e10;
1685// //cpl_msg_info(cpl_func, "median: %g mad: %g low: %g hi : %g", median, mad, low_cut, hi_cut);
1686// /* detect object (and other previously flagged values) as good pixels */
1687// object_mask = cpl_mask_threshold_image_create(std_med_ima, low_cut, hi_cut);
1688// //cpl_mask_not(object_mask);
1689
1690// //cpl_mask_save(object_mask,"object_mask.fits",NULL,CPL_IO_CREATE);
1691// /* to detect only object pixels as good one make or with previously flagged
1692// * NAN values
1693// */
1694
1695// cpl_mask_or(object_mask,nan_mask);
1696// //cpl_mask_save(object_mask,"or_mask.fits",NULL,CPL_IO_CREATE);
1697// /* increase the mask of a few pixels to be sure to include all object values
1698// * for this reason we use method dilation and we make a NOT on the image
1699// * because the DILATION enlarge the area of BAD pixels, not the GOOD ones
1700// * but we want to increase the area covered by the object (flagged as good)
1701// */
1702// cpl_filter_mode filter_mode = CPL_FILTER_EROSION ;
1703// cpl_msg_warning(cpl_func, "pfm = %s", pfm);
1704// if( strcmp(pfm,"erosion") == 0 ) {
1705// cpl_msg_info(cpl_func, "Filter erosion");
1706// filter_mode = CPL_FILTER_EROSION ;
1707// } else if( strcmp(pfm, "dilation") == 0 ) {
1708// cpl_msg_info(cpl_func, "Filter dilation");
1709// filter_mode = CPL_FILTER_DILATION ;
1710// } else if( strcmp(pfm,"closing") == 0 ) {
1711// cpl_msg_info(cpl_func, "Filter closing");
1712// filter_mode = CPL_FILTER_CLOSING ;
1713// }
1714// //cpl_filter_mode filter_mode = CPL_FILTER_DILATION ;
1715
1716// cpl_mask_not(object_mask);
1717// cpl_mask* obj_mask_filtered = eris_bpm_filter(object_mask, pfx, pfy, filter_mode);
1718// /* To have again a proper mask with object flagged as good we do a NOT */
1719// cpl_mask_not(obj_mask_filtered);
1720// //cpl_mask_save(obj_mask_filtered,"filtered_mask.fits",NULL,CPL_IO_CREATE);
1721// /* clean-up memory */
1722// cpl_mask_delete(object_mask);
1723// //cpl_mask_delete(nan_mask);
1724// eris_check_error_code("eris_object_mask");
1725// return obj_mask_filtered;
1726//}
1727
1751cpl_error_code
1753 cpl_frameset* frameset,
1754 const cpl_parameterlist* parlist,
1755 const char* recipe_name,
1756// const char* filenameSpec,
1757 cpl_boolean apply_flat,
1758 cpl_boolean is_pupil)
1759{
1760 hdrl_image* cube_collapsed = NULL;
1761 cpl_image* cube_cmap = NULL;
1762
1763 cpl_frame* frame = cpl_frameset_find(frameset, pro_catg);
1764 const char* cube_fname = cpl_frame_get_filename(frame);
1765 cpl_propertylist* pheader = cpl_propertylist_load(cube_fname, 0);
1766 cpl_propertylist* head_wcs_2d = eris_ifu_plist_extract_wcs2D(pheader);
1767 cpl_propertylist* dheader = cpl_propertylist_load(cube_fname, 1);
1768 cpl_propertylist* eheader = cpl_propertylist_load(cube_fname, 2);
1769 cpl_propertylist* qheader = cpl_propertylist_load(cube_fname, 3);
1770
1771 cpl_imagelist* iml_data = cpl_imagelist_load(cube_fname, CPL_TYPE_DOUBLE, 1);
1772 cpl_imagelist* iml_errs = cpl_imagelist_load(cube_fname, CPL_TYPE_DOUBLE, 2);
1773 cpl_imagelist* iml_qual = cpl_imagelist_load(cube_fname, CPL_TYPE_DOUBLE, 3);
1774
1775 cpl_propertylist_append(dheader, head_wcs_2d);
1776 cpl_propertylist_append(eheader, head_wcs_2d);
1777 cpl_propertylist_append(qheader, head_wcs_2d);
1778 cpl_propertylist_delete(head_wcs_2d);
1779 double dit = eris_pfits_get_dit(pheader);
1780 int ndit = eris_pfits_get_ndit(pheader);
1781 double exptime = dit * ndit;
1782
1783 cpl_size sz = cpl_imagelist_get_size(iml_data);
1784 cpl_image* data = NULL;
1785 cpl_image* errs = NULL;
1786 cpl_image* qual = NULL;
1787 cpl_mask* bpm_data = NULL;
1788 cpl_mask* bpm_errs = NULL;
1789 cpl_boolean edge_trim = CPL_TRUE;
1790 cpl_size k_center = 0.5 * sz;
1791 for(cpl_size k = 0; k < sz; k++) {
1792 data = cpl_imagelist_get(iml_data, k);
1793 errs = cpl_imagelist_get(iml_errs, k);
1794 qual = cpl_imagelist_get(iml_qual, k);
1795 cpl_image_reject_value(data, CPL_VALUE_NAN);
1796 cpl_image_reject_value(errs, CPL_VALUE_NAN);
1797 bpm_data = cpl_image_get_bpm(data);
1798 bpm_errs = cpl_image_get_bpm(errs);
1799 cpl_mask * bpm_mask = cpl_mask_threshold_image_create(qual, -0.5, 0.5) ;
1800 cpl_mask_not(bpm_mask) ;
1801 cpl_mask_or(bpm_mask,bpm_data);
1802 cpl_mask_or(bpm_mask,bpm_errs);
1803 cpl_image_reject_from_mask(data, bpm_mask);
1804 cpl_image_reject_from_mask(errs, bpm_mask);
1805
1806 //EXPOSURE MAP CREATION
1807
1808
1809 /* exposure map: use the plane at the center of the wavelength range
1810 * then replace intensity with EXPTIME value and error with 1 and
1811 * generate HDRL image. This will be re-sampled as the cube and give
1812 * the exposure map associated to the combined cube
1813 */
1814
1815 cpl_image* ima_time = NULL;
1816 cpl_image* err_time = NULL;
1817 hdrl_image *himg_exptime = NULL;
1818 if(apply_flat && k == k_center) {
1819 ima_time = cpl_image_new(cpl_image_get_size_x(data),
1820 cpl_image_get_size_y(data), CPL_TYPE_DOUBLE);
1821 err_time = cpl_image_new(cpl_image_get_size_x(data),
1822 cpl_image_get_size_y(data), CPL_TYPE_DOUBLE);
1823
1824
1825 cpl_image_add_scalar(ima_time, exptime);
1826 cpl_image_add_scalar(err_time, sqrt(exptime));
1827 himg_exptime = hdrl_image_create(ima_time, err_time);
1828 hdrl_image_reject_from_mask(himg_exptime, bpm_mask);
1829
1830 if (edge_trim > 0) {
1831 eris_ifu_resample_trim_edge(himg_exptime, edge_trim);
1832 }
1833 /*
1834 char* fname = cpl_sprintf("ima_time_%lld.fits",i);
1835 cpl_image_save(ima_time,fname,CPL_TYPE_DOUBLE,NULL,CPL_IO_DEFAULT);
1836 cpl_free(fname);
1837 fname = cpl_sprintf("err_time_%lld.fits",i);
1838 cpl_image_save(err_time,fname,CPL_TYPE_DOUBLE,NULL,CPL_IO_DEFAULT);
1839 cpl_free(fname);
1840 */
1841 cpl_image_delete(ima_time);
1842 cpl_image_delete(err_time);
1843
1844 /* SAVE EXPOSURE MAP */
1845
1846 cpl_propertylist* phead2D = cpl_propertylist_duplicate(dheader);
1847 char* fname = cpl_sprintf("%s_exposure_map.fits", recipe_name);
1848 eris_ifu_plist_erase_wcs3D(phead2D);
1849 eris_ifu_plist_erase_expmap_extra_keys(phead2D);
1850 cpl_propertylist_update_string(phead2D, CPL_DFS_PRO_CATG, ERIS_IFU_PRO_JITTER_EXPOSURE_MAP);
1851 cpl_propertylist_update_string(phead2D, "PRODCATG", PRODCATG_EXPOSUREMAP);
1852 eris_ifu_save_image_phase3(frameset, /*NULL, */parlist,frameset, /*NULL,*/
1853 recipe_name, phead2D, "RADECSYS", fname,
1854 hdrl_image_get_image(himg_exptime), CPL_TYPE_DOUBLE,
1855 UNIT_TIME);
1856 hdrl_image_delete(himg_exptime);
1857 cpl_propertylist_delete(phead2D);
1858 cpl_free(fname);
1859 }
1860
1861
1862
1863 cpl_mask_delete(bpm_mask) ;
1864 }
1865
1866 hdrl_imagelist* himlist = hdrl_imagelist_create(iml_data, iml_errs);
1867
1868 hdrl_imagelist_collapse_mean(himlist, &cube_collapsed, &cube_cmap);
1869 hdrl_imagelist_delete(himlist);
1870
1871 /* TODO: check header is correct */
1872 cpl_propertylist* phead2D = cpl_propertylist_duplicate(pheader);
1873 cpl_propertylist_erase_regexp(phead2D,NAXIS3,0);
1874 cpl_propertylist_erase_regexp(phead2D,CTYPE3,0);
1875 cpl_propertylist_erase_regexp(phead2D,CRVAL3,0);
1876 cpl_propertylist_erase_regexp(phead2D,CRPIX3,0);
1877 cpl_propertylist_erase_regexp(phead2D,CUNIT3, 0);
1878 cpl_propertylist_erase_regexp(phead2D,"CD3_*",0);
1879 cpl_propertylist_erase_regexp(phead2D,CD1_3,0);
1880 cpl_propertylist_erase_regexp(phead2D,CD2_3,0);
1881
1882 char* fname = NULL;
1883 char* prefix = NULL;
1884 char* suffix = NULL;
1885 if(apply_flat) {
1886 suffix = cpl_sprintf("%s","");
1887 } else {
1888 suffix = cpl_sprintf("%s","_no_flat");
1889 }
1890 if(strstr(recipe_name, "jitter") != NULL ) {
1891 prefix = cpl_sprintf("eris_ifu_jitter%s", suffix);
1892 } else {
1893 prefix = cpl_sprintf("eris_ifu_stdstar%s", suffix);
1894 }
1895 cpl_free(suffix);
1896 const char* pcatg = NULL;
1897 if(strstr(pro_catg, "OBJ") != NULL ) {
1898 fname = cpl_sprintf("%s_obj_cube_mean.fits", prefix);
1899 if(apply_flat) {
1900 pcatg = ERIS_IFU_PRO_JITTER_OBJ_CUBE_MEAN;
1901 } else {
1902 pcatg = ERIS_IFU_PRO_JITTER_OBJ_CUBE_NOFLAT_MEAN;
1903 }
1904 } else if(strstr(pro_catg, ERIS_IFU_PRO_JITTER_DAR_CUBE) != NULL ) {
1905 fname = cpl_sprintf("%s_dar_cube_mean.fits", prefix);
1906 if(apply_flat) {
1907 pcatg = ERIS_IFU_PRO_JITTER_OBJ_DAR_CUBE_MEAN;
1908 } else {
1909 pcatg = ERIS_IFU_PRO_JITTER_OBJ_CUBE_NOFLAT_MEAN;
1910 }
1911 } else if(strstr(pro_catg, ERIS_IFU_PRO_JITTER_TWK_CUBE) != NULL ) {
1912 fname = cpl_sprintf("%s_twk_cube_mean.fits", prefix);
1913 if(apply_flat) {
1914 pcatg = ERIS_IFU_PRO_JITTER_TWK_CUBE_MEAN;
1915 } else {
1916 pcatg = ERIS_IFU_PRO_JITTER_OBJ_CUBE_NOFLAT_MEAN;
1917 }
1918 } else if(strstr(pro_catg, "STD") != NULL ) {
1919 fname = cpl_sprintf("%s_std_cube_mean.fits", prefix);
1920 if(apply_flat) {
1921 pcatg = ERIS_IFU_PRO_JITTER_STD_CUBE_MEAN;
1922 } else {
1923 pcatg = ERIS_IFU_PRO_JITTER_STD_CUBE_NOFLAT_MEAN;
1924 }
1925 } else if(strstr(pro_catg, "PSF") != NULL ) {
1926 fname = cpl_sprintf("%s_psf_cube_mean.fits", prefix);
1927 if(apply_flat) {
1928 pcatg = ERIS_IFU_PRO_JITTER_PSF_CUBE_MEAN;
1929 } else {
1930 pcatg = ERIS_IFU_PRO_JITTER_PSF_CUBE_NOFLAT_MEAN;
1931 }
1932 }
1933 cpl_free(prefix);
1934 cpl_propertylist_update_string(phead2D, CPL_DFS_PRO_CATG, pcatg);
1935 /* remove WCS from primary header as data go into extensions */
1936// eris_ifu_plist_erase_wcs(phead2D);
1937
1938
1939 if( is_pupil) {
1940 /* compute PUPIL X/Y SHIFT for QC */
1941// hdrl_imagelist* himl = NULL;
1942
1943 cpl_table* qclog_tbl = eris_qclog_init();
1944 cpl_image* img = hdrl_image_get_image(cube_collapsed);
1945 cpl_size max_ima_x = 0;
1946 cpl_size max_ima_y = 0;
1947 cpl_image_get_maxpos(img, &max_ima_x, &max_ima_y);
1948 cpl_size sx = hdrl_image_get_size_x(cube_collapsed);
1949 cpl_size sy = hdrl_image_get_size_y(cube_collapsed);
1950 double max_ima_cx = cpl_image_get_centroid_x_window(img, 1, 1, sx, sy);
1951 double max_ima_cy = cpl_image_get_centroid_y_window(img, 1, 1, sx, sy);
1952 double xshift = max_ima_cx - (double) sx * 0.5;
1953 double yshift = max_ima_cy - (double) sy * 0.5;
1954 char* key_name = cpl_sprintf("QC PUPIL%d SHIFTX", 0);
1955
1956 eris_qclog_add_double_f(qclog_tbl, key_name, xshift,
1957 "[pix] X shift centroid - center image");
1958 cpl_free(key_name);
1959
1960 key_name = cpl_sprintf("QC PUPIL%d SHIFTY", 0);
1961 eris_qclog_add_double_f(qclog_tbl, key_name, yshift,
1962 "[pix] Y shift centroid - center image");
1963 cpl_msg_info(cpl_func, "xshift: %g yshift: %g", xshift, yshift);
1964
1965 cpl_free(key_name);
1966 eris_pfits_put_qc(phead2D, qclog_tbl);
1967 cpl_free(qclog_tbl);
1968 }
1969
1970 cpl_image* mask_ima = cpl_image_new_from_mask(hdrl_image_get_mask(cube_collapsed));
1971 eris_ifu_save_deq_image(frameset, NULL, parlist,frameset, NULL,
1972 recipe_name, phead2D, "RADECSYS", fname,
1973 hdrl_image_get_image(cube_collapsed),
1974 hdrl_image_get_error(cube_collapsed),
1975 rmse, mask_ima, flag16bit, UNIT_ADU);
1976
1977 cpl_image_delete(mask_ima);
1978 cpl_free(fname);
1980 cpl_propertylist_delete(dheader);
1981 cpl_propertylist_delete(eheader);
1982 cpl_propertylist_delete(qheader);
1983 cpl_propertylist_delete(pheader);
1984 cpl_imagelist_delete(iml_data);
1985 cpl_imagelist_delete(iml_errs);
1986 cpl_imagelist_delete(iml_qual);
1987 hdrl_image_delete(cube_collapsed);
1988 cpl_image_delete(cube_cmap);
1989 eris_check_error_code("eris_ifu_cube_collapse_mean_and_save");
1990
1991 return cpl_error_get_code();
1992}
1993
1994/* Helper function for computing maximum distance between cube centers */
1995static double
1996eris_ifu_compute_max_cubes_dist(cpl_frameset* cube_set)
1997{
1998 /* First check if cube size is too large NOT YET VERIFIED */
1999 double dist_max = -1;
2000 double dist_sqr = -1;
2001 cpl_size ncubes = cpl_frameset_get_size(cube_set);
2002 cpl_vector* vec_ra = cpl_vector_new(ncubes);
2003 cpl_vector* vec_dec = cpl_vector_new(ncubes);
2004 cpl_propertylist* plist = NULL;
2005 const char* name = NULL;
2006 double ra = 0, dec = 0;
2007 double ra1 = 0, dec1 = 0;
2008 double ra2 = 0, dec2 = 0;
2009 double dist = 0;
2010 //double cd1_1 = 0;
2011 //double cd2_2 = 0;
2012 double pix_size = 0;
2013 //const char* scale = NULL;
2014 cpl_propertylist* phead = NULL;
2015 //TODO: why cd3_3 is set and not used?
2016// double cd3_3 = 0;
2017 int extnum_raw = 1;
2018 cpl_frame* data_frame = NULL;
2019
2020 data_frame = cpl_frameset_get_position(cube_set, 0);
2021 name = cpl_frame_get_filename(data_frame);
2022 plist = cpl_propertylist_load(name, extnum_raw);
2023 phead = cpl_propertylist_load(name, 0);
2024// cd1_1 = cpl_propertylist_get_double(plist, "CD1_1");
2025// cd2_2 = cpl_propertylist_get_double(plist, "CD2_2");
2026// cd3_3 = cpl_propertylist_get_double(plist, "CD3_3");
2027
2028 ifsPreopticsScale scale_mas = eris_ifu_get_preopticsScale(phead);
2029
2030 cpl_propertylist_delete(plist);
2031 cpl_propertylist_delete(phead);
2032
2033 switch(scale_mas){
2034 case S250MAS: pix_size = 0.250; break;
2035 case S100MAS: pix_size = 0.100; break;
2036 case S25MAS: pix_size = 0.025; break;
2037 case PUPIL: pix_size = 0.025; break;
2038 case UNDEFINED_SCALE: pix_size = 0.100; break;
2039 }
2040
2041 pix_size /= 3600;
2042
2043 /*
2044 cpl_msg_warning(cpl_func, "cd1_1: %g cd2_2: %g cd3_3: %g",
2045 cd1_1, cd2_2, cd3_3);
2046 */
2047 for (cpl_size i = 0; i < ncubes ; i++) {
2048
2049 data_frame = cpl_frameset_get_position(cube_set, i);
2050 name = cpl_frame_get_filename(data_frame);
2051 plist = cpl_propertylist_load(name, extnum_raw);
2052
2053 ra = cpl_propertylist_get_double(plist, "RA");
2054 dec = cpl_propertylist_get_double(plist, "DEC");
2055 //cpl_msg_warning(cpl_func, "ra: %13.10g dec: %13.10g",ra,dec);
2056 cpl_vector_set(vec_ra, i, ra);
2057 cpl_vector_set(vec_dec, i, dec);
2058
2059 cpl_propertylist_delete(plist);
2060 }
2061 for (cpl_size i = 0; i < ncubes; i++) {
2062 ra1 = cpl_vector_get(vec_ra, i);
2063 dec1 = cpl_vector_get(vec_dec, i);
2064 for (cpl_size j = 0; j < ncubes && i != j; j++) {
2065 ra2 = cpl_vector_get(vec_ra, j);
2066 dec2 = cpl_vector_get(vec_dec, j);
2067 /*
2068 cpl_msg_warning(cpl_func, "ra1: %13.10g ra2: %13.10g",ra1, ra2);
2069 cpl_msg_warning(cpl_func, "dist_ra: %13.10g",(ra2 - ra1)/pix_size);
2070 cpl_msg_warning(cpl_func, "dec1: %13.10g dec2: %13.10g",dec1, dec2);
2071 cpl_msg_warning(cpl_func, "dist_dec: %13.10g",(dec2 - dec1)/pix_size);
2072 */
2073 dist = (ra2 - ra1) * (ra2 - ra1) / pix_size / pix_size +
2074 (dec2 - dec1) * (dec2 - dec1) / pix_size / pix_size;
2075 //cpl_msg_warning(cpl_func, "dist: %13.10g",dist);
2076 if(dist > dist_sqr) {
2077 dist_sqr = dist;
2078 }
2079 }
2080 }
2081 cpl_vector_delete(vec_ra);
2082 cpl_vector_delete(vec_dec);
2083 dist_max = sqrt(dist_sqr);
2084 cpl_msg_info(cpl_func, "Max distance between contributing cubes centers: %g",dist_max);
2085
2086 eris_check_error_code("eris_ifu_compute_max_cubes_dist");
2087 return dist_max;
2088
2089}
2090
2091/*
2092static cpl_error_code
2093eris_resample_compute_size(hdrl_parameter *aParams_outputgrid,
2094 int *aX, int *aY, int *aZ)
2095{
2096 cpl_ensure_code(aParams_outputgrid && aX && aY && aZ, CPL_ERROR_NULL_INPUT);
2097 const char func[] = "hdrl_resample_compute_size";
2098 double x1, y1, x2, y2;
2099
2100 double ramin = aParams_outputgrid->ra_min;
2101 double ramax = aParams_outputgrid->ra_max;
2102 double decmin = aParams_outputgrid->dec_min;
2103 double decmax = aParams_outputgrid->dec_max;
2104 hdrl_resample_wcs_projplane_from_celestial(aParams_outputgrid, ramin, decmin,
2105 &x1, &y1);
2106 hdrl_resample_wcs_projplane_from_celestial(aParams_outputgrid, ramax, decmax,
2107 &x2, &y2);
2108 *aX = lround(fabs(x2 - x1) / aParams_outputgrid->delta_ra) + 1;
2109 *aY = lround(fabs(y2 - y1) / aParams_outputgrid->delta_dec) + 1;
2110
2111 double lmin = aParams_outputgrid->lambda_min;
2112 double lmax = aParams_outputgrid->lambda_max;
2113
2114 *aZ = (int)ceil((lmax - lmin) / aParams_outputgrid->delta_lambda) + 1;
2115
2116 cpl_msg_debug(func, "Output cube size %d x %d x %d (fit to data)",
2117 *aX, *aY, *aZ);
2118 return CPL_ERROR_NONE;
2119}
2120*/
2141cpl_error_code
2142eris_ifu_combine(cubeType obj_type, cpl_frameset* frameset,
2143 const cpl_parameterlist * parlist,
2144 const char* recipe_name,
2145 const char* pipefile_prefix)
2146{
2147 cpl_size nframes = cpl_frameset_get_size(frameset);
2148
2149 cpl_frameset* cube_set = NULL;
2150 char *proCatg = NULL;
2151 char *filenameSpec = NULL;
2152 cpl_msg_info(cpl_func,"Combine cubes into a single one, compute mean, median along wavelength axis");
2153 eris_ifu_jitter_get_procatg_and_filename(obj_type, &proCatg, &filenameSpec);
2154
2155 if ((cube_set = eris_ifu_extract_frameset(frameset, proCatg )) == NULL){
2156
2157 cpl_msg_warning(cpl_func,"No %s files", proCatg);
2158 cpl_error_set_message(cpl_func, CPL_ERROR_FILE_NOT_FOUND, "Missing CUBE "
2159 "files");
2160 return cpl_error_get_code();
2161 } else {
2162 cpl_msg_info(cpl_func, "%02d file(s) of type %s",
2163 (int)cpl_frameset_get_size(cube_set), proCatg);
2164 }
2165
2166 /*convert the images in a table of predefined structure*/
2167 cpl_table* restable = NULL; /*Final table to resample*/
2168 cpl_table* tab_tmp = NULL;
2169 const cpl_parameter * par = NULL;
2170 cpl_frame * data_frame = NULL;
2171 cpl_frame * errs_frame = NULL;
2172 cpl_frame * qual_frame = NULL;
2173 int extnum_raw = 1;
2174 int extnum_err = 2;
2175 int extnum_bpm = 3;
2176 cpl_boolean is_variance = CPL_TRUE;
2177 cpl_boolean subtract_bkg = CPL_FALSE;
2178 int edge_trim = 0;
2179 /*Loops over all frames and adds the images/imagelists to the final table */
2180 cpl_size ncubes = cpl_frameset_get_size(cube_set);
2181 char* pname = cpl_sprintf("eris.%s.max-cubes-centres-dist",recipe_name);
2182 int max_cubes_dist = 20;
2183 //cpl_parameterlist_dump(parlist,stdout);
2184
2185 max_cubes_dist = cpl_parameter_get_int(
2186 cpl_parameterlist_find_const(parlist, pname));
2187
2188 if ( max_cubes_dist < eris_ifu_compute_max_cubes_dist(cube_set)) {
2189 cpl_msg_info(cpl_func,"max-cubes-centres-dist: %d",max_cubes_dist);
2190 cpl_msg_warning(cpl_func,"max distance between cube centres greater than threshold. Skip cube combination.");
2191 return CPL_ERROR_INCOMPATIBLE_INPUT;
2192 }
2193
2194 for (cpl_size i = 0; i < ncubes ; i++) {
2195
2196 data_frame = cpl_frameset_get_position(cube_set, i);
2197 errs_frame = data_frame;
2198 qual_frame = data_frame;
2199
2200 tab_tmp = eris_ifu_resample_get_table_from_frameset(data_frame,
2201 errs_frame,
2202 qual_frame, extnum_raw,
2203 extnum_err, extnum_bpm,
2204 is_variance, subtract_bkg, edge_trim);
2205
2206
2207 if(nframes == 1) {
2208 restable = tab_tmp;
2209 } else {
2210 eris_ifu_resample_update_table(tab_tmp, &restable);
2211 cpl_table_delete(tab_tmp);
2212 }
2213
2214
2215 }
2216 cpl_ensure_code(restable, CPL_ERROR_NULL_INPUT);
2217
2218 /* Save the final table if required */
2219 eris_ifu_resample_tablesave (recipe_name, par, parlist, frameset, restable,
2220 pipefile_prefix);
2221 /* Getting the input wcs needed for determine the scaling of output cube */
2222 cpl_wcs *wcs = eris_ifu_resample_get_wcs_from_frameset(cube_set, proCatg);
2223 cpl_ensure_code(wcs, CPL_ERROR_NULL_INPUT);
2224
2225 cpl_msg_info(cpl_func, "WCS used for the output grid definition and passed to "
2226 "hdrl_resample_compute():");
2227 cpl_msg_indent_more();
2228 eris_ifu_resample_wcs_print(wcs);
2229 cpl_msg_indent_less();
2230
2231 hdrl_parameter *aParams_method = NULL;
2232 /* set re-sampling method */
2233
2234 char *context = NULL;
2235 context = cpl_sprintf("eris.%s", recipe_name);
2236 aParams_method = eris_ifu_resample_set_method(parlist, context);
2237 cpl_free(context);
2238
2239 cpl_ensure_code(aParams_method, CPL_ERROR_NULL_INPUT);
2240
2241 /* set re-sampling outputgrid */
2242 hdrl_parameter *aParams_outputgrid = NULL;
2243
2244 aParams_outputgrid = eris_ifu_resample_set_outputgrid(recipe_name, parlist,
2245 restable, wcs);
2246 if(cpl_error_get_code() != CPL_ERROR_NONE) {
2247
2248 eris_check_error_code("eris_ifu_combine");
2249 return CPL_ERROR_INCOMPATIBLE_INPUT;
2250 }
2251
2252 cpl_ensure_code(aParams_outputgrid, CPL_ERROR_NULL_INPUT);
2253
2254 hdrl_resample_result *result = NULL;
2255
2256 /* !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! */
2257 /* !!! Calling the hdrl_resample_compute function doing the real work !!! */
2258 /* !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! */
2259 cpl_msg_severity level = cpl_msg_get_level();
2260 cpl_msg_set_level(CPL_MSG_INFO);
2261 result = hdrl_resample_compute(restable, aParams_method, aParams_outputgrid,
2262 wcs);
2263 cpl_msg_set_level(level);
2264 /* TO BE REMOVED: this was added to check position of object (verify DAR)
2265 hdrl_image* hima = NULL;
2266 double max_x = 0, max_y = 0, peak = 0;
2267 cpl_size margin = 100;
2268 for(cpl_size k= margin; k < hdrl_imagelist_get_size(result->himlist) - margin; k++) {
2269 hima = hdrl_imagelist_get(result->himlist, k);
2270 eris_gaussian_maxpos(hdrl_image_get_image(hima), 5, &max_x, &max_y, &peak);
2271 cpl_msg_info(cpl_func,"z=%lld x: %g y: %g",k,max_x,max_y);
2272 }
2273 */
2274
2275 cpl_ensure_code(result, CPL_ERROR_NULL_INPUT);
2276 char* pro_catg = NULL;
2277 if(strstr(pipefile_prefix,"no_flat") != NULL) {
2278 pro_catg = cpl_sprintf("%s_%s",proCatg,"COADD_NOFLAT");
2279 } else {
2280 pro_catg = cpl_sprintf("%s_%s",proCatg,"COADD");
2281 }
2282 char* fname = cpl_sprintf("%s_%s_coadd.fits", pipefile_prefix, filenameSpec);
2283
2284 eris_ifu_resample_save_cube(result, pro_catg, recipe_name, fname, parlist,
2285 frameset, CPL_FALSE);
2286
2287 eris_ifu_free_string(&fname);
2288 eris_ifu_free_string(&pro_catg);
2289 hdrl_image* cube_collapsed = NULL;
2290 cpl_image* cube_cmap = NULL;
2291 hdrl_imagelist_collapse_mean(result->himlist, &cube_collapsed, &cube_cmap);
2292
2293 /* TODO: check header is correct */
2294 pro_catg = cpl_sprintf("%s_%s",proCatg,"MEAN");
2295 fname = cpl_sprintf("%s_%s_mean.fits", pipefile_prefix, filenameSpec);
2296 cpl_propertylist* phead2D = cpl_propertylist_duplicate(result->header);
2297 /* remove WCS from primary header as data go into extensions */
2298// eris_ifu_plist_erase_wcs(phead2D);
2299 cpl_propertylist_append_string(phead2D, "PRODCATG",PRODCATG_WHITELIGHT);
2300
2301 eris_ifu_save_image(frameset, phead2D, parlist, recipe_name, pro_catg,
2302 fname, CPL_TYPE_DOUBLE, hdrl_image_get_image(cube_collapsed));
2303 cpl_image_save(hdrl_image_get_error(cube_collapsed), fname, CPL_TYPE_DOUBLE,
2304 NULL, CPL_IO_EXTEND);
2305 cpl_image_save(cpl_image_new_from_mask(hdrl_image_get_mask(cube_collapsed)),
2306 fname, CPL_TYPE_INT, NULL, CPL_IO_EXTEND);
2307
2308 eris_ifu_free_string(&fname);
2309 eris_ifu_free_string(&pro_catg);
2310 eris_ifu_free_string(&filenameSpec);
2311 eris_ifu_free_string(&proCatg);
2313 hdrl_image_delete(cube_collapsed);
2314 cpl_image_delete(cube_cmap);
2315
2316
2317 /* Printing the wcs after resampling: */
2318
2319 cpl_wcs_delete(wcs);
2320 wcs = cpl_wcs_new_from_propertylist(result->header);
2321 cpl_msg_info(cpl_func, "Final WCS after resampling: ");
2322 cpl_msg_indent_more();
2323 eris_ifu_resample_wcs_print(wcs);
2324 cpl_msg_indent_less();
2325
2326 /* Cleanup */
2327
2328 /* Delete the parameters */
2329 hdrl_parameter_delete(aParams_method);
2330 hdrl_parameter_delete(aParams_outputgrid);
2331
2333 cpl_table_delete(restable);
2334 cpl_wcs_delete(wcs);
2335 cpl_frameset_delete(cube_set);
2336 //cpl_frameset_delete(var_err_set);
2337 //cpl_frameset_delete(bpm_set);
2338 eris_check_error_code("eris_ifu_combine");
2339
2340 return cpl_error_get_code();
2341}
2342
2356static cpl_propertylist* eris_ifu_resample_get_header2D(cpl_propertylist *header3D){
2357 cpl_propertylist *phead2D = cpl_propertylist_duplicate(header3D);
2358 cpl_propertylist_set_int(phead2D,"NAXIS",2);
2359 cpl_propertylist_erase_regexp(phead2D,"NAXIS3",0);
2360 cpl_propertylist_erase_regexp(phead2D,"CTYPE3",0);
2361 cpl_propertylist_erase_regexp(phead2D,"CRVAL3",0);
2362 cpl_propertylist_erase_regexp(phead2D,"CRPIX3",0);
2363 cpl_propertylist_erase_regexp(phead2D,"CUNIT3",0);
2364 cpl_propertylist_erase_regexp(phead2D,"CDELT3",0);
2365 cpl_propertylist_erase_regexp(phead2D,"CD3_*",0);
2366 cpl_propertylist_erase_regexp(phead2D,"CD1_3",0);
2367 cpl_propertylist_erase_regexp(phead2D,"CD2_3",0);
2368 eris_check_error_code("eris_ifu_resample_get_header2D");
2369 return phead2D;
2370}
2371
2387cpl_error_code eris_ifu_resample_trim_edge(hdrl_image *himg, int edge_trim){
2388 if (edge_trim > 0) {
2389 //cpl_msg_info(cpl_func, "Trim input image edges of %d pixels", edge_trim);
2390 /* trim each product plane edge of edgetrim pixels */
2391 cpl_size sx = hdrl_image_get_size_x(himg);
2392 cpl_size sy = hdrl_image_get_size_y(himg);
2393 if (edge_trim >= 0.5* sx) {
2394 edge_trim = 0;
2395 cpl_msg_warning(cpl_func, "edge-trim must be smaller than half image "
2396 "size. Reset to 0");
2397 }
2398 if (edge_trim >= 0.5* sy) {
2399 edge_trim = 0;
2400 cpl_msg_warning(cpl_func, "edge-trim must be smaller than half image "
2401 "size. Reset to 0");
2402 }
2403
2404 /* Note FITS convention to loop over pixels */
2405 /* trim lower image border along X direction */
2406 for(cpl_size j = 1; j <= edge_trim; j++) {
2407 for(cpl_size i = 1; i <= sx; i++) {
2408 hdrl_image_reject(himg, i, j);
2409 }
2410 }
2411 /* trim upper image border along X direction */
2412 for(cpl_size j = sy - edge_trim + 1; j <= sy; j++) {
2413 for(cpl_size i = 1; i <= sx; i++) {
2414 hdrl_image_reject(himg, i, j);
2415 }
2416 }
2417 /* trim left image border along Y direction */
2418 for(cpl_size j = 1; j <= sy; j++) {
2419 for(cpl_size i = 1; i <= edge_trim; i++) {
2420 hdrl_image_reject(himg, i, j);
2421 }
2422 }
2423 /* trim right image border along Y direction */
2424 for(cpl_size j = 1; j <= sy; j++) {
2425 for(cpl_size i = sx - edge_trim + 1; i <= sx; i++) {
2426 hdrl_image_reject(himg, i, j);
2427 }
2428 }
2429 }
2430 else {
2431 cpl_msg_warning(cpl_func, "edge-trim is set to 0");
2432 }
2433 eris_check_error_code("eris_ifu_resample_trim_edge");
2434 return cpl_error_get_code();
2435}
2436
2461cpl_error_code
2462eris_ifu_combine_pbp(cpl_frameset* frameset,
2463 const cpl_parameterlist * parlist,
2464 const char *input_cube_pro_catg,
2465 const char *filenameSpec,
2466 float *offsetx,
2467 float *offsety,
2468 const char* offunit,
2469 const char* recipe_name,
2470 const char* pipefile_prefix){
2471
2472 cpl_frameset* cube_set = NULL;
2473 cpl_msg_info(cpl_func,"Combine cubes into a single one plane-by-plane, compute mean, median along wavelength axis");
2474
2475 if ((cube_set = eris_ifu_extract_frameset(frameset, input_cube_pro_catg )) == NULL){
2476
2477 cpl_msg_warning(cpl_func,"No %s files", input_cube_pro_catg);
2478 //cpl_frameset_dump(frameset, stdout);
2479 cpl_error_set_message(cpl_func, CPL_ERROR_FILE_NOT_FOUND,
2480 "Missing cube files tagged as %s", input_cube_pro_catg);
2481 return cpl_error_get_code();
2482 } else {
2483 cpl_msg_info(cpl_func, "%02d file(s) of type %s",
2484 (int)cpl_frameset_get_size(cube_set), input_cube_pro_catg);
2485 }
2486
2487 /*convert the images in a table of predefined structure*/
2488 //const cpl_parameter * par = NULL;
2489 const cpl_frame * data_frame = NULL;
2490 hdrl_parameter *aParams_outputgrid = NULL;
2491 cpl_propertylist *phead2D = NULL;
2492 int extnum_raw = 1;
2493 int extnum_err = 2;
2494 int extnum_bpm = 3;
2495
2496 //cpl_boolean is_variance = CPL_TRUE;
2497 cpl_boolean subtract_bkg = CPL_FALSE;
2498 int edge_trim = 0;
2499 /*Loops over all frames and adds the images/imagelists to the final table */
2500 cpl_size ncubes = cpl_frameset_get_size(cube_set);
2501 if (ncubes == 1){
2502 cpl_error_set_message(cpl_func, CPL_ERROR_FILE_NOT_FOUND, "Only one CUBE "
2503 "files");
2504 return cpl_error_get_code();
2505 }
2506
2507 char* pname = cpl_sprintf("eris.%s.max-cubes-centres-dist",recipe_name);
2508 int max_cubes_dist = 20;
2509 //cpl_parameterlist_dump(parlist,stdout);
2510 if (offsetx !=NULL && offsety !=NULL) // user offset mode
2511 max_cubes_dist = INT_MAX ;//set to the largest, i.e. ignore this param.
2512 else
2513 max_cubes_dist = cpl_parameter_get_int(
2514 cpl_parameterlist_find_const(parlist, pname));
2515 cpl_free(pname);
2516
2517 pname = cpl_sprintf("eris.%s.subtract-background",recipe_name);
2518 subtract_bkg = cpl_parameter_get_bool(
2519 cpl_parameterlist_find_const(parlist, pname));
2520 if(subtract_bkg == CPL_TRUE)
2521 cpl_msg_info(cpl_func, "Median of each plane will be subtracted.");
2522 cpl_free(pname);
2523
2524 pname = cpl_sprintf("eris.%s.edge-trim",recipe_name);
2525 edge_trim = cpl_parameter_get_int(
2526 cpl_parameterlist_find_const(parlist, pname));
2527 if (edge_trim > 0)
2528 cpl_msg_info(cpl_func, "%d pixel(s) on edges will be trimmed.", edge_trim);
2529
2530 cpl_free(pname);
2531
2532
2533
2534 if ( max_cubes_dist < eris_ifu_compute_max_cubes_dist(cube_set)) {
2535 cpl_msg_info(cpl_func,"max-cubes-centres-dist: %d",max_cubes_dist);
2536 cpl_msg_warning(cpl_func,"max distance between cube centres greater than threshold. Skip cube combination.");
2537 return CPL_ERROR_INCOMPATIBLE_INPUT;
2538 }
2539
2540 hdrl_parameter *aParams_method = NULL;
2541 /* set re-sampling method */
2542 char *context = NULL;
2543 context = cpl_sprintf("eris.%s", recipe_name);
2544 aParams_method = eris_ifu_resample_set_method(parlist, context);
2545 cpl_free(context);
2546 cpl_ensure_code(aParams_method, CPL_ERROR_NULL_INPUT);
2547
2548 hdrl_imagelist *res_cube = hdrl_imagelist_new();
2549 cpl_wcs *wcs = NULL;
2550 cpl_wcs **tmp_wcs_array = (cpl_wcs **) cpl_calloc(ncubes, sizeof(wcs));
2551 char const **name_array = (char const **) cpl_calloc(ncubes, sizeof(const char*));
2552
2553 cpl_wcs *wcs3d = NULL;
2554 double mjd_obs = 0;
2555 double crpix3 = 0;
2556 double crval3 = 0;
2557 // double cdelt3 = 0;
2558 double cd1_3 = 0;
2559 double cd2_3 = 0;
2560 double cd3_1 = 0;
2561 double cd3_2 = 0;
2562 double cd3_3 = 0;
2563 double cd1_1 = 0;
2564 double cd1_2 = 0;
2565 double cd2_1 = 0;
2566 double cd2_2 = 0;
2567 int naxis3 = 0;
2568 const char* ctype3 = "WAVE";
2569 // const char* cunit3 = "um";
2570 const cpl_matrix *cd = NULL;
2571 cpl_vector* exptime_vec = cpl_vector_new(ncubes);
2572 double exptime = 0;
2573 double dit = 0;
2574 int ndit = 0;
2575 for (cpl_size i = 0; i < ncubes ; i++) {
2576 data_frame = cpl_frameset_get_position_const(cube_set, i);
2577 const char *name = cpl_frame_get_filename(data_frame);
2578 cpl_propertylist *xheader_data = cpl_propertylist_load(name,
2579 extnum_raw);
2580 cpl_propertylist* phead = cpl_propertylist_load(name, 0);
2581 dit = eris_pfits_get_dit(phead);
2582 ndit = eris_pfits_get_ndit(phead);
2583
2584 cpl_propertylist_delete(phead);
2585 exptime = dit * ndit;
2586 cpl_vector_set(exptime_vec, i, exptime);
2587 cpl_propertylist *header2D = eris_ifu_resample_get_header2D(xheader_data);
2588 cpl_wcs *tmp_wcs = cpl_wcs_new_from_propertylist(header2D);
2589 cpl_propertylist_delete(header2D);
2590
2591 tmp_wcs_array[i] = tmp_wcs;
2592 name_array[i] = name;
2593
2594 if (i==0){
2595 wcs = tmp_wcs;
2596 cpl_ensure_code(wcs, CPL_ERROR_NULL_INPUT);
2597 cpl_msg_info(cpl_func, "WCS used for the output grid definition and passed to "
2598 "hdrl_resample_compute():");
2599 cpl_msg_indent_more();
2600 eris_ifu_resample_wcs_print(wcs);
2601 cpl_msg_indent_less();
2602 wcs3d = cpl_wcs_new_from_propertylist(xheader_data);
2603 mjd_obs = cpl_propertylist_get_double(xheader_data, "MJD-OBS");
2604 crpix3 = cpl_propertylist_get_double(xheader_data, "CRPIX3");
2605 crval3 = cpl_propertylist_get_double(xheader_data, "CRVAL3");
2606 naxis3 = cpl_propertylist_get_int(xheader_data, "NAXIS3");
2607 //cdelt3 = cpl_propertylist_get_double(xheader_data, "CDELT3");
2608 /*
2609 if(cpl_propertylist_has(xheader_data,"CTYPE3")) {
2610 ctype3 = cpl_propertylist_get_string(xheader_data, "CTYPE3");
2611 cpl_msg_warning(cpl_func,"check ctype3: >%s< >%s<",
2612 ctype3,cpl_propertylist_get_string(xheader_data, "CTYPE3"));
2613 } else {
2614 ctype3 = "WAVE";
2615 }
2616 if(cpl_propertylist_has(xheader_data,"CUNIT3")) {
2617 cunit3 = cpl_propertylist_get_string(xheader_data, "CUNIT3");
2618 } else {
2619 cunit3 = "um";
2620 }
2621 */
2622 cd3_3 = cpl_propertylist_get_double(xheader_data, "CD3_3");
2623 // cdelt3 = cd3_3;
2624
2625 cd = cpl_wcs_get_cd(wcs);
2626 cd1_1 = cpl_matrix_get(cd, 0, 0);
2627 cd1_2 = cpl_matrix_get(cd, 1, 0);
2628 cd2_1 = cpl_matrix_get(cd, 0, 1);
2629 cd2_2 = cpl_matrix_get(cd, 1, 1);
2630
2631 cpl_ensure_code(wcs3d, CPL_ERROR_NULL_INPUT);
2632 }
2633 cpl_propertylist_delete(xheader_data);
2634 }
2635
2636 /*Convert offset unit to degree*/
2637 double as2deg = 1.0/3600.0;
2638 double scalex = 1.0;
2639 double scaley = 1.0;
2640
2641 int offcode = -1;
2642 if (offunit != NULL && offsetx !=NULL && offsety !=NULL){
2643 if (strcmp(offunit, "PIXEL") == 0){
2644 offcode = 0;
2645 cpl_msg_info(cpl_func, "User-defined offset unit in PIXEL");
2646 }else if (strcmp(offunit, "DEGREE") == 0){
2647 offcode = 1;
2648 cpl_msg_info(cpl_func, "User-defined offset unit is DEGREE");
2649 }else if (strcmp(offunit, "ARCSEC") == 0){
2650 offcode = 2;
2651 cpl_msg_info(cpl_func, "User-defined offset unit is ARCSEC");
2652 scalex = as2deg;
2653 scaley = as2deg;
2654 }else
2655 cpl_msg_error(cpl_func, "User-defined offset unit or offset list is not correct.");
2656 }
2657
2658
2659 const cpl_array *dims = cpl_wcs_get_image_dims(wcs3d);
2660 int nchan = cpl_array_get_int(dims, 2, NULL);
2661 cpl_size ich_center = 0.5 * nchan;
2662 cpl_msg_info(cpl_func, "Number of planes: %4.4d", nchan);
2663 //cpl_table* restable_exptime = NULL;
2664// hdrl_image *res_himg_exptime = NULL;
2665 hdrl_imagelist* res_himlist_exptime = hdrl_imagelist_new();
2666
2667 cpl_msg_severity level = cpl_msg_get_level();
2668 cpl_msg_set_level(CPL_MSG_INFO);
2669 for (cpl_size ich = 0; ich < nchan; ich++){
2670 cpl_table* restable = NULL; /*Final table to resample*/
2671
2672 for (cpl_size i = 0; i < ncubes ; i++) {
2673 const char *name = name_array[i];
2674 cpl_image *im = cpl_image_load(name, CPL_TYPE_DOUBLE, ich, extnum_raw);
2675 cpl_image *err = cpl_image_load(name, CPL_TYPE_DOUBLE, ich, extnum_err);
2676 cpl_image* qual = cpl_image_load(name, CPL_TYPE_INT, ich, extnum_bpm);
2677 cpl_mask *mask = cpl_mask_threshold_image_create(qual, 0, INT_MAX);
2678 cpl_image_delete(qual);
2679 // EXPOSURE MAP
2680 /* exposure map: use the plane at the center of the wavelength range
2681 * then replace intensity with EXPTIME value and error with 1 and
2682 * generate HDRL image. This will be re-sampled as the cube and give
2683 * the exposure map associated to the combined cube
2684 */
2685
2686 cpl_image* ima_time = NULL;
2687 cpl_image* err_time = NULL;
2688 hdrl_image *himg_exptime = NULL;
2689 if(ich == ich_center) {
2690 ima_time = cpl_image_new(cpl_image_get_size_x(im),
2691 cpl_image_get_size_y(im), CPL_TYPE_DOUBLE);
2692 err_time = cpl_image_new(cpl_image_get_size_x(im),
2693 cpl_image_get_size_y(im), CPL_TYPE_DOUBLE);
2694
2695 exptime = cpl_vector_get(exptime_vec, i);
2696 cpl_image_add_scalar(ima_time, exptime);
2697 cpl_image_add_scalar(err_time, sqrt(exptime));
2698 himg_exptime = hdrl_image_create(ima_time, err_time);
2699 /* for EXPPOSURE MAP we assume all pixels are good */
2700 //hdrl_image_reject_from_mask(himg_exptime, mask);
2701
2702 if (edge_trim > 0) {
2703 eris_ifu_resample_trim_edge(himg_exptime, edge_trim);
2704 }
2705
2706 cpl_image_delete(ima_time);
2707 cpl_image_delete(err_time);
2708 }
2709
2710 hdrl_image *himg = hdrl_image_create(im, err);
2711 cpl_image_delete(im);
2712 cpl_image_delete(err);
2713 hdrl_image_reject_from_mask(himg, mask);
2714 cpl_mask_delete(mask);
2715
2716 /* Subtract the background on request */
2717 if(subtract_bkg == CPL_TRUE) {
2718 //cpl_msg_info(cpl_func, "Subtracting the median as requested ...");
2720 }
2721
2722 if (edge_trim > 0) {
2723 eris_ifu_resample_trim_edge(himg, edge_trim);
2724 }
2725
2726 cpl_table *tab_tmp = NULL;
2727 cpl_table *tab_exptime_tmp = NULL;
2728 if (offsetx !=NULL && offsety !=NULL){
2729 if (offcode==0){
2730 tab_tmp = hdrl_resample_image_to_table(himg, tmp_wcs_array[0]);
2731 cpl_table_subtract_scalar(tab_tmp, "ra", offsetx[i]*cd1_1+ offsety[i]*cd1_2);
2732 cpl_table_subtract_scalar(tab_tmp, "dec", offsetx[i]*cd2_1+ offsety[i]*cd2_2);
2733 if (ich ==0)
2734 cpl_msg_info(__func__, "User-defined offest in RA, DEC (ARCSEC) %g, %g",
2735 (offsetx[i]*cd1_1+ offsety[i]*cd1_2)*3600,
2736 (offsetx[i]*cd2_1+ offsety[i]*cd2_2)*3600);
2737
2738 if(ich == ich_center) {
2739 tab_exptime_tmp = hdrl_resample_image_to_table(himg_exptime,
2740 tmp_wcs_array[0]);
2741 cpl_table_subtract_scalar(tab_exptime_tmp, "ra", offsetx[i]*cd1_1+ offsety[i]*cd1_2);
2742 cpl_table_subtract_scalar(tab_exptime_tmp, "dec", offsetx[i]*cd2_1+ offsety[i]*cd2_2);
2743 }
2744 } else if (offcode==1 || offcode==2 ){
2745 tab_tmp = hdrl_resample_image_to_table(himg, tmp_wcs_array[0]);
2746 cpl_table_subtract_scalar(tab_tmp, "ra", offsetx[i]*scalex);
2747 cpl_table_subtract_scalar(tab_tmp, "dec", offsety[i]*scaley);
2748 if (ich ==0)
2749 cpl_msg_info(__func__, "User-defined offest in RA, DEC (%s) %g, %g",
2750 offunit,
2751 offsetx[i],
2752 offsety[i]);
2753
2754 if(ich == ich_center) {
2755 tab_exptime_tmp = hdrl_resample_image_to_table(himg_exptime, tmp_wcs_array[0]);
2756 cpl_table_subtract_scalar(tab_exptime_tmp, "ra", offsetx[i]*scalex);
2757 cpl_table_subtract_scalar(tab_exptime_tmp, "dec", offsety[i]*scaley);
2758 }
2759 }
2760 } else {
2761 tab_tmp = hdrl_resample_image_to_table(himg, tmp_wcs_array[i]);
2762 if(ich == ich_center) {
2763 tab_exptime_tmp = hdrl_resample_image_to_table(himg_exptime, tmp_wcs_array[i]);
2764 }
2765 }
2766 hdrl_image_delete(himg);
2767 eris_ifu_resample_update_table(tab_tmp, &restable);
2768
2769 if(ich == ich_center) {
2770
2771 hdrl_image_delete(himg_exptime);
2772 //eris_ifu_resample_update_table(tab_exptime_tmp, &restable_exptime);
2773 hdrl_resample_result *himl_exptime_tmp = NULL;
2774 himl_exptime_tmp = hdrl_resample_compute(tab_exptime_tmp,
2775 aParams_method, aParams_outputgrid, wcs);
2776 //const char* fname = cpl_sprintf("ima_time_%lld.fits",i);
2777 hdrl_image *himg_tmp = hdrl_image_duplicate(
2778 hdrl_imagelist_get_const(himl_exptime_tmp->himlist, 0));
2779 hdrl_imagelist_set(res_himlist_exptime, himg_tmp, i);
2780 hdrl_resample_result_delete(himl_exptime_tmp);
2781
2782 /*
2783 cpl_image_save(hdrl_image_get_image(himg_tmp), fname,
2784 CPL_TYPE_DOUBLE, NULL, CPL_IO_DEFAULT);
2785 */
2786 cpl_table_delete(tab_exptime_tmp);
2787
2788 }
2789
2790 cpl_table_delete(tab_tmp);
2791 if (ich%100 == 0)
2792 cpl_msg_debug(__func__, "Created table for channel %d cube %d", (int)ich, (int)i);
2793 } /* end loop over input data cubes */
2794 cpl_ensure_code(restable, CPL_ERROR_NULL_INPUT);
2795
2796 if (ich%100 == 0)
2797 cpl_msg_debug(__func__, "Table of the chanel %d size %lld x %lld", (int)ich,
2798 cpl_table_get_nrow(restable), cpl_table_get_ncol(restable));
2799
2800 if(ich==0){
2801 aParams_outputgrid = eris_ifu_resample_set_outputgrid2D(recipe_name, parlist,
2802 restable, wcs);
2803 if ( cpl_error_get_code() != CPL_ERROR_NONE) {
2804 eris_check_error_code("eris_ifu_combine_pbp");
2805 return CPL_ERROR_NONE;
2806 }
2807 cpl_ensure_code(aParams_outputgrid, CPL_ERROR_INCOMPATIBLE_INPUT);
2808
2809 //int err_code = cpl_table_save(restable, NULL, NULL, "test_table.fits", CPL_IO_CREATE);
2810 }
2811 cpl_ensure_code(aParams_outputgrid, CPL_ERROR_NULL_INPUT);
2812
2813 hdrl_resample_result *result = NULL;
2814 result = hdrl_resample_compute(restable, aParams_method, aParams_outputgrid, wcs);
2815
2816 /* OLD WAY
2817 if(ich == ich_center) {
2818 hdrl_resample_result *result_exptime = NULL;
2819 result_exptime = hdrl_resample_compute(restable_exptime, aParams_method, aParams_outputgrid, wcs);
2820 res_himg_exptime = hdrl_image_duplicate(hdrl_imagelist_get_const(result_exptime->himlist, 0));
2821
2822 }
2823 */
2824
2825 cpl_ensure_code(result, CPL_ERROR_NULL_INPUT);
2826 cpl_table_delete(restable);
2827
2828 if(ich==0){
2829 /* Printing the wcs after resampling: */
2830 cpl_wcs *res_wcs = cpl_wcs_new_from_propertylist(result->header);
2831 cpl_msg_info(cpl_func, "Final WCS after resampling: ");
2832 cpl_msg_indent_more();
2833 eris_ifu_resample_wcs_print(res_wcs);
2834 cpl_msg_indent_less();
2835 cpl_wcs_delete(res_wcs);
2836 phead2D = cpl_propertylist_duplicate(result->header);
2837 }
2838
2839 hdrl_image *res_himg = hdrl_image_duplicate(hdrl_imagelist_get_const(result->himlist, 0));
2840 hdrl_imagelist_set(res_cube, res_himg, ich);
2842
2843 if (ich==nchan-1){
2844 for (int i = 0; i < ncubes; i++){
2845 cpl_wcs_delete(tmp_wcs_array[i]);
2846 }
2847 cpl_free(tmp_wcs_array);
2848 cpl_free(name_array);
2849 }
2850
2851 if (ich%100 == 0)
2852 cpl_msg_info(__func__, "Coadding cubes: range [%4.4d,%4.4d] of %d",
2853 (int)ich,(int)ich+100, nchan);
2854 }
2855 cpl_msg_set_level(level);
2856
2857 hdrl_image* hima_exptime_mean = NULL;
2858 cpl_image* contrib_exptime_map = NULL;
2859 hdrl_imagelist_collapse_mean(res_himlist_exptime, &hima_exptime_mean,
2860 &contrib_exptime_map);
2861 hdrl_imagelist_delete(res_himlist_exptime);
2862 hdrl_image_delete(hima_exptime_mean);
2863 /*
2864 cpl_image_save(hdrl_image_get_image(hima_exptime_mean), "exptime_mean.fits",
2865 CPL_TYPE_DOUBLE, NULL, CPL_IO_DEFAULT);
2866 */
2867 cpl_image_multiply_scalar(contrib_exptime_map, cpl_vector_get_median(exptime_vec));
2868 //cpl_image_save(contrib, "contrib.fits", CPL_TYPE_DOUBLE, NULL, CPL_IO_DEFAULT);
2869 cpl_ensure_code(res_cube, CPL_ERROR_NULL_INPUT);
2870 cpl_vector_delete(exptime_vec);
2871
2872 hdrl_parameter_delete(aParams_method);
2873 hdrl_parameter_delete(aParams_outputgrid);
2874
2875 /* Save resample cube*/
2876
2877
2878 hdrl_resample_result *result = cpl_calloc(1, sizeof(hdrl_resample_result));
2879 result->himlist = res_cube;
2880 result->header = phead2D;
2881
2882
2883 cpl_propertylist_append_double(result->header, "CRPIX3", crpix3);
2884 cpl_propertylist_append_double(result->header, "CRVAL3", crval3);
2885 cpl_propertylist_append_int(result->header, "NAXIS3", naxis3);
2886 //cpl_propertylist_append_double(result->header, "CDELT3", cdelt3);
2887 cpl_propertylist_append_string(result->header, "CTYPE3", ctype3);
2888 cpl_propertylist_append_double(result->header, "CD1_3", cd1_3);
2889 cpl_propertylist_append_double(result->header, "CD2_3", cd2_3);
2890 cpl_propertylist_append_double(result->header, "CD3_3", cd3_3);
2891 cpl_propertylist_append_double(result->header, "CD3_1", cd3_1);
2892 cpl_propertylist_append_double(result->header, "CD3_2", cd3_2);
2893 cpl_propertylist_append_double(result->header, "MJD-OBS", mjd_obs);
2894 cpl_propertylist* wcs2D = eris_ifu_plist_extract_wcs2D(result->header);
2895
2896 hdrl_image* cube_collapsed = NULL;
2897 cpl_image* cube_cmap = NULL;
2898 hdrl_imagelist_collapse_mean(res_cube, &cube_collapsed, &cube_cmap);
2899 /* This makes 2D products not FITS standard
2900 cpl_propertylist* wcs2D = eris_ifu_plist_extract_wcs(result->header);
2901 cpl_propertylist_append(phead2D, wcs2D);
2902 cpl_propertylist_delete(wcs2D);
2903 */
2904
2905 char* pro_catg = NULL;
2906 if(strstr(pipefile_prefix,"no_flat") != NULL) {
2907 pro_catg = cpl_sprintf("%s_%s",input_cube_pro_catg,"COADD_NOFLAT");
2908 } else {
2909 pro_catg = cpl_sprintf("%s_%s",input_cube_pro_catg,"COADD");
2910 }
2911 char* fname = cpl_sprintf("%s_%s_coadd.fits", pipefile_prefix, filenameSpec);
2912
2913 if(strcmp(recipe_name,"eris_ifu_jitter") == 0) {
2914 eris_ifu_resample_save_cube(result, pro_catg, recipe_name, fname, parlist,
2915 frameset, CPL_FALSE);
2916 } else {
2917 eris_ifu_resample_save_cube(result, pro_catg, recipe_name, fname, parlist,
2918 frameset, CPL_FALSE);
2919 }
2920 eris_ifu_free_string(&fname);
2921 eris_ifu_free_string(&pro_catg);
2922
2923
2924 /* TODO: check header is correct */
2925 pro_catg = cpl_sprintf("%s_%s",input_cube_pro_catg,"MEAN");
2926 fname = cpl_sprintf("%s_%s_mean.fits", pipefile_prefix, filenameSpec);
2927 cpl_propertylist_append_string(phead2D, "PRODCATG", PRODCATG_WHITELIGHT);
2928 /* remove WCS from primary header as data go into extensions */
2929// eris_pfits_erase_wcs2D(phead2D);
2930 cpl_propertylist_update_string(phead2D, CPL_DFS_PRO_CATG, pro_catg);
2931 cpl_image* mask_ima = cpl_image_new_from_mask(hdrl_image_get_mask(cube_collapsed));
2932 eris_ifu_save_deq_image(frameset, NULL, parlist,frameset, NULL,
2933 recipe_name, phead2D, "RADECSYS", fname,
2934 hdrl_image_get_image(cube_collapsed),
2935 hdrl_image_get_error(cube_collapsed),
2936 rmse, mask_ima, flag16bit, UNIT_ADU);
2937 cpl_free(fname);
2938 fname = cpl_sprintf("%s_%s_exposure_map.fits", pipefile_prefix, filenameSpec);
2939
2940 cpl_propertylist_update_string(phead2D, CPL_DFS_PRO_CATG, ERIS_IFU_PRO_JITTER_EXPOSURE_MAP);
2941 cpl_propertylist_update_string(phead2D, "PRODCATG", PRODCATG_EXPOSUREMAP);
2942
2943 cpl_propertylist_append(phead2D, wcs2D);
2944 cpl_propertylist_delete(wcs2D);
2945 eris_ifu_save_image_phase3(frameset, /*NULL, */parlist,frameset, /*NULL,*/
2946 recipe_name, phead2D, "RADECSYS", fname, contrib_exptime_map,
2947 CPL_TYPE_DOUBLE, UNIT_TIME);
2948
2949
2950 cpl_image_delete(mask_ima);
2951 cpl_image_delete(contrib_exptime_map);
2952 eris_ifu_free_string(&fname);
2953 eris_ifu_free_string(&pro_catg);
2954 hdrl_image_delete(cube_collapsed);
2955 cpl_image_delete(cube_cmap);
2956 cpl_wcs_delete(wcs3d);
2958
2959 cpl_frameset_delete(cube_set);
2960 eris_check_error_code("eris_ifu_combine_pbp");
2961
2962 return cpl_error_get_code();
2963}
2964
ifsPreopticsScale eris_ifu_get_preopticsScale(cpl_propertylist *header)
Return the the pre-optics scaling.
cpl_error_code eris_ifu_heades_add_hduclass_qual(cpl_propertylist *plist, deqQualityType qualityType)
Add quality extension classification properties to FITS header.
Definition: eris_ifu_dfs.c:324
cpl_frameset * eris_ifu_extract_frameset(const cpl_frameset *in, const char *tag)
Extract frames with a specific tag from a frameset.
Definition: eris_ifu_dfs.c:201
cpl_error_code eris_ifu_save_image_phase3(cpl_frameset *allframes, const cpl_parameterlist *parlist, const cpl_frameset *usedframes, const char *recipe, const cpl_propertylist *applist, const char *remregexp, const char *filename, const cpl_image *image, cpl_type type, const char *unit)
Save a DFS-compliant Phase 3 image product.
Definition: eris_ifu_dfs.c:666
cpl_error_code eris_ifu_heades_add_hduclass_common(cpl_propertylist *plist, const char *deq_hduclas)
Add common HDU classification properties to FITS header.
Definition: eris_ifu_dfs.c:248
cpl_error_code eris_ifu_save_deq_image(cpl_frameset *allframes, cpl_propertylist *header, const cpl_parameterlist *parlist, const cpl_frameset *usedframes, const cpl_frame *inherit, const char *recipe, const cpl_propertylist *applist, const char *remregexp, const char *filename, const cpl_image *image, const cpl_image *error, deqErrorType errorType, const cpl_image *dataQualityMask, deqQualityType dataQualityType, const char *unit)
Save a DFS-compliant image product with data, error, and quality extensions.
Definition: eris_ifu_dfs.c:367
cpl_error_code eris_ifu_heades_add_hduclass_data(cpl_propertylist *plist)
Add data extension classification properties to FITS header.
Definition: eris_ifu_dfs.c:273
cpl_error_code eris_ifu_heades_add_hduclass_errs(cpl_propertylist *plist, deqErrorType errorType)
Add error extension classification properties to FITS header.
Definition: eris_ifu_dfs.c:297
eris_ifu_sdp_properties * eris_ifu_sdp_properties_collect(hdrl_resample_result *aCube, cpl_frameset *set, const cpl_parameterlist *parlist, const char *recipe_name)
Collect all SDP metadata from cube, frameset, and parameters.
void eris_ifu_sdp_properties_delete(eris_ifu_sdp_properties *aProperties)
Free SDP properties structure and all contained data.
Definition: eris_ifu_sdp.c:813
cpl_error_code eris_ifu_sdp_properties_update(cpl_propertylist *aHeader, const eris_ifu_sdp_properties *aProperties)
Update FITS header with ESO Science Data Product keywords.
cpl_error_code eris_ifu_jitter_get_procatg_and_filename(cubeType type, char **proCatg, char **filenamePrefix)
Get the value of the PRO.CATG and the filename as function of the type of cube.
cpl_error_code eris_ifu_combine(cubeType obj_type, cpl_frameset *frameset, const cpl_parameterlist *parlist, const char *recipe_name, const char *pipefile_prefix)
Resample and combine multiple IFU cubes into a single cube.
cpl_error_code eris_ifu_combine_pbp(cpl_frameset *frameset, const cpl_parameterlist *parlist, const char *input_cube_pro_catg, const char *filenameSpec, float *offsetx, float *offsety, const char *offunit, const char *recipe_name, const char *pipefile_prefix)
Resample and combine cubes plane-by-plane (2D spatial resampling per wavelength)
cpl_error_code eris_ifu_cube_collapse_mean_and_save(const char *pro_catg, cpl_frameset *frameset, const cpl_parameterlist *parlist, const char *recipe_name, cpl_boolean apply_flat, cpl_boolean is_pupil)
-------------------------------------------------------------------------—‍/
cpl_error_code eris_ifu_resample_save_cube(hdrl_resample_result *aCube, const char *procatg, const char *recipe, const char *filename, const cpl_parameterlist *parlist, cpl_frameset *frameset, cpl_boolean gen_phase3)
Save resampled cube to FITS file in DEQ (Data-Error-Quality) format.
cpl_error_code eris_ifu_resample_trim_edge(hdrl_image *himg, int edge_trim)
Trim edge pixels from hdrl_image by rejecting them.
void eris_ifu_free_propertylist(cpl_propertylist **item)
Free memory and set pointer to null.
void eris_ifu_free_string(char **item)
Free memory and set pointer to null.
cpl_table * eris_qclog_init(void)
Initialize QC table.
cpl_error_code eris_pfits_put_qc(cpl_propertylist *plist, cpl_table *qclog)
convert table with QC parameter information to a propertylist
cpl_error_code eris_qclog_add_double_f(cpl_table *table, const char *key_name, const double value, const char *key_help)
add QC float info to table
cpl_error_code eris_ifu_split3_hdrl_imagelist(hdrl_imagelist *cube, cpl_imagelist *dataCube, cpl_imagelist *errorCube, cpl_imagelist *qualCube)
extract from hdrl_imagelist the data and error information
cpl_error_code eris_ifu_save_image(cpl_frameset *fs, const cpl_propertylist *plist, const cpl_parameterlist *parlist, const char *recipe, const char *procatg, const char *filename, cpl_type type, const cpl_image *image)
Save image with DFS compliance.
int eris_pfits_get_ndit(const cpl_propertylist *plist)
find out the DIT value
Definition: eris_pfits.c:117
double eris_pfits_get_dit(const cpl_propertylist *plist)
find out the DIT value
Definition: eris_pfits.c:98
cpl_error_code eris_parameters_get_double(const cpl_parameterlist *parlist, const char *pname, double *pvalue)
get double parameter value if changed by the user
Definition: eris_utils.c:765
cpl_error_code eris_check_error_code(const char *func_id)
handle CPL errors
Definition: eris_utils.c:56
cpl_boolean eris_param_has_changed(const cpl_parameter *p)
verify if a parameter value has been changed (from command line or or rc file by a user)
Definition: eris_utils.c:809
cpl_error_code hdrl_image_reject_from_mask(hdrl_image *self, const cpl_mask *map)
set bpm of hdrl_image
Definition: hdrl_image.c:407
hdrl_value hdrl_image_get_median(const hdrl_image *self)
computes the median and associated error of an image.
hdrl_image * hdrl_image_duplicate(const hdrl_image *himg)
copy hdrl_image
Definition: hdrl_image.c:391
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
cpl_size hdrl_image_get_size_y(const hdrl_image *self)
return size of Y dimension of image
Definition: hdrl_image.c:540
cpl_size hdrl_image_get_size_x(const hdrl_image *self)
return size of X dimension of image
Definition: hdrl_image.c:525
hdrl_image * hdrl_image_create(const cpl_image *image, const cpl_image *error)
create a new hdrl_image from to existing images by copying them
Definition: hdrl_image.c:295
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
void hdrl_image_delete(hdrl_image *himg)
delete hdrl_image
Definition: hdrl_image.c:379
cpl_error_code hdrl_image_sub_scalar(hdrl_image *self, hdrl_value value)
Elementwise subtraction of a scalar from an image.
cpl_error_code hdrl_imagelist_collapse_mean(const hdrl_imagelist *himlist, hdrl_image **out, cpl_image **contrib)
Mean collapsing of image list.
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.
hdrl_imagelist * hdrl_imagelist_create(cpl_imagelist *imlist, cpl_imagelist *errlist)
Create an hdrl_imagelist out of 2 cpl_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.
void hdrl_parameter_delete(hdrl_parameter *obj)
shallow delete of a parameter
hdrl_parameter * hdrl_resample_parameter_create_nearest(void)
Creates a resample nearest neighbor hdrl parameter object. The algorithm does not use any weighting f...
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_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,...
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,...
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().