33 #include "muse_lingain_z.h" 36 typedef void (*muse_free_func)(
void *);
38 struct muse_gain_fit_t
44 typedef struct muse_gain_fit_t muse_gain_fit_t;
46 struct muse_linearity_fit_t
48 cpl_array *coefficients;
53 typedef struct muse_linearity_fit_t muse_linearity_fit_t;
61 muse_vfree(
void **array, cpl_size size, muse_free_func deallocator)
65 for (idx = 0; idx < size; ++idx) {
67 deallocator(array[idx]);
76 muse_linearity_fit_clear(muse_linearity_fit_t *aFit)
79 if (aFit->coefficients) {
80 cpl_array_delete(aFit->coefficients);
101 if (aParams->
ybox <= 0) {
102 cpl_msg_error(__func__,
"Invalid measurement window size: window must " 103 "be larger than 0!");
107 if (aParams->
xgap < 0) {
108 cpl_msg_error(__func__,
"Invalid tracing edge offset: offset is " 114 cpl_msg_error(__func__,
"Invalid offset from detector edge: offset is " 119 if (aParams->
order < 0) {
120 cpl_msg_error(__func__,
"Invalid polynomial fit order: polynomial " 121 "degree is less than 0!");
126 cpl_msg_error(__func__,
"Invalid exposure time offset: offset is less " 131 if (aParams->
sigma < 0.) {
132 cpl_msg_error(__func__,
"Invalid sigma clipping threshold used for " 133 "signal value cleaning: threshold is less than 0.!");
138 cpl_msg_error(__func__,
"Invalid signal level grid: minimum signal is " 144 cpl_msg_error(__func__,
"Invalid signal level grid: signal bin size " 145 "must larger than 0.!");
150 cpl_msg_error(__func__,
"Invalid signal level grid: maximum signal is " 151 "too small for chosen minimum value and bin size!");
156 cpl_msg_error(__func__,
"Invalid minimum signal used for gain fit: " 157 "minimum signal is less than 0.!");
162 cpl_msg_error(__func__,
"Invalid sigma clipping threshold used for " 163 "gain value cleaning: threshold is less than 0.!");
167 if (aParams->
ctsmin < 0.) {
168 cpl_msg_error(__func__,
"Invalid signal level grid used for " 169 "non-linearity analysis: minimum signal is " 174 if (aParams->
ctsbin <= 0.) {
175 cpl_msg_error(__func__,
"Invalid signal level grid used for " 176 "non-linearity analysis: signal bin size must larger " 182 cpl_msg_error(__func__,
"Invalid signal level grid used for " 183 "non-linearity analysis: maximum signal is " 184 "too small for chosen minimum value and bin size!");
189 cpl_msg_error(__func__,
"Invalid signal range used for residual " 190 "non-linearity fit: minimum signal is less than 0.!");
195 cpl_msg_error(__func__,
"Invalid signal range used for residual " 196 "non-linearity fit: maximum signal is too small for " 197 "the cosen minimum value!");
213 muse_lingain_load_images(
muse_processing *aProcessing,
unsigned short aIFU)
221 MUSE_TAG_MASTER_BIAS, aIFU);
225 cpl_propertylist *properties = cpl_propertylist_load(cpl_frame_get_filename(bias), 0);
226 cpl_frame_delete(bias);
228 cpl_propertylist_delete(properties);
250 cpl_ensure(aList, CPL_ERROR_NULL_INPUT, NULL);
252 cpl_errorstate status = cpl_errorstate_get();
254 cpl_table *exposures = cpl_table_new(aList->
size);
256 cpl_table_new_column(exposures,
"INDEX", CPL_TYPE_INT);
257 cpl_table_new_column(exposures,
"TYPE", CPL_TYPE_STRING);
258 cpl_table_new_column(exposures,
"MJDOBS", CPL_TYPE_DOUBLE);
259 cpl_table_new_column(exposures,
"EXPTIME", CPL_TYPE_DOUBLE);
261 unsigned int iexposure;
262 for (iexposure = 0; iexposure < aList->
size; ++iexposure) {
264 cpl_table_set_int(exposures,
"INDEX", iexposure, iexposure);
266 cpl_table_set_string(exposures,
"TYPE", iexposure,
267 MUSE_TAG_LINEARITY_BIAS);
269 cpl_table_set_string(exposures,
"TYPE", iexposure,
270 MUSE_TAG_LINEARITY_FLAT);
272 cpl_table_set_double(exposures,
"MJDOBS", iexposure,
274 cpl_table_set_double(exposures,
"EXPTIME", iexposure,
281 cpl_propertylist *order = cpl_propertylist_new();
282 cpl_propertylist_append_bool(order,
"EXPTIME", CPL_FALSE);
283 cpl_propertylist_append_bool(order,
"MJDOBS", CPL_FALSE);
285 cpl_table_sort(exposures, order);
286 cpl_propertylist_delete(order);
288 if (!cpl_errorstate_is_equal(status)) {
289 cpl_table_delete(exposures);
295 cpl_boolean invalid = CPL_FALSE;
296 cpl_size npairs = (cpl_size)(aList->
size / 2);
298 for (ipair = 0; ipair < npairs; ++ipair) {
299 register cpl_size jpair = 2 * ipair;
300 register cpl_size kpair = jpair + 1;
301 if (strcmp(cpl_table_get_string(exposures,
"TYPE", jpair),
302 cpl_table_get_string(exposures,
"TYPE", kpair)) != 0) {
303 cpl_msg_error(__func__,
"Invalid pair of exposures: Type of exposure " 304 "does not match: Got '%s', expected '%s'!",
305 cpl_table_get_string(exposures,
"TYPE", jpair),
306 cpl_table_get_string(exposures,
"TYPE", kpair));
310 if ((cpl_table_get_double(exposures,
"EXPTIME", jpair, NULL) -
311 cpl_table_get_double(exposures,
"EXPTIME", kpair, NULL)) > 0.1) {
312 cpl_msg_error(__func__,
"Invalid pair of exposures: exposure times " 313 "do not match: Got '%.4f', expected '%.4f'!",
314 cpl_table_get_double(exposures,
"EXPTIME", jpair, NULL),
315 cpl_table_get_double(exposures,
"EXPTIME", kpair, NULL));
322 cpl_table_delete(exposures);
323 cpl_error_set(__func__, CPL_ERROR_ILLEGAL_OUTPUT);
328 cpl_size nexposures = 2 * npairs;
329 if (aList->
size > nexposures) {
330 cpl_table_erase_window(exposures, nexposures, aList->
size - nexposures);
348 muse_lingain_create_windowlist(
const cpl_table *aTrace,
349 const cpl_table *aExposureList,
353 cpl_ensure(aTrace && aExposureList && aList && aParams,
354 CPL_ERROR_NULL_INPUT, NULL);
356 cpl_table *windowlist = cpl_table_new(0);
357 cpl_table_new_column(windowlist,
"Quadrant", CPL_TYPE_INT);
358 cpl_table_new_column(windowlist,
"SliceNo", CPL_TYPE_INT);
359 cpl_table_new_column(windowlist,
"Xmin", CPL_TYPE_INT);
360 cpl_table_new_column(windowlist,
"Ymin", CPL_TYPE_INT);
361 cpl_table_new_column(windowlist,
"Xmax", CPL_TYPE_INT);
362 cpl_table_new_column(windowlist,
"Ymax", CPL_TYPE_INT);
367 int idx = cpl_table_get_int(aExposureList,
"INDEX", 0, NULL);
369 cpl_errorstate status = cpl_errorstate_get();
372 cpl_size firstrow = 0;
374 unsigned char iquadrant;
375 for (iquadrant = 1; iquadrant <= 4; ++iquadrant) {
378 unsigned short islice;
379 for (islice = 1; islice <= kMuseSlicesPerCCD; ++islice) {
380 cpl_errorstate _status = cpl_errorstate_get();
386 cpl_msg_warning(__func__,
"slice %2d of IFU %hhu: tracing polynomials " 387 "missing!", islice, aParams->
nifu);
389 cpl_errorstate_set(_status);
393 double ledge = cpl_polynomial_get_coeff(traces[MUSE_TRACE_LEFT], &torder);
394 double redge = cpl_polynomial_get_coeff(traces[MUSE_TRACE_RIGHT], &torder);
395 if ((ledge < qwindow[0] + aParams->
xborder) ||
396 (redge > qwindow[1] - aParams->
xborder)) {
401 unsigned int nwindows = (qwindow[3] - qwindow[2] + 1) / aParams->
ybox;
404 cpl_table_set_size(windowlist, firstrow + nwindows);
406 unsigned int iwindow;
407 for (iwindow = 0; iwindow < nwindows; ++iwindow) {
408 cpl_size irow = firstrow + iwindow;
410 int ymin = iwindow * aParams->
ybox + qwindow[2];
411 int ymax = ymin + aParams->
ybox - 1;
413 double ymid = ymin + 0.5 * aParams->
ybox;
414 double xleft = cpl_polynomial_eval_1d(traces[MUSE_TRACE_LEFT], ymid, NULL);
415 double xright = cpl_polynomial_eval_1d(traces[MUSE_TRACE_RIGHT], ymid, NULL);
417 int xmin = (int)rint(xleft) + aParams->
xgap;
418 int xmax = (int)rint(xright) - aParams->
xgap;
420 cpl_table_set_int(windowlist,
"Quadrant", irow, iquadrant);
421 cpl_table_set_int(windowlist,
"SliceNo", irow, islice);
422 cpl_table_set_int(windowlist,
"Xmin", irow, xmin);
423 cpl_table_set_int(windowlist,
"Xmax", irow, xmax);
424 cpl_table_set_int(windowlist,
"Ymin", irow, ymin);
425 cpl_table_set_int(windowlist,
"Ymax", irow, ymax);
427 firstrow += nwindows;
436 if (!cpl_errorstate_is_equal(status)) {
437 cpl_table_delete(windowlist);
457 muse_lingain_quadrant_get_windows(cpl_table *aWindowList,
458 unsigned char aQuadrant)
460 cpl_table_select_all(aWindowList);
463 cpl_size nwindows = cpl_table_and_selected_int(aWindowList,
"Quadrant",
464 CPL_EQUAL_TO, aQuadrant);
469 cpl_array *windex = cpl_array_new(nwindows, CPL_TYPE_SIZE);
471 cpl_size iwindow = 0;
473 for (irow = 0; irow < cpl_table_get_nrow(aWindowList); ++irow) {
474 if (cpl_table_is_selected(aWindowList, irow)) {
475 cpl_array_set_cplsize(windex, iwindow, irow);
481 cpl_table_select_all(aWindowList);
500 muse_lingain_quadrant_extract_data(
const cpl_table *aGainTable,
501 const cpl_array *aWindows,
506 const cpl_array *src = cpl_table_get_array(aGainTable, aName, aRow);
512 cpl_size nwindows = cpl_array_get_size(aWindows);
513 cpl_array *data = cpl_array_new(nwindows, CPL_TYPE_DOUBLE);
516 for (iwindow = 0; iwindow < nwindows; ++iwindow) {
517 cpl_size jwindow = cpl_array_get_cplsize(aWindows, iwindow, NULL);
519 double value = cpl_array_get_double(src, jwindow, &invalid);
522 cpl_array_set_double(data, iwindow, value);
524 cpl_array_set_invalid(data, iwindow);
551 muse_lingain_compute_ron(
const cpl_table *aWindowList,
552 cpl_table *aExposureList,
556 cpl_table_select_all(aExposureList);
557 cpl_size nbias = cpl_table_and_selected_string(aExposureList,
"TYPE", CPL_EQUAL_TO,
558 MUSE_TAG_LINEARITY_BIAS);
560 cpl_msg_warning(__func__,
"Found more than 2 (%" CPL_SIZE_FORMAT
") " 561 "input images of type '%s' for IFU %u! Only the first 2 " 562 "images will be used!", nbias, MUSE_TAG_LINEARITY_BIAS,
568 for (iexposure = 0; iexposure < cpl_table_get_nrow(aExposureList); ++iexposure) {
569 if (cpl_table_is_selected(aExposureList, iexposure)) {
570 ibias = cpl_table_get_int(aExposureList,
"INDEX", iexposure, NULL);
574 cpl_table_select_all(aExposureList);
585 cpl_image *dbias = cpl_image_subtract_create(bias1->
data, bias2->
data);
587 cpl_array *ron = cpl_array_new(cpl_table_get_nrow(aWindowList), CPL_TYPE_DOUBLE);
590 for (iwindow = 0; iwindow < cpl_table_get_nrow(aWindowList); ++ iwindow) {
591 int xmin = cpl_table_get_int(aWindowList,
"Xmin", iwindow, NULL);
592 int ymin = cpl_table_get_int(aWindowList,
"Ymin", iwindow, NULL);
593 int xmax = cpl_table_get_int(aWindowList,
"Xmax", iwindow, NULL);
594 int ymax = cpl_table_get_int(aWindowList,
"Ymax", iwindow, NULL);
596 double _ron = cpl_image_get_stdev_window(dbias, xmin, ymin, xmax, ymax);
598 _ron = sqrt(0.5 * _ron * _ron);
599 cpl_array_set_double(ron, iwindow, _ron);
601 cpl_image_delete(dbias);
630 muse_lingain_window_get_gain(
const cpl_table *aWindowList,
631 const cpl_array *aRon,
632 cpl_table *aExposureList,
639 cpl_table_select_all(aExposureList);
640 cpl_size nflats = cpl_table_and_selected_string(aExposureList,
"TYPE", CPL_EQUAL_TO,
641 MUSE_TAG_LINEARITY_FLAT);
643 cpl_msg_warning(__func__,
"Found odd number of input images (%" CPL_SIZE_FORMAT
644 ") input images of type '%s' for IFU %u! Exposure series " 645 "may be incorrect!", nflats, MUSE_TAG_LINEARITY_FLAT,
649 cpl_table *flats = cpl_table_extract_selected(aExposureList);
650 cpl_table_select_all(aExposureList);
653 cpl_errorstate state = cpl_errorstate_get();
655 cpl_size npairs = nflats / 2;
656 cpl_size nwindows = cpl_table_get_nrow(aWindowList);
658 cpl_table *gain = cpl_table_new(npairs);
659 cpl_table_new_column(gain,
"ExpTime", CPL_TYPE_DOUBLE);
660 cpl_table_new_column_array(gain,
"Signal", CPL_TYPE_DOUBLE, nwindows);
661 cpl_table_new_column_array(gain,
"Variance", CPL_TYPE_DOUBLE, nwindows);
662 cpl_table_new_column_array(gain,
"Gain", CPL_TYPE_DOUBLE, nwindows);
664 if (!cpl_errorstate_is_equal(state)) {
665 cpl_table_delete(gain);
666 cpl_table_delete(flats);
672 double mron_value = cpl_array_get_median(aRon);
673 double mron_sdev = cpl_array_get_stdev(aRon);
679 for (ipair = 0; ipair < npairs; ++ipair) {
680 cpl_size iexposure = 2 * ipair;
681 cpl_size jexposure = iexposure + 1;
682 int iflat = cpl_table_get_int(flats,
"INDEX", iexposure, NULL);
683 int jflat = cpl_table_get_int(flats,
"INDEX", jexposure, NULL);
687 double exptime1 = cpl_table_get_double(flats,
"EXPTIME", iexposure, NULL);
688 double exptime2 = cpl_table_get_double(flats,
"EXPTIME", jexposure, NULL);
689 double flux1 = cpl_image_get_mean(flat1->
data);
690 double flux2 = cpl_image_get_mean(flat2->
data);
692 if (fabs((flux1 - flux2) / flux1) > aParams->
fluxtol) {
693 cpl_msg_warning(__func__,
"Inconsistent overall flux level (%.4f, " 694 "%.4f) detected for flat field pair (%d, %d) of IFU %u " 695 "with exposure time %.4f! Skipping this pair!",
696 flux1, flux2, iflat, jflat, aParams->
nifu, exptime1);
700 cpl_table_set_double(gain,
"ExpTime", kpair, 0.5 * (exptime1 + exptime2));
710 cpl_array *_signal = cpl_array_new(nwindows, CPL_TYPE_DOUBLE);
711 cpl_array *_variance = cpl_array_new(nwindows, CPL_TYPE_DOUBLE);
712 cpl_array *_gain = cpl_array_new(nwindows, CPL_TYPE_DOUBLE);
714 if (aParams->
sigma > 0.) {
717 for (iwindow = 0; iwindow < nwindows; ++iwindow) {
718 int xmin = cpl_table_get_int(aWindowList,
"Xmin", iwindow, NULL);
719 int ymin = cpl_table_get_int(aWindowList,
"Ymin", iwindow, NULL);
720 int xmax = cpl_table_get_int(aWindowList,
"Xmax", iwindow, NULL);
721 int ymax = cpl_table_get_int(aWindowList,
"Ymax", iwindow, NULL);
726 cpl_image *_flat1 = cpl_image_extract(flat1->
data,
727 xmin, ymin, xmax, ymax);
728 cpl_image *_flat2 = cpl_image_extract(flat2->
data,
729 xmin, ymin, xmax, ymax);
730 cpl_mask *mask = cpl_image_get_bpm(_flat1);
731 cpl_mask_or(mask, cpl_image_get_bpm(_flat2));
732 cpl_image_reject_from_mask(_flat2, mask);
734 double median1 = cpl_image_get_median(_flat1);
735 double median2 = cpl_image_get_median(_flat2);
736 double clip1 = aParams->
sigma * cpl_image_get_stdev(_flat1);
737 double clip2 = aParams->
sigma * cpl_image_get_stdev(_flat2);
744 cpl_mask_threshold_image_create(_flat1,
748 cpl_mask_threshold_image_create(_flat2,
757 cpl_mask_or(mask1, mask2);
758 cpl_mask_delete(mask2);
760 cpl_mask_or(mask, mask1);
761 cpl_mask_delete(mask1);
766 cpl_image_reject_from_mask(_flat2, mask);
767 cpl_image *dflat = cpl_image_subtract_create(_flat1, _flat2);
769 if ((median1 < kMuseSaturationLimit) && (median2 < kMuseSaturationLimit)) {
770 double s1 = cpl_image_get_mean(_flat1);
771 double s2 = cpl_image_get_mean(_flat2);
772 double s = 0.5 * (s1 + s2);
776 if ((cpl_array_get_double(aRon, iwindow, NULL) - mron_value) <= 5. * mron_sdev) {
777 double sdev = cpl_image_get_stdev(dflat);
778 v = 0.5 * sdev * sdev - mron_value * mron_value;
781 cpl_array_set_double(_signal, iwindow, s);
782 cpl_array_set_double(_variance, iwindow, v);
783 cpl_array_set_double(_gain, iwindow, g);
785 cpl_array_set_invalid(_signal, iwindow);
786 cpl_array_set_invalid(_variance, iwindow);
787 cpl_array_set_invalid(_gain, iwindow);
790 cpl_image_delete(dflat);
791 cpl_image_delete(_flat2);
792 cpl_image_delete(_flat1);
804 cpl_image *dflat = cpl_image_subtract_create(flat1->
data, flat2->
data);
807 cpl_image_reject_from_mask(flat1->
data, cpl_image_get_bpm(dflat));
808 cpl_image_reject_from_mask(flat2->
data, cpl_image_get_bpm(dflat));
811 for (iwindow = 0; iwindow < nwindows; ++iwindow) {
812 int xmin = cpl_table_get_int(aWindowList,
"Xmin", iwindow, NULL);
813 int ymin = cpl_table_get_int(aWindowList,
"Ymin", iwindow, NULL);
814 int xmax = cpl_table_get_int(aWindowList,
"Xmax", iwindow, NULL);
815 int ymax = cpl_table_get_int(aWindowList,
"Ymax", iwindow, NULL);
817 double median1 = cpl_image_get_median_window(flat1->
data,
818 xmin, ymin, xmax, ymax);
819 double median2 = cpl_image_get_median_window(flat2->
data,
820 xmin, ymin, xmax, ymax);
822 if ((median1 < kMuseSaturationLimit) && (median2 < kMuseSaturationLimit)) {
823 double s1 = cpl_image_get_mean_window(flat1->
data,
824 xmin, ymin, xmax, ymax);
825 double s2 = cpl_image_get_mean_window(flat2->
data,
826 xmin, ymin, xmax, ymax);
827 double s = 0.5 * (s1 + s2);
831 if ((cpl_array_get_double(aRon, iwindow, NULL) - mron_value) <= 5. * mron_sdev) {
832 double sdev = cpl_image_get_stdev_window(dflat,
833 xmin, ymin, xmax, ymax);
834 v = 0.5 * sdev * sdev - mron_value * mron_value;
837 cpl_array_set_double(_signal, iwindow, s);
838 cpl_array_set_double(_variance, iwindow, v);
839 cpl_array_set_double(_gain, iwindow, g);
841 cpl_array_set_invalid(_signal, iwindow);
842 cpl_array_set_invalid(_variance, iwindow);
843 cpl_array_set_invalid(_gain, iwindow);
847 cpl_image_delete(dflat);
851 cpl_table_set_array(gain,
"Signal", kpair, _signal);
852 cpl_table_set_array(gain,
"Variance", kpair, _variance);
853 cpl_table_set_array(gain,
"Gain", kpair, _gain);
856 cpl_array_delete(_signal);
857 cpl_array_delete(_variance);
858 cpl_array_delete(_gain);
861 cpl_table_delete(flats);
864 cpl_table_erase_invalid_rows(gain);
865 if (cpl_table_get_nrow(gain) == 0) {
866 cpl_table_delete(gain);
888 static cpl_error_code
889 muse_lingain_quadrant_get_gain(muse_gain_fit_t *aGain,
890 const cpl_table *aGainTable,
891 const cpl_array *aWindows,
898 cpl_size nrows = cpl_table_get_nrow(aGainTable);
899 cpl_size nwindows = cpl_array_get_size(aWindows);
900 cpl_size ndata = nrows * nwindows;
902 cpl_array *qsignal = cpl_array_new(ndata, CPL_TYPE_DOUBLE);
903 cpl_array *qgain = cpl_array_new(ndata, CPL_TYPE_DOUBLE);
906 for (irow = 0; irow < nrows; ++irow) {
908 cpl_size stride = irow * nwindows;
910 const cpl_array *_qsignal = cpl_table_get_array(aGainTable,
912 const cpl_array *_qgain = cpl_table_get_array(aGainTable,
915 for (iwindow = 0; iwindow < nwindows; ++iwindow) {
916 cpl_size jwindow = cpl_array_get_cplsize(aWindows, iwindow, NULL);
920 double s = cpl_array_get_double(_qsignal, jwindow, &snull);
921 double g = cpl_array_get_double(_qgain, jwindow, &gnull);
923 cpl_size idata = stride + iwindow;
924 if (((snull == 0) && (gnull == 0)) && (s > 0.)) {
925 cpl_array_set_double(qsignal, idata, log10(s));
926 cpl_array_set_double(qgain, idata, g);
928 cpl_array_set_invalid(qsignal, idata);
929 cpl_array_set_invalid(qgain, idata);
937 cpl_array *gx = cpl_array_new(nbins, CPL_TYPE_DOUBLE);
938 cpl_array *gy = cpl_array_new(nbins, CPL_TYPE_DOUBLE);
941 for (ibin = 0; ibin < nbins; ++ibin) {
944 double bmid = 0.5 * (bmin + bmax);
947 cpl_array_set_invalid(gx, ibin);
948 cpl_array_set_invalid(gy, ibin);
953 cpl_array_set(gx, ibin, pow(10., bmid));
955 cpl_array *gdata = cpl_array_duplicate(qgain);
958 for (idata = 0; idata < ndata; ++idata) {
962 double s = cpl_array_get_double(qsignal, idata, &invalid);
964 if (invalid || (s < bmin) || (s >= bmax)) {
965 cpl_array_set_invalid(gdata, idata);
970 if (cpl_array_count_invalid(gdata) == ndata) {
971 cpl_array_set_invalid(gx, ibin);
972 cpl_array_set_invalid(gy, ibin);
973 cpl_array_delete(gdata);
980 double median = cpl_array_get_median(gdata);
981 double stdev = cpl_array_get_stdev(gdata);
984 for (idata = 0; idata < ndata; ++idata) {
988 double g = cpl_array_get_double(gdata, idata, &invalid);
990 if (invalid || (g <= median - stdev) || (g >= median + stdev)) {
991 cpl_array_set_invalid(gdata, idata);
995 if (cpl_array_count_invalid(gdata) == ndata) {
996 cpl_array_set_invalid(gx, ibin);
997 cpl_array_set_invalid(gy, ibin);
1000 double gmean = cpl_array_get_mean(gdata);
1002 cpl_array_set_double(gy, ibin, gmean);
1005 cpl_array_delete(gdata);
1008 cpl_array_delete(qgain);
1009 cpl_array_delete(qsignal);
1016 cpl_size nx = cpl_array_get_size(gx);
1017 cpl_size ny = cpl_array_get_size(gy);
1019 if ((nx == 0) || (ny == 0)) {
1021 cpl_msg_debug(__func__,
"Fitting the gain relation failed. No data " 1023 cpl_array_delete(gy);
1024 cpl_array_delete(gx);
1026 return CPL_ERROR_DATA_NOT_FOUND;
1030 cpl_msg_debug(__func__,
"Fitting the gain relation failed. Gain " 1031 "and signal data sets do not match!");
1032 cpl_array_delete(gy);
1033 cpl_array_delete(gx);
1035 return CPL_ERROR_INCOMPATIBLE_INPUT;
1040 cpl_msg_debug(__func__,
"Fitting the gain relation failed. Insufficient " 1041 "data points; cannot fit a first order polynomial!");
1042 cpl_array_delete(gy);
1043 cpl_array_delete(gx);
1045 return CPL_ERROR_DATA_NOT_FOUND;
1049 double *_x = cpl_array_get_data_double(gx);
1050 double *_y = cpl_array_get_data_double(gy);
1052 cpl_matrix *x = cpl_matrix_wrap(1, nx, _x);
1053 cpl_vector *y = cpl_vector_wrap(ny, _y);
1054 cpl_boolean symetric = CPL_FALSE;
1055 const cpl_size mindeg = 0;
1056 const cpl_size maxdeg = 1;
1058 cpl_polynomial *gfit = cpl_polynomial_new(1);
1059 cpl_error_code ecode = cpl_polynomial_fit(gfit, x, &symetric, y, NULL,
1060 CPL_FALSE, &mindeg, &maxdeg);
1062 cpl_vector_unwrap(y);
1063 cpl_matrix_unwrap(x);
1065 if (ecode != CPL_ERROR_NONE) {
1066 cpl_polynomial_delete(gfit);
1067 cpl_array_delete(gy);
1068 cpl_array_delete(gx);
1070 cpl_msg_debug(__func__,
"Fitting the gain relation failed. Fitting " 1071 "first order polynomial to gain data failed!");
1074 return CPL_ERROR_ILLEGAL_OUTPUT;
1079 for (ix = 0; ix < nx; ++ix) {
1080 double rx = cpl_array_get_double(gx, ix, NULL);
1081 double ry = cpl_array_get_double(gy, ix, NULL);
1082 double ryf = cpl_polynomial_eval_1d(gfit, rx, NULL);
1083 rms += (ry - ryf) * (ry - ryf);
1087 cpl_array_delete(gy);
1088 cpl_array_delete(gx);
1090 const cpl_size order = 0;
1091 aGain->gain = cpl_polynomial_get_coeff(gfit, &order);
1094 cpl_polynomial_delete(gfit);
1096 return CPL_ERROR_NONE;
1112 static cpl_error_code
1113 muse_lingain_quadrant_get_nonlinearity(muse_linearity_fit_t *aFit,
1114 const cpl_table *aGainTable,
1115 const cpl_array *aWindows,
1126 cpl_size ntimes = cpl_table_get_nrow(aGainTable);
1128 cpl_error_code ecode = cpl_table_get_column_maxpos(aGainTable,
"ExpTime",
1130 if (ecode != CPL_ERROR_NONE) {
1131 cpl_msg_debug(__func__,
"Maximum exposure time cannot be determined from " 1132 "linearity data set of IFU %u!", aParams->
nifu);
1136 if (irow != ntimes - 1) {
1137 cpl_msg_debug(__func__,
"Linearity data set of IFU %u is not sorted in " 1138 "ascending order of the exposure time!", aParams->
nifu);
1139 return CPL_ERROR_ILLEGAL_INPUT;
1147 muse_lingain_quadrant_extract_data(aGainTable, aWindows,
"Signal", irow);
1148 cpl_array_logarithm(sigref, 10.);
1150 cpl_size nsignal = cpl_array_get_size(sigref);
1151 cpl_size nlevels = (cpl_size)((aParams->
ctsmax - aParams->
ctsmin) /
1154 cpl_array **counts = cpl_calloc(nlevels,
sizeof *counts);
1155 const double qgain = log10(aGain);
1156 cpl_size klevel = 0;
1158 for (ilevel = 0; ilevel < nlevels; ++ilevel) {
1159 double shigh = aParams->
ctsmax - ilevel * aParams->
ctsbin;
1160 double slow = shigh - aParams->
ctsbin;
1167 cpl_size *sindex = cpl_malloc(nsignal *
sizeof *sindex);
1170 for (idata = 0; idata < nsignal; ++idata) {
1172 double s = cpl_array_get_double(sigref, idata, &invalid);
1174 if (!invalid && ((s >= slow) && (s < shigh))) {
1175 sindex[kdata] = cpl_array_get_cplsize(aWindows, idata, NULL);
1182 cpl_msg_debug(__func__,
"No valid measurement found within signal " 1183 "range [%.6f, %.6f[ [log10(counts)] of the linearity " 1184 "data set of IFU %u!", slow + qgain, shigh + qgain,
1194 cpl_array *smean = cpl_array_new(ntimes, CPL_TYPE_DOUBLE);
1195 cpl_array *windex = cpl_array_wrap_cplsize(sindex, kdata);
1197 for (itime = 0; itime < ntimes; ++itime) {
1198 cpl_array *qsignal =
1199 muse_lingain_quadrant_extract_data(aGainTable, windex,
1201 if (qsignal && (cpl_array_count_invalid(qsignal) < kdata)) {
1202 double mean = cpl_array_get_mean(qsignal);
1203 cpl_array_set_double(smean, itime, mean);
1205 cpl_array_set_invalid(smean, itime);
1207 cpl_array_delete(qsignal);
1209 cpl_array_unwrap(windex);
1213 counts[klevel] = smean;
1219 cpl_array_delete(sigref);
1222 cpl_array *exptime = cpl_array_new(ntimes, CPL_TYPE_DOUBLE);
1225 for (itime = 0; itime < ntimes; ++itime) {
1226 double _exptime = cpl_table_get_double(aGainTable,
"ExpTime",
1228 cpl_array_set_double(exptime, itime, _exptime);
1230 cpl_array_add_scalar(exptime, aParams->
toffset);
1233 cpl_array **crate = cpl_calloc(klevel,
sizeof *crate);
1234 for (ilevel = 0; ilevel < klevel; ++ilevel) {
1235 crate[ilevel] = cpl_array_duplicate(counts[ilevel]);
1236 cpl_array_divide(crate[ilevel], exptime);
1246 double linmin = pow(10., aParams->
linearmin - qgain);
1247 double linmax = pow(10., aParams->
linearmax - qgain);
1249 cpl_array **cresidual = cpl_malloc(klevel *
sizeof *cresidual);
1251 for (ilevel = 0; ilevel < klevel; ++ilevel) {
1252 cpl_array *_exptime = cpl_array_duplicate(exptime);
1253 cpl_array *_crate = cpl_array_duplicate(crate[ilevel]);
1255 for (itime = 0; itime < ntimes; ++itime) {
1257 double cnts = cpl_array_get_double(counts[ilevel], itime, &invalid);
1259 if (invalid || ((cnts < linmin) || (cnts >= linmax))) {
1260 cpl_array_set_invalid(_exptime, itime);
1261 cpl_array_set_invalid(_crate, itime);
1265 if (cpl_array_get_size(_crate) == cpl_array_count_invalid(_crate)) {
1269 cresidual[ilevel] = NULL;
1271 cpl_array_delete(_crate);
1272 cpl_array_delete(_exptime);
1285 cpl_size nt = cpl_array_get_size(_exptime);
1286 cpl_size nc = cpl_array_get_size(_crate);
1291 cpl_msg_debug(__func__,
"Fitting a constant to the count rate data " 1292 "failed. Exposure time and count rate data sets do not " 1294 cpl_array_delete(_crate);
1295 cpl_array_delete(_exptime);
1297 muse_vfree((
void **)cresidual, klevel, (muse_free_func)cpl_array_delete);
1298 muse_vfree((
void **)crate, klevel, (muse_free_func)cpl_array_delete);
1299 muse_vfree((
void **)counts, nlevels, (muse_free_func)cpl_array_delete);
1301 cpl_array_delete(exptime);
1303 return CPL_ERROR_INCOMPATIBLE_INPUT;
1307 double *_et = cpl_array_get_data_double(_exptime);
1308 double *_cr = cpl_array_get_data_double(_crate);
1310 cpl_matrix *et = cpl_matrix_wrap(1, nt, _et);
1311 cpl_vector *cr = cpl_vector_wrap(nc, _cr);
1312 cpl_boolean symetric = CPL_FALSE;
1313 const cpl_size mindeg = 0;
1314 const cpl_size maxdeg = 0;
1316 cpl_polynomial *crfit = cpl_polynomial_new(1);
1317 cpl_error_code ecode = cpl_polynomial_fit(crfit, et, &symetric, cr, NULL,
1318 CPL_FALSE, &mindeg, &maxdeg);
1320 cpl_vector_unwrap(et);
1321 cpl_matrix_unwrap(cr);
1323 cpl_array_delete(_crate);
1324 cpl_array_delete(_exptime);
1326 if (ecode != CPL_ERROR_NONE) {
1327 cpl_polynomial_delete(crfit);
1328 cpl_msg_debug(__func__,
"Fitting a constant to the count rate level " 1329 "failed. Fitting a zero order polynomial to count rate " 1332 muse_vfree((
void **)cresidual, klevel, (muse_free_func)cpl_array_delete);
1333 muse_vfree((
void **)crate, klevel, (muse_free_func)cpl_array_delete);
1334 muse_vfree((
void **)counts, nlevels, (muse_free_func)cpl_array_delete);
1336 cpl_array_delete(exptime);
1338 return CPL_ERROR_ILLEGAL_OUTPUT;
1341 const cpl_size order = 0;
1342 double crlevel = cpl_polynomial_get_coeff(crfit, &order);
1343 cpl_polynomial_delete(crfit);
1345 double crlevel = cpl_array_get_mean(_crate);
1347 cpl_array_delete(_crate);
1348 cpl_array_delete(_exptime);
1351 cresidual[ilevel] = cpl_array_duplicate(crate[ilevel]);
1352 cpl_array_subtract_scalar(cresidual[ilevel], crlevel);
1353 cpl_array_divide(cresidual[ilevel], crate[ilevel]);
1357 cpl_array_delete(exptime);
1369 cpl_size nbins = (cpl_size)((signalmax - signalmin) / signalbin + 0.5);
1370 cpl_array *lx = cpl_array_new(nbins, CPL_TYPE_DOUBLE);
1371 cpl_array *ly = cpl_array_new(nbins, CPL_TYPE_DOUBLE);
1372 cpl_array *lyse = cpl_array_new(nbins, CPL_TYPE_DOUBLE);
1375 for (ibin = 0; ibin < nbins; ++ibin) {
1376 double bmin = signalmin + ibin * signalbin;
1377 double bmax = bmin + signalbin;
1378 double bmid = 0.5 * (bmin + bmax);
1380 cpl_array_set_double(lx, ibin, bmid);
1382 bmin = pow(10., bmin);
1383 bmax = pow(10., bmax);
1385 double *rdata = cpl_malloc((klevel * ntimes) *
sizeof *rdata);
1386 cpl_size kresidual = 0;
1387 for (ilevel = 0; ilevel < klevel; ++ilevel) {
1388 cpl_array *residual = cresidual[ilevel];
1390 if (residual != NULL) {
1391 for (itime = 0; itime < ntimes; ++itime) {
1394 double _counts = cpl_array_get_double(counts[ilevel], itime, &cinvalid);
1395 double _residual = cpl_array_get_double(residual, itime, &rinvalid);
1397 if (!(cinvalid || rinvalid) && ((_counts >= bmin) && (_counts < bmax))) {
1398 rdata[kresidual] = _residual;
1404 if (kresidual == 0) {
1405 cpl_array_set_invalid(lx, ibin);
1406 cpl_array_set_invalid(ly, ibin);
1407 cpl_array_set_invalid(lyse, ibin);
1409 cpl_array *r = cpl_array_wrap_double(rdata, kresidual);
1411 double rmean = cpl_array_get_mean(r);
1412 double rsdev = cpl_array_get_stdev(r);
1413 cpl_array_set_double(ly, ibin, rmean);
1414 cpl_array_set_double(lyse, ibin, rsdev);
1416 cpl_array_unwrap(r);
1421 muse_vfree((
void **)cresidual, klevel, (muse_free_func)cpl_array_delete);
1422 muse_vfree((
void **)crate, klevel, (muse_free_func)cpl_array_delete);
1423 muse_vfree((
void **)counts, nlevels, (muse_free_func)cpl_array_delete);
1430 cpl_size nx = cpl_array_get_size(lx);
1431 cpl_size ny = cpl_array_get_size(ly);
1433 if ((nx == 0) || (ny == 0)) {
1435 cpl_msg_debug(__func__,
"Fitting the non-linearity relation failed. " 1436 "No data available!");
1437 cpl_array_delete(lyse);
1438 cpl_array_delete(ly);
1439 cpl_array_delete(lx);
1441 return CPL_ERROR_DATA_NOT_FOUND;
1445 cpl_msg_debug(__func__,
"Fitting the non-linearity relation failed. " 1446 "Linearity residual and signal data sets do not match!");
1447 cpl_array_delete(lyse);
1448 cpl_array_delete(ly);
1449 cpl_array_delete(lx);
1451 return CPL_ERROR_INCOMPATIBLE_INPUT;
1455 if (nx < aParams->order + 2) {
1456 cpl_msg_debug(__func__,
"Fitting the non-linearity relation failed. " 1457 "Insufficient data points; cannot fit a %d order " 1458 "polynomial!", aParams->
order);
1459 cpl_array_delete(lyse);
1460 cpl_array_delete(ly);
1461 cpl_array_delete(lx);
1463 return CPL_ERROR_DATA_NOT_FOUND;
1467 double *_x = cpl_array_get_data_double(lx);
1468 double *_y = cpl_array_get_data_double(ly);
1469 double *_yse = cpl_array_get_data_double(lyse);
1471 cpl_matrix *x = cpl_matrix_wrap(1, nx, _x);
1472 cpl_vector *y = cpl_vector_wrap(ny, _y);
1473 cpl_vector *yse = cpl_vector_wrap(ny, _yse);
1474 cpl_boolean symetric = CPL_FALSE;
1475 const cpl_size mindeg = 0;
1476 const cpl_size maxdeg = aParams->
order;
1478 cpl_polynomial *lfit = cpl_polynomial_new(1);
1479 ecode = cpl_polynomial_fit(lfit, x, &symetric, y, NULL, CPL_FALSE,
1482 cpl_vector_unwrap(yse);
1483 cpl_vector_unwrap(y);
1484 cpl_matrix_unwrap(x);
1486 if (ecode != CPL_ERROR_NONE) {
1487 cpl_polynomial_delete(lfit);
1488 cpl_array_delete(lyse);
1489 cpl_array_delete(ly);
1490 cpl_array_delete(lx);
1492 cpl_msg_debug(__func__,
"Fitting the non-linearity relation failed. " 1493 "Fitting %d order polynomial to residual non-linearity " 1494 "data failed!", aParams->
order);
1495 aFit->coefficients = NULL;
1499 return CPL_ERROR_ILLEGAL_OUTPUT;
1503 cpl_size ncoeff = aParams->
order + 1;
1504 cpl_array *coefficients = cpl_array_new(ncoeff, CPL_TYPE_DOUBLE);
1507 for (icoeff = 0; icoeff < ncoeff; ++icoeff) {
1508 const cpl_size order = icoeff;
1509 double coeff = cpl_polynomial_get_coeff(lfit, &order);
1510 cpl_array_set_double(coefficients, icoeff, coeff);
1515 for (ix = 0; ix < nx; ++ix) {
1516 double rx = cpl_array_get_double(lx, ix, NULL);
1517 double ry = cpl_array_get_double(ly, ix, NULL);
1518 double ryf = cpl_polynomial_eval_1d(lfit, rx, NULL);
1519 rms += (ry - ryf) * (ry - ryf);
1522 cpl_polynomial_delete(lfit);
1524 aFit->coefficients = coefficients;
1525 aFit->range[0] = cpl_array_get_min(lx);
1526 aFit->range[1] = cpl_array_get_max(lx);
1529 cpl_array_delete(lyse);
1530 cpl_array_delete(ly);
1531 cpl_array_delete(lx);
1533 return CPL_ERROR_NONE;
1560 muse_lingain_compute_gain(
const cpl_table *aGainTable, cpl_table *aWindowList,
1563 const unsigned char nquadrants = 4;
1564 unsigned char iquadrant;
1566 cpl_table *gain = cpl_table_new(nquadrants);
1572 cpl_table_new_column(gain,
"Gain", CPL_TYPE_DOUBLE);
1573 cpl_table_new_column(gain,
"FitRms", CPL_TYPE_DOUBLE);
1575 for (iquadrant = 0; iquadrant < nquadrants; ++iquadrant) {
1576 unsigned int quadrant = iquadrant + 1;
1579 cpl_array *windex = muse_lingain_quadrant_get_windows(aWindowList, quadrant);
1582 cpl_msg_warning(__func__,
"No measurement windows defined for " 1583 "quadrant %hhu of IFU %u!", quadrant, aParams->
nifu);
1584 cpl_table_set_invalid(gain,
"Gain", iquadrant);
1585 cpl_table_set_invalid(gain,
"FitRms", iquadrant);
1592 muse_gain_fit_t qgain = {0., 0.};
1593 cpl_error_code ecode = muse_lingain_quadrant_get_gain(&qgain, aGainTable,
1595 if (ecode != CPL_ERROR_NONE) {
1596 cpl_msg_warning(__func__,
"Detector gain value for quadrant %hhu of " 1597 "IFU %u is invalid!", quadrant, aParams->
nifu);
1598 cpl_table_set_invalid(gain,
"Gain", iquadrant);
1599 cpl_table_set_invalid(gain,
"FitRms", iquadrant);
1601 cpl_msg_info(__func__,
"Detector gain value for quadrant %hhu of " 1602 "IFU %u: %.6f [counts/ADU]", quadrant, aParams->
nifu,
1606 cpl_table_set_double(gain,
"Gain", iquadrant, qgain.gain);
1607 cpl_table_set_double(gain,
"FitRms", iquadrant, qgain.rms);
1608 cpl_array_delete(windex);
1631 muse_lingain_compute_nonlinearity(
const cpl_table *aGainTable,
1632 cpl_table *aWindowList,
1633 const cpl_table *aGain,
1636 const unsigned char nquadrants = 4;
1637 unsigned char iquadrant;
1639 cpl_table *linearity = cpl_table_new(nquadrants);
1640 cpl_table_new_column(linearity,
"Gain", CPL_TYPE_DOUBLE);
1641 cpl_table_new_column(linearity,
"SignalMin", CPL_TYPE_DOUBLE);
1642 cpl_table_new_column(linearity,
"SignalMax", CPL_TYPE_DOUBLE);
1643 cpl_table_new_column_array(linearity,
"Coefficients", CPL_TYPE_DOUBLE,
1644 aParams->
order + 1);
1645 cpl_table_new_column(linearity,
"FitRms", CPL_TYPE_DOUBLE);
1647 cpl_table_set_column_unit(linearity,
"Gain",
"counts/ADU)");
1648 cpl_table_set_column_unit(linearity,
"SignalMin",
"log10(ADU)");
1649 cpl_table_set_column_unit(linearity,
"SignalMax",
"log10(ADU)");
1651 for (iquadrant = 0; iquadrant < nquadrants; ++iquadrant) {
1652 unsigned int quadrant = iquadrant + 1;
1655 cpl_array *windex = muse_lingain_quadrant_get_windows(aWindowList, quadrant);
1658 cpl_msg_warning(__func__,
"No measurement windows defined for " 1659 "quadrant %hhu of IFU %u!", quadrant, aParams->
nifu);
1667 double qgain = cpl_table_get_double(aGain,
"Gain", iquadrant, &invalid);
1670 cpl_msg_warning(__func__,
"Got invalid gain for quadrant %hhu of IFU " 1671 "%u! Skipping non-linearity computation!", quadrant,
1673 cpl_array_delete(windex);
1677 muse_linearity_fit_t qlinearity = {NULL, {0., 0.}, 0.};
1678 cpl_error_code ecode =
1679 muse_lingain_quadrant_get_nonlinearity(&qlinearity, aGainTable,
1680 windex, qgain, aParams);
1682 if (ecode != CPL_ERROR_NONE) {
1683 cpl_msg_warning(__func__,
"Computation of the detector non-linearity " 1684 "for quadrant %hhu of IFU %u failed!", quadrant,
1686 cpl_table_set_invalid(linearity,
"Gain", iquadrant);
1687 cpl_table_set_invalid(linearity,
"SignalMin", iquadrant);
1688 cpl_table_set_invalid(linearity,
"SignalMax", iquadrant);
1689 cpl_table_set_invalid(linearity,
"Coefficients", iquadrant);
1690 cpl_table_set_invalid(linearity,
"FitRms", iquadrant);
1692 cpl_msg_info(__func__,
"RMS of residual non-linearity model over " 1693 "signal range [%.4f, %.4f] [log10(ADU)] for quadrant " 1694 "%hhu of IFU %u: %.6e", qlinearity.range[0],
1695 qlinearity.range[1], quadrant, aParams->
nifu,
1698 cpl_table_set_double(linearity,
"Gain", iquadrant, qgain);
1699 cpl_table_set_double(linearity,
"SignalMin", iquadrant,
1700 qlinearity.range[0]);
1701 cpl_table_set_double(linearity,
"SignalMax", iquadrant,
1702 qlinearity.range[1]);
1703 cpl_table_set_array(linearity,
"Coefficients", iquadrant,
1704 qlinearity.coefficients);
1705 cpl_table_set_double(linearity,
"FitRms", iquadrant,
1708 muse_linearity_fit_clear(&qlinearity);
1709 cpl_array_delete(windex);
1729 if (!muse_lingain_validate_parameters(aParams)) {
1738 cpl_msg_error(__func__,
"%s could not be loaded!", MUSE_TAG_TRACE_TABLE);
1752 cpl_msg_error(__func__,
"Loading and basic processing of the raw input " 1753 "images failed for IFU %u!", aParams->
nifu);
1754 cpl_table_delete(trace);
1760 cpl_table *exposures = muse_lingain_sort_images(images);
1762 cpl_msg_error(__func__,
"Creating a sorted list of exposures failed " 1763 "for IFU %u!", aParams->
nifu);
1765 cpl_table_delete(trace);
1771 cpl_table *windowlist = muse_lingain_create_windowlist(trace, exposures,
1775 cpl_msg_error(__func__,
"Creating the list of measurement windows failed " 1776 "for IFU %u!", aParams->
nifu);
1777 cpl_table_delete(exposures);
1779 cpl_table_delete(trace);
1782 cpl_table_delete(trace);
1785 cpl_array *ron = muse_lingain_compute_ron(windowlist, exposures, images,
1788 cpl_msg_error(__func__,
"Calculating the RON for individual measurement " 1789 "windows failed for IFU %u!", aParams->
nifu);
1790 cpl_table_delete(windowlist);
1791 cpl_table_delete(exposures);
1798 cpl_table *gaindata = muse_lingain_window_get_gain(windowlist, ron,
1802 cpl_msg_error(__func__,
"Calculating the gain for individual measurement " 1803 "windows failed for IFU %u!", aParams->
nifu);
1804 cpl_array_delete(ron);
1805 cpl_table_delete(windowlist);
1806 cpl_table_delete(exposures);
1810 cpl_array_delete(ron);
1811 cpl_table_delete(exposures);
1815 cpl_table *gain = muse_lingain_compute_gain(gaindata, windowlist, aParams);
1821 if (cpl_table_count_invalid(gain,
"Gain") != 0) {
1822 cpl_msg_error(__func__,
"Calculating the detector gain failed for IFU %u!",
1824 cpl_table_delete(gain);
1825 cpl_table_delete(gaindata);
1826 cpl_table_delete(windowlist);
1832 cpl_table *linearity = muse_lingain_compute_nonlinearity(gaindata, windowlist,
1834 if (cpl_table_has_invalid(linearity,
"Coefficients")) {
1835 cpl_msg_error(__func__,
"Computing the detector residual non-linearity " 1836 "failed for IFU %u!", aParams->
nifu);
1837 cpl_table_delete(linearity);
1838 cpl_table_delete(gain);
1839 cpl_table_delete(gaindata);
1840 cpl_table_delete(windowlist);
1845 cpl_table_delete(gaindata);
1846 cpl_table_delete(windowlist);
1849 cpl_propertylist *header =
1854 cpl_propertylist_erase_regexp(header,
1855 "^SIMPLE$|^BITPIX$|^NAXIS|^EXTEND$|^XTENSION$|" 1856 "^DATASUM$|^DATAMIN$|^DATAMAX$|^DATAMD5$|" 1857 "^PCOUNT$|^GCOUNT$|^HDUVERS$|^BLANK$|" 1858 "^BZERO$|^BSCALE$|^BUNIT$|^CHECKSUM$|^INHERIT$|" 1859 "^PIPEFILE$|^ESO PRO ", 0);
1861 const unsigned char nquadrant = 4;
1862 unsigned char iquadrant;
1863 for (iquadrant = 0; iquadrant < nquadrant; ++iquadrant) {
1864 char keyword[KEYWORD_LENGTH];
1865 char comment[KEYWORD_LENGTH];
1866 unsigned char _iquadrant = iquadrant + 1;
1867 int ncoefficient = (int)cpl_table_get_column_depth(linearity,
1869 snprintf(keyword, KEYWORD_LENGTH,
"ESO DET OUT%d GAIN", _iquadrant);
1870 cpl_propertylist_update_double(header, keyword,
1871 cpl_table_get_double(gain,
"Gain",
1878 cpl_propertylist_set_comment(header, keyword,
1879 "[e-/ADU] Conversion ADUs to electrons");
1880 snprintf(keyword, KEYWORD_LENGTH,
"ESO DET OUT%d CONAD", _iquadrant);
1881 cpl_propertylist_set_comment(header, keyword,
1882 "[ADU/e-] Conversion electrons to ADUs");
1885 snprintf(keyword, KEYWORD_LENGTH, MUSE_HDR_NONLINn_LLO, _iquadrant);
1886 cpl_propertylist_append_double(header, keyword,
1887 cpl_table_get_double(linearity,
"SignalMin",
1889 cpl_propertylist_set_comment(header, keyword,
1890 "[log10(ADU)] Minimum signal used for fit");
1892 snprintf(keyword, KEYWORD_LENGTH, MUSE_HDR_NONLINn_LHI, _iquadrant);
1893 cpl_propertylist_append_double(header, keyword,
1894 cpl_table_get_double(linearity,
"SignalMax",
1896 cpl_propertylist_set_comment(header, keyword,
1897 "[log10(ADU)] Maximum signal used for fit");
1899 snprintf(keyword, KEYWORD_LENGTH, MUSE_HDR_NONLINn_ORDER, _iquadrant);
1900 cpl_propertylist_append_double(header, keyword, ncoefficient - 1);
1901 cpl_propertylist_set_comment(header, keyword,
1902 "Order of the polynomial fit");
1904 const cpl_array *coefficients =
1905 cpl_table_get_array(linearity,
"Coefficients", iquadrant);
1907 for (icoefficient = 0; icoefficient < ncoefficient; ++icoefficient) {
1908 double c = cpl_array_get_double(coefficients, icoefficient, NULL);
1909 snprintf(keyword, KEYWORD_LENGTH, MUSE_HDR_NONLINn_COEFFo,
1910 _iquadrant, (
unsigned char)icoefficient);
1911 cpl_propertylist_append_double(header, keyword, c);
1913 snprintf(comment, KEYWORD_LENGTH,
1914 "%d order coefficient of the polynomial fit", icoefficient);
1915 cpl_propertylist_set_comment(header, keyword, comment);
1918 snprintf(keyword, KEYWORD_LENGTH, QC_LINGAIN_GFITi_RMS, _iquadrant);
1919 cpl_propertylist_append_double(header, keyword,
1920 cpl_table_get_double(gain,
"FitRms",
1922 snprintf(keyword, KEYWORD_LENGTH, QC_LINGAIN_NLFITi_RMS, _iquadrant);
1923 cpl_propertylist_append_double(header, keyword,
1924 cpl_table_get_double(linearity,
"FitRms",
1927 cpl_table_erase_column(linearity,
"FitRms");
1932 cpl_propertylist_delete(header);
1933 cpl_table_delete(linearity);
1934 cpl_table_delete(gain);
const char * muse_pfits_get_dpr_type(const cpl_propertylist *aHeaders)
find out the DPR type
cpl_error_code muse_cplarray_erase_invalid(cpl_array *aArray)
Erase all invalid values from an array.
muse_imagelist * muse_basicproc_load(muse_processing *aProcessing, unsigned char aIFU, muse_basicproc_params *aBPars)
Load the raw input files from disk and do basic processing.
cpl_polynomial ** muse_trace_table_get_polys_for_slice(const cpl_table *aTable, const unsigned short aSlice)
construct polynomial from the trace table entry for the given slice
Structure definition for a collection of muse_images.
cpl_size * muse_quadrants_get_window(const muse_image *aImage, unsigned char aQuadrant)
Determine the data window of a given quadrant on the CCD.
cpl_image * data
the data extension
Structure to hold the parameters of the muse_lingain recipe.
void muse_imagelist_delete(muse_imagelist *aList)
Free the memory of the MUSE image list.
double signalmin
Minimum signal value in log(ADU) used for the gain analysis and the non-linearity polynomial model...
double linearmin
Lower limit of desired linear range in log10(counts).
double sigma
Sigma value used for signal value clipping.
int nifu
IFU to handle. If set to 0, all IFUs are processed serially. If set to -1, all IFUs are processed in ...
Structure definition of MUSE three extension FITS file.
void muse_basicproc_params_delete(muse_basicproc_params *aBPars)
Free a structure of basic processing parameters.
cpl_propertylist * header
the FITS header
double signalbin
Size of a signal bin in log10(ADU) used for the gain analysis and the non-linearity polynomial model...
void muse_trace_polys_delete(cpl_polynomial *aPolys[])
Delete the multi-polynomial array created in relation to tracing.
muse_image * muse_imagelist_get(muse_imagelist *aList, unsigned int aIdx)
Get the muse_image of given list index.
double signalmax
Maximum signal value in log(ADU) used for the gain analysis and the non-linearity polynomial model...
double fluxtol
Tolerance value for the overall flux consistency check of a pair of flat fields. The value is the max...
int order
Order of the polynomial used to fit the non-linearity residuals.
double toffset
Exposure time offset in seconds to apply to linearity flat fields.
double gainsigma
Sigma value for gain value clipping.
double linearmax
Upper limit of desired linear range in log10(counts).
cpl_table * muse_processing_load_ctable(muse_processing *aProcessing, const char *aTag, unsigned char aIFU)
Load a CPL table according to its tag and IFU/channel number.
double ctsmin
Minimum signal value in log(counts) to consider for the non-linearity analysis.
double ctsbin
Size of a signal bin in log10(counts) used for the non-linearity analysis.
double muse_pfits_get_mjdobs(const cpl_propertylist *aHeaders)
find out the Julian Date of the observation
cpl_error_code muse_image_reject_from_dq(muse_image *aImage)
Reject pixels of a muse_image depending on its DQ data.
int xborder
Extra offset from the detector edge used for the selection of slices.
double ctsmax
Maximum signal value in log(counts) to consider for the non-linearity analysis.
muse_basicproc_params * muse_basicproc_params_new_from_propertylist(const cpl_propertylist *aHeader)
Create a structure of basic processing parameters from a FITS header.
double muse_pfits_get_exptime(const cpl_propertylist *aHeaders)
find out the exposure time
int xgap
Extra offset from tracing edge.
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.
Structure of basic processing parameters.
double gainlimit
Minimum signal value [ADU] used for fitting the gain relation.
int ybox
Size of windows along the traces of the slices.
cpl_frame * muse_frameset_find_master(const cpl_frameset *aFrames, const char *aTag, unsigned char aIFU)
find the master frame according to its CCD number and tag