X-shooter Pipeline Reference Manual 3.8.15
xsh_localize_obj.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 * $Author: amodigli $
22 * $Date: 2012-12-18 14:15:45 $
23 * $Revision: 1.64 $
24 */
25
26#ifdef HAVE_CONFIG_H
27#include <config.h>
28#endif
29
30#define CPL_MODE 0
31/*----------------------------------------------------------------------------*/
38/*----------------------------------------------------------------------------*/
41/*-----------------------------------------------------------------------------
42 Includes
43 -----------------------------------------------------------------------------*/
44
45#include <math.h>
46#include <xsh_drl.h>
47
48#include <xsh_utils_table.h>
49#include <xsh_badpixelmap.h>
50#include <xsh_utils_wrappers.h>
51#include <xsh_data_pre.h>
52#include <xsh_dfs.h>
53#include <xsh_pfits.h>
54#include <xsh_error.h>
55#include <xsh_msg.h>
56#include <xsh_fit.h>
57#include <xsh_data_instrument.h>
59#include <xsh_data_spectrum.h>
60#include <xsh_data_rec.h>
61#include <xsh_utils_ifu.h>
62#include <xsh_ifu_defs.h>
63#include <xsh_parameters.h>
64#include <cpl.h>
65
66/*-----------------------------------------------------------------------------
67 Functions prototypes
68 -----------------------------------------------------------------------------*/
69static void chunk_coadd( double* data, double* flux, int* qual, int nx, int ny,
70 int ifirst, int ilast, const int decode_bp, int *skymask);
71
72static void tab_minmax_get_index( double* data, int size,
73 xsh_slit_limit_param * slit_limit_par,
74 int* min, int* max);
75
76static xsh_localization* xsh_localize_obj_auto( cpl_frame * rec_frame,
77 cpl_frame *skymask_frame,
79 xsh_slit_limit_param * slit_limit_par,const int decode_bp);
80
83/*****************************************************************************/
84
85
86/* not used
87static void mask_badpixels( double *errs, int *qual, int nx, int ny)
88{
89 int i,j;
90
91 XSH_ASSURE_NOT_NULL( errs);
92 XSH_ASSURE_NOT_NULL( qual);
93
94 for( j=0 ; j<ny ; j++) {
95 for(i=0; i<nx; i++){
96 if ( errs[i+j*nx] > 200){
97 qual[i+j*nx] = XSH_BAD_PIXEL;
98 }
99 }
100 }
101
102 cleanup:
103 return;
104}
105*/
106static void chunk_coadd( double* data, double *flux, int* qual, int nx, int ny,
107 int ifirst, int ilast, const int decode_bp,int *skymask)
108{
109 int i,j;
110 int* nbgoodpixels = NULL;
111
115 XSH_ASSURE_NOT_NULL(skymask);
116 XSH_ASSURE_NOT_ILLEGAL(ifirst >= 0 && ifirst <= ilast && ilast < nx);
117
118 XSH_CALLOC( nbgoodpixels, int, ny);
119 /* initialise data */
120 for(j=0; j<ny; j++){
121 data[j] = 0;
122 }
123 /* Note that here mask is a vector sampling wavelength direction,
124 * used to mask occurrence of sky lines on the re-sampled frame.
125 * For this reason we loop first along X (i) then on Y (j)
126 * */
127 int j_nx=0;
128 for(i=ifirst; i<=ilast; i++){
129 if (skymask[i] == 0){
130 for( j=0 ; j<ny ; j++) {
131 j_nx=i+j*nx;
132 if ( (qual[j_nx] & decode_bp) == 0 ){
133 nbgoodpixels[j]++;
134 data[j] +=flux[j_nx];
135 }
136 /*
137 else{
138 xsh_msg_dbg_high("bad pixel at %d %d",i,j);
139 }
140 */
141 }
142 }
143 else{
144 xsh_msg_dbg_medium("Mask by sky lines : column %d", i);
145 }
146 }
147
148 for(j=0; j<ny; j++){
149 if ( nbgoodpixels[j] != 0 ) {
150 data[j] = data[j] / (float)(nbgoodpixels[j]);
151 }
152 else data[j] = 0. ;
153 }
154
155 cleanup:
156 XSH_FREE(nbgoodpixels);
157 return;
158}
159
160static void tab_minmax_get_index( double* data, int size,
161 xsh_slit_limit_param * slit_limit_par,
162 int* min, int* max)
163{
164 int i, maxi, mini ;
165 int slit0, slit1 ;
166
170
171 if ( slit_limit_par == NULL ) {
172 slit0 = 0 ;
173 slit1 = size ;
174 }
175 else {
176 slit0 = slit_limit_par->min_slit_idx + 1 ;
177 slit1 = slit_limit_par->max_slit_idx ;
178 }
179 maxi = 0 ;
180 mini = 0 ;
181
182 for(i = slit0; i<slit1; i++) {
183 if (data[i]>data[maxi]){
184 maxi = i;
185 }
186 else if (data[i]<data[mini]){
187 mini = i;
188 }
189 }
190 *min = mini;
191 *max = maxi;
192
193 cleanup:
194 return;
195}
196
197/*-----------------------------------------------------------------------------
198 Implementation
199 -----------------------------------------------------------------------------*/
200
201/*****************************************************************************/
219/*****************************************************************************/
220cpl_frame* xsh_localize_obj (cpl_frame *rec_frame,
221 cpl_frame *skymask_frame,
223 xsh_localize_obj_param * loc_obj_par,
224 xsh_slit_limit_param * slitlimit_par,
225 const char * fname)
226{
227 cpl_frame *merge_frame = NULL;
228 cpl_frame *res_frame = NULL;
229 xsh_localization *loc_list = NULL;
230 int merge_par = 0;
231 const char* rec_prefix = "LOCALIZE";
232 char filename[256];
233
234 /* Check parameters */
235
236 XSH_ASSURE_NOT_NULL( loc_obj_par);
238
239 xsh_msg_dbg_medium( "Entering xsh_localize_obj") ;
240
241 if ( rec_frame != NULL){
242 const char* filename = NULL;
243 check( merge_frame = xsh_merge_ord( rec_frame, instrument, merge_par,
244 rec_prefix));
245 check( filename = cpl_frame_get_filename( merge_frame));
246 check( cpl_frame_set_level( merge_frame, CPL_FRAME_LEVEL_TEMPORARY));
247 xsh_add_temporary_file( filename);
248 }
249
250 /* parameters */
251 xsh_msg_dbg_medium("method %s",LOCALIZE_METHOD_PRINT( loc_obj_par->method));
252
253 if ( loc_obj_par->method == LOC_MANUAL_METHOD){
254 xsh_msg_dbg_medium("slit position %f slit half height %f",
255 loc_obj_par->slit_position, loc_obj_par->slit_hheight);
256 check( loc_list = xsh_localize_obj_manual( instrument, loc_obj_par));
257 }
258 else{
259 check( loc_list = xsh_localize_obj_auto( merge_frame, skymask_frame,
261 loc_obj_par, slitlimit_par,instrument->decode_bp));
262 }
263
264 xsh_msg_dbg_low( "Saving Localization Table" ) ;
265 if ( fname == NULL ) {
266 /*
267 fname = "LOCALIZATION_TABLE_" ||
268 xsh_instrument_arm_tostring(instrument) || ".fits";
269 */
270 sprintf(filename,"LOCALIZATION_TABLE_%s.fits",
272 xsh_add_temporary_file(filename);
273 } else {
274 sprintf(filename,fname);
275 }
276 check( res_frame = xsh_localization_save( loc_list, filename, instrument));
277
278 cleanup:
279 xsh_localization_free( &loc_list);
280 xsh_free_frame( &merge_frame);
281
282 return res_frame ;
283}
284/*****************************************************************************/
285
288{
289 xsh_localization *loc_list = NULL;
290 int deg_poly = 0;
291 double slit_cen, slit_up, slit_low;
292 cpl_size pows = 0;
293
294 /* Check parameters */
295 XSH_ASSURE_NOT_NULL( loc_obj_par);
297
298 /* Parameters */
299 xsh_msg_dbg_medium("slit position %f slit half height %f",
300 loc_obj_par->slit_position, loc_obj_par->slit_hheight);
301
302 check( loc_list = xsh_localization_create());
303
304 slit_cen = loc_obj_par->slit_position;
305 slit_up = slit_cen + loc_obj_par->slit_hheight;
306 slit_low = slit_cen - loc_obj_par->slit_hheight;
307
308 loc_list->pol_degree = deg_poly;
309 check( loc_list->cenpoly = cpl_polynomial_new( 1));
310 check( loc_list->edglopoly = cpl_polynomial_new( 1));
311 check( loc_list->edguppoly = cpl_polynomial_new( 1));
312 check( cpl_polynomial_set_coeff( loc_list->cenpoly,
313 &pows, slit_cen));
314 check( cpl_polynomial_set_coeff( loc_list->edglopoly,
315 &pows, slit_low));
316 check( cpl_polynomial_set_coeff( loc_list->edguppoly,
317 &pows, slit_up));
318
319 cleanup:
320 return loc_list;
321}
322
323/*****************************************************************************/
341/*****************************************************************************/
342static xsh_localization*
343xsh_localize_obj_auto( cpl_frame *merge_frame,
344 cpl_frame *skymask_frame,
346 xsh_localize_obj_param *loc_obj_par,
347 xsh_slit_limit_param *slitlimit_par,const int decode_bp)
348{
349 xsh_localization *loc_list = NULL;
350 xsh_spectrum *spectrum = NULL;
351 double* slit_center = NULL;
352 double* slit_upper = NULL;
353 double* slit_lower = NULL;
354 double* chunk_center = NULL;
355 cpl_vector* slit_center_v = NULL;
356 cpl_vector* slit_upper_v = NULL;
357 cpl_vector* slit_lower_v = NULL;
358 cpl_vector* chunk_center_v = NULL;
359
360#if 0
361 int i = 0;
362 int skip_order_nb = 0;
363 int first_order_index = 0;
364 int slit_nod, slit_lim ;
365#endif
366
367 int level;
368
369 int nslit, nlambda, chunk_size;
370 double* slit_tab = NULL;
371 double *flux = NULL;
372 //double *errs = NULL;
373 int *qual = NULL;
374 int ichunk, skip_chunk=0, nb_chunk;
375 cpl_polynomial *cen_poly = NULL;
376 int loc_degree;
377 double lambda_min, lambda_step;
378 double slit_min, slit_step, kappa=2.0;
379 int islit, iter, niter=5, nbrej=-1;
380 cpl_vector *dispslit = NULL;
381 cpl_vector *gausspos = NULL;
382 cpl_vector *gaussval = NULL;
383 int *sky_mask = NULL;
384 cpl_table *skymask_table = NULL;
385 const char* skymask_name = NULL;
386
387 /* Check parameters */
388 XSH_ASSURE_NOT_NULL( loc_obj_par);
390 XSH_ASSURE_NOT_NULL( merge_frame);
391
392 /* parameters */
393 kappa = loc_obj_par->kappa;
394 niter = loc_obj_par->niter;
395 loc_degree = loc_obj_par->loc_deg_poly;
396 xsh_msg_dbg_medium( "Entering xsh_localize_obj_auto") ;
397 xsh_msg_dbg_medium("Localize deg_poly %d chunk %d Thresh %f kappa %f niter %d",
398 loc_degree, loc_obj_par->loc_chunk_nb, loc_obj_par->loc_thresh,
399 kappa, niter);
400 if ( loc_obj_par->use_skymask){
401 XSH_ASSURE_NOT_NULL( skymask_frame);
402 check( skymask_name = cpl_frame_get_filename( skymask_frame));
403 xsh_msg_dbg_medium("Sky mask %s", skymask_name);
404 }
405
406 level = xsh_debug_level_get();
407
408 /* Load merge spectrum 2D */
409 check( spectrum = xsh_spectrum_load( merge_frame));
410 check( nslit = xsh_spectrum_get_size_slit( spectrum));
411 check( nlambda = xsh_spectrum_get_size_lambda( spectrum));
412 lambda_min = spectrum->lambda_min;
413 lambda_step = spectrum->lambda_step;
414 slit_min = spectrum->slit_min;
415 slit_step = spectrum->slit_step;
416
417 check( flux = cpl_image_get_data_double( spectrum->flux));
418 //check( errs = cpl_image_get_data_double( spectrum->errs));
419 check( qual = cpl_image_get_data_int( spectrum->qual));
420
421 /* allocate memory */
422 XSH_CALLOC( chunk_center, double, loc_obj_par->loc_chunk_nb);
423 XSH_CALLOC( slit_center, double, loc_obj_par->loc_chunk_nb);
424 XSH_CALLOC( slit_upper, double, loc_obj_par->loc_chunk_nb);
425 XSH_CALLOC( slit_lower, double, loc_obj_par->loc_chunk_nb);
426
427 /* create sky mask */
428 XSH_CALLOC( sky_mask, int, nlambda);
429
430 if ( loc_obj_par->use_skymask){
431 float *skymask_data = NULL;
432 int irow, nrow;
433 double fwhm =0.0, sky_min, sky_max;
434 int isky_min, isky_max, imask;
435 double width, resolution;
436
437 XSH_TABLE_LOAD( skymask_table, skymask_name);
438 /* sort wavelength */
439 check( xsh_sort_table_1( skymask_table, "WAVELENGTH", CPL_FALSE));
440 check( skymask_data = cpl_table_get_data_float( skymask_table,
441 "WAVELENGTH"));
442 check( nrow = cpl_table_get_nrow( skymask_table));
443
444 xsh_msg_dbg_low("lambda min %f, step %f", lambda_min, lambda_step);
445
446 for( irow=0; irow < nrow; irow++){
448 resolution = xsh_resolution_get( instrument, width);
449 fwhm = skymask_data[irow]/resolution;
450 sky_min = skymask_data[irow]-fwhm;
451 sky_max = skymask_data[irow]+fwhm;
452 isky_min = (int)xsh_round_double((sky_min-lambda_min)/lambda_step);
453 isky_max = (int)xsh_round_double((sky_max-lambda_min)/lambda_step);
454 for( imask=isky_min; imask <=isky_max; imask++){
455 sky_mask[imask] = 1;
456 }
457 }
458 if (level >= XSH_DEBUG_LEVEL_MEDIUM){
459 FILE *mask_file = NULL;
460 char mask_name[256];
461 int idbg=0;
462
463 sprintf( mask_name, "skymask.reg");
464 mask_file = fopen( mask_name, "w");
465
466 fprintf( mask_file,"# Region file format: DS9 version 4.1\n");
467 fprintf( mask_file,"global color=green dashlist=8 3 width=1 font=\"helvetica 10 normal\"\
468 select=1 highlite=1 dash=0 fixed=0 edit=1 move=1 delete=1 include=1 source=1\n");
469 fprintf(mask_file,"physical\n");
470
471 for(idbg=0; idbg< nlambda; idbg++){
472
473 if (sky_mask[idbg] == 1){
474 fprintf( mask_file, "line(%d,%d,%d,%d)\n",idbg+1, nslit, idbg+1, 1);
475 }
476 }
477 fclose( mask_file);
478 }
479 }
480 /* Create the result */
481 check( loc_list = xsh_localization_create());
482
483#if 0 /* REGDEBUG COMMENT for now don't what is it ?? */
484 /* Calculate slit_nod from arcsec to index */
485 {
486 float arcsec ;
488 case XSH_ARM_NIR:
489 arcsec = XSH_ARCSEC_NIR ;
490 break ;
491 case XSH_ARM_UVB:
492 arcsec = XSH_ARCSEC_UVB ;
493 break ;
494 case XSH_ARM_VIS:
495 arcsec = XSH_ARCSEC_VIS ;
496 break ;
497 default:
498 arcsec = 0.15 ;
499 break ;
500 }
501 slit_nod = loc_obj_par->nod_step/arcsec ;
502 xsh_msg_dbg_medium( "nod_step: %lf, slit_nod: %d",
503 loc_obj_par->nod_step, slit_nod ) ;
504 }
505#endif
506 XSH_CALLOC( slit_tab, double, nslit);
507 check( gaussval = cpl_vector_wrap( nslit, slit_tab));
508 check( gausspos = cpl_vector_new( nslit));
509 for( islit=0; islit < nslit; islit++){
510 cpl_vector_set( gausspos, islit, islit);
511 }
512
513 chunk_size = (nlambda-1) / (double)(loc_obj_par->loc_chunk_nb);
514
515 /* mask the bad pixel */
516 //check( mask_badpixels( errs, qual, nlambda, nslit));
517
518 for( ichunk=0; ichunk< loc_obj_par->loc_chunk_nb; ichunk++){
519 /*
520 If in nodding mode, yup and ylow must be relative to ymax:
521 ylow = ymax - nod_step/2
522 yup = ymax + nod_step/2
523 */
524 int ifirst=0,ilast=0,icenter=0;
525 double slit_cen=0, slit_low=0, slit_up=0;
526 double threshold = 0.0;
527 int is_valid_chunk = 1;
528
529 ifirst = ichunk*chunk_size;
530 ilast = (ichunk+1)*chunk_size;
531 icenter = (int)((ifirst+ilast)/2);
532
533 xsh_msg_dbg_medium("%d-%d", ifirst,ilast);
534
535 /* Sum by chunk the good pixels of flux image and normalize by */
536 chunk_coadd( slit_tab, flux, qual, nlambda, nslit, ifirst, ilast, decode_bp, sky_mask);
537
538 if (level >= XSH_DEBUG_LEVEL_MEDIUM){
539 FILE *coadd_file = NULL;
540 char coadd_name[256];
541 int icoadd=0;
542 const char* filename = cpl_frame_get_filename( merge_frame);
543
544 sprintf( coadd_name, "%s_%d_%d.dat",filename,ifirst,ilast);
545 coadd_file = fopen( coadd_name, "w");
546
547 for(icoadd=0; icoadd<nslit; icoadd++){
548 fprintf(coadd_file, "%d %f\n",icoadd,slit_tab[icoadd]);
549 }
550 fclose( coadd_file);
551 }
552
553 if ( loc_obj_par->method == LOC_GAUSSIAN_METHOD){
554 double cenpos=0, sigma=0,area=0, offset=0;
555 double kappa=3.0;
556
557 xsh_msg_dbg_medium("Using GAUSSIAN method to fit chunk %d_%d", ifirst,ilast);
558 cpl_vector_fit_gaussian( gausspos, NULL, gaussval, NULL, CPL_FIT_ALL, &cenpos,
559 &sigma,
560 &area, &offset, NULL,
561 NULL, NULL);
562 if( cpl_error_get_code() == CPL_ERROR_CONTINUE ){
563 xsh_msg_dbg_low("CONTINUE to fit %d_%d x0 %f sigma %f", ifirst,ilast,cenpos,sigma);
564 }
565 if( cpl_error_get_code() == CPL_ERROR_NONE ){
566 slit_cen = cenpos;
567 slit_low = cenpos-kappa*sigma;
568 slit_up = cenpos+kappa*sigma;
569 if (slit_low < 0) slit_low = 0;
570 if (slit_up >= nslit) slit_up = nslit-1;
571 }
572 else{
573 xsh_msg_dbg_low("FAILED to fit %d_%d", ifirst,ilast);
574 is_valid_chunk=0;
576 }
577 }
578 else
579 {
580 int ymax=0, ymin=0;
581 int ylow=0, yup=0;
582 tab_minmax_get_index( slit_tab, nslit, slitlimit_par, &ymin, &ymax);
583 if (ymin == ymax){
584 xsh_msg_warning( "No maximum found in chunk (skip it)");
585 is_valid_chunk=0;
586 }
587 else{
588
589
590 xsh_msg_dbg_medium("ymin [%d] = %f ymax [%d] = %f",ymin,
591 slit_tab[ymin], ymax, slit_tab[ymax]);
592
593 threshold = slit_tab[ymin]+(slit_tab[ymax]-slit_tab[ymin])*
594 loc_obj_par->loc_thresh;
595 xsh_msg_dbg_medium("Threshold at %f",threshold);
596
597 /* search upper edge */
598 yup = ymax+1;
599 while( yup < nslit && (slit_tab[yup] >= threshold)) {
600 yup++;
601 }
602
603 /* search lower edge */
604 ylow = ymax-1;
605 while( ylow > 0 && (slit_tab[ylow] >= threshold)) {
606 ylow--;
607 }
608
609 slit_cen = ymax;
610 slit_low = ylow;
611 slit_up = yup;
612 }
613 }
614 /* put result in tab */
615 if ( is_valid_chunk == 1){
616 chunk_center[ichunk-skip_chunk] = lambda_min+lambda_step*icenter;
617 slit_center[ichunk-skip_chunk] = slit_min+slit_step*slit_cen;
618 slit_lower[ichunk-skip_chunk] = slit_min+slit_step*slit_low;
619 slit_upper[ichunk-skip_chunk] = slit_min+slit_step*slit_up;
620 }
621 else{
622 skip_chunk++;
623 }
624 }
625
626 nb_chunk = loc_obj_par->loc_chunk_nb-skip_chunk;
627
628 XSH_ASSURE_NOT_ILLEGAL_MSG(loc_obj_par->loc_deg_poly < nb_chunk,
629 "To fit a polynomial of deg n (loc_obj_par->loc_deg_poly) "\
630 "you need at least n+1 points. This condition has failed "\
631 "after clipping out other chunks not matching predefined criteria. "\
632 "This can be due to the data (not well defined order trace) and/or "\
633 "to a too high degree (localize-deg-lambda) of the polynomial fit");
634
635 check( chunk_center_v = cpl_vector_wrap(
636 nb_chunk, chunk_center));
637 check( slit_center_v = cpl_vector_wrap(
638 nb_chunk, slit_center));
639
640 check( cen_poly =
641 xsh_polynomial_fit_1d_create( chunk_center_v, slit_center_v,
642 loc_degree, NULL));
643
644 /* sigma clipping */
645 if ( niter > 0){
646 xsh_msg_dbg_low("Doing sigma clipping");
647 iter = 0;
648
649 while ( (loc_degree < (nb_chunk-nbrej)) && iter < niter
650 && nbrej != 0){
651 double sigma_med;
652
653 nbrej = 0;
654 xsh_msg_dbg_medium( " *** NBITER = %d / %d ***", iter+1, niter);
655 dispslit = cpl_vector_new( nb_chunk);
656
657 for(ichunk = 0; ichunk < nb_chunk; ichunk++){
658 double slit_fit;
659 double lambda;
660 double slit_diff;
661
662 lambda = chunk_center[ichunk];
663
664 check( slit_fit = cpl_polynomial_eval_1d( cen_poly,
665 lambda, NULL));
666
667 slit_diff = slit_center[ichunk]-slit_fit;
668 xsh_msg_dbg_low("slit_center %f FIT %f, DIFF %d %f", slit_center[ichunk], slit_fit, ichunk, slit_diff);
669 check( cpl_vector_set( dispslit, ichunk, slit_diff));
670 }
671
672 check( sigma_med = cpl_vector_get_stdev( dispslit));
673 xsh_msg_dbg_medium(" kappa %f SIGMA MEDIAN = %f", kappa, sigma_med);
674
675 for(ichunk = 0; ichunk < nb_chunk; ichunk++){
676 if ( fabs(cpl_vector_get( dispslit, ichunk)) > (kappa * sigma_med) ){
677 nbrej++;
678 }
679 else{
680 chunk_center[ichunk-nbrej] = chunk_center[ichunk];
681 slit_center[ichunk-nbrej] = slit_center[ichunk];
682 slit_lower[ichunk-nbrej] = slit_lower[ichunk];
683 slit_upper[ichunk-nbrej] = slit_upper[ichunk];
684 }
685 }
686
687 xsh_msg_dbg_medium(" Removed %d points", nbrej);
688
689 nb_chunk -= nbrej;
690
691 /* Redo the fit */
692 xsh_unwrap_vector( &chunk_center_v);
693 xsh_unwrap_vector( &slit_center_v);
694 xsh_free_polynomial( &cen_poly);
695
696 check( chunk_center_v = cpl_vector_wrap(
697 nb_chunk, chunk_center));
698 check( slit_center_v = cpl_vector_wrap(
699 nb_chunk, slit_center));
700
701 check( cen_poly =
702 xsh_polynomial_fit_1d_create( chunk_center_v, slit_center_v,
703 loc_degree, NULL));
704 xsh_free_vector( &dispslit);
705 iter++;
706 }
707 }
708
709 if (level >= XSH_DEBUG_LEVEL_MEDIUM){
710 FILE *debug_file = NULL;
711 char debug_name[256];
712 int idebug=0;
713 const char* filename = cpl_frame_get_filename( merge_frame);
714
715 sprintf( debug_name, "%s_points.dat",filename);
716 debug_file = fopen( debug_name, "w");
717
718 fprintf( debug_file,"#chunk_pos slit_low slit_cen slit_up\n");
719
720 for(idebug=0; idebug<nb_chunk; idebug++){
721 fprintf( debug_file, "%f %f %f %f\n",chunk_center[idebug],
722 slit_lower[idebug], slit_center[idebug], slit_upper[idebug]);
723 }
724 fclose( debug_file);
725 }
726
727 /* create edges polynomials */
728 check( slit_lower_v = cpl_vector_wrap(
729 nb_chunk, slit_lower));
730 check( slit_upper_v = cpl_vector_wrap(
731 nb_chunk, slit_upper));
732 loc_list->pol_degree = loc_degree;
733
734 check(loc_list->cenpoly = cpl_polynomial_duplicate( cen_poly));
735 check(loc_list->edglopoly =
736 xsh_polynomial_fit_1d_create( chunk_center_v, slit_lower_v,
737 loc_degree, NULL));
738 check(loc_list->edguppoly =
739 xsh_polynomial_fit_1d_create( chunk_center_v, slit_upper_v,
740 loc_degree, NULL));
741
742 cleanup:
743 xsh_unwrap_vector( &chunk_center_v);
744 xsh_unwrap_vector( &slit_center_v);
745 xsh_unwrap_vector( &slit_upper_v);
746 xsh_unwrap_vector( &slit_lower_v);
747 xsh_unwrap_vector( &gaussval);
748 xsh_free_vector( &gausspos);
749 XSH_FREE( chunk_center);
750 XSH_FREE( slit_center);
751 XSH_FREE( slit_upper);
752 XSH_FREE( slit_lower);
753 XSH_FREE( slit_tab);
754 XSH_FREE( sky_mask);
755 xsh_free_table( &skymask_table);
756 xsh_free_polynomial( &cen_poly);
757 xsh_spectrum_free( &spectrum);
758 return loc_list;
759
760}
761/*****************************************************************************/
762
763/*
764static int comp_lambda( const void * one, const void * two )
765{
766 float * un = (float *)one, * deux = (float *)two ;
767
768 if ( *un < *deux ) return -1 ;
769 else if ( *un == *deux ) return 0 ;
770 else return 1 ;
771}
772
773
774static cpl_frame * concat_rec( cpl_frame * rect_frame,
775 xsh_instrument * instrument,
776 int slitlet, int * nb_orders )
777{
778 cpl_frame * result = NULL ;
779 xsh_rec_list * res_list = NULL, * rect_list = NULL ;
780 int norders, i, ns, nlambda, depth ;
781 float * res_slit = NULL, * rect_slit = NULL ;
782 double * res_lambda = NULL, * rect_lambda = NULL, * plambda ;
783 const char * tag = "ORDER2D_DOWN_IFU_VIS" ;
784 char * fname = NULL ;
785 float * res_data = NULL, * res_errs = NULL ;
786 int * res_qual = NULL ;
787 float * pdata, * perrs, * pdata1, * perrs1 ;
788 int * pqual, * pqual1 ;
789 float * rect_data = NULL, * rect_errs = NULL ;
790 int * rect_qual = NULL ;
791
792 XSH_ASSURE_NOT_NULL( rect_frame ) ;
793
794 check( rect_list = xsh_rec_list_load( rect_frame, instrument ) ) ;
795 // Create the list
796 check( res_list = xsh_rec_list_create_with_size( 1, instrument ) ) ;
797 // Populate with rect_frame
798 ns = xsh_rec_list_get_nslit( rect_list, 0);
799 norders = rect_list->size ;
800 *nb_orders = norders ;
801 for( i = 0, nlambda = 0 ; i<norders ; i++ ) {
802 nlambda += xsh_rec_list_get_nlambda( rect_list, i);
803 xsh_msg_dbg_medium( "Nlambda[%d] = %d, sum = %d", ns,
804 rect_list->list[i].nlambda, nlambda);
805 }
806 xsh_msg_dbg_medium( "Single Rectify list: ns = %d, Nlambda = %d",
807 ns, nlambda ) ;
808 // Now fill the list
809 // Prepare the entry
810 check( xsh_rec_list_set_data_size( res_list, 0, 16, nlambda, ns ) ) ;
811 // Fill slit array
812 res_slit = xsh_rec_list_get_slit( res_list, 0 ) ;
813 rect_slit = xsh_rec_list_get_slit( rect_list, 0 ) ;
814 memcpy( res_slit, rect_slit, ns*sizeof( float ) ) ;
815
816 // Fill Lambd array
817 res_lambda = xsh_rec_list_get_lambda( res_list, 0 ) ;
818 plambda = res_lambda ;
819 for( i = 0 ; i<norders ; i++ ) {
820 int nl ;
821
822 rect_lambda = xsh_rec_list_get_lambda( rect_list, i ) ;
823 nl = xsh_rec_list_get_nlambda( rect_list, i) ;
824 memcpy( plambda, rect_lambda, nl*sizeof(double) ) ;;
825 plambda += nl ;
826 }
827
828 // Now fill data, errs and qual
829 depth = ns*nlambda ;
830 res_data = xsh_rec_list_get_data1( res_list, 0 ) ;
831 pdata = res_data ;
832 pdata1 = res_data ;
833 res_errs = xsh_rec_list_get_errs1( res_list, 0 ) ;
834 perrs = res_errs ;
835 perrs1 = res_errs ;
836 res_qual = xsh_rec_list_get_qual1( res_list, 0 ) ;
837 pqual = res_qual ;
838 pqual1 = pqual ;
839 for( i = 0 ; i<norders ; i++ ) {
840 int nl ;
841 float * rdata, * pdata2 ;
842 float * rerrs, * perrs2 ;
843 int * rqual, * pqual2 ;
844 int ks, kl ;
845
846 check( rect_data = xsh_rec_list_get_data1( rect_list, i ) ) ;
847 check( rect_errs = xsh_rec_list_get_errs1( rect_list, i ) ) ;
848 check( rect_qual = xsh_rec_list_get_qual1( rect_list, i ) ) ;
849 nl = xsh_rec_list_get_nlambda( rect_list, i) ;
850 xsh_msg_dbg_medium( "Order %d, nlambdas: %d [%d]", i, nl, nlambda ) ;
851 rdata = rect_data ;
852 pdata2 = pdata1 ;
853 rerrs = rect_errs ;
854 perrs2 = perrs1 ;
855 rqual = rect_qual ;
856 pqual2 = pqual1 ;
857 for ( ks = 0 ; ks < ns ; ks++ ) {
858 for ( kl = 0 ; kl < nl ; kl++ ) {
859 *pdata++ = *rdata++ ;
860 *perrs++ = *rerrs++ ;
861 *pqual++ = *rqual++ ;
862 }
863 pdata2 += nlambda ;
864 pdata = pdata2 ;
865 perrs2 += nlambda ;
866 perrs = perrs2 ;
867 pqual2 += nlambda ;
868 pqual = pqual2 ;
869 }
870 pdata1 += nl ;
871 pdata = pdata1 ;
872 perrs1 += nl ;
873 perrs = perrs1 ;
874 pqual1 += nl ;
875 pqual = pqual1 ;
876 }
877
878 // Sort by Lambdas (lambda array, data, errs and qual !)
879 {
880 int * idx ;
881 int ks ;
882
883 idx = xsh_sort( res_lambda, nlambda, sizeof( double ), comp_lambda ) ;
884 pdata = res_data ;
885 perrs = res_errs ;
886 pqual = res_qual ;
887 for( ks = 0 ; ks < ns ; ks++ ) {
888 xsh_reindex_float( pdata, idx, nlambda ) ;
889 pdata += nlambda ;
890 xsh_reindex_float( perrs, idx, nlambda ) ;
891 perrs += nlambda ;
892 xsh_reindex_int( pqual, idx, nlambda ) ;
893 pqual += nlambda ;
894 }
895 }
896
897 fname = xsh_stringcat_any( "SINGLE_RECTIFY_", SlitletName[slitlet],
898 ".fits", (void*)NULL ) ;
899 check( result = xsh_rec_list_save( res_list, fname, tag, 1 ) ) ;
900 XSH_FREE( fname ) ;
901 fname = xsh_stringcat_any( "MULTIPLE_RECTIFY_", SlitletName[slitlet],
902 ".fits", (void*)NULL ) ;
903 check( rect_frame = xsh_rec_list_save( rect_list,fname, tag, 1 ) ) ;
904 cleanup:
905 xsh_rec_list_free( &res_list ) ;
906 xsh_rec_list_free( &rect_list ) ;
907 return result ;
908}
909*/
910
911/*****************************************************************************/
912/*****************************************************************************/
913cpl_frameset * xsh_localize_obj_ifu( cpl_frameset *rec_frameset,
914 cpl_frame *skymask_frame,
916 xsh_localize_obj_param * locobj_par,
917 xsh_slit_limit_param *slitlimit_par)
918{
919 int i, slitlet;
920 cpl_frameset *result_frameset = NULL;
921 char fname[256];
922
923 XSH_ASSURE_NOT_NULL( rec_frameset);
925 XSH_ASSURE_NOT_NULL( locobj_par);
926
927 check( result_frameset = cpl_frameset_new());
928
929 for( i = 0, slitlet = LOWER_IFU_SLITLET ; i < 3 ; i++, slitlet++ ) {
930 cpl_frame * loc_frame = NULL;
931 cpl_frame *rec_frame = NULL;
932
933 sprintf( fname ,"LOCALIZE_TABLE_%s_IFU_%s.fits", SlitletName[slitlet],
935
936 xsh_msg( "Localizing slitlet %s, frame '%s'", SlitletName[slitlet], fname);
937
938 check( rec_frame = cpl_frameset_get_frame( rec_frameset, i));
939
940 check( loc_frame = xsh_localize_obj( rec_frame, skymask_frame, instrument,
941 locobj_par,slitlimit_par, fname));
942 check( cpl_frameset_insert( result_frameset, loc_frame));
943 }
944
945 cleanup:
946 if ( cpl_error_get_code() != CPL_ERROR_NONE ) {
947 xsh_free_frameset( &result_frameset);
948 }
949 return result_frameset;
950}
951/*--------------------------------------------------------------------------*/
952
953/*---------------------------------------------------------------------------*/
959double xsh_convert_seeing( cpl_frame* frame)
960{
961 double mu=-1.0;
962 double avg_seeing, seeing_start, seeing_end;
963 double avg_airmass;
964 double k;
965 const char* filename = NULL;
966 cpl_propertylist *header = NULL;
967
968 XSH_ASSURE_NOT_NULL( frame);
969
970 check( filename = cpl_frame_get_filename( frame));
971 check( header = cpl_propertylist_load( filename, 0));
972 check( avg_airmass = xsh_pfits_get_airm_mean( header));
973 check( seeing_start = xsh_pfits_get_seeing_start( header));
974 check( seeing_end = xsh_pfits_get_seeing_end( header));
975 avg_seeing = 0.5*(seeing_start+seeing_end);
976
977 k = sqrt( 1.0-78.08*pow(XSH_LAMBDA_DIMM,0.4)*
978 pow(avg_airmass,-0.6)*pow(avg_seeing,-1.0/3.0));
979
980 xsh_msg("K %f", k);
981
982 mu = 1.0/CPL_MATH_FWHM_SIG*avg_seeing*pow(XSH_LAMBDA_DIMM,0.2)*
983 pow(avg_airmass,0.6)*k;
984
985 xsh_msg("Mu %f", mu);
986
987 cleanup:
988 if ( cpl_error_get_code() != CPL_ERROR_NONE ) {
989 mu =-1;
990 }
991 xsh_free_propertylist( &header);
992 return mu;
993}
994/*--------------------------------------------------------------------------*/
995
996/*---------------------------------------------------------------------------*/
1018cpl_frame* xsh_localize_ifu_slitlet( cpl_frame *merge2d_slitlet,
1019 cpl_frame *skymask_frame, int smooth_hsize,
1020 int nscales, int HF_skip, const char* resname, double cut_sigma_low,
1021 double cut_sigma_up, double cut_snr_low, double cut_snr_up,
1022 double slit_min, double slit_max, int deg, int box_hsize,
1024
1025{
1026 cpl_frame *result = NULL;
1027 double wmin, wstep;
1028 //double wmax;
1029 double smin, sstep;
1030 //double smax;
1031 int wsize, ssize;
1032 xsh_spectrum *spectrum2d = NULL;
1033 double *flux = NULL;
1034 double *errs = NULL;
1035 int *qual = NULL;
1036 int i, j, k;
1037 int level;
1038 double *slit_vect_data = NULL;
1039 double *sliterr_vect_data = NULL;
1040 double *slit_pos_data = NULL;
1041 int slit_vect_size;
1042 cpl_vector *slit_vect = NULL;
1043 cpl_vector *sliterr_vect = NULL;
1044 cpl_vector *slit_pos = NULL;
1045 cpl_vector *smooth_vect = NULL;
1046 double *spos_data = NULL;
1047 double *errpos_data = NULL;
1048 double *wpos_data = NULL;
1049 double *sigma_data = NULL;
1050 double *snr_data = NULL;
1051 int data_size=0, ndata_size=0;
1052 cpl_vector *spos_vect = NULL;
1053 cpl_vector *wpos_vect = NULL;
1054 cpl_vector *sigma_vect = NULL;
1055 cpl_vector *snr_vect = NULL;
1056 cpl_vector *sigma_sort_vect = NULL;
1057 cpl_vector *snr_sort_vect = NULL;
1058 double *sposg_data = NULL;
1059 double *wposg_data = NULL;
1060 cpl_matrix *decomp = NULL;
1061 int nb_scales;
1062 cpl_table *table = NULL;
1063 char tablename[256];
1064 cpl_propertylist *header = NULL;
1065 double n05, n95, snr_05, snr_95, sigma_05, sigma_95;
1066 int jmin, jmax;
1067 /* sky mask */
1068 int *sky_mask = NULL;
1069 cpl_table *skymask_table = NULL;
1070 const char* skymask_name = NULL;
1071
1072 XSH_ASSURE_NOT_NULL( merge2d_slitlet);
1073 XSH_ASSURE_NOT_ILLEGAL( smooth_hsize >= 0);
1074 XSH_ASSURE_NOT_ILLEGAL( nscales > 0);
1075 XSH_ASSURE_NOT_ILLEGAL( nscales > HF_skip);
1076 XSH_ASSURE_NOT_ILLEGAL( box_hsize >= 0);
1077 /*
1078 xsh_msg_dbg_medium("PARAMS merge2d_slitlet %s",cpl_frame_get_filename(merge2d_slitlet));
1079 if (skymask_frame != NULL){
1080 xsh_msg_dbg_medium("PARAMS skymask_frame %s",cpl_frame_get_filename(skymask_frame));
1081 }
1082
1083 xsh_msg_dbg_medium("PARAMS smooth_hsize %d", smooth_hsize);
1084 xsh_msg_dbg_medium("PARAMS nscales %d", nscales);
1085 xsh_msg_dbg_medium("PARAMS HF_skip %d", HF_skip);
1086 xsh_msg_dbg_medium("PARAMS res_name %s", resname);
1087 xsh_msg_dbg_medium("PARAMS cut sigma %f %f", cut_sigma_low, cut_sigma_up);
1088 xsh_msg_dbg_medium("PARAMS cut snr %f %f", cut_snr_low, cut_snr_up);
1089 xsh_msg_dbg_medium("PARAMS slit %f %f", slit_min, slit_max);
1090 xsh_msg_dbg_medium("PARAMS deg %d", deg);
1091 */
1092 level = xsh_debug_level_get();
1093 check( spectrum2d = xsh_spectrum_load( merge2d_slitlet));
1094 check( wmin = xsh_spectrum_get_lambda_min( spectrum2d));
1095 //check( wmax = xsh_spectrum_get_lambda_max( spectrum2d));
1096 check( wstep = xsh_spectrum_get_lambda_step( spectrum2d));
1097
1098if ( skymask_frame != NULL){
1099 check( skymask_name = cpl_frame_get_filename( skymask_frame));
1100 xsh_msg_dbg_medium("Sky mask %s", skymask_name);
1101 }
1102 wsize = spectrum2d->size_lambda;
1103
1104
1105 smin = spectrum2d->slit_min;
1106 //smax = spectrum2d->slit_max;
1107 sstep = spectrum2d->slit_step;
1108 ssize = spectrum2d->size_slit;
1109
1110 jmin = 0;
1111 jmax = ssize;
1112
1113
1114 while( (smin+jmin*sstep) < slit_min){
1115 jmin++;
1116 }
1117 while( (smin+jmax*sstep) > slit_max){
1118 jmax--;
1119 }
1120
1121 /* in case of wrong value put default */
1122 if ( jmin > jmax || jmax <= 0 || jmin >= ssize){
1123 jmin = 0;
1124 jmax = ssize;
1125 }
1126 xsh_msg_dbg_medium( "Use [%d-%d] from [%d %d] slitlet",
1127 jmin, jmax, 0, ssize);
1128
1129 check( flux = xsh_spectrum_get_flux( spectrum2d));
1130 check( errs = xsh_spectrum_get_errs( spectrum2d));
1131 check( qual = xsh_spectrum_get_qual( spectrum2d));
1132
1133 XSH_MALLOC( slit_vect_data, double , ssize);
1134 XSH_MALLOC( sliterr_vect_data, double, ssize);
1135 XSH_MALLOC( slit_pos_data, double , ssize);
1136
1137 XSH_MALLOC( spos_data, double , wsize);
1138 XSH_MALLOC( errpos_data, double, wsize);
1139 XSH_MALLOC( wpos_data, double , wsize);
1140
1141 XSH_MALLOC( sigma_data, double , wsize);
1142 XSH_MALLOC( snr_data, double , wsize);
1143
1144/* create sky mask */
1145 XSH_CALLOC( sky_mask, int, wsize);
1146
1147 if ( skymask_frame != NULL){
1148 float *skymask_data = NULL;
1149 int irow, nrow;
1150 double fwhm =0.0, sky_min, sky_max;
1151 int isky_min, isky_max, imask;
1152 double width, resolution;
1153
1154 XSH_TABLE_LOAD( skymask_table, skymask_name);
1155 /* sort wavelength */
1156 check( xsh_sort_table_1( skymask_table, "WAVELENGTH", CPL_FALSE));
1157 check( skymask_data = cpl_table_get_data_float( skymask_table,
1158 "WAVELENGTH"));
1159 check( nrow = cpl_table_get_nrow( skymask_table));
1160
1161 for( irow=0; irow < nrow; irow++){
1163 resolution = xsh_resolution_get( instrument, width);
1164 fwhm = skymask_data[irow]/resolution;
1165 sky_min = skymask_data[irow]-fwhm;
1166 sky_max = skymask_data[irow]+fwhm;
1167 isky_min = (int)xsh_round_double((sky_min-wmin)/wstep);
1168 if (isky_min < 0){
1169 isky_min=0;
1170 }
1171 isky_max = (int)xsh_round_double((sky_max-wmin)/wstep);
1172 if ( isky_max >= wsize){
1173 isky_max = wsize-1;
1174 }
1175 for( imask=isky_min; imask <=isky_max; imask++){
1176 sky_mask[imask] = 1;
1177 }
1178 }
1179
1180 if (level >= XSH_DEBUG_LEVEL_MEDIUM){
1181 FILE *mask_file = NULL;
1182 char mask_name[256];
1183 int idbg=0;
1184
1185 sprintf( mask_name, "skymask.reg");
1186 mask_file = fopen( mask_name, "w");
1187
1188 fprintf( mask_file,"# Region file format: DS9 version 4.1\n");
1189 fprintf( mask_file,"global color=green dashlist=8 3 width=1 font=\"helvetica 10 normal\"\
1190 select=1 highlite=1 dash=0 fixed=0 edit=1 move=1 delete=1 include=1 source=1\n");
1191 fprintf(mask_file,"physical\n");
1192
1193 for(idbg=0; idbg< wsize; idbg++){
1194
1195 if (sky_mask[idbg] == 1){
1196 fprintf( mask_file, "line(%d,%d,%d,%d)\n",idbg+1, ssize, idbg+1, 1);
1197 }
1198 }
1199 fclose( mask_file);
1200 }
1201 }
1202
1203
1204 data_size = 0;
1205
1206 for (i = 0; i < wsize; i++) {
1207 double slitcen = 0;
1208 double sigma = 0;
1209 double area = 0;
1210 double offset = 0;
1211 double a = 0, b = 0, c = 0;
1212 double frac_rej = 0;
1213
1214 cpl_matrix *covariance = NULL;
1215 double fit_err;
1216 int good_fit = 0;
1217
1218 slit_vect_size = 0;
1219
1220 /* get slit data */
1221 for (j = jmin; j < jmax; j++) {
1222 double val = 0, err = 0;
1223 int start, end;
1224 int ngood = 0, nbad = 0;
1225
1226 /* running box hsize */
1227 start = i - box_hsize;
1228 end = i + box_hsize;
1229
1230 if (start < 0) {
1231 start = 0;
1232 }
1233 if (end >= wsize) {
1234 end = wsize - 1;
1235 }
1236
1237 for (k = start; k <= end; k++) {
1238 if ( (qual[k + j * wsize] & instrument->decode_bp) == 0 ) {
1239 val += flux[k + j * wsize];
1240 ngood++;
1241 } else {
1242 nbad++;
1243 }
1244
1245 }
1246
1247 val += flux[i+j*wsize];
1248 err = errs[i+j*wsize];
1249
1250 if ( ngood > 0){
1251 val *= (double)(nbad+ngood)/(double)ngood;
1252
1253 slit_pos_data[slit_vect_size] = j;
1254 sliterr_vect_data[slit_vect_size] = err;
1255 slit_vect_data[slit_vect_size] = val;
1256
1257 slit_vect_size++;
1258 }
1259 }
1260
1261 frac_rej = 1-((double)slit_vect_size/(double)ssize);
1262 if ( frac_rej < 0.5){
1263 check( slit_pos = cpl_vector_wrap( slit_vect_size, slit_pos_data));
1264 check( slit_vect = cpl_vector_wrap( slit_vect_size, slit_vect_data));
1265 check( sliterr_vect = cpl_vector_wrap( slit_vect_size, sliterr_vect_data));
1266 /* smooth with median filter */
1267 check( smooth_vect = cpl_vector_filter_median_create( slit_vect,
1268 smooth_hsize));
1269
1270 /* gaussian fit */
1271#if CPL_MODE
1272 cpl_vector_fit_gaussian( slit_pos, NULL, slit_vect, sliterr_vect,
1273 CPL_FIT_ALL, &slitcen, &sigma, &area, &offset, &mse, NULL, &covariance);
1274 a = offset;
1275 fit_err = sqrt(cpl_matrix_get(covariance,0,0))*sstep;
1276 good_fit = (cpl_error_get_code() == CPL_ERROR_NONE);
1277#else
1278 {
1279 double init_par[6];
1280 double fit_errs[6];
1281 int status = 0;
1282 xsh_gsl_init_gaussian_fit( slit_pos, slit_vect, init_par);
1283 xsh_gsl_fit_gaussian( slit_pos, slit_vect, deg, init_par, fit_errs, &status);
1284 area = init_par[0];
1285 a = init_par[1];
1286 b = init_par[2];
1287 c = init_par[3];
1288 slitcen = init_par[4];
1289 sigma = init_par[5];
1290 offset = a+b*slitcen+c*slitcen*slitcen;
1291 fit_err = fit_errs[4]*sstep;
1292 good_fit= (status == 0);
1293 }
1294#endif
1295
1296 if ( good_fit){
1297 double height=0;
1298 double snr=0;
1299
1300 height = area / sqrt(2*M_PI*sigma*sigma);
1301 snr = (height+offset)/offset;
1302
1303 /* clean bad fit */
1304 if ( (area > 0) && (slitcen > 0) && (slitcen < ssize)
1305 && (sigma < ssize) && (snr > 0)){
1306 spos_data[data_size] = smin+slitcen*sstep;
1307 errpos_data[data_size] = fit_err;
1308 wpos_data[data_size] = wmin+i*wstep;
1309 sigma_data[data_size] = sigma*sstep;
1310 snr_data[data_size] = snr;
1311#if 0
1312 /* DEBUG only */
1313 if ( (wpos_data[data_size] >= (1500-wstep)) && (wpos_data[data_size] <= (1500+wstep))){
1314 char test_name[256];
1315 FILE* test_file = NULL;
1316 int itest;
1317 double dtest;
1318
1319 sprintf( test_name, "data_w%f_%s.dat", wpos_data[data_size],resname);
1320 XSH_REGDEBUG("Produce test file %s", test_name);
1321 test_file = fopen( test_name, "w+");
1322 fprintf( test_file, "# pos slit\n");
1323
1324 for( itest=0; itest < slit_vect_size; itest++){
1325 fprintf( test_file, "%f %f\n", cpl_vector_get( slit_pos, itest), cpl_vector_get( slit_vect, itest));
1326 }
1327
1328 fclose( test_file);
1329
1330 sprintf( test_name, "gauss_w%f_%s.dat", wpos_data[data_size],resname);
1331 XSH_REGDEBUG("Produce test file %s", test_name);
1332
1333 test_file = fopen( test_name, "w+");
1334 fprintf( test_file, "# pos gauss offset x0=%f_sig=%f_area=%f_offset=%f\n", slitcen, sigma, area, offset);
1335
1336 for( dtest=0; dtest < slit_vect_size; dtest+=0.1){
1337 double z, gauss;
1338 double off;
1339
1340 off = a+b*dtest+c*dtest*dtest;
1341 z = ( dtest-slitcen)/(sigma*XSH_MATH_SQRT_2);
1342 gauss = height*exp(-(z*z))+off;
1343
1344 fprintf( test_file, "%f %f %f\n", dtest, gauss, off);
1345 }
1346
1347 fclose( test_file);
1348 }
1349#endif
1350 data_size++;
1351 }
1352 }
1353 else{
1355 }
1356 xsh_free_matrix( &covariance);
1357 xsh_free_vector( &smooth_vect);
1358 xsh_unwrap_vector( &slit_pos);
1359 xsh_unwrap_vector( &slit_vect);
1360 xsh_unwrap_vector( &sliterr_vect);
1361 }
1362 }
1363
1364
1365 /* sort sigma and snr */
1366 check( sigma_vect = cpl_vector_wrap( data_size, sigma_data));
1367 check( snr_vect = cpl_vector_wrap( data_size, snr_data));
1368 check( sigma_sort_vect = cpl_vector_duplicate( sigma_vect));
1369 check( snr_sort_vect = cpl_vector_duplicate( snr_vect));
1370
1371 check( cpl_vector_sort( sigma_sort_vect, 1));
1372 check( cpl_vector_sort( snr_sort_vect, 1));
1373
1374 n05 = (data_size-1)*cut_sigma_low;
1375 i = ceil( n05);
1376 sigma_05 = cpl_vector_get( sigma_sort_vect, i);
1377 n05 = (data_size-1)*cut_snr_low;
1378 i = ceil( n05);
1379 snr_05 = cpl_vector_get( snr_sort_vect, i);
1380
1381
1382 n95 = (data_size-1)*cut_sigma_up;
1383 i = floor( n95);
1384 sigma_95 = cpl_vector_get( sigma_sort_vect, i);
1385 n95 = (data_size-1)*cut_snr_up;
1386 i = floor( n95);
1387 snr_95 = cpl_vector_get( snr_sort_vect, i);
1388
1389 /* filter on snr and sigma */
1390 for( i=0; i< data_size; i++){
1391 double sigma, snr;
1392
1393 sigma = sigma_data[i];
1394 snr = snr_data[i];
1395 if ( sigma_05 <= sigma && sigma <= sigma_95 && snr_05 <= snr && snr <=snr_95){
1396 wpos_data[ndata_size] = wpos_data[i];
1397 spos_data[ndata_size] = spos_data[i];
1398 errpos_data[ndata_size] = errpos_data[i];
1399 ndata_size++;
1400 }
1401 }
1402 xsh_msg_dbg_low( "Filtering by snr [%f,%f] and sigma [%f,%f] from %d lines to %d lines",
1403 sigma_05, sigma_95, snr_05, snr_95, data_size, ndata_size);
1404
1405 /* REGDEBUG */
1406 if (level >= XSH_DEBUG_LEVEL_MEDIUM){
1407 FILE *test_file = NULL;
1408 int itest;
1409 char test_name[256];
1410#if CPL_MODE
1411 sprintf( test_name, "cpl_gaussian_fit_%s.dat", resname);
1412#else
1413 sprintf( test_name, "gsl_gaussian_fit_%s.dat", resname);
1414#endif
1415 test_file = fopen( test_name, "w");
1416 fprintf( test_file, "# wavelength slit_fit fit_err sigma\n");
1417
1418 for(itest=0; itest<ndata_size; itest++){
1419 fprintf( test_file, "%f %f %f %f\n", wpos_data[itest], spos_data[itest], errpos_data[itest],
1420 sigma_data[itest]);
1421 }
1422 xsh_msg_dbg_medium( "Produce file %s", test_name);
1423 fclose( test_file);
1424 }
1425
1426 check( wpos_vect = cpl_vector_new( wsize));
1427 check( spos_vect = cpl_vector_new( wsize));
1428 check( sposg_data = cpl_vector_get_data( spos_vect));
1429 check( wposg_data = cpl_vector_get_data( wpos_vect));
1430 j=0;
1431
1432 for(i=0; i< wsize; i++){
1433 double wave, wkeep;
1434 double slit;
1435
1436 wave = wmin+i*wstep;
1437 wkeep = wpos_data[j];
1438 check( cpl_vector_set( wpos_vect, i, wave));
1439 if ( fabs(wave -wkeep) < 0.0000001){
1440 slit = spos_data[j];
1441 if (j < (ndata_size-1)){
1442 j++;
1443 }
1444 }
1445 else{
1446 slit = xsh_data_interpolate( wave, ndata_size, wpos_data, spos_data);
1447 }
1448 check( cpl_vector_set( spos_vect, i, slit));
1449 }
1450
1451 check( decomp = xsh_atrous( spos_vect, nscales));
1452
1453 nb_scales = nscales-HF_skip;
1454
1455 /* filter High frequency */
1456 for( i=0; i< wsize; i++){
1457 sposg_data[i] = 0;
1458 for(j=0; j<nb_scales; j++){
1459 sposg_data[i] += cpl_matrix_get( decomp, j, i);
1460 }
1461 }
1462
1463 /* Save result in a table */
1464 check( table = cpl_table_new( wsize));
1466 XSH_OBJPOS_UNIT_WAVELENGTH, CPL_TYPE_DOUBLE);
1468 XSH_OBJPOS_UNIT_SLIT, CPL_TYPE_DOUBLE);
1469
1470 for( i=0; i< wsize; i++){
1471 check( cpl_table_set_double( table, XSH_OBJPOS_COLNAME_WAVELENGTH,
1472 i, wposg_data[i]));
1473 check( cpl_table_set_double( table, XSH_OBJPOS_COLNAME_SLIT,
1474 i, sposg_data[i]));
1475 }
1476 sprintf( tablename, resname);
1477 header = cpl_propertylist_new();
1478 check( cpl_table_save( table, header, NULL, tablename, CPL_IO_DEFAULT));
1479 /* Create the frame */
1480 check(result=xsh_frame_product( tablename,
1481 "OBJPOS_TAB",
1482 CPL_FRAME_TYPE_TABLE,
1483 CPL_FRAME_GROUP_PRODUCT,
1484 CPL_FRAME_LEVEL_TEMPORARY));
1485
1486
1487 check (xsh_add_temporary_file( tablename));
1488
1489 cleanup:
1490 XSH_TABLE_FREE( table);
1491 xsh_free_propertylist( &header);
1492 xsh_free_matrix( &decomp);
1493 xsh_free_vector( &spos_vect);
1494 xsh_free_vector( &wpos_vect);
1495 xsh_free_vector( &sigma_sort_vect);
1496 xsh_free_vector( &snr_sort_vect);
1497 xsh_unwrap_vector( &sigma_vect);
1498 xsh_unwrap_vector( &snr_vect);
1499 XSH_FREE( sigma_data);
1500 XSH_FREE( snr_data);
1501 XSH_FREE( spos_data);
1502 XSH_FREE( errpos_data);
1503 XSH_FREE( wpos_data);
1504 XSH_FREE( slit_pos_data);
1505 XSH_FREE( slit_vect_data);
1506 XSH_FREE( sliterr_vect_data);
1507 xsh_free_vector( &smooth_vect);
1508 xsh_spectrum_free( &spectrum2d);
1509 XSH_FREE(sky_mask);
1510 xsh_free_table( &skymask_table);
1511 return result;
1512}
1513/*---------------------------------------------------------------------------*/
1514
1515/*---------------------------------------------------------------------------*/
1527cpl_frameset* xsh_localize_ifu( cpl_frameset *merge2d_frameset,
1528 cpl_frame *skymask_frame,
1530 const char* prefix)
1531{
1532 int i, slitlet;
1533 cpl_frameset *result_frameset = NULL;
1534 char fname[256];
1535 int smooth_hsize;
1536 int nscales;
1537 int HF_skip;
1538 double cut_sigma_low;
1539 double cut_sigma_up;
1540 double cut_snr_low;
1541 double cut_snr_up;
1542 double slit_min = -6.0;
1543 double slit_max = 6.0;
1544 cpl_frame *frame = NULL;
1545 const char *frame_name = NULL;
1546 cpl_propertylist *header = NULL;
1547 int deg = 2, skymask=0;
1548 int box_hsize;
1549 cpl_frame *mask_frame = NULL;
1550
1551 XSH_ASSURE_NOT_NULL( merge2d_frameset);
1553 XSH_ASSURE_NOT_NULL( locifu_par);
1554
1555 smooth_hsize = locifu_par->smooth_hsize;
1556 nscales = locifu_par->nscales;
1557 HF_skip = locifu_par->HF_skip;
1558 cut_sigma_low = locifu_par->cut_sigma_low;
1559 cut_sigma_up = locifu_par->cut_sigma_up;
1560 cut_snr_low = locifu_par->cut_snr_low;
1561 cut_snr_up = locifu_par->cut_snr_up;
1562 skymask = locifu_par->use_skymask;
1563 box_hsize = locifu_par->box_hsize;
1564 /* get slit borders */
1565 if (skymask){
1566 mask_frame = skymask_frame;
1567 }
1568
1569 check( frame = cpl_frameset_get_frame( merge2d_frameset, 0));
1570 check( frame_name = cpl_frame_get_filename( frame));
1571 check( header = cpl_propertylist_load( frame_name, 0));
1572 check( slit_min = xsh_pfits_get_rectify_space_min( header));
1573 xsh_free_propertylist( &header);
1574
1575 check( frame = cpl_frameset_get_frame( merge2d_frameset, 2));
1576 check( frame_name = cpl_frame_get_filename( frame));
1577 check( header = cpl_propertylist_load( frame_name, 0));
1578 check( slit_max = xsh_pfits_get_rectify_space_max( header));
1579 xsh_free_propertylist( &header);
1580
1581 slit_min += locifu_par->slitlow_edges_mask;
1582 slit_max -= locifu_par->slitup_edges_mask;
1583
1584 deg = locifu_par->bckg_deg;
1585
1586 check( result_frameset = cpl_frameset_new());
1587
1588 for( i = 0, slitlet = LOWER_IFU_SLITLET ; i < 3 ; i++, slitlet++ ) {
1589 cpl_frame * loc_frame = NULL;
1590 cpl_frame *merge2d_frame = NULL;
1591
1592 sprintf( fname ,"%s_LOCIFU_%s_%s.fits", prefix, SlitletName[slitlet],
1594
1595 xsh_msg( "Localizing IFU in [%f,%f] slitlet %s, frame '%s'", slit_min, slit_max,
1596 SlitletName[slitlet], fname);
1597
1598 check( merge2d_frame = cpl_frameset_get_frame( merge2d_frameset, i));
1599
1600
1601 check( loc_frame = xsh_localize_ifu_slitlet( merge2d_frame, mask_frame,
1602 smooth_hsize, nscales,
1603 HF_skip, fname, cut_sigma_low, cut_sigma_up, cut_snr_low, cut_snr_up,
1604 slit_min, slit_max, deg, box_hsize, instrument));
1605
1606 check( cpl_frameset_insert( result_frameset, loc_frame));
1607 }
1608
1609 cleanup:
1610 if ( cpl_error_get_code() != CPL_ERROR_NONE ) {
1611 xsh_free_frameset( &result_frameset);
1612 xsh_free_propertylist( &header);
1613 }
1614 return result_frameset;
1615}
1616/*---------------------------------------------------------------------------*/
static int width
static double sigma
static xsh_instrument * instrument
static float slit_step
static float lambda_step
cpl_frame * xsh_localization_save(xsh_localization *list, const char *filename, xsh_instrument *instrument)
save a localization to a frame
void xsh_localization_free(xsh_localization **list)
free memory associated to a localization_list
xsh_localization * xsh_localization_create(void)
Create an empty localization list.
int xsh_spectrum_get_size_lambda(xsh_spectrum *s)
Get lambda axis size of spectrum.
double xsh_spectrum_get_lambda_min(xsh_spectrum *s)
Get minimum lambda of spectrum.
xsh_spectrum * xsh_spectrum_load(cpl_frame *s1d_frame)
Load a 1D spectrum structure.
double * xsh_spectrum_get_errs(xsh_spectrum *s)
Get errs of spectrum.
double xsh_spectrum_get_lambda_step(xsh_spectrum *s)
Get bin in lambda of spectrum.
int * xsh_spectrum_get_qual(xsh_spectrum *s)
Get qual of spectrum.
double * xsh_spectrum_get_flux(xsh_spectrum *s)
Get flux of spectrum.
void xsh_spectrum_free(xsh_spectrum **s)
free memory associated to an 1D spectrum
int xsh_spectrum_get_size_slit(xsh_spectrum *s)
Get slit axis ize of spectrum.
#define XSH_REGDEBUG(...)
Definition: xsh_error.h:132
#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_ILLEGAL_MSG(cond, msg)
Definition: xsh_error.h:111
#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.
double xsh_resolution_get(xsh_instrument *instrument, double slit)
Get the resoltion.
XSH_ARM xsh_instrument_get_arm(xsh_instrument *i)
Get an arm on instrument structure.
double xsh_convert_seeing(cpl_frame *frame)
Convert seeing keywork in mu sigma.
static void chunk_coadd(double *data, double *flux, int *qual, int nx, int ny, int ifirst, int ilast, const int decode_bp, int *skymask)
cpl_frameset * xsh_localize_ifu(cpl_frameset *merge2d_frameset, cpl_frame *skymask_frame, xsh_localize_ifu_param *locifu_par, xsh_instrument *instrument, const char *prefix)
Localize center of object on a merge 2D IFU slitlet.
static void tab_minmax_get_index(double *data, int size, xsh_slit_limit_param *slit_limit_par, int *min, int *max)
static xsh_localization * xsh_localize_obj_auto(cpl_frame *rec_frame, cpl_frame *skymask_frame, xsh_instrument *instrument, xsh_localize_obj_param *loc_obj_par, xsh_slit_limit_param *slit_limit_par, const int decode_bp)
Build localization table frame from the rectified table frame.
cpl_frame * xsh_localize_ifu_slitlet(cpl_frame *merge2d_slitlet, cpl_frame *skymask_frame, int smooth_hsize, int nscales, int HF_skip, const char *resname, double cut_sigma_low, double cut_sigma_up, double cut_snr_low, double cut_snr_up, double slit_min, double slit_max, int deg, int box_hsize, xsh_instrument *instrument)
Localize center of object on a merge 2D IFU slitlet.
cpl_frame * xsh_localize_obj(cpl_frame *rec_frame, cpl_frame *skymask_frame, xsh_instrument *instrument, xsh_localize_obj_param *loc_obj_par, xsh_slit_limit_param *slitlimit_par, const char *fname)
Build the localization table.
cpl_frameset * xsh_localize_obj_ifu(cpl_frameset *rec_frameset, cpl_frame *skymask_frame, xsh_instrument *instrument, xsh_localize_obj_param *locobj_par, xsh_slit_limit_param *slitlimit_par)
static xsh_localization * xsh_localize_obj_manual(xsh_instrument *instrument, xsh_localize_obj_param *loc_obj_par)
cpl_frame * xsh_merge_ord(cpl_frame *sci_frame, xsh_instrument *instrument, int merge, const char *rec_prefix)
Merge orders of the rectified frame using merge parameters.
int size
#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(...)
Print a message on info level.
Definition: xsh_msg.h:121
#define xsh_msg_dbg_low(...)
Definition: xsh_msg.h:48
double xsh_pfits_get_rectify_space_max(cpl_propertylist *plist)
find out the rectify SPACE max
Definition: xsh_pfits.c:3370
double xsh_pfits_get_rectify_space_min(cpl_propertylist *plist)
find out the rectify space min
Definition: xsh_pfits.c:3353
double xsh_pfits_get_seeing_start(const cpl_propertylist *plist)
find out the TEL AMBI START value (Seeing)
Definition: xsh_pfits.c:544
double xsh_pfits_get_slit_width(const cpl_propertylist *plist, xsh_instrument *instrument)
find out the INS OPTIx NAME value (the width of the slit)
Definition: xsh_pfits.c:580
double xsh_pfits_get_airm_mean(const cpl_propertylist *plist)
find out the mean airmass value
Definition: xsh_pfits.c:511
double xsh_pfits_get_seeing_end(const cpl_propertylist *plist)
find out the TEL AMBI END value (Seeing)
Definition: xsh_pfits.c:561
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
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_free_frameset(cpl_frameset **f)
Deallocate a frame set and set the pointer to NULL.
Definition: xsh_utils.c:2254
cpl_error_code xsh_sort_table_1(cpl_table *t, const char *column, cpl_boolean reverse)
Sort a table by one column.
double xsh_data_interpolate(double wav, int nrow, double *pw, double *pe)
Interpolate data points.
void xsh_gsl_init_gaussian_fit(cpl_vector *xpos_vect, cpl_vector *ypos_vect, double *init_par)
Definition: xsh_utils.c:6862
long xsh_round_double(double x)
Computes round(x)
Definition: xsh_utils.c:4383
void xsh_free_matrix(cpl_matrix **m)
Deallocate a matrix and set the pointer to NULL.
Definition: xsh_utils.c:2209
void xsh_gsl_fit_gaussian(cpl_vector *xpos_vect, cpl_vector *ypos_vect, int deg, double *params, double *errs, int *status)
Definition: xsh_utils.c:6936
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_add_temporary_file(const char *name)
Add temporary file to temprary files list.
Definition: xsh_utils.c:1432
cpl_polynomial * cenpoly
cpl_polynomial * edguppoly
cpl_polynomial * edglopoly
enum localize_method method
cpl_propertylist * flux_header
cpl_image * flux
cpl_image * qual
#define XSH_ARCSEC_UVB
@ XSH_ARM_UVB
@ XSH_ARM_NIR
@ XSH_ARM_VIS
#define XSH_ARCSEC_NIR
#define XSH_ARCSEC_VIS
int nx
double kappa
Definition: xsh_detmon_lg.c:81
int niter
Definition: xsh_detmon_lg.c:82
int threshold
Definition: xsh_detmon_lg.c:90
int ny
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_MATH_SQRT_2
Definition: xsh_drl.h:567
#define XSH_OBJPOS_UNIT_SLIT
Definition: xsh_drl.h:499
#define XSH_OBJPOS_COLNAME_SLIT
Definition: xsh_drl.h:498
#define XSH_LAMBDA_DIMM
Definition: xsh_drl.h:500
#define XSH_OBJPOS_UNIT_WAVELENGTH
Definition: xsh_drl.h:497
#define XSH_OBJPOS_COLNAME_WAVELENGTH
Definition: xsh_drl.h:496
#define max(a, b)
@ LOWER_IFU_SLITLET
Definition: xsh_ifu_defs.h:42
static const char * SlitletName[]
Definition: xsh_ifu_defs.h:48
#define LOCALIZE_METHOD_PRINT(method)
@ LOC_GAUSSIAN_METHOD
@ LOC_MANUAL_METHOD
#define XSH_TABLE_NEW_COL(TABLE, NAME, UNIT, TYPE)
Definition: xsh_utils.h:113
#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
#define M_PI
Definition: xsh_utils.h:43
cpl_matrix * xsh_atrous(cpl_vector *spec, int nscales)
Do a wavelet decomposition using atrous from IDL.
#define XSH_TABLE_LOAD(TABLE, NAME)
#define XSH_TABLE_FREE(TABLE)
cpl_polynomial * xsh_polynomial_fit_1d_create(const cpl_vector *x_pos, const cpl_vector *values, int degree, double *mse)