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, int **s_index, double **xthe, double **ythe,
00078 xsh_instrument* instr, int nb_pinhole);
00079
00080 static void data_wavesol_fit_with_sigma( xsh_wavesol *wavesol,
00081 double *A, double *lambda, double *n, double *s, int size,
00082 int max_iter, double min_frac, double sigma, int* rejected);
00083
00084 static int lines_filter_by_sn( xsh_pre* pre, double sn_ref, double x,
00085 double y);
00086
00087
00088
00089
00090 static int lines_filter_by_sn( xsh_pre* pre, double sn_ref, double x,
00091 double y)
00092 {
00093 int res = 0;
00094 int xpix, ypix, rej;
00095 double flux, noise, sn;
00096
00097 XSH_ASSURE_NOT_NULL( pre);
00098
00099
00100 xpix =(int) rint(x);
00101 ypix = (int)rint(y);
00102
00103 check( flux = cpl_image_get( pre->data, xpix, ypix, &rej));
00104 check( noise = cpl_image_get( pre->errs, xpix, ypix, &rej));
00105 sn = flux / noise;
00106 res = (sn > sn_ref);
00107 cleanup:
00108 return res;
00109 }
00110
00111 static void theo_tab_model( xsh_xs_3* config_model, xsh_arclist *arclist,
00112 xsh_spectralformat_list *spectralformat_list, int* size, double **lambda,
00113 double **n, double **s, int **s_index, double **xthe, double **ythe,
00114 xsh_instrument* instr, int nb_pinhole)
00115 {
00116 int i,j;
00117 cpl_vector** spectral_tab = NULL;
00118 int global_size = 0;
00119 int arc_size;
00120 int index_loc = 0;
00121 int nb_slit = 1;
00122
00123 XSH_ASSURE_NOT_NULL( config_model);
00124 XSH_ASSURE_NOT_NULL( arclist);
00125 XSH_ASSURE_NOT_NULL( spectralformat_list);
00126
00127 check( arc_size = xsh_arclist_get_size( arclist));
00128 XSH_CALLOC( spectral_tab, cpl_vector*, arc_size);
00129
00130 nb_slit = nb_pinhole;
00131
00132 XSH_REGDEBUG("nb pinhole %d ", nb_slit);
00133
00134 for( i=0; i< arc_size; i++){
00135 cpl_vector* res = NULL;
00136 float lambdaARC;
00137 int res_size;
00138
00139 check( lambdaARC = xsh_arclist_get_wavelength(arclist, i));
00140 check( res = xsh_spectralformat_list_get_orders(spectralformat_list,
00141 lambdaARC));
00142 if (res != NULL){
00143 check( res_size = cpl_vector_get_size( res));
00144 res_size *= nb_slit;
00145 }
00146 else{
00147 res_size = 0;
00148 check(xsh_arclist_reject(arclist,i));
00149 }
00150 global_size += res_size;
00151 spectral_tab[i] = res;
00152 }
00153
00154 XSH_MALLOC( *lambda, double, global_size);
00155 XSH_MALLOC( *n, double, global_size);
00156 XSH_MALLOC( *s, double, global_size);
00157 XSH_MALLOC( *s_index, int, global_size);
00158 XSH_MALLOC( *xthe, double, global_size);
00159 XSH_MALLOC( *ythe, double, global_size);
00160
00161
00162 for( i=0; i< arc_size; i++){
00163 double model_x, model_y;
00164 cpl_vector *spectral_res = NULL;
00165 int spectral_res_size;
00166 float lambdaARC;
00167
00168 check( lambdaARC = xsh_arclist_get_wavelength(arclist, i));
00169 check( spectral_res = spectral_tab[i]);
00170
00171 if (spectral_res != NULL){
00172 check( spectral_res_size = cpl_vector_get_size( spectral_res));
00173 for( j=0; j< spectral_res_size; j++){
00174 int absorder = 0;
00175 int islit;
00176 double slit;
00177
00178 if (nb_slit > 1){
00179
00180 for(islit=0; islit < nb_slit; islit++){
00181
00182 slit = config_model->slit[islit];
00183 check( absorder = (int) cpl_vector_get( spectral_res, j));
00184
00185
00186
00187 check( xsh_model_get_xy( config_model, instr, lambdaARC, absorder, slit,
00188 &model_x, &model_y));
00189 (*lambda)[index_loc] = lambdaARC;
00190 (*n)[index_loc] = absorder;
00191 (*s_index)[index_loc] = islit;
00192 (*s)[index_loc] = slit;
00193 (*xthe)[index_loc] = model_x;
00194 (*ythe)[index_loc] = model_y;
00195 index_loc++;
00196 }
00197 }
00198 else{
00199 islit = 4;
00200 slit = config_model->slit[islit];
00201 check( absorder = (int) cpl_vector_get( spectral_res, j));
00202 check( xsh_model_get_xy( config_model, instr, lambdaARC, absorder, slit,
00203 &model_x, &model_y));
00204 (*lambda)[index_loc] = lambdaARC;
00205 (*n)[index_loc] = absorder;
00206 (*s_index)[index_loc] = islit;
00207 (*s)[index_loc] = slit;
00208 (*xthe)[index_loc] = model_x;
00209 (*ythe)[index_loc] = model_y;
00210 index_loc++;
00211 }
00212
00213
00214 }
00215 }
00216 }
00217 *size = global_size;
00218
00219 cleanup:
00220 if ( spectral_tab != NULL){
00221 for(i=0; i< arc_size; i++){
00222 xsh_free_vector( &spectral_tab[i]);
00223 }
00224 XSH_FREE( spectral_tab);
00225 }
00226 if ( cpl_error_get_code() != CPL_ERROR_NONE){
00227 XSH_FREE( *lambda);
00228 XSH_FREE( *n);
00229 XSH_FREE( *s);
00230 XSH_FREE( *s_index);
00231 XSH_FREE( *xthe);
00232 XSH_FREE( *ythe);
00233 }
00234 return;
00235 }
00236
00237
00238
00239
00240
00241
00242
00243
00244
00245
00246
00247
00248
00249
00250
00251
00252
00253
00254
00255
00256
00257 static void theo_tab_filter( xsh_the_map *the_tab, xsh_arclist *arclist,
00258 int* size, double **lambda, double **n, double **s, int **s_index,
00259 double **xthe, double **ythe, int nb_pinhole)
00260 {
00261 int the_size = 0;
00262 int arc_size = 0;
00263 int i=0, j=0;
00264 int sol_size = 0;
00265 int nb_slit = 1;
00266
00267 XSH_ASSURE_NOT_NULL( the_tab);
00268 XSH_ASSURE_NOT_NULL( arclist);
00269
00270 check( arc_size = xsh_arclist_get_size( arclist));
00271 check( the_size = xsh_the_map_get_size( the_tab));
00272
00273 XSH_MALLOC( *lambda, double, the_size);
00274 XSH_MALLOC( *n, double, the_size);
00275 XSH_MALLOC( *s, double, the_size);
00276 XSH_MALLOC( *s_index, int, the_size);
00277 XSH_MALLOC( *xthe, double, the_size);
00278 XSH_MALLOC( *ythe, double, the_size);
00279
00280 nb_slit = nb_pinhole;
00281 XSH_REGDEBUG("nb pinhole %d ", nb_slit);
00282
00283 for(i=0; i< arc_size; i++){
00284 float lambdaARC, lambdaTHE;
00285 int nb_match = 0, max_match = 0;
00286
00287 check(lambdaARC = xsh_arclist_get_wavelength(arclist, i));
00288 check(lambdaTHE = xsh_the_map_get_wavelength(the_tab, j));
00289
00290 xsh_msg_dbg_medium("LINETABLE : Line %d / %d Lambda %f",i+1, arc_size,
00291 lambdaARC);
00292
00293 while ( (j < the_size-1) &&
00294 ( (lambdaARC-lambdaTHE) > WAVELENGTH_PRECISION) ){
00295 j++;
00296 check(lambdaTHE = xsh_the_map_get_wavelength(the_tab, j));
00297 }
00298
00299 xsh_msg_dbg_medium("THETABLE : Line %d / %d Lambda %f",j+1, the_size, lambdaTHE);
00300
00301 max_match = the_size-j;
00302
00303 while ( nb_match < max_match && fabs(lambdaARC-lambdaTHE) <= WAVELENGTH_PRECISION ){
00304 double order, slit, xtheval, ytheval;
00305 int islit;
00306
00307 check( slit = (double) xsh_the_map_get_slit_position( the_tab, j));
00308 check( islit = xsh_the_map_get_slit_index( the_tab, j));
00309
00310 if ( (nb_slit > 1) || (islit == 4) ){
00311
00312 check( order = (double)xsh_the_map_get_order(the_tab, j));
00313 check( xtheval = xsh_the_map_get_detx(the_tab, j));
00314 check( ytheval = xsh_the_map_get_dety(the_tab, j));
00315 (*lambda)[sol_size] = lambdaTHE;
00316 (*n)[sol_size] = order;
00317 (*s)[sol_size] = slit;
00318 (*s_index)[sol_size] = islit;
00319 (*xthe)[sol_size] = xtheval;
00320 (*ythe)[sol_size] = ytheval;
00321
00322 XSH_REGDEBUG("sol_size %d order %f slit %f lambda %f",
00323 sol_size, (*n)[sol_size], (*s)[sol_size], (*lambda)[sol_size]);
00324
00325 sol_size++;
00326 nb_match++;
00327 }
00328 if ( j < (the_size-1) ){
00329 j++;
00330 check( lambdaTHE = (double) xsh_the_map_get_wavelength( the_tab, j));
00331 }
00332
00333 }
00334
00335 if (nb_match == 0) {
00336 check(xsh_arclist_reject(arclist, i));
00337 }
00338 }
00339 *size = sol_size;
00340
00341 cleanup:
00342 if ( cpl_error_get_code() != CPL_ERROR_NONE){
00343 XSH_FREE( lambda);
00344 XSH_FREE( n);
00345 XSH_FREE( s);
00346 XSH_FREE( s_index);
00347 XSH_FREE( xthe);
00348 XSH_FREE( ythe);
00349 }
00350 return;
00351 }
00352
00353
00354
00355
00356
00357
00358
00359
00360
00361
00362
00363
00364
00365
00366
00367
00368
00369
00370
00371
00372
00373
00374
00375
00376
00377 static void data_wavesol_fit_with_sigma( xsh_wavesol *wavesol,
00378 double *A, double *lambda, double *n, double *s, int size,
00379 int max_iter, double min_frac, double sigma, int* rejected)
00380 {
00381 int nbiter = 0;
00382 float frac = 1;
00383 int nbrejected = 0;
00384 int new_rejected = 1;
00385 cpl_polynomial *fitpoly = NULL;
00386 cpl_vector *dispy = NULL;
00387 int * idx = NULL;
00388 int index_size = 0;
00389 int i;
00390
00391 XSH_ASSURE_NOT_NULL( wavesol);
00392 XSH_ASSURE_NOT_NULL( A);
00393 XSH_ASSURE_NOT_NULL( lambda);
00394 XSH_ASSURE_NOT_NULL( n);
00395 XSH_ASSURE_NOT_NULL( s);
00396 XSH_ASSURE_NOT_ILLEGAL( size > 0);
00397
00398 XSH_CALLOC( idx, int, size);
00399 check(fitpoly = xsh_wavesol_get_poly( wavesol));
00400 check(xsh_wavesol_compute( wavesol, size, A,
00401 &(wavesol->min_y), &(wavesol->max_y), lambda, n, s, fitpoly));
00402
00403 index_size = size;
00404 for(i=0; i< index_size; i++){
00405 idx[i] = i;
00406 }
00407
00408 xsh_msg( "Fit wavesol with sigma clipping");
00409
00410 while (nbiter < max_iter && frac > min_frac && new_rejected > 0){
00411 double sigma_med;
00412
00413 new_rejected = 0;
00414 xsh_msg_dbg_high( " *** NBITER = %d / %d ***", nbiter+1, max_iter);
00415 dispy = cpl_vector_new( index_size);
00416
00417 for(i = 0; i < index_size; i ++){
00418 double soly;
00419 double diffy;
00420
00421 check( soly = xsh_wavesol_eval_poly(wavesol, lambda[i], n[i], s[i]));
00422 diffy = A[i]-soly;
00423 check( cpl_vector_set( dispy, i, diffy));
00424 }
00425 check( sigma_med = cpl_vector_get_stdev( dispy));
00426 xsh_msg_dbg_high(" sigma %f SIGMA MEDIAN = %f", sigma, sigma_med);
00427
00428 for(i = 0; i < index_size; i++){
00429 if ( fabs(cpl_vector_get( dispy, i)) > (sigma * sigma_med) ){
00430 new_rejected++;
00431 rejected[idx[i]] = 1;
00432 }
00433 else{
00434 idx[i-new_rejected] = idx[i];
00435 A[i-new_rejected] = A[i];
00436 lambda[i-new_rejected] = lambda[i];
00437 n[i-new_rejected] = n[i];
00438 s[i-new_rejected] = s[i];
00439 }
00440 }
00441 xsh_free_vector(&dispy);
00442 index_size -= new_rejected;
00443 nbrejected += new_rejected;
00444
00445 frac = 1 - ( (float) nbrejected / (float) size);
00446 xsh_msg_dbg_high(" NB_REJECT = %d", nbrejected);
00447 xsh_msg_dbg_high(" GOOD PIXEL FRACTION = %f", frac);
00448
00449 check( xsh_wavesol_compute(wavesol, index_size, A,
00450 &(wavesol->min_y), &(wavesol->max_y), lambda, n, s, fitpoly));
00451 nbiter++;
00452 }
00453
00454 cleanup:
00455 XSH_FREE( idx);
00456 xsh_free_vector(&dispy);
00457 return;
00458 }
00459
00460
00461
00462
00463
00501 void xsh_detect_arclines( cpl_frame *frame,
00502 cpl_frame *theo_tab_frame,
00503 cpl_frame *arc_lines_tab_frame,
00504 cpl_frame* wave_tab_guess_frame,
00505 cpl_frame *order_tab_recov_frame,
00506 cpl_frame *config_model_frame,
00507 cpl_frame *spectralformat_frame,
00508 cpl_frame **resid_tab_orders_frame,
00509 cpl_frame **arc_lines_clean_tab_frame,
00510 cpl_frame **wave_tab_frame,
00511 cpl_frame **resid_tab_frame,
00512 xsh_sol_wavelength solwave_type,
00513 xsh_detect_arclines_param *da,
00514 xsh_clipping_param *dac,
00515 xsh_instrument *instr,
00516 const char* rec_id)
00517 {
00518
00519
00520
00521 xsh_the_map *themap = NULL;
00522 xsh_resid_tab *resid = NULL;
00523 xsh_resid_tab *resid_orders = NULL;
00524 xsh_arclist *arclist = NULL;
00525 xsh_pre *pre = NULL;
00526 cpl_polynomial *fit2dx = NULL;
00527 cpl_propertylist* header = NULL;
00528 int * sort = NULL;
00529 double *vlambdadata = NULL, *vsdata = NULL, *vorderdata = NULL;
00530 int *vsindexdata = NULL;
00531 double *vxthedata = NULL, *vythedata = NULL;
00532 double *corr_x = NULL, *corr_y = NULL;
00533 double *gaussian_pos_x = NULL, *gaussian_pos_y = NULL,
00534 *gaussian_sigma_x = NULL, *gaussian_sigma_y = NULL;
00535 cpl_vector *gaussian_sigma_x_vect =NULL;
00536 cpl_vector *gaussian_sigma_y_vect =NULL;
00537
00538 double *diffx = NULL, *diffy = NULL;
00539 double *diffxmean = NULL, *diffymean = NULL, *diffxsig = NULL, *diffysig = NULL;
00540 xsh_order_list *order_tab_recov = NULL;
00541
00542 int nlinematched, nlinecat_clean=0, nlinecat;
00543 int i, sol_size = 0;
00544 xsh_wavesol* wave_table = NULL;
00545 char wave_table_name[80];
00546 const char* wave_table_tag = NULL;
00547 xsh_wavesol* wave_tab_guess = NULL;
00548 xsh_xs_3 config_model;
00549 xsh_xs_3* p_xs_3;
00550 xsh_spectralformat_list *spectralformat_list = NULL;
00551 int lines_not_in_image = 0;
00552 int lines_not_good_sn = 0;
00553 int lines_not_gauss_fit = 0;
00554 int lines_not_valid_pixels = 0;
00555 int lines_too_few_ph=0;
00556
00557 int* rejected = NULL;
00558 int nb_rejected = 0;
00559 int solution_type = 0;
00560 const char* solution_type_name[2] = { "POLY", "MODEL"};
00561 int detection_mode = 0;
00562 const char* detection_mode_name[3] = { "NORMAL", "CORRECTED", "RECOVER"};
00563 int nb_pinhole;
00564 char fname[256];
00565 char rname[256];
00566 char rtag[80];
00567 const char* tag=NULL;
00568
00569 char dpr_type[80];
00570 char type[80];
00571 cpl_table* wave_trace=NULL;
00572
00573
00574 cpl_table* cfg_tab=NULL;
00575 const char* cfg_name=NULL;
00576 char new_name[256];
00577 char basename[256];
00578
00579
00580 cpl_propertylist* plist=NULL;
00581
00582
00583
00584 int starti=0;
00585 int newi=0;
00586 int j,k;
00587 int min_slit_match=7;
00588 int found_temp=true;
00589
00590
00591 XSH_ASSURE_NOT_NULL( frame);
00592 XSH_ASSURE_NOT_NULL( arc_lines_tab_frame);
00593 XSH_ASSURE_NOT_NULL( spectralformat_frame);
00594 XSH_ASSURE_NOT_NULL( instr);
00595 XSH_ASSURE_NOT_NULL( da);
00596 XSH_ASSURE_NOT_NULL( dac);
00597
00598
00599 XSH_ASSURE_NOT_NULL( arc_lines_clean_tab_frame);
00600 XSH_ASSURE_NOT_NULL( resid_tab_frame);
00601
00602
00603
00604 xsh_msg_dbg_medium("---Detect arclines parameters");
00605 xsh_msg_dbg_medium("min_sn %f fit-window-half-size %d", da->min_sn, da->fit_window_hsize);
00606 xsh_msg_dbg_medium("Clipping sigma %f niter %d frac %f", dac->sigma, dac->niter, dac->frac);
00607
00608 check( pre = xsh_pre_load (frame, instr));
00609 check_msg( nb_pinhole = xsh_pfits_get_nb_pinhole( pre->data_header),
00610 "fail to find number of pinholes info in frame %s",
00611 cpl_frame_get_filename(frame) );
00612 check(strcpy(dpr_type,xsh_pfits_get_dpr_type(pre->data_header)));
00613
00614 if(strstr(dpr_type,"FMTCHK") != NULL) {
00615 strcpy(type,"FMTCHK_");
00616 } else if(strstr(dpr_type,"WAVE") != NULL) {
00617 strcpy(type,"WAVE_");
00618 } else if(strstr(dpr_type,"WAVE") != NULL) {
00619 strcpy(type,"ARC_");
00620 } else {
00621 strcpy(type,"");
00622 }
00623
00624 check( arclist = xsh_arclist_load( arc_lines_tab_frame));
00625 check( spectralformat_list = xsh_spectralformat_list_load(
00626 spectralformat_frame, instr));
00627 check(nlinecat = xsh_arclist_get_size( arclist));
00628
00629 check(xsh_arclist_lambda_sort( arclist));
00630
00631
00632
00633 if ( theo_tab_frame != NULL) {
00634 XSH_ASSURE_NOT_ILLEGAL( config_model_frame == NULL);
00635 solution_type = XSH_DETECT_ARCLINES_TYPE_POLY;
00636 check( themap = xsh_the_map_load( theo_tab_frame));
00637 check( xsh_the_map_lambda_order_slit_sort( themap));
00638 check( theo_tab_filter( themap, arclist, &sol_size, &vlambdadata,
00639 &vorderdata, &vsdata, &vsindexdata, &vxthedata, &vythedata, nb_pinhole));
00640
00641 }
00642 else if ( config_model_frame != NULL){
00643 solution_type = XSH_DETECT_ARCLINES_TYPE_MODEL;
00644 p_xs_3=&config_model;
00645
00646
00647
00648 check(cfg_name=cpl_frame_get_filename(config_model_frame));
00649 check(plist=cpl_propertylist_load(cfg_name,0));
00650 check(cfg_tab=cpl_table_load(cfg_name,1,0));
00651
00652 tag=XSH_GET_TAG_FROM_ARM(XSH_MOD_CFG_TAB,instr);
00653 sprintf(basename,"%s.fits",tag);
00654 sprintf(new_name,"local_%s",basename);
00655 check( xsh_pfits_set_pcatg(plist,tag));
00656 check(cpl_table_save(cfg_tab,plist,NULL,new_name,CPL_IO_DEFAULT));
00657 check(cpl_frame_set_filename(config_model_frame,new_name));
00658 xsh_free_table(&cfg_tab);
00659 xsh_free_propertylist(&plist);
00660
00661 check( xsh_model_config_load_best( config_model_frame, &config_model));
00662 XSH_REGDEBUG("load config model ok");
00663
00664 check(xsh_model_temperature_update_frame(&config_model_frame,frame,
00665 instr,&found_temp));
00666 check(xsh_model_temperature_update_structure(&config_model,frame,instr));
00667 check( theo_tab_model( &config_model, arclist, spectralformat_list,
00668 &sol_size, &vlambdadata, &vorderdata, &vsdata, &vsindexdata,
00669 &vxthedata, &vythedata, instr, nb_pinhole));
00670 }
00671 else{
00672 XSH_ASSURE_NOT_ILLEGAL_MSG(1==0,
00673 "Undefined solution type (POLY or MODEL). See your input sof");
00674 }
00675
00676 xsh_spectralformat_list_free(&spectralformat_list);
00677 xsh_msg_dbg_high("Solution type %s", solution_type_name[solution_type]);
00678
00679
00680 if ( wave_tab_guess_frame == NULL){
00681 if ( order_tab_recov_frame == NULL){
00682 detection_mode = XSH_DETECT_ARCLINES_MODE_NORMAL;
00683 }
00684 else{
00685 detection_mode = XSH_DETECT_ARCLINES_MODE_RECOVER;
00686 check( order_tab_recov = xsh_order_list_load( order_tab_recov_frame,
00687 instr));
00688 }
00689 }
00690 else {
00691 if ( order_tab_recov_frame == NULL){
00692 detection_mode = XSH_DETECT_ARCLINES_MODE_CORRECTED;
00693
00694 check( wave_tab_guess = xsh_wavesol_load( wave_tab_guess_frame,
00695 instr ));
00696 xsh_msg( "BinX,Y: %d, %d", wave_tab_guess->bin_x,
00697 wave_tab_guess->bin_y ) ;
00698 }
00699 else{
00700 detection_mode = -1;
00701 }
00702 }
00703
00704 XSH_ASSURE_NOT_ILLEGAL( detection_mode >=XSH_DETECT_ARCLINES_MODE_NORMAL
00705 && detection_mode <= XSH_DETECT_ARCLINES_MODE_RECOVER);
00706
00707 xsh_msg( "Detection mode : %s", detection_mode_name[detection_mode]);
00708
00709 xsh_msg( "Solution size %d", sol_size);
00710
00711
00712 XSH_MALLOC( corr_x, double, sol_size);
00713 XSH_MALLOC( corr_y, double, sol_size);
00714
00715 XSH_MALLOC( gaussian_pos_x, double, sol_size);
00716 XSH_MALLOC( gaussian_pos_y, double, sol_size);
00717 XSH_MALLOC( gaussian_sigma_x, double, sol_size);
00718 XSH_MALLOC( gaussian_sigma_y, double, sol_size);
00719
00720 XSH_MALLOC( diffx, double, sol_size);
00721 XSH_MALLOC( diffy, double, sol_size);
00722
00723 XSH_MALLOC( diffxmean, double, sol_size);
00724 XSH_MALLOC( diffymean, double, sol_size);
00725
00726 XSH_MALLOC( diffxsig, double, sol_size);
00727 XSH_MALLOC( diffysig, double, sol_size);
00728
00729
00730 nlinematched = 0;
00731
00732 if ( da->find_center_method == XSH_GAUSSIAN_METHOD){
00733 xsh_msg_dbg_medium( "USE method GAUSSIAN");
00734 }
00735 else{
00736 xsh_msg_dbg_medium( "USE method BARYCENTER");
00737 }
00738
00739 if (xsh_debug_level_get() >= XSH_DEBUG_LEVEL_HIGH) {
00740 FILE* regfile = NULL;
00741
00742 regfile = fopen( "FIT.reg", "w");
00743 fprintf( regfile, "# Region file format: DS9 version 4.0\n"\
00744 "global color=red font=\"helvetica 4 normal\""\
00745 "select=1 highlite=1 edit=1 move=1 delete=1 include=1 fixed=0 "\
00746 "source\nimage\n");
00747 fclose( regfile);
00748 regfile = fopen( "NOFIT.reg", "w");
00749 fprintf( regfile, "# Region file format: DS9 version 4.0\n"\
00750 "global color=red font=\"helvetica 4 normal\""\
00751 "select=1 highlite=1 edit=1 move=1 delete=1 include=1 fixed=0 "\
00752 "source\nimage\n");
00753 fclose( regfile);
00754 }
00755
00756 for(i=0; i< sol_size; i++){
00757 int xpos = 0, ypos = 0;
00758 int slit_index;
00759 double corrxv = 0.0, corryv = 0.0;
00760 double lambda, order, slit, xthe, ythe;
00761 double xcor, ycor, x, y, sig_x, sig_y;
00762
00763
00764 lambda = vlambdadata[i];
00765 order = vorderdata[i];
00766 slit = vsdata[i];
00767 slit_index = vsindexdata[i];
00768 xthe = vxthedata[i];
00769 ythe = vythedata[i];
00770 xsh_msg_dbg_high( "THE LAMBDA %f ORDER %f SLIT %f X %f Y %f", lambda,
00771 order, slit, xthe, ythe);
00772
00773
00774 xcor = xthe;
00775 ycor = ythe;
00776 xsh_msg_dbg_high("THE x%f y %f", xthe, ythe);
00777
00778
00779 if ( order_tab_recov != NULL){
00780 int iorder = 0;
00781 cpl_polynomial *center_poly_recov = NULL;
00782
00783 check (iorder = xsh_order_list_get_index_by_absorder(
00784 order_tab_recov, order));
00785 center_poly_recov = order_tab_recov->list[iorder].cenpoly;
00786 check( xcor = cpl_polynomial_eval_1d( center_poly_recov, ycor, NULL));
00787 corrxv = xcor-xthe;
00788 }
00789
00790 if ( wave_tab_guess != NULL){
00791 double diffxv, diffyv;
00792
00793 check( diffxv = xsh_wavesol_eval_polx( wave_tab_guess, lambda,
00794 order, 0));
00795 check( diffyv = xsh_wavesol_eval_poly( wave_tab_guess, lambda,
00796 order, 0));
00797
00798 xcor = xthe+diffxv;
00799 ycor = ythe+diffyv;
00800 corrxv = diffxv;
00801 corryv = diffyv;
00802 }
00803
00804 xsh_msg_dbg_high("x %f y %f CORR_x %f CORR_y %f", xcor, ycor,
00805 corrxv, corryv);
00806
00807 if ( xcor >=0.5 && ycor >=0.5 &&
00808 xcor < (pre->nx+0.5) && ycor < (pre->ny+0.5)){
00809 int best_med_res = 0;
00810 check ( best_med_res = xsh_pre_window_best_median_flux_pos( pre,
00811 (int)(xcor+0.5)-1, (int)(ycor+0.5)-1,
00812 da->search_window_hsize, da->running_median_hsize, &xpos, &ypos));
00813
00814 if (best_med_res != 0){
00815 lines_not_valid_pixels++;
00816 continue;
00817 }
00818
00819 xpos = xpos+1;
00820 ypos = ypos+1;
00821 xsh_msg_dbg_high("ADJ x %d y %d", xpos, ypos);
00822
00823 if ( da->find_center_method == XSH_GAUSSIAN_METHOD){
00824 cpl_image_fit_gaussian(pre->data, xpos, ypos, 1+2*da->fit_window_hsize,
00825 NULL, &x, &y, &sig_x, &sig_y, NULL, NULL);
00826 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);
00827 }
00828 else{
00829 xsh_image_find_barycenter( pre->data, xpos, ypos,
00830 1+2*da->fit_window_hsize,
00831 NULL, &x, &y, &sig_x, &sig_y, NULL, NULL);
00832 xsh_msg_dbg_high("BARY x %f y %f", x, y);
00833 }
00834
00835
00836
00837
00838
00839
00840
00841
00842
00843
00844
00845
00846
00847
00848 if (nb_pinhole==1 && solution_type==XSH_DETECT_ARCLINES_TYPE_MODEL) {
00849 if (p_xs_3->arm==2) {
00850 x=x+0.125;
00851 }
00852 else if (p_xs_3->arm==0) {
00853 x=x-0.51;
00854 }
00855 }
00856
00857 if (xsh_debug_level_get() >= XSH_DEBUG_LEVEL_HIGH) {
00858 const char* fit = NULL;
00859 const char* color = NULL;
00860 FILE* regfile = NULL;
00861
00862 if (cpl_error_get_code() == CPL_ERROR_NONE){
00863 fit = "FIT.reg";
00864 color = "green";
00865 }
00866 else{
00867 fit = "NOFIT.reg";
00868 color = "red";
00869 }
00870 regfile = fopen( fit, "a");
00871 fprintf( regfile, "point(%d,%d) #point=cross color=blue\n", (int)(xcor), (int)(ycor));
00872 fprintf( regfile, "box(%d,%d,%d,%d) #color=blue\n", (int)(xcor), (int)(ycor),
00873 da->search_window_hsize*2,da->search_window_hsize*2);
00874 fprintf( regfile, "point(%d,%d) #point=cross color=%s font=\"helvetica 10 normal\""\
00875 " text={%.3f}\n", xpos, ypos, color, lambda);
00876 fprintf( regfile, "box(%d,%d,%d,%d) #color=%s\n", (int)xpos, (int)ypos,
00877 1+2*da->fit_window_hsize, 1+2*da->fit_window_hsize, color);
00878 fclose( regfile);
00879 }
00880
00881 if( cpl_error_get_code() == CPL_ERROR_NONE ){
00882 if ( lines_filter_by_sn( pre, da->min_sn, x, y) ){
00883 vlambdadata[nlinematched] = lambda;
00884 vorderdata[nlinematched] = order;
00885 vsdata[nlinematched] = slit;
00886 vsindexdata[nlinematched] = slit_index;
00887 vxthedata[nlinematched] = xthe;
00888 vythedata[nlinematched] = ythe;
00889 corr_x[nlinematched] = corrxv;
00890 corr_y[nlinematched] = corryv;
00891 gaussian_pos_x[nlinematched] = x;
00892 gaussian_pos_y[nlinematched] = y;
00893 gaussian_sigma_x[nlinematched] = sig_x;
00894 gaussian_sigma_y[nlinematched] = sig_y;
00895 diffx[nlinematched] = x-xthe;
00896 diffy[nlinematched] = y-ythe;
00897
00898
00899
00900
00901
00902
00903 nlinematched++;
00904
00905 }
00906 else{
00907 lines_not_good_sn++;
00908 xsh_msg_dbg_medium("Not good s/n for this line");
00909 if (xsh_debug_level_get() >= XSH_DEBUG_LEVEL_HIGH) {
00910
00911 int xpix, ypix, rej;
00912 float flux, noise, sn;
00913 FILE* regfile = NULL;
00914 char sn_name[256];
00915
00916 xpix =(int) rint(x);
00917 ypix = (int)rint(y);
00918
00919 check( flux = cpl_image_get( pre->data, xpix, ypix, &rej));
00920 check( noise = cpl_image_get( pre->errs, xpix, ypix, &rej));
00921 sn = flux / noise;
00922 sprintf( sn_name, "bad_sn_%.3f.reg", da->min_sn);
00923 regfile = fopen( sn_name, "a");
00924 fprintf( regfile, "point(%d,%d) #point=cross color=red font=\"helvetica 10 normal\""\
00925 " text={%.3f [%.3f}\n", xpix, ypix, lambda, sn);
00926 fclose( regfile);
00927 }
00928 }
00929 }
00930 else{
00931 lines_not_gauss_fit++;
00932 xsh_msg_dbg_medium("No Fit Gaussian for this line");
00933 xsh_error_reset();
00934 }
00935 }
00936 else{
00937 lines_not_in_image++;
00938 xsh_msg_dbg_medium("Coordinates are not in the image");
00939 }
00940 }
00941
00942
00943
00944
00945
00946 if (nb_pinhole>1) {
00947 for (i=1; i<nlinematched; i++) {
00948 xsh_msg_dbg_high("filter 7 pinhole : line %d : lambda %f order %f slit %f",
00949 i, vlambdadata[i], vorderdata[i], vsdata[i]);
00950
00951 if (vlambdadata[i]!=vlambdadata[i-1] || vorderdata[i]!=vorderdata[i-1] || i==nlinematched-1) {
00952
00953
00954 if (i==nlinematched-1) {
00955 i++;
00956 }
00957 xsh_msg_dbg_high("filter 7 pinhole : from %d to %d find %d pinholes",
00958 starti, i, i-starti);
00959
00960 if (i-starti>=min_slit_match) {
00961 for (k=starti; k<=i-1; k++) {
00962 diffxmean[k]=0.0;
00963 diffymean[k]=0.0;
00964 for (j=starti; j<=i-1; j++) {
00965 if (j!=k) {
00966 diffxmean[k]+=diffx[j];
00967 diffymean[k]+=diffy[j];
00968 }
00969 }
00970 diffxmean[k]/=(float)(i-starti-1);
00971 diffymean[k]/=(float)(i-starti-1);
00972 }
00973 for (k=starti; k<=i-1; k++) {
00974 diffxsig[k]=0.0;
00975 diffysig[k]=0.0;
00976 for (j=starti; j<=i-1; j++) {
00977 if (j!=k) {
00978 diffysig[k]+=(diffy[j]-diffymean[k])*(diffy[j]-diffymean[k]);
00979 diffxsig[k]+=(diffx[j]-diffxmean[k])*(diffx[j]-diffxmean[k]);
00980 }
00981 }
00982 diffxsig[k]=sqrt(diffxsig[k]/(float)(i-starti-1));
00983 diffysig[k]=sqrt(diffysig[k]/(float)(i-starti-1));
00984 }
00985 for (j=starti; j<=i-1; j++) {
00986
00987
00988 vlambdadata[newi] = vlambdadata[j];
00989 vorderdata[newi] = vorderdata[j];
00990 vsdata[newi] =vsdata[j];
00991 vsindexdata[newi] = vsindexdata[j];
00992 vxthedata[newi] =vxthedata[j];
00993 vythedata[newi] =vythedata[j];
00994 corr_x[newi] = corr_x[j];
00995 corr_y[newi] = corr_y[j];
00996 gaussian_pos_x[newi] = gaussian_pos_x[j];
00997 gaussian_pos_y[newi] = gaussian_pos_y[j];
00998 gaussian_sigma_x[newi] = gaussian_sigma_x[j];
00999 gaussian_sigma_y[newi] = gaussian_sigma_y[j];
01000 diffx[newi] = diffx[j];
01001 diffy[newi] = diffy[j];
01002 newi++;
01003
01004
01005
01006
01007
01008 }
01009 }
01010 else {
01011 lines_too_few_ph+=i-starti;
01012 }
01013 starti=i;
01014 }
01015 }
01016 nlinematched=newi;
01017 }
01018
01019 xsh_msg("nlinematched / sol_size = %d / %d", nlinematched, sol_size);
01020 xsh_msg_dbg_medium(" %d lines not found in image", lines_not_in_image);
01021 xsh_msg_dbg_medium(" %d lines not good S/N", lines_not_good_sn);
01022 xsh_msg_dbg_medium(" %d lines no fit gaussian", lines_not_gauss_fit);
01023 xsh_msg_dbg_medium(" %d lines no valid pixels", lines_not_valid_pixels);
01024 xsh_msg_dbg_medium(" %d lines detected in less than %d/9 ph positions", lines_too_few_ph,min_slit_match);
01025 check( gaussian_sigma_x_vect = cpl_vector_wrap( nlinematched, gaussian_sigma_x));
01026 check( gaussian_sigma_y_vect = cpl_vector_wrap( nlinematched, gaussian_sigma_y));
01027 xsh_msg("sigma gaussian median in x %lg",
01028 cpl_vector_get_median_const( gaussian_sigma_x_vect));
01029 xsh_msg("sigma gaussian median in y %lg",
01030 cpl_vector_get_median_const( gaussian_sigma_y_vect));
01031 xsh_unwrap_vector( &gaussian_sigma_x_vect);
01032 xsh_unwrap_vector( &gaussian_sigma_y_vect);
01033
01034
01035
01036 if ( resid_tab_orders_frame != NULL){
01037 check( resid_orders = xsh_resid_tab_create( nlinematched, vlambdadata,
01038 vorderdata, vsdata, vsindexdata, vxthedata, vythedata, corr_x, corr_y,
01039 gaussian_pos_x, gaussian_pos_y, gaussian_sigma_x, gaussian_sigma_y,
01040 wave_table, solution_type));
01041
01042 sprintf(rtag,"%s%s%s",type,"RESID_TAB_ORDERS_",
01043 xsh_instrument_arm_tostring( instr ));
01044
01045 sprintf(rname,"%s%s",rtag,".fits");
01046
01047 check( *resid_tab_orders_frame = xsh_resid_tab_save( resid_orders,
01048 rname, instr,rtag));
01049 }
01050 xsh_msg_dbg_high("solution_type=%d poly=%d model=%d",solution_type,
01051 XSH_DETECT_ARCLINES_TYPE_POLY, XSH_DETECT_ARCLINES_TYPE_MODEL);
01052
01053
01054
01055
01056 if ( solution_type == XSH_DETECT_ARCLINES_TYPE_POLY &&
01057 (wave_tab_frame != NULL)){
01058 check( wave_table = xsh_wavesol_create( spectralformat_frame, da, instr));
01059 XSH_CALLOC( rejected, int, nlinematched);
01060 if ( solwave_type == XSH_SOLUTION_RELATIVE) {
01061
01062 check( xsh_wavesol_set_type( wave_table, XSH_WAVESOL_GUESS));
01063 check( data_wavesol_fit_with_sigma( wave_table, diffy, vlambdadata,
01064 vorderdata, vsdata, nlinematched, dac->niter, dac->frac,
01065 dac->sigma, rejected));
01066 nb_rejected = 0;
01067 for(i=0; i< nlinematched; i++){
01068 if (rejected[i] == 1){
01069 nb_rejected++;
01070 }
01071 else{
01072
01073 vxthedata[i-nb_rejected] = vxthedata[i];
01074 vythedata[i-nb_rejected] = vythedata[i];
01075 vsindexdata[i-nb_rejected] = vsindexdata[i];
01076 corr_x[i-nb_rejected] = corr_x[i];
01077 corr_y[i-nb_rejected] = corr_y[i];
01078 gaussian_pos_x[i-nb_rejected] = gaussian_pos_x[i];
01079 gaussian_pos_y[i-nb_rejected] = gaussian_pos_y[i];
01080 gaussian_sigma_x[i-nb_rejected] = gaussian_sigma_x[i];
01081 gaussian_sigma_y[i-nb_rejected] = gaussian_sigma_y[i];
01082 diffx[i-nb_rejected] = diffx[i];
01083 }
01084 }
01085 nlinematched = nlinematched-nb_rejected;
01086
01087 check( fit2dx = xsh_wavesol_get_polx( wave_table));
01088 check( xsh_wavesol_compute( wave_table,
01089 nlinematched, diffx,
01090 &(wave_table->min_x), &(wave_table->max_x),
01091 vlambdadata, vorderdata, vsdata, fit2dx));
01092
01093
01094 sprintf(wave_table_name,"%s%s.fits","WAVE_TAB_GUESS_",
01095 xsh_instrument_arm_tostring( instr )) ;
01096
01097 wave_table_tag = XSH_GET_TAG_FROM_ARM( XSH_WAVE_TAB_GUESS, instr);
01098
01099 }
01100 else{
01101 check( xsh_wavesol_set_type( wave_table, XSH_WAVESOL_2D));
01102
01103 check( data_wavesol_fit_with_sigma( wave_table, gaussian_pos_y,
01104 vlambdadata, vorderdata, vsdata, nlinematched, dac->niter, dac->frac,
01105 dac->sigma, rejected));
01106 nb_rejected = 0;
01107 for(i=0; i< nlinematched; i++){
01108 if (rejected[i] == 1){
01109 nb_rejected++;
01110 }
01111 else{
01112
01113 vxthedata[i-nb_rejected] = vxthedata[i];
01114 vythedata[i-nb_rejected] = vythedata[i];
01115 vsindexdata[i-nb_rejected] = vsindexdata[i];
01116 corr_x[i-nb_rejected] = corr_x[i];
01117 corr_y[i-nb_rejected] = corr_y[i];
01118 gaussian_pos_x[i-nb_rejected] = gaussian_pos_x[i];
01119 gaussian_sigma_x[i-nb_rejected] = gaussian_sigma_x[i];
01120 gaussian_sigma_y[i-nb_rejected] = gaussian_sigma_y[i];
01121 }
01122 }
01123 nlinematched = nlinematched-nb_rejected;
01124
01125 check(fit2dx = xsh_wavesol_get_polx(wave_table));
01126 check(xsh_wavesol_compute(wave_table, nlinematched, gaussian_pos_x,
01127 &wave_table->min_x, &wave_table->max_x, vlambdadata, vorderdata,
01128 vsdata, fit2dx));
01129
01130 sprintf(wave_table_name,"%s%s.fits","WAVE_TAB_2D_",
01131 xsh_instrument_arm_tostring( instr )) ;
01132
01133
01134 wave_table_tag = XSH_GET_TAG_FROM_ARM( XSH_WAVE_TAB_2D, instr);
01135 }
01136
01137
01138 check(header = xsh_wavesol_get_header(wave_table));
01139 check(xsh_pfits_set_qc_nlinecat(header,nlinecat));
01140 check(xsh_pfits_set_qc_nlinefound(header,sol_size));
01141
01142
01143
01144 check( wave_trace=xsh_wavesol_trace(wave_table,vlambdadata,vorderdata,
01145 vsdata,nlinematched));
01146
01147 check( *wave_tab_frame=xsh_wavesol_save(wave_table,wave_trace,
01148 wave_table_name,wave_table_tag));
01149 xsh_add_temporary_file( wave_table_name) ;
01150 check( cpl_frame_set_tag( *wave_tab_frame, wave_table_tag));
01151 }
01152
01153
01154 check( xsh_arclist_clean_from_list( arclist, vlambdadata, nlinematched));
01155 nlinecat_clean = arclist->size-arclist->nbrejected;
01156 check( header = xsh_arclist_get_header( arclist));
01157 check( xsh_pfits_set_qc_nlinecat( header,nlinecat));
01158 check( xsh_pfits_set_qc_nlinefound( header, sol_size));
01159 check( xsh_pfits_set_qc_nlinecat_clean( header,nlinecat_clean));
01160 check( xsh_pfits_set_qc_nlinefound_clean( header,nlinematched));
01161 if(strcmp(rec_id,"xsh_predict") == 0) {
01162 tag=XSH_GET_TAG_FROM_ARM( XSH_ARC_LINE_LIST_PREDICT, instr);
01163 } else {
01164 tag=XSH_GET_TAG_FROM_ARM( XSH_ARC_LINE_LIST_2DMAP, instr);
01165 }
01166 sprintf(fname,"%s%s",tag,".fits");
01167 check( *arc_lines_clean_tab_frame = xsh_arclist_save( arclist,fname,tag));
01168
01169
01170
01171 check( resid = xsh_resid_tab_create( nlinematched, vlambdadata, vorderdata,
01172 vsdata, vsindexdata, vxthedata, vythedata,
01173 corr_x, corr_y, gaussian_pos_x,
01174 gaussian_pos_y, gaussian_sigma_x, gaussian_sigma_y, wave_table,
01175 solution_type));
01176
01177 sprintf(rtag,"%s%s%s",type,"RESID_TAB_LINES_",
01178 xsh_instrument_arm_tostring( instr ));
01179 check( tag = cpl_frame_get_tag( frame));
01180
01181 XSH_REGDEBUG("TAG %s", tag);
01182
01183 if ( strstr( tag, XSH_AFC_CAL) ){
01184 sprintf(rname,"AFC_CAL_%s%s",rtag,".fits") ;
01185 }
01186 else if ( strstr( tag, XSH_AFC_ATT) ){
01187 sprintf(rname,"AFC_ATT_%s%s",rtag,".fits") ;
01188 }
01189 else {
01190 sprintf(rname,"%s%s",rtag,".fits") ;
01191 }
01192 check( *resid_tab_frame = xsh_resid_tab_save( resid,rname,instr,rtag));
01193
01194
01195 cleanup:
01196 XSH_FREE( vlambdadata);
01197 XSH_FREE( vorderdata);
01198 XSH_FREE( vsdata);
01199 XSH_FREE( vsindexdata);
01200 XSH_FREE( vxthedata);
01201 XSH_FREE( vythedata);
01202 XSH_FREE( corr_x);
01203 XSH_FREE( corr_y);
01204 XSH_FREE( gaussian_pos_x);
01205 XSH_FREE( gaussian_pos_y);
01206 XSH_FREE( gaussian_sigma_x);
01207 XSH_FREE( gaussian_sigma_y);
01208 XSH_FREE( sort);
01209 XSH_FREE( diffx);
01210 XSH_FREE( diffy);
01211 XSH_FREE( diffxmean);
01212 XSH_FREE( diffymean);
01213 XSH_FREE( diffxsig);
01214 XSH_FREE( diffysig);
01215 XSH_FREE( rejected);
01216
01217 xsh_free_table( &wave_trace);
01218 xsh_spectralformat_list_free(&spectralformat_list);
01219 xsh_pre_free( &pre);
01220 xsh_the_map_free( &themap);
01221 xsh_wavesol_free( &wave_table);
01222 xsh_wavesol_free( &wave_tab_guess);
01223 xsh_order_list_free( &order_tab_recov);
01224 xsh_resid_tab_free( &resid_orders);
01225 xsh_resid_tab_free( &resid);
01226 xsh_arclist_free( &arclist);
01227 return;
01228 }
01229
01230
01231
01232