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>
65#define MAX_COMBINE_CUBE_PLANE_DIAGONAL 1000
84eris_ifu_resample_update_table(
const cpl_table* append, cpl_table** storage) {
86 cpl_ensure_code(append, CPL_ERROR_NULL_INPUT);
88 if(*storage == NULL) {
89 *storage = cpl_table_new(0);
90 cpl_table_copy_structure(*storage, append);
91 cpl_table_insert(*storage, append, 0);
93 cpl_size nrow = cpl_table_get_nrow(*storage);
94 cpl_table_insert(*storage, append, nrow);
97 return cpl_error_get_code();
123eris_ifu_resample_tablesave (
const char* recipe_name,
const cpl_parameter *par,
124 const cpl_parameterlist *parlist, cpl_frameset *frameset,
125 cpl_table *tab_final,
const char* pipefile_prefix)
128 char* param_name = cpl_sprintf(
"eris.%s.save-table", recipe_name);
129 par = cpl_parameterlist_find_const (parlist, param_name);
130 cpl_free(param_name);
132 const cpl_boolean save_table = cpl_parameter_get_bool (par);
136 cpl_msg_info (cpl_func,
"Saving the table before resampling");
137 cpl_propertylist *applist = cpl_propertylist_new ();
138 cpl_propertylist_update_string (applist, CPL_DFS_PRO_CATG,
140 char* fname = cpl_sprintf(
"%s_resample_table.fits", pipefile_prefix);
141 cpl_dfs_save_table (frameset, NULL, parlist, frameset, NULL, tab_final,
142 NULL, recipe_name, applist, NULL, PACKAGE
"/" PACKAGE_VERSION,
145 cpl_propertylist_delete (applist);
148 return cpl_error_get_code();
173eris_ifu_resample_parameters_get_int(
const cpl_parameterlist* parlist,
174 const char* pname,
int *pvalue)
176 const cpl_parameter* p = cpl_parameterlist_find_const(parlist, pname);
179 if ( cpl_parameter_get_default_flag(p) && p_has_changed != CPL_FALSE) {
180 *pvalue = cpl_parameter_get_int(p);
182 *pvalue = cpl_parameter_get_default_int(p);
185 return cpl_error_get_code();
209set_double_from_parameters(
const cpl_parameterlist* parlist,
211 const double table_value,
216 if (temp_value != -1)
222 *value = table_value;
225 return cpl_error_get_code();
250set_int_from_parameters(
const cpl_parameterlist* parlist,
252 const int table_value,
256 eris_ifu_resample_parameters_get_int (parlist, pname, &temp_value);
257 if (temp_value != -1)
263 *value = table_value;
266 return cpl_error_get_code();
296static hdrl_parameter *
297eris_ifu_resample_set_method (
const cpl_parameterlist *parlist,
300 hdrl_parameter *aParams_method = NULL;
301 const cpl_parameter *p = NULL;
302 cpl_boolean p_has_changed = 0;
303 int loop_distance = 0;
304 double critical_radius_renka = 0.;
305 int kernel_size_lanczos = 0;
306 double pix_frac_drizzle_x = 0.;
307 double pix_frac_drizzle_y= 0.;
308 double pix_frac_drizzle_l= 0.;
309 cpl_boolean use_errorweights = FALSE;
311 const cpl_parameter *par = NULL;
312 char* param_name = NULL;
315 param_name = cpl_sprintf(
"%s.method.use-errorweights",context);
316 par = cpl_parameterlist_find_const (parlist, param_name);
317 use_errorweights = cpl_parameter_get_bool (par);
318 cpl_free(param_name);
321 param_name = cpl_sprintf(
"%s.method.loop-distance",context);
322 p = cpl_parameterlist_find_const (parlist, param_name);
325 if (cpl_parameter_get_default_flag (p) && p_has_changed != CPL_FALSE)
327 loop_distance = cpl_parameter_get_int (p);
331 loop_distance = LOOP_DISTANCE;
333 cpl_free(param_name);
336 param_name = cpl_sprintf(
"%s.method.renka.critical-radius",context);
337 set_double_from_parameters(parlist, param_name,
338 RENKA_CRITICAL_RADIUS, &critical_radius_renka);
339 cpl_free(param_name);
342 param_name = cpl_sprintf(
"%s.method.lanczos.kernel-size",context);
343 set_int_from_parameters(parlist, param_name,
344 LANCZOS_KERNEL_SIZE, &kernel_size_lanczos);
345 cpl_free(param_name);
352 param_name = cpl_sprintf(
"%s.method.drizzle.downscale-x",context);
353 set_double_from_parameters(parlist, param_name, DRIZZLE_DOWN_SCALING_FACTOR_X,
354 &pix_frac_drizzle_x);
355 cpl_free(param_name);
357 param_name = cpl_sprintf(
"%s.method.drizzle.downscale-y",context);
358 set_double_from_parameters(parlist, param_name, DRIZZLE_DOWN_SCALING_FACTOR_Y,
359 &pix_frac_drizzle_y);
360 cpl_free(param_name);
362 param_name = cpl_sprintf(
"%s.method.drizzle.downscale-z",context);
363 set_double_from_parameters(parlist, param_name, DRIZZLE_DOWN_SCALING_FACTOR_Z,
364 &pix_frac_drizzle_l);
365 cpl_free(param_name);
368 param_name = cpl_sprintf(
"%s.method",context);
369 par = cpl_parameterlist_find_const (parlist, param_name);
370 const char *method = cpl_parameter_get_string (par);
372 if (strcmp (method,
"NEAREST") == 0)
376 else if (strcmp (method,
"RENKA") == 0)
380 critical_radius_renka);
382 else if (strcmp (method,
"LINEAR") == 0)
387 else if (strcmp (method,
"QUADRATIC") == 0)
392 else if (strcmp (method,
"DRIZZLE") == 0)
400 else if (strcmp (method,
"LANCZOS") == 0)
404 kernel_size_lanczos);
410 kernel_size_lanczos);
411 cpl_msg_warning (cpl_func,
412 "%s is an unsupported method! Default to LANCZOS",
415 cpl_free(param_name);
418 return aParams_method;
444eris_ifu_resample_get_cdelt123(cpl_wcs *wcs,
double *cdelt1,
double *cdelt2,
447 cpl_ensure_code(wcs, CPL_ERROR_NULL_INPUT);
459 const cpl_matrix *cd = cpl_wcs_get_cd(wcs);
460 cd1_1 = cpl_matrix_get(cd, 0, 0);
461 cd1_2 = cpl_matrix_get(cd, 1, 0);
462 cd2_1 = cpl_matrix_get(cd, 0, 1);
463 cd2_2 = cpl_matrix_get(cd, 1, 1);
465 dx = sqrt(pow(cd1_1,2) + pow(cd2_1,2));
466 dy = sqrt(pow(cd1_2,2)+ pow(cd2_2,2));
478 cpl_size cdncol = cpl_matrix_get_ncol(cd);
480 cd3_3= fabs(cpl_matrix_get(cd, 2, 2));
489 cpl_msg_debug(cpl_func,
"cdelt1: %g, cdelt2: %g, cdelt3: %g",
490 *cdelt1, *cdelt2, *cdelt3);
492 return cpl_error_get_code();
520static hdrl_parameter *
521eris_ifu_resample_set_outputgrid (
const char* recipe_name,
522 const cpl_parameterlist *parlist, cpl_table *muse_table, cpl_wcs *wcs)
525 cpl_error_ensure(parlist != NULL, CPL_ERROR_NULL_INPUT,
526 return NULL,
"NULL Input Parameters");
527 cpl_error_ensure(muse_table != NULL, CPL_ERROR_NULL_INPUT,
528 return NULL,
"NULL Input table");
529 cpl_error_ensure(wcs != NULL, CPL_ERROR_NULL_INPUT,
530 return NULL,
"NULL Input wcs");
532 hdrl_parameter *aParams_outputgrid = NULL;
533 char* param_name = NULL;
534 char* context = NULL;
539 double lambda_min = 0.;
540 double lambda_max = 0.;
546 context = cpl_sprintf(
"eris.%s", recipe_name);
550 cpl_table_get_column_min (muse_table, HDRL_RESAMPLE_TABLE_RA);
552 cpl_table_get_column_max (muse_table, HDRL_RESAMPLE_TABLE_RA);
554 cpl_table_get_column_min (muse_table, HDRL_RESAMPLE_TABLE_DEC);
556 cpl_table_get_column_max (muse_table, HDRL_RESAMPLE_TABLE_DEC);
557 double lambda_min_tmp =
558 cpl_table_get_column_min (muse_table, HDRL_RESAMPLE_TABLE_LAMBDA);
559 double lambda_max_tmp =
560 cpl_table_get_column_max (muse_table, HDRL_RESAMPLE_TABLE_LAMBDA);
563 if(ra_max_tmp - ra_min_tmp > 180){
564 const double *ra = cpl_table_get_data_double_const(muse_table,
565 HDRL_RESAMPLE_TABLE_RA);
574 cpl_size nrow = cpl_table_get_nrow(muse_table);
576 for (cpl_size i = 0; i < nrow; i++) {
577 if (ra[i] > ra_min_tmp && ra[i] <= 180.) ra_min_tmp = ra[i];
578 if (ra[i] < ra_max_tmp && ra[i] > 180.) ra_max_tmp = ra[i];
583 cpl_msg_debug (cpl_func,
"min x %10.7f", ra_min_tmp);
584 cpl_msg_debug (cpl_func,
"max x %10.7f", ra_max_tmp);
585 cpl_msg_debug (cpl_func,
"min y %10.7f", dec_min_tmp);
586 cpl_msg_debug (cpl_func,
"max y %10.7f", dec_max_tmp);
587 cpl_msg_debug (cpl_func,
"min lambda %10.7f", lambda_min_tmp);
588 cpl_msg_debug (cpl_func,
"max lambda %10.7f", lambda_max_tmp);
607 dec_min = dec_min_tmp;
614 dec_max = dec_max_tmp;
625 lambda_max = lambda_max_tmp;
629 double cdelt1 = 0., cdelt2 = 0., cdelt3 = 0.;
630 eris_ifu_resample_get_cdelt123 (wcs, &cdelt1, &cdelt2, &cdelt3);
654 const cpl_parameter *par = NULL;
655 param_name = cpl_sprintf(
"%s.fieldmargin", context);
656 par = cpl_parameterlist_find_const (parlist, param_name);
657 cpl_free(param_name);
659 double fieldmargin = 0.;
660 fieldmargin = cpl_parameter_get_double(par);
663 int naxis = cpl_wcs_get_image_naxis(wcs);
679 lambda_min, lambda_max,
683 double dist_x = (ra_max-ra_min)/dx;
684 double dist_y = (dec_max-dec_min)/dy;
686 cpl_msg_debug(cpl_func,
"size1: %g size2: %g",dist_x, dist_y);
687 double dist = sqrt(dist_x * dist_x + dist_y * dist_y);
688 cpl_msg_info(cpl_func,
"Output cube diagonal size: %g", dist);
689 if ( dist > MAX_COMBINE_CUBE_PLANE_DIAGONAL) {
690 cpl_msg_info(cpl_func,
"Internal threshold on max cube plane diagonal size: %d",MAX_COMBINE_CUBE_PLANE_DIAGONAL);
691 cpl_msg_warning(cpl_func,
"Resampled cube plane diagonal greater than threshold. Skip cube combination.");
692 cpl_error_set(cpl_func,CPL_ERROR_INCOMPATIBLE_INPUT);
695 return aParams_outputgrid;
700 return aParams_outputgrid;
725static hdrl_parameter *
726eris_ifu_resample_set_outputgrid2D (
const char* recipe_name,
727 const cpl_parameterlist *parlist, cpl_table *muse_table, cpl_wcs *wcs)
730 cpl_error_ensure(parlist != NULL, CPL_ERROR_NULL_INPUT,
731 return NULL,
"NULL Input Parameters");
732 cpl_error_ensure(muse_table != NULL, CPL_ERROR_NULL_INPUT,
733 return NULL,
"NULL Input table");
734 cpl_error_ensure(wcs != NULL, CPL_ERROR_NULL_INPUT,
735 return NULL,
"NULL Input wcs");
737 hdrl_parameter *aParams_outputgrid = NULL;
738 char* param_name = NULL;
739 char* context = NULL;
748 context = cpl_sprintf(
"eris.%s", recipe_name);
752 cpl_table_get_column_min (muse_table, HDRL_RESAMPLE_TABLE_RA);
754 cpl_table_get_column_max (muse_table, HDRL_RESAMPLE_TABLE_RA);
756 cpl_table_get_column_min (muse_table, HDRL_RESAMPLE_TABLE_DEC);
758 cpl_table_get_column_max (muse_table, HDRL_RESAMPLE_TABLE_DEC);
763 if(ra_max_tmp - ra_min_tmp > 180){
764 const double *ra = cpl_table_get_data_double_const(muse_table,
765 HDRL_RESAMPLE_TABLE_RA);
774 cpl_size nrow = cpl_table_get_nrow(muse_table);
776 for (cpl_size i = 0; i < nrow; i++) {
777 if (ra[i] > ra_min_tmp && ra[i] <= 180.) ra_min_tmp = ra[i];
778 if (ra[i] < ra_max_tmp && ra[i] > 180.) ra_max_tmp = ra[i];
783 cpl_msg_debug (cpl_func,
"min x %10.7f", ra_min_tmp);
784 cpl_msg_debug (cpl_func,
"max x %10.7f", ra_max_tmp);
785 cpl_msg_debug (cpl_func,
"min y %10.7f", dec_min_tmp);
786 cpl_msg_debug (cpl_func,
"max y %10.7f", dec_max_tmp);
805 dec_min = dec_min_tmp;
812 dec_max = dec_max_tmp;
817 double cdelt1 = 0., cdelt2 = 0., cdelt3 = 0.;
818 eris_ifu_resample_get_cdelt123 (wcs, &cdelt1, &cdelt2, &cdelt3);
833 const cpl_parameter *par = NULL;
834 param_name = cpl_sprintf(
"%s.fieldmargin", context);
835 par = cpl_parameterlist_find_const (parlist, param_name);
836 cpl_free(param_name);
838 double fieldmargin = 0.;
839 fieldmargin = cpl_parameter_get_double(par);
851 double dist_x = (ra_max-ra_min)/dx;
852 double dist_y = (dec_max-dec_min)/dy;
854 cpl_msg_debug(cpl_func,
"size1: %g size2: %g",dist_x, dist_y);
855 double dist = sqrt(dist_x * dist_x + dist_y * dist_y);
856 cpl_msg_debug(cpl_func,
"Output cube diagonal size: %g", dist);
857 if ( dist > MAX_COMBINE_CUBE_PLANE_DIAGONAL) {
858 cpl_msg_info(cpl_func,
"Internal threshold on max cube plane diagonal size: %d",MAX_COMBINE_CUBE_PLANE_DIAGONAL);
859 cpl_msg_warning(cpl_func,
"Resampled cube plane diagonal %g greater than threshold. Skip cube combination.", dist);
860 cpl_error_set(cpl_func,CPL_ERROR_INCOMPATIBLE_INPUT);
862 return aParams_outputgrid;
867 return aParams_outputgrid;
903 const char *filename,
904 const cpl_parameterlist *parlist,
905 cpl_frameset *frameset,
906 cpl_boolean gen_phase3)
909 cpl_ensure_code(aCube && aCube->header, CPL_ERROR_NULL_INPUT);
913 cpl_propertylist* dataHeader;
914 cpl_propertylist* errsHeader;
915 cpl_propertylist* qualHeader;
918 deqErrorType errorType = rmse;
919 deqQualityType qualityType = flag16bit;
920 char* pipe_id = cpl_sprintf(
"%s%s%s", PACKAGE,
"/", PACKAGE_VERSION);
923 if(cpl_propertylist_has(aCube->header,
"CDELT3")) {
924 cpl_propertylist_erase(aCube->header,
"CDELT3");
927 if(!cpl_propertylist_has(aCube->header, CUNIT3)) {
928 cpl_propertylist_append_string(aCube->header,
"CUNIT3",
"um");
930 if(!cpl_propertylist_has(aCube->header,
"BUNIT")) {
931 cpl_propertylist_append_string(aCube->header,
"BUNIT", UNIT_ADU);
935 dataHeader = cpl_propertylist_duplicate(aCube->header);
936 errsHeader = cpl_propertylist_duplicate(aCube->header);
937 qualHeader = cpl_propertylist_duplicate(aCube->header);
938 cpl_propertylist_update_string(qualHeader,
"BUNIT",
"");
949 if(cpl_propertylist_has(dataHeader,
"CTYPE1")) {
950 const char* ctype1 = cpl_propertylist_get_string(dataHeader,
"CTYPE1");
951 if (!cpl_propertylist_has(errsHeader,
"CTYPE1"))
952 cpl_propertylist_append_string(errsHeader,
"CTYPE1", ctype1);
953 if (!cpl_propertylist_has(qualHeader,
"CTYPE1"))
954 cpl_propertylist_append_string(qualHeader,
"CTYPE1", ctype1);
956 if(cpl_propertylist_has(dataHeader,
"CTYPE2")) {
957 const char* ctype2 = cpl_propertylist_get_string(dataHeader,
"CTYPE2");
958 if (!cpl_propertylist_has(errsHeader,
"CTYPE2"))
959 cpl_propertylist_append_string(errsHeader,
"CTYPE2",ctype2);
960 if (!cpl_propertylist_has(errsHeader,
"CTYPE2"))
961 cpl_propertylist_append_string(qualHeader,
"CTYPE2", ctype2);
964 cpl_propertylist* applist = cpl_propertylist_duplicate(aCube->header);
973 frameset, parlist, recipe);
981 eris_ifu_plist_erase_wcs(applist);
982 eris_pfits_clean_header(dataHeader, CPL_TRUE);
983 eris_pfits_clean_header(errsHeader, CPL_TRUE);
984 eris_pfits_clean_header(qualHeader, CPL_TRUE);
985 cpl_propertylist_append_string(applist,
"PRODCATG",
"SCIENCE.CUBE.IFS");
987 cpl_propertylist_update_string(applist, CPL_DFS_PRO_CATG, procatg);
988 if(cpl_propertylist_has(applist,
"BUNIT")) {
989 cpl_propertylist_erase(applist,
"BUNIT");
991 if(cpl_propertylist_has(applist,
"COMMENT")) {
992 cpl_msg_warning(cpl_func,
"remove comments");
993 cpl_propertylist_erase(applist,
"COMMENT");
995 if(cpl_propertylist_has(dataHeader,
"COMMENT")) {
996 cpl_msg_warning(cpl_func,
"remove comments");
997 cpl_propertylist_erase(dataHeader,
"COMMENT");
999 cpl_dfs_save_propertylist(frameset, NULL, parlist, frameset,
1000 NULL, recipe, applist,
"RADECSYS", pipe_id, filename);
1002 cpl_imagelist* data = cpl_imagelist_new();
1003 cpl_imagelist* errs = cpl_imagelist_new();
1004 cpl_imagelist* qual = cpl_imagelist_new();
1009 cpl_type data_type = cpl_image_get_type(cpl_imagelist_get(data,0));
1011 cpl_imagelist_save(data, filename, data_type, dataHeader, CPL_IO_EXTEND);
1012 cpl_imagelist_save(errs, filename, data_type, errsHeader, CPL_IO_EXTEND);
1013 cpl_imagelist_save(qual, filename, CPL_TYPE_INT, qualHeader, CPL_IO_EXTEND);
1014 cpl_size ix_max = cpl_imagelist_get_size(data);
1016 for (cpl_size ix=ix_max-1; ix >= 0; ix--) {
1017 cpl_imagelist_unset(data, ix);
1018 cpl_imagelist_unset(errs, ix);
1019 q=cpl_imagelist_unset(qual, ix);
1020 cpl_image_delete(q);
1023 cpl_imagelist_delete(data);
1024 cpl_imagelist_delete(errs);
1025 cpl_imagelist_delete(qual);
1027 cpl_propertylist_delete(applist);
1047 return cpl_error_get_code();
1068eris_ifu_resample_get_wcs_from_frameset(cpl_frameset* frameset,
1069 const char* procatg) {
1071 cpl_ensure(frameset, CPL_ERROR_NULL_INPUT, NULL);
1074 cpl_frameset* in_set = NULL;
1077 cpl_error_set_message(cpl_func, CPL_ERROR_FILE_NOT_FOUND,
"Missing RAW "
1082 cpl_frame* frame = cpl_frameset_get_position(in_set, 0);
1083 const char* fname = cpl_frame_get_filename(frame);
1086 cpl_wcs *wcs = NULL;
1087 cpl_propertylist* head = NULL;
1089 cpl_errorstate prestate = cpl_errorstate_get();
1090 head = cpl_propertylist_load(fname, 0);
1091 wcs = cpl_wcs_new_from_propertylist(head);
1094 cpl_errorstate_set(prestate);
1095 cpl_propertylist_delete(head);
1097 cpl_propertylist_delete(head);
1098 cpl_frameset_delete(in_set);
1102 prestate = cpl_errorstate_get();
1103 head = cpl_propertylist_load(fname, 1);
1104 wcs = cpl_wcs_new_from_propertylist(head);
1105 cpl_propertylist_delete(head);
1106 cpl_frameset_delete(in_set);
1141eris_ifu_resample_frameset_to_table(
const cpl_frame *data_frm,
1142 const cpl_size data_ext_id,
1143 const cpl_frame *errs_frm,
1144 const cpl_size errs_ext_id,
1145 const cpl_frame *qual_frm,
1146 const cpl_size qual_ext_id,
1147 const cpl_boolean is_variance,
1148 const cpl_boolean subtract_bkg,
1152 cpl_ensure(data_frm, CPL_ERROR_NULL_INPUT, NULL);
1153 cpl_msg_severity level = cpl_msg_get_level();
1154 cpl_msg_set_level(CPL_MSG_INFO);
1155 const char *name = cpl_frame_get_filename(data_frm);
1156 cpl_imagelist *dlist = NULL;
1157 cpl_imagelist *elist = NULL;
1158 cpl_imagelist *qlist = NULL;
1160 cpl_table *tab_ext = NULL;
1162 cpl_errorstate prestate = cpl_errorstate_get();
1164 dlist = cpl_imagelist_load(name, CPL_TYPE_DOUBLE, data_ext_id);
1166 if (dlist == NULL) {
1168 if (!cpl_errorstate_is_equal(prestate)) {
1169 cpl_errorstate_set(prestate);
1174 cpl_propertylist *xheader_data = cpl_propertylist_load(name,
1176 cpl_wcs *wcs = cpl_wcs_new_from_propertylist(xheader_data);
1178 if (errs_frm != NULL) {
1179 name = cpl_frame_get_filename(errs_frm);
1180 elist = cpl_imagelist_load(name, CPL_TYPE_DOUBLE, errs_ext_id);
1183 cpl_imagelist_power(elist, 0.5);
1186 if (qual_frm != NULL) {
1187 name = cpl_frame_get_filename(qual_frm);
1188 qlist = cpl_imagelist_load(name, CPL_TYPE_INT, qual_ext_id);
1191 cpl_size size = cpl_imagelist_get_size(dlist);
1194 if (qual_frm != NULL && size > 0){
1195 for(cpl_size k = 0; k < size; k++) {
1196 cpl_image* data = cpl_imagelist_get(dlist, k);
1197 cpl_image* qual = cpl_imagelist_get(qlist, k);
1201 cpl_mask* mask = cpl_mask_threshold_image_create(qual, 0, INT_MAX);
1203 cpl_image_reject_from_mask(data, mask);
1204 cpl_mask_delete(mask);
1205 cpl_imagelist_set(dlist, data, k);
1212 if(elist != NULL && is_variance){
1213 for(cpl_size k = 0; k < size; k++) {
1214 cpl_image* data = cpl_imagelist_get(dlist, k);
1215 cpl_mask* data_mask = cpl_image_get_bpm(data);
1217 cpl_image* error = cpl_imagelist_get(elist, k);
1218 cpl_mask* error_mask = cpl_image_get_bpm(error);
1220 cpl_mask* merged = cpl_mask_duplicate(data_mask);
1223 cpl_mask_or(merged, error_mask);
1224 cpl_image_reject_from_mask(data, merged);
1225 cpl_image_reject_from_mask(error, merged);
1226 cpl_mask_delete(merged);
1231 hdrl_imagelist* hlist = NULL;
1235 int edge_trim = edgetrim;
1236 if (edge_trim > 0) {
1237 cpl_msg_info(cpl_func,
"Trim input image edges of %d pixels", edge_trim);
1242 if (edge_trim >= 0.5* sx) {
1244 cpl_msg_warning(cpl_func,
"edge-trim must be smaller than half image "
1245 "size. Reset to 0");
1247 if (edge_trim >= 0.5* sy) {
1249 cpl_msg_warning(cpl_func,
"edge-trim must be smaller than half image "
1250 "size. Reset to 0");
1252 for(cpl_size k = 0; k < sz; k++) {
1257 for(cpl_size j = 1; j <= edge_trim; j++) {
1258 for(cpl_size i = 1; i <= sx; i++) {
1263 for(cpl_size j = sy - edge_trim + 1; j <= sy; j++) {
1264 for(cpl_size i = 1; i <= sx; i++) {
1269 for(cpl_size j = 1; j <= sy; j++) {
1270 for(cpl_size i = 1; i <= edge_trim; i++) {
1275 for(cpl_size j = 1; j <= sy; j++) {
1276 for(cpl_size i = sx - edge_trim + 1; i <= sx; i++) {
1285 cpl_imagelist_delete(dlist);
1286 if (qual_frm != NULL)
1287 cpl_imagelist_delete(qlist);
1288 if (errs_frm != NULL)
1289 cpl_imagelist_delete(elist);
1294 cpl_msg_info(cpl_func,
"Reading the image ...");
1298 if(subtract_bkg == CPL_TRUE) {
1299 cpl_msg_info(cpl_func,
"Subtracting the median as requested ...");
1306 cpl_msg_info(cpl_func,
"Converting imagelist to table with hdrl function");
1313 cpl_wcs_delete(wcs);
1314 cpl_propertylist_delete(xheader_data);
1316 cpl_msg_set_level(level);
1338static cpl_error_code
1339eris_ifu_resample_wcs_print(cpl_wcs *wcs)
1341 cpl_ensure_code(wcs, CPL_ERROR_NULL_INPUT);
1343 const cpl_array *crval = cpl_wcs_get_crval(wcs);
1344 const cpl_array *crpix = cpl_wcs_get_crpix(wcs);
1345 const cpl_array *ctype = cpl_wcs_get_ctype(wcs);
1346 const cpl_array *cunit = cpl_wcs_get_cunit(wcs);
1348 const cpl_matrix *cd = cpl_wcs_get_cd(wcs);
1349 const cpl_array *dims = cpl_wcs_get_image_dims(wcs);
1350 cpl_size naxis = cpl_wcs_get_image_naxis(wcs);
1352 cpl_msg_info(cpl_func,
"NAXIS: %lld", naxis);
1355 cpl_msg_indent_more();
1357 for (cpl_size i = 0; i < naxis; i++) {
1358 cpl_msg_info(cpl_func,
"NAXIS%lld: %d", i + 1,
1359 cpl_array_get_int(dims, i, &testerr));
1361 cpl_msg_indent_less();
1363 double cd11 = cpl_matrix_get(cd, 0, 0);
1364 double cd12 = cpl_matrix_get(cd, 0, 1);
1365 double cd21 = cpl_matrix_get(cd, 1, 0);
1366 double cd22 = cpl_matrix_get(cd, 1, 1);
1367 double crpix1 = cpl_array_get_double(crpix, 0, &testerr);
1368 double crpix2 = cpl_array_get_double(crpix, 1, &testerr);
1369 double crval1 = cpl_array_get_double(crval, 0, &testerr);
1370 double crval2 = cpl_array_get_double(crval, 1, &testerr);
1372 cpl_msg_info(cpl_func,
"1st and 2nd dimension");
1373 cpl_msg_indent_more();
1374 cpl_msg_info(cpl_func,
"CD1_1: %g", cd11);
1375 cpl_msg_info(cpl_func,
"CD1_2: %g", cd12);
1376 cpl_msg_info(cpl_func,
"CD2_1: %g", cd21);
1377 cpl_msg_info(cpl_func,
"CD2_2: %g", cd22);
1379 cpl_msg_info(cpl_func,
"CRPIX1: %g", crpix1);
1380 cpl_msg_info(cpl_func,
"CRPIX2: %g", crpix2);
1381 cpl_msg_info(cpl_func,
"CRVAL1: %f", crval1);
1382 cpl_msg_info(cpl_func,
"CRVAL2: %f", crval2);
1384 cpl_msg_info(cpl_func,
"CTYPE1: %s", cpl_array_get_string(ctype, 0));
1385 cpl_msg_info(cpl_func,
"CTYPE2: %s", cpl_array_get_string(ctype, 1));
1389 cpl_msg_info(cpl_func,
"CUNIT1: %s", cpl_array_get_string(cunit, 0));
1390 cpl_msg_info(cpl_func,
"CUNIT2: %s", cpl_array_get_string(cunit, 1));
1392 cpl_msg_indent_less();
1395 cpl_size cdncol = cpl_matrix_get_ncol(cd);
1398 double cd13 = cpl_matrix_get(cd, 0, 2);
1399 double cd23 = cpl_matrix_get(cd, 1, 2);
1400 double cd31 = cpl_matrix_get(cd, 2, 0);
1401 double cd32 = cpl_matrix_get(cd, 2, 1);
1402 double cd33 = cpl_matrix_get(cd, 2, 2);
1403 double crval3 = cpl_array_get_double(crval, 2, &testerr);
1404 double crpix3 = cpl_array_get_double(crpix, 2, &testerr);
1406 cpl_msg_info(cpl_func,
"3rd dimension");
1407 cpl_msg_indent_more();
1408 cpl_msg_info(cpl_func,
"CD1_3: %g", cd13);
1409 cpl_msg_info(cpl_func,
"CD2_3: %g", cd23);
1410 cpl_msg_info(cpl_func,
"CD3_1: %g", cd31);
1411 cpl_msg_info(cpl_func,
"CD3_2: %g", cd32);
1412 cpl_msg_info(cpl_func,
"CD3_3: %g", cd33);
1414 cpl_msg_info(cpl_func,
"CRPIX3: %g", crpix3);
1415 cpl_msg_info(cpl_func,
"CRVAL3: %g", crval3);
1418 cpl_msg_info(cpl_func,
"CTYPE3: %s", cpl_array_get_string(ctype, 2));
1422 cpl_msg_info(cpl_func,
"CUNIT3: %s", cpl_array_get_string(cunit, 2));
1425 cpl_msg_indent_less();
1428 return cpl_error_get_code();
1456eris_ifu_resample_get_table_from_frameset(
const cpl_frame* data_frm,
1457 const cpl_frame* errs_frm,
1458 const cpl_frame* qual_frm,
1459 const cpl_size data_ext_id,
1460 const cpl_size errs_ext_id,
1461 const cpl_size qual_ext_id,
1462 const cpl_boolean is_variance,
1463 const cpl_boolean subtract_bkg,
1464 const int edge_trim) {
1466 cpl_ensure(data_frm, CPL_ERROR_NULL_INPUT, NULL);
1468 cpl_size next = cpl_frame_get_nextensions(data_frm);
1469 cpl_msg_info(cpl_func,
"Analysing and processing file %s",
1470 cpl_frame_get_filename(data_frm));
1478 cpl_msg_indent_more();
1479 cpl_table* table = NULL;
1480 if(data_ext_id != -1) {
1481 cpl_msg_info(cpl_func,
"Extension: %02lld", data_ext_id);
1483 table = eris_ifu_resample_frameset_to_table(data_frm, data_ext_id,
1484 errs_frm, errs_ext_id,
1485 qual_frm, qual_ext_id,
1486 is_variance, subtract_bkg, edge_trim);
1489 for (cpl_size i = 0; i < next + 1; i++ ) {
1490 cpl_table * table_local = NULL;
1491 cpl_msg_info(cpl_func,
"Extension: %02lld", i);
1492 table_local = eris_ifu_resample_frameset_to_table(data_frm, i,
1495 is_variance, subtract_bkg, edge_trim);
1498 if (table_local == NULL) {
1499 cpl_msg_info(cpl_func,
"No siutable data found - continuing");
1505 table = table_local;
1508 cpl_size nrow = cpl_table_get_nrow(table);
1509 cpl_table_insert(table, table_local, nrow);
1510 cpl_table_delete(table_local);
1515 cpl_msg_indent_less();
1753 cpl_frameset* frameset,
1754 const cpl_parameterlist* parlist,
1755 const char* recipe_name,
1757 cpl_boolean apply_flat,
1758 cpl_boolean is_pupil)
1760 hdrl_image* cube_collapsed = NULL;
1761 cpl_image* cube_cmap = NULL;
1763 cpl_frame* frame = cpl_frameset_find(frameset, pro_catg);
1764 const char* cube_fname = cpl_frame_get_filename(frame);
1765 cpl_propertylist* pheader = cpl_propertylist_load(cube_fname, 0);
1766 cpl_propertylist* head_wcs_2d = eris_ifu_plist_extract_wcs2D(pheader);
1767 cpl_propertylist* dheader = cpl_propertylist_load(cube_fname, 1);
1768 cpl_propertylist* eheader = cpl_propertylist_load(cube_fname, 2);
1769 cpl_propertylist* qheader = cpl_propertylist_load(cube_fname, 3);
1771 cpl_imagelist* iml_data = cpl_imagelist_load(cube_fname, CPL_TYPE_DOUBLE, 1);
1772 cpl_imagelist* iml_errs = cpl_imagelist_load(cube_fname, CPL_TYPE_DOUBLE, 2);
1773 cpl_imagelist* iml_qual = cpl_imagelist_load(cube_fname, CPL_TYPE_DOUBLE, 3);
1775 cpl_propertylist_append(dheader, head_wcs_2d);
1776 cpl_propertylist_append(eheader, head_wcs_2d);
1777 cpl_propertylist_append(qheader, head_wcs_2d);
1778 cpl_propertylist_delete(head_wcs_2d);
1781 double exptime = dit * ndit;
1783 cpl_size sz = cpl_imagelist_get_size(iml_data);
1784 cpl_image* data = NULL;
1785 cpl_image* errs = NULL;
1786 cpl_image* qual = NULL;
1787 cpl_mask* bpm_data = NULL;
1788 cpl_mask* bpm_errs = NULL;
1789 cpl_boolean edge_trim = CPL_TRUE;
1790 cpl_size k_center = 0.5 * sz;
1791 for(cpl_size k = 0; k < sz; k++) {
1792 data = cpl_imagelist_get(iml_data, k);
1793 errs = cpl_imagelist_get(iml_errs, k);
1794 qual = cpl_imagelist_get(iml_qual, k);
1795 cpl_image_reject_value(data, CPL_VALUE_NAN);
1796 cpl_image_reject_value(errs, CPL_VALUE_NAN);
1797 bpm_data = cpl_image_get_bpm(data);
1798 bpm_errs = cpl_image_get_bpm(errs);
1799 cpl_mask * bpm_mask = cpl_mask_threshold_image_create(qual, -0.5, 0.5) ;
1800 cpl_mask_not(bpm_mask) ;
1801 cpl_mask_or(bpm_mask,bpm_data);
1802 cpl_mask_or(bpm_mask,bpm_errs);
1803 cpl_image_reject_from_mask(data, bpm_mask);
1804 cpl_image_reject_from_mask(errs, bpm_mask);
1815 cpl_image* ima_time = NULL;
1816 cpl_image* err_time = NULL;
1817 hdrl_image *himg_exptime = NULL;
1818 if(apply_flat && k == k_center) {
1819 ima_time = cpl_image_new(cpl_image_get_size_x(data),
1820 cpl_image_get_size_y(data), CPL_TYPE_DOUBLE);
1821 err_time = cpl_image_new(cpl_image_get_size_x(data),
1822 cpl_image_get_size_y(data), CPL_TYPE_DOUBLE);
1825 cpl_image_add_scalar(ima_time, exptime);
1826 cpl_image_add_scalar(err_time, sqrt(exptime));
1830 if (edge_trim > 0) {
1841 cpl_image_delete(ima_time);
1842 cpl_image_delete(err_time);
1846 cpl_propertylist* phead2D = cpl_propertylist_duplicate(dheader);
1847 char* fname = cpl_sprintf(
"%s_exposure_map.fits", recipe_name);
1848 eris_ifu_plist_erase_wcs3D(phead2D);
1849 eris_ifu_plist_erase_expmap_extra_keys(phead2D);
1850 cpl_propertylist_update_string(phead2D, CPL_DFS_PRO_CATG, ERIS_IFU_PRO_JITTER_EXPOSURE_MAP);
1851 cpl_propertylist_update_string(phead2D,
"PRODCATG", PRODCATG_EXPOSUREMAP);
1853 recipe_name, phead2D,
"RADECSYS", fname,
1857 cpl_propertylist_delete(phead2D);
1863 cpl_mask_delete(bpm_mask) ;
1872 cpl_propertylist* phead2D = cpl_propertylist_duplicate(pheader);
1873 cpl_propertylist_erase_regexp(phead2D,NAXIS3,0);
1874 cpl_propertylist_erase_regexp(phead2D,CTYPE3,0);
1875 cpl_propertylist_erase_regexp(phead2D,CRVAL3,0);
1876 cpl_propertylist_erase_regexp(phead2D,CRPIX3,0);
1877 cpl_propertylist_erase_regexp(phead2D,CUNIT3, 0);
1878 cpl_propertylist_erase_regexp(phead2D,
"CD3_*",0);
1879 cpl_propertylist_erase_regexp(phead2D,CD1_3,0);
1880 cpl_propertylist_erase_regexp(phead2D,CD2_3,0);
1883 char* prefix = NULL;
1884 char* suffix = NULL;
1886 suffix = cpl_sprintf(
"%s",
"");
1888 suffix = cpl_sprintf(
"%s",
"_no_flat");
1890 if(strstr(recipe_name,
"jitter") != NULL ) {
1891 prefix = cpl_sprintf(
"eris_ifu_jitter%s", suffix);
1893 prefix = cpl_sprintf(
"eris_ifu_stdstar%s", suffix);
1896 const char* pcatg = NULL;
1897 if(strstr(pro_catg,
"OBJ") != NULL ) {
1898 fname = cpl_sprintf(
"%s_obj_cube_mean.fits", prefix);
1900 pcatg = ERIS_IFU_PRO_JITTER_OBJ_CUBE_MEAN;
1902 pcatg = ERIS_IFU_PRO_JITTER_OBJ_CUBE_NOFLAT_MEAN;
1904 }
else if(strstr(pro_catg, ERIS_IFU_PRO_JITTER_DAR_CUBE) != NULL ) {
1905 fname = cpl_sprintf(
"%s_dar_cube_mean.fits", prefix);
1907 pcatg = ERIS_IFU_PRO_JITTER_OBJ_DAR_CUBE_MEAN;
1909 pcatg = ERIS_IFU_PRO_JITTER_OBJ_CUBE_NOFLAT_MEAN;
1911 }
else if(strstr(pro_catg, ERIS_IFU_PRO_JITTER_TWK_CUBE) != NULL ) {
1912 fname = cpl_sprintf(
"%s_twk_cube_mean.fits", prefix);
1914 pcatg = ERIS_IFU_PRO_JITTER_TWK_CUBE_MEAN;
1916 pcatg = ERIS_IFU_PRO_JITTER_OBJ_CUBE_NOFLAT_MEAN;
1918 }
else if(strstr(pro_catg,
"STD") != NULL ) {
1919 fname = cpl_sprintf(
"%s_std_cube_mean.fits", prefix);
1921 pcatg = ERIS_IFU_PRO_JITTER_STD_CUBE_MEAN;
1923 pcatg = ERIS_IFU_PRO_JITTER_STD_CUBE_NOFLAT_MEAN;
1925 }
else if(strstr(pro_catg,
"PSF") != NULL ) {
1926 fname = cpl_sprintf(
"%s_psf_cube_mean.fits", prefix);
1928 pcatg = ERIS_IFU_PRO_JITTER_PSF_CUBE_MEAN;
1930 pcatg = ERIS_IFU_PRO_JITTER_PSF_CUBE_NOFLAT_MEAN;
1934 cpl_propertylist_update_string(phead2D, CPL_DFS_PRO_CATG, pcatg);
1945 cpl_size max_ima_x = 0;
1946 cpl_size max_ima_y = 0;
1947 cpl_image_get_maxpos(img, &max_ima_x, &max_ima_y);
1950 double max_ima_cx = cpl_image_get_centroid_x_window(img, 1, 1, sx, sy);
1951 double max_ima_cy = cpl_image_get_centroid_y_window(img, 1, 1, sx, sy);
1952 double xshift = max_ima_cx - (double) sx * 0.5;
1953 double yshift = max_ima_cy - (double) sy * 0.5;
1954 char* key_name = cpl_sprintf(
"QC PUPIL%d SHIFTX", 0);
1957 "[pix] X shift centroid - center image");
1960 key_name = cpl_sprintf(
"QC PUPIL%d SHIFTY", 0);
1962 "[pix] Y shift centroid - center image");
1963 cpl_msg_info(cpl_func,
"xshift: %g yshift: %g", xshift, yshift);
1967 cpl_free(qclog_tbl);
1972 recipe_name, phead2D,
"RADECSYS", fname,
1975 rmse, mask_ima, flag16bit, UNIT_ADU);
1977 cpl_image_delete(mask_ima);
1980 cpl_propertylist_delete(dheader);
1981 cpl_propertylist_delete(eheader);
1982 cpl_propertylist_delete(qheader);
1983 cpl_propertylist_delete(pheader);
1984 cpl_imagelist_delete(iml_data);
1985 cpl_imagelist_delete(iml_errs);
1986 cpl_imagelist_delete(iml_qual);
1988 cpl_image_delete(cube_cmap);
1991 return cpl_error_get_code();
1996eris_ifu_compute_max_cubes_dist(cpl_frameset* cube_set)
1999 double dist_max = -1;
2000 double dist_sqr = -1;
2001 cpl_size ncubes = cpl_frameset_get_size(cube_set);
2002 cpl_vector* vec_ra = cpl_vector_new(ncubes);
2003 cpl_vector* vec_dec = cpl_vector_new(ncubes);
2004 cpl_propertylist* plist = NULL;
2005 const char* name = NULL;
2006 double ra = 0, dec = 0;
2007 double ra1 = 0, dec1 = 0;
2008 double ra2 = 0, dec2 = 0;
2012 double pix_size = 0;
2014 cpl_propertylist* phead = NULL;
2018 cpl_frame* data_frame = NULL;
2020 data_frame = cpl_frameset_get_position(cube_set, 0);
2021 name = cpl_frame_get_filename(data_frame);
2022 plist = cpl_propertylist_load(name, extnum_raw);
2023 phead = cpl_propertylist_load(name, 0);
2030 cpl_propertylist_delete(plist);
2031 cpl_propertylist_delete(phead);
2034 case S250MAS: pix_size = 0.250;
break;
2035 case S100MAS: pix_size = 0.100;
break;
2036 case S25MAS: pix_size = 0.025;
break;
2037 case PUPIL: pix_size = 0.025;
break;
2038 case UNDEFINED_SCALE: pix_size = 0.100;
break;
2047 for (cpl_size i = 0; i < ncubes ; i++) {
2049 data_frame = cpl_frameset_get_position(cube_set, i);
2050 name = cpl_frame_get_filename(data_frame);
2051 plist = cpl_propertylist_load(name, extnum_raw);
2053 ra = cpl_propertylist_get_double(plist,
"RA");
2054 dec = cpl_propertylist_get_double(plist,
"DEC");
2056 cpl_vector_set(vec_ra, i, ra);
2057 cpl_vector_set(vec_dec, i, dec);
2059 cpl_propertylist_delete(plist);
2061 for (cpl_size i = 0; i < ncubes; i++) {
2062 ra1 = cpl_vector_get(vec_ra, i);
2063 dec1 = cpl_vector_get(vec_dec, i);
2064 for (cpl_size j = 0; j < ncubes && i != j; j++) {
2065 ra2 = cpl_vector_get(vec_ra, j);
2066 dec2 = cpl_vector_get(vec_dec, j);
2073 dist = (ra2 - ra1) * (ra2 - ra1) / pix_size / pix_size +
2074 (dec2 - dec1) * (dec2 - dec1) / pix_size / pix_size;
2076 if(dist > dist_sqr) {
2081 cpl_vector_delete(vec_ra);
2082 cpl_vector_delete(vec_dec);
2083 dist_max = sqrt(dist_sqr);
2084 cpl_msg_info(cpl_func,
"Max distance between contributing cubes centers: %g",dist_max);
2143 const cpl_parameterlist * parlist,
2144 const char* recipe_name,
2145 const char* pipefile_prefix)
2147 cpl_size nframes = cpl_frameset_get_size(frameset);
2149 cpl_frameset* cube_set = NULL;
2150 char *proCatg = NULL;
2151 char *filenameSpec = NULL;
2152 cpl_msg_info(cpl_func,
"Combine cubes into a single one, compute mean, median along wavelength axis");
2157 cpl_msg_warning(cpl_func,
"No %s files", proCatg);
2158 cpl_error_set_message(cpl_func, CPL_ERROR_FILE_NOT_FOUND,
"Missing CUBE "
2160 return cpl_error_get_code();
2162 cpl_msg_info(cpl_func,
"%02d file(s) of type %s",
2163 (
int)cpl_frameset_get_size(cube_set), proCatg);
2167 cpl_table* restable = NULL;
2168 cpl_table* tab_tmp = NULL;
2169 const cpl_parameter * par = NULL;
2170 cpl_frame * data_frame = NULL;
2171 cpl_frame * errs_frame = NULL;
2172 cpl_frame * qual_frame = NULL;
2176 cpl_boolean is_variance = CPL_TRUE;
2177 cpl_boolean subtract_bkg = CPL_FALSE;
2180 cpl_size ncubes = cpl_frameset_get_size(cube_set);
2181 char* pname = cpl_sprintf(
"eris.%s.max-cubes-centres-dist",recipe_name);
2182 int max_cubes_dist = 20;
2185 max_cubes_dist = cpl_parameter_get_int(
2186 cpl_parameterlist_find_const(parlist, pname));
2188 if ( max_cubes_dist < eris_ifu_compute_max_cubes_dist(cube_set)) {
2189 cpl_msg_info(cpl_func,
"max-cubes-centres-dist: %d",max_cubes_dist);
2190 cpl_msg_warning(cpl_func,
"max distance between cube centres greater than threshold. Skip cube combination.");
2191 return CPL_ERROR_INCOMPATIBLE_INPUT;
2194 for (cpl_size i = 0; i < ncubes ; i++) {
2196 data_frame = cpl_frameset_get_position(cube_set, i);
2197 errs_frame = data_frame;
2198 qual_frame = data_frame;
2200 tab_tmp = eris_ifu_resample_get_table_from_frameset(data_frame,
2202 qual_frame, extnum_raw,
2203 extnum_err, extnum_bpm,
2204 is_variance, subtract_bkg, edge_trim);
2210 eris_ifu_resample_update_table(tab_tmp, &restable);
2211 cpl_table_delete(tab_tmp);
2216 cpl_ensure_code(restable, CPL_ERROR_NULL_INPUT);
2219 eris_ifu_resample_tablesave (recipe_name, par, parlist, frameset, restable,
2222 cpl_wcs *wcs = eris_ifu_resample_get_wcs_from_frameset(cube_set, proCatg);
2223 cpl_ensure_code(wcs, CPL_ERROR_NULL_INPUT);
2225 cpl_msg_info(cpl_func,
"WCS used for the output grid definition and passed to "
2226 "hdrl_resample_compute():");
2227 cpl_msg_indent_more();
2228 eris_ifu_resample_wcs_print(wcs);
2229 cpl_msg_indent_less();
2231 hdrl_parameter *aParams_method = NULL;
2234 char *context = NULL;
2235 context = cpl_sprintf(
"eris.%s", recipe_name);
2236 aParams_method = eris_ifu_resample_set_method(parlist, context);
2239 cpl_ensure_code(aParams_method, CPL_ERROR_NULL_INPUT);
2242 hdrl_parameter *aParams_outputgrid = NULL;
2244 aParams_outputgrid = eris_ifu_resample_set_outputgrid(recipe_name, parlist,
2246 if(cpl_error_get_code() != CPL_ERROR_NONE) {
2249 return CPL_ERROR_INCOMPATIBLE_INPUT;
2252 cpl_ensure_code(aParams_outputgrid, CPL_ERROR_NULL_INPUT);
2254 hdrl_resample_result *result = NULL;
2259 cpl_msg_severity level = cpl_msg_get_level();
2260 cpl_msg_set_level(CPL_MSG_INFO);
2263 cpl_msg_set_level(level);
2275 cpl_ensure_code(result, CPL_ERROR_NULL_INPUT);
2276 char* pro_catg = NULL;
2277 if(strstr(pipefile_prefix,
"no_flat") != NULL) {
2278 pro_catg = cpl_sprintf(
"%s_%s",proCatg,
"COADD_NOFLAT");
2280 pro_catg = cpl_sprintf(
"%s_%s",proCatg,
"COADD");
2282 char* fname = cpl_sprintf(
"%s_%s_coadd.fits", pipefile_prefix, filenameSpec);
2285 frameset, CPL_FALSE);
2289 hdrl_image* cube_collapsed = NULL;
2290 cpl_image* cube_cmap = NULL;
2294 pro_catg = cpl_sprintf(
"%s_%s",proCatg,
"MEAN");
2295 fname = cpl_sprintf(
"%s_%s_mean.fits", pipefile_prefix, filenameSpec);
2296 cpl_propertylist* phead2D = cpl_propertylist_duplicate(result->header);
2299 cpl_propertylist_append_string(phead2D,
"PRODCATG",PRODCATG_WHITELIGHT);
2304 NULL, CPL_IO_EXTEND);
2306 fname, CPL_TYPE_INT, NULL, CPL_IO_EXTEND);
2314 cpl_image_delete(cube_cmap);
2319 cpl_wcs_delete(wcs);
2320 wcs = cpl_wcs_new_from_propertylist(result->header);
2321 cpl_msg_info(cpl_func,
"Final WCS after resampling: ");
2322 cpl_msg_indent_more();
2323 eris_ifu_resample_wcs_print(wcs);
2324 cpl_msg_indent_less();
2333 cpl_table_delete(restable);
2334 cpl_wcs_delete(wcs);
2335 cpl_frameset_delete(cube_set);
2340 return cpl_error_get_code();
2356static cpl_propertylist* eris_ifu_resample_get_header2D(cpl_propertylist *header3D){
2357 cpl_propertylist *phead2D = cpl_propertylist_duplicate(header3D);
2358 cpl_propertylist_set_int(phead2D,
"NAXIS",2);
2359 cpl_propertylist_erase_regexp(phead2D,
"NAXIS3",0);
2360 cpl_propertylist_erase_regexp(phead2D,
"CTYPE3",0);
2361 cpl_propertylist_erase_regexp(phead2D,
"CRVAL3",0);
2362 cpl_propertylist_erase_regexp(phead2D,
"CRPIX3",0);
2363 cpl_propertylist_erase_regexp(phead2D,
"CUNIT3",0);
2364 cpl_propertylist_erase_regexp(phead2D,
"CDELT3",0);
2365 cpl_propertylist_erase_regexp(phead2D,
"CD3_*",0);
2366 cpl_propertylist_erase_regexp(phead2D,
"CD1_3",0);
2367 cpl_propertylist_erase_regexp(phead2D,
"CD2_3",0);
2388 if (edge_trim > 0) {
2393 if (edge_trim >= 0.5* sx) {
2395 cpl_msg_warning(cpl_func,
"edge-trim must be smaller than half image "
2396 "size. Reset to 0");
2398 if (edge_trim >= 0.5* sy) {
2400 cpl_msg_warning(cpl_func,
"edge-trim must be smaller than half image "
2401 "size. Reset to 0");
2406 for(cpl_size j = 1; j <= edge_trim; j++) {
2407 for(cpl_size i = 1; i <= sx; i++) {
2412 for(cpl_size j = sy - edge_trim + 1; j <= sy; j++) {
2413 for(cpl_size i = 1; i <= sx; i++) {
2418 for(cpl_size j = 1; j <= sy; j++) {
2419 for(cpl_size i = 1; i <= edge_trim; i++) {
2424 for(cpl_size j = 1; j <= sy; j++) {
2425 for(cpl_size i = sx - edge_trim + 1; i <= sx; i++) {
2431 cpl_msg_warning(cpl_func,
"edge-trim is set to 0");
2434 return cpl_error_get_code();
2463 const cpl_parameterlist * parlist,
2464 const char *input_cube_pro_catg,
2465 const char *filenameSpec,
2468 const char* offunit,
2469 const char* recipe_name,
2470 const char* pipefile_prefix){
2472 cpl_frameset* cube_set = NULL;
2473 cpl_msg_info(cpl_func,
"Combine cubes into a single one plane-by-plane, compute mean, median along wavelength axis");
2477 cpl_msg_warning(cpl_func,
"No %s files", input_cube_pro_catg);
2479 cpl_error_set_message(cpl_func, CPL_ERROR_FILE_NOT_FOUND,
2480 "Missing cube files tagged as %s", input_cube_pro_catg);
2481 return cpl_error_get_code();
2483 cpl_msg_info(cpl_func,
"%02d file(s) of type %s",
2484 (
int)cpl_frameset_get_size(cube_set), input_cube_pro_catg);
2489 const cpl_frame * data_frame = NULL;
2490 hdrl_parameter *aParams_outputgrid = NULL;
2491 cpl_propertylist *phead2D = NULL;
2497 cpl_boolean subtract_bkg = CPL_FALSE;
2500 cpl_size ncubes = cpl_frameset_get_size(cube_set);
2502 cpl_error_set_message(cpl_func, CPL_ERROR_FILE_NOT_FOUND,
"Only one CUBE "
2504 return cpl_error_get_code();
2507 char* pname = cpl_sprintf(
"eris.%s.max-cubes-centres-dist",recipe_name);
2508 int max_cubes_dist = 20;
2510 if (offsetx !=NULL && offsety !=NULL)
2511 max_cubes_dist = INT_MAX ;
2513 max_cubes_dist = cpl_parameter_get_int(
2514 cpl_parameterlist_find_const(parlist, pname));
2517 pname = cpl_sprintf(
"eris.%s.subtract-background",recipe_name);
2518 subtract_bkg = cpl_parameter_get_bool(
2519 cpl_parameterlist_find_const(parlist, pname));
2520 if(subtract_bkg == CPL_TRUE)
2521 cpl_msg_info(cpl_func,
"Median of each plane will be subtracted.");
2524 pname = cpl_sprintf(
"eris.%s.edge-trim",recipe_name);
2525 edge_trim = cpl_parameter_get_int(
2526 cpl_parameterlist_find_const(parlist, pname));
2528 cpl_msg_info(cpl_func,
"%d pixel(s) on edges will be trimmed.", edge_trim);
2534 if ( max_cubes_dist < eris_ifu_compute_max_cubes_dist(cube_set)) {
2535 cpl_msg_info(cpl_func,
"max-cubes-centres-dist: %d",max_cubes_dist);
2536 cpl_msg_warning(cpl_func,
"max distance between cube centres greater than threshold. Skip cube combination.");
2537 return CPL_ERROR_INCOMPATIBLE_INPUT;
2540 hdrl_parameter *aParams_method = NULL;
2542 char *context = NULL;
2543 context = cpl_sprintf(
"eris.%s", recipe_name);
2544 aParams_method = eris_ifu_resample_set_method(parlist, context);
2546 cpl_ensure_code(aParams_method, CPL_ERROR_NULL_INPUT);
2549 cpl_wcs *wcs = NULL;
2550 cpl_wcs **tmp_wcs_array = (cpl_wcs **) cpl_calloc(ncubes,
sizeof(wcs));
2551 char const **name_array = (
char const **) cpl_calloc(ncubes,
sizeof(
const char*));
2553 cpl_wcs *wcs3d = NULL;
2568 const char* ctype3 =
"WAVE";
2570 const cpl_matrix *cd = NULL;
2571 cpl_vector* exptime_vec = cpl_vector_new(ncubes);
2575 for (cpl_size i = 0; i < ncubes ; i++) {
2576 data_frame = cpl_frameset_get_position_const(cube_set, i);
2577 const char *name = cpl_frame_get_filename(data_frame);
2578 cpl_propertylist *xheader_data = cpl_propertylist_load(name,
2580 cpl_propertylist* phead = cpl_propertylist_load(name, 0);
2584 cpl_propertylist_delete(phead);
2585 exptime = dit * ndit;
2586 cpl_vector_set(exptime_vec, i, exptime);
2587 cpl_propertylist *header2D = eris_ifu_resample_get_header2D(xheader_data);
2588 cpl_wcs *tmp_wcs = cpl_wcs_new_from_propertylist(header2D);
2589 cpl_propertylist_delete(header2D);
2591 tmp_wcs_array[i] = tmp_wcs;
2592 name_array[i] = name;
2596 cpl_ensure_code(wcs, CPL_ERROR_NULL_INPUT);
2597 cpl_msg_info(cpl_func,
"WCS used for the output grid definition and passed to "
2598 "hdrl_resample_compute():");
2599 cpl_msg_indent_more();
2600 eris_ifu_resample_wcs_print(wcs);
2601 cpl_msg_indent_less();
2602 wcs3d = cpl_wcs_new_from_propertylist(xheader_data);
2603 mjd_obs = cpl_propertylist_get_double(xheader_data,
"MJD-OBS");
2604 crpix3 = cpl_propertylist_get_double(xheader_data,
"CRPIX3");
2605 crval3 = cpl_propertylist_get_double(xheader_data,
"CRVAL3");
2606 naxis3 = cpl_propertylist_get_int(xheader_data,
"NAXIS3");
2622 cd3_3 = cpl_propertylist_get_double(xheader_data,
"CD3_3");
2625 cd = cpl_wcs_get_cd(wcs);
2626 cd1_1 = cpl_matrix_get(cd, 0, 0);
2627 cd1_2 = cpl_matrix_get(cd, 1, 0);
2628 cd2_1 = cpl_matrix_get(cd, 0, 1);
2629 cd2_2 = cpl_matrix_get(cd, 1, 1);
2631 cpl_ensure_code(wcs3d, CPL_ERROR_NULL_INPUT);
2633 cpl_propertylist_delete(xheader_data);
2637 double as2deg = 1.0/3600.0;
2638 double scalex = 1.0;
2639 double scaley = 1.0;
2642 if (offunit != NULL && offsetx !=NULL && offsety !=NULL){
2643 if (strcmp(offunit,
"PIXEL") == 0){
2645 cpl_msg_info(cpl_func,
"User-defined offset unit in PIXEL");
2646 }
else if (strcmp(offunit,
"DEGREE") == 0){
2648 cpl_msg_info(cpl_func,
"User-defined offset unit is DEGREE");
2649 }
else if (strcmp(offunit,
"ARCSEC") == 0){
2651 cpl_msg_info(cpl_func,
"User-defined offset unit is ARCSEC");
2655 cpl_msg_error(cpl_func,
"User-defined offset unit or offset list is not correct.");
2659 const cpl_array *dims = cpl_wcs_get_image_dims(wcs3d);
2660 int nchan = cpl_array_get_int(dims, 2, NULL);
2661 cpl_size ich_center = 0.5 * nchan;
2662 cpl_msg_info(cpl_func,
"Number of planes: %4.4d", nchan);
2667 cpl_msg_severity level = cpl_msg_get_level();
2668 cpl_msg_set_level(CPL_MSG_INFO);
2669 for (cpl_size ich = 0; ich < nchan; ich++){
2670 cpl_table* restable = NULL;
2672 for (cpl_size i = 0; i < ncubes ; i++) {
2673 const char *name = name_array[i];
2674 cpl_image *im = cpl_image_load(name, CPL_TYPE_DOUBLE, ich, extnum_raw);
2675 cpl_image *err = cpl_image_load(name, CPL_TYPE_DOUBLE, ich, extnum_err);
2676 cpl_image* qual = cpl_image_load(name, CPL_TYPE_INT, ich, extnum_bpm);
2677 cpl_mask *mask = cpl_mask_threshold_image_create(qual, 0, INT_MAX);
2678 cpl_image_delete(qual);
2686 cpl_image* ima_time = NULL;
2687 cpl_image* err_time = NULL;
2688 hdrl_image *himg_exptime = NULL;
2689 if(ich == ich_center) {
2690 ima_time = cpl_image_new(cpl_image_get_size_x(im),
2691 cpl_image_get_size_y(im), CPL_TYPE_DOUBLE);
2692 err_time = cpl_image_new(cpl_image_get_size_x(im),
2693 cpl_image_get_size_y(im), CPL_TYPE_DOUBLE);
2695 exptime = cpl_vector_get(exptime_vec, i);
2696 cpl_image_add_scalar(ima_time, exptime);
2697 cpl_image_add_scalar(err_time, sqrt(exptime));
2702 if (edge_trim > 0) {
2706 cpl_image_delete(ima_time);
2707 cpl_image_delete(err_time);
2711 cpl_image_delete(im);
2712 cpl_image_delete(err);
2714 cpl_mask_delete(mask);
2717 if(subtract_bkg == CPL_TRUE) {
2722 if (edge_trim > 0) {
2726 cpl_table *tab_tmp = NULL;
2727 cpl_table *tab_exptime_tmp = NULL;
2728 if (offsetx !=NULL && offsety !=NULL){
2731 cpl_table_subtract_scalar(tab_tmp,
"ra", offsetx[i]*cd1_1+ offsety[i]*cd1_2);
2732 cpl_table_subtract_scalar(tab_tmp,
"dec", offsetx[i]*cd2_1+ offsety[i]*cd2_2);
2734 cpl_msg_info(__func__,
"User-defined offest in RA, DEC (ARCSEC) %g, %g",
2735 (offsetx[i]*cd1_1+ offsety[i]*cd1_2)*3600,
2736 (offsetx[i]*cd2_1+ offsety[i]*cd2_2)*3600);
2738 if(ich == ich_center) {
2741 cpl_table_subtract_scalar(tab_exptime_tmp,
"ra", offsetx[i]*cd1_1+ offsety[i]*cd1_2);
2742 cpl_table_subtract_scalar(tab_exptime_tmp,
"dec", offsetx[i]*cd2_1+ offsety[i]*cd2_2);
2744 }
else if (offcode==1 || offcode==2 ){
2746 cpl_table_subtract_scalar(tab_tmp,
"ra", offsetx[i]*scalex);
2747 cpl_table_subtract_scalar(tab_tmp,
"dec", offsety[i]*scaley);
2749 cpl_msg_info(__func__,
"User-defined offest in RA, DEC (%s) %g, %g",
2754 if(ich == ich_center) {
2756 cpl_table_subtract_scalar(tab_exptime_tmp,
"ra", offsetx[i]*scalex);
2757 cpl_table_subtract_scalar(tab_exptime_tmp,
"dec", offsety[i]*scaley);
2762 if(ich == ich_center) {
2767 eris_ifu_resample_update_table(tab_tmp, &restable);
2769 if(ich == ich_center) {
2773 hdrl_resample_result *himl_exptime_tmp = NULL;
2775 aParams_method, aParams_outputgrid, wcs);
2786 cpl_table_delete(tab_exptime_tmp);
2790 cpl_table_delete(tab_tmp);
2792 cpl_msg_debug(__func__,
"Created table for channel %d cube %d", (
int)ich, (
int)i);
2794 cpl_ensure_code(restable, CPL_ERROR_NULL_INPUT);
2797 cpl_msg_debug(__func__,
"Table of the chanel %d size %lld x %lld", (
int)ich,
2798 cpl_table_get_nrow(restable), cpl_table_get_ncol(restable));
2801 aParams_outputgrid = eris_ifu_resample_set_outputgrid2D(recipe_name, parlist,
2803 if ( cpl_error_get_code() != CPL_ERROR_NONE) {
2805 return CPL_ERROR_NONE;
2807 cpl_ensure_code(aParams_outputgrid, CPL_ERROR_INCOMPATIBLE_INPUT);
2811 cpl_ensure_code(aParams_outputgrid, CPL_ERROR_NULL_INPUT);
2813 hdrl_resample_result *result = NULL;
2825 cpl_ensure_code(result, CPL_ERROR_NULL_INPUT);
2826 cpl_table_delete(restable);
2830 cpl_wcs *res_wcs = cpl_wcs_new_from_propertylist(result->header);
2831 cpl_msg_info(cpl_func,
"Final WCS after resampling: ");
2832 cpl_msg_indent_more();
2833 eris_ifu_resample_wcs_print(res_wcs);
2834 cpl_msg_indent_less();
2835 cpl_wcs_delete(res_wcs);
2836 phead2D = cpl_propertylist_duplicate(result->header);
2844 for (
int i = 0; i < ncubes; i++){
2845 cpl_wcs_delete(tmp_wcs_array[i]);
2847 cpl_free(tmp_wcs_array);
2848 cpl_free(name_array);
2852 cpl_msg_info(__func__,
"Coadding cubes: range [%4.4d,%4.4d] of %d",
2853 (
int)ich,(
int)ich+100, nchan);
2855 cpl_msg_set_level(level);
2857 hdrl_image* hima_exptime_mean = NULL;
2858 cpl_image* contrib_exptime_map = NULL;
2860 &contrib_exptime_map);
2867 cpl_image_multiply_scalar(contrib_exptime_map, cpl_vector_get_median(exptime_vec));
2869 cpl_ensure_code(res_cube, CPL_ERROR_NULL_INPUT);
2870 cpl_vector_delete(exptime_vec);
2878 hdrl_resample_result *result = cpl_calloc(1,
sizeof(hdrl_resample_result));
2879 result->himlist = res_cube;
2880 result->header = phead2D;
2883 cpl_propertylist_append_double(result->header,
"CRPIX3", crpix3);
2884 cpl_propertylist_append_double(result->header,
"CRVAL3", crval3);
2885 cpl_propertylist_append_int(result->header,
"NAXIS3", naxis3);
2887 cpl_propertylist_append_string(result->header,
"CTYPE3", ctype3);
2888 cpl_propertylist_append_double(result->header,
"CD1_3", cd1_3);
2889 cpl_propertylist_append_double(result->header,
"CD2_3", cd2_3);
2890 cpl_propertylist_append_double(result->header,
"CD3_3", cd3_3);
2891 cpl_propertylist_append_double(result->header,
"CD3_1", cd3_1);
2892 cpl_propertylist_append_double(result->header,
"CD3_2", cd3_2);
2893 cpl_propertylist_append_double(result->header,
"MJD-OBS", mjd_obs);
2894 cpl_propertylist* wcs2D = eris_ifu_plist_extract_wcs2D(result->header);
2896 hdrl_image* cube_collapsed = NULL;
2897 cpl_image* cube_cmap = NULL;
2905 char* pro_catg = NULL;
2906 if(strstr(pipefile_prefix,
"no_flat") != NULL) {
2907 pro_catg = cpl_sprintf(
"%s_%s",input_cube_pro_catg,
"COADD_NOFLAT");
2909 pro_catg = cpl_sprintf(
"%s_%s",input_cube_pro_catg,
"COADD");
2911 char* fname = cpl_sprintf(
"%s_%s_coadd.fits", pipefile_prefix, filenameSpec);
2913 if(strcmp(recipe_name,
"eris_ifu_jitter") == 0) {
2915 frameset, CPL_FALSE);
2918 frameset, CPL_FALSE);
2925 pro_catg = cpl_sprintf(
"%s_%s",input_cube_pro_catg,
"MEAN");
2926 fname = cpl_sprintf(
"%s_%s_mean.fits", pipefile_prefix, filenameSpec);
2927 cpl_propertylist_append_string(phead2D,
"PRODCATG", PRODCATG_WHITELIGHT);
2930 cpl_propertylist_update_string(phead2D, CPL_DFS_PRO_CATG, pro_catg);
2933 recipe_name, phead2D,
"RADECSYS", fname,
2936 rmse, mask_ima, flag16bit, UNIT_ADU);
2938 fname = cpl_sprintf(
"%s_%s_exposure_map.fits", pipefile_prefix, filenameSpec);
2940 cpl_propertylist_update_string(phead2D, CPL_DFS_PRO_CATG, ERIS_IFU_PRO_JITTER_EXPOSURE_MAP);
2941 cpl_propertylist_update_string(phead2D,
"PRODCATG", PRODCATG_EXPOSUREMAP);
2943 cpl_propertylist_append(phead2D, wcs2D);
2944 cpl_propertylist_delete(wcs2D);
2946 recipe_name, phead2D,
"RADECSYS", fname, contrib_exptime_map,
2947 CPL_TYPE_DOUBLE, UNIT_TIME);
2950 cpl_image_delete(mask_ima);
2951 cpl_image_delete(contrib_exptime_map);
2955 cpl_image_delete(cube_cmap);
2956 cpl_wcs_delete(wcs3d);
2959 cpl_frameset_delete(cube_set);
2962 return cpl_error_get_code();
ifsPreopticsScale eris_ifu_get_preopticsScale(cpl_propertylist *header)
Return the the pre-optics scaling.
cpl_error_code eris_ifu_heades_add_hduclass_qual(cpl_propertylist *plist, deqQualityType qualityType)
Add quality extension classification properties to FITS header.
cpl_frameset * eris_ifu_extract_frameset(const cpl_frameset *in, const char *tag)
Extract frames with a specific tag from a frameset.
cpl_error_code eris_ifu_save_image_phase3(cpl_frameset *allframes, const cpl_parameterlist *parlist, const cpl_frameset *usedframes, const char *recipe, const cpl_propertylist *applist, const char *remregexp, const char *filename, const cpl_image *image, cpl_type type, const char *unit)
Save a DFS-compliant Phase 3 image product.
cpl_error_code eris_ifu_heades_add_hduclass_common(cpl_propertylist *plist, const char *deq_hduclas)
Add common HDU classification properties to FITS header.
cpl_error_code eris_ifu_save_deq_image(cpl_frameset *allframes, cpl_propertylist *header, const cpl_parameterlist *parlist, const cpl_frameset *usedframes, const cpl_frame *inherit, const char *recipe, const cpl_propertylist *applist, const char *remregexp, const char *filename, const cpl_image *image, const cpl_image *error, deqErrorType errorType, const cpl_image *dataQualityMask, deqQualityType dataQualityType, const char *unit)
Save a DFS-compliant image product with data, error, and quality extensions.
cpl_error_code eris_ifu_heades_add_hduclass_data(cpl_propertylist *plist)
Add data extension classification properties to FITS header.
cpl_error_code eris_ifu_heades_add_hduclass_errs(cpl_propertylist *plist, deqErrorType errorType)
Add error extension classification properties to FITS header.
eris_ifu_sdp_properties * eris_ifu_sdp_properties_collect(hdrl_resample_result *aCube, cpl_frameset *set, const cpl_parameterlist *parlist, const char *recipe_name)
Collect all SDP metadata from cube, frameset, and parameters.
void eris_ifu_sdp_properties_delete(eris_ifu_sdp_properties *aProperties)
Free SDP properties structure and all contained data.
cpl_error_code eris_ifu_sdp_properties_update(cpl_propertylist *aHeader, const eris_ifu_sdp_properties *aProperties)
Update FITS header with ESO Science Data Product keywords.
cpl_error_code eris_ifu_jitter_get_procatg_and_filename(cubeType type, char **proCatg, char **filenamePrefix)
Get the value of the PRO.CATG and the filename as function of the type of cube.
cpl_error_code eris_ifu_combine(cubeType obj_type, cpl_frameset *frameset, const cpl_parameterlist *parlist, const char *recipe_name, const char *pipefile_prefix)
Resample and combine multiple IFU cubes into a single cube.
cpl_error_code eris_ifu_combine_pbp(cpl_frameset *frameset, const cpl_parameterlist *parlist, const char *input_cube_pro_catg, const char *filenameSpec, float *offsetx, float *offsety, const char *offunit, const char *recipe_name, const char *pipefile_prefix)
Resample and combine cubes plane-by-plane (2D spatial resampling per wavelength)
cpl_error_code eris_ifu_cube_collapse_mean_and_save(const char *pro_catg, cpl_frameset *frameset, const cpl_parameterlist *parlist, const char *recipe_name, cpl_boolean apply_flat, cpl_boolean is_pupil)
-------------------------------------------------------------------------—/
cpl_error_code eris_ifu_resample_save_cube(hdrl_resample_result *aCube, const char *procatg, const char *recipe, const char *filename, const cpl_parameterlist *parlist, cpl_frameset *frameset, cpl_boolean gen_phase3)
Save resampled cube to FITS file in DEQ (Data-Error-Quality) format.
cpl_error_code eris_ifu_resample_trim_edge(hdrl_image *himg, int edge_trim)
Trim edge pixels from hdrl_image by rejecting them.
void eris_ifu_free_propertylist(cpl_propertylist **item)
Free memory and set pointer to null.
void eris_ifu_free_string(char **item)
Free memory and set pointer to null.
cpl_table * eris_qclog_init(void)
Initialize QC table.
cpl_error_code eris_pfits_put_qc(cpl_propertylist *plist, cpl_table *qclog)
convert table with QC parameter information to a propertylist
cpl_error_code eris_qclog_add_double_f(cpl_table *table, const char *key_name, const double value, const char *key_help)
add QC float info to table
cpl_error_code eris_ifu_split3_hdrl_imagelist(hdrl_imagelist *cube, cpl_imagelist *dataCube, cpl_imagelist *errorCube, cpl_imagelist *qualCube)
extract from hdrl_imagelist the data and error information
cpl_error_code eris_ifu_save_image(cpl_frameset *fs, const cpl_propertylist *plist, const cpl_parameterlist *parlist, const char *recipe, const char *procatg, const char *filename, cpl_type type, const cpl_image *image)
Save image with DFS compliance.
int eris_pfits_get_ndit(const cpl_propertylist *plist)
find out the DIT value
double eris_pfits_get_dit(const cpl_propertylist *plist)
find out the DIT value
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
cpl_error_code eris_check_error_code(const char *func_id)
handle CPL errors
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)
cpl_error_code hdrl_image_reject_from_mask(hdrl_image *self, const cpl_mask *map)
set bpm of hdrl_image
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
cpl_mask * hdrl_image_get_mask(hdrl_image *himg)
get cpl bad pixel mask from image
cpl_image * hdrl_image_get_error(hdrl_image *himg)
get error as cpl image
cpl_size hdrl_image_get_size_y(const hdrl_image *self)
return size of Y dimension of image
cpl_size hdrl_image_get_size_x(const hdrl_image *self)
return size of X dimension of image
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
cpl_image * hdrl_image_get_image(hdrl_image *himg)
get data as cpl image
cpl_error_code hdrl_image_reject(hdrl_image *self, cpl_size xpos, cpl_size ypos)
mark pixel as bad
void hdrl_image_delete(hdrl_image *himg)
delete hdrl_image
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().