30 #include <cpl_image.h>
31 #include <cpl_vector.h>
32 #include <cpl_matrix.h>
34 #include <cpl_parameterlist.h>
47 #include "gimessages.h"
49 #include "gilocalize.h"
66 static const cxchar* _task =
"giraffe_localize_spectra";
75 GILOCALIZE_HALF_WIDTH,
79 typedef enum GiLocalizeMethod GiLocalizeMethod;
86 enum GiThresholdMethod
88 GILOCALIZE_THRESHOLD_GLOBAL,
89 GILOCALIZE_THRESHOLD_LOCAL,
90 GILOCALIZE_THRESHOLD_ROW
93 typedef enum GiThresholdMethod GiThresholdMethod;
130 _giraffe_validate_pixel(cxint *pixels, cxint xsize, cxint ysize,
131 cxint xpos, cxint ypos, cxint xwidth, cxint ywidth,
136 cxint xstart = xpos - xwidth;
137 cxint ystart = ypos - ywidth;
138 cxint xend = xpos + xwidth;
139 cxint yend = ypos + ywidth;
149 xstart = CX_MAX(0, xstart);
150 ystart = CX_MAX(0, ystart);
152 xend = CX_MIN(xsize - 1, xend);
153 yend = CX_MIN(ysize - 1, yend);
155 xwidth = CX_MAX(xwidth,1 );
156 ywidth = CX_MAX(ywidth,1 );
164 for (i = ystart; i <= yend; i++) {
183 for (j = xstart; j <= xend; j++) {
184 if (pixels[row + j]) {
188 if (_count >= count) {
227 inline static cpl_matrix*
228 _giraffe_fit_border(cpl_matrix* mborder, cpl_matrix* mbase,
229 cpl_matrix* mxok, cxint nspectra, cxdouble sigma,
230 cxint niter, cxdouble mfrac, cpl_matrix* mcoeff)
233 const cxchar*
const fctid =
"_giraffe_fit_border";
235 register cxint x = 0;
236 register cxint naccept = 0;
237 register cxint ntotal = 0;
238 register cxint iteration = 0;
239 register cxint nx = cpl_matrix_get_ncol(mbase);
240 register cxint yorder = cpl_matrix_get_nrow(mbase);
241 register cxint nxok = cpl_matrix_get_nrow(mxok);
243 register cxdouble ratio = 1.0;
245 cpl_matrix* mtmp = NULL;
246 cpl_matrix* yraw = NULL;
247 cpl_matrix* ydiff = NULL;
248 cpl_matrix* mfit = NULL;
249 cpl_matrix* coeffs = NULL;
254 cpl_error_set(fctid, CPL_ERROR_INCOMPATIBLE_INPUT);
256 GIDEBUG(
gi_warning(
"%s: not enough points mxok[%d] for %d order fit",
257 fctid, nxok, yorder));
267 yraw = cpl_matrix_new(1, nxok);
268 ydiff = cpl_matrix_new(nxok, 1);
270 mtmp = cpl_matrix_duplicate(mxok);
276 for (x = 0; x < nxok; x++) {
277 cxdouble data = cpl_matrix_get(mborder, x, nspectra);
278 cpl_matrix_set(yraw, 0, x, data);
289 while (naccept > 0 && iteration < niter && ratio > mfrac) {
291 register cxint k = 0;
292 register cxint l = 0;
294 register cxdouble ysigma = 0.;
296 cpl_matrix* rawbase = giraffe_chebyshev_base1d(0., nx, yorder, mtmp);
297 cx_assert(rawbase != NULL);
299 if (coeffs != NULL) {
300 cpl_matrix_delete(coeffs);
304 if (coeffs == NULL) {
305 gi_warning(
"%s: error in giraffe_matrix_leastsq(), spectrum %d",
310 cpl_matrix_delete(rawbase);
314 cpl_matrix_delete(mfit);
317 mfit = cpl_matrix_product_create(coeffs, mbase);
319 for (x = 0; x < cpl_matrix_get_nrow(ydiff); x++) {
321 cxint xok = (cxint) cpl_matrix_get(mtmp, x, 0);
324 cpl_matrix_get(yraw, 0, x) - cpl_matrix_get(mfit, 0, xok);
327 cpl_matrix_set(ydiff, x , 0, diff);
339 for (l = 0; l < cpl_matrix_get_nrow(ydiff); l++) {
341 if (fabs(cpl_matrix_get(ydiff, l, 0)) <= ysigma) {
343 cxint xok = cpl_matrix_get(mtmp, l, 0);
344 cxdouble data = cpl_matrix_get(yraw, 0, l);
346 cpl_matrix_set(mtmp, k, 0, xok);
347 cpl_matrix_set(yraw, 0, k, data);
369 ratio = (cxdouble) naccept / (cxdouble) ntotal;
371 GIDEBUG(
gi_message(
"Iteration %d: Sigma %f, accepted bins: %d, "
372 "rejected %d\n", iteration, ysigma, naccept,
379 cpl_matrix_resize(mtmp, 0,
380 naccept - cpl_matrix_get_nrow(mtmp), 0, 0);
381 cpl_matrix_resize(yraw, 0,
382 0, 0, naccept - cpl_matrix_get_ncol(yraw));
383 cpl_matrix_resize(ydiff, 0,
384 naccept - cpl_matrix_get_nrow(ydiff), 0, 0);
389 if (coeffs != NULL) {
392 for (l = 0; l < cpl_matrix_get_nrow(mcoeff); l++) {
393 cpl_matrix_set(mcoeff, l, 0, cpl_matrix_get(coeffs, 0, l));
402 cpl_matrix_delete(coeffs);
403 cpl_matrix_delete(ydiff);
404 cpl_matrix_delete(yraw);
405 cpl_matrix_delete(mtmp);
412 inline static cpl_image*
413 _giraffe_filter_gauss1d(
const cpl_image* image, cxint radius, cxdouble width)
416 cxdouble w2 = width * width;
420 cpl_matrix* kernel = cpl_matrix_new(1, 2 * radius + 1);
422 cpl_image* fimage = NULL;
425 if (kernel == NULL) {
429 for (i = -radius; i <= radius; ++i) {
431 cxdouble y = exp(-x2 / (2. * w2));
433 cpl_matrix_set(kernel, 0, i + radius, y);
437 fimage = cpl_image_new(cpl_image_get_size_x(image),
438 cpl_image_get_size_y(image),
439 cpl_image_get_type(image));
441 if (fimage == NULL) {
442 cpl_matrix_delete(kernel);
446 cpl_image_filter(fimage, image, kernel, CPL_FILTER_LINEAR,
448 cpl_matrix_delete(kernel);
455 inline static cpl_image*
456 _giraffe_filter_sobel(
const cpl_image* image, cxbool vertical)
458 cpl_matrix* kernel = cpl_matrix_new(3, 3);
460 cpl_image* fimage = NULL;
463 if (kernel == NULL) {
470 cpl_matrix_set(kernel, 0, 0, -1);
471 cpl_matrix_set(kernel, 1, 0, -2);
472 cpl_matrix_set(kernel, 2, 0, -1);
474 cpl_matrix_set(kernel, 0, 2, 1);
475 cpl_matrix_set(kernel, 1, 2, 2);
476 cpl_matrix_set(kernel, 2, 2, 1);
478 cpl_matrix_set(kernel, 0, 0, 0);
479 cpl_matrix_set(kernel, 1, 0, -0.5);
480 cpl_matrix_set(kernel, 2, 0, 0);
482 cpl_matrix_set(kernel, 0, 2, 0);
483 cpl_matrix_set(kernel, 1, 2, 0.5);
484 cpl_matrix_set(kernel, 2, 2, 0);
489 cpl_matrix_set(kernel, 0, 0, 1);
490 cpl_matrix_set(kernel, 0, 1, 2);
491 cpl_matrix_set(kernel, 0, 2, 1);
493 cpl_matrix_set(kernel, 2, 0, -1);
494 cpl_matrix_set(kernel, 2, 1, -2);
495 cpl_matrix_set(kernel, 2, 2, -1);
499 fimage = cpl_image_new(cpl_image_get_size_x(image),
500 cpl_image_get_size_y(image),
501 cpl_image_get_type(image));
503 if (fimage == NULL) {
504 cpl_matrix_delete(kernel);
508 cpl_image_filter(fimage, image, kernel, CPL_FILTER_LINEAR,
510 cpl_matrix_delete(kernel);
518 _giraffe_build_edge_mask(cpl_image* raw, cpl_image* bpixel, cxint nspectra,
519 cxdouble noise, GiMaskParameters* config,
520 cxint* ndetect, cpl_matrix* mxok, cpl_matrix* myup,
524 const cxint margin = 5;
532 cxint nrows = cpl_image_get_size_y(raw);
533 cxint ncols = cpl_image_get_size_x(raw);
537 cxdouble* buffer = NULL;
539 cpl_mask* kernel = NULL;
541 cpl_image* fraw = NULL;
542 cpl_image* sraw = NULL;
543 cpl_image* vertical1 = NULL;
544 cpl_image* vertical2 = NULL;
545 cpl_image* center = NULL;
558 kernel = cpl_mask_new(1, 15);
560 if (kernel != NULL) {
562 cpl_mask_not(kernel);
564 fraw = cpl_image_new(ncols, nrows, cpl_image_get_type(raw));
567 cpl_mask_delete(kernel);
573 cpl_image_filter_mask(fraw, raw, kernel, CPL_FILTER_MEDIAN,
578 cpl_mask_delete(kernel);
582 sraw = _giraffe_filter_gauss1d(fraw, 6, 1.);
586 cpl_image_delete(fraw);
593 vertical1 = _giraffe_filter_sobel(sraw, TRUE);
594 vertical2 = _giraffe_filter_sobel(vertical1, TRUE);
596 cpl_image_save(sraw,
"master_flat_smooth.fits", -32, 0, CPL_IO_DEFAULT);
597 cpl_image_save(vertical1,
"vertical.fits", -32, 0, CPL_IO_DEFAULT);
598 cpl_image_save(vertical2,
"vertical2.fits", -32, 0, CPL_IO_DEFAULT);
605 center = cpl_image_new(ncols, nrows, CPL_TYPE_INT);
607 flags = cx_calloc(ncols,
sizeof(cxint));
608 buffer = cx_calloc(ncols,
sizeof(cxdouble));
610 if ((center == NULL) || (flags ==NULL) || (buffer == NULL)) {
618 cpl_image_delete(center);
621 cpl_image_delete(vertical2);
624 cpl_image_delete(vertical1);
627 cpl_image_delete(sraw);
630 cpl_image_delete(fraw);
638 for (m = 0; m < nrows; ++m) {
640 register cxint irow = m * ncols;
641 register cxint n = 0;
646 cxint* _center = cpl_image_get_data_int(center) + irow;
648 const cxdouble* _vt1 = cpl_image_get_data_double_const(vertical1) +
650 const cxdouble* _vt2 = cpl_image_get_data_double_const(vertical2) +
652 const cxdouble* _fraw = cpl_image_get_data_double_const(fraw) +
656 memset(buffer, 0, ncols *
sizeof(cxdouble));
657 memset(flags, 0, ncols *
sizeof(cxint));
660 for (n = 0; n < ncols; ++n) {
670 if ((n - 1 >= 0) && (_vt2[n - 1] > 0.)) {
671 buffer[n - 1] = _vt1[n - 1];
673 if ((n + 1 < ncols) && (_vt2[n + 1] > 0.)) {
674 buffer[n + 1] = _vt1[n + 1];
680 while (iteration < ncols) {
684 cxdouble dx = 3. * 2. * noise;
687 for (n = 0; n < ncols; ++n) {
689 if (!flags[n] && (buffer[n] > dx)) {
699 register cxint k = 0;
706 cxdouble signal = 0.;
712 while ((k >= 0) && (buffer[k] > 0.)) {
719 while ((k < ncols) && (buffer[k] > 0.)) {
725 while ((k < ncols) && (buffer[k] < 0.)) {
730 width = end - start + 1;
738 signal = (_fraw[pos] > 0.) ? _fraw[pos] : 0.;
739 sigma = sqrt((noise * noise + signal) / config->xbin);
741 if ((signal / sigma > 10.) && (width > 1)) {
743 start = (start == pos) ? start - 1 : start;
744 end = (end == pos) ? end + 1 : end;
747 _center[start] += -1;
758 for (n = 0; n < ncols; ++n) {
760 if (_center[n] == 1) {
766 if (scount >= smax) {
782 cx_print(
"scount: %d (%d) at %d\n", smax, nspectra, mmax);
790 const cxint limit = 0.85 * nrows;
795 const cxdouble hwf = sqrt(2. * log(2.));
797 cxint* xtrace = cx_calloc(nrows,
sizeof(cxint));
798 cxint* ytrace = cx_calloc(nrows,
sizeof(cxint));
800 cpl_image* mask = cpl_image_new(ncols, nrows, CPL_TYPE_INT);
802 for (m = 0; m < ncols; ++m) {
804 const cxint* _center = cpl_image_get_data_int(center);
805 const cxint* _reference = _center + mmax * ncols;
807 cxbool out_of_bounds = FALSE;
812 if (_reference[m] == 1) {
814 register cxint j = mmax;
815 register cxint pos = m;
820 xtrace[connected] = pos;
821 ytrace[connected] = j;
827 register cxint k = 0;
828 register cxint l = j * ncols;
829 register cxint kmin = (pos - 1 >= 0) ? pos - 1 : 0;
830 register cxint kmax = (pos + 1 < ncols) ? pos + 1 : ncols - 1;
832 for (k = kmin; k <= kmax; ++k) {
834 if (_center[l + k] == 1) {
836 if ((pos <= margin) || (pos >= ncols - margin)) {
837 out_of_bounds = TRUE;
841 xtrace[connected] = k;
842 ytrace[connected] = j;
859 register cxint k = 0;
860 register cxint l = j * ncols;
861 register cxint kmin = (pos - 1 >= 0) ? pos - 1 : 0;
862 register cxint kmax = (pos + 1 < ncols) ? pos + 1 : ncols - 1;
864 for (k = kmin; k <= kmax; ++k) {
866 if (_center[l + k] == 1) {
868 if ((pos <= margin) || (pos >= ncols - margin)) {
869 out_of_bounds = TRUE;
873 xtrace[connected] = k;
874 ytrace[connected] = j;
886 if ((connected < limit) || (out_of_bounds == TRUE)) {
888 memset(xtrace, 0, nrows *
sizeof(cxint));
889 memset(ytrace, 0, nrows *
sizeof(cxint));
891 if (out_of_bounds == TRUE) {
892 cx_print(
"discarded candidate %d, going out of detector "
893 "boundaries.\n", itrace);
897 cx_print(
"discarded candidate %d, not enough connected "
898 "centers (%d, required: %d)\n", itrace, connected,
905 cxint* _mask = cpl_image_get_data_int(mask);
907 for (j = 0; j < connected; ++j) {
909 register cxint x = xtrace[j];
910 register cxint y = ytrace[j] * ncols;
911 register cxint ix = x;
915 while ((_center[y + ix] != -1) && (ix > 0)) {
921 while ((_center[y + ix] != -2) && (ix < ncols - 1)) {
936 cx_print(
"scount: %d (expected: %d)\n", ispectra, nspectra);
944 for (m = 0; m < nrows; ++m) {
946 register cxint j = 0;
947 register cxint ns = 0;
949 const cxint* _mask = cpl_image_get_data_int(mask) + m * ncols;
950 const cxint* _center = cpl_image_get_data_int(center) + m * ncols;
953 for (j = 0; j < ncols; ++j) {
957 register cxint x = j;
958 register cxint ix = x;
961 while ((_center[ix] != -1) && (ix > 0)) {
964 cpl_matrix_set(mylo, naccepted, ns, x - hwf * fabs(x - ix));
967 while ((_center[ix] != -2) && (ix < ncols - 1)) {
970 cpl_matrix_set(myup, naccepted, ns, x + hwf * fabs(ix - x));
977 if (ns == ispectra) {
978 cpl_matrix_set(mxok, naccepted, 0, m);
987 cpl_image_save(center,
"center.fits", -32, 0, CPL_IO_DEFAULT);
988 cpl_image_save(mask,
"mask.fits", -32, 0, CPL_IO_DEFAULT);
990 cpl_image_delete(mask);
991 cpl_image_delete(center);
992 cpl_image_delete(vertical2);
993 cpl_image_delete(vertical1);
994 cpl_image_delete(sraw);
995 cpl_image_delete(fraw);
1036 _giraffe_build_raw_mask(cpl_image *raw, cpl_image *bpixel, cxint nspectra,
1037 cxdouble noise, GiMaskParameters *config,
1038 cxint *ndetect, cpl_matrix *mxok, cpl_matrix *myup,
1042 register cxint x = 0;
1043 register cxint y = 0;
1044 register cxint xretry = 0;
1045 register cxint xok = 0;
1050 cxint *yabove = NULL;
1051 cxint *ybelow = NULL;
1052 cxint *good_pixels = NULL;
1053 cxint ywidth = config->ywidth > 1 ? config->ywidth : 2;
1054 cxint ckwidth = config->ckdata.width;
1055 cxint ckheight = config->ckdata.height;
1056 cxint ckcount = config->ckdata.count;
1059 cxdouble* pixels = NULL;
1061 cpl_mask* med = NULL;
1063 cpl_image* img = raw;
1068 med = cpl_mask_new(1, 15);
1074 img = cpl_image_new(cpl_image_get_size_x(raw),
1075 cpl_image_get_size_y(raw),
1076 cpl_image_get_type(raw));
1078 cpl_image_filter_mask(img, raw, med, CPL_FILTER_MEDIAN,
1083 cpl_mask_delete(med);
1088 GIDEBUG(
gi_message(
"noise = %g start = %d tries = %d xbin = %d "
1089 "ywidth = %d", noise, config->start, config->retry,
1090 config->xbin, ywidth));
1092 pixels = cpl_image_get_data_double(img);
1094 nrows = cpl_image_get_size_y(img);
1095 ncols = cpl_image_get_size_x(img);
1098 if (config->xbin > 1) {
1102 cxdouble* _pixels = NULL;
1105 nrows = (cxint) ceil(nrows / config->xbin);
1106 config->start = (cxint) ceil(config->start / config->xbin);
1108 _pixels = cx_calloc(ncols * nrows,
sizeof(cxdouble));
1110 for (y = 0; y < ncols; ++y) {
1112 for (x = 0; x < nrows; ++x) {
1114 register cxint xx = 0;
1115 register cxint zx = x * ncols;
1116 register cxint xr = x * config->xbin;
1117 register cxint zr = xr * ncols;
1120 _pixels[zx + y] = 0.;
1122 for (xx = 0; xx < config->xbin && xr < nx; ++xx) {
1123 _pixels[zx + y] += pixels[zr + y];
1126 _pixels[zx + y] /= config->xbin;
1136 good_pixels = cx_calloc(nrows * ncols,
sizeof(cxint));
1138 switch (config->method) {
1140 case GILOCALIZE_THRESHOLD_LOCAL:
1143 cxint ywidth2 = ywidth / 2;
1144 cxint sz = 2 * ywidth2 + 1;
1146 cpl_vector* ymins = cpl_vector_new(sz);
1155 for (x = 0; x < nrows; x++) {
1157 cpl_vector_fill(ymins, 0.);
1159 for (y = 0; y < ncols; y++) {
1161 register cxint k = 0;
1162 register cxint kk = 0;
1164 cxdouble value = 0.;
1166 cxdouble threshold = 0.;
1169 for (kk = 0, k = -ywidth2; k <= ywidth2; k++) {
1171 register cxint ky = y + k;
1173 if (ky < 0 || ky >= ncols) {
1177 cpl_vector_set(ymins, kk, pixels[x * ncols + ky]);
1185 if (config->threshold > 0.) {
1187 const cxint count = 2;
1202 for (i = 0; i < count; i++) {
1203 bkg += fabs(cpl_vector_get(ymins, i));
1205 bkg /= (cxdouble)count;
1207 threshold = sqrt((2. * noise * noise +
1208 fabs(pixels[x * ncols + y]) + bkg / count) / config->xbin);
1214 register cxdouble mean = 0.;
1217 for (i = 0; i < kk; i++) {
1218 mean += cpl_vector_get(ymins, i);
1224 bkg = (cpl_vector_get(ymins, 0) +
1225 cpl_vector_get(ymins, 1)) / 2.0;
1226 threshold = mean - bkg;
1235 value = pixels[x * ncols + y] - bkg;
1241 if (value > fabs(config->threshold) * threshold) {
1242 good_pixels[x * ncols + y] = 1;
1247 cpl_vector_delete(ymins);
1254 case GILOCALIZE_THRESHOLD_ROW:
1257 cpl_image* snr = cpl_image_abs_create(raw);
1259 cxint sx = cpl_image_get_size_x(snr);
1262 cpl_image_power(snr, 0.5);
1264 for (x = 0; x < nrows; ++x) {
1266 const cxdouble* _snr = cpl_image_get_data_double_const(snr);
1268 cxdouble avsnr = giraffe_array_median(_snr + x * sx, sx);
1271 for (y = 0; y < ncols; ++y) {
1273 if (pixels[x * ncols + y] <= 0.) {
1277 if (_snr[x * ncols + y] > avsnr * fabs(config->threshold)) {
1278 good_pixels[x * ncols + y] = 1;
1285 cpl_image_delete(snr);
1295 cxdouble threshold = 0.;
1302 if (config->threshold > 0.) {
1303 threshold = config->threshold * noise;
1307 cxdouble mean = cpl_image_get_mean(raw);
1309 threshold = -config->threshold * mean *
1310 (nspectra * config->wavg / ncols);
1314 for (x = 0; x < nrows; x++) {
1316 for (y = 0; y < ncols; y++) {
1318 if (pixels[x * ncols + y] > threshold) {
1319 good_pixels[x * ncols + y] = 1;
1332 GIDEBUG(cxint *data = cx_calloc(nrows * ncols,
sizeof(cxint));
1333 memcpy(data, good_pixels, nrows * ncols *
sizeof(cxint));
1334 cpl_image *gp = cpl_image_wrap_int(ncols, nrows, data);
1335 cpl_image_save(gp,
"locmask.fits", 32, NULL, CPL_IO_DEFAULT);
1336 cpl_image_unwrap(gp);
1344 yabove = cx_calloc(nspectra + 1,
sizeof(cxint));
1345 ybelow = cx_calloc(nspectra + 1,
sizeof(cxint));
1357 for (x = config->start; (x >= 0) && (xretry <= config->retry); x--) {
1359 register cxint zx = x * ncols;
1360 register cxint nborders = 0;
1361 register cxint nbelow = 0;
1362 register cxint nabove = 0;
1363 register cxint in_spectrum = 0;
1366 for (y = 1; y < ny; y++) {
1368 register cxint tmp = 2 * good_pixels[zx + y];
1374 nborders = CX_MAX(CX_MAX(nborders, nbelow), nabove);
1376 if (nborders > nspectra) {
1386 if (good_pixels[zx + y + 1]) {
1393 if ((tmp - good_pixels[zx + y - 1]) == 2) {
1407 ybelow[nbelow++] = y;
1416 if (good_pixels[zx + y - 1]) {
1423 if ((tmp - good_pixels[zx + y + 1]) == 2) {
1437 yabove[nabove++] = y;
1449 !good_pixels[zx + y - 1] && !good_pixels[zx + y + 1]) {
1451 if (_giraffe_validate_pixel(good_pixels, ncols, nrows, y, x,
1452 ckwidth, ckheight, ckcount)) {
1454 yabove[nabove++] = y;
1455 ybelow[nbelow++] = y;
1468 *ndetect = nborders;
1470 if (!in_spectrum && (nbelow == nspectra) && (nbelow == nabove)) {
1479 for (y = 0; y < nspectra; y++) {
1480 cpl_matrix_set(mylo, xok, y, (cxdouble) ybelow[y]);
1481 cpl_matrix_set(myup, xok, y, (cxdouble) yabove[y]);
1482 cpl_matrix_set(mxok, xok, 0, (config->xbin > 1) ?
1483 (cxdouble) (x + 0.5) * config->xbin :
1489 else if (xretry++ < config->retry) {
1519 for (x = config->start + 1; (x < nrows) &&
1520 (xretry <= config->retry); x++) {
1522 register cxint zx = x * ncols;
1523 register cxint nborders = 0;
1524 register cxint nbelow = 0;
1525 register cxint nabove = 0;
1526 register cxint in_spectrum = 0;
1529 for (y = 1; y < ny; y++) {
1531 register cxint tmp = 2 * good_pixels[zx + y];
1533 nborders = CX_MAX(CX_MAX(nborders, nbelow), nabove);
1535 if (nborders > nspectra) {
1539 if (good_pixels[zx + y + 1]) {
1540 if ((tmp - good_pixels[zx + y - 1]) == 2) {
1542 ybelow[nbelow++] = y;
1548 if (good_pixels[zx + y - 1]) {
1549 if ((tmp - good_pixels[zx + y + 1]) == 2) {
1551 yabove[nabove++] = y;
1560 !good_pixels[zx + y - 1] && !good_pixels[zx + y + 1]) {
1562 if (_giraffe_validate_pixel(good_pixels, ncols, nrows, y, x,
1563 ckwidth, ckheight, ckcount)) {
1565 yabove[nabove++] = y;
1566 ybelow[nbelow++] = y;
1579 *ndetect = nborders;
1581 if (!in_spectrum && (nbelow == nspectra) && (nbelow == nabove)) {
1583 for (y = 0; y < nspectra; y++) {
1584 cpl_matrix_set(mylo, xok, y, (cxdouble) ybelow[y]);
1585 cpl_matrix_set(myup, xok, y, (cxdouble) yabove[y]);
1586 cpl_matrix_set(mxok, xok, 0, (config->xbin > 1) ?
1587 (cxdouble) (x + 0.5) * config->xbin :
1593 else if (xretry++ < config->retry) {
1604 cx_free(good_pixels);
1606 if (pixels != cpl_image_get_data_double(img)) {
1612 cpl_image_delete(img);
1617 if (*ndetect < nspectra) {
1620 else if (*ndetect > nspectra) {
1628 *ndetect = nspectra;
1662 _giraffe_fit_raw_mask(cpl_matrix *mxok, cpl_matrix *myup, cpl_matrix *mylo,
1663 cpl_table *fibers, GiMaskParameters *config,
1664 GiMaskPosition *position)
1667 register cxint nn, x, nspectra;
1668 register cxint nx = cpl_matrix_get_nrow(position->my);
1669 register cxint ns = cpl_table_get_nrow(fibers);
1677 mxraw = cpl_matrix_new(nx, 1);
1678 mcoeff = cpl_matrix_new(config->ydeg + 1, 1);
1685 for (x = 0; x < nx; x++) {
1686 cpl_matrix_set(mxraw, x, 0, x);
1693 base = giraffe_chebyshev_base1d(0., nx, config->ydeg + 1, mxraw);
1694 cpl_matrix_delete(mxraw);
1697 for (nn = 0; nn < ns; nn++) {
1698 cpl_matrix *ylofit = NULL;
1699 cpl_matrix *yupfit = NULL;
1712 ylofit = _giraffe_fit_border(mylo, base, mxok, nspectra,
1713 config->sigma, config->niter,
1714 config->mfrac, mcoeff);
1715 if (ylofit == NULL) {
1716 cpl_msg_warning(_task,
"Could not compute low border for "
1723 yupfit = _giraffe_fit_border(myup, base, mxok, nspectra,
1724 config->sigma, config->niter,
1725 config->mfrac, mcoeff);
1726 if (yupfit == NULL) {
1727 cpl_msg_warning(_task,
"Could not compute up border for "
1739 for (x = 0; x < nx; x++) {
1741 cpl_matrix_set(position->my, x, nn, 0.5 *
1742 (cpl_matrix_get(yupfit, x, 0) +
1743 cpl_matrix_get(ylofit, x, 0)));
1745 cpl_matrix_set(position->my, x, nn, 0.5 *
1746 (cpl_matrix_get(yupfit, x, 0) -
1747 cpl_matrix_get(ylofit, x, 0)) + config->ewid);
1750 cpl_matrix_delete(ylofit);
1751 cpl_matrix_delete(yupfit);
1756 cpl_msg_info(_task,
"%03d spectrum positions fitted", nspectra);
1758 cpl_matrix_delete(base);
1759 cpl_matrix_delete(mcoeff);
1761 if (nspectra == 0) {
1762 cpl_msg_warning(_task,
"could not fit any spectra, check number "
1763 "of good wavelength bins");
1807 _giraffe_fit_raw_centroid(cpl_image* mz, cpl_matrix* mxok, cpl_matrix* myup,
1808 cpl_matrix* mylo, cpl_table* fibers,
1809 GiMaskParameters* config, GiMaskPosition* position,
1810 GiMaskPosition* coeffs)
1813 const cxchar*
const fctid =
"_giraffe_fit_raw_centroid";
1815 register cxint nn = 0;
1816 register cxint x = 0;
1817 register cxint y = 0;
1818 register cxint nspectra = 0;
1819 register cxint nx = cpl_image_get_size_y(mz);
1820 register cxint ny = cpl_image_get_size_x(mz);
1821 register cxint ns = cpl_table_get_nrow(fibers);
1823 cxint yorder = config->ydeg + 1;
1824 cxint worder = config->wdeg + 1;
1826 cpl_matrix* mxraw = NULL;
1827 cpl_matrix* base = NULL;
1828 cpl_matrix* mycenter = NULL;
1829 cpl_matrix* mywidth = NULL;
1830 cpl_matrix* mx = NULL;
1831 cpl_matrix* my = NULL;
1832 cpl_matrix* mw = NULL;
1833 cpl_matrix* chebcoeff = NULL;
1834 cpl_matrix* mfitlocw = NULL;
1835 cpl_matrix* ycenfit = NULL;
1836 cpl_matrix* ycencoeff = NULL;
1840 if (cpl_matrix_get_nrow(position->my) != nx ||
1841 cpl_matrix_get_ncol(position->my) != ns) {
1842 gi_error(
"%s: invalid size for position->my[%" CPL_SIZE_FORMAT
",%"
1843 CPL_SIZE_FORMAT
"], expected [%d,%d]", fctid,
1844 cpl_matrix_get_nrow(position->my),
1845 cpl_matrix_get_ncol(position->my), nx, ns);
1849 if (cpl_matrix_get_nrow(position->mw) != nx ||
1850 cpl_matrix_get_ncol(position->mw) != ns) {
1851 gi_error(
"%s: invalid size for position->mw[%" CPL_SIZE_FORMAT
",%"
1852 CPL_SIZE_FORMAT
"], expected [%d,%d]", fctid,
1853 cpl_matrix_get_nrow(position->my),
1854 cpl_matrix_get_ncol(position->my), nx, ns);
1863 mxraw = cpl_matrix_new(nx, 1);
1865 for (x = 0; x < nx; x++) {
1866 cpl_matrix_set(mxraw, x, 0, x);
1874 base = giraffe_chebyshev_base1d(0., nx, yorder, mxraw);
1875 cpl_matrix_delete(mxraw);
1877 mycenter = cpl_matrix_new(cpl_matrix_get_nrow(mxok), ns);
1878 mywidth = cpl_matrix_new(1, cpl_matrix_get_nrow(mxok) * ns);
1880 ycencoeff = cpl_matrix_new(yorder, 1);
1882 for (nn = 0; nn < ns; nn++) {
1899 cxdouble* pixels = cpl_image_get_data_double(mz);
1901 for (x = 0; x < cpl_matrix_get_nrow(mxok); x++) {
1903 register cxint zx = (cxint) cpl_matrix_get(mxok, x, 0);
1905 register cxdouble zz = 0.;
1906 register cxdouble yy = 0.;
1908 cxdouble lower = cpl_matrix_get(mylo, x, nspectra);
1909 cxdouble upper = cpl_matrix_get(myup, x, nspectra);
1912 for (y = (cxint) lower; y <= (cxint) upper; y++) {
1913 yy += pixels[zx * ny + y] * y;
1914 zz += pixels[zx * ny + y];
1917 cpl_matrix_set(mycenter, x, nspectra, yy / zz);
1918 cpl_matrix_set(mywidth, 0, x * ns + nspectra, config->ewid +
1919 (upper - lower) / 2.0);
1927 cpl_matrix_fill(ycencoeff, 0.);
1928 ycenfit = _giraffe_fit_border(mycenter, base, mxok, nspectra,
1929 config->sigma, config->niter,
1930 config->mfrac, ycencoeff);
1931 if (ycenfit == NULL) {
1932 cpl_msg_warning(_task,
"Could not fit centroid for spectrum %d",
1942 for (x = 0; x < yorder; x++) {
1943 cpl_matrix_set(coeffs->my, x, nn,
1944 cpl_matrix_get(ycencoeff, x, 0));
1952 for (x = 0; x < nx; x++) {
1953 cpl_matrix_set(position->my, x, nn,
1954 cpl_matrix_get(ycenfit, 0, x));
1957 cpl_matrix_delete(ycenfit);
1963 cpl_image_save(lycenter,
"lycenter.fits", -32, NULL,
1965 cpl_image_delete(lycenter);
1968 cpl_image_save(lycenter,
"lycenterfit.fits", -32, NULL,
1970 cpl_image_delete(lycenter);
1973 cpl_image_save(lyxok,
"lyxok.fits", -32, NULL,
1975 cpl_image_delete(lyxok));
1978 cpl_msg_info(_task,
"%03d spectrum positions fitted", nspectra);
1980 cpl_matrix_delete(base);
1981 cpl_matrix_delete(mycenter);
1982 cpl_matrix_delete(ycencoeff);
1984 if (nspectra == 0) {
1985 cpl_msg_warning(_task,
"Could not fit any spectra, check number of "
1986 "good wavelength bins");
1988 cpl_matrix_delete(mywidth);
1996 cpl_msg_info(_task,
"2D fit (order %dx%d) of mask width", worder,
2003 mx = cpl_matrix_new(cpl_matrix_get_nrow(mxok) * nspectra, 1);
2004 my = cpl_matrix_new(cpl_matrix_get_nrow(mxok) * nspectra, 1);
2005 mw = cpl_matrix_new(1, cpl_matrix_get_nrow(mxok) * nspectra);
2007 for (y = 0, nn = 0; nn < nspectra; nn++) {
2019 for (x = 0; x < cpl_matrix_get_nrow(mxok); x++) {
2021 register cxint zx = (cxint) cpl_matrix_get(mxok, x, 0);
2022 register cxint lx = x * nspectra + y;
2025 cpl_matrix_set(mx, lx, 0, cpl_matrix_get(mxok, x, 0));
2026 cpl_matrix_set(my, lx, 0, cpl_matrix_get(position->my, zx, nn));
2027 cpl_matrix_set(mw, 0, lx, cpl_matrix_get(mywidth, 0, x * ns + y));
2032 base = giraffe_chebyshev_base2d(0., 0., nx, ny, worder, worder, mx, my);
2034 cpl_matrix_delete(my);
2035 cpl_matrix_delete(mx);
2038 cpl_matrix_delete(base);
2039 cpl_matrix_delete(mw);
2041 cpl_matrix_delete(mywidth);
2043 if (chebcoeff == NULL) {
2044 gi_warning(
"%s: error in giraffe_matrix_leastsq() for width 2D fit",
2053 for (nn = 0; nn < cpl_matrix_get_ncol(chebcoeff); nn++) {
2054 cpl_matrix_set(coeffs->mw, 0, nn, cpl_matrix_get(chebcoeff, 0, nn));
2061 mx = cpl_matrix_new(nx * nspectra, 1);
2062 my = cpl_matrix_new(nx * nspectra, 1);
2064 for (y = 0, nn = 0; nn < nspectra; nn++) {
2076 for (x = 0; x < nx; x++) {
2078 register cxint lx = x * nspectra + y;
2080 cpl_matrix_set(mx, lx, 0, x);
2081 cpl_matrix_set(my, lx, 0, cpl_matrix_get(position->my, x, nn));
2087 cpl_matrix_set_size(chebcoeff, worder, worder);
2089 mfitlocw = giraffe_chebyshev_fit2d(0., 0., nx, ny, chebcoeff, mx, my);
2090 cpl_matrix_delete(chebcoeff);
2092 cpl_matrix_delete(my);
2093 cpl_matrix_delete(mx);
2095 for (y = 0, nn = 0; nn < nspectra; nn++) {
2107 for (x = 0; x < nx; x++) {
2109 register cxint lx = x * nspectra + y;
2111 cpl_matrix_set(position->mw, x, nn,
2112 cpl_matrix_get(mfitlocw, lx, 0));
2118 cpl_matrix_delete(mfitlocw);
2151 _giraffe_localize_spectra(cpl_image *mzraw, cpl_image *bpixel,
2152 cpl_table *fibers, GiLocalizeMethod method,
2153 cxbool normalize, cxdouble noise,
2154 GiMaskParameters *config, GiMaskPosition *position,
2155 GiMaskPosition *coeffs)
2160 cxint ndetect, nspectra;
2163 cxdouble uplost = 0.;
2164 cxdouble lolost = 0.;
2165 cxdouble avglost = 0.;
2166 cxdouble avgmask = 0.;
2167 cxdouble sigmask = 0.;
2168 cxdouble sigmean = 0.;
2169 cxdouble avgborders = 0.;
2178 cpl_image *mz = NULL;
2179 cpl_image *mznorm = NULL;
2183 nx = cpl_image_get_size_y(mzraw);
2184 ny = cpl_image_get_size_x(mzraw);
2185 _mzraw = cpl_image_get_data_double(mzraw);
2188 if (normalize == TRUE) {
2190 cxdouble zxmax = 0.0;
2191 cxdouble *_mzx = NULL;
2192 cxdouble *_mznorm = NULL;
2194 cpl_image *mzx = NULL;
2197 cpl_msg_info(_task,
"Using normalized spectra for localization");
2206 mznorm = cpl_image_new(ny, nx, CPL_TYPE_DOUBLE);
2207 _mznorm = cpl_image_get_data_double(mznorm);
2209 mzx = cpl_image_new(1, nx, CPL_TYPE_DOUBLE);
2210 _mzx = cpl_image_get_data_double(mzx);
2217 for (x = 0 ; x < nx; x++) {
2218 for (y = 0 ; y < ny; y++) {
2219 _mzx[x] += _mzraw[x * ny + y];
2226 if (_mzx[x] > zxmax) {
2231 GIDEBUG(cpl_image_save(mzx,
"mzx.fits", -32, NULL, CPL_IO_DEFAULT));
2233 for (x = 0 ; x < nx; x++) {
2235 register cxdouble zxnorm = zxmax / _mzx[x];
2237 for (y = 0 ; y < ny; y++) {
2238 _mznorm[x * ny + y] = _mzraw[x * ny + y] * zxnorm;
2243 cpl_image_delete(mzx);
2252 cpl_msg_info(_task,
"Using raw spectra for localization");
2261 nspectra = cpl_table_get_nrow(fibers);
2263 mxok = cpl_matrix_new(nx, 1);
2264 myup = cpl_matrix_new(nx, nspectra);
2265 mylo = cpl_matrix_new(nx, nspectra);
2272 config->xbin = (config->xbin > 1) ? 2 * (config->xbin / 2) : 1;
2274 GIDEBUG(cpl_image_save(mz,
"mz.fits", -32, NULL, CPL_IO_DEFAULT));
2281 cpl_msg_info(_task,
"Generating mask (%d spectra expected) ...",
2287 nxok = _giraffe_build_edge_mask(mz, bpixel, nspectra, noise, config,
2288 &ndetect, mxok, myup, mylo);
2293 nxok = _giraffe_build_raw_mask(mz, bpixel, nspectra, noise, config,
2294 &ndetect, mxok, myup, mylo);
2298 cpl_matrix_delete(mxok);
2299 cpl_matrix_delete(myup);
2300 cpl_matrix_delete(mylo);
2304 cpl_msg_warning(_task,
"Invalid number of spectra detected: "
2305 "%d != %d", ndetect, nspectra);
2309 cpl_msg_warning(_task,
"No abcissa with good number "
2314 cpl_msg_warning(_task,
"Error while searching for spectra");
2322 cpl_msg_info(_task,
"%d spectra detected in %d wavelength bins",
2331 cpl_matrix_resize(mxok, 0, nxok - cpl_matrix_get_nrow(mxok), 0, 0);
2332 cpl_matrix_resize(myup, 0, nxok - cpl_matrix_get_nrow(myup), 0, 0);
2333 cpl_matrix_resize(mylo, 0, nxok - cpl_matrix_get_nrow(mylo), 0, 0);
2335 GIDEBUG(
gi_message(
"%s: mxok[0-%d]=[%g-%g]", __func__,
2336 cpl_matrix_get_nrow(mxok) - 1,
2337 cpl_matrix_get_min(mxok),
2338 cpl_matrix_get_max(mxok)));
2341 cpl_msg_info(_task,
"Computing spectrum positions and widths in "
2342 "pixel range [%g,%g]", cpl_matrix_get_min(mxok),
2343 cpl_matrix_get_max(mxok));
2345 if (cpl_matrix_get_nrow(mxok) <= config->ydeg) {
2346 cpl_msg_info(_task,
"Not enough data points %" CPL_SIZE_FORMAT
2347 " for %d order fit", cpl_matrix_get_nrow(mxok),
2354 case GILOCALIZE_HALF_WIDTH:
2355 cpl_msg_info(_task,
"Using half-width for localization");
2356 _giraffe_fit_raw_mask(mxok, myup, mylo, fibers, config,
2360 case GILOCALIZE_BARYCENTER:
2362 cpl_msg_info(_task,
"Using barycenter for localization");
2363 _giraffe_fit_raw_centroid(mz, mxok, myup, mylo, fibers, config,
2368 if (normalize == 1) {
2369 cpl_image_delete(mznorm);
2382 mwid = cpl_matrix_new(nxok, nspectra);
2384 for (n = 0, nn = 0; nn < cpl_table_get_nrow(fibers); nn++) {
2386 for (x = 0; x < nxok; x++) {
2387 register cxint lx = (cxint) cpl_matrix_get(mxok, x, 0);
2389 cxdouble lower = cpl_matrix_get(mylo, x, n);
2390 cxdouble upper = cpl_matrix_get(myup, x, n);
2391 cxdouble width = cpl_matrix_get(position->mw, lx, nn);
2393 uplost += cpl_matrix_get(position->my, lx, nn) + width - upper;
2394 lolost += cpl_matrix_get(position->my, lx, nn) - width - lower;
2396 avgborders += upper - lower;
2399 cpl_matrix_set(mwid, x, n, 2. * width);
2404 sigmean = cpl_matrix_get_mean(mwid);
2406 avglost = (lolost + uplost) / (nspectra * nxok);
2407 avgmask = 2.0 * avgmask / nspectra;
2409 cpl_msg_info(_task,
"Mask was computed using %d of %d wavelength bins",
2411 cpl_msg_info(_task,
"Average # of pixels per spectra: %.4g",
2413 cpl_msg_info(_task,
"Average # of in-borders pixels per spectra: %.4g",
2414 avgborders / nspectra);
2415 cpl_msg_info(_task,
"Average lost pixels per spectra: %.4g",
2417 cpl_msg_info(_task,
"Average lost pixels at upper border: %.4g",
2418 uplost / (nspectra * nxok));
2419 cpl_msg_info(_task,
"Average lost pixels at lower border: %.4g",
2420 lolost / (nspectra * nxok));
2421 cpl_msg_info(_task,
"Average spectrum width: %.4g +/- %.4g, "
2422 "(min, max) = (%.4g, %.4g)", sigmean, sigmask,
2423 cpl_matrix_get_min(mwid), cpl_matrix_get_max(mwid));
2425 cpl_matrix_delete(mwid);
2427 cpl_matrix_delete(mylo);
2428 cpl_matrix_delete(myup);
2429 cpl_matrix_delete(mxok);
2437 _giraffe_finalize_fibers(cpl_table *fibers, cpl_matrix *locy, GiImage *mlocy,
2438 cxdouble maxoffset, cxdouble* maxshift)
2450 cxdouble max_shift = 0.;
2451 cxdouble *positions = NULL;
2453 cpl_image *_mlocy = NULL;
2456 if (fibers == NULL || locy == NULL || mlocy == NULL) {
2460 if (cpl_table_has_column(fibers,
"RINDEX") == FALSE) {
2464 nx = cpl_matrix_get_ncol(locy);
2465 ny = cpl_matrix_get_nrow(locy);
2467 nfibers = cpl_table_get_nrow(fibers);
2470 _nx = cpl_image_get_size_x(_mlocy);
2471 _ny = cpl_image_get_size_y(_mlocy);
2477 if (nfibers > _nx) {
2481 cpl_table_select_all(fibers);
2488 irow = (_ny - 1) / 2;
2489 positions = (cxdouble *)cpl_image_get_data(_mlocy) + irow * _nx;
2500 for (i = 0; i < nfibers; i++) {
2504 cxint pos = cpl_table_get_int(fibers,
"RINDEX", i, NULL) - 1;
2506 cxdouble yc = cpl_matrix_get(locy, irow, j);
2507 cxdouble shift = fabs(yc - positions[pos]);
2509 if (shift <= maxoffset) {
2510 cpl_table_unselect_row(fibers, i);
2514 max_shift = CX_MAX(max_shift, shift);
2520 cpl_table_erase_selected(fibers);
2522 if (maxshift != NULL) {
2523 *maxshift = max_shift;
2561 GiTable *fibers, GiLocalization *master,
2562 GiImage *badpixels, GiLocalizeConfig *config)
2565 const cxchar *fctid =
"giraffe_localize_spectra";
2577 cxdouble conad = 0.;
2578 cxdouble bias_ron = 0.;
2579 cxdouble mask_sigma = 0.;
2583 cpl_propertylist *properties;
2587 cpl_image *_result = NULL;
2591 cpl_table *_fibers = NULL;
2592 cpl_table *fiber_setup = NULL;
2595 GiLocalizeMethod method;
2597 GiInstrumentMode mode;
2599 GiMaskParameters mask_config;
2601 GiMaskPosition mask_position;
2602 GiMaskPosition mask_coeffs;
2610 if (result == NULL || image == NULL || fibers == NULL || config == NULL) {
2611 cpl_error_set(fctid, CPL_ERROR_NULL_INPUT);
2615 if (badpixels != NULL) {
2616 cpl_msg_debug(fctid,
"Bad pixel correction is not available. Bad "
2617 "pixel map will be ignored.");
2622 if (_fibers == NULL) {
2623 cpl_error_set(fctid, CPL_ERROR_DATA_NOT_FOUND);
2627 fiber_setup = _fibers;
2637 nfibers = cpl_table_get_nrow(_fibers);
2639 cpl_msg_info(fctid,
"Setting number of fibers (%s) to %d",
2640 GIALIAS_NFIBERS, nfibers);
2642 cpl_propertylist_update_int(properties, GIALIAS_NFIBERS, nfibers);
2643 cpl_propertylist_set_comment(properties, GIALIAS_NFIBERS,
2644 "Number of fibres");
2647 giraffe_error_push();
2651 if (cpl_error_get_code() != CPL_ERROR_NONE) {
2655 giraffe_error_pop();
2662 if (config->ron > 0.) {
2663 cpl_msg_info(fctid,
"Setting bias sigma value (%s) to %.5g",
2664 GIALIAS_BIASSIGMA, config->ron);
2665 cpl_propertylist_update_double(properties, GIALIAS_BIASSIGMA,
2670 cpl_msg_info(fctid,
"Bias sigma value: %.3g e-", bias_ron);
2673 if (cpl_propertylist_has(properties, GIALIAS_DATANCOM)) {
2674 nframes = cpl_propertylist_get_int(properties, GIALIAS_DATANCOM);
2678 if (config->noise > 0.) {
2679 cpl_msg_info(fctid,
"Noise multiplier: %.3g",
2683 cpl_msg_info(fctid,
"Threshold multiplier: %.3g",
2684 fabs(config->noise));
2692 nrows = cpl_image_get_size_y(_image);
2694 if (config->start < 0) {
2695 config->start = nrows / 2;
2705 if (config->ywidth < 1) {
2707 cpl_msg_info(fctid,
"Configuring equilizing filter width from "
2712 config->ywidth = 16;
2724 cpl_msg_error(fctid,
"Invalid instrument mode!");
2730 if (!cpl_propertylist_has(properties, GIALIAS_SLITNAME)) {
2731 cpl_msg_error(fctid,
"Property (%s) not found in raw image",
2736 const cxchar *slit =
2737 cpl_propertylist_get_string(properties, GIALIAS_SLITNAME);
2739 cpl_msg_info(fctid,
"Setting equilizing filter to %d [pxl] "
2740 "for slit configuration `%s'", config->ywidth,
2753 mwidth = GISPECTRUM_MWIDTH_MEDUSA;
2762 mwidth = GISPECTRUM_MWIDTH_IFU;
2771 mwidth = GISPECTRUM_MWIDTH_IFU;
2780 cpl_msg_error(fctid,
"Invalid instrument mode!");
2790 if (config->centroid == TRUE) {
2791 method = GILOCALIZE_BARYCENTER;
2794 method = GILOCALIZE_HALF_WIDTH;
2803 mask_config.ywidth = config->ywidth;
2804 mask_config.method = config->threshold;
2805 mask_config.threshold = config->noise;
2806 mask_config.ydeg = config->yorder;
2807 mask_config.wdeg = config->worder;
2808 mask_config.ewid = config->ewidth;
2809 mask_config.wavg = mwidth;
2810 mask_config.ckdata.width = ckwidth;
2811 mask_config.ckdata.height = ckheight;
2812 mask_config.ckdata.count = ckcount;
2813 mask_config.sigma = config->sigma;
2814 mask_config.niter = config->iterations;
2815 mask_config.mfrac = config->fraction;
2816 mask_config.start = config->start;
2817 mask_config.retry = config->retries;
2818 mask_config.xbin = config->binsize;
2828 if (config->noise > 0.) {
2829 mask_config.threshold *= sqrt(nframes * conad);
2843 if (config->full != TRUE) {
2845 cpl_msg_info(fctid,
"Computing spectrum localization using SIWC "
2848 if (!master || !master->locy || !master->locy) {
2849 cpl_msg_error(fctid,
"Required full master localization is "
2860 cpl_table_unselect_all(_fibers);
2861 cpl_table_or_selected_int(_fibers,
"RP", CPL_EQUAL_TO, -1);
2863 fiber_setup = cpl_table_extract_selected(_fibers);
2864 nfibers = cpl_table_get_nrow(fiber_setup);
2874 mask_position.type = GIMASK_FITTED_DATA;
2875 mask_position.my = cpl_matrix_new(nrows, nfibers);
2876 mask_position.mw = cpl_matrix_new(nrows, nfibers);
2878 mask_coeffs.type = GIMASK_FIT_COEFFS;
2879 mask_coeffs.my = cpl_matrix_new(mask_config.ydeg + 1, nfibers);
2880 mask_coeffs.mw = cpl_matrix_new(1, (mask_config.wdeg + 1) *
2881 (mask_config.wdeg + 1));
2889 _image = cpl_image_multiply_scalar_create(_image, nframes * conad);
2891 mask_sigma = sqrt(nframes) * bias_ron;
2898 status = _giraffe_localize_spectra(_image, _bpixel, fiber_setup,
2899 method, config->normalize,
2901 &mask_config, &mask_position,
2904 cpl_image_delete(_image);
2908 result->locy = NULL;
2909 result->locw = NULL;
2910 result->locc = NULL;
2913 cpl_matrix_delete(mask_position.my);
2914 cpl_matrix_delete(mask_position.mw);
2916 cpl_matrix_delete(mask_coeffs.my);
2917 cpl_matrix_delete(mask_coeffs.mw);
2919 if (config->full != TRUE) {
2920 cpl_table_delete(fiber_setup);
2923 cpl_msg_error(fctid,
"Spectrum localization computation failed!");
2933 if (config->full != TRUE) {
2940 cpl_table_delete(fiber_setup);
2945 if (master != NULL && master->locy != NULL) {
2947 cxint nf = cpl_table_get_nrow(_fibers);
2949 cxdouble maxoffset = 0.5 * mask_config.wavg;
2950 cxdouble maxshift = 0.;
2953 cpl_msg_info(fctid,
"Comparing detected and expected fiber "
2956 status = _giraffe_finalize_fibers(_fibers, mask_position.my,
2957 master->locy, maxoffset,
2965 cxint _nf = cpl_image_get_size_x(mlocy);
2967 cpl_msg_error(fctid,
"More fibers (%d) than expected "
2968 "(%d) were found!", nf, _nf);
2972 result->locy = NULL;
2973 result->locw = NULL;
2974 result->locc = NULL;
2977 cpl_matrix_delete(mask_position.my);
2978 cpl_matrix_delete(mask_position.mw);
2980 cpl_matrix_delete(mask_coeffs.my);
2981 cpl_matrix_delete(mask_coeffs.mw);
2983 if (config->full != TRUE) {
2984 cpl_table_delete(fiber_setup);
2987 cpl_msg_error(fctid,
"Comparison of fiber positions "
2993 cx_assert(cpl_table_get_nrow(_fibers) <= nf);
2995 cpl_msg_info(fctid,
"%" CPL_SIZE_FORMAT
" of %d expected fibers "
2996 "were detected.", cpl_table_get_nrow(_fibers), nf);
2998 if (cpl_table_get_nrow(_fibers) < nf) {
2999 cpl_msg_debug(fctid,
"Maximum offset from the expected "
3000 "position is %.2f, maximum allowed offset is %.2f",
3001 maxshift, maxoffset);
3002 cpl_msg_warning(fctid,
"%" CPL_SIZE_FORMAT
" fibers are "
3003 "missing!", nf - cpl_table_get_nrow(_fibers));
3020 cpl_matrix_get_ncol(mask_position.my),
3021 cpl_matrix_get_nrow(mask_position.my));
3024 cpl_matrix_delete(mask_position.my);
3031 cpl_propertylist_set_int(properties, GIALIAS_NAXIS1,
3032 cpl_image_get_size_x(_result));
3033 cpl_propertylist_set_int(properties, GIALIAS_NAXIS2,
3034 cpl_image_get_size_y(_result));
3035 cpl_propertylist_set_int(properties, GIALIAS_BITPIX, -32);
3036 cpl_propertylist_set_double(properties, GIALIAS_BZERO, 0.);
3037 cpl_propertylist_set_double(properties, GIALIAS_BSCALE, 1.);
3039 cpl_propertylist_append_int(properties, GIALIAS_LOCNX,
3040 cpl_image_get_size_y(_result));
3041 cpl_propertylist_append_int(properties, GIALIAS_LOCNS,
3042 cpl_image_get_size_x(_result));
3044 if (config->centroid) {
3045 cpl_propertylist_append_string(properties, GIALIAS_LMETHOD,
3049 cpl_propertylist_append_string(properties, GIALIAS_LMETHOD,
3053 if (config->normalize) {
3054 cpl_propertylist_append_int(properties, GIALIAS_LNORMALIZE,
3058 cpl_propertylist_append_int(properties, GIALIAS_LNORMALIZE,
3062 cpl_propertylist_append_bool(properties, GIALIAS_LFULLLOC, config->full);
3063 cpl_propertylist_append_int(properties, GIALIAS_LOCYDEG, config->yorder);
3064 cpl_propertylist_append_int(properties, GIALIAS_LOCWDEG, config->worder);
3065 cpl_propertylist_append_double(properties, GIALIAS_LEXTRAWID,
3067 cpl_propertylist_append_double(properties, GIALIAS_LNOISEMULT,
3070 cpl_propertylist_append_double(properties, GIALIAS_LCLIPSIGMA,
3072 cpl_propertylist_append_int(properties, GIALIAS_LCLIPNITER,
3073 config->iterations);
3074 cpl_propertylist_append_double(properties, GIALIAS_LCLIPMFRAC,
3078 if (cpl_propertylist_has(properties, GIALIAS_GIRFTYPE)) {
3079 cpl_propertylist_set_string(properties, GIALIAS_GIRFTYPE,
"LOCY");
3082 cpl_propertylist_append_string(properties, GIALIAS_GIRFTYPE,
"LOCY");
3084 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
"GIRAFFE "
3085 "localization centroid");
3092 cpl_matrix_get_ncol(mask_position.mw),
3093 cpl_matrix_get_nrow(mask_position.mw));
3096 cpl_matrix_delete(mask_position.mw);
3103 cpl_propertylist_set_int(properties, GIALIAS_NAXIS1,
3104 cpl_image_get_size_x(_result));
3105 cpl_propertylist_set_int(properties, GIALIAS_NAXIS2,
3106 cpl_image_get_size_y(_result));
3108 if (cpl_propertylist_has(properties, GIALIAS_GIRFTYPE)) {
3109 cpl_propertylist_set_string(properties, GIALIAS_GIRFTYPE,
3113 cpl_propertylist_append_string(properties, GIALIAS_GIRFTYPE,
3116 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
"GIRAFFE "
3117 "localization half-width");
3122 locc = cpl_table_new(cpl_matrix_get_ncol(mask_coeffs.my));
3124 cpl_table_new_column(locc,
"BUTTON", CPL_TYPE_INT);
3125 for (i = 0; i < cpl_table_get_nrow(locc); i++) {
3126 cpl_table_set_int(locc,
"BUTTON", i, i);
3129 for (i = 0; i < cpl_matrix_get_nrow(mask_coeffs.my); i++) {
3130 cxchar *label = NULL;
3132 cx_asprintf(&label,
"YC%d", i);
3133 cpl_table_new_column(locc, label, CPL_TYPE_DOUBLE);
3139 cpl_table_delete(locc);
3141 _my = cpl_matrix_transpose_create(mask_coeffs.my);
3143 cpl_matrix_delete(_my);
3144 cpl_matrix_delete(mask_coeffs.my);
3151 pname = cx_string_new();
3153 for (i = 0; i < cpl_matrix_get_ncol(mask_coeffs.mw); i++) {
3154 cx_string_sprintf(pname,
"%s%d", GIALIAS_LOCWIDCOEF, i);
3155 cpl_propertylist_append_double(properties, cx_string_get(pname),
3156 cpl_matrix_get(mask_coeffs.mw, 0, i));
3159 cx_string_delete(pname);
3160 cpl_matrix_delete(mask_coeffs.mw);
3162 cpl_propertylist_update_string(properties, GIALIAS_GIRFTYPE,
3164 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
"GIRAFFE "
3165 "localization fit coefficients");
3194 GiLocalizeConfig *config = NULL;
3201 config = cx_calloc(1,
sizeof *config);
3208 config->full = TRUE;
3209 config->centroid = TRUE;
3210 config->threshold = GILOCALIZE_THRESHOLD_LOCAL;
3213 p = cpl_parameterlist_find(list,
"giraffe.localization.mode");
3214 s = cpl_parameter_get_string(p);
3215 if (strcmp(s,
"siwc") == 0) {
3216 config->full = FALSE;
3219 p = cpl_parameterlist_find(list,
"giraffe.localization.start");
3220 config->start = cpl_parameter_get_int(p);
3222 p = cpl_parameterlist_find(list,
"giraffe.localization.retries");
3223 config->retries = cpl_parameter_get_int(p);
3225 p = cpl_parameterlist_find(list,
"giraffe.localization.binsize");
3226 config->binsize = cpl_parameter_get_int(p);
3228 p = cpl_parameterlist_find(list,
"giraffe.localization.ewidth");
3229 config->ewidth = cpl_parameter_get_double(p);
3231 p = cpl_parameterlist_find(list,
"giraffe.localization.ywidth");
3232 config->ywidth = cpl_parameter_get_int(p);
3234 p = cpl_parameterlist_find(list,
"giraffe.localization.center");
3235 s = cpl_parameter_get_string(p);
3236 if (!strcmp(s,
"hwidth")) {
3237 config->centroid = FALSE;
3240 p = cpl_parameterlist_find(list,
"giraffe.localization.normalize");
3241 config->normalize = cpl_parameter_get_bool(p);
3243 p = cpl_parameterlist_find(list,
"giraffe.localization.threshold");
3244 s = cpl_parameter_get_string(p);
3246 if (strncmp(s,
"global", 6) == 0) {
3247 config->threshold = GILOCALIZE_THRESHOLD_GLOBAL;
3249 else if (strncmp(s,
"row", 3) == 0) {
3250 config->threshold = GILOCALIZE_THRESHOLD_ROW;
3253 config->threshold = GILOCALIZE_THRESHOLD_LOCAL;
3256 p = cpl_parameterlist_find(list,
"giraffe.localization.noise");
3257 config->noise = cpl_parameter_get_double(p);
3259 p = cpl_parameterlist_find(list,
"giraffe.localization.ron");
3260 config->ron = cpl_parameter_get_double(p);
3262 p = cpl_parameterlist_find(list,
"giraffe.localization.yorder");
3263 config->yorder = cpl_parameter_get_int(p);
3265 p = cpl_parameterlist_find(list,
"giraffe.localization.worder");
3266 config->worder = cpl_parameter_get_int(p);
3268 p = cpl_parameterlist_find(list,
"giraffe.localization.sigma");
3269 config->sigma = cpl_parameter_get_double(p);
3271 p = cpl_parameterlist_find(list,
"giraffe.localization.iterations");
3272 config->iterations = cpl_parameter_get_int(p);
3274 p = cpl_parameterlist_find(list,
"giraffe.localization.fraction");
3275 config->fraction = cpl_parameter_get_double(p);
3329 p = cpl_parameter_new_enum(
"giraffe.localization.mode",
3331 "Localization mode: Use all spectra "
3332 "or the 5 SIWC spectra",
3333 "giraffe.localization",
3334 "all", 2,
"all",
"siwc");
3335 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"sloc-mode");
3336 cpl_parameterlist_append(list, p);
3339 p = cpl_parameter_new_value(
"giraffe.localization.start",
3342 "giraffe.localization",
3344 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"sloc-start");
3345 cpl_parameterlist_append(list, p);
3348 p = cpl_parameter_new_value(
"giraffe.localization.retries",
3350 "Initial localization detection "
3352 "giraffe.localization",
3354 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"sloc-retries");
3355 cpl_parameterlist_append(list, p);
3358 p = cpl_parameter_new_value(
"giraffe.localization.binsize",
3360 "Initial localization detection "
3362 "giraffe.localization",
3364 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"sloc-binsize");
3365 cpl_parameterlist_append(list, p);
3368 p = cpl_parameter_new_value(
"giraffe.localization.ewidth",
3370 "Localization detection extra width.",
3371 "giraffe.localization",
3373 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"sloc-ewidth");
3374 cpl_parameterlist_append(list, p);
3377 p = cpl_parameter_new_value(
"giraffe.localization.ywidth",
3379 "Full width [pxl] of the equilizing "
3380 "filter (distance between two "
3381 "adjacent fibers).",
3382 "giraffe.localization",
3384 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"sloc-ywidth");
3385 cpl_parameterlist_append(list, p);
3388 p = cpl_parameter_new_enum(
"giraffe.localization.center",
3390 "Method used for mask center "
3392 "giraffe.localization",
3393 "centroid", 2,
"centroid",
3395 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"sloc-center");
3396 cpl_parameterlist_append(list, p);
3399 p = cpl_parameter_new_value(
"giraffe.localization.normalize",
3401 "Enable spectrum normalization along "
3402 "the dispersion axis.",
3403 "giraffe.localization",
3405 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"sloc-norm");
3406 cpl_parameterlist_append(list, p);
3409 p = cpl_parameter_new_value(
"giraffe.localization.noise",
3411 "Threshold multiplier.",
3412 "giraffe.localization",
3414 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"sloc-noise");
3415 cpl_parameterlist_append(list, p);
3418 p = cpl_parameter_new_enum(
"giraffe.localization.threshold",
3420 "Selects thresholding algorithm: local, "
3422 "giraffe.localization",
3423 "local", 3,
"local",
"row",
"global");
3424 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"sloc-threshold");
3425 cpl_parameterlist_append(list, p);
3428 p = cpl_parameter_new_value(
"giraffe.localization.ron",
3430 "New bias sigma (RON) value for dark "
3432 "giraffe.localization",
3434 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"sloc-ron");
3435 cpl_parameterlist_append(list, p);
3438 p = cpl_parameter_new_value(
"giraffe.localization.yorder",
3440 "Order of Chebyshev polynomial fit.",
3441 "giraffe.localization",
3443 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"sloc-yorder");
3444 cpl_parameterlist_append(list, p);
3447 p = cpl_parameter_new_value(
"giraffe.localization.worder",
3449 "Order of Chebyshev 2D polynomial fit.",
3450 "giraffe.localization",
3452 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"sloc-worder");
3453 cpl_parameterlist_append(list, p);
3456 p = cpl_parameter_new_value(
"giraffe.localization.sigma",
3458 "Localization clipping: sigma threshold "
3460 "giraffe.localization",
3462 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"sloc-sigma");
3463 cpl_parameterlist_append(list, p);
3466 p = cpl_parameter_new_value(
"giraffe.localization.iterations",
3468 "Localization clipping: number of "
3470 "giraffe.localization",
3472 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"sloc-niter");
3473 cpl_parameterlist_append(list, p);
3476 p = cpl_parameter_new_range(
"giraffe.localization.fraction",
3478 "Localization clipping: minimum fraction "
3479 "of points accepted/total.",
3480 "giraffe.localization",
3482 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"sloc-mfrac");
3483 cpl_parameterlist_append(list, p);
cxint giraffe_array_sort(cxdouble *array, cxsize size)
Sorts an array in ascending order.
cpl_propertylist * giraffe_image_get_properties(const GiImage *self)
Get the properties of an image.
cxint giraffe_image_copy_matrix(GiImage *self, cpl_matrix *matrix)
Copies matrix elements into an image.
cpl_image * giraffe_image_get(const GiImage *self)
Gets the image data.
GiImage * giraffe_image_create(cpl_type type, cxint nx, cxint ny)
Creates an image container of a given type.
cxint giraffe_image_set_properties(GiImage *self, cpl_propertylist *properties)
Attaches a property list to an image.
GiLocalizeConfig * giraffe_localize_config_create(cpl_parameterlist *list)
Creates a setup structure for the spectrum localization.
void giraffe_localize_config_destroy(GiLocalizeConfig *config)
Destroys a spectrum localization setup structure.
void giraffe_localize_config_add(cpl_parameterlist *list)
Adds parameters for the spectrum localization.
cxint giraffe_localize_spectra(GiLocalization *result, GiImage *image, GiTable *fibers, GiLocalization *master, GiImage *badpixels, GiLocalizeConfig *config)
Finds the location of spectra in a Giraffe observation.
cpl_image * giraffe_matrix_create_image(const cpl_matrix *matrix)
Converts a matrix into an image.
cpl_matrix * giraffe_matrix_leastsq(const cpl_matrix *mA, const cpl_matrix *mB)
Computes the solution of an equation using a pseudo-inverse.
cxdouble giraffe_matrix_sigma_mean(const cpl_matrix *matrix, cxdouble mean)
Compute sigma of matrix elements, with a given mean value.
void gi_warning(const cxchar *format,...)
Log a warning.
void gi_error(const cxchar *format,...)
Log an error message.
void gi_message(const cxchar *format,...)
Log a normal message.
cxint giraffe_table_copy_matrix(GiTable *table, const cxchar *name, cpl_matrix *matrix)
Copies matrix elements into a table.
GiTable * giraffe_table_create(cpl_table *table, cpl_propertylist *properties)
Creates a Giraffe table from a table and a property list.
cpl_propertylist * giraffe_table_get_properties(const GiTable *self)
Gets the table properties.
cpl_table * giraffe_table_get(const GiTable *self)
Get the table data from a Giraffe table.
cxdouble giraffe_propertylist_get_conad(const cpl_propertylist *properties)
Retrieve the ADU to electrons conversion factor from the given properties.
cxdouble giraffe_propertylist_get_ron(const cpl_propertylist *properties)
Retrieve the read-out noise from the given properties.
GiInstrumentMode giraffe_get_mode(cpl_propertylist *properties)
Determines the instrument mode from a property list.