X-shooter Pipeline Reference Manual 3.8.15
xsh_follow_arclines.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 16:08:02 $
23 * $Revision: 1.169 $
24 */
25
26#ifdef HAVE_CONFIG_H
27# include <config.h>
28#endif
29#define REGDEBUG_WAVECAL 0
30/*----------------------------------------------------------------------------*/
43/*----------------------------------------------------------------------------*/
46/*-----------------------------------------------------------------------------
47 Includes
48 -----------------------------------------------------------------------------*/
49#include <math.h>
50#include <xsh_drl.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_utils_wrappers.h>
57#include <xsh_badpixelmap.h>
58#include <xsh_data_order.h>
59#include <xsh_data_the_map.h>
60#include <xsh_data_arclist.h>
61#include <xsh_data_wavesol.h>
62#include <xsh_data_wavemap.h>
63#include <xsh_data_resid_tab.h>
64#include <xsh_data_linetilt.h>
65#include <xsh_utils.h>
66#include <xsh_ifu_defs.h>
67#include <xsh_model_io.h>
68#include <xsh_model_kernel.h>
69#include <xsh_model_utils.h>
71#include <xsh_data_shift_tab.h>
72#include <xsh_data_dispersol.h>
73
74#define XSH_SPECRES_CLIP_KAPPA 3.
75#define XSH_SPECRES_CLIP_NITER 2
76#define XSH_SPECRES_CLIP_FRAC 0.5
77
78/*-----------------------------------------------------------------------------
79 Functions prototypes
80 -----------------------------------------------------------------------------*/
81
82//static const char* fwhm_debug_mode = NULL;
83
84enum {
86} ;
87
88enum {
91} ;
92
93typedef struct {
94 double xpos ;
95 double ypos ;
96 double fwhm ;
97 double area;
98 int good ;
99} CENTROIDS ;
100
101static int detect_centroid( xsh_pre * pre,
102 float lambda, int ordnum, double xpix,
103 double ypix,
104 xsh_follow_arclines_param * follow_param,
105 int is_center,
106 XSH_GAUSSIAN_FIT * fit_res );
107
108/*-----------------------------------------------------------------------------
109 Implementation
110 -----------------------------------------------------------------------------*/
111
127// FIXME is_center parameter NOT USED!!
128static int
130 float lambda,
131 int ordnum,
132 double xpix,
133 double y,
134 xsh_follow_arclines_param * follow_param,
135 int is_center,
136 XSH_GAUSSIAN_FIT * fit_res )
137{
138 cpl_vector * pixpos = NULL ;
139 cpl_vector * pixval = NULL ;
140 int first, last, i, nelem, ny, ix ;
141 double dy, ypix;
142 int ret = 1 ;
143
144 if ( xpix <= 0 || xpix > pre->nx || y <= 0 || y > pre->ny )
145 return ret ;
146
147 XSH_ASSURE_NOT_ILLEGAL( xpix > 0 && xpix <= pre->nx &&
148 y > 0 && y <= pre->ny );
149
150 ix = xsh_round_double(xpix);
151 ypix = xsh_round_double( y);
152
153 first = ypix - follow_param->range;
154
155 if ( first < 0 ) {
156 xsh_msg_dbg_high( "%lf - %d : Y pixel < 0 (%d)", lambda, ordnum, first ) ;
157 first = 1 ;
158 }
159 last = ypix + follow_param->range ;
160 if( last > pre->ny ) {
161 xsh_msg_dbg_low( "Y = %d too high, ignore", last ) ;
162 return ret ; /* Y too high, ignore */
163 }
164
165 dy = ypix - follow_param->range ;
166 nelem = (follow_param->range*2) + 1 ;
167
168 check( pixpos = cpl_vector_new( nelem ) ) ;
169 check( pixval = cpl_vector_new( nelem ) ) ;
170
171 for( i = 0, ny = first ; ny <= last ; ny++, i++, dy += 1. ) {
172 int rej ;
173 double value ;
174
175 cpl_error_reset() ;
176 value = cpl_image_get( pre->data, ix, ny, &rej ) ;
177 if ( cpl_error_get_code() != CPL_ERROR_NONE ) {
178 xsh_msg_dbg_high( " *** X,Y out of range %d,%d", ix, ny ) ;
179 cpl_error_reset() ;
180 continue ;
181 }
182 cpl_vector_set( pixval, i, value ) ;
183 cpl_vector_set( pixpos, i, dy ) ;
184 }
185
186 cpl_error_reset() ;
187 xsh_vector_fit_gaussian( pixpos, pixval, fit_res ) ;
188 if ( cpl_error_get_code() != CPL_ERROR_NONE ) {
189 xsh_msg_dbg_high("Failed to fit Y centroid for x %d y [%d, %d], wavelength %f",
190 ix, first, last, lambda);
191 cpl_error_reset() ;
192 }
193 else {
195 FILE *debug_file = NULL;
196 char name[256];
197 int idebug = 0;
198
199
200 sprintf( name, "profil_%d_%f_%f.dat", ordnum, lambda, xpix);
201
202 debug_file = fopen( name, "w");
203 fprintf( debug_file, "# ypos yflux yfit\n");
204 for( idebug=0; idebug< nelem; idebug++){
205 double gauss;
206 double x;
207
208 x = cpl_vector_get( pixpos, idebug);
209 gauss = fit_res->area / sqrt(2 *M_PI*fit_res->sigma*fit_res->sigma)*
210 exp( -(x-fit_res->peakpos)*(x-fit_res->peakpos)/(2*fit_res->sigma*fit_res->sigma)) + fit_res->offset;
211
212 fprintf( debug_file, "%f %f %f\n", x, cpl_vector_get( pixval, idebug), gauss);
213 }
214 fclose( debug_file);
215 }
216 xsh_msg_dbg_medium( "For x %d, y [%d, %d] centroid fit in Y = %lf",
217 ix, first, last, fit_res->peakpos);
218 ret = 0;
219 }
220
221 cleanup:
222 xsh_free_vector( &pixpos ) ;
223 xsh_free_vector( &pixval ) ;
224 return ret ;
225}
226
227static cpl_polynomial *
228get_slit_ifu_lo_poly( xsh_order * porder,int ifu_flag )
229{
230 switch ( ifu_flag ) {
231 case CENTER_SLIT:
232 return porder->edglopoly ;
234 return porder->edglopoly ;
236 return porder->sliclopoly ;
238 return porder->slicuppoly ;
239 default:
240 return NULL ;
241 }
242}
243
244static cpl_polynomial *
245get_slit_ifu_up_poly( xsh_order * porder,int ifu_flag )
246{
247 switch ( ifu_flag ) {
248 case CENTER_SLIT:
249 return porder->edguppoly ;
251 return porder->sliclopoly ;
253 return porder->slicuppoly ;
255 return porder->edguppoly ;
256 default:
257 return NULL ;
258 }
259}
260
261/*****************************************************************************/
262
263/*****************************************************************************/
264/*****************************************************************************/
265static int
266find_tilt( double yp0, double xc, float lambda, int ordnum,
267 xsh_follow_arclines_param * follow_param,
268 xsh_pre * pre, xsh_order_list * orders,
269 xsh_instrument * instrument, double * slope,
270 double * chisq, double * minx, double * maxx,
271 int * nt, int * ng, double * fwhm_center,
272 double * good_center, int ifu_flag )
273{
274 int result=0;
275 cpl_polynomial * lo_poly = NULL;
276 cpl_polynomial * up_poly = NULL;
277
278 CENTROIDS * centers = NULL ;
279 double x0, x, ycenter ;
280 XSH_GAUSSIAN_FIT fit_res ;
281 cpl_vector *pos_vect =NULL, *val_vect = NULL;
282 int *flags = NULL;
283 cpl_polynomial * poly_tilt = NULL ;
284 double xmin, xmax ;
285 int iorder ;
286 int i, j, k, ngood;
287 int csize ;
288 cpl_vector *fwhm_median_vect = NULL;
289 double *median_data = NULL;
290 double fwhm_median = 0.0;
291 int tilt_niter;
292 double tilt_kappa, tilt_frac, good_frac;
293 int ntot, no_centroid =0;
294 cpl_size deg1 = 1;
295
296 XSH_ASSURE_NOT_NULL( follow_param);
297 XSH_ASSURE_NOT_NULL( follow_param->tilt_clipping);
298
299 tilt_niter = follow_param->tilt_clipping->niter;
300 tilt_kappa = follow_param->tilt_clipping->sigma;
301 tilt_frac = follow_param->tilt_clipping->frac;
302
303 /* Extract from order table lower and upper edge of the arc line */
304 check( iorder = xsh_order_list_get_order( orders, ordnum));
305
306 if ( iorder < 0 ) {
307 xsh_msg( "Unknown Absolute Order %d", ordnum);
309 }
310
311 lo_poly = get_slit_ifu_lo_poly( &orders->list[iorder], ifu_flag ) ;
312 up_poly = get_slit_ifu_up_poly( &orders->list[iorder], ifu_flag ) ;
313
314 check( xmin = xsh_order_list_eval( orders, lo_poly, yp0));
315 check( xmax = xsh_order_list_eval( orders, up_poly, yp0));
316
317 /* Suppress points too close to the order edges */
318 xmin = *minx;
319 xmax = *maxx;
320
321 xmin += follow_param->margin ;
322 xmax -= follow_param->margin ;
323
324 *minx = xmin ;
325 *maxx = xmax ;
326
327 xsh_msg_dbg_medium( " Find_tilt - order %d, idx %d, xmin = %.2lf, xmax = %.2lf, yc = %.2lf",
328 ordnum, iorder, xmin, xmax, yp0 ) ;
329
330 if ( xmax <= xmin ) {
331 xsh_msg( "******* Something Wrong in Edges!" ) ;
332 return FIND_TILT_BAD_EDGES ;
333 }
334
335 if ( xc <= xmin || xc >= xmax ) {
336 xsh_msg( "lambda %f, order %d: Center (from wavesol) %lf not in ] %lf, %lf [",
337 lambda, ordnum, xc, xmin, xmax ) ;
338 return FIND_TILT_BAD_CENTER ;
339 }
340
341
342 /* Beware !!! (int)(xmax - xmin + 1 ) may lead to problems ! */
343 csize = ceil(xmax) - floor(xmin);
344 XSH_CALLOC( centers, CENTROIDS, csize);
345
346 x0 = xsh_round_double(xc);
347 ngood = 0 ;
348 ntot = 0;
349
350 ycenter = xsh_round_double(yp0);
351
352 /* First go from center to lower X edge */
353 for( x = x0 ; x >= xmin ; x--) {
354 int ok;
355
356 centers[ntot].xpos = x;
357 ok = detect_centroid( pre, lambda, ordnum, x, ycenter,
358 follow_param, 0, &fit_res);
359 if ( ok == 0 && fit_res.peakpos > 0. ) {
360 centers[ntot].ypos = fit_res.peakpos ;
361 centers[ntot].fwhm = CPL_MATH_FWHM_SIG*fit_res.sigma;
362 centers[ntot].area = fit_res.area;
363 centers[ntot].good = 0 ;
364 ycenter = fit_res.peakpos ;
365 ngood++ ;
366 }
367 else {
368 no_centroid++ ;
369 centers[ntot].good = 1 ;
370 }
371 ntot++;
372 }
373 ycenter = xsh_round_double( yp0);
374
375 /* Last go from center(+1) to upper X edge */
376 for( x = x0+1 ; x <= xmax ; x++) {
377 int ok;
378
379 centers[ntot].xpos = x;
380 ok = detect_centroid( pre, lambda, ordnum, x, ycenter,
381 follow_param, 0, &fit_res);
382 if ( ok == 0 && fit_res.peakpos > 0. ) {
383 centers[ntot].ypos = fit_res.peakpos ;
384 centers[ntot].fwhm = CPL_MATH_FWHM_SIG*fit_res.sigma ;
385 centers[ntot].area = fit_res.area;
386 centers[ntot].good = 0 ;
387 ycenter = fit_res.peakpos ;
388 ngood++ ;
389 }
390 else {
391 no_centroid++ ;
392 centers[ntot].good = 1 ;
393 }
394 ntot++;
395 }
396 xsh_msg_dbg_medium(" For lambda %f in order %d : Fitted %d/%d (%d no centroid)",
397 lambda, ordnum, ngood, ntot, no_centroid);
398
399 good_frac = (double)ngood/(double) ntot;
400
401 if ( ngood >= 2 && good_frac >= tilt_frac){
402 /* Try to fit tilt */
403 xsh_msg_dbg_medium("---Clipping on tilt");
404
405 check( pos_vect = cpl_vector_new( ngood));
406 check( val_vect = cpl_vector_new( ngood));
407
408 j = 0;
409 for( i = 0 ; i< ntot; i++ ) {
410 if ( centers[i].good == 0){
411 cpl_vector_set( pos_vect, j, centers[i].xpos);
412 cpl_vector_set( val_vect, j, centers[i].ypos);
413 j++;
414 }
415 }
416
417 check( xsh_array_clip_poly1d( pos_vect, val_vect, tilt_kappa, tilt_niter,
418 tilt_frac, 1, &poly_tilt, chisq, &flags));
419
420 /* compute fwhm median on good lines */
421 j = 0;
422 k = 0;
423 XSH_CALLOC( median_data, double, ngood);
424
425 for( i = 0 ; i< ntot; i++ ) {
426 if ( centers[i].good == 0){
427 if (flags[j] == 0){
428 median_data[k] = centers[i].fwhm;
429 k++;
430 }
431 else{
432 centers[i].good = 2;
433 }
434 j++;
435 }
436 }
437
438 check( fwhm_median_vect = cpl_vector_wrap( k, median_data));
439
440 check( fwhm_median = cpl_vector_get_median( fwhm_median_vect));
441 *fwhm_center = fwhm_median;
442
443 /* get slope and center */
444 check( *good_center = cpl_polynomial_eval_1d( poly_tilt, x0, NULL));
445 check( *slope = cpl_polynomial_get_coeff( poly_tilt, &deg1));
446
447 xsh_msg_dbg_medium(" find_tilt : slope %f fwhm_median %f tilt at center %f %f ",
448 *slope, fwhm_median, x0, *good_center);
449
450 *nt = ntot;
451 *ng = ngood;
452
453#if REGDEBUG_WAVECAL
454 /* REGDEBUG */
455 {
456 FILE *fwhm_debug_file = NULL;
457 char name[256];
458 int idebug = 0;
459
460 sprintf( name, "fwhm_%d_%f.dat", ordnum, lambda);
461
462 if ( fwhm_debug_mode == NULL){
463 fwhm_debug_file = fopen( name, "w");
464 fwhm_debug_mode = "a+";
465 fprintf( fwhm_debug_file,
466 "# wavelength order xcen xcen-x0 ygauss ytilt fwhm area good\n");
467 }
468 else{
469 fwhm_debug_file = fopen( name, fwhm_debug_mode);
470 }
471 for( idebug=0; idebug< ntot; idebug++){
472 double fit;
473 double xpos;
474
475 xpos = xsh_round_double( centers[idebug].xpos );
476
477 fit = cpl_polynomial_eval_1d( poly_tilt, xpos, NULL);
478 fprintf( fwhm_debug_file, "%f %d %f %f %f %f %f %f %d\n",
479 lambda, ordnum,
480 xpos, xpos-x0,
481 centers[idebug].ypos, fit, centers[idebug].fwhm,
482 centers[idebug].area, centers[idebug].good);
483 }
484 fclose( fwhm_debug_file);
485 }
486#endif
487 }
488 else{
489 result = FIND_TILT_BAD_FIT;
490 xsh_msg( "Not enough points to do the fit: Good fraction points (%f < %f)",
491 good_frac, tilt_frac);
492 }
493
494 cleanup:
495 xsh_free_vector( &pos_vect);
496 xsh_free_vector(&val_vect);
497 xsh_free_polynomial( &poly_tilt);
498 XSH_FREE( flags);
499 xsh_unwrap_vector( &fwhm_median_vect);
500 XSH_FREE( median_data);
501 XSH_FREE( centers);
502 return result;
503}
504/*****************************************************************************/
505
506static float linear_interpol( float xa, float ya, float xb, float yb, float x )
507{
508 float xb_xa ;
509
510 xb_xa = xb-xa ;
511 return (ya*(xb-x))/(xb_xa) + (yb*(x-xa))/(xb_xa) ;
512
513}
514
515static float get_lambda( float * data, float x, float y, int nx, int ny )
516{
517 /* Should calculate lambda with interpolation, because in the wavemap
518 image there are only integer pixels
519 */
520 float lambda, lm, lp;
521 int ym, yp;
522 int xpos;
523
524 if(y>ny) return 0;
525
526 /* take pixels around value */
527 yp = (int) ceil( y);
528 ym = (int) floor( y);
529
530 xpos = (int) xsh_round_double( x);
531 lm = data[(ym-1)*nx+xpos];
532 lp = data[(yp-1)*nx+xpos];
533
534 // Now interpolate
535 lambda = linear_interpol( ym, lm, yp, lp, y ) ;
536 xsh_msg_dbg_high( " ym: %d, lm: %f, yp: %d, lp: %f, y: %f ==> %f",
537 ym, lm, yp, lp, y, lambda ) ;
538
539 return lambda ;
540}
541
542/*****************************************************************************/
556/*****************************************************************************/
557static void
558compute_specres( cpl_frame *wavemap_frame,
559 cpl_frame *disptab_frame,
560 xsh_instrument* instr,
561 xsh_linetilt_list *tilt_list,
562 int niter,
563 double kappa,
564 double frac,
565 double *specres_med,
566 double *specres_stdev)
567{
568 cpl_array *a_specres = NULL ;
569 int i=0, j=0, status=0 ;
570 xsh_dispersol_list* displist = NULL;
571 cpl_vector* positions = NULL;
572 const char *wavemap_name = NULL;
573 cpl_image *wavemap_img = NULL;
574 float *wavemap_data = NULL ;
575 int nx=0, ny=0;
576 int tilt_list_size =0;
577
578 XSH_ASSURE_NOT_NULL( tilt_list);
579
580 /* Load Wavemap image */
581 if (disptab_frame != NULL){
582 xsh_msg("Use the DISP_TAB");
583 check( displist = xsh_dispersol_list_load( disptab_frame, instr));
584 check( positions = cpl_vector_new(2));
585 }
586 else{
587 xsh_msg("Use the WAVE_MAP");
588 XSH_ASSURE_NOT_NULL( wavemap_frame);
589 check( wavemap_name = cpl_frame_get_filename( wavemap_frame));
590 check( wavemap_img = cpl_image_load( wavemap_name, CPL_TYPE_FLOAT, 0, 0));
591 check( wavemap_data = cpl_image_get_data_float( wavemap_img));
592 nx = cpl_image_get_size_x( wavemap_img);
593 ny = cpl_image_get_size_y( wavemap_img);
594 }
595
596 /* Create vector to house line's specres */
597 tilt_list_size = tilt_list->size;
598 check( a_specres = cpl_array_new( tilt_list_size, CPL_TYPE_DOUBLE));
599
600 /* Loop over Lines */
601 for( i=0 ; i< tilt_list_size; i++) {
602 double xc=0, yc=0, yc0=0, yc1=0;
603 double pix_width=0 ;
604 float lambda0=0, lambda1=0, lambdac=0;
605 double tilt =0.0;
606
607 /* get X, Y position of the center */
608 xc = tilt_list->list[i]->cenposx;
609 yc = tilt_list->list[i]->tilt_y;
610
611 pix_width = tilt_list->list[i]->deltay;
612 lambdac = tilt_list->list[i]->wavelength;
613 tilt = tilt_list->list[i]->tilt;
614
615 /* line is init to invalid data */
616 tilt_list->list[i]->specres = 0;
617 tilt_list->list[i]->flag = 1;
618
619 /* calculate lower/upper positions of mid width in pixels */
620 yc0 = yc - (pix_width/2.);
621 yc1 = yc + (pix_width/2.);
622
623 if ( displist != NULL){
624 int absorder = tilt_list->list[i]->order;
625 int iorder=0;
626
627 while( absorder != displist->list[iorder].absorder){
628 iorder++;
629 }
630 cpl_vector_set( positions, 0, xc);
631 cpl_vector_set( positions, 1, yc0);
632
633 check( lambda0 = xsh_dispersol_list_eval( displist,
634 displist->list[iorder].lambda_poly, positions));
635
636 cpl_vector_set( positions, 1, yc1);
637 check( lambda1 = xsh_dispersol_list_eval( displist,
638 displist->list[iorder].lambda_poly, positions));
639 }
640 else{
641 if( (yc0 >= 1) && (yc1 <= ny-1) ) {
642 // convert to lambda ==> get lambda from wavemap
643 lambda0 = get_lambda( wavemap_data, xc, yc0, nx,ny);
644 lambda1 = get_lambda( wavemap_data, xc, yc1, nx,ny);
645 }
646 }
647 xsh_msg_dbg_high("Lambda1: %f, Lambda0: %f, Lambdac: %f, Pix_width: %lf",
648 lambda0, lambda1, lambdac, pix_width);
649
650 if ( lambda0 != 0. && lambda1 != 0.) {
651 tilt_list->list[i]->flag = 0;
652 /* tilt = tan(alpha) cos(atan(alpha)) = 1/sqrt(1+tilt²) */
653 tilt_list->list[i]->specres = lambdac/fabs(lambda1 - lambda0)*
654 sqrt( 1+tilt*tilt);
655 }
656 check( cpl_array_set( a_specres, i, tilt_list->list[i]->specres));
657 }
658
659 check( *specres_med = cpl_array_get_median( a_specres));
660 check( *specres_stdev = cpl_array_get_stdev( a_specres));
661
662 xsh_msg("---Clipping on spectral resolution");
664 frac, specres_med, specres_stdev));
665
666 for( j=0; j < tilt_list_size; j++){
667 //double val=0.0;
668
669
670 //val = cpl_array_get_double( a_specres, j, &status);
671 if (status == 1){
672 tilt_list->list[j]->flag = 2;
673 }
674 }
675 cleanup:
676 xsh_free_array( &a_specres);
677 xsh_free_image( &wavemap_img);
678 xsh_free_vector( &positions);
679 xsh_dispersol_list_free( &displist);
680 return ;
681}
682/*****************************************************************************/
683
684/*****************************************************************************/
685/*****************************************************************************/
686static void
687set_qc_parameters( cpl_propertylist *tilt_header,
688 cpl_propertylist *shift_header,
689 double * ypos,
690 double * width,
691 double * intens,
692 xsh_linetilt_list *tilt_list,
694 int nlinecat,
695 int nb_lines,
696 double specres_med,
697 double specres_stdev,
698 float exptime)
699{
700 cpl_array * yarr = NULL, * warr = NULL, * iarr = NULL ;
701 double avg, med, std, wavg, wstd, iavg;//, iavg_norm;
702 int n=tilt_list->size;
703 int i=0;
704 int k=0;
705 check( yarr = cpl_array_wrap_double( ypos, n ) ) ;
706 check( warr = cpl_array_wrap_double( width, n ) ) ;
707 check( iarr = cpl_array_wrap_double( intens, n ) ) ;
708 /* remove from vector flagged values */
709 for(i=0;i<n;i++) {
710 if(tilt_list->list[i]->flag >0 ) {
711 cpl_array_set_invalid(yarr,i);
712 cpl_array_set_invalid(warr,i);
713 cpl_array_set_invalid(iarr,i);
714 k++;
715 }
716 }
717 check( avg = cpl_array_get_mean( yarr ) ) ;
718 check( med = cpl_array_get_median( yarr ) ) ;
719 check( std = cpl_array_get_stdev( yarr ) ) ;
720
721 check( wstd = cpl_array_get_stdev( warr ) ) ;
722 check( wavg = cpl_array_get_mean( warr ) ) ;
723
724 check( iavg = cpl_array_get_mean( iarr ) ) ;
725
726 xsh_msg("diffy avg=%16.14g med=%16.14g std=%16.14g",avg,med,std);
727 xsh_msg("FWHM avg=%16.14g std=%16.14g",wavg,wstd);
728 xsh_msg("nlineint iavg=%16.14g",iavg);
729
730 xsh_msg_dbg_low( "FWHM Avg: %.2lf, RMS: %.2lf", wavg, wstd ) ;
731
732 check( xsh_pfits_set_qc( tilt_header, &avg, QC_WAVECAL_DIFFYAVG,
733 instrument ) ) ;
734 check( xsh_pfits_set_qc( tilt_header, &med, QC_WAVECAL_DIFFYMED,
735 instrument ) ) ;
736 check( xsh_pfits_set_qc( tilt_header, &std, QC_WAVECAL_DIFFYSTD,
737 instrument ) ) ;
738 check( xsh_pfits_set_qc( tilt_header, &iavg, QC_WAVECAL_NLININT,
739 instrument));
740 /*if(exptime > 0){
741 check(iavg_norm = iavg/exptime);
742 check( xsh_pfits_set_qc( tilt_header, &iavg_norm,QC_WAVECAL_LININT_NORM,
743 instrument));
744 }*/
745 check( xsh_pfits_set_qc( tilt_header, &wavg, QC_WAVECAL_FWHMAVG,
746 instrument ) ) ;
747 check( xsh_pfits_set_qc( tilt_header, &wstd, QC_WAVECAL_FWHMRMS,
748 instrument ) ) ;
749 check( xsh_pfits_set_qc( tilt_header, &nlinecat, QC_WAVECAL_CATLINE,
750 instrument ) ) ;
751 check( xsh_pfits_set_qc( tilt_header, &nb_lines, QC_WAVECAL_MATCHLINE,
752 instrument ) ) ;
754 instrument));
755 check( xsh_pfits_set_qc( tilt_header, &specres_med, QC_RESOLMED, instrument));
756 if(isnan(specres_stdev)) {
757 /* dummy value to allow saving of product */
758 specres_stdev = -999;
759 }
760 check( xsh_pfits_set_qc( tilt_header, &specres_stdev, QC_RESOLRMS, instrument));
761
762
763
764 check( xsh_pfits_set_qc( shift_header, &avg, QC_WAVECAL_DIFFYAVG,
765 instrument ) ) ;
766 check( xsh_pfits_set_qc( shift_header, &med, QC_WAVECAL_DIFFYMED,
767 instrument ) ) ;
768 check( xsh_pfits_set_qc( shift_header, &std, QC_WAVECAL_DIFFYSTD,
769 instrument ) ) ;
770 check( xsh_pfits_set_qc( shift_header, &iavg, QC_WAVECAL_NLININT,
771 instrument));
772
773
774 cleanup:
775 xsh_unwrap_array( &yarr ) ;
776 xsh_unwrap_array( &warr ) ;
777 xsh_unwrap_array( &iarr ) ;
778 return ;
779}
780/*****************************************************************************/
781
782/*****************************************************************************/
783/*****************************************************************************/
784static void
785clean_arclist_data( cpl_frame *wavesol_frame,
786 cpl_frame *arclines_frame,
787 xsh_order_list *orders,
788 cpl_frame *config_model_frame,
789 cpl_frame *pre_frame,
790 cpl_frame *spectralformat_frame,
791 double **lambda,
792 double **n,
793 double **x,
794 double **y,
795 double **xmin,
796 double **xmax,
797 int *size,
798 double slit,
799 double slit_min,
800 double slit_max,
802{
803 xsh_arclist *arclist = NULL;
804 int arclist_size =0, global_size=0;
805 cpl_vector** spectral_tab = NULL;
806 xsh_wavesol * wsol = NULL;
807 xsh_xs_3 config_model;
808 int found_temp=true;
809 xsh_spectralformat_list *spectralformat_list = NULL;
810 int bad_positions=0;
811 int i, j, index_loc = 0;
812
813 check( arclist = xsh_arclist_load( arclines_frame));
814 check( spectralformat_list =
815 xsh_spectralformat_list_load( spectralformat_frame, instrument));
816 check( arclist_size = xsh_arclist_get_size( arclist));
817
818 if ( wavesol_frame != NULL) {
819 check( wsol = xsh_wavesol_load( wavesol_frame, instrument));
820 }
821 else if ( config_model_frame != NULL ) {
822 /* update cfg table for temperature */
823 check(xsh_model_temperature_update_frame( &config_model_frame, pre_frame,
824 instrument, &found_temp));
825 check( xsh_model_config_load_best( config_model_frame, &config_model));
827 }
828 else{
830 "Undefined solution type (POLY or MODEL).");
831 }
832
833 XSH_CALLOC( spectral_tab, cpl_vector*, arclist_size);
834
835 for( i=0; i< arclist_size; i++){
836 cpl_vector* res = NULL;
837 float lambdaARC;
838 int res_size;
839
840 check( lambdaARC = xsh_arclist_get_wavelength( arclist, i));
841 check( res = xsh_spectralformat_list_get_orders( spectralformat_list,
842 lambdaARC));
843
844 if ( res != NULL){
845 check( res_size = cpl_vector_get_size( res));
846 }
847 else{
848 res_size = 0;
849 check( xsh_arclist_reject(arclist,i));
850 }
851 global_size += res_size;
852 spectral_tab[i] = res;
853 }
854 xsh_msg (" Arc List x Order size = %d", global_size);
855
856 XSH_MALLOC( *lambda, double, global_size);
857 XSH_MALLOC( *n, double, global_size);
858 XSH_MALLOC( *x, double, global_size);
859 XSH_MALLOC( *xmin, double, global_size);
860 XSH_MALLOC( *xmax, double, global_size);
861 XSH_MALLOC( *y, double, global_size);
862
863 for( i=0; i< arclist_size; i++){
864 double comp_x=0.0, comp_y=0.0, comp_xmin=0.0, comp_ymin=0.0;
865 double comp_xmax=0.0, comp_ymax=0.0;
866 cpl_vector *spectral_res = NULL;
867 int spectral_res_size=0;
868 float lambdaARC;
869
870 check( lambdaARC = xsh_arclist_get_wavelength(arclist, i));
871 check( spectral_res = spectral_tab[i]);
872
873 if (spectral_res != NULL){
874 check( spectral_res_size = cpl_vector_get_size( spectral_res));
875 for( j=0; j< spectral_res_size; j++){
876 int absorder = 0;
877 int oidx;
878
879 check( absorder = (int) cpl_vector_get( spectral_res, j));
880 check( oidx = xsh_order_list_get_index_by_absorder( orders, absorder));
881 if (wsol != NULL){
882 check( comp_x = xsh_wavesol_eval_polx( wsol, lambdaARC, absorder,
883 slit));
884 check( comp_y = xsh_wavesol_eval_poly( wsol, lambdaARC, absorder,
885 slit));
886 check( comp_xmin = xsh_wavesol_eval_polx( wsol, lambdaARC, absorder,
887 slit_min));
888 check( comp_xmax = xsh_wavesol_eval_polx( wsol, lambdaARC, absorder,
889 slit_max));
890 }
891 else{
892 check( xsh_model_get_xy( &config_model, instrument,
893 lambdaARC, absorder, slit, &comp_x, &comp_y));
894 check( xsh_model_get_xy( &config_model, instrument,
895 lambdaARC, absorder, slit_min, &comp_xmin, &comp_ymin));
896 check( xsh_model_get_xy( &config_model, instrument,
897 lambdaARC, absorder, slit_max, &comp_xmax, &comp_ymax));
898 }
899
900 if ( (comp_y >= xsh_order_list_get_starty( orders, oidx )) &&
901 (comp_y <= xsh_order_list_get_endy( orders, oidx))) {
902 (*lambda)[index_loc] = lambdaARC;
903 (*n)[index_loc] = absorder;
904 (*x)[index_loc] = comp_x;
905 (*y)[index_loc] = comp_y;
906 if ( comp_xmin < comp_xmax){
907 (*xmin)[index_loc] = comp_xmin;
908 (*xmax)[index_loc] = comp_xmax;
909 xsh_msg_dbg_medium("The relation f(s)->X is increasing, X and s axes have same direction.");
910 }
911 else{
912 (*xmin)[index_loc] = comp_xmax;
913 (*xmax)[index_loc] = comp_xmin;
914 xsh_msg_dbg_medium("The relation f(s)->X is decreasing, X and s axes have opposite directions");
915 }
917 "lambda %f order %d slit %f smin %f smax %f x %f y %f xmin %f xmax %f",
918 lambdaARC, absorder, slit, slit_min, slit_max, comp_x, comp_y, comp_xmin, comp_xmax);
919 index_loc++;
920 }
921 else{
922 xsh_msg_dbg_low( "Ypix out of Order" ) ;
923 bad_positions++ ;
924 }
925 }
926 }
927 }
928 xsh_msg("Clean size %d (%d bad positions)", index_loc, bad_positions);
929 *size = index_loc;
930
931 cleanup:
932 xsh_arclist_free( &arclist);
933 xsh_wavesol_free( &wsol);
934 xsh_spectralformat_list_free( &spectralformat_list);
935
936 if ( spectral_tab != NULL){
937 for(i=0; i< arclist_size; i++){
938 xsh_free_vector( &spectral_tab[i]);
939 }
940 XSH_FREE( spectral_tab);
941 }
942 if ( cpl_error_get_code() != CPL_ERROR_NONE){
943 XSH_FREE( *lambda);
944 XSH_FREE( *n);
945 XSH_FREE( *x);
946 XSH_FREE( *xmin);
947 XSH_FREE( *xmax);
948 XSH_FREE( *y);
949 }
950 return;
951}
952/*****************************************************************************/
953
954
955/*****************************************************************************/
983/*****************************************************************************/
984static void
985xsh_follow_arclines( cpl_frame *pre_frame,
986 cpl_frame *arclines_frame,
987 cpl_frame *wavesol_frame,
988 cpl_frame *order_frame,
989 cpl_frame *spectralformat_frame,
990 cpl_frame * config_model_frame,
991 cpl_frame *wavemap_frame,
992 cpl_frame *disptab_frame,
993 xsh_follow_arclines_param *follow_param,
994 double slit,
995 double slit_min,
996 double slit_max,
997 const char* tag_id,
998 int ifu_flag,
1000 cpl_frame **tilt_frame,
1001 cpl_frame **shift_frame, const int clean_tmp)
1002{
1003 xsh_pre * pre = NULL;
1004 xsh_order_list *orders = NULL;
1005 double * vlambdadata = NULL, *vorderdata = NULL;
1006 double * vxthedata = NULL, * vythedata = NULL;
1007 double *vxmindata = NULL, * vxmaxdata = NULL;
1008 int sol_size = 0;
1009 const char* arclines_name = NULL;
1010 int s_n_low = 0, no_centroid = 0, bad_position = 0 , no_tilt=0;
1011
1012 XSH_GAUSSIAN_FIT fit_res ;
1013 xsh_shift_tab *shift_tab = NULL ;
1014 xsh_linetilt_list *tilt_list = NULL ;
1015 cpl_propertylist *tilt_header = NULL ;
1016 double *ydelta = NULL, *width = NULL, *intensity = NULL;
1017 int goodlines = 0 ;
1018
1019 int i;
1020
1021 const char * shift_tag = NULL ;
1022 char tag[256];
1023 char fname[256];
1024 int specres_niter=0;
1025 double specres_kappa = 0, specres_frac = 1;
1026 double specres_med = 0., specres_stdev = 0 ;
1027 cpl_vector *svect = NULL;
1028 int cat_lines, match_lines, found_lines;
1029
1030 XSH_ASSURE_NOT_NULL( pre_frame);
1031 XSH_ASSURE_NOT_NULL( arclines_frame);
1032 XSH_ASSURE_NOT_NULL( order_frame);
1033
1034 XSH_ASSURE_NOT_NULL( follow_param);
1036 XSH_ASSURE_NOT_NULL( spectralformat_frame);
1037 XSH_ASSURE_NOT_NULL( tag_id);
1038
1039 /* parameters */
1040 XSH_ASSURE_NOT_NULL( follow_param->specres_clipping);
1041 specres_niter = follow_param->specres_clipping->niter;
1042 specres_kappa = follow_param->specres_clipping->sigma;
1043 specres_frac = follow_param->specres_clipping->frac;
1044
1045 check( pre = xsh_pre_load( pre_frame, instrument));
1046 check( orders = xsh_order_list_load( order_frame, instrument));
1047 /* Load data */
1048
1050
1051 /* clean arc list */
1052 check( clean_arclist_data( wavesol_frame, arclines_frame, orders,
1053 config_model_frame, pre_frame, spectralformat_frame,
1054 &vlambdadata, &vorderdata,
1055 &vxthedata, &vythedata, &vxmindata, &vxmaxdata,
1056 &sol_size, slit, slit_min, slit_max, instrument));
1057
1058 /* Create shifttab */
1059 check( shift_tab = xsh_shift_tab_create( instrument));
1060
1061 XSH_CALLOC( ydelta, double, sol_size);
1062 XSH_CALLOC( width, double, sol_size);
1063 XSH_CALLOC( intensity, double, sol_size);
1064
1065 check( arclines_name = cpl_frame_get_filename( arclines_frame));
1066 check( tilt_header = cpl_propertylist_load( arclines_name, 0));
1067 check( tilt_list = xsh_linetilt_list_new( sol_size, tilt_header));
1068
1069 /* Loop over all clean arc lines:
1070 Loop over arclines
1071 get a lambda and slit
1072 while ( get the order from Theoretical map )
1073 compute X,Y position from wavesol
1074 compute center
1075 fit gaussian
1076 if fit OK then
1077 calculate tilt
1078 save
1079 fi
1080 endwhile
1081 */
1082
1083 for ( i = 0 ; i < sol_size ; i++ ) {
1084 double lambda, order;
1085 double xpix, ypix;
1086 int ok, ret, rej;
1087 int nt = 0, ng = 0 ;
1088 double s_n =0;
1089
1090 lambda = vlambdadata[i];
1091 order = vorderdata[i];
1092 xpix = vxthedata[i];
1093 ypix = vythedata[i];
1094
1095 xsh_msg_dbg_low( "LAMBDA %f, ORDER %f, SLIT %f, X %lf, Y %lf", lambda,
1096 order, slit, xpix, ypix );
1097 /* Detect Centroid */
1098 check( ok = detect_centroid( pre, lambda, order, xpix, ypix,
1099 follow_param, 1, &fit_res));
1100
1101 if ( ok != 0 ) {
1102 xsh_msg_dbg_low( " ******* Could NOT Fit Centroid" ) ;
1103 xsh_msg_dbg_low( "Lambda: %f, Order: %f - Line NOT Fitted",
1104 lambda, order);
1105 no_centroid++ ;
1106 }
1107 else{
1108 /* IF S/N is < limit, dont keep this line */
1109 if( (xpix<=pre->nx) && (xpix>0) && (ypix<=pre->ny) && (ypix>0)) {
1110 double flux_val, err_val;
1111
1112 check( flux_val = cpl_image_get( pre->data, xpix, ypix, &rej ));
1113 check( err_val = cpl_image_get( pre->data, xpix, ypix, &rej ));
1114 check(s_n = cpl_image_get( pre->data, xpix, ypix, &rej )/
1115 cpl_image_get( pre->errs, xpix, ypix, &rej));
1116
1117 if ( s_n < follow_param->s_n_min ) {
1118 xsh_msg_dbg_low( " %f/%f < %f => S/N too low, line skipped",
1119 flux_val, err_val, follow_param->s_n_min);
1120 s_n_low++;
1121 }
1122 else{
1123 double slope = 0.0, chisq = 0.0, xmin = 0.0, xmax = 0.0;
1124 double fwhm_center = 0.0;
1125 double good_center = 0;
1126
1127 /* Then fit tilt:
1128 Along the X axis, [xc-length,xc+length], detect
1129 all centroids.
1130 Then fit a line (ax+b) between all centroids.
1131 The 'a' should be 0 if the line is // x-axis
1132 */
1133 xmin = vxmindata[i];
1134 xmax = vxmaxdata[i];
1135 if ( (ret = find_tilt( fit_res.peakpos, xpix, lambda, order,
1136 follow_param, pre, orders,
1137 instrument, &slope, &chisq,
1138 &xmin, &xmax, &nt, &ng,
1139 &fwhm_center, &good_center,
1140 ifu_flag )) == 0 ) {
1141 int rej, ix, iy;
1142 xsh_linetilt * tilt_line = NULL ;
1143
1144 ix = xsh_round_double( xpix);
1145 iy = xsh_round_double( ypix);
1146
1147 /* Prepare for QC param */
1148 ydelta[goodlines] = good_center-ypix;
1149 width[goodlines] = CPL_MATH_FWHM_SIG*fit_res.sigma;
1150 intensity[goodlines] = cpl_image_get( pre->data, ix, iy, &rej);
1151 check( tilt_line = xsh_linetilt_new());
1152 tilt_line->wavelength = lambda;
1153 tilt_line->slit = slit;
1154 tilt_line->tilt = slope ;
1155 tilt_line->chisq = chisq;
1156 tilt_line->xmin = xmin ;
1157 tilt_line->xmax = xmax ;
1158 tilt_line->ntot = nt ;
1159 tilt_line->ngood = ng ;
1160 tilt_line->deltay = fwhm_center;
1161 tilt_line->name = NULL ;
1162 tilt_line->order = order ;
1163 tilt_line->cenposx = xpix;
1164 tilt_line->pre_pos_y = ypix;
1165 tilt_line->cenposy = fit_res.peakpos ;
1166 tilt_line->tilt_y = good_center;
1167 tilt_line->shift_y = good_center-ypix;
1168 tilt_line->specres = fwhm_center; /* in pixels ! */
1169 tilt_line->area = fit_res.area;
1170 tilt_line->intensity = cpl_image_get( pre->data, ix, iy, &rej);
1171 check( xsh_linetilt_list_add( tilt_list, tilt_line, goodlines));
1172 xsh_msg_dbg_low( "Lambda: %f, Order: %f - Line Fitted: %lf,%lf",
1173 lambda, order, xpix, tilt_line->cenposy);
1174 goodlines++;
1175 }
1176 else{
1177 xsh_msg_dbg_low( " Can't compute tilt, line skipped" );
1178 no_tilt++;
1179 }
1180 }
1181 }
1182 else{
1183 xsh_msg_dbg_low( " Bad position, line skipped" ) ;
1184 bad_position++;
1185 }
1186 }
1187 }
1188
1189 xsh_msg( "No centroid %d Bad Position: %d, S/N too low: %d No tilt %d",
1190 no_centroid, bad_position, s_n_low, no_tilt);
1191 if ( goodlines == 0 ){
1192 xsh_msg( "***** NO FITTED LINE !!!!" );
1193 }
1194 XSH_ASSURE_NOT_ILLEGAL( goodlines > 0 ) ;
1195 xsh_msg( " Good Fitted lines : %d/%d", goodlines, sol_size);
1196
1197 /* If we have produced a wavemap, we can compute the spectral resolution */
1198 if ( wavemap_frame != NULL || disptab_frame != NULL) {
1199 compute_specres( wavemap_frame, disptab_frame, instrument,
1200 tilt_list, specres_niter,
1201 specres_kappa, specres_frac,
1202 &specres_med, &specres_stdev);
1203 if ( cpl_error_get_code() != CPL_ERROR_NONE ) {
1204 xsh_msg( "ERROR while computing Spectral Resolution" ) ;
1205 cpl_error_reset() ;
1206 specres_med = 0. ;
1207 specres_stdev = 0. ;
1208 }
1209 }
1210 else {
1211 xsh_msg_error( "NO Wavemap (WAVE_MAP) or dispersion table (DISP_TAB), can't compute Spectral Resolution");
1212 }
1213
1214 cat_lines = sol_size;
1215 found_lines = goodlines;
1216 match_lines = tilt_list->size;
1217
1218 set_qc_parameters( tilt_header, shift_tab->header, ydelta, width, intensity,
1219 tilt_list,instrument, sol_size, goodlines,
1220 specres_med, specres_stdev,pre->exptime);
1222
1223
1224 sprintf( tag, "FOLLOW_LINETILT_%s_%s", tag_id,
1226 sprintf(fname,"%s.fits",tag);
1227
1228 check( *tilt_frame = xsh_linetilt_list_save( tilt_list, instrument, fname, tag,
1229 specres_kappa, specres_niter,pre->exptime));
1230
1231 sprintf( tag, "TILT_TAB_%s_%s", tag_id,
1233 check( cpl_frame_set_tag( *tilt_frame, tag));
1234
1235 xsh_add_temporary_file(fname) ;
1236 xsh_msg( "Tilt Frame saved" ) ;
1237
1238 /* Save the shift table */
1239 xsh_msg( "Saving shift table SLIT, %d lines", tilt_list->size);
1240
1241 check( svect = cpl_vector_wrap( tilt_list->size, ydelta));
1242
1243 if ( ifu_flag == CENTER_SLIT){
1244 check( shift_tab->shift_y = cpl_vector_get_median( svect));
1246 }
1247 else{
1248 check( shift_tab->shift_y_cen = cpl_vector_get_median( svect));
1249 sprintf( tag, "SHIFT_TAB_%s_%s", tag_id,
1251 shift_tag = tag;
1252 }
1253
1254 check( xsh_pfits_set_qc( shift_tab->header, &cat_lines, QC_WAVECAL_CATLINE,
1255 instrument ) ) ;
1256 check( xsh_pfits_set_qc( shift_tab->header, &match_lines, QC_WAVECAL_MATCHLINE,
1257 instrument ) ) ;
1258 check( xsh_pfits_set_qc( shift_tab->header, &found_lines, QC_WAVECAL_FOUNDLINE,
1259 instrument ) ) ;
1260 check( *shift_frame = xsh_shift_tab_save( shift_tab, shift_tag,clean_tmp));
1261
1262 cleanup:
1263 xsh_pre_free( &pre);
1264 xsh_order_list_free( &orders);
1265 xsh_linetilt_list_free( &tilt_list);
1266 xsh_shift_tab_free( &shift_tab);
1267 xsh_unwrap_vector( &svect);
1268
1269 XSH_FREE( vlambdadata);
1270 XSH_FREE( vorderdata);
1271 XSH_FREE( vxthedata);
1272 XSH_FREE( vxmindata);
1273 XSH_FREE( vxmaxdata);
1274 XSH_FREE( vythedata);
1275 XSH_FREE( ydelta);
1276 XSH_FREE( width);
1277 XSH_FREE( intensity);
1278 return ;
1279}
1280
1281void
1282xsh_follow_arclines_slit( cpl_frame *pre_frame,
1283 cpl_frame *arclines_frame,
1284 cpl_frame *wavesol_frame,
1285 cpl_frame *order_frame,
1286 cpl_frame *spectralformat_frame,
1287 cpl_frame * config_model_frame,
1288 cpl_frame *wavemap_frame,
1289 cpl_frame *slitmap_frame,
1290 cpl_frame *disptab_frame,
1291 xsh_follow_arclines_param *follow_param,
1293 cpl_frame **tilt_frame,
1294 cpl_frame **shift_frame)
1295{
1296 double smin=-6.0, smax=6.0;
1297
1298/*
1299 check( xsh_get_slit_edges( slitmap_frame, &smin, &smax,
1300 NULL, NULL, instrument));
1301*/
1302
1303 check( xsh_follow_arclines( pre_frame, arclines_frame, wavesol_frame,
1304 order_frame, spectralformat_frame,
1305 config_model_frame, wavemap_frame, disptab_frame,
1306 follow_param, 0.0, smin, smax, "SLIT",
1308 tilt_frame, shift_frame,0));
1309
1310 cleanup:
1311 return;
1312}
1313
1314/*****************************************************************************/
1315
1316/*****************************************************************************/
1339/*****************************************************************************/
1340void
1341xsh_follow_arclines_ifu( cpl_frame *pre_frame,
1342 cpl_frame *arclines_frame,
1343 cpl_frame *wavesol_frame,
1344 cpl_frame *order_frame,
1345 cpl_frame *spectralformat_frame,
1346 cpl_frame *config_model_frame,
1347 cpl_frame *wavemap_frame,
1348 cpl_frame *slitmap_frame,
1349 cpl_frame *disptab_frame,
1350 xsh_follow_arclines_param *follow_param,
1352 cpl_frameset *tilt_set,
1353 cpl_frame **shift_frame)
1354{
1355
1356
1357 int ifu, spos=0;
1358 double slit_min=-6, slitcen_lo=-2.0, slitcen_up=2.0, slit_max=6;
1359 char tag_id[256];
1360 cpl_frameset *shift_ifu_set = NULL;
1361 cpl_frame *shift_tab_frame = NULL;
1362 xsh_shift_tab *shift_tab_low = NULL;
1363 xsh_shift_tab *shift_tab_cen = NULL;
1364 xsh_shift_tab *shift_tab_up = NULL;
1365 const char* shift_tag = NULL;
1366 double slit_pos[4];
1367 double shift_center=0, shift_lower=0, shift_upper=0;
1368
1369 XSH_ASSURE_NOT_NULL( tilt_set);
1370 XSH_ASSURE_NOT_NULL( spectralformat_frame);
1371
1372 shift_ifu_set = cpl_frameset_new();
1373
1374 check( xsh_get_slit_edges( slitmap_frame, &slit_min, &slit_max,
1375 &slitcen_lo, &slitcen_up, instrument));
1376 slit_pos[0] = slit_min;
1377 slit_pos[1] = slitcen_lo;
1378 slit_pos[2] = slitcen_up;
1379 slit_pos[3] = slit_max;
1380
1381 /* Loop over 3 slitlets */
1382
1383 for( ifu = LOWER_IFU_SLITLET ; ifu <= UPPER_IFU_SLITLET; ifu++) {
1384 cpl_frame *tmp_tilt_frame = NULL;
1385 cpl_frame *tmp_shift_frame = NULL;
1386 double slitlet_min, slitlet_max, slitlet_center;
1387
1388 slitlet_min = slit_pos[spos];
1389 slitlet_max = slit_pos[spos+1];
1390 slitlet_center = (slitlet_min+slitlet_max)/2.0;
1391
1392 xsh_msg( "IFU Slitlet %s: center = %f arcsec; edges = [%f, %f] arcsec",
1393 SlitletName[ifu], slitlet_center, slitlet_min, slitlet_max);
1394
1395 sprintf( tag_id, "%s_IFU", SlitletName[ifu]) ;
1396
1397 check( xsh_follow_arclines( pre_frame, arclines_frame, wavesol_frame,
1398 order_frame, spectralformat_frame, config_model_frame, wavemap_frame,
1399 disptab_frame,
1400 follow_param, slitlet_center, slitlet_min, slitlet_max, tag_id, ifu, instrument,
1401 &tmp_tilt_frame, &tmp_shift_frame,1));
1402
1403 spos++;
1404 check( cpl_frameset_insert( tilt_set, tmp_tilt_frame));
1405 check( cpl_frameset_insert( shift_ifu_set, tmp_shift_frame));
1406 }
1407
1408 check( shift_tab_frame = cpl_frameset_get_frame( shift_ifu_set, 0));
1409 xsh_msg_dbg_medium( "Lower Shift tab frame %s",
1410 cpl_frame_get_filename( shift_tab_frame));
1411 check( shift_tab_low = xsh_shift_tab_load( shift_tab_frame, instrument));
1412
1413 check( shift_tab_frame = cpl_frameset_get_frame( shift_ifu_set, 1));
1414 xsh_msg_dbg_medium( "Center Shift tab frame %s",
1415 cpl_frame_get_filename( shift_tab_frame));
1416 check( shift_tab_cen = xsh_shift_tab_load( shift_tab_frame, instrument));
1417
1418 check( shift_tab_frame = cpl_frameset_get_frame( shift_ifu_set, 2));
1419 xsh_msg_dbg_medium( "Upper Shift tab frame %s",
1420 cpl_frame_get_filename( shift_tab_frame));
1421 check( shift_tab_up = xsh_shift_tab_load( shift_tab_frame, instrument));
1422
1423 /* Merge all the shift tab in the center shift tab */
1424 shift_lower = shift_tab_low->shift_y_cen;
1425 shift_center = shift_tab_cen->shift_y_cen;
1426 shift_upper = shift_tab_up->shift_y_cen;
1427
1428 xsh_msg_dbg_medium( "Measured shift for slitlet: LOWER %f, CENTER %f, UPPER %f",
1429 shift_lower, shift_center, shift_upper);
1430
1431 shift_tab_cen->shift_y_down = shift_lower;
1432 shift_tab_cen->shift_y_up = shift_upper;
1433
1435 check( *shift_frame = xsh_shift_tab_save( shift_tab_cen, shift_tag,0));
1436
1437
1438 cleanup:
1439 xsh_shift_tab_free( &shift_tab_low);
1440 xsh_shift_tab_free( &shift_tab_cen);
1441 xsh_shift_tab_free( &shift_tab_up);
1442 xsh_free_frameset( &shift_ifu_set);
1443 return;
1444}
1445
static double exptime
static int width
static xsh_instrument * instrument
void xsh_set_image_cpl_bpmap(cpl_image *image, cpl_image *bpmap, const int decode_bp)
float xsh_arclist_get_wavelength(xsh_arclist *list, int idx)
get wavelength of a line in the arcline list
void xsh_arclist_free(xsh_arclist **list)
free memory associated to a arclist
int xsh_arclist_get_size(xsh_arclist *list)
get size of arcline list
void xsh_arclist_reject(xsh_arclist *list, int idx)
reject a line from the list
xsh_arclist * xsh_arclist_load(cpl_frame *frame)
load an arcline list frame in arclist structure
xsh_dispersol_list * xsh_dispersol_list_load(cpl_frame *frame, xsh_instrument *instrument)
Load a dispersion list from a frame.
void xsh_dispersol_list_free(xsh_dispersol_list **list)
Free the dispersion list.
double xsh_dispersol_list_eval(xsh_dispersol_list *list, cpl_polynomial *poly, cpl_vector *pos)
Evaluate the polynomial according the binning.
void xsh_linetilt_list_free(xsh_linetilt_list **list)
free memory associated to a arclist
xsh_linetilt * xsh_linetilt_new(void)
void xsh_linetilt_list_add(xsh_linetilt_list *list, xsh_linetilt *line, int idx)
xsh_linetilt_list * xsh_linetilt_list_new(int size, cpl_propertylist *header)
cpl_frame * xsh_linetilt_list_save(xsh_linetilt_list *list, xsh_instrument *instr, const char *filename, const char *tag, const double kappa, const int niter, float exptime)
save a (ks clip clean) linetilt list to a frame
int xsh_order_list_get_index_by_absorder(xsh_order_list *list, double absorder)
xsh_order_list * xsh_order_list_load(cpl_frame *frame, xsh_instrument *instr)
load an order list from a frame
int xsh_order_list_get_starty(xsh_order_list *list, int i)
get position on Y axis of first pixel detected on order
int xsh_order_list_get_endy(xsh_order_list *list, int i)
get position on Y axis of last pixel detected on order
int xsh_order_list_get_order(xsh_order_list *list, int absorder)
void xsh_order_list_free(xsh_order_list **list)
free memory associated to an order_list
double xsh_order_list_eval(xsh_order_list *list, cpl_polynomial *poly, double y)
Evaluate an order list poly.
xsh_pre * xsh_pre_load(cpl_frame *frame, xsh_instrument *instr)
Load a xsh_pre structure from a frame.
Definition: xsh_data_pre.c:849
void xsh_pre_free(xsh_pre **pre)
Free a xsh_pre structure.
Definition: xsh_data_pre.c:823
cpl_frame * xsh_shift_tab_save(xsh_shift_tab *tab, const char *tag, const int clean_tmp)
xsh_shift_tab * xsh_shift_tab_create(xsh_instrument *instrument)
void xsh_shift_tab_free(xsh_shift_tab **tab)
Free memory associated to a the_arcline.
xsh_shift_tab * xsh_shift_tab_load(cpl_frame *frame, xsh_instrument *instr)
Load a shift table.
cpl_vector * xsh_spectralformat_list_get_orders(xsh_spectralformat_list *list, float lambda)
Returns list of absolute orders containing lambda.
xsh_spectralformat_list * xsh_spectralformat_list_load(cpl_frame *frame, xsh_instrument *instr)
Load a spectralformat list from a frame.
void xsh_spectralformat_list_free(xsh_spectralformat_list **list)
Free memory associated to an spactralformat_list.
xsh_wavesol * xsh_wavesol_load(cpl_frame *frame, xsh_instrument *instrument)
load a wavelength solution
double xsh_wavesol_eval_poly(xsh_wavesol *sol, double lambda, double order, double slit)
eval the polynomial solution in Y
void xsh_wavesol_free(xsh_wavesol **w)
free wavelength solution structure
double xsh_wavesol_eval_polx(xsh_wavesol *sol, double lambda, double order, double slit)
eval the polynomial solution in X
#define XSH_ASSURE_NOT_ILLEGAL(cond)
Definition: xsh_error.h:107
#define check(COMMAND)
Definition: xsh_error.h:71
#define XSH_ASSURE_NOT_ILLEGAL_MSG(cond, msg)
Definition: xsh_error.h:111
#define XSH_ASSURE_NOT_NULL(pointer)
Definition: xsh_error.h:99
static int detect_centroid(xsh_pre *pre, float lambda, int ordnum, double xpix, double ypix, xsh_follow_arclines_param *follow_param, int is_center, XSH_GAUSSIAN_FIT *fit_res)
static cpl_polynomial * get_slit_ifu_lo_poly(xsh_order *porder, int ifu_flag)
static cpl_polynomial * get_slit_ifu_up_poly(xsh_order *porder, int ifu_flag)
static void compute_specres(cpl_frame *wavemap_frame, cpl_frame *disptab_frame, xsh_instrument *instr, xsh_linetilt_list *tilt_list, int niter, double kappa, double frac, double *specres_med, double *specres_stdev)
static float get_lambda(float *data, float x, float y, int nx, int ny)
static float linear_interpol(float xa, float ya, float xb, float yb, float x)
static void xsh_follow_arclines(cpl_frame *pre_frame, cpl_frame *arclines_frame, cpl_frame *wavesol_frame, cpl_frame *order_frame, cpl_frame *spectralformat_frame, cpl_frame *config_model_frame, cpl_frame *wavemap_frame, cpl_frame *disptab_frame, xsh_follow_arclines_param *follow_param, double slit, double slit_min, double slit_max, const char *tag_id, int ifu_flag, xsh_instrument *instrument, cpl_frame **tilt_frame, cpl_frame **shift_frame, const int clean_tmp)
Detect and follow arc lines. Computes center, width and tilt parameters. The position of the center o...
void xsh_follow_arclines_slit(cpl_frame *pre_frame, cpl_frame *arclines_frame, cpl_frame *wavesol_frame, cpl_frame *order_frame, cpl_frame *spectralformat_frame, cpl_frame *config_model_frame, cpl_frame *wavemap_frame, cpl_frame *slitmap_frame, cpl_frame *disptab_frame, xsh_follow_arclines_param *follow_param, xsh_instrument *instrument, cpl_frame **tilt_frame, cpl_frame **shift_frame)
static int find_tilt(double yp0, double xc, float lambda, int ordnum, xsh_follow_arclines_param *follow_param, xsh_pre *pre, xsh_order_list *orders, xsh_instrument *instrument, double *slope, double *chisq, double *minx, double *maxx, int *nt, int *ng, double *fwhm_center, double *good_center, int ifu_flag)
static void set_qc_parameters(cpl_propertylist *tilt_header, cpl_propertylist *shift_header, double *ypos, double *width, double *intens, xsh_linetilt_list *tilt_list, xsh_instrument *instrument, int nlinecat, int nb_lines, double specres_med, double specres_stdev, float exptime)
static void clean_arclist_data(cpl_frame *wavesol_frame, cpl_frame *arclines_frame, xsh_order_list *orders, cpl_frame *config_model_frame, cpl_frame *pre_frame, cpl_frame *spectralformat_frame, double **lambda, double **n, double **x, double **y, double **xmin, double **xmax, int *size, double slit, double slit_min, double slit_max, xsh_instrument *instrument)
void xsh_follow_arclines_ifu(cpl_frame *pre_frame, cpl_frame *arclines_frame, cpl_frame *wavesol_frame, cpl_frame *order_frame, cpl_frame *spectralformat_frame, cpl_frame *config_model_frame, cpl_frame *wavemap_frame, cpl_frame *slitmap_frame, cpl_frame *disptab_frame, xsh_follow_arclines_param *follow_param, xsh_instrument *instrument, cpl_frameset *tilt_set, cpl_frame **shift_frame)
Detect and follow arc lines. Computes center, width and tilt parameters. The position of the center o...
@ LAMBDA_TOO_SMALL
@ LAMBDA_NOT_FOUND
@ LAMBDA_FOUND
@ FIND_TILT_BAD_EDGES
@ FIND_TILT_BAD_FIT
@ FIND_TILT_BAD_CENTER
@ FIND_TILT_UNKNOW_ORDER
@ FIND_TILT_CLIPPED
const char * xsh_instrument_arm_tostring(xsh_instrument *i)
Get the string associated with an arm.
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
#define xsh_msg_dbg_medium(...)
Definition: xsh_msg.h:44
#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
#define xsh_msg_dbg_low(...)
Definition: xsh_msg.h:48
#define xsh_msg_dbg_high(...)
Definition: xsh_msg.h:40
void xsh_pfits_set_qc_slit_width(cpl_propertylist *source_list, cpl_propertylist *set_list, xsh_instrument *instrument)
void xsh_pfits_set_qc(cpl_propertylist *plist, void *value, const char *kw, xsh_instrument *instrument)
void xsh_get_slit_edges(cpl_frame *slitmap_frame, double *sdown, double *sup, double *sldown, double *slup, xsh_instrument *instrument)
Trace slit edges in a master flat.
Definition: xsh_rectify.c:671
void xsh_unwrap_vector(cpl_vector **v)
Unwrap a vector and set the pointer to NULL.
Definition: xsh_utils.c:2345
void xsh_vector_fit_gaussian(cpl_vector *x, cpl_vector *y, XSH_GAUSSIAN_FIT *result)
set debug level
Definition: xsh_utils.c:3101
void xsh_free_polynomial(cpl_polynomial **p)
Deallocate a polynomial and set the pointer to NULL.
Definition: xsh_utils.c:2194
void xsh_free_vector(cpl_vector **v)
Deallocate a vector and set the pointer to NULL.
Definition: xsh_utils.c:2284
void xsh_free_image(cpl_image **i)
Deallocate an image and set the pointer to NULL.
Definition: xsh_utils.c:2116
int xsh_debug_level_get(void)
get debug level
Definition: xsh_utils.c:3142
void xsh_free_frameset(cpl_frameset **f)
Deallocate a frame set and set the pointer to NULL.
Definition: xsh_utils.c:2254
void xsh_free_array(cpl_array **m)
Deallocate an array and set the pointer to NULL.
Definition: xsh_utils.c:2299
void xsh_array_clip_poly1d(cpl_vector *pos_vect, cpl_vector *val_vect, double kappa, int niter, double frac_min, int deg, cpl_polynomial **polyp, double *chisq, int **flagsp)
clip outliers from a 1D poly fit
Definition: xsh_utils.c:6568
long xsh_round_double(double x)
Computes round(x)
Definition: xsh_utils.c:4383
void xsh_array_clip_median(cpl_array *array, double kappa, int niter, double frac_min, double *median, double *stdev)
median clip of an array
Definition: xsh_utils.c:6493
void xsh_unwrap_array(cpl_array **a)
Unwrap an array and set the pointer to NULL.
Definition: xsh_utils.c:2361
void xsh_add_temporary_file(const char *name)
Add temporary file to temprary files list.
Definition: xsh_utils.c:1432
unsigned int first
Definition: irplib_error.c:88
unsigned int last
Definition: irplib_error.c:89
xsh_dispersol * list
cpl_polynomial * lambda_poly
xsh_clipping_param * tilt_clipping
xsh_clipping_param * specres_clipping
xsh_linetilt ** list
xsh_order * list
cpl_polynomial * edguppoly
cpl_polynomial * edglopoly
cpl_polynomial * slicuppoly
cpl_polynomial * sliclopoly
cpl_image * qual
Definition: xsh_data_pre.h:71
float exptime
Definition: xsh_data_pre.h:92
cpl_image * data
Definition: xsh_data_pre.h:65
cpl_propertylist * data_header
Definition: xsh_data_pre.h:66
cpl_image * errs
Definition: xsh_data_pre.h:68
cpl_propertylist * header
int nx
double kappa
Definition: xsh_detmon_lg.c:81
int niter
Definition: xsh_detmon_lg.c:82
int n
Definition: xsh_detmon_lg.c:92
int ny
int order
Definition: xsh_detmon_lg.c:80
#define XSH_SHIFT_TAB_SLIT
Definition: xsh_dfs.h:577
#define XSH_SHIFT_TAB_IFU
Definition: xsh_dfs.h:578
#define XSH_GET_TAG_FROM_ARM(TAG, instr)
Definition: xsh_dfs.h:1548
@ CENTER_IFU_SLITLET
Definition: xsh_ifu_defs.h:42
@ LOWER_IFU_SLITLET
Definition: xsh_ifu_defs.h:42
@ UPPER_IFU_SLITLET
Definition: xsh_ifu_defs.h:42
@ CENTER_SLIT
Definition: xsh_ifu_defs.h:41
static const char * SlitletName[]
Definition: xsh_ifu_defs.h:48
cpl_error_code xsh_model_temperature_update_frame(cpl_frame **model_config_frame, cpl_frame *ref_frame, xsh_instrument *instrument, int *found_temp)
#define QC_WAVECAL_NLININT
Definition: xsh_pfits_qc.h:198
#define QC_RESOLRMS
Definition: xsh_pfits_qc.h:210
#define QC_WAVECAL_DIFFYSTD
Definition: xsh_pfits_qc.h:193
#define QC_WAVECAL_CATLINE
Definition: xsh_pfits_qc.h:204
#define QC_WAVECAL_FOUNDLINE
Definition: xsh_pfits_qc.h:205
#define QC_WAVECAL_DIFFYMED
Definition: xsh_pfits_qc.h:192
#define QC_WAVECAL_MATCHLINE
Definition: xsh_pfits_qc.h:206
#define QC_RESOLMED
Definition: xsh_pfits_qc.h:209
#define QC_WAVECAL_FWHMAVG
Definition: xsh_pfits_qc.h:197
#define QC_WAVECAL_FWHMRMS
Definition: xsh_pfits_qc.h:196
#define QC_WAVECAL_DIFFYAVG
Definition: xsh_pfits_qc.h:191
#define XSH_FREE(POINTER)
Definition: xsh_utils.h:92
@ XSH_DEBUG_LEVEL_HIGH
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