X-shooter Pipeline Reference Manual 3.8.15
xsh_opt_extract.c
Go to the documentation of this file.
1/* *
2 * This file is part of the ESO X-shooter Pipeline *
3 * Copyright (C) 2006 European Southern Observatory *
4 * *
5 * This library 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, 51 Franklin St, Fifth Floor, Boston, MA 02111-1307 USA *
18 * */
19
20
21/*
22 * $Author: amodigli $
23 * $Date: 2012-12-18 16:10:42 $
24 * $Revision: 1.110 $
25*/
26
27#ifdef HAVE_CONFIG_H
28#include <config.h>
29#endif
30/*----------------------------------------------------------------------------*/
35/*----------------------------------------------------------------------------*/
38/*-----------------------------------------------------------------------------
39 Includes
40 -----------------------------------------------------------------------------*/
41
42#include <math.h>
43#include <xsh_drl.h>
44
45#include <xsh_badpixelmap.h>
46#include <xsh_data_rec.h>
47#include <xsh_data_pre.h>
48#include <xsh_utils_wrappers.h>
49#include <xsh_data_order.h>
50#include <xsh_data_wavesol.h>
53#include <xsh_dfs.h>
54#include <xsh_pfits.h>
55#include <xsh_error.h>
56#include <xsh_msg.h>
57#include <xsh_fit.h>
58#include <xsh_model_io.h>
59#include <xsh_model_kernel.h>
60#include <xsh_rectify.h>
61#include <cpl.h>
62#include <gsl/gsl_sf_erf.h>
63#include <xsh_blaze.h>
64/*-----------------------------------------------------------------------------
65 Functions prototypes
66 -----------------------------------------------------------------------------*/
67#define SLIT_USEFUL_WINDOW_FACTOR 0.80
68#define FIT_FWHM_LIMIT 30
69#define OPT_EXTRACT_SLIT_SIZE 80
70
71#define FLUX_MODE 0
72#define WAVEMAP_MODE 1
73
74#define REGDEBUG_PIXELSIZE 0
75#define REGDEBUG_INTEGRATE 0
76/*-----------------------------------------------------------------------------
77 *
78 Prototype
79 -----------------------------------------------------------------------------*/
80
81static void xsh_image_gaussian_fit_y( cpl_image* img, int chunk_size,
82 int deg_poly, int oversample,
83 cpl_polynomial** center, cpl_polynomial** height,
84 cpl_polynomial** width, cpl_polynomial** offset);
85
86static cpl_image* xsh_image_create_gaussian_image( cpl_image *src,
87 cpl_polynomial* centerp, cpl_polynomial* heightp, cpl_polynomial* widthp,
88 cpl_polynomial* offsetp);
89
90static cpl_image* xsh_image_create_model_image( cpl_image *src_img,
91 double *x_data, double *y_data, double kappa, int niter, double frac_min);
92
93static cpl_vector* xsh_vector_integrate( int biny, int absorder, int oversample,
94 cpl_vector *ref_pos, double step,
95 cpl_vector *ref_values, cpl_vector *err_values, cpl_vector *qual_values,
96 cpl_vector *new_pos,
97 cpl_vector **spectrum_err, cpl_vector **spectrum_qual);
98
99/*-----------------------------------------------------------------------------
100 *
101 Implementation
102 -----------------------------------------------------------------------------*/
103
104static void xsh_interpolate_spectrum( int biny, int oversample,
105 int absorder, double lambda_step,
106 cpl_vector *init_pos,
107 cpl_vector *std_flux, cpl_vector *std_err, cpl_vector *std_qual,
108 cpl_vector *opt_flux, cpl_vector *opt_err, cpl_vector *opt_qual,
109 cpl_vector **res_pos,
110 cpl_vector **res_std_flux, cpl_vector **res_std_err, cpl_vector **res_std_qual,
111 cpl_vector **res_opt_flux, cpl_vector **res_opt_err, cpl_vector **res_opt_qual)
112{
113
114 int init_size;
115 double lambda_min, lambda_max;
116 double binlambda_min, binlambda_max;
117 int mult_min, mult_max;
118 int res_size, i;
119
120 XSH_ASSURE_NOT_NULL( init_pos);
121 XSH_ASSURE_NOT_NULL( std_flux);
122 XSH_ASSURE_NOT_NULL( std_err);
123 XSH_ASSURE_NOT_NULL( std_qual);
124 XSH_ASSURE_NOT_NULL( opt_flux);
125 XSH_ASSURE_NOT_NULL( opt_err);
126 XSH_ASSURE_NOT_NULL( opt_qual);
127 XSH_ASSURE_NOT_NULL( res_pos);
128 XSH_ASSURE_NOT_NULL( res_std_flux);
129 XSH_ASSURE_NOT_NULL( res_std_err);
130 XSH_ASSURE_NOT_NULL( res_std_qual);
131 XSH_ASSURE_NOT_NULL( res_opt_flux);
132 XSH_ASSURE_NOT_NULL( res_opt_err);
133 XSH_ASSURE_NOT_NULL( res_opt_qual);
134
135 check( init_size = cpl_vector_get_size( init_pos));
136
137 check( lambda_min = cpl_vector_get( init_pos, 0));
138 check( lambda_max = cpl_vector_get( init_pos, init_size-1));
139
140 /* round value to be multiple of lambda_step */
141 mult_min = (int) ceil(lambda_min/lambda_step)+1;
142 binlambda_min = mult_min*lambda_step;
143
144 mult_max = (int) floor(lambda_max/lambda_step)-1;
145 binlambda_max = mult_max*lambda_step;
146
147 xsh_msg("lambda_min %f lambda_max %f TEST %f %f : step %f", lambda_min, lambda_max, binlambda_min, binlambda_max, lambda_step);
148
149 res_size = mult_max-mult_min+1;
150
151 check( *res_pos = cpl_vector_new( res_size));
152
153 for( i=0; i< res_size; i++){
154 check( cpl_vector_set( *res_pos, i,binlambda_min+i*lambda_step));
155 }
156
157 check( *res_std_flux = xsh_vector_integrate(
158 biny, absorder, oversample, init_pos, lambda_step,
159 std_flux, std_err, std_qual, *res_pos, res_std_err, res_std_qual));
160 check( *res_opt_flux = xsh_vector_integrate(
161 biny, absorder, oversample, init_pos, lambda_step,
162 opt_flux, opt_err, opt_qual, *res_pos, res_opt_err, res_opt_qual));
163
164 cleanup:
165 if (cpl_error_get_code() != CPL_ERROR_NONE){
166 xsh_free_vector( res_pos);
167 xsh_free_vector( res_std_flux);
168 xsh_free_vector( res_std_err);
169 xsh_free_vector( res_std_qual);
170 xsh_free_vector( res_opt_flux);
171 xsh_free_vector( res_opt_err);
172 xsh_free_vector( res_opt_qual);
173 }
174 return;
175}
176
177
178cpl_image* xsh_optextract_produce_model( cpl_image* s2Dby1D_img,
179 int method, int chunk_ovsamp_size, int deg_poly, int oversample,
180 double *extract_x_data, double *extract_y_data, int abs_order,
181 double kappa, int niter, double frac_min){
182
183 cpl_polynomial *center_gauss_poly = NULL;
184 cpl_polynomial *height_gauss_poly = NULL;
185 cpl_polynomial *offset_gauss_poly = NULL;
186 cpl_polynomial *width_gauss_poly = NULL;
187 cpl_image *model_img = NULL;
188 int model_x, model_y, ix, iy;
189 double *model_data = NULL;
190
191 XSH_ASSURE_NOT_NULL( s2Dby1D_img);
192 XSH_ASSURE_NOT_NULL( extract_x_data);
193 XSH_ASSURE_NOT_NULL( extract_y_data);
194
195 xsh_msg_dbg_high( "Produce model for order %d", abs_order);
196
197 if ( method == GAUSS_METHOD){
198 check( xsh_image_gaussian_fit_y( s2Dby1D_img, chunk_ovsamp_size,
199 deg_poly, oversample,
200 &center_gauss_poly, &height_gauss_poly, &width_gauss_poly,
201 &offset_gauss_poly));
202
203 check( model_img = xsh_image_create_gaussian_image( s2Dby1D_img,
204 center_gauss_poly,
205 height_gauss_poly, width_gauss_poly, offset_gauss_poly));
206
207#if REGDEBUG_GAUSSIAN
209 char poly_name[256];
210 FILE* poly_file = NULL;
211 int nx = cpl_image_get_size_x( model_img);
212 int i;
213
214 sprintf( poly_name, "gaussian_center_poly_%d.dat", abs_order);
215 poly_file = fopen( poly_name, "w");
216
217 for(i=0; i< nx; i++){
218 fprintf( poly_file, "%d %d\n", i+1,
219 (int) cpl_polynomial_eval_1d( center_gauss_poly, i, NULL)+1);
220 }
221 fclose( poly_file);
222
223 sprintf( poly_name, "gaussian_sigma_poly_%d.dat", abs_order);
224 poly_file = fopen( poly_name, "w");
225
226 for(i=0; i< nx; i++){
227 fprintf( poly_file, "%d %d\n", i+1,
228 (int) cpl_polynomial_eval_1d( width_gauss_poly, i, NULL)+1);
229 }
230 fclose( poly_file);
231 }
232#endif
233 }
234 else{
235 check( model_img = xsh_image_create_model_image( s2Dby1D_img,
236 extract_x_data, extract_y_data, kappa, niter, frac_min));
237 }
238
239 /* force positivity */
240 model_x = cpl_image_get_size_x( model_img);
241 model_y = cpl_image_get_size_y( model_img);
242 model_data = cpl_image_get_data_double( model_img);
243
244 for( ix=0; ix < model_x; ix++){
245 double psum =0.0;
246
247 for( iy=0; iy < model_y; iy++){
248 int ipos = iy*model_x+ix;
249
250 if ( model_data[ipos] < 0){
251 model_data[ipos] = 0;
252 }
253 psum += model_data[ipos];
254 }
255 /* normalize */
256 for( iy=0; iy < model_y; iy++){
257 int ipos = iy*model_x+ix;
258
259 model_data[ipos] /= psum;
260 }
261 }
262
263 cleanup:
264 if (cpl_error_get_code() != CPL_ERROR_NONE){
265 xsh_free_image( &model_img);
266 }
267 xsh_free_polynomial( &center_gauss_poly);
268 xsh_free_polynomial( &width_gauss_poly);
269 xsh_free_polynomial( &height_gauss_poly);
270 xsh_free_polynomial( &offset_gauss_poly);
271 return model_img;
272}
273
274
275
276
277
278 /* SLIT aproximative position */
279static void xsh_object_localize( cpl_frame *slitmap_frame,
280 cpl_frame *loc_frame, int oversample, int box_hsize,
281 int nlambdas, int ny_extract, double *extract_x_data,
282 double *extract_y_data, int *ymin, int *ymax)
283{
284 int u;
285 double scen=0.0;
286 cpl_image *slitmap = NULL;
287 float *slit_data = NULL;
288 int slit_nx;
289 int ny_extract_min =-1, ny_extract_max=-1, ny_extract_cen=-1;
290 int mid_lambda;
291 double diff_s_min = 50.0;
292 xsh_localization *loc = NULL;
293 int ref_ov=0;
294
295 XSH_ASSURE_NOT_NULL( slitmap_frame);
296 XSH_ASSURE_NOT_NULL( loc_frame);
297 XSH_ASSURE_NOT_NULL( extract_x_data);
298 XSH_ASSURE_NOT_NULL( extract_y_data);
299 XSH_ASSURE_NOT_NULL( ymin);
300 XSH_ASSURE_NOT_NULL( ymax);
301
302 check( slitmap = cpl_image_load( cpl_frame_get_filename( slitmap_frame),
303 CPL_TYPE_FLOAT, 0, 0));
304 check( loc = xsh_localization_load( loc_frame));
305
306 check( slit_nx = cpl_image_get_size_x( slitmap));
307 check( slit_data = cpl_image_get_data_float( slitmap));
308
309 check( scen = cpl_polynomial_eval_1d( loc->cenpoly, 0, NULL));
310
311 mid_lambda = 0.5*nlambdas;
312
313 for( u=1; u < ny_extract; u++){
314 int x,y;
315 double s;
316 double diff;
317
318 x = (int) (extract_x_data[u*nlambdas+mid_lambda]+0.5);
319 y = (int) (extract_y_data[u*nlambdas+mid_lambda]+0.5);
320 s = slit_data[x+y*slit_nx];
321
322 diff = fabs(s-scen);
323
324 if (diff < diff_s_min){
325 ny_extract_cen= u;
326 diff_s_min = diff;
327 }
328 }
329
330 /* and the over-sample of reference pixel */
331 ref_ov += floor(0.5*(oversample-1));
332
333 ny_extract_min= ny_extract_cen -box_hsize-ref_ov;
334 ny_extract_max= ny_extract_cen +box_hsize+ref_ov;
335
336 if ( ny_extract_min < 0){
337 ny_extract_min = 0;
338 }
339 if ( ny_extract_max >= ny_extract){
340 ny_extract_max = ny_extract-1;
341 }
342 *ymin = ny_extract_min;
343 *ymax = ny_extract_max;
344
345 cleanup:
346 xsh_free_image( &slitmap);
348 return;
349}
350
351
352static int xsh_interpolate_linear(float *fluxtab, float *errtab, int *qualtab,
353 int nx, int ny, float pos_x,
354 float pos_y, double *flux, double *err, int* qual, int mode){
355
356 float f00=0.0, f10=0.0, f01=0.0, f11=0.0;
357 float e00=0.0, e10=0.0, e01=0.0, e11=0.0;
358 int intx, inty;
359 float fx=0.0, fy=0.0;
360 //int nbval = 0;
361 int qual_comb = -1;
362 double A1, A2, A3, A4;
363 int ret=0;
364
365 intx = (int) pos_x;
366 inty = (int) pos_y;
367
368 XSH_ASSURE_NOT_ILLEGAL( intx >= 0 && intx <nx);
369 XSH_ASSURE_NOT_ILLEGAL( inty >= 0 && inty <ny);
370 XSH_ASSURE_NOT_NULL( flux);
372 int pix=intx+inty*nx;
373 int pix_plus_1=pix+1;
374 /* combine the bad pixels but flag them */
375 f00 = fluxtab[pix];
376 e00 = errtab[pix];
377 qual_comb = qualtab[pix];
378
379 if ( (intx +1) < nx){
380 f10 = fluxtab[pix_plus_1];
381 e10 = errtab[pix_plus_1];
382 fx = pos_x-intx;
383 qual_comb |= qualtab[pix_plus_1];
384 }
385 if ( (inty +1) < ny){
386 f01 = fluxtab[pix+nx];
387 e01 = errtab[pix+nx];
388 fy = pos_y-inty;
389 qual_comb |= qualtab[pix+nx];
390 }
391 if ( ((intx +1) < nx) && ((inty +1) < ny)){
392 f11 = fluxtab[pix_plus_1+nx];
393 e11 = errtab[pix_plus_1+nx];
394 qual_comb |= qualtab[pix_plus_1+nx];
395 }
396
397 if ( mode == WAVEMAP_MODE &&
398 ( f00 == 0.0 || f10 == 0.0 || f01 == 0 || f11 == 0)){
399 xsh_msg_dbg_medium("pixel %f, %f at zero, interpolate with "\
400 "(%d,%d)%f, (%d,%d)%f (%d,%d)%f, (%d,%d)%f", pos_x, pos_y,
401 intx, inty, f00, intx+1, inty, f10, intx, inty+1, f01,
402 intx+1, inty+1, f11);
403 ret = 1;
404 }
405 A1 = (1-fx)*(1-fy);
406 A2 = fx*(1-fy);
407 A3 = (1-fx)*fy;
408 A4 = fx*fy;
409
410 *flux = f00*A1 + f10*A2 + f01*A3 + f11*A4;
411
412 *err = sqrt( A1*A1*e00*e00+A2*A2*e10*e10+A3*A3*e01*e01+A4*A4*e11*e11);
413 *qual = qual_comb;
414
415 cleanup:
416 return ret;
417}
418/*****************************************************************************/
419
420/*****************************************************************************/
455static void xsh_wavemap_lambda_range( cpl_frame *wavemap_frame,
456 cpl_frame *slitmap_frame,
457 int starty,
458 int endy, int oversample, xsh_order_list *order_list, int iorder,double *xtab,
459 double *ytab, int order, xsh_spectralformat_list *spectralformat,
460 double *lambdastab, double *slitstab, int* sizetab, xsh_instrument *instr)
461{
462 cpl_image* wavemap = NULL;
463 cpl_image *slitmap = NULL;
464 int* qual = NULL;
465 float *wavemap_data = NULL;
466 float *slitmap_data = NULL;
467 int i, nx=0, ny=0, qual_val=0;
468 double lambda=0.0, lambda_next=0.0, lambda_old=0, err=0;
469 double lambda_sf_max=0;
470 double slit= -1;
471 float y_ovs=0, x=0, y=0;
472 int size=0;
473 int ycen=0;
474 double xcen=0;
475 double lambda_min=0, lambda_max=0;
476 const char* wave_filename = NULL;
477 const char* slitmap_name = NULL;
478
479 XSH_ASSURE_NOT_NULL( wavemap_frame);
480 XSH_ASSURE_NOT_NULL( sizetab);
481 XSH_ASSURE_NOT_NULL( order_list);
482 XSH_ASSURE_NOT_NULL( xtab);
483 XSH_ASSURE_NOT_NULL( ytab);
484 XSH_ASSURE_NOT_NULL( spectralformat);
485 XSH_ASSURE_NOT_NULL( lambdastab);
486
487 check( wave_filename = cpl_frame_get_filename( wavemap_frame));
488 check( slitmap_name = cpl_frame_get_filename( slitmap_frame));
489 xsh_msg_dbg_medium("Wave map name %s", wave_filename);
490 xsh_msg_dbg_medium("Slit map name %s", slitmap_name);
491
492 check( wavemap = cpl_image_load( wave_filename,
493 CPL_TYPE_FLOAT, 0, 0));
494 check( slitmap = cpl_image_load( slitmap_name,
495 CPL_TYPE_FLOAT, 0, 0));
496
497 check( nx = cpl_image_get_size_x( wavemap));
498 check( ny = cpl_image_get_size_y( wavemap));
499
500 XSH_CALLOC(qual, int, nx*ny);
501 check( wavemap_data = cpl_image_get_data_float( wavemap));
502 check( slitmap_data = cpl_image_get_data_float( slitmap));
503
505 spectralformat, order));
507 spectralformat, order));
508
509 xsh_msg_dbg_high( "SPECTRAL FORMAT lambda min %f max %f", lambda_min, lambda_sf_max);
510
511 lambda_max = lambda_sf_max;
512
513
514 if ( xsh_instrument_get_arm( instr) !=XSH_ARM_NIR){
515
516 ycen = starty-1;
517
518 xsh_msg_dbg_high("DBG ycen %d", ycen);
519
520 while( !((lambda >= lambda_min) && (lambda_next > lambda)) ){
521 double start_x_next;
522
523 ycen++;
524 check( xcen = xsh_order_list_eval( order_list, order_list->list[iorder].cenpoly,
525 ycen));
526 xsh_msg_dbg_high("DBG xcen %f ycen %d", xcen, ycen);
527 check( xsh_interpolate_linear( wavemap_data, wavemap_data, qual,
528 nx, ny, xcen-1, ycen-1, &lambda, &err, &qual_val, WAVEMAP_MODE));
529 xsh_msg_dbg_high("interLambda %f", lambda);
530
531 check( start_x_next= xsh_order_list_eval( order_list, order_list->list[iorder].cenpoly,
532 ycen+1));
533 check( xsh_interpolate_linear( wavemap_data, wavemap_data, qual,
534 nx, ny, start_x_next-1, ycen, &lambda_next, &err, &qual_val, WAVEMAP_MODE));
535 }
536 y_ovs=ycen*oversample+1;
537 }
538 else{
539
540 ycen = endy+1;
541
542// xsh_msg_dbg_high("DBG ycen %d", ycen);
543
544 while( !((lambda >= lambda_min) && (lambda_next > lambda)) ){
545 double start_x_next;
546
547 ycen--;
548 check( xcen = xsh_order_list_eval( order_list, order_list->list[iorder].cenpoly,
549 ycen));
550// xsh_msg_dbg_high("DBG xcen %f ycen %d", xcen, ycen);
551 check( xsh_interpolate_linear( wavemap_data, wavemap_data, qual,
552 nx, ny, xcen-1, ycen-1, &lambda, &err, &qual_val, WAVEMAP_MODE));
553// xsh_msg_dbg_high("interLambda %f", lambda);
554
555 check( start_x_next= xsh_order_list_eval( order_list, order_list->list[iorder].cenpoly,
556 ycen-1));
557 check( xsh_interpolate_linear( wavemap_data, wavemap_data, qual,
558 nx, ny, start_x_next-1, ycen-2, &lambda_next, &err, &qual_val, WAVEMAP_MODE));
559 }
560 y_ovs=ycen*oversample-1;
561 }
562
563 /* initialize */
564 size=0;
565 x= (float) xcen;
566 y =(float) ycen;
567 lambda_old = -1;
568
569 xsh_msg_dbg_high("first X,Y %f,%f lambda %f in [%f,%f]", x,y, lambda, lambda_min, lambda_max);
570
571 if ( xsh_instrument_get_arm( instr) !=XSH_ARM_NIR){
572
573 while( (lambda >= lambda_min) && (lambda > lambda_old) &&
574 (lambda <= lambda_max) &&
575 ( y_ovs <= (endy * oversample)) ){
576
577 xsh_msg_dbg_high("size %d lambda %f ,x %f y %f", size, lambda, x,y);
578 /* take the X,Y pixel */
579 lambdastab[size] = lambda;
580 xtab[size] = x;
581 ytab[size] = y;
582 size++;
583 lambda_old = lambda;
584 y = (float) (y_ovs) / (float) oversample;
585
586 check( x = xsh_order_list_eval( order_list, order_list->list[iorder].cenpoly,
587 y));
588
589 check( xsh_interpolate_linear( wavemap_data, wavemap_data, qual,
590 nx, ny, x-1, y-1, &lambda, &err, &qual_val, WAVEMAP_MODE));
591 y_ovs++;
592 }
593 }
594 else{
595 while( (lambda >= lambda_min) &&
596 (lambda <= lambda_max) &&
597 ( y_ovs >= (starty * oversample)) ){
598
599 xsh_msg_dbg_high("size %d lambda %f ,x %f y %f", size, lambda, x,y);
600 /* take the X,Y pixel */
601 lambdastab[size] = lambda;
602 xtab[size] = x;
603 ytab[size] = y;
604 size++;
605 lambda_old = lambda;
606 y = (float) (y_ovs) / (float) oversample;
607
608 check( x = xsh_order_list_eval( order_list, order_list->list[iorder].cenpoly,
609 y));
610
611 check( xsh_interpolate_linear( wavemap_data, wavemap_data, qual,
612 nx, ny, x-1, y-1, &lambda, &err, &qual_val, WAVEMAP_MODE));
613 y_ovs--;
614// xsh_msg_dbg_high("X Y %f, %f", x, y);
615// xsh_msg_dbg_high(" END : lambda %f [%f,%f] (lambda %f > old %f) (Y %f >= %d)", lambda, lambda_min, lambda_max, lambda, lambda_old, y_ovs, starty * oversample);
616 }
617 }
618
619 *sizetab = size;
620
621 for (i=0; i< size; i++){
622 x = xtab[i];
623 y = ytab[i];
624 check( xsh_interpolate_linear( slitmap_data, slitmap_data, qual,
625 nx, ny, x-1, y-1, &slit, &err, &qual_val, WAVEMAP_MODE));
626 slitstab[i] = slit;
627 }
628
629/* REGDEBUG to remove after */
630/*
631{
632 int itest;
633
634 for(itest=0; itest < size; itest++){
635 xsh_msg("CMP order, X,Y, lambda slit : %d %f %f %f %f", order, xtab[itest], ytab[itest], lambdastab[itest],
636 slitstab[itest]);
637 }
638}
639*/
640
641 cleanup:
642 XSH_FREE( qual);
643 xsh_free_image( &wavemap);
644 return;
645}
646/*****************************************************************************/
647
648
649/*****************************************************************************/
664/*****************************************************************************/
665static void xsh_vector_divide_poly( cpl_vector *vector, double oversample,
666 cpl_polynomial *poly, int shift, xsh_instrument *instr)
667{
668 int i, size;
669
670 XSH_ASSURE_NOT_NULL( vector);
671 XSH_ASSURE_NOT_NULL( poly);
672
673 check( size = cpl_vector_get_size( vector));
674
675 for( i=0; i< size; i++){
676 double flux, val;
677
678 check( val = cpl_polynomial_eval_1d( poly,
679 (double)i/oversample+shift, NULL));
680 if ( xsh_instrument_get_arm( instr) == XSH_ARM_NIR){
681 check( val = cpl_polynomial_eval_1d( poly,
682 (double)(size-1-i)/oversample+shift, NULL));
683 }
684 else{
685 check( val = cpl_polynomial_eval_1d( poly,
686 (double)i/oversample+shift, NULL));
687 }
688 check( flux = cpl_vector_get( vector, i));
689 check( cpl_vector_set( vector, i, flux/val));
690 }
691
692 cleanup:
693 return;
694}
695
696
697/*****************************************************************************/
719/*****************************************************************************/
720static void xsh_image_gaussian_fit_y( cpl_image* img, int chunk_size,
721 int deg_poly, int oversample,
722 cpl_polynomial** center, cpl_polynomial** height,
723 cpl_polynomial** width, cpl_polynomial** offset)
724{
725
726 int nx, ny, nb_chunk, i;
727 cpl_vector *center_gauss = NULL;
728 cpl_vector *height_gauss = NULL;
729 cpl_vector *width_gauss = NULL;
730 cpl_vector *offset_gauss = NULL;
731 cpl_vector *xpos_gauss = NULL;
732 cpl_vector *chunk_vector = NULL;
733 cpl_vector *chunk_pos_vector = NULL;
734 cpl_vector *median_vect = NULL;
735 double *data = NULL;
736 double *center_gauss_data = NULL;
737 double *height_gauss_data = NULL;
738 double *width_gauss_data = NULL;
739 double *offset_gauss_data = NULL;
740 double *xpos_gauss_data = NULL;
741 double *median_data = NULL;
742
743 int nbfit = 0;
744
746 XSH_ASSURE_NOT_ILLEGAL( chunk_size >=0);
747 XSH_ASSURE_NOT_ILLEGAL( deg_poly >=0);
748 XSH_ASSURE_NOT_NULL( center);
749 XSH_ASSURE_NOT_NULL( height);
751 XSH_ASSURE_NOT_NULL( offset);
752
753 check( nx = cpl_image_get_size_x( img));
754 check( ny = cpl_image_get_size_y( img));
755 check( data = cpl_image_get_data_double( img));
756
757 if ( (chunk_size >= nx) || (chunk_size == 0)){
758 chunk_size = nx;
759 }
760
761 nb_chunk = nx / chunk_size;
762
763 XSH_MALLOC( center_gauss_data, double, nb_chunk*2);
764 XSH_MALLOC( height_gauss_data, double, nb_chunk*2);
765 XSH_MALLOC( width_gauss_data, double, nb_chunk*2);
766 XSH_MALLOC( offset_gauss_data, double, nb_chunk*2);
767 XSH_MALLOC( xpos_gauss_data, double, nb_chunk*2);
768 XSH_CALLOC( median_data, double, chunk_size);
769
770 check( chunk_vector = cpl_vector_new( ny));
771 check( chunk_pos_vector = cpl_vector_new( ny));
772
773 xsh_msg_dbg_medium( "nx %d nb_chunks %d of size %d", nx, nb_chunk,
774 chunk_size);
775
776 for( i=0; i< nb_chunk; i++){
777 double x0, sigma, area, offset_val, height_val, fwhm;
778 double med;
779 int j;
780 int ideb, ilast;
781
782 ideb = i*chunk_size;
783 ilast = ideb+chunk_size;
784
785 for( j=0; j< ny; j++){
786 int l, val_nb=0;
787
788 for( l=ideb; l< ilast; l++){
789 val_nb++;
790 median_data[l-ideb] = data[j*nx+l];
791 }
792 check( median_vect = cpl_vector_wrap( val_nb, median_data));
793 check( med = cpl_vector_get_median( median_vect));
794 check( cpl_vector_set( chunk_vector, j, med));
795 check( cpl_vector_set( chunk_pos_vector, j, j));
796 }
797#if REGDEBUG_OPTEXTRACT
798#if REGDEBUG_GAUSSIAN
799 /* REGDEBUG */
800 {
801 FILE* regdebug = NULL;
802 char filename[256];
803
804 sprintf(filename, "gaussian_points_%d.dat",i);
805 regdebug = fopen(filename,"w");
806 for( j=0; j< ny; j++){
807 double pos, val;
808
809 pos = cpl_vector_get( chunk_pos_vector,j);
810 val = cpl_vector_get( chunk_vector,j);
811
812 fprintf(regdebug, "%f %f\n", pos, val);
813 }
814 fclose(regdebug);
815 }
816#endif
817#endif
818 cpl_vector_fit_gaussian( chunk_pos_vector, NULL, chunk_vector, NULL,
819 CPL_FIT_ALL, &x0, &sigma, &area, &offset_val, NULL, NULL, NULL);
820
821 if ( cpl_error_get_code() == CPL_ERROR_NONE){
822 height_val = area / sqrt(2*CPL_MATH_PI*sigma*sigma);
823 fwhm = CPL_MATH_FWHM_SIG*sigma;
824 if ( fwhm < (oversample*FIT_FWHM_LIMIT)){
825 center_gauss_data[nbfit] = x0;
826 height_gauss_data[nbfit] = height_val;
827 width_gauss_data[nbfit] = sigma;
828 xpos_gauss_data[nbfit] = ideb;
829 offset_gauss_data[nbfit] = offset_val;
830 nbfit++;
831 center_gauss_data[nbfit] = x0;
832 height_gauss_data[nbfit] = height_val;
833 width_gauss_data[nbfit] = sigma;
834 xpos_gauss_data[nbfit] = ilast;
835 offset_gauss_data[nbfit] = offset_val;
836 nbfit++;
837 }
838 else{
839 xsh_msg("WARNING chunk %d-%d so big FWHM %f expected < %d",
840 ideb, ilast, fwhm, oversample*FIT_FWHM_LIMIT);
841 }
842 }
843 else {
844 xsh_msg("WARNING chunk %d-%d don't converge with gaussian fit",
845 ideb, ilast);
847 }
848 }
849 /* fit with order 3 polynomials */
850
851 if ( nbfit == 0){
852 xsh_msg("WARNING nbfit == 0 force poly degree to 0");
853 nbfit = 1;
854 deg_poly = 0;
855 xpos_gauss_data[0] = 0;
856 center_gauss_data[0] = ny/2.0;
857 height_gauss_data[0] = 1;
858 width_gauss_data[0] = ny;
859 offset_gauss_data[0] = 0.0;
860 }
861 else if ( nbfit <= deg_poly){
862 xsh_msg("WARNING %d (nbfit) < %d (poly degree) force poly degree to %d",
863 nbfit, deg_poly, nbfit-1);
864 deg_poly = nbfit-1;
865 }
866
867 check(xpos_gauss = cpl_vector_wrap( nbfit, xpos_gauss_data));
868 check( center_gauss = cpl_vector_wrap( nbfit, center_gauss_data));
869 check( height_gauss = cpl_vector_wrap( nbfit, height_gauss_data));
870 check( width_gauss = cpl_vector_wrap( nbfit, width_gauss_data));
871 check( offset_gauss = cpl_vector_wrap( nbfit, offset_gauss_data));
872
873 check( *center = xsh_polynomial_fit_1d_create( xpos_gauss,
874 center_gauss, deg_poly, NULL));
875 check( *height = xsh_polynomial_fit_1d_create( xpos_gauss,
876 height_gauss, deg_poly, NULL));
878 width_gauss, deg_poly, NULL));
879 check( *offset = xsh_polynomial_fit_1d_create( xpos_gauss,
880 offset_gauss, deg_poly, NULL));
881
882 cleanup:
883 if ( cpl_error_get_code() != CPL_ERROR_NONE){
884 xsh_free_polynomial( center);
885 xsh_free_polynomial( height);
887 xsh_free_polynomial( offset);
888 }
889 xsh_unwrap_vector( &center_gauss);
890 xsh_unwrap_vector( &width_gauss);
891 xsh_unwrap_vector( &height_gauss);
892 xsh_unwrap_vector( &offset_gauss);
893 xsh_unwrap_vector( &xpos_gauss);
894 xsh_unwrap_vector( &median_vect);
895 XSH_FREE( center_gauss_data);
896 XSH_FREE( width_gauss_data);
897 XSH_FREE( height_gauss_data);
898 XSH_FREE( offset_gauss_data);
899 XSH_FREE( xpos_gauss_data);
900 XSH_FREE( median_data);
901 xsh_free_vector( &chunk_vector);
902 xsh_free_vector( &chunk_pos_vector);
903 return;
904}
905
906
907static cpl_image* xsh_image_create_model_image( cpl_image *src_img,
908 double *x_data, double *y_data, double kappa, int niter, double frac_min)
909{
910 cpl_image *result = NULL;
911 int nx,ny;
912 int i,j;
913 double *data = NULL;
914 double *src_data = NULL;
915 double *pos = NULL;
916 double *line = NULL;
917 cpl_vector *pos_vect = NULL;
918 cpl_vector *line_vect = NULL;
919 int size;
920 cpl_polynomial *fit = NULL;
921
922 /* Check input parameters */
923 XSH_ASSURE_NOT_NULL( src_img);
924 check( nx = cpl_image_get_size_x( src_img));
925 check( ny = cpl_image_get_size_y( src_img));
926 check( result = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE));
927 check( data = cpl_image_get_data_double( result));
928 check( src_data = cpl_image_get_data_double( src_img));
929 XSH_CALLOC( pos, double, nx);
930 XSH_CALLOC( line, double, nx);
931 double before=cpl_test_get_walltime();
932 xsh_msg_debug("setting timing");
933
934 for(j=0; j<ny; j++){
935 double val;
936 double chisq;
937 int *flags;
938 int j_nx=j*nx;
939 size = 0;
940 for(i=0; i<nx; i++){
941 double diffx, diffy;
942 int pix=i+j_nx;
943 diffx = x_data[pix]-(int)x_data[pix];
944 if(diffx<0.2) {
945 diffy = y_data[pix]-(int)y_data[pix];
946
947 if ( diffy < 0.2){
948 pos[size] = i;
949 line[size] = src_data[pix];
950 size++;
951 }
952 }
953 }
954 check( pos_vect = cpl_vector_wrap( size, pos));
955 check( line_vect = cpl_vector_wrap( size, line));
956
957 check( xsh_array_clip_poly1d( pos_vect, line_vect, kappa, niter,
958 frac_min, 3, &fit, &chisq, &flags));
959
960 for(i=0; i<nx; i++){
961 check( val = cpl_polynomial_eval_1d( fit, i, NULL));
962 data[i+j_nx] = fabs(val);
963 }
964 xsh_unwrap_vector( &pos_vect);
965 xsh_unwrap_vector( &line_vect);
966 xsh_free_polynomial( &fit);
967 XSH_FREE( flags);
968 }
969 double after=cpl_test_get_walltime();
970 xsh_msg_debug("setting timing");
971 xsh_msg("Time spent %g %g %g",before,after,after-before);
972
973 cleanup:
974 if ( cpl_error_get_code() != CPL_ERROR_NONE){
975 xsh_free_image( &result);
976 xsh_unwrap_vector( &pos_vect);
977 xsh_unwrap_vector( &line_vect);
978 xsh_free_polynomial( &fit);
979 }
980 XSH_FREE( pos);
981 XSH_FREE( line);
982 return result;
983}
984
985static cpl_image* xsh_image_create_gaussian_image( cpl_image *src,
986 cpl_polynomial* centerp, cpl_polynomial* heightp, cpl_polynomial* widthp,
987 cpl_polynomial* offsetp)
988{
989 cpl_image *result = NULL;
990 int nx,ny;
991 int i,j;
992 double *data = NULL;
993
994 /* Check input parameters */
996 XSH_ASSURE_NOT_NULL( centerp);
997 XSH_ASSURE_NOT_NULL( heightp);
998 XSH_ASSURE_NOT_NULL( widthp);
999 XSH_ASSURE_NOT_NULL( offsetp);
1000
1001 check( nx = cpl_image_get_size_x( src));
1002 check( ny = cpl_image_get_size_y( src));
1003 check( result = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE));
1004 check( data = cpl_image_get_data_double( result));
1005
1006 for(i=0; i<nx; i++){
1007 double x0, height, width, offset, z;
1008 //double zlow, zup;
1009 /* not needed: weight; */
1010
1011 check( x0 = cpl_polynomial_eval_1d( centerp, i, NULL));
1012 check( height = cpl_polynomial_eval_1d( heightp, i, NULL));
1013 check( width = cpl_polynomial_eval_1d( widthp, i, NULL));
1014 check( offset = cpl_polynomial_eval_1d( offsetp, i, NULL));
1015
1016 double fct=width*XSH_MATH_SQRT_2;
1017 double inv_fct=1.0/fct;
1018 //double x0_plus=x0+0.5;
1019 //double x0_minus=x0-0.5;
1020
1021 for(j=0; j<ny; j++){
1022 z = (j-x0)*inv_fct;
1023 //zlow = (j-x0_plus)*inv_fct;
1024 //zup = (j-x0_minus)*inv_fct;
1025 data[i+j*nx] = height*exp(-(z*z))+offset;
1026/* not needed
1027
1028 weight = fct*(gsl_sf_erf(zup)-gsl_sf_erf(zlow));
1029 */
1030 }
1031
1032 }
1033
1034
1035 cleanup:
1036 if ( cpl_error_get_code() != CPL_ERROR_NONE){
1037 xsh_free_image( &result);
1038 }
1039 return result;
1040}
1041/*****************************************************************************/
1058/*****************************************************************************/
1059
1060static cpl_vector* xsh_image_extract_standard( cpl_image* img,
1061 cpl_image *err_img, cpl_image *qual_img,
1062 cpl_vector **err_v, cpl_vector **qual_v)
1063{
1064 cpl_vector *result = NULL;
1065 int nx, ny, i;
1066 double *data = NULL;
1067 double *err = NULL;
1068 int *qual = NULL;
1069 double* pres=NULL;
1070 double* perr=NULL;
1071 double* pqua=NULL;
1072
1073
1074 XSH_ASSURE_NOT_NULL( img);
1075 XSH_ASSURE_NOT_NULL( err_img);
1076 XSH_ASSURE_NOT_NULL( qual_img);
1077 XSH_ASSURE_NOT_NULL( err_v);
1078 XSH_ASSURE_NOT_NULL( qual_v);
1079
1080 /* get data */
1081 check (nx = cpl_image_get_size_x( img));
1082 check (ny = cpl_image_get_size_y( img));
1083 check( data = cpl_image_get_data_double( img));
1084 check( err = cpl_image_get_data_double( err_img));
1085 check( qual = cpl_image_get_data_int( qual_img));
1086
1087 /* allocate memory */
1088 check( result = cpl_vector_new( nx));
1089 check( *err_v = cpl_vector_new( nx));
1090 check( *qual_v = cpl_vector_new( nx));
1091 pres=cpl_vector_get_data(result);
1092 perr=cpl_vector_get_data(*err_v);
1093 pqua=cpl_vector_get_data(*qual_v);
1094 /* standard extraction */
1095 for( i=0; i< nx; i++){
1096 int j;
1097 double pix_val = 0.0, err_val = 0.0;
1098 int qual_val = 0;
1099
1100 for( j=0; j< ny; j++){
1101 int pix=j*nx+i;
1102 pix_val += data[pix];
1103 err_val += err[pix]*err[pix];
1104 qual_val |= qual[pix];
1105 }
1106
1107 pres[i]=pix_val;
1108 perr[i]=sqrt(err_val);
1109 pqua[i]=(double)qual_val;
1110
1111 }
1112
1113 cleanup:
1114 if ( cpl_error_get_code() != CPL_ERROR_NONE){
1115 xsh_free_vector( &result);
1116 xsh_free_vector( err_v);
1117 xsh_free_vector( qual_v);
1118 }
1119 return result;
1120}
1121/*****************************************************************************/
1122
1123/*****************************************************************************/
1140/*****************************************************************************/
1141static cpl_image* xsh_image_divide_1D( cpl_image *image, cpl_image *err_img,
1142 cpl_vector *s1D, cpl_vector *err_s1D, cpl_image **err_2dby1D_img)
1143{
1144 cpl_image *result = NULL;
1145 int i,j, nx, ny, size;
1146 double *data = NULL;
1147 double *err = NULL;
1148 double *data_res= NULL;
1149 double *err_res = NULL;
1150
1151 /* Check parameters */
1152 XSH_ASSURE_NOT_NULL( image);
1153 XSH_ASSURE_NOT_NULL( err_img);
1154 XSH_ASSURE_NOT_NULL( s1D);
1155 XSH_ASSURE_NOT_NULL( err_s1D);
1156 XSH_ASSURE_NOT_NULL( err_2dby1D_img);
1157
1158 check( size = cpl_vector_get_size( s1D));
1159 check( nx = cpl_image_get_size_x( image));
1160 check( ny = cpl_image_get_size_y( image));
1162
1163 check( data = cpl_image_get_data_double( image));
1164 check( err = cpl_image_get_data_double( err_img));
1165 check( result = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE));
1166 check( *err_2dby1D_img = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE));
1167 check( data_res = cpl_image_get_data_double( result));
1168 check( err_res = cpl_image_get_data_double( *err_2dby1D_img));
1169
1170 for(i=0; i< nx; i++){
1171 double d, d2, e2, d1, e1;
1172
1173 check( d2 = cpl_vector_get( s1D, i));
1174 check( e2 = cpl_vector_get( err_s1D, i));
1175
1176 if ( d2 == 0.0){
1177 if (i > 0 && i < (nx-1)){
1178 d2 = 0.5*(cpl_vector_get( s1D, i-1)+cpl_vector_get( s1D, i+1));
1179 e2 = 0.5*sqrt(cpl_vector_get( err_s1D, i-1)+cpl_vector_get( err_s1D, i+1));
1180 }
1181 }
1182 for( j=0; j <ny; j++){
1183 int pix=i+j*nx;
1184 d1 = data[pix];
1185 e1 = err[pix];
1186 d = d1/d2;
1187
1188 data_res[pix] = d;
1189 err_res[pix] = d * sqrt( (e1*e1)/(d1*d1) + (e2*e2)/(d2*d2) );
1190 }
1191 }
1192
1193 cleanup:
1194 if ( cpl_error_get_code() != CPL_ERROR_NONE){
1195 xsh_free_image( &result);
1196 xsh_free_image( err_2dby1D_img);
1197 }
1198 return result;
1199}
1200/*****************************************************************************/
1201
1202
1203/*****************************************************************************/
1228/*****************************************************************************/
1229
1230static cpl_vector* xsh_image_extract_optimal( cpl_image* img,
1231 cpl_image *errs_img, cpl_image *qual_img, xsh_opt_extract_param* opt_par,
1232 cpl_image *model_img,const int decode_bp,
1233 cpl_image** corr_img, cpl_vector **err_v, cpl_vector **qual_v)
1234{
1235 int nb_rejected =-1, iter=1, clip_niter = 1, nbtotal_rejected=0;
1236 double clip_kappa = 3, clip_frac=0, clip_frac_rej =0.4;
1237 int *rejected = NULL;
1238 double *ab = NULL;
1239 cpl_vector *median_vect = NULL;
1240
1241 cpl_vector *result = NULL;
1242 int nx, ny, i;
1243 double *data = NULL;
1244 double *errs = NULL;
1245 int *qual_data = NULL;
1246 cpl_image *corr_image = NULL;
1247 double *corr_data = NULL;
1248 double *model_data = NULL;
1249
1250
1251 XSH_ASSURE_NOT_NULL( img);
1252 XSH_ASSURE_NOT_NULL( errs_img);
1253 XSH_ASSURE_NOT_NULL( qual_img);
1254 XSH_ASSURE_NOT_NULL( model_img);
1255 XSH_ASSURE_NOT_NULL( corr_img);
1256 XSH_ASSURE_NOT_NULL( err_v);
1257 XSH_ASSURE_NOT_NULL( qual_v);
1258 XSH_ASSURE_NOT_NULL( opt_par);
1259
1260 /* get parameters */
1261 clip_kappa = opt_par->clip_kappa;
1262 clip_niter = opt_par->clip_niter;
1263 clip_frac = opt_par->clip_frac;
1264
1265 /* get data */
1266 check (nx = cpl_image_get_size_x( img));
1267 check (ny = cpl_image_get_size_y( img));
1268 check( data = cpl_image_get_data_double( img));
1269 check( errs = cpl_image_get_data_double( errs_img));
1270 check( qual_data = cpl_image_get_data_int( qual_img));
1271
1272 check( corr_image = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE));
1273 check( corr_data = cpl_image_get_data_double( corr_image));
1274 check( model_data = cpl_image_get_data_double( model_img));
1275
1276 /* allocate memory */
1277 check( result = cpl_vector_new( nx));
1278 check( *err_v = cpl_vector_new( nx));
1279 check( *qual_v = cpl_vector_new( nx));
1280
1281 XSH_CALLOC( rejected, int, ny);
1282 XSH_CALLOC( ab, double, nx);
1283
1284 /* optimal extraction */
1285 for( i=0; i< nx; i++){
1286 int j;
1287 double pix_val = 0.0;
1288 double f_err = 0.0;
1289 double f_num = 0.0;
1290 double f_denom = 0.0;
1291 double sum_MP = 0.0;
1292 //double f_denom_mod = 0.0;
1293 double err_val=0.0;
1294 int qual_val = 0;
1295 //double s2dby1d_variance, crh_threshold;
1296 int count = 0;
1297 double med_flux;
1298 double sum, sum_err;
1299 int sum_qual = 0;
1300
1301 for( j=0; j< ny; j++){
1302 rejected[j] = 1;
1303 }
1304
1305 /* First detect bad points */
1306
1307 count = 0;
1308 nbtotal_rejected=0;
1309 sum = 0;
1310 sum_err = 0;
1311
1312 for( j=0; j< ny; j++){
1313 int qual;
1314 double flux = 0.0, err_val = 0.0;
1315 double model_flux = 0.0;
1316
1317 flux = data[i+j*nx];
1318 err_val = errs[i+j*nx];
1319 qual = qual_data[i+j*nx];
1320 model_flux = model_data[i+j*nx];
1321
1322 if ( ((qual & decode_bp ) == 0) && model_flux != 0){
1323 /* good pix */
1324 ab[count] = flux/model_flux;
1325 count++;
1326 qual_val |= qual;
1327 }
1328 else{
1329 rejected[j] =0;
1330 nbtotal_rejected++;
1331 sum += flux;
1332 sum_err += err_val*err_val;
1333 sum_qual |= qual;
1334 }
1335 }
1336 clip_frac = nbtotal_rejected/(double)ny;
1337 xsh_msg_dbg_medium("Bad pixel fraction %f for %d", clip_frac, i);
1338
1339 if (clip_frac < 1){
1340 check( median_vect = cpl_vector_wrap( count, ab));
1341 check( med_flux = cpl_vector_get_median( median_vect));
1342 xsh_unwrap_vector( &median_vect);
1343 }
1344 else{
1345 med_flux = sum;
1346 sum_err = sqrt( sum_err);
1347 }
1348
1349 iter=1;
1350 nb_rejected =-1;
1351
1352 /* sigma clipping */
1353 while( iter <= clip_niter && nb_rejected != 0 && count != 0
1354 && clip_frac <= clip_frac_rej){
1355 nb_rejected = 0;
1356 int nbdiff =0;
1357 double sigma;
1358
1359 xsh_msg_dbg_medium("Sigma clipping : Iteration %d/%d", iter, clip_niter);
1360
1361 /* compute sigma */
1362 for( j=0; j< ny; j++){
1363 double flux = 0.0;
1364 double model_flux = 0.0;
1365 double diff = 0.0;
1366
1367 model_flux = model_data[i+j*nx];
1368 flux = data[i+j*nx];
1369
1370 /* Good pixel only */
1371 if ( rejected[j] == 1){
1372 diff = flux-med_flux*model_flux;
1373 ab[nbdiff] = diff;
1374 nbdiff++;
1375 }
1376 }
1377 check( median_vect = cpl_vector_wrap( nbdiff, ab));
1378 check( sigma = cpl_vector_get_stdev( median_vect));
1379 xsh_unwrap_vector( &median_vect);
1380 xsh_msg_dbg_medium( "Sigma %f", sigma);
1381
1382 /* rejected point */
1383 for( j=0; j< ny; j++){
1384 double flux = 0.0;
1385 double model_flux = 0.0;
1386 double diff = 0.0;
1387
1388 model_flux = model_data[i+j*nx];
1389 flux = data[i+j*nx];
1390 diff = flux-med_flux*model_flux;
1391
1392 if ( rejected[j] == 1 && fabs( diff) > clip_kappa*sigma){
1393 nb_rejected++;
1394 rejected[j] = 0;
1395 qual_val = QFLAG_COSMIC_RAY_REMOVED;
1396 }
1397 }
1398 nbtotal_rejected += nb_rejected;
1399 /* recompute median */
1400 count = 0;
1401
1402 for( j=0; j< ny; j++){
1403 double flux = 0.0;
1404 double model_flux = 0.0;
1405
1406 flux = data[i+j*nx];
1407 model_flux = model_data[i+j*nx];
1408
1409 if ( rejected[j] == 1 && model_flux != 0){
1410 ab[count] = flux/model_flux;
1411 count++;
1412 }
1413 }
1414 if ( count != 0){
1415 check( median_vect = cpl_vector_wrap( count, ab));
1416 check( med_flux = cpl_vector_get_median( median_vect));
1417 xsh_unwrap_vector( &median_vect);
1418 }
1419 clip_frac = nbtotal_rejected/(double)ny;
1420 iter++;
1421 }
1422 xsh_msg_dbg_medium("Fraction rejected %f", clip_frac);
1423
1424 /* correct bad points */
1425 for( j=0; j< ny; j++){
1426 int qual;
1427 //double flux = 0.0;
1428 qual = qual_data[i+j*nx];
1429 if ( ( (qual & decode_bp) == 0 ) && rejected[j] == 1 ){
1430 corr_data[i+j*nx] = data[i+j*nx];
1431 }
1432 else{
1433 corr_data[i+j*nx] = 0;
1434 }
1435 }
1436
1437 /* optimal extraction */
1438 if ( clip_frac < 1){
1439 for( j=0; j< ny; j++){
1440 double p, p2;
1441 double variance;
1442
1443 p = model_data[i+j*nx];
1444 p2 = p*p;
1445 variance = errs[i+j*nx]*errs[i+j*nx];
1446 // xsh_msg(" Rejected %d %f %f %f\n", rejected[j], f_num, f_denom,corr_data[i+j*nx]);
1447 sum_MP += rejected[j]*p;
1448 //f_err += rejected[j]*p2;
1449 f_err += rejected[j]*p/variance;
1450 //xsh_msg("f_err %f", f_err);
1451 f_num += (double)rejected[j]*p*corr_data[i+j*nx]/variance;
1452 f_denom += (double)rejected[j]*p2/variance;
1453 }
1454 pix_val = f_num / f_denom;
1455 err_val = sqrt( sum_MP / f_denom);
1456 }
1457 else{
1458 pix_val = sum;
1459 err_val = sum_err;
1460 }
1461 check( cpl_vector_set( result, i, pix_val));
1462 check( cpl_vector_set( *err_v, i, err_val));
1463 check( cpl_vector_set( *qual_v, i, (double)qual_val));
1464 }
1465
1466
1467#if REGDEBUG_N
1469 char std_spectrum_name[256];
1470 FILE* std_spectrum_file = NULL;
1471
1472 sprintf( std_spectrum_name, "n.dat");
1473 std_spectrum_file = fopen( std_spectrum_name, "w");
1474
1475 fprintf( std_spectrum_file, "#index n\n");
1476
1477 for(i=0; i< nx; i++){
1478 fprintf(std_spectrum_file, "%d %f\n", i, n[i]);
1479 }
1480 fclose( std_spectrum_file);
1481
1482 sprintf( std_spectrum_name, "rejfrac.dat");
1483 std_spectrum_file = fopen( std_spectrum_name, "w");
1484
1485 fprintf( std_spectrum_file, "#index rejfrac\n");
1486
1487 for(i=0; i< nx; i++){
1488 fprintf(std_spectrum_file, "%d %f\n", i, rej[i]);
1489 }
1490 fclose( std_spectrum_file);
1491 }
1492#endif
1493 *corr_img = corr_image;
1494
1495 cleanup:
1496 if ( cpl_error_get_code() != CPL_ERROR_NONE){
1497 xsh_free_vector( &result);
1498 xsh_free_vector( err_v);
1499 xsh_free_vector( qual_v);
1500 }
1501 XSH_FREE( rejected);
1502 XSH_FREE( ab);
1503 return result;
1504}
1505
1506/*****************************************************************************/
1533/*****************************************************************************/
1534
1535static cpl_vector* xsh_vector_integrate( int biny,
1536 int oversample, int absorder,
1537 cpl_vector *ref_pos, double step,
1538 cpl_vector *ref_values, cpl_vector *err_values, cpl_vector *qual_values,
1539 cpl_vector *new_pos,
1540 cpl_vector **spectrum_err, cpl_vector **spectrum_qual)
1541{
1542
1543 int new_pos_size, ref_pos_size, pixel_pos_size;
1544 int i, j;
1545 cpl_vector *result=NULL;
1546 cpl_vector *pixel_pos_vect = NULL, *pixel_size_vect = NULL;
1547 cpl_polynomial *pixel_size_poly = NULL;
1548 double hstep;
1549 double lambda_bin_factor =0.0;
1550
1551 XSH_ASSURE_NOT_NULL( ref_pos);
1552 XSH_ASSURE_NOT_NULL( ref_values);
1553 XSH_ASSURE_NOT_NULL( err_values);
1554 XSH_ASSURE_NOT_NULL( qual_values);
1555 XSH_ASSURE_NOT_NULL( new_pos);
1556 XSH_ASSURE_NOT_NULL( spectrum_err);
1557 XSH_ASSURE_NOT_NULL( spectrum_qual);
1558
1559 check( new_pos_size = cpl_vector_get_size( new_pos));
1560 check( ref_pos_size = cpl_vector_get_size( ref_pos));
1561 check( result = cpl_vector_new( new_pos_size));
1562
1563 check( *spectrum_err = cpl_vector_new( new_pos_size));
1564 check( *spectrum_qual = cpl_vector_new( new_pos_size));
1565
1566 /* normalize by pixels size */
1567 pixel_pos_size = ref_pos_size-2*oversample;
1568 check( pixel_pos_vect = cpl_vector_new( pixel_pos_size));
1569 check( pixel_size_vect = cpl_vector_new( pixel_pos_size));
1570
1571 for( j=oversample; j< ref_pos_size-oversample; j++){
1572 double xA, xB, xC;
1573 double sizeAB, sizeBC, size;
1574
1575 check( xA = cpl_vector_get( ref_pos, j-oversample));
1576 check( xB = cpl_vector_get( ref_pos, j));
1577 check( xC = cpl_vector_get( ref_pos, j+oversample));
1578 sizeAB = xB-xA;
1579 sizeBC = xC-xB;
1580 size = (sizeAB+sizeBC)/2.0;
1581 check( cpl_vector_set( pixel_pos_vect, j-oversample, xB));
1582 check( cpl_vector_set( pixel_size_vect, j-oversample, size));
1583 }
1584
1585 check( pixel_size_poly = xsh_polynomial_fit_1d_create( pixel_pos_vect,
1586 pixel_size_vect, 1, NULL));
1587
1588 j =1;
1589 hstep = step/2.0;
1590
1591 for( i=0; i< new_pos_size; i++){
1592 double xA=0, xB=0, xC=0;
1593 double deltaA1=0, deltaA2=0;
1594
1595 //double y_min, y_max;
1596 double size_min=0, pos_min=0, frac_min=0, sigy_min=0, flux_min=0;
1597 double size_max=0, pos_max=0, frac_max=0, sigy_max=0, flux_max=0;
1598
1599 double sig_min=0, sig_max=0;
1600 int qual_min=0, qual_max=0, qual=0;
1601 //double yA, yB, sig_yA, sig_yB;
1602
1603 double y=0, ymin=0, ymax=0;
1604 //double border_xmin, border_xmax, delta_xmin, delta_xmax;
1605 double xmin=0, xmax=0;
1606 /* TODO: why sig_y has been set to 26? */
1607 double x=0, sig_y = 26;
1608 int j_start=0, j_end=0, k=0;
1609 //double A, B, min_cont, max_cont;
1610
1611 double flux=0, err=0, pixel_size_comp=0;
1612 double nb_x=0;
1613
1614 check( x= cpl_vector_get( new_pos, i));
1615 xmin = x-hstep;
1616 xmax = x+hstep;
1617 check( xA = cpl_vector_get( ref_pos, j-1));
1618 check( xB = cpl_vector_get( ref_pos, j));
1619 deltaA2 = (xB-xA)/2.0;
1620
1621 if ( j>1){
1622 check( xC = cpl_vector_get( ref_pos, j-2));
1623 deltaA1 = (xA-xC)/2.0;
1624 }
1625 else{
1626 deltaA1 = deltaA2;
1627 }
1628
1629 while( !(xmin >= (xA-deltaA1) && xmin <= (xA+deltaA2)) ){
1630 j++;
1631 check( xA = cpl_vector_get( ref_pos, j-1));
1632 check( xB = cpl_vector_get( ref_pos, j));
1633 check( xC = cpl_vector_get( ref_pos, j-2));
1634 deltaA1 = (xA-xC)/2.0;
1635 deltaA2 = (xB-xA)/2.0;
1636 }
1637
1638 j_start = j-1;
1639 size_min = deltaA1+deltaA2;
1640 pos_min = xA+deltaA2;
1641 frac_min = (pos_min-xmin)/size_min;
1642 check( ymin = cpl_vector_get( ref_values, j_start));
1643 check( sigy_min = cpl_vector_get( err_values, j_start));
1644 check( qual_min = cpl_vector_get( qual_values, j_start));
1645
1646 flux_min = frac_min*ymin;
1647 sig_min = frac_min*frac_min*sigy_min*sigy_min;
1648
1649 while( !(xmax >= (xA-deltaA1) && xmax <= (xA+deltaA2)) ){
1650 j++;
1651 check( xA = cpl_vector_get( ref_pos, j-1));
1652 check( xC = cpl_vector_get( ref_pos, j-2));
1653
1654 deltaA1 = (xA-xC)/2.0;
1655
1656 if ( j < ref_pos_size){
1657 check( xB = cpl_vector_get( ref_pos, j));
1658 deltaA2 = (xB-xA)/2.0;
1659 }
1660 else{
1661 deltaA2 = deltaA1;
1662 }
1663 }
1664
1665 j_end = j-1;
1666 j = j_end;
1667
1668 if ( j_start == j_end){
1669 flux_min =0;
1670 frac_min = 0;
1671 frac_max = (xmax-xmin)/size_min;
1672 }
1673 else{
1674 size_max = deltaA1+deltaA2;
1675 pos_max = xA-deltaA1;
1676 frac_max = (xmax-pos_max)/size_max;
1677 }
1678 check( ymax = cpl_vector_get( ref_values, j_end));
1679 check( sigy_max = cpl_vector_get( err_values, j_end));
1680 check( qual_max = cpl_vector_get( qual_values, j_end));
1681
1682 flux_max = frac_max*ymax;
1683 sig_max = frac_max*frac_max*sigy_max*sigy_max;
1684
1685 xsh_msg_dbg_high( "xmin %f pos_min %f diff %f size %f frac*size %f", xmin , pos_min, pos_min-xmin, size_min, frac_min*size_min);
1686 xsh_msg_dbg_high( "xmax %f pos_max %f diff %f size %f frac*size %f", xmax , pos_max, xmax-pos_max, size_max, frac_max*size_max);
1687
1688 nb_x = j_end-j_start-1;
1689 if ( nb_x > 0){
1690 nb_x += frac_min+frac_max;
1691 }
1692 else{
1693 nb_x = frac_min+frac_max;;
1694 }
1695
1696 qual = qual_min;
1697 qual |= qual_max;
1698 flux = flux_min+flux_max;
1699 err = sig_min+sig_max;
1700
1701 xsh_msg_dbg_high("nb_x %f size %f %f cont_min %f max_cont %f", nb_x, size_min, size_max, frac_min, frac_max);
1702
1703 for( k=j_start+1; k <= j_end-1; k++){
1704 check( y = cpl_vector_get( ref_values, k));
1705 check( sig_y = cpl_vector_get( err_values, k));
1706 flux +=y;
1707 err += sig_y*sig_y;
1708 qual |= (int)xsh_round_double(cpl_vector_get( qual_values, k));
1709 }
1710 err = sqrt( err);
1711
1712 /* correct pixel size conversion */
1713 check( pixel_size_comp = cpl_polynomial_eval_1d( pixel_size_poly, x, NULL));
1714
1715 lambda_bin_factor = 1.0/(step*biny);
1716
1717 flux *= pixel_size_comp*lambda_bin_factor;
1718 err *= pixel_size_comp*lambda_bin_factor;
1719
1720 check( cpl_vector_set( result, i, flux));
1721 check( cpl_vector_set( *spectrum_err, i, err));
1722 check( cpl_vector_set( *spectrum_qual, i, qual));
1723 }
1724
1725#if REGDEBUG_PIXELSIZE
1727 FILE* debug_file = NULL;
1728 char debug_name[256];
1729
1730 sprintf( debug_name, "pixel_size_%d.dat", absorder);
1731 xsh_msg("Create %s", debug_name);
1732
1733 debug_file = fopen( debug_name, "w");
1734
1735 fprintf( debug_file,"# wavelength size fit\n");
1736
1737 for( j=0; j< pixel_pos_size; j++){
1738 double pos, size, pol;
1739
1740 check( pos = cpl_vector_get( pixel_pos_vect, j));
1741 check( size = cpl_vector_get( pixel_size_vect, j));
1742 check( pol = cpl_polynomial_eval_1d( pixel_size_poly, pos, NULL));
1743 fprintf( debug_file ,"%f %f %f\n", pos, size, pol);
1744 }
1745 fclose( debug_file);
1746 }
1747#endif
1748
1749#if REGDEBUG_INTEGRATE
1751 FILE* debug_file = NULL;
1752 char debug_name[256];
1753
1754 sprintf( debug_name, "integrate_flux_%d.dat", absorder);
1755 xsh_msg("Create %s", debug_name);
1756
1757 debug_file = fopen( debug_name, "w");
1758
1759 fprintf( debug_file,"# wavelength flux err qual\n");
1760
1761 for( j=0; j< new_pos_size; j++){
1762 double xA, yA, sig_yA, qual_yA;
1763
1764 check( xA = cpl_vector_get( new_pos, j));
1765 check( yA = cpl_vector_get( result, j));
1766 check( sig_yA = cpl_vector_get( *spectrum_err, j));
1767 check( qual_yA = cpl_vector_get( *spectrum_qual, j));
1768 fprintf( debug_file ,"%f %f %f %f\n", xA, yA, sig_yA, qual_yA);
1769 }
1770
1771 fclose( debug_file);
1772 }
1773#endif
1774/* ENDREGDEBUG */
1775
1776 cleanup:
1777 xsh_free_vector( &pixel_pos_vect);
1778 xsh_free_vector( &pixel_size_vect);
1779 xsh_free_polynomial( &pixel_size_poly);
1780 if( cpl_error_get_code() != CPL_ERROR_NONE){
1781 xsh_free_vector( &result);
1782 xsh_free_vector( spectrum_err);
1783 xsh_free_vector( spectrum_qual);
1784 }
1785 return result;
1786}
1787/*****************************************************************************/
1788
1789
1790/*****************************************************************************/
1831/*****************************************************************************/
1832void xsh_opt_extract(cpl_frame *sci_frame, cpl_frame *orderlist_frame,
1833 cpl_frame *wavesol_frame, cpl_frame *model_frame, cpl_frame *wavemap_frame,
1834 cpl_frame *slitmap_frame, cpl_frame *loc_frame,
1835 cpl_frame *spectralformat_frame, cpl_frame *masterflat_frame,
1837 const char* rec_prefix, cpl_frame** orderext1d_frame,
1838 cpl_frame** orderoxt1d_frame, cpl_frame** orderoxt1d_eso_frame,
1839 cpl_frame** qc_subextract_frame, cpl_frame** qc_s2ddiv1d_frame,
1840 cpl_frame** qc_model_frame, cpl_frame** qc_weight_frame)
1841
1842{
1843
1844 check(
1845 xsh_opt_extract_orders( sci_frame, orderlist_frame, wavesol_frame,
1846 model_frame, wavemap_frame, slitmap_frame, loc_frame,spectralformat_frame,
1847 masterflat_frame, instrument, opt_extract_par, 0, 100, rec_prefix,
1848 orderext1d_frame, orderoxt1d_frame, orderoxt1d_eso_frame, qc_subextract_frame,
1849 qc_s2ddiv1d_frame, qc_model_frame, qc_weight_frame));
1850
1851 cleanup: return;
1852}
1853
1854/*****************************************************************************/
1897/*****************************************************************************/
1898void
1899xsh_opt_extract_orders( cpl_frame *sci_frame,
1900 cpl_frame *orderlist_frame,
1901 cpl_frame *wavesol_frame,
1902 cpl_frame *model_frame,
1903 cpl_frame *wavemap_frame,
1904 cpl_frame *slitmap_frame,
1905 cpl_frame *loc_frame,
1906 cpl_frame *spectralformat_frame,
1907 cpl_frame *masterflat_frame,
1909 xsh_opt_extract_param* opt_extract_par,
1910 int min_index, int max_index,
1911 const char* rec_prefix,
1912 cpl_frame** orderext1d_frame,
1913 cpl_frame** orderoxt1d_frame,
1914 cpl_frame** res_frame_ext,
1915 cpl_frame** qc_subextract_frame,
1916 cpl_frame** qc_s2ddiv1d_frame,
1917 cpl_frame** qc_model_frame,
1918 cpl_frame** qc_weight_frame)
1919{
1920 xsh_pre *sci_pre = NULL;
1921 xsh_order_list *order_list = NULL;
1922 cpl_image *blaze_img = NULL;
1923 float *sci_pre_data = NULL;
1924 float *sci_pre_errs = NULL;
1925 int *sci_pre_qual = NULL;
1926 xsh_wavesol *wavesol = NULL;
1927 xsh_spectralformat_list *spectralformat_list= NULL;
1928
1929 int iorder, i;
1930 double *slits_cen = NULL;
1931 double *lambdas = NULL;
1932 double *pos_cen_x = NULL;
1933 double *pos_cen_y = NULL;
1934 int nlambdas = 0;
1935 cpl_vector *extract_lambda_vect = NULL;
1936 cpl_vector *result_lambda_vect = NULL;
1937 int result_lambda_vect_size;
1938
1939 cpl_vector *div_flux_vect = NULL;
1940 cpl_vector *div_err_vect = NULL;
1941
1942 cpl_vector *std_flux_vect = NULL;
1943 cpl_vector *std_err_vect = NULL;
1944 cpl_vector *std_qual_vect = NULL;
1945
1946 cpl_vector *opt_flux_vect = NULL;
1947 cpl_vector *opt_err_vect = NULL;
1948 cpl_vector *opt_qual_vect = NULL;
1949
1950 cpl_image *extract_img = NULL;
1951 cpl_image *extract_errs_img = NULL;
1952 cpl_image *extract_qual_img = NULL;
1953
1954 cpl_image *x_img = NULL;
1955 cpl_image *y_img = NULL;
1956
1957 cpl_image *sub_extract_img = NULL;
1958 cpl_image *sub_extract_errs_img = NULL;
1959 cpl_image *sub_extract_qual_img = NULL;
1960
1961 cpl_image *sub_x_img = NULL;
1962 cpl_image *sub_y_img = NULL;
1963
1964 cpl_image *s2Dby1D_img = NULL;
1965 cpl_image *s2Dby1D_err_img = NULL;
1966 cpl_image *weight_img = NULL;
1967 double *extract_data = NULL;
1968 double *extract_errs = NULL;
1969 int *extract_qual = NULL;
1970 double *extract_x_data = NULL;
1971 double *extract_y_data = NULL;
1972
1973 int nx_extract, ny_extract;
1974 int box_hsize;
1975 int niter=1, iiter;
1976 int chunk_ovsamp_size;
1977 int extract_slit_hsize=0;
1978 int deg_poly = 3;
1979 cpl_vector* std_extract_vect = NULL;
1980 cpl_vector* std_err_extract_vect = NULL;
1981 cpl_vector* std_qual_extract_vect = NULL;
1982
1983 cpl_vector *opt_extract_vect = NULL;
1984 cpl_vector *opt_err_extract_vect = NULL;
1985 cpl_vector *opt_qual_extract_vect = NULL;
1986
1987 double lambda_step=0.0;
1988 xsh_rec_list *orderext1d = NULL;
1989 xsh_rec_list *orderoxt1d = NULL;
1990 char tag_drl[256];
1991 char tag[256];
1992
1993 char filename[256];
1994 char filename_drl[256];
1995
1996 cpl_image *model_img = NULL;
1997 int rec_min_index, rec_max_index;
1998 double slit_down, slit_up;
1999 double lambda_min, lambda_max;
2000 int blaze_shift;
2001 xsh_xs_3 config_model;
2002 double conad = 0.0;
2003 double square_oversample;
2004 //int binx;
2005 int biny;
2006
2007 /* qc products */
2008 char qc_extname[256];
2009 cpl_propertylist *qc_header = NULL;
2010 char qc_subextract_name[256];
2011 char qc_s2ddiv1d_name[256];
2012 char qc_model_name[256];
2013 char qc_weight_name[256];
2014 int mode = CPL_IO_DEFAULT;
2015
2016
2017 int decode_bp=instrument->decode_bp;
2018 /* Checking input parameters */
2019 XSH_ASSURE_NOT_NULL( sci_frame);
2020 XSH_ASSURE_NOT_NULL( orderlist_frame);
2021 XSH_ASSURE_NOT_NULL( wavemap_frame);
2022 XSH_ASSURE_NOT_NULL( slitmap_frame);
2023 XSH_ASSURE_NOT_NULL( loc_frame);
2024 XSH_ASSURE_NOT_NULL( spectralformat_frame);
2025 XSH_ASSURE_NOT_NULL( masterflat_frame);
2027 XSH_ASSURE_NOT_NULL( opt_extract_par);
2028 XSH_ASSURE_NOT_NULL( orderext1d_frame);
2029 XSH_ASSURE_NOT_NULL( orderoxt1d_frame);
2030
2031 /* Loading data and initialize binning*/
2032 check( sci_pre = xsh_pre_load( sci_frame, instrument));
2033
2034 /* Square oversample */
2035 square_oversample = opt_extract_par->oversample*opt_extract_par->oversample;
2036 /* ADU to e- */
2037 conad = sci_pre->conad;
2038
2039 //binx = instrument->binx;
2040 biny = instrument->biny;
2041 if (model_frame == NULL){
2042 XSH_ASSURE_NOT_NULL( wavesol_frame);
2043 check( wavesol = xsh_wavesol_load( wavesol_frame, instrument));
2044 }
2045 else{
2046 check( xsh_model_config_load_best( model_frame, &config_model));
2047 check( xsh_model_binxy( &config_model, instrument->binx, instrument->biny));
2048 }
2049
2050 lambda_step = opt_extract_par->lambda_step;
2051 niter = opt_extract_par->niter;
2052 /* Loading data */
2053 check( sci_pre_data = cpl_image_get_data_float( sci_pre->data));
2054 check( sci_pre_errs = cpl_image_get_data_float( sci_pre->errs));
2055 check( sci_pre_qual = cpl_image_get_data_int( sci_pre->qual));
2056
2057 check( order_list = xsh_order_list_load( orderlist_frame, instrument));
2058
2059 check( spectralformat_list = xsh_spectralformat_list_load(
2060 spectralformat_frame, instrument));
2061
2062 nx_extract = sci_pre->ny*opt_extract_par->oversample;
2063 ny_extract = (OPT_EXTRACT_SLIT_SIZE * opt_extract_par->oversample)+1;
2064
2065 XSH_MALLOC( lambdas, double, nx_extract);
2066 XSH_MALLOC( slits_cen, double, nx_extract);
2067 XSH_MALLOC( pos_cen_x, double, nx_extract);
2068 XSH_MALLOC( pos_cen_y, double, nx_extract);
2069 XSH_MALLOC( extract_data, double, nx_extract*ny_extract);
2070 XSH_MALLOC( extract_errs, double, nx_extract*ny_extract);
2071 XSH_MALLOC( extract_qual, int, nx_extract*ny_extract);
2072
2073 XSH_MALLOC( extract_x_data, double, nx_extract*ny_extract);
2074 XSH_MALLOC( extract_y_data, double, nx_extract*ny_extract);
2075
2076 chunk_ovsamp_size = opt_extract_par->chunk_size*opt_extract_par->oversample;
2077 box_hsize = opt_extract_par->box_hsize*opt_extract_par->oversample;
2078
2079 xsh_msg( "Create and multiply by blaze function");
2080
2081 xsh_msg_dbg_low("---opt_extract : create blaze");
2082 check( blaze_img = xsh_create_blaze( masterflat_frame,
2083 order_list, instrument));
2084
2085 xsh_msg_dbg_low("---opt_extract : multiply by blaze");
2086 check( xsh_pre_multiply_image( sci_pre, blaze_img));
2087
2088#if REGDEBUG_BLAZE
2090 cpl_frame *norm_frame = NULL;
2091
2092 check( norm_frame = xsh_pre_save( sci_pre,"OPTEXT_MULT_BLAZE.fits",
2093 "OPTEXT_MULT_BLAZE",1));
2094 xsh_free_frame( &norm_frame);
2095 }
2096#endif
2097
2098 /* Create result */
2099 check( orderext1d = xsh_rec_list_create( instrument));
2100 check( orderoxt1d = xsh_rec_list_create( instrument));
2101
2102 if (min_index >= 0){
2103 rec_min_index = min_index;
2104 }
2105 else{
2106 rec_min_index = 0;
2107 }
2108 if ( max_index < orderext1d->size){
2109 rec_max_index = max_index;
2110 }
2111 else{
2112 rec_max_index = orderext1d->size-1;
2113 }
2114
2115 check( xsh_get_slit_edges( slitmap_frame, &slit_down, &slit_up,
2116 NULL, NULL, instrument));
2117
2118
2119 XSH_NEW_PROPERTYLIST( qc_header);
2120 /* Loop on order */
2121 for( iorder= rec_min_index; iorder <= rec_max_index; iorder++) {
2122 int abs_order,start_y, end_y;
2123 double cos_tilt=0, sin_tilt=0;
2124 int ny_extract_min=-1, ny_extract_max = -1;
2125
2126
2127
2128 abs_order = order_list->list[iorder].absorder;
2129 check( start_y = xsh_order_list_get_starty( order_list, iorder));
2130 check( end_y = xsh_order_list_get_endy( order_list, iorder));
2131 xsh_msg("Extracting order %d", abs_order);
2132 xsh_msg_dbg_low("---opt_extract : order %d [%d --> %d]", abs_order, start_y, end_y);
2133
2134 check( xsh_wavemap_lambda_range( wavemap_frame, slitmap_frame, start_y, end_y,
2135 opt_extract_par->oversample, order_list, iorder,
2136 pos_cen_x, pos_cen_y, abs_order, spectralformat_list, lambdas, slits_cen,
2137 &nlambdas, instrument));
2138
2139#if REGDEBUG_OPTEXTRACT
2140/* REGDEBUG */
2141{
2142 FILE *test_file;
2143 char test_name[256];
2144 float x, y, lambda, slit;
2145 double x_slit_min=0, y_slit_min=0;
2146 double x_slit_max=0, y_slit_max=0;
2147 double x_slit_cen=0, y_slit_cen=0;
2148 double diffx, diffy;
2149
2150 if (model_frame == NULL){
2151 sprintf( test_name, "poly%dx%d_lambda_range_%d.dat", binx, biny, abs_order);
2152 }
2153 else{
2154 sprintf( test_name, "model%dx%d_lambda_range_%d.dat", binx, biny, abs_order);
2155 }
2156 xsh_msg("Write %s", test_name);
2157 test_file = fopen( test_name, "w");
2158 fprintf( test_file, "# x y lambda slit slit_min slit_max x_slit_min y_slit_min x_slit_cen y_slit_cen x_slit_max y_slit_max diffx diffy\n");
2159
2160 for(i=0; i< nlambdas; i++){
2161 x = pos_cen_x[i];
2162 y = pos_cen_y[i];
2163 lambda = lambdas[i];
2164 slit = slits_cen[i];
2165
2166 if (model_frame == NULL){
2167 check( x_slit_min = xsh_wavesol_eval_polx( wavesol, lambda,(double) abs_order,
2168 slit_down));
2169 check( x_slit_max = xsh_wavesol_eval_polx( wavesol, lambda, (double)abs_order,
2170 slit_up));
2171 check( y_slit_min = xsh_wavesol_eval_poly( wavesol, lambda, (double)abs_order,
2172 slit_down));
2173 check( y_slit_max = xsh_wavesol_eval_poly( wavesol, lambda, (double)abs_order,
2174 slit_up));
2175 check( x_slit_cen = xsh_wavesol_eval_polx( wavesol, lambda,(double) abs_order,
2176 slit));
2177 check( y_slit_cen = xsh_wavesol_eval_poly( wavesol, lambda, (double)abs_order,
2178 slit));
2179 }
2180 else{
2181 check( xsh_model_get_xy( &config_model, instrument, lambda, (double) abs_order,
2182 slit_down, &x_slit_min, &y_slit_min));
2183 check( xsh_model_get_xy( &config_model, instrument, lambda, (double) abs_order,
2184 slit_up, &x_slit_max, &y_slit_max));
2185 check( xsh_model_get_xy( &config_model, instrument, lambda, (double) abs_order,
2186 slit, &x_slit_cen, &y_slit_cen));
2187 }
2188 diffx = x-x_slit_cen;
2189 diffy = y-y_slit_cen;
2190
2191 fprintf( test_file, "%f %f %f %f %f %f %f %f %f %f %f %f %f %f\n", x, y, lambda, slit, slit_down, slit_up,
2192 x_slit_min, y_slit_min, x_slit_cen, y_slit_cen, x_slit_max, y_slit_max, diffx, diffy);
2193 }
2194 fclose( test_file);
2195}
2196#endif
2197/* END OF REGDEBUG */
2198
2199
2200 xsh_msg_dbg_low("---opt_extract : (%f --> %f) in %d", lambdas[0], lambdas[nlambdas-1],
2201 nlambdas);
2202
2203 xsh_msg_dbg_low("---opt_extract : s (%f --> %f)", slit_down, slit_up);
2204
2205 for(i=0; i< nlambdas; i++){
2206 float x, y, lambda, slit;
2207 double x_slit_min=0, y_slit_min=0;
2208 double x_slit_max=0, y_slit_max=0;
2209 float x_min, x_max;
2210 float y_min, y_max;
2211 int l;
2212
2213 x = pos_cen_x[i];
2214 y = pos_cen_y[i];
2215 lambda = lambdas[i];
2216 slit = slits_cen[i];
2217 xsh_msg_dbg_high(" x y lambda slit => %f %f %f %f", x, y, lambda, slit);
2218
2219 if ( ( i%opt_extract_par->oversample) == 0){
2220 xsh_msg_dbg_high("evaluate tilt lambda order %f %f", lambda, (double) abs_order);
2221 /* compute the tilt */
2222 if (model_frame == NULL){
2223 check( x_slit_min = xsh_wavesol_eval_polx( wavesol, lambda,(double) abs_order,
2224 slit_down));
2225 check( x_slit_max = xsh_wavesol_eval_polx( wavesol, lambda, (double)abs_order,
2226 slit_up));
2227 check( y_slit_min = xsh_wavesol_eval_poly( wavesol, lambda, (double)abs_order,
2228 slit_down));
2229 check( y_slit_max = xsh_wavesol_eval_poly( wavesol, lambda, (double)abs_order,
2230 slit_up));
2231 }
2232 else{
2233 /* evaluate at slit_max slit_min */
2234 check( xsh_model_get_xy( &config_model, instrument, lambda, (double) abs_order,
2235 slit_down, &x_slit_min, &y_slit_min));
2236 check( xsh_model_get_xy( &config_model, instrument, lambda, (double) abs_order,
2237 slit_up, &x_slit_max, &y_slit_max));
2238 }
2239 if ( x_slit_min < x_slit_max){
2240 x_min = x_slit_min;
2241 x_max = x_slit_max;
2242 }
2243 else{
2244 x_min = x_slit_max;
2245 x_max = x_slit_min;
2246 }
2247 if ( y_slit_min < y_slit_max){
2248 y_min = y_slit_min;
2249 y_max = y_slit_max;
2250 }
2251 else{
2252 y_min = y_slit_max;
2253 y_max = y_slit_min;
2254 }
2255 xsh_msg_dbg_high("min(%f,%f) max(%f,%f)", x_min, y_min, x_max, y_max);
2256 /* Dont check solution in Y */
2257 if ( x_min >= 1 && x_max <= sci_pre->nx && (x_min <= x) &&
2258 ( x <= x_max) &&
2259 y_min >= 1 && y_max <= sci_pre->ny){
2260 double tilt_x, tilt_y;
2261 double hypotenuse_tilt;
2262
2263 if ((y_min >= y) || (y >= y_max)){
2264 xsh_msg_warning("y out of bounds : %f [%f,%f]", y, y_min,y_max);
2265 }
2266 tilt_x = x_slit_max-x_slit_min;
2267 tilt_y = y_slit_max-y_slit_min;
2268 xsh_msg_dbg_high("tilt_x %f tilt_y %f", tilt_x, tilt_y);
2269
2270 if ( extract_slit_hsize <= 0){
2271 extract_slit_hsize = (int)(fabs(tilt_x)/2.0*SLIT_USEFUL_WINDOW_FACTOR)+1;
2272 ny_extract = extract_slit_hsize*opt_extract_par->oversample*2+1;
2273 xsh_msg_dbg_high("SLIT hsize %d global %d", extract_slit_hsize,ny_extract);
2274 }
2275 hypotenuse_tilt = sqrt( tilt_y*tilt_y+tilt_x*tilt_x);
2276 cos_tilt = tilt_x/hypotenuse_tilt;
2277 sin_tilt = tilt_y/hypotenuse_tilt;
2278 }
2279 }
2280
2281 for (l=0; l< ny_extract; l++){
2282 float xs, fy, fx;
2283 double flux=-1, err=-1;
2284 int qual_val=0;
2285
2286 xs = (float)l/(float)opt_extract_par->oversample-
2287 extract_slit_hsize;
2288
2289 /* FITS to C convention */
2290 fy = xs*sin_tilt+y-1;
2291
2292 if ( (fy >= 0) && (fy < sci_pre->ny) ){
2293 int pixel=i+l*nlambdas;
2294 /* FITS to C convention */
2295 fx = xs*cos_tilt+x-1;
2296
2297 check(xsh_interpolate_linear( sci_pre_data, sci_pre_errs, sci_pre_qual, sci_pre->nx,
2298 sci_pre->ny, fx, fy, &flux, &err, &qual_val, FLUX_MODE));
2299
2300 extract_x_data[pixel] = fx;
2301 extract_y_data[pixel] = fy;
2302
2303 /* Flux in e- and divided by the factor square oversample */
2304 flux *= conad/square_oversample;
2305 extract_data[pixel] = flux;
2306 extract_errs[pixel] = err * conad /
2307 opt_extract_par->oversample;
2308 extract_qual[pixel] = qual_val;
2309 }
2310 else{
2311 break;
2312 }
2313 }
2314 }/* end extraction */
2315
2316 /* SLIT approximate position */
2317
2318 check( xsh_object_localize( slitmap_frame, loc_frame, opt_extract_par->oversample,
2319 box_hsize, nlambdas,
2320 ny_extract, extract_x_data,
2321 extract_y_data, &ny_extract_min, &ny_extract_max));
2322
2323 check( extract_img = cpl_image_wrap_double( nlambdas,
2324 ny_extract, extract_data));
2325 check( extract_errs_img = cpl_image_wrap_double( nlambdas,
2326 ny_extract, extract_errs));
2327 check( extract_qual_img = cpl_image_wrap_int( nlambdas,
2328 ny_extract, extract_qual));
2329
2330#if REGDEBUG_EXTRACT
2331 { char name[256];
2332 sprintf( name, "extract_%d.fits", abs_order);
2333 cpl_image_save( extract_img, name, CPL_BPP_IEEE_DOUBLE,
2334 NULL,CPL_IO_DEFAULT);
2335 xsh_msg("EXTRACT save %s",name);
2336 sprintf( name, "extract_errs_%d.fits", abs_order);
2337 cpl_image_save( extract_errs_img, name, CPL_BPP_IEEE_DOUBLE,
2338 NULL,CPL_IO_DEFAULT);
2339 sprintf( name, "extract_qual_%d.fits", abs_order);
2340 cpl_image_save( extract_qual_img, name, XSH_PRE_QUAL_BPP,
2341 NULL,CPL_IO_DEFAULT);
2342 }
2343#endif
2344
2345
2346 check( x_img = cpl_image_wrap_double( nlambdas,
2347 ny_extract, extract_x_data));
2348 check( y_img = cpl_image_wrap_double( nlambdas,
2349 ny_extract, extract_y_data));
2350
2351#if REGDEBUG_EXTRACT_XY
2352 { char name[256];
2353
2354 sprintf( name, "extract_x_%d.fits", abs_order);
2355 xsh_msg("EXTRACT_XY save %s",name);
2356 cpl_image_save( x_img, name, CPL_BPP_IEEE_DOUBLE, NULL,CPL_IO_DEFAULT);
2357 sprintf( name, "extract_y_%d.fits", abs_order);
2358 xsh_msg("EXTRACT_XY save %s",name);
2359 cpl_image_save( y_img, name, CPL_BPP_IEEE_DOUBLE, NULL,CPL_IO_DEFAULT);
2360 }
2361#endif
2362
2363
2364 check( sub_extract_img = cpl_image_extract( extract_img, 1,
2365 ny_extract_min+1, nlambdas, ny_extract_max+1));
2366 check( sub_extract_errs_img = cpl_image_extract( extract_errs_img, 1,
2367 ny_extract_min+1, nlambdas, ny_extract_max+1));
2368 check( sub_extract_qual_img = cpl_image_extract( extract_qual_img, 1,
2369 ny_extract_min+1, nlambdas, ny_extract_max+1));
2370
2371
2372 /* Write QC SUB EXTRACT PROD */
2373 sprintf( qc_subextract_name, "sub_extract.fits");
2374
2375 sprintf( qc_extname, "ORD%d_FLUX", abs_order);
2376 check( xsh_pfits_set_extname ( qc_header, qc_extname));
2377 cpl_image_save( sub_extract_img, qc_subextract_name, CPL_BPP_IEEE_DOUBLE,
2378 qc_header,mode);
2379
2380 sprintf( qc_extname, "ORD%d_ERRS", abs_order);
2381 check( xsh_pfits_set_extname ( qc_header, qc_extname));
2382 cpl_image_save( sub_extract_errs_img, qc_subextract_name, CPL_BPP_IEEE_DOUBLE,
2383 qc_header, CPL_IO_EXTEND);
2384
2385 sprintf( qc_extname, "ORD%d_QUAL", abs_order);
2386 check( xsh_pfits_set_extname ( qc_header, qc_extname));
2387 cpl_image_save( sub_extract_qual_img, qc_subextract_name, XSH_PRE_QUAL_BPP,
2388 qc_header, CPL_IO_EXTEND);
2389
2390
2391
2392 check( sub_x_img = cpl_image_extract( x_img, 1,
2393 ny_extract_min+1, nlambdas, ny_extract_max+1));
2394 check( sub_y_img = cpl_image_extract( y_img, 1,
2395 ny_extract_min+1, nlambdas, ny_extract_max+1));
2396
2397
2398#if REGDEBUG_SUBEXTRACT_XY
2399 { char name[256];
2400
2401 sprintf( name, "sub_x_%d.fits", abs_order);
2402 xsh_msg("EXTRACT_XY save %s",name);
2403 cpl_image_save( sub_x_img, name, CPL_BPP_IEEE_DOUBLE, NULL,CPL_IO_DEFAULT);
2404 sprintf( name, "sub_y_%d.fits", abs_order);
2405 xsh_msg("EXTRACT_XY save %s",name);
2406 cpl_image_save( sub_y_img, name, CPL_BPP_IEEE_DOUBLE, NULL,CPL_IO_DEFAULT);
2407 }
2408#endif
2409
2410 check( extract_lambda_vect = cpl_vector_wrap( nlambdas, lambdas));
2411
2412 /* do standard extraction */
2413 check( std_extract_vect = xsh_image_extract_standard( sub_extract_img,
2414 sub_extract_errs_img, sub_extract_qual_img, &std_err_extract_vect,
2415 &std_qual_extract_vect));
2416
2417#if REGDEBUG_EXTRACT_STD
2418 {
2419 char std_spectrum_name[256];
2420 FILE* std_spectrum_file = NULL;
2421 int nx = cpl_vector_get_size( std_extract_vect);
2422
2423 sprintf( std_spectrum_name, "extract_std_%d.dat", abs_order);
2424 std_spectrum_file = fopen( std_spectrum_name, "w");
2425
2426 fprintf(std_spectrum_file, "#index lambda flux err qual\n");
2427
2428 for(i=0; i< nx; i++){
2429 double lambdav;
2430 double fluxv;
2431 double lambda_err;
2432 double qual_v;
2433
2434 lambdav = cpl_vector_get( extract_lambda_vect, i);
2435 fluxv = cpl_vector_get( std_extract_vect, i);
2436 lambda_err = cpl_vector_get( std_err_extract_vect, i);
2437 qual_v = cpl_vector_get( std_qual_extract_vect, i);
2438 fprintf(std_spectrum_file, "%d %f %f %f %f\n", i, lambdav, fluxv,
2439 lambda_err, qual_v);
2440 }
2441 fclose( std_spectrum_file);
2442 }
2443#endif
2444
2445
2446 div_flux_vect = std_extract_vect;
2447 div_err_vect = std_err_extract_vect;
2448
2449 for( iiter=1; iiter <= niter; iiter++){
2450 double kappa = 3.0, frac_min=0.7;
2451 int niter = 2;
2452
2453 xsh_free_image( &s2Dby1D_img);
2454 xsh_free_image( &s2Dby1D_err_img);
2455 xsh_free_image( &model_img);
2456 xsh_free_image( &weight_img);
2457
2458 xsh_msg_dbg_low("---opt_extract : Do optimal extraction %d / %d", iiter, niter);
2459 /* Create 2D/1D image */
2460 check( s2Dby1D_img = xsh_image_divide_1D( sub_extract_img,
2461 sub_extract_errs_img, div_flux_vect, div_err_vect,
2462 &s2Dby1D_err_img));
2463
2464 sprintf( qc_s2ddiv1d_name, "s2Dby1D.fits");
2465 sprintf( qc_extname, "ORD%d_FLUX", abs_order);
2466 check( xsh_pfits_set_extname ( qc_header, qc_extname));
2467 cpl_image_save( s2Dby1D_img, qc_s2ddiv1d_name, CPL_BPP_IEEE_DOUBLE, qc_header,
2468 mode);
2469 sprintf( qc_extname, "ORD%d_ERRS", abs_order);
2470 cpl_image_save( s2Dby1D_err_img, qc_s2ddiv1d_name, CPL_BPP_IEEE_DOUBLE, qc_header,
2471 CPL_IO_EXTEND);
2472
2473 check( model_img = xsh_optextract_produce_model( s2Dby1D_img,
2474 opt_extract_par->method, chunk_ovsamp_size, deg_poly,
2475 opt_extract_par->oversample, extract_x_data, extract_y_data,
2476 abs_order, kappa, niter, frac_min));
2477
2478 sprintf( qc_model_name, "model.fits");
2479 sprintf( qc_extname, "ORD%d_FLUX", abs_order);
2480 check( xsh_pfits_set_extname ( qc_header, qc_extname));
2481 cpl_image_save( model_img, qc_model_name, CPL_BPP_IEEE_DOUBLE, qc_header, mode);
2482
2483 /* optimal extraction */
2484 xsh_free_vector( &opt_extract_vect);
2485 xsh_free_vector( &opt_err_extract_vect);
2486 xsh_free_vector( &opt_qual_extract_vect);
2487 check( opt_extract_vect = xsh_image_extract_optimal( sub_extract_img,
2488 sub_extract_errs_img, sub_extract_qual_img, opt_extract_par,
2489 model_img, decode_bp,&weight_img, &opt_err_extract_vect, &opt_qual_extract_vect));
2490 div_flux_vect = opt_extract_vect;
2491 div_err_vect = opt_err_extract_vect;
2492 }
2493
2494 sprintf( qc_weight_name, "weight.fits");
2495 sprintf( qc_extname, "ORD%d_FLUX", abs_order);
2496 check( xsh_pfits_set_extname ( qc_header, qc_extname));
2497 cpl_image_save( weight_img, qc_weight_name, CPL_BPP_IEEE_DOUBLE, qc_header,
2498 mode);
2499
2500 /* correct blaze */
2502 blaze_shift = (int)pos_cen_y[0];
2503 }
2504 else{
2505 int nx = cpl_vector_get_size( std_extract_vect);
2506
2507 blaze_shift = (int)pos_cen_y[nx-1];
2508 }
2509
2510 check( xsh_vector_divide_poly( std_extract_vect,
2511 opt_extract_par->oversample,
2512 order_list->list[iorder].blazepoly, blaze_shift, instrument));
2513
2514 check( xsh_vector_divide_poly( std_err_extract_vect,
2515 opt_extract_par->oversample,
2516 order_list->list[iorder].blazepoly, blaze_shift, instrument));
2517
2518#if REGDEBUG_BLAZE
2520 char std_spectrum_name[256];
2521 FILE* std_spectrum_file = NULL;
2522 int nx = cpl_vector_get_size( std_extract_vect);
2523
2524 sprintf( std_spectrum_name, "spectrum_divbyblaze_std_%d.dat", abs_order);
2525 std_spectrum_file = fopen( std_spectrum_name, "w");
2526
2527 fprintf( std_spectrum_file, "#index lambda flux err qual\n");
2528 for(i=0; i< nx; i++){
2529 double lambdav;
2530 double fluxv, flux_err, qual_v;
2531
2532 lambdav = cpl_vector_get( extract_lambda_vect, i);
2533 fluxv = cpl_vector_get( std_extract_vect, i);
2534 flux_err = cpl_vector_get( std_err_extract_vect, i);
2535 qual_v = cpl_vector_get( std_qual_extract_vect, i);
2536 fprintf(std_spectrum_file, "%d %f %f %f %f\n", i, lambdav, fluxv,
2537 flux_err, qual_v);
2538 }
2539 fclose( std_spectrum_file);
2540 }
2541#endif
2542
2543 check( xsh_vector_divide_poly( opt_extract_vect,
2544 opt_extract_par->oversample,
2545 order_list->list[iorder].blazepoly, blaze_shift, instrument));
2546
2547 check( xsh_vector_divide_poly( opt_err_extract_vect,
2548 opt_extract_par->oversample,
2549 order_list->list[iorder].blazepoly, blaze_shift, instrument));
2550
2551 /* convert in lambda spaces */
2552 check( xsh_interpolate_spectrum( biny, abs_order, opt_extract_par->oversample,
2553 opt_extract_par->lambda_step,
2554 extract_lambda_vect,
2555 std_extract_vect, std_err_extract_vect, std_qual_extract_vect,
2556 opt_extract_vect, opt_err_extract_vect, opt_qual_extract_vect,
2557 &result_lambda_vect,
2558 &std_flux_vect, &std_err_vect, &std_qual_vect,
2559 &opt_flux_vect, &opt_err_vect, &opt_qual_vect));
2560
2561 result_lambda_vect_size = cpl_vector_get_size( result_lambda_vect);
2562
2563 /* REGDEBUG */
2565 char std_spectrum_name[256];
2566 FILE* std_spectrum_file = NULL;
2567
2568 sprintf( std_spectrum_name, "spectrum_std_%d.dat", abs_order);
2569 std_spectrum_file = fopen( std_spectrum_name, "w");
2570
2571 fprintf(std_spectrum_file, "#lambda flux err qual sn %d\n",
2572 opt_extract_par->oversample);
2573 for(i=0; i< result_lambda_vect_size; i++){
2574 double lambdav, flux, err;
2575 double qual;
2576
2577 lambdav = cpl_vector_get( result_lambda_vect, i);
2578 flux = cpl_vector_get( std_flux_vect, i);
2579 err = cpl_vector_get( std_err_vect, i);
2580 qual = cpl_vector_get( std_qual_vect, i);
2581 fprintf(std_spectrum_file, "%f %f %f %f %f\n", lambdav, flux / conad,
2582 err / conad, qual, flux/err);
2583 }
2584 fclose( std_spectrum_file);
2585
2586 sprintf( std_spectrum_name, "spectrum_opt_%d.dat", abs_order);
2587 std_spectrum_file = fopen( std_spectrum_name, "w");
2588
2589 fprintf(std_spectrum_file, "#lambda flux err qual sn %d\n",
2590 opt_extract_par->oversample);
2591 for(i=0; i< result_lambda_vect_size; i++){
2592 double lambdav, flux, err;
2593 double qual;
2594
2595 lambdav = cpl_vector_get( result_lambda_vect, i);
2596 flux = cpl_vector_get( opt_flux_vect, i);
2597 err = cpl_vector_get( opt_err_vect, i);
2598 qual = cpl_vector_get( opt_qual_vect, i);
2599
2600 fprintf(std_spectrum_file, "%f %f %f %f %f\n", lambdav, flux / conad,
2601 err / conad, qual, flux/err);
2602 }
2603 fclose( std_spectrum_file);
2604 }
2605 /* save extracted order in rectify list */
2606 xsh_rec_list_set_data_size( orderext1d, iorder, abs_order,
2607 result_lambda_vect_size, 1);
2608 xsh_rec_list_set_data_size( orderoxt1d, iorder, abs_order,
2609 result_lambda_vect_size, 1);
2610 for(i=0; i< result_lambda_vect_size; i++){
2611 double lambdav, extflux,oxtflux;
2612 double exterr, oxterr;
2613 int extqual, oxtqual;
2614
2615 lambdav = cpl_vector_get( result_lambda_vect, i);
2616 extflux = cpl_vector_get( std_flux_vect, i);
2617 oxtflux = cpl_vector_get( opt_flux_vect, i);
2618
2619 exterr = cpl_vector_get( std_err_vect, i);
2620 oxterr = cpl_vector_get( opt_err_vect, i);
2621
2622 extqual = (int) cpl_vector_get( std_qual_vect, i);
2623 oxtqual = (int) cpl_vector_get( opt_qual_vect, i);
2624
2625 orderext1d->list[iorder].lambda[i] = lambdav;
2626 orderoxt1d->list[iorder].lambda[i] = lambdav;
2627
2628 /* e- to ADU */
2629 orderext1d->list[iorder].data1[i] = extflux / conad;
2630 orderoxt1d->list[iorder].data1[i] = oxtflux / conad;
2631
2632 orderext1d->list[iorder].errs1[i] = exterr / conad;
2633 orderoxt1d->list[iorder].errs1[i] = oxterr / conad;
2634
2635 orderext1d->list[iorder].qual1[i] = extqual;
2636 orderoxt1d->list[iorder].qual1[i] = oxtqual;
2637 mode = CPL_IO_EXTEND;
2638
2639 }
2640 /* free memory */
2641 xsh_free_image( &s2Dby1D_img);
2642 xsh_free_image( &s2Dby1D_err_img);
2643 xsh_free_image( &model_img);
2644 xsh_free_image( &weight_img);
2645 xsh_unwrap_image( &extract_img);
2646 xsh_unwrap_image( &extract_errs_img);
2647 xsh_unwrap_image( &extract_qual_img);
2648 xsh_unwrap_image( &x_img);
2649 xsh_unwrap_image( &y_img);
2650 xsh_free_image( &sub_extract_img);
2651 xsh_free_image( &sub_extract_errs_img);
2652 xsh_free_image( &sub_extract_qual_img);
2653 xsh_free_image( &sub_x_img);
2654 xsh_free_image( &sub_y_img);
2655
2656 xsh_free_vector( &std_extract_vect);
2657 xsh_free_vector( &std_err_extract_vect);
2658 xsh_free_vector( &std_qual_extract_vect);
2659
2660 xsh_free_vector( &opt_extract_vect);
2661 xsh_free_vector( &opt_err_extract_vect);
2662 xsh_free_vector( &opt_qual_extract_vect);
2663
2664 xsh_unwrap_vector( &extract_lambda_vect);
2665 xsh_free_vector( &result_lambda_vect);
2666 xsh_free_vector( &std_flux_vect);
2667 xsh_free_vector( &std_err_vect);
2668 xsh_free_vector( &std_qual_vect);
2669 xsh_free_vector( &opt_flux_vect);
2670 xsh_free_vector( &opt_err_vect);
2671 xsh_free_vector( &opt_qual_vect);
2672 }
2673
2674 /* Save orderext1d rectify list */
2675 /* Set header from science header */
2676 check( cpl_propertylist_append( orderext1d->header, sci_pre->data_header));
2677 /* Add useful keywords */
2679 check( xsh_pfits_set_rectify_bin_space( orderext1d->header, 0.0));
2680
2681 check( lambda_min = xsh_rec_list_get_lambda_min( orderext1d));
2682 check( xsh_pfits_set_rectify_lambda_min( orderext1d->header, lambda_min));
2683 check( lambda_max = xsh_rec_list_get_lambda_max( orderext1d));
2684 check( xsh_pfits_set_rectify_lambda_max( orderext1d->header, lambda_max));
2685 sprintf(tag,"%s_%s",rec_prefix,XSH_GET_TAG_FROM_ARM(XSH_ORDER_EXT1D,instrument));
2686
2687 check( xsh_pfits_set_pcatg( orderext1d->header, tag));
2688 check( xsh_pfits_set_rectify_space_min( orderext1d->header, 0.0)) ;
2689 check( xsh_pfits_set_rectify_space_max( orderext1d->header, 0.0));
2690 sprintf(filename,"%s.fits",tag);
2691 check( *orderext1d_frame = xsh_rec_list_save( orderext1d,
2692 filename, tag, CPL_TRUE));
2693
2694
2695 /* Save orderoxt1d rectify list */
2696 /* Set header from science header */
2697 check( cpl_propertylist_append( orderoxt1d->header, sci_pre->data_header));
2698 /* Add useful keywords */
2700 check( xsh_pfits_set_rectify_bin_space( orderoxt1d->header, 0.0));
2701 check( lambda_min = xsh_rec_list_get_lambda_min( orderoxt1d));
2702 check( xsh_pfits_set_rectify_lambda_min(orderoxt1d->header, lambda_min));
2703 check( lambda_max = xsh_rec_list_get_lambda_max( orderoxt1d));
2704 check( xsh_pfits_set_rectify_lambda_max( orderoxt1d->header, lambda_max));
2705 sprintf(tag,"%s_%s",rec_prefix,XSH_GET_TAG_FROM_ARM(XSH_ORDER_OXT1D,instrument));
2706 sprintf(tag_drl,"%s_DRL",tag);
2707
2708 check( xsh_pfits_set_pcatg( orderoxt1d->header, tag_drl));
2709 check( xsh_pfits_set_rectify_space_min( orderoxt1d->header, 0.0)) ;
2710 check( xsh_pfits_set_rectify_space_max( orderoxt1d->header, 0.0));
2711 sprintf(filename,"%s.fits",tag);
2712 sprintf(filename_drl,"DRL_%s.fits",filename);
2713 xsh_add_temporary_file(filename);
2714 xsh_add_temporary_file(filename_drl);
2715
2716 check( *orderoxt1d_frame = xsh_rec_list_save( orderoxt1d, filename_drl,
2717 tag_drl,CPL_TRUE));
2718
2719 check( *res_frame_ext = xsh_rec_list_save2( orderoxt1d,filename,tag));
2720
2721 /* qc frames */
2722 XSH_ASSURE_NOT_NULL( qc_subextract_frame);
2723 XSH_ASSURE_NOT_NULL( qc_s2ddiv1d_frame);
2724 XSH_ASSURE_NOT_NULL( qc_model_frame);
2725 XSH_ASSURE_NOT_NULL( qc_weight_frame);
2726
2727 *qc_subextract_frame = cpl_frame_new();
2728
2729 sprintf(tag,"%s_OXT_SUBEXTRACT",rec_prefix);
2730 check ((cpl_frame_set_filename ( *qc_subextract_frame, qc_subextract_name),
2731 cpl_frame_set_tag( *qc_subextract_frame, tag),
2732 cpl_frame_set_type ( *qc_subextract_frame, CPL_FRAME_TYPE_IMAGE) )) ;
2733
2734 *qc_s2ddiv1d_frame = cpl_frame_new();
2735 sprintf(tag,"%s_OXT_S2DDIV1D",rec_prefix);
2736 check ((cpl_frame_set_filename ( *qc_s2ddiv1d_frame, qc_s2ddiv1d_name),
2737 cpl_frame_set_tag( *qc_s2ddiv1d_frame, tag),
2738 cpl_frame_set_type ( *qc_s2ddiv1d_frame, CPL_FRAME_TYPE_IMAGE) )) ;
2739
2740 *qc_model_frame = cpl_frame_new();
2741 sprintf(tag,"%s_OXT_MODEL",rec_prefix);
2742 check ((cpl_frame_set_filename ( *qc_model_frame, qc_model_name),
2743 cpl_frame_set_tag( *qc_model_frame, tag),
2744 cpl_frame_set_type ( *qc_model_frame, CPL_FRAME_TYPE_IMAGE) ));
2745
2746 *qc_weight_frame = cpl_frame_new();
2747 sprintf(tag,"%s_OXT_WEIGHT",rec_prefix);
2748 check ((cpl_frame_set_filename ( *qc_weight_frame, qc_weight_name),
2749 cpl_frame_set_tag( *qc_weight_frame, tag),
2750 cpl_frame_set_type ( *qc_weight_frame, CPL_FRAME_TYPE_IMAGE) ));
2751
2752 check( xsh_add_temporary_file( qc_subextract_name));
2753 check( xsh_add_temporary_file( qc_s2ddiv1d_name));
2754 check( xsh_add_temporary_file( qc_model_name));
2755 check( xsh_add_temporary_file( qc_weight_name));
2756
2757 cleanup:
2758 if ( cpl_error_get_code() != CPL_ERROR_NONE) {
2759 xsh_unwrap_image( &sub_extract_img);
2760 xsh_unwrap_image( &sub_extract_errs_img);
2761 xsh_unwrap_image( &sub_extract_qual_img);
2762
2763 xsh_free_vector( &std_extract_vect);
2764 xsh_free_vector( &std_err_extract_vect);
2765 xsh_free_vector( &std_qual_extract_vect);
2766
2767 xsh_free_vector( &opt_extract_vect);
2768 xsh_free_vector( &opt_err_extract_vect);
2769 xsh_free_vector( &opt_qual_extract_vect);
2770
2771 xsh_unwrap_vector( &extract_lambda_vect);
2772
2773 xsh_free_vector( &result_lambda_vect);
2774 xsh_free_vector( &std_flux_vect);
2775 xsh_free_vector( &std_err_vect);
2776 xsh_free_vector( &std_qual_vect);
2777 xsh_free_vector( &opt_flux_vect);
2778 xsh_free_vector( &opt_err_vect);
2779 xsh_free_vector( &opt_qual_vect);
2780 }
2781 xsh_pre_free( &sci_pre);
2782 xsh_order_list_free( &order_list);
2783 xsh_spectralformat_list_free( &spectralformat_list);
2784 xsh_wavesol_free( &wavesol);
2785 XSH_FREE( lambdas);
2786 XSH_FREE( pos_cen_x);
2787 XSH_FREE( pos_cen_y);
2788 XSH_FREE( extract_data);
2789 XSH_FREE( extract_x_data);
2790 XSH_FREE( extract_y_data);
2791 XSH_FREE( extract_errs);
2792 XSH_FREE( extract_qual);
2793 xsh_free_image( &blaze_img);
2794 xsh_rec_list_free( &orderext1d);
2795 xsh_rec_list_free( &orderoxt1d);
2796 xsh_free_propertylist( &qc_header);
2797 return;
2798}
2799/*----------------------------------------------------------------------------*/
2800/* END xsh_opt_extract */
2801/*----------------------------------------------------------------------------*/
2802
static int starty
static int width
static int endy
static const double step
static char mode[32]
static double sigma
static xsh_instrument * instrument
int binx
int biny
static float lambda_step
cpl_image * xsh_create_blaze(cpl_frame *masterflat_frame, xsh_order_list *order_list, xsh_instrument *instrument)
Normalize a master flat frame order by order.
Definition: xsh_blaze.c:197
xsh_localization * xsh_localization_load(cpl_frame *frame)
Load a localization list from a frame.
void xsh_localization_free(xsh_localization **list)
free memory associated to a localization_list
xsh_order_list * xsh_order_list_load(cpl_frame *frame, xsh_instrument *instr)
load an order list from a frame
int xsh_order_list_get_starty(xsh_order_list *list, int i)
get position on Y axis of first pixel detected on order
int xsh_order_list_get_endy(xsh_order_list *list, int i)
get position on Y axis of last pixel detected on order
void xsh_order_list_free(xsh_order_list **list)
free memory associated to an order_list
double xsh_order_list_eval(xsh_order_list *list, cpl_polynomial *poly, double y)
Evaluate an order list poly.
xsh_pre * xsh_pre_load(cpl_frame *frame, xsh_instrument *instr)
Load a xsh_pre structure from a frame.
Definition: xsh_data_pre.c:849
void xsh_pre_multiply_image(const xsh_pre *pre, cpl_image *img)
multiply a frame in PRE format by an image
void xsh_pre_free(xsh_pre **pre)
Free a xsh_pre structure.
Definition: xsh_data_pre.c:823
cpl_frame * xsh_pre_save(const xsh_pre *pre, const char *filename, const char *tag, int temp)
Save PRE on disk.
cpl_frame * xsh_rec_list_save(xsh_rec_list *list, const char *filename, const char *tag, int is_temp)
Save a rec list in a frame.
cpl_frame * xsh_rec_list_save2(xsh_rec_list *list, const char *filename, const char *tag)
save an rec list to a frame
void xsh_rec_list_set_data_size(xsh_rec_list *list, int idx, int absorder, int nlambda, int ns)
Allocate memory for the order idx of the rectify list.
Definition: xsh_data_rec.c:191
double xsh_rec_list_get_lambda_max(xsh_rec_list *list)
Definition: xsh_data_rec.c:957
double xsh_rec_list_get_lambda_min(xsh_rec_list *list)
Definition: xsh_data_rec.c:935
void xsh_rec_list_free(xsh_rec_list **list)
free memory associated to a rec_list
Definition: xsh_data_rec.c:760
xsh_rec_list * xsh_rec_list_create(xsh_instrument *instr)
Create an empty order list.
Definition: xsh_data_rec.c:100
xsh_spectralformat_list * xsh_spectralformat_list_load(cpl_frame *frame, xsh_instrument *instr)
Load a spectralformat list from a frame.
float xsh_spectralformat_list_get_lambda_max(xsh_spectralformat_list *list, int absorder)
float xsh_spectralformat_list_get_lambda_min(xsh_spectralformat_list *list, int absorder)
Returns lambda min for a given absolute order.
void xsh_spectralformat_list_free(xsh_spectralformat_list **list)
Free memory associated to an spactralformat_list.
xsh_wavesol * xsh_wavesol_load(cpl_frame *frame, xsh_instrument *instrument)
load a wavelength solution
double xsh_wavesol_eval_poly(xsh_wavesol *sol, double lambda, double order, double slit)
eval the polynomial solution in Y
void xsh_wavesol_free(xsh_wavesol **w)
free wavelength solution structure
double xsh_wavesol_eval_polx(xsh_wavesol *sol, double lambda, double order, double slit)
eval the polynomial solution in X
#define XSH_ASSURE_NOT_ILLEGAL(cond)
Definition: xsh_error.h:107
#define check(COMMAND)
Definition: xsh_error.h:71
#define xsh_error_reset()
Definition: xsh_error.h:87
#define XSH_ASSURE_NOT_NULL(pointer)
Definition: xsh_error.h:99
XSH_ARM xsh_instrument_get_arm(xsh_instrument *i)
Get an arm on instrument structure.
cpl_error_code xsh_model_config_load_best(cpl_frame *config_frame, xsh_xs_3 *p_xs_3)
Load the config model table and fill the struct.
Definition: xsh_model_io.c:174
void xsh_model_get_xy(xsh_xs_3 *p_xs_3, xsh_instrument *instr, double lambda_nm, int morder, double ent_slit_pos, double *x, double *y)
Compute the detector location (floating point pixels) of a given wavelength/entrance slit position.
void xsh_model_binxy(xsh_xs_3 *p_xs_3, int bin_X, int bin_Y)
corrects model for detector's binning
int size
int * y
int * x
static SimAnneal s
Definition: xsh_model_sa.c:99
#define xsh_msg_warning(...)
Print an warning message.
Definition: xsh_msg.h:88
#define xsh_msg_dbg_medium(...)
Definition: xsh_msg.h:44
#define xsh_msg_debug(...)
Print a debug message.
Definition: xsh_msg.h:99
#define xsh_msg(...)
Print a message on info level.
Definition: xsh_msg.h:121
#define xsh_msg_dbg_low(...)
Definition: xsh_msg.h:48
#define xsh_msg_dbg_high(...)
Definition: xsh_msg.h:40
static void xsh_wavemap_lambda_range(cpl_frame *wavemap_frame, cpl_frame *slitmap_frame, int starty, int endy, int oversample, xsh_order_list *order_list, int iorder, double *xtab, double *ytab, int order, xsh_spectralformat_list *spectralformat, double *lambdastab, double *slitstab, int *sizetab, xsh_instrument *instr)
Give the value of lambda and x central position for a given y.
static cpl_vector * xsh_image_extract_standard(cpl_image *img, cpl_image *err_img, cpl_image *qual_img, cpl_vector **err_v, cpl_vector **qual_v)
Do a standard extraction on image.
#define OPT_EXTRACT_SLIT_SIZE
static cpl_vector * xsh_image_extract_optimal(cpl_image *img, cpl_image *errs_img, cpl_image *qual_img, xsh_opt_extract_param *opt_par, cpl_image *model_img, const int decode_bp, cpl_image **corr_img, cpl_vector **err_v, cpl_vector **qual_v)
Do an optimal extraction on image.
static int xsh_interpolate_linear(float *fluxtab, float *errtab, int *qualtab, int nx, int ny, float pos_x, float pos_y, double *flux, double *err, int *qual, int mode)
static cpl_image * xsh_image_create_model_image(cpl_image *src_img, double *x_data, double *y_data, double kappa, int niter, double frac_min)
#define FLUX_MODE
static cpl_image * xsh_image_create_gaussian_image(cpl_image *src, cpl_polynomial *centerp, cpl_polynomial *heightp, cpl_polynomial *widthp, cpl_polynomial *offsetp)
static void xsh_interpolate_spectrum(int biny, int oversample, int absorder, double lambda_step, cpl_vector *init_pos, cpl_vector *std_flux, cpl_vector *std_err, cpl_vector *std_qual, cpl_vector *opt_flux, cpl_vector *opt_err, cpl_vector *opt_qual, cpl_vector **res_pos, cpl_vector **res_std_flux, cpl_vector **res_std_err, cpl_vector **res_std_qual, cpl_vector **res_opt_flux, cpl_vector **res_opt_err, cpl_vector **res_opt_qual)
static cpl_vector * xsh_vector_integrate(int biny, int absorder, int oversample, cpl_vector *ref_pos, double step, cpl_vector *ref_values, cpl_vector *err_values, cpl_vector *qual_values, cpl_vector *new_pos, cpl_vector **spectrum_err, cpl_vector **spectrum_qual)
Interpolate values following given positions.
static void xsh_image_gaussian_fit_y(cpl_image *img, int chunk_size, int deg_poly, int oversample, cpl_polynomial **center, cpl_polynomial **height, cpl_polynomial **width, cpl_polynomial **offset)
Do a gaussian fit of Y columns by chunk and fit the position by polynomials in X.
cpl_image * xsh_optextract_produce_model(cpl_image *s2Dby1D_img, int method, int chunk_ovsamp_size, int deg_poly, int oversample, double *extract_x_data, double *extract_y_data, int abs_order, double kappa, int niter, double frac_min)
static cpl_image * xsh_image_divide_1D(cpl_image *image, cpl_image *err_img, cpl_vector *s1D, cpl_vector *err_s1D, cpl_image **err_2dby1D_img)
Divide an image by a 1D spectrum (same scale)
void xsh_opt_extract(cpl_frame *sci_frame, cpl_frame *orderlist_frame, cpl_frame *wavesol_frame, cpl_frame *model_frame, cpl_frame *wavemap_frame, cpl_frame *slitmap_frame, cpl_frame *loc_frame, cpl_frame *spectralformat_frame, cpl_frame *masterflat_frame, xsh_instrument *instrument, xsh_opt_extract_param *opt_extract_par, const char *rec_prefix, cpl_frame **orderext1d_frame, cpl_frame **orderoxt1d_frame, cpl_frame **orderoxt1d_eso_frame, cpl_frame **qc_subextract_frame, cpl_frame **qc_s2ddiv1d_frame, cpl_frame **qc_model_frame, cpl_frame **qc_weight_frame)
Create a SPECTRUM 1D from a Science frame.
static void xsh_object_localize(cpl_frame *slitmap_frame, cpl_frame *loc_frame, int oversample, int box_hsize, int nlambdas, int ny_extract, double *extract_x_data, double *extract_y_data, int *ymin, int *ymax)
#define WAVEMAP_MODE
#define SLIT_USEFUL_WINDOW_FACTOR
void xsh_opt_extract_orders(cpl_frame *sci_frame, cpl_frame *orderlist_frame, cpl_frame *wavesol_frame, cpl_frame *model_frame, cpl_frame *wavemap_frame, cpl_frame *slitmap_frame, cpl_frame *loc_frame, cpl_frame *spectralformat_frame, cpl_frame *masterflat_frame, xsh_instrument *instrument, xsh_opt_extract_param *opt_extract_par, int min_index, int max_index, const char *rec_prefix, cpl_frame **orderext1d_frame, cpl_frame **orderoxt1d_frame, cpl_frame **res_frame_ext, cpl_frame **qc_subextract_frame, cpl_frame **qc_s2ddiv1d_frame, cpl_frame **qc_model_frame, cpl_frame **qc_weight_frame)
Create a SPECTRUM 1D from a Science frame.
static void xsh_vector_divide_poly(cpl_vector *vector, double oversample, cpl_polynomial *poly, int shift, xsh_instrument *instr)
Divide vector values by a polynomial.
#define FIT_FWHM_LIMIT
void xsh_pfits_set_pcatg(cpl_propertylist *plist, const char *value)
Write the PCATG value.
Definition: xsh_pfits.c:1008
void xsh_pfits_set_rectify_bin_lambda(cpl_propertylist *plist, double value)
WRITE the lambda binning.
Definition: xsh_pfits.c:3188
void xsh_pfits_set_rectify_lambda_max(cpl_propertylist *plist, double value)
WRITE the lambda max value.
Definition: xsh_pfits.c:3236
void xsh_pfits_set_rectify_space_max(cpl_propertylist *plist, double value)
WRITE the space (slit) max value.
Definition: xsh_pfits.c:3268
void xsh_pfits_set_extname(cpl_propertylist *plist, const char *value)
Write the EXTNAME value.
Definition: xsh_pfits.c:979
void xsh_pfits_set_rectify_space_min(cpl_propertylist *plist, double value)
WRITE the space (slit) min value.
Definition: xsh_pfits.c:3252
void xsh_pfits_set_rectify_bin_space(cpl_propertylist *plist, double value)
WRITE the space (slit) binning.
Definition: xsh_pfits.c:3204
void xsh_pfits_set_rectify_lambda_min(cpl_propertylist *plist, double value)
WRITE the lambda min value.
Definition: xsh_pfits.c:3220
void xsh_get_slit_edges(cpl_frame *slitmap_frame, double *sdown, double *sup, double *sldown, double *slup, xsh_instrument *instrument)
Trace slit edges in a master flat.
Definition: xsh_rectify.c:671
void xsh_unwrap_vector(cpl_vector **v)
Unwrap a vector and set the pointer to NULL.
Definition: xsh_utils.c:2345
void xsh_free_polynomial(cpl_polynomial **p)
Deallocate a polynomial and set the pointer to NULL.
Definition: xsh_utils.c:2194
void xsh_free_vector(cpl_vector **v)
Deallocate a vector and set the pointer to NULL.
Definition: xsh_utils.c:2284
void xsh_free_image(cpl_image **i)
Deallocate an image and set the pointer to NULL.
Definition: xsh_utils.c:2116
int xsh_debug_level_get(void)
get debug level
Definition: xsh_utils.c:3142
void xsh_free_frame(cpl_frame **f)
Deallocate a frame and set the pointer to NULL.
Definition: xsh_utils.c:2269
void xsh_array_clip_poly1d(cpl_vector *pos_vect, cpl_vector *val_vect, double kappa, int niter, double frac_min, int deg, cpl_polynomial **polyp, double *chisq, int **flagsp)
clip outliers from a 1D poly fit
Definition: xsh_utils.c:6568
long xsh_round_double(double x)
Computes round(x)
Definition: xsh_utils.c:4383
void xsh_unwrap_image(cpl_image **i)
Unwrap an image and set the pointer to NULL.
Definition: xsh_utils.c:2330
void xsh_free_propertylist(cpl_propertylist **p)
Deallocate a property list and set the pointer to NULL.
Definition: xsh_utils.c:2179
void xsh_add_temporary_file(const char *name)
Add temporary file to temprary files list.
Definition: xsh_utils.c:1432
cpl_polynomial * cenpoly
enum optextract_method method
xsh_order * list
cpl_polynomial * blazepoly
cpl_polynomial * cenpoly
cpl_image * qual
Definition: xsh_data_pre.h:71
double conad
Definition: xsh_data_pre.h:96
cpl_image * data
Definition: xsh_data_pre.h:65
cpl_propertylist * data_header
Definition: xsh_data_pre.h:66
cpl_image * errs
Definition: xsh_data_pre.h:68
cpl_propertylist * header
Definition: xsh_data_rec.h:89
xsh_rec * list
Definition: xsh_data_rec.h:87
float * errs1
Definition: xsh_data_rec.h:71
float * data1
Definition: xsh_data_rec.h:68
int * qual1
Definition: xsh_data_rec.h:75
double * lambda
Definition: xsh_data_rec.h:67
#define QFLAG_COSMIC_RAY_REMOVED
@ XSH_ARM_NIR
#define XSH_PRE_QUAL_BPP
Definition: xsh_data_pre.h:47
int nx
double kappa
Definition: xsh_detmon_lg.c:81
int niter
Definition: xsh_detmon_lg.c:82
int n
Definition: xsh_detmon_lg.c:92
const char * method
Definition: xsh_detmon_lg.c:78
int ny
int order
Definition: xsh_detmon_lg.c:80
#define XSH_ORDER_EXT1D
Definition: xsh_dfs.h:1029
#define XSH_GET_TAG_FROM_ARM(TAG, instr)
Definition: xsh_dfs.h:1548
#define XSH_ORDER_OXT1D
Definition: xsh_dfs.h:1034
#define XSH_MATH_SQRT_2
Definition: xsh_drl.h:567
@ GAUSS_METHOD
#define XSH_NEW_PROPERTYLIST(POINTER)
Definition: xsh_utils.h:70
#define XSH_FREE(POINTER)
Definition: xsh_utils.h:92
@ XSH_DEBUG_LEVEL_MEDIUM
Definition: xsh_utils.h:138
#define XSH_MALLOC(POINTER, TYPE, SIZE)
Definition: xsh_utils.h:49
#define XSH_CALLOC(POINTER, TYPE, SIZE)
Definition: xsh_utils.h:56
cpl_polynomial * xsh_polynomial_fit_1d_create(const cpl_vector *x_pos, const cpl_vector *values, int degree, double *mse)