29 #include <cxstrutils.h>
31 #include <cpl_parameterlist.h>
32 #include <cpl_matrix.h>
33 #include <cpl_table.h>
44 #include "gipsfdata.h"
47 #include "gilocalization.h"
48 #include "gimessages.h"
49 #include "gifiberutils.h"
51 #include "giextract.h"
63 PROFILE_PSFEXP = 1 << 1,
64 PROFILE_PSFEXP2 = 1 << 2,
65 PROFILE_GAUSSIAN = 1 << 3
68 typedef enum GiProfileId GiProfileId;
75 struct GiExtractOptimalConfig {
89 typedef struct GiExtractOptimalConfig GiExtractOptimalConfig;
96 struct GiExtractHorneConfig {
105 typedef struct GiExtractHorneConfig GiExtractHorneConfig;
108 struct GiExtractionData {
115 typedef struct GiExtractionData GiExtractionData;
118 struct GiExtractionSlice {
126 cpl_matrix* variance;
130 typedef struct GiExtractionSlice GiExtractionSlice;
133 struct GiExtractionPsfLimits {
140 typedef struct GiExtractionPsfLimits GiExtractionPsfLimits;
143 struct GiExtractionWorkspace {
151 typedef struct GiExtractionWorkspace GiExtractionWorkspace;
154 struct GiVirtualSlit {
158 cxdouble extra_width;
169 typedef struct GiVirtualSlit GiVirtualSlit;
176 inline static GiExtractionSlice*
177 _giraffe_extractionslice_new(cxint nflx, cxint ndata, cxint nbkg)
180 GiExtractionSlice*
self = cx_malloc(
sizeof *
self);
185 self->fsize = nflx + nbkg;
188 self->flux = cpl_matrix_new(self->fsize, 1);
189 self->variance = cpl_matrix_new(self->fsize, 1);
190 self->model = cpl_matrix_new(self->msize, 1);
198 _giraffe_extractionslice_delete(GiExtractionSlice*
self)
202 if (self->model != NULL) {
203 cpl_matrix_delete(self->model);
207 if (self->variance != NULL) {
208 cpl_matrix_delete(self->variance);
209 self->variance = NULL;
212 if (self->flux != NULL) {
213 cpl_matrix_delete(self->flux);
225 inline static GiExtractionPsfLimits*
226 _giraffe_extraction_psflimits_new(cxint size)
229 GiExtractionPsfLimits*
self = cx_malloc(
sizeof *
self);
233 self->ymin = cx_calloc(self->size,
sizeof(cxint));
234 self->ymax = cx_calloc(self->size,
sizeof(cxint));
242 _giraffe_extraction_psflimits_delete(GiExtractionPsfLimits*
self)
246 if (self->ymin != NULL) {
250 if (self->ymax != NULL) {
262 inline static GiExtractionWorkspace*
263 _giraffe_optimal_workspace_new(cxint m, cxint n)
266 GiExtractionWorkspace*
self = cx_malloc(
sizeof *
self);
269 self->atw = cpl_matrix_new(m, n);
270 self->atwa = cpl_matrix_new(m, m);
271 self->c = cpl_matrix_new(m, m);
272 self->atws = cpl_matrix_new(m, 1);
274 self->tmp = cpl_matrix_new(m, m);
282 _giraffe_optimal_workspace_delete(GiExtractionWorkspace*
self)
286 if (self->atws != NULL) {
287 cpl_matrix_delete(self->atws);
290 if (self->atwa != NULL) {
291 cpl_matrix_delete(self->atwa);
294 if (self->c != NULL) {
295 cpl_matrix_delete(self->c);
298 if (self->atw != NULL) {
299 cpl_matrix_delete(self->atw);
302 if (self->tmp != NULL) {
303 cpl_matrix_delete(self->tmp);
320 _giraffe_virtualslit_allocate(GiVirtualSlit*
self)
323 if ((
self != NULL) && (self->width > 0)) {
325 self->position = cx_calloc(self->width,
sizeof(cxdouble));
326 self->signal = cx_calloc(self->width,
sizeof(cxdouble));
327 self->variance = cx_calloc(self->width,
sizeof(cxdouble));
328 self->fraction = cx_calloc(self->width,
sizeof(cxdouble));
330 self->mask = cx_calloc(self->width,
sizeof(cxdouble));
331 self->offset = cx_calloc(self->width,
sizeof(cxdouble));
340 inline static GiVirtualSlit*
341 _giraffe_virtualslit_new(cxdouble extra_width)
344 GiVirtualSlit*
self = cx_calloc(1,
sizeof *
self);
348 self->extra_width = extra_width;
350 self->position = NULL;
352 self->variance = NULL;
353 self->fraction = NULL;
363 _giraffe_virtualslit_clear(GiVirtualSlit*
self)
368 if (self->position != NULL) {
369 cx_free(self->position);
370 self->position = NULL;
373 if (self->signal != NULL) {
374 cx_free(self->signal);
378 if (self->variance != NULL) {
379 cx_free(self->variance);
380 self->variance = NULL;
383 if (self->fraction != NULL) {
384 cx_free(self->fraction);
385 self->fraction = NULL;
388 if (self->mask != NULL) {
393 if (self->offset != NULL) {
394 cx_free(self->offset);
398 self->extra_width = 0.;
410 _giraffe_virtualslit_delete(GiVirtualSlit*
self)
414 _giraffe_virtualslit_clear(
self);
425 _giraffe_virtualslit_setup(GiVirtualSlit*
self, cxint bin,
426 cxdouble center, cxdouble width,
427 const cpl_image* signal,
const cpl_image* variance,
428 const cpl_image* bpixel)
431 register cxint ny = cpl_image_get_size_x(signal);
432 register cxint offset = bin * cpl_image_get_size_x(signal);
434 register cxdouble lower = center - (width +
self->extra_width);
435 register cxdouble upper = center + (width +
self->extra_width);
437 register cxint first = (cxint) floor(lower);
438 register cxint last = (cxint) ceil(upper);
440 const cxdouble* s = cpl_image_get_data_double_const(signal);
441 const cxdouble* v = cpl_image_get_data_double_const(variance);
448 lower = CX_MAX(0., lower);
449 upper = CX_MIN(ny, upper);
451 first = CX_MAX(0, first);
452 last = CX_MIN(ny, last);
454 self->center = center;
455 self->width = last - first + 1;
462 _giraffe_virtualslit_allocate(
self);
464 if (bpixel != NULL) {
466 register cxint k = 0;
467 register cxint y = 0;
469 const cxint* _bpixel = cpl_image_get_data_int_const(bpixel);
472 for (y = first; y <= last; y++) {
474 register cxint ypos = offset + y;
476 cxint ok = (_bpixel[ypos] & GIR_M_PIX_SET) == 0 ? 1 : 0;
479 self->position[k] = y - center;
480 self->fraction[k] = 1.;
482 self->signal[k] = s[ypos];
483 self->variance[k] = v[ypos];
486 self->offset[k] = ypos;
494 register cxint k = 0;
495 register cxint y = 0;
498 for (y = first; y <= last; y++) {
500 register cxint ypos = offset + y;
505 self->position[k] = y - center;
506 self->fraction[k] = 1.;
508 self->signal[k] = s[ypos];
509 self->variance[k] = v[ypos];
512 self->offset[k] = ypos;
526 self->fraction[0] = ((cxdouble)first + 1.) - lower;
527 self->fraction[
self->width - 1] = upper - ((cxdouble)last - 1.);
539 _giraffe_matrix_invert(cpl_matrix* m_inv,
const cpl_matrix* m, cpl_matrix* lu)
544 cxint n = cpl_matrix_get_ncol(m);
546 register cxint sz = n * n *
sizeof(cxdouble);
548 const cxdouble* _m = cpl_matrix_get_data_const(m);
550 cxdouble* _m_inv = cpl_matrix_get_data(m_inv);
551 cxdouble* _m_lu = cpl_matrix_get_data(lu);
553 cpl_array* perm = cpl_array_new(n, CPL_TYPE_INT);
555 register cxint* perm_data = cpl_array_get_data_int(perm);
558 memset(_m_inv, 0, sz);
559 memcpy(_m_lu, _m, sz);
561 if (cpl_matrix_decomp_lu(lu, perm, &i) != 0) {
562 cpl_array_delete(perm);
571 for (i = 0; i < n; ++i) {
572 _m_inv[i * n + perm_data[i]] = 1.;
575 cpl_array_delete(perm);
578 status = cpl_matrix_solve_lu(lu, m_inv, NULL);
581 cpl_matrix_delete(m_inv);
594 inline static cpl_matrix*
595 _giraffe_compute_psf(GiModel* psf,
const cpl_matrix* x)
598 register cxint i = 0;
599 register cxint n = 0;
603 const cxdouble* _x = NULL;
607 cpl_matrix* y = NULL;
609 cx_assert(psf != NULL);
610 cx_assert(x != NULL);
611 cx_assert(cpl_matrix_get_ncol(x) == 1);
613 n = cpl_matrix_get_nrow(x);
615 y = cpl_matrix_new(n, 1);
617 _x = cpl_matrix_get_data_const(x);
618 _y = cpl_matrix_get_data(y);
620 for (i = 0; i < n; i++) {
621 giraffe_model_set_argument(psf,
"x", _x[i]);
622 giraffe_model_evaluate(psf, &_y[i], &status);
625 cpl_matrix_delete(y);
641 _giraffe_horne_extract_slit(GiExtractionData* result,
642 const GiVirtualSlit* vslit, GiModel* psf,
643 const GiExtractHorneConfig* config)
653 cxdouble* tdata = NULL;
654 cxdouble* _mnpsf = NULL;
656 cpl_matrix* mnpsf = NULL;
657 cpl_matrix* mvslit = NULL;
665 mvslit = cpl_matrix_wrap(vslit->width, 1, vslit->position);
666 mnpsf = _giraffe_compute_psf(psf, mvslit);
668 cpl_matrix_unwrap(mvslit);
680 _mnpsf = cpl_matrix_get_data(mnpsf);
684 for (i = 0; i < vslit->width; ++i) {
685 _mnpsf[i] = CX_MAX(_mnpsf[i], 0.);
689 for (i = 0; i < vslit->width; ++i) {
698 tdata = cx_malloc(vslit->width *
sizeof(cxdouble));
703 while (i < vslit->width) {
704 if (vslit->mask[i] > 0) {
705 tdata[ngood] = CX_MAX(vslit->signal[i], 0.);
713 bkg = 0.5 * (tdata[0] + tdata[1]);
730 cxint niter = config->clip.iterations;
731 cxint nmin = (cxint)config->clip.fraction;
733 cxdouble sigma = config->clip.level * config->clip.level;
734 cxdouble* variance = NULL;
745 for (i = 0; i < vslit->width; ++i) {
746 if (vslit->mask[i] != 0) {
747 flx += (vslit->signal[i] - bkg) * vslit->fraction[i];
748 norm += vslit->fraction[i] * _mnpsf[i];
760 variance = cx_calloc(vslit->width,
sizeof(cxdouble));
762 for (i = 0; i < vslit->width; ++i) {
764 register cxdouble ve = flx * _mnpsf[i] + bkg;
766 variance[i] = vslit->variance[i] + fabs(vslit->fraction[i] * ve);
777 while ((iteration < niter) && (ngood > nmin) && (nreject != 0)) {
794 for (i = 0; i < vslit->width; ++i) {
796 if (vslit->mask[i] != 0) {
798 cxdouble m = vslit->signal[i] - bkg - flx * _mnpsf[i];
800 m *= vslit->fraction[i];
801 m *= m / variance[i] ;
812 if ((sigma > 0.) && (mmax > sigma)) {
813 vslit->mask[imax] = 0;
823 for (i = 0; i < vslit->width; ++i) {
825 if (vslit->mask[i] != 0) {
827 register cxdouble data = vslit->signal[i] - bkg;
828 register cxdouble p = _mnpsf[i];
830 data *= vslit->fraction[i];
831 p *= vslit->fraction[i];
833 norm += p * p / variance[i];
834 _flx += p * data / variance[i];
849 for (i = 0; i < vslit->width; ++i) {
851 register cxdouble ve = flx * _mnpsf[i] + bkg;
853 variance[i] = vslit->variance[i] + fabs(vslit->fraction[i] * ve);
866 cpl_matrix_delete(mnpsf);
870 result->error = sqrt(var);
871 result->position = vslit->center;
872 result->npixels = ngood;
874 return ngood == 0 ? 1 : 0;
923 _giraffe_optimal_extract_slice(GiExtractionSlice* slice,
924 const cpl_matrix* AT,
927 GiExtractionPsfLimits* limits,
928 GiExtractionWorkspace* ws)
931 register cxint i = 0;
932 register cxint n = cpl_matrix_get_ncol(AT);
933 register cxint m = cpl_matrix_get_nrow(AT);
937 const cxdouble* at = cpl_matrix_get_data_const(AT);
938 const cxdouble* w = cpl_matrix_get_data_const(W);
939 const cxdouble* s = cpl_matrix_get_data_const(S);
940 const cxdouble* c = cpl_matrix_get_data_const(ws->c);
942 cxdouble* atw = cpl_matrix_get_data(ws->atw);
943 cxdouble* atwa = cpl_matrix_get_data(ws->atwa);
944 cxdouble* atws = cpl_matrix_get_data(ws->atws);
945 cxdouble* sf = cpl_matrix_get_data(slice->flux);
946 cxdouble* sv = cpl_matrix_get_data(slice->variance);
947 cxdouble* sm = cpl_matrix_get_data(slice->model);
950 for (i = 0; i < m; ++i) {
952 register cxint j = 0;
953 register cxint im = i * m;
954 register cxint in = i * n;
955 register cxint ymin = limits->ymin[i];
956 register cxint ymax = limits->ymax[i];
961 for (j = 0; j < n; ++j) {
963 register cxint k = in + j;
966 atw[k] = w[j] * at[k];
967 atws[i] += atw[k] * s[j];
971 for (j = 0; j < i; ++j) {
973 register cxint k = 0;
974 register cxint l = im + j;
977 for (k = ymin; k < ymax; ++k) {
978 atwa[l] += atw[in + k] * at[j * n + k];
981 atwa[j * m + i] = atwa[l];
987 for (j = ymin; j < ymax; ++j) {
988 atwa[im + i] += atw[in + j] * at[in + j];
994 status = _giraffe_matrix_invert(ws->c, ws->atwa, ws->tmp);
1000 for (i = 0; i < m; ++i) {
1002 register cxint j = 0;
1003 register cxint im = i * m;
1009 for (j = 0; j < m; ++j) {
1010 sf[i] += c[im + j] * atws[j];
1015 for (i = 0; i < n; ++i) {
1017 register cxint j = 0;
1022 for (j = 0; j < m; ++j) {
1023 sm[i] += at[j * n + i] * sf[j];
1057 _giraffe_extract_summation(
const cpl_image* mz,
const cpl_image* mvarz,
1058 const cpl_table* fibers,
const cpl_image* my,
1059 const cpl_image* mw, cpl_image* mbpx,
1060 cpl_image* ms, cpl_image* mse,
1061 cpl_image* msn, cpl_image* msy)
1066 const cxchar* idx = NULL;
1068 cxint ny = cpl_image_get_size_x(mz);
1069 cxint nfibers = cpl_table_get_nrow(fibers);
1070 cxint nspectra = cpl_image_get_size_x(my);
1071 cxint nbins = cpl_image_get_size_y(my);
1073 const cxdouble* pixels = cpl_image_get_data_double_const(mz);
1074 const cxdouble* variances = cpl_image_get_data_double_const(mvarz);
1075 const cxdouble* locy = cpl_image_get_data_double_const(my);
1076 const cxdouble* locw = cpl_image_get_data_double_const(mw);
1078 cxdouble* flux = cpl_image_get_data_double(ms);
1079 cxdouble* flux_error = cpl_image_get_data_double(mse);
1080 cxdouble* flux_npixels = cpl_image_get_data_double(msn);
1081 cxdouble* flux_ypos = cpl_image_get_data_double(msy);
1089 cx_assert(nfibers <= nspectra);
1093 cx_assert(cpl_table_has_column(fibers, idx) != 0);
1097 const cxint* bpx = cpl_image_get_data_int(mbpx);
1099 for (nn = 0; nn < nfibers; nn++) {
1102 register cxint ns = cpl_table_get_int(fibers, idx, nn, NULL) - 1;
1105 for (x = 0; x < cpl_image_get_size_y(mz) && x < nbins; x++) {
1109 cxint lx = x * nspectra + ns;
1110 cxint sx = x * nfibers + nn;
1112 cxdouble ylower = locy[lx] - locw[lx];
1113 cxdouble yupper = locy[lx] + locw[lx];
1116 cxdouble error2 = 0.;
1120 flux_npixels[sx] = 0.;
1121 flux_error[sx] = 0.;
1129 if (locw[lx] <= 0.0) {
1141 ylo = (cxint) ceil(ylower);
1142 yup = (cxint) floor(yupper);
1145 if (yup < 0. || ylo - 1 >= ny) {
1165 if (!(bpx[x * ny + y] & GIR_M_PIX_SET)) {
1167 cxdouble extcoeff = (cxdouble)ylo - ylower;
1168 cxdouble extcoeff2 = extcoeff * extcoeff;
1169 cxdouble px = CX_MAX(pixels[x * ny + y], 0.);
1171 flux[sx] = pixels[x * ny + y] * extcoeff;
1172 flux_npixels[sx] = extcoeff;
1173 error2 = variances[x * ny + y] * extcoeff2;
1175 zsum = px * extcoeff;
1176 ysum = y * px * extcoeff;
1187 for (y = CX_MAX(ylo, 0); (y < yup) && (y < ny); y++) {
1189 if (!(bpx[x * ny + y] & GIR_M_PIX_SET)) {
1191 cxdouble px = CX_MAX(pixels[x * ny + y], 0.);
1193 flux[sx] += pixels[x * ny + y];
1194 flux_npixels[sx] += 1.0;
1195 error2 += variances[x * ny + y];
1213 if (!(bpx[x * ny + y] & GIR_M_PIX_SET)) {
1215 cxdouble extcoeff = yupper - (cxdouble)yup;
1216 cxdouble extcoeff2 = extcoeff * extcoeff;
1217 cxdouble px = CX_MAX(pixels[x * ny + y], 0.);
1219 flux[sx] += pixels[x * ny + y] * extcoeff;
1220 flux_npixels[sx] += extcoeff;
1221 error2 += variances[x * ny + y] * extcoeff2;
1223 zsum += px * extcoeff;
1224 ysum += y * px * extcoeff;
1230 flux_error[sx] = sqrt(error2);
1235 if (fabs(ysum) < DBL_EPSILON || fabs(zsum) < DBL_EPSILON) {
1236 flux_ypos[sx] = 0.5 * (yupper + ylower);
1239 flux_ypos[sx] = ysum / zsum;
1249 for (nn = 0; nn < nfibers; nn++) {
1252 register cxint ns = cpl_table_get_int(fibers, idx,
1256 for (x = 0; x < cpl_image_get_size_y(mz) && x < nbins; x++) {
1260 cxint lx = x * nspectra + ns;
1261 cxint sx = x * nfibers + nn;
1263 cxdouble yupper, ylower;
1266 cxdouble error2 = 0.;
1270 flux_npixels[sx] = 0.;
1271 flux_error[sx] = 0.;
1279 if (locw[lx] <= 0.0) {
1291 yupper = locy[lx] + locw[lx];
1292 ylower = locy[lx] - locw[lx];
1294 ylo = (cxint) ceil(ylower);
1295 yup = (cxint) floor(yupper);
1298 if (yup < 0. || ylo - 1 >= ny) {
1318 cxdouble extcoeff = (cxdouble)ylo - ylower;
1319 cxdouble extcoeff2 = extcoeff * extcoeff;
1320 cxdouble px = CX_MAX(pixels[x * ny + y], 0.);
1322 flux[sx] = pixels[x * ny + y] * extcoeff;
1323 flux_npixels[sx] = extcoeff;
1324 error2 = variances[x * ny + y] * extcoeff2;
1326 zsum = px * extcoeff;
1327 ysum = y * px * extcoeff;
1336 for (y = CX_MAX(ylo, 0); (y < yup) && (y < ny); y++) {
1338 cxdouble px = CX_MAX(pixels[x * ny + y], 0.);
1340 flux[sx] += pixels[x * ny + y];
1341 flux_npixels[sx] += 1.0;
1342 error2 += variances[x * ny + y];
1357 cxdouble extcoeff = yupper - (cxdouble)yup;
1358 cxdouble extcoeff2 = extcoeff * extcoeff;
1359 cxdouble px = CX_MAX(pixels[x * ny + y], 0.);
1361 flux[sx] += pixels[x * ny + y] * extcoeff;
1362 flux_npixels[sx] += extcoeff;
1363 error2 += variances[x * ny + y] * extcoeff2;
1365 zsum += px * extcoeff;
1366 ysum += y * px * extcoeff;
1370 flux_error[sx] = sqrt(error2);
1375 if (fabs(ysum) < DBL_EPSILON || fabs(zsum) < DBL_EPSILON) {
1376 flux_ypos[sx] = 0.5 * (yupper + ylower);
1379 flux_ypos[sx] = ysum / zsum;
1415 _giraffe_extract_horne(
const cpl_image* mz,
const cpl_image* mzvar,
1416 const cpl_table* fibers,
const cpl_image* my,
1417 const cpl_image* mw,
const GiPsfData* psfdata,
1418 cpl_image* mbpx, cpl_image* ms, cpl_image* mse,
1419 cpl_image* msn, cpl_image* msy,
1420 const GiExtractHorneConfig* config)
1423 const cxchar* idx = NULL;
1430 const cxdouble* locy = NULL;
1431 const cxdouble* locw = NULL;
1432 const cxdouble* width = NULL;
1433 const cxdouble* exponent = NULL;
1435 GiModel* psfmodel = NULL;
1438 cx_assert(mz != NULL);
1439 cx_assert(mzvar != NULL);
1441 cx_assert(fibers != NULL);
1443 cx_assert(my != NULL);
1444 cx_assert(mw != NULL);
1446 cx_assert(psfdata != NULL);
1448 cx_assert(ms != NULL);
1449 cx_assert(mse != NULL);
1450 cx_assert(msn != NULL);
1451 cx_assert(msy != NULL);
1453 cx_assert(config != NULL);
1455 ny = cpl_image_get_size_x(mz);
1456 nx = cpl_image_get_size_y(mz);
1457 nfibers = cpl_table_get_nrow(fibers);
1459 locy = cpl_image_get_data_double_const(my);
1460 locw = cpl_image_get_data_double_const(mw);
1462 cx_assert((ny == cpl_image_get_size_x(mzvar)) &&
1463 (nx == cpl_image_get_size_y(mzvar)));
1465 cx_assert(cpl_image_get_size_x(my) == cpl_image_get_size_x(mw));
1466 cx_assert(cpl_image_get_size_y(my) == cpl_image_get_size_y(mw));
1468 cx_assert(giraffe_psfdata_fibers(psfdata) ==
1469 (cxsize)cpl_image_get_size_x(my));
1470 cx_assert(giraffe_psfdata_bins(psfdata) ==
1471 (cxsize)cpl_image_get_size_y(my));
1473 cx_assert((nfibers == cpl_image_get_size_x(ms)) &&
1474 (nx == cpl_image_get_size_y(ms)));
1475 cx_assert((nfibers == cpl_image_get_size_x(mse)) &&
1476 (nx == cpl_image_get_size_y(mse)));
1477 cx_assert((nfibers == cpl_image_get_size_x(msn)) &&
1478 (nx == cpl_image_get_size_y(msn)));
1479 cx_assert((nfibers == cpl_image_get_size_x(msy)) &&
1480 (nx == cpl_image_get_size_y(msy)));
1482 cx_assert((mbpx == NULL) || ((ny == cpl_image_get_size_x(mbpx)) &&
1483 (nx == cpl_image_get_size_y(mbpx))));
1493 cx_assert(cpl_table_has_column(fibers, idx) != 0);
1501 if (giraffe_psfdata_contains(psfdata,
"Center") == FALSE) {
1505 if (giraffe_psfdata_contains(psfdata,
"Width2") == TRUE) {
1506 exponent = cpl_image_get_data_const(giraffe_psfdata_get_data(psfdata,
1510 width = cpl_image_get_data_const(giraffe_psfdata_get_data(psfdata,
1518 psfmodel = giraffe_model_new(giraffe_psfdata_get_model(psfdata));
1520 if (psfmodel == NULL) {
1524 giraffe_model_set_parameter(psfmodel,
"Center", 0.);
1525 giraffe_model_set_parameter(psfmodel,
"Amplitude", 1.);
1526 giraffe_model_set_parameter(psfmodel,
"Background", 0.);
1533 for (fiber = 0; fiber < nfibers; ++fiber) {
1535 register cxint bin = 0;
1536 register cxint fidx = cpl_table_get_int(fibers, idx, fiber, NULL) - 1;
1538 cxint nbins = CX_MIN(nx, cpl_image_get_size_y(my));
1540 cxdouble* _ms = cpl_image_get_data_double(ms);
1541 cxdouble* _mse = cpl_image_get_data_double(mse);
1542 cxdouble* _msy = cpl_image_get_data_double(msy);
1543 cxdouble* _msn = cpl_image_get_data_double(msn);
1546 for (bin = 0; bin < nbins; bin++) {
1548 register cxint lpos = bin * cpl_image_get_size_x(my) + fidx;
1549 register cxint spos = bin * nfibers + fiber;
1554 register cxdouble lcenter = locy[lpos];
1555 register cxdouble lwidth = locw[lpos];
1557 register cxdouble ylower = lcenter - lwidth;
1558 register cxdouble yupper = lcenter + lwidth;
1560 GiVirtualSlit* vslit = NULL;
1562 GiExtractionData result = {0., 0., 0., 0.};
1569 if ((lwidth <= 0.) || (yupper < 0.) || (ylower > ny)) {
1577 vslit = _giraffe_virtualslit_new(config->ewidth);
1579 vwidth = _giraffe_virtualslit_setup(vslit, bin, lcenter, lwidth,
1583 _giraffe_virtualslit_delete(vslit);
1594 giraffe_model_set_parameter(psfmodel,
"Width1", width[lpos]);
1596 if (exponent != NULL) {
1597 giraffe_model_set_parameter(psfmodel,
"Width2",
1607 status = _giraffe_horne_extract_slit(&result, vslit, psfmodel,
1610 _giraffe_virtualslit_delete(vslit);
1615 giraffe_model_delete(psfmodel);
1621 _ms[spos] = result.value;
1622 _mse[spos] = result.error;
1623 _msy[spos] = result.position;
1624 _msn[spos] = result.npixels;
1631 giraffe_model_delete(psfmodel);
1645 _giraffe_optimal_build_profiles(cpl_matrix* profiles,
1646 GiExtractionPsfLimits* limits,
1647 const cpl_image* my,
const cpl_image* mw,
1648 const cpl_table* fibers, cxint bin,
1649 GiModel* psf,
const cxdouble* width,
1650 const cxdouble* exponent, cxdouble wfactor)
1656 cxint nfibers = cpl_table_get_nrow(fibers);
1657 cxint ny = cpl_matrix_get_ncol(profiles);
1659 const cxdouble* locy = cpl_image_get_data_double_const(my);
1660 const cxdouble* locw = cpl_image_get_data_double_const(mw);
1662 cxdouble* _profiles = cpl_matrix_get_data(profiles);
1664 cxdouble* ypos = NULL;
1667 cx_assert(cpl_table_has_column(fibers, idx) != 0);
1668 cx_assert((limits == NULL) ||
1669 (cpl_matrix_get_nrow(profiles) == limits->size));
1671 ypos = cx_calloc(ny,
sizeof(cxdouble));
1673 for (fiber = 0; fiber < nfibers; ++fiber) {
1675 register cxint i = 0;
1676 register cxint y = 0;
1677 register cxint k = 0;
1679 cxint fidx = cpl_table_get_int(fibers, idx, fiber, NULL) - 1;
1680 cxint lpos = bin * cpl_image_get_size_x(my) + fidx;
1682 register cxdouble lcenter = locy[lpos];
1683 register cxdouble lwidth = locw[lpos];
1685 register cxdouble ylower = lcenter - fabs(wfactor) * lwidth;
1686 register cxdouble yupper = lcenter + fabs(wfactor) * lwidth;
1688 register cxint first = (cxint) floor(ylower);
1689 register cxint last = (cxint) ceil(yupper);
1691 register cxint vwidth = 0;
1694 cxdouble* _mnpsf = NULL;
1696 cpl_matrix* positions = NULL;
1697 cpl_matrix* mnpsf = NULL;
1704 ylower = CX_MAX(0., ylower);
1705 yupper = CX_MIN(ny - 1., yupper);
1707 first = CX_MAX(0, first);
1708 last = CX_MIN(ny - 1, last);
1710 vwidth = last - first + 1;
1712 if (limits != NULL) {
1713 limits->ymin[fiber] = first;
1714 limits->ymax[fiber] = last + 1;
1722 giraffe_model_set_parameter(psf,
"Width1", width[lpos]);
1724 if (exponent != NULL) {
1725 giraffe_model_set_parameter(psf,
"Width2", exponent[lpos]);
1734 for (y = first; y <= last; ++y) {
1735 ypos[k] = y - lcenter;
1739 positions = cpl_matrix_wrap(vwidth, 1, ypos);
1740 mnpsf = _giraffe_compute_psf(psf, positions);
1742 cpl_matrix_unwrap(positions);
1745 if (mnpsf == NULL) {
1752 _mnpsf = cpl_matrix_get_data(mnpsf);
1754 for (i = 0; i < vwidth; ++i) {
1755 _mnpsf[i] = CX_MAX(_mnpsf[i], 0.);
1759 for (i = 0; i < vwidth; ++i) {
1763 k = fiber * ny + first;
1764 for (y = 0; y < vwidth; ++y) {
1765 _profiles[k + y] = _mnpsf[y];
1768 cpl_matrix_delete(mnpsf);
1782 _giraffe_extract_optimal(
const cpl_image* mz,
const cpl_image* mzvar,
1783 const cpl_table* fibers,
const cpl_image* my,
1784 const cpl_image* mw,
const GiPsfData* psfdata,
1785 cpl_image* mbpx, cpl_image* ms, cpl_image* mse,
1786 cpl_image* msm, cpl_image* msy,
1787 const GiExtractOptimalConfig* config)
1790 const cxbool nolimits = (config->limits == TRUE) ? FALSE : TRUE;
1792 const cxint bkg_nc = config->bkgorder + 1;
1793 const cxint niter = config->clip.iterations;
1795 register cxint i = 0;
1803 const cxdouble wfactor = config->ewidth;
1804 const cxdouble sigma = config->clip.level * config->clip.level;
1805 const cxdouble fraction = config->clip.fraction;
1807 const cxdouble* width = NULL;
1808 const cxdouble* exponent = NULL;
1810 cxdouble* _ypos = NULL;
1811 cxdouble* _bkg_base = NULL;
1812 cxdouble* _profiles = NULL;
1813 cxdouble* _signal = NULL;
1814 cxdouble* _variance = NULL;
1815 cxdouble* _mask = NULL;
1816 cxdouble* _weights = NULL;
1818 cpl_matrix* ypos = NULL;
1819 cpl_matrix* bkg_base = NULL;
1820 cpl_matrix* profiles = NULL;
1821 cpl_matrix* weights = NULL;
1822 cpl_matrix* signal = NULL;
1823 cpl_matrix* variance = NULL;
1824 cpl_matrix* mask = NULL;
1826 GiModel* psfmodel = NULL;
1828 GiExtractionPsfLimits* limits = NULL;
1830 GiExtractionSlice* slice = NULL;
1832 GiExtractionWorkspace* workspace;
1835 cx_assert(mz != NULL);
1836 cx_assert(mzvar != NULL);
1838 cx_assert(fibers != NULL);
1840 cx_assert(my != NULL);
1841 cx_assert(mw != NULL);
1843 cx_assert(psfdata != NULL);
1845 cx_assert(ms != NULL);
1846 cx_assert(mse != NULL);
1847 cx_assert(msm != NULL);
1848 cx_assert(msy != NULL);
1850 ny = cpl_image_get_size_x(mz);
1851 nx = cpl_image_get_size_y(mz);
1853 nfibers = cpl_table_get_nrow(fibers);
1854 nbins = CX_MIN(nx, cpl_image_get_size_y(my));
1856 cx_assert((ny == cpl_image_get_size_x(mzvar)) &&
1857 (nx == cpl_image_get_size_y(mzvar)));
1859 cx_assert(cpl_image_get_size_x(my) == cpl_image_get_size_x(mw));
1860 cx_assert(cpl_image_get_size_y(my) == cpl_image_get_size_y(mw));
1862 cx_assert(giraffe_psfdata_fibers(psfdata) ==
1863 (cxsize)cpl_image_get_size_x(my));
1864 cx_assert(giraffe_psfdata_bins(psfdata) ==
1865 (cxsize)cpl_image_get_size_y(my));
1867 cx_assert((nfibers == cpl_image_get_size_x(ms)) &&
1868 (nx == cpl_image_get_size_y(ms)));
1869 cx_assert((nfibers == cpl_image_get_size_x(mse)) &&
1870 (nx == cpl_image_get_size_y(mse)));
1871 cx_assert((nfibers == cpl_image_get_size_x(msy)) &&
1872 (nx == cpl_image_get_size_y(msy)));
1873 cx_assert((ny == cpl_image_get_size_x(msm)) &&
1874 (nx == cpl_image_get_size_y(msm)));
1876 cx_assert((mbpx == NULL) || ((ny == cpl_image_get_size_x(mbpx)) &&
1877 (nx == cpl_image_get_size_y(mbpx))));
1885 if (giraffe_psfdata_contains(psfdata,
"Center") == FALSE) {
1889 if (giraffe_psfdata_contains(psfdata,
"Width2") == TRUE) {
1890 exponent = cpl_image_get_data_const(giraffe_psfdata_get_data(psfdata,
1894 width = cpl_image_get_data_const(giraffe_psfdata_get_data(psfdata,
1902 psfmodel = giraffe_model_new(giraffe_psfdata_get_model(psfdata));
1904 if (psfmodel == NULL) {
1908 giraffe_model_set_parameter(psfmodel,
"Amplitude", 1.);
1909 giraffe_model_set_parameter(psfmodel,
"Background", 0.);
1910 giraffe_model_set_parameter(psfmodel,
"Center", 0.);
1917 ypos = cpl_matrix_new(ny, 1);
1920 giraffe_model_delete(psfmodel);
1926 _ypos = cpl_matrix_get_data(ypos);
1928 for (i = 0; i < ny; ++i) {
1938 profiles = cpl_matrix_new(nfibers + bkg_nc, ny);
1940 if (profiles == NULL) {
1941 cpl_matrix_delete(ypos);
1944 giraffe_model_delete(psfmodel);
1950 _profiles = cpl_matrix_get_data(profiles);
1953 signal = cpl_matrix_new(ny, 1);
1955 if (signal == NULL) {
1956 cpl_matrix_delete(profiles);
1959 cpl_matrix_delete(ypos);
1962 giraffe_model_delete(psfmodel);
1968 _signal = cpl_matrix_get_data(signal);
1971 variance = cpl_matrix_new(ny, 1);
1973 if (variance == NULL) {
1974 cpl_matrix_delete(signal);
1977 cpl_matrix_delete(profiles);
1980 cpl_matrix_delete(ypos);
1983 giraffe_model_delete(psfmodel);
1989 _variance = cpl_matrix_get_data(variance);
1992 mask = cpl_matrix_new(ny, 1);
1995 cpl_matrix_delete(variance);
1998 cpl_matrix_delete(signal);
2001 cpl_matrix_delete(profiles);
2004 cpl_matrix_delete(ypos);
2007 giraffe_model_delete(psfmodel);
2013 _mask = cpl_matrix_get_data(mask);
2016 weights = cpl_matrix_new(ny, 1);
2018 if (weights == NULL) {
2019 cpl_matrix_delete(mask);
2022 cpl_matrix_delete(variance);
2025 cpl_matrix_delete(signal);
2028 cpl_matrix_delete(profiles);
2031 cpl_matrix_delete(ypos);
2034 giraffe_model_delete(psfmodel);
2040 _weights = cpl_matrix_get_data(weights);
2048 bkg_base = giraffe_chebyshev_base1d(0., ny, bkg_nc, ypos);
2050 cpl_matrix_delete(ypos);
2053 if (bkg_base == NULL) {
2054 cpl_matrix_delete(weights);
2057 cpl_matrix_delete(mask);
2060 cpl_matrix_delete(variance);
2063 cpl_matrix_delete(signal);
2066 cpl_matrix_delete(profiles);
2069 cpl_matrix_delete(ypos);
2072 giraffe_model_delete(psfmodel);
2078 _bkg_base = cpl_matrix_get_data(bkg_base);
2080 for (i = 0; i < bkg_nc; ++i) {
2082 register cxint j = 0;
2083 register cxint offset = nfibers * ny;
2085 for (j = 0; j < ny; ++j) {
2086 _profiles[i * ny + j + offset] = _bkg_base[i * ny + j];
2093 cpl_matrix_delete(bkg_base);
2101 slice = _giraffe_extractionslice_new(nfibers, ny, bkg_nc);
2103 if (slice == NULL) {
2104 cpl_matrix_delete(weights);
2107 cpl_matrix_delete(mask);
2110 cpl_matrix_delete(variance);
2113 cpl_matrix_delete(signal);
2116 cpl_matrix_delete(profiles);
2119 cpl_matrix_delete(ypos);
2122 giraffe_model_delete(psfmodel);
2129 limits = _giraffe_extraction_psflimits_new(nfibers + bkg_nc);
2131 if (limits == NULL) {
2133 _giraffe_extractionslice_delete(slice);
2136 cpl_matrix_delete(weights);
2139 cpl_matrix_delete(mask);
2142 cpl_matrix_delete(variance);
2145 cpl_matrix_delete(signal);
2148 cpl_matrix_delete(profiles);
2151 cpl_matrix_delete(ypos);
2154 giraffe_model_delete(psfmodel);
2161 for (i = 0; i < limits->size; ++i) {
2162 limits->ymin[i] = 0;
2163 limits->ymax[i] = ny;
2171 workspace = _giraffe_optimal_workspace_new(nfibers + bkg_nc, ny);
2173 for (bin = 0; bin < nbins; ++bin) {
2175 cxbool stop = FALSE;
2181 const cxdouble* _my = cpl_image_get_data_double_const(my);
2182 const cxdouble* _mz = cpl_image_get_data_double_const(mz);
2183 const cxdouble* _mzvar = cpl_image_get_data_double_const(mzvar);
2185 cxdouble* _ms = cpl_image_get_data_double(ms);
2186 cxdouble* _mse = cpl_image_get_data_double(mse);
2187 cxdouble* _msy = cpl_image_get_data_double(msy);
2188 cxdouble* _msm = cpl_image_get_data_double(msm);
2192 GiExtractionPsfLimits* _limits = (nolimits == FALSE) ? limits : NULL;
2194 cx_assert(_mz != NULL);
2195 cx_assert(_mzvar != NULL);
2203 status = _giraffe_optimal_build_profiles(profiles, _limits, my, mw,
2204 fibers, bin, psfmodel, width,
2208 _giraffe_optimal_workspace_delete(workspace);
2211 _giraffe_extraction_psflimits_delete(limits);
2214 _giraffe_extractionslice_delete(slice);
2217 cpl_matrix_delete(weights);
2220 cpl_matrix_delete(mask);
2223 cpl_matrix_delete(variance);
2226 cpl_matrix_delete(signal);
2229 cpl_matrix_delete(profiles);
2232 cpl_matrix_delete(ypos);
2235 giraffe_model_delete(psfmodel);
2249 const cxint* _mbpx = cpl_image_get_data_int_const(mbpx);
2252 cx_assert(_mbpx != NULL);
2254 for (i = 0; i < ny; ++i) {
2256 cxbool bad = (_mbpx[bin * ny + i] & GIR_M_PIX_SET) ||
2257 (_mz[bin * ny + i] < 0.);
2259 _signal[i] = _mz[bin * ny + i];
2260 _variance[i] = _signal[i] + _mzvar[bin * ny + i];
2268 _weights[i] = _mask[i] / _variance[i];
2275 for (i = 0; i < ny; ++i) {
2277 cxbool bad = (_mz[bin * ny + i] < 0.);
2279 _signal[i] = _mz[bin * ny + i];
2280 _variance[i] = _signal[i] + _mzvar[bin * ny + i];
2288 _weights[i] = _mask[i] / _variance[i];
2300 nmin = (cxint)(fraction * ngood);
2302 while ((iter < niter) && (stop == FALSE)) {
2306 const cxdouble* _model = NULL;
2309 status = _giraffe_optimal_extract_slice(slice, profiles,
2310 signal, weights, limits, workspace);
2313 _giraffe_optimal_workspace_delete(workspace);
2316 _giraffe_extraction_psflimits_delete(limits);
2319 _giraffe_extractionslice_delete(slice);
2322 cpl_matrix_delete(weights);
2325 cpl_matrix_delete(mask);
2328 cpl_matrix_delete(variance);
2331 cpl_matrix_delete(signal);
2334 cpl_matrix_delete(profiles);
2337 cpl_matrix_delete(ypos);
2340 giraffe_model_delete(psfmodel);
2351 _model = cpl_matrix_get_data(slice->model);
2353 for (i = 0; i < ny; ++i) {
2355 if (_mask[i] > 0.) {
2358 cxdouble residual = _signal[i] - _model[i];
2361 _variance[i] = _model[i] + _mzvar[bin * ny + i];
2363 bad = (residual * residual) > (sigma * _variance[i]) ?
2372 _weights[i] = _mask[i] / _variance[i];
2378 if ((nreject == 0) || (ngood <= nmin)) {
2392 memcpy(&_ms[bin * nfibers], cpl_matrix_get_data(slice->flux),
2393 slice->nflx *
sizeof(cxdouble));
2394 memcpy(&_mse[bin * nfibers], cpl_matrix_get_data(slice->variance),
2395 slice->nflx *
sizeof(cxdouble));
2396 memcpy(&_msm[bin * ny], cpl_matrix_get_data(slice->model),
2397 slice->msize *
sizeof(cxdouble));
2399 memcpy(&_msy[bin * nfibers], &_my[bin * nfibers],
2400 nfibers *
sizeof(cxdouble));
2407 cpl_matrix_fill_window(profiles, 0., 0, 0, nfibers, ny);
2416 cpl_image_power(mse, 0.5);
2418 _giraffe_optimal_workspace_delete(workspace);
2421 _giraffe_extraction_psflimits_delete(limits);
2424 _giraffe_extractionslice_delete(slice);
2427 cpl_matrix_delete(weights);
2430 cpl_matrix_delete(mask);
2433 cpl_matrix_delete(variance);
2436 cpl_matrix_delete(signal);
2439 cpl_matrix_delete(profiles);
2442 giraffe_model_delete(psfmodel);
2476 GiTable* fibers, GiLocalization* sloc,
2477 GiImage* bpixel, GiImage* slight,
2478 GiExtractConfig* config)
2481 const cxchar *fctid =
"giraffe_extract_spectra";
2490 cxdouble bias_ron = 0.;
2491 cxdouble bias_sigma = 0.;
2492 cxdouble dark_value = 0.;
2493 cxdouble exptime = 0.;
2494 cxdouble conad = 1.;
2496 cpl_propertylist *properties;
2498 cpl_image* _image = NULL;
2499 cpl_image* _locy = NULL;
2500 cpl_image* _locw = NULL;
2501 cpl_image* _spectra = NULL;
2502 cpl_image* _error = NULL;
2503 cpl_image* _npixels = NULL;
2504 cpl_image* _centroid = NULL;
2505 cpl_image* _model = NULL;
2507 cpl_table* _fibers = NULL;
2514 if (!result || !image || !fibers || !sloc || !config) {
2515 cpl_error_set(fctid, CPL_ERROR_NULL_INPUT);
2520 if ((sloc->locy == NULL) || (sloc->locw == NULL)) {
2521 cpl_error_set(fctid, CPL_ERROR_NULL_INPUT);
2526 if (result->spectra != NULL || result->error != NULL ||
2527 result->npixels != NULL || result->centroid != NULL ||
2528 result->model != NULL) {
2529 gi_warning(
"%s: Results structure at %p is not empty! Contents "
2530 "might be lost.", fctid, result);
2536 if (_fibers == NULL) {
2537 cpl_error_set(fctid, CPL_ERROR_DATA_NOT_FOUND);
2542 if ((config->emethod == GIEXTRACT_OPTIMAL) && (sloc->psf == NULL)) {
2543 cpl_msg_error(fctid,
"Missing data: PSF profile data is required "
2544 "for optimal spectrum extraction!");
2545 cpl_error_set(fctid, CPL_ERROR_NULL_INPUT);
2553 if (properties == NULL) {
2554 cpl_error_set(fctid, CPL_ERROR_DATA_NOT_FOUND);
2559 giraffe_error_push();
2563 if (cpl_error_get_code() != CPL_ERROR_NONE) {
2567 giraffe_error_pop();
2570 if (!cpl_propertylist_has(properties, GIALIAS_BIASERROR)) {
2571 cpl_msg_warning(fctid,
"Missing bias error property (%s)! Setting "
2572 "bias error to 0.", GIALIAS_BIASERROR);
2576 bias_sigma = cpl_propertylist_get_double(properties, GIALIAS_BIASERROR);
2580 if (config->ron > 0.) {
2582 cpl_msg_info(fctid,
"Setting bias RMS property (%s) to %.4g ADU",
2583 GIALIAS_BIASSIGMA, config->ron);
2585 cpl_propertylist_update_double(properties, GIALIAS_BIASSIGMA,
2592 if (!cpl_propertylist_has(properties, GIALIAS_DARKVALUE)) {
2596 cpl_msg_warning(fctid,
"Missing dark value property (%s), will be "
2597 "set to 0.!", GIALIAS_DARKVALUE);
2598 cpl_propertylist_append_double(properties, GIALIAS_DARKVALUE,
2603 dark_value = cpl_propertylist_get_double(properties,
2608 if (!cpl_propertylist_has(properties, GIALIAS_EXPTIME)) {
2609 cpl_msg_error(fctid,
"Missing exposure time property (%s)!",
2614 exptime = cpl_propertylist_get_double(properties, GIALIAS_EXPTIME);
2618 if (cpl_propertylist_has(properties, GIALIAS_DATANCOM)) {
2619 nframes = cpl_propertylist_get_int(properties, GIALIAS_DATANCOM);
2631 bias_sigma *= conad;
2632 dark_value *= conad;
2652 ny = cpl_image_get_size_x(_image);
2653 nx = cpl_image_get_size_y(_locw);
2654 ns = cpl_table_get_nrow(_fibers);
2657 switch (config->emethod) {
2661 cxint xsize = cpl_image_get_size_x(_image);
2662 cxint ysize = cpl_image_get_size_y(_image);
2664 cxdouble ron_variance = bias_ron * bias_ron;
2665 cxdouble bias_variance = bias_sigma * bias_sigma;
2666 cxdouble dark_variance = dark_value * exptime;
2668 cpl_image* bpixmap = NULL;
2669 cpl_image* variance = NULL;
2676 result->model = NULL;
2683 if (bpixel != NULL) {
2687 if (cpl_image_get_size_x(bpixmap) != xsize ||
2688 cpl_image_get_size_y(bpixmap) != ysize) {
2690 cxbool crop = FALSE;
2692 cpl_propertylist *p =
2695 GiWindow w = {1, 1, 0, 0};
2698 w.x1 = cpl_image_get_size_x(bpixmap);
2699 w.y1 = cpl_image_get_size_y(bpixmap);
2701 if (cpl_propertylist_has(p, GIALIAS_PRSCX)) {
2702 w.x0 += cpl_propertylist_get_int(p, GIALIAS_PRSCX);
2706 if (cpl_propertylist_has(p, GIALIAS_OVSCX)) {
2707 w.x1 -= cpl_propertylist_get_int(p, GIALIAS_OVSCX);
2711 if (cpl_propertylist_has(p, GIALIAS_PRSCY)) {
2712 w.y0 += cpl_propertylist_get_int(p, GIALIAS_PRSCY);
2716 if (cpl_propertylist_has(p, GIALIAS_OVSCY)) {
2717 w.y1 -= cpl_propertylist_get_int(p, GIALIAS_OVSCY);
2721 if ((w.x1 - w.x0 + 1) != xsize ||
2722 (w.y1 - w.y0 + 1) != ysize) {
2723 cpl_msg_error(fctid,
"Invalid bad pixel map! Image "
2724 "sizes do not match!");
2727 result->spectra = NULL;
2730 result->error = NULL;
2733 result->npixels = NULL;
2736 result->centroid = NULL;
2739 result->model = NULL;
2741 cpl_image_delete(_image);
2748 bpixmap = cpl_image_extract(bpixmap, w.x0, w.y0,
2756 if (slight != NULL) {
2757 cpl_msg_warning(fctid,
"Scattered light model will be "
2758 "ignored for extraction method `SUM'");
2761 variance = cpl_image_abs_create(_image);
2769 cpl_image_add_scalar(variance, nframes * (ron_variance + nframes *
2770 (bias_variance + dark_variance)));
2772 status = _giraffe_extract_summation(_image, variance, _fibers,
2773 _locy, _locw, bpixmap,
2774 _spectra, _error, _npixels,
2777 cpl_image_delete(variance);
2779 cpl_image_delete(bpixmap);
2787 case GIEXTRACT_OPTIMAL:
2790 cxint xsize = cpl_image_get_size_x(_image);
2791 cxint ysize = cpl_image_get_size_y(_image);
2794 cxdouble ron_variance = bias_ron * bias_ron;
2795 cxdouble bias_variance = bias_sigma * bias_sigma;
2796 cxdouble dark_variance = dark_value * exptime;
2798 cpl_image* variance = NULL;
2799 cpl_image* bpixmap = NULL;
2801 GiExtractOptimalConfig setup;
2806 result->npixels = NULL;
2815 setup.clip.iterations = config->psf.iterations;
2816 setup.clip.level = config->psf.sigma;
2817 setup.clip.fraction = config->optimal.fraction;
2818 setup.limits = config->optimal.wfactor < 0. ? FALSE : TRUE;
2819 setup.ewidth = CX_MAX(1., fabs(config->optimal.wfactor));
2820 setup.bkgorder = config->optimal.bkgorder;
2821 setup.exptime = exptime;
2822 setup.ron = bias_sigma;
2823 setup.dark = dark_value;
2826 if (bpixel != NULL) {
2830 if (cpl_image_get_size_x(bpixmap) != xsize ||
2831 cpl_image_get_size_y(bpixmap) != ysize) {
2833 cxbool crop = FALSE;
2835 cpl_propertylist *p =
2838 GiWindow w = {1, 1, 0, 0};
2841 w.x1 = cpl_image_get_size_x(bpixmap);
2842 w.y1 = cpl_image_get_size_y(bpixmap);
2844 if (cpl_propertylist_has(p, GIALIAS_PRSCX)) {
2845 w.x0 += cpl_propertylist_get_int(p, GIALIAS_PRSCX);
2849 if (cpl_propertylist_has(p, GIALIAS_OVSCX)) {
2850 w.x1 -= cpl_propertylist_get_int(p, GIALIAS_OVSCX);
2854 if (cpl_propertylist_has(p, GIALIAS_PRSCY)) {
2855 w.y0 += cpl_propertylist_get_int(p, GIALIAS_PRSCY);
2859 if (cpl_propertylist_has(p, GIALIAS_OVSCY)) {
2860 w.y1 -= cpl_propertylist_get_int(p, GIALIAS_OVSCY);
2864 if ((w.x1 - w.x0 + 1) != xsize ||
2865 (w.y1 - w.y0 + 1) != ysize) {
2867 cpl_msg_error(fctid,
"Invalid bad pixel map! "
2868 "Image sizes do not match!");
2871 result->spectra = NULL;
2874 result->error = NULL;
2877 result->npixels = NULL;
2880 result->centroid = NULL;
2883 result->model = NULL;
2885 cpl_image_delete(_image);
2893 bpixmap = cpl_image_extract(bpixmap, w.x0, w.y0,
2901 variance = cpl_image_new(xsize, ysize, CPL_TYPE_DOUBLE);
2909 v0 = nframes * (ron_variance + nframes *
2910 (bias_variance + dark_variance));
2919 if (slight != NULL) {
2921 register cxsize i = 0;
2922 register cxsize npixels = xsize * ysize;
2924 const cxdouble* _slight =
2927 cxdouble* _variance = cpl_image_get_data_double(variance);
2929 for (i = 0; i < npixels; i++) {
2930 _variance[i] = v0 + fabs(_slight[i]) * conad * nframes;
2936 register cxsize i = 0;
2937 register cxsize npixels = xsize * ysize;
2939 cxdouble* _variance = cpl_image_get_data_double(variance);
2941 for (i = 0; i < npixels; i++) {
2948 status = _giraffe_extract_optimal(_image, variance, _fibers,
2949 _locy, _locw, sloc->psf,
2950 bpixmap, _spectra, _error,
2951 _model, _centroid, &setup);
2953 cpl_image_delete(variance);
2957 cpl_image_delete(bpixmap);
2965 case GIEXTRACT_HORNE:
2968 cxint xsize = cpl_image_get_size_x(_image);
2969 cxint ysize = cpl_image_get_size_y(_image);
2972 cxdouble ron_variance = bias_ron * bias_ron;
2973 cxdouble bias_variance = bias_sigma * bias_sigma;
2974 cxdouble dark_variance = dark_value * exptime;
2976 cpl_image* variance = NULL;
2977 cpl_image* bpixmap = NULL;
2979 GiExtractHorneConfig setup;
2986 result->model = NULL;
2993 setup.clip.iterations = config->psf.iterations;
2994 setup.clip.level = config->psf.sigma;
2995 setup.clip.fraction = config->horne.mingood;
2996 setup.ewidth = config->horne.ewidth;
2997 setup.exptime = exptime;
2998 setup.ron = bias_sigma;
2999 setup.dark = dark_value;
3001 if (bpixel != NULL) {
3005 if (cpl_image_get_size_x(bpixmap) != xsize ||
3006 cpl_image_get_size_y(bpixmap) != ysize) {
3008 cxbool crop = FALSE;
3010 cpl_propertylist *p =
3013 GiWindow w = {1, 1, 0, 0};
3016 w.x1 = cpl_image_get_size_x(bpixmap);
3017 w.y1 = cpl_image_get_size_y(bpixmap);
3019 if (cpl_propertylist_has(p, GIALIAS_PRSCX)) {
3020 w.x0 += cpl_propertylist_get_int(p, GIALIAS_PRSCX);
3024 if (cpl_propertylist_has(p, GIALIAS_OVSCX)) {
3025 w.x1 -= cpl_propertylist_get_int(p, GIALIAS_OVSCX);
3029 if (cpl_propertylist_has(p, GIALIAS_PRSCY)) {
3030 w.y0 += cpl_propertylist_get_int(p, GIALIAS_PRSCY);
3034 if (cpl_propertylist_has(p, GIALIAS_OVSCY)) {
3035 w.y1 -= cpl_propertylist_get_int(p, GIALIAS_OVSCY);
3039 if ((w.x1 - w.x0 + 1) != xsize ||
3040 (w.y1 - w.y0 + 1) != ysize) {
3042 cpl_msg_error(fctid,
"Invalid bad pixel map! "
3043 "Image sizes do not match!");
3046 result->spectra = NULL;
3049 result->error = NULL;
3052 result->npixels = NULL;
3055 result->centroid = NULL;
3058 result->model = NULL;
3060 cpl_image_delete(_image);
3068 bpixmap = cpl_image_extract(bpixmap, w.x0, w.y0,
3076 variance = cpl_image_new(xsize, ysize, CPL_TYPE_DOUBLE);
3084 v0 = nframes * (ron_variance + nframes *
3085 (bias_variance + dark_variance));
3095 if (slight != NULL) {
3097 register cxsize i = 0;
3098 register cxsize npixels = xsize * ysize;
3100 const cxdouble* _slight =
3103 cxdouble* _variance = cpl_image_get_data_double(variance);
3105 for (i = 0; i < npixels; i++) {
3106 _variance[i] = v0 + fabs(_slight[i]) * nframes * conad;
3112 register cxsize i = 0;
3113 register cxsize npixels = xsize * ysize;
3115 cxdouble* _variance = cpl_image_get_data_double(variance);
3117 for (i = 0; i < npixels; i++) {
3124 status = _giraffe_extract_horne(_image, variance, _fibers,
3125 _locy, _locw, sloc->psf,
3126 bpixmap, _spectra, _error,
3127 _npixels, _centroid, &setup);
3129 cpl_image_delete(variance);
3133 cpl_image_delete(bpixmap);
3142 gi_message(
"%s: Method %d selected for spectrum extraction.",
3143 fctid, config->emethod);
3144 cpl_msg_error(fctid,
"Invalid extraction method!");
3150 cpl_image_delete(_image);
3156 result->spectra = NULL;
3159 result->error = NULL;
3162 result->npixels = NULL;
3165 result->centroid = NULL;
3168 result->model = NULL;
3170 cpl_msg_error(fctid,
"Spectrum extraction (method %d) failed!",
3173 cpl_image_delete(_image);
3191 if (result->spectra) {
3196 if (result->model) {
3201 if (result->error) {
3226 cpl_propertylist_set_int(properties, GIALIAS_NAXIS1,
3227 cpl_image_get_size_x(_spectra));
3228 cpl_propertylist_set_int(properties, GIALIAS_NAXIS2,
3229 cpl_image_get_size_y(_spectra));
3231 cpl_propertylist_set_int(properties, GIALIAS_BITPIX, -32);
3232 cpl_propertylist_set_double(properties, GIALIAS_BZERO, 0.);
3233 cpl_propertylist_set_double(properties, GIALIAS_BSCALE, 1.);
3235 cpl_propertylist_update_int(properties, GIALIAS_NFIBERS,
3236 cpl_image_get_size_x(_spectra));
3238 cpl_propertylist_append_int(properties, GIALIAS_EXT_NX,
3239 cpl_image_get_size_y(_spectra));
3240 cpl_propertylist_append_int(properties, GIALIAS_EXT_NS,
3241 cpl_image_get_size_x(_spectra));
3243 switch (config->emethod) {
3245 cpl_propertylist_append_string(properties, GIALIAS_EXT_METHOD,
3247 cpl_propertylist_set_comment(properties, GIALIAS_EXT_METHOD,
3248 "Spectrum extraction method");
3251 case GIEXTRACT_HORNE:
3254 cpl_propertylist_append_string(properties, GIALIAS_EXT_METHOD,
3256 cpl_propertylist_set_comment(properties, GIALIAS_EXT_METHOD,
3257 "Spectrum extraction method");
3259 cpl_propertylist_append_string(properties, GIALIAS_EXTPSF_MODEL,
3261 cpl_propertylist_set_comment(properties, GIALIAS_EXTPSF_MODEL,
3263 cpl_propertylist_append_double(properties, GIALIAS_EXTPSF_SIGMA,
3265 cpl_propertylist_set_comment(properties, GIALIAS_EXTPSF_SIGMA,
3266 "PSF fit sigma clipping threshold");
3267 cpl_propertylist_append_int(properties, GIALIAS_EXTPSF_NITER,
3268 config->psf.iterations);
3269 cpl_propertylist_set_comment(properties, GIALIAS_EXTPSF_NITER,
3270 "PSF fit maximum number of "
3273 cpl_propertylist_append_int(properties, GIALIAS_EXTHRN_EWIDTH,
3274 config->horne.ewidth);
3275 cpl_propertylist_set_comment(properties, GIALIAS_EXTHRN_EWIDTH,
3276 "Number of extra pixels used");
3277 cpl_propertylist_append_int(properties, GIALIAS_EXTHRN_MINGOOD,
3278 config->horne.mingood);
3279 cpl_propertylist_set_comment(properties, GIALIAS_EXTHRN_MINGOOD,
3280 "Minimum number of pixels to keep");
3286 case GIEXTRACT_OPTIMAL:
3287 cpl_propertylist_append_string(properties, GIALIAS_EXT_METHOD,
3289 cpl_propertylist_set_comment(properties, GIALIAS_EXT_METHOD,
3290 "Spectrum extraction method");
3292 cpl_propertylist_append_string(properties, GIALIAS_EXTPSF_MODEL,
3294 cpl_propertylist_set_comment(properties, GIALIAS_EXTPSF_MODEL,
3296 cpl_propertylist_append_double(properties, GIALIAS_EXTPSF_SIGMA,
3298 cpl_propertylist_set_comment(properties, GIALIAS_EXTPSF_SIGMA,
3299 "PSF fit sigma clipping threshold");
3300 cpl_propertylist_append_int(properties, GIALIAS_EXTPSF_NITER,
3301 config->psf.iterations);
3302 cpl_propertylist_set_comment(properties, GIALIAS_EXTPSF_NITER,
3303 "PSF fit maximum number of "
3306 cpl_propertylist_append_double(properties, GIALIAS_EXTOPT_FRACTION,
3307 config->optimal.fraction);
3308 cpl_propertylist_set_comment(properties, GIALIAS_EXTOPT_FRACTION,
3309 "Minimum fraction of pixels used.");
3310 cpl_propertylist_append_double(properties, GIALIAS_EXTOPT_WFACTOR,
3311 config->optimal.wfactor);
3312 cpl_propertylist_set_comment(properties, GIALIAS_EXTOPT_WFACTOR,
3313 "Multiple of the fiber PSF half "
3314 "width used for spectrum "
3316 cpl_propertylist_append_int(properties, GIALIAS_EXTOPT_BGORDER,
3317 config->optimal.bkgorder);
3318 cpl_propertylist_set_comment(properties, GIALIAS_EXTOPT_BGORDER,
3319 "Order of the background polynomial "
3320 "model along the spatial direction.");
3328 cpl_propertylist_update_string(properties, GIALIAS_GIRFTYPE,
"EXTSP");
3329 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
3330 "Extracted spectra");
3340 cpl_propertylist_set_string(properties, GIALIAS_GIRFTYPE,
"EXTERRS");
3341 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
3342 "Extracted spectra errors");
3352 cpl_propertylist_set_string(properties, GIALIAS_GIRFTYPE,
"EXTYCEN");
3353 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
3354 "Extracted spectra centroids");
3361 if (result->npixels != NULL) {
3365 cpl_propertylist_set_string(properties, GIALIAS_GIRFTYPE,
"EXTNPIX");
3366 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
3367 "Extracted spectra npixels");
3375 if (result->model != NULL) {
3379 cpl_propertylist_set_string(properties, GIALIAS_GIRFTYPE,
"EXTMODEL");
3380 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
3381 "Model spectra used for extraction");
3406 GiExtractConfig* config = NULL;
3413 config = cx_calloc(1,
sizeof *config);
3415 p = cpl_parameterlist_find(list,
"giraffe.extraction.method");
3416 s = cpl_parameter_get_string(p);
3417 if (!strcmp(s,
"OPTIMAL")) {
3418 config->emethod = GIEXTRACT_OPTIMAL;
3420 else if (!strcmp(s,
"HORNE")) {
3421 config->emethod = GIEXTRACT_HORNE;
3424 config->emethod = GIEXTRACT_SUM;
3427 p = cpl_parameterlist_find(list,
"giraffe.extraction.ron");
3428 config->ron = cpl_parameter_get_double(p);
3430 p = cpl_parameterlist_find(list,
"giraffe.extraction.psf.model");
3431 config->psf.model = cx_strdup(cpl_parameter_get_string(p));
3433 p = cpl_parameterlist_find(list,
"giraffe.extraction.psf.sigma");
3434 config->psf.sigma = cpl_parameter_get_double(p);
3436 p = cpl_parameterlist_find(list,
"giraffe.extraction.psf.iterations");
3437 config->psf.iterations = cpl_parameter_get_int(p);
3440 p = cpl_parameterlist_find(list,
"giraffe.extraction.horne.extrawidth");
3441 config->horne.ewidth = cpl_parameter_get_int(p);
3443 p = cpl_parameterlist_find(list,
"giraffe.extraction.horne.mingood");
3444 config->horne.mingood = cpl_parameter_get_double(p);
3447 p = cpl_parameterlist_find(list,
"giraffe.extraction.optimal.fraction");
3448 config->optimal.fraction = cpl_parameter_get_double(p);
3450 p = cpl_parameterlist_find(list,
"giraffe.extraction.optimal.wfactor");
3451 config->optimal.wfactor = cpl_parameter_get_double(p);
3453 p = cpl_parameterlist_find(list,
"giraffe.extraction.optimal.bkgorder");
3454 config->optimal.bkgorder = cpl_parameter_get_int(p);
3479 if (config->psf.model) {
3480 cx_free(config->psf.model);
3507 cpl_parameter* p = NULL;
3514 p = cpl_parameter_new_enum(
"giraffe.extraction.method",
3516 "Extraction method: 'SUM', 'HORNE' or "
3518 "giraffe.extraction",
3519 "SUM", 3,
"SUM",
"OPTIMAL",
"HORNE");
3520 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"extr-method");
3521 cpl_parameterlist_append(list, p);
3524 p = cpl_parameter_new_value(
"giraffe.extraction.ron",
3526 "New bias sigma (RON) value for "
3529 "giraffe.extraction",
3531 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"extr-ron");
3532 cpl_parameterlist_append(list, p);
3535 p = cpl_parameter_new_enum(
"giraffe.extraction.psf.model",
3537 "PSF profile model: `psfexp', `psfexp2'",
3538 "giraffe.extraction.psf",
3539 "psfexp2", 2,
"psfexp",
"psfexp2");
3540 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"extr-psfmodel");
3541 cpl_parameterlist_append(list, p);
3544 p = cpl_parameter_new_value(
"giraffe.extraction.psf.sigma",
3546 "Sigma clippging threshold used for "
3547 "rejecting data points during PSF fitting "
3548 "(Horne's sigma). It is used to reject bad "
3549 "pixels and cosmics.",
3550 "giraffe.extraction.psf",
3552 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"extr-psfsigma");
3553 cpl_parameterlist_append(list, p);
3556 p = cpl_parameter_new_value(
"giraffe.extraction.psf.iterations",
3558 "Maximum number of iterations used for "
3559 "fitting the PSF profile.",
3560 "giraffe.extraction.psf",
3562 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"extr-psfniter");
3563 cpl_parameterlist_append(list, p);
3566 p = cpl_parameter_new_value(
"giraffe.extraction.horne.extrawidth",
3568 "Horne extraction method: Number of "
3569 "extra pixels added to the fiber "
3571 "giraffe.extraction.horne",
3573 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"extr-hewidth");
3574 cpl_parameterlist_append(list, p);
3577 p = cpl_parameter_new_value(
"giraffe.extraction.horne.mingood",
3579 "Horne extraction method: Minimum number of "
3580 "points used for the profile fit. It sets "
3581 "the lower limit of data points for the "
3583 "giraffe.extraction.horne",
3585 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"extr-hmingood");
3586 cpl_parameterlist_append(list, p);
3589 p = cpl_parameter_new_range(
"giraffe.extraction.optimal.fraction",
3591 "Optimal extraction method: Minimum fraction "
3592 "of the data points used for fitting the "
3593 "fiber profiles. It sets the lower limit "
3594 "for the pixel rejection.",
3595 "giraffe.extraction.optimal",
3597 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"extr-omfrac");
3598 cpl_parameterlist_append(list, p);
3601 p = cpl_parameter_new_value(
"giraffe.extraction.optimal.wfactor",
3603 "Optimal extraction method: Factor by which "
3604 "the fiber PSF half width is multiplied. "
3605 "Adjacent spectra within this area are "
3606 "assumed to affect the spectrum being "
3608 "giraffe.extraction.optimal",
3610 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"extr-owfactor");
3611 cpl_parameterlist_append(list, p);
3614 p = cpl_parameter_new_value(
"giraffe.extraction.optimal.bkgorder",
3616 "Optimal extraction method: Order of the "
3617 "polynomial background model, which is "
3618 "fitted for each wavelength bin along the "
3619 "spatial direction.",
3620 "giraffe.extraction.optimal",
3622 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"extr-obkgorder");
3623 cpl_parameterlist_append(list, p);
cxint giraffe_array_sort(cxdouble *array, cxsize size)
Sorts an array in ascending order.
const cxchar * giraffe_fiberlist_query_index(const cpl_table *fibers)
Query a fiber list for the name of the fiber reference index column.
void giraffe_image_delete(GiImage *self)
Destroys an image.
cpl_propertylist * giraffe_image_get_properties(const GiImage *self)
Get the properties of 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.
void gi_warning(const cxchar *format,...)
Log a warning.
void gi_message(const cxchar *format,...)
Log a normal message.
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.
cxint giraffe_propertylist_update(cpl_propertylist *self, cpl_propertylist *properties, const cxchar *regexp)
Update a property list.
cxdouble giraffe_propertylist_get_ron(const cpl_propertylist *properties)
Retrieve the read-out noise from the given properties.