X-shooter Pipeline Reference Manual 3.8.15
xsh_detect_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 $
23 * $Revision: 1.166 $
24*/
25
26#ifdef HAVE_CONFIG_H
27#include <config.h>
28#endif
29
30/*---------------------------------------------------------------------------*/
38/*---------------------------------------------------------------------------*/
41/*---------------------------------------------------------------------------
42 Includes
43----------------------------------------------------------------------------*/
44
45#include <math.h>
46#include <xsh_drl.h>
47#include <xsh_pfits.h>
48#include <xsh_utils.h>
49#include <xsh_model_utils.h>
50#include <xsh_data_order.h>
51#include <xsh_error.h>
52#include <xsh_utils.h>
53#include <xsh_msg.h>
54#include <xsh_data_pre.h>
55#include <xsh_data_order.h>
56#include <xsh_data_the_map.h>
57#include <xsh_data_arclist.h>
58#include <xsh_data_wavesol.h>
59#include <xsh_data_resid_tab.h>
60#include <xsh_data_wavemap.h>
62#include <xsh_model_io.h>
63#include <xsh_model_kernel.h>
64#include <xsh_fit.h>
65#include <gsl/gsl_multifit.h>
66#include <cpl.h>
67
68/*---------------------------------------------------------------------------
69 Typedefs
70 ---------------------------------------------------------------------------*/
71static void theo_tab_filter( xsh_the_map *the_tab, xsh_arclist *arclist,
72 int* size, double **lambda, double **n, double **s, int **s_index,
73 double **xthe, double **ythe, int nb_pinhole);
74
75static void theo_tab_model( xsh_xs_3* config_model, xsh_arclist *arclist,
76 xsh_spectralformat_list *spectralformat_list, int* size, double **lambda,
77 double **n, double **s, double **sn,
78 int **s_index, double **xthe, double **ythe,
79 xsh_instrument* instr, int nb_pinhole);
80
81static void data_wavesol_fit_with_sigma( xsh_wavesol *wavesol,
82 double *A, double *lambda, double *n, double *s, int size,
83 int max_iter, double min_frac, double sigma, int* rejected);
84
85static int lines_filter_by_sn( xsh_pre* pre, double sn_ref, double x,
86 double y, double* sn);
87/*---------------------------------------------------------------------------
88 Functions prototypes
89 ---------------------------------------------------------------------------*/
90
91static int lines_filter_by_sn( xsh_pre* pre, double sn_ref, double x,
92 double y, double* sn)
93{
94 int res = 0;
95 int xpix, ypix, rej;
96 double flux, noise;
97
99
100 /* compute the S/N for the center pixel */
101 xpix =(int) rint(x);
102 ypix = (int)rint(y);
103
104 check( flux = cpl_image_get( pre->data, xpix, ypix, &rej));
105 check( noise = cpl_image_get( pre->errs, xpix, ypix, &rej));
106 *sn = flux / noise;
107 res = (*sn > sn_ref);
108 cleanup:
109 return res;
110}
111
112static void theo_tab_model( xsh_xs_3* config_model, xsh_arclist *arclist,
113 xsh_spectralformat_list *spectralformat_list, int* size, double **lambda,
114 double **n, double **s, double** sn,
115 int **s_index, double **xthe, double **ythe,
116 xsh_instrument* instr, int nb_pinhole)
117{
118 int i,j;
119 cpl_vector** spectral_tab = NULL;
120 int global_size = 0;
121 int arc_size=0;
122 int index_loc = 0;
123 int nb_slit = 1;
124
125 XSH_ASSURE_NOT_NULL( config_model);
126 XSH_ASSURE_NOT_NULL( arclist);
127 XSH_ASSURE_NOT_NULL( spectralformat_list);
128
129 check( arc_size = xsh_arclist_get_size( arclist));
130 XSH_CALLOC( spectral_tab, cpl_vector*, arc_size);
131
132 nb_slit = nb_pinhole;
133
134 XSH_REGDEBUG("nb pinhole %d ", nb_slit);
135
136 for( i=0; i< arc_size; i++){
137 cpl_vector* res = NULL;
138 float lambdaARC;
139 int res_size;
140
141 check( lambdaARC = xsh_arclist_get_wavelength(arclist, i));
142 check( res = xsh_spectralformat_list_get_orders(spectralformat_list,
143 lambdaARC));
144 if (res != NULL){
145 check( res_size = cpl_vector_get_size( res));
146 res_size *= nb_slit;
147 }
148 else{
149 res_size = 0;
150 check(xsh_arclist_reject(arclist,i));
151 }
152 global_size += res_size;
153 spectral_tab[i] = res;
154 }
155
156 XSH_MALLOC( *lambda, double, global_size);
157 XSH_MALLOC( *n, double, global_size);
158 XSH_MALLOC( *s, double, global_size);
159 XSH_MALLOC( *sn, double, global_size);
160 XSH_MALLOC( *s_index, int, global_size);
161 XSH_MALLOC( *xthe, double, global_size);
162 XSH_MALLOC( *ythe, double, global_size);
163
164
165 for( i=0; i< arc_size; i++){
166 double model_x, model_y;
167 cpl_vector *spectral_res = NULL;
168 int spectral_res_size;
169 float lambdaARC;
170
171 check( lambdaARC = xsh_arclist_get_wavelength(arclist, i));
172 check( spectral_res = spectral_tab[i]);
173
174 if (spectral_res != NULL){
175 check( spectral_res_size = cpl_vector_get_size( spectral_res));
176 for( j=0; j< spectral_res_size; j++){
177 int absorder = 0;
178 int islit;
179 double slit;
180
181 if (nb_slit > 1){
182
183 for(islit=0; islit < nb_slit; islit++){
184
185 slit = config_model->slit[islit];
186 check( absorder = (int) cpl_vector_get( spectral_res, j));
187
188/* check( xsh_model_get_xy( config_model, instr, 1.000274*lambdaARC, absorder, slit, */
189/* &model_x, &model_y)); */
190 check( xsh_model_get_xy( config_model, instr, lambdaARC, absorder, slit,
191 &model_x, &model_y));
192 (*lambda)[index_loc] = lambdaARC;
193 (*n)[index_loc] = absorder;
194 (*s_index)[index_loc] = islit;
195 (*s)[index_loc] = slit;
196 (*sn)[index_loc] = 0;
197 (*xthe)[index_loc] = model_x;
198 (*ythe)[index_loc] = model_y;
199 index_loc++;
200 }
201 }
202 else{
203 islit = 4;
204 slit = config_model->slit[islit];
205 check( absorder = (int) cpl_vector_get( spectral_res, j));
206 check( xsh_model_get_xy( config_model, instr, lambdaARC, absorder, slit,
207 &model_x, &model_y));
208 (*lambda)[index_loc] = lambdaARC;
209 (*n)[index_loc] = absorder;
210 (*s_index)[index_loc] = islit;
211 (*s)[index_loc] = slit;
212 (*sn)[index_loc] = 0;
213 (*xthe)[index_loc] = model_x;
214 (*ythe)[index_loc] = model_y;
215 index_loc++;
216 }
217/* xsh_msg_dbg_medium( "lambda %f order %d slit %f x %f y %f",lambdaARC, absorder, slit,
218 model_x, model_y); */
219 }
220 }
221 }
222 *size = global_size;
223
224 cleanup:
225 if ( spectral_tab != NULL){
226 for(i=0; i< arc_size; i++){
227 xsh_free_vector( &spectral_tab[i]);
228 }
229 XSH_FREE( spectral_tab);
230 }
231 if ( cpl_error_get_code() != CPL_ERROR_NONE){
232 XSH_FREE( *lambda);
233 XSH_FREE( *n);
234 XSH_FREE( *s);
235 XSH_FREE( *s_index);
236 XSH_FREE( *xthe);
237 XSH_FREE( *ythe);
238 }
239 return;
240}
241
242/*
243 @brief
244 Filters data from theoretical tab with arcline list and clean this list
245 @param[in] the_tab
246 The theoretical tab
247 @param[in] arclist
248 The arcline list
249 @param[out] size
250 Size of output data arrays
251 @param[out] lambda
252 Lambda array
253 @param[out] n
254 Order array
255 @param[out] s
256 Slit array
257 @param[out] xthe
258 Theoretical positions in X array
259 @param[out] ythe
260 Theoretical positions in Y array
261*/
262static void theo_tab_filter( xsh_the_map *the_tab, xsh_arclist *arclist,
263 int* size, double **lambda, double **n, double **s, int **s_index,
264 double **xthe, double **ythe, int nb_pinhole)
265{
266 int the_size = 0;
267 int arc_size = 0;
268 int i=0, j=0;
269 int sol_size = 0;
270 int nb_slit = 1;
271
272 XSH_ASSURE_NOT_NULL( the_tab);
273 XSH_ASSURE_NOT_NULL( arclist);
274
275 check( arc_size = xsh_arclist_get_size( arclist));
276 check( the_size = xsh_the_map_get_size( the_tab));
277
278 XSH_MALLOC( *lambda, double, the_size);
279 XSH_MALLOC( *n, double, the_size);
280 XSH_MALLOC( *s, double, the_size);
281 XSH_MALLOC( *s_index, int, the_size);
282 XSH_MALLOC( *xthe, double, the_size);
283 XSH_MALLOC( *ythe, double, the_size);
284
285 nb_slit = nb_pinhole;
286 XSH_REGDEBUG("nb pinhole %d ", nb_slit);
287
288 for(i=0; i< arc_size; i++){
289 float lambdaARC, lambdaTHE;
290 int nb_match = 0, max_match = 0;
291
292 check(lambdaARC = xsh_arclist_get_wavelength(arclist, i));
293 check(lambdaTHE = xsh_the_map_get_wavelength(the_tab, j));
294
295 xsh_msg_dbg_medium("LINETABLE : Line %d / %d Lambda %f",i+1, arc_size,
296 lambdaARC);
297
298 while ( (j < the_size-1) &&
299 ( (lambdaARC-lambdaTHE) > WAVELENGTH_PRECISION) ){
300 j++;
301 check(lambdaTHE = xsh_the_map_get_wavelength(the_tab, j));
302 }
303
304 xsh_msg_dbg_medium("THETABLE : Line %d / %d Lambda %f",j+1, the_size, lambdaTHE);
305
306 max_match = the_size-j;
307
308 while ( nb_match < max_match && fabs(lambdaARC-lambdaTHE) <= WAVELENGTH_PRECISION ){
309 double order, slit, xtheval, ytheval;
310 int islit;
311
312 check( slit = (double) xsh_the_map_get_slit_position( the_tab, j));
313 check( islit = xsh_the_map_get_slit_index( the_tab, j));
314
315 if ( (nb_slit > 1) || (islit == 4) ){
316
317 check( order = (double)xsh_the_map_get_order(the_tab, j));
318 check( xtheval = xsh_the_map_get_detx(the_tab, j));
319 check( ytheval = xsh_the_map_get_dety(the_tab, j));
320 (*lambda)[sol_size] = lambdaTHE;
321 (*n)[sol_size] = order;
322 (*s)[sol_size] = slit;
323 (*s_index)[sol_size] = islit;
324 (*xthe)[sol_size] = xtheval;
325 (*ythe)[sol_size] = ytheval;
326
327 XSH_REGDEBUG("sol_size %d order %f slit %f lambda %f",
328 sol_size, (*n)[sol_size], (*s)[sol_size], (*lambda)[sol_size]);
329
330 sol_size++;
331 nb_match++;
332 }
333 if ( j < (the_size-1) ){
334 j++;
335 check( lambdaTHE = (double) xsh_the_map_get_wavelength( the_tab, j));
336 }
337
338 }
339
340 if (nb_match == 0) {
341 check(xsh_arclist_reject(arclist, i));
342 }
343 }
344 *size = sol_size;
345
346 cleanup:
347 if ( cpl_error_get_code() != CPL_ERROR_NONE){
348 XSH_FREE( *lambda);
349 XSH_FREE( *n);
350 XSH_FREE( *s);
351 XSH_FREE( *s_index);
352 XSH_FREE( *xthe);
353 XSH_FREE( *ythe);
354 }
355 return;
356}
357
358/*
359 @brief
360 Fit a polynomial like A = f(lambda,n,s) using sigma clipping
361 @param[in] wavesol
362 the wave solution structure
363 @param[in] A
364 The A array data
365 @param[in] lambda
366 The lambda array data
367 @param[in] n
368 The order array data
369 @param[in] s
370 The slit array data
371 @param[in] size
372 The size of A, n and s arrays
373 @param[in] max_iter
374 The maximum number of iterations for sigma clipping
375 @param[in] min_frac
376 The minimum fraction allowed for sigma clipping
377 @param[in] sigma
378 The sigma for sigma clipping
379 @param[out] rejected
380 A pre allocated array which contains rejected lines
381*/
383 double *A, double *lambda, double *n, double *s, int size,
384 int max_iter, double min_frac, double sigma, int* rejected)
385{
386 int nbiter = 0;
387 float frac = 1;
388 int nbrejected = 0;
389 int new_rejected = 1;
390 cpl_polynomial *fitpoly = NULL;
391 cpl_vector *dispy = NULL;
392 int * idx = NULL;
393 int index_size = 0;
394 int i;
395
396 XSH_ASSURE_NOT_NULL( wavesol);
398 XSH_ASSURE_NOT_NULL( lambda);
402
403 XSH_CALLOC( idx, int, size);
404 check(fitpoly = xsh_wavesol_get_poly( wavesol));
405 check(xsh_wavesol_compute( wavesol, size, A,
406 &(wavesol->min_y), &(wavesol->max_y), lambda, n, s, fitpoly));
407
408 index_size = size;
409 for(i=0; i< index_size; i++){
410 idx[i] = i;
411 }
412
413 xsh_msg( "Fit wavesol with sigma clipping");
414
415 while (nbiter < max_iter && frac > min_frac && new_rejected > 0){
416 double sigma_med;
417
418 new_rejected = 0;
419 xsh_msg_dbg_high( " *** NBITER = %d / %d ***", nbiter+1, max_iter);
420 dispy = cpl_vector_new( index_size);
421
422 for(i = 0; i < index_size; i ++){
423 double soly;
424 double diffy;
425 check( soly = xsh_wavesol_eval_poly(wavesol, lambda[i], n[i], s[i]));
426 diffy = A[i]-soly;
427 check( cpl_vector_set( dispy, i, diffy));
428
429 }
430 check( sigma_med = cpl_vector_get_stdev( dispy));
431 xsh_msg_dbg_high(" sigma %f SIGMA MEDIAN = %f", sigma, sigma_med);
432
433 for(i = 0; i < index_size; i++){
434 if ( fabs(cpl_vector_get( dispy, i)) > (sigma * sigma_med) ){
435 new_rejected++;
436 rejected[idx[i]] = 1;
437 }
438 else{
439 idx[i-new_rejected] = idx[i];
440 A[i-new_rejected] = A[i];
441 lambda[i-new_rejected] = lambda[i];
442 n[i-new_rejected] = n[i];
443 s[i-new_rejected] = s[i];
444 }
445 }
446 xsh_free_vector(&dispy);
447 index_size -= new_rejected;
448 nbrejected += new_rejected;
449
450 frac = 1 - ( (float) nbrejected / (float) size);
451 xsh_msg_dbg_high(" NB_REJECT = %d", nbrejected);
452 xsh_msg_dbg_high(" GOOD PIXEL FRACTION = %f", frac);
453
454 check( xsh_wavesol_compute(wavesol, index_size, A,
455 &(wavesol->min_y), &(wavesol->max_y), lambda, n, s, fitpoly));
456 nbiter++;
457 }
458
459 cleanup:
460 XSH_FREE( idx);
461 xsh_free_vector(&dispy);
462 return;
463}
464
465/*---------------------------------------------------------------------------
466 Implementation
467 ---------------------------------------------------------------------------*/
468
506void
507xsh_detect_arclines_dan( cpl_frame *frame,
508 cpl_frame *theo_tab_frame,
509 cpl_frame *arc_lines_tab_frame,
510 cpl_frame* wave_tab_guess_frame,
511 cpl_frame *order_tab_recov_frame,
512 cpl_frame *config_model_frame,
513 cpl_frame *spectralformat_frame,
514 cpl_frame **resid_tab_orders_frame,
515 cpl_frame **arc_lines_clean_tab_frame,
516 cpl_frame **wave_tab_frame,
517 cpl_frame **resid_tab_frame,
518 xsh_sol_wavelength solwave_type,
521 xsh_instrument *instr,
522 const char* rec_id,
523 const int clean_tmp,
524 const int resid_tab_name_sw)
525{
526
527 xsh_the_map *themap = NULL;
528 xsh_resid_tab *resid = NULL;
529 xsh_resid_tab *resid_orders = NULL;
530 xsh_arclist *arclist = NULL;
531 xsh_pre *pre = NULL;
532 cpl_polynomial *fit2dx = NULL;
533 cpl_propertylist* header = NULL;
534 int * sort = NULL;
535 double *vlambdadata = NULL, *vsdata = NULL, *vorderdata = NULL, *vsndata=NULL;
536 int *vsindexdata = NULL;
537 double *vxthedata = NULL, *vythedata = NULL;
538 double *corr_x = NULL, *corr_y = NULL;
539 double *gaussian_pos_x = NULL, *gaussian_pos_y = NULL,
540 *gaussian_sigma_x = NULL, *gaussian_sigma_y = NULL,
541 *gaussian_fwhm_x = NULL, *gaussian_fwhm_y = NULL, *gaussian_norm=NULL;
542 int* flag=NULL;
543 cpl_array *gaussian_sigma_x_arr =NULL;
544 cpl_array *gaussian_sigma_y_arr =NULL;
545 cpl_array *gaussian_fwhm_x_arr =NULL;
546 cpl_array *gaussian_fwhm_y_arr =NULL;
547
548
549 double *diffx = NULL, *diffy = NULL;
550 double *diffxmean = NULL, *diffymean = NULL, *diffxsig = NULL, *diffysig = NULL;
551 xsh_order_list *order_tab_recov = NULL;
552 /* Others */
553 int nlinematched, nlinecat_clean=0, nlinecat;
554 int i, sol_size = 0;
555 xsh_wavesol* wave_table = NULL;
556 char wave_table_name[256];
557 const char* wave_table_tag = NULL;
558 xsh_wavesol* wave_tab_guess = NULL;
559 xsh_xs_3 config_model;
560 xsh_xs_3* p_xs_3;
561 xsh_spectralformat_list *spectralformat_list = NULL;
562 int lines_not_in_image = 0;
563 int lines_not_good_sn = 0;
564 int lines_not_gauss_fit = 0;
565 int lines_not_valid_pixels = 0;
566 int lines_too_few_ph=0;
567
568 int* rejected = NULL;
569 int nb_rejected = 0;
570 int solution_type = 0;
571 const char* solution_type_name[2] = { "POLY", "MODEL"};
572 int detection_mode = 0;
573 const char* detection_mode_name[3] = { "NORMAL", "CORRECTED", "RECOVER"};
574 int nb_pinhole;
575 char fname[256];
576 char rname[256];
577 char rtag[256];
578 const char* tag=NULL;
579
580 char dpr_type[256];
581 char type[256];
582 cpl_table* wave_trace=NULL;
583
584
585 cpl_table* cfg_tab=NULL;
586 const char* cfg_name=NULL;
587 char new_name[256];
588 char basename[256];
589 cpl_propertylist* plist=NULL;
590
591 int min_slit_match=4;
592 double line_devs=2.5;
593 int found_temp=true;
594 int starti=0;
595 int nbflag=0;
596 int newi=0;
597 int j,k;
598
599
600 /* check input parameters */
601 XSH_ASSURE_NOT_NULL( frame);
602 XSH_ASSURE_NOT_NULL( arc_lines_tab_frame);
603 XSH_ASSURE_NOT_NULL( spectralformat_frame);
604 XSH_ASSURE_NOT_NULL( instr);
607
608 /* Ouput parameters */
609 XSH_ASSURE_NOT_NULL( arc_lines_clean_tab_frame);
610 XSH_ASSURE_NOT_NULL( resid_tab_frame);
611
612
613 /* Debug parameters */
614 xsh_msg_dbg_medium("---Detect arclines parameters");
615 xsh_msg_dbg_medium("min_sn %f fit-window-half-size %d", da->min_sn, da->fit_window_hsize);
616 xsh_msg_dbg_medium("Clipping sigma %f niter %d frac %f", dac->sigma, dac->niter, dac->frac);
617
618 check( pre = xsh_pre_load (frame, instr));
620 "fail to find number of pinholes info in frame %s",
621 cpl_frame_get_filename(frame) );
622 check(strcpy(dpr_type,xsh_pfits_get_dpr_type(pre->data_header)));
623
624 if(strstr(dpr_type,"FMTCHK") != NULL) {
625 strcpy(type,"FMTCHK_");
626 } else if(strstr(dpr_type,"WAVE") != NULL) {
627 strcpy(type,"WAVE_");
628 } else if(strstr(dpr_type,"WAVE") != NULL) {
629 strcpy(type,"ARC_");
630 } else {
631 strcpy(type,"");
632 }
633
634
635
636 check( arclist = xsh_arclist_load( arc_lines_tab_frame));
637 check( spectralformat_list = xsh_spectralformat_list_load(
638 spectralformat_frame, instr));
639 check(nlinecat = xsh_arclist_get_size( arclist));
640
641 /* sort data by wavelength */
643
644
645 /* First define the solution type of detect_arclines */
646 if ( theo_tab_frame != NULL) {
647 /* poly mode */
648 XSH_ASSURE_NOT_ILLEGAL( config_model_frame == NULL);
649 solution_type = XSH_DETECT_ARCLINES_TYPE_POLY;
650 check( themap = xsh_the_map_load( theo_tab_frame));
652 check( theo_tab_filter( themap, arclist, &sol_size, &vlambdadata,
653 &vorderdata, &vsdata, &vsindexdata, &vxthedata, &vythedata, nb_pinhole));
654
655 }
656 else if ( config_model_frame != NULL){
657 /* physical model mode */
658 solution_type = XSH_DETECT_ARCLINES_TYPE_MODEL;
659 p_xs_3=&config_model;
660
661
662 //Make sure to modify only current copy
663 check(cfg_name=cpl_frame_get_filename(config_model_frame));
664 check(plist=cpl_propertylist_load(cfg_name,0));
665 check(cfg_tab=cpl_table_load(cfg_name,1,0));
666 //check(basename=xsh_get_basename(cfg_name));
668 sprintf(basename,"%s.fits",tag);
669 sprintf(new_name,"local_%s",basename);
670 check( xsh_pfits_set_pcatg(plist,tag));
671 check(cpl_table_save(cfg_tab,plist,NULL,new_name,CPL_IO_DEFAULT));
672 xsh_add_temporary_file(new_name);
673 check(cpl_frame_set_filename(config_model_frame,new_name));
674 xsh_free_table(&cfg_tab);
675 xsh_free_propertylist(&plist);
676
677 check( xsh_model_config_load_best( config_model_frame, &config_model));
678 XSH_REGDEBUG("load config model ok");
679 /* update cfg table frame and corresponding model structure for temperature */
680 check(xsh_model_temperature_update_frame(&config_model_frame,frame,
681 instr,&found_temp));
682 check(xsh_model_temperature_update_structure(&config_model,frame,instr));
683 /* predict line positions on detector */
684 check( theo_tab_model( &config_model, arclist, spectralformat_list,
685 &sol_size, &vlambdadata, &vorderdata, &vsdata,
686 &vsndata, &vsindexdata, &vxthedata, &vythedata,
687 instr, nb_pinhole));
688
689 }
690 else{
692 "Undefined solution type (POLY or MODEL). See your input sof");
693 }
694
695
696 xsh_spectralformat_list_free(&spectralformat_list);
697 xsh_msg_dbg_high("Solution type %s", solution_type_name[solution_type]);
698
699 /* Define the mode of use */
700 if ( wave_tab_guess_frame == NULL){
701 if ( order_tab_recov_frame == NULL){
702 detection_mode = XSH_DETECT_ARCLINES_MODE_NORMAL;
703 }
704 else{
705 detection_mode = XSH_DETECT_ARCLINES_MODE_RECOVER;
706 check( order_tab_recov = xsh_order_list_load( order_tab_recov_frame,
707 instr));
708 }
709 }
710 else {
711 if ( order_tab_recov_frame == NULL){
712 detection_mode = XSH_DETECT_ARCLINES_MODE_CORRECTED;
713
714 check( wave_tab_guess = xsh_wavesol_load( wave_tab_guess_frame,
715 instr ));
716 xsh_msg( "BinX,Y: %d, %d", wave_tab_guess->bin_x,
717 wave_tab_guess->bin_y ) ;
718 }
719 else{
720 detection_mode = -1;
721 }
722 }
723
725 && detection_mode <= XSH_DETECT_ARCLINES_MODE_RECOVER);
726
727 xsh_msg( "Detection mode : %s", detection_mode_name[detection_mode]);
728
729 /* allocate some memory */
730 XSH_MALLOC( corr_x, double, sol_size);
731 XSH_MALLOC( corr_y, double, sol_size);
732
733 XSH_MALLOC( gaussian_pos_x, double, sol_size);
734 XSH_MALLOC( gaussian_pos_y, double, sol_size);
735 XSH_MALLOC( gaussian_sigma_x, double, sol_size);
736 XSH_MALLOC( gaussian_sigma_y, double, sol_size);
737
738 XSH_MALLOC( gaussian_fwhm_x, double, sol_size);
739 XSH_MALLOC( gaussian_fwhm_y, double, sol_size);
740 XSH_MALLOC( gaussian_norm, double, sol_size);
741
742 XSH_MALLOC( diffx, double, sol_size);
743 XSH_MALLOC( diffy, double, sol_size);
744
745 XSH_MALLOC( diffxmean, double, sol_size);
746 XSH_MALLOC( diffymean, double, sol_size);
747
748 XSH_MALLOC( diffxsig, double, sol_size);
749 XSH_MALLOC( diffysig, double, sol_size);
750 XSH_MALLOC( flag, int, sol_size);
751
752 /* initialize array value to default not to have garbage on resid tab*/
753
754 memset(corr_x,0,sizeof(double)*sol_size);
755 memset(corr_y,0,sizeof(double)*sol_size);
756
757 memset(gaussian_pos_x,0,sizeof(double)*sol_size);
758 memset(gaussian_pos_y,0,sizeof(double)*sol_size);
759 memset(gaussian_sigma_x,0,sizeof(double)*sol_size);
760 memset(gaussian_sigma_y,0,sizeof(double)*sol_size);
761 memset(gaussian_fwhm_x,0,sizeof(double)*sol_size);
762 memset(gaussian_fwhm_y,0,sizeof(double)*sol_size);
763 memset(gaussian_norm,0,sizeof(double)*sol_size);
764 memset(diffx,0,sizeof(double)*sol_size);
765 memset(diffy,0,sizeof(double)*sol_size);
766 memset(diffxmean,0,sizeof(double)*sol_size);
767 memset(diffymean,0,sizeof(double)*sol_size);
768
769 memset(diffxsig,0,sizeof(double)*sol_size);
770 memset(diffysig,0,sizeof(double)*sol_size);
771
772 memset(flag,0,sizeof(int)*sol_size);
773
774
775 /* init */
776 nlinematched = 0;
777
779 xsh_msg_dbg_medium( "USE method GAUSSIAN");
780 }
781 else{
782 xsh_msg_dbg_medium( "USE method BARYCENTER");
783 }
784
785
786 /* perform line detection and sn thresholding */
787 for(i=0; i< sol_size; i++){
788 int xpos = 0, ypos = 0;
789 int slit_index;
790 double corrxv = 0.0, corryv = 0.0;
791 double lambda, order, slit, xthe, ythe,sn=0;
792 double xcor, ycor, x, y, sig_x, sig_y, norm, fwhm_x, fwhm_y;
793 /* double max_off=1.0; */
794
795 lambda = vlambdadata[i];
796 order = vorderdata[i];
797 slit = vsdata[i];
798
799 if(vsndata != NULL) {
800 sn = vsndata[i];
801 }
802
803 slit_index = vsindexdata[i];
804 xthe = vxthedata[i];
805 ythe = vythedata[i];
806 xsh_msg_dbg_high( "THE LAMBDA %f ORDER %f SLIT %f X %f Y %f", lambda,
807 order, slit, xthe, ythe);
808
809 /* Theoretical using pixel FITS convention */
810 xcor = xthe;
811 ycor = ythe;
812 xsh_msg_dbg_high("THE x%f y %f", xthe, ythe);
813
814 /* Find X position from RECOVER ORDER TAB */
815 if ( order_tab_recov != NULL){
816 int iorder = 0;
817 cpl_polynomial *center_poly_recov = NULL;
818
820 order_tab_recov, order));
821 center_poly_recov = order_tab_recov->list[iorder].cenpoly;
822 check( xcor = cpl_polynomial_eval_1d( center_poly_recov, ycor, NULL));
823 corrxv = xcor-xthe;
824 }
825 /* apply x and y translation from predict data to the model */
826 if ( wave_tab_guess != NULL){
827 double diffxv, diffyv;
828
829 check( diffxv = xsh_wavesol_eval_polx( wave_tab_guess, lambda,
830 order, 0));
831 check( diffyv = xsh_wavesol_eval_poly( wave_tab_guess, lambda,
832 order, 0));
833
834 xcor = xthe+diffxv;
835 ycor = ythe+diffyv;
836 corrxv = diffxv;
837 corryv = diffyv;
838 }
839 /* Corrected position using pixel FITS convention */
840 xsh_msg_dbg_high("x %f y %f CORR_x %f CORR_y %f", xcor, ycor,
841 corrxv, corryv);
842
843 if ( xcor >=0.5 && ycor >=0.5 &&
844 xcor < (pre->nx+0.5) && ycor < (pre->ny+0.5)){
845 int best_med_res = 0;
846 check ( best_med_res = xsh_pre_window_best_median_flux_pos( pre,
847 (int)(xcor+0.5)-1, (int)(ycor+0.5)-1,
848 da->search_window_hsize, da->running_median_hsize, &xpos, &ypos));
849
850 if (best_med_res != 0){
851 lines_not_valid_pixels++;
852 flag[i]=1;
853 xsh_msg_dbg_medium("Not valid pixels for this line %d",i);
854 continue;
855 }
856 /* Put positions in CPL coordinates */
857 xpos = xpos+1;
858 ypos = ypos+1;
859 xsh_msg_dbg_high("ADJ x %d y %d", xpos, ypos);
860
862 cpl_image_fit_gaussian(pre->data, xpos, ypos, 1+2*da->fit_window_hsize,
863 &norm, &x, &y, &sig_x, &sig_y, &fwhm_x, &fwhm_y);
864 xsh_msg_dbg_high("GAUS x %f y %f %d %d %d", x, y,da->fit_window_hsize, da->search_window_hsize, da->running_median_hsize);
865
866
867 /* IN CASE OF usage of cpl_fit_image_gauss we need this code
868
869 cpl_array* parameters=cpl_array_new(7, CPL_TYPE_DOUBLE);
870 cpl_array* err_params=cpl_array_new(7, CPL_TYPE_DOUBLE);
871 cpl_array* fit_params=cpl_array_new(7, CPL_TYPE_INT);
872 int kk=0;
873 //All parameter should be fitted
874 for (kk = 0; kk < 7; kk++)
875 cpl_array_set(fit_params, kk, 1);
876 double rms=0;
877 double red_chisq=0;
878 cpl_matrix* covariance=NULL;
879 cpl_matrix* phys_cov=NULL;
880 double major=0;
881 double minor=0;
882 double angle=0;
883 int size=1+2*da->fit_window_hsize;
884 cpl_fit_image_gaussian(pre->data, pre->errs, xpos, ypos,size,size,
885 parameters,err_params,fit_params,
886 &rms,&red_chisq,&covariance,
887 &major,&minor,&angle,&phys_cov);
888
889 double rho=0;
890 double back=cpl_array_get(parameters,0,NULL);
891 norm=cpl_array_get(parameters,1,NULL);
892 rho=cpl_array_get(parameters,2,NULL);
893
894 x=cpl_array_get(parameters,3,NULL);
895 y=cpl_array_get(parameters,4,NULL);
896 sig_x=cpl_array_get(parameters,5,NULL);
897 sig_y=cpl_array_get(parameters,6,NULL);
898 fwhm_x=sig_x*CPL_MATH_FWHM_SIG;
899 fwhm_y=sig_y*CPL_MATH_FWHM_SIG;
900
901 xsh_free_array(&parameters);
902 xsh_free_array(&err_params);
903 xsh_free_array(&fit_params);
904 */
905
906
907
908 }
909 else{
910 xsh_image_find_barycenter( pre->data, xpos, ypos,
911 1+2*da->fit_window_hsize,
912 &norm, &x, &y, &sig_x, &sig_y, &fwhm_x, &fwhm_y);
913 xsh_msg_dbg_high("BARY x %f y %f", x, y);
914 }
915
916 /* (PBR) Adjust 1ph centroids to be consistent with 9ph
917 - This is cleaner than adjusting the ph position in the cfg because:
918 a. we do it just once here -> it gets used in the matching, then
919 propogates to the model optimisation
920 b. We don't need to know if the input config was itself optimised
921 to 9ph or 1ph (i.e. do we need the shift or not)
922 c. Potentially this can be used for poly-mode too, though I
923 restrict it to phys mod below (feel free to change it)
924 - The NIR value of 0.125pix for x I got by using a config optimised
925 to 9ph data (i.e. 2dmap) in xsh_predict and checking the input
926 offsets. For y the offset was less than the s.d.
927 - Same procedure for UVB. Also same for VIS, but value was ~0
928 - A correction from 9ph to full slit centre is still required*/
929 if (nb_pinhole==1 && solution_type==XSH_DETECT_ARCLINES_TYPE_MODEL) {
930 if (p_xs_3->arm==2) {
931 x=x+0.125;
932 }
933 else if (p_xs_3->arm==0) {
934 x=x-0.51;
935 }
936 }
937
938 if( cpl_error_get_code() == CPL_ERROR_NONE ){
939 int code_sn=lines_filter_by_sn( pre, da->min_sn, x, y, &sn);
940 vlambdadata[i] = lambda;
941 vorderdata[i] = order;
942 vsdata[i] = slit;
943 if(vsndata!=NULL) {
944 vsndata[i] = sn;
945 }
946
947 vsindexdata[i] = slit_index;
948 vxthedata[i] = xthe;
949 vythedata[i] = ythe;
950 corr_x[i] = corrxv;
951 corr_y[i] = corryv;
952 gaussian_pos_x[i] = x;
953 gaussian_pos_y[i] = y;
954 gaussian_sigma_x[i] = sig_x;
955 gaussian_sigma_y[i] = sig_y;
956 gaussian_fwhm_x[i] = fwhm_x;
957 gaussian_fwhm_y[i] = fwhm_y;
958 gaussian_norm[i] = norm;
959 diffx[i] = x-xthe;
960 diffy[i] = y-ythe;
961
962
963 if ( code_sn ){
964 nlinematched++;
965 }
966 else{
967 lines_not_good_sn++;
968 flag[i]=2;
969 }
970 }
971 else{
972 flag[i]=3;
973 lines_not_gauss_fit++;
974 xsh_msg_dbg_medium("No Fit Gaussian for this line %d",i);
976 }
977 }
978 else{
979 flag[i]=4;
980 lines_not_in_image++;
981 xsh_msg_dbg_medium("Coordinates are not in the image");
982 }
983 }
984
985 /*In case of multi-pinhole arc lamp frame, remove all matches where there
986 are less than 7 out of 9 slit positions matched
987 the following assumes that the vectors are ordered so that, within
988 each order, all matches for a single wavelength are grouped together
989 */
990
991 if (nb_pinhole>1) {
992 for (i=1; i<sol_size; i++) {
993 xsh_msg_dbg_high("filter 7 pinhole : line %d : lambda %f order %f slit %f",
994 i, vlambdadata[i], vorderdata[i], vsdata[i]);
995
996 if (vlambdadata[i]!=vlambdadata[i-1] ||
997 vorderdata[i]!=vorderdata[i-1] ||
998 i==sol_size-1) {
999
1000 //special case for last line: TODO AMO: Why this???
1001 if (i==sol_size-1) {
1002 i++;
1003 }
1004 xsh_msg_dbg_high("filter n pinhole : from %d to %d find %d pinholes",
1005 starti, i, i-starti);
1006
1007 if (i-starti-nbflag>=min_slit_match) {
1008 for (k=starti; k<=i-1; k++) {
1009 diffxmean[k]=0.0;
1010 diffymean[k]=0.0;
1011 for (j=starti; j<=i-1; j++) {
1012 if (j!=k && flag[k]==0 && flag[j]==0) {
1013 diffxmean[k]+=diffx[j];
1014 diffymean[k]+=diffy[j];
1015 }
1016 }
1017 diffxmean[k]/=(float)(i-starti-1);
1018 diffymean[k]/=(float)(i-starti-1);
1019 }
1020 for (k=starti; k<=i-1; k++) {
1021 diffxsig[k]=0.0;
1022 diffysig[k]=0.0;
1023 for (j=starti; j<=i-1; j++) {
1024 if (j!=k && flag[k]==0 && flag[j]==0) {
1025 diffysig[k]+=(diffy[j]-diffymean[k])*(diffy[j]-diffymean[k]);
1026 diffxsig[k]+=(diffx[j]-diffxmean[k])*(diffx[j]-diffxmean[k]);
1027 }
1028 }
1029 diffxsig[k]=sqrt(diffxsig[k]/(float)(i-starti-1));
1030 diffysig[k]=sqrt(diffysig[k]/(float)(i-starti-1));
1031 }
1032 for (j=starti; j<=i-1; j++) {
1033 //Below removes individual ph with very different resids than the mean for this line
1034 if (fabs(diffx[j]-diffxmean[j])>line_devs*diffxsig[j] && fabs(diffy[j]-diffymean[j])>line_devs*diffysig[j] && flag[j]==0) {
1035 flag[j]=6;
1036 nbflag++;
1037 }
1038 }
1039 newi+=i-starti-nbflag;
1040 } // end min_slit_match
1041 // check min_slit again (rather than else) because nbflag may have increased
1042 if (i-starti-nbflag<min_slit_match) {
1043 lines_too_few_ph+=i-starti;
1044 for (j=starti; j<=i-1; j++) {
1045 if (flag[j]==0) {
1046 flag[j]=5;
1047 xsh_msg_dbg_medium("Too few pin-holes for this line %d",i);
1048 }
1049 }
1050 }
1051 starti=i;
1052 nbflag=0;
1053 }
1054
1055 if (i<sol_size && flag[i]!=0) {
1056 nbflag+=1;
1057 }
1058 }// endfor
1059 nlinematched=newi;
1060 } // if multi-pinhole (xsh_2dmap)
1061
1062
1063
1064
1065 xsh_msg("nlinematched / sol_size = %d / %d", nlinematched, sol_size);
1066 assure(nlinematched > 0, CPL_ERROR_ILLEGAL_INPUT,
1067 "No line matched, value of parameter "
1068 "detectarclines-search-win-hsize may be too large or too small "
1069 "detectarclines-fit-win-hsize=%d may be too large or too small "
1070 "detectarclines-running-median-hsize may be too large or too small "
1071 "detectarclines-min-sn may be too small "
1072 "or value of parameter detectarclines-min-sn may be too large",
1073 da->fit_window_hsize);
1074
1075
1076 xsh_msg_dbg_medium(" %d lines not found in image", lines_not_in_image);
1077 xsh_msg_dbg_medium(" %d lines not good S/N", lines_not_good_sn);
1078 xsh_msg_dbg_medium(" %d lines no fit gaussian", lines_not_gauss_fit);
1079 xsh_msg_dbg_medium(" %d lines no valid pixels", lines_not_valid_pixels);
1080 xsh_msg_dbg_medium(" %d lines detected in less than %d/9 ph positions", lines_too_few_ph,min_slit_match);
1081
1082 /* wrap variables into arrays in order to set to invalid flag==0 values
1083 and compute proper statistics */
1084 gaussian_sigma_x_arr = cpl_array_wrap_double(gaussian_sigma_x, nlinematched);
1085 gaussian_sigma_y_arr = cpl_array_wrap_double(gaussian_sigma_y, nlinematched);
1086 gaussian_fwhm_x_arr = cpl_array_wrap_double( gaussian_fwhm_x,nlinematched);
1087 gaussian_fwhm_y_arr = cpl_array_wrap_double( gaussian_fwhm_y,nlinematched);
1088
1089 /* Because we have initialised variabbles that hold results to 0, and not
1090 all elements are later filled by actual values from the Gauss fit,
1091 to get proper statistics we need to set to invalid all array elements
1092 corresponding to flag value 0
1093 */
1094
1095 for(i=0;i<nlinematched;i++){
1096 if(flag[i] > 0) {
1097 cpl_array_set_invalid(gaussian_sigma_x_arr,i);
1098 cpl_array_set_invalid(gaussian_sigma_y_arr,i);
1099 cpl_array_set_invalid(gaussian_fwhm_x_arr,i);
1100 cpl_array_set_invalid(gaussian_fwhm_y_arr,i);
1101 }
1102 }
1103
1104 xsh_msg("sigma gaussian median in x %lg",
1105 cpl_array_get_median( gaussian_sigma_x_arr));
1106 xsh_msg("sigma gaussian median in y %lg",
1107 cpl_array_get_median( gaussian_sigma_y_arr));
1108
1109 xsh_msg("FWHM gaussian median in x %lg",
1110 cpl_array_get_median( gaussian_fwhm_x_arr));
1111 xsh_msg("FWHM gaussian median in y %lg",
1112 cpl_array_get_median( gaussian_fwhm_y_arr));
1113
1114 xsh_unwrap_array( &gaussian_sigma_x_arr);
1115 xsh_unwrap_array( &gaussian_sigma_y_arr);
1116
1117 xsh_unwrap_array( &gaussian_fwhm_x_arr);
1118 xsh_unwrap_array( &gaussian_fwhm_y_arr);
1119
1120
1121 /* Produce the residual orders tab */
1122 /* make sure result is with flag==0 only */
1123 if ( resid_tab_orders_frame != NULL){
1124 check(resid_orders=xsh_resid_tab_create_not_flagged(sol_size,
1125 vlambdadata,
1126 vorderdata,
1127 vsdata,vsndata,
1128 vsindexdata,
1129 vxthedata, vythedata,
1130 corr_x, corr_y,
1131 gaussian_norm,
1132 gaussian_pos_x,
1133 gaussian_pos_y,
1134 gaussian_sigma_x,
1135 gaussian_sigma_y,
1136 gaussian_fwhm_x,
1137 gaussian_fwhm_y,flag,
1138 wave_table,
1139 solution_type));
1140
1141 if(resid_tab_name_sw) {
1142 sprintf(rtag,"%s%s%s",type,"RESID_TAB_ORDERS_",
1144 } else {
1145 sprintf(rtag,"%s%s%s",type,"RESID_ALL_TAB_ORDERS_",
1147 }
1148
1149 sprintf(rname,"%s%s",rtag,".fits");
1150
1151 check( *resid_tab_orders_frame = xsh_resid_tab_save( resid_orders,
1152 rname, instr,rtag));
1154
1155 }
1156
1157
1158 xsh_msg_dbg_high("solution_type=%d poly=%d model=%d",solution_type,
1160
1161 /* Produce the wave tab if need */
1162 /* TODO: CHECK
1163 Here poly mode results may have changed due to having stored
1164 results also with flag!=0 */
1165 if ( solution_type == XSH_DETECT_ARCLINES_TYPE_POLY &&
1166 (wave_tab_frame != NULL)){
1167 check( wave_table = xsh_wavesol_create( spectralformat_frame, da, instr));
1168 XSH_CALLOC( rejected, int, nlinematched);
1169 if ( solwave_type == XSH_SOLUTION_RELATIVE) {
1170 /* compute dY = f(lambda,n) */
1172 check( data_wavesol_fit_with_sigma( wave_table, diffy, vlambdadata,
1173 vorderdata, vsdata, nlinematched, dac->niter, dac->frac,
1174 dac->sigma, rejected));
1175 nb_rejected = 0;
1176 for(i=0; i< nlinematched; i++){
1177 if (rejected[i] == 1){
1178 nb_rejected++;
1179 }
1180 else{
1181 /* reindex data */
1182 vxthedata[i-nb_rejected] = vxthedata[i];
1183 vythedata[i-nb_rejected] = vythedata[i];
1184 vsindexdata[i-nb_rejected] = vsindexdata[i];
1185 corr_x[i-nb_rejected] = corr_x[i];
1186 corr_y[i-nb_rejected] = corr_y[i];
1187 gaussian_pos_x[i-nb_rejected] = gaussian_pos_x[i];
1188 gaussian_pos_y[i-nb_rejected] = gaussian_pos_y[i];
1189 gaussian_sigma_x[i-nb_rejected] = gaussian_sigma_x[i];
1190 gaussian_sigma_y[i-nb_rejected] = gaussian_sigma_y[i];
1191 gaussian_fwhm_x[i-nb_rejected] = gaussian_fwhm_x[i];
1192 gaussian_fwhm_y[i-nb_rejected] = gaussian_fwhm_y[i];
1193 gaussian_norm[i-nb_rejected] = gaussian_norm[i];
1194 diffx[i-nb_rejected] = diffx[i];
1195 }
1196 }
1197 nlinematched = nlinematched-nb_rejected;
1198 /* compute dX = f(lambda,n) */
1199 check( fit2dx = xsh_wavesol_get_polx( wave_table));
1200 check( xsh_wavesol_compute( wave_table,
1201 nlinematched, diffx,
1202 &(wave_table->min_x), &(wave_table->max_x),
1203 vlambdadata, vorderdata, vsdata, fit2dx));
1204
1205
1206 sprintf(wave_table_name,"%s%s.fits","WAVE_TAB_GUESS_",
1207 xsh_instrument_arm_tostring( instr )) ;
1208
1209 wave_table_tag = XSH_GET_TAG_FROM_ARM( XSH_WAVE_TAB_GUESS, instr);
1210
1211 }
1212 else{
1214 /* Y = f(lambda,n,s) with sigma clipping */
1215 check( data_wavesol_fit_with_sigma( wave_table, gaussian_pos_y,
1216 vlambdadata, vorderdata, vsdata, nlinematched, dac->niter, dac->frac,
1217 dac->sigma, rejected));
1218 nb_rejected = 0;
1219 for(i=0; i< nlinematched; i++){
1220 if (rejected[i] == 1){
1221 nb_rejected++;
1222 }
1223 else{
1224 /* reindex data */
1225 vxthedata[i-nb_rejected] = vxthedata[i];
1226 vythedata[i-nb_rejected] = vythedata[i];
1227 vsindexdata[i-nb_rejected] = vsindexdata[i];
1228 corr_x[i-nb_rejected] = corr_x[i];
1229 corr_y[i-nb_rejected] = corr_y[i];
1230 gaussian_pos_x[i-nb_rejected] = gaussian_pos_x[i];
1231 gaussian_sigma_x[i-nb_rejected] = gaussian_sigma_x[i];
1232 gaussian_sigma_y[i-nb_rejected] = gaussian_sigma_y[i];
1233 gaussian_fwhm_x[i-nb_rejected] = gaussian_fwhm_x[i];
1234 gaussian_fwhm_y[i-nb_rejected] = gaussian_fwhm_y[i];
1235 }
1236 }
1237 nlinematched = nlinematched-nb_rejected;
1238 /* X = f(lambda,n,s) */
1239 check(fit2dx = xsh_wavesol_get_polx(wave_table));
1240 check(xsh_wavesol_compute(wave_table, nlinematched, gaussian_pos_x,
1241 &wave_table->min_x, &wave_table->max_x, vlambdadata, vorderdata,
1242 vsdata, fit2dx));
1243
1244 sprintf(wave_table_name,"%s%s.fits","WAVE_TAB_2D_",
1245 xsh_instrument_arm_tostring( instr )) ;
1246
1247
1248 wave_table_tag = XSH_GET_TAG_FROM_ARM( XSH_WAVE_TAB_2D, instr);
1249 }
1250
1251
1252
1253 /* save the wavelength solution */
1254 check(header = xsh_wavesol_get_header(wave_table));
1255 check(xsh_pfits_set_qc_nlinecat(header,nlinecat));
1256 xsh_msg("poly nlinecat=%d",nlinecat);
1257
1258 check(xsh_pfits_set_qc_nlinefound(header,sol_size));
1259 // check(xsh_pfits_set_qc_nlinecat_clean(header,nlineclean));
1260 /* check(xsh_pfits_set_qc(header, (void*)model_date, XSH_QC_MODEL_FMTCHK_DATE,
1261 instr));*/
1262 check( wave_trace=xsh_wavesol_trace(wave_table,vlambdadata,vorderdata,
1263 vsdata,nlinematched));
1264
1265 check( *wave_tab_frame=xsh_wavesol_save(wave_table,wave_trace,
1266 wave_table_name,wave_table_tag));
1267 //xsh_add_temporary_file( wave_table_name) ;
1268 check( cpl_frame_set_tag( *wave_tab_frame, wave_table_tag));
1269 } /* end of poly mode case (if wavetab is provided in input) */
1270
1271
1272 /* Produce the clean arcline list using only not flagged points (flag==0) */
1273 xsh_arclist_clean_from_list_not_flagged( arclist, vlambdadata, flag,sol_size);
1274
1275 nlinecat_clean = arclist->size-arclist->nbrejected;
1276 check( header = xsh_arclist_get_header( arclist));
1277 check( xsh_pfits_set_qc_nlinecat( header,nlinecat));
1278
1279 check( xsh_pfits_set_qc_nlinefound( header, sol_size));
1280 check( xsh_pfits_set_qc_nlinecat_clean( header,nlinecat_clean));
1281 check( xsh_pfits_set_qc_nlinefound_clean( header,nlinematched));
1282 if(strcmp(rec_id,"xsh_predict") == 0) {
1284 } else {
1286 }
1287
1288 sprintf(fname,"%s%s",tag,".fits");
1289 check( *arc_lines_clean_tab_frame = xsh_arclist_save( arclist,fname,tag));
1291
1292
1293 /* Produce the residual tab */
1294 check( resid = xsh_resid_tab_create( sol_size, vlambdadata, vorderdata,
1295 vsdata,vsndata,vsindexdata,vxthedata,vythedata,
1296 corr_x,corr_y,gaussian_norm,
1297 gaussian_pos_x,gaussian_pos_y,
1298 gaussian_sigma_x, gaussian_sigma_y,
1299 gaussian_fwhm_x, gaussian_fwhm_y,flag,
1300 wave_table,solution_type));
1301
1302 if(resid_tab_name_sw) {
1303 sprintf(rtag,"%s%s%s",type,"RESID_TAB_LINES_",
1305 } else {
1306 sprintf(rtag,"%s%s%s",type,"RESID_TAB_ALL_LINES_",
1308 }
1309 check( tag = cpl_frame_get_tag( frame));
1310
1311 XSH_REGDEBUG("TAG %s", tag);
1312
1313 if ( strstr( tag, XSH_AFC_CAL) ){
1314 sprintf(rname,"AFC_CAL_%s%s",rtag,".fits") ;
1315 }
1316 else if ( strstr( tag, XSH_AFC_ATT) ){
1317 sprintf(rname,"AFC_ATT_%s%s",rtag,".fits") ;
1318 }
1319 else {
1320 sprintf(rname,"%s%s",rtag,".fits") ;
1321 }
1322 check( *resid_tab_frame = xsh_resid_tab_save( resid,rname,instr,rtag));
1323 if(clean_tmp) {
1325 }
1326
1327 cleanup:
1328 XSH_FREE( vlambdadata);
1329 XSH_FREE( vorderdata);
1330 XSH_FREE( vsdata);
1331 XSH_FREE( vsndata);
1332 XSH_FREE( vsindexdata);
1333 XSH_FREE( vxthedata);
1334 XSH_FREE( vythedata);
1335 XSH_FREE( corr_x);
1336 XSH_FREE( corr_y);
1337 XSH_FREE( gaussian_pos_x);
1338 XSH_FREE( gaussian_pos_y);
1339 XSH_FREE( gaussian_sigma_x);
1340 XSH_FREE( gaussian_sigma_y);
1341 XSH_FREE( gaussian_fwhm_x);
1342 XSH_FREE( gaussian_fwhm_y);
1343 XSH_FREE( gaussian_norm);
1344 XSH_FREE( sort);
1345 XSH_FREE( diffx);
1346 XSH_FREE( diffy);
1347 XSH_FREE( diffxmean);
1348 XSH_FREE( diffymean);
1349 XSH_FREE( diffxsig);
1350 XSH_FREE( diffysig);
1351 XSH_FREE( rejected);
1352 XSH_FREE( flag);
1353
1354 xsh_free_table( &wave_trace);
1355 xsh_spectralformat_list_free(&spectralformat_list);
1356 xsh_pre_free( &pre);
1357 xsh_the_map_free( &themap);
1358 xsh_wavesol_free( &wave_table);
1359 xsh_wavesol_free( &wave_tab_guess);
1360 xsh_order_list_free( &order_tab_recov);
1361 xsh_resid_tab_free( &resid_orders);
1362 xsh_resid_tab_free( &resid);
1363 xsh_arclist_free( &arclist);
1364 return;
1365}
1366/*---------------------------------------------------------------------------*/
1405void
1406xsh_detect_arclines( cpl_frame *frame,
1407 cpl_frame *theo_tab_frame,
1408 cpl_frame *arc_lines_tab_frame,
1409 cpl_frame* wave_tab_guess_frame,
1410 cpl_frame *order_tab_recov_frame,
1411 cpl_frame *config_model_frame,
1412 cpl_frame *spectralformat_frame,
1413 cpl_frame **resid_tab_orders_frame,
1414 cpl_frame **arc_lines_clean_tab_frame,
1415 cpl_frame **wave_tab_frame,
1416 cpl_frame **resid_tab_frame,
1417 xsh_sol_wavelength solwave_type,
1419 xsh_clipping_param *dac,
1420 xsh_instrument *instr,
1421 const char* rec_id,
1422 const int clean_tmp,const int resid_tab_name_sw)
1423{
1424 /* MUST BE DEALLOCATED in caller */
1425
1426 /* ALLOCATED locally */
1427 xsh_the_map *themap = NULL;
1428 xsh_resid_tab *resid = NULL;
1429 xsh_resid_tab *resid_orders = NULL;
1430 xsh_arclist *arclist = NULL;
1431 xsh_pre *pre = NULL;
1432 cpl_polynomial *fit2dx = NULL;
1433 cpl_propertylist* header = NULL;
1434 int * sort = NULL;
1435 double *vlambdadata = NULL, *vsdata = NULL, *vorderdata = NULL, *vsndata=NULL;
1436 int *vsindexdata = NULL;
1437 double *vxthedata = NULL, *vythedata = NULL;
1438 double *corr_x = NULL, *corr_y = NULL;
1439 double *gaussian_pos_x = NULL, *gaussian_pos_y = NULL,
1440 *gaussian_sigma_x = NULL, *gaussian_sigma_y = NULL,
1441 *gaussian_fwhm_x = NULL, *gaussian_fwhm_y = NULL, *gaussian_norm=NULL;
1442 cpl_vector *gaussian_sigma_x_vect =NULL;
1443 cpl_vector *gaussian_sigma_y_vect =NULL;
1444 cpl_vector *gaussian_fwhm_x_vect =NULL;
1445 cpl_vector *gaussian_fwhm_y_vect =NULL;
1446 //cpl_vector *gaussian_norm_vect =NULL;
1447
1448 double *diffx = NULL, *diffy = NULL;
1449 double *diffxmean = NULL, *diffymean = NULL, *diffxsig = NULL, *diffysig = NULL;
1450 xsh_order_list *order_tab_recov = NULL;
1451 /* Others */
1452 int nlinematched, nlinecat_clean=0, nlinecat;
1453 int i, sol_size = 0;
1454 xsh_wavesol* wave_table = NULL;
1455 char wave_table_name[256];
1456 const char* wave_table_tag = NULL;
1457 xsh_wavesol* wave_tab_guess = NULL;
1458 xsh_xs_3 config_model;
1459 xsh_xs_3* p_xs_3;
1460 xsh_spectralformat_list *spectralformat_list = NULL;
1461 int lines_not_in_image = 0;
1462 int lines_not_good_sn = 0;
1463 int lines_not_gauss_fit = 0;
1464 int lines_not_valid_pixels = 0;
1465 int lines_too_few_ph=0;
1466 //const char* model_date = NULL;
1467 int* rejected = NULL;
1468 int nb_rejected = 0;
1469 int solution_type = 0;
1470 const char* solution_type_name[2] = { "POLY", "MODEL"};
1471 int detection_mode = 0;
1472 const char* detection_mode_name[3] = { "NORMAL", "CORRECTED", "RECOVER"};
1473 int nb_pinhole;
1474 char fname[256];
1475 char rname[256];
1476 char rtag[256];
1477 const char* tag=NULL;
1478
1479 char dpr_type[256];
1480 char type[256];
1481 cpl_table* wave_trace=NULL;
1482 int* flag=NULL;
1483
1484 cpl_table* cfg_tab=NULL;
1485 const char* cfg_name=NULL;
1486 char new_name[256];
1487 char basename[256];
1488
1489
1490 cpl_propertylist* plist=NULL;
1491 //double* p_cfg_val=NULL;
1492 //double zeroK=-273.15;
1493
1494 int starti=0;
1495 int newi=0;
1496 int j,k;
1497 int min_slit_match=7;
1498 int found_temp=true;
1499
1500 /* check input parameters */
1501 XSH_ASSURE_NOT_NULL( frame);
1502 XSH_ASSURE_NOT_NULL( arc_lines_tab_frame);
1503 XSH_ASSURE_NOT_NULL( spectralformat_frame);
1504 XSH_ASSURE_NOT_NULL( instr);
1506 XSH_ASSURE_NOT_NULL( dac);
1507
1508 /* Ouput parameters */
1509 XSH_ASSURE_NOT_NULL( arc_lines_clean_tab_frame);
1510 XSH_ASSURE_NOT_NULL( resid_tab_frame);
1511
1512
1513 /* Debug parameters */
1514 xsh_msg_dbg_medium("---Detect arclines parameters");
1515 xsh_msg_dbg_medium("min_sn %f fit-window-half-size %d", da->min_sn, da->fit_window_hsize);
1516 xsh_msg_dbg_medium("Clipping sigma %f niter %d frac %f", dac->sigma, dac->niter, dac->frac);
1517
1518 check( pre = xsh_pre_load (frame, instr));
1519 check_msg( nb_pinhole = xsh_pfits_get_nb_pinhole( pre->data_header),
1520 "fail to find number of pinholes info in frame %s",
1521 cpl_frame_get_filename(frame) );
1522 check(strcpy(dpr_type,xsh_pfits_get_dpr_type(pre->data_header)));
1523
1524 if(strstr(dpr_type,"FMTCHK") != NULL) {
1525 strcpy(type,"FMTCHK_");
1526 } else if(strstr(dpr_type,"WAVE") != NULL) {
1527 strcpy(type,"WAVE_");
1528 } else if(strstr(dpr_type,"WAVE") != NULL) {
1529 strcpy(type,"ARC_");
1530 } else {
1531 strcpy(type,"");
1532 }
1533
1534 check( arclist = xsh_arclist_load( arc_lines_tab_frame));
1535 check( spectralformat_list = xsh_spectralformat_list_load(
1536 spectralformat_frame, instr));
1537 check(nlinecat = xsh_arclist_get_size( arclist));
1538 /* sort data by wavelength */
1539 check(xsh_arclist_lambda_sort( arclist));
1540
1541
1542 /* First define the solution type of detect_arclines */
1543 if ( theo_tab_frame != NULL) {
1544 XSH_ASSURE_NOT_ILLEGAL( config_model_frame == NULL);
1545 solution_type = XSH_DETECT_ARCLINES_TYPE_POLY;
1546 check( themap = xsh_the_map_load( theo_tab_frame));
1548 check( theo_tab_filter( themap, arclist, &sol_size, &vlambdadata,
1549 &vorderdata, &vsdata, &vsindexdata, &vxthedata, &vythedata, nb_pinhole));
1550
1551 }
1552 else if ( config_model_frame != NULL){
1553 solution_type = XSH_DETECT_ARCLINES_TYPE_MODEL;
1554 p_xs_3=&config_model;
1555
1556
1557 //Make sure to modify only current copy
1558 check(cfg_name=cpl_frame_get_filename(config_model_frame));
1559 check(plist=cpl_propertylist_load(cfg_name,0));
1560 check(cfg_tab=cpl_table_load(cfg_name,1,0));
1561 //check(basename=xsh_get_basename(cfg_name));
1563 sprintf(basename,"%s.fits",tag);
1564 sprintf(new_name,"local_%s",basename);
1565 check( xsh_pfits_set_pcatg(plist,tag));
1566 check(cpl_table_save(cfg_tab,plist,NULL,new_name,CPL_IO_DEFAULT));
1567 xsh_add_temporary_file(new_name);
1568 check(cpl_frame_set_filename(config_model_frame,new_name));
1569 xsh_free_table(&cfg_tab);
1570 xsh_free_propertylist(&plist);
1571
1572 check( xsh_model_config_load_best( config_model_frame, &config_model));
1573 XSH_REGDEBUG("load config model ok");
1574 /* update cfg table frame and corresponding model structure for temperature */
1575 check(xsh_model_temperature_update_frame(&config_model_frame,frame,
1576 instr,&found_temp));
1577 check(xsh_model_temperature_update_structure(&config_model,frame,instr));
1578 check( theo_tab_model( &config_model, arclist, spectralformat_list,
1579 &sol_size, &vlambdadata, &vorderdata, &vsdata,
1580 &vsndata, &vsindexdata, &vxthedata, &vythedata,
1581 instr, nb_pinhole));
1582 }
1583 else{
1585 "Undefined solution type (POLY or MODEL). See your input sof");
1586 }
1587
1588 xsh_spectralformat_list_free(&spectralformat_list);
1589 xsh_msg_dbg_high("Solution type %s", solution_type_name[solution_type]);
1590
1591 /* TODO detection mode checks and correction by order_tab_recov_frame is
1592 * not needed now that XSH is operational */
1593
1594 /* Define the mode of use */
1595 if ( wave_tab_guess_frame == NULL){
1596 if ( order_tab_recov_frame == NULL){
1597 detection_mode = XSH_DETECT_ARCLINES_MODE_NORMAL;
1598 }
1599 else{
1600 detection_mode = XSH_DETECT_ARCLINES_MODE_RECOVER;
1601 check( order_tab_recov = xsh_order_list_load( order_tab_recov_frame,
1602 instr));
1603 }
1604 }
1605 else {
1606 if ( order_tab_recov_frame == NULL){
1607 detection_mode = XSH_DETECT_ARCLINES_MODE_CORRECTED;
1608
1609 check( wave_tab_guess = xsh_wavesol_load( wave_tab_guess_frame,
1610 instr ));
1611 xsh_msg( "BinX,Y: %d, %d", wave_tab_guess->bin_x,
1612 wave_tab_guess->bin_y ) ;
1613 }
1614 else{
1615 detection_mode = -1;
1616 }
1617 }
1618
1620 && detection_mode <= XSH_DETECT_ARCLINES_MODE_RECOVER);
1621
1622 xsh_msg( "Detection mode : %s", detection_mode_name[detection_mode]);
1623
1624 xsh_msg( "Solution size %d", sol_size);
1625
1626 /* allocate some memory */
1627 XSH_MALLOC( corr_x, double, sol_size);
1628 XSH_MALLOC( corr_y, double, sol_size);
1629
1630 XSH_MALLOC( gaussian_pos_x, double, sol_size);
1631 XSH_MALLOC( gaussian_pos_y, double, sol_size);
1632 XSH_MALLOC( gaussian_sigma_x, double, sol_size);
1633 XSH_MALLOC( gaussian_sigma_y, double, sol_size);
1634
1635 XSH_MALLOC( gaussian_fwhm_x, double, sol_size);
1636 XSH_MALLOC( gaussian_fwhm_y, double, sol_size);
1637 XSH_MALLOC( gaussian_norm, double, sol_size);
1638
1639 XSH_MALLOC( diffx, double, sol_size);
1640 XSH_MALLOC( diffy, double, sol_size);
1641
1642 XSH_MALLOC( diffxmean, double, sol_size);
1643 XSH_MALLOC( diffymean, double, sol_size);
1644
1645 XSH_MALLOC( diffxsig, double, sol_size);
1646 XSH_MALLOC( diffysig, double, sol_size);
1647
1648 /* init */
1649 nlinematched = 0;
1650
1652 xsh_msg_dbg_medium( "USE method GAUSSIAN");
1653 }
1654 else{
1655 xsh_msg_dbg_medium( "USE method BARYCENTER");
1656 }
1657
1658 /* TODO: this ASCII input should be removed */
1660 FILE* regfile = NULL;
1661
1662 regfile = fopen( "FIT.reg", "w");
1663 fprintf( regfile, "# Region file format: DS9 version 4.0\n"\
1664 "global color=red font=\"helvetica 4 normal\""\
1665 "select=1 highlite=1 edit=1 move=1 delete=1 include=1 fixed=0 "\
1666 "source\nimage\n");
1667 fclose( regfile);
1668 regfile = fopen( "NOFIT.reg", "w");
1669 fprintf( regfile, "# Region file format: DS9 version 4.0\n"\
1670 "global color=red font=\"helvetica 4 normal\""\
1671 "select=1 highlite=1 edit=1 move=1 delete=1 include=1 fixed=0 "\
1672 "source\nimage\n");
1673 fclose( regfile);
1674 }
1675
1676 for(i=0; i< sol_size; i++){
1677 int xpos = 0, ypos = 0;
1678 int slit_index;
1679 double corrxv = 0.0, corryv = 0.0;
1680 double lambda, order, slit, xthe, ythe,sn=0;
1681 double xcor, ycor, x, y, sig_x, sig_y, norm, fwhm_x, fwhm_y;
1682 /* double max_off=1.0; */
1683
1684 lambda = vlambdadata[i];
1685 order = vorderdata[i];
1686 slit = vsdata[i];
1687
1688 if(vsndata != NULL) {
1689 sn = vsndata[i];
1690 }
1691
1692 slit_index = vsindexdata[i];
1693 xthe = vxthedata[i];
1694 ythe = vythedata[i];
1695 xsh_msg_dbg_high( "THE LAMBDA %f ORDER %f SLIT %f X %f Y %f", lambda,
1696 order, slit, xthe, ythe);
1697
1698 /* Theoretical using pixel FITS convention */
1699 xcor = xthe;
1700 ycor = ythe;
1701 xsh_msg_dbg_high("THE x%f y %f", xthe, ythe);
1702
1703 /* Find X position from RECOVER ORDER TAB */
1704 if ( order_tab_recov != NULL){
1705 int iorder = 0;
1706 cpl_polynomial *center_poly_recov = NULL;
1707
1709 order_tab_recov, order));
1710 center_poly_recov = order_tab_recov->list[iorder].cenpoly;
1711 check( xcor = cpl_polynomial_eval_1d( center_poly_recov, ycor, NULL));
1712 corrxv = xcor-xthe;
1713 }
1714 /* apply x and y translation from predict data to the model */
1715 if ( wave_tab_guess != NULL){
1716 double diffxv, diffyv;
1717
1718 check( diffxv = xsh_wavesol_eval_polx( wave_tab_guess, lambda,
1719 order, 0));
1720 check( diffyv = xsh_wavesol_eval_poly( wave_tab_guess, lambda,
1721 order, 0));
1722
1723 xcor = xthe+diffxv;
1724 ycor = ythe+diffyv;
1725 corrxv = diffxv;
1726 corryv = diffyv;
1727 }
1728 /* Corrected position using pixel FITS convention */
1729 xsh_msg_dbg_high("x %f y %f CORR_x %f CORR_y %f", xcor, ycor,
1730 corrxv, corryv);
1731
1732 if ( xcor >=0.5 && ycor >=0.5 &&
1733 xcor < (pre->nx+0.5) && ycor < (pre->ny+0.5)){
1734 int best_med_res = 0;
1735 check ( best_med_res = xsh_pre_window_best_median_flux_pos( pre,
1736 (int)(xcor+0.5)-1, (int)(ycor+0.5)-1,
1737 da->search_window_hsize, da->running_median_hsize, &xpos, &ypos));
1738
1739 if (best_med_res != 0){
1740 lines_not_valid_pixels++;
1741 continue;
1742 }
1743 /* Put positions in CPL coordinates */
1744 xpos = xpos+1;
1745 ypos = ypos+1;
1746 xsh_msg_dbg_high("ADJ x %d y %d", xpos, ypos);
1747
1749 cpl_image_fit_gaussian(pre->data, xpos, ypos, 1+2*da->fit_window_hsize,
1750 &norm, &x, &y, &sig_x, &sig_y, &fwhm_x, &fwhm_y);
1751 xsh_msg_dbg_high("GAUS x %f y %f %d %d %d", x, y,da->fit_window_hsize, da->search_window_hsize, da->running_median_hsize);
1752
1753
1754 /*
1755 cpl_array* parameters=cpl_array_new(7, CPL_TYPE_DOUBLE);
1756 cpl_array* err_params=cpl_array_new(7, CPL_TYPE_DOUBLE);
1757 cpl_array* fit_params=cpl_array_new(7, CPL_TYPE_INT);
1758 int kk=0;
1759 //All parameter should be fitted
1760 for (kk = 0; kk < 7; kk++)
1761 cpl_array_set(fit_params, kk, 1);
1762 double rms=0;
1763 double red_chisq=0;
1764 cpl_matrix* covariance=NULL;
1765 cpl_matrix* phys_cov=NULL;
1766 double major=0;
1767 double minor=0;
1768 double angle=0;
1769 int size=1+2*da->fit_window_hsize;
1770 cpl_fit_image_gaussian(pre->data, pre->errs, xpos, ypos,size,size,
1771 parameters,err_params,fit_params,
1772 &rms,&red_chisq,&covariance,
1773 &major,&minor,&angle,&phys_cov);
1774
1775 double rho=0;
1776 double back=cpl_array_get(parameters,0,NULL);
1777 norm=cpl_array_get(parameters,1,NULL);
1778 rho=cpl_array_get(parameters,2,NULL);
1779
1780 x=cpl_array_get(parameters,3,NULL);
1781 y=cpl_array_get(parameters,4,NULL);
1782 sig_x=cpl_array_get(parameters,5,NULL);
1783 sig_y=cpl_array_get(parameters,6,NULL);
1784 fwhm_x=sig_x*CPL_MATH_FWHM_SIG;
1785 fwhm_y=sig_y*CPL_MATH_FWHM_SIG;
1786
1787 xsh_free_array(&parameters);
1788 xsh_free_array(&err_params);
1789 xsh_free_array(&fit_params);
1790 */
1791
1792
1793
1794
1795 }
1796 else{
1797 xsh_image_find_barycenter( pre->data, xpos, ypos,
1798 1+2*da->fit_window_hsize,
1799 &norm, &x, &y, &sig_x, &sig_y, &fwhm_x, &fwhm_y);
1800 xsh_msg_dbg_high("BARY x %f y %f", x, y);
1801 }
1802
1803 /* (PBR) Adjust 1ph centroids to be consistent with 9ph
1804 - This is cleaner than adjusting the ph position in the cfg because:
1805 a. we do it just once here -> it gets used in the matching, then
1806 propogates to the model optimisation
1807 b. We don't need to know if the input config was itself optimised
1808 to 9ph or 1ph (i.e. do we need the shift or not)
1809 c. Potentially this can be used for poly-mode too, though I
1810 restrict it to phys mod below (feel free to change it)
1811 - The NIR value of 0.125pix for x I got by using a config optimised
1812 to 9ph data (i.e. 2dmap) in xsh_predict and checking the input
1813 offsets. For y the offset was less than the s.d.
1814 - Same procedure for UVB. Also same for VIS, but value was ~0
1815 - A correction from 9ph to full slit centre is still required*/
1816 if (nb_pinhole==1 && solution_type==XSH_DETECT_ARCLINES_TYPE_MODEL) {
1817 if (p_xs_3->arm==2) {
1818 x=x+0.125;
1819 }
1820 else if (p_xs_3->arm==0) {
1821 x=x-0.51;
1822 }
1823 }
1824
1826 const char* fit = NULL;
1827 const char* color = NULL;
1828 FILE* regfile = NULL;
1829
1830 if (cpl_error_get_code() == CPL_ERROR_NONE){
1831 fit = "FIT.reg";
1832 color = "green";
1833 }
1834 else{
1835 fit = "NOFIT.reg";
1836 color = "red";
1837 }
1838 regfile = fopen( fit, "a");
1839 fprintf( regfile, "point(%d,%d) #point=cross color=blue\n", (int)(xcor), (int)(ycor));
1840 fprintf( regfile, "box(%d,%d,%d,%d) #color=blue\n", (int)(xcor), (int)(ycor),
1842 fprintf( regfile, "point(%d,%d) #point=cross color=%s font=\"helvetica 10 normal\""\
1843 " text={%.3f}\n", xpos, ypos, color, lambda);
1844 fprintf( regfile, "box(%d,%d,%d,%d) #color=%s\n", (int)xpos, (int)ypos,
1845 1+2*da->fit_window_hsize, 1+2*da->fit_window_hsize, color);
1846 fclose( regfile);
1847 }
1848
1849 if( cpl_error_get_code() == CPL_ERROR_NONE ){
1850 if ( lines_filter_by_sn( pre, da->min_sn, x, y, &sn) ){
1851 vlambdadata[nlinematched] = lambda;
1852 vorderdata[nlinematched] = order;
1853 vsdata[nlinematched] = slit;
1854 if(vsndata!=NULL) {
1855 vsndata[nlinematched] = sn;
1856 }
1857
1858 vsindexdata[nlinematched] = slit_index;
1859 vxthedata[nlinematched] = xthe;
1860 vythedata[nlinematched] = ythe;
1861 corr_x[nlinematched] = corrxv;
1862 corr_y[nlinematched] = corryv;
1863 gaussian_pos_x[nlinematched] = x;
1864 gaussian_pos_y[nlinematched] = y;
1865 gaussian_sigma_x[nlinematched] = sig_x;
1866 gaussian_sigma_y[nlinematched] = sig_y;
1867 gaussian_fwhm_x[nlinematched] = fwhm_x;
1868 gaussian_fwhm_y[nlinematched] = fwhm_y;
1869 gaussian_norm[nlinematched] = norm;
1870 diffx[nlinematched] = x-xthe;
1871 diffy[nlinematched] = y-ythe;
1872/* if ((fabs(diffx[nlinematched])>max_off || */
1873/* fabs(diffy[nlinematched])>max_off) */
1874/* && nb_pinhole>1) { */
1875/* xsh_msg_dbg_medium("Offset too large for this line"); */
1876/* } */
1877/* else { */
1878 nlinematched++;
1879/* } */
1880 }
1881 else{
1882 lines_not_good_sn++;
1883 xsh_msg_dbg_medium("Not good s/n for this line");
1885
1886 int xpix, ypix, rej;
1887 float flux, noise, sn;
1888 FILE* regfile = NULL;
1889 char sn_name[256];
1890
1891 xpix =(int) rint(x);
1892 ypix = (int)rint(y);
1893
1894 check( flux = cpl_image_get( pre->data, xpix, ypix, &rej));
1895 check( noise = cpl_image_get( pre->errs, xpix, ypix, &rej));
1896 sn = flux / noise;
1897 sprintf( sn_name, "bad_sn_%.3f.reg", da->min_sn);
1898 regfile = fopen( sn_name, "a");
1899 fprintf( regfile, "point(%d,%d) #point=cross color=red font=\"helvetica 10 normal\""\
1900 " text={%.3f [%.3f}\n", xpix, ypix, lambda, sn);
1901 fclose( regfile);
1902 }
1903 }
1904 }
1905 else{
1906 lines_not_gauss_fit++;
1907 xsh_msg_dbg_medium("No Fit Gaussian for this line");
1909 }
1910 }
1911 else{
1912 lines_not_in_image++;
1913 xsh_msg_dbg_medium("Coordinates are not in the image");
1914 }
1915 }
1916
1917
1918 /*Remove all matches where there are less than 7 out of 9 slit positions matched*/
1919 /*The following assumes that the vectors are ordered so that, within each order,
1920 all matches for a single wavelength are grouped together*/
1921 if (nb_pinhole>1) {
1922 for (i=1; i<nlinematched; i++) {
1923 xsh_msg_dbg_high("filter 7 pinhole : line %d : lambda %f order %f slit %f",
1924 i, vlambdadata[i], vorderdata[i], vsdata[i]);
1925
1926 if (vlambdadata[i]!=vlambdadata[i-1] || vorderdata[i]!=vorderdata[i-1] || i==nlinematched-1) {
1927
1928 /*special case for last line*/
1929 if (i==nlinematched-1) {
1930 i++;
1931 }
1932 xsh_msg_dbg_high("filter 7 pinhole : from %d to %d find %d pinholes",
1933 starti, i, i-starti);
1934
1935 if (i-starti>=min_slit_match) {
1936 for (k=starti; k<=i-1; k++) {
1937 diffxmean[k]=0.0;
1938 diffymean[k]=0.0;
1939 for (j=starti; j<=i-1; j++) {
1940 if (j!=k) {
1941 diffxmean[k]+=diffx[j];
1942 diffymean[k]+=diffy[j];
1943 }
1944 }
1945 diffxmean[k]/=(float)(i-starti-1);
1946 diffymean[k]/=(float)(i-starti-1);
1947 }
1948 for (k=starti; k<=i-1; k++) {
1949 diffxsig[k]=0.0;
1950 diffysig[k]=0.0;
1951 for (j=starti; j<=i-1; j++) {
1952 if (j!=k) {
1953 diffysig[k]+=(diffy[j]-diffymean[k])*(diffy[j]-diffymean[k]);
1954 diffxsig[k]+=(diffx[j]-diffxmean[k])*(diffx[j]-diffxmean[k]);
1955 }
1956 }
1957 diffxsig[k]=sqrt(diffxsig[k]/(float)(i-starti-1));
1958 diffysig[k]=sqrt(diffysig[k]/(float)(i-starti-1));
1959 }
1960 for (j=starti; j<=i-1; j++) {
1961 /*Below would remove individual ph with very different resids than the mean for this line*/
1962 //if (fabs(diffx[j]-diffxmean)<1.5*diffxsig && fabs(diffy[j]-diffymean)<1.5*diffysig) {
1963 vlambdadata[newi] = vlambdadata[j];
1964 vorderdata[newi] = vorderdata[j];
1965 vsdata[newi] =vsdata[j];
1966 vsindexdata[newi] = vsindexdata[j];
1967 vxthedata[newi] =vxthedata[j];
1968 vythedata[newi] =vythedata[j];
1969 corr_x[newi] = corr_x[j];
1970 corr_y[newi] = corr_y[j];
1971 gaussian_pos_x[newi] = gaussian_pos_x[j];
1972 gaussian_pos_y[newi] = gaussian_pos_y[j];
1973 gaussian_sigma_x[newi] = gaussian_sigma_x[j];
1974 gaussian_sigma_y[newi] = gaussian_sigma_y[j];
1975 gaussian_fwhm_x[newi] = gaussian_fwhm_x[j];
1976 gaussian_fwhm_y[newi] = gaussian_fwhm_y[j];
1977 gaussian_norm[newi] = gaussian_norm[j];
1978 diffx[newi] = diffx[j];
1979 diffy[newi] = diffy[j];
1980 newi++;
1981 //}
1982/* else { */
1983/* printf("no %lf %lf %lf %lf %lf\n",vlambdadata[starti], fabs(diffx[j]),diffxmean,fabs(diffx[j]-diffxmean),diffxsig); */
1984/* printf("no %lf %lf %lf %lf %lf\n",vlambdadata[starti], fabs(diffx[j]),diffxmean,fabs(diffx[j]-diffxmean),diffxsig); */
1985/* } */
1986 }
1987 } /*end min_slit_match*/
1988 else {
1989 lines_too_few_ph+=i-starti;
1990 }
1991 starti=i;
1992 }
1993 }/*endfor*/
1994 nlinematched=newi;
1995 }
1996
1997 xsh_msg("nlinematched / sol_size = %d / %d", nlinematched, sol_size);
1998 assure(nlinematched > 0, CPL_ERROR_ILLEGAL_INPUT,
1999 "No line matched, "
2000 "detectarclines-search-win-hsize, "
2001 "detectarclines-fit-win-hsize=%d, "
2002 "detectarclines-running-median-hsize too large or too small "
2003 "or detectarclines-min-sn may be too large",
2004 da->fit_window_hsize);
2005
2006 xsh_msg_dbg_medium(" %d lines not found in image", lines_not_in_image);
2007 xsh_msg_dbg_medium(" %d lines not good S/N", lines_not_good_sn);
2008 xsh_msg_dbg_medium(" %d lines no fit gaussian", lines_not_gauss_fit);
2009 xsh_msg_dbg_medium(" %d lines no valid pixels", lines_not_valid_pixels);
2010 xsh_msg_dbg_medium(" %d lines detected in less than %d/9 ph positions", lines_too_few_ph,min_slit_match);
2011 check( gaussian_sigma_x_vect = cpl_vector_wrap( nlinematched, gaussian_sigma_x));
2012 check( gaussian_sigma_y_vect = cpl_vector_wrap( nlinematched, gaussian_sigma_y));
2013 check( gaussian_fwhm_x_vect = cpl_vector_wrap( nlinematched, gaussian_fwhm_x));
2014 check( gaussian_fwhm_y_vect = cpl_vector_wrap( nlinematched, gaussian_fwhm_y));
2015 xsh_msg("sigma gaussian median in x %lg",
2016 cpl_vector_get_median_const( gaussian_sigma_x_vect));
2017 xsh_msg("sigma gaussian median in y %lg",
2018 cpl_vector_get_median_const( gaussian_sigma_y_vect));
2019
2020 xsh_msg("FWHM gaussian median in x %lg",
2021 cpl_vector_get_median_const( gaussian_fwhm_x_vect));
2022 xsh_msg("FWHM gaussian median in y %lg",
2023 cpl_vector_get_median_const( gaussian_fwhm_y_vect));
2024
2025
2026 xsh_unwrap_vector( &gaussian_sigma_x_vect);
2027 xsh_unwrap_vector( &gaussian_sigma_y_vect);
2028
2029 xsh_unwrap_vector( &gaussian_fwhm_x_vect);
2030 xsh_unwrap_vector( &gaussian_fwhm_y_vect);
2031
2032
2033 /* Produce the residual orders tab */
2034
2035 if ( resid_tab_orders_frame != NULL){
2036 check(resid_orders=xsh_resid_tab_create(nlinematched,vlambdadata,vorderdata,
2037 vsdata, vsndata,vsindexdata,
2038 vxthedata, vythedata,
2039 corr_x, corr_y,
2040 gaussian_norm,
2041 gaussian_pos_x, gaussian_pos_y,
2042 gaussian_sigma_x, gaussian_sigma_y,
2043 gaussian_fwhm_x, gaussian_fwhm_y,flag,
2044 wave_table, solution_type));
2045
2046
2047 sprintf(rtag,"%s%s%s",type,"RESID_TAB_ORDERS_",
2049
2050 sprintf(rname,"%s%s",rtag,".fits");
2051
2052 check( *resid_tab_orders_frame = xsh_resid_tab_save( resid_orders,
2053 rname, instr,rtag));
2055
2056 }
2057 xsh_msg_dbg_high("solution_type=%d poly=%d model=%d",solution_type,
2059
2060
2061 /* Produce the wave tab if need */
2062 if ( solution_type == XSH_DETECT_ARCLINES_TYPE_POLY &&
2063 (wave_tab_frame != NULL)){
2064 check( wave_table = xsh_wavesol_create( spectralformat_frame, da, instr));
2065 XSH_CALLOC( rejected, int, nlinematched);
2066 if ( solwave_type == XSH_SOLUTION_RELATIVE) {
2067 /* compute dY = f(lambda,n) */
2069 check( data_wavesol_fit_with_sigma( wave_table, diffy, vlambdadata,
2070 vorderdata, vsdata, nlinematched, dac->niter, dac->frac,
2071 dac->sigma, rejected));
2072 nb_rejected = 0;
2073 for(i=0; i< nlinematched; i++){
2074 if (rejected[i] == 1){
2075 nb_rejected++;
2076 }
2077 else{
2078 /* reindex data */
2079 vxthedata[i-nb_rejected] = vxthedata[i];
2080 vythedata[i-nb_rejected] = vythedata[i];
2081 vsindexdata[i-nb_rejected] = vsindexdata[i];
2082 corr_x[i-nb_rejected] = corr_x[i];
2083 corr_y[i-nb_rejected] = corr_y[i];
2084 gaussian_pos_x[i-nb_rejected] = gaussian_pos_x[i];
2085 gaussian_pos_y[i-nb_rejected] = gaussian_pos_y[i];
2086 gaussian_sigma_x[i-nb_rejected] = gaussian_sigma_x[i];
2087 gaussian_sigma_y[i-nb_rejected] = gaussian_sigma_y[i];
2088 gaussian_fwhm_x[i-nb_rejected] = gaussian_fwhm_x[i];
2089 gaussian_fwhm_y[i-nb_rejected] = gaussian_fwhm_y[i];
2090 gaussian_norm[i-nb_rejected] = gaussian_norm[i];
2091 diffx[i-nb_rejected] = diffx[i];
2092 }
2093 }
2094 nlinematched = nlinematched-nb_rejected;
2095 /* compute dX = f(lambda,n) */
2096 check( fit2dx = xsh_wavesol_get_polx( wave_table));
2097 check( xsh_wavesol_compute( wave_table,
2098 nlinematched, diffx,
2099 &(wave_table->min_x), &(wave_table->max_x),
2100 vlambdadata, vorderdata, vsdata, fit2dx));
2101
2102
2103 sprintf(wave_table_name,"%s%s.fits","WAVE_TAB_GUESS_",
2104 xsh_instrument_arm_tostring( instr )) ;
2105
2106 wave_table_tag = XSH_GET_TAG_FROM_ARM( XSH_WAVE_TAB_GUESS, instr);
2107
2108 }
2109 else{
2111 /* Y = f(lambda,n,s) with sigma clipping */
2112 check( data_wavesol_fit_with_sigma( wave_table, gaussian_pos_y,
2113 vlambdadata, vorderdata, vsdata, nlinematched, dac->niter, dac->frac,
2114 dac->sigma, rejected));
2115 nb_rejected = 0;
2116 for(i=0; i< nlinematched; i++){
2117 if (rejected[i] == 1){
2118 nb_rejected++;
2119 }
2120 else{
2121 /* reindex data */
2122 vxthedata[i-nb_rejected] = vxthedata[i];
2123 vythedata[i-nb_rejected] = vythedata[i];
2124 vsindexdata[i-nb_rejected] = vsindexdata[i];
2125 corr_x[i-nb_rejected] = corr_x[i];
2126 corr_y[i-nb_rejected] = corr_y[i];
2127 gaussian_pos_x[i-nb_rejected] = gaussian_pos_x[i];
2128 gaussian_sigma_x[i-nb_rejected] = gaussian_sigma_x[i];
2129 gaussian_sigma_y[i-nb_rejected] = gaussian_sigma_y[i];
2130 gaussian_fwhm_x[i-nb_rejected] = gaussian_fwhm_x[i];
2131 gaussian_fwhm_y[i-nb_rejected] = gaussian_fwhm_y[i];
2132 }
2133 }
2134 nlinematched = nlinematched-nb_rejected;
2135 /* X = f(lambda,n,s) */
2136 check(fit2dx = xsh_wavesol_get_polx(wave_table));
2137 check(xsh_wavesol_compute(wave_table, nlinematched, gaussian_pos_x,
2138 &wave_table->min_x, &wave_table->max_x, vlambdadata, vorderdata,
2139 vsdata, fit2dx));
2140
2141 sprintf(wave_table_name,"%s%s.fits","WAVE_TAB_2D_",
2142 xsh_instrument_arm_tostring( instr )) ;
2143
2144
2145 wave_table_tag = XSH_GET_TAG_FROM_ARM( XSH_WAVE_TAB_2D, instr);
2146 }
2147
2148 /* save the wavelength solution */
2149 check(header = xsh_wavesol_get_header(wave_table));
2150 check(xsh_pfits_set_qc_nlinecat(header,nlinecat));
2151 xsh_msg("nlinecat=%d",nlinecat);
2152
2153 check(xsh_pfits_set_qc_nlinefound(header,sol_size));
2154// check(xsh_pfits_set_qc_nlinecat_clean(header,nlineclean));
2155/* check(xsh_pfits_set_qc(header, (void*)model_date, XSH_QC_MODEL_FMTCHK_DATE,
2156 instr));*/
2157 check( wave_trace=xsh_wavesol_trace(wave_table,vlambdadata,vorderdata,
2158 vsdata,nlinematched));
2159
2160 check( *wave_tab_frame=xsh_wavesol_save(wave_table,wave_trace,
2161 wave_table_name,wave_table_tag));
2162 //xsh_add_temporary_file( wave_table_name) ;
2163 check( cpl_frame_set_tag( *wave_tab_frame, wave_table_tag));
2164 }
2165
2166 /* Produce the clean arcline list */
2167 check( xsh_arclist_clean_from_list( arclist, vlambdadata, nlinematched));
2168 nlinecat_clean = arclist->size-arclist->nbrejected;
2169 check( header = xsh_arclist_get_header( arclist));
2170 check( xsh_pfits_set_qc_nlinecat( header,nlinecat));
2171 xsh_msg("nlinecat=%d",nlinecat);
2172
2173 check( xsh_pfits_set_qc_nlinefound( header, sol_size));
2174 check( xsh_pfits_set_qc_nlinecat_clean( header,nlinecat_clean));
2175 check( xsh_pfits_set_qc_nlinefound_clean( header,nlinematched));
2176 if(strcmp(rec_id,"xsh_predict") == 0) {
2178 } else {
2180 }
2181 sprintf(fname,"%s%s",tag,".fits");
2182 check( *arc_lines_clean_tab_frame = xsh_arclist_save( arclist,fname,tag));
2183 if(clean_tmp) {
2185 }
2186 /* Produce the residual tab */
2187 check( resid = xsh_resid_tab_create( nlinematched, vlambdadata, vorderdata,
2188 vsdata,vsndata,vsindexdata,vxthedata,vythedata,
2189 corr_x,corr_y,gaussian_norm,
2190 gaussian_pos_x,gaussian_pos_y,
2191 gaussian_sigma_x, gaussian_sigma_y,
2192 gaussian_fwhm_x, gaussian_fwhm_y,flag,
2193 wave_table,solution_type));
2194
2195 if(resid_tab_name_sw) {
2196 sprintf(rtag,"%s%s%s",type,"RESID_TAB_DRL_LINES_",
2198 } else {
2199 sprintf(rtag,"%s%s%s",type,"RESID_TAB_LINES_",
2201 }
2202 check( tag = cpl_frame_get_tag( frame));
2203
2204 XSH_REGDEBUG("TAG %s", tag);
2205
2206 if ( strstr( tag, XSH_AFC_CAL) ){
2207 sprintf(rname,"AFC_CAL_%s%s",rtag,".fits") ;
2208 }
2209 else if ( strstr( tag, XSH_AFC_ATT) ){
2210 sprintf(rname,"AFC_ATT_%s%s",rtag,".fits") ;
2211 }
2212 else {
2213 sprintf(rname,"%s%s",rtag,".fits") ;
2214 }
2215 check( *resid_tab_frame = xsh_resid_tab_save( resid,rname,instr,rtag));
2216 if(clean_tmp || resid_tab_name_sw) {
2218 }
2219
2220 cleanup:
2221 XSH_FREE( vlambdadata);
2222 XSH_FREE( vorderdata);
2223 XSH_FREE( vsdata);
2224 XSH_FREE( vsndata);
2225 XSH_FREE( vsindexdata);
2226 XSH_FREE( vxthedata);
2227 XSH_FREE( vythedata);
2228 XSH_FREE( corr_x);
2229 XSH_FREE( corr_y);
2230 XSH_FREE( gaussian_pos_x);
2231 XSH_FREE( gaussian_pos_y);
2232 XSH_FREE( gaussian_sigma_x);
2233 XSH_FREE( gaussian_sigma_y);
2234 XSH_FREE( gaussian_fwhm_x);
2235 XSH_FREE( gaussian_fwhm_y);
2236 XSH_FREE( gaussian_norm);
2237 XSH_FREE( sort);
2238 XSH_FREE( diffx);
2239 XSH_FREE( diffy);
2240 XSH_FREE( diffxmean);
2241 XSH_FREE( diffymean);
2242 XSH_FREE( diffxsig);
2243 XSH_FREE( diffysig);
2244 XSH_FREE( rejected);
2245
2246 xsh_free_table( &wave_trace);
2247 xsh_spectralformat_list_free(&spectralformat_list);
2248 xsh_pre_free( &pre);
2249 xsh_the_map_free( &themap);
2250 xsh_wavesol_free( &wave_table);
2251 xsh_wavesol_free( &wave_tab_guess);
2252 xsh_order_list_free( &order_tab_recov);
2253 xsh_resid_tab_free( &resid_orders);
2254 xsh_resid_tab_free( &resid);
2255 xsh_arclist_free( &arclist);
2256 return;
2257}
2258
2259/*---------------------------------------------------------------------------*/
static double sigma
void xsh_arclist_lambda_sort(xsh_arclist *list)
sort arcline list by increasing lambda
void xsh_arclist_clean_from_list_not_flagged(xsh_arclist *list, double *lambda, int *flag, int size)
float xsh_arclist_get_wavelength(xsh_arclist *list, int idx)
get wavelength of a line in the arcline list
void xsh_arclist_clean_from_list(xsh_arclist *list, double *lambda, int size)
Clean an arclist according to a list of valid lambda.
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
cpl_propertylist * xsh_arclist_get_header(xsh_arclist *list)
get header of the table
cpl_frame * xsh_arclist_save(xsh_arclist *list, const char *filename, const char *tag)
save a arclist 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
void xsh_order_list_free(xsh_order_list **list)
free memory associated to an order_list
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
int xsh_pre_window_best_median_flux_pos(xsh_pre *pre, int xcen, int ycen, int search_window_hsize, int running_median_hsize, int *xadj, int *yadj)
Search pixel position of best running median flux in the search window.
void xsh_pre_free(xsh_pre **pre)
Free a xsh_pre structure.
Definition: xsh_data_pre.c:823
void xsh_resid_tab_free(xsh_resid_tab **resid)
Free memory associated to a resid_tab.
xsh_resid_tab * xsh_resid_tab_create(int size, double *lambda, double *order, double *slit, double *sn, int *slit_index, double *thpre_x, double *thpre_y, double *corr_x, double *corr_y, double *gaussian_norm, double *gaussian_fit_x, double *gaussian_fit_y, double *gaussian_sigma_x, double *gaussian_sigma_y, double *gaussian_fwhm_x, double *gaussian_fwhm_y, int *flag, xsh_wavesol *wavesol, int wavesol_type)
Create a residual tab structure.
xsh_resid_tab * xsh_resid_tab_create_not_flagged(int size, double *lambda, double *order, double *slit, double *sn, int *slit_index, double *thpre_x, double *thpre_y, double *corr_x, double *corr_y, double *gaussian_norm, double *gaussian_fit_x, double *gaussian_fit_y, double *gaussian_sigma_x, double *gaussian_sigma_y, double *gaussian_fwhm_x, double *gaussian_fwhm_y, int *flag, xsh_wavesol *wavesol, int wavesol_type)
cpl_frame * xsh_resid_tab_save(xsh_resid_tab *resid, const char *filename, xsh_instrument *instr, const char *tag)
Save a residual tab to a frame.
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.
void xsh_the_map_free(xsh_the_map **list)
free memory associated to a the_map
double xsh_the_map_get_detx(xsh_the_map *list, int idx)
get detx of the map list
float xsh_the_map_get_slit_position(xsh_the_map *list, int idx)
get slit position of the map list
int xsh_the_map_get_slit_index(xsh_the_map *list, int idx)
get slit position of the map list
double xsh_the_map_get_dety(xsh_the_map *list, int idx)
get dety of the map list
int xsh_the_map_get_order(xsh_the_map *list, int idx)
get order of the map list
int xsh_the_map_get_size(xsh_the_map *list)
get size of the map list
void xsh_the_map_lambda_order_slit_sort(xsh_the_map *list)
float xsh_the_map_get_wavelength(xsh_the_map *list, int idx)
get wavelength of the map list
xsh_the_map * xsh_the_map_load(cpl_frame *frame)
load a theoretical map frame in the_map structure. Suppress spurious entries in the THE MAP (marked w...
cpl_table * xsh_wavesol_trace(xsh_wavesol *wsol, double *lambda, double *order, double *slit, int size)
void xsh_wavesol_compute(xsh_wavesol *sol, int size, double *pos, double *posmin, double *posmax, double *lambda, double *order, double *slit, cpl_polynomial *result)
compute a wavelength solution
cpl_polynomial * xsh_wavesol_get_polx(xsh_wavesol *sol)
get the solution in X
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
xsh_wavesol * xsh_wavesol_create(cpl_frame *spectral_format_frame, xsh_detect_arclines_param *params, xsh_instrument *instrument)
Create a new wavelength solution structure.
cpl_polynomial * xsh_wavesol_get_poly(xsh_wavesol *sol)
get the 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
void xsh_wavesol_set_type(xsh_wavesol *wsol, enum wavesol_type type)
set the type of the wave table
cpl_propertylist * xsh_wavesol_get_header(xsh_wavesol *sol)
get header of the table
cpl_frame * xsh_wavesol_save(xsh_wavesol *w, cpl_table *trace, const char *filename, const char *tag)
save a wavelength solution
void xsh_detect_arclines_dan(cpl_frame *frame, cpl_frame *theo_tab_frame, cpl_frame *arc_lines_tab_frame, cpl_frame *wave_tab_guess_frame, cpl_frame *order_tab_recov_frame, cpl_frame *config_model_frame, cpl_frame *spectralformat_frame, cpl_frame **resid_tab_orders_frame, cpl_frame **arc_lines_clean_tab_frame, cpl_frame **wave_tab_frame, cpl_frame **resid_tab_frame, xsh_sol_wavelength solwave_type, xsh_detect_arclines_param *da, xsh_clipping_param *dac, xsh_instrument *instr, const char *rec_id, const int clean_tmp, const int resid_tab_name_sw)
detect the position on the detector of emission lines listed in a catalogue, from expected position v...
static int lines_filter_by_sn(xsh_pre *pre, double sn_ref, double x, double y, double *sn)
static void theo_tab_model(xsh_xs_3 *config_model, xsh_arclist *arclist, xsh_spectralformat_list *spectralformat_list, int *size, double **lambda, double **n, double **s, double **sn, int **s_index, double **xthe, double **ythe, xsh_instrument *instr, int nb_pinhole)
static void theo_tab_filter(xsh_the_map *the_tab, xsh_arclist *arclist, int *size, double **lambda, double **n, double **s, int **s_index, double **xthe, double **ythe, int nb_pinhole)
void xsh_detect_arclines(cpl_frame *frame, cpl_frame *theo_tab_frame, cpl_frame *arc_lines_tab_frame, cpl_frame *wave_tab_guess_frame, cpl_frame *order_tab_recov_frame, cpl_frame *config_model_frame, cpl_frame *spectralformat_frame, cpl_frame **resid_tab_orders_frame, cpl_frame **arc_lines_clean_tab_frame, cpl_frame **wave_tab_frame, cpl_frame **resid_tab_frame, xsh_sol_wavelength solwave_type, xsh_detect_arclines_param *da, xsh_clipping_param *dac, xsh_instrument *instr, const char *rec_id, const int clean_tmp, const int resid_tab_name_sw)
detect the position on the detector of emission lines listed in a catalogue, from expected position v...
static void data_wavesol_fit_with_sigma(xsh_wavesol *wavesol, double *A, double *lambda, double *n, double *s, int size, int max_iter, double min_frac, double sigma, int *rejected)
#define XSH_REGDEBUG(...)
Definition: xsh_error.h:132
#define XSH_ASSURE_NOT_ILLEGAL(cond)
Definition: xsh_error.h:107
#define assure(CONDITION, ERROR_CODE,...)
Definition: xsh_error.h:54
#define check(COMMAND)
Definition: xsh_error.h:71
#define check_msg(COMMAND,...)
Definition: xsh_error.h:62
#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
cpl_error_code xsh_image_find_barycenter(const cpl_image *im, int xpos, int ypos, int size, double *norm, double *xcen, double *ycen, double *sig_x, double *sig_y, double *fwhm_x, double *fwhm_y)
Apply a gaussian fit on an image sub window.
Definition: xsh_fit.c:696
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.
int size
int * y
int * x
static SimAnneal s
Definition: xsh_model_sa.c:99
#define xsh_msg_dbg_medium(...)
Definition: xsh_msg.h:44
#define xsh_msg(...)
Print a message on info level.
Definition: xsh_msg.h:121
#define xsh_msg_dbg_high(...)
Definition: xsh_msg.h:40
const char * xsh_pfits_get_dpr_type(const cpl_propertylist *plist)
find out the DPR TECH
Definition: xsh_pfits.c:1667
void xsh_pfits_set_pcatg(cpl_propertylist *plist, const char *value)
Write the PCATG value.
Definition: xsh_pfits.c:1008
int xsh_pfits_get_nb_pinhole(const cpl_propertylist *plist)
Get the number of pinhole.
Definition: xsh_pfits.c:3409
void xsh_pfits_set_qc_nlinecat(cpl_propertylist *plist, double value)
Definition: xsh_pfits_qc.c:587
void xsh_pfits_set_qc_nlinecat_clean(cpl_propertylist *plist, double value)
Definition: xsh_pfits_qc.c:597
void xsh_pfits_set_qc_nlinefound_clean(cpl_propertylist *plist, double value)
Definition: xsh_pfits_qc.c:635
void xsh_pfits_set_qc_nlinefound(cpl_propertylist *plist, double value)
Definition: xsh_pfits_qc.c:626
void xsh_unwrap_vector(cpl_vector **v)
Unwrap a vector and set the pointer to NULL.
Definition: xsh_utils.c:2345
void xsh_free_vector(cpl_vector **v)
Deallocate a vector and set the pointer to NULL.
Definition: xsh_utils.c:2284
int xsh_debug_level_get(void)
get debug level
Definition: xsh_utils.c:3142
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_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
DOUBLE slit[10]
xsh_find_center_method find_center_method
xsh_order * list
cpl_polynomial * cenpoly
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
@ XSH_WAVESOL_GUESS
@ XSH_WAVESOL_2D
int n
Definition: xsh_detmon_lg.c:92
int order
Definition: xsh_detmon_lg.c:80
#define XSH_MOD_CFG_TAB
Definition: xsh_dfs.h:1251
#define XSH_ARC_LINE_LIST_2DMAP
Definition: xsh_dfs.h:965
#define XSH_WAVE_TAB_GUESS
Definition: xsh_dfs.h:555
#define XSH_GET_TAG_FROM_ARM(TAG, instr)
Definition: xsh_dfs.h:1548
#define XSH_ARC_LINE_LIST_PREDICT
Definition: xsh_dfs.h:960
#define XSH_WAVE_TAB_2D
Definition: xsh_dfs.h:556
#define XSH_AFC_CAL
Definition: xsh_dfs.h:539
#define XSH_AFC_ATT
Definition: xsh_dfs.h:543
#define WAVELENGTH_PRECISION
Definition: xsh_drl.h:190
xsh_sol_wavelength
Definition: xsh_drl.h:231
@ XSH_SOLUTION_RELATIVE
Definition: xsh_drl.h:233
#define XSH_DETECT_ARCLINES_MODE_RECOVER
Definition: xsh_drl.h:229
#define XSH_DETECT_ARCLINES_MODE_NORMAL
Definition: xsh_drl.h:227
#define XSH_DETECT_ARCLINES_TYPE_POLY
Definition: xsh_drl.h:225
#define XSH_DETECT_ARCLINES_MODE_CORRECTED
Definition: xsh_drl.h:228
#define XSH_DETECT_ARCLINES_TYPE_MODEL
Definition: xsh_drl.h:226
cpl_error_code xsh_model_temperature_update_frame(cpl_frame **model_config_frame, cpl_frame *ref_frame, xsh_instrument *instrument, int *found_temp)
cpl_error_code xsh_model_temperature_update_structure(struct xs_3 *p_xs_3_config, cpl_frame *ref_frame, xsh_instrument *instrument)
@ XSH_GAUSSIAN_METHOD
#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