ERIS Pipeline Reference Manual 1.8.10
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/*----------------------------------------------------------------------------*/
46/*----------------------------------------------------------------------------*/
47#define MAX_COMBINE_CUBE_PLANE_DIAGONAL 1000
57static cpl_error_code
58eris_ifu_resample_update_table(const cpl_table* append, cpl_table** storage) {
59
60 cpl_ensure_code(append, CPL_ERROR_NULL_INPUT);
61
62 if(*storage == NULL) {
63 *storage = cpl_table_new(0);
64 cpl_table_copy_structure(*storage, append);
65 cpl_table_insert(*storage, append, 0);
66 } else {
67 cpl_size nrow = cpl_table_get_nrow(*storage);
68 cpl_table_insert(*storage, append, nrow);
69 }
70 eris_check_error_code("eris_ifu_resample_update_table");
71 return cpl_error_get_code();
72
73}
74
75
86static cpl_error_code
87eris_ifu_resample_tablesave (const char* recipe_name, const cpl_parameter *par,
88 const cpl_parameterlist *parlist, cpl_frameset *frameset,
89 cpl_table *tab_final, const char* pipefile_prefix)
90{
91 /* Save the final table if required */
92 char* param_name = cpl_sprintf("eris.%s.save-table", recipe_name);
93 par = cpl_parameterlist_find_const (parlist, param_name);
94 cpl_free(param_name);
95
96 const cpl_boolean save_table = cpl_parameter_get_bool (par);
97
98 if (save_table)
99 {
100 cpl_msg_info (cpl_func, "Saving the table before resampling");
101 cpl_propertylist *applist = cpl_propertylist_new ();
102 cpl_propertylist_update_string (applist, CPL_DFS_PRO_CATG,
103 "RESAMPLE_TABLE");
104 char* fname = cpl_sprintf("%s_resample_table.fits", pipefile_prefix);
105 cpl_dfs_save_table (frameset, NULL, parlist, frameset, NULL, tab_final,
106 NULL, recipe_name, applist, NULL, PACKAGE "/" PACKAGE_VERSION,
107 fname);
108 cpl_free(fname);
109 cpl_propertylist_delete (applist);
110 }
111 eris_check_error_code("eris_ifu_resample_tablesave");
112 return cpl_error_get_code();
113}
114
115
116
117/*----------------------------------------------------------------------------*/
125/*----------------------------------------------------------------------------*/
126
127/* set propertylist (int) value */
128static cpl_error_code
129eris_ifu_resample_parameters_get_int(const cpl_parameterlist* parlist,
130 const char* pname, int *pvalue)
131{
132 const cpl_parameter* p = cpl_parameterlist_find_const(parlist, pname);
133 cpl_boolean p_has_changed = eris_param_has_changed(p);
134
135 if ( cpl_parameter_get_default_flag(p) && p_has_changed != CPL_FALSE) {
136 *pvalue = cpl_parameter_get_int(p);
137 } else {
138 *pvalue = cpl_parameter_get_default_int(p);
139 }
140 eris_check_error_code("eris_ifu_resample_parameters_get_int");
141 return cpl_error_get_code();
142}
143
144/*----------------------------------------------------------------------------*/
154/*----------------------------------------------------------------------------*/
155
156
157static cpl_error_code
158set_double_from_parameters(const cpl_parameterlist* parlist,
159 const char* pname,
160 const double table_value,
161 double *value) {
162
163 double temp_value;
164 eris_parameters_get_double (parlist, pname, &temp_value);
165 if (temp_value != -1)
166 {
167 *value = temp_value;
168 }
169 else
170 {
171 *value = table_value;
172 }
173 eris_check_error_code("set_double_from_parameters");
174 return cpl_error_get_code();
175}
176
177
178/*----------------------------------------------------------------------------*/
188/*----------------------------------------------------------------------------*/
189
190
191static cpl_error_code
192set_int_from_parameters(const cpl_parameterlist* parlist,
193 const char* pname,
194 const int table_value,
195 int *value) {
196
197 int temp_value;
198 eris_ifu_resample_parameters_get_int (parlist, pname, &temp_value);
199 if (temp_value != -1)
200 {
201 *value = temp_value;
202 }
203 else
204 {
205 *value = table_value;
206 }
207 eris_check_error_code("set_int_from_parameters");
208 return cpl_error_get_code();
209}
210
211
212
213/*----------------------------------------------------------------------------*/
219/*----------------------------------------------------------------------------*/
220
221static hdrl_parameter *
222eris_ifu_resample_set_method (const cpl_parameterlist *parlist,
223 const char* context)
224{
225 hdrl_parameter *aParams_method = NULL;
226 const cpl_parameter *p = NULL;
227 cpl_boolean p_has_changed = 0;
228 int loop_distance = 0;
229 double critical_radius_renka = 0.;
230 int kernel_size_lanczos = 0;
231 double pix_frac_drizzle_x = 0.;
232 double pix_frac_drizzle_y= 0.;
233 double pix_frac_drizzle_l= 0.;
234 cpl_boolean use_errorweights = FALSE;
235
236 const cpl_parameter *par = NULL;
237 char* param_name = NULL;
238 /* Apply additional weights based on the error if required */
239
240 param_name = cpl_sprintf("%s.method.use-errorweights",context);
241 par = cpl_parameterlist_find_const (parlist, param_name);
242 use_errorweights = cpl_parameter_get_bool (par);
243 cpl_free(param_name);
244
246 param_name = cpl_sprintf("%s.method.loop-distance",context);
247 p = cpl_parameterlist_find_const (parlist, param_name);
248
249 p_has_changed = eris_param_has_changed (p);
250 if (cpl_parameter_get_default_flag (p) && p_has_changed != CPL_FALSE)
251 {
252 loop_distance = cpl_parameter_get_int (p);
253 }
254 else
255 {
256 loop_distance = LOOP_DISTANCE;
257 }
258 cpl_free(param_name);
259
261 param_name = cpl_sprintf("%s.method.renka.critical-radius",context);
262 set_double_from_parameters(parlist, param_name,
263 RENKA_CRITICAL_RADIUS, &critical_radius_renka);
264 cpl_free(param_name);
265
267 param_name = cpl_sprintf("%s.method.lanczos.kernel-size",context);
268 set_int_from_parameters(parlist, param_name,
269 LANCZOS_KERNEL_SIZE, &kernel_size_lanczos);
270 cpl_free(param_name);
271
277 param_name = cpl_sprintf("%s.method.drizzle.downscale-x",context);
278 set_double_from_parameters(parlist, param_name, DRIZZLE_DOWN_SCALING_FACTOR_X,
279 &pix_frac_drizzle_x);
280 cpl_free(param_name);
281
282 param_name = cpl_sprintf("%s.method.drizzle.downscale-y",context);
283 set_double_from_parameters(parlist, param_name, DRIZZLE_DOWN_SCALING_FACTOR_Y,
284 &pix_frac_drizzle_y);
285 cpl_free(param_name);
286
287 param_name = cpl_sprintf("%s.method.drizzle.downscale-z",context);
288 set_double_from_parameters(parlist, param_name, DRIZZLE_DOWN_SCALING_FACTOR_Z,
289 &pix_frac_drizzle_l);
290 cpl_free(param_name);
291
292 /* Create the right re-sampling parameter */
293 param_name = cpl_sprintf("%s.method",context);
294 par = cpl_parameterlist_find_const (parlist, param_name);
295 const char *method = cpl_parameter_get_string (par);
296
297 if (strcmp (method, "NEAREST") == 0)
298 {
300 }
301 else if (strcmp (method, "RENKA") == 0)
302 {
303 aParams_method = hdrl_resample_parameter_create_renka(loop_distance,
304 use_errorweights,
305 critical_radius_renka);
306 }
307 else if (strcmp (method, "LINEAR") == 0)
308 {
309 aParams_method = hdrl_resample_parameter_create_linear(loop_distance,
310 use_errorweights);
311 }
312 else if (strcmp (method, "QUADRATIC") == 0)
313 {
314 aParams_method = hdrl_resample_parameter_create_quadratic(loop_distance,
315 use_errorweights);
316 }
317 else if (strcmp (method, "DRIZZLE") == 0)
318 {
319 aParams_method = hdrl_resample_parameter_create_drizzle(loop_distance,
320 use_errorweights,
321 pix_frac_drizzle_x,
322 pix_frac_drizzle_y,
323 pix_frac_drizzle_l);
324 }
325 else if (strcmp (method, "LANCZOS") == 0)
326 {
327 aParams_method = hdrl_resample_parameter_create_lanczos(loop_distance,
328 use_errorweights,
329 kernel_size_lanczos);
330 }
331 else
332 {
333 aParams_method = hdrl_resample_parameter_create_lanczos(loop_distance,
334 use_errorweights,
335 kernel_size_lanczos);
336 cpl_msg_warning (cpl_func,
337 "%s is an unsupported method! Default to LANCZOS",
338 method);
339 }
340 cpl_free(param_name);
341 //eris_ifu_resample_print_method_params(aParams_method);
342 eris_check_error_code("eris_ifu_resample_set_method");
343 return aParams_method;
344}
345
346/*----------------------------------------------------------------------------*/
356/*----------------------------------------------------------------------------*/
357
358static cpl_error_code
359eris_ifu_resample_get_cdelt123(cpl_wcs *wcs, double *cdelt1, double *cdelt2,
360 double *cdelt3)
361{
362 cpl_ensure_code(wcs, CPL_ERROR_NULL_INPUT);
363 cpl_ensure_code(cdelt1, CPL_ERROR_NULL_INPUT);
364 cpl_ensure_code(cdelt2, CPL_ERROR_NULL_INPUT);
365
366 double cd1_1 = 0;
367 double cd1_2 = 0;
368 double cd2_1 = 0;
369 double cd2_2 = 0;
370 double cd3_3 = 0;
371
372 double dx=0;
373 double dy=0;
374
375 const cpl_matrix *cd = cpl_wcs_get_cd(wcs);
376 cd1_1 = cpl_matrix_get(cd, 0, 0);
377 cd1_2 = cpl_matrix_get(cd, 1, 0);
378 cd2_1 = cpl_matrix_get(cd, 0, 1);
379 cd2_2 = cpl_matrix_get(cd, 1, 1);
380
381 dx = sqrt(pow(cd1_1,2) + pow(cd2_1,2));
382 dy = sqrt(pow(cd1_2,2)+ pow(cd2_2,2));
383
384// double det = cd1_1 * cd2_2 - cd1_2 * cd2_1;
385// if (det < 0)
386// dx *= -1.0;
387
388// *cdelt1 = fabs(cpl_matrix_get(cd, 0, 0)); /* CD1_1 */
389// *cdelt2 = fabs(cpl_matrix_get(cd, 1, 1)); /* CD2_2 */
390 *cdelt1 = dx;
391 *cdelt2 = dy;
392
393 /* Changes for Cube and image */
394 cpl_size cdncol = cpl_matrix_get_ncol(cd);
395 if (cdncol == 3) {
396 cd3_3= fabs(cpl_matrix_get(cd, 2, 2)); /* CD3_3 */
397 }
398 else {
399 cd3_3 = 1.; /* CD3_3 */
400 }
401
402 *cdelt3 = cd3_3;
403
404 /* Check keys (for debug) */
405 cpl_msg_debug(cpl_func, "cdelt1: %g, cdelt2: %g, cdelt3: %g",
406 *cdelt1, *cdelt2, *cdelt3);
407 eris_check_error_code("eris_ifu_resample_get_cdelt123");
408 return cpl_error_get_code();
409}
410
411
412/*----------------------------------------------------------------------------*/
420/*----------------------------------------------------------------------------*/
421
422static hdrl_parameter *
423eris_ifu_resample_set_outputgrid (const char* recipe_name,
424 const cpl_parameterlist *parlist, cpl_table *muse_table, cpl_wcs *wcs)
425{
426 /* Should be done in a dedicated _new function */
427 cpl_error_ensure(parlist != NULL, CPL_ERROR_NULL_INPUT,
428 return NULL, "NULL Input Parameters");
429 cpl_error_ensure(muse_table != NULL, CPL_ERROR_NULL_INPUT,
430 return NULL, "NULL Input table");
431 cpl_error_ensure(wcs != NULL, CPL_ERROR_NULL_INPUT,
432 return NULL, "NULL Input wcs");
433
434 hdrl_parameter *aParams_outputgrid = NULL;
435 char* param_name = NULL;
436 char* context = NULL;
437 double ra_min = 0.;
438 double ra_max = 0.;
439 double dec_min = 0.;
440 double dec_max = 0.;
441 double lambda_min = 0.;
442 double lambda_max = 0.;
443
444 double dx = 0.;
445 double dy = 0.;
446 double dlambda = 0.;
447
448 context = cpl_sprintf("eris.%s", recipe_name);
449
450 /* init relevant parameters */
451 double ra_min_tmp =
452 cpl_table_get_column_min (muse_table, HDRL_RESAMPLE_TABLE_RA);
453 double ra_max_tmp =
454 cpl_table_get_column_max (muse_table, HDRL_RESAMPLE_TABLE_RA);
455 double dec_min_tmp =
456 cpl_table_get_column_min (muse_table, HDRL_RESAMPLE_TABLE_DEC);
457 double dec_max_tmp =
458 cpl_table_get_column_max (muse_table, HDRL_RESAMPLE_TABLE_DEC);
459 double lambda_min_tmp =
460 cpl_table_get_column_min (muse_table, HDRL_RESAMPLE_TABLE_LAMBDA);
461 double lambda_max_tmp =
462 cpl_table_get_column_max (muse_table, HDRL_RESAMPLE_TABLE_LAMBDA);
463
464 /* We have the rare case that the image spans over ra = 0.*/
465 if(ra_max_tmp - ra_min_tmp > 180){
466 const double *ra = cpl_table_get_data_double_const(muse_table,
467 HDRL_RESAMPLE_TABLE_RA);
468 /* Should we also take the bpm into account ?
469 const int *bpm = cpl_table_get_data_int_const(ResTable,
470 HDRL_RESAMPLE_TABLE_BPM);
471 */
472
473 /* set to extreme values for a start */
474 ra_min_tmp = 0.;
475 ra_max_tmp = 360.;
476 cpl_size nrow = cpl_table_get_nrow(muse_table);
477
478 for (cpl_size i = 0; i < nrow; i++) {
479 if (ra[i] > ra_min_tmp && ra[i] <= 180.) ra_min_tmp = ra[i]; /* get the maximum */
480 if (ra[i] < ra_max_tmp && ra[i] > 180.) ra_max_tmp = ra[i]; /* get the minimum */
481 }
482 }
483
484 /* Check output (for debug) */
485 cpl_msg_debug (cpl_func, "min x %10.7f", ra_min_tmp);
486 cpl_msg_debug (cpl_func, "max x %10.7f", ra_max_tmp);
487 cpl_msg_debug (cpl_func, "min y %10.7f", dec_min_tmp);
488 cpl_msg_debug (cpl_func, "max y %10.7f", dec_max_tmp);
489 cpl_msg_debug (cpl_func, "min lambda %10.7f", lambda_min_tmp);
490 cpl_msg_debug (cpl_func, "max lambda %10.7f", lambda_max_tmp);
491 /*
492 param_name = cpl_sprintf("%s.outgrid.ra-min", context);
493 set_double_from_parameters(parlist, param_name, ra_min_tmp, &ra_min);
494 cpl_free(param_name);
495 */
496 ra_min = ra_min_tmp;
497 /*
498 param_name = cpl_sprintf("%s.outgrid.ra-max", context);
499 set_double_from_parameters(parlist, param_name, ra_max_tmp, &ra_max);
500 cpl_free(param_name);
501 */
502 ra_max = ra_max_tmp;
503
504 /*
505 param_name = cpl_sprintf("%s.outgrid.dec-min", context);
506 set_double_from_parameters(parlist, param_name, dec_min_tmp, &dec_min);
507 cpl_free(param_name);
508 */
509 dec_min = dec_min_tmp;
510
511 /*
512 param_name = cpl_sprintf("%s.outgrid.dec-max", context);
513 set_double_from_parameters(parlist, param_name, dec_max_tmp, &dec_max);
514 cpl_free(param_name);
515 */
516 dec_max = dec_max_tmp;
517
518 /*
519 param_name = cpl_sprintf("%s.outgrid.lambda-min", context);
520 set_double_from_parameters(parlist, param_name, lambda_min_tmp, &lambda_min);
521 cpl_free(param_name);
522
523 param_name = cpl_sprintf("%s.outgrid.lambda-max", context);
524 set_double_from_parameters(parlist, param_name, lambda_max_tmp, &lambda_max);
525 cpl_free(param_name);
526 */
527 lambda_max = lambda_max_tmp;
528
529
530 /* Reset all delta values */
531 double cdelt1 = 0., cdelt2 = 0., cdelt3 = 0.;
532 eris_ifu_resample_get_cdelt123 (wcs, &cdelt1, &cdelt2, &cdelt3);
533 /*
534 param_name = cpl_sprintf("%s.outgrid.delta-ra", context);
535 set_double_from_parameters(parlist, param_name, cdelt1, &dx);
536 cpl_free(param_name);
537 */
538 dx = cdelt1;
539
540 /*
541 param_name = cpl_sprintf("%s.outgrid.delta-dec", context);
542 set_double_from_parameters(parlist, param_name, cdelt2, &dy);
543 cpl_free(param_name);
544 */
545 dy = cdelt2;
546
547 /*
548 param_name = cpl_sprintf("%s.outgrid.delta-lambda", context);
549 set_double_from_parameters(parlist, param_name, cdelt3, &dlambda);
550 cpl_free(param_name);
551 */
552 dlambda = cdelt3;
553
554
555 /* Assign the field margin */
556 const cpl_parameter *par = NULL;
557 param_name = cpl_sprintf("%s.fieldmargin", context);
558 par = cpl_parameterlist_find_const (parlist, param_name);
559 cpl_free(param_name);
560
561 double fieldmargin = 0.;
562 fieldmargin = cpl_parameter_get_double(par);
563
564 /* create final outgrid parameter structure */
565 int naxis = cpl_wcs_get_image_naxis(wcs);
566
567 if(naxis == 2) {
568 aParams_outputgrid =
570 ra_min, ra_max,
571 dec_min, dec_max,
572 fieldmargin);
573 } else {
574 //cpl_msg_info(cpl_func,"dx: %g dy: %g dlambda: %g", dx, dy, dlambda);
575 //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",
576 // ra_min, ra_max, dec_min, dec_max, lambda_min, lambda_max, fieldmargin);
577 aParams_outputgrid =
579 ra_min, ra_max,
580 dec_min, dec_max,
581 lambda_min, lambda_max,
582 fieldmargin);
583 }
584
585 double dist_x = (ra_max-ra_min)/dx;
586 double dist_y = (dec_max-dec_min)/dy;
587
588 cpl_msg_debug(cpl_func,"size1: %g size2: %g",dist_x, dist_y);
589 double dist = sqrt(dist_x * dist_x + dist_y * dist_y);
590 cpl_msg_info(cpl_func,"Output cube diagonal size: %g", dist);
591 if ( dist > MAX_COMBINE_CUBE_PLANE_DIAGONAL) {
592 cpl_msg_info(cpl_func,"Internal threshold on max cube plane diagonal size: %d",MAX_COMBINE_CUBE_PLANE_DIAGONAL);
593 cpl_msg_warning(cpl_func,"Resampled cube plane diagonal greater than threshold. Skip cube combination.");
594 cpl_error_set(cpl_func,CPL_ERROR_INCOMPATIBLE_INPUT);
595 eris_check_error_code("eris_ifu_resample_set_outputgrid");
596 return aParams_outputgrid;
597 }
598 cpl_free(context);
599 eris_check_error_code("eris_ifu_resample_set_outputgrid");
600
601 return aParams_outputgrid;
602}
603
604/*----------------------------------------------------------------------------*/
612/*----------------------------------------------------------------------------*/
613
614static hdrl_parameter *
615eris_ifu_resample_set_outputgrid2D (const char* recipe_name,
616 const cpl_parameterlist *parlist, cpl_table *muse_table, cpl_wcs *wcs)
617{
618 /* Should be done in a dedicated _new function */
619 cpl_error_ensure(parlist != NULL, CPL_ERROR_NULL_INPUT,
620 return NULL, "NULL Input Parameters");
621 cpl_error_ensure(muse_table != NULL, CPL_ERROR_NULL_INPUT,
622 return NULL, "NULL Input table");
623 cpl_error_ensure(wcs != NULL, CPL_ERROR_NULL_INPUT,
624 return NULL, "NULL Input wcs");
625
626 hdrl_parameter *aParams_outputgrid = NULL;
627 char* param_name = NULL;
628 char* context = NULL;
629 double ra_min = 0.;
630 double ra_max = 0.;
631 double dec_min = 0.;
632 double dec_max = 0.;
633
634 double dx = 0.;
635 double dy = 0.;
636
637 context = cpl_sprintf("eris.%s", recipe_name);
638 //cpl_table_save (muse_table, NULL, NULL, "muse_table.fits", CPL_IO_CREATE);
639 /* init relevant parameters */
640 double ra_min_tmp =
641 cpl_table_get_column_min (muse_table, HDRL_RESAMPLE_TABLE_RA);
642 double ra_max_tmp =
643 cpl_table_get_column_max (muse_table, HDRL_RESAMPLE_TABLE_RA);
644 double dec_min_tmp =
645 cpl_table_get_column_min (muse_table, HDRL_RESAMPLE_TABLE_DEC);
646 double dec_max_tmp =
647 cpl_table_get_column_max (muse_table, HDRL_RESAMPLE_TABLE_DEC);
648
649 //cpl_msg_info(cpl_func,"ra_max_tmp: %g ra_min_tmp: %g",ra_max_tmp, ra_min_tmp);
650 //cpl_msg_info(cpl_func,"dec_max_tmp: %g dec_min_tmp: %g",dec_max_tmp, dec_min_tmp);
651 /* We have the rare case that the image spans over ra = 0.*/
652 if(ra_max_tmp - ra_min_tmp > 180){
653 const double *ra = cpl_table_get_data_double_const(muse_table,
654 HDRL_RESAMPLE_TABLE_RA);
655 /* Should we also take the bpm into account ?
656 const int *bpm = cpl_table_get_data_int_const(ResTable,
657 HDRL_RESAMPLE_TABLE_BPM);
658 */
659
660 /* set to extreme values for a start */
661 ra_min_tmp = 0.;
662 ra_max_tmp = 360.;
663 cpl_size nrow = cpl_table_get_nrow(muse_table);
664
665 for (cpl_size i = 0; i < nrow; i++) {
666 if (ra[i] > ra_min_tmp && ra[i] <= 180.) ra_min_tmp = ra[i]; /* get the maximum */
667 if (ra[i] < ra_max_tmp && ra[i] > 180.) ra_max_tmp = ra[i]; /* get the minimum */
668 }
669 }
670
671 /* Check output (for debug) */
672 cpl_msg_debug (cpl_func, "min x %10.7f", ra_min_tmp);
673 cpl_msg_debug (cpl_func, "max x %10.7f", ra_max_tmp);
674 cpl_msg_debug (cpl_func, "min y %10.7f", dec_min_tmp);
675 cpl_msg_debug (cpl_func, "max y %10.7f", dec_max_tmp);
676
677 /*
678 param_name = cpl_sprintf("%s.outgrid.ra-min", context);
679 set_double_from_parameters(parlist, param_name, ra_min_tmp, &ra_min);
680 cpl_free(param_name);
681 */
682 ra_min = ra_min_tmp;
683 /*
684 param_name = cpl_sprintf("%s.outgrid.ra-max", context);
685 set_double_from_parameters(parlist, param_name, ra_max_tmp, &ra_max);
686 cpl_free(param_name);
687 */
688 ra_max = ra_max_tmp;
689 /*
690 param_name = cpl_sprintf("%s.outgrid.dec-min", context);
691 set_double_from_parameters(parlist, param_name, dec_min_tmp, &dec_min);
692 cpl_free(param_name);
693 */
694 dec_min = dec_min_tmp;
695
696 /*
697 param_name = cpl_sprintf("%s.outgrid.dec-max", context);
698 set_double_from_parameters(parlist, param_name, dec_max_tmp, &dec_max);
699 cpl_free(param_name);
700 */
701 dec_max = dec_max_tmp;
702
703
704
705 /* Reset all delta values */
706 double cdelt1 = 0., cdelt2 = 0., cdelt3 = 0.;
707 eris_ifu_resample_get_cdelt123 (wcs, &cdelt1, &cdelt2, &cdelt3);
708 /*
709 param_name = cpl_sprintf("%s.outgrid.delta-ra", context);
710 set_double_from_parameters(parlist, param_name, cdelt1, &dx);
711 cpl_free(param_name);
712 */
713 dx = cdelt1;
714 /*
715 param_name = cpl_sprintf("%s.outgrid.delta-dec", context);
716 set_double_from_parameters(parlist, param_name, cdelt2, &dy);
717 cpl_free(param_name);
718 */
719 dy = cdelt2;
720
721 /* Assign the field margin */
722 const cpl_parameter *par = NULL;
723 param_name = cpl_sprintf("%s.fieldmargin", context);
724 par = cpl_parameterlist_find_const (parlist, param_name);
725 cpl_free(param_name);
726
727 double fieldmargin = 0.;
728 fieldmargin = cpl_parameter_get_double(par);
729
730 /* create final outgrid parameter structure */
731 //int naxis = cpl_wcs_get_image_naxis(wcs);
732
733 aParams_outputgrid =
735 ra_min, ra_max,
736 dec_min, dec_max,
737 fieldmargin);
738 //cpl_msg_info(cpl_func,"ra_max: %g ra_min: %g dx: %g",ra_max, ra_min, dx);
739 //cpl_msg_info(cpl_func,"dec_max: %g dec_min: %g dy: %g",dec_max, dec_min, dy);
740 double dist_x = (ra_max-ra_min)/dx;
741 double dist_y = (dec_max-dec_min)/dy;
742
743 cpl_msg_debug(cpl_func,"size1: %g size2: %g",dist_x, dist_y);
744 double dist = sqrt(dist_x * dist_x + dist_y * dist_y);
745 cpl_msg_debug(cpl_func,"Output cube diagonal size: %g", dist);
746 if ( dist > MAX_COMBINE_CUBE_PLANE_DIAGONAL) {
747 cpl_msg_info(cpl_func,"Internal threshold on max cube plane diagonal size: %d",MAX_COMBINE_CUBE_PLANE_DIAGONAL);
748 cpl_msg_warning(cpl_func,"Resampled cube plane diagonal %g greater than threshold. Skip cube combination.", dist);
749 cpl_error_set(cpl_func,CPL_ERROR_INCOMPATIBLE_INPUT);
750 eris_check_error_code("eris_ifu_resample_set_outputgrid2D");
751 return aParams_outputgrid;
752 }
753 cpl_free(context);
754 eris_check_error_code("eris_ifu_resample_set_outputgrid2D");
755
756 return aParams_outputgrid;
757}
758
759
760/*----------------------------------------------------------------------------*/
775cpl_error_code
776eris_ifu_resample_save_cube(hdrl_resample_result *aCube,
777 const char *procatg,
778 const char *recipe,
779 const char *filename,
780 const cpl_parameterlist *parlist,
781 cpl_frameset *frameset,
782 cpl_boolean gen_phase3)
783{
784
785 cpl_ensure_code(aCube && aCube->header, CPL_ERROR_NULL_INPUT);
786
787// cpl_error_code rc;
788 //cpl_propertylist *header = NULL;
789 cpl_propertylist* dataHeader;
790 cpl_propertylist* errsHeader;
791 cpl_propertylist* qualHeader;
792
793
794 deqErrorType errorType = rmse;
795 deqQualityType qualityType = flag16bit;
796 char* pipe_id = cpl_sprintf("%s%s%s", PACKAGE, "/", PACKAGE_VERSION);
797
798 /* Create FITS headers for the extensions */
799 if(cpl_propertylist_has(aCube->header, "CDELT3")) {
800 cpl_propertylist_erase(aCube->header, "CDELT3");
801 }
802
803 if(!cpl_propertylist_has(aCube->header, CUNIT3)) {
804 cpl_propertylist_append_string(aCube->header,"CUNIT3","um");
805 }
806 if(!cpl_propertylist_has(aCube->header, "BUNIT")) {
807 cpl_propertylist_append_string(aCube->header,"BUNIT", UNIT_ADU);
808 }
809
810
811 dataHeader = cpl_propertylist_duplicate(aCube->header);
812 errsHeader = cpl_propertylist_duplicate(aCube->header);
813 qualHeader = cpl_propertylist_duplicate(aCube->header);
814 cpl_propertylist_update_string(qualHeader,"BUNIT", "");
815
816 eris_ifu_heades_add_hduclass_common(dataHeader, "IMAGE");
818
819 eris_ifu_heades_add_hduclass_common(errsHeader, "IMAGE");
820 eris_ifu_heades_add_hduclass_errs(errsHeader, errorType);
821
822 eris_ifu_heades_add_hduclass_common(qualHeader, "IMAGE");
823 eris_ifu_heades_add_hduclass_qual(qualHeader, qualityType);
824
825 if(cpl_propertylist_has(dataHeader,"CTYPE1")) {
826 const char* ctype1 = cpl_propertylist_get_string(dataHeader,"CTYPE1");
827 if (!cpl_propertylist_has(errsHeader,"CTYPE1"))
828 cpl_propertylist_append_string(errsHeader, "CTYPE1", ctype1);
829 if (!cpl_propertylist_has(qualHeader,"CTYPE1"))
830 cpl_propertylist_append_string(qualHeader, "CTYPE1", ctype1);
831 }
832 if(cpl_propertylist_has(dataHeader,"CTYPE2")) {
833 const char* ctype2 = cpl_propertylist_get_string(dataHeader,"CTYPE2");
834 if (!cpl_propertylist_has(errsHeader,"CTYPE2"))
835 cpl_propertylist_append_string(errsHeader, "CTYPE2",ctype2);
836 if (!cpl_propertylist_has(errsHeader,"CTYPE2"))
837 cpl_propertylist_append_string(qualHeader, "CTYPE2", ctype2);
838 }
839
840 cpl_propertylist* applist = cpl_propertylist_duplicate(aCube->header);
841 if(gen_phase3) {
842 /*
843 if(cpl_propertylist_has(aCube->header,"MJD-OBS")) {
844 double mjd_obs = cpl_propertylist_get_double(aCube->header, "MJD-OBS");
845 cpl_propertylist_append_double(header, "MJD-OBS",mjd_obs);
846 }
847 */
848 eris_ifu_sdp_properties * prop = eris_ifu_sdp_properties_collect(aCube,
849 frameset, parlist, recipe);
850
851 eris_ifu_sdp_properties_update(applist, prop);
852 //cpl_propertylist_dump(applist,stdout);
853 eris_ifu_sdp_properties_delete(prop);
854 }
855
856 /* Create an empty primary HDU with only the header information*/
857 eris_ifu_plist_erase_wcs(applist);
858 eris_pfits_clean_header(dataHeader, CPL_TRUE);
859 eris_pfits_clean_header(errsHeader, CPL_TRUE);
860 eris_pfits_clean_header(qualHeader, CPL_TRUE);
861 cpl_propertylist_append_string(applist, "PRODCATG", "SCIENCE.CUBE.IFS");
862
863 cpl_propertylist_update_string(applist, CPL_DFS_PRO_CATG, procatg);
864 if(cpl_propertylist_has(applist,"BUNIT")) {
865 cpl_propertylist_erase(applist,"BUNIT");
866 }
867 if(cpl_propertylist_has(applist,"COMMENT")) {
868 cpl_msg_warning(cpl_func,"remove comments");
869 cpl_propertylist_erase(applist,"COMMENT");
870 }
871 if(cpl_propertylist_has(dataHeader,"COMMENT")) {
872 cpl_msg_warning(cpl_func,"remove comments");
873 cpl_propertylist_erase(dataHeader,"COMMENT");
874 }
875 cpl_dfs_save_propertylist(frameset, NULL, parlist, frameset,
876 NULL, recipe, applist, "RADECSYS", pipe_id, filename);
877 cpl_free(pipe_id);
878 cpl_imagelist* data = cpl_imagelist_new();
879 cpl_imagelist* errs = cpl_imagelist_new();
880 cpl_imagelist* qual = cpl_imagelist_new();
881
882 eris_ifu_split3_hdrl_imagelist(aCube->himlist, data, errs, qual);
883
884 /* Write the extension units */
885 cpl_type data_type = cpl_image_get_type(cpl_imagelist_get(data,0));
886
887 cpl_imagelist_save(data, filename, data_type, dataHeader, CPL_IO_EXTEND);
888 cpl_imagelist_save(errs, filename, data_type, errsHeader, CPL_IO_EXTEND);
889 cpl_imagelist_save(qual, filename, CPL_TYPE_INT, qualHeader, CPL_IO_EXTEND);
890 cpl_size ix_max = cpl_imagelist_get_size(data);
891 cpl_image * q;
892 for (cpl_size ix=ix_max-1; ix >= 0; ix--) {
893 cpl_imagelist_unset(data, ix);
894 cpl_imagelist_unset(errs, ix);
895 q=cpl_imagelist_unset(qual, ix);
896 cpl_image_delete(q);
897 }
898
899 cpl_imagelist_delete(data);
900 cpl_imagelist_delete(errs);
901 cpl_imagelist_delete(qual);
902
903 cpl_propertylist_delete(applist);
904 /* changed up to here */
905
906 /* this creates a problem */
907// eris_ifu_plist_erase_wcs(header);
908
909
910// if (rc != CPL_ERROR_NONE) {
911// return cpl_error_get_code();
912// }
913
914
915// cpl_size planes = hdrl_imagelist_get_size(aCube->himlist);
916
917
918 //cpl_propertylist_delete(header);
919 eris_ifu_free_propertylist(&dataHeader);
920 eris_ifu_free_propertylist(&errsHeader);
921 eris_ifu_free_propertylist(&qualHeader);
922 eris_check_error_code("eris_ifu_resample_save_cube");
923 return cpl_error_get_code();
924}
925
926
940static cpl_wcs *
941eris_ifu_resample_get_wcs_from_frameset(cpl_frameset* frameset,
942 const char* procatg) {
943
944 cpl_ensure(frameset, CPL_ERROR_NULL_INPUT, NULL);
945
946 /* Read wcs from firts raw image */
947 cpl_frameset* in_set = NULL;
948
949 if ((in_set = eris_ifu_extract_frameset(frameset, procatg)) == NULL){
950 cpl_error_set_message(cpl_func, CPL_ERROR_FILE_NOT_FOUND, "Missing RAW "
951 "files");
952 return NULL;
953 }
954
955 cpl_frame* frame = cpl_frameset_get_position(in_set, 0);
956 const char* fname = cpl_frame_get_filename(frame);
957
958
959 cpl_wcs *wcs = NULL;
960 cpl_propertylist* head = NULL;
961
962 cpl_errorstate prestate = cpl_errorstate_get();
963 head = cpl_propertylist_load(fname, 0);
964 wcs = cpl_wcs_new_from_propertylist(head);
965 if (wcs == NULL) {
966 /* Not possible to read wcs - trying from extension */
967 cpl_errorstate_set(prestate);
968 cpl_propertylist_delete(head);
969 } else {
970 cpl_propertylist_delete(head);
971 cpl_frameset_delete(in_set);
972 return wcs;
973 }
974
975 prestate = cpl_errorstate_get();
976 head = cpl_propertylist_load(fname, 1);
977 wcs = cpl_wcs_new_from_propertylist(head);
978 cpl_propertylist_delete(head);
979 cpl_frameset_delete(in_set);
980 eris_check_error_code("eris_ifu_resample_get_wcs_from_frameset");
981 return wcs ;
982
983}
984
985
1002static cpl_table *
1003eris_ifu_resample_frameset_to_table(const cpl_frame *data_frm,
1004 const cpl_size data_ext_id,
1005 const cpl_frame *errs_frm,
1006 const cpl_size errs_ext_id,
1007 const cpl_frame *qual_frm,
1008 const cpl_size qual_ext_id,
1009 const cpl_boolean is_variance,
1010 const cpl_boolean subtract_bkg,
1011 const int edgetrim)
1012{
1013
1014 cpl_ensure(data_frm, CPL_ERROR_NULL_INPUT, NULL);
1015 cpl_msg_severity level = cpl_msg_get_level();
1016 cpl_msg_set_level(CPL_MSG_INFO);
1017 const char *name = cpl_frame_get_filename(data_frm);
1018 cpl_imagelist *dlist = NULL;
1019 cpl_imagelist *elist = NULL;
1020 cpl_imagelist *qlist = NULL;
1021
1022 cpl_table *tab_ext = NULL;
1023
1024 cpl_errorstate prestate = cpl_errorstate_get();
1025
1026 dlist = cpl_imagelist_load(name, CPL_TYPE_DOUBLE, data_ext_id);
1027
1028 if (dlist == NULL) {
1029 /*It was not an image nor a imagelist - reset cpl error and return NULL*/
1030 if (!cpl_errorstate_is_equal(prestate)) {
1031 cpl_errorstate_set(prestate);
1032 }
1033 return NULL;
1034 }
1035
1036 cpl_propertylist *xheader_data = cpl_propertylist_load(name,
1037 data_ext_id);
1038 cpl_wcs *wcs = cpl_wcs_new_from_propertylist(xheader_data);
1039
1040 if (errs_frm != NULL) {
1041 name = cpl_frame_get_filename(errs_frm);
1042 elist = cpl_imagelist_load(name, CPL_TYPE_DOUBLE, errs_ext_id);
1043 if(is_variance) {
1044 /* transform variance into error */
1045 cpl_imagelist_power(elist, 0.5);
1046 }
1047 }
1048 if (qual_frm != NULL) {
1049 name = cpl_frame_get_filename(qual_frm);
1050 qlist = cpl_imagelist_load(name, CPL_TYPE_INT, qual_ext_id);
1051 }
1052
1053 cpl_size size = cpl_imagelist_get_size(dlist);
1054 /* ingest pixel quality in data */
1055
1056 if (qual_frm != NULL && size > 0){
1057 for(cpl_size k = 0; k < size; k++) {
1058 cpl_image* data = cpl_imagelist_get(dlist, k);
1059 cpl_image* qual = cpl_imagelist_get(qlist, k);
1060
1061 /*we use INT_MAX instead of 1.1 as some pipeline
1062 * may use pixel codes as qualifier */
1063 cpl_mask* mask = cpl_mask_threshold_image_create(qual, 0, INT_MAX);
1064
1065 cpl_image_reject_from_mask(data, mask);
1066 cpl_mask_delete(mask);
1067 cpl_imagelist_set(dlist, data, k);
1068 }
1069 }
1070
1071 /* In case the error is passed as variance we take the square root. This could
1072 * cause on some data to add a bad pixel mask - thus we have to synchronize it
1073 * */
1074 if(elist != NULL && is_variance){
1075 for(cpl_size k = 0; k < size; k++) {
1076 cpl_image* data = cpl_imagelist_get(dlist, k);
1077 cpl_mask* data_mask = cpl_image_get_bpm(data);
1078
1079 cpl_image* error = cpl_imagelist_get(elist, k);
1080 cpl_mask* error_mask = cpl_image_get_bpm(error);
1081
1082 cpl_mask* merged = cpl_mask_duplicate(data_mask);
1083
1084 /*Add original bad pixels to previous iteration*/
1085 cpl_mask_or(merged, error_mask);
1086 cpl_image_reject_from_mask(data, merged);
1087 cpl_image_reject_from_mask(error, merged);
1088 cpl_mask_delete(merged);
1089
1090 }
1091 }
1092
1093 hdrl_imagelist* hlist = NULL;
1094 if (size > 0) {
1095 hlist = hdrl_imagelist_create(dlist, elist);
1096 }
1097 int edge_trim = edgetrim;
1098 if (edge_trim > 0) {
1099 cpl_msg_info(cpl_func, "Trim input image edges of %d pixels", edge_trim);
1100 /* trim each product plane edge of edgetrim pixels */
1101 cpl_size sx = hdrl_imagelist_get_size_x(hlist);
1102 cpl_size sy = hdrl_imagelist_get_size_y(hlist);
1103 cpl_size sz = hdrl_imagelist_get_size(hlist);
1104 if (edge_trim >= 0.5* sx) {
1105 edge_trim = 0;
1106 cpl_msg_warning(cpl_func, "edge-trim must be smaller than half image "
1107 "size. Reset to 0");
1108 }
1109 if (edge_trim >= 0.5* sy) {
1110 edge_trim = 0;
1111 cpl_msg_warning(cpl_func, "edge-trim must be smaller than half image "
1112 "size. Reset to 0");
1113 }
1114 for(cpl_size k = 0; k < sz; k++) {
1115
1116 hdrl_image* hima = hdrl_imagelist_get(hlist, k);
1117 /* Note FITS convention to loop over pixels */
1118 /* trim lower image border along X direction */
1119 for(cpl_size j = 1; j <= edge_trim; j++) {
1120 for(cpl_size i = 1; i <= sx; i++) {
1121 hdrl_image_reject(hima, i, j);
1122 }
1123 }
1124 /* trim upper image border along X direction */
1125 for(cpl_size j = sy - edge_trim + 1; j <= sy; j++) {
1126 for(cpl_size i = 1; i <= sx; i++) {
1127 hdrl_image_reject(hima, i, j);
1128 }
1129 }
1130 /* trim left image border along Y direction */
1131 for(cpl_size j = 1; j <= sy; j++) {
1132 for(cpl_size i = 1; i <= edge_trim; i++) {
1133 hdrl_image_reject(hima, i, j);
1134 }
1135 }
1136 /* trim right image border along Y direction */
1137 for(cpl_size j = 1; j <= sy; j++) {
1138 for(cpl_size i = sx - edge_trim + 1; i <= sx; i++) {
1139 hdrl_image_reject(hima, i, j);
1140 }
1141 }
1142 }
1143 }
1144
1145 /* single list are not anymore needed as data are copied to the hdrl_imagelist
1146 * therefore we free them */
1147 cpl_imagelist_delete(dlist);
1148 if (qual_frm != NULL)
1149 cpl_imagelist_delete(qlist);
1150 if (errs_frm != NULL)
1151 cpl_imagelist_delete(elist);
1152
1153 /* Do the calculation */
1154 if (size == 1){
1155 /* Subtract the background if we have an image and not a cube */
1156 cpl_msg_info(cpl_func, "Reading the image ...");
1157 hdrl_image* hima = hdrl_imagelist_get(hlist, 0);
1158
1159 /* Subtract the background on request */
1160 if(subtract_bkg == CPL_TRUE) {
1161 cpl_msg_info(cpl_func, "Subtracting the median as requested ...");
1163 }
1164
1166 tab_ext = hdrl_resample_image_to_table(hima, wcs);
1167 } else {
1168 cpl_msg_info(cpl_func, "Converting imagelist to table with hdrl function");
1169 tab_ext = hdrl_resample_imagelist_to_table(hlist, wcs);
1170 }
1171
1172 /*Cleanup the memory */
1173 if (hlist != NULL)
1174 hdrl_imagelist_delete(hlist);
1175 cpl_wcs_delete(wcs);
1176 cpl_propertylist_delete(xheader_data);
1177 eris_check_error_code("eris_ifu_resample_frameset_to_table");
1178 cpl_msg_set_level(level);
1179 return tab_ext;
1180
1181}
1182
1183
1184/*----------------------------------------------------------------------------*/
1190/*----------------------------------------------------------------------------*/
1191
1192static cpl_error_code
1193eris_ifu_resample_wcs_print(cpl_wcs *wcs)
1194{
1195 cpl_ensure_code(wcs, CPL_ERROR_NULL_INPUT);
1196
1197 const cpl_array *crval = cpl_wcs_get_crval(wcs);
1198 const cpl_array *crpix = cpl_wcs_get_crpix(wcs);
1199 const cpl_array *ctype = cpl_wcs_get_ctype(wcs);
1200 const cpl_array *cunit = cpl_wcs_get_cunit(wcs);
1201
1202 const cpl_matrix *cd = cpl_wcs_get_cd(wcs);
1203 const cpl_array *dims = cpl_wcs_get_image_dims(wcs);
1204 cpl_size naxis = cpl_wcs_get_image_naxis(wcs);
1205
1206 cpl_msg_info(cpl_func, "NAXIS: %lld", naxis);
1207 int testerr = 0;
1208
1209 cpl_msg_indent_more();
1210 /* Check NAXIS */
1211 for (cpl_size i = 0; i < naxis; i++) {
1212 cpl_msg_info(cpl_func, "NAXIS%lld: %d", i + 1,
1213 cpl_array_get_int(dims, i, &testerr));
1214 }
1215 cpl_msg_indent_less();
1216
1217 double cd11 = cpl_matrix_get(cd, 0, 0);
1218 double cd12 = cpl_matrix_get(cd, 0, 1);
1219 double cd21 = cpl_matrix_get(cd, 1, 0);
1220 double cd22 = cpl_matrix_get(cd, 1, 1);
1221 double crpix1 = cpl_array_get_double(crpix, 0, &testerr);
1222 double crpix2 = cpl_array_get_double(crpix, 1, &testerr);
1223 double crval1 = cpl_array_get_double(crval, 0, &testerr);
1224 double crval2 = cpl_array_get_double(crval, 1, &testerr);
1225
1226 cpl_msg_info(cpl_func, "1st and 2nd dimension");
1227 cpl_msg_indent_more();
1228 cpl_msg_info(cpl_func, "CD1_1: %g", cd11);
1229 cpl_msg_info(cpl_func, "CD1_2: %g", cd12);
1230 cpl_msg_info(cpl_func, "CD2_1: %g", cd21);
1231 cpl_msg_info(cpl_func, "CD2_2: %g", cd22);
1232
1233 cpl_msg_info(cpl_func, "CRPIX1: %g", crpix1);
1234 cpl_msg_info(cpl_func, "CRPIX2: %g", crpix2);
1235 cpl_msg_info(cpl_func, "CRVAL1: %f", crval1);
1236 cpl_msg_info(cpl_func, "CRVAL2: %f", crval2);
1237 if (ctype) {
1238 cpl_msg_info(cpl_func, "CTYPE1: %s", cpl_array_get_string(ctype, 0));
1239 cpl_msg_info(cpl_func, "CTYPE2: %s", cpl_array_get_string(ctype, 1));
1240 }
1241
1242 if (cunit) {
1243 cpl_msg_info(cpl_func, "CUNIT1: %s", cpl_array_get_string(cunit, 0));
1244 cpl_msg_info(cpl_func, "CUNIT2: %s", cpl_array_get_string(cunit, 1));
1245 }
1246 cpl_msg_indent_less();
1247
1248 /* Is it a 3D cube or a 2D image */
1249 cpl_size cdncol = cpl_matrix_get_ncol(cd);
1250 if (cdncol == 3) {
1251
1252 double cd13 = cpl_matrix_get(cd, 0, 2);
1253 double cd23 = cpl_matrix_get(cd, 1, 2);
1254 double cd31 = cpl_matrix_get(cd, 2, 0);
1255 double cd32 = cpl_matrix_get(cd, 2, 1);
1256 double cd33 = cpl_matrix_get(cd, 2, 2);
1257 double crval3 = cpl_array_get_double(crval, 2, &testerr);
1258 double crpix3 = cpl_array_get_double(crpix, 2, &testerr);
1259
1260 cpl_msg_info(cpl_func, "3rd dimension");
1261 cpl_msg_indent_more();
1262 cpl_msg_info(cpl_func, "CD1_3: %g", cd13);
1263 cpl_msg_info(cpl_func, "CD2_3: %g", cd23);
1264 cpl_msg_info(cpl_func, "CD3_1: %g", cd31);
1265 cpl_msg_info(cpl_func, "CD3_2: %g", cd32);
1266 cpl_msg_info(cpl_func, "CD3_3: %g", cd33);
1267
1268 cpl_msg_info(cpl_func, "CRPIX3: %g", crpix3);
1269 cpl_msg_info(cpl_func, "CRVAL3: %g", crval3);
1270
1271 if (ctype) {
1272 cpl_msg_info(cpl_func, "CTYPE3: %s", cpl_array_get_string(ctype, 2));
1273 }
1274
1275 if (cunit) {
1276 cpl_msg_info(cpl_func, "CUNIT3: %s", cpl_array_get_string(cunit, 2));
1277 }
1278
1279 cpl_msg_indent_less();
1280 }
1281 eris_check_error_code("eris_ifu_resample_wcs_print");
1282 return cpl_error_get_code();
1283}
1284
1285
1302static cpl_table*
1303eris_ifu_resample_get_table_from_frameset(const cpl_frame* data_frm,
1304 const cpl_frame* errs_frm,
1305 const cpl_frame* qual_frm,
1306 const cpl_size data_ext_id,
1307 const cpl_size errs_ext_id,
1308 const cpl_size qual_ext_id,
1309 const cpl_boolean is_variance,
1310 const cpl_boolean subtract_bkg,
1311 const int edge_trim) {
1312
1313 cpl_ensure(data_frm, CPL_ERROR_NULL_INPUT, NULL);
1314 /* Assumption - all have the same number of extensions */
1315 cpl_size next = cpl_frame_get_nextensions(data_frm);
1316 cpl_msg_info(cpl_func, "Analysing and processing file %s",
1317 cpl_frame_get_filename(data_frm));
1318
1319 /* Assumption
1320 * if data_ext_id == -1 the loop can be done from 0 to next for data, error,
1321 * and quality; only the filename changes
1322 * if data_ext_id != -1 no loop is needed and is difficult to do in the loop
1323 * */
1324
1325 cpl_msg_indent_more();
1326 cpl_table* table = NULL;
1327 if(data_ext_id != -1) {
1328 cpl_msg_info(cpl_func, "Extension: %02lld", data_ext_id);
1329 /* No Loop, only use the extension given by the user */
1330 table = eris_ifu_resample_frameset_to_table(data_frm, data_ext_id,
1331 errs_frm, errs_ext_id,
1332 qual_frm, qual_ext_id,
1333 is_variance, subtract_bkg, edge_trim);
1334 } else {
1335 /* Loop */
1336 for (cpl_size i = 0; i < next + 1; i++ ) {
1337 cpl_table * table_local = NULL;
1338 cpl_msg_info(cpl_func,"Extension: %02lld", i);
1339 table_local = eris_ifu_resample_frameset_to_table(data_frm, i,
1340 errs_frm, i,
1341 qual_frm, i,
1342 is_variance, subtract_bkg, edge_trim);
1343
1344 /* If table_local == NULL go to the next extension */
1345 if (table_local == NULL) {
1346 cpl_msg_info(cpl_func, "No siutable data found - continuing");
1347 continue;
1348 }
1349
1350 /* Now table_local is always != NULL */
1351 if (table == NULL){
1352 table = table_local;
1353 continue;
1354 } else {
1355 cpl_size nrow = cpl_table_get_nrow(table);
1356 cpl_table_insert(table, table_local, nrow);
1357 cpl_table_delete(table_local);
1358 table_local = NULL;
1359 }
1360 }
1361 }
1362 cpl_msg_indent_less();
1363 eris_check_error_code("eris_ifu_resample_get_table_from_frameset");
1364
1365 return table;
1366}
1367
1370// * @internal
1371// * @brief Find the aperture(s) with the greatest flux
1372// * @param self The aperture object
1373// * @param ind The aperture-indices in order of decreasing flux
1374// * @param nfind Number of indices to find
1375// * @return CPL_ERROR_NONE or the relevant _cpl_error_code_ on error
1376// *
1377// * nfind must be at least 1 and at most the size of the aperture object.
1378// *
1379// * The ind array must be able to hold (at least) nfind integers.
1380// * On success the first nfind elements of ind point to indices of the
1381// * aperture object.
1382// *
1383// * To find the single ind of the aperture with the maximum flux use simply:
1384// * int ind;
1385// * eris_apertures_find_max_flux(self, &ind, 1);
1386// *
1387// */
1389//static cpl_error_code eris_apertures_find_max_flux(const cpl_apertures * self,
1390// int * ind, int nfind)
1391//{
1392// const int nsize = cpl_apertures_get_size(self);
1393// int ifind;
1394
1395
1396// cpl_ensure_code(nsize > 0, cpl_error_get_code());
1397// cpl_ensure_code(ind, CPL_ERROR_NULL_INPUT);
1398// cpl_ensure_code(nfind > 0, CPL_ERROR_ILLEGAL_INPUT);
1399// cpl_ensure_code(nfind <= nsize, CPL_ERROR_ILLEGAL_INPUT);
1400
1401// for (ifind=0; ifind < nfind; ifind++) {
1402// double maxflux = -1;
1403// int maxind = -1;
1404// int i;
1405// for (i=1; i <= nsize; i++) {
1406// int k;
1407
1408// /* The flux has to be the highest among those not already found */
1409// for (k=0; k < ifind; k++) if (ind[k] == i) break;
1410
1411// if (k == ifind) {
1412// /* i has not been inserted into ind */
1413// const double flux = cpl_apertures_get_flux(self, i);
1414
1415// if (maxind < 0 || flux > maxflux) {
1416// maxind = i;
1417// maxflux = flux;
1418// }
1419// }
1420// }
1421// ind[ifind] = maxind;
1422// }
1423// eris_check_error_code("eris_apertures_find_max_flux");
1424// return CPL_ERROR_NONE;
1425
1426//}
1427
1429// * @brief filter input mask by kernel of given sizes and mode
1430// *
1431// * @param input_mask input mask
1432// * @param kernel_nx kernel X size
1433// * @param kernel_ny kernel X size
1434// * @param filter filter mode
1435// * @return mask corresponding to what filtered
1436// * @doc
1437// *
1438// * */
1439//static cpl_mask *
1440//eris_bpm_filter(
1441// const cpl_mask * input_mask,
1442// cpl_size kernel_nx,
1443// cpl_size kernel_ny,
1444// cpl_filter_mode filter)
1445//{
1446// cpl_mask * kernel = NULL;
1447// cpl_mask * filtered_mask = NULL;
1448// cpl_mask * expanded_mask = NULL;
1449// cpl_mask * expanded_filtered_mask = NULL;
1450
1451// /* Check Entries */
1452// cpl_ensure(input_mask != NULL, CPL_ERROR_NULL_INPUT, NULL);
1453// cpl_ensure(kernel_nx >= 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
1454// cpl_ensure(kernel_ny >= 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
1455// cpl_ensure(filter == CPL_FILTER_EROSION || filter == CPL_FILTER_DILATION ||
1456// filter == CPL_FILTER_OPENING || filter == CPL_FILTER_CLOSING,
1457// CPL_ERROR_ILLEGAL_INPUT, NULL);
1458
1459// /* Only odd-sized masks allowed */
1460// cpl_ensure((kernel_nx&1) == 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
1461// cpl_ensure((kernel_ny&1) == 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
1462
1463// kernel = cpl_mask_new(kernel_nx, kernel_ny);
1464// cpl_mask_not(kernel); /* All values set to unity*/
1465
1466// /* Enlarge the original mask with the kernel size and assume that outside
1467// * all pixels are good */
1468// expanded_mask = cpl_mask_new(
1469// cpl_mask_get_size_x(input_mask) + 2 * kernel_nx,
1470// cpl_mask_get_size_y(input_mask) + 2 * kernel_ny);
1471
1472// cpl_mask_copy(expanded_mask, input_mask, kernel_nx + 1, kernel_ny +1 );
1473
1474// expanded_filtered_mask = cpl_mask_new(cpl_mask_get_size_x(expanded_mask),
1475// cpl_mask_get_size_y(expanded_mask));
1476
1477// if(cpl_mask_filter(expanded_filtered_mask, expanded_mask, kernel, filter,
1478// CPL_BORDER_ZERO) != CPL_ERROR_NONE) {
1479
1480// cpl_mask_delete(kernel);
1481// cpl_mask_delete(expanded_filtered_mask);
1482// cpl_mask_delete(expanded_mask);
1483// return NULL;
1484// }
1485
1486// /* Extract the original mask from the expanded mask */
1487// filtered_mask = cpl_mask_extract(expanded_filtered_mask,
1488// kernel_nx+1, kernel_ny + 1,
1489// cpl_mask_get_size_x(input_mask) + kernel_nx,
1490// cpl_mask_get_size_y(input_mask) + kernel_ny);
1491
1492
1493// /* Free memory */
1494// cpl_mask_delete(kernel);
1495// cpl_mask_delete(expanded_filtered_mask);
1496// cpl_mask_delete(expanded_mask);
1497// eris_check_error_code("eris_bpm_filter");
1498// return filtered_mask;
1499//}
1500
1502// * @brief generate mask corresponding to post-filtered image according to method
1503// *
1504// * @param std_med_ima input image
1505// * @param pfm post filter mask method
1506// * @param pfx post filer x size
1507// * @param pfy post filer x size
1508// * @param kappa value used in kappa-sigma clipping
1509// * @return cpl_error_code
1510// * @doc
1511// *
1512// * */
1513//static cpl_mask*
1514//eris_object_mask(cpl_image* std_med_ima, const char* pfm, const int pfx,
1515// const int pfy, const double kappa)
1516//{
1517// cpl_mask* object_mask;
1518// double mad;
1519// double median;
1520// cpl_mask* nan_mask;
1521// /* flag NAN in image (due to SPIFFIER BP conventions) before computing mad */
1522// cpl_image_reject_value(std_med_ima,CPL_VALUE_NAN);
1523// /* extract mask of NANs to be used later */
1524// nan_mask = cpl_image_get_bpm(std_med_ima);
1525// //cpl_mask_save(nan_mask,"nan_mask.fits",NULL,CPL_IO_CREATE);
1526// median = cpl_image_get_mad(std_med_ima, &mad);
1527// double low_cut = median - 10 * kappa * mad * CPL_MATH_STD_MAD;
1528// double hi_cut = median + kappa * mad * CPL_MATH_STD_MAD;
1529// /* takes as low cut a very small value to be sure to get all negative pixels
1530// * as good */
1531// low_cut = -1.e10;
1532// //cpl_msg_info(cpl_func, "median: %g mad: %g low: %g hi : %g", median, mad, low_cut, hi_cut);
1533// /* detect object (and other previously flagged values) as good pixels */
1534// object_mask = cpl_mask_threshold_image_create(std_med_ima, low_cut, hi_cut);
1535// //cpl_mask_not(object_mask);
1536
1537// //cpl_mask_save(object_mask,"object_mask.fits",NULL,CPL_IO_CREATE);
1538// /* to detect only object pixels as good one make or with previously flagged
1539// * NAN values
1540// */
1541
1542// cpl_mask_or(object_mask,nan_mask);
1543// //cpl_mask_save(object_mask,"or_mask.fits",NULL,CPL_IO_CREATE);
1544// /* increase the mask of a few pixels to be sure to include all object values
1545// * for this reason we use method dilation and we make a NOT on the image
1546// * because the DILATION enlarge the area of BAD pixels, not the GOOD ones
1547// * but we want to increase the area covered by the object (flagged as good)
1548// */
1549// cpl_filter_mode filter_mode = CPL_FILTER_EROSION ;
1550// cpl_msg_warning(cpl_func, "pfm = %s", pfm);
1551// if( strcmp(pfm,"erosion") == 0 ) {
1552// cpl_msg_info(cpl_func, "Filter erosion");
1553// filter_mode = CPL_FILTER_EROSION ;
1554// } else if( strcmp(pfm, "dilation") == 0 ) {
1555// cpl_msg_info(cpl_func, "Filter dilation");
1556// filter_mode = CPL_FILTER_DILATION ;
1557// } else if( strcmp(pfm,"closing") == 0 ) {
1558// cpl_msg_info(cpl_func, "Filter closing");
1559// filter_mode = CPL_FILTER_CLOSING ;
1560// }
1561// //cpl_filter_mode filter_mode = CPL_FILTER_DILATION ;
1562
1563// cpl_mask_not(object_mask);
1564// cpl_mask* obj_mask_filtered = eris_bpm_filter(object_mask, pfx, pfy, filter_mode);
1565// /* To have again a proper mask with object flagged as good we do a NOT */
1566// cpl_mask_not(obj_mask_filtered);
1567// //cpl_mask_save(obj_mask_filtered,"filtered_mask.fits",NULL,CPL_IO_CREATE);
1568// /* clean-up memory */
1569// cpl_mask_delete(object_mask);
1570// //cpl_mask_delete(nan_mask);
1571// eris_check_error_code("eris_object_mask");
1572// return obj_mask_filtered;
1573//}
1574
1586cpl_error_code
1588 cpl_frameset* frameset,
1589 const cpl_parameterlist* parlist,
1590 const char* recipe_name,
1591// const char* filenameSpec,
1592 cpl_boolean apply_flat,
1593 cpl_boolean is_pupil)
1594{
1595 hdrl_image* cube_collapsed = NULL;
1596 cpl_image* cube_cmap = NULL;
1597
1598 cpl_frame* frame = cpl_frameset_find(frameset, pro_catg);
1599 const char* cube_fname = cpl_frame_get_filename(frame);
1600 cpl_propertylist* pheader = cpl_propertylist_load(cube_fname, 0);
1601 cpl_propertylist* head_wcs_2d = eris_ifu_plist_extract_wcs2D(pheader);
1602 cpl_propertylist* dheader = cpl_propertylist_load(cube_fname, 1);
1603 cpl_propertylist* eheader = cpl_propertylist_load(cube_fname, 2);
1604 cpl_propertylist* qheader = cpl_propertylist_load(cube_fname, 3);
1605
1606 cpl_imagelist* iml_data = cpl_imagelist_load(cube_fname, CPL_TYPE_DOUBLE, 1);
1607 cpl_imagelist* iml_errs = cpl_imagelist_load(cube_fname, CPL_TYPE_DOUBLE, 2);
1608 cpl_imagelist* iml_qual = cpl_imagelist_load(cube_fname, CPL_TYPE_DOUBLE, 3);
1609
1610 cpl_propertylist_append(dheader, head_wcs_2d);
1611 cpl_propertylist_append(eheader, head_wcs_2d);
1612 cpl_propertylist_append(qheader, head_wcs_2d);
1613 cpl_propertylist_delete(head_wcs_2d);
1614 double dit = eris_pfits_get_dit(pheader);
1615 int ndit = eris_pfits_get_ndit(pheader);
1616 double exptime = dit * ndit;
1617
1618 cpl_size sz = cpl_imagelist_get_size(iml_data);
1619 cpl_image* data = NULL;
1620 cpl_image* errs = NULL;
1621 cpl_image* qual = NULL;
1622 cpl_mask* bpm_data = NULL;
1623 cpl_mask* bpm_errs = NULL;
1624 cpl_boolean edge_trim = CPL_TRUE;
1625 cpl_size k_center = 0.5 * sz;
1626 for(cpl_size k = 0; k < sz; k++) {
1627 data = cpl_imagelist_get(iml_data, k);
1628 errs = cpl_imagelist_get(iml_errs, k);
1629 qual = cpl_imagelist_get(iml_qual, k);
1630 cpl_image_reject_value(data, CPL_VALUE_NAN);
1631 cpl_image_reject_value(errs, CPL_VALUE_NAN);
1632 bpm_data = cpl_image_get_bpm(data);
1633 bpm_errs = cpl_image_get_bpm(errs);
1634 cpl_mask * bpm_mask = cpl_mask_threshold_image_create(qual, -0.5, 0.5) ;
1635 cpl_mask_not(bpm_mask) ;
1636 cpl_mask_or(bpm_mask,bpm_data);
1637 cpl_mask_or(bpm_mask,bpm_errs);
1638 cpl_image_reject_from_mask(data, bpm_mask);
1639 cpl_image_reject_from_mask(errs, bpm_mask);
1640
1641 //EXPOSURE MAP CREATION
1642
1643
1644 /* exposure map: use the plane at the center of the wavelength range
1645 * then replace intensity with EXPTIME value and error with 1 and
1646 * generate HDRL image. This will be re-sampled as the cube and give
1647 * the exposure map associated to the combined cube
1648 */
1649
1650 cpl_image* ima_time = NULL;
1651 cpl_image* err_time = NULL;
1652 hdrl_image *himg_exptime = NULL;
1653 if(apply_flat && k == k_center) {
1654 ima_time = cpl_image_new(cpl_image_get_size_x(data),
1655 cpl_image_get_size_y(data), CPL_TYPE_DOUBLE);
1656 err_time = cpl_image_new(cpl_image_get_size_x(data),
1657 cpl_image_get_size_y(data), CPL_TYPE_DOUBLE);
1658
1659
1660 cpl_image_add_scalar(ima_time, exptime);
1661 cpl_image_add_scalar(err_time, sqrt(exptime));
1662 himg_exptime = hdrl_image_create(ima_time, err_time);
1663 hdrl_image_reject_from_mask(himg_exptime, bpm_mask);
1664
1665 if (edge_trim > 0) {
1666 eris_ifu_resample_trim_edge(himg_exptime, edge_trim);
1667 }
1668 /*
1669 char* fname = cpl_sprintf("ima_time_%lld.fits",i);
1670 cpl_image_save(ima_time,fname,CPL_TYPE_DOUBLE,NULL,CPL_IO_DEFAULT);
1671 cpl_free(fname);
1672 fname = cpl_sprintf("err_time_%lld.fits",i);
1673 cpl_image_save(err_time,fname,CPL_TYPE_DOUBLE,NULL,CPL_IO_DEFAULT);
1674 cpl_free(fname);
1675 */
1676 cpl_image_delete(ima_time);
1677 cpl_image_delete(err_time);
1678
1679 /* SAVE EXPOSURE MAP */
1680
1681 cpl_propertylist* phead2D = cpl_propertylist_duplicate(dheader);
1682 char* fname = cpl_sprintf("%s_exposure_map.fits", recipe_name);
1683 eris_ifu_plist_erase_wcs3D(phead2D);
1684 eris_ifu_plist_erase_expmap_extra_keys(phead2D);
1685 cpl_propertylist_update_string(phead2D, CPL_DFS_PRO_CATG, ERIS_IFU_PRO_JITTER_EXPOSURE_MAP);
1686 cpl_propertylist_update_string(phead2D, "PRODCATG", PRODCATG_EXPOSUREMAP);
1687 eris_ifu_save_image_phase3(frameset, /*NULL, */parlist,frameset, /*NULL,*/
1688 recipe_name, phead2D, "RADECSYS", fname,
1689 hdrl_image_get_image(himg_exptime), CPL_TYPE_DOUBLE,
1690 UNIT_TIME);
1691 hdrl_image_delete(himg_exptime);
1692 cpl_propertylist_delete(phead2D);
1693 cpl_free(fname);
1694 }
1695
1696
1697
1698 cpl_mask_delete(bpm_mask) ;
1699 }
1700
1701 hdrl_imagelist* himlist = hdrl_imagelist_create(iml_data, iml_errs);
1702
1703 hdrl_imagelist_collapse_mean(himlist, &cube_collapsed, &cube_cmap);
1704 hdrl_imagelist_delete(himlist);
1705
1706 /* TODO: check header is correct */
1707 cpl_propertylist* phead2D = cpl_propertylist_duplicate(pheader);
1708 cpl_propertylist_erase_regexp(phead2D,NAXIS3,0);
1709 cpl_propertylist_erase_regexp(phead2D,CTYPE3,0);
1710 cpl_propertylist_erase_regexp(phead2D,CRVAL3,0);
1711 cpl_propertylist_erase_regexp(phead2D,CRPIX3,0);
1712 cpl_propertylist_erase_regexp(phead2D,CUNIT3, 0);
1713 cpl_propertylist_erase_regexp(phead2D,"CD3_*",0);
1714 cpl_propertylist_erase_regexp(phead2D,CD1_3,0);
1715 cpl_propertylist_erase_regexp(phead2D,CD2_3,0);
1716
1717 char* fname = NULL;
1718 char* prefix = NULL;
1719 char* suffix = NULL;
1720 if(apply_flat) {
1721 suffix = cpl_sprintf("%s","");
1722 } else {
1723 suffix = cpl_sprintf("%s","_no_flat");
1724 }
1725 if(strstr(recipe_name, "jitter") != NULL ) {
1726 prefix = cpl_sprintf("eris_ifu_jitter%s", suffix);
1727 } else {
1728 prefix = cpl_sprintf("eris_ifu_stdstar%s", suffix);
1729 }
1730 cpl_free(suffix);
1731 const char* pcatg = NULL;
1732 if(strstr(pro_catg, "OBJ") != NULL ) {
1733 fname = cpl_sprintf("%s_obj_cube_mean.fits", prefix);
1734 if(apply_flat) {
1735 pcatg = ERIS_IFU_PRO_JITTER_OBJ_CUBE_MEAN;
1736 } else {
1737 pcatg = ERIS_IFU_PRO_JITTER_OBJ_CUBE_NOFLAT_MEAN;
1738 }
1739 } else if(strstr(pro_catg, ERIS_IFU_PRO_JITTER_DAR_CUBE) != NULL ) {
1740 fname = cpl_sprintf("%s_dar_cube_mean.fits", prefix);
1741 if(apply_flat) {
1742 pcatg = ERIS_IFU_PRO_JITTER_OBJ_DAR_CUBE_MEAN;
1743 } else {
1744 pcatg = ERIS_IFU_PRO_JITTER_OBJ_CUBE_NOFLAT_MEAN;
1745 }
1746 } else if(strstr(pro_catg, ERIS_IFU_PRO_JITTER_TWK_CUBE) != NULL ) {
1747 fname = cpl_sprintf("%s_twk_cube_mean.fits", prefix);
1748 if(apply_flat) {
1749 pcatg = ERIS_IFU_PRO_JITTER_TWK_CUBE_MEAN;
1750 } else {
1751 pcatg = ERIS_IFU_PRO_JITTER_OBJ_CUBE_NOFLAT_MEAN;
1752 }
1753 } else if(strstr(pro_catg, "STD") != NULL ) {
1754 fname = cpl_sprintf("%s_std_cube_mean.fits", prefix);
1755 if(apply_flat) {
1756 pcatg = ERIS_IFU_PRO_JITTER_STD_CUBE_MEAN;
1757 } else {
1758 pcatg = ERIS_IFU_PRO_JITTER_STD_CUBE_NOFLAT_MEAN;
1759 }
1760 } else if(strstr(pro_catg, "PSF") != NULL ) {
1761 fname = cpl_sprintf("%s_psf_cube_mean.fits", prefix);
1762 if(apply_flat) {
1763 pcatg = ERIS_IFU_PRO_JITTER_PSF_CUBE_MEAN;
1764 } else {
1765 pcatg = ERIS_IFU_PRO_JITTER_PSF_CUBE_NOFLAT_MEAN;
1766 }
1767 }
1768 cpl_free(prefix);
1769 cpl_propertylist_update_string(phead2D, CPL_DFS_PRO_CATG, pcatg);
1770 /* remove WCS from primary header as data go into extensions */
1771// eris_ifu_plist_erase_wcs(phead2D);
1772
1773
1774 if( is_pupil) {
1775 /* compute PUPIL X/Y SHIFT for QC */
1776// hdrl_imagelist* himl = NULL;
1777
1778 cpl_table* qclog_tbl = eris_qclog_init();
1779 cpl_image* img = hdrl_image_get_image(cube_collapsed);
1780 cpl_size max_ima_x = 0;
1781 cpl_size max_ima_y = 0;
1782 cpl_image_get_maxpos(img, &max_ima_x, &max_ima_y);
1783 cpl_size sx = hdrl_image_get_size_x(cube_collapsed);
1784 cpl_size sy = hdrl_image_get_size_y(cube_collapsed);
1785 double max_ima_cx = cpl_image_get_centroid_x_window(img, 1, 1, sx, sy);
1786 double max_ima_cy = cpl_image_get_centroid_y_window(img, 1, 1, sx, sy);
1787 double xshift = max_ima_cx - (double) sx * 0.5;
1788 double yshift = max_ima_cy - (double) sy * 0.5;
1789 char* key_name = cpl_sprintf("QC PUPIL%d SHIFTX", 0);
1790
1791 eris_qclog_add_double_f(qclog_tbl, key_name, xshift,
1792 "[pix] X shift centroid - center image");
1793 cpl_free(key_name);
1794
1795 key_name = cpl_sprintf("QC PUPIL%d SHIFTY", 0);
1796 eris_qclog_add_double_f(qclog_tbl, key_name, yshift,
1797 "[pix] Y shift centroid - center image");
1798 cpl_msg_info(cpl_func, "xshift: %g yshift: %g", xshift, yshift);
1799
1800 cpl_free(key_name);
1801 eris_pfits_put_qc(phead2D, qclog_tbl);
1802 cpl_free(qclog_tbl);
1803 }
1804
1805 cpl_image* mask_ima = cpl_image_new_from_mask(hdrl_image_get_mask(cube_collapsed));
1806 eris_ifu_save_deq_image(frameset, NULL, parlist,frameset, NULL,
1807 recipe_name, phead2D, "RADECSYS", fname,
1808 hdrl_image_get_image(cube_collapsed),
1809 hdrl_image_get_error(cube_collapsed),
1810 rmse, mask_ima, flag16bit, UNIT_ADU);
1811
1812 cpl_image_delete(mask_ima);
1813 cpl_free(fname);
1815 cpl_propertylist_delete(dheader);
1816 cpl_propertylist_delete(eheader);
1817 cpl_propertylist_delete(qheader);
1818 cpl_propertylist_delete(pheader);
1819 cpl_imagelist_delete(iml_data);
1820 cpl_imagelist_delete(iml_errs);
1821 cpl_imagelist_delete(iml_qual);
1822 hdrl_image_delete(cube_collapsed);
1823 cpl_image_delete(cube_cmap);
1824 eris_check_error_code("eris_ifu_cube_collapse_mean_and_save");
1825
1826 return cpl_error_get_code();
1827}
1828
1829static double
1830eris_ifu_compute_max_cubes_dist(cpl_frameset* cube_set)
1831{
1832 /* First check if cube size is too large NOT YET VERIFIED */
1833 double dist_max = -1;
1834 double dist_sqr = -1;
1835 cpl_size ncubes = cpl_frameset_get_size(cube_set);
1836 cpl_vector* vec_ra = cpl_vector_new(ncubes);
1837 cpl_vector* vec_dec = cpl_vector_new(ncubes);
1838 cpl_propertylist* plist = NULL;
1839 const char* name = NULL;
1840 double ra = 0, dec = 0;
1841 double ra1 = 0, dec1 = 0;
1842 double ra2 = 0, dec2 = 0;
1843 double dist = 0;
1844 //double cd1_1 = 0;
1845 //double cd2_2 = 0;
1846 double pix_size = 0;
1847 //const char* scale = NULL;
1848 cpl_propertylist* phead = NULL;
1849 //TODO: why cd3_3 is set and not used?
1850// double cd3_3 = 0;
1851 int extnum_raw = 1;
1852 cpl_frame* data_frame = NULL;
1853
1854 data_frame = cpl_frameset_get_position(cube_set, 0);
1855 name = cpl_frame_get_filename(data_frame);
1856 plist = cpl_propertylist_load(name, extnum_raw);
1857 phead = cpl_propertylist_load(name, 0);
1858// cd1_1 = cpl_propertylist_get_double(plist, "CD1_1");
1859// cd2_2 = cpl_propertylist_get_double(plist, "CD2_2");
1860// cd3_3 = cpl_propertylist_get_double(plist, "CD3_3");
1861
1862 ifsPreopticsScale scale_mas = eris_ifu_get_preopticsScale(phead);
1863
1864 cpl_propertylist_delete(plist);
1865 cpl_propertylist_delete(phead);
1866
1867 switch(scale_mas){
1868 case S250MAS: pix_size = 0.250; break;
1869 case S100MAS: pix_size = 0.100; break;
1870 case S25MAS: pix_size = 0.025; break;
1871 case PUPIL: pix_size = 0.025; break;
1872 case UNDEFINED_SCALE: pix_size = 0.100; break;
1873 }
1874
1875 pix_size /= 3600;
1876
1877 /*
1878 cpl_msg_warning(cpl_func, "cd1_1: %g cd2_2: %g cd3_3: %g",
1879 cd1_1, cd2_2, cd3_3);
1880 */
1881 for (cpl_size i = 0; i < ncubes ; i++) {
1882
1883 data_frame = cpl_frameset_get_position(cube_set, i);
1884 name = cpl_frame_get_filename(data_frame);
1885 plist = cpl_propertylist_load(name, extnum_raw);
1886
1887 ra = cpl_propertylist_get_double(plist, "RA");
1888 dec = cpl_propertylist_get_double(plist, "DEC");
1889 //cpl_msg_warning(cpl_func, "ra: %13.10g dec: %13.10g",ra,dec);
1890 cpl_vector_set(vec_ra, i, ra);
1891 cpl_vector_set(vec_dec, i, dec);
1892
1893 cpl_propertylist_delete(plist);
1894 }
1895 for (cpl_size i = 0; i < ncubes; i++) {
1896 ra1 = cpl_vector_get(vec_ra, i);
1897 dec1 = cpl_vector_get(vec_dec, i);
1898 for (cpl_size j = 0; j < ncubes && i != j; j++) {
1899 ra2 = cpl_vector_get(vec_ra, j);
1900 dec2 = cpl_vector_get(vec_dec, j);
1901 /*
1902 cpl_msg_warning(cpl_func, "ra1: %13.10g ra2: %13.10g",ra1, ra2);
1903 cpl_msg_warning(cpl_func, "dist_ra: %13.10g",(ra2 - ra1)/pix_size);
1904 cpl_msg_warning(cpl_func, "dec1: %13.10g dec2: %13.10g",dec1, dec2);
1905 cpl_msg_warning(cpl_func, "dist_dec: %13.10g",(dec2 - dec1)/pix_size);
1906 */
1907 dist = (ra2 - ra1) * (ra2 - ra1) / pix_size / pix_size +
1908 (dec2 - dec1) * (dec2 - dec1) / pix_size / pix_size;
1909 //cpl_msg_warning(cpl_func, "dist: %13.10g",dist);
1910 if(dist > dist_sqr) {
1911 dist_sqr = dist;
1912 }
1913 }
1914 }
1915 cpl_vector_delete(vec_ra);
1916 cpl_vector_delete(vec_dec);
1917 dist_max = sqrt(dist_sqr);
1918 cpl_msg_info(cpl_func, "Max distance between contributing cubes centers: %g",dist_max);
1919
1920 eris_check_error_code("eris_ifu_compute_max_cubes_dist");
1921 return dist_max;
1922
1923}
1924
1925/*
1926static cpl_error_code
1927eris_resample_compute_size(hdrl_parameter *aParams_outputgrid,
1928 int *aX, int *aY, int *aZ)
1929{
1930 cpl_ensure_code(aParams_outputgrid && aX && aY && aZ, CPL_ERROR_NULL_INPUT);
1931 const char func[] = "hdrl_resample_compute_size";
1932 double x1, y1, x2, y2;
1933
1934 double ramin = aParams_outputgrid->ra_min;
1935 double ramax = aParams_outputgrid->ra_max;
1936 double decmin = aParams_outputgrid->dec_min;
1937 double decmax = aParams_outputgrid->dec_max;
1938 hdrl_resample_wcs_projplane_from_celestial(aParams_outputgrid, ramin, decmin,
1939 &x1, &y1);
1940 hdrl_resample_wcs_projplane_from_celestial(aParams_outputgrid, ramax, decmax,
1941 &x2, &y2);
1942 *aX = lround(fabs(x2 - x1) / aParams_outputgrid->delta_ra) + 1;
1943 *aY = lround(fabs(y2 - y1) / aParams_outputgrid->delta_dec) + 1;
1944
1945 double lmin = aParams_outputgrid->lambda_min;
1946 double lmax = aParams_outputgrid->lambda_max;
1947
1948 *aZ = (int)ceil((lmax - lmin) / aParams_outputgrid->delta_lambda) + 1;
1949
1950 cpl_msg_debug(func, "Output cube size %d x %d x %d (fit to data)",
1951 *aX, *aY, *aZ);
1952 return CPL_ERROR_NONE;
1953}
1954*/
1967cpl_error_code
1968eris_ifu_combine(cubeType obj_type, cpl_frameset* frameset,
1969 const cpl_parameterlist * parlist,
1970 const char* recipe_name,
1971 const char* pipefile_prefix)
1972{
1973 cpl_size nframes = cpl_frameset_get_size(frameset);
1974
1975 cpl_frameset* cube_set = NULL;
1976 char *proCatg = NULL;
1977 char *filenameSpec = NULL;
1978 cpl_msg_info(cpl_func,"Combine cubes into a single one, compute mean, median along wavelength axis");
1979 eris_ifu_jitter_get_procatg_and_filename(obj_type, &proCatg, &filenameSpec);
1980
1981 if ((cube_set = eris_ifu_extract_frameset(frameset, proCatg )) == NULL){
1982
1983 cpl_msg_warning(cpl_func,"No %s files", proCatg);
1984 cpl_error_set_message(cpl_func, CPL_ERROR_FILE_NOT_FOUND, "Missing CUBE "
1985 "files");
1986 return cpl_error_get_code();
1987 } else {
1988 cpl_msg_info(cpl_func, "%02d file(s) of type %s",
1989 (int)cpl_frameset_get_size(cube_set), proCatg);
1990 }
1991
1992 /*convert the images in a table of predefined structure*/
1993 cpl_table* restable = NULL; /*Final table to resample*/
1994 cpl_table* tab_tmp = NULL;
1995 const cpl_parameter * par = NULL;
1996 cpl_frame * data_frame = NULL;
1997 cpl_frame * errs_frame = NULL;
1998 cpl_frame * qual_frame = NULL;
1999 int extnum_raw = 1;
2000 int extnum_err = 2;
2001 int extnum_bpm = 3;
2002 cpl_boolean is_variance = CPL_TRUE;
2003 cpl_boolean subtract_bkg = CPL_FALSE;
2004 int edge_trim = 0;
2005 /*Loops over all frames and adds the images/imagelists to the final table */
2006 cpl_size ncubes = cpl_frameset_get_size(cube_set);
2007 char* pname = cpl_sprintf("eris.%s.max-cubes-centres-dist",recipe_name);
2008 int max_cubes_dist = 20;
2009 //cpl_parameterlist_dump(parlist,stdout);
2010
2011 max_cubes_dist = cpl_parameter_get_int(
2012 cpl_parameterlist_find_const(parlist, pname));
2013
2014 if ( max_cubes_dist < eris_ifu_compute_max_cubes_dist(cube_set)) {
2015 cpl_msg_info(cpl_func,"max-cubes-centres-dist: %d",max_cubes_dist);
2016 cpl_msg_warning(cpl_func,"max distance between cube centres greater than threshold. Skip cube combination.");
2017 return CPL_ERROR_INCOMPATIBLE_INPUT;
2018 }
2019
2020 for (cpl_size i = 0; i < ncubes ; i++) {
2021
2022 data_frame = cpl_frameset_get_position(cube_set, i);
2023 errs_frame = data_frame;
2024 qual_frame = data_frame;
2025
2026 tab_tmp = eris_ifu_resample_get_table_from_frameset(data_frame,
2027 errs_frame,
2028 qual_frame, extnum_raw,
2029 extnum_err, extnum_bpm,
2030 is_variance, subtract_bkg, edge_trim);
2031
2032
2033 if(nframes == 1) {
2034 restable = tab_tmp;
2035 } else {
2036 eris_ifu_resample_update_table(tab_tmp, &restable);
2037 cpl_table_delete(tab_tmp);
2038 }
2039
2040
2041 }
2042 cpl_ensure_code(restable, CPL_ERROR_NULL_INPUT);
2043
2044 /* Save the final table if required */
2045 eris_ifu_resample_tablesave (recipe_name, par, parlist, frameset, restable,
2046 pipefile_prefix);
2047 /* Getting the input wcs needed for determine the scaling of output cube */
2048 cpl_wcs *wcs = eris_ifu_resample_get_wcs_from_frameset(cube_set, proCatg);
2049 cpl_ensure_code(wcs, CPL_ERROR_NULL_INPUT);
2050
2051 cpl_msg_info(cpl_func, "WCS used for the output grid definition and passed to "
2052 "hdrl_resample_compute():");
2053 cpl_msg_indent_more();
2054 eris_ifu_resample_wcs_print(wcs);
2055 cpl_msg_indent_less();
2056
2057 hdrl_parameter *aParams_method = NULL;
2058 /* set re-sampling method */
2059
2060 char *context = NULL;
2061 context = cpl_sprintf("eris.%s", recipe_name);
2062 aParams_method = eris_ifu_resample_set_method(parlist, context);
2063 cpl_free(context);
2064
2065 cpl_ensure_code(aParams_method, CPL_ERROR_NULL_INPUT);
2066
2067 /* set re-sampling outputgrid */
2068 hdrl_parameter *aParams_outputgrid = NULL;
2069
2070 aParams_outputgrid = eris_ifu_resample_set_outputgrid(recipe_name, parlist,
2071 restable, wcs);
2072 if(cpl_error_get_code() != CPL_ERROR_NONE) {
2073
2074 eris_check_error_code("eris_ifu_combine");
2075 return CPL_ERROR_INCOMPATIBLE_INPUT;
2076 }
2077
2078 cpl_ensure_code(aParams_outputgrid, CPL_ERROR_NULL_INPUT);
2079
2080 hdrl_resample_result *result = NULL;
2081
2082 /* !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! */
2083 /* !!! Calling the hdrl_resample_compute function doing the real work !!! */
2084 /* !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! */
2085 cpl_msg_severity level = cpl_msg_get_level();
2086 cpl_msg_set_level(CPL_MSG_INFO);
2087 result = hdrl_resample_compute(restable, aParams_method, aParams_outputgrid,
2088 wcs);
2089 cpl_msg_set_level(level);
2090 /* TO BE REMOVED: this was added to check position of object (verify DAR)
2091 hdrl_image* hima = NULL;
2092 double max_x = 0, max_y = 0, peak = 0;
2093 cpl_size margin = 100;
2094 for(cpl_size k= margin; k < hdrl_imagelist_get_size(result->himlist) - margin; k++) {
2095 hima = hdrl_imagelist_get(result->himlist, k);
2096 eris_gaussian_maxpos(hdrl_image_get_image(hima), 5, &max_x, &max_y, &peak);
2097 cpl_msg_info(cpl_func,"z=%lld x: %g y: %g",k,max_x,max_y);
2098 }
2099 */
2100
2101 cpl_ensure_code(result, CPL_ERROR_NULL_INPUT);
2102 char* pro_catg = NULL;
2103 if(strstr(pipefile_prefix,"no_flat") != NULL) {
2104 pro_catg = cpl_sprintf("%s_%s",proCatg,"COADD_NOFLAT");
2105 } else {
2106 pro_catg = cpl_sprintf("%s_%s",proCatg,"COADD");
2107 }
2108 char* fname = cpl_sprintf("%s_%s_coadd.fits", pipefile_prefix, filenameSpec);
2109
2110 eris_ifu_resample_save_cube(result, pro_catg, recipe_name, fname, parlist,
2111 frameset, CPL_FALSE);
2112
2113 eris_ifu_free_string(&fname);
2114 eris_ifu_free_string(&pro_catg);
2115 hdrl_image* cube_collapsed = NULL;
2116 cpl_image* cube_cmap = NULL;
2117 hdrl_imagelist_collapse_mean(result->himlist, &cube_collapsed, &cube_cmap);
2118
2119 /* TODO: check header is correct */
2120 pro_catg = cpl_sprintf("%s_%s",proCatg,"MEAN");
2121 fname = cpl_sprintf("%s_%s_mean.fits", pipefile_prefix, filenameSpec);
2122 cpl_propertylist* phead2D = cpl_propertylist_duplicate(result->header);
2123 /* remove WCS from primary header as data go into extensions */
2124// eris_ifu_plist_erase_wcs(phead2D);
2125 cpl_propertylist_append_string(phead2D, "PRODCATG",PRODCATG_WHITELIGHT);
2126
2127 eris_ifu_save_image(frameset, phead2D, parlist, recipe_name, pro_catg,
2128 fname, CPL_TYPE_DOUBLE, hdrl_image_get_image(cube_collapsed));
2129 cpl_image_save(hdrl_image_get_error(cube_collapsed), fname, CPL_TYPE_DOUBLE,
2130 NULL, CPL_IO_EXTEND);
2131 cpl_image_save(cpl_image_new_from_mask(hdrl_image_get_mask(cube_collapsed)),
2132 fname, CPL_TYPE_INT, NULL, CPL_IO_EXTEND);
2133
2134 eris_ifu_free_string(&fname);
2135 eris_ifu_free_string(&pro_catg);
2136 eris_ifu_free_string(&filenameSpec);
2137 eris_ifu_free_string(&proCatg);
2139 hdrl_image_delete(cube_collapsed);
2140 cpl_image_delete(cube_cmap);
2141
2142
2143 /* Printing the wcs after resampling: */
2144
2145 cpl_wcs_delete(wcs);
2146 wcs = cpl_wcs_new_from_propertylist(result->header);
2147 cpl_msg_info(cpl_func, "Final WCS after resampling: ");
2148 cpl_msg_indent_more();
2149 eris_ifu_resample_wcs_print(wcs);
2150 cpl_msg_indent_less();
2151
2152 /* Cleanup */
2153
2154 /* Delete the parameters */
2155 hdrl_parameter_delete(aParams_method);
2156 hdrl_parameter_delete(aParams_outputgrid);
2157
2159 cpl_table_delete(restable);
2160 cpl_wcs_delete(wcs);
2161 cpl_frameset_delete(cube_set);
2162 //cpl_frameset_delete(var_err_set);
2163 //cpl_frameset_delete(bpm_set);
2164 eris_check_error_code("eris_ifu_combine");
2165
2166 return cpl_error_get_code();
2167}
2168
2176static cpl_propertylist* eris_ifu_resample_get_header2D(cpl_propertylist *header3D){
2177 cpl_propertylist *phead2D = cpl_propertylist_duplicate(header3D);
2178 cpl_propertylist_set_int(phead2D,"NAXIS",2);
2179 cpl_propertylist_erase_regexp(phead2D,"NAXIS3",0);
2180 cpl_propertylist_erase_regexp(phead2D,"CTYPE3",0);
2181 cpl_propertylist_erase_regexp(phead2D,"CRVAL3",0);
2182 cpl_propertylist_erase_regexp(phead2D,"CRPIX3",0);
2183 cpl_propertylist_erase_regexp(phead2D,"CUNIT3",0);
2184 cpl_propertylist_erase_regexp(phead2D,"CDELT3",0);
2185 cpl_propertylist_erase_regexp(phead2D,"CD3_*",0);
2186 cpl_propertylist_erase_regexp(phead2D,"CD1_3",0);
2187 cpl_propertylist_erase_regexp(phead2D,"CD2_3",0);
2188 eris_check_error_code("eris_ifu_resample_get_header2D");
2189 return phead2D;
2190}
2191
2200cpl_error_code eris_ifu_resample_trim_edge(hdrl_image *himg, int edge_trim){
2201 if (edge_trim > 0) {
2202 //cpl_msg_info(cpl_func, "Trim input image edges of %d pixels", edge_trim);
2203 /* trim each product plane edge of edgetrim pixels */
2204 cpl_size sx = hdrl_image_get_size_x(himg);
2205 cpl_size sy = hdrl_image_get_size_y(himg);
2206 if (edge_trim >= 0.5* sx) {
2207 edge_trim = 0;
2208 cpl_msg_warning(cpl_func, "edge-trim must be smaller than half image "
2209 "size. Reset to 0");
2210 }
2211 if (edge_trim >= 0.5* sy) {
2212 edge_trim = 0;
2213 cpl_msg_warning(cpl_func, "edge-trim must be smaller than half image "
2214 "size. Reset to 0");
2215 }
2216
2217 /* Note FITS convention to loop over pixels */
2218 /* trim lower image border along X direction */
2219 for(cpl_size j = 1; j <= edge_trim; j++) {
2220 for(cpl_size i = 1; i <= sx; i++) {
2221 hdrl_image_reject(himg, i, j);
2222 }
2223 }
2224 /* trim upper image border along X direction */
2225 for(cpl_size j = sy - edge_trim + 1; j <= sy; j++) {
2226 for(cpl_size i = 1; i <= sx; i++) {
2227 hdrl_image_reject(himg, i, j);
2228 }
2229 }
2230 /* trim left image border along Y direction */
2231 for(cpl_size j = 1; j <= sy; j++) {
2232 for(cpl_size i = 1; i <= edge_trim; i++) {
2233 hdrl_image_reject(himg, i, j);
2234 }
2235 }
2236 /* trim right image border along Y direction */
2237 for(cpl_size j = 1; j <= sy; j++) {
2238 for(cpl_size i = sx - edge_trim + 1; i <= sx; i++) {
2239 hdrl_image_reject(himg, i, j);
2240 }
2241 }
2242 }
2243 else {
2244 cpl_msg_warning(cpl_func, "edge-trim is set to 0");
2245 }
2246 eris_check_error_code("eris_ifu_resample_trim_edge");
2247 return cpl_error_get_code();
2248}
2249
2264cpl_error_code
2265eris_ifu_combine_pbp(cpl_frameset* frameset,
2266 const cpl_parameterlist * parlist,
2267 const char *input_cube_pro_catg,
2268 const char *filenameSpec,
2269 float *offsetx,
2270 float *offsety,
2271 const char* offunit,
2272 const char* recipe_name,
2273 const char* pipefile_prefix){
2274
2275 cpl_frameset* cube_set = NULL;
2276 cpl_msg_info(cpl_func,"Combine cubes into a single one plane-by-plane, compute mean, median along wavelength axis");
2277
2278 if ((cube_set = eris_ifu_extract_frameset(frameset, input_cube_pro_catg )) == NULL){
2279
2280 cpl_msg_warning(cpl_func,"No %s files", input_cube_pro_catg);
2281 //cpl_frameset_dump(frameset, stdout);
2282 cpl_error_set_message(cpl_func, CPL_ERROR_FILE_NOT_FOUND,
2283 "Missing cube files tagged as %s", input_cube_pro_catg);
2284 return cpl_error_get_code();
2285 } else {
2286 cpl_msg_info(cpl_func, "%02d file(s) of type %s",
2287 (int)cpl_frameset_get_size(cube_set), input_cube_pro_catg);
2288 }
2289
2290 /*convert the images in a table of predefined structure*/
2291 //const cpl_parameter * par = NULL;
2292 const cpl_frame * data_frame = NULL;
2293 hdrl_parameter *aParams_outputgrid = NULL;
2294 cpl_propertylist *phead2D = NULL;
2295 int extnum_raw = 1;
2296 int extnum_err = 2;
2297 int extnum_bpm = 3;
2298
2299 //cpl_boolean is_variance = CPL_TRUE;
2300 cpl_boolean subtract_bkg = CPL_FALSE;
2301 int edge_trim = 0;
2302 /*Loops over all frames and adds the images/imagelists to the final table */
2303 cpl_size ncubes = cpl_frameset_get_size(cube_set);
2304 if (ncubes == 1){
2305 cpl_error_set_message(cpl_func, CPL_ERROR_FILE_NOT_FOUND, "Only one CUBE "
2306 "files");
2307 return cpl_error_get_code();
2308 }
2309
2310 char* pname = cpl_sprintf("eris.%s.max-cubes-centres-dist",recipe_name);
2311 int max_cubes_dist = 20;
2312 //cpl_parameterlist_dump(parlist,stdout);
2313 if (offsetx !=NULL && offsety !=NULL) // user offset mode
2314 max_cubes_dist = INT_MAX ;//set to the largest, i.e. ignore this param.
2315 else
2316 max_cubes_dist = cpl_parameter_get_int(
2317 cpl_parameterlist_find_const(parlist, pname));
2318 cpl_free(pname);
2319
2320 pname = cpl_sprintf("eris.%s.subtract-background",recipe_name);
2321 subtract_bkg = cpl_parameter_get_bool(
2322 cpl_parameterlist_find_const(parlist, pname));
2323 if(subtract_bkg == CPL_TRUE)
2324 cpl_msg_info(cpl_func, "Median of each plane will be subtracted.");
2325 cpl_free(pname);
2326
2327 pname = cpl_sprintf("eris.%s.edge-trim",recipe_name);
2328 edge_trim = cpl_parameter_get_int(
2329 cpl_parameterlist_find_const(parlist, pname));
2330 if (edge_trim > 0)
2331 cpl_msg_info(cpl_func, "%d pixel(s) on edges will be trimmed.", edge_trim);
2332
2333 cpl_free(pname);
2334
2335
2336
2337 if ( max_cubes_dist < eris_ifu_compute_max_cubes_dist(cube_set)) {
2338 cpl_msg_info(cpl_func,"max-cubes-centres-dist: %d",max_cubes_dist);
2339 cpl_msg_warning(cpl_func,"max distance between cube centres greater than threshold. Skip cube combination.");
2340 return CPL_ERROR_INCOMPATIBLE_INPUT;
2341 }
2342
2343 hdrl_parameter *aParams_method = NULL;
2344 /* set re-sampling method */
2345 char *context = NULL;
2346 context = cpl_sprintf("eris.%s", recipe_name);
2347 aParams_method = eris_ifu_resample_set_method(parlist, context);
2348 cpl_free(context);
2349 cpl_ensure_code(aParams_method, CPL_ERROR_NULL_INPUT);
2350
2351 hdrl_imagelist *res_cube = hdrl_imagelist_new();
2352 cpl_wcs *wcs = NULL;
2353 cpl_wcs **tmp_wcs_array = (cpl_wcs **) cpl_calloc(ncubes, sizeof(wcs));
2354 char const **name_array = (char const **) cpl_calloc(ncubes, sizeof(const char*));
2355
2356 cpl_wcs *wcs3d = NULL;
2357 double mjd_obs = 0;
2358 double crpix3 = 0;
2359 double crval3 = 0;
2360 // double cdelt3 = 0;
2361 double cd1_3 = 0;
2362 double cd2_3 = 0;
2363 double cd3_1 = 0;
2364 double cd3_2 = 0;
2365 double cd3_3 = 0;
2366 double cd1_1 = 0;
2367 double cd1_2 = 0;
2368 double cd2_1 = 0;
2369 double cd2_2 = 0;
2370 int naxis3 = 0;
2371 const char* ctype3 = "WAVE";
2372 // const char* cunit3 = "um";
2373 const cpl_matrix *cd = NULL;
2374 cpl_vector* exptime_vec = cpl_vector_new(ncubes);
2375 double exptime = 0;
2376 double dit = 0;
2377 int ndit = 0;
2378 for (cpl_size i = 0; i < ncubes ; i++) {
2379 data_frame = cpl_frameset_get_position_const(cube_set, i);
2380 const char *name = cpl_frame_get_filename(data_frame);
2381 cpl_propertylist *xheader_data = cpl_propertylist_load(name,
2382 extnum_raw);
2383 cpl_propertylist* phead = cpl_propertylist_load(name, 0);
2384 dit = eris_pfits_get_dit(phead);
2385 ndit = eris_pfits_get_ndit(phead);
2386
2387 cpl_propertylist_delete(phead);
2388 exptime = dit * ndit;
2389 cpl_vector_set(exptime_vec, i, exptime);
2390 cpl_propertylist *header2D = eris_ifu_resample_get_header2D(xheader_data);
2391 cpl_wcs *tmp_wcs = cpl_wcs_new_from_propertylist(header2D);
2392 cpl_propertylist_delete(header2D);
2393
2394 tmp_wcs_array[i] = tmp_wcs;
2395 name_array[i] = name;
2396
2397 if (i==0){
2398 wcs = tmp_wcs;
2399 cpl_ensure_code(wcs, CPL_ERROR_NULL_INPUT);
2400 cpl_msg_info(cpl_func, "WCS used for the output grid definition and passed to "
2401 "hdrl_resample_compute():");
2402 cpl_msg_indent_more();
2403 eris_ifu_resample_wcs_print(wcs);
2404 cpl_msg_indent_less();
2405 wcs3d = cpl_wcs_new_from_propertylist(xheader_data);
2406 mjd_obs = cpl_propertylist_get_double(xheader_data, "MJD-OBS");
2407 crpix3 = cpl_propertylist_get_double(xheader_data, "CRPIX3");
2408 crval3 = cpl_propertylist_get_double(xheader_data, "CRVAL3");
2409 naxis3 = cpl_propertylist_get_int(xheader_data, "NAXIS3");
2410 //cdelt3 = cpl_propertylist_get_double(xheader_data, "CDELT3");
2411 /*
2412 if(cpl_propertylist_has(xheader_data,"CTYPE3")) {
2413 ctype3 = cpl_propertylist_get_string(xheader_data, "CTYPE3");
2414 cpl_msg_warning(cpl_func,"check ctype3: >%s< >%s<",
2415 ctype3,cpl_propertylist_get_string(xheader_data, "CTYPE3"));
2416 } else {
2417 ctype3 = "WAVE";
2418 }
2419 if(cpl_propertylist_has(xheader_data,"CUNIT3")) {
2420 cunit3 = cpl_propertylist_get_string(xheader_data, "CUNIT3");
2421 } else {
2422 cunit3 = "um";
2423 }
2424 */
2425 cd3_3 = cpl_propertylist_get_double(xheader_data, "CD3_3");
2426 // cdelt3 = cd3_3;
2427
2428 cd = cpl_wcs_get_cd(wcs);
2429 cd1_1 = cpl_matrix_get(cd, 0, 0);
2430 cd1_2 = cpl_matrix_get(cd, 1, 0);
2431 cd2_1 = cpl_matrix_get(cd, 0, 1);
2432 cd2_2 = cpl_matrix_get(cd, 1, 1);
2433
2434 cpl_ensure_code(wcs3d, CPL_ERROR_NULL_INPUT);
2435 }
2436 cpl_propertylist_delete(xheader_data);
2437 }
2438
2439 /*Convert offset unit to degree*/
2440 double as2deg = 1.0/3600.0;
2441 double scalex = 1.0;
2442 double scaley = 1.0;
2443
2444 int offcode = -1;
2445 if (offunit != NULL && offsetx !=NULL && offsety !=NULL){
2446 if (strcmp(offunit, "PIXEL") == 0){
2447 offcode = 0;
2448 cpl_msg_info(cpl_func, "User-defined offset unit in PIXEL");
2449 }else if (strcmp(offunit, "DEGREE") == 0){
2450 offcode = 1;
2451 cpl_msg_info(cpl_func, "User-defined offset unit is DEGREE");
2452 }else if (strcmp(offunit, "ARCSEC") == 0){
2453 offcode = 2;
2454 cpl_msg_info(cpl_func, "User-defined offset unit is ARCSEC");
2455 scalex = as2deg;
2456 scaley = as2deg;
2457 }else
2458 cpl_msg_error(cpl_func, "User-defined offset unit or offset list is not correct.");
2459 }
2460
2461
2462 const cpl_array *dims = cpl_wcs_get_image_dims(wcs3d);
2463 int nchan = cpl_array_get_int(dims, 2, NULL);
2464 cpl_size ich_center = 0.5 * nchan;
2465 cpl_msg_info(cpl_func, "Number of planes: %4.4d", nchan);
2466 //cpl_table* restable_exptime = NULL;
2467// hdrl_image *res_himg_exptime = NULL;
2468 hdrl_imagelist* res_himlist_exptime = hdrl_imagelist_new();
2469
2470 cpl_msg_severity level = cpl_msg_get_level();
2471 cpl_msg_set_level(CPL_MSG_INFO);
2472 for (cpl_size ich = 0; ich < nchan; ich++){
2473 cpl_table* restable = NULL; /*Final table to resample*/
2474
2475 for (cpl_size i = 0; i < ncubes ; i++) {
2476 const char *name = name_array[i];
2477 cpl_image *im = cpl_image_load(name, CPL_TYPE_DOUBLE, ich, extnum_raw);
2478 cpl_image *err = cpl_image_load(name, CPL_TYPE_DOUBLE, ich, extnum_err);
2479 cpl_image* qual = cpl_image_load(name, CPL_TYPE_INT, ich, extnum_bpm);
2480 cpl_mask *mask = cpl_mask_threshold_image_create(qual, 0, INT_MAX);
2481 cpl_image_delete(qual);
2482 // EXPOSURE MAP
2483 /* exposure map: use the plane at the center of the wavelength range
2484 * then replace intensity with EXPTIME value and error with 1 and
2485 * generate HDRL image. This will be re-sampled as the cube and give
2486 * the exposure map associated to the combined cube
2487 */
2488
2489 cpl_image* ima_time = NULL;
2490 cpl_image* err_time = NULL;
2491 hdrl_image *himg_exptime = NULL;
2492 if(ich == ich_center) {
2493 ima_time = cpl_image_new(cpl_image_get_size_x(im),
2494 cpl_image_get_size_y(im), CPL_TYPE_DOUBLE);
2495 err_time = cpl_image_new(cpl_image_get_size_x(im),
2496 cpl_image_get_size_y(im), CPL_TYPE_DOUBLE);
2497
2498 exptime = cpl_vector_get(exptime_vec, i);
2499 cpl_image_add_scalar(ima_time, exptime);
2500 cpl_image_add_scalar(err_time, sqrt(exptime));
2501 himg_exptime = hdrl_image_create(ima_time, err_time);
2502 /* for EXPPOSURE MAP we assume all pixels are good */
2503 //hdrl_image_reject_from_mask(himg_exptime, mask);
2504
2505 if (edge_trim > 0) {
2506 eris_ifu_resample_trim_edge(himg_exptime, edge_trim);
2507 }
2508
2509 cpl_image_delete(ima_time);
2510 cpl_image_delete(err_time);
2511 }
2512
2513 hdrl_image *himg = hdrl_image_create(im, err);
2514 cpl_image_delete(im);
2515 cpl_image_delete(err);
2516 hdrl_image_reject_from_mask(himg, mask);
2517 cpl_mask_delete(mask);
2518
2519 /* Subtract the background on request */
2520 if(subtract_bkg == CPL_TRUE) {
2521 //cpl_msg_info(cpl_func, "Subtracting the median as requested ...");
2523 }
2524
2525 if (edge_trim > 0) {
2526 eris_ifu_resample_trim_edge(himg, edge_trim);
2527 }
2528
2529 cpl_table *tab_tmp = NULL;
2530 cpl_table *tab_exptime_tmp = NULL;
2531 if (offsetx !=NULL && offsety !=NULL){
2532 if (offcode==0){
2533 tab_tmp = hdrl_resample_image_to_table(himg, tmp_wcs_array[0]);
2534 cpl_table_subtract_scalar(tab_tmp, "ra", offsetx[i]*cd1_1+ offsety[i]*cd1_2);
2535 cpl_table_subtract_scalar(tab_tmp, "dec", offsetx[i]*cd2_1+ offsety[i]*cd2_2);
2536 if (ich ==0)
2537 cpl_msg_info(__func__, "User-defined offest in RA, DEC (ARCSEC) %g, %g",
2538 (offsetx[i]*cd1_1+ offsety[i]*cd1_2)*3600,
2539 (offsetx[i]*cd2_1+ offsety[i]*cd2_2)*3600);
2540
2541 if(ich == ich_center) {
2542 tab_exptime_tmp = hdrl_resample_image_to_table(himg_exptime,
2543 tmp_wcs_array[0]);
2544 cpl_table_subtract_scalar(tab_exptime_tmp, "ra", offsetx[i]*cd1_1+ offsety[i]*cd1_2);
2545 cpl_table_subtract_scalar(tab_exptime_tmp, "dec", offsetx[i]*cd2_1+ offsety[i]*cd2_2);
2546 }
2547 } else if (offcode==1 || offcode==2 ){
2548 tab_tmp = hdrl_resample_image_to_table(himg, tmp_wcs_array[0]);
2549 cpl_table_subtract_scalar(tab_tmp, "ra", offsetx[i]*scalex);
2550 cpl_table_subtract_scalar(tab_tmp, "dec", offsety[i]*scaley);
2551 if (ich ==0)
2552 cpl_msg_info(__func__, "User-defined offest in RA, DEC (%s) %g, %g",
2553 offunit,
2554 offsetx[i],
2555 offsety[i]);
2556
2557 if(ich == ich_center) {
2558 tab_exptime_tmp = hdrl_resample_image_to_table(himg_exptime, tmp_wcs_array[0]);
2559 cpl_table_subtract_scalar(tab_exptime_tmp, "ra", offsetx[i]*scalex);
2560 cpl_table_subtract_scalar(tab_exptime_tmp, "dec", offsety[i]*scaley);
2561 }
2562 }
2563 } else {
2564 tab_tmp = hdrl_resample_image_to_table(himg, tmp_wcs_array[i]);
2565 if(ich == ich_center) {
2566 tab_exptime_tmp = hdrl_resample_image_to_table(himg_exptime, tmp_wcs_array[i]);
2567 }
2568 }
2569 hdrl_image_delete(himg);
2570 eris_ifu_resample_update_table(tab_tmp, &restable);
2571
2572 if(ich == ich_center) {
2573
2574 hdrl_image_delete(himg_exptime);
2575 //eris_ifu_resample_update_table(tab_exptime_tmp, &restable_exptime);
2576 hdrl_resample_result *himl_exptime_tmp = NULL;
2577 himl_exptime_tmp = hdrl_resample_compute(tab_exptime_tmp,
2578 aParams_method, aParams_outputgrid, wcs);
2579 //const char* fname = cpl_sprintf("ima_time_%lld.fits",i);
2580 hdrl_image *himg_tmp = hdrl_image_duplicate(
2581 hdrl_imagelist_get_const(himl_exptime_tmp->himlist, 0));
2582 hdrl_imagelist_set(res_himlist_exptime, himg_tmp, i);
2583 hdrl_resample_result_delete(himl_exptime_tmp);
2584
2585 /*
2586 cpl_image_save(hdrl_image_get_image(himg_tmp), fname,
2587 CPL_TYPE_DOUBLE, NULL, CPL_IO_DEFAULT);
2588 */
2589 cpl_table_delete(tab_exptime_tmp);
2590
2591 }
2592
2593 cpl_table_delete(tab_tmp);
2594 if (ich%100 == 0)
2595 cpl_msg_debug(__func__, "Created table for channel %d cube %d", (int)ich, (int)i);
2596 } /* end loop over input data cubes */
2597 cpl_ensure_code(restable, CPL_ERROR_NULL_INPUT);
2598
2599 if (ich%100 == 0)
2600 cpl_msg_debug(__func__, "Table of the chanel %d size %lld x %lld", (int)ich,
2601 cpl_table_get_nrow(restable), cpl_table_get_ncol(restable));
2602
2603 if(ich==0){
2604 aParams_outputgrid = eris_ifu_resample_set_outputgrid2D(recipe_name, parlist,
2605 restable, wcs);
2606 if ( cpl_error_get_code() != CPL_ERROR_NONE) {
2607 eris_check_error_code("eris_ifu_combine_pbp");
2608 return CPL_ERROR_NONE;
2609 }
2610 cpl_ensure_code(aParams_outputgrid, CPL_ERROR_INCOMPATIBLE_INPUT);
2611
2612 //int err_code = cpl_table_save(restable, NULL, NULL, "test_table.fits", CPL_IO_CREATE);
2613 }
2614 cpl_ensure_code(aParams_outputgrid, CPL_ERROR_NULL_INPUT);
2615
2616 hdrl_resample_result *result = NULL;
2617 result = hdrl_resample_compute(restable, aParams_method, aParams_outputgrid, wcs);
2618
2619 /* OLD WAY
2620 if(ich == ich_center) {
2621 hdrl_resample_result *result_exptime = NULL;
2622 result_exptime = hdrl_resample_compute(restable_exptime, aParams_method, aParams_outputgrid, wcs);
2623 res_himg_exptime = hdrl_image_duplicate(hdrl_imagelist_get_const(result_exptime->himlist, 0));
2624
2625 }
2626 */
2627
2628 cpl_ensure_code(result, CPL_ERROR_NULL_INPUT);
2629 cpl_table_delete(restable);
2630
2631 if(ich==0){
2632 /* Printing the wcs after resampling: */
2633 cpl_wcs *res_wcs = cpl_wcs_new_from_propertylist(result->header);
2634 cpl_msg_info(cpl_func, "Final WCS after resampling: ");
2635 cpl_msg_indent_more();
2636 eris_ifu_resample_wcs_print(res_wcs);
2637 cpl_msg_indent_less();
2638 cpl_wcs_delete(res_wcs);
2639 phead2D = cpl_propertylist_duplicate(result->header);
2640 }
2641
2642 hdrl_image *res_himg = hdrl_image_duplicate(hdrl_imagelist_get_const(result->himlist, 0));
2643 hdrl_imagelist_set(res_cube, res_himg, ich);
2645
2646 if (ich==nchan-1){
2647 for (int i = 0; i < ncubes; i++){
2648 cpl_wcs_delete(tmp_wcs_array[i]);
2649 }
2650 cpl_free(tmp_wcs_array);
2651 cpl_free(name_array);
2652 }
2653
2654 if (ich%100 == 0)
2655 cpl_msg_info(__func__, "Coadding cubes: range [%4.4d,%4.4d] of %d",
2656 (int)ich,(int)ich+100, nchan);
2657 }
2658 cpl_msg_set_level(level);
2659
2660 hdrl_image* hima_exptime_mean = NULL;
2661 cpl_image* contrib_exptime_map = NULL;
2662 hdrl_imagelist_collapse_mean(res_himlist_exptime, &hima_exptime_mean,
2663 &contrib_exptime_map);
2664 hdrl_imagelist_delete(res_himlist_exptime);
2665 hdrl_image_delete(hima_exptime_mean);
2666 /*
2667 cpl_image_save(hdrl_image_get_image(hima_exptime_mean), "exptime_mean.fits",
2668 CPL_TYPE_DOUBLE, NULL, CPL_IO_DEFAULT);
2669 */
2670 cpl_image_multiply_scalar(contrib_exptime_map, cpl_vector_get_median(exptime_vec));
2671 //cpl_image_save(contrib, "contrib.fits", CPL_TYPE_DOUBLE, NULL, CPL_IO_DEFAULT);
2672 cpl_ensure_code(res_cube, CPL_ERROR_NULL_INPUT);
2673 cpl_vector_delete(exptime_vec);
2674
2675 hdrl_parameter_delete(aParams_method);
2676 hdrl_parameter_delete(aParams_outputgrid);
2677
2678 /* Save resample cube*/
2679
2680
2681 hdrl_resample_result *result = cpl_calloc(1, sizeof(hdrl_resample_result));
2682 result->himlist = res_cube;
2683 result->header = phead2D;
2684
2685
2686 cpl_propertylist_append_double(result->header, "CRPIX3", crpix3);
2687 cpl_propertylist_append_double(result->header, "CRVAL3", crval3);
2688 cpl_propertylist_append_int(result->header, "NAXIS3", naxis3);
2689 //cpl_propertylist_append_double(result->header, "CDELT3", cdelt3);
2690 cpl_propertylist_append_string(result->header, "CTYPE3", ctype3);
2691 cpl_propertylist_append_double(result->header, "CD1_3", cd1_3);
2692 cpl_propertylist_append_double(result->header, "CD2_3", cd2_3);
2693 cpl_propertylist_append_double(result->header, "CD3_3", cd3_3);
2694 cpl_propertylist_append_double(result->header, "CD3_1", cd3_1);
2695 cpl_propertylist_append_double(result->header, "CD3_2", cd3_2);
2696 cpl_propertylist_append_double(result->header, "MJD-OBS", mjd_obs);
2697 cpl_propertylist* wcs2D = eris_ifu_plist_extract_wcs2D(result->header);
2698
2699 hdrl_image* cube_collapsed = NULL;
2700 cpl_image* cube_cmap = NULL;
2701 hdrl_imagelist_collapse_mean(res_cube, &cube_collapsed, &cube_cmap);
2702 /* This makes 2D products not FITS standard
2703 cpl_propertylist* wcs2D = eris_ifu_plist_extract_wcs(result->header);
2704 cpl_propertylist_append(phead2D, wcs2D);
2705 cpl_propertylist_delete(wcs2D);
2706 */
2707
2708 char* pro_catg = NULL;
2709 if(strstr(pipefile_prefix,"no_flat") != NULL) {
2710 pro_catg = cpl_sprintf("%s_%s",input_cube_pro_catg,"COADD_NOFLAT");
2711 } else {
2712 pro_catg = cpl_sprintf("%s_%s",input_cube_pro_catg,"COADD");
2713 }
2714 char* fname = cpl_sprintf("%s_%s_coadd.fits", pipefile_prefix, filenameSpec);
2715
2716 if(strcmp(recipe_name,"eris_ifu_jitter") == 0) {
2717 eris_ifu_resample_save_cube(result, pro_catg, recipe_name, fname, parlist,
2718 frameset, CPL_FALSE);
2719 } else {
2720 eris_ifu_resample_save_cube(result, pro_catg, recipe_name, fname, parlist,
2721 frameset, CPL_FALSE);
2722 }
2723 eris_ifu_free_string(&fname);
2724 eris_ifu_free_string(&pro_catg);
2725
2726
2727 /* TODO: check header is correct */
2728 pro_catg = cpl_sprintf("%s_%s",input_cube_pro_catg,"MEAN");
2729 fname = cpl_sprintf("%s_%s_mean.fits", pipefile_prefix, filenameSpec);
2730 cpl_propertylist_append_string(phead2D, "PRODCATG", PRODCATG_WHITELIGHT);
2731 /* remove WCS from primary header as data go into extensions */
2732// eris_pfits_erase_wcs2D(phead2D);
2733 cpl_propertylist_update_string(phead2D, CPL_DFS_PRO_CATG, pro_catg);
2734 cpl_image* mask_ima = cpl_image_new_from_mask(hdrl_image_get_mask(cube_collapsed));
2735 eris_ifu_save_deq_image(frameset, NULL, parlist,frameset, NULL,
2736 recipe_name, phead2D, "RADECSYS", fname,
2737 hdrl_image_get_image(cube_collapsed),
2738 hdrl_image_get_error(cube_collapsed),
2739 rmse, mask_ima, flag16bit, UNIT_ADU);
2740 cpl_free(fname);
2741 fname = cpl_sprintf("%s_%s_exposure_map.fits", pipefile_prefix, filenameSpec);
2742
2743 cpl_propertylist_update_string(phead2D, CPL_DFS_PRO_CATG, ERIS_IFU_PRO_JITTER_EXPOSURE_MAP);
2744 cpl_propertylist_update_string(phead2D, "PRODCATG", PRODCATG_EXPOSUREMAP);
2745
2746 cpl_propertylist_append(phead2D, wcs2D);
2747 cpl_propertylist_delete(wcs2D);
2748 eris_ifu_save_image_phase3(frameset, /*NULL, */parlist,frameset, /*NULL,*/
2749 recipe_name, phead2D, "RADECSYS", fname, contrib_exptime_map,
2750 CPL_TYPE_DOUBLE, UNIT_TIME);
2751
2752
2753 cpl_image_delete(mask_ima);
2754 cpl_image_delete(contrib_exptime_map);
2755 eris_ifu_free_string(&fname);
2756 eris_ifu_free_string(&pro_catg);
2757 hdrl_image_delete(cube_collapsed);
2758 cpl_image_delete(cube_cmap);
2759 cpl_wcs_delete(wcs3d);
2761
2762 cpl_frameset_delete(cube_set);
2763 eris_check_error_code("eris_ifu_combine_pbp");
2764
2765 return cpl_error_get_code();
2766}
2767
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 common quality info properties on errors
Definition: eris_ifu_dfs.c:265
cpl_frameset * eris_ifu_extract_frameset(const cpl_frameset *in, const char *tag)
Extract the frames with the given tag from a frameset.
Definition: eris_ifu_dfs.c:171
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)
Write DFS pipeline product with image, error and data qual.
Definition: eris_ifu_dfs.c:586
cpl_error_code eris_ifu_heades_add_hduclass_common(cpl_propertylist *plist, const char *deq_hduclas)
add common properties
Definition: eris_ifu_dfs.c:211
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)
Write DFS pipeline product with image, error and data qual.
Definition: eris_ifu_dfs.c:300
cpl_error_code eris_ifu_heades_add_hduclass_data(cpl_propertylist *plist)
add common data info properties
Definition: eris_ifu_dfs.c:229
cpl_error_code eris_ifu_heades_add_hduclass_errs(cpl_propertylist *plist, deqErrorType errorType)
add common errors info properties
Definition: eris_ifu_dfs.c:246
cpl_error_code eris_ifu_combine(cubeType obj_type, cpl_frameset *frameset, const cpl_parameterlist *parlist, const char *recipe_name, const char *pipefile_prefix)
This function resample input cubes using HDRL.
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)
This function resample input cubes plane by plane using HDRL.
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 an cube (imagelist, product of pixel resampling) on a FITS file with proper FITS header
cpl_error_code eris_ifu_resample_trim_edge(hdrl_image *himg, int edge_trim)
Duplicate a 2D header from a 3D header.
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
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().