35#include "eris_nix_utils.h"
36#include "eris_utils.h"
37#include "eris_nix_casu_utils.h"
38#include "eris_nix_dfs.h"
39#include "eris_nix_match.h"
42#include <casu_stats.h>
100 const int read_offsets,
101 const cpl_table * refine_wcs,
102 const master_dark * mdark,
103 const gain_linearity * gain_lin,
104 const master_flat * flatfield_1,
105 const master_flat * flatfield_2,
106 const master_bpm * mbad_pix_map,
108 const char * fill_rejected,
109 const double fill_value,
110 const cpl_size x_probe,
111 const cpl_size y_probe) {
113 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
115 if (limage == NULL)
return CPL_ERROR_NONE;
124 if (limage->himage) {
141 enu_check_error_code(
"failure in basic calibration");
143 }
else if (limage->himagelist) {
145 cpl_error_code didfail = CPL_ERROR_NONE;
148 cpl_errorstate cleanstate = cpl_errorstate_get();
150#pragma omp parallel for
152 for (cpl_size plane = 0; plane < nplanes; plane++) {
154 cpl_error_code errori = cpl_error_get_code();
162 if (didfail)
continue;
166 cpl_msg_info(cpl_func,
"basic calibrate cube: plane %d",
189 errori = cpl_error_set_message(cpl_func, cpl_error_get_code(),
190 "Cannot process plane %d/%d",
191 (
int)plane+1, (
int)nplanes);
200 cpl_errorstate_dump(cleanstate, CPL_FALSE, NULL);
201 cpl_errorstate_set(cleanstate);
203#pragma omp critical(enu_basic_calibrate)
208 if (didfail) (void)cpl_error_set_message(cpl_func, didfail,
209 "Cannot process %d plane(s)",
215 return cpl_error_get_code();
258 cpl_image * confidence,
259 cpl_propertylist * plist,
260 const cpl_frame * frame,
261 const int read_offsets,
262 const master_dark * mdark,
263 const gain_linearity * gain_lin,
264 const master_flat * flatfield_1,
265 const master_flat * flatfield_2,
266 const master_bpm * mbad_pix_map,
267 const int set_confidence,
269 const char * fill_rejected,
270 const double fill_value,
271 const cpl_size x_probe,
272 const cpl_size y_probe) {
274 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
276 cpl_ensure_code(himage, CPL_ERROR_NULL_INPUT);
285 cpl_size nx_chip = 0;
286 cpl_size ny_chip = 0;
287 cpl_boolean windowed = 0;
297 enu_check_error_code(
"failed to read detector window information");
299 cpl_msg_info(cpl_func,
"detector is windowed!");
304 const char * instrume = cpl_propertylist_get_string(plist,
"INSTRUME");
305 const cpl_boolean eris = !strcmp(instrume,
"ERIS");
309 int probe = CPL_FALSE;
310 if ((x_probe >= 0) && (y_probe >= 0)) {
311 if (x_probe < nx && y_probe < ny) {
322 cpl_msg_info(cpl_func,
323 "..basic_calibrate probe (%d,%d) entry(v,e,r)=%5.3e %5.3e %d",
324 (
int) x_probe, (
int)y_probe, val.data, val.error,
327 cpl_msg_info(cpl_func,
328 "..basic_calibrate probe (median) entry(v,e)=%5.3e %5.3e",
329 val.data, val.error);
337 cpl_msg_info(cpl_func,
"basic calibrate: correcting wcs for windowed data");
338 const double crpix1 = cpl_propertylist_get_double(plist,
"CRPIX1");
339 cpl_propertylist_update_double(plist,
"CRPIX1",
340 crpix1 - (
double)strx + 1.0);
341 const double crpix2 = cpl_propertylist_get_double(plist,
"CRPIX2");
342 cpl_propertylist_update_double(plist,
"CRPIX2",
343 crpix2 - (
double)stry + 1.0);
344 cpl_msg_info(cpl_func,
"..crpix %f %f -> %f %f", crpix1, crpix2,
345 crpix1 - (
double)strx + 1.0,
346 crpix2 - (
double)stry + 1.0);
351 if (read_offsets && eris) {
352 cpl_msg_info(cpl_func,
"basic calibrate: removing read-offsets and "
355 enu_check_error_code(
"error removing read-offsets");
363 cpl_msg_info(cpl_func,
364 "..basic_calibrate probe (%d,%d) after "
365 "read_offsets(v,e,r)=%5.3e %5.3e %d",
366 (
int) x_probe, (
int)y_probe,
367 val.data, val.error, reject);
369 cpl_msg_info(cpl_func,
370 "..basic_calibrate probe (median) after "
371 "read_offsets(v,e)=%5.3e %5.3e", val.data, val.error);
376 cpl_msg_info(cpl_func,
"basic calibrate: subtracting master dark");
384 "^ESO DET SEQ1 DIT$|"
385 "^ESO DET READ CURNAME$|"
386 "^ESO DET SEQ1 WIN NX$|"
387 "^ESO DET SEQ1 WIN NY$|"
388 "^ESO DET SEQ1 WIN ROT$|"
389 "^ESO DET SEQ1 WIN STRX$|"
390 "^ESO DET SEQ1 WIN STRY$|"
392 "^ESO DET NCORRS NAME$|"
395 "^ESO DET WIN STARTX$|"
396 "^ESO DET WIN STARTY$");
397 const char * filename = NULL;
399 filename = cpl_frame_get_filename(frame);
401 if (filename == NULL) {
402 enu_check(compatible, CPL_ERROR_INCOMPATIBLE_INPUT,
403 "MASTER_DARK does not match frame");
405 enu_check(compatible, CPL_ERROR_INCOMPATIBLE_INPUT,
406 "MASTER_DARK does not match %s", filename);
412 double gain = cpl_propertylist_get_double(mdark->plist,
414 cpl_propertylist_update_double(plist,
"ESO DARK GAIN", gain);
415 cpl_propertylist_set_comment(plist,
"ESO DARK GAIN",
416 "GAIN from master_dark [e-/ADU]");
417 enu_check_error_code(
"error subtracting DARK");
430 cpl_msg_info(cpl_func,
431 "..basic_calibrate probe (%d,%d) after dark_subtract "
432 "dark(v,e,r)=%5.3e %5.3e %d "
433 "data(v,e,r)=%5.3e %5.3e %d",
434 (
int) x_probe, (
int)y_probe,
435 dark.data, dark.error, reject2,
436 val.data, val.error, reject);
439 cpl_msg_info(cpl_func,
440 "..basic_calibrate probe (median) after dark_subtract "
441 "dark(v,e)=%5.3e %5.3e data(v,e)=%5.3e %5.3e",
442 dark.data, dark.error,
443 val.data, val.error);
448 hdrl_image * himage_copy = NULL;
452 cpl_msg_info(cpl_func,
"basic calibrate: linearizing and estimating "
455 himage_copy = end_linearize_and_variance_detmon(gain_lin,
463 enu_check_error_code(
"error estimating variance");
468 cpl_propertylist_update_double(plist,
469 "ESO DETMON SATURATION",
470 gain_lin->saturation_limit);
471 cpl_propertylist_set_comment(plist,
472 "ESO DETMON SATURATION",
473 "saturation level from "
474 "linearisation [ADU]");
479 cpl_propertylist_update_string(plist,
"BUNIT",
"adu/s");
487 cpl_msg_info(cpl_func,
488 "..basic_calibrate probe (%d,%d) "
489 "after linearize/calc_variance "
490 "data(v,e,r)=%5.3e %5.3e %d",
491 (
int) x_probe, (
int)y_probe, val.data, val.error,
494 cpl_msg_info(cpl_func,
495 "..basic_calibrate probe (median) "
496 "after linearize/calc_variance "
497 "data(v,e)=%5.3e %5.3e", val.data, val.error);
502 cpl_msg_info(cpl_func,
"basic calibrate: dividing by flatfield 1");
509 hdrl_image * flat_window = !windowed ? NULL :
514 const hdrl_image * flat_use = !windowed ? flatfield_1->flat
517 enu_check_error_code(
"error extracting window from flatfield 1");
525 if (set_confidence) {
532 cpl_image * flat_conf_window = !windowed ? NULL :
533 cpl_image_extract(flatfield_1->confidence,
537 const cpl_image * flat_conf_use = !windowed ? flatfield_1->confidence
540 enu_check_error_code(
"error extracting window from flatfield 1 "
542 cpl_image_copy(confidence, flat_conf_use, 1, 1);
543 cpl_image_delete(flat_conf_window);
547 cpl_mask * old_mask =
548 cpl_image_set_bpm(confidence,
550 cpl_mask_delete(old_mask);
551 cpl_image_fill_rejected(confidence, 0.0);
552 (void)cpl_image_unset_bpm(confidence);
554 cpl_image_accept_all(confidence);
558 enu_check_error_code(
"error dividing by flatfield 1");
571 cpl_msg_info(cpl_func,
572 "..basic_calibrate probe (%d,%d) after flatfield1 "
573 "flat(v,e,r)=%5.3e %5.3e %d "
574 "data(v,e,r)=%5.3e %5.3e %d",
575 (
int) x_probe, (
int)y_probe,
576 flat.data, flat.error, reject2,
577 val.data, val.error, reject);
580 cpl_msg_info(cpl_func,
581 "..basic_calibrate probe (median) after flatfield1 "
582 "flat(v,e)=%5.3e %5.3e data(v,e)=%5.3e %5.3e",
583 flat.data, flat.error,
584 val.data, val.error);
591 cpl_msg_info(cpl_func,
"basic calibrate: dividing by flatfield 2");
593 hdrl_image * flat_window = !windowed ? NULL :
598 const hdrl_image * flat_use = !windowed ? flatfield_2->flat
601 enu_check_error_code(
"error extracting window from flatfield 2");
605 if (set_confidence) {
606 cpl_image * flat_conf_window = !windowed ? NULL :
607 cpl_image_extract(flatfield_2->confidence,
611 const cpl_image * flat_conf_use = !windowed ? flatfield_2->confidence
614 enu_check_error_code(
"error extracting window from flatfield 1 "
617 cpl_image_copy(confidence, flat_conf_use, 1, 1);
618 cpl_image_delete(flat_conf_window);
622 cpl_mask * old_mask =
623 cpl_image_set_bpm(confidence,
625 cpl_mask_delete(old_mask);
626 cpl_image_fill_rejected(confidence, 0.0);
627 (void)cpl_image_unset_bpm(confidence);
629 cpl_image_accept_all(confidence);
633 enu_check_error_code(
"error dividing by flatfield 2");
646 cpl_msg_info(cpl_func,
647 "..basic_calibrate probe (%d,%d) after flatfield2 "
648 "flat(v,e,r)=%5.3e %5.3e %d "
649 "data(v,e,r)=%5.3e %5.3e %d",
650 (
int) x_probe, (
int)y_probe,
651 flat.data, flat.error, reject2,
652 val.data, val.error, reject);
655 cpl_msg_info(cpl_func,
656 "..basic_calibrate probe (median) after flatfield2 "
657 "flat(v,e)=%5.3e %5.3e data(v,e)=%5.3e %5.3e",
658 flat.data, flat.error,
659 val.data, val.error);
665 cpl_msg_info(cpl_func,
"basic calibrate: associating with master_bpm");
666 cpl_mask * new_mask = en_master_bpm_get_mask(mbad_pix_map, flag_mask);
670 cpl_mask * new_mask_window = !windowed ? NULL :
671 cpl_mask_extract(new_mask,
675 const cpl_mask * new_mask_use = !windowed ? new_mask : new_mask_window;
677 enu_check_error_code(
"error extracting window from bpm mask");
682 cpl_mask_or(current_mask, new_mask_use);
683 enu_check_error_code(
"error setting new mask");
685 if (set_confidence) {
689 double * conf_data = cpl_image_get_data(confidence);
690 cpl_binary * mask_data = cpl_mask_get_data(current_mask);
691 cpl_size nx_mask = cpl_mask_get_size_x(current_mask);
692 cpl_size ny_mask = cpl_mask_get_size_y(current_mask);
693 for (cpl_size i = 0; i < nx_mask * ny_mask; i++) {
694 if (mask_data[i] == CPL_BINARY_1) {
707 double conf = cpl_image_get(confidence, x_probe, y_probe, &ignore);
708 cpl_msg_info(cpl_func,
"..basic_calibrate probe (%d,%d) "
709 "after associate_mask(v,e,r,c)=%5.3e %5.3e %d %5.3e",
710 (
int) x_probe, (
int)y_probe, val.data, val.error,
713 cpl_mask_delete(new_mask);
714 cpl_mask_delete(new_mask_window);
717 if (!strcmp(fill_rejected,
"noop")) {
718 cpl_msg_info(cpl_func,
"not modifying rejected pixels");
721 if (!strcmp(fill_rejected,
"set_NaN")) {
723 cpl_msg_info(cpl_func,
"basic calibrate: setting rejected pixels "
725 }
else if (!strcmp(fill_rejected,
"set_value")) {
727 cpl_msg_info(cpl_func,
"basic calibrate: setting rejected pixels "
730 enu_check(CPL_FALSE, CPL_ERROR_UNSPECIFIED,
731 "bad fill-rejected value: programming error");
735 enu_check_error_code(
"error setting rejected pixels to %f",
745 double conf = cpl_image_get(confidence, x_probe, y_probe, &ignore);
746 cpl_msg_info(cpl_func,
747 "..basic_calibrate probe (%d,%d) after "
748 "fill-rejected(v,e,r,c)=%5.3e %5.3e %d %5.3e",
749 (
int) x_probe, (
int)y_probe, val.data, val.error,
756 return cpl_error_get_code();;
774 const int min_coadds,
775 const hdrl_parameter * collapse_params,
776 const cpl_size filter_size_x,
777 const cpl_size filter_size_y,
778 const hdrl_flat_method method) {
782 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
783 cpl_ensure(himlist, CPL_ERROR_NULL_INPUT, NULL);
784 cpl_ensure(collapse_params, CPL_ERROR_NULL_INPUT, NULL);
786 cpl_image * contribution = NULL;
787 hdrl_parameter * flat_params = NULL;
788 hdrl_image * result = NULL;
793 enu_check(size >= min_coadds, CPL_ERROR_ILLEGAL_INPUT,
794 "insufficient frames for flat (%d<min:%d)",
795 (
int) size, min_coadds);
802 if (method == HDRL_FLAT_FREQ_HIGH) {
803 cpl_msg_info(cpl_func,
"Calculating high freq flatfield");
804 }
else if (method == HDRL_FLAT_FREQ_LOW) {
805 cpl_msg_info(cpl_func,
"Calculating low freq flatfield");
807 cpl_msg_info(cpl_func,
"..smoothing filter sizes %d %d",
808 (
int)filter_size_x, (
int)filter_size_y);
810 flat_params, &result, &contribution);
811 enu_check_error_code(
"error computing flat-field");
814 cpl_image_delete(contribution);
816 if (cpl_error_get_code() != CPL_ERROR_NONE) {
840 const double fwhm_pix,
843 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
844 cpl_ensure_code(limage, CPL_ERROR_NULL_INPUT);
845 cpl_ensure_code(abmaglim, CPL_ERROR_NULL_INPUT);
848 if (fwhm_pix <= 0.0) {
850 return CPL_ERROR_NONE;
856 double histo_min = 10.;
857 double histo_max = 1.;
860 double bin_size = 0.;
861 cpl_size error_niter = 0;
862 hdrl_mode_type mode_method = HDRL_MODE_MEDIAN;
863 hdrl_parameter * mode_parameter = NULL;
871 double* pconf = cpl_image_get_data_double(limage->confidence);
876 cpl_binary* pbpm_data = cpl_mask_get_data(mask_data);
877 cpl_binary* pbpm_errs = cpl_mask_get_data(mask_errs);
878 for(cpl_size indexj = 0; indexj < sy; indexj++) {
879 for(cpl_size indexi = 0; indexi < sx; indexi++) {
880 if(!isfinite(pdata[indexi+indexj*sx]) ) {
881 pbpm_data[indexi+indexj*sx] = CPL_BINARY_1;
882 pbpm_errs[indexi+indexj*sx] = CPL_BINARY_1;
883 pconf[indexi+indexj*sx] = CPL_BINARY_1;
893 ((
int) (3.0 * fwhm_pix) / 2) * 2 + 1,
894 ((
int) (3.0 * fwhm_pix) / 2) * 2 + 1,
895 HDRL_IMAGE_EXTEND_MIRROR,
901 return cpl_error_get_code();
919 const cpl_propertylist * wcs_plist) {
921 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
922 cpl_ensure_code(catalogue, CPL_ERROR_NULL_INPUT);
923 cpl_ensure_code(wcs_plist, CPL_ERROR_NULL_INPUT);
925 cpl_array * conv_status = NULL;
926 cpl_matrix * ref_phys = NULL;
927 cpl_matrix * ref_world = NULL;
928 cpl_wcs * wcs = NULL;
940 cpl_size nobj = cpl_table_get_nrow(catalogue);
941 const double * ra = cpl_table_get_data_double_const(catalogue,
"RA");
942 const double * dec = NULL;
946 if (cpl_table_has_column(catalogue,
"DEC")) {
947 dec = cpl_table_get_data_double_const(catalogue,
"DEC");
949 ref_world = cpl_matrix_new(nobj, 2);
950 for (cpl_size iobj = 0; iobj < nobj; iobj++) {
951 cpl_matrix_set(ref_world, iobj, 0, ra[iobj]);
952 cpl_matrix_set(ref_world, iobj, 1, dec[iobj]);
955 wcs = cpl_wcs_new_from_propertylist(wcs_plist);
957 cpl_wcs_convert(wcs, ref_world, &ref_phys, &conv_status,
959 enu_check_error_code(
"failure calculating predicted positions");
964 if (!cpl_table_has_column(catalogue,
"X_coordinate")) {
965 cpl_table_new_column(catalogue,
"X_coordinate", CPL_TYPE_DOUBLE);
967 if (!cpl_table_has_column(catalogue,
"Y_coordinate")) {
968 cpl_table_new_column(catalogue,
"Y_coordinate", CPL_TYPE_DOUBLE);
971 for (cpl_size iobj = 0; iobj < nobj; iobj++) {
972 cpl_table_set_double(catalogue,
"X_coordinate", iobj,
973 cpl_matrix_get(ref_phys, iobj, 0));
974 cpl_table_set_double(catalogue,
"Y_coordinate", iobj,
975 cpl_matrix_get(ref_phys, iobj, 1));
985 cpl_array_delete(conv_status);
986 cpl_matrix_delete(ref_phys);
987 cpl_matrix_delete(ref_world);
990 return cpl_error_get_code();
1012 const hdrl_image * himage,
1013 const cpl_image * confidence,
1014 const cpl_wcs * wcs,
1015 hdrl_parameter * params) {
1017 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
1019 cpl_ensure(himage, CPL_ERROR_NULL_INPUT, NULL);
1023 cpl_ensure(wcs, CPL_ERROR_NULL_INPUT, NULL);
1024 cpl_ensure(params, CPL_ERROR_NULL_INPUT, NULL);
1026 hdrl_catalogue_result * result = NULL;
1034 cpl_ensure(cpl_image_get_max(confidence) > 0.0, CPL_ERROR_ILLEGAL_INPUT,
1039 confidence, wcs, params);
1057 if (cpl_error_get_code() != CPL_ERROR_NONE) {
1067 const char * message = cpl_error_get_message();
1068 cpl_msg_info(cpl_func,
"catalogue %s", message);
1070 if (strstr(message,
"No objects found in image")) {
1071 cpl_msg_warning(cpl_func,
"No objects found in image");
1092 const double timerange,
1093 const located_imagelist * pool,
1096 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
1098 cpl_ensure(target, CPL_ERROR_NULL_INPUT, NULL);
1099 cpl_ensure(timerange > 0.0, CPL_ERROR_ILLEGAL_INPUT, NULL);
1100 cpl_ensure(pool, CPL_ERROR_NULL_INPUT, NULL);
1102 cpl_vector * result = NULL;
1103 cpl_vector * result_diff = NULL;
1105 double target_mjd = cpl_propertylist_get_double(target->plist,
"MJD-OBS");
1107 for (cpl_size j = 0; j < pool->size; j++) {
1108 double mjd_j = cpl_propertylist_get_double((pool->limages[j])->
1110 double diff = (target_mjd - mjd_j) * 3600.0 * 24.0;
1112 if (fabs(diff) < timerange/2.0) {
1113 if (result == NULL) {
1114 result = cpl_vector_new(1);
1115 result_diff = cpl_vector_new(1);
1116 cpl_vector_set(result, 0, (
double)j);
1117 cpl_vector_set(result_diff, 0, diff);
1119 cpl_size size = cpl_vector_get_size(result);
1120 cpl_vector_set_size(result, size + 1);
1121 cpl_vector_set_size(result_diff, size + 1);
1122 cpl_vector_set(result, size, (
double)j);
1123 cpl_vector_set(result_diff, size, diff);
1127 enu_check_error_code(
"error in enu_bracket_skys: %s",
1128 cpl_error_get_message());
1129 enu_check(result && cpl_vector_get_size(result) > 0,
1130 CPL_ERROR_INCOMPATIBLE_INPUT,
1131 "no sky frames found within %f sec of target frame", timerange/2.0);
1134 cpl_size size = cpl_vector_get_size(result);
1136 cpl_msg_info(
" ",
".. enu_bracket_skys: timerange %5.1f", timerange);
1137 cpl_msg_info(
" ",
".. j# diff(sec)");
1139 for (cpl_size i = 0; i < size; i++) {
1140 int j = cpl_vector_get(result, i);
1141 double diff = cpl_vector_get(result_diff, i);
1142 cpl_msg_info(
" ",
".. %3d %7.1f", j, diff);
1147 if (cpl_error_get_code() != CPL_ERROR_NONE) {
1148 cpl_vector_delete(result);
1151 cpl_vector_delete(result_diff);
1168 hdrl_parameter * params) {
1170 cpl_ensure_code(cpl_error_get_code() == CPL_ERROR_NONE,
1171 cpl_error_get_code());
1172 cpl_ensure_code(limlist, CPL_ERROR_NULL_INPUT);
1173 cpl_ensure_code(params, CPL_ERROR_NULL_INPUT);
1175 cpl_size njitters = limlist->size;
1176 cpl_msg_info(cpl_func,
"cataloguing %d images", (
int) njitters);
1178 for (cpl_size i = 0; i < njitters; i++) {
1179 cpl_wcs * wcs = cpl_wcs_new_from_propertylist
1180 (limlist->limages[i]->plist);
1181 cpl_msg_info(cpl_func,
"cataloguing image %d", (
int) i);
1183 double* pdata = cpl_image_get_data_double(
hdrl_image_get_image((limlist->limages[i])->himage));
1189 double* pconf = cpl_image_get_data_double(limlist->limages[i]->confidence);
1190 cpl_binary* pbpm_data = cpl_mask_get_data(mask);
1191 cpl_binary* pbpm_errs = cpl_mask_get_data(mskerr);
1192 for(cpl_size indexj = 0; indexj < sy; indexj++) {
1193 for(cpl_size indexi = 0; indexi < sx; indexi++) {
1194 if(!isfinite(pdata[indexi+indexj*sx]) ) {
1195 pbpm_data[indexi+indexj*sx] = CPL_BINARY_1;
1196 pbpm_errs[indexi+indexj*sx] = CPL_BINARY_1;
1197 pconf[indexi+indexj*sx] = 0;
1204 limlist->limages[i]->himage,
1205 limlist->limages[i]->confidence,
1208 cpl_wcs_delete(wcs);
1210 if (cpl_error_get_code() != CPL_ERROR_NONE)
break;
1213 return cpl_error_get_code();
1235 const cpl_propertylist * plist2,
1236 const char * regexp) {
1238 if (cpl_error_get_code() != CPL_ERROR_NONE)
return 0;
1242 cpl_ensure(plist1 && plist2 && regexp, CPL_ERROR_NULL_INPUT, 0);
1244 int result = CPL_TRUE;
1245 char * reason = NULL;
1249 cpl_propertylist * subset1 = cpl_propertylist_new();
1250 cpl_propertylist_copy_property_regexp(subset1, plist1, regexp, 0);
1251 cpl_propertylist * subset2 = cpl_propertylist_new();
1252 cpl_propertylist_copy_property_regexp(subset2, plist2, regexp, 0);
1253 enu_check_error_code(
"error deriving subset of propertylist");
1257 cpl_size size1 = cpl_propertylist_get_size(subset1);
1258 cpl_size size2 = cpl_propertylist_get_size(subset2);
1259 if (size1 != size2) {
1261 reason = cpl_sprintf(
"propertylist sizes do not match %s",
" ");
1263 for (cpl_size i=0; i<size1; i++) {
1264 const cpl_property * prop1 = cpl_propertylist_get(subset1, i);
1265 const char * name = cpl_property_get_name(prop1);
1269 if (!cpl_propertylist_has(subset2, name)) {
1271 reason = cpl_sprintf(
"name mismatch: %s", name);
1277 const cpl_property * prop2 = cpl_propertylist_get_property_const(
1280 const cpl_type type1 = cpl_property_get_type(prop1);
1281 const cpl_type type2 = cpl_property_get_type(prop2);
1282 if (type1 != type2) {
1284 reason = cpl_sprintf(
"type mismatch: name %s: types %s and %s",
1286 cpl_type_get_name(cpl_property_get_type(prop1)),
1287 cpl_type_get_name(cpl_property_get_type(prop2)));
1295 const char * sval1 = NULL;
1296 const char * sval2 = NULL;
1300 ival1 = cpl_property_get_int(prop1);
1301 ival2 = cpl_property_get_int(prop2);
1302 if (ival1 != ival2) {
1304 reason = cpl_sprintf(
"value mismatch: name %s: values "
1306 name, ival1, ival2);
1310 case CPL_TYPE_DOUBLE:
1311 dval1 = cpl_property_get_double(prop1);
1312 dval2 = cpl_property_get_double(prop2);
1313 if (dval1 != dval2) {
1315 reason = cpl_sprintf(
"value mismatch: name %s: values "
1317 name, dval1, dval2);
1321 case CPL_TYPE_STRING:
1322 sval1 = cpl_property_get_string(prop1);
1323 sval2 = cpl_property_get_string(prop2);
1324 if (strcmp(sval1, sval2)) {
1326 reason = cpl_sprintf(
"value mismatch: name %s: values "
1327 "%s, %s", name, sval1, sval2);
1331 cpl_error_set_message(cpl_func, CPL_ERROR_INVALID_TYPE,
1332 "name %s, type not handled %s",
1333 name, cpl_type_get_name(
1334 cpl_property_get_type(prop1)));
1337 if (result == CPL_FALSE)
break;
1344 cpl_propertylist_delete(subset1);
1345 cpl_propertylist_delete(subset2);
1350 if (cpl_error_get_code() != CPL_ERROR_NONE) {
1353 if (result == CPL_FALSE) {
1354 cpl_msg_info(cpl_func,
"propertylist match failed: %s", reason);
1379 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
1380 cpl_ensure_code(limage, CPL_ERROR_NULL_INPUT);
1382 cpl_boolean bad = CPL_FALSE;
1384 if (!cpl_propertylist_has(limage->plist,
"CTYPE1")) {
1388 const char * ctype1 = cpl_propertylist_get_string(
1389 limage->plist,
"CTYPE1");
1390 if (strstr(ctype1,
"PIXEL")) {
1396 if (!cpl_propertylist_has(limage->plist,
"CTYPE2")) {
1401 const char * ctype2 = cpl_propertylist_get_string(limage->plist,
1403 if (strstr(ctype2,
"PIXEL")) {
1408 cpl_error_code result = cpl_error_get_code();
1409 if (result == CPL_ERROR_NONE && bad) {
1410 result = cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
1411 "missing or bad WCS keywords in "
1413 cpl_frame_get_filename(limage->frame));
1436 const located_imagelist * limlist,
1437 const char * nameroot,
1438 const char * recipename,
1439 cpl_frameset * frameset,
1440 const cpl_parameterlist * parlist,
1441 const cpl_frameset * used) {
1443 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
1444 if (!debug)
return CPL_ERROR_NONE;
1445 cpl_ensure_code(limlist, CPL_ERROR_NULL_INPUT);
1446 cpl_ensure_code(nameroot, CPL_ERROR_NULL_INPUT);
1447 cpl_ensure_code(recipename, CPL_ERROR_NULL_INPUT);
1448 cpl_ensure_code(frameset, CPL_ERROR_NULL_INPUT);
1449 cpl_ensure_code(parlist, CPL_ERROR_NULL_INPUT);
1450 cpl_ensure_code(used, CPL_ERROR_NULL_INPUT);
1452 char * fname_copy = NULL;
1453 char * out_fname = NULL;
1457 cpl_size nimages = limlist->size;
1458 for (cpl_size j = 0; j < nimages; j++) {
1462 if (limlist->limages[j]->frame) {
1463 const char * fname = cpl_frame_get_filename(
1464 limlist->limages[j]->frame);
1465 fname_copy = cpl_strdup(fname);
1467 const char * fname =
"debug.fits";
1468 fname_copy = cpl_strdup(fname);
1470 out_fname = cpl_sprintf(
"%s_%s", nameroot, basename(fname_copy));
1472 cpl_msg_info(cpl_func,
"..(debug) writing %s", out_fname);
1476 cpl_propertylist * applist = cpl_propertylist_new();
1477 cpl_propertylist_update_string(applist,
1480 cpl_propertylist_update_string(applist,
"PRODCATG",
"ANCILLARY.IMAGE");
1482 cpl_frameset * provenance = cpl_frameset_new();
1483 cpl_frameset_insert(provenance, cpl_frame_duplicate(
1484 limlist->limages[j]->frame));
1490 limlist->limages[j],
1492 limlist->limages[j]->frame,
1497 cpl_frameset_delete(provenance);
1498 cpl_propertylist_delete(applist);
1503 cpl_free(fname_copy);
1504 cpl_free(out_fname);
1506 return cpl_error_get_code();
1526 const hdrl_image * flat,
1527 const cpl_image * confidence,
1528 const cpl_mask * cold_bpm,
1529 cpl_frameset * frameset,
1530 const cpl_parameterlist * parlist,
1531 const char * filename,
1532 const char * recipe_name,
1533 const cpl_propertylist* qclog) {
1537 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
1538 cpl_ensure_code(pro_catg, CPL_ERROR_NULL_INPUT);
1539 cpl_ensure_code(flat, CPL_ERROR_NULL_INPUT);
1540 cpl_ensure_code(frameset, CPL_ERROR_NULL_INPUT);
1541 cpl_ensure_code(parlist, CPL_ERROR_NULL_INPUT);
1542 cpl_ensure_code(filename, CPL_ERROR_NULL_INPUT);
1544 mef_extension_list * mefs = NULL;
1545 cpl_propertylist * plist = NULL;
1549 plist = cpl_propertylist_new();
1550 cpl_propertylist_append_string(plist, CPL_DFS_PRO_CATG, pro_catg);
1551 cpl_propertylist_update_string(plist,
"PRODCATG",
"ANCILLARY.IMAGE");
1556 cpl_propertylist_append_double(plist,
"ESO QC FLAT MED",
1557 (
double) flat_median.data);
1558 cpl_propertylist_set_comment(plist,
"ESO QC FLAT MED",
1559 "[ADU] Median flat value");
1561 cpl_propertylist_append_double(plist,
"ESO QC FLAT MEAN",
1562 (
double) flat_mean.data);
1563 cpl_propertylist_set_comment(plist,
"ESO QC FLAT MEAN",
1564 "[ADU] Mean flat value");
1566 cpl_propertylist_append_double(plist,
"ESO QC FLAT RMS", flat_rms);
1567 cpl_propertylist_set_comment(plist,
"ESO QC FLAT RMS",
1568 "[ADU] RMS flat value");
1572 cpl_msg_warning(cpl_func,
"TBD: calculate normalisation value");
1574 cpl_propertylist_append_double(plist,
"ESO QC FLAT NORM", norm);
1575 cpl_propertylist_set_comment(plist,
"ESO QC FLAT NORM",
1576 "Flat normalisation value");
1578 cpl_size ncold = cpl_mask_count(cold_bpm);
1579 cpl_propertylist_append_double(plist,
"ESO QC NUMBER COLD PIXELS",
1581 cpl_propertylist_set_comment(plist,
"ESO QC NUMBER COLD PIXELS",
1582 "Number of cold pixels");
1584 enu_check_error_code(
"error constructing output propertylist");
1597 cpl_propertylist_append(plist, qclog);
1610 PACKAGE
"/" PACKAGE_VERSION,
1615 cpl_propertylist_delete(plist);
1617 return cpl_error_get_code();
1638 if (cpl_error_get_code() != CPL_ERROR_NONE)
return 0;
1640 double airmass = 0.0;
1641 if (cpl_propertylist_has(plist,
"AIRMASS")) {
1643 airmass = cpl_propertylist_get_double(plist,
"AIRMASS");
1644 }
else if (cpl_propertylist_has(plist,
"ESO TEL AIRM START")) {
1646 const double astart = cpl_propertylist_get_double(plist,
1647 "ESO TEL AIRM START");
1648 const double aend = cpl_propertylist_get_double(plist,
1649 "ESO TEL AIRM END");
1650 airmass = (astart + aend) / 2.0;
1653 const double alt = cpl_propertylist_get_double(plist,
1655 double airmass_alt = 1.0 /cos ((90.0 - alt) * CPL_MATH_PI / 180.0);
1656 if (fabs(airmass - airmass_alt) > 1.0e-5) {
1657 cpl_msg_warning(cpl_func,
"discrepency - airmass %f -> %f", airmass,
1659 airmass = airmass_alt;
1662 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
1663 "unable to find airmass information");
1692 const cpl_parameterlist * parlist,
1693 const cpl_propertylist * plist,
1694 double * obj_core_radius,
1695 int * bkg_mesh_size) {
1697 char * param_name = NULL;
1698 const cpl_parameter * p = NULL;
1699 cpl_parameter * par = NULL;
1701 cpl_error_code error_code = CPL_ERROR_NONE;
1703 if (cpl_error_get_code() != CPL_ERROR_NONE)
return CPL_ERROR_NONE;
1704 cpl_ensure_code(context, CPL_ERROR_NULL_INPUT);
1705 cpl_ensure_code(parlist, CPL_ERROR_NULL_INPUT);
1706 cpl_ensure_code(plist, CPL_ERROR_NULL_INPUT);
1707 cpl_ensure_code(obj_core_radius, CPL_ERROR_NULL_INPUT);
1708 cpl_ensure_code(bkg_mesh_size, CPL_ERROR_NULL_INPUT);
1710 *obj_core_radius = -1.0;
1711 *bkg_mesh_size = -1;
1716 param_name = cpl_sprintf(
"%s.catalogue.ao-params", context);
1717 cpl_msg_info(cpl_func,
"%s", param_name);
1718 p = cpl_parameterlist_find_const(parlist, param_name);
1719 const char * ao_params = cpl_parameter_get_string(p);
1720 cpl_free(param_name);
1721 enu_check_error_code(
"failed to read ao-params");
1723 if (!strcmp(ao_params,
"auto") &&
1724 cpl_propertylist_has(plist,
"ESO OBS AOMODE")) {
1725 const char * aomode = cpl_propertylist_get_string(plist,
1727 cpl_msg_info(cpl_func,
"AOMODE: %s", aomode);
1729 if (!strcmp(aomode,
"FULL_AO")) {
1730 *obj_core_radius = 10.0;
1731 *bkg_mesh_size = 64;
1732 if(!strstr(context,
"skysub")) {
1734 param_name = cpl_sprintf(
"%s.catalogue.obj.core-radius", context);
1735 par = (cpl_parameter*) cpl_parameterlist_find_const(parlist, param_name);
1736 cpl_parameter_set_double(par, *obj_core_radius);
1737 cpl_free(param_name);
1741 param_name = cpl_sprintf(
"%s.catalogue.bkg.mesh-size", context);
1742 par = (cpl_parameter*) cpl_parameterlist_find_const(parlist, param_name);
1743 cpl_parameter_set_int(par, *bkg_mesh_size);
1744 cpl_free(param_name);
1746 }
else if (!strcmp(aomode,
"NO_AO")) {
1747 *obj_core_radius = 25.0;
1748 *bkg_mesh_size = 128;
1749 if(!strstr(context,
"skysub")) {
1751 param_name = cpl_sprintf(
"%s.catalogue.obj.core-radius", context);
1752 par = (cpl_parameter*) cpl_parameterlist_find_const(parlist, param_name);
1753 cpl_parameter_set_double(par, *obj_core_radius);
1754 cpl_free(param_name);
1757 param_name = cpl_sprintf(
"%s.catalogue.bkg.mesh-size", context);
1758 par = (cpl_parameter*) cpl_parameterlist_find_const(parlist, param_name);
1759 cpl_parameter_set_int(par, *bkg_mesh_size);
1760 cpl_free(param_name);
1762 error_code = CPL_ERROR_DATA_NOT_FOUND;
1766 if(!strstr(context,
"skysub")) {
1767 param_name = cpl_sprintf(
"%s.catalogue.obj.core-radius", context);
1768 p = cpl_parameterlist_find_const(parlist, param_name);
1769 *obj_core_radius = cpl_parameter_get_double(p);
1770 cpl_free(param_name);
1773 param_name = cpl_sprintf(
"%s.catalogue.bkg.mesh-size", context);
1774 p = cpl_parameterlist_find_const(parlist, param_name);
1775 *bkg_mesh_size = cpl_parameter_get_int(p);
1776 cpl_free(param_name);
1778 cpl_msg_info(cpl_func,
"catalogue AO-related params: %s %f %d",
1779 ao_params, *obj_core_radius, *bkg_mesh_size);
1801 if (cpl_error_get_code() != CPL_ERROR_NONE)
return 0;
1803 double tel_alt = 0.0;
1804 if (cpl_propertylist_has(plist,
"ESO TEL ALT")) {
1806 tel_alt = cpl_propertylist_get_double(plist,
"ESO TEL ALT");
1808 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
1809 "unable to find ESO TEL ALT information");
1831 if (cpl_error_get_code() != CPL_ERROR_NONE)
return 0;
1834 if (cpl_propertylist_has(plist,
"ESO DET DIT")) {
1835 dit = cpl_propertylist_get_double(plist,
"ESO DET DIT");
1836 }
else if (cpl_propertylist_has(plist,
"ESO DET SEQ1 DIT")) {
1837 dit = cpl_propertylist_get_double(plist,
"ESO DET SEQ1 DIT");
1839 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
1840 "no keyword ESO.DET.DIT or ESO.DET.SEQ1.DIT");
1862 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
1864 const char * result = NULL;
1865 if (cpl_propertylist_has(plist,
"ESO DET NCORRS NAME")) {
1866 result = cpl_propertylist_get_string(plist,
"ESO DET NCORRS NAME");
1867 }
else if (cpl_propertylist_has(plist,
"ESO DET READ CURNAME")) {
1868 result = cpl_propertylist_get_string(plist,
"ESO DET READ CURNAME");
1870 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
1871 "no keyword ESO.DET.NCORRS.NAME or "
1872 "ESO.DET.READ.CURNAME");
1898 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
1899 cpl_ensure(plist, CPL_ERROR_NULL_INPUT, NULL);
1901 const char * filter = NULL;
1903 if (cpl_propertylist_has(plist,
"FILTER")) {
1904 filter = cpl_propertylist_get_string(plist,
"FILTER");
1905 }
else if (cpl_propertylist_has(plist,
"INSTRUME")) {
1906 const char * instrume = cpl_propertylist_get_string(plist,
"INSTRUME");
1908 if (!strcmp(instrume,
"NAOS+CONICA")) {
1909 filter = cpl_propertylist_get_string(plist,
"ESO INS OPTI4 NAME");
1910 cpl_msg_info(cpl_func,
"OPTI4 %s", filter);
1911 if (!strcmp(filter,
"clear") || !strcmp(filter,
"empty")) {
1912 filter = cpl_propertylist_get_string(plist,
"ESO INS OPTI5 NAME");
1913 cpl_msg_info(cpl_func,
"OPTI5 %s", filter);
1915 if (!strcmp(filter,
"clear") || !strcmp(filter,
"empty")) {
1916 filter = cpl_propertylist_get_string(plist,
"ESO INS OPTI6 NAME");
1917 cpl_msg_info(cpl_func,
"OPTI6 %s", filter);
1919 }
else if (!strcmp(instrume,
"ERIS")) {
1920 filter = cpl_propertylist_get_string(plist,
"ESO INS2 NXFW NAME");
1922 cpl_error_set(cpl_func, CPL_ERROR_INCOMPATIBLE_INPUT);
1949 if (cpl_error_get_code() != CPL_ERROR_NONE)
return -1.0;
1950 cpl_ensure(filter, CPL_ERROR_NULL_INPUT, -1.0);
1952 double lambda = 1.0;
1953 if (!strcmp(filter,
"J")) {
1955 }
else if (!strcmp(filter,
"H")) {
1957 }
else if (!strcmp(filter,
"Ks")) {
1959 }
else if (!strcmp(filter,
"Short-Lp")) {
1961 }
else if (!strcmp(filter,
"L-Broad")) {
1963 }
else if (!strcmp(filter,
"Lp")) {
1965 }
else if (!strcmp(filter,
"Mp")) {
1967 }
else if (!strcmp(filter,
"Pa-b")) {
1969 }
else if (!strcmp(filter,
"Fe-II")) {
1971 }
else if (!strcmp(filter,
"H2-cont")) {
1973 }
else if (!strcmp(filter,
"H2-1-0S")) {
1975 }
else if (!strcmp(filter,
"Br-g")) {
1977 }
else if (!strcmp(filter,
"K-peak")) {
1979 }
else if (!strcmp(filter,
"IB-2.42")) {
1981 }
else if (!strcmp(filter,
"IB-2.48")) {
1983 }
else if (!strcmp(filter,
"Br-a-cont")) {
1985 }
else if (!strcmp(filter,
"Br-a")) {
1988 cpl_msg_warning(cpl_func,
"filter %s not recognized", filter);
2006 const char * eris_nix_license =
2007 "This file is part of the ERIS/NIX Instrument Pipeline\n"
2008 "Copyright (C) 2017 European Southern Observatory\n"
2010 "This program is free software; you can redistribute it and/or modify\n"
2011 "it under the terms of the GNU General Public License as published by\n"
2012 "the Free Software Foundation; either version 2 of the License, or\n"
2013 "(at your option) any later version.\n"
2015 "This program is distributed in the hope that it will be useful,\n"
2016 "but WITHOUT ANY WARRANTY; without even the implied warranty of\n"
2017 "MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the\n"
2018 "GNU General Public License for more details.\n"
2020 "You should have received a copy of the GNU General Public License\n"
2021 "along with this program; if not, write to the Free Software\n"
2022 "Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, \n"
2023 "MA 02110-1301 USA";
2024 return eris_nix_license;
2039 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
2040 cpl_ensure(wcs, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
2041 cpl_ensure(ra, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
2042 cpl_ensure(dec, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
2046 const cpl_array * dims = cpl_wcs_get_image_dims(wcs);
2047 cpl_matrix * from = cpl_matrix_new(1, 2);
2049 double dim = cpl_array_get(dims, 0, &ignore);
2050 cpl_matrix_set(from, 0, 0, dim / 2.0);
2051 dim = cpl_array_get(dims, 1, &ignore);
2052 cpl_matrix_set(from, 0, 1, dim / 2.0);
2056 cpl_matrix * to = NULL;
2057 cpl_array * status = NULL;
2058 cpl_wcs_convert(wcs, from, &to, &status, CPL_WCS_PHYS2WORLD);
2060 *ra = cpl_matrix_get(to, 0, 0);
2061 *dec = cpl_matrix_get(to, 0, 1);
2065 cpl_matrix_delete(from);
2066 cpl_matrix_delete(to);
2067 cpl_array_delete(status);
2069 return cpl_error_get_code();
2072#define MJD_OBJ_REF_WINDOW_STRY_POS_FAST_UNCORR 61014.7
2093 cpl_size * strx, cpl_size * stry,
2094 cpl_size * nx_chip, cpl_size * ny_chip,
2095 cpl_boolean * windowed,
2096 const cpl_propertylist * plist) {
2098 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
2100 cpl_ensure_code(plist, CPL_ERROR_NULL_INPUT);
2101 enu_check(cpl_propertylist_has(plist,
"INSTRUME") == 1,
2102 CPL_ERROR_DATA_NOT_FOUND,
2103 "propertylist does not contain INSTRUME keyword");
2107 const char * instrume = cpl_propertylist_get_string(plist,
"INSTRUME");
2115 if (!strcmp(instrume,
"ERIS")) {
2116 nxw = cpl_propertylist_get_int(plist,
"ESO DET SEQ1 WIN NX");
2117 nyw = cpl_propertylist_get_int(plist,
"ESO DET SEQ1 WIN NY");
2118 rotw = cpl_propertylist_get_int(plist,
"ESO DET SEQ1 WIN ROT");
2119 strxw = cpl_propertylist_get_int(plist,
"ESO DET SEQ1 WIN STRX");
2120 stryw = cpl_propertylist_get_int(plist,
"ESO DET SEQ1 WIN STRY");
2124 nxw = cpl_propertylist_get_int(plist,
"ESO DET WIN NX");
2125 nyw = cpl_propertylist_get_int(plist,
"ESO DET WIN NY");
2126 if (cpl_propertylist_get_type(plist,
"ESO DET WIN STARTX") ==
2128 strxw = (int) cpl_propertylist_get_double(plist,
2129 "ESO DET WIN STARTX");
2130 stryw = (int) cpl_propertylist_get_double(plist,
2131 "ESO DET WIN STARTY");
2133 strxw = cpl_propertylist_get_int(plist,
"ESO DET WIN STARTX");
2134 stryw = cpl_propertylist_get_int(plist,
"ESO DET WIN STARTY");
2138 cpl_size nx_chipw = 2048;
2139 cpl_size ny_chipw = 2048;
2141 if(cpl_propertylist_has(plist,
"ESO DET CHIP NX")) {
2142 nx_chipw = cpl_propertylist_get_int(plist,
"ESO DET CHIP NX");
2144 nx_chipw = cpl_propertylist_get_int(plist,
"ESO DET CHIP1 NX");
2147 if(cpl_propertylist_has(plist,
"ESO DET CHIP NY")) {
2148 ny_chipw = cpl_propertylist_get_int(plist,
"ESO DET CHIP NY");
2150 ny_chipw = cpl_propertylist_get_int(plist,
"ESO DET CHIP1 NY");
2153 enu_check_error_code(
"failed to read detector mode information");
2154 enu_check(rotw == 0, CPL_ERROR_UNSUPPORTED_MODE,
2155 "detector window rot must be 0");
2157 cpl_boolean windowedw = (rotw != 0) ||
2160 (nxw != nx_chipw) ||
2162 cpl_msg_info(cpl_func,
"flatfield rot=%d strx=%d stry=%d nx=%d ny=%d windowedw: %d",
2163 (
int)rotw, (
int)strxw, (
int)stryw, (
int)nxw, (
int)nyw, windowedw);
2165 if (!strcmp(instrume,
"ERIS")) {
2167 double mjd_obs = MJD_OBJ_REF_WINDOW_STRY_POS_FAST_UNCORR;
2168 if(cpl_propertylist_has(plist,
"ESO DET SEQ1 WIN NAME")) {
2169 mjd_obs = cpl_propertylist_get_double(plist,
"MJD-OBS");
2172 const char * det_mode =
"SLOW_GR_UTR";
2173 if(cpl_propertylist_has(plist,
"ESO DET READ CURNAME")) {
2177 const char * win_name =
"windowF";
2178 if(cpl_propertylist_has(plist,
"ESO DET SEQ1 WIN NAME")) {
2179 win_name = cpl_propertylist_get_string(plist,
"ESO DET SEQ1 WIN NAME");
2182 if ( windowedw && (strcmp(win_name,
"windowF") != 0) &&
2183 (strcmp(det_mode,
"FAST_UNCORR") == 0) && ((stryw % 2) == 0) && (mjd_obs < MJD_OBJ_REF_WINDOW_STRY_POS_FAST_UNCORR) ){
2185 cpl_msg_warning(cpl_func,
"windowing special case windowedw: %d win_name: %s", windowedw, win_name);
2187 cpl_msg_warning(cpl_func,
"windowing normal case windowedw: %d win_name: %s", windowedw, win_name);
2200 *nx_chip = nx_chipw;
2201 *ny_chip = ny_chipw;
2202 *windowed = windowedw;
2206 return cpl_error_get_code();
2218 const hdrl_catalogue_result * target) {
2219 cpl_ensure(target, CPL_ERROR_NULL_INPUT, NULL);
2221 hdrl_catalogue_result * result = cpl_malloc(
sizeof(hdrl_catalogue_result));
2223 result->catalogue = cpl_table_duplicate(target->catalogue);
2224 result->background = cpl_image_duplicate(target->background);
2225 result->segmentation_map = cpl_image_duplicate(target->segmentation_map);
2226 result->qclist = cpl_propertylist_duplicate(target->qclist);
2246 hdrl_image ** result,
2247 mef_extension_list ** mef_extensions,
2248 cpl_propertylist ** plist) {
2250 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
2252 cpl_propertylist * plist1 = NULL;
2253 cpl_propertylist * plist2 = NULL;
2254 cpl_propertylist * plist3 = NULL;
2255 cpl_image * data = NULL;
2256 cpl_image * error = NULL;
2257 cpl_mask * bpm = NULL;
2258 cpl_mask * old_bpm = NULL;
2262 cpl_ensure_code(filename, CPL_ERROR_NULL_INPUT);
2266 cpl_size size = cpl_fits_count_extensions(filename);
2267 cpl_ensure_code(size >= 3, CPL_ERROR_BAD_FILE_FORMAT);
2271 *plist = cpl_propertylist_load(filename, 0);
2272 enu_check_error_code(
"accessing HDU 0 header");
2276 plist1 = cpl_propertylist_load(filename, 1);
2277 enu_check_error_code(
"accessing HDU 1 header");
2279 const char * extname = cpl_propertylist_get_string(plist1,
"EXTNAME");
2280 const char * hduclass = cpl_propertylist_get_string(plist1,
"HDUCLASS");
2281 const char * hdudoc = cpl_propertylist_get_string(plist1,
"HDUDOC");
2283 const char * hduclas1 = cpl_propertylist_get_string(plist1,
"HDUCLAS1");
2284 const char * hduclas2 = cpl_propertylist_get_string(plist1,
"HDUCLAS2");
2285 enu_check_error_code(
"accessing HDU 1 propertylist");
2286 enu_check(strstr(extname,
"DATA") != NULL &&
2287 strstr(hduclass,
"ESO") != NULL &&
2288 strstr(hdudoc,
"SDP") != NULL &&
2289 strstr(hduclas1,
"IMAGE") != NULL &&
2290 strstr(hduclas2,
"DATA") != NULL,
2291 CPL_ERROR_BAD_FILE_FORMAT,
"bad HD1 format");
2292 data = cpl_image_load(filename, HDRL_TYPE_DATA, 0, 1);
2293 enu_check_error_code(
"accessing HDU 1 data");
2300 const char * wcs_items =
"NAXIS WCSAXES CRPIX1 CRPIX2 "
2301 "CD1_1 CD1_2 CD2_1 CD2_2 CUNIT1 CUNIT2 "
2302 "CTYPE1 CTYPE2 CRVAL1 CRVAL2 LONPOLE "
2303 "LATPOLE CSYER1 CSYER2 EQUINOX MJD-OBS "
2304 "DATE-OBS NAXIS2 NAXIS1";
2305 cpl_size plist1_size = cpl_propertylist_get_size(plist1);
2306 for (cpl_size ip=plist1_size; ip>0; ip--) {
2307 const cpl_property * newprop = cpl_propertylist_get_const(plist1, ip-1);
2308 const char * name = cpl_property_get_name(newprop);
2309 if (strstr(wcs_items, name) != NULL) {
2310 if (cpl_propertylist_has(*plist, name)) {
2311 cpl_propertylist_copy_property(*plist, plist1, name);
2313 cpl_propertylist_insert_after_property(*plist,
"NAXIS",
2321 plist2 = cpl_propertylist_load(filename, 2);
2322 enu_check_error_code(
"accessing HDU 2");
2323 extname = cpl_propertylist_get_string(plist2,
"EXTNAME");
2324 hduclass = cpl_propertylist_get_string(plist2,
"HDUCLASS");
2325 hdudoc = cpl_propertylist_get_string(plist2,
"HDUDOC");
2327 hduclas1 = cpl_propertylist_get_string(plist2,
"HDUCLAS1");
2328 hduclas2 = cpl_propertylist_get_string(plist2,
"HDUCLAS2");
2329 enu_check_error_code(
"accessing HDU 2 propertylist");
2330 enu_check(strstr(extname,
"ERR") != NULL &&
2331 strstr(hduclass,
"ESO") != NULL &&
2332 strstr(hdudoc,
"SDP") != NULL &&
2333 strstr(hduclas1,
"IMAGE") != NULL &&
2334 strstr(hduclas2,
"ERROR") != NULL,
2335 CPL_ERROR_BAD_FILE_FORMAT,
"bad HD2 format");
2336 error = cpl_image_load(filename, HDRL_TYPE_ERROR, 0, 2);
2337 enu_check_error_code(
"accessing HDU 1 data error");
2341 plist3 = cpl_propertylist_load(filename, 3);
2342 enu_check_error_code(
"accessing HDU 3");
2343 extname = cpl_propertylist_get_string(plist3,
"EXTNAME");
2344 hduclass = cpl_propertylist_get_string(plist3,
"HDUCLASS");
2345 hdudoc = cpl_propertylist_get_string(plist3,
"HDUDOC");
2347 hduclas1 = cpl_propertylist_get_string(plist3,
"HDUCLAS1");
2348 hduclas2 = cpl_propertylist_get_string(plist3,
"HDUCLAS2");
2349 enu_check_error_code(
"accessing HDU 3 propertylist");
2350 enu_check(strstr(extname,
"DQ") != NULL &&
2351 strstr(hduclass,
"ESO") != NULL &&
2352 strstr(hdudoc,
"SDP") != NULL &&
2353 strstr(hduclas1,
"IMAGE") != NULL &&
2354 strstr(hduclas2,
"QUALITY") != NULL,
2355 CPL_ERROR_BAD_FILE_FORMAT,
"bad HD3 format");
2359 bpm = cpl_mask_load(filename, 0, 3);
2360 enu_check_error_code(
"accessing HDU 3 data mask");
2361 old_bpm = cpl_image_set_bpm(data, bpm);
2362 cpl_mask_delete(old_bpm); old_bpm = NULL;
2367 int n_ext = size - 3;
2369 for (
int ext=0; ext<n_ext; ext++) {
2373 cpl_propertylist * ext_plist = cpl_propertylist_load(filename, ext+4);
2374 const char * name = cpl_propertylist_get_string(ext_plist,
"EXTNAME");
2375 const char * mef_type = cpl_propertylist_get_string(ext_plist,
2376 "ERIS_NIX_MEF_TYPE");
2380 if (!strcmp(mef_type, MEF_EXTENSION_CONTAINING_MASK)) {
2381 cpl_mask * mask = cpl_mask_load(filename, 0, ext+4);
2383 (*mef_extensions)->mef[ext] =
new;
2384 cpl_mask_delete(mask);
2385 }
else if (!strcmp(mef_type, MEF_EXTENSION_CONTAINING_IMAGE)) {
2386 cpl_image * image = cpl_image_load(filename, HDRL_TYPE_DATA, 0,
2389 (*mef_extensions)->mef[ext] =
new;
2390 cpl_image_delete(image);
2391 }
else if (!strcmp(mef_type, MEF_EXTENSION_CONTAINING_TABLE)) {
2392 cpl_table * table = cpl_table_load(filename, ext+4, 0);
2394 (*mef_extensions)->mef[ext] =
new;
2395 cpl_table_delete(table);
2397 cpl_error_set_message(
"enu_himage_load_from_fits",
2398 CPL_ERROR_UNSUPPORTED_MODE,
2399 "unsupported extension name: %s", name);
2402 cpl_propertylist_delete(ext_plist);
2406 cpl_propertylist_delete(plist1);
2407 cpl_propertylist_delete(plist2);
2408 cpl_propertylist_delete(plist3);
2409 cpl_image_delete(data);
2410 cpl_image_delete(error);
2411 cpl_mask_delete(old_bpm);
2412 if (cpl_error_get_code() != CPL_ERROR_NONE) {
2413 cpl_propertylist_delete(*plist);
2418 *mef_extensions = NULL;
2421 return cpl_error_get_code();
2443 cpl_image ** pcopyconf,
2444 const cpl_boolean collapse_cube) {
2446 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
2448 cpl_ensure(frame, CPL_ERROR_NULL_INPUT, NULL);
2449 cpl_ensure(cpl_frame_get_filename(frame), CPL_ERROR_NULL_INPUT, NULL);
2451 located_image * result = NULL;
2452 mef_extension_list * mefs = NULL;
2453 cpl_image * data = NULL;
2454 cpl_image * err = NULL;
2455 cpl_image * qual = NULL;
2456 cpl_image * conf = NULL;
2457 cpl_image * bkg = NULL;
2458 cpl_image * bkg_err = NULL;
2459 cpl_image * bkg_conf = NULL;
2461 const char * filename = cpl_frame_get_filename(frame);
2462 cpl_msg_info(cpl_func,
"..reading %s", filename);
2468 cpl_size next = cpl_fits_count_extensions(filename);
2470 cpl_propertylist * plist = cpl_propertylist_load(filename, 0);
2473 const char * instrume = cpl_propertylist_get_string(plist,
"INSTRUME");
2474 const cpl_boolean eris = !strcmp(instrume,
"ERIS");
2476 char * first_ext = NULL;
2477 if (eris && (next > 0)) {
2478 cpl_propertylist * plist1 = cpl_propertylist_load(filename, 1);
2479 if (cpl_propertylist_has(plist1,
"EXTNAME")) {
2480 first_ext = cpl_sprintf(
"%s", cpl_propertylist_get_string(plist1,
"EXTNAME"));
2482 cpl_propertylist_delete(plist1);
2486 char * frame_format = NULL;
2487 if (cpl_propertylist_has(plist,
"ESO DET FRAM FORMAT")) {
2488 frame_format = cpl_sprintf(
"%s", cpl_propertylist_get_string(plist,
"ESO DET FRAM FORMAT"));
2489 cpl_msg_info(cpl_func,
"..format %s", frame_format);
2495 cpl_msg_info(cpl_func,
"NACO data");
2500 data = cpl_image_load(filename, HDRL_TYPE_DATA, 0, 0);
2501 enu_check_error_code(
"failed to read file: %s", filename);
2503 cpl_image_delete(data); data = NULL;
2509 cpl_image * confidence;
2510 if (pcopyconf != NULL && *pcopyconf != NULL &&
2511 cpl_image_get_size_x(*pcopyconf) == nx &&
2512 cpl_image_get_size_y(*pcopyconf) == ny) {
2513 confidence = *pcopyconf;
2516 confidence = cpl_image_new(nx, ny, HDRL_TYPE_DATA);
2517 cpl_image_fill_window(confidence, 1, 1, nx, ny, 100.0);
2520 cpl_frame_dump(frame, NULL);
2531 cpl_frame_duplicate(frame));
2534 }
else if ((next == 0) ||
2537 !strcmp(first_ext,
"ASM_DATA"))) {
2544 cpl_table * asm_table = cpl_table_load(filename, 1, CPL_TRUE);
2545 cpl_size nrow = cpl_table_get_nrow(asm_table);
2546 if(cpl_table_has_column(asm_table,
"IA_FWHM")) {
2547 double * ia_fwhm_data = cpl_table_get_data_double(asm_table,
2549 cpl_vector * ia_fwhm = cpl_vector_wrap(nrow, ia_fwhm_data);
2550 double ia_fwhm_mean = cpl_vector_get_mean(ia_fwhm);
2558 double ia_fwhm_corr = ia_fwhm_mean * pow(0.5/lambda, 1./5.);
2560 cpl_propertylist_update_double(plist,
"PSF_FWHM", ia_fwhm_corr);
2561 cpl_propertylist_set_comment(plist,
"PSF_FWHM",
"PSF fwhm derived "
2562 "from AO ASM table [arcsec]");
2564 cpl_vector_unwrap(ia_fwhm);
2566 cpl_msg_warning(cpl_func,
"ASM table miss IA_FWHM column. PSF_FWHM not computed");
2568 cpl_table_delete(asm_table);
2575 if (cpl_propertylist_has(plist,
"ESO TEL AMBI FWHM START") &&
2576 cpl_propertylist_has(plist,
"ESO TEL AMBI FWHM END") &&
2577 cpl_propertylist_has(plist,
"ESO TEL AIRM START") &&
2578 cpl_propertylist_has(plist,
"ESO TEL AIRM END")) {
2580 double dimm1 = cpl_propertylist_get_double(plist,
2581 "ESO TEL AMBI FWHM START");
2582 double dimm2 = cpl_propertylist_get_double(plist,
2583 "ESO TEL AMBI FWHM END");
2584 double airmass1 = cpl_propertylist_get_double(plist,
2585 "ESO TEL AIRM START");
2586 double airmass2 = cpl_propertylist_get_double(plist,
2587 "ESO TEL AIRM END");
2595 pow(0.5/lambda, 1./5.) *
2596 pow(airmass1, 3./5.) *
2598 (pow(lambda/1.e6, 2./5.) *
2599 pow(airmass1, -1./5.) /
2604 pow(0.5/lambda, 1./5.) *
2605 pow(airmass2, 3./5.) *
2607 (pow(lambda/1.e6, 2./5.) *
2608 pow(airmass2, -1./5.) /
2613 cpl_propertylist_update_double(plist,
"PSF_FWHM",
2614 (dimm1 + dimm2) / 2.0);
2615 cpl_propertylist_set_comment(plist,
"PSF_FWHM",
2617 "from DIMM [arcsec]");
2621 if (!strcmp(frame_format,
"single")) {
2625 data = cpl_image_load(filename, HDRL_TYPE_DATA, 0, 0);
2626 enu_check_error_code(
"failed to read file: %s", filename);
2628 cpl_image_delete(data); data = NULL;
2634 cpl_image * confidence;
2635 if (pcopyconf != NULL && *pcopyconf != NULL &&
2636 cpl_image_get_size_x(*pcopyconf) == nx &&
2637 cpl_image_get_size_y(*pcopyconf) == ny) {
2638 confidence = *pcopyconf;
2641 confidence = cpl_image_new(nx, ny, HDRL_TYPE_DATA);
2642 cpl_image_fill_window(confidence, 1, 1, nx, ny, 100.0);
2655 cpl_frame_duplicate(frame));
2658 }
else if (!strcmp(frame_format,
"cube")) {
2662 cpl_imagelist * datalist = cpl_imagelist_load(filename,
2664 enu_check_error_code(
"failed to read file: %s", filename);
2667 cpl_imagelist_delete(datalist); datalist = NULL;
2673 cpl_image * confidence;
2674 if (pcopyconf != NULL && *pcopyconf != NULL &&
2675 cpl_image_get_size_x(*pcopyconf) == nx &&
2676 cpl_image_get_size_y(*pcopyconf) == ny) {
2677 confidence = *pcopyconf;
2680 confidence = cpl_image_new(nx, ny, HDRL_TYPE_DATA);
2681 cpl_image_fill_window(confidence, 1, 1, nx, ny, 100.0);
2685 hdrl_image * himage = NULL;
2686 cpl_image* conf_ima = NULL;
2687 if(collapse_cube == 1) {
2689 }
else if(collapse_cube == 2) {
2691 }
else if(collapse_cube == 3) {
2694 cpl_image_delete(conf_ima);
2705 cpl_frame_duplicate(frame));
2717 cpl_frame_duplicate(frame));
2725 cpl_propertylist_delete(plist);
2728 enu_check_error_code(
"Unable to load components for file");
2732 const char * hduclass = cpl_propertylist_get_string(
2733 mefs->mef[0]->plist,
"HDUCLASS");
2734 const char * hdudoc = cpl_propertylist_get_string(
2735 mefs->mef[0]->plist,
"HDUDOC");
2736 const char * hduvers = cpl_propertylist_get_string(
2737 mefs->mef[0]->plist,
"HDUVERS");
2738 enu_check_error_code(
"Unable to read ESO DPS keywords from header");
2740 enu_check(strstr(hduclass,
"ESO") != NULL &&
2741 strstr(hdudoc,
"SDP") != NULL &&
2742 strstr(hduvers,
"SDP version") != NULL,
2743 CPL_ERROR_BAD_FILE_FORMAT,
"file not in ESO DPS format");
2748 const char * scidata = NULL;
2749 if (mefs->size > 0 &&
2750 cpl_propertylist_has(mefs->mef[1]->plist,
"SCIDATA")){
2751 scidata = cpl_propertylist_get_string(mefs->mef[1]->plist,
2754 const char * errdata = NULL;
2755 if (cpl_propertylist_has(mefs->mef[0]->plist,
"ERRDATA")){
2756 errdata = cpl_propertylist_get_string(mefs->mef[0]->plist,
2759 const char * qualdata = NULL;
2760 if (cpl_propertylist_has(mefs->mef[0]->plist,
"QUALDATA")){
2761 qualdata = cpl_propertylist_get_string(mefs->mef[0]->plist,
2764 const char * confdata = NULL;
2765 if (cpl_propertylist_has(mefs->mef[0]->plist,
"CONFDATA")){
2766 confdata = cpl_propertylist_get_string(mefs->mef[0]->plist,
2769 const char * bkgdata = NULL;
2770 if (cpl_propertylist_has(mefs->mef[0]->plist,
"BKGDATA")){
2771 bkgdata = cpl_propertylist_get_string(mefs->mef[0]->plist,
2774 const char * bkgerr = NULL;
2775 if (cpl_propertylist_has(mefs->mef[0]->plist,
"BKGERR")){
2776 bkgerr = cpl_propertylist_get_string(mefs->mef[0]->plist,
2779 const char * bkgconf = NULL;
2780 if (cpl_propertylist_has(mefs->mef[0]->plist,
"BKGCONF")){
2781 bkgconf = cpl_propertylist_get_string(mefs->mef[0]->plist,
2784 enu_check_error_code(
"accessing DPS standard extensions");
2788 for (cpl_size i=0; i<next; i++) {
2789 const char * extname = cpl_propertylist_get_string(
2790 mefs->mef[i]->plist,
"EXTNAME");
2795 if (scidata && !strcmp(scidata, extname)) {
2798 enu_check_error_code(
"accessing SCIDATA extension");
2803 cpl_propertylist * primary_plist = plist;
2810 plist = cpl_propertylist_duplicate(mefs->mef[i]->plist);
2811 const char * sdp7_items =
" EXTNAME HDUCLAS1 HDUCLAS2"
2813 cpl_size psize = cpl_propertylist_get_size(plist);
2814 for (cpl_size ip = psize; ip > 0; ip--) {
2815 const cpl_property * prop = cpl_propertylist_get_const(
2817 const char * name = cpl_property_get_name(prop);
2823 char * padded_name = cpl_sprintf(
" %s ", name);
2824 if (strstr(sdp7_items, padded_name) != NULL) {
2825 cpl_propertylist_erase(plist, name);
2827 cpl_free(padded_name);
2832 psize = cpl_propertylist_get_size(primary_plist);
2833 for (cpl_size ip = psize; ip > 0; ip--) {
2834 const cpl_property * prop = cpl_propertylist_get_const(
2835 primary_plist, ip-1);
2836 const char * name = cpl_property_get_name(prop);
2837 if (!cpl_propertylist_has(plist, name)) {
2838 cpl_propertylist_copy_property(plist, primary_plist,
2842 cpl_propertylist_delete(primary_plist);
2844 }
else if (errdata && !strcmp(errdata, extname)) {
2847 enu_check_error_code(
"accessing SCIERR extension");
2848 }
else if (qualdata && !strcmp(qualdata, extname)) {
2850 "IMAGE",
"QUALITY");
2851 enu_check_error_code(
"accessing QUALDATA extension");
2852 }
else if (confdata && !strcmp(confdata, extname)) {
2855 enu_check_error_code(
"accessing CONFDATA extension");
2856 }
else if (bkgdata && !strcmp(bkgdata, extname)) {
2858 "IMAGE",
"BKG_DATA");
2859 enu_check_error_code(
"accessing BKGDATA extension");
2860 }
else if (bkgerr && !strcmp(bkgerr, extname)) {
2862 "IMAGE",
"BKG_ERR");
2863 enu_check_error_code(
"accessing BKGERR extension");
2864 }
else if (bkgconf && !strcmp(bkgconf, extname)) {
2866 "IMAGE",
"BKG_CONF");
2867 enu_check_error_code(
"accessing BKGCONF extension");
2873 enu_check(data != NULL, CPL_ERROR_BAD_FILE_FORMAT,
"no data in file");
2874 cpl_mask * bpm = cpl_image_get_bpm(data);
2875 cpl_mask_threshold_image(bpm, qual, 0.9, 1.1, CPL_BINARY_1);
2876 bpm = cpl_image_get_bpm(err);
2877 cpl_mask_threshold_image(bpm, qual, 0.9, 1.1, CPL_BINARY_1);
2879 hdrl_image * hbkg = NULL;
2894 cpl_frame_duplicate(frame));
2897 cpl_image_delete(data);
2898 cpl_image_delete(err);
2899 cpl_image_delete(qual);
2900 cpl_image_delete(bkg);
2901 cpl_image_delete(bkg_err);
2912 cpl_free(frame_format);
2913 cpl_free(first_ext);
2914 cpl_propertylist_delete(plist);
2915 cpl_image_delete(data);
2916 cpl_image_delete(err);
2917 cpl_image_delete(qual);
2918 cpl_image_delete(conf);
2919 cpl_image_delete(bkg);
2920 cpl_image_delete(bkg_err);
2921 cpl_image_delete(bkg_conf);
2923 if (cpl_error_get_code() != CPL_ERROR_NONE) {
2946 cpl_propertylist ** plist) {
2948 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
2950 cpl_ensure(filename, CPL_ERROR_NULL_INPUT, NULL);
2955 cpl_size next = cpl_fits_count_extensions(filename);
2957 mef_extension *
new = NULL;
2958 cpl_propertylist * eplist = NULL;
2962 *plist = cpl_propertylist_load(filename, 0);
2963 enu_check_error_code(
"accessing HDU 0 header");
2967 for (
int i=0; i<next; i++) {
2972 cpl_propertylist_delete(eplist);
2973 eplist = cpl_propertylist_load(filename, hdu);
2974 const char * name = cpl_propertylist_get_string(eplist,
"EXTNAME");
2975 const char * hduclas1 = cpl_propertylist_get_string(eplist,
2977 enu_check_error_code(
"error accessing keywords in HDU %d", hdu);
2981 if (!strcmp(hduclas1, MEF_EXTENSION_CONTAINING_MASK)) {
2982 cpl_mask * mask = cpl_mask_load(filename, 0, hdu);
2984 cpl_mask_delete(mask);
2985 }
else if (!strcmp(hduclas1,
"IMAGE")) {
2986 cpl_image * image = cpl_image_load(filename, HDRL_TYPE_DATA, 0, hdu);
2988 cpl_image_delete(image);
2989 }
else if (!strcmp(hduclas1, MEF_EXTENSION_CONTAINING_TABLE)) {
2990 cpl_table * table = cpl_table_load(filename, hdu, 0);
2992 cpl_table_delete(table);
2994 cpl_error_set_message(
"enu_load_mef_components",
2995 CPL_ERROR_UNSUPPORTED_MODE,
2996 "extension has unsupported HDUCLAS1: %s %s", name, hduclas1);
2998 enu_check_error_code(
"error accessing HDU %d", hdu);
3000 result->mef[i] =
new;
3006 cpl_propertylist_delete(eplist);
3007 if (cpl_error_get_code() != CPL_ERROR_NONE) {
3010 cpl_propertylist_delete(*plist);
3032 cpl_image_delete(limage->confidence);
3034 cpl_image_delete(limage->bkg_confidence);
3035 cpl_propertylist_delete(limage->plist);
3037 cpl_mask_delete(limage->object_mask);
3040 cpl_frame_delete(limage->frame);
3058 cpl_ensure(limage, CPL_ERROR_NULL_INPUT, NULL);
3062 hdrl_image * himage = NULL;
3063 if (limage->himage) {
3066 hdrl_imagelist * himagelist = NULL;
3067 if (limage->himagelist) {
3070 cpl_image * confidence = NULL;
3071 if (limage->confidence) {
3072 confidence = cpl_image_duplicate(limage->confidence);
3074 hdrl_image * bkg = NULL;
3078 cpl_image * bkg_confidence = NULL;
3079 if (limage->bkg_confidence) {
3080 bkg_confidence = cpl_image_duplicate(limage->bkg_confidence);
3082 cpl_propertylist * plist = NULL;
3083 if (limage->plist) {
3084 plist = cpl_propertylist_duplicate(limage->plist);
3086 hdrl_catalogue_result * objects = NULL;
3087 if (limage->objects) {
3090 cpl_mask * object_mask = NULL;
3091 if (limage->object_mask) {
3092 object_mask = cpl_mask_duplicate(limage->object_mask);
3094 hdrl_catalogue_result * matchstd_wcs = NULL;
3095 if (limage->matchstd_wcs) {
3097 limage->matchstd_wcs);
3099 hdrl_catalogue_result * matchstd_phot = NULL;
3100 if (limage->matchstd_phot) {
3102 limage->matchstd_phot);
3104 cpl_frame * frame = NULL;
3105 if (limage->frame) {
3106 frame = cpl_frame_duplicate(limage->frame);
3149 hdrl_imagelist * himagelist,
3150 cpl_image * confidence,
3152 cpl_image * bkg_confidence,
3153 cpl_propertylist * plist,
3154 hdrl_catalogue_result * objects,
3155 cpl_mask * object_mask,
3156 hdrl_catalogue_result * wcs,
3157 hdrl_catalogue_result * photom,
3158 cpl_frame * frame) {
3160 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
3163 cpl_ensure(!(himage && himagelist), CPL_ERROR_ILLEGAL_INPUT, NULL);
3165 located_image * result = cpl_malloc(
sizeof(located_image));
3166 result->himage = himage;
3167 result->himagelist = himagelist;
3168 result->confidence = confidence;
3170 result->bkg_confidence = bkg_confidence;
3171 result->plist = plist;
3172 result->objects = objects;
3173 result->object_mask = object_mask;
3174 result->matchstd_wcs = wcs;
3175 result->matchstd_phot = photom;
3176 result->frame = frame;
3196 cpl_frameset * frameset,
3198 cpl_frameset * used) {
3200 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
3201 cpl_ensure(frameset, CPL_ERROR_NULL_INPUT, NULL);
3202 cpl_ensure(tag, CPL_ERROR_NULL_INPUT, NULL);
3203 cpl_ensure(used, CPL_ERROR_NULL_INPUT, NULL);
3205 located_imagelist * result = NULL;
3206 cpl_frameset_iterator * frameset_iter = NULL;
3207 mef_extension_list * mefs = NULL;
3208 cpl_image * data = NULL;
3209 cpl_image * err = NULL;
3210 cpl_image * qual = NULL;
3211 cpl_image * conf = NULL;
3212 cpl_image * bkg = NULL;
3213 cpl_image * bkg_err = NULL;
3214 cpl_image * bkg_conf = NULL;
3215 cpl_propertylist * plist = NULL;
3219 frameset_iter = cpl_frameset_iterator_new(frameset);
3222 for (cpl_frame * frame = NULL;
3223 (frame = cpl_frameset_iterator_get(frameset_iter)) &&
3224 (cpl_error_get_code() == CPL_ERROR_NONE);
3225 cpl_frameset_iterator_advance(frameset_iter, 1)) {
3227 const char * frametag = cpl_frame_get_tag(frame);
3228 if (!strcmp(tag, frametag)) {
3236 cpl_frameset_iterator_reset(frameset_iter);
3237 cpl_size position = 0;
3238 for (cpl_frame * frame = NULL;
3239 (frame = cpl_frameset_iterator_get(frameset_iter)) &&
3240 (cpl_error_get_code() == CPL_ERROR_NONE);
3241 cpl_frameset_iterator_advance(frameset_iter, 1)) {
3246 const char * frametag = cpl_frame_get_tag(frame);
3248 if (!strcmp(tag, frametag)) {
3491 cpl_frame * dup_frame = cpl_frame_duplicate(frame);
3492 cpl_frameset_insert(used, dup_frame);
3500 cpl_propertylist_delete(plist);
3501 cpl_image_delete(data);
3502 cpl_image_delete(err);
3503 cpl_image_delete(qual);
3504 cpl_image_delete(conf);
3505 cpl_image_delete(bkg);
3506 cpl_image_delete(bkg_err);
3507 cpl_image_delete(bkg_conf);
3509 cpl_frameset_iterator_delete(frameset_iter);
3510 if (cpl_error_get_code() != CPL_ERROR_NONE) {
3534 const char * dps_name,
3535 const char * hduclas1,
3536 const char * hduclas2) {
3538 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
3539 cpl_ensure(mef, CPL_ERROR_NULL_INPUT, NULL);
3540 cpl_ensure(dps_name, CPL_ERROR_NULL_INPUT, NULL);
3541 cpl_ensure(hduclas1, CPL_ERROR_NULL_INPUT, NULL);
3542 cpl_ensure(hduclas2, CPL_ERROR_NULL_INPUT, NULL);
3544 cpl_image * result = NULL;
3546 cpl_propertylist * plist = mef->plist;
3547 const char * hduclass = cpl_propertylist_get_string(plist,
"HDUCLASS");
3548 const char * hdudoc = cpl_propertylist_get_string(plist,
"HDUDOC");
3550 const char * my_hduclas1 = cpl_propertylist_get_string(plist,
"HDUCLAS1");
3551 const char * my_hduclas2 = cpl_propertylist_get_string(plist,
"HDUCLAS2");
3553 enu_check(strstr(hduclass,
"ESO") != NULL &&
3554 strstr(hdudoc,
"SDP") != NULL &&
3556 strstr(my_hduclas1, hduclas1) != NULL &&
3557 strstr(my_hduclas2, hduclas2) != NULL,
3558 CPL_ERROR_BAD_FILE_FORMAT,
"bad %s format", dps_name);
3563 if (cpl_error_get_code() != CPL_ERROR_NONE) {
3566 if (result != NULL) mef->data = NULL;
3585 for (cpl_size i = 0; i < limlist->size; i++) {
3588 cpl_free(limlist->limages);
3605 const located_imagelist * limlist) {
3607 cpl_ensure(limlist, CPL_ERROR_NULL_INPUT, NULL);
3610 for (cpl_size i = 0; i < limlist->size; i++) {
3613 limlist->limages[i]), i);
3634 located_image * limage,
3635 cpl_size position) {
3637 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
3641 cpl_ensure_code(position >= 0 && position < limlist->size,
3642 CPL_ERROR_ACCESS_OUT_OF_RANGE);
3643 limlist->limages[position] = limage;
3645 return cpl_error_get_code();
3662 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
3664 located_imagelist * limlist = cpl_malloc(
sizeof(located_imagelist));
3665 limlist->size = size;
3666 limlist->limages = cpl_calloc(size,
sizeof(located_image));
3667 for (cpl_size i = 0; i < limlist->size; i++) {
3668 limlist->limages[i] = NULL;
3687 if (!strcmp(mef->data_type, MEF_EXTENSION_CONTAINING_IMAGE)) {
3688 cpl_image_delete((cpl_image *) (mef->data));
3689 }
else if (!strcmp(mef->data_type, MEF_EXTENSION_CONTAINING_MASK)) {
3690 cpl_mask_delete((cpl_mask *) mef->data);
3691 }
else if (!strcmp(mef->data_type, MEF_EXTENSION_CONTAINING_TABLE)) {
3692 cpl_table_delete((cpl_table *) mef->data);
3693 }
else if (!strcmp(mef->data_type, MEF_EXTENSION_CONTAINING_VECTOR)) {
3694 cpl_vector_delete((cpl_vector *) mef->data);
3696 cpl_error_set(
"enu_mef_extension_delete",
3697 CPL_ERROR_UNSUPPORTED_MODE);
3699 cpl_propertylist_delete(mef->plist);
3700 cpl_free(mef->name);
3719 for (cpl_size i = 0; i < list->size; i++) {
3722 cpl_free(list->mef);
3741 const char * target) {
3743 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
3745 cpl_ensure(meflist, CPL_ERROR_NULL_INPUT, NULL);
3746 cpl_ensure(target, CPL_ERROR_NULL_INPUT, NULL);
3750 cpl_mask * result = NULL;
3751 for (
int i=0; i<meflist->size; i++) {
3752 if (!strcmp(meflist->mef[i]->name, target)) {
3753 cpl_ensure(!strcmp(meflist->mef[i]->data_type,
3754 MEF_EXTENSION_CONTAINING_MASK),
3755 CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
3756 result = cpl_mask_duplicate((
const cpl_mask *) meflist->mef[i]->data);
3763 cpl_ensure(result, CPL_ERROR_DATA_NOT_FOUND, NULL);
3782 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
3784 mef_extension_list * mef_list = cpl_malloc(
sizeof(mef_extension_list));
3785 mef_list->size = size;
3786 mef_list->mef = cpl_calloc(size,
sizeof(mef_extension));
3787 for (cpl_size i = 0; i < size; i++) {
3788 mef_list->mef[i] = NULL;
3812 const cpl_image * data,
3813 const cpl_propertylist * plist) {
3815 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
3819 cpl_ensure(name && data, CPL_ERROR_NULL_INPUT, NULL);
3821 mef_extension * result = cpl_malloc(
sizeof(mef_extension));
3822 result->name = cpl_strdup(name);
3823 result->data = (
void *) cpl_image_duplicate(data);
3824 result->data_type = MEF_EXTENSION_CONTAINING_IMAGE;
3826 result->plist = cpl_propertylist_duplicate(plist);
3828 result->plist = NULL;
3851 const cpl_mask * data,
3852 const cpl_propertylist * plist) {
3854 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
3858 cpl_ensure(name && data, CPL_ERROR_NULL_INPUT, NULL);
3860 mef_extension * result = cpl_malloc(
sizeof(mef_extension));
3861 result->name = cpl_strdup(name);
3862 result->data = (
void *) cpl_mask_duplicate(data);
3863 result->data_type = MEF_EXTENSION_CONTAINING_MASK;
3865 result->plist = cpl_propertylist_duplicate(plist);
3867 result->plist = NULL;
3890 const cpl_table * table,
3891 const cpl_propertylist * plist) {
3893 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
3894 cpl_ensure(name, CPL_ERROR_NULL_INPUT, NULL);
3895 cpl_ensure(table, CPL_ERROR_NULL_INPUT, NULL);
3897 mef_extension * result = cpl_malloc(
sizeof(mef_extension));
3898 result->name = cpl_strdup(name);
3899 result->data = (
void *) cpl_table_duplicate(table);
3900 result->plist = NULL;
3902 result->plist = cpl_propertylist_duplicate(plist);
3904 result->plist = NULL;
3906 result->data_type = MEF_EXTENSION_CONTAINING_TABLE;
3929 const cpl_vector * vector,
3930 const cpl_propertylist * plist) {
3932 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
3933 cpl_ensure(name, CPL_ERROR_NULL_INPUT, NULL);
3934 cpl_ensure(vector, CPL_ERROR_NULL_INPUT, NULL);
3936 mef_extension * result = cpl_malloc(
sizeof(mef_extension));
3937 result->name = cpl_strdup(name);
3938 result->data = (
void *) cpl_vector_duplicate(vector);
3939 result->plist = NULL;
3941 result->plist = cpl_propertylist_duplicate(plist);
3943 result->plist = NULL;
3945 result->data_type = MEF_EXTENSION_CONTAINING_VECTOR;
3966 const char * filename,
3967 const cpl_propertylist * plist,
3970 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
3972 cpl_msg_debug(cpl_func,
"enu_mef_extension_save entry");
3974 cpl_propertylist * plist_copy = NULL;
3978 cpl_ensure_code(mef, CPL_ERROR_NULL_INPUT);
3979 cpl_ensure_code(filename, CPL_ERROR_NULL_INPUT);
3980 cpl_ensure_code(plist, CPL_ERROR_NULL_INPUT);
3982 plist_copy = cpl_propertylist_duplicate(plist);
3983 cpl_propertylist_update_string(plist_copy,
"ERIS_NIX_MEF_TYPE",
3986 if (!strcmp(mef->data_type, MEF_EXTENSION_CONTAINING_IMAGE)) {
3987 if (cpl_propertylist_has(plist,
"CD3_3")) {
3998 cpl_msg_info(cpl_func,
"..saving data with NDIM=3");
4000 cpl_imagelist * imlist = cpl_imagelist_new();
4001 cpl_imagelist_set(imlist, mef->data, 0);
4002 cpl_imagelist_save(imlist,
4004 CPL_TYPE_UNSPECIFIED,
4008 cpl_imagelist_unset(imlist, 0);
4009 cpl_imagelist_delete(imlist);
4012 cpl_image_save((
const cpl_image *) (mef->data),
4014 CPL_TYPE_UNSPECIFIED,
4018 }
else if (!strcmp(mef->data_type, MEF_EXTENSION_CONTAINING_MASK)) {
4019 cpl_mask_save((
const cpl_mask *) (mef->data),
4023 }
else if (!strcmp(mef->data_type, MEF_EXTENSION_CONTAINING_TABLE)) {
4024 cpl_table_save((
const cpl_table *) (mef->data),
4029 }
else if (!strcmp(mef->data_type, MEF_EXTENSION_CONTAINING_VECTOR)) {
4030 cpl_vector_save((
const cpl_vector *) (mef->data),
4036 cpl_error_set_message(
"enu_mef_extension_save",
4037 CPL_ERROR_UNSUPPORTED_MODE,
4038 "unsupported extension type: %s", mef->data_type);
4041 cpl_propertylist_delete(plist_copy);
4043 return cpl_error_get_code();
4064 const double obj_threshold,
4065 const int bkg_mesh_size,
4066 const double bkg_smooth_fwhm,
4067 located_imagelist * limlist) {
4069 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
4073 cpl_ensure_code(limlist, CPL_ERROR_NULL_INPUT);
4075 casu_fits ** indata = NULL;
4076 casu_fits ** inconf = NULL;
4077 casu_fits ** invar = NULL;
4086 int casu_status = CASU_OK;
4089 cpl_msg_info(cpl_func,
".. calculating object masks");
4090 for (cpl_size j = 0; j < limlist->size; j++) {
4091 casu_opm(indata[j], inconf[j], obj_min_pixels, obj_threshold,
4092 bkg_mesh_size, bkg_smooth_fwhm, niter, &casu_status);
4093 enu_check(casu_status == CASU_OK, CPL_ERROR_UNSPECIFIED,
4094 "CASU error calculating object mask");
4095 limlist->limages[j]->object_mask = cpl_mask_duplicate(
4096 cpl_image_get_bpm(indata[j]->image));
4101 for (cpl_size j = 0; j < limlist->size; j++) {
4102 casu_fits_delete(indata[j]);
4103 casu_fits_delete(inconf[j]);
4104 casu_fits_delete(invar[j]);
4110 return cpl_error_get_code();
4133 const int nsigma_cut) {
4135 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
4136 cpl_ensure_code(limlist, CPL_ERROR_NULL_INPUT);
4138 cpl_msg_info(cpl_func,
"..calculating LSS object pixel mask");
4143 for (cpl_size j = 0; j < limlist->size; j++) {
4146 cpl_image * collapsed = cpl_image_collapse_median_create(
4148 limlist->limages[j]->himage), 0, 0, 0);
4152 double csigma = 0.0;
4153 double cmedian = cpl_image_get_median_dev(collapsed, &csigma);
4154 double maxval = cpl_image_get_max(collapsed);
4155 cpl_mask * collapsed_opm = cpl_mask_threshold_image_create(
4156 collapsed, cmedian + nsigma_cut * csigma,
4161 cpl_mask * kernel = cpl_mask_new(11, 1);
4162 cpl_mask_not(kernel);
4163 cpl_mask * dilated_collapsed_opm = cpl_mask_new(nx, 1);
4164 cpl_mask_filter(dilated_collapsed_opm, collapsed_opm, kernel,
4165 CPL_FILTER_DILATION, CPL_BORDER_NOP);
4166 cpl_mask_delete(kernel);
4170 const cpl_binary * collapsed_opm_data = cpl_mask_get_data_const(
4171 dilated_collapsed_opm);
4172 cpl_mask * opm = cpl_mask_new(nx, ny);
4173 cpl_binary * opm_data = cpl_mask_get_data(opm);
4174 for (cpl_size iy=0; iy < ny; iy++) {
4175 memcpy(opm_data + iy * (nx *
sizeof(*opm_data)),
4177 nx *
sizeof(*opm_data));
4182 cpl_mask_delete(limlist->limages[j]->object_mask);
4183 limlist->limages[j]->object_mask = opm;
4185 cpl_image_delete(collapsed);
4186 cpl_mask_delete(collapsed_opm);
4187 cpl_mask_delete(dilated_collapsed_opm);
4190 return cpl_error_get_code();
4210 const char * wcs_method,
4211 const char * catalogue,
4212 located_image * limage,
4213 const double match_rad,
4214 cpl_table ** matched_stds,
4215 cpl_matrix ** xy_shift) {
4217 cpl_matrix * celestial = NULL;
4218 cpl_array * colnames = NULL;
4219 cpl_array * conv_status = NULL;
4220 cpl_table * refcat_copy = NULL;
4221 cpl_matrix * ref_phys = NULL;
4222 cpl_matrix * ref_world = NULL;
4223 cpl_matrix * xy = NULL;
4224 cpl_wcs * wcs = NULL;
4226 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
4228 cpl_msg_info(cpl_func,
"correcting %s", cpl_frame_get_filename(limage->frame));
4232 cpl_msg_info(cpl_func,
" ..%s", cpl_frame_get_filename(limage->frame));
4236 cpl_size nref = cpl_table_get_nrow(refcat);
4237 const double * ra = cpl_table_get_data_double_const(refcat,
"RA");
4238 const double * dec = NULL;
4239 if (cpl_table_has_column(refcat,
"DEC")) {
4240 dec = cpl_table_get_data_double_const(refcat,
"DEC");
4242 ref_world = cpl_matrix_new(nref, 2);
4243 for (cpl_size iobj = 0; iobj < nref; iobj++) {
4244 cpl_matrix_set(ref_world, iobj, 0, ra[iobj]);
4245 cpl_matrix_set(ref_world, iobj, 1, dec[iobj]);
4247 wcs = cpl_wcs_new_from_propertylist(limage->plist);
4248 cpl_wcs_convert(wcs, ref_world, &ref_phys, &conv_status,
4249 CPL_WCS_WORLD2PHYS);
4250 enu_check_error_code(
"failure calculating predicted positions");
4256 refcat_copy = cpl_table_duplicate(refcat);
4257 colnames = cpl_table_get_column_names(refcat_copy);
4258 cpl_size ncol = cpl_array_get_size(colnames);
4259 for (cpl_size ic = 0; ic < ncol; ic++) {
4260 const char * colname = cpl_array_get_string(colnames, ic);
4261 if (strcmp(colname,
"RA") &&
4262 strcmp(colname,
"DEC") &&
4263 strcmp(colname,
"Aper_flux_3") &&
4264 strcmp(colname,
"Std index")) {
4265 cpl_table_erase_column(refcat_copy, colname);
4271 cpl_table_new_column(refcat_copy,
"xpredict", CPL_TYPE_DOUBLE);
4272 cpl_table_new_column(refcat_copy,
"ypredict", CPL_TYPE_DOUBLE);
4273 for (cpl_size iobj = 0; iobj < nref; iobj++) {
4274 cpl_table_set_double(refcat_copy,
"xpredict", iobj,
4275 cpl_matrix_get(ref_phys, iobj, 0));
4276 cpl_table_set_double(refcat_copy,
"ypredict", iobj,
4277 cpl_matrix_get(ref_phys, iobj, 1));
4279 enu_check_error_code(
"failure while creating predicted positions");
4283 cpl_msg_info(cpl_func,
" ..matching sources");
4285 if (limage->objects && limage->objects->catalogue) {
4286 nobj = cpl_table_get_nrow(limage->objects->catalogue);
4289 if (min(nref, nobj) > 0) {
4291 cpl_msg_info(cpl_func,
" ..calling enm_try_all_associations");
4297 enu_check_error_code(
"enu_correct_wcs failure at "
4298 "try_all_associations");
4305 enu_check_error_code(
"failure to correct wcs");
4308 cpl_msg_info(cpl_func,
" ..no catalogued objects, cannot correct wcs");
4312 cpl_matrix_delete(celestial);
4313 cpl_array_delete(colnames);
4314 cpl_array_delete(conv_status);
4315 cpl_table_delete(refcat_copy);
4316 cpl_matrix_delete(ref_phys);
4317 cpl_matrix_delete(ref_world);
4318 cpl_matrix_delete(xy);
4319 cpl_wcs_delete(wcs);
4321 return cpl_error_get_code();
4338 const cpl_table * refine_wcs) {
4340 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
4341 cpl_ensure_code(limage, CPL_ERROR_NULL_INPUT);
4343 if (!refine_wcs)
return CPL_ERROR_NONE;
4345 const char * camera = cpl_propertylist_get_string(limage->plist,
4346 "ESO INS2 NXCW NAME");
4347 enu_check(camera != NULL,
4348 CPL_ERROR_ILLEGAL_INPUT,
4349 "failed to read NIX camera from ESO.INS2.NXCW.NAME");
4351 cpl_size nrow = cpl_table_get_nrow(refine_wcs);
4352 cpl_boolean match = CPL_FALSE;
4359 for (cpl_size row = 0; row < nrow; row++) {
4360 const char * row_cam = cpl_table_get_string(refine_wcs,
"camera", row);
4362 if (!strcmp(camera, row_cam)) {
4366 cd1_1 = cpl_table_get_double(refine_wcs,
"CD1_1", row, &ignore);
4367 cd1_2 = cpl_table_get_double(refine_wcs,
"CD1_2", row, &ignore);
4368 cd2_1 = cpl_table_get_double(refine_wcs,
"CD2_1", row, &ignore);
4369 cd2_2 = cpl_table_get_double(refine_wcs,
"CD2_2", row, &ignore);
4375 cpl_propertylist_update_double(limage->plist,
"CD1_1", cd1_1);
4376 cpl_propertylist_update_double(limage->plist,
"CD1_2", cd1_2);
4377 cpl_propertylist_update_double(limage->plist,
"CD2_1", cd2_1);
4378 cpl_propertylist_update_double(limage->plist,
"CD2_2", cd2_2);
4383 cpl_error_set_message(cpl_func,
4384 CPL_ERROR_ILLEGAL_INPUT,
4385 "camera name not matched: %s",
4388 cd1_1 = cpl_propertylist_get_double(limage->plist,
"CD1_1");
4389 cd1_2 = cpl_propertylist_get_double(limage->plist,
"CD1_2");
4390 cd2_1 = cpl_propertylist_get_double(limage->plist,
"CD2_1");
4391 cd2_2 = cpl_propertylist_get_double(limage->plist,
"CD2_2");
4392 cpl_msg_info(cpl_func,
"..set CD matrix cd1_1=%10.8e cd1_2=%10.8e "
4393 "cd_2_1=%10.8e cd2_2=%10.8e",
4394 cd1_1, cd1_2, cd2_1, cd2_2);
4399 return cpl_error_get_code();
4417 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
4418 cpl_ensure_code(confidence, CPL_ERROR_NULL_INPUT);
4420 const cpl_size nx = cpl_image_get_size_x(confidence);
4421 const cpl_size ny = cpl_image_get_size_y(confidence);
4422 double * data = cpl_image_get_data_double(confidence);
4424 const cpl_mask * bpm = cpl_image_get_bpm_const(confidence);
4425 cpl_mask * mask = bpm ? cpl_image_get_bpm(confidence) : NULL;
4426 cpl_binary * mask_data = bpm ? cpl_mask_get_data(mask) : NULL;
4431 for (cpl_size i = 0; i < nx*ny; i++) {
4432 if (mask_data != NULL && mask_data[i] == CPL_TRUE) {
4433 mask_data[i] = CPL_FALSE;
4436 if (data[i] < 0.0) data[i] = 0.0;
4438 if (data[i] > 0.0) nonzero++;
4442 cpl_image_multiply_scalar(confidence, 100.0 * (
double)nonzero / sum);
4445 return cpl_error_get_code();
4487 const cpl_propertylist * plist,
4488 cpl_image * confidence,
4489 const int set_confidence) {
4493 double * hdata = cpl_image_get_data_double(himg);
4494 double * edata = cpl_image_get_data_double(eimg);
4495 cpl_binary* hbpm = cpl_mask_get_data(cpl_image_get_bpm(himg));
4496 cpl_binary* ebpm = cpl_mask_get_data(cpl_image_get_bpm(eimg));
4497 const int mx = cpl_image_get_size_x(himg);
4498 cpl_boolean is_fast_uncorr = CPL_FALSE;
4499 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
4502 cpl_ensure_code(hdata != NULL, CPL_ERROR_UNSUPPORTED_MODE);
4503 cpl_ensure_code(edata != NULL, CPL_ERROR_UNSUPPORTED_MODE);
4504 cpl_ensure_code(hbpm != NULL, CPL_ERROR_UNSUPPORTED_MODE);
4505 cpl_ensure_code(ebpm != NULL, CPL_ERROR_UNSUPPORTED_MODE);
4507 cpl_ensure_code(plist, CPL_ERROR_NULL_INPUT);
4516 cpl_size nx_chip = 0;
4517 cpl_size ny_chip = 0;
4518 cpl_boolean windowed = 0;
4528 enu_check_error_code(
"failed to read detector window information");
4534 double mjd = cpl_propertylist_get_double(plist,
"MJD-OBS");
4535 enu_check_error_code(
"failed to read MJD-OBS");
4537 cpl_msg_info(cpl_func,
"PIPE-10417 DITDELAY bug, "
4538 "remove_read_offsets will ignore bottom 4 rows "
4541 const char * curr_name = cpl_propertylist_get_string(plist,
"ESO DET READ CURNAME");
4542 if (strcmp(curr_name,
"FAST_UNCORR") == 0) {
4543 is_fast_uncorr = CPL_TRUE;
4545 cpl_size top_width = 4;
4546 cpl_size bottom_width = 4;
4547 cpl_msg_debug(cpl_func,
4548 "..offsets from %d top rows and %d bottom rows of chip",
4549 (
int) top_width, (
int) bottom_width);
4551 cpl_size nread_chan = nx_chip / 64;
4555 cpl_vector* voff = cpl_vector_new(nread_chan);
4557 double* pvoff = cpl_vector_get_data(voff);
4558 for (cpl_size read_chan = 0; read_chan < nread_chan; read_chan++) {
4562 cpl_array * offset_data = cpl_array_new(0, CPL_TYPE_DOUBLE);
4565 for (cpl_size x = 64 * read_chan; x < 64 * (read_chan+1); x++) {
4570 if (x+1 >= strx && x+1 < strx+nx) {
4572 for (cpl_size y = 0; y < bottom_width; y++) {
4576 if (y+1 >= stry && y+1 < stry+ny) {
4577 const int ipos = x-strx+1 + mx * (y-stry+1);
4581 const double pixval = hdata[ipos];
4585 cpl_array_set_size(offset_data, ndata);
4586 cpl_array_set(offset_data, ndata-1, pixval);
4589 hbpm[ipos] = CPL_BINARY_1;
4590 ebpm[ipos] = CPL_BINARY_1;
4592 if (set_confidence && confidence) {
4593 cpl_image_set(confidence,
4602 for (cpl_size y = ny_chip-top_width; y < ny_chip; y++) {
4603 if (y+1 >= stry && y+1 < stry+ny) {
4604 const int ipos = x-strx+1 + mx * (y-stry+1);
4607 const double pixval = hdata[ipos];
4609 cpl_array_set_size(offset_data, ndata);
4610 cpl_array_set(offset_data, ndata-1, pixval);
4611 hbpm[ipos] = CPL_BINARY_1;
4612 ebpm[ipos] = CPL_BINARY_1;
4613 if (set_confidence && confidence) {
4614 cpl_image_set(confidence,
4627 hdrl_value offset = {0.0, 0.0};
4629 cpl_msg_debug(cpl_func,
"read_chan %d", (
int)read_chan);
4630 cpl_msg_debug(cpl_func,
"ndata %d", (
int)ndata);
4631 cpl_msg_debug(cpl_func,
"median %f", cpl_array_get_median(offset_data));
4632 cpl_msg_debug(cpl_func,
"mean %f", cpl_array_get_mean(offset_data));
4633 cpl_msg_debug(cpl_func,
"std %f", cpl_array_get_stdev(offset_data));
4638 double mean = cpl_array_get_mean(offset_data);
4639 double stdev = cpl_array_get_stdev(offset_data);
4641 for (cpl_size i = 0; i < cpl_array_get_size(offset_data); i++) {
4643 double datum = cpl_array_get_double(offset_data, i, &invalid);
4644 if (!invalid && fabs(datum - mean) > 5 * stdev) {
4645 cpl_array_set_invalid(offset_data, i);
4650 offset.data = cpl_array_get_mean(offset_data);
4651 offset.error = cpl_array_get_stdev(offset_data);
4652 pvoff[read_chan] = offset.data;
4654 cpl_msg_debug(cpl_func,
"..removing offset %f (%f) for read "
4655 "channel %d", offset.data, offset.error,
4658 for (cpl_size x = 64 * read_chan; x < 64 * (read_chan+1); x++) {
4659 for (cpl_size y = 0; y < ny_chip; y++) {
4661 if (x+1 >= strx && x+1 < strx+nx &&
4662 y+1 >= stry && y+1 < stry+ny) {
4663 const int ipos = x-strx+1 + mx * (y-stry+1);
4666 hdata[ipos] -= offset.data;
4667 edata[ipos] = hypot(edata[ipos], offset.error);
4673 cpl_array_delete(offset_data);
4678 if(is_fast_uncorr && ndata > 0) {
4679 cpl_image_add_scalar(himg,cpl_vector_get_mean(voff));
4684 for (cpl_size y = 0; y < ny_chip; y++) {
4685 for (cpl_size x = 0; x < 4; x++) {
4687 if (x+1 >= strx && x+1 < strx+nx &&
4688 y+1 >= stry && y+1 < stry+ny) {
4689 const int ipos = x-strx+1 + mx * (y-stry+1);
4691 hbpm[ipos] = CPL_BINARY_1;
4692 ebpm[ipos] = CPL_BINARY_1;
4694 if (set_confidence && confidence) {
4695 cpl_image_set(confidence,
4704 for (cpl_size y = 0; y < ny_chip; y++) {
4705 for (cpl_size x = nx_chip-4; x < nx_chip; x++) {
4707 if (x+1 >= strx && x+1 < strx+nx &&
4708 y+1 >= stry && y+1 < stry+ny) {
4713 if (set_confidence && confidence) {
4714 cpl_image_set(confidence,
4723 cpl_vector_delete(voff);
4726 if (set_confidence && confidence) {
4732 return cpl_error_get_code();
4746 const char * preface) {
4748 if (cpl_error_get_code() != CPL_ERROR_NONE)
return NULL;
4750 char * in_name = cpl_strdup(filename);
4753 char * pch = strstr(filename,
"ERIS.20");
4755 pch = strstr(filename,
"NACO.20");
4758 char * outname = cpl_sprintf(
"%s.%s", preface, pch);
4784 const char * select_method,
4785 const double timerange,
4786 const located_imagelist * target_list,
4787 const located_imagelist * sky_list,
4788 const cpl_size x_probe,
4789 const cpl_size y_probe,
4790 hdrl_imagelist ** sky_himagelist,
4791 cpl_imagelist ** sky_confidence_list) {
4793 if (cpl_error_get_code() != CPL_ERROR_NONE)
return cpl_error_get_code();
4795 cpl_image * contrib = NULL;
4796 hdrl_imagelist * skylist = NULL;
4797 cpl_image * sky_confidence = NULL;
4798 hdrl_image * sky_himage = NULL;
4799 cpl_vector * sky_indeces = NULL;
4801 cpl_msg_debug(cpl_func,
"enu_sky_backgrounds entry");
4805 cpl_ensure_code(method, CPL_ERROR_NULL_INPUT);
4806 cpl_ensure_code(select_method, CPL_ERROR_NULL_INPUT);
4807 cpl_ensure_code(target_list, CPL_ERROR_NULL_INPUT);
4808 cpl_ensure_code(sky_list, CPL_ERROR_NULL_INPUT);
4809 cpl_ensure_code(!strcmp(select_method,
"bracket") &&
4810 (timerange > 0.0), CPL_ERROR_ILLEGAL_INPUT);
4813 *sky_confidence_list = cpl_imagelist_new();
4817 int debug = CPL_FALSE;
4818 cpl_vector * data = NULL;
4819 cpl_vector * error = NULL;
4820 cpl_vector * reject = NULL;
4821 cpl_vector * confidence = NULL;
4823 cpl_msg_info(cpl_func,
"..calculating skys");
4827 cpl_size njitters = target_list->size;
4828 for (cpl_size j = 0; j < njitters; j++) {
4832 if (strstr(select_method,
"bracket") != NULL) {
4834 sky_list, x_probe > -1);
4836 enu_check_error_code(
"failure in enu_bracket_skys");
4838 cpl_size skysize = cpl_vector_get_size(sky_indeces);
4840 if (x_probe > 0 && skysize > 0) {
4842 sky_list->limages[0]->himage);
4844 sky_list->limages[0]->himage);
4845 if (x_probe > 0 && x_probe <= x_size &&
4846 y_probe > 0 && y_probe <= y_size) {
4849 data = cpl_vector_new(skysize);
4850 error = cpl_vector_new(skysize);
4851 reject = cpl_vector_new(skysize);
4852 confidence = cpl_vector_new(skysize);
4857 if (strstr(method,
"collapse-median") != NULL) {
4864 for (
int k = 0; k < (int)skysize; k++) {
4865 int i = (int)cpl_vector_get(sky_indeces, k);
4867 limages[i]->himage);
4868 cpl_image * confcopy = cpl_image_duplicate(sky_list->
4869 limages[i]->confidence);
4870 if (sky_list->limages[i]->object_mask != NULL) {
4872 sky_list->limages[i]->object_mask);
4873 cpl_mask_or(cpl_image_get_bpm(confcopy),
4874 sky_list->limages[i]->object_mask);
4875 cpl_image_fill_rejected(confcopy, 0.0);
4876 cpl_image_accept_all(confcopy);
4882 x_probe, y_probe, &rejected);
4884 double conf = cpl_image_get(confcopy,
4887 cpl_vector_set(data, k, himval.data);
4888 cpl_vector_set(error, k, himval.error);
4889 cpl_vector_set(reject, k, (
double)rejected);
4890 cpl_vector_set(confidence, k, conf);
4896 if (sky_confidence == NULL) {
4897 sky_confidence = confcopy;
4899 cpl_image_add(sky_confidence, confcopy);
4900 cpl_image_delete(confcopy);
4903 cpl_msg_debug(cpl_func,
"after skylist %s",
4904 cpl_error_get_message());
4910 cpl_image_delete(contrib);
4912 cpl_msg_debug(cpl_func,
"after collapse %s",
4913 cpl_error_get_message());
4921 ".. enu_sky_backgrounds: probe pixel %d,%d",
4922 (
int)x_probe, (
int)y_probe);
4923 cpl_msg_info(
" ",
".. jitter %d" , (
int)j);
4924 cpl_msg_info(
" ",
"... j# q c d (e)");
4926 for (cpl_size k = 0; k < skysize; k++) {
4927 cpl_msg_info(
" ",
"... %3d %3d %3d %8.2f(%4.2f)",
4928 (
int)cpl_vector_get(sky_indeces, k),
4929 (
int)cpl_vector_get(reject, k),
4930 (
int)cpl_vector_get(confidence, k),
4931 cpl_vector_get(data, k),
4932 cpl_vector_get(error, k));
4937 x_probe, y_probe, &rejected);
4938 double cval = cpl_image_get(sky_confidence, x_probe,
4940 cpl_msg_info(
" ",
"... -> %3d %3d %8.2f(%4.2f)",
4941 rejected, (
int)cval, hval.data, hval.error);
4943 cpl_vector_delete(data);
4944 cpl_vector_delete(error);
4945 cpl_vector_delete(reject);
4946 cpl_vector_delete(confidence);
4950 enu_check_error_code(
"failure in collapse-median");
4952 }
else if (strstr(method,
"median-median") != NULL) {
4957 cpl_array * sky_medians = cpl_array_new(0, CPL_TYPE_DOUBLE);
4959 for (
int k = 0; k < (int)skysize; k++) {
4960 int i = (int)cpl_vector_get(sky_indeces, k);
4965 sky_list->limages[i]->himage);
4966 if (sky_list->limages[i]->object_mask != NULL) {
4968 sky_list->limages[i]->object_mask);
4973 cpl_size asize = cpl_array_get_size(sky_medians);
4974 cpl_array_set_size(sky_medians, asize+1);
4975 cpl_array_set_double(sky_medians, asize, bkg.data);
4978 cpl_vector_set(data, k, bkg.data);
4979 cpl_vector_set(error, k, bkg.error);
4980 cpl_vector_set(confidence, k, 100.0);
4981 cpl_vector_set(reject, k, 0.0);
4988 double sky_med_med = cpl_array_get_median(sky_medians);
4989 double sky_med_stdev = cpl_array_get_stdev(sky_medians);
4992 cpl_image * sky_data = cpl_image_new(nx, ny, HDRL_TYPE_DATA);
4993 cpl_image_fill_window(sky_data, 1, 1, nx, ny, sky_med_med);
4994 cpl_image * sky_error = cpl_image_new(nx, ny, HDRL_TYPE_ERROR);
4995 cpl_image_fill_window(sky_error, 1, 1, nx, ny, sky_med_stdev);
4997 cpl_image_delete(sky_data);
4998 cpl_image_delete(sky_error);
5000 sky_confidence = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
5001 cpl_image_fill_window(sky_confidence, 1, 1, nx, ny, 100.0);
5003 cpl_array_delete(sky_medians);
5007 ".. enu_sky_backgrounds: median of images");
5008 cpl_msg_info(
" ",
".. jitter %d" , (
int)j);
5009 cpl_msg_info(
" ",
"... j# c d (e)");
5011 for (cpl_size k = 0; k < skysize; k++) {
5012 cpl_msg_info(
" ",
"... %3d %3d %8.2f(%4.2f)",
5013 (
int)cpl_vector_get(sky_indeces, k),
5014 (
int)cpl_vector_get(confidence, k),
5015 cpl_vector_get(data, k),
5016 cpl_vector_get(error, k));
5022 double cval = cpl_image_get(sky_confidence, 1, 1, &ignore);
5023 cpl_msg_info(
" ",
"... -> %3d %3d %8.2f(%4.2f)",
5024 rejected, (
int)cval, hval.data, hval.error);
5026 cpl_vector_delete(data);
5027 cpl_vector_delete(error);
5028 cpl_vector_delete(confidence);
5030 enu_check_error_code(
"failure in median-median");
5037 cpl_imagelist_set(*sky_confidence_list, sky_confidence,
5038 cpl_imagelist_get_size(*sky_confidence_list));
5039 cpl_vector_delete(sky_indeces); sky_indeces = NULL;
5044 cpl_vector_delete(sky_indeces);
5045 if (cpl_error_get_code() != CPL_ERROR_NONE) {
5047 *sky_himagelist = NULL;
5048 cpl_imagelist_delete(*sky_confidence_list);
5049 *sky_confidence_list = NULL;
5052 return cpl_error_get_code();
5071 const double stk_lthr,
const double stk_hthr,
5072 const char * stk_method,
const int stk_fast,
5073 located_image ** result) {
5075 cpl_ensure_code(cpl_error_get_code() == CPL_ERROR_NONE,
5076 cpl_error_get_code());
5078 casu_fits ** indata = NULL;
5079 casu_fits ** inconf = NULL;
5080 casu_fits ** invar = NULL;
5081 casu_fits * out = NULL;
5082 casu_fits * outc = NULL;
5083 casu_fits * outv = NULL;
5087 cpl_ensure_code(limages, CPL_ERROR_NULL_INPUT);
5088 cpl_ensure_code((stk_fast == 0) || (stk_fast == 1),
5089 CPL_ERROR_ILLEGAL_INPUT);
5093 const int casu_stk_method = (!strcmp(stk_method,
"nearest") ? 0 : 1);
5101 int njitters = (int) limages->size;
5102 int casu_status = CASU_OK;
5103 casu_imstack(indata, inconf, invar, NULL, (
int)njitters, (int)njitters,
5104 stk_lthr, stk_hthr, casu_stk_method, 0, stk_fast, 0,
5105 "ESO DET DIT", &out, &outc, &outv, &casu_status);
5108 cpl_image_duplicate(outc->image),
5111 NULL, NULL, NULL, NULL, NULL);
5116 return cpl_error_get_code();
5136 const char * select_method,
5137 const double timerange,
5138 const located_imagelist * sky_data,
5139 const cpl_size x_probe,
5140 const cpl_size y_probe,
5141 located_imagelist * target_data) {
5143 cpl_ensure_code(cpl_error_get_code() == CPL_ERROR_NONE,
5144 cpl_error_get_code());
5145 cpl_ensure_code(method, CPL_ERROR_NULL_INPUT);
5146 cpl_ensure_code(select_method, CPL_ERROR_NULL_INPUT);
5147 cpl_ensure_code(target_data, CPL_ERROR_NULL_INPUT);
5149 cpl_imagelist * sky_conf_list = NULL;
5150 hdrl_imagelist * sky_himlist = NULL;
5152 cpl_msg_info(cpl_func,
5153 "enu_sky_subtract called: method=%s, sky_selector=%s, "
5155 method, select_method, timerange);
5156 cpl_size njitters = target_data->size;
5157 cpl_msg_info(cpl_func,
"number of jitters to be reduced %d",
5163 sky_data, x_probe, y_probe, &sky_himlist,
5165 enu_check_error_code(
"failure in enu_sky_backgrounds");
5169 hdrl_value before = {0.0, 0.0};
5171 hdrl_value bkg = {0.0, 0.0};
5172 double bkgconf = 0.0;
5173 hdrl_value after = {0.0, 0.0};
5176 int bkgrejected = 1;
5179 int header_printed = CPL_FALSE;
5181 for (cpl_size i = 0; i < njitters; i++) {
5186 target_data->limages[i]->himage);
5188 target_data->limages[i]->himage);
5189 if (x_probe > 0 && x_probe <= x_size &&
5190 y_probe > 0 && y_probe <= y_size) {
5192 x_probe, y_probe, &brejected);
5193 bconf = cpl_image_get(target_data->limages[i]->confidence,
5194 x_probe, y_probe, &ignore);
5203 target_data->limages[i]->bkg_confidence = cpl_image_duplicate(
5204 cpl_imagelist_get(sky_conf_list, i));
5223 cpl_image_fill_rejected(target_data->limages[i]->confidence, 0.0);
5224 cpl_image_accept_all(target_data->limages[i]->confidence);
5229 if (x_probe > 0 && x_probe <= x_size &&
5230 y_probe > 0 && y_probe <= y_size) {
5232 x_probe, y_probe, &arejected);
5233 aconf = cpl_image_get(target_data->limages[i]->confidence,
5234 x_probe, y_probe, &ignore);
5236 x_probe, y_probe, &bkgrejected);
5237 bkgconf = cpl_image_get(target_data->limages[i]->bkg_confidence,
5238 x_probe, y_probe, &ignore);
5240 if (!header_printed) {
5241 header_printed = CPL_TRUE;
5242 cpl_msg_info(
" ",
"..Subtracting sky backgrounds: probe pixel %d,%d",
5243 (
int)x_probe, (
int)y_probe);
5244 cpl_msg_info(
" ",
" j# before sky after");
5245 cpl_msg_info(
" ",
" q c d e q c d e q c d e");
5248 "..%3d %d %3.0f %4.2e(%4.2e)) %d %3.0f %4.2e(%4.2e) %d %3.0f %10.2e(%4.2e)",
5250 brejected, bconf, before.data, before.error,
5251 bkgrejected, bkgconf, bkg.data, bkg.error,
5252 arejected, aconf, after.data, after.error);
5257 cpl_imagelist_delete(sky_conf_list);
5260 return cpl_error_get_code();
5277 cpl_propertylist* qc_list,
5280 cpl_ensure_code(mbad_pix_map != NULL, CPL_ERROR_NULL_INPUT);
5281 cpl_ensure_code(qc_list != NULL, CPL_ERROR_NULL_INPUT);
5282 cpl_ensure_code(prefix != NULL, CPL_ERROR_NULL_INPUT);
5285 cpl_mask * bpm_mask = NULL;
5287 flag_mask = ~flag_mask;
5288 bpm_mask = en_master_bpm_get_mask(mbad_pix_map, flag_mask);
5289 int nbad = cpl_mask_count(bpm_mask);
5290 cpl_propertylist_append_int(qc_list,
"ESO QC NUMBER BAD PIXELS", nbad);
5291 cpl_size sx = cpl_mask_get_size_x(bpm_mask);
5292 cpl_size sy = cpl_mask_get_size_y(bpm_mask);
5293 double fraction = (double) nbad / (sx * sy);
5294 cpl_propertylist_append_double(qc_list,
"ESO QC FRACTION BAD PIXELS",
fraction);
5295 return cpl_error_get_code();
5315 cpl_mask* bp_map_nl_mask,
5316 const cpl_parameterlist* parlist,
5317 const char* context,
5318 const double threshold_pos,
5319 const cpl_boolean verbose,
5320 const cpl_boolean rescale_by_dit){
5322 cpl_ensure(lamp_on_limlist != NULL, CPL_ERROR_NULL_INPUT, NULL);
5323 cpl_ensure(bp_map_nl_mask != NULL, CPL_ERROR_NULL_INPUT, NULL);
5324 cpl_ensure(parlist != NULL, CPL_ERROR_NULL_INPUT, NULL);
5325 cpl_ensure(context != NULL, CPL_ERROR_NULL_INPUT, NULL);
5328 cpl_propertylist* qc;
5337 qc = cpl_propertylist_new();
5338 double threshold_neg = - 4.3e7;
5339 threshold_neg = 300.;
5340 long nthresh_pos = 0;
5341 long nthresh_neg = 0;
5343 const cpl_parameter* param;
5347 pname = cpl_sprintf(
"%s%s",context,
".saturation_neg");
5348 cpl_msg_warning(cpl_func,
"pname: %s",pname);
5349 param = cpl_parameterlist_find_const(parlist,pname);
5350 threshold_neg = cpl_parameter_get_double(param);
5356 double mean_min = DBL_MAX;
5357 double median_min = DBL_MAX;
5359 double mean_max = -DBL_MIN;
5360 double median_max = -DBL_MIN;
5361 double nthresh_pos_max = 0;
5362 double nthresh_neg_max = 0;
5370 cpl_propertylist* plist = NULL;
5371 for(cpl_size i = 0; i < lamp_on_limlist->size; i++) {
5379 cpl_mask_or(mask, bp_map_nl_mask);
5383 plist = lamp_on_limlist->limages[i]->plist;
5384 if(cpl_propertylist_has(plist,
"ESO DET SEQ1 DIT")) {
5385 dit = cpl_propertylist_get_double(plist,
"ESO DET SEQ1 DIT");
5387 nthresh_neg = eris_image_get_threshpix(ima, threshold_neg, CPL_FALSE);
5388 nthresh_pos = eris_image_get_threshpix(ima, threshold_pos, CPL_TRUE);
5390 eris_image_flag_threshpix(&ima, threshold_neg, CPL_FALSE);
5391 mean = cpl_image_get_mean(ima);
5392 median = cpl_image_get_median(ima);
5393 rms = cpl_image_get_stdev(ima);
5398 if(rescale_by_dit) {
5404 key_name = cpl_sprintf(
"ESO QC FLAT_ON%lld NPOSSAT",i);
5405 key_comm = cpl_sprintf(
"Number of flat pixels above threshold");
5406 cpl_propertylist_append_int(qc, key_name, nthresh_pos);
5407 cpl_propertylist_set_comment(qc, key_name, key_comm);
5411 key_name = cpl_sprintf(
"ESO QC FLAT_ON%lld NNEGSAT",i);
5412 key_comm = cpl_sprintf(
"Number of flat pixels below neg threshold");
5413 cpl_propertylist_append_int(qc, key_name, nthresh_neg);
5414 cpl_propertylist_set_comment(qc, key_name, key_comm);
5418 key_name = cpl_sprintf(
"ESO QC FLAT_ON%lld MEAN",i);
5419 key_comm = cpl_sprintf(
"[ADU] Mean of flat");
5420 cpl_propertylist_append_double(qc, key_name, mean);
5421 cpl_propertylist_set_comment(qc,key_name,key_comm);
5425 key_name = cpl_sprintf(
"ESO QC FLAT_ON%lld MEDIAN",i);
5426 key_comm = cpl_sprintf(
"[ADU] Median of flat");
5427 cpl_propertylist_append_double(qc, key_name, median);
5428 cpl_propertylist_set_comment(qc,key_name,key_comm);
5432 key_name = cpl_sprintf(
"ESO QC FLAT_ON%lld STDEV",i);
5433 key_comm = cpl_sprintf(
"[ADU] Stdev of flat");
5434 cpl_propertylist_append_double(qc, key_name, rms);
5435 cpl_propertylist_set_comment(qc,key_name,key_comm);
5440 if(mean > mean_max) {
5443 if(median > median_max) {
5444 median_max = median;
5448 if(mean < mean_min) {
5451 if(median < median_min) {
5452 median_min = median;
5456 if(nthresh_pos > nthresh_pos_max) {
5457 nthresh_pos_max = nthresh_pos;
5459 if(nthresh_neg > nthresh_neg_max) {
5460 nthresh_neg_max = nthresh_neg;
5463 cpl_msg_info(cpl_func,
"mean_min: %g",mean_min);
5464 cpl_msg_info(cpl_func,
"mean_max: %g",mean_max);
5465 cpl_msg_info(cpl_func,
"median_min: %g",median_min);
5466 cpl_msg_info(cpl_func,
"median_max: %g",median_max);
5467 cpl_msg_info(cpl_func,
"f1: %g",mean_min / mean_max);
5468 cpl_msg_info(cpl_func,
"f2: %g",median_min / median_max);
5470 key_name = cpl_sprintf(
"ESO QC FLAT_ON MEAN MIN");
5471 key_comm = cpl_sprintf(
"[ADU] Min of Means of flat");
5472 cpl_propertylist_append_double(qc, key_name, mean_min);
5473 cpl_propertylist_set_comment(qc,key_name,key_comm);
5477 key_name = cpl_sprintf(
"ESO QC FLAT_ON MEDIAN MIN");
5478 key_comm = cpl_sprintf(
"[ADU] Min of Medians of flat");
5479 cpl_propertylist_append_double(qc, key_name, median_min);
5480 cpl_propertylist_set_comment(qc,key_name,key_comm);
5484 key_name = cpl_sprintf(
"ESO QC FLAT_ON MEAN MAX");
5485 key_comm = cpl_sprintf(
"[ADU] Max of Means of flat");
5486 cpl_propertylist_append_double(qc, key_name, mean_max);
5487 cpl_propertylist_set_comment(qc,key_name,key_comm);
5491 key_name = cpl_sprintf(
"ESO QC FLAT_ON MEDIAN MAX");
5492 key_comm = cpl_sprintf(
"[ADU] Max of Medians of flat");
5493 cpl_propertylist_append_double(qc, key_name, median_max);
5494 cpl_propertylist_set_comment(qc,key_name,key_comm);
5499 key_name = cpl_sprintf(
"ESO QC FLAT_ON MEAN FRAC");
5500 key_comm = cpl_sprintf(
"Min/Max of Means of flat");
5501 cpl_propertylist_append_double(qc, key_name, mean_min / mean_max);
5502 cpl_propertylist_set_comment(qc,key_name,key_comm);
5507 key_name = cpl_sprintf(
"ESO QC FLAT_ON MEDIAN FRAC");
5508 key_comm = cpl_sprintf(
"Min/Max of Medians of flat");
5509 cpl_propertylist_append_double(qc, key_name, median_min / median_max);
5510 cpl_propertylist_set_comment(qc,key_name,key_comm);
5516 key_name = cpl_sprintf(
"ESO QC FLAT_ON NPOSSAT MAX");
5517 key_comm = cpl_sprintf(
"Max of Number of flat pixels above threshold");
5518 cpl_propertylist_append_int(qc, key_name, nthresh_pos_max);
5519 cpl_propertylist_set_comment(qc, key_name, key_comm);
5523 key_name = cpl_sprintf(
"ESO QC FLAT_ON NNEGSAT MAX");
5524 key_comm = cpl_sprintf(
"MAX Number of flat pixels below neg threshold");
5525 cpl_propertylist_append_int(qc, key_name, nthresh_neg_max);
5526 cpl_propertylist_set_comment(qc, key_name, key_comm);
cpl_error_code encu_limlist_to_casu_fits(located_imagelist *limlist, casu_fits ***indata, casu_fits ***inconf, casu_fits ***invar)
Translate a located_imagelist to arrays of casu_fits structs.
cpl_error_code enu_dfs_save_limage(cpl_frameset *allframes, const cpl_parameterlist *parlist, const cpl_frameset *provenance, const cpl_boolean prov_raw, const located_image *limage, const char *recipe, const cpl_frame *inherit, cpl_propertylist *applist, const char *pipe_id, const char *filename)
Save a located image structure to a MEF.
cpl_error_code enu_dfs_save_himage(cpl_frameset *allframes, const cpl_parameterlist *parlist, const cpl_frameset *provenance, const cpl_boolean prov_raw, const hdrl_image *image, const hdrl_imagelist *imagelist, const mef_extension_list *mefs, const char *recipe, const cpl_frame *inherit_frame, const cpl_propertylist *applist, const cpl_propertylist *wcs_plist, const char *pipe_id, const char *filename)
Save an hdrl_image/imagelist as a DFS-compliant MEF pipeline product.
cpl_error_code enm_correct_crpix(const char *wcs_method, const char *catalogue, cpl_table *matched_stds, located_image *limage, cpl_matrix **xy_shift)
Apply median offsets from matched standards table to image.
cpl_error_code enm_try_all_associations(cpl_table *objtab, cpl_table *stdtab, const double match_rad, cpl_table **matchtab)
Match standards with objects by trying all possibilities and selecting best.
located_image * enu_located_image_duplicate(const located_image *limage)
Make a deep copy of a located_image and its contents.
cpl_mask * enu_mef_extension_list_get_mask(mef_extension_list *meflist, const char *target)
Get a cpl_mask from a named mef_extension in a list.
cpl_vector * enu_bracket_skys(const located_image *target, const double timerange, const located_imagelist *pool, const int debug)
Find images taken within a given time of a target image.
cpl_error_code enu_modify_CD_matrix(located_image *limage, const cpl_table *refine_wcs)
Update the CD matrix to reduce image distortion and rotation.
void enu_located_imagelist_delete(located_imagelist *limlist)
Delete a located_imagelist and its contents.
cpl_error_code enu_located_imagelist_insert(located_imagelist *limlist, located_image *limage, cpl_size position)
Insert a located_image at a specified point in a located_imagelist.
void enu_mef_extension_delete(mef_extension *mef)
Delete a mef_extension and its contents.
cpl_error_code eris_nix_get_badpix_qc_from_ima(const master_bpm *mbad_pix_map, cpl_propertylist *qc_list, const char *prefix)
Compute QC on input master bpm.
mef_extension * enu_mef_new_mask(const char *name, const cpl_mask *data, const cpl_propertylist *plist)
Create a mef_extension struct holding a cpl_mask.
cpl_propertylist * enu_raw_flats_qc(located_imagelist *lamp_on_limlist, cpl_mask *bp_map_nl_mask, const cpl_parameterlist *parlist, const char *context, const double threshold_pos, const cpl_boolean verbose, const cpl_boolean rescale_by_dit)
Compute QC on input raw flats.
cpl_error_code enu_flat_save(const char *pro_catg, const hdrl_image *flat, const cpl_image *confidence, const cpl_mask *cold_bpm, cpl_frameset *frameset, const cpl_parameterlist *parlist, const char *filename, const char *recipe_name, const cpl_propertylist *qclog)
Save a flatfield result.
int enu_check_conformance(const cpl_propertylist *plist1, const cpl_propertylist *plist2, const char *regexp)
Check that the specified subset of two propertylists match.
located_image * enu_load_limage_from_frame(const cpl_frame *frame, cpl_image **pcopyconf, const cpl_boolean collapse_cube)
Load components of a located_image from a frame.
hdrl_catalogue_result * enu_hdrl_catalogue_result_duplicate(const hdrl_catalogue_result *target)
Return a deep copy of the given hdrl_catalogue_result object.
located_imagelist * enu_located_imagelist_duplicate(const located_imagelist *limlist)
Make a deep copy of a located_imagelist and its contents.
cpl_error_code enu_calc_maglim(const located_image *limage, const double photzp, const double fwhm_pix, double *abmaglim)
Calculate magnitude limit of image.
cpl_error_code enu_opm_lss_limlist(located_imagelist *limlist, const int nsigma_cut)
Calculate object masks for LSS images in a located_imagelist.
double enu_get_tel_alt(const cpl_propertylist *plist)
Get telescope altitude of an observation.
mef_extension * enu_mef_new_image(const char *name, const cpl_image *data, const cpl_propertylist *plist)
Create a mef_extension to hold a cpl_image.
cpl_error_code enu_get_ra_dec(const cpl_wcs *wcs, double *ra, double *dec)
Get RA and Dec at centre of image with given wcs.
hdrl_image * enu_calc_flat(hdrl_imagelist *himlist, const int min_coadds, const hdrl_parameter *collapse_params, const cpl_size filter_size_x, const cpl_size filter_size_y, const hdrl_flat_method method)
Calculate a flatfield result.
cpl_error_code enu_get_window_info(cpl_size *nx, cpl_size *ny, int *rot, cpl_size *strx, cpl_size *stry, cpl_size *nx_chip, cpl_size *ny_chip, cpl_boolean *windowed, const cpl_propertylist *plist)
Get the detector 'window' information.
cpl_error_code enu_basic_calibrate(located_image *limage, const int read_offsets, const cpl_table *refine_wcs, const master_dark *mdark, const gain_linearity *gain_lin, const master_flat *flatfield_1, const master_flat *flatfield_2, const master_bpm *mbad_pix_map, const int flag_mask, const char *fill_rejected, const double fill_value, const cpl_size x_probe, const cpl_size y_probe)
Do basic calibration of located_image (single or cube)
cpl_error_code enu_catalogue_limlist(located_imagelist *limlist, hdrl_parameter *params)
Calculate object catalogues for a list of images.
cpl_error_code enu_sky_subtract_limlist(const char *method, const char *select_method, const double timerange, const located_imagelist *sky_data, const cpl_size x_probe, const cpl_size y_probe, located_imagelist *target_data)
Estimate and subtract sky backgrounds for a list of target images.
cpl_error_code enu_himage_load_from_fits(const char *filename, hdrl_image **result, mef_extension_list **mef_extensions, cpl_propertylist **plist)
Load an hdrl_image from a multi-extension FITS file.
cpl_error_code enu_debug_limlist_save(const int debug, const located_imagelist *limlist, const char *nameroot, const char *recipename, cpl_frameset *frameset, const cpl_parameterlist *parlist, const cpl_frameset *used)
Save a list of intermediate image results for use in debugging.
mef_extension * enu_mef_new_table(const char *name, const cpl_table *table, const cpl_propertylist *plist)
Create a mef_extension struct holding a cpl_table.
mef_extension_list * enu_mef_extension_list_new(cpl_size size)
Construct a new mef_extension_list.
located_image * enu_located_image_new(hdrl_image *himage, hdrl_imagelist *himagelist, cpl_image *confidence, hdrl_image *bkg, cpl_image *bkg_confidence, cpl_propertylist *plist, hdrl_catalogue_result *objects, cpl_mask *object_mask, hdrl_catalogue_result *wcs, hdrl_catalogue_result *photom, cpl_frame *frame)
Create a located_image structure and initialise the contents.
double enu_get_airmass(const cpl_propertylist *plist)
Get the mean airmass of an observation.
cpl_error_code enu_calc_pixel_coords(cpl_table *catalogue, const cpl_propertylist *wcs_plist)
Calculate predicted positions of catalogue objects for given wcs.
located_imagelist * enu_limlist_load_from_frameset(cpl_frameset *frameset, const char *tag, cpl_frameset *used)
Load tagged data from a frameset into a located_imagelist.
void enu_located_image_delete(located_image *limage)
Delete a located_image and its contents.
mef_extension * enu_mef_new_vector(const char *name, const cpl_vector *vector, const cpl_propertylist *plist)
Create a mef_extension struct holding a cpl_vector.
cpl_error_code enu_remove_read_offsets(hdrl_image *himage, const cpl_propertylist *plist, cpl_image *confidence, const int set_confidence)
Function to remove read offsets from an himage and set blank pixels
cpl_error_code enu_normalise_confidence(cpl_image *confidence)
Normalise confidence array so that mean of good pixels is 100.
cpl_error_code enu_stack(located_imagelist *limages, const double stk_lthr, const double stk_hthr, const char *stk_method, const int stk_fast, located_image **result)
Stack a set of images.
cpl_error_code enu_correct_wcs(const cpl_table *refcat, const char *wcs_method, const char *catalogue, located_image *limage, const double match_rad, cpl_table **matched_stds, cpl_matrix **xy_shift)
Correct the wcs of an image.
cpl_error_code enu_basic_calibrate_himage(hdrl_image *himage, cpl_image *confidence, cpl_propertylist *plist, const cpl_frame *frame, const int read_offsets, const master_dark *mdark, const gain_linearity *gain_lin, const master_flat *flatfield_1, const master_flat *flatfield_2, const master_bpm *mbad_pix_map, const int set_confidence, const int flag_mask, const char *fill_rejected, const double fill_value, const cpl_size x_probe, const cpl_size y_probe)
Do basic calibration of an hdrl_image.
cpl_error_code enu_mef_extension_save(const mef_extension *mef, const char *filename, const cpl_propertylist *plist, unsigned mode)
Save a mef_extension struct to FITS file.
char * enu_repreface(const char *filename, const char *preface)
Preface a raw filename with a string.
const char * enu_get_license(void)
Get the pipeline copyright and license.
void enu_mef_extension_list_delete(mef_extension_list *list)
Delete a mef_extension_list and its contents.
cpl_image * enu_load_component(mef_extension *mef, const char *dps_name, const char *hduclas1, const char *hduclas2)
Load data from a MEF component of specified type.
const char * enu_get_filter(const cpl_propertylist *plist)
Get the filter used in an observation.
cpl_error_code enu_check_wcs(const located_image *limage)
Check that the image has valid WCS keywords.
cpl_error_code enu_get_rcore_and_mesh_size(const char *context, const cpl_parameterlist *parlist, const cpl_propertylist *plist, double *obj_core_radius, int *bkg_mesh_size)
Get catalogue core-radius and mesh-size appropriate to AO mode.
located_imagelist * enu_located_imagelist_new(cpl_size size)
Return a pointer to a new located_imagelist.
hdrl_catalogue_result * enu_catalogue_compute(const hdrl_image *himage, const cpl_image *confidence, const cpl_wcs *wcs, hdrl_parameter *params)
Wrapper for hdrl_catalogue_compute.
cpl_error_code enu_opm_limlist(const int obj_min_pixels, const double obj_threshold, const int bkg_mesh_size, const double bkg_smooth_fwhm, located_imagelist *limlist)
Calculate object masks for images in a located_imagelist.
const char * enu_get_det_mode(const cpl_propertylist *plist)
Get the detector mode of an integration.
cpl_error_code enu_sky_backgrounds(const char *method, const char *select_method, const double timerange, const located_imagelist *target_list, const located_imagelist *sky_list, const cpl_size x_probe, const cpl_size y_probe, hdrl_imagelist **sky_himagelist, cpl_imagelist **sky_confidence_list)
Get sky backgrounds for a list of target images.
double enu_get_dit(const cpl_propertylist *plist)
Get the DIT of an integration.
double enu_get_filter_wavelength(const char *filter)
Get the effective wavelength of a filter.
mef_extension_list * enu_load_mef_components(const char *filename, cpl_propertylist **plist)
Load components of a multi-extension FITS file.
cpl_error_code eris_check_error_code(const char *func_id)
handle CPL errors
hdrl_catalogue_result * hdrl_catalogue_compute(const cpl_image *image_, const cpl_image *confidence_map, const cpl_wcs *wcs, hdrl_parameter *param_)
build object catalog
void hdrl_catalogue_result_delete(hdrl_catalogue_result *result)
delete hdrl parameter result object
hdrl_parameter * hdrl_collapse_mode_parameter_create(double histo_min, double histo_max, double bin_size, hdrl_mode_type mode_method, cpl_size error_niter)
create a parameter object for the mode
hdrl_parameter * hdrl_flat_parameter_create(cpl_size filter_size_x, cpl_size filter_size_y, hdrl_flat_method method)
Creates FLAT Parameters object.
cpl_error_code hdrl_flat_compute(hdrl_imagelist *hdrl_data, const cpl_mask *stat_mask, const hdrl_parameter *collapse_params, hdrl_parameter *flat_params, hdrl_image **master, cpl_image **contrib_map)
compute high or low frequency master flat with median filtering
hdrl_value hdrl_image_get_pixel(const hdrl_image *self, cpl_size xpos, cpl_size ypos, int *pis_rejected)
get pixel values of hdrl_image
cpl_error_code hdrl_image_sub_image(hdrl_image *self, const hdrl_image *other)
Subtract two images, store the result in the first image.
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.
cpl_error_code hdrl_image_div_image(hdrl_image *self, const hdrl_image *other)
Divide two images, store the result in the first image.
hdrl_image * hdrl_image_duplicate(const hdrl_image *himg)
copy hdrl_image
double hdrl_image_get_stdev(const hdrl_image *self)
computes the standard deviation of the data of an image
cpl_error_code hdrl_image_copy(hdrl_image *dst, const hdrl_image *src, cpl_size xpos, cpl_size ypos)
Copy one image into another.
cpl_mask * hdrl_image_get_mask(hdrl_image *himg)
get cpl bad pixel mask from image
hdrl_image * hdrl_image_extract(const hdrl_image *self, cpl_size llx, cpl_size lly, cpl_size urx, cpl_size ury)
extract copy of window from image
cpl_image * hdrl_image_get_error(hdrl_image *himg)
get error as cpl image
hdrl_value hdrl_image_get_mean(const hdrl_image *self)
computes mean pixel value and associated error of an 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
const cpl_image * hdrl_image_get_image_const(const hdrl_image *himg)
get data as cpl image
void hdrl_image_delete(hdrl_image *himg)
delete hdrl_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.
cpl_error_code hdrl_imagelist_collapse(const hdrl_imagelist *himlist, const hdrl_parameter *param, hdrl_image **out, cpl_image **contrib)
collapsing of image list
hdrl_imagelist * hdrl_imagelist_new(void)
Create an empty imagelist.
cpl_error_code hdrl_imagelist_collapse_median(const hdrl_imagelist *himlist, hdrl_image **out, cpl_image **contrib)
Median collapsing of image list.
hdrl_imagelist * hdrl_imagelist_duplicate(const hdrl_imagelist *himlist)
Duplicate an image list.
cpl_size hdrl_imagelist_get_size_x(const hdrl_imagelist *himlist)
Get number of colums of images in the imagelist.
cpl_error_code hdrl_imagelist_collapse_weighted_mean(const hdrl_imagelist *himlist, hdrl_image **out, cpl_image **contrib)
Weighted Mean collapsing of image list.
hdrl_image * hdrl_imagelist_get(const hdrl_imagelist *himlist, cpl_size inum)
Get an image from a list of images.
cpl_error_code hdrl_maglim_compute(const cpl_image *image, const double zeropoint, const double fwhm, const cpl_size kernel_size_x, const cpl_size kernel_size_y, const hdrl_image_extend_method image_extend_method, const hdrl_parameter *mode_parameter, double *limiting_magnitude)
Computes the limiting magnitude of an image.
void hdrl_parameter_delete(hdrl_parameter *obj)
shallow delete of a parameter
double fraction(double x, double y, double r_out)
Fraction of pixel bounded.