X-shooter Pipeline Reference Manual 3.8.15
xsh_utils_ifu.c
Go to the documentation of this file.
1/* a
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 14:15:45 $
24 * $Revision: 1.42 $
25 * $Name: not supported by cvs2svn $
26 */
27/*----------------------------------------------------------------------------
28 Includes
29 ----------------------------------------------------------------------------*/
30#include <math.h>
31#include <string.h>
32#include <xsh_error.h>
33#include <xsh_utils_ifu.h>
34#include <xsh_utils_image.h>
35#include <xsh_utils_imagelist.h>
36#include <xsh_utils_wrappers.h>
37#include <xsh_utils_table.h>
38#include <xsh_pfits.h>
39#include <xsh_model_io.h>
40#include <xsh_data_rec.h>
41#include <xsh_data_spectrum.h>
42#include <xsh_dfs.h>
43#include <xsh_drl.h>
44#include <xsh_globals.h>
45
46#define PI_NUMB (3.1415926535897932384626433832795)
47
48static double
49xsh_interpol(double x,double x1,double x2,double y1,double y2)
50{
51 return y1+(y2-y1)/(x2-x1)*(x-x1);
52}
53
54
55static void
56xsh_plist_set_coord1(cpl_propertylist** plist,
57 const double crpix1,
58 const double crval1,
59 const double cdelt1);
60static void
61xsh_plist_set_coord2(cpl_propertylist** plist,
62 const double crpix2,
63 const double crval2,
64 const double cdelt2);
65static void
66xsh_plist_set_coord3(cpl_propertylist** plist,
67 const int crpix3,
68 const double crval3,
69 const double cdelt3);
70
71
72static void
73xsh_plist_set_cd_matrix2(cpl_propertylist** plist,
74 const double cd1_1,
75 const double cd1_2,
76 const double cd2_1,
77 const double cd2_2);
78
79
80static void
81xsh_plist_set_cd_matrix3(cpl_propertylist** plist,
82 const double cd1_3,
83 const double cd2_3,
84 const double cd3_1,
85 const double cd3_2,
86 const double cd3_3);
87
88
89
98static void
99xsh_plist_set_coord1(cpl_propertylist** plist,
100 const double crpix1,
101 const double crval1,
102 const double cdelt1)
103{
104 cpl_propertylist_erase_regexp(*plist, "^CTYPE1",0);
105 cpl_propertylist_insert_after_string(*plist,"EXPTIME","CTYPE1","RA---TAN");
106 cpl_propertylist_set_comment(*plist, "CTYPE1", "Projected Rectascension");
107 cpl_propertylist_erase_regexp(*plist, "^CRPIX1",0);
108 cpl_propertylist_insert_after_double(*plist,"CTYPE1","CRPIX1", crpix1) ;
109 cpl_propertylist_set_comment(*plist, "CRPIX1","Reference pixel in RA" ) ;
110
111 cpl_propertylist_erase_regexp(*plist, "^CRVAL1",0);
112 cpl_propertylist_insert_after_double(*plist, "CRPIX1", "CRVAL1", crval1 ) ;
113 cpl_propertylist_set_comment(*plist, "CRVAL1","Reference RA" ) ;
114
115 cpl_propertylist_erase_regexp(*plist, "^CDELT1",0);
116 cpl_propertylist_insert_after_double(*plist,"CRVAL1","CDELT1",cdelt1 ) ;
117 cpl_propertylist_set_comment(*plist, "CDELT1","pixel scale" ) ;
118
119 cpl_propertylist_erase_regexp(*plist, "^CUNIT1",0);
120 cpl_propertylist_insert_after_string(*plist, "CDELT1", "CUNIT1", "deg" ) ;
121 cpl_propertylist_set_comment(*plist, "CUNIT1","RA-UNIT" ) ;
122
123 return;
124}
125
126
135static void
136xsh_plist_set_coord2(cpl_propertylist** plist,
137 const double crpix2,
138 const double crval2,
139 const double cdelt2)
140{
141 cpl_propertylist_erase_regexp(*plist, "^CTYPE2",0);
142 cpl_propertylist_insert_after_string(*plist,"CUNIT1","CTYPE2","DEC--TAN");
143 cpl_propertylist_set_comment(*plist, "CTYPE2", "Projected Declination") ;
144
145 cpl_propertylist_erase_regexp(*plist, "^CRPIX2",0);
146 cpl_propertylist_insert_after_double(*plist,"CTYPE2","CRPIX2",crpix2 ) ;
147 cpl_propertylist_set_comment(*plist, "CRPIX2", "Reference pixel in DEC") ;
148
149 cpl_propertylist_erase_regexp(*plist,"^CRVAL2",0);
150 cpl_propertylist_insert_after_double(*plist,"CRPIX2","CRVAL2",crval2) ;
151 cpl_propertylist_set_comment(*plist,"CRVAL2","Reference DEC") ;
152
153 cpl_propertylist_erase_regexp(*plist,"^CDELT2",0);
154 cpl_propertylist_insert_after_double(*plist,"CRVAL2","CDELT2",cdelt2 ) ;
155 cpl_propertylist_set_comment(*plist,"CDELT2","pixel scale") ;
156
157 cpl_propertylist_erase_regexp(*plist,"^CUNIT2",0);
158 cpl_propertylist_insert_after_string(*plist,"CDELT2","CUNIT2", "deg" ) ;
159 cpl_propertylist_set_comment(*plist,"CUNIT2","DEC-UNIT") ;
160
161
162}
163
164
165
174static void
175xsh_plist_set_coord3(cpl_propertylist** plist,
176 const int crpix3,
177 const double crval3,
178 const double cdelt3)
179{
180 cpl_propertylist_erase_regexp(*plist, "^CTYPE3",0);
181 cpl_propertylist_insert_after_string(*plist,"EXPTIME", "CTYPE3", "WAVE" ) ;
182 cpl_propertylist_set_comment(*plist,"CTYPE3","wavelength axis in nm") ;
183
184
185 cpl_propertylist_erase_regexp(*plist, "^CRPIX3",0);
186 cpl_propertylist_insert_after_double(*plist, "CTYPE3", "CRPIX3", (double)crpix3 ) ;
187 cpl_propertylist_set_comment(*plist, "CRPIX3", "Reference pixel in z") ;
188
189 cpl_propertylist_erase_regexp(*plist, "^CRVAL3",0);
190 cpl_propertylist_insert_after_double(*plist,"CRPIX3", "CRVAL3", crval3) ;
191 cpl_propertylist_set_comment(*plist, "CRVAL3", "central wavelength") ;
192
193 cpl_propertylist_erase_regexp(*plist, "^CDELT3",0);
194
195 cpl_propertylist_insert_after_double(*plist,"CRVAL3","CDELT3",cdelt3) ;
196 cpl_propertylist_set_comment(*plist, "CDELT3", "nm per pixel") ;
197
198 cpl_propertylist_erase_regexp(*plist, "^CUNIT3",0);
199 cpl_propertylist_insert_after_string(*plist,"CDELT3", "CUNIT3", "nm" ) ;
200 cpl_propertylist_set_comment(*plist, "CUNIT3", "spectral unit" ) ;
201
202 cpl_propertylist_erase_regexp(*plist, "^SPECSYS",0);
203 cpl_propertylist_insert_after_string(*plist,"CUNIT3", "SPECSYS", "TOPOCENT" ) ;
204 cpl_propertylist_set_comment(*plist, "SPECSYS", "Coordinate reference frame" ) ;
205
206}
207
208
219static void
220xsh_plist_set_cd_matrix2(cpl_propertylist** plist,
221 const double cd1_1,
222 const double cd1_2,
223 const double cd2_1,
224 const double cd2_2)
225{
226
227 check(cpl_propertylist_erase_regexp(*plist, "^CD1_1",0));
228 check(cpl_propertylist_insert_after_double(*plist,"EXPTIME",
229 "CD1_1", cd1_1 )) ;
230 check(cpl_propertylist_set_comment(*plist, "CD1_1",
231 "CD rotation matrix" )) ;
232
233 check(cpl_propertylist_erase_regexp(*plist, "^CD1_2",0));
234 check(cpl_propertylist_insert_after_double(*plist, "CD1_1",
235 "CD1_2", cd1_2 )) ;
236 check(cpl_propertylist_set_comment(*plist, "CD1_2",
237 "CD rotation matrix" )) ;
238
239 check(cpl_propertylist_erase_regexp(*plist, "^CD2_1",0));
240 check(cpl_propertylist_insert_after_double(*plist, "CD1_2",
241 "CD2_1", cd2_1 )) ;
242 check(cpl_propertylist_set_comment(*plist, "CD2_1",
243 "CD rotation matrix" )) ;
244
245 check(cpl_propertylist_erase_regexp(*plist, "^CD2_2",0));
246 check(cpl_propertylist_insert_after_double(*plist, "CD2_1",
247 "CD2_2", cd2_2 )) ;
248 check(cpl_propertylist_set_comment(*plist, "CD2_2",
249 "CD rotation matrix" )) ;
250
251 cleanup:
252 return;
253
254
255}
256
257
270static void
271xsh_plist_set_cd_matrix3(cpl_propertylist** plist,
272 const double cd1_3,
273 const double cd2_3,
274 const double cd3_1,
275 const double cd3_2,
276 const double cd3_3)
277{
278
279
280 check(cpl_propertylist_erase_regexp(*plist, "^CD1_3",0));
281 check(cpl_propertylist_insert_after_double(*plist,"EXPTIME",
282 "CD1_3", cd1_3 )) ;
283 check(cpl_propertylist_set_comment(*plist, "CD1_3",
284 "CD rotation matrix" )) ;
285
286
287 check(cpl_propertylist_erase_regexp(*plist, "^CD2_3",0));
288 check(cpl_propertylist_insert_after_double(*plist,"CD1_3",
289 "CD2_3", cd2_3 )) ;
290 check(cpl_propertylist_set_comment(*plist, "CD2_3",
291 "CD rotation matrix" )) ;
292
293
294
295 check(cpl_propertylist_erase_regexp(*plist, "^CD3_1",0));
296 check(cpl_propertylist_insert_after_double(*plist,"CD2_3",
297 "CD3_1", cd3_1 )) ;
298 check(cpl_propertylist_set_comment(*plist, "CD3_1",
299 "CD rotation matrix" )) ;
300
301 check(cpl_propertylist_erase_regexp(*plist, "^CD3_2",0));
302 check(cpl_propertylist_insert_after_double(*plist, "CD3_1",
303 "CD3_2", cd3_2 )) ;
304 check(cpl_propertylist_set_comment(*plist, "CD3_2",
305 "CD rotation matrix" )) ;
306
307 check(cpl_propertylist_erase_regexp(*plist, "^CD3_3",0));
308 check(cpl_propertylist_insert_after_double(*plist, "CD3_2",
309 "CD3_3", cd3_3 )) ;
310 check(cpl_propertylist_set_comment(*plist, "CD3_3",
311 "CD rotation matrix" )) ;
312
313 cleanup:
314 return;
315
316
317}
318
330cpl_error_code
331xsh_cube_set_wcs(cpl_propertylist * plist,
332 float cenLambda,
333 float disp_x,
334 float disp_y,
335 float disp_z,
336 float center_x,
337 float center_y,
338 int center_z)
339{
340
341
342 double ra ;
343 double dec ;
344 double angle ;
345 float radangle ;
346 double cd1_1, cd1_2, cd2_1, cd2_2 ;
347 int sign_swap = -1;
348
349
350 double cdelt1=disp_x;
351 double cdelt2=disp_y;
352 double cdelt3=disp_z;
353
354 double crpix1=center_x;
355 double crpix2=center_y;
356 int crpix3=center_z;
357
358 double crval1=0;
359 double crval2=0;
360 double crval3=cenLambda;
361
362 ra = xsh_pfits_get_ra(plist) ;
363 dec = xsh_pfits_get_dec(plist) ;
364
365 //get better coordinate values
368
369 //xsh_msg("ra=%f",ra);
370 //xsh_msg("dec=%f",dec);
371 ra=xsh_hms2deg(ra);
372 dec=xsh_sess2deg(dec);
373 //xsh_msg("ra=%f",ra);
374 //xsh_msg("dec=%f",dec);
375
376 crval1=ra;
377 crval2=dec;
378
379 angle = xsh_pfits_get_posangle(plist) ;
380 /* in PUPIL data there is not posangle info: we reset the error */
381 if(cpl_error_get_code() != CPL_ERROR_NONE) {
382 cpl_error_reset();
383 }
384 cdelt1=sign_swap*disp_x / 3600.;
385 cdelt2= +disp_y / 3600.;
386
387
388 radangle = angle * PI_NUMB / 180. ;
389 cd1_1 = +cdelt1*cos(radangle);
390 cd1_2 = -cdelt2*sin(radangle);
391 cd2_1 = +cdelt1*sin(radangle);
392 cd2_2 = +cdelt2*cos(radangle);
393
394
395 xsh_plist_set_coord1(&plist,crpix1,crval1,cdelt1);
396 xsh_plist_set_coord2(&plist,crpix2,crval2,cdelt2);
397
398 xsh_plist_set_coord3(&plist,crpix3,crval3,cdelt3);
399 xsh_plist_set_cd_matrix2(&plist,cd1_1,cd1_2,cd2_1,cd2_2);
400 xsh_plist_set_cd_matrix3(&plist,0,0,0,0,disp_z);
401
402 return cpl_error_get_code();
403}
404
405
406
407
408
409
410
413/*---------------------------------------------------------------------------*/
424/*---------------------------------------------------------------------------*/
425
426
427void
428xsh_edge_check(const int px,const int nx,const int rad_x,
429 int* llx,int* urx)
430{
431
432
433 *llx=px-rad_x;
434 *urx=px+rad_x;
435 *llx=(*llx>1) ? *llx : 1;
436 *urx=(*urx<nx) ? *urx : nx;
437 return;
438
439}
440
441
442/*---------------------------------------------------------------------------*/
457/*---------------------------------------------------------------------------*/
458
459void
460xsh_convert_xy_to_ws(double x_centroid,
461 double* p_obj_cen,
462 double* p_slit,
463 double* p_wave,
464 const int lly,
465 const int nx,
466 const int row,
467 double* p_obj_cen_s,
468 double* p_obj_cen_w)
469{
470 int x_ceil=0;
471 int x_floor=0;
472 double s_ceil=0;
473 double s_floor=0;
474 double w_ceil=0;
475 double w_floor=0;
476
477 p_obj_cen[row]=x_centroid;
478 x_ceil=ceil(x_centroid);
479 x_floor=floor(x_centroid);
480
481 s_ceil =p_slit[lly*nx+x_ceil];
482 s_floor=p_slit[lly*nx+x_floor];
483
484 w_ceil =p_wave[lly*nx+x_ceil];
485 w_floor=p_wave[lly*nx+x_floor];
486
487 p_obj_cen_s[row]=xsh_interpol(x_centroid,x_floor,x_ceil,s_floor,s_ceil);
488 p_obj_cen_w[row]=xsh_interpol(x_centroid,x_floor,x_ceil,w_floor,w_ceil);
489
490 return;
491}
492
493
494/*---------------------------------------------------------------------------*/
501/*---------------------------------------------------------------------------*/
502
503cpl_error_code
505{
506 /* swap LOW/UPP in case of VIS */
507 cpl_table_duplicate_column(*tab,"OBJ_LOW_S_TMP",*tab,"OBJ_LOW_S");
508 cpl_table_duplicate_column(*tab,"OBJ_UPP_S_TMP",*tab,"OBJ_UPP_S");
509
510 cpl_table_erase_column(*tab,"OBJ_LOW_S");
511 cpl_table_erase_column(*tab,"OBJ_UPP_S");
512
513 cpl_table_duplicate_column(*tab,"OBJ_UPP_S",*tab,"OBJ_LOW_S_TMP");
514 cpl_table_duplicate_column(*tab,"OBJ_LOW_S",*tab,"OBJ_UPP_S_TMP");
515
516 cpl_table_erase_column(*tab,"OBJ_LOW_S_TMP");
517 cpl_table_erase_column(*tab,"OBJ_UPP_S_TMP");
518
519 cpl_table_duplicate_column(*tab,"OBJ_LOW_W_TMP",*tab,"OBJ_LOW_W");
520 cpl_table_duplicate_column(*tab,"OBJ_UPP_W_TMP",*tab,"OBJ_UPP_W");
521
522 cpl_table_erase_column(*tab,"OBJ_LOW_W");
523 cpl_table_erase_column(*tab,"OBJ_UPP_W");
524
525 cpl_table_duplicate_column(*tab,"OBJ_UPP_W",*tab,"OBJ_LOW_W_TMP");
526 cpl_table_duplicate_column(*tab,"OBJ_LOW_W",*tab,"OBJ_UPP_W_TMP");
527
528 cpl_table_erase_column(*tab,"OBJ_LOW_W_TMP");
529 cpl_table_erase_column(*tab,"OBJ_UPP_W_TMP");
530
531
532 return cpl_error_get_code();
533}
534
535
536
537/*---------------------------------------------------------------------------*/
544/*---------------------------------------------------------------------------*/
545
546cpl_table*
547xsh_table_edge_prepare(const char* name)
548{
549 cpl_table* tab=NULL;
550 int nrow=0;
551
552 check(tab=cpl_table_load(name,2,0));
553
554 nrow=cpl_table_get_nrow(tab);
555
556 cpl_table_new_column(tab,"OBJ_LOW_X",CPL_TYPE_DOUBLE);
557 cpl_table_new_column(tab,"OBJ_CEN_X",CPL_TYPE_DOUBLE);
558 cpl_table_new_column(tab,"OBJ_UPP_X",CPL_TYPE_DOUBLE);
559
560 cpl_table_fill_column_window(tab,"OBJ_LOW_X",0,nrow,-1);
561 cpl_table_fill_column_window(tab,"OBJ_CEN_X",0,nrow,-1);
562 cpl_table_fill_column_window(tab,"OBJ_UPP_X",0,nrow,-1);
563
564 cpl_table_new_column(tab,"OBJ_LOW_S",CPL_TYPE_DOUBLE);
565 cpl_table_new_column(tab,"OBJ_LOW_W",CPL_TYPE_DOUBLE);
566 cpl_table_new_column(tab,"OBJ_CEN_S",CPL_TYPE_DOUBLE);
567 cpl_table_new_column(tab,"OBJ_CEN_W",CPL_TYPE_DOUBLE);
568 cpl_table_new_column(tab,"OBJ_UPP_S",CPL_TYPE_DOUBLE);
569 cpl_table_new_column(tab,"OBJ_UPP_W",CPL_TYPE_DOUBLE);
570
571 cpl_table_fill_column_window(tab,"OBJ_LOW_S",0,nrow,-1);
572 cpl_table_fill_column_window(tab,"OBJ_LOW_W",0,nrow,-1);
573 cpl_table_fill_column_window(tab,"OBJ_CEN_S",0,nrow,-1);
574 cpl_table_fill_column_window(tab,"OBJ_CEN_W",0,nrow,-1);
575 cpl_table_fill_column_window(tab,"OBJ_UPP_S",0,nrow,-1);
576 cpl_table_fill_column_window(tab,"OBJ_UPP_W",0,nrow,-1);
577
578 cleanup:
579 return tab;
580
581}
582
583/* Not used
584static cpl_error_code
585xsh_tab_clean_traces(cpl_table** tab,
586 const double fit_wmin,
587 const double fit_wmax,
588 const char* wave_col_name,
589 const char* trace_col_name,
590 const char* fit_col_name,
591 const char* res_col_name,
592 const char* flag_col_name,
593 const int niter,
594 const double kappa)
595{
596
597 double* pdata=NULL;
598 double* pwave=NULL;
599 double* pres=NULL;
600 double* pfit=NULL;
601 int* pflag=NULL;
602
603 cpl_polynomial* p_fit=NULL;
604 cpl_vector* vdata=NULL;
605 cpl_vector* vwave=NULL;
606 int fit_pows[2];
607 double fit_off=0;
608 double fit_slope=0;
609 int row=0;
610 int nrow=0;
611 double mse_fit=0;
612 double res_stdev;
613 double res_max=1;
614 int i=0;
615 cpl_table* ext=NULL;
616 int next=0;
617
618 check(nrow=cpl_table_get_nrow(*tab));
619 cpl_table_and_selected_double(*tab,wave_col_name,CPL_GREATER_THAN,fit_wmin);
620 cpl_table_and_selected_double(*tab,wave_col_name,CPL_LESS_THAN,fit_wmax);
621 ext=cpl_table_extract_selected(*tab);
622 check(next=cpl_table_get_nrow(ext));
623
624
625 pdata=cpl_table_get_data_double(ext,trace_col_name);
626 pwave=cpl_table_get_data_double(ext,wave_col_name);
627 vdata=cpl_vector_wrap(next,pdata);
628 vwave=cpl_vector_wrap(next,pwave);
629 p_fit=xsh_polynomial_fit_1d_create(vwave,vdata,1,&mse_fit);
630
631 cpl_polynomial_dump(p_fit,stdout);
632 fit_pows[0]=0;
633 check(fit_off=cpl_polynomial_get_coeff(p_fit,fit_pows));
634 fit_pows[0]=1;
635 check(fit_slope=cpl_polynomial_get_coeff(p_fit,fit_pows));
636
637
638
639 cpl_table_new_column(*tab,fit_col_name,CPL_TYPE_DOUBLE);
640 cpl_table_fill_column_window(*tab,fit_col_name,0,nrow,-1);
641 pfit=cpl_table_get_data_double(*tab,fit_col_name);
642
643 cpl_table_new_column(*tab,res_col_name,CPL_TYPE_DOUBLE);
644 cpl_table_fill_column_window(*tab,res_col_name,0,nrow,-1);
645 pres=cpl_table_get_data_double(*tab,res_col_name);
646
647 cpl_table_new_column(*tab,flag_col_name,CPL_TYPE_INT);
648 cpl_table_fill_column_window(*tab,flag_col_name,0,nrow,-1);
649 pflag=cpl_table_get_data_int(*tab,flag_col_name);
650
651
652
653 pdata=cpl_table_get_data_double(*tab,trace_col_name);
654 pwave=cpl_table_get_data_double(*tab,wave_col_name);
655
656 for(row=0;row<nrow;row++) {
657 pfit[row]=fit_off+pwave[row]*fit_slope;
658 pres[row]=pfit[row]-pdata[row];
659 }
660
661
662 // kappa-sigma-clip upper trace
663 for(i=0;i<niter;i++) {
664 res_stdev=cpl_table_get_column_stdev(*tab,res_col_name);
665 xsh_msg("res_stdev=%g",res_stdev);
666 for(row=0;row<nrow;row++) {
667 if((fabs(pres[row])>kappa*res_stdev) ||
668 (fabs(pres[row])>res_max)) {
669 //xsh_msg("set invalid %d res=%g",row,pres[row]);
670 cpl_table_set_invalid(*tab,res_col_name,row);
671 } else {
672 pflag[row]=1;
673 }
674 }
675 }
676
677 cpl_table_erase_invalid(*tab);
678
679
680 cleanup:
681 return cpl_error_get_code();
682
683}
684 */
685
686
687
688
689/*---------------------------------------------------------------------------*/
698/*---------------------------------------------------------------------------*/
699
700cpl_error_code
701xsh_ifu_trace_object_calibrate(const char* ifu_object_ff_name,
702 const char* order_tab_edges_ifu_name,
703 const char* slit_map_name,
704 const char* wave_map_name)
705{
706 cpl_image* ifu_object_ff_ima=NULL;
707 cpl_image* wave_map_ima=NULL;
708 cpl_image* slit_map_ima=NULL;
709 cpl_image* ratio_ima=NULL;
710 cpl_table* tab=NULL;
711
712 //int sx=0;
713 //int sy=0;
714
715 int nrow=0;
716
717 double* p_edge_lo_x=NULL;
718 double* p_edge_up_x=NULL;
719
720 double* p_slice_lo_x=NULL;
721 double* p_slice_up_x=NULL;
722
723 double* p_center_y=NULL;
724
725 double* p_obj_low=NULL;
726 double* p_obj_cen=NULL;
727 double* p_obj_upp=NULL;
728
729 double* p_obj_low_s=NULL;
730 double* p_obj_low_w=NULL;
731 double* p_obj_cen_s=NULL;
732 double* p_obj_cen_w=NULL;
733 double* p_obj_upp_s=NULL;
734 double* p_obj_upp_w=NULL;
735
736 double* p_slit=NULL;
737 double* p_wave=NULL;
738
739 double x_centroid=0;
740
741 int llx=0;
742 int lly=0;
743 int urx=0;
744 int ury=0;
745
746 int row=0;
747
748 //double range1_w_min=0;
749 //double range1_w_max=0;
750 //double range2_w_min=0;
751 //double range2_w_max=0;
752
753 //double range_ks_w_min=0;
754 //double range_ks_w_max=0;
755
756 int nx=0;
757 int ny=0;
758 cpl_propertylist* plist=NULL;
760 const char* pcatg=NULL;
761 char tag[25];
762 char name[256];
763 cpl_size px=0;
764 cpl_size py=0;
765
766 int rad_x=10;
767 /*
768 int niter=2;
769 int kappa=2.0;
770 */
771 int biny=1;
772
773 typedef enum {centroid, gaussian} xsh_fit_method;
774 xsh_fit_method fit_method = centroid;
775
776
777 check(ifu_object_ff_ima=cpl_image_load(ifu_object_ff_name,CPL_TYPE_DOUBLE,0,0));
778 plist=cpl_propertylist_load(ifu_object_ff_name,0);
779 pcatg=xsh_pfits_get_pcatg(plist);
780
781 xsh_msg("pcatg=%s",pcatg);
782
783 if(strstr(pcatg,"UVB") != NULL) {
784 arm=XSH_ARM_UVB;
785 //range1_w_min=300;
786 //range1_w_max=500;
787 //range_ks_w_min=380;
788 //range_ks_w_max=550;
790
791
792 } else if(strstr(pcatg,"VIS") != NULL) {
793 arm=XSH_ARM_VIS;
794 //range1_w_min=600;
795 //range1_w_max=950;
796
797 //range_ks_w_min=600;
798 //range_ks_w_max=900;
800
801 } else if(strstr(pcatg,"NIR") != NULL) {
802 arm=XSH_ARM_NIR;
803 //range1_w_min=1050;
804 //range1_w_max=1250;
805 //range2_w_min=2000;
806 //range2_w_max=2200;
807
808 //range_ks_w_min=1100;
809 //range_ks_w_max=2200;
810
811
812 }
813
814 check(slit_map_ima=cpl_image_load(slit_map_name,CPL_TYPE_DOUBLE,0,0));
815 check(wave_map_ima=cpl_image_load(wave_map_name,CPL_TYPE_DOUBLE,0,0));
816 //sx=cpl_image_get_size_x(ifu_object_ff_ima);
817 //sy=cpl_image_get_size_y(ifu_object_ff_ima);
818
819 check(tab=xsh_table_edge_prepare(order_tab_edges_ifu_name));
820 nrow=cpl_table_get_nrow(tab);
821
822 p_edge_lo_x=cpl_table_get_data_double(tab,"EDG_LO_X");
823 p_slice_lo_x=cpl_table_get_data_double(tab,"SLIC_LO_X");
824 p_obj_low=cpl_table_get_data_double(tab,"OBJ_LOW_X");
825 p_obj_low_s=cpl_table_get_data_double(tab,"OBJ_LOW_S");
826 p_obj_low_w=cpl_table_get_data_double(tab,"OBJ_LOW_W");
827
828 p_edge_up_x=cpl_table_get_data_double(tab,"EDG_UP_X");
829 p_slice_up_x=cpl_table_get_data_double(tab,"SLIC_UP_X");
830 p_obj_upp=cpl_table_get_data_double(tab,"OBJ_UPP_X");
831 p_obj_upp_s=cpl_table_get_data_double(tab,"OBJ_UPP_S");
832 p_obj_upp_w=cpl_table_get_data_double(tab,"OBJ_UPP_W");
833
834 p_center_y=cpl_table_get_data_double(tab,"CENTER_Y");
835 p_obj_cen=cpl_table_get_data_double(tab,"OBJ_CEN_X");
836 p_obj_cen_s=cpl_table_get_data_double(tab,"OBJ_CEN_S");
837 p_obj_cen_w=cpl_table_get_data_double(tab,"OBJ_CEN_W");
838
839 p_slit=cpl_image_get_data_double(slit_map_ima);
840 p_wave=cpl_image_get_data_double(wave_map_ima);
841 nx=cpl_image_get_size_x(wave_map_ima);
842 ny=cpl_image_get_size_y(wave_map_ima);
843
844 for(row=0;row<nrow;row++) {
845
846 lly=(int)p_center_y[row]/biny;
847 ury=(int)p_center_y[row]/biny;
848
849 //xsh_msg("lly=%d ury=%d nx=%d ny=%d",lly,ury,nx,ny);
850 if((llx<nx) && (ury<ny) ) {
851
852 /* CEN IFU slice obj position */
853 llx=floor(p_slice_lo_x[row]);
854 urx=ceil(p_slice_up_x[row]);
855
856 check(cpl_image_get_maxpos_window(ifu_object_ff_ima,
857 llx,lly,urx,ury,&px,&py));
858
859 xsh_edge_check(px,nx,rad_x,&llx,&urx);
860 if(fit_method==centroid) {
861 x_centroid=cpl_image_get_centroid_x_window(ifu_object_ff_ima,
862 llx,lly,urx,ury);
863
864 } else {
865 check(x_centroid=xsh_image_fit_gaussian_max_pos_x_window(ifu_object_ff_ima,
866 llx,urx,lly));
867 }
868
869
870 p_obj_cen[row]=x_centroid;
871 /* now measure position in arcsecs and wavelength */
872 if(urx<nx) {
873 check(xsh_convert_xy_to_ws(x_centroid,p_obj_cen,p_slit,p_wave,
874 lly,nx,row,p_obj_cen_s,p_obj_cen_w));
875 }
876
877 /* LOW IFU slice obj position */
878 llx=floor(p_edge_lo_x[row]);
879 urx=ceil(p_slice_lo_x[row]);
880 cpl_image_get_maxpos_window(ifu_object_ff_ima,
881 llx,lly,urx,ury,&px,&py);
882
883 xsh_edge_check(px,nx,rad_x,&llx,&urx);
884 if(fit_method==centroid) {
885 check(x_centroid=cpl_image_get_centroid_x_window(ifu_object_ff_ima,
886 llx,lly,urx,ury));
887 } else {
888 x_centroid=xsh_image_fit_gaussian_max_pos_x_window(ifu_object_ff_ima,
889 llx,urx,lly);
890 }
891
892
893 /* now measure position in arcsecs and wavelength */
894 if(urx<nx) {
895 check(xsh_convert_xy_to_ws(x_centroid,p_obj_low,p_slit,p_wave,
896 lly,nx,row,p_obj_low_s,p_obj_low_w));
897 }
898
899 /* UPP IFU slice obj position */
900 llx=floor(p_slice_up_x[row]);
901 urx=ceil(p_edge_up_x[row]);
902 cpl_image_get_maxpos_window(ifu_object_ff_ima,
903 llx,lly,urx,ury,&px,&py);
904
905 xsh_edge_check(px,nx,rad_x,&llx,&urx);
906 if(fit_method==centroid) {
907 x_centroid=cpl_image_get_centroid_x_window(ifu_object_ff_ima,
908 llx,lly,urx,ury);
909 } else {
910 x_centroid=xsh_image_fit_gaussian_max_pos_x_window(ifu_object_ff_ima,
911 llx,urx,lly);
912 }
913
914 p_obj_upp[row]=x_centroid;
915 /* now measure position in arcsecs and wavelength */
916 if(urx<nx) {
917 check(xsh_convert_xy_to_ws(x_centroid,p_obj_upp,p_slit,p_wave,
918 lly,nx,row,p_obj_upp_s,p_obj_upp_w));
919 }
920 }
921 }
922
923 if(arm==XSH_ARM_VIS) {
925 }
926 cpl_table_duplicate_column(tab,"OBJ_CEN_PLUS_UPP_S",tab,"OBJ_CEN_S");
927
928 cpl_table_duplicate_column(tab,"OBJ_CEN_PLUS_LOW_S",tab,"OBJ_CEN_S");
929
930 cpl_table_add_columns(tab,"OBJ_CEN_PLUS_UPP_S","OBJ_UPP_S");
931 cpl_table_add_columns(tab,"OBJ_CEN_PLUS_LOW_S","OBJ_LOW_S");
932
933
934
935 /* we now clean results:
936 -fit linear to OBJ_CEN_PLUS_UPP/LOW_S
937 -compute residuals fit-value
938 -compute stdev residuals
939 -apply kappa-sigma clip (2 times) flagging as null clipped values
940 */
941
942 /*
943 check(xsh_tab_clean_traces(&tab,range_ks_w_min,range_ks_w_max,
944 "OBJ_CEN_W","OBJ_CEN_PLUS_UPP_S",
945 "CEN_UPP_S_FIT","CEN_UPP_S_RES",
946 "CEN_UPP_S_FLAG",niter,kappa));
947
948 check(xsh_tab_clean_traces(&tab,range_ks_w_min,range_ks_w_max,
949 "OBJ_CEN_W","OBJ_CEN_PLUS_LOW_S",
950 "CEN_LOW_S_FIT","CEN_LOW_S_RES",
951 "CEN_LOW_S_FLAG",niter,kappa));
952
953 */
954 /* save result */
955 sprintf(tag,"TRACE_OBJ_%s",xsh_arm_tostring(arm));
956 sprintf(name,"%s.fits",tag);
957 xsh_pfits_set_pcatg(plist,tag);
958 check(cpl_table_save(tab,plist,NULL,name,CPL_IO_DEFAULT));
959
960
961
962 cleanup:
963 xsh_free_image(&ifu_object_ff_ima);
964 xsh_free_image(&slit_map_ima);
965 xsh_free_image(&wave_map_ima);
966 xsh_free_image(&ratio_ima);
967 xsh_free_propertylist(&plist);
968 return cpl_error_get_code();
969
970}
971
972
973/*---------------------------------------------------------------------------*/
981/*---------------------------------------------------------------------------*/
982
983cpl_error_code
985 cpl_frame* sci_frame,
987{
988
989 cpl_propertylist* afc_plist=NULL;
990 cpl_propertylist* sci_plist=NULL;
991 const char* name=NULL;
992 const char* sci_obs_targ_name=NULL;
993 const char* afc_obs_targ_name=NULL;
994 const char* afc_slit_value=NULL;
995 //const char* sci_slit_value=NULL;
996 int sci_obs_id=0;
997 int afc_obs_id=0;
998
999 check(name=cpl_frame_get_filename(model_config_frame));
1000 check(afc_plist=cpl_propertylist_load(name,0));
1001
1002 check(name=cpl_frame_get_filename(sci_frame));
1003 check(sci_plist=cpl_propertylist_load(name,0));
1004
1005 check(afc_slit_value=xsh_pfits_get_slit_value(afc_plist,instrument ));
1006
1007 //check(sci_slit_value=xsh_pfits_get_slit_value(sci_plist,instrument ));
1008
1009 if(strstr(afc_slit_value,"Pin_0.5") == NULL) {
1010 xsh_msg_error("You have used uncorrect AFC corrected model cfg frame");
1011 xsh_msg_error("IFU AFC corrected model CFG must have INS.OPTIi.NAME='Pin_0.5'");
1012 cpl_error_set(cpl_func,CPL_ERROR_ILLEGAL_INPUT);
1013 }
1014
1015
1016 check(sci_obs_targ_name=xsh_pfits_get_obs_targ_name(sci_plist));
1017 check(afc_obs_targ_name=xsh_pfits_get_obs_targ_name(afc_plist));
1018 if(strcmp(sci_obs_targ_name,afc_obs_targ_name) != 0) {
1019 xsh_msg_error("Improper AFC corrected model cfg frame to reduce sci frame");
1020 xsh_msg_error("Their OBS.TARG.NAME values must match");
1021 cpl_error_set(cpl_func,CPL_ERROR_ILLEGAL_INPUT);
1022
1023 }
1024
1025 check(sci_obs_id=xsh_pfits_get_obs_id(sci_plist));
1026 check(afc_obs_id=xsh_pfits_get_obs_id(afc_plist));
1027
1028 if(sci_obs_id != afc_obs_id) {
1029 xsh_msg_error("Improper AFC corrected model cfg frame to reduce sci frame");
1030 xsh_msg_error("Their OBS.ID values must match");
1031 cpl_error_set(cpl_func,CPL_ERROR_ILLEGAL_INPUT);
1032 }
1033
1034 cleanup:
1035
1036 xsh_free_propertylist(&sci_plist);
1037 xsh_free_propertylist(&afc_plist);
1038
1039 return cpl_error_get_code();
1040}
1041
1042
1043/*---------------------------------------------------------------------------*/
1050/*---------------------------------------------------------------------------*/
1051
1052cpl_error_code
1053xsh_frame_check_model_cfg_is_afc_corrected(cpl_frame* model_config_frame){
1054
1055 cpl_propertylist* plist=NULL;
1056 const char* name=NULL;
1057 const char* raw1_catg=NULL;
1058
1059 check(name=cpl_frame_get_filename(model_config_frame));
1060 check(plist=cpl_propertylist_load(name,0));
1061 check(raw1_catg=xsh_pfits_get_raw1catg(plist));
1062 if(strstr(raw1_catg,"AFC_ATT") == NULL) {
1063 xsh_msg_error("model cfg frame seems not to be AFC corrected");
1064 xsh_msg_error("Their PRO.REC1.RAW1..NAME values must contain AFC_ATT");
1065 cpl_error_set(cpl_func,CPL_ERROR_ILLEGAL_INPUT);
1066 }
1067
1068 cleanup:
1069
1070 xsh_free_propertylist(&plist);
1071
1072 return cpl_error_get_code();
1073
1074}
1075
1076static cpl_error_code
1077xsh_add_correct_coeff(cpl_table** table, cpl_propertylist* plist,
1078 const char* prefix, const int row,
1079 const char* col12, const char* col32)
1080{
1081 double coeff_t1=0;
1082 double coeff_t2=0;
1083 double coeff_t3=0;
1084 double diff_c12=0;
1085 double diff_c32=0;
1086 char key_name[40];
1087
1088
1089 sprintf(key_name,"%s_%s",prefix,"T1");
1090 check(coeff_t1=cpl_propertylist_get_double(plist,key_name));
1091
1092 sprintf(key_name,"%s_%s",prefix,"T2");
1093 coeff_t2=cpl_propertylist_get_double(plist,key_name);
1094
1095 sprintf(key_name,"%s_%s",prefix,"T3");
1096 coeff_t3=cpl_propertylist_get_double(plist,key_name);
1097 diff_c12=coeff_t1-coeff_t2;
1098 diff_c32=coeff_t3-coeff_t2;
1099
1100 cpl_table_set_double(*table,col12 ,row, diff_c12);
1101 cpl_table_set_double(*table,col32 ,row, diff_c32);
1102
1103 cleanup:
1104
1105 return cpl_error_get_code();
1106
1107}
1108
1109static cpl_frame*
1110xsh_crea_correct_coeff(cpl_frame* qc_trace_merged_frame,xsh_instrument* inst)
1111{
1112
1113 cpl_propertylist* plist=NULL;
1114 const char* fname=NULL;
1115 cpl_table* table=NULL;
1116 cpl_frame* result=NULL;
1117 char pname[256];
1118 char ptag[40];
1119
1120 fname=cpl_frame_get_filename(qc_trace_merged_frame);
1121 plist=cpl_propertylist_load(fname,0);
1122
1123 table=cpl_table_new(3);
1124 check(cpl_table_new_column(table,"DIFF_T12", CPL_TYPE_DOUBLE));
1125 check(cpl_table_new_column(table,"DIFF_T32", CPL_TYPE_DOUBLE));
1126
1127
1129 "DIFF_T12","DIFF_T32"));
1130
1132 "DIFF_T12","DIFF_T32"));
1133
1135 "DIFF_T12","DIFF_T32"));
1136
1137 sprintf(ptag,"IFU_CFG_COR_%s",xsh_instrument_arm_tostring(inst));
1138 sprintf(pname,"%s.fits",ptag);
1139 xsh_msg("tag=%s name=%s",ptag,pname);
1140
1141 check(cpl_table_save(table,plist,NULL,pname,CPL_IO_DEFAULT));
1142 result=xsh_frame_product(pname,ptag,CPL_FRAME_TYPE_TABLE,
1143 CPL_FRAME_GROUP_PRODUCT,CPL_FRAME_LEVEL_FINAL);
1144
1145
1146 cleanup:
1147
1148 xsh_free_propertylist(&plist);
1149 xsh_free_table(&table);
1150
1151 return result;
1152}
1153
1154static cpl_frame*
1155xsh_frame_build_sky_area(cpl_frame* slitmap_frame, const char* tag) {
1156
1157 char name_o[256];
1158 cpl_frame* result = NULL;
1159
1160 cpl_image* ima_slit = NULL;
1161 cpl_image* ima_sky = NULL;
1162 float* pslit = NULL;
1163 float* psky = NULL;
1164
1165 const char* name = NULL;
1166 cpl_propertylist* plist = NULL;
1167 int sx = 0;
1168 int sy = 0;
1169 int i = 0;
1170 int j = 0;
1171
1172 check(name = cpl_frame_get_filename(slitmap_frame));
1173 cpl_frame_dump(slitmap_frame, stdout);
1174 check(ima_slit = cpl_image_load(name, CPL_TYPE_FLOAT, 0, 0));
1175 check(plist = cpl_propertylist_load(name, 0));
1176 pslit = cpl_image_get_data_float(ima_slit);
1177
1178 sx = cpl_image_get_size_x(ima_slit);
1179 sy = cpl_image_get_size_y(ima_slit);
1180 ima_sky = cpl_image_new(sx, sy, CPL_TYPE_FLOAT);
1181 psky = cpl_image_get_data_float(ima_sky);
1182 for (j = 1; j < sy - 1; j++) {
1183 for (i = 1; i < sx - 1; i++) {
1184
1185 psky[j * sx + i] = 0.25 * (pslit[j * sx + i + 1] - pslit[j * sx + i - 1])
1186 * (pslit[(j + 1) * sx + i] - pslit[(j - 1) * sx + i]);
1187
1188 }
1189 }
1190
1191 sprintf(name_o, "%s.fits", tag);
1192 check(cpl_image_save(ima_sky, name_o, XSH_PRE_DATA_BPP, plist, CPL_IO_DEFAULT));
1193
1194 result = cpl_frame_duplicate(slitmap_frame);
1195 cpl_frame_set_filename(result, name_o);
1196 cpl_frame_set_tag(result, tag);
1197
1198 cleanup: xsh_free_propertylist(&plist);
1199 xsh_free_image(&ima_slit);
1200 xsh_free_image(&ima_sky);
1201 return result;
1202
1203}
1204
1205static cpl_frame*
1206xsh_frame_build_sky_map(cpl_frame* slitmap_frame,const double value,const char* tag)
1207{
1208
1209 char name_o[256];
1210 cpl_frame* result=NULL;
1211
1212 cpl_image* ima = NULL;
1213 const char* name = NULL;
1214 cpl_propertylist* plist = NULL;
1215
1216 check(name = cpl_frame_get_filename(slitmap_frame));
1217 check(ima = cpl_image_load(name, XSH_PRE_DATA_TYPE, 0, 0));
1218 check(plist = cpl_propertylist_load(name, 0));
1219
1220 check(cpl_image_add_scalar(ima,value));
1221 sprintf(name_o,"%s.fits",tag);
1222 check(cpl_image_save(ima, name_o, XSH_PRE_DATA_BPP, plist, CPL_IO_DEFAULT));
1223
1224 result = cpl_frame_duplicate(slitmap_frame);
1225 cpl_frame_set_filename(result, name_o);
1226 cpl_frame_set_tag(result, tag);
1227 xsh_add_temporary_file(name_o);
1228
1229 cleanup:
1230 xsh_free_propertylist(&plist);
1231 xsh_free_image(&ima);
1232
1233 return result;
1234
1235}
1236/*---------------------------------------------------------------------------*/
1246/*---------------------------------------------------------------------------*/
1247cpl_frame*
1248xsh_build_ifu_map(cpl_frame* div_frame,
1249 cpl_frame* wavemap_frame,
1250 cpl_frame* slitmap_frame,
1252{
1253
1254 cpl_frame* map = NULL;
1255 cpl_frame* ra_map = NULL;
1256 cpl_frame* dec_map = NULL;
1257 cpl_frame* sky_area = NULL;
1258
1259 char name_o[256];
1260 char tag_o[256];
1261
1262 cpl_image* ima = NULL;
1263 const char* name = NULL;
1264 cpl_propertylist* plist = NULL;
1265 double RA=0;
1266 double DEC=0;
1267
1268 name = cpl_frame_get_filename(div_frame);
1269 plist = cpl_propertylist_load(name, 0);
1270
1272 sprintf(name_o,"%s.fits",tag_o);
1273 /*
1274 map=cpl_frame_duplicate(div_frame);
1275 cpl_frame_set_filename(map,name_o);
1276 cpl_frame_set_tag(map,tag_o);
1277 */
1278 xsh_frame_image_save2ext(div_frame,name_o,0,0);
1279 xsh_frame_image_save2ext(wavemap_frame,name_o,0,1);
1280 xsh_frame_image_save2ext(slitmap_frame,name_o,0,2);
1281 xsh_frame_image_save2ext(slitmap_frame,name_o,1,3);
1282
1283 /*
1284 RA=xsh_pfits_get_ra(plist);
1285 DEC=xsh_pfits_get_dec(plist);
1286 get better coordinate values
1287 */
1290
1291 /* converts from hours-min-sec to degree */
1292 RA=xsh_hms2deg(RA);
1293 /* converts from sessagesimal angles to degree */
1294 DEC=xsh_sess2deg(DEC);
1295
1296 /* converts to arcsecs */
1297 RA*=3600;
1298 DEC*=3600;
1299
1300 /* prepare RA map */
1301 sprintf(tag_o,"%s_%s","IFU_MAP_SKY_RA",xsh_instrument_arm_tostring(instrument));
1302 check(ra_map=xsh_frame_build_sky_map(slitmap_frame,RA,tag_o));
1303 check(xsh_frame_image_save2ext(ra_map,name_o,0,4));
1304
1305 /* prepare DEC map */
1306 sprintf(tag_o,"%s_%s","IFU_MAP_SKY_DEC",xsh_instrument_arm_tostring(instrument));
1307 check(dec_map=xsh_frame_build_sky_map(slitmap_frame,DEC,tag_o));
1308 check(xsh_frame_image_save2ext(dec_map,name_o,0,5));
1309
1310 /* prepare DEC map */
1311 sprintf(tag_o,"%s_%s","IFU_MAP_SKY_AREA",xsh_instrument_arm_tostring(instrument));
1312 check(sky_area=xsh_frame_build_sky_area(slitmap_frame,tag_o));
1313 check(xsh_frame_image_save2ext(sky_area,name_o,0,6));
1314 xsh_add_temporary_file(name_o);
1315 /*
1316 xsh_frame_image_save2ext(slice_map,name_o,6);
1317 */
1318 //sprintf(tag_o,"%s_%s",XSH_IFU_MAP_SKY,xsh_instrument_arm_tostring(instrument));
1319
1320 map = xsh_frame_product(name_o, tag_o, CPL_FRAME_TYPE_IMAGE,
1321 CPL_FRAME_GROUP_PRODUCT, CPL_FRAME_LEVEL_FINAL);
1322
1323 cleanup:
1324 xsh_free_propertylist(&plist);
1325 xsh_free_image(&ima);
1326
1327 return map;
1328
1329}
1330
1331
1332
1333
1334static cpl_error_code
1336 cpl_imagelist* errs,
1337 cpl_imagelist* qual,
1338 cpl_propertylist* phead,
1339 const double wave_s,
1340 const double wave_e,
1341 const int index)
1342{
1343
1344 int zstart=0;
1345 int zend=0;
1346
1347 double flux=0;
1348 double err=0;
1349 double sn=0;
1350 double ws=wave_s;
1351 double we=wave_e;
1352
1353 char comment[40];
1354 char qc_key[20];
1355
1356 double cdelt3=0;
1357 double crval3=0;
1358 int naxis3=0;
1359 int j=0;
1360 cpl_imagelist* data_swap=NULL;
1361 cpl_imagelist* errs_swap=NULL;
1362 //cpl_imagelist* qual_swap=NULL;
1363 cpl_image* data_mean_ima=NULL;
1364 cpl_image* errs_mean_ima=NULL;
1365 //cpl_image* qual_mean_ima=NULL;
1366
1367 /* for the cube case we need to get size, wave start, wave step */
1368 check(crval3=xsh_pfits_get_crval3(phead));
1369 check(cdelt3=xsh_pfits_get_cdelt3(phead));
1370 //check(naxis3=xsh_pfits_get_naxis3(phead));
1371 naxis3=(int)cpl_imagelist_get_size(data);
1372
1373 /* special case: index=0 is used to use full range for statistics */
1374 if (index != 0) {
1375 zstart=((ws - crval3)/cdelt3+0.5);
1376 zend= ((we - crval3)/cdelt3-0.5);
1377 } else {
1378 zstart=1;
1379 zend=naxis3;
1380 ws=crval3;
1381 we=crval3+(naxis3-1.)*cdelt3;
1382 }
1383
1384
1385 /* make consistency check */
1386 zstart=(zstart < naxis3) ? zstart: naxis3;
1387 zend=(zend < naxis3) ? zend: naxis3;
1388 zstart=(zstart < zend) ? zstart: zend-1;
1389
1390 /* swap Y with Z axis to be able to collapse over Y (spatial direction) */
1391 data_swap=cpl_imagelist_swap_axis_create(data,CPL_SWAP_AXIS_YZ);
1392 errs_swap=cpl_imagelist_swap_axis_create(errs,CPL_SWAP_AXIS_YZ);
1393
1394 /* average extraction of IFU spectrum (average over original Y direction) */
1395 data_mean_ima=cpl_imagelist_collapse_create(data_swap);
1396 errs_mean_ima=cpl_imagelist_collapse_create(errs_swap);
1397 xsh_free_imagelist(&data_swap);
1398 xsh_free_imagelist(&errs_swap);
1399 /* loop over each IFU slice */
1400 for(j=1;j<=3;j++) {
1401
1402 /* collapse the cube along direction cross to the ifu slices j */
1403 check(flux=cpl_image_get_mean_window(data_mean_ima,j,zstart,j,zend));
1404 if (index != 0) {
1405 sprintf(qc_key,"%s%d%d VAL",XSH_QC_FLUX,j,index);
1406 } else {
1407 sprintf(qc_key,"%s%d VAL",XSH_QC_FLUX,j);
1408 }
1409 sprintf(comment,"Flux in slic %d, %4.0f-%4.0f nm",j,ws,we);
1410
1411 check(cpl_propertylist_append_double(phead,qc_key,flux));
1412 check(cpl_propertylist_set_comment(phead,qc_key,comment));
1413 err=cpl_image_get_mean_window(errs_mean_ima,j,zstart,j,zend);
1414 if (index != 0) {
1415 sprintf(qc_key,"%s%d%d ERR",XSH_QC_FLUX,j,index);
1416 } else {
1417 sprintf(qc_key,"%s%d ERR",XSH_QC_FLUX,j);
1418 }
1419 sprintf(comment,"Err Flux in slic %d, %4.0f-%4.0f nm",j,ws,we);
1420 cpl_propertylist_update_double(phead,qc_key,err);
1421 cpl_propertylist_set_comment(phead,qc_key,comment);
1422 if (index != 0) {
1423 sprintf(qc_key,"%s%d%d SN",XSH_QC_FLUX,j,index);
1424 } else {
1425 sprintf(qc_key,"%s%d SN",XSH_QC_FLUX,j);
1426 }
1427 sprintf(comment,"SNR in slic %d, %4.0f-%4.0f nm",j,ws,we);
1428
1429 sn = (fabs(err) > 1.e-37) ? flux/err : -999;
1430
1431 cpl_propertylist_append_double(phead,qc_key,sn);
1432 cpl_propertylist_set_comment(phead,qc_key,comment);
1433
1434 }
1435 xsh_free_image(&data_mean_ima);
1436 xsh_free_image(&errs_mean_ima);
1437 cleanup:
1438 return cpl_error_get_code();
1439
1440}
1441
1442/*---------------------------------------------------------------------------*/
1455/*---------------------------------------------------------------------------*/
1456
1457static cpl_error_code
1459 cpl_imagelist* errs,
1460 cpl_imagelist* qual,
1462 cpl_propertylist* phead)
1463{
1464
1465
1466 cpl_ensure(data != NULL, CPL_ERROR_NULL_INPUT,CPL_ERROR_NULL_INPUT);
1467 cpl_ensure(errs != NULL, CPL_ERROR_NULL_INPUT,CPL_ERROR_NULL_INPUT);
1468 cpl_ensure(qual != NULL, CPL_ERROR_NULL_INPUT,CPL_ERROR_NULL_INPUT);
1469 cpl_ensure(instrument != NULL, CPL_ERROR_NULL_INPUT,CPL_ERROR_NULL_INPUT);
1470 cpl_ensure(phead != NULL, CPL_ERROR_NULL_INPUT,CPL_ERROR_NULL_INPUT);
1471
1473 check( xsh_monitor_spectrum3D_flux_qc(data,errs,qual,phead,450,550,0));
1474 check( xsh_monitor_spectrum3D_flux_qc(data,errs,qual,phead,450,470,1));
1475 check( xsh_monitor_spectrum3D_flux_qc(data,errs,qual,phead,510,530,2));
1476 }
1478 xsh_monitor_spectrum3D_flux_qc(data,errs,qual,phead,600,999,0);
1479 xsh_monitor_spectrum3D_flux_qc(data,errs,qual,phead,672,680,1);
1480 xsh_monitor_spectrum3D_flux_qc(data,errs,qual,phead,745,756,2);
1481 xsh_monitor_spectrum3D_flux_qc(data,errs,qual,phead,992,999,3);
1482
1483 }
1485 xsh_monitor_spectrum3D_flux_qc(data,errs,qual,phead,1100,2450,0);
1486 check(xsh_monitor_spectrum3D_flux_qc(data,errs,qual,phead,1514,1548,1));
1487 /* TODO support JH
1488 if(!xsh_instrument_nir_is_JH(in_frm,instrument)) {
1489 check(xsh_monitor_spectrum3D_flux_qc(data,errs,qual,phead,2214,2243,2));
1490 }
1491 */
1492 }
1493
1494 cleanup:
1495
1496 return cpl_error_get_code();
1497}
1498
1499
1500
1501
1502
1503/*---------------------------------------------------------------------------*/
1521/*---------------------------------------------------------------------------*/
1522cpl_error_code
1523xsh_build_ifu_cube(cpl_frame* div_frame,
1524 cpl_frame* ifu_cfg_tab_frame,
1525 cpl_frame* ifu_cfg_cor_frame,
1526 cpl_frame* spectral_format_frame,
1527 cpl_frame* model_config_frame,
1528 cpl_frame* wavesol_frame,
1530 cpl_frameset* frameset,
1531 cpl_parameterlist* parameters,
1532 xsh_rectify_param * rectify_par,
1533 const char* recipe_id,
1534 const char* rec_prefix,
1535 const int frame_is_object)
1536{
1537
1538
1539 /* Parameters */
1540
1541 cpl_image* data_img=NULL;
1542 cpl_image* errs_img=NULL;
1543 cpl_image* qual_img=NULL;
1544
1545 cpl_image* data_tmp=NULL;
1546 cpl_image* errs_tmp=NULL;
1547 cpl_image* qual_tmp=NULL;
1548 cpl_image* mask_tmp=NULL;
1549 cpl_frame* qc_trace_frame=NULL;
1550 cpl_frame* qc_trace_merged_frame=NULL;
1551 cpl_frame* frame_data_cube=NULL;
1552 cpl_frame* frame_errs_cube=NULL;
1553 cpl_frame* frame_qual_cube=NULL;
1554 cpl_frame* frame_merged_cube=NULL;
1555
1556
1557 /* Intermediate frames */
1558
1559 /* variables */
1560 XSH_ARM arm;
1561 xsh_xs_3 model_config ;
1562 xsh_rec_list * rec_list = NULL ;
1563
1564 cpl_table* ifu_cfg_tab=NULL;
1565 cpl_table* sp_fmt_tab=NULL;
1566 cpl_vector* profile=NULL;
1567 cpl_imagelist* data_cube=NULL;
1568 cpl_imagelist* errs_cube=NULL;
1569 cpl_imagelist* qual_cube=NULL;
1570
1571 cpl_imagelist* data_cube_merge=NULL;
1572 cpl_imagelist* errs_cube_merge=NULL;
1573 cpl_imagelist* qual_cube_merge=NULL;
1574 cpl_imagelist* mask_data_merge=NULL;
1575 cpl_imagelist* mask_errs_merge=NULL;
1576 cpl_imagelist* mask_qual_merge=NULL;
1577
1578
1579 float* pima=NULL;
1580 float* perr=NULL;
1581 int* pqua=NULL;
1582 int* pmsk=NULL;
1583 int nord=0;
1584 int save_size_uvb=4;
1585 int save_size_vis=2;
1586 int save_size_nir=2;
1587 int save_size=0;
1588 const int peack_search_hsize=5;
1589 int method=0; /* Gaussian */
1590
1591
1592
1593 cpl_propertylist* data_plist=NULL;
1594 cpl_propertylist* errs_plist=NULL;
1595 cpl_propertylist* qual_plist=NULL;
1596
1597
1598 const char* sci_name=NULL;
1599 const char* sp_fmt_name=NULL;
1600 const char* ifu_cfg_name=NULL;
1601
1602 char name[256];
1603 char data_extid[40];
1604 char errs_extid[40];
1605 char qual_extid[40];
1606
1607 char qualifier[10];
1608 int binx=1;
1609 int biny=1;
1610
1611 double flux_upp=0;
1612 double flux_cen=0;
1613 double flux_low=0;
1614
1615
1616 double errs_upp=0;
1617 double errs_cen=0;
1618 double errs_low=0;
1619 double cube_wave_min=0;
1620 double cube_wave_max=0;
1621
1622 //double cube_wmin=550;
1623 //double cube_wmax=1000;
1624 double cube_wstep=0;
1625 double cube_sstep=0;
1626
1627 int is=0;
1628
1629 int ord=0;
1630 int ord_min=0;
1631 int ord_max=0;
1632
1633 double wave=0;
1634 double wave_min=0;
1635 double wave_max=0;
1636 double wave_step=0;
1637
1638 double wave_min_old=0;
1639 //double wave_max_old=0;
1640
1641 int nstep_off=0;
1642
1643 double s=0;
1644 double s_min=-2;
1645 double s_max=2;
1646 double s_step=0;
1647
1648 double x=0;
1649 double y=0;
1650 double s_upp=0;
1651 double s_low=0;
1652
1653
1654 int radius=rectify_par->rectif_radius;
1655
1656 double confidence=0;
1657 int naxis1=0;
1658 int naxis2=0;
1659 int naxis3=0;
1660 int ik=0;
1661 //double wave_cen=0;
1662
1663 /*
1664 double stupid_offset_uvb=0.145;
1665 double stupid_offset_vis=0.09;
1666 double stupid_offset_nir=0.16;
1667 */
1668 double w_low_coeff1=0;
1669 double w_low_coeff2=0;
1670
1671 double w_upp_coeff1=0;
1672 double w_upp_coeff2=0;
1673
1674 double w_low_coeff1_uvb=-0.002972; /*-0.00085*/
1675 double w_low_coeff2_uvb=2.26497e-6+0.2e-6;/*-2.26e-6;*/
1676 double w_upp_coeff1_uvb=8.3355331e-5; /*-0.00085*/
1677 double w_upp_coeff2_uvb=-1.0682e-6;/*-2.26e-6;*/
1678
1679
1680
1681 double w_low_coeff1_vis=-0.0016549569;
1682 double w_low_coeff2_vis=1.183805e-6+0.04e-6;
1683 double w_upp_coeff1_vis=-0.0016610719;
1684 double w_upp_coeff2_vis=1.0823013e-6+0.02e-6;
1685
1686 double w_low_coeff1_nir=-0.00010;
1687 double w_low_coeff2_nir=-5.0e-9;
1688
1689 double w_upp_coeff1_nir=0;
1690 double w_upp_coeff2_nir=0;
1691
1692
1693 double s_upp_off_uvb=4.0514429+0.16;/*1.91319135+stupid_offset_uvb;*/
1694 double s_upp_off_vis=5.2504895+0.16-0.703;/*+0.05; */
1695 double s_upp_off_nir=3.5452+0.31;
1696 double s_upp_off=0;
1697
1698 double s_low_off_uvb=-3.4636662+0.16;/*-2.1704691+stupid_offset_uvb;*/
1699
1700 double s_low_off_vis=-2.9428071+0.16-0.645;/*-0.07; */
1701 double s_low_off_nir=-4.451682+0.21+0.31;
1702
1703 double s_low_off=0;
1704
1705 int w_step_fct=1;
1706
1707 int status=0;
1708 int mk=0;
1709 char tag[256];
1710 xsh_wavesol *wavesol=NULL;
1711 cpl_frame* trace_corr_tab=NULL;
1712 double s_upp_d0=0;
1713 double s_low_d0=0;
1714 double s_upp_d1=0;
1715 double s_low_d1=0;
1716 double s_upp_d2=0;
1717 double s_low_d2=0;
1718 int cut_uvb_spectrum=0;
1719
1720
1721
1722 if(frame_is_object) {
1723 sprintf(qualifier,"OBJ");
1724 } else {
1725 sprintf(qualifier,"SKY");
1726 }
1727 /* NEW */
1728 /* Find STD star centroids (for each IFU slice):
1729 input STD star frame
1730 input order centre traces
1731 input order edges traces
1732 input slit map
1733 output STD star centroids traces table+ QC on centroid measure
1734 */
1735
1736 /* Find STD star centroids traces slices interdistance
1737 (each measured at the same wavelength)
1738 Find STD star centroids distances to IFU edges traces
1739 (each measured at the same wavelength)
1740 input STD star centroid traces table
1741 input order edges/IFU slices traces table
1742 input slit map (orders are spatially over-sized)
1743 input wave map (orders are spatially over-sized)
1744 output table with STD star centroids traces slices interdistance
1745 output corresponding QC(statistics)
1746 output table with STD star distances to IFU edges traces
1747 output corresponding QC(statistics)
1748 ? output "measured" slit map
1749 ? output "measured" wave map
1750 */
1751
1752
1753 check(sci_name=cpl_frame_get_filename(div_frame));
1754 check(data_plist=cpl_propertylist_load(sci_name,0));
1755 errs_plist=cpl_propertylist_duplicate(data_plist);
1756 qual_plist=cpl_propertylist_duplicate(data_plist);
1757
1758 arm=instrument->arm;
1759
1760
1761 cube_wstep=rectify_par->rectif_bin_lambda;
1762 cube_sstep=rectify_par->rectif_bin_space;
1763
1764 if(arm == XSH_ARM_UVB) {
1765 //cube_wmin=300;
1766 //cube_wmax=550;
1767
1768 binx=xsh_pfits_get_binx(data_plist);
1769 biny=xsh_pfits_get_biny(data_plist);
1770
1771 s_low_off=s_low_off_uvb;
1772 s_upp_off=s_upp_off_uvb;
1773 w_low_coeff1=w_low_coeff1_uvb;
1774 w_low_coeff2=w_low_coeff2_uvb;
1775
1776 w_upp_coeff1=w_upp_coeff1_uvb;
1777 w_upp_coeff2=w_upp_coeff2_uvb;
1778 save_size=save_size_uvb;
1779 cut_uvb_spectrum=xsh_parameters_cut_uvb_spectrum_get(recipe_id,
1780 parameters);
1781
1782 //model_config.temper= xsh_pfits_get_temp2(data_plist);
1783 } else if(arm == XSH_ARM_VIS) {
1784 //cube_wmin=550;
1785 //cube_wmax=1000;
1786
1787 //model_config.temper= xsh_pfits_get_temp5(data_plist);
1788 binx=xsh_pfits_get_binx(data_plist);
1789 biny=xsh_pfits_get_biny(data_plist);
1790
1791 s_low_off=s_low_off_vis;
1792 s_upp_off=s_upp_off_vis;
1793
1794 w_low_coeff1=w_low_coeff1_vis;
1795 w_low_coeff2=w_low_coeff2_vis;
1796
1797 w_upp_coeff1=w_upp_coeff1_vis;
1798 w_upp_coeff2=w_upp_coeff2_vis;
1799 save_size=save_size_vis;
1800
1801
1802 } else if(arm == XSH_ARM_NIR) {
1803 //cube_wmin=1000;
1804 //cube_wmax=2500;
1805
1806 //model_config.temper= xsh_pfits_get_temp82(data_plist);
1807
1808 s_low_off=s_low_off_nir;
1809 s_upp_off=s_upp_off_nir;
1810
1811
1812 w_low_coeff1=w_low_coeff1_nir;
1813 w_low_coeff2=w_low_coeff2_nir;
1814
1815 w_upp_coeff1=w_upp_coeff1_nir;
1816 w_upp_coeff2=w_upp_coeff2_nir;
1817 save_size=save_size_nir;
1818
1819 }
1820
1821
1822
1823 if(ifu_cfg_tab_frame) {
1824
1825 ifu_cfg_name=cpl_frame_get_filename(ifu_cfg_tab_frame);
1826 ifu_cfg_tab=cpl_table_load(ifu_cfg_name,1,0);
1827 s_upp_off=cpl_table_get_double(ifu_cfg_tab,"S_UPP_OFF",arm,&status);
1828 s_low_off=cpl_table_get_double(ifu_cfg_tab,"S_LOW_OFF",arm,&status);
1829
1830 w_upp_coeff1=cpl_table_get_double(ifu_cfg_tab,"W_UPP_COEF1",arm,&status);
1831 w_low_coeff1=cpl_table_get_double(ifu_cfg_tab,"W_LOW_COEF1",arm,&status);
1832 w_upp_coeff2=cpl_table_get_double(ifu_cfg_tab,"W_UPP_COEF2",arm,&status);
1833 w_low_coeff2=cpl_table_get_double(ifu_cfg_tab,"W_LOW_COEF2",arm,&status);
1834
1835 w_step_fct=cpl_table_get_int(ifu_cfg_tab,"W_STEP_FCT",arm,&status);
1836 xsh_msg("s_upp_off=%10.8g",s_upp_off);
1837 xsh_msg("s_low_off=%10.8g",s_low_off);
1838
1839
1840 xsh_msg("w_upp_coeff1=%10.8g",w_upp_coeff1);
1841 xsh_msg("w_upp_coeff2=%10.8g",w_upp_coeff2);
1842 xsh_msg("w_low_coeff1=%10.8g",w_low_coeff1);
1843 xsh_msg("w_low_coeff2=%10.8g",w_low_coeff2);
1844
1845 xsh_msg("w_step_fct=%d",w_step_fct);
1846
1847 }
1848
1849 if(ifu_cfg_cor_frame) {
1850 ifu_cfg_name=cpl_frame_get_filename(ifu_cfg_cor_frame);
1851 ifu_cfg_tab=cpl_table_load(ifu_cfg_name,1,0);
1852 s_upp_d0=cpl_table_get_double(ifu_cfg_tab,"DIFF_T12",0,&status);
1853 s_low_d0=cpl_table_get_double(ifu_cfg_tab,"DIFF_T32",0,&status);
1854 s_upp_d1=cpl_table_get_double(ifu_cfg_tab,"DIFF_T12",1,&status);
1855 s_low_d1=cpl_table_get_double(ifu_cfg_tab,"DIFF_T32",1,&status);
1856 s_upp_d2=cpl_table_get_double(ifu_cfg_tab,"DIFF_T12",2,&status);
1857 s_low_d2=cpl_table_get_double(ifu_cfg_tab,"DIFF_T32",2,&status);
1858 xsh_msg("upp_d0=%g upp_d1=%g upp_d2=%g low_d0=%g low_d1=%g low_d2=%g",
1859 s_upp_d0,s_upp_d1,s_upp_d2,s_low_d0,s_low_d1,s_low_d2);
1860
1861 }
1862
1863 cube_wstep*=w_step_fct;
1864 if(model_config_frame!=NULL) {
1865 check( xsh_model_config_load_best( model_config_frame, &model_config));
1866 xsh_model_binxy(&model_config,binx,biny);
1867 } else {
1868 check( wavesol = xsh_wavesol_load( wavesol_frame, instrument));
1869 }
1870 check( rec_list = xsh_rec_list_create( instrument));
1871 check( profile = cpl_vector_new( CPL_KERNEL_DEF_SAMPLES));
1872 check(cpl_vector_fill_kernel_profile( profile,rectify_par->kernel_type,
1873 rectify_par->rectif_radius));
1874
1875 check(sci_name=cpl_frame_get_filename(div_frame));
1876 check(data_img=cpl_image_load(sci_name,XSH_PRE_DATA_TYPE,0,0));
1877 check(errs_img=cpl_image_load(sci_name,XSH_PRE_ERRS_TYPE,0,1));
1878 check(qual_img=cpl_image_load(sci_name,XSH_PRE_QUAL_TYPE,0,2));
1879 int nx=cpl_image_get_size_x(data_img);
1880 int ny=cpl_image_get_size_y(data_img);
1881 int* pqual_img=cpl_image_get_data_int(qual_img);
1882 cpl_image* var_img=cpl_image_duplicate(data_img);
1883 cpl_image_power(var_img,2.);
1884 cpl_vector* profile2=cpl_vector_duplicate(profile);
1885 cpl_vector_power(profile2,2.);
1886 check(sp_fmt_name=cpl_frame_get_filename(spectral_format_frame));
1887 check(sp_fmt_tab=cpl_table_load(sp_fmt_name,1,0));
1888 check(ord_min=cpl_table_get_column_min(sp_fmt_tab,"ORDER"));
1889 check(ord_max=cpl_table_get_column_max(sp_fmt_tab,"ORDER"));
1890 nord=ord_max-ord_min+1;
1891
1892 check(cube_wave_min=cpl_table_get(sp_fmt_tab,"WLMIN",nord-1,&status));
1893 check(cube_wave_max=cpl_table_get(sp_fmt_tab,"WLMAX",0,&status));
1894
1895
1896 xsh_msg_debug("cube_wave_min=%10.8g cube_wave_max=%10.8g",cube_wave_min,cube_wave_max);
1897
1898 /* on X we have the different IFU slices */
1899 naxis1=3;
1900 xsh_msg_debug("ord_min=%d ord_max=%d",ord_min,ord_max);
1901 wave_step= cube_wstep;
1902 s_step= cube_sstep;
1903 naxis2=(int)((s_max-s_min)/s_step+0.5)+1;
1904
1905
1906 xsh_free_imagelist(&data_cube_merge);
1907 xsh_free_imagelist(&errs_cube_merge);
1908 xsh_free_imagelist(&qual_cube_merge);
1909 xsh_free_imagelist(&mask_data_merge);
1910 xsh_free_imagelist(&mask_errs_merge);
1911 xsh_free_imagelist(&mask_qual_merge);
1912
1913 data_cube_merge=cpl_imagelist_new();
1914 errs_cube_merge=cpl_imagelist_new();
1915 qual_cube_merge=cpl_imagelist_new();
1916 mask_data_merge=cpl_imagelist_new();
1917 mask_errs_merge=cpl_imagelist_new();
1918 mask_qual_merge=cpl_imagelist_new();
1919
1920
1921 xsh_free_image(&data_tmp);
1922 xsh_free_image(&errs_tmp);
1923 xsh_free_image(&qual_tmp);
1924 xsh_free_image(&mask_tmp);
1925 //double var_low=0, var_cen=0, var_upp=0;
1926 mk=0;
1927
1928 for( ord = ord_max; ord >= ord_min; ord-- ) {
1929
1930 xsh_free_imagelist(&data_cube);
1931 xsh_free_imagelist(&errs_cube);
1932 xsh_free_imagelist(&qual_cube);
1933
1934 data_cube=cpl_imagelist_new();
1935 errs_cube=cpl_imagelist_new();
1936 qual_cube=cpl_imagelist_new();
1937
1938
1939 data_tmp=cpl_image_new(naxis1,naxis2,XSH_PRE_DATA_TYPE);
1940 errs_tmp=cpl_image_new(naxis1,naxis2,XSH_PRE_ERRS_TYPE);
1941 qual_tmp=cpl_image_new(naxis1,naxis2,XSH_PRE_QUAL_TYPE);
1942 mask_tmp=cpl_image_new(naxis1,naxis2,XSH_PRE_QUAL_TYPE);
1943
1944 check(pima=cpl_image_get_data_float(data_tmp));
1945 check(perr=cpl_image_get_data_float(errs_tmp));
1946 check(pqua=cpl_image_get_data_int(qual_tmp));
1947 check(pmsk=cpl_image_get_data_int(mask_tmp));
1948
1949
1950
1951 check(wave_min=cpl_table_get(sp_fmt_tab,"WLMIN",ord-ord_min,&status));
1952 check(wave_max=cpl_table_get(sp_fmt_tab,"WLMAX",ord-ord_min,&status));
1953 //cube_wave_min=wave_min;
1954 naxis3=(wave_max-wave_min)/wave_step+1;
1955 //wave_cen=wave_min+naxis3/2*wave_step;
1956 xsh_msg_debug("order=%d naxis1=%d,naxis2=%d naxis3=%d, wave_min=%g wave_max=%g",
1957 ord,naxis1,naxis2,naxis3,wave_min,wave_max);
1958
1959 /* in case we go from ord_max to ord_min */
1960 if(ord<ord_max) {
1961 nstep_off=(int)((wave_min-wave_min_old)/wave_step+0.5);
1962 wave_min=wave_min_old+nstep_off*wave_step;
1963 }
1964
1965 ik=0;
1966 for( wave = wave_min; wave <= wave_max; wave+=wave_step) {
1967 is=0;
1968 for( s = s_min; s <= s_max; s+=s_step ) {
1969 int qual_upp=0;
1970 int qual_cen=0;
1971 int qual_low=0;
1972 /* upper slice */
1973 s_upp=-s+s_upp_off+w_upp_coeff1*wave+w_upp_coeff2*wave*wave;
1974 s_upp-=(s_upp_d0+s_upp_d1*wave+s_upp_d2*wave*wave)*s_step;
1975 if(wavesol!=NULL) {
1976 check( x= xsh_wavesol_eval_polx( wavesol,wave, ord,s_upp));
1977 check( y= xsh_wavesol_eval_poly( wavesol,wave, ord,s_upp));
1978 } else {
1979
1980 check(xsh_model_get_xy(&model_config,instrument,wave,ord,s_upp,
1981 &x,&y));
1982 }
1983 //y+=s_upp_d0+s_upp_d1*wave+s_upp_d2*wave*wave;
1984
1985 //xsh_msg("is=%d s_min=%g s_max=%g s=%g",is,s_min,s_max,s);
1986 check(flux_upp=cpl_image_get_interpolated(data_img,x,y,
1987 profile,radius,
1988 profile,radius,
1989 &confidence));
1990
1991 if (confidence <= 0) {
1992 pima[is*naxis1+0]=0;
1993 perr[is*naxis1+0]=1;
1994 pqua[is*naxis1+0]=QFLAG_MISSING_DATA;
1995 pmsk[is*naxis1+0]=1;
1996 /* If flux interpolation ok */
1997 } else {
1998
1999 check(errs_upp=cpl_image_get_interpolated(errs_img,x,y,
2000 profile,radius,
2001 profile,radius,
2002 &confidence));
2003 /*
2004 if(confidence > 0 && var_upp>0) {
2005 errs_upp=sqrt(var_upp);
2006 } else {
2007 errs_upp=0;
2008 }
2009 */
2010 /* Now propagate flags */
2011 int ix=(int)(x+0.5);
2012 int iy=(int)(x+0.5);
2013 int jj_nx=0;
2014 int jj_min= ( iy - radius >= 0 ) ? iy - radius: 0;
2015 int jj_max= ( iy + radius < ny ) ? iy + radius: 0;
2016 //int ii_min= ( ix - radius >= 0 ) ? ix - radius: 0;
2017 //int ii_max= ( ix + radius < nx ) ? ix + radius: 0;
2018 for(int jj=jj_min;jj<jj_max;jj++) {
2019 jj_nx=jj*nx;
2020 for(int ii=ix-radius;ii<ix+radius-1;ii++) {
2021 qual_upp |= pqual_img[jj_nx+ii]; //UMR
2022 }
2023 }
2024
2025 pima[is*naxis1+0]=flux_upp;
2026 perr[is*naxis1+0]=errs_upp;
2027 pqua[is*naxis1+0]=qual_upp;
2028 pmsk[is*naxis1+0]=1;
2029 }
2030 /* central slice */
2031 if(wavesol!=NULL) {
2032 check( x= xsh_wavesol_eval_polx( wavesol,wave, ord,s));
2033 check( y= xsh_wavesol_eval_poly( wavesol,wave, ord,s));
2034 } else {
2035 check(xsh_model_get_xy(&model_config,instrument,wave,ord,s,
2036 &x,&y));
2037 }
2038
2039
2040 check(flux_cen=cpl_image_get_interpolated(data_img,x,y,
2041 profile,radius,
2042 profile,radius,
2043 &confidence));
2044
2045 if (confidence <= 0) {
2046 pima[is*naxis1+1]=0;
2047 perr[is*naxis1+1]=1;
2048 pqua[is*naxis1+1]=QFLAG_MISSING_DATA;
2049 pmsk[is*naxis1+1]=1;
2050 /* If flux interpolation ok */
2051 } else {
2052
2053 check(errs_cen=cpl_image_get_interpolated(errs_img,x,y,
2054 profile,radius,
2055 profile,radius,
2056 &confidence));
2057 /*
2058 if(confidence > 0 && var_cen>0) {
2059 errs_cen=sqrt(var_cen);
2060 } else {
2061 errs_cen=0;
2062 }
2063 */
2064
2065 /* Now propagate flags */
2066 int ix=(int)(x+0.5);
2067 int iy=(int)(x+0.5);
2068 int jj_nx=0;
2069 int jj_min= ( iy - radius >= 0 ) ? iy - radius: 0;
2070 int jj_max= ( iy + radius < ny ) ? iy + radius: 0;
2071 //int ii_min= ( ix - radius >= 0 ) ? ix - radius: 0;
2072 //int ii_max= ( ix + radius < nx ) ? ix + radius: 0;
2073 for(int jj=jj_min;jj<jj_max;jj++) {
2074 jj_nx=jj*nx;
2075 for(int ii=ix-radius;ii<ix+radius-1;ii++) {
2076 qual_cen |= pqual_img[jj_nx+ii]; //UMR
2077 }
2078 }
2079
2080 pima[is*naxis1+1]=flux_cen;
2081 perr[is*naxis1+1]=errs_cen;
2082 pqua[is*naxis1+1]=qual_cen;
2083 pmsk[is*naxis1+1]=1;
2084 }
2085 /* lower slice */
2086 s_low=-s+s_low_off+w_low_coeff1*wave+w_low_coeff2*wave*wave;
2087 s_low-=(s_low_d0+s_low_d1*wave+s_low_d2*wave*wave)*s_step;
2088 if(wavesol!=NULL) {
2089 check( x= xsh_wavesol_eval_polx( wavesol,wave, ord,s_low));
2090 check( y= xsh_wavesol_eval_poly( wavesol,wave, ord,s_low));
2091 } else {
2092
2093 check(xsh_model_get_xy(&model_config,instrument,wave,ord,s_low,
2094 &x,&y));
2095 }
2096
2097 check(flux_low=cpl_image_get_interpolated(data_img,x,y,
2098 profile,radius,
2099 profile,radius,
2100 &confidence));
2101
2102 if (confidence <= 0 ) {
2103 pima[is*naxis1+2]=0;
2104 perr[is*naxis1+2]=1;
2105 pqua[is*naxis1+2]=QFLAG_MISSING_DATA;
2106 pmsk[is*naxis1+2]=1;
2107 /* If flux interpolation ok */
2108 } else {
2109
2110 check(errs_low=cpl_image_get_interpolated(errs_img,x,y,
2111 profile,radius,
2112 profile,radius,
2113 &confidence));
2114 /*
2115 if(var_low>0 ) {
2116 errs_low=sqrt(confidence > 0 && var_low);
2117 } else {
2118 errs_low=0;
2119 }
2120 */
2121
2122 /* Now propagate flags */
2123 int ix=(int)(x+0.5);
2124 int iy=(int)(x+0.5);
2125 int jj_nx=0;
2126 int jj_min= ( iy - radius >= 0 ) ? iy - radius: 0;
2127 int jj_max= ( iy + radius < ny ) ? iy + radius: 0;
2128 //int ii_min= ( ix - radius >= 0 ) ? ix - radius: 0;
2129 //int ii_max= ( ix + radius < nx ) ? ix + radius: 0;
2130 for(int jj=jj_min;jj<jj_max;jj++) {
2131 jj_nx=jj*nx;
2132 for(int ii=ix-radius;ii<ix+radius-1;ii++) {
2133 qual_low |= pqual_img[jj_nx+ii]; //UMR
2134 }
2135 }
2136
2137 pima[is*naxis1+2]=flux_low;
2138 perr[is*naxis1+2]=errs_low;
2139 pqua[is*naxis1+2]=qual_low;
2140 pmsk[is*naxis1+2]=1;
2141 }
2142
2143 is++;
2144 } /* end loop on s */
2145 //wave_max_old=wave_max;
2146 wave_min_old=wave_min;
2147
2148 check(cpl_imagelist_set(data_cube,cpl_image_duplicate(data_tmp),ik));
2149 check(cpl_imagelist_set(errs_cube,cpl_image_duplicate(errs_tmp),ik));
2150 check(cpl_imagelist_set(qual_cube,cpl_image_duplicate(qual_tmp),ik));
2151
2152 /*
2153 check(xsh_iml_merge_avg(&data_cube_merge,&mask_data_merge,
2154 data_tmp,mask_tmp,mk));
2155 check(xsh_iml_merge_avg(&errs_cube_merge,&mask_errs_merge,
2156 errs_tmp,mask_tmp,mk));
2157 check(xsh_iml_merge_avg(&qual_cube_merge,&mask_qual_merge,
2158 qual_tmp,mask_tmp,mk));
2159 */
2160 check(xsh_iml_merge_wgt(&data_cube_merge, &errs_cube_merge, &qual_cube_merge,
2161 data_tmp, errs_tmp, qual_tmp, mk,
2163 ik++;
2164 if(ord==ord_max) {
2165 mk++;
2166 } else {
2167 mk=(int)((wave-cube_wave_min)/wave_step+0.5);
2168 }
2169 } /* end loop on wave */
2170
2171 sprintf(tag,"%s_ORDER3D_DATA_%s_%s",rec_prefix,qualifier,
2172 xsh_arm_tostring(arm));
2173 sprintf(name,"%s.fits",tag);
2174 if(frame_is_object==0) {
2176 }
2177 xsh_cube_set_wcs(data_plist,wave_min,
2178 0.6,s_step,wave_step,
2179 1,naxis2/2.,0);
2180
2181 xsh_cube_set_wcs(errs_plist,wave_min,
2182 0.6,s_step,wave_step,
2183 1,naxis2/2.,0);
2184
2185 xsh_cube_set_wcs(qual_plist,wave_min,
2186 0.6,s_step,wave_step,
2187 1,naxis2/2.,0);
2188
2189 sprintf(data_extid,"ORD%2.2d_FLUX",ord);
2190 xsh_pfits_set_extname (data_plist, data_extid);
2191 sprintf(errs_extid,"ORD%2.2d_ERRS",ord);
2192 xsh_pfits_set_extname (errs_plist, errs_extid);
2193 sprintf(qual_extid,"ORD%2.2d_QUAL",ord);
2194 xsh_pfits_set_extname (qual_plist, qual_extid);
2195
2196 xsh_plist_set_extra_keys(data_plist,"IMAGE","DATA","RMSE",
2197 data_extid,errs_extid,qual_extid,0);
2198
2199 xsh_plist_set_extra_keys(errs_plist,"IMAGE","ERROR","RMSE",
2200 data_extid,errs_extid,qual_extid,1);
2201
2202 xsh_plist_set_extra_keys(qual_plist,"IMAGE","QUALITY","FLAG32BIT",
2203 data_extid,errs_extid,qual_extid,2);
2204
2205
2206
2207 if(ord==ord_max) {
2208 xsh_pfits_set_pcatg(data_plist,tag);
2209 frame_data_cube=xsh_frame_product(name,tag,CPL_FRAME_TYPE_IMAGE,CPL_FRAME_GROUP_PRODUCT,CPL_FRAME_LEVEL_FINAL);
2210 cpl_imagelist_save(data_cube,name,XSH_SPECTRUM_DATA_BPP,data_plist,CPL_IO_DEFAULT);
2211 cpl_imagelist_save(errs_cube,name,XSH_SPECTRUM_ERRS_BPP,errs_plist,CPL_IO_EXTEND);
2212 cpl_imagelist_save(qual_cube,name,XSH_SPECTRUM_QUAL_BPP,qual_plist,CPL_IO_EXTEND);
2213 } else {
2214 cpl_imagelist_save(data_cube,name,XSH_SPECTRUM_DATA_BPP,data_plist,CPL_IO_EXTEND);
2215 cpl_imagelist_save(errs_cube,name,XSH_SPECTRUM_ERRS_BPP,errs_plist,CPL_IO_EXTEND);
2216 cpl_imagelist_save(qual_cube,name,XSH_SPECTRUM_QUAL_BPP,qual_plist,CPL_IO_EXTEND);
2217 }
2218
2219
2220 sprintf(name,"QC%s_%2.2d.fits",tag,ord);
2221
2222 xsh_free_image(&data_tmp);
2223 xsh_free_image(&errs_tmp);
2224 xsh_free_image(&qual_tmp);
2225 xsh_free_image(&mask_tmp);
2226
2227 } /* end loop on ord */
2228
2229 sprintf(tag,"%s_MERGE3D_DATA_%s_%s",rec_prefix,qualifier,
2230 xsh_arm_tostring(arm));
2231 sprintf(name,"%s.fits",tag);
2232 xsh_cube_set_wcs(data_plist,cube_wave_min,
2233 0.6,s_step,wave_step,
2234 1,naxis2/2.,0);
2235
2236 xsh_cube_set_wcs(errs_plist,cube_wave_min,
2237 0.6,s_step,wave_step,
2238 1,naxis2/2.,0);
2239
2240 xsh_cube_set_wcs(qual_plist,cube_wave_min,
2241 0.6,s_step,wave_step,
2242 1,naxis2/2.,0);
2243
2244 sprintf(data_extid,"FLUX");
2245 xsh_pfits_set_extname (data_plist, data_extid);
2246 sprintf(errs_extid,"ERRS");
2247 xsh_pfits_set_extname (errs_plist, errs_extid);
2248 sprintf(qual_extid,"QUAL");
2249 xsh_pfits_set_extname (qual_plist, qual_extid);
2250
2251 xsh_plist_set_extra_keys(data_plist,"IMAGE","DATA","RMSE",
2252 data_extid,errs_extid,qual_extid,0);
2253 xsh_plist_set_extra_keys(errs_plist,"IMAGE","ERROR","RMSE",
2254 data_extid,errs_extid,qual_extid,1);
2255 xsh_plist_set_extra_keys(qual_plist,"IMAGE","QUALITY","FLAG32BIT",
2256 data_extid,errs_extid,qual_extid,2);
2257 xsh_pfits_set_pcatg(data_plist,tag);
2258
2259 if(cut_uvb_spectrum) {
2260 check(xsh_imagelist_cut_dichroic_uvb(data_cube_merge,errs_cube_merge,
2261 qual_cube_merge,data_plist));
2262 }
2263
2264 if( model_config_frame != NULL) {
2265
2266 xsh_monitor_spectrum3D_flux(data_cube_merge,errs_cube_merge,
2267 qual_cube_merge,instrument,data_plist);
2268 }
2269
2270 cpl_imagelist_save(data_cube_merge,name,XSH_SPECTRUM_DATA_BPP,data_plist,CPL_IO_DEFAULT);
2271
2272 cpl_imagelist_save(errs_cube_merge,name,XSH_SPECTRUM_ERRS_BPP,errs_plist,CPL_IO_EXTEND);
2273
2274 cpl_imagelist_save(qual_cube_merge,name,XSH_SPECTRUM_QUAL_BPP,qual_plist,CPL_IO_EXTEND);
2275
2276 xsh_msg_debug("merge cube size=%" CPL_SIZE_FORMAT "",cpl_imagelist_get_size(data_cube_merge));
2277
2278 frame_merged_cube=xsh_frame_product(name,tag,CPL_FRAME_TYPE_IMAGE,
2279 CPL_FRAME_GROUP_PRODUCT,CPL_FRAME_LEVEL_FINAL);
2280
2281
2282 check(xsh_add_product_imagelist(frame_merged_cube,frameset,parameters,
2283 recipe_id,instrument,NULL));
2284
2285
2286 if(frame_is_object) {
2287
2288 sprintf(tag,"%s_MERGE3D_%s",rec_prefix,qualifier);
2289
2290 check(qc_trace_merged_frame=xsh_cube_qc_trace_window(frame_merged_cube,
2291 instrument,
2292 "MERGE3D",rec_prefix,
2293 save_size+1,
2294 naxis2-save_size,
2295 peack_search_hsize,
2296 method,1));
2297
2298 trace_corr_tab=xsh_crea_correct_coeff(qc_trace_merged_frame,instrument);
2299
2300 check( xsh_add_product_table( qc_trace_merged_frame, frameset,
2301 parameters, recipe_id, instrument,
2302 NULL));
2303
2304 sprintf(tag,"%s_ORDER3D_DATA_%s",rec_prefix,qualifier);
2305
2306 check(xsh_add_product_imagelist(frame_data_cube,frameset,parameters,
2307 recipe_id,instrument,tag));
2308 /*
2309 check(xsh_add_product_imagelist(frame_errs_cube,frameset,parameters,
2310 recipe_id,instrument,NULL));
2311
2312 check(xsh_add_product_imagelist(frame_qual_cube,frameset,parameters,
2313 recipe_id,instrument,NULL));
2314 */
2315 }
2316
2317 if(trace_corr_tab) {
2318 xsh_add_product_table( trace_corr_tab, frameset,
2319 parameters, recipe_id, instrument,NULL);
2320 }
2321
2322 cleanup:
2323
2324 xsh_free_propertylist(&data_plist);
2325 xsh_free_propertylist(&errs_plist);
2326 xsh_free_propertylist(&qual_plist);
2327 xsh_free_vector(&profile);
2328
2329 xsh_free_image(&data_tmp);
2330 xsh_free_image(&errs_tmp);
2331 xsh_free_image(&qual_tmp);
2332 xsh_free_image(&mask_tmp);
2333
2334 xsh_free_image(&data_img);
2335 xsh_free_image(&errs_img);
2336 xsh_free_image(&qual_img);
2337
2338 xsh_free_image(&var_img);
2339 xsh_free_vector(&profile2);
2340 xsh_free_table(&sp_fmt_tab);
2341
2342 xsh_free_imagelist(&data_cube);
2343 xsh_free_imagelist(&errs_cube);
2344 xsh_free_imagelist(&qual_cube);
2345
2346
2347 xsh_free_imagelist(&data_cube_merge);
2348 xsh_free_imagelist(&errs_cube_merge);
2349 xsh_free_imagelist(&qual_cube_merge);
2350 xsh_free_imagelist(&mask_data_merge);
2351 xsh_free_imagelist(&mask_errs_merge);
2352 xsh_free_imagelist(&mask_qual_merge);
2353
2354 xsh_rec_list_free(&rec_list);
2355 xsh_free_frame(&qc_trace_frame);
2356 xsh_free_frame(&qc_trace_merged_frame);
2357 xsh_free_frame(&trace_corr_tab);
2358 xsh_free_frame(&frame_data_cube);
2359 xsh_free_frame(&frame_errs_cube);
2360 xsh_free_frame(&frame_qual_cube);
2361 xsh_free_frame(&frame_merged_cube);
2362 /*
2363 xsh_free_frame(&spectral_format_frame);
2364 xsh_free_frame(&model_config_frame);
2365 */
2366
2367 return cpl_error_get_code();
2368
2369}
2370
2371
2372/*---------------------------------------------------------------------------*/
2379cpl_matrix * xsh_atrous( cpl_vector *spec, int nscales){
2380
2381 cpl_matrix *decomp = NULL;
2382 cpl_vector *filter = NULL;
2383 double orig_filter[3] = {3./8., 1./4., 1./16.};
2384 double *new_filter = NULL;
2385 int filter_size = 3;
2386
2387 int spec_size;
2388 cpl_vector *sp = NULL;
2389 cpl_vector *smooth = NULL;
2390 int i,k;
2391 cpl_vector *tmp_filter = NULL;
2392
2393 /* Check parameters */
2394 XSH_ASSURE_NOT_NULL( spec);
2395 XSH_ASSURE_NOT_ILLEGAL( nscales >= 1);
2396
2397 /* Init */
2398 check( spec_size = cpl_vector_get_size( spec));
2399 check( decomp = cpl_matrix_new( nscales+1, spec_size));
2400 sp = cpl_vector_duplicate( spec);
2401 smooth = cpl_vector_duplicate( spec);
2402 check( filter = cpl_vector_wrap( filter_size, orig_filter));
2403
2404
2405 for(k=0; k< nscales; k++){
2406 xsh_msg_dbg_medium("scale %d use filter with size %d", k, filter_size);
2407 /* smooth by convolution filter */
2408 check( cpl_wlcalib_xc_convolve( smooth, filter));
2409
2410 /* fill result matrix */
2411 for(i=0; i<spec_size; i++){
2412 double sp_val, smooth_val;
2413
2414 check( sp_val = cpl_vector_get( sp, i));
2415 check( smooth_val = cpl_vector_get( smooth, i));
2416 check( cpl_matrix_set( decomp, nscales-k, i, sp_val-smooth_val));
2417 }
2418
2419 /* sp = smooth */
2420 xsh_free_vector( &sp);
2421 check( sp = cpl_vector_duplicate(smooth));
2422
2423 /* generate new filter */
2424 check( tmp_filter = cpl_vector_duplicate( filter));
2425 XSH_FREE( new_filter);
2426 filter_size = cpl_vector_get_size( tmp_filter);
2427 XSH_CALLOC( new_filter, double , filter_size*2-1);
2428
2429 for(i=0; i<filter_size; i++){
2430 new_filter[2*i] = cpl_vector_get( tmp_filter, i);
2431 }
2432 xsh_free_vector( &tmp_filter);
2434 check( filter = cpl_vector_wrap( filter_size*2-1, new_filter));
2435 filter_size = filter_size*2-1;
2436 }
2437
2438
2439 /* decomp[0,*] = smooth */
2440 for(i=0; i<spec_size; i++){
2441 double val;
2442
2443 check( val = cpl_vector_get( sp, i));
2444 check( cpl_matrix_set( decomp, 0, i, val));
2445 }
2446
2447 cleanup:
2448 if ( cpl_error_get_code() != CPL_ERROR_NONE){
2449 xsh_free_matrix( &decomp);
2450 }
2452 xsh_free_vector( &sp);
2454 XSH_FREE( new_filter);
2455 return decomp;
2456}
2457/*---------------------------------------------------------------------------*/
2458
2459
2460/*---------------------------------------------------------------------------*/
2467cpl_frameset * xsh_shift_offsettab( cpl_frameset *shiftifu_frameset,
2468 double offset_low, double offset_up)
2469{
2470 cpl_frameset *result = NULL;
2471 cpl_frame *lo_frame = NULL;
2472 const char *lo_name = NULL;
2473 cpl_table *lo_table = NULL;
2474 double *lo_data = NULL;
2475 cpl_frame *up_frame = NULL;
2476 const char *up_name = NULL;
2477 cpl_table *up_table = NULL;
2478 double *up_data = NULL;
2479 int i, size;
2480 cpl_propertylist *lo_header = NULL;
2481 cpl_propertylist *up_header = NULL;
2482 cpl_frame *prod_frame = NULL;
2483
2484 lo_frame = cpl_frameset_get_frame(shiftifu_frameset,0);
2485 lo_name = cpl_frame_get_filename(lo_frame);
2486 xsh_msg("Name %s", lo_name);
2487 XSH_TABLE_LOAD( lo_table, lo_name);
2488 check( lo_data = cpl_table_get_data_double( lo_table,
2490
2491 up_frame = cpl_frameset_get_frame(shiftifu_frameset,2);
2492 up_name = cpl_frame_get_filename(up_frame);
2493 xsh_msg("Name %s", up_name);
2494 XSH_TABLE_LOAD( up_table, up_name);
2495 check( up_data = cpl_table_get_data_double( up_table,
2497
2498 size = cpl_table_get_nrow( lo_table);
2499
2500 for( i=0; i<size; i++){
2501 lo_data[i] += offset_low;
2502 up_data[i] += offset_up;
2503 }
2504
2505 lo_header = cpl_propertylist_load( lo_name,0);
2506 check( cpl_table_save( lo_table, lo_header, NULL, "tmp_OFFSET_TAB_LOW.fits", CPL_IO_DEFAULT));
2507 up_header = cpl_propertylist_load( up_name,0);
2508 check( cpl_table_save( up_table, up_header, NULL, "tmp_OFFSET_TAB_UP.fits", CPL_IO_DEFAULT));
2509
2510 result = cpl_frameset_new();
2511 check(prod_frame = xsh_frame_product( "tmp_OFFSET_TAB_LOW.fits",
2512 "OFFSET_TAB",
2513 CPL_FRAME_TYPE_TABLE,
2514 CPL_FRAME_GROUP_PRODUCT,
2515 CPL_FRAME_LEVEL_TEMPORARY));
2516 cpl_frameset_insert( result, prod_frame);
2517
2518 check(prod_frame = cpl_frame_duplicate( cpl_frameset_get_frame(shiftifu_frameset,1)));
2519 cpl_frameset_insert( result, prod_frame);
2520
2521 check(prod_frame = xsh_frame_product( "tmp_OFFSET_TAB_UP.fits",
2522 "OFFSET_TAB",
2523 CPL_FRAME_TYPE_TABLE,
2524 CPL_FRAME_GROUP_PRODUCT,
2525 CPL_FRAME_LEVEL_TEMPORARY));
2526 cpl_frameset_insert( result, prod_frame);
2527
2528 cleanup:
2529 if ( cpl_error_get_code() != CPL_ERROR_NONE){
2530 xsh_free_frameset( &result);
2531 }
2532 XSH_TABLE_FREE( lo_table);
2533 XSH_TABLE_FREE( up_table);
2534 xsh_free_propertylist( &lo_header);
2535 xsh_free_propertylist( &up_header);
2536 return result;
2537}
2538/*---------------------------------------------------------------------------*/
static xsh_instrument * instrument
int binx
int biny
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_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
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_ASSURE_NOT_NULL(pointer)
Definition: xsh_error.h:99
const char * xsh_instrument_arm_tostring(xsh_instrument *i)
Get the string associated with an arm.
const char * xsh_arm_tostring(XSH_ARM arm)
Get the string associated with an arm.
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_dbg_medium(...)
Definition: xsh_msg.h:44
#define xsh_msg_debug(...)
Print a debug message.
Definition: xsh_msg.h:99
#define xsh_msg_error(...)
Print an error message.
Definition: xsh_msg.h:62
#define xsh_msg(...)
Print a message on info level.
Definition: xsh_msg.h:121
const char * xsh_pfits_get_raw1catg(const cpl_propertylist *plist)
find out the ESO.RAW1.CATG
Definition: xsh_pfits.c:1548
cpl_error_code xsh_plist_set_extra_keys(cpl_propertylist *plist, const char *hduclas1, const char *hduclas2, const char *hduclas3, const char *scidata, const char *errdata, const char *qualdata, const int type)
set hdu keys
Definition: xsh_pfits.c:4250
void xsh_pfits_set_pcatg(cpl_propertylist *plist, const char *value)
Write the PCATG value.
Definition: xsh_pfits.c:1008
double xsh_pfits_get_cdelt3(const cpl_propertylist *plist)
find out the cdelt3
Definition: xsh_pfits.c:2235
int xsh_pfits_get_binx(const cpl_propertylist *plist)
find out the BINX value
Definition: xsh_pfits.c:289
void xsh_pfits_set_extname(cpl_propertylist *plist, const char *value)
Write the EXTNAME value.
Definition: xsh_pfits.c:979
double xsh_pfits_get_tel_targ_delta(const cpl_propertylist *plist)
Get the TEL TARG DELTA.
Definition: xsh_pfits.c:3459
const char * xsh_pfits_get_pcatg(const cpl_propertylist *plist)
find out the pcatg
Definition: xsh_pfits.c:1627
double xsh_pfits_get_dec(const cpl_propertylist *plist)
Get the Right Ascension.
Definition: xsh_pfits.c:3500
double xsh_pfits_get_tel_targ_alpha(const cpl_propertylist *plist)
Get the TEL TARG ALPHA.
Definition: xsh_pfits.c:3439
const char * xsh_pfits_get_obs_targ_name(const cpl_propertylist *plist)
find out the ESO.OBS.TARG.NAME
Definition: xsh_pfits.c:1569
int xsh_pfits_get_biny(const cpl_propertylist *plist)
find out the BINY value
Definition: xsh_pfits.c:306
int xsh_pfits_get_obs_id(cpl_propertylist *plist)
find out the OBS ID
Definition: xsh_pfits.c:3388
double xsh_pfits_get_ra(const cpl_propertylist *plist)
Get the Right Ascension.
Definition: xsh_pfits.c:3479
char * xsh_pfits_get_slit_value(const cpl_propertylist *plist, xsh_instrument *instrument)
find out the INS OPTIx NAME value (the width of the slit)
Definition: xsh_pfits.c:619
double xsh_pfits_get_posangle(const cpl_propertylist *plist)
find out the value of the CUMOFFSETX keyword in a header
Definition: xsh_pfits.c:172
double xsh_pfits_get_crval3(const cpl_propertylist *plist)
find out the crval3
Definition: xsh_pfits.c:1946
void xsh_unwrap_vector(cpl_vector **v)
Unwrap a vector and set the pointer to NULL.
Definition: xsh_utils.c:2345
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
void xsh_free_frame(cpl_frame **f)
Deallocate a frame and set the pointer to NULL.
Definition: xsh_utils.c:2269
void xsh_free_frameset(cpl_frameset **f)
Deallocate a frame set and set the pointer to NULL.
Definition: xsh_utils.c:2254
double xsh_hms2deg(const double hms)
Convert a double from hours minute seconds to deg:
Definition: xsh_utils.c:312
double xsh_sess2deg(const double sess)
Convert a double from ssessagesimal to deg: 203049.197= 20:30:49.197 = 20.5136658333.
Definition: xsh_utils.c:345
void xsh_free_matrix(cpl_matrix **m)
Deallocate a matrix and set the pointer to NULL.
Definition: xsh_utils.c:2209
void xsh_free_table(cpl_table **t)
Deallocate a table and set the pointer to NULL.
Definition: xsh_utils.c:2133
void xsh_free_propertylist(cpl_propertylist **p)
Deallocate a property list and set the pointer to NULL.
Definition: xsh_utils.c:2179
void xsh_free_imagelist(cpl_imagelist **i)
Deallocate an image list and set the pointer to NULL.
Definition: xsh_utils.c:2164
void xsh_add_temporary_file(const char *name)
Add temporary file to temprary files list.
Definition: xsh_utils.c:1432
cpl_kernel kernel_type
#define QFLAG_MISSING_DATA
@ XSH_ARM_UNDEFINED
@ XSH_ARM_UVB
@ XSH_ARM_NIR
@ XSH_ARM_VIS
#define XSH_PRE_DATA_TYPE
Definition: xsh_data_pre.h:42
#define XSH_PRE_QUAL_TYPE
Definition: xsh_data_pre.h:46
#define XSH_PRE_ERRS_TYPE
Definition: xsh_data_pre.h:44
#define XSH_PRE_DATA_BPP
Definition: xsh_data_pre.h:43
#define XSH_SPECTRUM_DATA_BPP
#define XSH_SPECTRUM_QUAL_BPP
#define XSH_SPECTRUM_ERRS_BPP
int nx
int lly
Definition: xsh_detmon_lg.c:86
int llx
Definition: xsh_detmon_lg.c:85
const char * method
Definition: xsh_detmon_lg.c:78
int urx
Definition: xsh_detmon_lg.c:87
int ury
Definition: xsh_detmon_lg.c:88
int ny
int filter
void xsh_add_product_imagelist(cpl_frame *frame, cpl_frameset *frameset, const cpl_parameterlist *parameters, const char *recipe_id, xsh_instrument *instrument, const char *final_prefix)
Definition: xsh_dfs.c:2825
void xsh_add_product_table(cpl_frame *frame, cpl_frameset *frameset, const cpl_parameterlist *parameters, const char *recipe_id, xsh_instrument *instrument, const char *final_prefix)
Save Table product (input frame has several extensions, 1 table per extension)
Definition: xsh_dfs.c:3146
cpl_frame * xsh_frame_product(const char *fname, const char *tag, cpl_frame_type type, cpl_frame_group group, cpl_frame_level level)
Creates a frame with given characteristics.
Definition: xsh_dfs.c:930
#define XSH_IFU_MAP_SKY
Definition: xsh_dfs.h:58
#define XSH_SHIFTIFU_COLNAME_SHIFTSLIT
Definition: xsh_drl.h:518
void smooth(double vec[], long n, int w, double svec[])
int xsh_parameters_cut_uvb_spectrum_get(const char *recipe_id, const cpl_parameterlist *list)
#define XSH_QC_TRACE_FIT_C0
Definition: xsh_pfits_qc.h:53
#define XSH_QC_FLUX
Definition: xsh_pfits_qc.h:112
#define XSH_QC_TRACE_FIT_C2
Definition: xsh_pfits_qc.h:55
#define XSH_QC_TRACE_FIT_C1
Definition: xsh_pfits_qc.h:54
#define XSH_FREE(POINTER)
Definition: xsh_utils.h:92
#define XSH_CALLOC(POINTER, TYPE, SIZE)
Definition: xsh_utils.h:56
static void xsh_plist_set_cd_matrix3(cpl_propertylist **plist, const double cd1_3, const double cd2_3, const double cd3_1, const double cd3_2, const double cd3_3)
set world coordinate system: CD matrix
static cpl_frame * xsh_frame_build_sky_area(cpl_frame *slitmap_frame, const char *tag)
static cpl_error_code xsh_monitor_spectrum3D_flux(cpl_imagelist *data, cpl_imagelist *errs, cpl_imagelist *qual, xsh_instrument *instrument, cpl_propertylist *phead)
Computes flux, err, SNR on IFU cube over full range and sub intervals.
cpl_error_code xsh_cube_set_wcs(cpl_propertylist *plist, float cenLambda, float disp_x, float disp_y, float disp_z, float center_x, float center_y, int center_z)
static double xsh_interpol(double x, double x1, double x2, double y1, double y2)
Definition: xsh_utils_ifu.c:49
static cpl_error_code xsh_add_correct_coeff(cpl_table **table, cpl_propertylist *plist, const char *prefix, const int row, const char *col12, const char *col32)
cpl_error_code xsh_table_edges_swap_low_upp(cpl_table **tab)
swap lower and upper edge columns
static cpl_error_code xsh_monitor_spectrum3D_flux_qc(cpl_imagelist *data, cpl_imagelist *errs, cpl_imagelist *qual, cpl_propertylist *phead, const double wave_s, const double wave_e, const int index)
cpl_error_code xsh_build_ifu_cube(cpl_frame *div_frame, cpl_frame *ifu_cfg_tab_frame, cpl_frame *ifu_cfg_cor_frame, cpl_frame *spectral_format_frame, cpl_frame *model_config_frame, cpl_frame *wavesol_frame, xsh_instrument *instrument, cpl_frameset *frameset, cpl_parameterlist *parameters, xsh_rectify_param *rectify_par, const char *recipe_id, const char *rec_prefix, const int frame_is_object)
Reconstruct IFU cube.
void xsh_edge_check(const int px, const int nx, const int rad_x, int *llx, int *urx)
check edge values (llx,urx) fit in image size (1,nx)
cpl_matrix * xsh_atrous(cpl_vector *spec, int nscales)
Do a wavelet decomposition using atrous from IDL.
static void xsh_plist_set_coord3(cpl_propertylist **plist, const int crpix3, const double crval3, const double cdelt3)
set world coordinate system
static void xsh_plist_set_coord2(cpl_propertylist **plist, const double crpix2, const double crval2, const double cdelt2)
set world coordinate system
static void xsh_plist_set_cd_matrix2(cpl_propertylist **plist, const double cd1_1, const double cd1_2, const double cd2_1, const double cd2_2)
set world coordinate system: CD matrix
cpl_table * xsh_table_edge_prepare(const char *name)
Add to edges table columns to store LOW/CEN/UPP object positions.
cpl_frameset * xsh_shift_offsettab(cpl_frameset *shiftifu_frameset, double offset_low, double offset_up)
Do a wavelet decomposition using atrous from IDL.
cpl_error_code xsh_frame_check_model_cfg_is_proper_for_sci(cpl_frame *model_config_frame, cpl_frame *sci_frame, xsh_instrument *instrument)
Check if a model configuration frame has been corrected for flexures.
static void xsh_plist_set_coord1(cpl_propertylist **plist, const double crpix1, const double crval1, const double cdelt1)
set world coordinate system
Definition: xsh_utils_ifu.c:99
#define PI_NUMB
Definition: xsh_utils_ifu.c:46
cpl_error_code xsh_frame_check_model_cfg_is_afc_corrected(cpl_frame *model_config_frame)
utility to check if a frame has been corrected for flexures
cpl_frame * xsh_build_ifu_map(cpl_frame *div_frame, cpl_frame *wavemap_frame, cpl_frame *slitmap_frame, xsh_instrument *instrument)
Reconstruct IFU cube.
static cpl_frame * xsh_frame_build_sky_map(cpl_frame *slitmap_frame, const double value, const char *tag)
cpl_error_code xsh_ifu_trace_object_calibrate(const char *ifu_object_ff_name, const char *order_tab_edges_ifu_name, const char *slit_map_name, const char *wave_map_name)
Function to calibrate object traces in IFU mode.
void xsh_convert_xy_to_ws(double x_centroid, double *p_obj_cen, double *p_slit, double *p_wave, const int lly, const int nx, const int row, double *p_obj_cen_s, double *p_obj_cen_w)
convert X-Y coordinates to Wavelengt Slit space
static cpl_frame * xsh_crea_correct_coeff(cpl_frame *qc_trace_merged_frame, xsh_instrument *inst)
double xsh_image_fit_gaussian_max_pos_x_window(const cpl_image *ima, const int llx, const int urx, const int ypos)
cpl_error_code xsh_frame_image_save2ext(cpl_frame *frm, const char *name_o, const int ext_i, const int ext_o)
cpl_frame * xsh_cube_qc_trace_window(cpl_frame *frm_cube, xsh_instrument *instrument, const char *suffix, const char *rec_prefix, const int win_min, const int win_max, const int hsize, const int method, const int compute_qc)
Trace object position in a cube.
cpl_error_code xsh_iml_merge_wgt(cpl_imagelist **data, cpl_imagelist **errs, cpl_imagelist **qual, const cpl_image *flux_b, const cpl_image *errs_b, const cpl_image *qual_b, const int mk, const int decode_bp)
merge imagelist via average
cpl_error_code xsh_imagelist_cut_dichroic_uvb(cpl_imagelist *org_data_iml, cpl_imagelist *org_errs_iml, cpl_imagelist *org_qual_iml, cpl_propertylist *phead)
#define XSH_TABLE_LOAD(TABLE, NAME)
#define XSH_TABLE_FREE(TABLE)