36 #include "muse_exp_align_z.h" 39 typedef struct muse_indexpair muse_indexpair;
41 struct muse_indexpair {
46 typedef void (*muse_free_func)(
void *);
50 static const double deg2as = 3600.;
58 muse_vfree(
void **array, cpl_size size, muse_free_func deallocator)
62 for (idx = 0; idx < size; ++idx) {
64 deallocator(array[idx]);
73 _muse_condition_less_than(
double aValue1,
double aValue2)
75 return (aValue1 < aValue2) ? TRUE : FALSE;
94 muse_utils_replace_nan(
muse_image *aImage,
float aValue)
96 cpl_ensure(aImage, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
97 cpl_ensure(aImage->
data, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_ILLEGAL_INPUT);
99 cpl_size npixel = cpl_image_get_size_x(aImage->
data) *
100 cpl_image_get_size_y(aImage->
data);
102 float *data = cpl_image_get_data_float(aImage->
data);
103 if (cpl_error_get_code() == CPL_ERROR_TYPE_MISMATCH) {
104 cpl_error_set_message(__func__, CPL_ERROR_TYPE_MISMATCH,
105 "MUSE image data buffer has invalid type!");
106 return CPL_ERROR_TYPE_MISMATCH;
110 for (ipixel = 0; ipixel < npixel; ++ipixel) {
111 if (isnan(data[ipixel])) {
112 data[ipixel] = aValue;
117 npixel = cpl_image_get_size_x(aImage->
stat) *
118 cpl_image_get_size_y(aImage->
stat);
120 data = cpl_image_get_data_float(aImage->
stat);
121 if (cpl_error_get_code() == CPL_ERROR_TYPE_MISMATCH) {
122 cpl_error_set_message(__func__, CPL_ERROR_TYPE_MISMATCH,
123 "MUSE image data buffer has invalid type!");
124 return CPL_ERROR_TYPE_MISMATCH;
127 for (ipixel = 0; ipixel < npixel; ++ipixel) {
128 if (isnan(data[ipixel])) {
129 data[ipixel] = aValue;
134 return CPL_ERROR_NONE;
161 image->
header = cpl_propertylist_load(aFilename, 0);
163 cpl_error_set_message(__func__, cpl_error_get_code(),
"Loading primary FITS " 164 "header of \"%s\" did not succeed", aFilename);
170 int extension = cpl_fits_find_extension(aFilename, EXTNAME_DATA);
171 cpl_propertylist *hdata = cpl_propertylist_load(aFilename, extension);
175 cpl_msg_debug(__func__,
"Skipping extension %d [%s]", extension,
177 cpl_propertylist_delete(hdata);
179 hdata = cpl_propertylist_load(aFilename, extension);
181 cpl_msg_debug(__func__,
"Taking extension %d [%s]", extension,
186 image->
data = cpl_image_load(aFilename, CPL_TYPE_FLOAT, 0, extension);
188 cpl_error_set_message(__func__, MUSE_ERROR_READ_DATA,
"Could not load " 189 "extension %s from \"%s\"", extname, aFilename);
195 if (cpl_propertylist_has(hdata,
"BUNIT")) {
196 cpl_propertylist_append_string(image->
header,
"BUNIT",
198 cpl_propertylist_set_comment(image->
header,
"BUNIT",
199 cpl_propertylist_get_comment(hdata,
"BUNIT"));
201 cpl_msg_warning(__func__,
"No BUNIT given in extension %d [%s] of \"%s\"!",
202 extension, extname, aFilename);
205 if (!cpl_propertylist_has(hdata,
"CUNIT1") ||
206 !cpl_propertylist_has(hdata,
"CUNIT2")) {
207 cpl_msg_warning(__func__,
"No WCS found in extension %d [%s] of \"%s\"!",
208 extension, extname, aFilename);
213 cpl_propertylist_erase_regexp(hdata,
"(^ESO |" MUSE_WCS_KEYS ")", 1);
214 cpl_propertylist_append(image->
header, hdata);
215 cpl_propertylist_delete(hdata);
219 if (extname && !strncmp(extname, EXTNAME_DATA, strlen(EXTNAME_DATA)+1)) {
221 extension = cpl_fits_find_extension(aFilename, EXTNAME_STAT);
224 char *statname = cpl_sprintf(
"%s_STAT", extname);
225 extension = cpl_fits_find_extension(aFilename, statname);
229 cpl_errorstate status = cpl_errorstate_get();
230 image->
stat = cpl_image_load(aFilename, CPL_TYPE_INT, 0, extension);
231 if (!cpl_errorstate_is_equal(status)) {
232 if (cpl_error_get_code() == CPL_ERROR_DATA_NOT_FOUND) {
233 cpl_errorstate_set(status);
234 cpl_msg_debug(__func__,
"Ignoring empty extension %s in \"%s\"",
235 EXTNAME_STAT, aFilename);
237 cpl_error_set_message(__func__, MUSE_ERROR_READ_STAT,
"Could not load " 238 "extension %s from \"%s\"", EXTNAME_STAT, aFilename);
248 if (extname && !strncmp(extname, EXTNAME_DATA, strlen(EXTNAME_DATA)+1)) {
250 extension = cpl_fits_find_extension(aFilename, EXTNAME_DQ);
253 char *dqname = cpl_sprintf(
"%s_DQ", extname);
254 extension = cpl_fits_find_extension(aFilename, dqname);
258 cpl_errorstate status = cpl_errorstate_get();
259 image->
dq = cpl_image_load(aFilename, CPL_TYPE_INT, 0, extension);
260 if (!cpl_errorstate_is_equal(status)) {
261 cpl_errorstate_set(status);
262 cpl_error_set_message(__func__, MUSE_ERROR_READ_DQ,
"Could not load " 263 "extension %s from \"%s\"", EXTNAME_DQ, aFilename);
287 cpl_ensure(aProcessing, CPL_ERROR_NULL_INPUT, NULL);
288 cpl_size nframes = cpl_frameset_get_size(aProcessing->
inframes);
289 cpl_ensure(nframes, CPL_ERROR_DATA_NOT_FOUND, NULL);
293 cpl_size iexposure = 0;
295 for (iframe = 0; iframe < nframes; ++iframe) {
296 const cpl_frame *frame = cpl_frameset_get_position(aProcessing->
inframes,
298 const char *tag = cpl_frame_get_tag(frame);
301 const char *filename = cpl_frame_get_filename(frame);
303 cpl_msg_debug(__func__,
"Loading FOV image '%s' as exposure %" 304 CPL_SIZE_FORMAT, filename, iexposure + 1);
308 cpl_msg_error(__func__,
"Could not load FOV image '%s'", filename);
315 CPL_FRAME_GROUP_RAW, 1);
334 cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
335 "Detection threshold must be greater than zero!");
338 if (aParams->
fwhm < 0.5) {
339 cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
340 "FWHM must be greater than 0.5 pixels!");
344 cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
345 "Maximum number of sources must be grater than the " 346 "minimum number of sources!");
350 cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
351 "Upper limit of the point-source roundness must be " 352 "greater than the lower limit!");
356 cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
357 "Low limit of the point-source sharpness must be " 358 "greater than zero!");
362 cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
363 "Upper limit of the point-source sharpness must be " 364 "greater than the lower limit!");
380 muse_align_parse_valuelist(
const char *aValuelist)
382 cpl_ensure(aValuelist, CPL_ERROR_NULL_INPUT, NULL);
383 if (strlen(aValuelist) == 0) {
384 cpl_error_set_message(__func__, CPL_ERROR_DATA_NOT_FOUND,
385 "List of values is empty!");
390 cpl_size nvalues = cpl_array_get_size(strings);
392 cpl_array *values = cpl_array_new(nvalues, CPL_TYPE_DOUBLE);
394 for (ivalue = 0; ivalue < nvalues; ++ivalue) {
395 const char *svalue = cpl_array_get_string(strings, ivalue);
396 if ((strncasecmp(svalue,
"inf", 3) == 0) ||
397 (strncasecmp(svalue,
"nan", 3) == 0)) {
398 cpl_error_set_message(__func__, CPL_ERROR_ILLEGAL_INPUT,
399 "Illegal value \"%s\" encountered!", svalue);
400 cpl_array_delete(values);
401 cpl_array_delete(strings);
406 double value = strtod(svalue, &last);
407 if( (errno == ERANGE) || ((value == 0.) && (last == svalue))) {
408 cpl_error_set_message(__func__, CPL_ERROR_TYPE_MISMATCH,
409 "Could not convert \"%s\" to a decimal number!", svalue);
410 cpl_array_delete(values);
411 cpl_array_delete(strings);
414 cpl_array_set_double(values, ivalue, value);
417 cpl_array_delete(strings);
432 muse_align_celestial_from_pixel(cpl_table *aTable,
433 cpl_propertylist *aWCSHeader)
440 if (unit1 && unit2) {
441 if ((strncmp(unit1, unit2, 4) == 0) && (strncmp(unit1,
"deg", 3) == 0)) {
448 cpl_errorstate status = cpl_errorstate_get();
450 cpl_table_new_column(aTable,
"RA", CPL_TYPE_DOUBLE);
451 cpl_table_new_column(aTable,
"DEC", CPL_TYPE_DOUBLE);
454 for (irow = 0; irow < cpl_table_get_nrow(aTable); ++irow) {
455 double x = cpl_table_get_double(aTable,
"X", irow, NULL);
456 double y = cpl_table_get_double(aTable,
"Y", irow, NULL);
461 cpl_table_set_double(aTable,
"RA", irow, ra);
462 cpl_table_set_double(aTable,
"DEC", irow, dec);
466 if (!cpl_errorstate_is_equal(status)) {
475 #define USE_SVD_SOLVER 476 #ifdef USE_SVD_SOLVER 478 #if CPL_VERSION_CODE < 459008 479 #define cpl_matrix_solve_svd muse_cplmatrix_solve_svd 522 static cpl_error_code
523 muse_cplmatrix_decomp_sv_jacobi(cpl_matrix *U, cpl_vector *S, cpl_matrix *V)
525 cpl_ensure_code((U != NULL) && (V != NULL) && (S != NULL),
526 CPL_ERROR_NULL_INPUT);
527 cpl_ensure_code(cpl_matrix_get_nrow(U) >= cpl_matrix_get_ncol(U),
528 CPL_ERROR_ILLEGAL_INPUT);
529 cpl_ensure_code(cpl_matrix_get_nrow(V) == cpl_matrix_get_ncol(V),
530 CPL_ERROR_ILLEGAL_INPUT);
531 cpl_ensure_code(cpl_matrix_get_nrow(V) == cpl_matrix_get_ncol(U),
532 CPL_ERROR_INCOMPATIBLE_INPUT);
533 cpl_ensure_code(cpl_vector_get_size(S) == cpl_matrix_get_ncol(U),
534 CPL_ERROR_INCOMPATIBLE_INPUT);
536 const cpl_size m = cpl_matrix_get_nrow(U);
537 const cpl_size n = cpl_matrix_get_ncol(U);
543 cpl_size sweepmax = CPL_MAX(5 * n, 12);
545 double tolerance = 10. * m * DBL_EPSILON;
548 cpl_matrix_fill(V, 0.);
549 cpl_matrix_fill_diagonal(V, 1., 0);
556 cpl_vector *cv1 = cpl_vector_new(m);
557 double *
const _cv1 = cpl_vector_get_data(cv1);
559 for (j = 0; j < n; ++j) {
560 for (i = 0; i < m; ++i) {
561 _cv1[i] = cpl_matrix_get(U, i, j);
564 double s = sqrt(cpl_vector_product(cv1, cv1));
565 cpl_vector_set(S, j, DBL_EPSILON * s);
569 cpl_vector *cv2 = cpl_vector_new(m);
570 double *
const _cv2 = cpl_vector_get_data(cv2);
572 while ((count > 0) && (sweep <= sweepmax)) {
575 count = n * (n - 1) / 2;
577 for (j = 0; j < n - 1; ++j) {
579 for (k = j + 1; k < n; ++k) {
580 for (i = 0; i < m; ++i) {
581 _cv1[i] = cpl_matrix_get(U, i, j);
583 for (i = 0; i < m; ++i) {
584 _cv2[i] = cpl_matrix_get(U, i, k);
589 double a = sqrt(cpl_vector_product(cv1, cv1));
590 double b = sqrt(cpl_vector_product(cv2, cv2));
591 double c = 2. * cpl_vector_product(cv1, cv2);
594 double abserr_a = cpl_vector_get(S, j);
595 double abserr_b = cpl_vector_get(S, k);
596 cpl_boolean orthogonal = (fabs(c) <= tolerance * a * b);
597 cpl_boolean sorted = (a >= b);
598 cpl_boolean noisy_a = (a < abserr_a);
599 cpl_boolean noisy_b = (b < abserr_b);
601 if (sorted && (orthogonal || noisy_a || noisy_b)) {
608 double q = a * a - b * b;
609 double v = sqrt(c * c + q * q);
613 if (v == 0. || !sorted) {
617 cs = sqrt((v + q) / (2. * v));
618 sn = c / (2. * v * cs);
622 for (i = 0; i < m; ++i) {
623 const double u_ij = cpl_matrix_get(U, i, j);
624 const double u_ik = cpl_matrix_get(U, i, k);
625 cpl_matrix_set(U, i, j, u_ij * cs + u_ik * sn);
626 cpl_matrix_set(U, i, k, -u_ij * sn + u_ik * cs);
629 for (i = 0; i < n; ++i) {
630 const double v_ij = cpl_matrix_get(V, i, j);
631 const double v_ik = cpl_matrix_get(V, i, k);
632 cpl_matrix_set(V, i, j, v_ij * cs + v_ik * sn);
633 cpl_matrix_set(V, i, k, -v_ij * sn + v_ik * cs);
636 cpl_vector_set(S, j, fabs(cs) * abserr_a +
637 fabs(sn) * abserr_b);
638 cpl_vector_set(S, k, fabs(sn) * abserr_a +
639 fabs(cs) * abserr_b);
645 cpl_vector_delete(cv2);
648 double prev_norm = -1.;
649 for (j = 0; j < n; ++j) {
650 for (i = 0; i < m; ++i) {
651 _cv1[i] = cpl_matrix_get(U, i, j);
654 double norm = sqrt(cpl_vector_product(cv1, cv1));
659 if ((norm == 0.) || (prev_norm == 0.) ||
660 ((j > 0) && (norm <= tolerance * prev_norm))) {
661 cpl_vector_set(S, j, 0.);
662 cpl_matrix_fill_column(U, 0., j);
665 cpl_vector_set(S, j, norm);
668 for (i = 0; i < m; ++i) {
669 cpl_matrix_set(U, i, j, _cv1[i] / norm);
674 cpl_vector_delete(cv1);
677 return CPL_ERROR_CONTINUE;
679 return CPL_ERROR_NONE;
693 muse_cplmatrix_solve_svd(
const cpl_matrix *coeff,
const cpl_matrix *rhs)
695 cpl_ensure((coeff != NULL) && (rhs != NULL), CPL_ERROR_NULL_INPUT, NULL);
696 cpl_ensure(cpl_matrix_get_nrow(coeff) == cpl_matrix_get_nrow(rhs),
697 CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
698 cpl_ensure(cpl_matrix_get_ncol(rhs) == 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
701 cpl_matrix *U = cpl_matrix_duplicate(coeff);
702 cpl_size nc = cpl_matrix_get_ncol(U);
704 cpl_matrix *V = cpl_matrix_new(nc, nc);
705 cpl_vector *S = cpl_vector_new(nc);
707 cpl_error_code status = muse_cplmatrix_decomp_sv_jacobi(U, S, V);
708 if (status != CPL_ERROR_NONE) {
709 cpl_vector_delete(S);
710 cpl_matrix_delete(V);
711 cpl_matrix_delete(U);
712 cpl_error_set_where(cpl_func);
717 cpl_matrix *Ut = cpl_matrix_transpose_create(U);
718 cpl_matrix *w = cpl_matrix_product_create(Ut, rhs);
719 cpl_matrix_delete(Ut);
721 const double *_s = cpl_vector_get_data_const(S);
722 double *_w = cpl_matrix_get_data(w);
725 for (ic = 0; ic < nc; ++ic) {
733 cpl_matrix *solution = cpl_matrix_product_create(V, w);
734 cpl_matrix_delete(w);
735 cpl_vector_delete(S);
736 cpl_matrix_delete(V);
737 cpl_matrix_delete(U);
756 muse_cplmatrix_solve_least_square(
const cpl_matrix *aMatrix1,
757 const cpl_matrix *aMatrix2)
759 cpl_ensure(aMatrix1 && aMatrix2, CPL_ERROR_NULL_INPUT, NULL);
761 cpl_size nc = cpl_matrix_get_ncol(aMatrix1);
762 cpl_size nr = cpl_matrix_get_nrow(aMatrix1);
763 cpl_ensure(cpl_matrix_get_nrow(aMatrix2) == nr,
764 CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
766 cpl_matrix *result = NULL;
767 cpl_errorstate status = cpl_errorstate_get();
769 #ifndef USE_SVD_SOLVER 771 result = cpl_matrix_solve(aMatrix1, aMatrix2);
775 result = cpl_matrix_solve_normal(aMatrix1, aMatrix2);
777 cpl_matrix *mt = cpl_matrix_transpose_create(aMatrix1);
778 cpl_matrix *mgram = cpl_matrix_product_create(aMatrix1, mt);
779 cpl_matrix *mw = cpl_matrix_solve(mgram, aMatrix2);
781 result = cpl_matrix_product_create(mt, mw);
782 cpl_matrix_delete(mw);
783 cpl_matrix_delete(mgram);
784 cpl_matrix_delete(mt);
789 result = cpl_matrix_solve_svd(aMatrix1, aMatrix2);
791 cpl_matrix *mt = cpl_matrix_transpose_create(aMatrix1);
792 cpl_matrix *mgram = cpl_matrix_product_create(aMatrix1, mt);
793 cpl_matrix *mw = cpl_matrix_solve(mgram, aMatrix2);
795 result = cpl_matrix_product_create(mt, mw);
796 cpl_matrix_delete(mw);
797 cpl_matrix_delete(mgram);
798 cpl_matrix_delete(mt);
802 if (status != cpl_errorstate_get()) {
803 cpl_matrix_delete(result);
820 static cpl_error_code
821 muse_cplmatrix_cosine(cpl_matrix *aMatrix)
823 cpl_ensure(aMatrix, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
825 cpl_size size = cpl_matrix_get_nrow(aMatrix) * cpl_matrix_get_ncol(aMatrix);
826 double *mdata = cpl_matrix_get_data(aMatrix);
829 for (i = 0; i < size; ++i) {
830 mdata[i] = cos(mdata[i]);
832 return CPL_ERROR_NONE;
853 muse_align_compute_distances(
const cpl_matrix *aDeltaRA,
const cpl_matrix *aDeltaDEC,
854 const cpl_matrix *aMeanDEC,
const double *aOffsetRA,
855 const double *aOffsetDEC)
858 cpl_errorstate status = cpl_errorstate_get();
865 cpl_matrix *tmp = cpl_matrix_duplicate(aDeltaDEC);
866 cpl_matrix_subtract_scalar(tmp, *aOffsetDEC);
869 cpl_matrix_delete(tmp);
874 cpl_matrix *dec = cpl_matrix_duplicate(aMeanDEC);
875 muse_cplmatrix_cosine(dec);
879 dRA = cpl_matrix_duplicate(aDeltaRA);
880 cpl_matrix_subtract_scalar(dRA, *aOffsetRA);
881 cpl_matrix_multiply(dRA, dec);
885 cpl_matrix_delete(dec);
888 cpl_matrix_delete(dRA);
890 cpl_matrix_add(distance, dDEC);
891 cpl_matrix_power(distance, 0.5);
892 cpl_matrix_delete(dDEC);
894 if (status != cpl_errorstate_get()) {
895 cpl_matrix_delete(distance);
924 muse_align_estimate_offsets(
double *aOffsetRA,
double *aOffsetDEC,
double *aWeight,
925 const cpl_matrix *aDeltaRA,
const cpl_matrix *aDeltaDEC,
926 const cpl_matrix *aMeanDEC,
double aRadius,
int aNbins)
928 cpl_ensure((aOffsetRA && aOffsetDEC) && aWeight, CPL_ERROR_NULL_INPUT, 0);
929 cpl_ensure((aDeltaRA && aDeltaDEC) && aMeanDEC, CPL_ERROR_NULL_INPUT, 0);
930 cpl_ensure(cpl_matrix_get_nrow(aDeltaRA) == cpl_matrix_get_nrow(aDeltaDEC),
931 CPL_ERROR_INCOMPATIBLE_INPUT, 0);
932 cpl_ensure(cpl_matrix_get_nrow(aDeltaRA) == cpl_matrix_get_nrow(aMeanDEC),
933 CPL_ERROR_INCOMPATIBLE_INPUT, 0);
934 cpl_ensure(aRadius > 0., CPL_ERROR_ILLEGAL_INPUT, 0);
935 cpl_ensure(aNbins > 0, CPL_ERROR_ILLEGAL_INPUT, 0);
937 cpl_matrix *distance = muse_align_compute_distances(aDeltaRA, aDeltaDEC,
938 aMeanDEC, NULL, NULL);
939 if (cpl_matrix_get_min(distance) >= aRadius) {
940 cpl_matrix_delete(distance);
948 _muse_condition_less_than);
949 cpl_matrix_delete(distance);
952 cpl_array_delete(selection);
955 double bstep = 2. * aRadius / (double)aNbins;
968 cpl_matrix *bins = cpl_matrix_new(1, aNbins);
970 for (ibin = 0; ibin < aNbins; ++ibin) {
971 cpl_matrix_set(bins, 0, ibin, -aRadius + ibin * bstep);
973 cpl_matrix *histogram = cpl_matrix_new(aNbins, aNbins);
974 double *hdata = cpl_matrix_get_data(histogram);
975 cpl_size npos = cpl_matrix_get_ncol(dRA);
978 for (ipos = 0; ipos < npos; ++ipos) {
979 cpl_size ix = (cpl_size)((cpl_matrix_get(dRA, 0, ipos) + aRadius) / bstep);
980 cpl_size iy = (cpl_size)((cpl_matrix_get(dDEC, 0, ipos) + aRadius) / bstep);
981 if (((ix >= 0) && (ix < aNbins)) && ((iy >= 0) && (iy < aNbins))) {
982 hdata[iy * aNbins + ix] += 1.;
986 cpl_matrix_delete(dRA);
987 cpl_matrix_delete(dDEC);
990 cpl_matrix_delete(histogram);
991 cpl_matrix_delete(bins);
997 cpl_matrix_get_maxpos(histogram, &row, &column);
999 *aOffsetRA = cpl_matrix_get(bins, 0, column) + 0.5 * bstep;
1000 *aOffsetDEC = cpl_matrix_get(bins, 0, row) + 0.5 * bstep;
1001 *aWeight = cpl_matrix_get_max(histogram);
1002 cpl_matrix_delete(histogram);
1003 cpl_matrix_delete(bins);
1030 muse_align_update_offsets(
double *aOffsetRA,
double *aOffsetDEC,
double *aWeight,
1031 const cpl_matrix *aDeltaRA,
const cpl_matrix *aDeltaDEC,
1032 const cpl_matrix *aMeanDEC,
double aRadius)
1034 cpl_ensure((aOffsetRA && aOffsetDEC) && aWeight, CPL_ERROR_NULL_INPUT,
1036 cpl_ensure((aDeltaRA && aDeltaDEC) && aMeanDEC, CPL_ERROR_NULL_INPUT,
1038 cpl_ensure(cpl_matrix_get_nrow(aDeltaRA) == cpl_matrix_get_nrow(aDeltaDEC),
1039 CPL_ERROR_INCOMPATIBLE_INPUT, CPL_FALSE);
1040 cpl_ensure(cpl_matrix_get_nrow(aDeltaRA) == cpl_matrix_get_nrow(aMeanDEC),
1041 CPL_ERROR_INCOMPATIBLE_INPUT, CPL_FALSE);
1042 cpl_ensure(aRadius > 0., CPL_ERROR_ILLEGAL_INPUT, CPL_FALSE);
1044 cpl_matrix *distance = muse_align_compute_distances(aDeltaRA, aDeltaDEC,
1045 aMeanDEC, aOffsetRA, aOffsetDEC);
1046 if (cpl_matrix_get_min(distance) >= aRadius) {
1047 cpl_matrix_delete(distance);
1055 _muse_condition_less_than);
1056 cpl_matrix_delete(distance);
1059 cpl_array_delete(selection);
1061 cpl_size nselected = cpl_matrix_get_ncol(dRA);
1062 double variance = 2. * aRadius * aRadius;
1063 if (nselected > 1) {
1064 double sdevRA = cpl_matrix_get_stdev(dRA);
1065 double sdevDEC = cpl_matrix_get_stdev(dDEC);
1066 double _variance = sdevRA * sdevRA + sdevDEC * sdevDEC;
1067 if (_variance > nselected / DBL_MAX) {
1068 variance = _variance;
1072 *aOffsetRA = cpl_matrix_get_median(dRA);
1073 *aOffsetDEC = cpl_matrix_get_median(dDEC);
1074 *aWeight = nselected / variance;
1076 cpl_matrix_delete(dRA);
1077 cpl_matrix_delete(dDEC);
1100 cpl_table **aCataloglist,
1101 const cpl_array *aSearchradius,
1102 int aNbins, cpl_boolean aWeight,
1103 cpl_propertylist *aHeader)
1105 cpl_ensure(aImagelist && aCataloglist, CPL_ERROR_NULL_INPUT, NULL);
1106 cpl_ensure(aSearchradius, CPL_ERROR_NULL_INPUT, NULL);
1107 cpl_ensure(aNbins > 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
1110 cpl_ensure(nfields >= 2, CPL_ERROR_DATA_NOT_FOUND, NULL);
1114 cpl_size npairs = nfields * (nfields - 1) / 2;
1116 muse_indexpair *combinations = cpl_calloc(npairs,
sizeof *combinations);
1119 for (ifield = 0; ifield < nfields; ++ifield) {
1121 for (jfield = ifield + 1; jfield < nfields; ++jfield) {
1122 combinations[ipair].first = ifield;
1123 combinations[ipair].second = jfield;
1129 for (ifield = 0; ifield < nfields; ++ifield) {
1130 char *kw = cpl_sprintf(QC_ALIGN_NDETi, (
int)ifield + 1);
1131 cpl_propertylist_append_int(aHeader, kw,
1132 cpl_table_get_nrow(aCataloglist[ifield]));
1140 cpl_array **aoverlaps = cpl_calloc(nfields,
sizeof(cpl_array *));
1141 for (ifield = 0; ifield < nfields; ++ifield) {
1142 aoverlaps[ifield] = cpl_array_new(nfields, CPL_TYPE_INT);
1146 cpl_matrix *ra_offsets = NULL;
1147 cpl_matrix *dec_offsets = NULL;
1148 cpl_size *has_neighbors = cpl_malloc(nfields *
sizeof *has_neighbors);
1150 for (isearch = 0; isearch < cpl_array_get_size(aSearchradius); ++isearch) {
1151 double radius = cpl_array_get_double(aSearchradius, isearch, NULL);
1159 for (ifield = 0; ifield < nfields; ++ifield) {
1160 cpl_array_fill_window_int(aoverlaps[ifield], 0, nfields, 0);
1164 cpl_array_set_invalid(aoverlaps[ifield], ifield);
1166 cpl_array_fill_window_invalid(aoverlaps[ifield], 0, nfields);
1172 cpl_matrix *weights = cpl_matrix_new(npairs + 1, nfields);
1173 cpl_matrix *offsets_ra = cpl_matrix_new(npairs + 1, 1);
1174 cpl_matrix *offsets_dec = cpl_matrix_new(npairs + 1, 1);
1176 cpl_matrix_fill_row(weights, 1., 0);
1177 cpl_size kpairs = 0;
1179 memset(has_neighbors, 0, nfields *
sizeof *has_neighbors);
1181 for (ipair = 0; ipair < npairs; ++ipair) {
1183 const cpl_table *catalog1 = aCataloglist[combinations[ipair].first];
1184 const cpl_table *catalog2 = aCataloglist[combinations[ipair].second];
1186 cpl_size nstars1 = cpl_table_get_nrow(catalog1);
1187 cpl_size nstars2 = cpl_table_get_nrow(catalog2);
1189 const double *ra1 = cpl_table_get_data_double_const(catalog1,
"RA");
1190 const double *ra2 = cpl_table_get_data_double_const(catalog2,
"RA");
1191 const double *dec1 = cpl_table_get_data_double_const(catalog1,
"DEC");
1192 const double *dec2 = cpl_table_get_data_double_const(catalog2,
"DEC");
1194 cpl_matrix *delta_ra = cpl_matrix_new(nstars1, nstars2);
1195 cpl_matrix *delta_dec = cpl_matrix_new(nstars1, nstars2);
1196 cpl_matrix *dec_mean = cpl_matrix_new(nstars1, nstars2);
1199 for (irow = 0; irow < nstars1; ++irow) {
1201 for (icol = 0; icol < nstars2; ++icol) {
1202 cpl_matrix_set(delta_ra, irow, icol, (ra1[irow] - ra2[icol]) * deg2as);
1203 cpl_matrix_set(delta_dec, irow, icol, (dec1[irow] - dec2[icol]) * deg2as);
1204 cpl_matrix_set(dec_mean, irow, icol,
1205 0.5 *(dec1[irow] + dec2[icol]) * CPL_MATH_RAD_DEG);
1209 double offset_ra = 0.;
1210 double offset_dec = 0.;
1212 cpl_size noverlap = 0;
1214 noverlap = muse_align_estimate_offsets(&offset_ra, &offset_dec, &weight,
1215 delta_ra, delta_dec, dec_mean,
1218 offset_ra = cpl_matrix_get(ra_offsets, combinations[ipair].first, 0) -
1219 cpl_matrix_get(ra_offsets, combinations[ipair].second, 0);
1220 offset_dec = cpl_matrix_get(dec_offsets, combinations[ipair].first, 0) -
1221 cpl_matrix_get(dec_offsets, combinations[ipair].second, 0);
1223 noverlap = muse_align_update_offsets(&offset_ra, &offset_dec, &weight,
1224 delta_ra, delta_dec, dec_mean,
1227 cpl_matrix_delete(dec_mean);
1228 cpl_matrix_delete(delta_ra);
1229 cpl_matrix_delete(delta_dec);
1237 cpl_matrix_set(offsets_ra, kpairs, 0, weight * offset_ra);
1238 cpl_matrix_set(offsets_dec, kpairs, 0, weight * offset_dec);
1239 cpl_matrix_set(weights, kpairs, combinations[ipair].first, weight);
1240 cpl_matrix_set(weights, kpairs, combinations[ipair].second, -weight);
1241 has_neighbors[combinations[ipair].first] += 1;
1242 has_neighbors[combinations[ipair].second] += 1;
1245 cpl_array_set_int(aoverlaps[combinations[ipair].first],
1246 combinations[ipair].second, noverlap);
1247 cpl_array_set_int(aoverlaps[combinations[ipair].second],
1248 combinations[ipair].first, noverlap);
1253 cpl_matrix_delete(ra_offsets);
1254 cpl_matrix_delete(dec_offsets);
1259 cpl_matrix_delete(offsets_dec);
1260 cpl_matrix_delete(offsets_ra);
1261 cpl_matrix_delete(weights);
1263 cpl_msg_warning(__func__,
"No overlapping fields of view were found " 1264 "for search radius %.4f!", radius);
1266 ra_offsets = cpl_matrix_new(nfields, 1);
1267 dec_offsets = cpl_matrix_new(nfields, 1);
1271 cpl_matrix *_weights = cpl_matrix_wrap(kpairs + 1, nfields,
1272 cpl_matrix_get_data(weights));
1273 cpl_matrix *_offsets_ra = cpl_matrix_wrap(kpairs + 1, 1,
1274 cpl_matrix_get_data(offsets_ra));
1275 cpl_matrix *_offsets_dec = cpl_matrix_wrap(kpairs + 1, 1,
1276 cpl_matrix_get_data(offsets_dec));
1278 ra_offsets = muse_cplmatrix_solve_least_square(_weights, _offsets_ra);
1279 dec_offsets = muse_cplmatrix_solve_least_square(_weights, _offsets_dec);
1281 cpl_matrix_unwrap(_offsets_dec);
1282 cpl_matrix_unwrap(_offsets_ra);
1283 cpl_matrix_unwrap(_weights);
1285 cpl_matrix_delete(offsets_dec);
1286 cpl_matrix_delete(offsets_ra);
1287 cpl_matrix_delete(weights);
1290 cpl_free(combinations);
1299 cpl_size minmatch = LLONG_MAX,
1301 for (ifield = 0; ifield < nfields; ++ifield) {
1305 cpl_table_set_string(offsets, MUSE_OFFSETS_DATEOBS, ifield, timestamp);
1306 cpl_table_set_double(offsets, MUSE_OFFSETS_MJDOBS, ifield, mjd);
1308 if (has_neighbors[ifield]) {
1309 double ra_offset = cpl_matrix_get(ra_offsets, ifield, 0);
1310 double dec_offset = cpl_matrix_get(dec_offsets, ifield, 0);
1311 cpl_table_set_double(offsets, MUSE_OFFSETS_DRA, ifield, ra_offset / deg2as);
1312 cpl_table_set_double(offsets, MUSE_OFFSETS_DDEC, ifield, dec_offset / deg2as);
1314 int nmedian = cpl_array_get_median(aoverlaps[ifield]);
1315 char *kw = cpl_sprintf(QC_ALIGN_NMATCHi, (
int)ifield + 1);
1316 cpl_propertylist_append_int(aHeader, kw, nmedian);
1318 if (nmedian < minmatch) {
1322 cpl_table_set_double(offsets, MUSE_OFFSETS_DRA, ifield, 0.);
1323 cpl_table_set_double(offsets, MUSE_OFFSETS_DDEC, ifield, 0.);
1324 char *kw = cpl_sprintf(QC_ALIGN_NMATCHi, (
int)ifield + 1);
1325 cpl_propertylist_append_int(aHeader, kw, 0);
1330 cpl_propertylist_append_int(aHeader, QC_ALIGN_NMATCH_MIN, minmatch);
1331 cpl_propertylist_append_int(aHeader, QC_ALIGN_NOMATCH, nomatch);
1332 muse_vfree((
void **)aoverlaps, nfields, (muse_free_func)cpl_array_delete);
1333 cpl_free(has_neighbors);
1334 cpl_matrix_delete(dec_offsets);
1335 cpl_matrix_delete(ra_offsets);
1340 double ra_shift = cpl_table_get(offsets, MUSE_OFFSETS_DRA, 0, NULL);
1341 double dec_shift = cpl_table_get(offsets, MUSE_OFFSETS_DDEC, 0, NULL);
1342 cpl_table_subtract_scalar(offsets, MUSE_OFFSETS_DRA, ra_shift);
1343 cpl_table_subtract_scalar(offsets, MUSE_OFFSETS_DDEC, dec_shift);
1344 cpl_table_set_column_invalid(offsets, MUSE_OFFSETS_FSCALE, 0, nfields);
1347 double qcmin = cpl_table_get_column_min(offsets, MUSE_OFFSETS_DRA) * 3600.;
1348 double qcmax = cpl_table_get_column_max(offsets, MUSE_OFFSETS_DRA) * 3600.;
1349 double qcmean = cpl_table_get_column_mean(offsets, MUSE_OFFSETS_DRA) * 3600.;
1350 double qcsdev = cpl_table_get_column_stdev(offsets, MUSE_OFFSETS_DRA) * 3600.;
1352 cpl_propertylist_append_double(aHeader, QC_ALIGN_DRA_MIN, qcmin);
1353 cpl_propertylist_append_double(aHeader, QC_ALIGN_DRA_MAX, qcmax);
1354 cpl_propertylist_append_double(aHeader, QC_ALIGN_DRA_MEAN, qcmean);
1355 cpl_propertylist_append_double(aHeader, QC_ALIGN_DRA_STDEV, qcsdev);
1357 qcmin = cpl_table_get_column_min(offsets, MUSE_OFFSETS_DDEC) * 3600.;
1358 qcmax = cpl_table_get_column_max(offsets, MUSE_OFFSETS_DDEC) * 3600.;
1359 qcmean = cpl_table_get_column_mean(offsets, MUSE_OFFSETS_DDEC) * 3600.;
1360 qcsdev = cpl_table_get_column_stdev(offsets, MUSE_OFFSETS_DDEC) * 3600.;
1362 cpl_propertylist_append_double(aHeader, QC_ALIGN_DDEC_MIN, qcmin);
1363 cpl_propertylist_append_double(aHeader, QC_ALIGN_DDEC_MAX, qcmax);
1364 cpl_propertylist_append_double(aHeader, QC_ALIGN_DDEC_MEAN, qcmean);
1365 cpl_propertylist_append_double(aHeader, QC_ALIGN_DDEC_STDEV, qcsdev);
1385 cpl_size nexposures = cpl_frameset_count_tags(aProcessing->
inframes,
1386 MUSE_TAG_IMAGE_FOV);
1387 if (nexposures < 2) {
1388 cpl_msg_error(__func__,
"This recipe requires at least 2 input FOV " 1394 cpl_msg_debug(__func__,
"Validating source detection parameters");
1395 if (!muse_align_check_detection_params(aParams)) {
1396 cpl_error_set_where(__func__);
1401 cpl_msg_info(__func__,
"Loading FOV images");
1403 muse_imagelist *fovimages = muse_processing_fov_load_all(aProcessing);
1405 cpl_msg_error(__func__,
"Could not load FOV input images!");
1411 for (iexposure = 0; iexposure < nexposures; ++iexposure) {
1413 cpl_errorstate es = cpl_errorstate_get();
1416 if (!cpl_errorstate_is_equal(es) &&
1417 cpl_error_get_code() == CPL_ERROR_DATA_NOT_FOUND) {
1418 cpl_msg_error(__func__,
"Missing \"MJD-OBS\" in FOV image %" 1419 CPL_SIZE_FORMAT, iexposure + 1);
1424 for (jexposure = iexposure + 1; jexposure < nexposures; ++jexposure) {
1427 if (_mjdobs == mjdobs) {
1428 cpl_msg_warning(__func__,
"Found FOV images with non-unique " 1429 "\"MJD-OBS\" value (image %" CPL_SIZE_FORMAT
1430 " and %" CPL_SIZE_FORMAT
")!", iexposure + 1,
1437 cpl_msg_info(__func__,
"Computing pointing corrections for %" CPL_SIZE_FORMAT
1438 " FOV images", nexposures);
1442 for (iexposure = 0; iexposure < nexposures; ++iexposure) {
1444 muse_utils_replace_nan(image, 0.);
1454 cpl_table **srclists = cpl_calloc(nexposures,
sizeof *srclists);
1456 for (iexposure = 0; iexposure < nexposures; ++iexposure) {
1457 cpl_msg_debug(__func__,
"Detecting sources on image %" CPL_SIZE_FORMAT
1458 " of %" CPL_SIZE_FORMAT, iexposure + 1, nexposures);
1462 cpl_boolean iterate = CPL_TRUE;
1464 cpl_size ndetections = 0;
1465 cpl_table *detections;
1466 cpl_image *image = cpl_image_cast(fov->
data, CPL_TYPE_DOUBLE);
1474 cpl_size nx = cpl_image_get_size_x(image);
1475 cpl_size ny = cpl_image_get_size_y(image);
1477 if ((nx < 6) || (ny < 6)) {
1478 cpl_image_delete(image);
1479 muse_vfree((
void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
1481 cpl_msg_error(__func__,
"Image %" CPL_SIZE_FORMAT
" of %" 1482 CPL_SIZE_FORMAT
"is too small!", iexposure + 1, nexposures);
1486 if ((nx % 2 != 0) || (ny % 2 != 0)) {
1487 cpl_msg_debug(__func__,
"Trimming image %" CPL_SIZE_FORMAT
" to an " 1488 "even number of pixels along the axes.", iexposure + 1);
1495 cpl_image *_image = cpl_image_extract(image, 1, 1, nx, ny);
1497 cpl_image_delete(image);
1498 muse_vfree((
void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
1500 cpl_msg_error(__func__,
"Trimming image %" CPL_SIZE_FORMAT
" failed",
1504 cpl_image_delete(image);
1509 cpl_errorstate status = cpl_errorstate_get();
1512 roundness, sharpness);
1514 if (cpl_error_get_code() == CPL_ERROR_DATA_NOT_FOUND) {
1515 cpl_errorstate_set(status);
1518 cpl_image_delete(image);
1519 muse_vfree((
void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
1521 cpl_msg_error(__func__,
"Source detection failed on image %" 1522 CPL_SIZE_FORMAT
" of %" CPL_SIZE_FORMAT, iexposure + 1,
1527 ndetections = cpl_table_get_nrow(detections);
1530 if ((ndetections > 0) &&
1531 ((ndetections >= nobjects[0]) && (ndetections <= nobjects[1]))) {
1532 iterate = CPL_FALSE;
1534 cpl_table_delete(detections);
1536 if (++iteration < aParams->iterations) {
1537 if (ndetections > nobjects[1]) {
1538 threshold += aParams->
step;
1540 threshold -= aParams->
step;
1541 if (threshold <= 0.) {
1542 iterate = CPL_FALSE;
1546 iterate = CPL_FALSE;
1550 cpl_image_delete(image);
1553 muse_vfree((
void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
1555 if (ndetections == 0) {
1556 cpl_msg_error(__func__,
"No point-sources were detected in image %" 1557 CPL_SIZE_FORMAT, iexposure + 1);
1558 }
else if (ndetections < nobjects[0]) {
1559 cpl_msg_error(__func__,
"Too few (%" CPL_SIZE_FORMAT
") point-sources " 1560 "were detected in image %" CPL_SIZE_FORMAT, ndetections,
1562 }
else if (ndetections > nobjects[1]) {
1563 cpl_msg_error(__func__,
"Too many (%" CPL_SIZE_FORMAT
") point-sources " 1564 "were detected in image %" CPL_SIZE_FORMAT, ndetections,
1567 cpl_msg_error(__func__,
"Invalid number (%" CPL_SIZE_FORMAT
") of " 1568 "point-sources were detected in image %" CPL_SIZE_FORMAT,
1569 ndetections, iexposure + 1);
1574 cpl_msg_info(__func__,
"%" CPL_SIZE_FORMAT
" point-sources detected in " 1575 "image %" CPL_SIZE_FORMAT, ndetections, iexposure + 1);
1579 if (muse_align_celestial_from_pixel(detections, fov->
header)) {
1580 muse_vfree((
void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
1582 cpl_msg_error(__func__,
"Computing WCS coordinates from pixel positions " 1583 "failed for image %" CPL_SIZE_FORMAT
" of %" 1584 CPL_SIZE_FORMAT, iexposure + 1, nexposures);
1587 srclists[iexposure] = detections;
1592 cpl_array *radii = muse_align_parse_valuelist(aParams->
rsearch);
1594 muse_vfree((
void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
1596 cpl_msg_error(__func__,
"Cannot compute field offsets! No valid search " 1597 "radius was given!");
1601 cpl_msg_info(__func__,
"Calculating field offsets...");
1602 cpl_propertylist *header = cpl_propertylist_new();
1603 cpl_table *offsets = muse_align_compute_field_offsets(fovimages, srclists, radii,
1606 cpl_array_delete(radii);
1607 muse_vfree((
void **)srclists, nexposures, (muse_free_func)cpl_table_delete);
1610 cpl_msg_error(__func__,
"Could not compute FOV offsets for %" 1611 CPL_SIZE_FORMAT
" images", nexposures);
1612 cpl_propertylist_delete(header);
1618 header, MUSE_TAG_OFFSET_LIST,
1620 cpl_propertylist_delete(header);
1621 cpl_table_delete(offsets);
1623 if (rc != CPL_ERROR_NONE) {
int srcmin
Minimum number of sources which must be found.
Structure definition for a collection of muse_images.
void muse_image_delete(muse_image *aImage)
Deallocate memory associated to a muse_image object.
const char * muse_pfits_get_extname(const cpl_propertylist *aHeaders)
find out the extension name
double fwhm
FWHM in pixels of the convolution filter.
A structure containing a spatial two-axis WCS.
cpl_size muse_pfits_get_naxis(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the size of a given axis
cpl_image * data
the data extension
double roundmax
Upper limit of the allowed point-source roundness.
cpl_image * stat
the statistics extension
const char * muse_pfits_get_dateobs(const cpl_propertylist *aHeaders)
find out the date of observations
void muse_imagelist_delete(muse_imagelist *aList)
Free the memory of the MUSE image list.
muse_image * muse_fov_load(const char *aFilename)
Load a FOV image into a MUSE image.
const char * muse_pfits_get_bunit(const cpl_propertylist *aHeaders)
find out the unit string
Structure definition of MUSE three extension FITS file.
const char * muse_pfits_get_cunit(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS axis unit string
cpl_propertylist * header
the FITS header
cpl_boolean muse_processing_check_intags(muse_processing *aProcessing, const char *aTag, int aNChars)
Check that a tag is part of the input tags of a processing structure.
double step
Increment/decrement of the threshold value in subsequent iterations.
unsigned int muse_imagelist_get_size(muse_imagelist *aList)
Return the number of stored images.
cpl_matrix * muse_cplmatrix_extract_selected(const cpl_matrix *aMatrix, const cpl_array *aIndices)
Extract the elements given by the index array.
cpl_image * dq
the data quality extension
cpl_table * muse_cpltable_new(const muse_cpltable_def *aDef, cpl_size aLength)
Create an empty table according to the specified definition.
muse_wcs * muse_wcs_new(cpl_propertylist *aHeader)
Create a new WCS structure from a given FITS header.
cpl_array * muse_cplarray_new_from_delimited_string(const char *aString, const char *aDelim)
Convert a delimited string into an array of strings.
double sharpmax
Upper limit of the allowed point-source sharpness.
Structure to hold the parameters of the muse_exp_align recipe.
int srcmax
Maximum number of sources which may be found.
#define MUSE_WCS_KEYS
regular expression for WCS properties
muse_image * muse_imagelist_get(muse_imagelist *aList, unsigned int aIdx)
Get the muse_image of given list index.
double threshold
Initial intensity threshold for detecting point sources in sigma above background RMS...
static void muse_wcs_celestial_from_pixel_fast(muse_wcs *aWCS, double aX, double aY, double *aRA, double *aDEC)
Convert from pixel coordinates to celestial spherical coordinates.
int nbins
Number of bins to use for 2D histogram on the first iteration of the offset computation.
void muse_processing_append_used(muse_processing *aProcessing, cpl_frame *aFrame, cpl_frame_group aGroup, int aDuplicate)
Add a frame to the set of used frames.
const char * rsearch
Search radius (in arcsec) for each iteration of the offset computation.
double muse_pfits_get_mjdobs(const cpl_propertylist *aHeaders)
find out the Julian Date of the observation
double roundmin
Lower limit of the allowed point-source roundness.
muse_imagelist * muse_imagelist_new(void)
Create a new (empty) MUSE image list.
cpl_matrix * muse_cplmatrix_multiply_create(const cpl_matrix *aMatrix1, const cpl_matrix *aMatrix2)
Compute the element-wise product of two matrices.
cpl_array * muse_cplmatrix_where(const cpl_matrix *aMatrix, double aValue, muse_cplmatrix_element_compare_func aCondition)
Select matrix elements according to a condition.
cpl_table * muse_find_stars(const cpl_image *aImage, double aHmin, double aFwhm, const double *aRoundLimits, const double *aSharpLimits)
Find positive brightness perturbations (i.e stars) in an image.
cpl_error_code muse_processing_save_table(muse_processing *aProcessing, int aIFU, void *aTable, cpl_propertylist *aHeader, const char *aTag, muse_table_type aType)
Save a computed table to disk.
muse_image * muse_image_new(void)
Allocate memory for a new muse_image object.
cpl_error_code muse_imagelist_set(muse_imagelist *aList, muse_image *aImage, unsigned int aIdx)
Set the muse_image of given list index.
cpl_error_code muse_image_dq_to_nan(muse_image *aImage)
Convert pixels flagged in the DQ extension to NANs in DATA (and STAT, if present).
double sharpmin
Lower limit of the allowed point-source sharpness.