00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #ifdef HAVE_CONFIG_H
00027 #include <config.h>
00028 #endif
00029
00030
00038
00041
00042
00043
00044
00045 #include <math.h>
00046 #include <xsh_drl.h>
00047 #include <xsh_pfits.h>
00048 #include <xsh_utils.h>
00049 #include <xsh_model_utils.h>
00050 #include <xsh_data_order.h>
00051 #include <xsh_error.h>
00052 #include <xsh_utils.h>
00053 #include <xsh_msg.h>
00054 #include <xsh_data_pre.h>
00055 #include <xsh_data_order.h>
00056 #include <xsh_data_the_map.h>
00057 #include <xsh_data_arclist.h>
00058 #include <xsh_data_wavesol.h>
00059 #include <xsh_data_resid_tab.h>
00060 #include <xsh_data_wavemap.h>
00061 #include <xsh_data_spectralformat.h>
00062 #include <xsh_model_io.h>
00063 #include <xsh_model_kernel.h>
00064 #include <xsh_fit.h>
00065 #include <gsl/gsl_multifit.h>
00066 #include <cpl.h>
00067
00068
00069
00070
00071 static void theo_tab_filter( xsh_the_map *the_tab, xsh_arclist *arclist,
00072 int* size, double **lambda, double **n, double **s, int **s_index,
00073 double **xthe, double **ythe, int nb_pinhole);
00074
00075 static void theo_tab_model( xsh_xs_3* config_model, xsh_arclist *arclist,
00076 xsh_spectralformat_list *spectralformat_list, int* size, double **lambda,
00077 double **n, double **s, double **sn,
00078 int **s_index, double **xthe, double **ythe,
00079 xsh_instrument* instr, int nb_pinhole);
00080
00081 static void data_wavesol_fit_with_sigma( xsh_wavesol *wavesol,
00082 double *A, double *lambda, double *n, double *s, int size,
00083 int max_iter, double min_frac, double sigma, int* rejected);
00084
00085 static int lines_filter_by_sn( xsh_pre* pre, double sn_ref, double x,
00086 double y, double* sn);
00087
00088
00089
00090
00091 static int lines_filter_by_sn( xsh_pre* pre, double sn_ref, double x,
00092 double y, double* sn)
00093 {
00094 int res = 0;
00095 int xpix, ypix, rej;
00096 double flux, noise;
00097
00098 XSH_ASSURE_NOT_NULL( pre);
00099
00100
00101 xpix =(int) rint(x);
00102 ypix = (int)rint(y);
00103
00104 check( flux = cpl_image_get( pre->data, xpix, ypix, &rej));
00105 check( noise = cpl_image_get( pre->errs, xpix, ypix, &rej));
00106 *sn = flux / noise;
00107 res = (*sn > sn_ref);
00108 cleanup:
00109 return res;
00110 }
00111
00112 static void theo_tab_model( xsh_xs_3* config_model, xsh_arclist *arclist,
00113 xsh_spectralformat_list *spectralformat_list, int* size, double **lambda,
00114 double **n, double **s, double** sn,
00115 int **s_index, double **xthe, double **ythe,
00116 xsh_instrument* instr, int nb_pinhole)
00117 {
00118 int i,j;
00119 cpl_vector** spectral_tab = NULL;
00120 int global_size = 0;
00121 int arc_size;
00122 int index_loc = 0;
00123 int nb_slit = 1;
00124
00125 XSH_ASSURE_NOT_NULL( config_model);
00126 XSH_ASSURE_NOT_NULL( arclist);
00127 XSH_ASSURE_NOT_NULL( spectralformat_list);
00128
00129 check( arc_size = xsh_arclist_get_size( arclist));
00130 XSH_CALLOC( spectral_tab, cpl_vector*, arc_size);
00131
00132 nb_slit = nb_pinhole;
00133
00134 XSH_REGDEBUG("nb pinhole %d ", nb_slit);
00135
00136 for( i=0; i< arc_size; i++){
00137 cpl_vector* res = NULL;
00138 float lambdaARC;
00139 int res_size;
00140
00141 check( lambdaARC = xsh_arclist_get_wavelength(arclist, i));
00142 check( res = xsh_spectralformat_list_get_orders(spectralformat_list,
00143 lambdaARC));
00144 if (res != NULL){
00145 check( res_size = cpl_vector_get_size( res));
00146 res_size *= nb_slit;
00147 }
00148 else{
00149 res_size = 0;
00150 check(xsh_arclist_reject(arclist,i));
00151 }
00152 global_size += res_size;
00153 spectral_tab[i] = res;
00154 }
00155
00156 XSH_MALLOC( *lambda, double, global_size);
00157 XSH_MALLOC( *n, double, global_size);
00158 XSH_MALLOC( *s, double, global_size);
00159 XSH_MALLOC( *sn, double, global_size);
00160 XSH_MALLOC( *s_index, int, global_size);
00161 XSH_MALLOC( *xthe, double, global_size);
00162 XSH_MALLOC( *ythe, double, global_size);
00163
00164
00165 for( i=0; i< arc_size; i++){
00166 double model_x, model_y;
00167 cpl_vector *spectral_res = NULL;
00168 int spectral_res_size;
00169 float lambdaARC;
00170
00171 check( lambdaARC = xsh_arclist_get_wavelength(arclist, i));
00172 check( spectral_res = spectral_tab[i]);
00173
00174 if (spectral_res != NULL){
00175 check( spectral_res_size = cpl_vector_get_size( spectral_res));
00176 for( j=0; j< spectral_res_size; j++){
00177 int absorder = 0;
00178 int islit;
00179 double slit;
00180
00181 if (nb_slit > 1){
00182
00183 for(islit=0; islit < nb_slit; islit++){
00184
00185 slit = config_model->slit[islit];
00186 check( absorder = (int) cpl_vector_get( spectral_res, j));
00187
00188
00189
00190 check( xsh_model_get_xy( config_model, instr, lambdaARC, absorder, slit,
00191 &model_x, &model_y));
00192 (*lambda)[index_loc] = lambdaARC;
00193 (*n)[index_loc] = absorder;
00194 (*s_index)[index_loc] = islit;
00195 (*s)[index_loc] = slit;
00196 (*sn)[index_loc] = 0;
00197 (*xthe)[index_loc] = model_x;
00198 (*ythe)[index_loc] = model_y;
00199 index_loc++;
00200 }
00201 }
00202 else{
00203 islit = 4;
00204 slit = config_model->slit[islit];
00205 check( absorder = (int) cpl_vector_get( spectral_res, j));
00206 check( xsh_model_get_xy( config_model, instr, lambdaARC, absorder, slit,
00207 &model_x, &model_y));
00208 (*lambda)[index_loc] = lambdaARC;
00209 (*n)[index_loc] = absorder;
00210 (*s_index)[index_loc] = islit;
00211 (*s)[index_loc] = slit;
00212 (*sn)[index_loc] = 0;
00213 (*xthe)[index_loc] = model_x;
00214 (*ythe)[index_loc] = model_y;
00215 index_loc++;
00216 }
00217
00218
00219 }
00220 }
00221 }
00222 *size = global_size;
00223
00224 cleanup:
00225 if ( spectral_tab != NULL){
00226 for(i=0; i< arc_size; i++){
00227 xsh_free_vector( &spectral_tab[i]);
00228 }
00229 XSH_FREE( spectral_tab);
00230 }
00231 if ( cpl_error_get_code() != CPL_ERROR_NONE){
00232 XSH_FREE( *lambda);
00233 XSH_FREE( *n);
00234 XSH_FREE( *s);
00235 XSH_FREE( *s_index);
00236 XSH_FREE( *xthe);
00237 XSH_FREE( *ythe);
00238 }
00239 return;
00240 }
00241
00242
00243
00244
00245
00246
00247
00248
00249
00250
00251
00252
00253
00254
00255
00256
00257
00258
00259
00260
00261
00262 static void theo_tab_filter( xsh_the_map *the_tab, xsh_arclist *arclist,
00263 int* size, double **lambda, double **n, double **s, int **s_index,
00264 double **xthe, double **ythe, int nb_pinhole)
00265 {
00266 int the_size = 0;
00267 int arc_size = 0;
00268 int i=0, j=0;
00269 int sol_size = 0;
00270 int nb_slit = 1;
00271
00272 XSH_ASSURE_NOT_NULL( the_tab);
00273 XSH_ASSURE_NOT_NULL( arclist);
00274
00275 check( arc_size = xsh_arclist_get_size( arclist));
00276 check( the_size = xsh_the_map_get_size( the_tab));
00277
00278 XSH_MALLOC( *lambda, double, the_size);
00279 XSH_MALLOC( *n, double, the_size);
00280 XSH_MALLOC( *s, double, the_size);
00281 XSH_MALLOC( *s_index, int, the_size);
00282 XSH_MALLOC( *xthe, double, the_size);
00283 XSH_MALLOC( *ythe, double, the_size);
00284
00285 nb_slit = nb_pinhole;
00286 XSH_REGDEBUG("nb pinhole %d ", nb_slit);
00287
00288 for(i=0; i< arc_size; i++){
00289 float lambdaARC, lambdaTHE;
00290 int nb_match = 0, max_match = 0;
00291
00292 check(lambdaARC = xsh_arclist_get_wavelength(arclist, i));
00293 check(lambdaTHE = xsh_the_map_get_wavelength(the_tab, j));
00294
00295 xsh_msg_dbg_medium("LINETABLE : Line %d / %d Lambda %f",i+1, arc_size,
00296 lambdaARC);
00297
00298 while ( (j < the_size-1) &&
00299 ( (lambdaARC-lambdaTHE) > WAVELENGTH_PRECISION) ){
00300 j++;
00301 check(lambdaTHE = xsh_the_map_get_wavelength(the_tab, j));
00302 }
00303
00304 xsh_msg_dbg_medium("THETABLE : Line %d / %d Lambda %f",j+1, the_size, lambdaTHE);
00305
00306 max_match = the_size-j;
00307
00308 while ( nb_match < max_match && fabs(lambdaARC-lambdaTHE) <= WAVELENGTH_PRECISION ){
00309 double order, slit, xtheval, ytheval;
00310 int islit;
00311
00312 check( slit = (double) xsh_the_map_get_slit_position( the_tab, j));
00313 check( islit = xsh_the_map_get_slit_index( the_tab, j));
00314
00315 if ( (nb_slit > 1) || (islit == 4) ){
00316
00317 check( order = (double)xsh_the_map_get_order(the_tab, j));
00318 check( xtheval = xsh_the_map_get_detx(the_tab, j));
00319 check( ytheval = xsh_the_map_get_dety(the_tab, j));
00320 (*lambda)[sol_size] = lambdaTHE;
00321 (*n)[sol_size] = order;
00322 (*s)[sol_size] = slit;
00323 (*s_index)[sol_size] = islit;
00324 (*xthe)[sol_size] = xtheval;
00325 (*ythe)[sol_size] = ytheval;
00326
00327 XSH_REGDEBUG("sol_size %d order %f slit %f lambda %f",
00328 sol_size, (*n)[sol_size], (*s)[sol_size], (*lambda)[sol_size]);
00329
00330 sol_size++;
00331 nb_match++;
00332 }
00333 if ( j < (the_size-1) ){
00334 j++;
00335 check( lambdaTHE = (double) xsh_the_map_get_wavelength( the_tab, j));
00336 }
00337
00338 }
00339
00340 if (nb_match == 0) {
00341 check(xsh_arclist_reject(arclist, i));
00342 }
00343 }
00344 *size = sol_size;
00345
00346 cleanup:
00347 if ( cpl_error_get_code() != CPL_ERROR_NONE){
00348 XSH_FREE( lambda);
00349 XSH_FREE( n);
00350 XSH_FREE( s);
00351 XSH_FREE( s_index);
00352 XSH_FREE( xthe);
00353 XSH_FREE( ythe);
00354 }
00355 return;
00356 }
00357
00358
00359
00360
00361
00362
00363
00364
00365
00366
00367
00368
00369
00370
00371
00372
00373
00374
00375
00376
00377
00378
00379
00380
00381
00382 static void data_wavesol_fit_with_sigma( xsh_wavesol *wavesol,
00383 double *A, double *lambda, double *n, double *s, int size,
00384 int max_iter, double min_frac, double sigma, int* rejected)
00385 {
00386 int nbiter = 0;
00387 float frac = 1;
00388 int nbrejected = 0;
00389 int new_rejected = 1;
00390 cpl_polynomial *fitpoly = NULL;
00391 cpl_vector *dispy = NULL;
00392 int * idx = NULL;
00393 int index_size = 0;
00394 int i;
00395
00396 XSH_ASSURE_NOT_NULL( wavesol);
00397 XSH_ASSURE_NOT_NULL( A);
00398 XSH_ASSURE_NOT_NULL( lambda);
00399 XSH_ASSURE_NOT_NULL( n);
00400 XSH_ASSURE_NOT_NULL( s);
00401 XSH_ASSURE_NOT_ILLEGAL( size > 0);
00402
00403 XSH_CALLOC( idx, int, size);
00404 check(fitpoly = xsh_wavesol_get_poly( wavesol));
00405 check(xsh_wavesol_compute( wavesol, size, A,
00406 &(wavesol->min_y), &(wavesol->max_y), lambda, n, s, fitpoly));
00407
00408 index_size = size;
00409 for(i=0; i< index_size; i++){
00410 idx[i] = i;
00411 }
00412
00413 xsh_msg( "Fit wavesol with sigma clipping");
00414
00415 while (nbiter < max_iter && frac > min_frac && new_rejected > 0){
00416 double sigma_med;
00417
00418 new_rejected = 0;
00419 xsh_msg_dbg_high( " *** NBITER = %d / %d ***", nbiter+1, max_iter);
00420 dispy = cpl_vector_new( index_size);
00421
00422 for(i = 0; i < index_size; i ++){
00423 double soly;
00424 double diffy;
00425 check( soly = xsh_wavesol_eval_poly(wavesol, lambda[i], n[i], s[i]));
00426 diffy = A[i]-soly;
00427 check( cpl_vector_set( dispy, i, diffy));
00428
00429 }
00430 check( sigma_med = cpl_vector_get_stdev( dispy));
00431 xsh_msg_dbg_high(" sigma %f SIGMA MEDIAN = %f", sigma, sigma_med);
00432
00433 for(i = 0; i < index_size; i++){
00434 if ( fabs(cpl_vector_get( dispy, i)) > (sigma * sigma_med) ){
00435 new_rejected++;
00436 rejected[idx[i]] = 1;
00437 }
00438 else{
00439 idx[i-new_rejected] = idx[i];
00440 A[i-new_rejected] = A[i];
00441 lambda[i-new_rejected] = lambda[i];
00442 n[i-new_rejected] = n[i];
00443 s[i-new_rejected] = s[i];
00444 }
00445 }
00446 xsh_free_vector(&dispy);
00447 index_size -= new_rejected;
00448 nbrejected += new_rejected;
00449
00450 frac = 1 - ( (float) nbrejected / (float) size);
00451 xsh_msg_dbg_high(" NB_REJECT = %d", nbrejected);
00452 xsh_msg_dbg_high(" GOOD PIXEL FRACTION = %f", frac);
00453
00454 check( xsh_wavesol_compute(wavesol, index_size, A,
00455 &(wavesol->min_y), &(wavesol->max_y), lambda, n, s, fitpoly));
00456 nbiter++;
00457 }
00458
00459 cleanup:
00460 XSH_FREE( idx);
00461 xsh_free_vector(&dispy);
00462 return;
00463 }
00464
00465
00466
00467
00468
00506 void
00507 xsh_detect_arclines_dan( cpl_frame *frame,
00508 cpl_frame *theo_tab_frame,
00509 cpl_frame *arc_lines_tab_frame,
00510 cpl_frame* wave_tab_guess_frame,
00511 cpl_frame *order_tab_recov_frame,
00512 cpl_frame *config_model_frame,
00513 cpl_frame *spectralformat_frame,
00514 cpl_frame **resid_tab_orders_frame,
00515 cpl_frame **arc_lines_clean_tab_frame,
00516 cpl_frame **wave_tab_frame,
00517 cpl_frame **resid_tab_frame,
00518 xsh_sol_wavelength solwave_type,
00519 xsh_detect_arclines_param *da,
00520 xsh_clipping_param *dac,
00521 xsh_instrument *instr,
00522 const char* rec_id,
00523 const int clean_tmp,
00524 const int resid_tab_name_sw)
00525 {
00526
00527 xsh_the_map *themap = NULL;
00528 xsh_resid_tab *resid = NULL;
00529 xsh_resid_tab *resid_orders = NULL;
00530 xsh_arclist *arclist = NULL;
00531 xsh_pre *pre = NULL;
00532 cpl_polynomial *fit2dx = NULL;
00533 cpl_propertylist* header = NULL;
00534 int * sort = NULL;
00535 double *vlambdadata = NULL, *vsdata = NULL, *vorderdata = NULL, *vsndata=NULL;
00536 int *vsindexdata = NULL;
00537 double *vxthedata = NULL, *vythedata = NULL;
00538 double *corr_x = NULL, *corr_y = NULL;
00539 double *gaussian_pos_x = NULL, *gaussian_pos_y = NULL,
00540 *gaussian_sigma_x = NULL, *gaussian_sigma_y = NULL,
00541 *gaussian_fwhm_x = NULL, *gaussian_fwhm_y = NULL, *gaussian_norm=NULL;
00542 int* flag=NULL;
00543 cpl_array *gaussian_sigma_x_arr =NULL;
00544 cpl_array *gaussian_sigma_y_arr =NULL;
00545 cpl_array *gaussian_fwhm_x_arr =NULL;
00546 cpl_array *gaussian_fwhm_y_arr =NULL;
00547
00548
00549 double *diffx = NULL, *diffy = NULL;
00550 double *diffxmean = NULL, *diffymean = NULL, *diffxsig = NULL, *diffysig = NULL;
00551 xsh_order_list *order_tab_recov = NULL;
00552
00553 int nlinematched, nlinecat_clean=0, nlinecat;
00554 int i, sol_size = 0;
00555 xsh_wavesol* wave_table = NULL;
00556 char wave_table_name[256];
00557 const char* wave_table_tag = NULL;
00558 xsh_wavesol* wave_tab_guess = NULL;
00559 xsh_xs_3 config_model;
00560 xsh_xs_3* p_xs_3;
00561 xsh_spectralformat_list *spectralformat_list = NULL;
00562 int lines_not_in_image = 0;
00563 int lines_not_good_sn = 0;
00564 int lines_not_gauss_fit = 0;
00565 int lines_not_valid_pixels = 0;
00566 int lines_too_few_ph=0;
00567
00568 int* rejected = NULL;
00569 int nb_rejected = 0;
00570 int solution_type = 0;
00571 const char* solution_type_name[2] = { "POLY", "MODEL"};
00572 int detection_mode = 0;
00573 const char* detection_mode_name[3] = { "NORMAL", "CORRECTED", "RECOVER"};
00574 int nb_pinhole;
00575 char fname[256];
00576 char rname[256];
00577 char rtag[256];
00578 const char* tag=NULL;
00579
00580 char dpr_type[256];
00581 char type[256];
00582 cpl_table* wave_trace=NULL;
00583
00584
00585 cpl_table* cfg_tab=NULL;
00586 const char* cfg_name=NULL;
00587 char new_name[256];
00588 char basename[256];
00589 cpl_propertylist* plist=NULL;
00590
00591 int min_slit_match=4;
00592 double line_devs=2.5;
00593 int found_temp=true;
00594 int starti=0;
00595 int nbflag=0;
00596 int newi=0;
00597 int j,k;
00598
00599
00600
00601 XSH_ASSURE_NOT_NULL( frame);
00602 XSH_ASSURE_NOT_NULL( arc_lines_tab_frame);
00603 XSH_ASSURE_NOT_NULL( spectralformat_frame);
00604 XSH_ASSURE_NOT_NULL( instr);
00605 XSH_ASSURE_NOT_NULL( da);
00606 XSH_ASSURE_NOT_NULL( dac);
00607
00608
00609 XSH_ASSURE_NOT_NULL( arc_lines_clean_tab_frame);
00610 XSH_ASSURE_NOT_NULL( resid_tab_frame);
00611
00612
00613
00614 xsh_msg_dbg_medium("---Detect arclines parameters");
00615 xsh_msg_dbg_medium("min_sn %f fit-window-half-size %d", da->min_sn, da->fit_window_hsize);
00616 xsh_msg_dbg_medium("Clipping sigma %f niter %d frac %f", dac->sigma, dac->niter, dac->frac);
00617
00618 check( pre = xsh_pre_load (frame, instr));
00619 check_msg( nb_pinhole = xsh_pfits_get_nb_pinhole( pre->data_header),
00620 "fail to find number of pinholes info in frame %s",
00621 cpl_frame_get_filename(frame) );
00622 check(strcpy(dpr_type,xsh_pfits_get_dpr_type(pre->data_header)));
00623
00624 if(strstr(dpr_type,"FMTCHK") != NULL) {
00625 strcpy(type,"FMTCHK_");
00626 } else if(strstr(dpr_type,"WAVE") != NULL) {
00627 strcpy(type,"WAVE_");
00628 } else if(strstr(dpr_type,"WAVE") != NULL) {
00629 strcpy(type,"ARC_");
00630 } else {
00631 strcpy(type,"");
00632 }
00633
00634
00635
00636 check( arclist = xsh_arclist_load( arc_lines_tab_frame));
00637 check( spectralformat_list = xsh_spectralformat_list_load(
00638 spectralformat_frame, instr));
00639 check(nlinecat = xsh_arclist_get_size( arclist));
00640
00641
00642 check(xsh_arclist_lambda_sort( arclist));
00643
00644
00645
00646 if ( theo_tab_frame != NULL) {
00647
00648 XSH_ASSURE_NOT_ILLEGAL( config_model_frame == NULL);
00649 solution_type = XSH_DETECT_ARCLINES_TYPE_POLY;
00650 check( themap = xsh_the_map_load( theo_tab_frame));
00651 check( xsh_the_map_lambda_order_slit_sort( themap));
00652 check( theo_tab_filter( themap, arclist, &sol_size, &vlambdadata,
00653 &vorderdata, &vsdata, &vsindexdata, &vxthedata, &vythedata, nb_pinhole));
00654
00655 }
00656 else if ( config_model_frame != NULL){
00657
00658 solution_type = XSH_DETECT_ARCLINES_TYPE_MODEL;
00659 p_xs_3=&config_model;
00660
00661
00662
00663 check(cfg_name=cpl_frame_get_filename(config_model_frame));
00664 check(plist=cpl_propertylist_load(cfg_name,0));
00665 check(cfg_tab=cpl_table_load(cfg_name,1,0));
00666
00667 tag=XSH_GET_TAG_FROM_ARM(XSH_MOD_CFG_TAB,instr);
00668 sprintf(basename,"%s.fits",tag);
00669 sprintf(new_name,"local_%s",basename);
00670 check( xsh_pfits_set_pcatg(plist,tag));
00671 check(cpl_table_save(cfg_tab,plist,NULL,new_name,CPL_IO_DEFAULT));
00672 xsh_add_temporary_file(new_name);
00673 check(cpl_frame_set_filename(config_model_frame,new_name));
00674 xsh_free_table(&cfg_tab);
00675 xsh_free_propertylist(&plist);
00676
00677 check( xsh_model_config_load_best( config_model_frame, &config_model));
00678 XSH_REGDEBUG("load config model ok");
00679
00680 check(xsh_model_temperature_update_frame(&config_model_frame,frame,
00681 instr,&found_temp));
00682 check(xsh_model_temperature_update_structure(&config_model,frame,instr));
00683
00684 check( theo_tab_model( &config_model, arclist, spectralformat_list,
00685 &sol_size, &vlambdadata, &vorderdata, &vsdata,
00686 &vsndata, &vsindexdata, &vxthedata, &vythedata,
00687 instr, nb_pinhole));
00688
00689 }
00690 else{
00691 XSH_ASSURE_NOT_ILLEGAL_MSG(1==0,
00692 "Undefined solution type (POLY or MODEL). See your input sof");
00693 }
00694
00695
00696 xsh_spectralformat_list_free(&spectralformat_list);
00697 xsh_msg_dbg_high("Solution type %s", solution_type_name[solution_type]);
00698
00699
00700 if ( wave_tab_guess_frame == NULL){
00701 if ( order_tab_recov_frame == NULL){
00702 detection_mode = XSH_DETECT_ARCLINES_MODE_NORMAL;
00703 }
00704 else{
00705 detection_mode = XSH_DETECT_ARCLINES_MODE_RECOVER;
00706 check( order_tab_recov = xsh_order_list_load( order_tab_recov_frame,
00707 instr));
00708 }
00709 }
00710 else {
00711 if ( order_tab_recov_frame == NULL){
00712 detection_mode = XSH_DETECT_ARCLINES_MODE_CORRECTED;
00713
00714 check( wave_tab_guess = xsh_wavesol_load( wave_tab_guess_frame,
00715 instr ));
00716 xsh_msg( "BinX,Y: %d, %d", wave_tab_guess->bin_x,
00717 wave_tab_guess->bin_y ) ;
00718 }
00719 else{
00720 detection_mode = -1;
00721 }
00722 }
00723
00724 XSH_ASSURE_NOT_ILLEGAL( detection_mode >=XSH_DETECT_ARCLINES_MODE_NORMAL
00725 && detection_mode <= XSH_DETECT_ARCLINES_MODE_RECOVER);
00726
00727 xsh_msg( "Detection mode : %s", detection_mode_name[detection_mode]);
00728
00729
00730 XSH_MALLOC( corr_x, double, sol_size);
00731 XSH_MALLOC( corr_y, double, sol_size);
00732
00733 XSH_MALLOC( gaussian_pos_x, double, sol_size);
00734 XSH_MALLOC( gaussian_pos_y, double, sol_size);
00735 XSH_MALLOC( gaussian_sigma_x, double, sol_size);
00736 XSH_MALLOC( gaussian_sigma_y, double, sol_size);
00737
00738 XSH_MALLOC( gaussian_fwhm_x, double, sol_size);
00739 XSH_MALLOC( gaussian_fwhm_y, double, sol_size);
00740 XSH_MALLOC( gaussian_norm, double, sol_size);
00741
00742 XSH_MALLOC( diffx, double, sol_size);
00743 XSH_MALLOC( diffy, double, sol_size);
00744
00745 XSH_MALLOC( diffxmean, double, sol_size);
00746 XSH_MALLOC( diffymean, double, sol_size);
00747
00748 XSH_MALLOC( diffxsig, double, sol_size);
00749 XSH_MALLOC( diffysig, double, sol_size);
00750 XSH_MALLOC( flag, int, sol_size);
00751
00752
00753
00754 memset(corr_x,0,sizeof(double)*sol_size);
00755 memset(corr_y,0,sizeof(double)*sol_size);
00756
00757 memset(gaussian_pos_x,0,sizeof(double)*sol_size);
00758 memset(gaussian_pos_y,0,sizeof(double)*sol_size);
00759 memset(gaussian_sigma_x,0,sizeof(double)*sol_size);
00760 memset(gaussian_sigma_y,0,sizeof(double)*sol_size);
00761 memset(gaussian_fwhm_x,0,sizeof(double)*sol_size);
00762 memset(gaussian_fwhm_y,0,sizeof(double)*sol_size);
00763 memset(gaussian_norm,0,sizeof(double)*sol_size);
00764 memset(diffx,0,sizeof(double)*sol_size);
00765 memset(diffy,0,sizeof(double)*sol_size);
00766 memset(diffxmean,0,sizeof(double)*sol_size);
00767 memset(diffymean,0,sizeof(double)*sol_size);
00768
00769 memset(diffxsig,0,sizeof(double)*sol_size);
00770 memset(diffysig,0,sizeof(double)*sol_size);
00771
00772 memset(flag,0,sizeof(int)*sol_size);
00773
00774
00775
00776 nlinematched = 0;
00777
00778 if ( da->find_center_method == XSH_GAUSSIAN_METHOD){
00779 xsh_msg_dbg_medium( "USE method GAUSSIAN");
00780 }
00781 else{
00782 xsh_msg_dbg_medium( "USE method BARYCENTER");
00783 }
00784
00785
00786
00787 for(i=0; i< sol_size; i++){
00788 int xpos = 0, ypos = 0;
00789 int slit_index;
00790 double corrxv = 0.0, corryv = 0.0;
00791 double lambda, order, slit, xthe, ythe,sn=0;
00792 double xcor, ycor, x, y, sig_x, sig_y, norm, fwhm_x, fwhm_y;
00793
00794
00795 lambda = vlambdadata[i];
00796 order = vorderdata[i];
00797 slit = vsdata[i];
00798
00799 if(vsndata != NULL) {
00800 sn = vsndata[i];
00801 }
00802
00803 slit_index = vsindexdata[i];
00804 xthe = vxthedata[i];
00805 ythe = vythedata[i];
00806 xsh_msg_dbg_high( "THE LAMBDA %f ORDER %f SLIT %f X %f Y %f", lambda,
00807 order, slit, xthe, ythe);
00808
00809
00810 xcor = xthe;
00811 ycor = ythe;
00812 xsh_msg_dbg_high("THE x%f y %f", xthe, ythe);
00813
00814
00815 if ( order_tab_recov != NULL){
00816 int iorder = 0;
00817 cpl_polynomial *center_poly_recov = NULL;
00818
00819 check (iorder = xsh_order_list_get_index_by_absorder(
00820 order_tab_recov, order));
00821 center_poly_recov = order_tab_recov->list[iorder].cenpoly;
00822 check( xcor = cpl_polynomial_eval_1d( center_poly_recov, ycor, NULL));
00823 corrxv = xcor-xthe;
00824 }
00825
00826 if ( wave_tab_guess != NULL){
00827 double diffxv, diffyv;
00828
00829 check( diffxv = xsh_wavesol_eval_polx( wave_tab_guess, lambda,
00830 order, 0));
00831 check( diffyv = xsh_wavesol_eval_poly( wave_tab_guess, lambda,
00832 order, 0));
00833
00834 xcor = xthe+diffxv;
00835 ycor = ythe+diffyv;
00836 corrxv = diffxv;
00837 corryv = diffyv;
00838 }
00839
00840 xsh_msg_dbg_high("x %f y %f CORR_x %f CORR_y %f", xcor, ycor,
00841 corrxv, corryv);
00842
00843 if ( xcor >=0.5 && ycor >=0.5 &&
00844 xcor < (pre->nx+0.5) && ycor < (pre->ny+0.5)){
00845 int best_med_res = 0;
00846 check ( best_med_res = xsh_pre_window_best_median_flux_pos( pre,
00847 (int)(xcor+0.5)-1, (int)(ycor+0.5)-1,
00848 da->search_window_hsize, da->running_median_hsize, &xpos, &ypos));
00849
00850 if (best_med_res != 0){
00851 lines_not_valid_pixels++;
00852 flag[i]=1;
00853 xsh_msg_dbg_medium("Not valid pixels for this line %d",i);
00854 continue;
00855 }
00856
00857 xpos = xpos+1;
00858 ypos = ypos+1;
00859 xsh_msg_dbg_high("ADJ x %d y %d", xpos, ypos);
00860
00861 if ( da->find_center_method == XSH_GAUSSIAN_METHOD){
00862 cpl_image_fit_gaussian(pre->data, xpos, ypos, 1+2*da->fit_window_hsize,
00863 &norm, &x, &y, &sig_x, &sig_y, &fwhm_x, &fwhm_y);
00864 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);
00865
00866
00867
00868
00869
00870
00871
00872
00873
00874
00875
00876
00877
00878
00879
00880
00881
00882
00883
00884
00885
00886
00887
00888
00889
00890
00891
00892
00893
00894
00895
00896
00897
00898
00899
00900
00901
00902
00903
00904
00905
00906
00907
00908 }
00909 else{
00910 xsh_image_find_barycenter( pre->data, xpos, ypos,
00911 1+2*da->fit_window_hsize,
00912 &norm, &x, &y, &sig_x, &sig_y, &fwhm_x, &fwhm_y);
00913 xsh_msg_dbg_high("BARY x %f y %f", x, y);
00914 }
00915
00916
00917
00918
00919
00920
00921
00922
00923
00924
00925
00926
00927
00928
00929 if (nb_pinhole==1 && solution_type==XSH_DETECT_ARCLINES_TYPE_MODEL) {
00930 if (p_xs_3->arm==2) {
00931 x=x+0.125;
00932 }
00933 else if (p_xs_3->arm==0) {
00934 x=x-0.51;
00935 }
00936 }
00937
00938 if( cpl_error_get_code() == CPL_ERROR_NONE ){
00939 int code_sn=lines_filter_by_sn( pre, da->min_sn, x, y, &sn);
00940 vlambdadata[i] = lambda;
00941 vorderdata[i] = order;
00942 vsdata[i] = slit;
00943 if(vsndata!=NULL) {
00944 vsndata[i] = sn;
00945 }
00946
00947 vsindexdata[i] = slit_index;
00948 vxthedata[i] = xthe;
00949 vythedata[i] = ythe;
00950 corr_x[i] = corrxv;
00951 corr_y[i] = corryv;
00952 gaussian_pos_x[i] = x;
00953 gaussian_pos_y[i] = y;
00954 gaussian_sigma_x[i] = sig_x;
00955 gaussian_sigma_y[i] = sig_y;
00956 gaussian_fwhm_x[i] = fwhm_x;
00957 gaussian_fwhm_y[i] = fwhm_y;
00958 gaussian_norm[i] = norm;
00959 diffx[i] = x-xthe;
00960 diffy[i] = y-ythe;
00961
00962
00963 if ( code_sn ){
00964 nlinematched++;
00965 }
00966 else{
00967 lines_not_good_sn++;
00968 flag[i]=2;
00969 }
00970 }
00971 else{
00972 flag[i]=3;
00973 lines_not_gauss_fit++;
00974 xsh_msg_dbg_medium("No Fit Gaussian for this line %d",i);
00975 xsh_error_reset();
00976 }
00977 }
00978 else{
00979 flag[i]=4;
00980 lines_not_in_image++;
00981 xsh_msg_dbg_medium("Coordinates are not in the image");
00982 }
00983 }
00984
00985
00986
00987
00988
00989
00990
00991 if (nb_pinhole>1) {
00992 for (i=1; i<sol_size; i++) {
00993 xsh_msg_dbg_high("filter 7 pinhole : line %d : lambda %f order %f slit %f",
00994 i, vlambdadata[i], vorderdata[i], vsdata[i]);
00995
00996 if (vlambdadata[i]!=vlambdadata[i-1] ||
00997 vorderdata[i]!=vorderdata[i-1] ||
00998 i==sol_size-1) {
00999
01000
01001 if (i==sol_size-1) {
01002 i++;
01003 }
01004 xsh_msg_dbg_high("filter n pinhole : from %d to %d find %d pinholes",
01005 starti, i, i-starti);
01006
01007 if (i-starti-nbflag>=min_slit_match) {
01008 for (k=starti; k<=i-1; k++) {
01009 diffxmean[k]=0.0;
01010 diffymean[k]=0.0;
01011 for (j=starti; j<=i-1; j++) {
01012 if (j!=k && flag[k]==0 && flag[j]==0) {
01013 diffxmean[k]+=diffx[j];
01014 diffymean[k]+=diffy[j];
01015 }
01016 }
01017 diffxmean[k]/=(float)(i-starti-1);
01018 diffymean[k]/=(float)(i-starti-1);
01019 }
01020 for (k=starti; k<=i-1; k++) {
01021 diffxsig[k]=0.0;
01022 diffysig[k]=0.0;
01023 for (j=starti; j<=i-1; j++) {
01024 if (j!=k && flag[k]==0 && flag[j]==0) {
01025 diffysig[k]+=(diffy[j]-diffymean[k])*(diffy[j]-diffymean[k]);
01026 diffxsig[k]+=(diffx[j]-diffxmean[k])*(diffx[j]-diffxmean[k]);
01027 }
01028 }
01029 diffxsig[k]=sqrt(diffxsig[k]/(float)(i-starti-1));
01030 diffysig[k]=sqrt(diffysig[k]/(float)(i-starti-1));
01031 }
01032 for (j=starti; j<=i-1; j++) {
01033
01034 if (fabs(diffx[j]-diffxmean[j])>line_devs*diffxsig[j] && fabs(diffy[j]-diffymean[j])>line_devs*diffysig[j] && flag[j]==0) {
01035 flag[j]=6;
01036 nbflag++;
01037 }
01038 }
01039 newi+=i-starti-nbflag;
01040 }
01041
01042 if (i-starti-nbflag<min_slit_match) {
01043 lines_too_few_ph+=i-starti;
01044 for (j=starti; j<=i-1; j++) {
01045 if (flag[j]==0) {
01046 flag[j]=5;
01047 xsh_msg_dbg_medium("Too few pin-holes for this line %d",i);
01048 }
01049 }
01050 }
01051 starti=i;
01052 nbflag=0;
01053 }
01054 if (flag[i]!=0) {
01055 nbflag+=1;
01056 }
01057 }
01058 nlinematched=newi;
01059 }
01060
01061
01062
01063
01064 xsh_msg("nlinematched / sol_size = %d / %d", nlinematched, sol_size);
01065 assure(nlinematched > 0, CPL_ERROR_ILLEGAL_INPUT,
01066 "No line matched, value of parameter "
01067 "detectarclines-search-window-half-size may be too large or too small "
01068 "detectarclines-fit-window-half-size=%d may be too large or too small "
01069 "detectarclines-running-median-half-size may be too large or too small "
01070 "detectarclines-min-sn may be too small "
01071 "or value of parameter detectarclines-min-sn may be too large",
01072 da->fit_window_hsize);
01073
01074
01075 xsh_msg_dbg_medium(" %d lines not found in image", lines_not_in_image);
01076 xsh_msg_dbg_medium(" %d lines not good S/N", lines_not_good_sn);
01077 xsh_msg_dbg_medium(" %d lines no fit gaussian", lines_not_gauss_fit);
01078 xsh_msg_dbg_medium(" %d lines no valid pixels", lines_not_valid_pixels);
01079 xsh_msg_dbg_medium(" %d lines detected in less than %d/9 ph positions", lines_too_few_ph,min_slit_match);
01080
01081
01082
01083 gaussian_sigma_x_arr = cpl_array_wrap_double(gaussian_sigma_x, nlinematched);
01084 gaussian_sigma_y_arr = cpl_array_wrap_double(gaussian_sigma_y, nlinematched);
01085 gaussian_fwhm_x_arr = cpl_array_wrap_double( gaussian_fwhm_x,nlinematched);
01086 gaussian_fwhm_y_arr = cpl_array_wrap_double( gaussian_fwhm_y,nlinematched);
01087
01088
01089
01090
01091
01092
01093
01094 for(i=0;i<nlinematched;i++){
01095 if(flag[i] > 0) {
01096 cpl_array_set_invalid(gaussian_sigma_x_arr,i);
01097 cpl_array_set_invalid(gaussian_sigma_y_arr,i);
01098 cpl_array_set_invalid(gaussian_fwhm_x_arr,i);
01099 cpl_array_set_invalid(gaussian_fwhm_y_arr,i);
01100 }
01101 }
01102
01103 xsh_msg("sigma gaussian median in x %lg",
01104 cpl_array_get_median( gaussian_sigma_x_arr));
01105 xsh_msg("sigma gaussian median in y %lg",
01106 cpl_array_get_median( gaussian_sigma_y_arr));
01107
01108 xsh_msg("FWHM gaussian median in x %lg",
01109 cpl_array_get_median( gaussian_fwhm_x_arr));
01110 xsh_msg("FWHM gaussian median in y %lg",
01111 cpl_array_get_median( gaussian_fwhm_y_arr));
01112
01113 xsh_unwrap_array( &gaussian_sigma_x_arr);
01114 xsh_unwrap_array( &gaussian_sigma_y_arr);
01115
01116 xsh_unwrap_array( &gaussian_fwhm_x_arr);
01117 xsh_unwrap_array( &gaussian_fwhm_y_arr);
01118
01119
01120
01121
01122 if ( resid_tab_orders_frame != NULL){
01123 check(resid_orders=xsh_resid_tab_create_not_flagged(sol_size,
01124 vlambdadata,
01125 vorderdata,
01126 vsdata,vsndata,
01127 vsindexdata,
01128 vxthedata, vythedata,
01129 corr_x, corr_y,
01130 gaussian_norm,
01131 gaussian_pos_x,
01132 gaussian_pos_y,
01133 gaussian_sigma_x,
01134 gaussian_sigma_y,
01135 gaussian_fwhm_x,
01136 gaussian_fwhm_y,flag,
01137 wave_table,
01138 solution_type));
01139
01140 if(resid_tab_name_sw) {
01141 sprintf(rtag,"%s%s%s",type,"RESID_TAB_ORDERS_",
01142 xsh_instrument_arm_tostring( instr ));
01143 } else {
01144 sprintf(rtag,"%s%s%s",type,"RESID_ALL_TAB_ORDERS_",
01145 xsh_instrument_arm_tostring( instr ));
01146 }
01147
01148 sprintf(rname,"%s%s",rtag,".fits");
01149
01150 check( *resid_tab_orders_frame = xsh_resid_tab_save( resid_orders,
01151 rname, instr,rtag));
01152 xsh_add_temporary_file(rname);
01153
01154 }
01155
01156
01157 xsh_msg_dbg_high("solution_type=%d poly=%d model=%d",solution_type,
01158 XSH_DETECT_ARCLINES_TYPE_POLY, XSH_DETECT_ARCLINES_TYPE_MODEL);
01159
01160
01161
01162
01163
01164 if ( solution_type == XSH_DETECT_ARCLINES_TYPE_POLY &&
01165 (wave_tab_frame != NULL)){
01166 check( wave_table = xsh_wavesol_create( spectralformat_frame, da, instr));
01167 XSH_CALLOC( rejected, int, nlinematched);
01168 if ( solwave_type == XSH_SOLUTION_RELATIVE) {
01169
01170 check( xsh_wavesol_set_type( wave_table, XSH_WAVESOL_GUESS));
01171 check( data_wavesol_fit_with_sigma( wave_table, diffy, vlambdadata,
01172 vorderdata, vsdata, nlinematched, dac->niter, dac->frac,
01173 dac->sigma, rejected));
01174 nb_rejected = 0;
01175 for(i=0; i< nlinematched; i++){
01176 if (rejected[i] == 1){
01177 nb_rejected++;
01178 }
01179 else{
01180
01181 vxthedata[i-nb_rejected] = vxthedata[i];
01182 vythedata[i-nb_rejected] = vythedata[i];
01183 vsindexdata[i-nb_rejected] = vsindexdata[i];
01184 corr_x[i-nb_rejected] = corr_x[i];
01185 corr_y[i-nb_rejected] = corr_y[i];
01186 gaussian_pos_x[i-nb_rejected] = gaussian_pos_x[i];
01187 gaussian_pos_y[i-nb_rejected] = gaussian_pos_y[i];
01188 gaussian_sigma_x[i-nb_rejected] = gaussian_sigma_x[i];
01189 gaussian_sigma_y[i-nb_rejected] = gaussian_sigma_y[i];
01190 gaussian_fwhm_x[i-nb_rejected] = gaussian_fwhm_x[i];
01191 gaussian_fwhm_y[i-nb_rejected] = gaussian_fwhm_y[i];
01192 gaussian_norm[i-nb_rejected] = gaussian_norm[i];
01193 diffx[i-nb_rejected] = diffx[i];
01194 }
01195 }
01196 nlinematched = nlinematched-nb_rejected;
01197
01198 check( fit2dx = xsh_wavesol_get_polx( wave_table));
01199 check( xsh_wavesol_compute( wave_table,
01200 nlinematched, diffx,
01201 &(wave_table->min_x), &(wave_table->max_x),
01202 vlambdadata, vorderdata, vsdata, fit2dx));
01203
01204
01205 sprintf(wave_table_name,"%s%s.fits","WAVE_TAB_GUESS_",
01206 xsh_instrument_arm_tostring( instr )) ;
01207
01208 wave_table_tag = XSH_GET_TAG_FROM_ARM( XSH_WAVE_TAB_GUESS, instr);
01209
01210 }
01211 else{
01212 check( xsh_wavesol_set_type( wave_table, XSH_WAVESOL_2D));
01213
01214 check( data_wavesol_fit_with_sigma( wave_table, gaussian_pos_y,
01215 vlambdadata, vorderdata, vsdata, nlinematched, dac->niter, dac->frac,
01216 dac->sigma, rejected));
01217 nb_rejected = 0;
01218 for(i=0; i< nlinematched; i++){
01219 if (rejected[i] == 1){
01220 nb_rejected++;
01221 }
01222 else{
01223
01224 vxthedata[i-nb_rejected] = vxthedata[i];
01225 vythedata[i-nb_rejected] = vythedata[i];
01226 vsindexdata[i-nb_rejected] = vsindexdata[i];
01227 corr_x[i-nb_rejected] = corr_x[i];
01228 corr_y[i-nb_rejected] = corr_y[i];
01229 gaussian_pos_x[i-nb_rejected] = gaussian_pos_x[i];
01230 gaussian_sigma_x[i-nb_rejected] = gaussian_sigma_x[i];
01231 gaussian_sigma_y[i-nb_rejected] = gaussian_sigma_y[i];
01232 gaussian_fwhm_x[i-nb_rejected] = gaussian_fwhm_x[i];
01233 gaussian_fwhm_y[i-nb_rejected] = gaussian_fwhm_y[i];
01234 }
01235 }
01236 nlinematched = nlinematched-nb_rejected;
01237
01238 check(fit2dx = xsh_wavesol_get_polx(wave_table));
01239 check(xsh_wavesol_compute(wave_table, nlinematched, gaussian_pos_x,
01240 &wave_table->min_x, &wave_table->max_x, vlambdadata, vorderdata,
01241 vsdata, fit2dx));
01242
01243 sprintf(wave_table_name,"%s%s.fits","WAVE_TAB_2D_",
01244 xsh_instrument_arm_tostring( instr )) ;
01245
01246
01247 wave_table_tag = XSH_GET_TAG_FROM_ARM( XSH_WAVE_TAB_2D, instr);
01248 }
01249
01250
01251
01252
01253 check(header = xsh_wavesol_get_header(wave_table));
01254 check(xsh_pfits_set_qc_nlinecat(header,nlinecat));
01255 xsh_msg("poly nlinecat=%d",nlinecat);
01256
01257 check(xsh_pfits_set_qc_nlinefound(header,sol_size));
01258
01259
01260
01261 check( wave_trace=xsh_wavesol_trace(wave_table,vlambdadata,vorderdata,
01262 vsdata,nlinematched));
01263
01264 check( *wave_tab_frame=xsh_wavesol_save(wave_table,wave_trace,
01265 wave_table_name,wave_table_tag));
01266
01267 check( cpl_frame_set_tag( *wave_tab_frame, wave_table_tag));
01268 }
01269
01270
01271
01272 xsh_arclist_clean_from_list_not_flagged( arclist, vlambdadata, flag,sol_size);
01273
01274 nlinecat_clean = arclist->size-arclist->nbrejected;
01275 check( header = xsh_arclist_get_header( arclist));
01276 check( xsh_pfits_set_qc_nlinecat( header,nlinecat));
01277
01278 check( xsh_pfits_set_qc_nlinefound( header, sol_size));
01279 check( xsh_pfits_set_qc_nlinecat_clean( header,nlinecat_clean));
01280 check( xsh_pfits_set_qc_nlinefound_clean( header,nlinematched));
01281 if(strcmp(rec_id,"xsh_predict") == 0) {
01282 tag=XSH_GET_TAG_FROM_ARM( XSH_ARC_LINE_LIST_PREDICT, instr);
01283 } else {
01284 tag=XSH_GET_TAG_FROM_ARM( XSH_ARC_LINE_LIST_2DMAP, instr);
01285 }
01286
01287 sprintf(fname,"%s%s",tag,".fits");
01288 check( *arc_lines_clean_tab_frame = xsh_arclist_save( arclist,fname,tag));
01289 xsh_add_temporary_file(fname);
01290
01291
01292
01293 check( resid = xsh_resid_tab_create( sol_size, vlambdadata, vorderdata,
01294 vsdata,vsndata,vsindexdata,vxthedata,vythedata,
01295 corr_x,corr_y,gaussian_norm,
01296 gaussian_pos_x,gaussian_pos_y,
01297 gaussian_sigma_x, gaussian_sigma_y,
01298 gaussian_fwhm_x, gaussian_fwhm_y,flag,
01299 wave_table,solution_type));
01300
01301 if(resid_tab_name_sw) {
01302 sprintf(rtag,"%s%s%s",type,"RESID_TAB_LINES_",
01303 xsh_instrument_arm_tostring( instr ));
01304 } else {
01305 sprintf(rtag,"%s%s%s",type,"RESID_TAB_ALL_LINES_",
01306 xsh_instrument_arm_tostring( instr ));
01307 }
01308 check( tag = cpl_frame_get_tag( frame));
01309
01310 XSH_REGDEBUG("TAG %s", tag);
01311
01312 if ( strstr( tag, XSH_AFC_CAL) ){
01313 sprintf(rname,"AFC_CAL_%s%s",rtag,".fits") ;
01314 }
01315 else if ( strstr( tag, XSH_AFC_ATT) ){
01316 sprintf(rname,"AFC_ATT_%s%s",rtag,".fits") ;
01317 }
01318 else {
01319 sprintf(rname,"%s%s",rtag,".fits") ;
01320 }
01321 check( *resid_tab_frame = xsh_resid_tab_save( resid,rname,instr,rtag));
01322 if(clean_tmp) {
01323 xsh_add_temporary_file(rname);
01324 }
01325
01326 cleanup:
01327 XSH_FREE( vlambdadata);
01328 XSH_FREE( vorderdata);
01329 XSH_FREE( vsdata);
01330 XSH_FREE( vsndata);
01331 XSH_FREE( vsindexdata);
01332 XSH_FREE( vxthedata);
01333 XSH_FREE( vythedata);
01334 XSH_FREE( corr_x);
01335 XSH_FREE( corr_y);
01336 XSH_FREE( gaussian_pos_x);
01337 XSH_FREE( gaussian_pos_y);
01338 XSH_FREE( gaussian_sigma_x);
01339 XSH_FREE( gaussian_sigma_y);
01340 XSH_FREE( gaussian_fwhm_x);
01341 XSH_FREE( gaussian_fwhm_y);
01342 XSH_FREE( gaussian_norm);
01343 XSH_FREE( sort);
01344 XSH_FREE( diffx);
01345 XSH_FREE( diffy);
01346 XSH_FREE( diffxmean);
01347 XSH_FREE( diffymean);
01348 XSH_FREE( diffxsig);
01349 XSH_FREE( diffysig);
01350 XSH_FREE( rejected);
01351 XSH_FREE( flag);
01352
01353 xsh_free_table( &wave_trace);
01354 xsh_spectralformat_list_free(&spectralformat_list);
01355 xsh_pre_free( &pre);
01356 xsh_the_map_free( &themap);
01357 xsh_wavesol_free( &wave_table);
01358 xsh_wavesol_free( &wave_tab_guess);
01359 xsh_order_list_free( &order_tab_recov);
01360 xsh_resid_tab_free( &resid_orders);
01361 xsh_resid_tab_free( &resid);
01362 xsh_arclist_free( &arclist);
01363 return;
01364 }
01365
01404 void
01405 xsh_detect_arclines( cpl_frame *frame,
01406 cpl_frame *theo_tab_frame,
01407 cpl_frame *arc_lines_tab_frame,
01408 cpl_frame* wave_tab_guess_frame,
01409 cpl_frame *order_tab_recov_frame,
01410 cpl_frame *config_model_frame,
01411 cpl_frame *spectralformat_frame,
01412 cpl_frame **resid_tab_orders_frame,
01413 cpl_frame **arc_lines_clean_tab_frame,
01414 cpl_frame **wave_tab_frame,
01415 cpl_frame **resid_tab_frame,
01416 xsh_sol_wavelength solwave_type,
01417 xsh_detect_arclines_param *da,
01418 xsh_clipping_param *dac,
01419 xsh_instrument *instr,
01420 const char* rec_id,
01421 const int clean_tmp,const int resid_tab_name_sw)
01422 {
01423
01424
01425
01426 xsh_the_map *themap = NULL;
01427 xsh_resid_tab *resid = NULL;
01428 xsh_resid_tab *resid_orders = NULL;
01429 xsh_arclist *arclist = NULL;
01430 xsh_pre *pre = NULL;
01431 cpl_polynomial *fit2dx = NULL;
01432 cpl_propertylist* header = NULL;
01433 int * sort = NULL;
01434 double *vlambdadata = NULL, *vsdata = NULL, *vorderdata = NULL, *vsndata=NULL;
01435 int *vsindexdata = NULL;
01436 double *vxthedata = NULL, *vythedata = NULL;
01437 double *corr_x = NULL, *corr_y = NULL;
01438 double *gaussian_pos_x = NULL, *gaussian_pos_y = NULL,
01439 *gaussian_sigma_x = NULL, *gaussian_sigma_y = NULL,
01440 *gaussian_fwhm_x = NULL, *gaussian_fwhm_y = NULL, *gaussian_norm=NULL;
01441 cpl_vector *gaussian_sigma_x_vect =NULL;
01442 cpl_vector *gaussian_sigma_y_vect =NULL;
01443 cpl_vector *gaussian_fwhm_x_vect =NULL;
01444 cpl_vector *gaussian_fwhm_y_vect =NULL;
01445
01446
01447 double *diffx = NULL, *diffy = NULL;
01448 double *diffxmean = NULL, *diffymean = NULL, *diffxsig = NULL, *diffysig = NULL;
01449 xsh_order_list *order_tab_recov = NULL;
01450
01451 int nlinematched, nlinecat_clean=0, nlinecat;
01452 int i, sol_size = 0;
01453 xsh_wavesol* wave_table = NULL;
01454 char wave_table_name[256];
01455 const char* wave_table_tag = NULL;
01456 xsh_wavesol* wave_tab_guess = NULL;
01457 xsh_xs_3 config_model;
01458 xsh_xs_3* p_xs_3;
01459 xsh_spectralformat_list *spectralformat_list = NULL;
01460 int lines_not_in_image = 0;
01461 int lines_not_good_sn = 0;
01462 int lines_not_gauss_fit = 0;
01463 int lines_not_valid_pixels = 0;
01464 int lines_too_few_ph=0;
01465
01466 int* rejected = NULL;
01467 int nb_rejected = 0;
01468 int solution_type = 0;
01469 const char* solution_type_name[2] = { "POLY", "MODEL"};
01470 int detection_mode = 0;
01471 const char* detection_mode_name[3] = { "NORMAL", "CORRECTED", "RECOVER"};
01472 int nb_pinhole;
01473 char fname[256];
01474 char rname[256];
01475 char rtag[256];
01476 const char* tag=NULL;
01477
01478 char dpr_type[256];
01479 char type[256];
01480 cpl_table* wave_trace=NULL;
01481 int* flag=NULL;
01482
01483 cpl_table* cfg_tab=NULL;
01484 const char* cfg_name=NULL;
01485 char new_name[256];
01486 char basename[256];
01487
01488
01489 cpl_propertylist* plist=NULL;
01490
01491
01492
01493 int starti=0;
01494 int newi=0;
01495 int j,k;
01496 int min_slit_match=7;
01497 int found_temp=true;
01498
01499
01500 XSH_ASSURE_NOT_NULL( frame);
01501 XSH_ASSURE_NOT_NULL( arc_lines_tab_frame);
01502 XSH_ASSURE_NOT_NULL( spectralformat_frame);
01503 XSH_ASSURE_NOT_NULL( instr);
01504 XSH_ASSURE_NOT_NULL( da);
01505 XSH_ASSURE_NOT_NULL( dac);
01506
01507
01508 XSH_ASSURE_NOT_NULL( arc_lines_clean_tab_frame);
01509 XSH_ASSURE_NOT_NULL( resid_tab_frame);
01510
01511
01512
01513 xsh_msg_dbg_medium("---Detect arclines parameters");
01514 xsh_msg_dbg_medium("min_sn %f fit-window-half-size %d", da->min_sn, da->fit_window_hsize);
01515 xsh_msg_dbg_medium("Clipping sigma %f niter %d frac %f", dac->sigma, dac->niter, dac->frac);
01516
01517 check( pre = xsh_pre_load (frame, instr));
01518 check_msg( nb_pinhole = xsh_pfits_get_nb_pinhole( pre->data_header),
01519 "fail to find number of pinholes info in frame %s",
01520 cpl_frame_get_filename(frame) );
01521 check(strcpy(dpr_type,xsh_pfits_get_dpr_type(pre->data_header)));
01522
01523 if(strstr(dpr_type,"FMTCHK") != NULL) {
01524 strcpy(type,"FMTCHK_");
01525 } else if(strstr(dpr_type,"WAVE") != NULL) {
01526 strcpy(type,"WAVE_");
01527 } else if(strstr(dpr_type,"WAVE") != NULL) {
01528 strcpy(type,"ARC_");
01529 } else {
01530 strcpy(type,"");
01531 }
01532
01533 check( arclist = xsh_arclist_load( arc_lines_tab_frame));
01534 check( spectralformat_list = xsh_spectralformat_list_load(
01535 spectralformat_frame, instr));
01536 check(nlinecat = xsh_arclist_get_size( arclist));
01537
01538 check(xsh_arclist_lambda_sort( arclist));
01539
01540
01541
01542 if ( theo_tab_frame != NULL) {
01543 XSH_ASSURE_NOT_ILLEGAL( config_model_frame == NULL);
01544 solution_type = XSH_DETECT_ARCLINES_TYPE_POLY;
01545 check( themap = xsh_the_map_load( theo_tab_frame));
01546 check( xsh_the_map_lambda_order_slit_sort( themap));
01547 check( theo_tab_filter( themap, arclist, &sol_size, &vlambdadata,
01548 &vorderdata, &vsdata, &vsindexdata, &vxthedata, &vythedata, nb_pinhole));
01549
01550 }
01551 else if ( config_model_frame != NULL){
01552 solution_type = XSH_DETECT_ARCLINES_TYPE_MODEL;
01553 p_xs_3=&config_model;
01554
01555
01556
01557 check(cfg_name=cpl_frame_get_filename(config_model_frame));
01558 check(plist=cpl_propertylist_load(cfg_name,0));
01559 check(cfg_tab=cpl_table_load(cfg_name,1,0));
01560
01561 tag=XSH_GET_TAG_FROM_ARM(XSH_MOD_CFG_TAB,instr);
01562 sprintf(basename,"%s.fits",tag);
01563 sprintf(new_name,"local_%s",basename);
01564 check( xsh_pfits_set_pcatg(plist,tag));
01565 check(cpl_table_save(cfg_tab,plist,NULL,new_name,CPL_IO_DEFAULT));
01566 xsh_add_temporary_file(new_name);
01567 check(cpl_frame_set_filename(config_model_frame,new_name));
01568 xsh_free_table(&cfg_tab);
01569 xsh_free_propertylist(&plist);
01570
01571 check( xsh_model_config_load_best( config_model_frame, &config_model));
01572 XSH_REGDEBUG("load config model ok");
01573
01574 check(xsh_model_temperature_update_frame(&config_model_frame,frame,
01575 instr,&found_temp));
01576 check(xsh_model_temperature_update_structure(&config_model,frame,instr));
01577 check( theo_tab_model( &config_model, arclist, spectralformat_list,
01578 &sol_size, &vlambdadata, &vorderdata, &vsdata,
01579 &vsndata, &vsindexdata, &vxthedata, &vythedata,
01580 instr, nb_pinhole));
01581 }
01582 else{
01583 XSH_ASSURE_NOT_ILLEGAL_MSG(1==0,
01584 "Undefined solution type (POLY or MODEL). See your input sof");
01585 }
01586
01587 xsh_spectralformat_list_free(&spectralformat_list);
01588 xsh_msg_dbg_high("Solution type %s", solution_type_name[solution_type]);
01589
01590
01591
01592
01593
01594 if ( wave_tab_guess_frame == NULL){
01595 if ( order_tab_recov_frame == NULL){
01596 detection_mode = XSH_DETECT_ARCLINES_MODE_NORMAL;
01597 }
01598 else{
01599 detection_mode = XSH_DETECT_ARCLINES_MODE_RECOVER;
01600 check( order_tab_recov = xsh_order_list_load( order_tab_recov_frame,
01601 instr));
01602 }
01603 }
01604 else {
01605 if ( order_tab_recov_frame == NULL){
01606 detection_mode = XSH_DETECT_ARCLINES_MODE_CORRECTED;
01607
01608 check( wave_tab_guess = xsh_wavesol_load( wave_tab_guess_frame,
01609 instr ));
01610 xsh_msg( "BinX,Y: %d, %d", wave_tab_guess->bin_x,
01611 wave_tab_guess->bin_y ) ;
01612 }
01613 else{
01614 detection_mode = -1;
01615 }
01616 }
01617
01618 XSH_ASSURE_NOT_ILLEGAL( detection_mode >=XSH_DETECT_ARCLINES_MODE_NORMAL
01619 && detection_mode <= XSH_DETECT_ARCLINES_MODE_RECOVER);
01620
01621 xsh_msg( "Detection mode : %s", detection_mode_name[detection_mode]);
01622
01623 xsh_msg( "Solution size %d", sol_size);
01624
01625
01626 XSH_MALLOC( corr_x, double, sol_size);
01627 XSH_MALLOC( corr_y, double, sol_size);
01628
01629 XSH_MALLOC( gaussian_pos_x, double, sol_size);
01630 XSH_MALLOC( gaussian_pos_y, double, sol_size);
01631 XSH_MALLOC( gaussian_sigma_x, double, sol_size);
01632 XSH_MALLOC( gaussian_sigma_y, double, sol_size);
01633
01634 XSH_MALLOC( gaussian_fwhm_x, double, sol_size);
01635 XSH_MALLOC( gaussian_fwhm_y, double, sol_size);
01636 XSH_MALLOC( gaussian_norm, double, sol_size);
01637
01638 XSH_MALLOC( diffx, double, sol_size);
01639 XSH_MALLOC( diffy, double, sol_size);
01640
01641 XSH_MALLOC( diffxmean, double, sol_size);
01642 XSH_MALLOC( diffymean, double, sol_size);
01643
01644 XSH_MALLOC( diffxsig, double, sol_size);
01645 XSH_MALLOC( diffysig, double, sol_size);
01646
01647
01648 nlinematched = 0;
01649
01650 if ( da->find_center_method == XSH_GAUSSIAN_METHOD){
01651 xsh_msg_dbg_medium( "USE method GAUSSIAN");
01652 }
01653 else{
01654 xsh_msg_dbg_medium( "USE method BARYCENTER");
01655 }
01656
01657
01658 if (xsh_debug_level_get() >= XSH_DEBUG_LEVEL_HIGH) {
01659 FILE* regfile = NULL;
01660
01661 regfile = fopen( "FIT.reg", "w");
01662 fprintf( regfile, "# Region file format: DS9 version 4.0\n"\
01663 "global color=red font=\"helvetica 4 normal\""\
01664 "select=1 highlite=1 edit=1 move=1 delete=1 include=1 fixed=0 "\
01665 "source\nimage\n");
01666 fclose( regfile);
01667 regfile = fopen( "NOFIT.reg", "w");
01668 fprintf( regfile, "# Region file format: DS9 version 4.0\n"\
01669 "global color=red font=\"helvetica 4 normal\""\
01670 "select=1 highlite=1 edit=1 move=1 delete=1 include=1 fixed=0 "\
01671 "source\nimage\n");
01672 fclose( regfile);
01673 }
01674
01675 for(i=0; i< sol_size; i++){
01676 int xpos = 0, ypos = 0;
01677 int slit_index;
01678 double corrxv = 0.0, corryv = 0.0;
01679 double lambda, order, slit, xthe, ythe,sn=0;
01680 double xcor, ycor, x, y, sig_x, sig_y, norm, fwhm_x, fwhm_y;
01681
01682
01683 lambda = vlambdadata[i];
01684 order = vorderdata[i];
01685 slit = vsdata[i];
01686
01687 if(vsndata != NULL) {
01688 sn = vsndata[i];
01689 }
01690
01691 slit_index = vsindexdata[i];
01692 xthe = vxthedata[i];
01693 ythe = vythedata[i];
01694 xsh_msg_dbg_high( "THE LAMBDA %f ORDER %f SLIT %f X %f Y %f", lambda,
01695 order, slit, xthe, ythe);
01696
01697
01698 xcor = xthe;
01699 ycor = ythe;
01700 xsh_msg_dbg_high("THE x%f y %f", xthe, ythe);
01701
01702
01703 if ( order_tab_recov != NULL){
01704 int iorder = 0;
01705 cpl_polynomial *center_poly_recov = NULL;
01706
01707 check (iorder = xsh_order_list_get_index_by_absorder(
01708 order_tab_recov, order));
01709 center_poly_recov = order_tab_recov->list[iorder].cenpoly;
01710 check( xcor = cpl_polynomial_eval_1d( center_poly_recov, ycor, NULL));
01711 corrxv = xcor-xthe;
01712 }
01713
01714 if ( wave_tab_guess != NULL){
01715 double diffxv, diffyv;
01716
01717 check( diffxv = xsh_wavesol_eval_polx( wave_tab_guess, lambda,
01718 order, 0));
01719 check( diffyv = xsh_wavesol_eval_poly( wave_tab_guess, lambda,
01720 order, 0));
01721
01722 xcor = xthe+diffxv;
01723 ycor = ythe+diffyv;
01724 corrxv = diffxv;
01725 corryv = diffyv;
01726 }
01727
01728 xsh_msg_dbg_high("x %f y %f CORR_x %f CORR_y %f", xcor, ycor,
01729 corrxv, corryv);
01730
01731 if ( xcor >=0.5 && ycor >=0.5 &&
01732 xcor < (pre->nx+0.5) && ycor < (pre->ny+0.5)){
01733 int best_med_res = 0;
01734 check ( best_med_res = xsh_pre_window_best_median_flux_pos( pre,
01735 (int)(xcor+0.5)-1, (int)(ycor+0.5)-1,
01736 da->search_window_hsize, da->running_median_hsize, &xpos, &ypos));
01737
01738 if (best_med_res != 0){
01739 lines_not_valid_pixels++;
01740 continue;
01741 }
01742
01743 xpos = xpos+1;
01744 ypos = ypos+1;
01745 xsh_msg_dbg_high("ADJ x %d y %d", xpos, ypos);
01746
01747 if ( da->find_center_method == XSH_GAUSSIAN_METHOD){
01748 cpl_image_fit_gaussian(pre->data, xpos, ypos, 1+2*da->fit_window_hsize,
01749 &norm, &x, &y, &sig_x, &sig_y, &fwhm_x, &fwhm_y);
01750 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);
01751
01752
01753
01754
01755
01756
01757
01758
01759
01760
01761
01762
01763
01764
01765
01766
01767
01768
01769
01770
01771
01772
01773
01774
01775
01776
01777
01778
01779
01780
01781
01782
01783
01784
01785
01786
01787
01788
01789
01790
01791
01792
01793
01794 }
01795 else{
01796 xsh_image_find_barycenter( pre->data, xpos, ypos,
01797 1+2*da->fit_window_hsize,
01798 &norm, &x, &y, &sig_x, &sig_y, &fwhm_x, &fwhm_y);
01799 xsh_msg_dbg_high("BARY x %f y %f", x, y);
01800 }
01801
01802
01803
01804
01805
01806
01807
01808
01809
01810
01811
01812
01813
01814
01815 if (nb_pinhole==1 && solution_type==XSH_DETECT_ARCLINES_TYPE_MODEL) {
01816 if (p_xs_3->arm==2) {
01817 x=x+0.125;
01818 }
01819 else if (p_xs_3->arm==0) {
01820 x=x-0.51;
01821 }
01822 }
01823
01824 if (xsh_debug_level_get() >= XSH_DEBUG_LEVEL_HIGH) {
01825 const char* fit = NULL;
01826 const char* color = NULL;
01827 FILE* regfile = NULL;
01828
01829 if (cpl_error_get_code() == CPL_ERROR_NONE){
01830 fit = "FIT.reg";
01831 color = "green";
01832 }
01833 else{
01834 fit = "NOFIT.reg";
01835 color = "red";
01836 }
01837 regfile = fopen( fit, "a");
01838 fprintf( regfile, "point(%d,%d) #point=cross color=blue\n", (int)(xcor), (int)(ycor));
01839 fprintf( regfile, "box(%d,%d,%d,%d) #color=blue\n", (int)(xcor), (int)(ycor),
01840 da->search_window_hsize*2,da->search_window_hsize*2);
01841 fprintf( regfile, "point(%d,%d) #point=cross color=%s font=\"helvetica 10 normal\""\
01842 " text={%.3f}\n", xpos, ypos, color, lambda);
01843 fprintf( regfile, "box(%d,%d,%d,%d) #color=%s\n", (int)xpos, (int)ypos,
01844 1+2*da->fit_window_hsize, 1+2*da->fit_window_hsize, color);
01845 fclose( regfile);
01846 }
01847
01848 if( cpl_error_get_code() == CPL_ERROR_NONE ){
01849 if ( lines_filter_by_sn( pre, da->min_sn, x, y, &sn) ){
01850 vlambdadata[nlinematched] = lambda;
01851 vorderdata[nlinematched] = order;
01852 vsdata[nlinematched] = slit;
01853 if(vsndata!=NULL) {
01854 vsndata[nlinematched] = sn;
01855 }
01856
01857 vsindexdata[nlinematched] = slit_index;
01858 vxthedata[nlinematched] = xthe;
01859 vythedata[nlinematched] = ythe;
01860 corr_x[nlinematched] = corrxv;
01861 corr_y[nlinematched] = corryv;
01862 gaussian_pos_x[nlinematched] = x;
01863 gaussian_pos_y[nlinematched] = y;
01864 gaussian_sigma_x[nlinematched] = sig_x;
01865 gaussian_sigma_y[nlinematched] = sig_y;
01866 gaussian_fwhm_x[nlinematched] = fwhm_x;
01867 gaussian_fwhm_y[nlinematched] = fwhm_y;
01868 gaussian_norm[nlinematched] = norm;
01869 diffx[nlinematched] = x-xthe;
01870 diffy[nlinematched] = y-ythe;
01871
01872
01873
01874
01875
01876
01877 nlinematched++;
01878
01879 }
01880 else{
01881 lines_not_good_sn++;
01882 xsh_msg_dbg_medium("Not good s/n for this line");
01883 if (xsh_debug_level_get() >= XSH_DEBUG_LEVEL_HIGH) {
01884
01885 int xpix, ypix, rej;
01886 float flux, noise, sn;
01887 FILE* regfile = NULL;
01888 char sn_name[256];
01889
01890 xpix =(int) rint(x);
01891 ypix = (int)rint(y);
01892
01893 check( flux = cpl_image_get( pre->data, xpix, ypix, &rej));
01894 check( noise = cpl_image_get( pre->errs, xpix, ypix, &rej));
01895 sn = flux / noise;
01896 sprintf( sn_name, "bad_sn_%.3f.reg", da->min_sn);
01897 regfile = fopen( sn_name, "a");
01898 fprintf( regfile, "point(%d,%d) #point=cross color=red font=\"helvetica 10 normal\""\
01899 " text={%.3f [%.3f}\n", xpix, ypix, lambda, sn);
01900 fclose( regfile);
01901 }
01902 }
01903 }
01904 else{
01905 lines_not_gauss_fit++;
01906 xsh_msg_dbg_medium("No Fit Gaussian for this line");
01907 xsh_error_reset();
01908 }
01909 }
01910 else{
01911 lines_not_in_image++;
01912 xsh_msg_dbg_medium("Coordinates are not in the image");
01913 }
01914 }
01915
01916
01917
01918
01919
01920 if (nb_pinhole>1) {
01921 for (i=1; i<nlinematched; i++) {
01922 xsh_msg_dbg_high("filter 7 pinhole : line %d : lambda %f order %f slit %f",
01923 i, vlambdadata[i], vorderdata[i], vsdata[i]);
01924
01925 if (vlambdadata[i]!=vlambdadata[i-1] || vorderdata[i]!=vorderdata[i-1] || i==nlinematched-1) {
01926
01927
01928 if (i==nlinematched-1) {
01929 i++;
01930 }
01931 xsh_msg_dbg_high("filter 7 pinhole : from %d to %d find %d pinholes",
01932 starti, i, i-starti);
01933
01934 if (i-starti>=min_slit_match) {
01935 for (k=starti; k<=i-1; k++) {
01936 diffxmean[k]=0.0;
01937 diffymean[k]=0.0;
01938 for (j=starti; j<=i-1; j++) {
01939 if (j!=k) {
01940 diffxmean[k]+=diffx[j];
01941 diffymean[k]+=diffy[j];
01942 }
01943 }
01944 diffxmean[k]/=(float)(i-starti-1);
01945 diffymean[k]/=(float)(i-starti-1);
01946 }
01947 for (k=starti; k<=i-1; k++) {
01948 diffxsig[k]=0.0;
01949 diffysig[k]=0.0;
01950 for (j=starti; j<=i-1; j++) {
01951 if (j!=k) {
01952 diffysig[k]+=(diffy[j]-diffymean[k])*(diffy[j]-diffymean[k]);
01953 diffxsig[k]+=(diffx[j]-diffxmean[k])*(diffx[j]-diffxmean[k]);
01954 }
01955 }
01956 diffxsig[k]=sqrt(diffxsig[k]/(float)(i-starti-1));
01957 diffysig[k]=sqrt(diffysig[k]/(float)(i-starti-1));
01958 }
01959 for (j=starti; j<=i-1; j++) {
01960
01961
01962 vlambdadata[newi] = vlambdadata[j];
01963 vorderdata[newi] = vorderdata[j];
01964 vsdata[newi] =vsdata[j];
01965 vsindexdata[newi] = vsindexdata[j];
01966 vxthedata[newi] =vxthedata[j];
01967 vythedata[newi] =vythedata[j];
01968 corr_x[newi] = corr_x[j];
01969 corr_y[newi] = corr_y[j];
01970 gaussian_pos_x[newi] = gaussian_pos_x[j];
01971 gaussian_pos_y[newi] = gaussian_pos_y[j];
01972 gaussian_sigma_x[newi] = gaussian_sigma_x[j];
01973 gaussian_sigma_y[newi] = gaussian_sigma_y[j];
01974 gaussian_fwhm_x[newi] = gaussian_fwhm_x[j];
01975 gaussian_fwhm_y[newi] = gaussian_fwhm_y[j];
01976 gaussian_norm[newi] = gaussian_norm[j];
01977 diffx[newi] = diffx[j];
01978 diffy[newi] = diffy[j];
01979 newi++;
01980
01981
01982
01983
01984
01985 }
01986 }
01987 else {
01988 lines_too_few_ph+=i-starti;
01989 }
01990 starti=i;
01991 }
01992 }
01993 nlinematched=newi;
01994 }
01995
01996 xsh_msg("nlinematched / sol_size = %d / %d", nlinematched, sol_size);
01997 assure(nlinematched > 0, CPL_ERROR_ILLEGAL_INPUT,
01998 "No line matched, "
01999 "detectarclines-search-window-half-size, "
02000 "detectarclines-fit-window-half-size=%d, "
02001 "detectarclines-running-median-half-size too large or too small "
02002 "or detectarclines-min-sn may be too large",
02003 da->fit_window_hsize);
02004
02005 xsh_msg_dbg_medium(" %d lines not found in image", lines_not_in_image);
02006 xsh_msg_dbg_medium(" %d lines not good S/N", lines_not_good_sn);
02007 xsh_msg_dbg_medium(" %d lines no fit gaussian", lines_not_gauss_fit);
02008 xsh_msg_dbg_medium(" %d lines no valid pixels", lines_not_valid_pixels);
02009 xsh_msg_dbg_medium(" %d lines detected in less than %d/9 ph positions", lines_too_few_ph,min_slit_match);
02010 check( gaussian_sigma_x_vect = cpl_vector_wrap( nlinematched, gaussian_sigma_x));
02011 check( gaussian_sigma_y_vect = cpl_vector_wrap( nlinematched, gaussian_sigma_y));
02012 check( gaussian_fwhm_x_vect = cpl_vector_wrap( nlinematched, gaussian_fwhm_x));
02013 check( gaussian_fwhm_y_vect = cpl_vector_wrap( nlinematched, gaussian_fwhm_y));
02014 xsh_msg("sigma gaussian median in x %lg",
02015 cpl_vector_get_median_const( gaussian_sigma_x_vect));
02016 xsh_msg("sigma gaussian median in y %lg",
02017 cpl_vector_get_median_const( gaussian_sigma_y_vect));
02018
02019 xsh_msg("FWHM gaussian median in x %lg",
02020 cpl_vector_get_median_const( gaussian_fwhm_x_vect));
02021 xsh_msg("FWHM gaussian median in y %lg",
02022 cpl_vector_get_median_const( gaussian_fwhm_y_vect));
02023
02024
02025 xsh_unwrap_vector( &gaussian_sigma_x_vect);
02026 xsh_unwrap_vector( &gaussian_sigma_y_vect);
02027
02028 xsh_unwrap_vector( &gaussian_fwhm_x_vect);
02029 xsh_unwrap_vector( &gaussian_fwhm_y_vect);
02030
02031
02032
02033
02034 if ( resid_tab_orders_frame != NULL){
02035 check(resid_orders=xsh_resid_tab_create(nlinematched,vlambdadata,vorderdata,
02036 vsdata, vsndata,vsindexdata,
02037 vxthedata, vythedata,
02038 corr_x, corr_y,
02039 gaussian_norm,
02040 gaussian_pos_x, gaussian_pos_y,
02041 gaussian_sigma_x, gaussian_sigma_y,
02042 gaussian_fwhm_x, gaussian_fwhm_y,flag,
02043 wave_table, solution_type));
02044
02045
02046 sprintf(rtag,"%s%s%s",type,"RESID_TAB_ORDERS_",
02047 xsh_instrument_arm_tostring( instr ));
02048
02049 sprintf(rname,"%s%s",rtag,".fits");
02050
02051 check( *resid_tab_orders_frame = xsh_resid_tab_save( resid_orders,
02052 rname, instr,rtag));
02053 xsh_add_temporary_file(rname);
02054
02055 }
02056 xsh_msg_dbg_high("solution_type=%d poly=%d model=%d",solution_type,
02057 XSH_DETECT_ARCLINES_TYPE_POLY, XSH_DETECT_ARCLINES_TYPE_MODEL);
02058
02059
02060
02061 if ( solution_type == XSH_DETECT_ARCLINES_TYPE_POLY &&
02062 (wave_tab_frame != NULL)){
02063 check( wave_table = xsh_wavesol_create( spectralformat_frame, da, instr));
02064 XSH_CALLOC( rejected, int, nlinematched);
02065 if ( solwave_type == XSH_SOLUTION_RELATIVE) {
02066
02067 check( xsh_wavesol_set_type( wave_table, XSH_WAVESOL_GUESS));
02068 check( data_wavesol_fit_with_sigma( wave_table, diffy, vlambdadata,
02069 vorderdata, vsdata, nlinematched, dac->niter, dac->frac,
02070 dac->sigma, rejected));
02071 nb_rejected = 0;
02072 for(i=0; i< nlinematched; i++){
02073 if (rejected[i] == 1){
02074 nb_rejected++;
02075 }
02076 else{
02077
02078 vxthedata[i-nb_rejected] = vxthedata[i];
02079 vythedata[i-nb_rejected] = vythedata[i];
02080 vsindexdata[i-nb_rejected] = vsindexdata[i];
02081 corr_x[i-nb_rejected] = corr_x[i];
02082 corr_y[i-nb_rejected] = corr_y[i];
02083 gaussian_pos_x[i-nb_rejected] = gaussian_pos_x[i];
02084 gaussian_pos_y[i-nb_rejected] = gaussian_pos_y[i];
02085 gaussian_sigma_x[i-nb_rejected] = gaussian_sigma_x[i];
02086 gaussian_sigma_y[i-nb_rejected] = gaussian_sigma_y[i];
02087 gaussian_fwhm_x[i-nb_rejected] = gaussian_fwhm_x[i];
02088 gaussian_fwhm_y[i-nb_rejected] = gaussian_fwhm_y[i];
02089 gaussian_norm[i-nb_rejected] = gaussian_norm[i];
02090 diffx[i-nb_rejected] = diffx[i];
02091 }
02092 }
02093 nlinematched = nlinematched-nb_rejected;
02094
02095 check( fit2dx = xsh_wavesol_get_polx( wave_table));
02096 check( xsh_wavesol_compute( wave_table,
02097 nlinematched, diffx,
02098 &(wave_table->min_x), &(wave_table->max_x),
02099 vlambdadata, vorderdata, vsdata, fit2dx));
02100
02101
02102 sprintf(wave_table_name,"%s%s.fits","WAVE_TAB_GUESS_",
02103 xsh_instrument_arm_tostring( instr )) ;
02104
02105 wave_table_tag = XSH_GET_TAG_FROM_ARM( XSH_WAVE_TAB_GUESS, instr);
02106
02107 }
02108 else{
02109 check( xsh_wavesol_set_type( wave_table, XSH_WAVESOL_2D));
02110
02111 check( data_wavesol_fit_with_sigma( wave_table, gaussian_pos_y,
02112 vlambdadata, vorderdata, vsdata, nlinematched, dac->niter, dac->frac,
02113 dac->sigma, rejected));
02114 nb_rejected = 0;
02115 for(i=0; i< nlinematched; i++){
02116 if (rejected[i] == 1){
02117 nb_rejected++;
02118 }
02119 else{
02120
02121 vxthedata[i-nb_rejected] = vxthedata[i];
02122 vythedata[i-nb_rejected] = vythedata[i];
02123 vsindexdata[i-nb_rejected] = vsindexdata[i];
02124 corr_x[i-nb_rejected] = corr_x[i];
02125 corr_y[i-nb_rejected] = corr_y[i];
02126 gaussian_pos_x[i-nb_rejected] = gaussian_pos_x[i];
02127 gaussian_sigma_x[i-nb_rejected] = gaussian_sigma_x[i];
02128 gaussian_sigma_y[i-nb_rejected] = gaussian_sigma_y[i];
02129 gaussian_fwhm_x[i-nb_rejected] = gaussian_fwhm_x[i];
02130 gaussian_fwhm_y[i-nb_rejected] = gaussian_fwhm_y[i];
02131 }
02132 }
02133 nlinematched = nlinematched-nb_rejected;
02134
02135 check(fit2dx = xsh_wavesol_get_polx(wave_table));
02136 check(xsh_wavesol_compute(wave_table, nlinematched, gaussian_pos_x,
02137 &wave_table->min_x, &wave_table->max_x, vlambdadata, vorderdata,
02138 vsdata, fit2dx));
02139
02140 sprintf(wave_table_name,"%s%s.fits","WAVE_TAB_2D_",
02141 xsh_instrument_arm_tostring( instr )) ;
02142
02143
02144 wave_table_tag = XSH_GET_TAG_FROM_ARM( XSH_WAVE_TAB_2D, instr);
02145 }
02146
02147
02148 check(header = xsh_wavesol_get_header(wave_table));
02149 check(xsh_pfits_set_qc_nlinecat(header,nlinecat));
02150 xsh_msg("nlinecat=%d",nlinecat);
02151
02152 check(xsh_pfits_set_qc_nlinefound(header,sol_size));
02153
02154
02155
02156 check( wave_trace=xsh_wavesol_trace(wave_table,vlambdadata,vorderdata,
02157 vsdata,nlinematched));
02158
02159 check( *wave_tab_frame=xsh_wavesol_save(wave_table,wave_trace,
02160 wave_table_name,wave_table_tag));
02161
02162 check( cpl_frame_set_tag( *wave_tab_frame, wave_table_tag));
02163 }
02164
02165
02166 check( xsh_arclist_clean_from_list( arclist, vlambdadata, nlinematched));
02167 nlinecat_clean = arclist->size-arclist->nbrejected;
02168 check( header = xsh_arclist_get_header( arclist));
02169 check( xsh_pfits_set_qc_nlinecat( header,nlinecat));
02170 xsh_msg("nlinecat=%d",nlinecat);
02171
02172 check( xsh_pfits_set_qc_nlinefound( header, sol_size));
02173 check( xsh_pfits_set_qc_nlinecat_clean( header,nlinecat_clean));
02174 check( xsh_pfits_set_qc_nlinefound_clean( header,nlinematched));
02175 if(strcmp(rec_id,"xsh_predict") == 0) {
02176 tag=XSH_GET_TAG_FROM_ARM( XSH_ARC_LINE_LIST_PREDICT, instr);
02177 } else {
02178 tag=XSH_GET_TAG_FROM_ARM( XSH_ARC_LINE_LIST_2DMAP, instr);
02179 }
02180 sprintf(fname,"%s%s",tag,".fits");
02181 check( *arc_lines_clean_tab_frame = xsh_arclist_save( arclist,fname,tag));
02182 if(clean_tmp) {
02183 xsh_add_temporary_file(fname);
02184 }
02185
02186 check( resid = xsh_resid_tab_create( nlinematched, vlambdadata, vorderdata,
02187 vsdata,vsndata,vsindexdata,vxthedata,vythedata,
02188 corr_x,corr_y,gaussian_norm,
02189 gaussian_pos_x,gaussian_pos_y,
02190 gaussian_sigma_x, gaussian_sigma_y,
02191 gaussian_fwhm_x, gaussian_fwhm_y,flag,
02192 wave_table,solution_type));
02193
02194 if(resid_tab_name_sw) {
02195 sprintf(rtag,"%s%s%s",type,"RESID_TAB_DRL_LINES_",
02196 xsh_instrument_arm_tostring( instr ));
02197 } else {
02198 sprintf(rtag,"%s%s%s",type,"RESID_TAB_LINES_",
02199 xsh_instrument_arm_tostring( instr ));
02200 }
02201 check( tag = cpl_frame_get_tag( frame));
02202
02203 XSH_REGDEBUG("TAG %s", tag);
02204
02205 if ( strstr( tag, XSH_AFC_CAL) ){
02206 sprintf(rname,"AFC_CAL_%s%s",rtag,".fits") ;
02207 }
02208 else if ( strstr( tag, XSH_AFC_ATT) ){
02209 sprintf(rname,"AFC_ATT_%s%s",rtag,".fits") ;
02210 }
02211 else {
02212 sprintf(rname,"%s%s",rtag,".fits") ;
02213 }
02214 check( *resid_tab_frame = xsh_resid_tab_save( resid,rname,instr,rtag));
02215 if(clean_tmp || resid_tab_name_sw) {
02216 xsh_add_temporary_file(rname);
02217 }
02218
02219 cleanup:
02220 XSH_FREE( vlambdadata);
02221 XSH_FREE( vorderdata);
02222 XSH_FREE( vsdata);
02223 XSH_FREE( vsndata);
02224 XSH_FREE( vsindexdata);
02225 XSH_FREE( vxthedata);
02226 XSH_FREE( vythedata);
02227 XSH_FREE( corr_x);
02228 XSH_FREE( corr_y);
02229 XSH_FREE( gaussian_pos_x);
02230 XSH_FREE( gaussian_pos_y);
02231 XSH_FREE( gaussian_sigma_x);
02232 XSH_FREE( gaussian_sigma_y);
02233 XSH_FREE( gaussian_fwhm_x);
02234 XSH_FREE( gaussian_fwhm_y);
02235 XSH_FREE( gaussian_norm);
02236 XSH_FREE( sort);
02237 XSH_FREE( diffx);
02238 XSH_FREE( diffy);
02239 XSH_FREE( diffxmean);
02240 XSH_FREE( diffymean);
02241 XSH_FREE( diffxsig);
02242 XSH_FREE( diffysig);
02243 XSH_FREE( rejected);
02244
02245 xsh_free_table( &wave_trace);
02246 xsh_spectralformat_list_free(&spectralformat_list);
02247 xsh_pre_free( &pre);
02248 xsh_the_map_free( &themap);
02249 xsh_wavesol_free( &wave_table);
02250 xsh_wavesol_free( &wave_tab_guess);
02251 xsh_order_list_free( &order_tab_recov);
02252 xsh_resid_tab_free( &resid_orders);
02253 xsh_resid_tab_free( &resid);
02254 xsh_arclist_free( &arclist);
02255 return;
02256 }
02257
02258