GIRAFFE Pipeline Reference Manual

giextract.c
1/*
2 * This file is part of the GIRAFFE Pipeline
3 * Copyright (C) 2002-2019 European Southern Observatory
4 *
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
18 */
19
20#ifdef HAVE_CONFIG_H
21# include <config.h>
22#endif
23
24#include <math.h>
25#include <float.h>
26
27#include <cxmemory.h>
28#include <cxstring.h>
29#include <cxstrutils.h>
30
31#include <cpl_parameterlist.h>
32#include <cpl_matrix.h>
33#include <cpl_table.h>
34#include <cpl_msg.h>
35
36#include "gimacros.h"
37#include "gierror.h"
38#include "gialias.h"
39#include "giclip.h"
40#include "giarray.h"
41#include "giimage.h"
42#include "gimatrix.h"
43#include "giwindow.h"
44#include "gipsfdata.h"
45#include "gimodel.h"
46#include "gimath.h"
47#include "gilocalization.h"
48#include "gimessages.h"
49#include "gifiberutils.h"
50#include "giutils.h"
51#include "giextract.h"
52
53
62enum GiProfileId {
63 PROFILE_PSFEXP = 1 << 1,
64 PROFILE_PSFEXP2 = 1 << 2,
65 PROFILE_GAUSSIAN = 1 << 3
66};
67
68typedef enum GiProfileId GiProfileId;
69
70
71/*
72 * Optimal spectrum extraction algorithm configuration data
73 */
74
75struct GiExtractOptimalConfig {
76
77 GiClipParams clip;
78
79 cxbool limits;
80
81 cxint bkgorder;
82
83 cxdouble exptime;
84 cxdouble ron;
85 cxdouble dark;
86 cxdouble ewidth;
87};
88
89typedef struct GiExtractOptimalConfig GiExtractOptimalConfig;
90
91
92/*
93 * Original Horne spectrum extraction algorithm configuration data
94 */
95
96struct GiExtractHorneConfig {
97 GiClipParams clip;
98
99 cxdouble exptime;
100 cxdouble ron;
101 cxdouble dark;
102 cxdouble ewidth;
103};
104
105typedef struct GiExtractHorneConfig GiExtractHorneConfig;
106
107
108struct GiExtractionData {
109 cxdouble value;
110 cxdouble error;
111 cxdouble position;
112 cxdouble npixels;
113};
114
115typedef struct GiExtractionData GiExtractionData;
116
117
118struct GiExtractionSlice {
119 cxint fsize;
120 cxint msize;
121
122 cxint nflx;
123 cxint nbkg;
124
125 cpl_matrix* flux;
126 cpl_matrix* variance;
127 cpl_matrix* model;
128};
129
130typedef struct GiExtractionSlice GiExtractionSlice;
131
132
133struct GiExtractionPsfLimits {
134 cxint size;
135
136 cxint* ymin;
137 cxint* ymax;
138};
139
140typedef struct GiExtractionPsfLimits GiExtractionPsfLimits;
141
142
143struct GiExtractionWorkspace {
144 cpl_matrix* atw;
145 cpl_matrix* atwa;
146 cpl_matrix* atws;
147 cpl_matrix* c;
148 cpl_matrix* tmp;
149};
150
151typedef struct GiExtractionWorkspace GiExtractionWorkspace;
152
153
154struct GiVirtualSlit {
155 cxint width;
156
157 cxdouble center;
158 cxdouble extra_width;
159
160 cxdouble* position;
161 cxdouble* signal;
162 cxdouble* variance;
163 cxdouble* fraction;
164
165 cxint* mask;
166 cxint* offset;
167};
168
169typedef struct GiVirtualSlit GiVirtualSlit;
170
171
172/*
173 * Extraction slice implementation
174 */
175
176inline static GiExtractionSlice*
177_giraffe_extractionslice_new(cxint nflx, cxint ndata, cxint nbkg)
178{
179
180 GiExtractionSlice* self = cx_malloc(sizeof *self);
181
182 self->nflx = nflx;
183 self->nbkg = nbkg;
184
185 self->fsize = nflx + nbkg;
186 self->msize = ndata;
187
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);
191
192 return self;
193
194}
195
196
197inline static void
198_giraffe_extractionslice_delete(GiExtractionSlice* self)
199{
200
201 if (self != NULL) {
202 if (self->model != NULL) {
203 cpl_matrix_delete(self->model);
204 self->model = NULL;
205 }
206
207 if (self->variance != NULL) {
208 cpl_matrix_delete(self->variance);
209 self->variance = NULL;
210 }
211
212 if (self->flux != NULL) {
213 cpl_matrix_delete(self->flux);
214 self->flux = NULL;
215 }
216
217 cx_free(self);
218 }
219
220 return;
221
222}
223
224
225inline static GiExtractionPsfLimits*
226_giraffe_extraction_psflimits_new(cxint size)
227{
228
229 GiExtractionPsfLimits* self = cx_malloc(sizeof *self);
230
231 self->size = size;
232
233 self->ymin = cx_calloc(self->size, sizeof(cxint));
234 self->ymax = cx_calloc(self->size, sizeof(cxint));
235
236 return self;
237
238}
239
240
241inline static void
242_giraffe_extraction_psflimits_delete(GiExtractionPsfLimits* self)
243{
244
245 if (self != NULL) {
246 if (self->ymin != NULL) {
247 cx_free(self->ymin);
248 }
249
250 if (self->ymax != NULL) {
251 cx_free(self->ymax);
252 }
253
254 cx_free(self);
255 }
256
257 return;
258
259}
260
261
262inline static GiExtractionWorkspace*
263_giraffe_optimal_workspace_new(cxint m, cxint n)
264{
265
266 GiExtractionWorkspace* self = cx_malloc(sizeof *self);
267
268
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);
273
274 self->tmp = cpl_matrix_new(m, m);
275
276 return self;
277
278}
279
280
281inline static void
282_giraffe_optimal_workspace_delete(GiExtractionWorkspace* self)
283{
284
285 if (self != NULL) {
286 if (self->atws != NULL) {
287 cpl_matrix_delete(self->atws);
288 }
289
290 if (self->atwa != NULL) {
291 cpl_matrix_delete(self->atwa);
292 }
293
294 if (self->c != NULL) {
295 cpl_matrix_delete(self->c);
296 }
297
298 if (self->atw != NULL) {
299 cpl_matrix_delete(self->atw);
300 }
301
302 if (self->tmp != NULL) {
303 cpl_matrix_delete(self->tmp);
304 }
305
306 cx_free(self);
307
308 }
309
310 return;
311
312}
313
314
315/*
316 * Virtual slit implementation
317 */
318
319inline static void
320_giraffe_virtualslit_allocate(GiVirtualSlit* self)
321{
322
323 if ((self != NULL) && (self->width > 0)) {
324
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));
329
330 self->mask = cx_calloc(self->width, sizeof(cxdouble));
331 self->offset = cx_calloc(self->width, sizeof(cxdouble));
332
333 }
334
335 return;
336
337}
338
339
340inline static GiVirtualSlit*
341_giraffe_virtualslit_new(cxdouble extra_width)
342{
343
344 GiVirtualSlit* self = cx_calloc(1, sizeof *self);
345
346 self->width = 0;
347 self->center = 0.;
348 self->extra_width = extra_width;
349
350 self->position = NULL;
351 self->signal = NULL;
352 self->variance = NULL;
353 self->fraction = NULL;
354 self->mask = NULL;
355 self->offset = NULL;
356
357 return self;
358
359}
360
361
362inline static void
363_giraffe_virtualslit_clear(GiVirtualSlit* self)
364{
365
366 if (self != NULL) {
367
368 if (self->position != NULL) {
369 cx_free(self->position);
370 self->position = NULL;
371 }
372
373 if (self->signal != NULL) {
374 cx_free(self->signal);
375 self->signal = NULL;
376 }
377
378 if (self->variance != NULL) {
379 cx_free(self->variance);
380 self->variance = NULL;
381 }
382
383 if (self->fraction != NULL) {
384 cx_free(self->fraction);
385 self->fraction = NULL;
386 }
387
388 if (self->mask != NULL) {
389 cx_free(self->mask);
390 self->mask = NULL;
391 }
392
393 if (self->offset != NULL) {
394 cx_free(self->offset);
395 self->offset = NULL;
396 }
397
398 self->extra_width = 0.;
399 self->center = 0.;
400 self->width = 0;
401
402 }
403
404 return;
405
406}
407
408
409inline static void
410_giraffe_virtualslit_delete(GiVirtualSlit* self)
411{
412
413 if (self != NULL) {
414 _giraffe_virtualslit_clear(self);
415
416 cx_free(self);
417 }
418
419 return;
420
421}
422
423
424inline static cxint
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)
429{
430
431 register cxint ny = cpl_image_get_size_x(signal);
432 register cxint offset = bin * cpl_image_get_size_x(signal);
433
434 register cxdouble lower = center - (width + self->extra_width);
435 register cxdouble upper = center + (width + self->extra_width);
436
437 register cxint first = (cxint) floor(lower);
438 register cxint last = (cxint) ceil(upper);
439
440 const cxdouble* s = cpl_image_get_data_double_const(signal);
441 const cxdouble* v = cpl_image_get_data_double_const(variance);
442
443
444 /*
445 * Upper, lower border and width of the virtual slit
446 */
447
448 lower = CX_MAX(0., lower);
449 upper = CX_MIN(ny, upper);
450
451 first = CX_MAX(0, first);
452 last = CX_MIN(ny, last);
453
454 self->center = center;
455 self->width = last - first + 1;
456
457
458 /*
459 * Create and fill the buffers
460 */
461
462 _giraffe_virtualslit_allocate(self);
463
464 if (bpixel != NULL) {
465
466 register cxint k = 0;
467 register cxint y = 0;
468
469 const cxint* _bpixel = cpl_image_get_data_int_const(bpixel);
470
471
472 for (y = first; y <= last; y++) {
473
474 register cxint ypos = offset + y;
475
476 cxint ok = (_bpixel[ypos] & GIR_M_PIX_SET) == 0 ? 1 : 0;
477
478
479 self->position[k] = y - center;
480 self->fraction[k] = 1.;
481
482 self->signal[k] = s[ypos];
483 self->variance[k] = v[ypos];
484
485 self->mask[k] = ok;
486 self->offset[k] = ypos;
487 ++k;
488
489 }
490
491 }
492 else {
493
494 register cxint k = 0;
495 register cxint y = 0;
496
497
498 for (y = first; y <= last; y++) {
499
500 register cxint ypos = offset + y;
501
502 cxint ok = 1;
503
504
505 self->position[k] = y - center;
506 self->fraction[k] = 1.;
507
508 self->signal[k] = s[ypos];
509 self->variance[k] = v[ypos];
510
511 self->mask[k] = ok;
512 self->offset[k] = ypos;
513 ++k;
514
515 }
516
517 }
518
519
520 /*
521 * Correct for pixel fractions at the borders of the
522 * virtual slit, since they have been set to the full
523 * pixel in the above loop.
524 */
525
526 self->fraction[0] = ((cxdouble)first + 1.) - lower;
527 self->fraction[self->width - 1] = upper - ((cxdouble)last - 1.);
528
529 return self->width;
530
531}
532
533
534/*
535 * Compute the inverse of a square matrix
536 */
537
538inline static cxint
539_giraffe_matrix_invert(cpl_matrix* m_inv, const cpl_matrix* m, cpl_matrix* lu)
540{
541
542 cxint i = 0;
543 cxint status = 0;
544 cxint n = cpl_matrix_get_ncol(m);
545
546 register cxint sz = n * n * sizeof(cxdouble);
547
548 const cxdouble* _m = cpl_matrix_get_data_const(m);
549
550 cxdouble* _m_inv = cpl_matrix_get_data(m_inv);
551 cxdouble* _m_lu = cpl_matrix_get_data(lu);
552
553 cpl_array* perm = cpl_array_new(n, CPL_TYPE_INT);
554
555 register cxint* perm_data = cpl_array_get_data_int(perm);
556
557
558 memset(_m_inv, 0, sz);
559 memcpy(_m_lu, _m, sz);
560
561 if (cpl_matrix_decomp_lu(lu, perm, &i) != 0) {
562 cpl_array_delete(perm);
563 return 1;
564 }
565
566
567 /*
568 * Create an identity matrix with the rows permuted
569 */
570
571 for (i = 0; i < n; ++i) {
572 _m_inv[i * n + perm_data[i]] = 1.;
573 }
574
575 cpl_array_delete(perm);
576
577
578 status = cpl_matrix_solve_lu(lu, m_inv, NULL);
579
580 if (status != 0) {
581 cpl_matrix_delete(m_inv);
582 return 2;
583 }
584
585 return 0;
586
587}
588
589
590/*
591 * Compute the PSF profile for a set of abscissa values.
592 */
593
594inline static cpl_matrix*
595_giraffe_compute_psf(GiModel* psf, const cpl_matrix* x)
596{
597
598 register cxint i = 0;
599 register cxint n = 0;
600
601 cxint status = 0;
602
603 const cxdouble* _x = NULL;
604
605 cxdouble* _y = NULL;
606
607 cpl_matrix* y = NULL;
608
609 cx_assert(psf != NULL);
610 cx_assert(x != NULL);
611 cx_assert(cpl_matrix_get_ncol(x) == 1);
612
613 n = cpl_matrix_get_nrow(x);
614
615 y = cpl_matrix_new(n, 1);
616
617 _x = cpl_matrix_get_data_const(x);
618 _y = cpl_matrix_get_data(y);
619
620 for (i = 0; i < n; i++) {
621 giraffe_model_set_argument(psf, "x", _x[i]);
622 giraffe_model_evaluate(psf, &_y[i], &status);
623
624 if (status != 0) {
625 cpl_matrix_delete(y);
626 return NULL;
627 }
628 }
629
630 return y;
631
632}
633
634
635/*
636 * Horne extraction of a single wavelength bin for the given virtual
637 * slit.
638 */
639
640inline static cxint
641_giraffe_horne_extract_slit(GiExtractionData* result,
642 const GiVirtualSlit* vslit, GiModel* psf,
643 const GiExtractHorneConfig* config)
644{
645
646 cxint i = 0;
647 cxint ngood = 0;
648
649 cxdouble var = 0.;
650 cxdouble bkg = 0.;
651 cxdouble flx = 0.;
652 cxdouble norm = 0.;
653 cxdouble* tdata = NULL;
654 cxdouble* _mnpsf = NULL;
655
656 cpl_matrix* mnpsf = NULL;
657 cpl_matrix* mvslit = NULL;
658
659
660
661 /*
662 * Compute the PSF model.
663 */
664
665 mvslit = cpl_matrix_wrap(vslit->width, 1, vslit->position);
666 mnpsf = _giraffe_compute_psf(psf, mvslit);
667
668 cpl_matrix_unwrap(mvslit);
669 mvslit = NULL;
670
671 if (mnpsf == NULL) {
672 return -1;
673 }
674
675
676 /*
677 * Enforce positivity and normalization of the profile model.
678 */
679
680 _mnpsf = cpl_matrix_get_data(mnpsf);
681
682 norm = 0.;
683
684 for (i = 0; i < vslit->width; ++i) {
685 _mnpsf[i] = CX_MAX(_mnpsf[i], 0.);
686 norm += _mnpsf[i];
687 }
688
689 for (i = 0; i < vslit->width; ++i) {
690 _mnpsf[i] /= norm;
691 }
692
693
694 /*
695 * Estimate background and determine the number of valid pixels
696 */
697
698 tdata = cx_malloc(vslit->width * sizeof(cxdouble));
699
700 i = 0;
701 ngood = 0;
702
703 while (i < vslit->width) {
704 if (vslit->mask[i] > 0) {
705 tdata[ngood] = CX_MAX(vslit->signal[i], 0.);
706 ++ngood;
707 }
708 ++i;
709 }
710
711 if (ngood > 1) {
712 giraffe_array_sort(tdata, ngood);
713 bkg = 0.5 * (tdata[0] + tdata[1]);
714 }
715
716 cx_free(tdata);
717 tdata = NULL;
718
719
720 /*
721 * Try extraction only if there are good pixels available. If no good
722 * pixels are left skip this spectral bin and set the flux and the variance
723 * to zero.
724 */
725
726 if (ngood > 0) {
727
728 cxint iteration = 0;
729 cxint nreject = -1;
730 cxint niter = config->clip.iterations;
731 cxint nmin = (cxint)config->clip.fraction;
732
733 cxdouble sigma = config->clip.level * config->clip.level;
734 cxdouble* variance = NULL;
735
736
737 /*
738 * Compute standard extraction flux and rescale it to account for
739 * bad pixels.
740 */
741
742 norm = 0.;
743 flx = 0.;
744
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];
749 }
750 }
751
752 flx /= norm;
753
754
755 /*
756 * Allocate buffer for the variance estimates and compute the initial
757 * variances from the expected profile.
758 */
759
760 variance = cx_calloc(vslit->width, sizeof(cxdouble));
761
762 for (i = 0; i < vslit->width; ++i) {
763
764 register cxdouble ve = flx * _mnpsf[i] + bkg;
765
766 variance[i] = vslit->variance[i] + fabs(vslit->fraction[i] * ve);
767
768 }
769
770
771 /*
772 * Reject cosmics and extract spectrum
773 */
774
775 nreject = -1;
776
777 while ((iteration < niter) && (ngood > nmin) && (nreject != 0)) {
778
779 cxint imax = 0;
780
781 cxdouble _flx = 0.;
782 cxdouble mmax = 0.;
783
784
785 norm = 0.;
786 var = 0.;
787 nreject = 0;
788
789
790 /*
791 * Reject cosmics
792 */
793
794 for (i = 0; i < vslit->width; ++i) {
795
796 if (vslit->mask[i] != 0) {
797
798 cxdouble m = vslit->signal[i] - bkg - flx * _mnpsf[i];
799
800 m *= vslit->fraction[i];
801 m *= m / variance[i] ;
802
803 if (m > mmax) {
804 mmax = m;
805 imax = i;
806 }
807
808 }
809
810 }
811
812 if ((sigma > 0.) && (mmax > sigma)) {
813 vslit->mask[imax] = 0;
814 ++nreject;
815 --ngood;
816 }
817
818
819 /*
820 * Compute flux and variance estimates.
821 */
822
823 for (i = 0; i < vslit->width; ++i) {
824
825 if (vslit->mask[i] != 0) {
826
827 register cxdouble data = vslit->signal[i] - bkg;
828 register cxdouble p = _mnpsf[i];
829
830 data *= vslit->fraction[i];
831 p *= vslit->fraction[i];
832
833 norm += p * p / variance[i];
834 _flx += p * data / variance[i];
835 var += p;
836
837 }
838
839 }
840
841 flx = _flx / norm;
842 var /= norm;
843
844
845 /*
846 * Update variance estimates
847 */
848
849 for (i = 0; i < vslit->width; ++i) {
850
851 register cxdouble ve = flx * _mnpsf[i] + bkg;
852
853 variance[i] = vslit->variance[i] + fabs(vslit->fraction[i] * ve);
854
855 }
856
857 ++iteration;
858
859 }
860
861 cx_free(variance);
862 variance = NULL;
863
864 }
865
866 cpl_matrix_delete(mnpsf);
867 mnpsf = NULL;
868
869 result->value = flx;
870 result->error = sqrt(var);
871 result->position = vslit->center;
872 result->npixels = ngood;
873
874 return ngood == 0 ? 1 : 0;
875
876}
877
878
879/*
880 * @brief
881 * Compute the optimal extracted flux and its variance for a single
882 * wavelength bin.
883 *
884 * @param slice The results container to store the flux, variance and
885 * extraction model.
886 * @param AT The transposed design matrix of the linear system.
887 * @param S Column vector of the measured signal.
888 * @param W Matrix of weights of the measured signals.
889 * @param limits Cutoff parameters.
890 * @param ws Workspace for the matrix operations.
891 *
892 * @return
893 * The function returns 0 on success, and a non-zero value if an error
894 * occurred.
895 *
896 * The functions computes the optimal extracted fluxes for a single wavelength
897 * bin by solving the linear system:
898 * @f[
899 * \mathbf{f}\left(x\right) =
900 * \left(\mathbf{A}^\mathrm{T}\mathbf{W}\mathbf{A}\right)^{-1}
901 * \mathbf{A}^\mathrm{T}\mathbf{W}\mathbf{s}
902 * @f]
903 * where the @f$\mathbf{s}@f$ is the column vector of the measured fluxes
904 * written as an @f$\left(n \times 1\right)@f$ matrix, @f$\mathbf{W}@f$ is the
905 * diagonal @f$\left(n \times n\right)@f$ matrix of the weights, and
906 * @f$\mathbf{A}^\mathrm{T}@f$ is the transposed of the
907 * @f$\left(n \times m\right)@f$ design matrix @f$\mathbf{A}@f$.
908 *
909 * Defining the matrix @f$\mathbf{C} = \left(c_\mathrm{ij}\right) \equiv
910 * \left(\mathbf{A}^\mathrm{T}\mathbf{W}\mathbf{A}\right)^{-1}@f$, and using
911 * @f$a_\mathrm{ij}@f$ and @f$w_\mathrm{ij}@f$ to denote the elements of the
912 * transposed design matrix and the weight matrix respectively, the extracted
913 * flux can be written as the following sum:
914 * @f[
915 * f_\mathrm{i} = \sum\limits_{\mathrm{l} = 0}^{\mathrm{m} - 1} c_\mathrm{il}
916 * \sum\limits_{\mathrm{k} = 0}^{\mathrm{n} - 1} a_\mathrm{lk}
917 * w_\mathrm{kk} s_\mathrm{k}
918 * @f]
919 *
920 */
921
922inline static cxint
923_giraffe_optimal_extract_slice(GiExtractionSlice* slice,
924 const cpl_matrix* AT,
925 const cpl_matrix* S,
926 const cpl_matrix* W,
927 GiExtractionPsfLimits* limits,
928 GiExtractionWorkspace* ws)
929{
930
931 register cxint i = 0;
932 register cxint n = cpl_matrix_get_ncol(AT);
933 register cxint m = cpl_matrix_get_nrow(AT);
934
935 cxint status = 0;
936
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);
941
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);
948
949
950 for (i = 0; i < m; ++i) {
951
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];
957
958
959 atws[i] = 0.;
960
961 for (j = 0; j < n; ++j) {
962
963 register cxint k = in + j;
964
965
966 atw[k] = w[j] * at[k];
967 atws[i] += atw[k] * s[j];
968
969 }
970
971 for (j = 0; j < i; ++j) {
972
973 register cxint k = 0;
974 register cxint l = im + j;
975
976 atwa[l] = 0.;
977 for (k = ymin; k < ymax; ++k) {
978 atwa[l] += atw[in + k] * at[j * n + k];
979 }
980
981 atwa[j * m + i] = atwa[l];
982
983 }
984
985 atwa[im + i] = 0.;
986
987 for (j = ymin; j < ymax; ++j) {
988 atwa[im + i] += atw[in + j] * at[in + j];
989 }
990
991 }
992
993
994 status = _giraffe_matrix_invert(ws->c, ws->atwa, ws->tmp);
995
996 if (status != 0) {
997 return 1;
998 }
999
1000 for (i = 0; i < m; ++i) {
1001
1002 register cxint j = 0;
1003 register cxint im = i * m;
1004
1005
1006 sf[i] = 0.;
1007 sv[i] = c[im + i];
1008
1009 for (j = 0; j < m; ++j) {
1010 sf[i] += c[im + j] * atws[j];
1011 }
1012
1013 }
1014
1015 for (i = 0; i < n; ++i) {
1016
1017 register cxint j = 0;
1018
1019
1020 sm[i] = 0.;
1021
1022 for (j = 0; j < m; ++j) {
1023 sm[i] += at[j * n + i] * sf[j];
1024 }
1025
1026 }
1027
1028 return 0;
1029
1030}
1031
1032
1033/*
1034 * @brief
1035 * Extract spectra by simple summation.
1036 *
1037 * @param mz Pixels values [nx, ny]
1038 * @param mslz Scattered light model pixel values [nx, ny]
1039 * @param fibers List of fibers to extract
1040 * @param my Fiber centroid positions [ns, ny]
1041 * @param mw Fiber widths [ns, ny]
1042 * @param ms Extracted flux [ns, ny]
1043 * @param mse Extracted flux error [ns, ny]
1044 * @param msn Number of extracted pixels [ns, ny]
1045 * @param msy Extracted flux centroid position [ns, ny]
1046 *
1047 * For each X bin and each spectrum the flux is computed as the sum of
1048 * pixels value from @em mz[nx,ny] along the virtual slit computed using
1049 * the localization mask @em my[ns, ny] and @em mw[ns, ny].
1050 * The table @em fibers specifies the spectra to extract.
1051 *
1052 * The images @em ms[ns, ny], @em mse[ns, ny], @em msn[ns, ny] and
1053 * @em msy[ns, ny] must be allocated by the caller.
1054 */
1055
1056inline static cxint
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)
1062{
1063
1064 register cxint nn;
1065
1066 const cxchar* idx = NULL;
1067
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);
1072
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);
1077
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);
1082
1083
1084 /*
1085 * The number of fibers to be process must be less or equal to the
1086 * number of spectra available in the localization.
1087 */
1088
1089 cx_assert(nfibers <= nspectra);
1090
1091 idx = giraffe_fiberlist_query_index(fibers);
1092
1093 cx_assert(cpl_table_has_column(fibers, idx) != 0);
1094
1095 if (mbpx != NULL) {
1096
1097 const cxint* bpx = cpl_image_get_data_int(mbpx);
1098
1099 for (nn = 0; nn < nfibers; nn++) {
1100
1101 register cxint x;
1102 register cxint ns = cpl_table_get_int(fibers, idx, nn, NULL) - 1;
1103
1104
1105 for (x = 0; x < cpl_image_get_size_y(mz) && x < nbins; x++) {
1106
1107 cxint y;
1108 cxint yup, ylo;
1109 cxint lx = x * nspectra + ns;
1110 cxint sx = x * nfibers + nn;
1111
1112 cxdouble ylower = locy[lx] - locw[lx];
1113 cxdouble yupper = locy[lx] + locw[lx];
1114 cxdouble zsum = 0.;
1115 cxdouble ysum = 0.;
1116 cxdouble error2 = 0.;
1117
1118
1119 flux[sx] = 0.;
1120 flux_npixels[sx] = 0.;
1121 flux_error[sx] = 0.;
1122 flux_ypos[sx] = 0.;
1123
1124
1125 /*
1126 * Skip zero-width (invalid) spectra
1127 */
1128
1129 if (locw[lx] <= 0.0) {
1130 continue;
1131 }
1132
1133
1134 /*
1135 * Upper and lower border of the virtual slit. The real ones
1136 * and the borders corrected for pixel fractions. If we are
1137 * out of the the image boundaries we skip the extraction
1138 * for this bin and fiber.
1139 */
1140
1141 ylo = (cxint) ceil(ylower);
1142 yup = (cxint) floor(yupper);
1143
1144
1145 if (yup < 0. || ylo - 1 >= ny) {
1146 continue;
1147 }
1148
1149
1150 /*
1151 * Summation along the virtual slit. Check that y is always
1152 * in the range of valid pixels [0, ny[. Take into account
1153 * pixel fractions at the beginning and the end of the
1154 * virtual slit.
1155 */
1156
1157 /*
1158 * Lower ordinate pixel fraction
1159 */
1160
1161 y = ylo - 1;
1162
1163 if (y >= 0) {
1164
1165 if (!(bpx[x * ny + y] & GIR_M_PIX_SET)) {
1166
1167 cxdouble extcoeff = (cxdouble)ylo - ylower;
1168 cxdouble extcoeff2 = extcoeff * extcoeff;
1169 cxdouble px = CX_MAX(pixels[x * ny + y], 0.);
1170
1171 flux[sx] = pixels[x * ny + y] * extcoeff;
1172 flux_npixels[sx] = extcoeff;
1173 error2 = variances[x * ny + y] * extcoeff2;
1174
1175 zsum = px * extcoeff;
1176 ysum = y * px * extcoeff;
1177
1178 }
1179
1180 }
1181
1182
1183 /*
1184 * Sum pixel values along virtual slit.
1185 */
1186
1187 for (y = CX_MAX(ylo, 0); (y < yup) && (y < ny); y++) {
1188
1189 if (!(bpx[x * ny + y] & GIR_M_PIX_SET)) {
1190
1191 cxdouble px = CX_MAX(pixels[x * ny + y], 0.);
1192
1193 flux[sx] += pixels[x * ny + y];
1194 flux_npixels[sx] += 1.0;
1195 error2 += variances[x * ny + y];
1196
1197 zsum += px;
1198 ysum += y * px;
1199
1200 }
1201
1202 }
1203
1204
1205 /*
1206 * Upper ordinate pixel fraction
1207 */
1208
1209 y = yup;
1210
1211 if (y < ny) {
1212
1213 if (!(bpx[x * ny + y] & GIR_M_PIX_SET)) {
1214
1215 cxdouble extcoeff = yupper - (cxdouble)yup;
1216 cxdouble extcoeff2 = extcoeff * extcoeff;
1217 cxdouble px = CX_MAX(pixels[x * ny + y], 0.);
1218
1219 flux[sx] += pixels[x * ny + y] * extcoeff;
1220 flux_npixels[sx] += extcoeff;
1221 error2 += variances[x * ny + y] * extcoeff2;
1222
1223 zsum += px * extcoeff;
1224 ysum += y * px * extcoeff;
1225
1226 }
1227
1228 }
1229
1230 flux_error[sx] = sqrt(error2);
1231
1232 // FIXME: Check this protection from division by zero. Also
1233 // the minimum condition for the pixel values above.
1234
1235 if (fabs(ysum) < DBL_EPSILON || fabs(zsum) < DBL_EPSILON) {
1236 flux_ypos[sx] = 0.5 * (yupper + ylower);
1237 }
1238 else {
1239 flux_ypos[sx] = ysum / zsum;
1240 }
1241
1242 }
1243
1244 }
1245
1246 }
1247 else {
1248
1249 for (nn = 0; nn < nfibers; nn++) {
1250
1251 register cxint x;
1252 register cxint ns = cpl_table_get_int(fibers, idx,
1253 nn, NULL) - 1;
1254
1255
1256 for (x = 0; x < cpl_image_get_size_y(mz) && x < nbins; x++) {
1257
1258 cxint y;
1259 cxint yup, ylo;
1260 cxint lx = x * nspectra + ns;
1261 cxint sx = x * nfibers + nn;
1262
1263 cxdouble yupper, ylower;
1264 cxdouble zsum = 0.;
1265 cxdouble ysum = 0.;
1266 cxdouble error2 = 0.;
1267
1268
1269 flux[sx] = 0.;
1270 flux_npixels[sx] = 0.;
1271 flux_error[sx] = 0.;
1272 flux_ypos[sx] = 0.;
1273
1274
1275 /*
1276 * Skip zero-width (invalid) spectra
1277 */
1278
1279 if (locw[lx] <= 0.0) {
1280 continue;
1281 }
1282
1283
1284 /*
1285 * Upper and lower border of the virtual slit. The real ones
1286 * and the borders corrected for pixel fractions. If we are
1287 * out of the the image boundaries we skip the extraction
1288 * for this bin and fiber.
1289 */
1290
1291 yupper = locy[lx] + locw[lx];
1292 ylower = locy[lx] - locw[lx];
1293
1294 ylo = (cxint) ceil(ylower);
1295 yup = (cxint) floor(yupper);
1296
1297
1298 if (yup < 0. || ylo - 1 >= ny) {
1299 continue;
1300 }
1301
1302
1303 /*
1304 * Summation along the virtual slit. Check that y is always
1305 * in the range of valid pixels [0, ny[. Take into account
1306 * pixel fractions at the beginning and the end of the
1307 * virtual slit.
1308 */
1309
1310 /*
1311 * Lower ordinate pixel fraction
1312 */
1313
1314 y = ylo - 1;
1315
1316 if (y >= 0) {
1317
1318 cxdouble extcoeff = (cxdouble)ylo - ylower;
1319 cxdouble extcoeff2 = extcoeff * extcoeff;
1320 cxdouble px = CX_MAX(pixels[x * ny + y], 0.);
1321
1322 flux[sx] = pixels[x * ny + y] * extcoeff;
1323 flux_npixels[sx] = extcoeff;
1324 error2 = variances[x * ny + y] * extcoeff2;
1325
1326 zsum = px * extcoeff;
1327 ysum = y * px * extcoeff;
1328
1329 }
1330
1331
1332 /*
1333 * Sum pixel values along virtual slit.
1334 */
1335
1336 for (y = CX_MAX(ylo, 0); (y < yup) && (y < ny); y++) {
1337
1338 cxdouble px = CX_MAX(pixels[x * ny + y], 0.);
1339
1340 flux[sx] += pixels[x * ny + y];
1341 flux_npixels[sx] += 1.0;
1342 error2 += variances[x * ny + y];
1343
1344 zsum += px;
1345 ysum += y * px;
1346 }
1347
1348
1349 /*
1350 * Upper ordinate pixel fraction
1351 */
1352
1353 y = yup;
1354
1355 if (y < ny) {
1356
1357 cxdouble extcoeff = yupper - (cxdouble)yup;
1358 cxdouble extcoeff2 = extcoeff * extcoeff;
1359 cxdouble px = CX_MAX(pixels[x * ny + y], 0.);
1360
1361 flux[sx] += pixels[x * ny + y] * extcoeff;
1362 flux_npixels[sx] += extcoeff;
1363 error2 += variances[x * ny + y] * extcoeff2;
1364
1365 zsum += px * extcoeff;
1366 ysum += y * px * extcoeff;
1367
1368 }
1369
1370 flux_error[sx] = sqrt(error2);
1371
1372 // FIXME: Check this protection from division by zero. Also
1373 // the minimum condition for the pixel values above.
1374
1375 if (fabs(ysum) < DBL_EPSILON || fabs(zsum) < DBL_EPSILON) {
1376 flux_ypos[sx] = 0.5 * (yupper + ylower);
1377 }
1378 else {
1379 flux_ypos[sx] = ysum / zsum;
1380 }
1381
1382 }
1383
1384 }
1385
1386 }
1387
1388 return 0;
1389
1390}
1391
1392
1393/*
1394 * @brief
1395 * Extract spectra using the optimal extraction method.
1396 *
1397 * @param mz Pixels values [nx, ny]
1398 * @param mvar Initial variance [nx, ny]
1399 * @param fibers List of fibers to extract
1400 * @param my Fiber centroid positions [ns, ny]
1401 * @param mw Fiber widths [ns, ny]
1402 * @param ms Extracted flux [ns, ny]
1403 * @param mse Extracted flux error [ns, ny]
1404 * @param msn Number of extracted pixels [ns, ny]
1405 * @param msy Extracted flux centroid position [ns, ny]
1406 * @param config Optimal extraction method setup
1407 *
1408 * TBD
1409 *
1410 * The images @em ms[ns, ny], @em mse[ns, ny], @em msn[ns, ny] and
1411 * @em msy[ns, ny] must be allocated by the caller.
1412 */
1413
1414inline static cxint
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)
1421{
1422
1423 const cxchar* idx = NULL;
1424
1425 cxint nx = 0;
1426 cxint ny = 0;
1427 cxint fiber = 0;
1428 cxint nfibers = 0;
1429
1430 const cxdouble* locy = NULL;
1431 const cxdouble* locw = NULL;
1432 const cxdouble* width = NULL;
1433 const cxdouble* exponent = NULL;
1434
1435 GiModel* psfmodel = NULL;
1436
1437
1438 cx_assert(mz != NULL);
1439 cx_assert(mzvar != NULL);
1440
1441 cx_assert(fibers != NULL);
1442
1443 cx_assert(my != NULL);
1444 cx_assert(mw != NULL);
1445
1446 cx_assert(psfdata != NULL);
1447
1448 cx_assert(ms != NULL);
1449 cx_assert(mse != NULL);
1450 cx_assert(msn != NULL);
1451 cx_assert(msy != NULL);
1452
1453 cx_assert(config != NULL);
1454
1455 ny = cpl_image_get_size_x(mz);
1456 nx = cpl_image_get_size_y(mz);
1457 nfibers = cpl_table_get_nrow(fibers);
1458
1459 locy = cpl_image_get_data_double_const(my);
1460 locw = cpl_image_get_data_double_const(mw);
1461
1462 cx_assert((ny == cpl_image_get_size_x(mzvar)) &&
1463 (nx == cpl_image_get_size_y(mzvar)));
1464
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));
1467
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));
1472
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)));
1481
1482 cx_assert((mbpx == NULL) || ((ny == cpl_image_get_size_x(mbpx)) &&
1483 (nx == cpl_image_get_size_y(mbpx))));
1484
1485
1486 /*
1487 * Get the index column mapping the current spectum number to the
1488 * corresponding reference localization spectrum number.
1489 */
1490
1491 idx = giraffe_fiberlist_query_index(fibers);
1492
1493 cx_assert(cpl_table_has_column(fibers, idx) != 0);
1494
1495
1496 /*
1497 * Get the PSF profile data arrays for efficency reasons in the
1498 * following loops.
1499 */
1500
1501 if (giraffe_psfdata_contains(psfdata, "Center") == FALSE) {
1502 return -1;
1503 }
1504
1505 if (giraffe_psfdata_contains(psfdata, "Width2") == TRUE) {
1506 exponent = cpl_image_get_data_const(giraffe_psfdata_get_data(psfdata,
1507 "Width2"));
1508 }
1509
1510 width = cpl_image_get_data_const(giraffe_psfdata_get_data(psfdata,
1511 "Width1"));
1512
1513
1514 /*
1515 * Create the PSF profile model from the PSF data object.
1516 */
1517
1518 psfmodel = giraffe_model_new(giraffe_psfdata_get_model(psfdata));
1519
1520 if (psfmodel == NULL) {
1521 return -2;
1522 }
1523
1524 giraffe_model_set_parameter(psfmodel, "Center", 0.);
1525 giraffe_model_set_parameter(psfmodel, "Amplitude", 1.);
1526 giraffe_model_set_parameter(psfmodel, "Background", 0.);
1527
1528
1529 /*
1530 * Extract each fiber spectrum
1531 */
1532
1533 for (fiber = 0; fiber < nfibers; ++fiber) {
1534
1535 register cxint bin = 0;
1536 register cxint fidx = cpl_table_get_int(fibers, idx, fiber, NULL) - 1;
1537
1538 cxint nbins = CX_MIN(nx, cpl_image_get_size_y(my));
1539
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);
1544
1545
1546 for (bin = 0; bin < nbins; bin++) {
1547
1548 register cxint lpos = bin * cpl_image_get_size_x(my) + fidx;
1549 register cxint spos = bin * nfibers + fiber;
1550
1551 cxint status = 0;
1552 cxint vwidth = 0;
1553
1554 register cxdouble lcenter = locy[lpos];
1555 register cxdouble lwidth = locw[lpos];
1556
1557 register cxdouble ylower = lcenter - lwidth;
1558 register cxdouble yupper = lcenter + lwidth;
1559
1560 GiVirtualSlit* vslit = NULL;
1561
1562 GiExtractionData result = {0., 0., 0., 0.};
1563
1564
1565 /*
1566 * Skip zero-width, invalid spectra
1567 */
1568
1569 if ((lwidth <= 0.) || (yupper < 0.) || (ylower > ny)) {
1570 continue;
1571 }
1572
1573 /*
1574 * Fill the virtual slit with data
1575 */
1576
1577 vslit = _giraffe_virtualslit_new(config->ewidth);
1578
1579 vwidth = _giraffe_virtualslit_setup(vslit, bin, lcenter, lwidth,
1580 mz, mzvar, mbpx);
1581
1582 if (vwidth == 0) {
1583 _giraffe_virtualslit_delete(vslit);
1584 vslit = NULL;
1585
1586 continue;
1587 }
1588
1589
1590 /*
1591 * Update PSF profile model width and exponent
1592 */
1593
1594 giraffe_model_set_parameter(psfmodel, "Width1", width[lpos]);
1595
1596 if (exponent != NULL) {
1597 giraffe_model_set_parameter(psfmodel, "Width2",
1598 exponent[lpos]);
1599 }
1600
1601
1602 /*
1603 * Compute flux from the virtual slit using Horne's optimal
1604 * extraction algorithm.
1605 */
1606
1607 status = _giraffe_horne_extract_slit(&result, vslit, psfmodel,
1608 config);
1609
1610 _giraffe_virtualslit_delete(vslit);
1611 vslit = NULL;
1612
1613 if (status < 0) {
1614
1615 giraffe_model_delete(psfmodel);
1616 psfmodel = NULL;
1617
1618 return 1;
1619 }
1620
1621 _ms[spos] = result.value;
1622 _mse[spos] = result.error;
1623 _msy[spos] = result.position;
1624 _msn[spos] = result.npixels;
1625
1626 }
1627
1628 }
1629
1630
1631 giraffe_model_delete(psfmodel);
1632 psfmodel = NULL;
1633
1634 return 0;
1635
1636}
1637
1638
1639/*
1640 * Fill extraction matrix with the fiber profiles and the coefficients of
1641 * the Chebyshev polynomial model of the background.
1642 */
1643
1644inline static cxint
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)
1651{
1652
1653 const cxchar* idx = giraffe_fiberlist_query_index(fibers);
1654
1655 cxint fiber = 0;
1656 cxint nfibers = cpl_table_get_nrow(fibers);
1657 cxint ny = cpl_matrix_get_ncol(profiles);
1658
1659 const cxdouble* locy = cpl_image_get_data_double_const(my);
1660 const cxdouble* locw = cpl_image_get_data_double_const(mw);
1661
1662 cxdouble* _profiles = cpl_matrix_get_data(profiles);
1663
1664 cxdouble* ypos = NULL;
1665
1666
1667 cx_assert(cpl_table_has_column(fibers, idx) != 0);
1668 cx_assert((limits == NULL) ||
1669 (cpl_matrix_get_nrow(profiles) == limits->size));
1670
1671 ypos = cx_calloc(ny, sizeof(cxdouble));
1672
1673 for (fiber = 0; fiber < nfibers; ++fiber) {
1674
1675 register cxint i = 0;
1676 register cxint y = 0;
1677 register cxint k = 0;
1678
1679 cxint fidx = cpl_table_get_int(fibers, idx, fiber, NULL) - 1;
1680 cxint lpos = bin * cpl_image_get_size_x(my) + fidx;
1681
1682 register cxdouble lcenter = locy[lpos];
1683 register cxdouble lwidth = locw[lpos];
1684
1685 register cxdouble ylower = lcenter - fabs(wfactor) * lwidth;
1686 register cxdouble yupper = lcenter + fabs(wfactor) * lwidth;
1687
1688 register cxint first = (cxint) floor(ylower);
1689 register cxint last = (cxint) ceil(yupper);
1690
1691 register cxint vwidth = 0;
1692
1693 cxdouble norm = 0.;
1694 cxdouble* _mnpsf = NULL;
1695
1696 cpl_matrix* positions = NULL;
1697 cpl_matrix* mnpsf = NULL;
1698
1699
1700 /*
1701 * Upper, lower border and width of the virtual slit
1702 */
1703
1704 ylower = CX_MAX(0., ylower);
1705 yupper = CX_MIN(ny - 1., yupper);
1706
1707 first = CX_MAX(0, first);
1708 last = CX_MIN(ny - 1, last);
1709
1710 vwidth = last - first + 1;
1711
1712 if (vwidth < 1) {
1713 continue;
1714 }
1715
1716 if (limits != NULL) {
1717 limits->ymin[fiber] = first;
1718 limits->ymax[fiber] = last + 1;
1719 }
1720
1721
1722 /*
1723 * Update PSF profile model width and exponent
1724 */
1725
1726 giraffe_model_set_parameter(psf, "Width1", width[lpos]);
1727
1728 if (exponent != NULL) {
1729 giraffe_model_set_parameter(psf, "Width2", exponent[lpos]);
1730 }
1731
1732
1733 /*
1734 * Compute normalized psf model
1735 */
1736
1737 k = 0;
1738 for (y = first; y <= last; ++y) {
1739 ypos[k] = y - lcenter;
1740 ++k;
1741 }
1742
1743 positions = cpl_matrix_wrap(vwidth, 1, ypos);
1744 mnpsf = _giraffe_compute_psf(psf, positions);
1745
1746 cpl_matrix_unwrap(positions);
1747 positions = NULL;
1748
1749 if (mnpsf == NULL) {
1750 cx_free(ypos);
1751 ypos = NULL;
1752
1753 return 1;
1754 }
1755
1756 _mnpsf = cpl_matrix_get_data(mnpsf);
1757
1758 for (i = 0; i < vwidth; ++i) {
1759 _mnpsf[i] = CX_MAX(_mnpsf[i], 0.);
1760 norm += _mnpsf[i];
1761 }
1762
1763 for (i = 0; i < vwidth; ++i) {
1764 _mnpsf[i] /= norm;
1765 }
1766
1767 k = fiber * ny + first;
1768 for (y = 0; y < vwidth; ++y) {
1769 _profiles[k + y] = _mnpsf[y];
1770 }
1771
1772 cpl_matrix_delete(mnpsf);
1773 mnpsf = NULL;
1774
1775 }
1776
1777 cx_free(ypos);
1778 ypos = NULL;
1779
1780 return 0;
1781
1782}
1783
1784
1785inline static cxint
1786_giraffe_extract_optimal(const cpl_image* mz, const cpl_image* mzvar,
1787 const cpl_table* fibers, const cpl_image* my,
1788 const cpl_image* mw, const GiPsfData* psfdata,
1789 cpl_image* mbpx, cpl_image* ms, cpl_image* mse,
1790 cpl_image* msm, cpl_image* msy,
1791 const GiExtractOptimalConfig* config)
1792{
1793
1794 const cxbool nolimits = (config->limits == TRUE) ? FALSE : TRUE;
1795
1796 const cxint bkg_nc = config->bkgorder + 1;
1797 const cxint niter = config->clip.iterations;
1798
1799 register cxint i = 0;
1800
1801 cxint nx = 0;
1802 cxint ny = 0;
1803 cxint bin = 0;
1804 cxint nbins = 0;
1805 cxint nfibers = 0;
1806
1807 const cxdouble wfactor = config->ewidth;
1808 const cxdouble sigma = config->clip.level * config->clip.level;
1809 const cxdouble fraction = config->clip.fraction;
1810
1811 const cxdouble* width = NULL;
1812 const cxdouble* exponent = NULL;
1813
1814 cxdouble* _ypos = NULL;
1815 cxdouble* _bkg_base = NULL;
1816 cxdouble* _profiles = NULL;
1817 cxdouble* _signal = NULL;
1818 cxdouble* _variance = NULL;
1819 cxdouble* _mask = NULL;
1820 cxdouble* _weights = NULL;
1821
1822 cpl_matrix* ypos = NULL;
1823 cpl_matrix* bkg_base = NULL;
1824 cpl_matrix* profiles = NULL;
1825 cpl_matrix* weights = NULL;
1826 cpl_matrix* signal = NULL;
1827 cpl_matrix* variance = NULL;
1828 cpl_matrix* mask = NULL;
1829
1830 GiModel* psfmodel = NULL;
1831
1832 GiExtractionPsfLimits* limits = NULL;
1833
1834 GiExtractionSlice* slice = NULL;
1835
1836 GiExtractionWorkspace* workspace;
1837
1838
1839 cx_assert(mz != NULL);
1840 cx_assert(mzvar != NULL);
1841
1842 cx_assert(fibers != NULL);
1843
1844 cx_assert(my != NULL);
1845 cx_assert(mw != NULL);
1846
1847 cx_assert(psfdata != NULL);
1848
1849 cx_assert(ms != NULL);
1850 cx_assert(mse != NULL);
1851 cx_assert(msm != NULL);
1852 cx_assert(msy != NULL);
1853
1854 ny = cpl_image_get_size_x(mz);
1855 nx = cpl_image_get_size_y(mz);
1856
1857 nfibers = cpl_table_get_nrow(fibers);
1858 nbins = CX_MIN(nx, cpl_image_get_size_y(my));
1859
1860 cx_assert((ny == cpl_image_get_size_x(mzvar)) &&
1861 (nx == cpl_image_get_size_y(mzvar)));
1862
1863 cx_assert(cpl_image_get_size_x(my) == cpl_image_get_size_x(mw));
1864 cx_assert(cpl_image_get_size_y(my) == cpl_image_get_size_y(mw));
1865
1866 cx_assert(giraffe_psfdata_fibers(psfdata) ==
1867 (cxsize)cpl_image_get_size_x(my));
1868 cx_assert(giraffe_psfdata_bins(psfdata) ==
1869 (cxsize)cpl_image_get_size_y(my));
1870
1871 cx_assert((nfibers == cpl_image_get_size_x(ms)) &&
1872 (nx == cpl_image_get_size_y(ms)));
1873 cx_assert((nfibers == cpl_image_get_size_x(mse)) &&
1874 (nx == cpl_image_get_size_y(mse)));
1875 cx_assert((nfibers == cpl_image_get_size_x(msy)) &&
1876 (nx == cpl_image_get_size_y(msy)));
1877 cx_assert((ny == cpl_image_get_size_x(msm)) &&
1878 (nx == cpl_image_get_size_y(msm)));
1879
1880 cx_assert((mbpx == NULL) || ((ny == cpl_image_get_size_x(mbpx)) &&
1881 (nx == cpl_image_get_size_y(mbpx))));
1882
1883
1884 /*
1885 * Get the PSF profile data arrays for efficiency reasons in the
1886 * following loops.
1887 */
1888
1889 if (giraffe_psfdata_contains(psfdata, "Center") == FALSE) {
1890 return -1;
1891 }
1892
1893 if (giraffe_psfdata_contains(psfdata, "Width2") == TRUE) {
1894 exponent = cpl_image_get_data_const(giraffe_psfdata_get_data(psfdata,
1895 "Width2"));
1896 }
1897
1898 width = cpl_image_get_data_const(giraffe_psfdata_get_data(psfdata,
1899 "Width1"));
1900
1901
1902 /*
1903 * Create the PSF profile model from the PSF data object.
1904 */
1905
1906 psfmodel = giraffe_model_new(giraffe_psfdata_get_model(psfdata));
1907
1908 if (psfmodel == NULL) {
1909 return -2;
1910 }
1911
1912 giraffe_model_set_parameter(psfmodel, "Amplitude", 1.);
1913 giraffe_model_set_parameter(psfmodel, "Background", 0.);
1914 giraffe_model_set_parameter(psfmodel, "Center", 0.);
1915
1916
1917 /*
1918 * Set up the vector of pixel positions
1919 */
1920
1921 ypos = cpl_matrix_new(ny, 1);
1922
1923 if (ypos == NULL) {
1924 giraffe_model_delete(psfmodel);
1925 psfmodel = NULL;
1926
1927 return -3;
1928 }
1929
1930 _ypos = cpl_matrix_get_data(ypos);
1931
1932 for (i = 0; i < ny; ++i) {
1933 _ypos[i] = i;
1934 }
1935
1936
1937 /*
1938 * Create profile matrix and the matrices for the signal, bad
1939 * pixel mask, variance and the weights.
1940 */
1941
1942 profiles = cpl_matrix_new(nfibers + bkg_nc, ny);
1943
1944 if (profiles == NULL) {
1945 cpl_matrix_delete(ypos);
1946 ypos = NULL;
1947
1948 giraffe_model_delete(psfmodel);
1949 psfmodel = NULL;
1950
1951 return -3;
1952 }
1953
1954 _profiles = cpl_matrix_get_data(profiles);
1955
1956
1957 signal = cpl_matrix_new(ny, 1);
1958
1959 if (signal == NULL) {
1960 cpl_matrix_delete(profiles);
1961 profiles = NULL;
1962
1963 cpl_matrix_delete(ypos);
1964 ypos = NULL;
1965
1966 giraffe_model_delete(psfmodel);
1967 psfmodel = NULL;
1968
1969 return -3;
1970 }
1971
1972 _signal = cpl_matrix_get_data(signal);
1973
1974
1975 variance = cpl_matrix_new(ny, 1);
1976
1977 if (variance == NULL) {
1978 cpl_matrix_delete(signal);
1979 signal = NULL;
1980
1981 cpl_matrix_delete(profiles);
1982 profiles = NULL;
1983
1984 cpl_matrix_delete(ypos);
1985 ypos = NULL;
1986
1987 giraffe_model_delete(psfmodel);
1988 psfmodel = NULL;
1989
1990 return -3;
1991 }
1992
1993 _variance = cpl_matrix_get_data(variance);
1994
1995
1996 mask = cpl_matrix_new(ny, 1);
1997
1998 if (mask == NULL) {
1999 cpl_matrix_delete(variance);
2000 variance = NULL;
2001
2002 cpl_matrix_delete(signal);
2003 signal = NULL;
2004
2005 cpl_matrix_delete(profiles);
2006 profiles = NULL;
2007
2008 cpl_matrix_delete(ypos);
2009 ypos = NULL;
2010
2011 giraffe_model_delete(psfmodel);
2012 psfmodel = NULL;
2013
2014 return -3;
2015 }
2016
2017 _mask = cpl_matrix_get_data(mask);
2018
2019
2020 weights = cpl_matrix_new(ny, 1);
2021
2022 if (weights == NULL) {
2023 cpl_matrix_delete(mask);
2024 mask = NULL;
2025
2026 cpl_matrix_delete(variance);
2027 variance = NULL;
2028
2029 cpl_matrix_delete(signal);
2030 signal = NULL;
2031
2032 cpl_matrix_delete(profiles);
2033 profiles = NULL;
2034
2035 cpl_matrix_delete(ypos);
2036 ypos = NULL;
2037
2038 giraffe_model_delete(psfmodel);
2039 psfmodel = NULL;
2040
2041 return -3;
2042 }
2043
2044 _weights = cpl_matrix_get_data(weights);
2045
2046
2047 /*
2048 * Fill design matrix with the basis functions of the
2049 * background polynomial model.
2050 */
2051
2052 bkg_base = giraffe_chebyshev_base1d(0., ny, bkg_nc, ypos);
2053
2054 cpl_matrix_delete(ypos);
2055 ypos = NULL;
2056
2057 if (bkg_base == NULL) {
2058 cpl_matrix_delete(weights);
2059 weights = NULL;
2060
2061 cpl_matrix_delete(mask);
2062 mask = NULL;
2063
2064 cpl_matrix_delete(variance);
2065 variance = NULL;
2066
2067 cpl_matrix_delete(signal);
2068 signal = NULL;
2069
2070 cpl_matrix_delete(profiles);
2071 profiles = NULL;
2072
2073 cpl_matrix_delete(ypos);
2074 ypos = NULL;
2075
2076 giraffe_model_delete(psfmodel);
2077 psfmodel = NULL;
2078
2079 return -3;
2080 }
2081
2082 _bkg_base = cpl_matrix_get_data(bkg_base);
2083
2084 for (i = 0; i < bkg_nc; ++i) {
2085
2086 register cxint j = 0;
2087 register cxint offset = nfibers * ny;
2088
2089 for (j = 0; j < ny; ++j) {
2090 _profiles[i * ny + j + offset] = _bkg_base[i * ny + j];
2091 }
2092
2093 }
2094
2095 _bkg_base = NULL;
2096
2097 cpl_matrix_delete(bkg_base);
2098 bkg_base = NULL;
2099
2100
2101 /*
2102 * Extract all fiber spectra simultaneously for each wavelength bin
2103 */
2104
2105 slice = _giraffe_extractionslice_new(nfibers, ny, bkg_nc);
2106
2107 if (slice == NULL) {
2108 cpl_matrix_delete(weights);
2109 weights = NULL;
2110
2111 cpl_matrix_delete(mask);
2112 mask = NULL;
2113
2114 cpl_matrix_delete(variance);
2115 variance = NULL;
2116
2117 cpl_matrix_delete(signal);
2118 signal = NULL;
2119
2120 cpl_matrix_delete(profiles);
2121 profiles = NULL;
2122
2123 cpl_matrix_delete(ypos);
2124 ypos = NULL;
2125
2126 giraffe_model_delete(psfmodel);
2127 psfmodel = NULL;
2128
2129 return -3;
2130 }
2131
2132
2133 limits = _giraffe_extraction_psflimits_new(nfibers + bkg_nc);
2134
2135 if (limits == NULL) {
2136
2137 _giraffe_extractionslice_delete(slice);
2138 slice = NULL;
2139
2140 cpl_matrix_delete(weights);
2141 weights = NULL;
2142
2143 cpl_matrix_delete(mask);
2144 mask = NULL;
2145
2146 cpl_matrix_delete(variance);
2147 variance = NULL;
2148
2149 cpl_matrix_delete(signal);
2150 signal = NULL;
2151
2152 cpl_matrix_delete(profiles);
2153 profiles = NULL;
2154
2155 cpl_matrix_delete(ypos);
2156 ypos = NULL;
2157
2158 giraffe_model_delete(psfmodel);
2159 psfmodel = NULL;
2160
2161 return -3;
2162
2163 }
2164
2165 for (i = 0; i < limits->size; ++i) {
2166 limits->ymin[i] = 0;
2167 limits->ymax[i] = ny;
2168 }
2169
2170
2171 /*
2172 * Allocate workspace for matrix multiplications
2173 */
2174
2175 workspace = _giraffe_optimal_workspace_new(nfibers + bkg_nc, ny);
2176
2177 for (bin = 0; bin < nbins; ++bin) {
2178
2179 cxbool stop = FALSE;
2180
2181 cxint iter = 0;
2182 cxint nmin = 0;
2183 cxint ngood = ny;
2184
2185 const cxdouble* _my = cpl_image_get_data_double_const(my);
2186 const cxdouble* _mz = cpl_image_get_data_double_const(mz);
2187 const cxdouble* _mzvar = cpl_image_get_data_double_const(mzvar);
2188
2189 cxdouble* _ms = cpl_image_get_data_double(ms);
2190 cxdouble* _mse = cpl_image_get_data_double(mse);
2191 cxdouble* _msy = cpl_image_get_data_double(msy);
2192 cxdouble* _msm = cpl_image_get_data_double(msm);
2193
2194 cxint status = 0;
2195
2196 GiExtractionPsfLimits* _limits = (nolimits == FALSE) ? limits : NULL;
2197
2198 cx_assert(_mz != NULL);
2199 cx_assert(_mzvar != NULL);
2200
2201
2202 /*
2203 * Fill the design matrix with the fiber profiles for the
2204 * current wavelength bin
2205 */
2206
2207 status = _giraffe_optimal_build_profiles(profiles, _limits, my, mw,
2208 fibers, bin, psfmodel, width,
2209 exponent, wfactor);
2210
2211 if (status != 0) {
2212 _giraffe_optimal_workspace_delete(workspace);
2213 workspace = NULL;
2214
2215 _giraffe_extraction_psflimits_delete(limits);
2216 limits = NULL;
2217
2218 _giraffe_extractionslice_delete(slice);
2219 slice = NULL;
2220
2221 cpl_matrix_delete(weights);
2222 weights = NULL;
2223
2224 cpl_matrix_delete(mask);
2225 mask = NULL;
2226
2227 cpl_matrix_delete(variance);
2228 variance = NULL;
2229
2230 cpl_matrix_delete(signal);
2231 signal = NULL;
2232
2233 cpl_matrix_delete(profiles);
2234 profiles = NULL;
2235
2236 cpl_matrix_delete(ypos);
2237 ypos = NULL;
2238
2239 giraffe_model_delete(psfmodel);
2240 psfmodel = NULL;
2241
2242 return -4;
2243 }
2244
2245
2246 /*
2247 * Fill the signal, variance, mask and weight matrices
2248 */
2249
2250
2251 if (mbpx != NULL) {
2252
2253 const cxint* _mbpx = cpl_image_get_data_int_const(mbpx);
2254
2255
2256 cx_assert(_mbpx != NULL);
2257
2258 for (i = 0; i < ny; ++i) {
2259
2260 cxbool bad = (_mbpx[bin * ny + i] & GIR_M_PIX_SET) ||
2261 (_mz[bin * ny + i] < 0.);
2262
2263 _signal[i] = _mz[bin * ny + i];
2264 _variance[i] = _signal[i] + _mzvar[bin * ny + i];
2265 _mask[i] = 1.;
2266
2267 if (bad == TRUE) {
2268 _mask[i] = 0.;
2269 --ngood;
2270 }
2271
2272 _weights[i] = _mask[i] / _variance[i];
2273
2274 }
2275
2276 }
2277 else {
2278
2279 for (i = 0; i < ny; ++i) {
2280
2281 cxbool bad = (_mz[bin * ny + i] < 0.);
2282
2283 _signal[i] = _mz[bin * ny + i];
2284 _variance[i] = _signal[i] + _mzvar[bin * ny + i];
2285 _mask[i] = 1.;
2286
2287 if (bad == TRUE) {
2288 _mask[i] = 0.;
2289 --ngood;
2290 }
2291
2292 _weights[i] = _mask[i] / _variance[i];
2293
2294 }
2295
2296 }
2297
2298
2299 /*
2300 * Extract simultaneously the fluxes of all fibers for the current
2301 * wavelength bin
2302 */
2303
2304 nmin = (cxint)(fraction * ngood);
2305
2306 while ((iter < niter) && (stop == FALSE)) {
2307
2308 cxint nreject = 0;
2309
2310 const cxdouble* _model = NULL;
2311
2312
2313 status = _giraffe_optimal_extract_slice(slice, profiles,
2314 signal, weights, limits, workspace);
2315
2316 if (status != 0) {
2317 _giraffe_optimal_workspace_delete(workspace);
2318 workspace = NULL;
2319
2320 _giraffe_extraction_psflimits_delete(limits);
2321 limits = NULL;
2322
2323 _giraffe_extractionslice_delete(slice);
2324 slice = NULL;
2325
2326 cpl_matrix_delete(weights);
2327 weights = NULL;
2328
2329 cpl_matrix_delete(mask);
2330 mask = NULL;
2331
2332 cpl_matrix_delete(variance);
2333 variance = NULL;
2334
2335 cpl_matrix_delete(signal);
2336 signal = NULL;
2337
2338 cpl_matrix_delete(profiles);
2339 profiles = NULL;
2340
2341 cpl_matrix_delete(ypos);
2342 ypos = NULL;
2343
2344 giraffe_model_delete(psfmodel);
2345 psfmodel = NULL;
2346
2347 return -5;
2348 }
2349
2350
2351 /*
2352 * Update weighting factors
2353 */
2354
2355 _model = cpl_matrix_get_data(slice->model);
2356
2357 for (i = 0; i < ny; ++i) {
2358
2359 if (_mask[i] > 0.) {
2360
2361 cxbool bad = FALSE;
2362 cxdouble residual = _signal[i] - _model[i];
2363
2364
2365 _variance[i] = _model[i] + _mzvar[bin * ny + i];
2366
2367 bad = (residual * residual) > (sigma * _variance[i]) ?
2368 TRUE : FALSE;
2369
2370 if (bad == TRUE) {
2371 _mask[i] = 0.;
2372 ++nreject;
2373 --ngood;
2374 }
2375
2376 _weights[i] = _mask[i] / _variance[i];
2377
2378 }
2379
2380 }
2381
2382 if ((nreject == 0) || (ngood <= nmin)) {
2383 stop = TRUE;
2384 }
2385
2386 ++iter;
2387
2388 }
2389
2390
2391 /*
2392 * Copy the extracted fluxes, their variance and the modeled signal
2393 * to the result images.
2394 */
2395
2396 memcpy(&_ms[bin * nfibers], cpl_matrix_get_data(slice->flux),
2397 slice->nflx * sizeof(cxdouble));
2398 memcpy(&_mse[bin * nfibers], cpl_matrix_get_data(slice->variance),
2399 slice->nflx * sizeof(cxdouble));
2400 memcpy(&_msm[bin * ny], cpl_matrix_get_data(slice->model),
2401 slice->msize * sizeof(cxdouble));
2402
2403 memcpy(&_msy[bin * nfibers], &_my[bin * nfibers],
2404 nfibers * sizeof(cxdouble));
2405
2406
2407 /*
2408 * Reset the profile part of the design matrix
2409 */
2410
2411 cpl_matrix_fill_window(profiles, 0., 0, 0, nfibers, ny);
2412
2413 }
2414
2415
2416 /*
2417 * Compute errors of the extracted spectra from the variance
2418 */
2419
2420 cpl_image_power(mse, 0.5);
2421
2422 _giraffe_optimal_workspace_delete(workspace);
2423 workspace = NULL;
2424
2425 _giraffe_extraction_psflimits_delete(limits);
2426 limits = NULL;
2427
2428 _giraffe_extractionslice_delete(slice);
2429 slice = NULL;
2430
2431 cpl_matrix_delete(weights);
2432 weights = NULL;
2433
2434 cpl_matrix_delete(mask);
2435 mask = NULL;
2436
2437 cpl_matrix_delete(variance);
2438 variance = NULL;
2439
2440 cpl_matrix_delete(signal);
2441 signal = NULL;
2442
2443 cpl_matrix_delete(profiles);
2444 profiles = NULL;
2445
2446 giraffe_model_delete(psfmodel);
2447 psfmodel = NULL;
2448
2449 return 0;
2450
2451}
2452
2453
2478cxint
2479giraffe_extract_spectra(GiExtraction* result, GiImage* image,
2480 GiTable* fibers, GiLocalization* sloc,
2481 GiImage* bpixel, GiImage* slight,
2482 GiExtractConfig* config)
2483{
2484
2485 const cxchar *fctid = "giraffe_extract_spectra";
2486
2487
2488 cxint ns = 0;
2489 cxint nx = 0;
2490 cxint ny = 0;
2491 cxint status = 0;
2492 cxint nframes = 1;
2493
2494 cxdouble bias_ron = 0.;
2495 cxdouble bias_sigma = 0.;
2496 cxdouble dark_value = 0.;
2497 cxdouble exptime = 0.;
2498 cxdouble conad = 1.;
2499
2500 cpl_propertylist *properties;
2501
2502 cpl_image* _image = NULL;
2503 cpl_image* _locy = NULL;
2504 cpl_image* _locw = NULL;
2505 cpl_image* _spectra = NULL;
2506 cpl_image* _error = NULL;
2507 cpl_image* _npixels = NULL;
2508 cpl_image* _centroid = NULL;
2509 cpl_image* _model = NULL;
2510
2511 cpl_table* _fibers = NULL;
2512
2513
2514 /*
2515 * Preprocessing
2516 */
2517
2518 if (!result || !image || !fibers || !sloc || !config) {
2519 cpl_error_set(fctid, CPL_ERROR_NULL_INPUT);
2520 return 1;
2521 }
2522
2523
2524 if ((sloc->locy == NULL) || (sloc->locw == NULL)) {
2525 cpl_error_set(fctid, CPL_ERROR_NULL_INPUT);
2526 return 1;
2527 }
2528
2529
2530 if (result->spectra != NULL || result->error != NULL ||
2531 result->npixels != NULL || result->centroid != NULL ||
2532 result->model != NULL) {
2533 gi_warning("%s: Results structure at %p is not empty! Contents "
2534 "might be lost.", fctid, result);
2535 }
2536
2537
2538 _fibers = giraffe_table_get(fibers);
2539
2540 if (_fibers == NULL) {
2541 cpl_error_set(fctid, CPL_ERROR_DATA_NOT_FOUND);
2542 return 1;
2543 }
2544
2545
2546 if ((config->emethod == GIEXTRACT_OPTIMAL) && (sloc->psf == NULL)) {
2547 cpl_msg_error(fctid, "Missing data: PSF profile data is required "
2548 "for optimal spectrum extraction!");
2549 cpl_error_set(fctid, CPL_ERROR_NULL_INPUT);
2550
2551 return 1;
2552 }
2553
2554
2555 properties = giraffe_image_get_properties(image);
2556
2557 if (properties == NULL) {
2558 cpl_error_set(fctid, CPL_ERROR_DATA_NOT_FOUND);
2559 return 1;
2560 }
2561
2562
2563 giraffe_error_push();
2564
2565 conad = giraffe_propertylist_get_conad(properties);
2566
2567 if (cpl_error_get_code() != CPL_ERROR_NONE) {
2568 return 1;
2569 }
2570
2571 giraffe_error_pop();
2572
2573
2574 if (!cpl_propertylist_has(properties, GIALIAS_BIASERROR)) {
2575 cpl_msg_warning(fctid, "Missing bias error property (%s)! Setting "
2576 "bias error to 0.", GIALIAS_BIASERROR);
2577 bias_sigma = 0.;
2578 }
2579 else {
2580 bias_sigma = cpl_propertylist_get_double(properties, GIALIAS_BIASERROR);
2581 }
2582
2583
2584 if (config->ron > 0.) {
2585
2586 cpl_msg_info(fctid, "Setting bias RMS property (%s) to %.4g ADU",
2587 GIALIAS_BIASSIGMA, config->ron);
2588
2589 cpl_propertylist_update_double(properties, GIALIAS_BIASSIGMA,
2590 config->ron);
2591 }
2592
2593 bias_ron = giraffe_propertylist_get_ron(properties);
2594
2595
2596 if (!cpl_propertylist_has(properties, GIALIAS_DARKVALUE)) {
2597
2598 dark_value = 0.;
2599
2600 cpl_msg_warning(fctid, "Missing dark value property (%s), will be "
2601 "set to 0.!", GIALIAS_DARKVALUE);
2602 cpl_propertylist_append_double(properties, GIALIAS_DARKVALUE,
2603 dark_value);
2604
2605 }
2606 else {
2607 dark_value = cpl_propertylist_get_double(properties,
2608 GIALIAS_DARKVALUE);
2609 }
2610
2611
2612 if (!cpl_propertylist_has(properties, GIALIAS_EXPTIME)) {
2613 cpl_msg_error(fctid, "Missing exposure time property (%s)!",
2614 GIALIAS_EXPTIME);
2615 return 1;
2616 }
2617 else {
2618 exptime = cpl_propertylist_get_double(properties, GIALIAS_EXPTIME);
2619 }
2620
2621
2622 if (cpl_propertylist_has(properties, GIALIAS_DATANCOM)) {
2623 nframes = cpl_propertylist_get_int(properties, GIALIAS_DATANCOM);
2624 }
2625
2626
2627 /*
2628 * Processing
2629 */
2630
2631 /*
2632 * Convert the bias and dark errors from ADU to electrons.
2633 */
2634
2635 bias_sigma *= conad;
2636 dark_value *= conad;
2637
2638 /*
2639 * For extracting the spectra, the bias and dark corrected raw image is
2640 * converted from ADU to electrons, and, in case it is an averaged frame,
2641 * it is scaled by the number of frames which were used. This turns the
2642 * raw frame into an image of the total number of the recorded
2643 * photoelectrons.
2644 *
2645 * To compensate for that, the extracted spectra, their errors, and,
2646 * possibly the spectrum model are rescaled after the extraction step
2647 * is completed.
2648 */
2649
2650 _image = cpl_image_multiply_scalar_create(giraffe_image_get(image),
2651 nframes * conad);
2652
2653 _locy = giraffe_image_get(sloc->locy);
2654 _locw = giraffe_image_get(sloc->locw);
2655
2656 ny = cpl_image_get_size_x(_image);
2657 nx = cpl_image_get_size_y(_locw);
2658 ns = cpl_table_get_nrow(_fibers);
2659
2660
2661 switch (config->emethod) {
2662 case GIEXTRACT_SUM:
2663 {
2664
2665 cxint xsize = cpl_image_get_size_x(_image);
2666 cxint ysize = cpl_image_get_size_y(_image);
2667
2668 cxdouble ron_variance = bias_ron * bias_ron;
2669 cxdouble bias_variance = bias_sigma * bias_sigma;
2670 cxdouble dark_variance = dark_value * exptime;
2671
2672 cpl_image* bpixmap = NULL;
2673 cpl_image* variance = NULL;
2674
2675
2676 result->spectra = giraffe_image_create(CPL_TYPE_DOUBLE, ns, nx);
2677 result->error = giraffe_image_create(CPL_TYPE_DOUBLE, ns, nx);
2678 result->npixels = giraffe_image_create(CPL_TYPE_DOUBLE, ns, nx);
2679 result->centroid = giraffe_image_create(CPL_TYPE_DOUBLE, ns, nx);
2680 result->model = NULL;
2681
2682 _spectra = giraffe_image_get(result->spectra);
2683 _error = giraffe_image_get(result->error);
2684 _npixels = giraffe_image_get(result->npixels);
2685 _centroid = giraffe_image_get(result->centroid);
2686
2687 if (bpixel != NULL) {
2688
2689 bpixmap = giraffe_image_get(bpixel);
2690
2691 if (cpl_image_get_size_x(bpixmap) != xsize ||
2692 cpl_image_get_size_y(bpixmap) != ysize) {
2693
2694 cxbool crop = FALSE;
2695
2696 cpl_propertylist *p =
2698
2699 GiWindow w = {1, 1, 0, 0};
2700
2701
2702 w.x1 = cpl_image_get_size_x(bpixmap);
2703 w.y1 = cpl_image_get_size_y(bpixmap);
2704
2705 if (cpl_propertylist_has(p, GIALIAS_PRSCX)) {
2706 w.x0 += cpl_propertylist_get_int(p, GIALIAS_PRSCX);
2707 crop = TRUE;
2708 }
2709
2710 if (cpl_propertylist_has(p, GIALIAS_OVSCX)) {
2711 w.x1 -= cpl_propertylist_get_int(p, GIALIAS_OVSCX);
2712 crop = TRUE;
2713 }
2714
2715 if (cpl_propertylist_has(p, GIALIAS_PRSCY)) {
2716 w.y0 += cpl_propertylist_get_int(p, GIALIAS_PRSCY);
2717 crop = TRUE;
2718 }
2719
2720 if (cpl_propertylist_has(p, GIALIAS_OVSCY)) {
2721 w.y1 -= cpl_propertylist_get_int(p, GIALIAS_OVSCY);
2722 crop = TRUE;
2723 }
2724
2725 if ((w.x1 - w.x0 + 1) != xsize ||
2726 (w.y1 - w.y0 + 1) != ysize) {
2727 cpl_msg_error(fctid, "Invalid bad pixel map! Image "
2728 "sizes do not match!");
2729
2730 giraffe_image_delete(result->spectra);
2731 result->spectra = NULL;
2732
2733 giraffe_image_delete(result->error);
2734 result->error = NULL;
2735
2736 giraffe_image_delete(result->npixels);
2737 result->npixels = NULL;
2738
2739 giraffe_image_delete(result->centroid);
2740 result->centroid = NULL;
2741
2742 giraffe_image_delete(result->model);
2743 result->model = NULL;
2744
2745 cpl_image_delete(_image);
2746 _image = NULL;
2747
2748 return 1;
2749 }
2750
2751 if (crop == TRUE) {
2752 bpixmap = cpl_image_extract(bpixmap, w.x0, w.y0,
2753 w.x1, w.y1);
2754 }
2755
2756 }
2757
2758 }
2759
2760 if (slight != NULL) {
2761 cpl_msg_warning(fctid, "Scattered light model will be "
2762 "ignored for extraction method `SUM'");
2763 }
2764
2765 variance = cpl_image_abs_create(_image);
2766
2767 /*
2768 * Add readout noise for the raw frame, and the errors due
2769 * to bias and dark subtraction, rescaled to the number of
2770 * frames used to create the input frame.
2771 */
2772
2773 cpl_image_add_scalar(variance, nframes * (ron_variance + nframes *
2774 (bias_variance + dark_variance)));
2775
2776 status = _giraffe_extract_summation(_image, variance, _fibers,
2777 _locy, _locw, bpixmap,
2778 _spectra, _error, _npixels,
2779 _centroid);
2780
2781 cpl_image_delete(variance);
2782 if (bpixmap != giraffe_image_get(bpixel)) {
2783 cpl_image_delete(bpixmap);
2784 }
2785 bpixmap = NULL;
2786
2787 break;
2788
2789 }
2790
2791 case GIEXTRACT_OPTIMAL:
2792 {
2793
2794 cxint xsize = cpl_image_get_size_x(_image);
2795 cxint ysize = cpl_image_get_size_y(_image);
2796
2797 cxdouble v0 = 0.;
2798 cxdouble ron_variance = bias_ron * bias_ron;
2799 cxdouble bias_variance = bias_sigma * bias_sigma;
2800 cxdouble dark_variance = dark_value * exptime;
2801
2802 cpl_image* variance = NULL;
2803 cpl_image* bpixmap = NULL;
2804
2805 GiExtractOptimalConfig setup;
2806
2807
2808 result->spectra = giraffe_image_create(CPL_TYPE_DOUBLE, ns, nx);
2809 result->error = giraffe_image_create(CPL_TYPE_DOUBLE, ns, nx);
2810 result->npixels = NULL;
2811 result->model = giraffe_image_create(CPL_TYPE_DOUBLE, ny, nx);
2812 result->centroid = giraffe_image_create(CPL_TYPE_DOUBLE, ns, nx);
2813
2814 _spectra = giraffe_image_get(result->spectra);
2815 _error = giraffe_image_get(result->error);
2816 _model = giraffe_image_get(result->model);
2817 _centroid = giraffe_image_get(result->centroid);
2818
2819 setup.clip.iterations = config->psf.iterations;
2820 setup.clip.level = config->psf.sigma;
2821 setup.clip.fraction = config->optimal.fraction;
2822 setup.limits = config->optimal.wfactor < 0. ? FALSE : TRUE;
2823 setup.ewidth = CX_MAX(1., fabs(config->optimal.wfactor));
2824 setup.bkgorder = config->optimal.bkgorder;
2825 setup.exptime = exptime;
2826 setup.ron = bias_sigma;
2827 setup.dark = dark_value;
2828
2829
2830 if (bpixel != NULL) {
2831
2832 bpixmap = giraffe_image_get(bpixel);
2833
2834 if (cpl_image_get_size_x(bpixmap) != xsize ||
2835 cpl_image_get_size_y(bpixmap) != ysize) {
2836
2837 cxbool crop = FALSE;
2838
2839 cpl_propertylist *p =
2841
2842 GiWindow w = {1, 1, 0, 0};
2843
2844
2845 w.x1 = cpl_image_get_size_x(bpixmap);
2846 w.y1 = cpl_image_get_size_y(bpixmap);
2847
2848 if (cpl_propertylist_has(p, GIALIAS_PRSCX)) {
2849 w.x0 += cpl_propertylist_get_int(p, GIALIAS_PRSCX);
2850 crop = TRUE;
2851 }
2852
2853 if (cpl_propertylist_has(p, GIALIAS_OVSCX)) {
2854 w.x1 -= cpl_propertylist_get_int(p, GIALIAS_OVSCX);
2855 crop = TRUE;
2856 }
2857
2858 if (cpl_propertylist_has(p, GIALIAS_PRSCY)) {
2859 w.y0 += cpl_propertylist_get_int(p, GIALIAS_PRSCY);
2860 crop = TRUE;
2861 }
2862
2863 if (cpl_propertylist_has(p, GIALIAS_OVSCY)) {
2864 w.y1 -= cpl_propertylist_get_int(p, GIALIAS_OVSCY);
2865 crop = TRUE;
2866 }
2867
2868 if ((w.x1 - w.x0 + 1) != xsize ||
2869 (w.y1 - w.y0 + 1) != ysize) {
2870
2871 cpl_msg_error(fctid, "Invalid bad pixel map! "
2872 "Image sizes do not match!");
2873
2874 giraffe_image_delete(result->spectra);
2875 result->spectra = NULL;
2876
2877 giraffe_image_delete(result->error);
2878 result->error = NULL;
2879
2880 giraffe_image_delete(result->npixels);
2881 result->npixels = NULL;
2882
2883 giraffe_image_delete(result->centroid);
2884 result->centroid = NULL;
2885
2886 giraffe_image_delete(result->model);
2887 result->model = NULL;
2888
2889 cpl_image_delete(_image);
2890 _image = NULL;
2891
2892 return 1;
2893
2894 }
2895
2896 if (crop == TRUE) {
2897 bpixmap = cpl_image_extract(bpixmap, w.x0, w.y0,
2898 w.x1, w.y1);
2899 }
2900
2901 }
2902
2903 }
2904
2905 variance = cpl_image_new(xsize, ysize, CPL_TYPE_DOUBLE);
2906
2907 /*
2908 * Add readout noise for the raw frame, and the errors due
2909 * to bias and dark subtraction, rescaled to the number of
2910 * frames used to create the input frame.
2911 */
2912
2913 v0 = nframes * (ron_variance + nframes *
2914 (bias_variance + dark_variance));
2915
2916
2917 /*
2918 * If a scattered light map has been used, add its contribution
2919 * to the variance, rescaled to the number of raw frames used, and
2920 * converted to photoelectrons.
2921 */
2922
2923 if (slight != NULL) {
2924
2925 register cxsize i = 0;
2926 register cxsize npixels = xsize * ysize;
2927
2928 const cxdouble* _slight =
2929 cpl_image_get_data_double(giraffe_image_get(slight));
2930
2931 cxdouble* _variance = cpl_image_get_data_double(variance);
2932
2933 for (i = 0; i < npixels; i++) {
2934 _variance[i] = v0 + fabs(_slight[i]) * conad * nframes;
2935 }
2936
2937 }
2938 else {
2939
2940 register cxsize i = 0;
2941 register cxsize npixels = xsize * ysize;
2942
2943 cxdouble* _variance = cpl_image_get_data_double(variance);
2944
2945 for (i = 0; i < npixels; i++) {
2946 _variance[i] = v0;
2947 }
2948
2949 }
2950
2951
2952 status = _giraffe_extract_optimal(_image, variance, _fibers,
2953 _locy, _locw, sloc->psf,
2954 bpixmap, _spectra, _error,
2955 _model, _centroid, &setup);
2956
2957 cpl_image_delete(variance);
2958 variance = NULL;
2959
2960 if (bpixmap != giraffe_image_get(bpixel)) {
2961 cpl_image_delete(bpixmap);
2962 }
2963 bpixmap = NULL;
2964
2965 break;
2966
2967 }
2968
2969 case GIEXTRACT_HORNE:
2970 {
2971
2972 cxint xsize = cpl_image_get_size_x(_image);
2973 cxint ysize = cpl_image_get_size_y(_image);
2974
2975 cxdouble v0 = 0.;
2976 cxdouble ron_variance = bias_ron * bias_ron;
2977 cxdouble bias_variance = bias_sigma * bias_sigma;
2978 cxdouble dark_variance = dark_value * exptime;
2979
2980 cpl_image* variance = NULL;
2981 cpl_image* bpixmap = NULL;
2982
2983 GiExtractHorneConfig setup;
2984
2985
2986 result->spectra = giraffe_image_create(CPL_TYPE_DOUBLE, ns, nx);
2987 result->error = giraffe_image_create(CPL_TYPE_DOUBLE, ns, nx);
2988 result->npixels = giraffe_image_create(CPL_TYPE_DOUBLE, ns, nx);
2989 result->centroid = giraffe_image_create(CPL_TYPE_DOUBLE, ns, nx);
2990 result->model = NULL;
2991
2992 _spectra = giraffe_image_get(result->spectra);
2993 _error = giraffe_image_get(result->error);
2994 _npixels = giraffe_image_get(result->npixels);
2995 _centroid = giraffe_image_get(result->centroid);
2996
2997 setup.clip.iterations = config->psf.iterations;
2998 setup.clip.level = config->psf.sigma;
2999 setup.clip.fraction = config->horne.mingood;
3000 setup.ewidth = config->horne.ewidth;
3001 setup.exptime = exptime;
3002 setup.ron = bias_sigma;
3003 setup.dark = dark_value;
3004
3005 if (bpixel != NULL) {
3006
3007 bpixmap = giraffe_image_get(bpixel);
3008
3009 if (cpl_image_get_size_x(bpixmap) != xsize ||
3010 cpl_image_get_size_y(bpixmap) != ysize) {
3011
3012 cxbool crop = FALSE;
3013
3014 cpl_propertylist *p =
3016
3017 GiWindow w = {1, 1, 0, 0};
3018
3019
3020 w.x1 = cpl_image_get_size_x(bpixmap);
3021 w.y1 = cpl_image_get_size_y(bpixmap);
3022
3023 if (cpl_propertylist_has(p, GIALIAS_PRSCX)) {
3024 w.x0 += cpl_propertylist_get_int(p, GIALIAS_PRSCX);
3025 crop = TRUE;
3026 }
3027
3028 if (cpl_propertylist_has(p, GIALIAS_OVSCX)) {
3029 w.x1 -= cpl_propertylist_get_int(p, GIALIAS_OVSCX);
3030 crop = TRUE;
3031 }
3032
3033 if (cpl_propertylist_has(p, GIALIAS_PRSCY)) {
3034 w.y0 += cpl_propertylist_get_int(p, GIALIAS_PRSCY);
3035 crop = TRUE;
3036 }
3037
3038 if (cpl_propertylist_has(p, GIALIAS_OVSCY)) {
3039 w.y1 -= cpl_propertylist_get_int(p, GIALIAS_OVSCY);
3040 crop = TRUE;
3041 }
3042
3043 if ((w.x1 - w.x0 + 1) != xsize ||
3044 (w.y1 - w.y0 + 1) != ysize) {
3045
3046 cpl_msg_error(fctid, "Invalid bad pixel map! "
3047 "Image sizes do not match!");
3048
3049 giraffe_image_delete(result->spectra);
3050 result->spectra = NULL;
3051
3052 giraffe_image_delete(result->error);
3053 result->error = NULL;
3054
3055 giraffe_image_delete(result->npixels);
3056 result->npixels = NULL;
3057
3058 giraffe_image_delete(result->centroid);
3059 result->centroid = NULL;
3060
3061 giraffe_image_delete(result->model);
3062 result->model = NULL;
3063
3064 cpl_image_delete(_image);
3065 _image = NULL;
3066
3067 return 1;
3068
3069 }
3070
3071 if (crop == TRUE) {
3072 bpixmap = cpl_image_extract(bpixmap, w.x0, w.y0,
3073 w.x1, w.y1);
3074 }
3075
3076 }
3077
3078 }
3079
3080 variance = cpl_image_new(xsize, ysize, CPL_TYPE_DOUBLE);
3081
3082 /*
3083 * Add readout noise for the raw frame, and the errors due
3084 * to bias and dark subtraction, rescaled to the number of
3085 * frames used to create the input frame.
3086 */
3087
3088 v0 = nframes * (ron_variance + nframes *
3089 (bias_variance + dark_variance));
3090
3091
3092 /*
3093 * If a scattered light map has been used, add its contribution
3094 * to the variance, rescaled to the number of raw frames used, and
3095 * converted to photoelectrons.
3096 */
3097
3098
3099 if (slight != NULL) {
3100
3101 register cxsize i = 0;
3102 register cxsize npixels = xsize * ysize;
3103
3104 const cxdouble* _slight =
3105 cpl_image_get_data_double(giraffe_image_get(slight));
3106
3107 cxdouble* _variance = cpl_image_get_data_double(variance);
3108
3109 for (i = 0; i < npixels; i++) {
3110 _variance[i] = v0 + fabs(_slight[i]) * nframes * conad;
3111 }
3112
3113 }
3114 else {
3115
3116 register cxsize i = 0;
3117 register cxsize npixels = xsize * ysize;
3118
3119 cxdouble* _variance = cpl_image_get_data_double(variance);
3120
3121 for (i = 0; i < npixels; i++) {
3122 _variance[i] = v0;
3123 }
3124
3125 }
3126
3127
3128 status = _giraffe_extract_horne(_image, variance, _fibers,
3129 _locy, _locw, sloc->psf,
3130 bpixmap, _spectra, _error,
3131 _npixels, _centroid, &setup);
3132
3133 cpl_image_delete(variance);
3134 variance = NULL;
3135
3136 if (bpixmap != giraffe_image_get(bpixel)) {
3137 cpl_image_delete(bpixmap);
3138 }
3139 bpixmap = NULL;
3140
3141 break;
3142
3143 }
3144
3145 default:
3146 gi_message("%s: Method %d selected for spectrum extraction.",
3147 fctid, config->emethod);
3148 cpl_msg_error(fctid, "Invalid extraction method!");
3149
3150 status = 1;
3151 break;
3152 }
3153
3154 cpl_image_delete(_image);
3155 _image = NULL;
3156
3157 if (status) {
3158
3159 giraffe_image_delete(result->spectra);
3160 result->spectra = NULL;
3161
3162 giraffe_image_delete(result->error);
3163 result->error = NULL;
3164
3165 giraffe_image_delete(result->npixels);
3166 result->npixels = NULL;
3167
3168 giraffe_image_delete(result->centroid);
3169 result->centroid = NULL;
3170
3171 giraffe_image_delete(result->model);
3172 result->model = NULL;
3173
3174 cpl_msg_error(fctid, "Spectrum extraction (method %d) failed!",
3175 config->emethod);
3176
3177 cpl_image_delete(_image);
3178 _image = NULL;
3179
3180 return 1;
3181
3182 }
3183
3184
3185 /*
3186 * Postprocessing
3187 */
3188
3189
3190 /*
3191 * Rescale the spectrum extraction products to the original, averaged
3192 * input raw frame.
3193 */
3194
3195 if (result->spectra) {
3196 cpl_image_divide_scalar(giraffe_image_get(result->spectra),
3197 nframes * conad);
3198 }
3199
3200 if (result->model) {
3201 cpl_image_divide_scalar(giraffe_image_get(result->model),
3202 nframes * conad);
3203 }
3204
3205 if (result->error) {
3206 cpl_image_divide_scalar(giraffe_image_get(result->error),
3207 nframes * conad);
3208 }
3209
3210
3211 /*
3212 * Extracted spectra frame
3213 */
3214
3215 properties = giraffe_image_get_properties(image);
3216 giraffe_image_set_properties(result->spectra, properties);
3217
3218 properties = giraffe_image_get_properties(result->spectra);
3219
3220 /*
3221 * Copy some properties from the localization frame.
3222 */
3223
3224 // FIXME: Is this really needed? (RP)
3225
3226 giraffe_propertylist_update(properties,
3227 giraffe_image_get_properties(sloc->locy),
3228 "^ESO PRO LOC.*");
3229
3230 cpl_propertylist_set_int(properties, GIALIAS_NAXIS1,
3231 cpl_image_get_size_x(_spectra));
3232 cpl_propertylist_set_int(properties, GIALIAS_NAXIS2,
3233 cpl_image_get_size_y(_spectra));
3234
3235 cpl_propertylist_set_int(properties, GIALIAS_BITPIX, -32);
3236 cpl_propertylist_set_double(properties, GIALIAS_BZERO, 0.);
3237 cpl_propertylist_set_double(properties, GIALIAS_BSCALE, 1.);
3238
3239 cpl_propertylist_update_int(properties, GIALIAS_NFIBERS,
3240 cpl_image_get_size_x(_spectra));
3241
3242 cpl_propertylist_append_int(properties, GIALIAS_EXT_NX,
3243 cpl_image_get_size_y(_spectra));
3244 cpl_propertylist_append_int(properties, GIALIAS_EXT_NS,
3245 cpl_image_get_size_x(_spectra));
3246
3247 switch (config->emethod) {
3248 case GIEXTRACT_SUM:
3249 cpl_propertylist_append_string(properties, GIALIAS_EXT_METHOD,
3250 "SUM");
3251 cpl_propertylist_set_comment(properties, GIALIAS_EXT_METHOD,
3252 "Spectrum extraction method");
3253 break;
3254
3255 case GIEXTRACT_HORNE:
3256 {
3257
3258 cpl_propertylist_append_string(properties, GIALIAS_EXT_METHOD,
3259 "HORNE");
3260 cpl_propertylist_set_comment(properties, GIALIAS_EXT_METHOD,
3261 "Spectrum extraction method");
3262
3263 cpl_propertylist_append_string(properties, GIALIAS_EXTPSF_MODEL,
3264 config->psf.model);
3265 cpl_propertylist_set_comment(properties, GIALIAS_EXTPSF_MODEL,
3266 "PSF model used");
3267 cpl_propertylist_append_double(properties, GIALIAS_EXTPSF_SIGMA,
3268 config->psf.sigma);
3269 cpl_propertylist_set_comment(properties, GIALIAS_EXTPSF_SIGMA,
3270 "PSF fit sigma clipping threshold");
3271 cpl_propertylist_append_int(properties, GIALIAS_EXTPSF_NITER,
3272 config->psf.iterations);
3273 cpl_propertylist_set_comment(properties, GIALIAS_EXTPSF_NITER,
3274 "PSF fit maximum number of "
3275 "iterations");
3276
3277 cpl_propertylist_append_int(properties, GIALIAS_EXTHRN_EWIDTH,
3278 config->horne.ewidth);
3279 cpl_propertylist_set_comment(properties, GIALIAS_EXTHRN_EWIDTH,
3280 "Number of extra pixels used");
3281 cpl_propertylist_append_int(properties, GIALIAS_EXTHRN_MINGOOD,
3282 config->horne.mingood);
3283 cpl_propertylist_set_comment(properties, GIALIAS_EXTHRN_MINGOOD,
3284 "Minimum number of pixels to keep");
3285
3286
3287 break;
3288 }
3289
3290 case GIEXTRACT_OPTIMAL:
3291 cpl_propertylist_append_string(properties, GIALIAS_EXT_METHOD,
3292 "OPTIMAL");
3293 cpl_propertylist_set_comment(properties, GIALIAS_EXT_METHOD,
3294 "Spectrum extraction method");
3295
3296 cpl_propertylist_append_string(properties, GIALIAS_EXTPSF_MODEL,
3297 config->psf.model);
3298 cpl_propertylist_set_comment(properties, GIALIAS_EXTPSF_MODEL,
3299 "PSF model used");
3300 cpl_propertylist_append_double(properties, GIALIAS_EXTPSF_SIGMA,
3301 config->psf.sigma);
3302 cpl_propertylist_set_comment(properties, GIALIAS_EXTPSF_SIGMA,
3303 "PSF fit sigma clipping threshold");
3304 cpl_propertylist_append_int(properties, GIALIAS_EXTPSF_NITER,
3305 config->psf.iterations);
3306 cpl_propertylist_set_comment(properties, GIALIAS_EXTPSF_NITER,
3307 "PSF fit maximum number of "
3308 "iterations");
3309
3310 cpl_propertylist_append_double(properties, GIALIAS_EXTOPT_FRACTION,
3311 config->optimal.fraction);
3312 cpl_propertylist_set_comment(properties, GIALIAS_EXTOPT_FRACTION,
3313 "Minimum fraction of pixels used.");
3314 cpl_propertylist_append_double(properties, GIALIAS_EXTOPT_WFACTOR,
3315 config->optimal.wfactor);
3316 cpl_propertylist_set_comment(properties, GIALIAS_EXTOPT_WFACTOR,
3317 "Multiple of the fiber PSF half "
3318 "width used for spectrum "
3319 "extraction.");
3320 cpl_propertylist_append_int(properties, GIALIAS_EXTOPT_BGORDER,
3321 config->optimal.bkgorder);
3322 cpl_propertylist_set_comment(properties, GIALIAS_EXTOPT_BGORDER,
3323 "Order of the background polynomial "
3324 "model along the spatial direction.");
3325
3326 break;
3327
3328 default:
3329 break;
3330 }
3331
3332 cpl_propertylist_update_string(properties, GIALIAS_GIRFTYPE, "EXTSP");
3333 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
3334 "Extracted spectra");
3335
3336
3337 /*
3338 * Extracted spectra errors frame
3339 */
3340
3341 giraffe_image_set_properties(result->error, properties);
3342 properties = giraffe_image_get_properties(result->error);
3343
3344 cpl_propertylist_set_string(properties, GIALIAS_GIRFTYPE, "EXTERRS");
3345 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
3346 "Extracted spectra errors");
3347
3348
3349 /*
3350 * Extracted spectra centroids frame
3351 */
3352
3353 giraffe_image_set_properties(result->centroid, properties);
3354 properties = giraffe_image_get_properties(result->centroid);
3355
3356 cpl_propertylist_set_string(properties, GIALIAS_GIRFTYPE, "EXTYCEN");
3357 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
3358 "Extracted spectra centroids");
3359
3360
3361 /*
3362 * Extracted spectra npixels frame
3363 */
3364
3365 if (result->npixels != NULL) {
3366 giraffe_image_set_properties(result->npixels, properties);
3367 properties = giraffe_image_get_properties(result->npixels);
3368
3369 cpl_propertylist_set_string(properties, GIALIAS_GIRFTYPE, "EXTNPIX");
3370 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
3371 "Extracted spectra npixels");
3372 }
3373
3374
3375 /*
3376 * Model spectra frame
3377 */
3378
3379 if (result->model != NULL) {
3380 giraffe_image_set_properties(result->model, properties);
3381 properties = giraffe_image_get_properties(result->model);
3382
3383 cpl_propertylist_set_string(properties, GIALIAS_GIRFTYPE, "EXTMODEL");
3384 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
3385 "Model spectra used for extraction");
3386 }
3387
3388 return 0;
3389
3390}
3391
3392
3403GiExtractConfig*
3404giraffe_extract_config_create(cpl_parameterlist* list)
3405{
3406
3407 const cxchar* s;
3408 cpl_parameter* p;
3409
3410 GiExtractConfig* config = NULL;
3411
3412
3413 if (!list) {
3414 return NULL;
3415 }
3416
3417 config = cx_calloc(1, sizeof *config);
3418
3419 p = cpl_parameterlist_find(list, "giraffe.extraction.method");
3420 s = cpl_parameter_get_string(p);
3421 if (!strcmp(s, "OPTIMAL")) {
3422 config->emethod = GIEXTRACT_OPTIMAL;
3423 }
3424 else if (!strcmp(s, "HORNE")) {
3425 config->emethod = GIEXTRACT_HORNE;
3426 }
3427 else {
3428 config->emethod = GIEXTRACT_SUM;
3429 }
3430
3431 p = cpl_parameterlist_find(list, "giraffe.extraction.ron");
3432 config->ron = cpl_parameter_get_double(p);
3433
3434 p = cpl_parameterlist_find(list, "giraffe.extraction.psf.model");
3435 config->psf.model = cx_strdup(cpl_parameter_get_string(p));
3436
3437 p = cpl_parameterlist_find(list, "giraffe.extraction.psf.sigma");
3438 config->psf.sigma = cpl_parameter_get_double(p);
3439
3440 p = cpl_parameterlist_find(list, "giraffe.extraction.psf.iterations");
3441 config->psf.iterations = cpl_parameter_get_int(p);
3442
3443
3444 p = cpl_parameterlist_find(list, "giraffe.extraction.horne.extrawidth");
3445 config->horne.ewidth = cpl_parameter_get_int(p);
3446
3447 p = cpl_parameterlist_find(list, "giraffe.extraction.horne.mingood");
3448 config->horne.mingood = cpl_parameter_get_double(p);
3449
3450
3451 p = cpl_parameterlist_find(list, "giraffe.extraction.optimal.fraction");
3452 config->optimal.fraction = cpl_parameter_get_double(p);
3453
3454 p = cpl_parameterlist_find(list, "giraffe.extraction.optimal.wfactor");
3455 config->optimal.wfactor = cpl_parameter_get_double(p);
3456
3457 p = cpl_parameterlist_find(list, "giraffe.extraction.optimal.bkgorder");
3458 config->optimal.bkgorder = cpl_parameter_get_int(p);
3459
3460 return config;
3461
3462}
3463
3464
3477void
3478giraffe_extract_config_destroy(GiExtractConfig* config)
3479{
3480
3481 if (config) {
3482
3483 if (config->psf.model) {
3484 cx_free(config->psf.model);
3485 }
3486
3487 cx_free(config);
3488
3489 }
3490
3491 return;
3492
3493}
3494
3495
3507void
3508giraffe_extract_config_add(cpl_parameterlist* list)
3509{
3510
3511 cpl_parameter* p = NULL;
3512
3513
3514 if (list == NULL) {
3515 return;
3516 }
3517
3518 p = cpl_parameter_new_enum("giraffe.extraction.method",
3519 CPL_TYPE_STRING,
3520 "Extraction method: 'SUM', 'HORNE' or "
3521 "'OPTIMAL'",
3522 "giraffe.extraction",
3523 "SUM", 3, "SUM", "OPTIMAL", "HORNE");
3524 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "extr-method");
3525 cpl_parameterlist_append(list, p);
3526
3527
3528 p = cpl_parameter_new_value("giraffe.extraction.ron",
3529 CPL_TYPE_DOUBLE,
3530 "New bias sigma (RON) value for "
3531 "bias and dark "
3532 "corrected image",
3533 "giraffe.extraction",
3534 -1.);
3535 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "extr-ron");
3536 cpl_parameterlist_append(list, p);
3537
3538
3539 p = cpl_parameter_new_enum("giraffe.extraction.psf.model",
3540 CPL_TYPE_STRING,
3541 "PSF profile model: `psfexp', `psfexp2'",
3542 "giraffe.extraction.psf",
3543 "psfexp2", 2, "psfexp", "psfexp2");
3544 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "extr-psfmodel");
3545 cpl_parameterlist_append(list, p);
3546
3547
3548 p = cpl_parameter_new_value("giraffe.extraction.psf.sigma",
3549 CPL_TYPE_DOUBLE,
3550 "Sigma clippging threshold used for "
3551 "rejecting data points during PSF fitting "
3552 "(Horne's sigma). It is used to reject bad "
3553 "pixels and cosmics.",
3554 "giraffe.extraction.psf",
3555 7.);
3556 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "extr-psfsigma");
3557 cpl_parameterlist_append(list, p);
3558
3559
3560 p = cpl_parameter_new_value("giraffe.extraction.psf.iterations",
3561 CPL_TYPE_INT,
3562 "Maximum number of iterations used for "
3563 "fitting the PSF profile.",
3564 "giraffe.extraction.psf",
3565 2);
3566 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "extr-psfniter");
3567 cpl_parameterlist_append(list, p);
3568
3569
3570 p = cpl_parameter_new_value("giraffe.extraction.horne.extrawidth",
3571 CPL_TYPE_INT,
3572 "Horne extraction method: Number of "
3573 "extra pixels added to the fiber "
3574 "half-width.",
3575 "giraffe.extraction.horne",
3576 2);
3577 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "extr-hewidth");
3578 cpl_parameterlist_append(list, p);
3579
3580
3581 p = cpl_parameter_new_value("giraffe.extraction.horne.mingood",
3582 CPL_TYPE_INT,
3583 "Horne extraction method: Minimum number of "
3584 "points used for the profile fit. It sets "
3585 "the lower limit of data points for the "
3586 "pixel rejection.",
3587 "giraffe.extraction.horne",
3588 3);
3589 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "extr-hmingood");
3590 cpl_parameterlist_append(list, p);
3591
3592
3593 p = cpl_parameter_new_range("giraffe.extraction.optimal.fraction",
3594 CPL_TYPE_DOUBLE,
3595 "Optimal extraction method: Minimum fraction "
3596 "of the data points used for fitting the "
3597 "fiber profiles. It sets the lower limit "
3598 "for the pixel rejection.",
3599 "giraffe.extraction.optimal",
3600 0.9, 0.0, 1.0);
3601 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "extr-omfrac");
3602 cpl_parameterlist_append(list, p);
3603
3604
3605 p = cpl_parameter_new_value("giraffe.extraction.optimal.wfactor",
3606 CPL_TYPE_DOUBLE,
3607 "Optimal extraction method: Factor by which "
3608 "the fiber PSF half width is multiplied. "
3609 "Adjacent spectra within this area are "
3610 "assumed to affect the spectrum being "
3611 "extracted.",
3612 "giraffe.extraction.optimal",
3613 3.);
3614 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "extr-owfactor");
3615 cpl_parameterlist_append(list, p);
3616
3617
3618 p = cpl_parameter_new_value("giraffe.extraction.optimal.bkgorder",
3619 CPL_TYPE_INT,
3620 "Optimal extraction method: Order of the "
3621 "polynomial background model, which is "
3622 "fitted for each wavelength bin along the "
3623 "spatial direction.",
3624 "giraffe.extraction.optimal",
3625 2);
3626 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "extr-obkgorder");
3627 cpl_parameterlist_append(list, p);
3628
3629
3630 return;
3631
3632}
cxint giraffe_array_sort(cxdouble *array, cxsize size)
Sorts an array in ascending order.
Definition: giarray.c:169
void giraffe_extract_config_add(cpl_parameterlist *list)
Adds parameters for the spectrum extraction.
Definition: giextract.c:3508
cxint giraffe_extract_spectra(GiExtraction *result, GiImage *image, GiTable *fibers, GiLocalization *sloc, GiImage *bpixel, GiImage *slight, GiExtractConfig *config)
Extracts the spectra from a preprocessed frame.
Definition: giextract.c:2479
GiExtractConfig * giraffe_extract_config_create(cpl_parameterlist *list)
Creates a setup structure for the spectrum extraction.
Definition: giextract.c:3404
void giraffe_extract_config_destroy(GiExtractConfig *config)
Destroys a spectrum extraction setup structure.
Definition: giextract.c:3478
const cxchar * giraffe_fiberlist_query_index(const cpl_table *fibers)
Query a fiber list for the name of the fiber reference index column.
cpl_image * giraffe_image_get(const GiImage *self)
Gets the image data.
Definition: giimage.c:218
cpl_propertylist * giraffe_image_get_properties(const GiImage *self)
Get the properties of an image.
Definition: giimage.c:282
void giraffe_image_delete(GiImage *self)
Destroys an image.
Definition: giimage.c:181
GiImage * giraffe_image_create(cpl_type type, cxint nx, cxint ny)
Creates an image container of a given type.
Definition: giimage.c:95
cxint giraffe_image_set_properties(GiImage *self, cpl_propertylist *properties)
Attaches a property list to an image.
Definition: giimage.c:312
void gi_warning(const cxchar *format,...)
Log a warning.
Definition: gimessages.c:119
void gi_message(const cxchar *format,...)
Log a normal message.
Definition: gimessages.c:146
cpl_table * giraffe_table_get(const GiTable *self)
Get the table data from a Giraffe table.
Definition: gitable.c:433
cxdouble giraffe_propertylist_get_conad(const cpl_propertylist *properties)
Retrieve the ADU to electrons conversion factor from the given properties.
Definition: giutils.c:1471
cxint giraffe_propertylist_update(cpl_propertylist *self, cpl_propertylist *properties, const cxchar *regexp)
Update a property list.
Definition: giutils.c:1015
cxdouble giraffe_propertylist_get_ron(const cpl_propertylist *properties)
Retrieve the read-out noise from the given properties.
Definition: giutils.c:1555

This file is part of the GIRAFFE Pipeline Reference Manual 2.18.3.
Documentation copyright © 2002-2006 European Southern Observatory.
Generated on Fri May 23 2025 16:13:45 by doxygen 1.9.6 written by Dimitri van Heesch, © 1997-2004