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 #define CPL_MODE 0
00031
00038
00041
00042
00043
00044
00045 #include <math.h>
00046 #include <xsh_drl.h>
00047
00048 #include <xsh_utils_table.h>
00049 #include <xsh_badpixelmap.h>
00050 #include <xsh_utils_wrappers.h>
00051 #include <xsh_data_pre.h>
00052 #include <xsh_dfs.h>
00053 #include <xsh_pfits.h>
00054 #include <xsh_error.h>
00055 #include <xsh_msg.h>
00056 #include <xsh_fit.h>
00057 #include <xsh_data_instrument.h>
00058 #include <xsh_data_localization.h>
00059 #include <xsh_data_spectrum.h>
00060 #include <xsh_data_rec.h>
00061 #include <xsh_utils_ifu.h>
00062 #include <xsh_ifu_defs.h>
00063 #include <xsh_parameters.h>
00064 #include <cpl.h>
00065
00066
00067
00068
00069 static void chunk_coadd( double* data, double* flux, int* qual, int nx, int ny,
00070 int ifirst, int ilast, const int decode_bp, int *skymask);
00071
00072 static void tab_minmax_get_index( double* data, int size,
00073 xsh_slit_limit_param * slit_limit_par,
00074 int* min, int* max);
00075
00076 static xsh_localization* xsh_localize_obj_auto( cpl_frame * rec_frame,
00077 cpl_frame *skymask_frame,
00078 xsh_instrument* instrument, xsh_localize_obj_param * loc_obj_par,
00079 xsh_slit_limit_param * slit_limit_par,const int decode_bp);
00080
00081 static xsh_localization* xsh_localize_obj_manual(
00082 xsh_instrument *instrument, xsh_localize_obj_param *loc_obj_par);
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106 static void chunk_coadd( double* data, double *flux, int* qual, int nx, int ny,
00107 int ifirst, int ilast, const int decode_bp,int *skymask)
00108 {
00109 int i,j;
00110 int* nbgoodpixels = NULL;
00111
00112 XSH_ASSURE_NOT_NULL(data);
00113 XSH_ASSURE_NOT_NULL(flux);
00114 XSH_ASSURE_NOT_NULL(qual);
00115 XSH_ASSURE_NOT_NULL(skymask);
00116 XSH_ASSURE_NOT_ILLEGAL(ifirst >= 0 && ifirst <= ilast && ilast < nx);
00117
00118 XSH_CALLOC( nbgoodpixels, int, ny);
00119
00120 for(j=0; j<ny; j++){
00121 data[j] = 0;
00122 }
00123
00124
00125
00126
00127 int j_nx=0;
00128 for(i=ifirst; i<=ilast; i++){
00129 if (skymask[i] == 0){
00130 for( j=0 ; j<ny ; j++) {
00131 j_nx=i+j*nx;
00132 if ( (qual[j_nx] & decode_bp) == 0 ){
00133 nbgoodpixels[j]++;
00134 data[j] +=flux[j_nx];
00135 }
00136
00137
00138
00139
00140
00141 }
00142 }
00143 else{
00144 xsh_msg_dbg_medium("Mask by sky lines : column %d", i);
00145 }
00146 }
00147
00148 for(j=0; j<ny; j++){
00149 if ( nbgoodpixels[j] != 0 ) {
00150 data[j] = data[j] / (float)(nbgoodpixels[j]);
00151 }
00152 else data[j] = 0. ;
00153 }
00154
00155 cleanup:
00156 XSH_FREE(nbgoodpixels);
00157 return;
00158 }
00159
00160 static void tab_minmax_get_index( double* data, int size,
00161 xsh_slit_limit_param * slit_limit_par,
00162 int* min, int* max)
00163 {
00164 int i, maxi, mini ;
00165 int slit0, slit1 ;
00166
00167 XSH_ASSURE_NOT_NULL(data);
00168 XSH_ASSURE_NOT_NULL(min);
00169 XSH_ASSURE_NOT_NULL(max);
00170
00171 if ( slit_limit_par == NULL ) {
00172 slit0 = 0 ;
00173 slit1 = size ;
00174 }
00175 else {
00176 slit0 = slit_limit_par->min_slit_idx + 1 ;
00177 slit1 = slit_limit_par->max_slit_idx ;
00178 }
00179 maxi = 0 ;
00180 mini = 0 ;
00181
00182 for(i = slit0; i<slit1; i++) {
00183 if (data[i]>data[maxi]){
00184 maxi = i;
00185 }
00186 else if (data[i]<data[mini]){
00187 mini = i;
00188 }
00189 }
00190 *min = mini;
00191 *max = maxi;
00192
00193 cleanup:
00194 return;
00195 }
00196
00197
00198
00199
00200
00201
00219
00220 cpl_frame* xsh_localize_obj (cpl_frame *rec_frame,
00221 cpl_frame *skymask_frame,
00222 xsh_instrument* instrument,
00223 xsh_localize_obj_param * loc_obj_par,
00224 xsh_slit_limit_param * slitlimit_par,
00225 const char * fname)
00226 {
00227 cpl_frame *merge_frame = NULL;
00228 cpl_frame *res_frame = NULL;
00229 xsh_localization *loc_list = NULL;
00230 int merge_par = 0;
00231 const char* rec_prefix = "LOCALIZE";
00232 char filename[256];
00233
00234
00235
00236 XSH_ASSURE_NOT_NULL( loc_obj_par);
00237 XSH_ASSURE_NOT_NULL( instrument);
00238
00239 xsh_msg_dbg_medium( "Entering xsh_localize_obj") ;
00240
00241 if ( rec_frame != NULL){
00242 const char* filename = NULL;
00243 check( merge_frame = xsh_merge_ord( rec_frame, instrument, merge_par,
00244 rec_prefix));
00245 check( filename = cpl_frame_get_filename( merge_frame));
00246 check( cpl_frame_set_level( merge_frame, CPL_FRAME_LEVEL_TEMPORARY));
00247 xsh_add_temporary_file( filename);
00248 }
00249
00250
00251 xsh_msg_dbg_medium("method %s",LOCALIZE_METHOD_PRINT( loc_obj_par->method));
00252
00253 if ( loc_obj_par->method == LOC_MANUAL_METHOD){
00254 xsh_msg_dbg_medium("slit position %f slit half height %f",
00255 loc_obj_par->slit_position, loc_obj_par->slit_hheight);
00256 check( loc_list = xsh_localize_obj_manual( instrument, loc_obj_par));
00257 }
00258 else{
00259 check( loc_list = xsh_localize_obj_auto( merge_frame, skymask_frame,
00260 instrument,
00261 loc_obj_par, slitlimit_par,instrument->decode_bp));
00262 }
00263
00264 xsh_msg_dbg_low( "Saving Localization Table" ) ;
00265 if ( fname == NULL ) {
00266
00267
00268
00269
00270 sprintf(filename,"LOCALIZATION_TABLE_%s.fits",
00271 xsh_instrument_arm_tostring(instrument));
00272 xsh_add_temporary_file(filename);
00273 } else {
00274 sprintf(filename,fname);
00275 }
00276 check( res_frame = xsh_localization_save( loc_list, filename, instrument));
00277
00278 cleanup:
00279 xsh_localization_free( &loc_list);
00280 xsh_free_frame( &merge_frame);
00281
00282 return res_frame ;
00283 }
00284
00285
00286 static xsh_localization* xsh_localize_obj_manual(
00287 xsh_instrument *instrument, xsh_localize_obj_param *loc_obj_par)
00288 {
00289 xsh_localization *loc_list = NULL;
00290 int deg_poly = 0;
00291 double slit_cen, slit_up, slit_low;
00292 cpl_size pows = 0;
00293
00294
00295 XSH_ASSURE_NOT_NULL( loc_obj_par);
00296 XSH_ASSURE_NOT_NULL( instrument);
00297
00298
00299 xsh_msg_dbg_medium("slit position %f slit half height %f",
00300 loc_obj_par->slit_position, loc_obj_par->slit_hheight);
00301
00302 check( loc_list = xsh_localization_create());
00303
00304 slit_cen = loc_obj_par->slit_position;
00305 slit_up = slit_cen + loc_obj_par->slit_hheight;
00306 slit_low = slit_cen - loc_obj_par->slit_hheight;
00307
00308 loc_list->pol_degree = deg_poly;
00309 check( loc_list->cenpoly = cpl_polynomial_new( 1));
00310 check( loc_list->edglopoly = cpl_polynomial_new( 1));
00311 check( loc_list->edguppoly = cpl_polynomial_new( 1));
00312 check( cpl_polynomial_set_coeff( loc_list->cenpoly,
00313 &pows, slit_cen));
00314 check( cpl_polynomial_set_coeff( loc_list->edglopoly,
00315 &pows, slit_low));
00316 check( cpl_polynomial_set_coeff( loc_list->edguppoly,
00317 &pows, slit_up));
00318
00319 cleanup:
00320 return loc_list;
00321 }
00322
00323
00341
00342 static xsh_localization*
00343 xsh_localize_obj_auto( cpl_frame *merge_frame,
00344 cpl_frame *skymask_frame,
00345 xsh_instrument *instrument,
00346 xsh_localize_obj_param *loc_obj_par,
00347 xsh_slit_limit_param *slitlimit_par,const int decode_bp)
00348 {
00349 xsh_localization *loc_list = NULL;
00350 xsh_spectrum *spectrum = NULL;
00351 double* slit_center = NULL;
00352 double* slit_upper = NULL;
00353 double* slit_lower = NULL;
00354 double* chunk_center = NULL;
00355 cpl_vector* slit_center_v = NULL;
00356 cpl_vector* slit_upper_v = NULL;
00357 cpl_vector* slit_lower_v = NULL;
00358 cpl_vector* chunk_center_v = NULL;
00359
00360 #if 0
00361 int i = 0;
00362 int skip_order_nb = 0;
00363 int first_order_index = 0;
00364 int slit_nod, slit_lim ;
00365 #endif
00366
00367 int level;
00368
00369 int nslit, nlambda, chunk_size;
00370 double* slit_tab = NULL;
00371 double *flux = NULL;
00372
00373 int *qual = NULL;
00374 int ichunk, skip_chunk=0, nb_chunk;
00375 cpl_polynomial *cen_poly = NULL;
00376 int loc_degree;
00377 double lambda_min, lambda_step;
00378 double slit_min, slit_step, kappa=2.0;
00379 int islit, iter, niter=5, nbrej=-1;
00380 cpl_vector *dispslit = NULL;
00381 cpl_vector *gausspos = NULL;
00382 cpl_vector *gaussval = NULL;
00383 int *sky_mask = NULL;
00384 cpl_table *skymask_table = NULL;
00385 const char* skymask_name = NULL;
00386
00387
00388 XSH_ASSURE_NOT_NULL( loc_obj_par);
00389 XSH_ASSURE_NOT_NULL( instrument);
00390 XSH_ASSURE_NOT_NULL( merge_frame);
00391
00392
00393 kappa = loc_obj_par->kappa;
00394 niter = loc_obj_par->niter;
00395 loc_degree = loc_obj_par->loc_deg_poly;
00396 xsh_msg_dbg_medium( "Entering xsh_localize_obj_auto") ;
00397 xsh_msg_dbg_medium("Localize deg_poly %d chunk %d Thresh %f kappa %f niter %d",
00398 loc_degree, loc_obj_par->loc_chunk_nb, loc_obj_par->loc_thresh,
00399 kappa, niter);
00400 if ( loc_obj_par->use_skymask){
00401 XSH_ASSURE_NOT_NULL( skymask_frame);
00402 check( skymask_name = cpl_frame_get_filename( skymask_frame));
00403 xsh_msg_dbg_medium("Sky mask %s", skymask_name);
00404 }
00405
00406 level = xsh_debug_level_get();
00407
00408
00409 check( spectrum = xsh_spectrum_load( merge_frame));
00410 check( nslit = xsh_spectrum_get_size_slit( spectrum));
00411 check( nlambda = xsh_spectrum_get_size_lambda( spectrum));
00412 lambda_min = spectrum->lambda_min;
00413 lambda_step = spectrum->lambda_step;
00414 slit_min = spectrum->slit_min;
00415 slit_step = spectrum->slit_step;
00416
00417 check( flux = cpl_image_get_data_double( spectrum->flux));
00418
00419 check( qual = cpl_image_get_data_int( spectrum->qual));
00420
00421
00422 XSH_CALLOC( chunk_center, double, loc_obj_par->loc_chunk_nb);
00423 XSH_CALLOC( slit_center, double, loc_obj_par->loc_chunk_nb);
00424 XSH_CALLOC( slit_upper, double, loc_obj_par->loc_chunk_nb);
00425 XSH_CALLOC( slit_lower, double, loc_obj_par->loc_chunk_nb);
00426
00427
00428 XSH_CALLOC( sky_mask, int, nlambda);
00429
00430 if ( loc_obj_par->use_skymask){
00431 float *skymask_data = NULL;
00432 int irow, nrow;
00433 double fwhm =0.0, sky_min, sky_max;
00434 int isky_min, isky_max, imask;
00435 double width, resolution;
00436
00437 XSH_TABLE_LOAD( skymask_table, skymask_name);
00438
00439 check( xsh_sort_table_1( skymask_table, "WAVELENGTH", CPL_FALSE));
00440 check( skymask_data = cpl_table_get_data_float( skymask_table,
00441 "WAVELENGTH"));
00442 check( nrow = cpl_table_get_nrow( skymask_table));
00443
00444 xsh_msg_dbg_low("lambda min %f, step %f", lambda_min, lambda_step);
00445
00446 for( irow=0; irow < nrow; irow++){
00447 check( width = xsh_pfits_get_slit_width( spectrum->flux_header, instrument));
00448 resolution = xsh_resolution_get( instrument, width);
00449 fwhm = skymask_data[irow]/resolution;
00450 sky_min = skymask_data[irow]-fwhm;
00451 sky_max = skymask_data[irow]+fwhm;
00452 isky_min = (int)xsh_round_double((sky_min-lambda_min)/lambda_step);
00453 isky_max = (int)xsh_round_double((sky_max-lambda_min)/lambda_step);
00454 for( imask=isky_min; imask <=isky_max; imask++){
00455 sky_mask[imask] = 1;
00456 }
00457 }
00458 if (level >= XSH_DEBUG_LEVEL_MEDIUM){
00459 FILE *mask_file = NULL;
00460 char mask_name[256];
00461 int idbg=0;
00462
00463 sprintf( mask_name, "skymask.reg");
00464 mask_file = fopen( mask_name, "w");
00465
00466 fprintf( mask_file,"# Region file format: DS9 version 4.1\n");
00467 fprintf( mask_file,"global color=green dashlist=8 3 width=1 font=\"helvetica 10 normal\"\
00468 select=1 highlite=1 dash=0 fixed=0 edit=1 move=1 delete=1 include=1 source=1\n");
00469 fprintf(mask_file,"physical\n");
00470
00471 for(idbg=0; idbg< nlambda; idbg++){
00472
00473 if (sky_mask[idbg] == 1){
00474 fprintf( mask_file, "line(%d,%d,%d,%d)\n",idbg+1, nslit, idbg+1, 1);
00475 }
00476 }
00477 fclose( mask_file);
00478 }
00479 }
00480
00481 check( loc_list = xsh_localization_create());
00482
00483 #if 0
00484
00485 {
00486 float arcsec ;
00487 switch( xsh_instrument_get_arm( instrument ) ) {
00488 case XSH_ARM_NIR:
00489 arcsec = XSH_ARCSEC_NIR ;
00490 break ;
00491 case XSH_ARM_UVB:
00492 arcsec = XSH_ARCSEC_UVB ;
00493 break ;
00494 case XSH_ARM_VIS:
00495 arcsec = XSH_ARCSEC_VIS ;
00496 break ;
00497 default:
00498 arcsec = 0.15 ;
00499 break ;
00500 }
00501 slit_nod = loc_obj_par->nod_step/arcsec ;
00502 xsh_msg_dbg_medium( "nod_step: %lf, slit_nod: %d",
00503 loc_obj_par->nod_step, slit_nod ) ;
00504 }
00505 #endif
00506 XSH_CALLOC( slit_tab, double, nslit);
00507 check( gaussval = cpl_vector_wrap( nslit, slit_tab));
00508 check( gausspos = cpl_vector_new( nslit));
00509 for( islit=0; islit < nslit; islit++){
00510 cpl_vector_set( gausspos, islit, islit);
00511 }
00512
00513 chunk_size = (nlambda-1) / (double)(loc_obj_par->loc_chunk_nb);
00514
00515
00516
00517
00518 for( ichunk=0; ichunk< loc_obj_par->loc_chunk_nb; ichunk++){
00519
00520
00521
00522
00523
00524 int ifirst=0,ilast=0,icenter=0;
00525 double slit_cen=0, slit_low=0, slit_up=0;
00526 double threshold = 0.0;
00527 int is_valid_chunk = 1;
00528
00529 ifirst = ichunk*chunk_size;
00530 ilast = (ichunk+1)*chunk_size;
00531 icenter = (int)((ifirst+ilast)/2);
00532
00533 xsh_msg_dbg_medium("%d-%d", ifirst,ilast);
00534
00535
00536 chunk_coadd( slit_tab, flux, qual, nlambda, nslit, ifirst, ilast, decode_bp, sky_mask);
00537
00538 if (level >= XSH_DEBUG_LEVEL_MEDIUM){
00539 FILE *coadd_file = NULL;
00540 char coadd_name[256];
00541 int icoadd=0;
00542 const char* filename = cpl_frame_get_filename( merge_frame);
00543
00544 sprintf( coadd_name, "%s_%d_%d.dat",filename,ifirst,ilast);
00545 coadd_file = fopen( coadd_name, "w");
00546
00547 for(icoadd=0; icoadd<nslit; icoadd++){
00548 fprintf(coadd_file, "%d %f\n",icoadd,slit_tab[icoadd]);
00549 }
00550 fclose( coadd_file);
00551 }
00552
00553 if ( loc_obj_par->method == LOC_GAUSSIAN_METHOD){
00554 double cenpos=0, sigma=0,area=0, offset=0;
00555 double kappa=3.0;
00556
00557 xsh_msg_dbg_medium("Using GAUSSIAN method to fit chunk %d_%d", ifirst,ilast);
00558 cpl_vector_fit_gaussian( gausspos, NULL, gaussval, NULL, CPL_FIT_ALL, &cenpos,
00559 &sigma,
00560 &area, &offset, NULL,
00561 NULL, NULL);
00562 if( cpl_error_get_code() == CPL_ERROR_CONTINUE ){
00563 xsh_msg_dbg_low("CONTINUE to fit %d_%d x0 %f sigma %f", ifirst,ilast,cenpos,sigma);
00564 }
00565 if( cpl_error_get_code() == CPL_ERROR_NONE ){
00566 slit_cen = cenpos;
00567 slit_low = cenpos-kappa*sigma;
00568 slit_up = cenpos+kappa*sigma;
00569 if (slit_low < 0) slit_low = 0;
00570 if (slit_up >= nslit) slit_up = nslit-1;
00571 }
00572 else{
00573 xsh_msg_dbg_low("FAILED to fit %d_%d", ifirst,ilast);
00574 is_valid_chunk=0;
00575 xsh_error_reset();
00576 }
00577 }
00578 else
00579 {
00580 int ymax=0, ymin=0;
00581 int ylow=0, yup=0;
00582 tab_minmax_get_index( slit_tab, nslit, slitlimit_par, &ymin, &ymax);
00583 if (ymin == ymax){
00584 xsh_msg_warning( "No maximum found in chunk (skip it)");
00585 is_valid_chunk=0;
00586 }
00587 else{
00588
00589
00590 xsh_msg_dbg_medium("ymin [%d] = %f ymax [%d] = %f",ymin,
00591 slit_tab[ymin], ymax, slit_tab[ymax]);
00592
00593 threshold = slit_tab[ymin]+(slit_tab[ymax]-slit_tab[ymin])*
00594 loc_obj_par->loc_thresh;
00595 xsh_msg_dbg_medium("Threshold at %f",threshold);
00596
00597
00598 yup = ymax+1;
00599 while( yup < nslit && (slit_tab[yup] >= threshold)) {
00600 yup++;
00601 }
00602
00603
00604 ylow = ymax-1;
00605 while( ylow > 0 && (slit_tab[ylow] >= threshold)) {
00606 ylow--;
00607 }
00608
00609 slit_cen = ymax;
00610 slit_low = ylow;
00611 slit_up = yup;
00612 }
00613 }
00614
00615 if ( is_valid_chunk == 1){
00616 chunk_center[ichunk-skip_chunk] = lambda_min+lambda_step*icenter;
00617 slit_center[ichunk-skip_chunk] = slit_min+slit_step*slit_cen;
00618 slit_lower[ichunk-skip_chunk] = slit_min+slit_step*slit_low;
00619 slit_upper[ichunk-skip_chunk] = slit_min+slit_step*slit_up;
00620 }
00621 else{
00622 skip_chunk++;
00623 }
00624 }
00625
00626 nb_chunk = loc_obj_par->loc_chunk_nb-skip_chunk;
00627
00628 XSH_ASSURE_NOT_ILLEGAL(loc_obj_par->loc_deg_poly < nb_chunk);
00629
00630 check( chunk_center_v = cpl_vector_wrap(
00631 nb_chunk, chunk_center));
00632 check( slit_center_v = cpl_vector_wrap(
00633 nb_chunk, slit_center));
00634
00635 check( cen_poly =
00636 xsh_polynomial_fit_1d_create( chunk_center_v, slit_center_v,
00637 loc_degree, NULL));
00638
00639
00640 if ( niter > 0){
00641 xsh_msg_dbg_low("Doing sigma clipping");
00642 iter = 0;
00643
00644 while ( (loc_degree < (nb_chunk-nbrej)) && iter < niter
00645 && nbrej != 0){
00646 double sigma_med;
00647
00648 nbrej = 0;
00649 xsh_msg_dbg_medium( " *** NBITER = %d / %d ***", iter+1, niter);
00650 dispslit = cpl_vector_new( nb_chunk);
00651
00652 for(ichunk = 0; ichunk < nb_chunk; ichunk++){
00653 double slit_fit;
00654 double lambda;
00655 double slit_diff;
00656
00657 lambda = chunk_center[ichunk];
00658
00659 check( slit_fit = cpl_polynomial_eval_1d( cen_poly,
00660 lambda, NULL));
00661
00662 slit_diff = slit_center[ichunk]-slit_fit;
00663 xsh_msg_dbg_low("slit_center %f FIT %f, DIFF %d %f", slit_center[ichunk], slit_fit, ichunk, slit_diff);
00664 check( cpl_vector_set( dispslit, ichunk, slit_diff));
00665 }
00666
00667 check( sigma_med = cpl_vector_get_stdev( dispslit));
00668 xsh_msg_dbg_medium(" kappa %f SIGMA MEDIAN = %f", kappa, sigma_med);
00669
00670 for(ichunk = 0; ichunk < nb_chunk; ichunk++){
00671 if ( fabs(cpl_vector_get( dispslit, ichunk)) > (kappa * sigma_med) ){
00672 nbrej++;
00673 }
00674 else{
00675 chunk_center[ichunk-nbrej] = chunk_center[ichunk];
00676 slit_center[ichunk-nbrej] = slit_center[ichunk];
00677 slit_lower[ichunk-nbrej] = slit_lower[ichunk];
00678 slit_upper[ichunk-nbrej] = slit_upper[ichunk];
00679 }
00680 }
00681
00682 xsh_msg_dbg_medium(" Removed %d points", nbrej);
00683
00684 nb_chunk -= nbrej;
00685
00686
00687 xsh_unwrap_vector( &chunk_center_v);
00688 xsh_unwrap_vector( &slit_center_v);
00689 xsh_free_polynomial( &cen_poly);
00690
00691 check( chunk_center_v = cpl_vector_wrap(
00692 nb_chunk, chunk_center));
00693 check( slit_center_v = cpl_vector_wrap(
00694 nb_chunk, slit_center));
00695
00696 check( cen_poly =
00697 xsh_polynomial_fit_1d_create( chunk_center_v, slit_center_v,
00698 loc_degree, NULL));
00699 xsh_free_vector( &dispslit);
00700 iter++;
00701 }
00702 }
00703
00704 if (level >= XSH_DEBUG_LEVEL_MEDIUM){
00705 FILE *debug_file = NULL;
00706 char debug_name[256];
00707 int idebug=0;
00708 const char* filename = cpl_frame_get_filename( merge_frame);
00709
00710 sprintf( debug_name, "%s_points.dat",filename);
00711 debug_file = fopen( debug_name, "w");
00712
00713 fprintf( debug_file,"#chunk_pos slit_low slit_cen slit_up\n");
00714
00715 for(idebug=0; idebug<nb_chunk; idebug++){
00716 fprintf( debug_file, "%f %f %f %f\n",chunk_center[idebug],
00717 slit_lower[idebug], slit_center[idebug], slit_upper[idebug]);
00718 }
00719 fclose( debug_file);
00720 }
00721
00722
00723 check( slit_lower_v = cpl_vector_wrap(
00724 nb_chunk, slit_lower));
00725 check( slit_upper_v = cpl_vector_wrap(
00726 nb_chunk, slit_upper));
00727 loc_list->pol_degree = loc_degree;
00728
00729 check(loc_list->cenpoly = cpl_polynomial_duplicate( cen_poly));
00730 check(loc_list->edglopoly =
00731 xsh_polynomial_fit_1d_create( chunk_center_v, slit_lower_v,
00732 loc_degree, NULL));
00733 check(loc_list->edguppoly =
00734 xsh_polynomial_fit_1d_create( chunk_center_v, slit_upper_v,
00735 loc_degree, NULL));
00736
00737 cleanup:
00738 xsh_unwrap_vector( &chunk_center_v);
00739 xsh_unwrap_vector( &slit_center_v);
00740 xsh_unwrap_vector( &slit_upper_v);
00741 xsh_unwrap_vector( &slit_lower_v);
00742 xsh_unwrap_vector( &gaussval);
00743 xsh_free_vector( &gausspos);
00744 XSH_FREE( chunk_center);
00745 XSH_FREE( slit_center);
00746 XSH_FREE( slit_upper);
00747 XSH_FREE( slit_lower);
00748 XSH_FREE( slit_tab);
00749 XSH_FREE( sky_mask);
00750 xsh_free_table( &skymask_table);
00751 xsh_free_polynomial( &cen_poly);
00752 xsh_spectrum_free( &spectrum);
00753 return loc_list;
00754
00755 }
00756
00757
00758
00759
00760
00761
00762
00763
00764
00765
00766
00767
00768
00769
00770
00771
00772
00773
00774
00775
00776
00777
00778
00779
00780
00781
00782
00783
00784
00785
00786
00787
00788
00789
00790
00791
00792
00793
00794
00795
00796
00797
00798
00799
00800
00801
00802
00803
00804
00805
00806
00807
00808
00809
00810
00811
00812
00813
00814
00815
00816
00817
00818
00819
00820
00821
00822
00823
00824
00825
00826
00827
00828
00829
00830
00831
00832
00833
00834
00835
00836
00837
00838
00839
00840
00841
00842
00843
00844
00845
00846
00847
00848
00849
00850
00851
00852
00853
00854
00855
00856
00857
00858
00859
00860
00861
00862
00863
00864
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 cpl_frameset * xsh_localize_obj_ifu( cpl_frameset *rec_frameset,
00909 cpl_frame *skymask_frame,
00910 xsh_instrument * instrument,
00911 xsh_localize_obj_param * locobj_par,
00912 xsh_slit_limit_param *slitlimit_par)
00913 {
00914 int i, slitlet;
00915 cpl_frameset *result_frameset = NULL;
00916 char fname[256];
00917
00918 XSH_ASSURE_NOT_NULL( rec_frameset);
00919 XSH_ASSURE_NOT_NULL( instrument);
00920 XSH_ASSURE_NOT_NULL( locobj_par);
00921
00922 check( result_frameset = cpl_frameset_new());
00923
00924 for( i = 0, slitlet = LOWER_IFU_SLITLET ; i < 3 ; i++, slitlet++ ) {
00925 cpl_frame * loc_frame = NULL;
00926 cpl_frame *rec_frame = NULL;
00927
00928 sprintf( fname ,"LOCALIZE_TABLE_%s_IFU_%s.fits", SlitletName[slitlet],
00929 xsh_instrument_arm_tostring( instrument));
00930
00931 xsh_msg( "Localizing slitlet %s, frame '%s'", SlitletName[slitlet], fname);
00932
00933 check( rec_frame = cpl_frameset_get_frame( rec_frameset, i));
00934
00935 check( loc_frame = xsh_localize_obj( rec_frame, skymask_frame, instrument,
00936 locobj_par,slitlimit_par, fname));
00937 check( cpl_frameset_insert( result_frameset, loc_frame));
00938 }
00939
00940 cleanup:
00941 if ( cpl_error_get_code() != CPL_ERROR_NONE ) {
00942 xsh_free_frameset( &result_frameset);
00943 }
00944 return result_frameset;
00945 }
00946
00947
00948
00954 double xsh_convert_seeing( cpl_frame* frame)
00955 {
00956 double mu=-1.0;
00957 double avg_seeing, seeing_start, seeing_end;
00958 double avg_airmass;
00959 double k;
00960 const char* filename = NULL;
00961 cpl_propertylist *header = NULL;
00962
00963 XSH_ASSURE_NOT_NULL( frame);
00964
00965 check( filename = cpl_frame_get_filename( frame));
00966 check( header = cpl_propertylist_load( filename, 0));
00967 check( avg_airmass = xsh_pfits_get_airm_mean( header));
00968 check( seeing_start = xsh_pfits_get_seeing_start( header));
00969 check( seeing_end = xsh_pfits_get_seeing_end( header));
00970 avg_seeing = 0.5*(seeing_start+seeing_end);
00971
00972 k = sqrt( 1.0-78.08*pow(XSH_LAMBDA_DIMM,0.4)*
00973 pow(avg_airmass,-0.6)*pow(avg_seeing,-1.0/3.0));
00974
00975 xsh_msg("K %f", k);
00976
00977 mu = 1.0/CPL_MATH_FWHM_SIG*avg_seeing*pow(XSH_LAMBDA_DIMM,0.2)*
00978 pow(avg_airmass,0.6)*k;
00979
00980 xsh_msg("Mu %f", mu);
00981
00982 cleanup:
00983 if ( cpl_error_get_code() != CPL_ERROR_NONE ) {
00984 mu =-1;
00985 }
00986 xsh_free_propertylist( &header);
00987 return mu;
00988 }
00989
00990
00991
01013 cpl_frame* xsh_localize_ifu_slitlet( cpl_frame *merge2d_slitlet,
01014 cpl_frame *skymask_frame, int smooth_hsize,
01015 int nscales, int HF_skip, const char* resname, double cut_sigma_low,
01016 double cut_sigma_up, double cut_snr_low, double cut_snr_up,
01017 double slit_min, double slit_max, int deg, int box_hsize,
01018 xsh_instrument *instrument)
01019
01020 {
01021 cpl_frame *result = NULL;
01022 double wmin, wstep;
01023
01024 double smin, sstep;
01025
01026 int wsize, ssize;
01027 xsh_spectrum *spectrum2d = NULL;
01028 double *flux = NULL;
01029 double *errs = NULL;
01030 int *qual = NULL;
01031 int i, j, k;
01032 int level;
01033 double *slit_vect_data = NULL;
01034 double *sliterr_vect_data = NULL;
01035 double *slit_pos_data = NULL;
01036 int slit_vect_size;
01037 cpl_vector *slit_vect = NULL;
01038 cpl_vector *sliterr_vect = NULL;
01039 cpl_vector *slit_pos = NULL;
01040 cpl_vector *smooth_vect = NULL;
01041 double *spos_data = NULL;
01042 double *errpos_data = NULL;
01043 double *wpos_data = NULL;
01044 double *sigma_data = NULL;
01045 double *snr_data = NULL;
01046 int data_size=0, ndata_size=0;
01047 cpl_vector *spos_vect = NULL;
01048 cpl_vector *wpos_vect = NULL;
01049 cpl_vector *sigma_vect = NULL;
01050 cpl_vector *snr_vect = NULL;
01051 cpl_vector *sigma_sort_vect = NULL;
01052 cpl_vector *snr_sort_vect = NULL;
01053 double *sposg_data = NULL;
01054 double *wposg_data = NULL;
01055 cpl_matrix *decomp = NULL;
01056 int nb_scales;
01057 cpl_table *table = NULL;
01058 char tablename[256];
01059 cpl_propertylist *header = NULL;
01060 double n05, n95, snr_05, snr_95, sigma_05, sigma_95;
01061 int jmin, jmax;
01062
01063 int *sky_mask = NULL;
01064 cpl_table *skymask_table = NULL;
01065 const char* skymask_name = NULL;
01066
01067 XSH_ASSURE_NOT_NULL( merge2d_slitlet);
01068 XSH_ASSURE_NOT_ILLEGAL( smooth_hsize >= 0);
01069 XSH_ASSURE_NOT_ILLEGAL( nscales > 0);
01070 XSH_ASSURE_NOT_ILLEGAL( nscales > HF_skip);
01071 XSH_ASSURE_NOT_ILLEGAL( box_hsize >= 0);
01072
01073
01074
01075
01076
01077
01078
01079
01080
01081
01082
01083
01084
01085
01086
01087 level = xsh_debug_level_get();
01088 check( spectrum2d = xsh_spectrum_load( merge2d_slitlet));
01089 check( wmin = xsh_spectrum_get_lambda_min( spectrum2d));
01090
01091 check( wstep = xsh_spectrum_get_lambda_step( spectrum2d));
01092
01093 if ( skymask_frame != NULL){
01094 check( skymask_name = cpl_frame_get_filename( skymask_frame));
01095 xsh_msg_dbg_medium("Sky mask %s", skymask_name);
01096 }
01097 wsize = spectrum2d->size_lambda;
01098
01099
01100 smin = spectrum2d->slit_min;
01101
01102 sstep = spectrum2d->slit_step;
01103 ssize = spectrum2d->size_slit;
01104
01105 jmin = 0;
01106 jmax = ssize;
01107
01108
01109 while( (smin+jmin*sstep) < slit_min){
01110 jmin++;
01111 }
01112 while( (smin+jmax*sstep) > slit_max){
01113 jmax--;
01114 }
01115
01116
01117 if ( jmin > jmax || jmax <= 0 || jmin >= ssize){
01118 jmin = 0;
01119 jmax = ssize;
01120 }
01121 xsh_msg_dbg_medium( "Use [%d-%d] from [%d %d] slitlet",
01122 jmin, jmax, 0, ssize);
01123
01124 check( flux = xsh_spectrum_get_flux( spectrum2d));
01125 check( errs = xsh_spectrum_get_errs( spectrum2d));
01126 check( qual = xsh_spectrum_get_qual( spectrum2d));
01127
01128 XSH_MALLOC( slit_vect_data, double , ssize);
01129 XSH_MALLOC( sliterr_vect_data, double, ssize);
01130 XSH_MALLOC( slit_pos_data, double , ssize);
01131
01132 XSH_MALLOC( spos_data, double , wsize);
01133 XSH_MALLOC( errpos_data, double, wsize);
01134 XSH_MALLOC( wpos_data, double , wsize);
01135
01136 XSH_MALLOC( sigma_data, double , wsize);
01137 XSH_MALLOC( snr_data, double , wsize);
01138
01139
01140 XSH_CALLOC( sky_mask, int, wsize);
01141
01142 if ( skymask_frame != NULL){
01143 float *skymask_data = NULL;
01144 int irow, nrow;
01145 double fwhm =0.0, sky_min, sky_max;
01146 int isky_min, isky_max, imask;
01147 double width, resolution;
01148
01149 XSH_TABLE_LOAD( skymask_table, skymask_name);
01150
01151 check( xsh_sort_table_1( skymask_table, "WAVELENGTH", CPL_FALSE));
01152 check( skymask_data = cpl_table_get_data_float( skymask_table,
01153 "WAVELENGTH"));
01154 check( nrow = cpl_table_get_nrow( skymask_table));
01155
01156 for( irow=0; irow < nrow; irow++){
01157 check( width = xsh_pfits_get_slit_width( spectrum2d->flux_header, instrument));
01158 resolution = xsh_resolution_get( instrument, width);
01159 fwhm = skymask_data[irow]/resolution;
01160 sky_min = skymask_data[irow]-fwhm;
01161 sky_max = skymask_data[irow]+fwhm;
01162 isky_min = (int)xsh_round_double((sky_min-wmin)/wstep);
01163 if (isky_min < 0){
01164 isky_min=0;
01165 }
01166 isky_max = (int)xsh_round_double((sky_max-wmin)/wstep);
01167 if ( isky_max >= wsize){
01168 isky_max = wsize-1;
01169 }
01170 for( imask=isky_min; imask <=isky_max; imask++){
01171 sky_mask[imask] = 1;
01172 }
01173 }
01174
01175 if (level >= XSH_DEBUG_LEVEL_MEDIUM){
01176 FILE *mask_file = NULL;
01177 char mask_name[256];
01178 int idbg=0;
01179
01180 sprintf( mask_name, "skymask.reg");
01181 mask_file = fopen( mask_name, "w");
01182
01183 fprintf( mask_file,"# Region file format: DS9 version 4.1\n");
01184 fprintf( mask_file,"global color=green dashlist=8 3 width=1 font=\"helvetica 10 normal\"\
01185 select=1 highlite=1 dash=0 fixed=0 edit=1 move=1 delete=1 include=1 source=1\n");
01186 fprintf(mask_file,"physical\n");
01187
01188 for(idbg=0; idbg< wsize; idbg++){
01189
01190 if (sky_mask[idbg] == 1){
01191 fprintf( mask_file, "line(%d,%d,%d,%d)\n",idbg+1, ssize, idbg+1, 1);
01192 }
01193 }
01194 fclose( mask_file);
01195 }
01196 }
01197
01198
01199 data_size = 0;
01200
01201 for (i = 0; i < wsize; i++) {
01202 double slitcen = 0;
01203 double sigma = 0;
01204 double area = 0;
01205 double offset = 0;
01206 double a = 0, b = 0, c = 0;
01207 double frac_rej = 0;
01208
01209 cpl_matrix *covariance = NULL;
01210 double fit_err;
01211 int good_fit = 0;
01212
01213 slit_vect_size = 0;
01214
01215
01216 for (j = jmin; j < jmax; j++) {
01217 double val = 0, err = 0;
01218 int start, end;
01219 int ngood = 0, nbad = 0;
01220
01221
01222 start = i - box_hsize;
01223 end = i + box_hsize;
01224
01225 if (start < 0) {
01226 start = 0;
01227 }
01228 if (end >= wsize) {
01229 end = wsize - 1;
01230 }
01231
01232 for (k = start; k <= end; k++) {
01233 if ( (qual[k + j * wsize] & instrument->decode_bp) == 0 ) {
01234 val += flux[k + j * wsize];
01235 ngood++;
01236 } else {
01237 nbad++;
01238 }
01239
01240 }
01241
01242 val += flux[i+j*wsize];
01243 err = errs[i+j*wsize];
01244
01245 if ( ngood > 0){
01246 val *= (double)(nbad+ngood)/(double)ngood;
01247
01248 slit_pos_data[slit_vect_size] = j;
01249 sliterr_vect_data[slit_vect_size] = err;
01250 slit_vect_data[slit_vect_size] = val;
01251
01252 slit_vect_size++;
01253 }
01254 }
01255
01256 frac_rej = 1-((double)slit_vect_size/(double)ssize);
01257 if ( frac_rej < 0.5){
01258 check( slit_pos = cpl_vector_wrap( slit_vect_size, slit_pos_data));
01259 check( slit_vect = cpl_vector_wrap( slit_vect_size, slit_vect_data));
01260 check( sliterr_vect = cpl_vector_wrap( slit_vect_size, sliterr_vect_data));
01261
01262 check( smooth_vect = cpl_vector_filter_median_create( slit_vect,
01263 smooth_hsize));
01264
01265
01266 #if CPL_MODE
01267 cpl_vector_fit_gaussian( slit_pos, NULL, slit_vect, sliterr_vect,
01268 CPL_FIT_ALL, &slitcen, &sigma, &area, &offset, &mse, NULL, &covariance);
01269 a = offset;
01270 fit_err = sqrt(cpl_matrix_get(covariance,0,0))*sstep;
01271 good_fit = (cpl_error_get_code() == CPL_ERROR_NONE);
01272 #else
01273 {
01274 double init_par[6];
01275 double fit_errs[6];
01276 int status = 0;
01277 xsh_gsl_init_gaussian_fit( slit_pos, slit_vect, init_par);
01278 xsh_gsl_fit_gaussian( slit_pos, slit_vect, deg, init_par, fit_errs, &status);
01279 area = init_par[0];
01280 a = init_par[1];
01281 b = init_par[2];
01282 c = init_par[3];
01283 slitcen = init_par[4];
01284 sigma = init_par[5];
01285 offset = a+b*slitcen+c*slitcen*slitcen;
01286 fit_err = fit_errs[4]*sstep;
01287 good_fit= (status == 0);
01288 }
01289 #endif
01290
01291 if ( good_fit){
01292 double height=0;
01293 double snr=0;
01294
01295 height = area / sqrt(2*M_PI*sigma*sigma);
01296 snr = (height+offset)/offset;
01297
01298
01299 if ( (area > 0) && (slitcen > 0) && (slitcen < ssize)
01300 && (sigma < ssize) && (snr > 0)){
01301 spos_data[data_size] = smin+slitcen*sstep;
01302 errpos_data[data_size] = fit_err;
01303 wpos_data[data_size] = wmin+i*wstep;
01304 sigma_data[data_size] = sigma*sstep;
01305 snr_data[data_size] = snr;
01306 #if 0
01307
01308 if ( (wpos_data[data_size] >= (1500-wstep)) && (wpos_data[data_size] <= (1500+wstep))){
01309 char test_name[256];
01310 FILE* test_file = NULL;
01311 int itest;
01312 double dtest;
01313
01314 sprintf( test_name, "data_w%f_%s.dat", wpos_data[data_size],resname);
01315 XSH_REGDEBUG("Produce test file %s", test_name);
01316 test_file = fopen( test_name, "w+");
01317 fprintf( test_file, "# pos slit\n");
01318
01319 for( itest=0; itest < slit_vect_size; itest++){
01320 fprintf( test_file, "%f %f\n", cpl_vector_get( slit_pos, itest), cpl_vector_get( slit_vect, itest));
01321 }
01322
01323 fclose( test_file);
01324
01325 sprintf( test_name, "gauss_w%f_%s.dat", wpos_data[data_size],resname);
01326 XSH_REGDEBUG("Produce test file %s", test_name);
01327
01328 test_file = fopen( test_name, "w+");
01329 fprintf( test_file, "# pos gauss offset x0=%f_sig=%f_area=%f_offset=%f\n", slitcen, sigma, area, offset);
01330
01331 for( dtest=0; dtest < slit_vect_size; dtest+=0.1){
01332 double z, gauss;
01333 double off;
01334
01335 off = a+b*dtest+c*dtest*dtest;
01336 z = ( dtest-slitcen)/(sigma*XSH_MATH_SQRT_2);
01337 gauss = height*exp(-(z*z))+off;
01338
01339 fprintf( test_file, "%f %f %f\n", dtest, gauss, off);
01340 }
01341
01342 fclose( test_file);
01343 }
01344 #endif
01345 data_size++;
01346 }
01347 }
01348 else{
01349 xsh_error_reset();
01350 }
01351 xsh_free_matrix( &covariance);
01352 xsh_free_vector( &smooth_vect);
01353 xsh_unwrap_vector( &slit_pos);
01354 xsh_unwrap_vector( &slit_vect);
01355 xsh_unwrap_vector( &sliterr_vect);
01356 }
01357 }
01358
01359
01360
01361 check( sigma_vect = cpl_vector_wrap( data_size, sigma_data));
01362 check( snr_vect = cpl_vector_wrap( data_size, snr_data));
01363 check( sigma_sort_vect = cpl_vector_duplicate( sigma_vect));
01364 check( snr_sort_vect = cpl_vector_duplicate( snr_vect));
01365
01366 check( cpl_vector_sort( sigma_sort_vect, 1));
01367 check( cpl_vector_sort( snr_sort_vect, 1));
01368
01369 n05 = (data_size-1)*cut_sigma_low;
01370 i = ceil( n05);
01371 sigma_05 = cpl_vector_get( sigma_sort_vect, i);
01372 n05 = (data_size-1)*cut_snr_low;
01373 i = ceil( n05);
01374 snr_05 = cpl_vector_get( snr_sort_vect, i);
01375
01376
01377 n95 = (data_size-1)*cut_sigma_up;
01378 i = floor( n95);
01379 sigma_95 = cpl_vector_get( sigma_sort_vect, i);
01380 n95 = (data_size-1)*cut_snr_up;
01381 i = floor( n95);
01382 snr_95 = cpl_vector_get( snr_sort_vect, i);
01383
01384
01385 for( i=0; i< data_size; i++){
01386 double sigma, snr;
01387
01388 sigma = sigma_data[i];
01389 snr = snr_data[i];
01390 if ( sigma_05 <= sigma && sigma <= sigma_95 && snr_05 <= snr && snr <=snr_95){
01391 wpos_data[ndata_size] = wpos_data[i];
01392 spos_data[ndata_size] = spos_data[i];
01393 errpos_data[ndata_size] = errpos_data[i];
01394 ndata_size++;
01395 }
01396 }
01397 xsh_msg_dbg_low( "Filtering by snr [%f,%f] and sigma [%f,%f] from %d lines to %d lines",
01398 sigma_05, sigma_95, snr_05, snr_95, data_size, ndata_size);
01399
01400
01401 if (level >= XSH_DEBUG_LEVEL_MEDIUM){
01402 FILE *test_file = NULL;
01403 int itest;
01404 char test_name[256];
01405 #if CPL_MODE
01406 sprintf( test_name, "cpl_gaussian_fit_%s.dat", resname);
01407 #else
01408 sprintf( test_name, "gsl_gaussian_fit_%s.dat", resname);
01409 #endif
01410 test_file = fopen( test_name, "w");
01411 fprintf( test_file, "# wavelength slit_fit fit_err sigma\n");
01412
01413 for(itest=0; itest<ndata_size; itest++){
01414 fprintf( test_file, "%f %f %f %f\n", wpos_data[itest], spos_data[itest], errpos_data[itest],
01415 sigma_data[itest]);
01416 }
01417 xsh_msg_dbg_medium( "Produce file %s", test_name);
01418 fclose( test_file);
01419 }
01420
01421 check( wpos_vect = cpl_vector_new( wsize));
01422 check( spos_vect = cpl_vector_new( wsize));
01423 check( sposg_data = cpl_vector_get_data( spos_vect));
01424 check( wposg_data = cpl_vector_get_data( wpos_vect));
01425 j=0;
01426
01427 for(i=0; i< wsize; i++){
01428 double wave, wkeep;
01429 double slit;
01430
01431 wave = wmin+i*wstep;
01432 wkeep = wpos_data[j];
01433 check( cpl_vector_set( wpos_vect, i, wave));
01434 if ( fabs(wave -wkeep) < 0.0000001){
01435 slit = spos_data[j];
01436 if (j < (ndata_size-1)){
01437 j++;
01438 }
01439 }
01440 else{
01441 slit = xsh_data_interpolate( wave, ndata_size, wpos_data, spos_data);
01442 }
01443 check( cpl_vector_set( spos_vect, i, slit));
01444 }
01445
01446 check( decomp = xsh_atrous( spos_vect, nscales));
01447
01448 nb_scales = nscales-HF_skip;
01449
01450
01451 for( i=0; i< wsize; i++){
01452 sposg_data[i] = 0;
01453 for(j=0; j<nb_scales; j++){
01454 sposg_data[i] += cpl_matrix_get( decomp, j, i);
01455 }
01456 }
01457
01458
01459 check( table = cpl_table_new( wsize));
01460 XSH_TABLE_NEW_COL(table, XSH_OBJPOS_COLNAME_WAVELENGTH,
01461 XSH_OBJPOS_UNIT_WAVELENGTH, CPL_TYPE_DOUBLE);
01462 XSH_TABLE_NEW_COL(table, XSH_OBJPOS_COLNAME_SLIT,
01463 XSH_OBJPOS_UNIT_SLIT, CPL_TYPE_DOUBLE);
01464
01465 for( i=0; i< wsize; i++){
01466 check( cpl_table_set_double( table, XSH_OBJPOS_COLNAME_WAVELENGTH,
01467 i, wposg_data[i]));
01468 check( cpl_table_set_double( table, XSH_OBJPOS_COLNAME_SLIT,
01469 i, sposg_data[i]));
01470 }
01471 sprintf( tablename, resname);
01472 header = cpl_propertylist_new();
01473 check( cpl_table_save( table, header, NULL, tablename, CPL_IO_DEFAULT));
01474
01475 check(result=xsh_frame_product( tablename,
01476 "OBJPOS_TAB",
01477 CPL_FRAME_TYPE_TABLE,
01478 CPL_FRAME_GROUP_PRODUCT,
01479 CPL_FRAME_LEVEL_TEMPORARY));
01480
01481
01482 check (xsh_add_temporary_file( tablename));
01483
01484 cleanup:
01485 XSH_TABLE_FREE( table);
01486 xsh_free_propertylist( &header);
01487 xsh_free_matrix( &decomp);
01488 xsh_free_vector( &spos_vect);
01489 xsh_free_vector( &wpos_vect);
01490 xsh_free_vector( &sigma_sort_vect);
01491 xsh_free_vector( &snr_sort_vect);
01492 xsh_unwrap_vector( &sigma_vect);
01493 xsh_unwrap_vector( &snr_vect);
01494 XSH_FREE( sigma_data);
01495 XSH_FREE( snr_data);
01496 XSH_FREE( spos_data);
01497 XSH_FREE( errpos_data);
01498 XSH_FREE( wpos_data);
01499 XSH_FREE( slit_pos_data);
01500 XSH_FREE( slit_vect_data);
01501 XSH_FREE( sliterr_vect_data);
01502 xsh_free_vector( &smooth_vect);
01503 xsh_spectrum_free( &spectrum2d);
01504 XSH_FREE(sky_mask);
01505 xsh_free_table( &skymask_table);
01506 return result;
01507 }
01508
01509
01510
01522 cpl_frameset* xsh_localize_ifu( cpl_frameset *merge2d_frameset,
01523 cpl_frame *skymask_frame,
01524 xsh_localize_ifu_param * locifu_par, xsh_instrument *instrument,
01525 const char* prefix)
01526 {
01527 int i, slitlet;
01528 cpl_frameset *result_frameset = NULL;
01529 char fname[256];
01530 int smooth_hsize;
01531 int nscales;
01532 int HF_skip;
01533 double cut_sigma_low;
01534 double cut_sigma_up;
01535 double cut_snr_low;
01536 double cut_snr_up;
01537 double slit_min = -6.0;
01538 double slit_max = 6.0;
01539 cpl_frame *frame = NULL;
01540 const char *frame_name = NULL;
01541 cpl_propertylist *header = NULL;
01542 int deg = 2, skymask=0;
01543 int box_hsize;
01544 cpl_frame *mask_frame = NULL;
01545
01546 XSH_ASSURE_NOT_NULL( merge2d_frameset);
01547 XSH_ASSURE_NOT_NULL( instrument);
01548 XSH_ASSURE_NOT_NULL( locifu_par);
01549
01550 smooth_hsize = locifu_par->smooth_hsize;
01551 nscales = locifu_par->nscales;
01552 HF_skip = locifu_par->HF_skip;
01553 cut_sigma_low = locifu_par->cut_sigma_low;
01554 cut_sigma_up = locifu_par->cut_sigma_up;
01555 cut_snr_low = locifu_par->cut_snr_low;
01556 cut_snr_up = locifu_par->cut_snr_up;
01557 skymask = locifu_par->use_skymask;
01558 box_hsize = locifu_par->box_hsize;
01559
01560 if (skymask){
01561 mask_frame = skymask_frame;
01562 }
01563
01564 check( frame = cpl_frameset_get_frame( merge2d_frameset, 0));
01565 check( frame_name = cpl_frame_get_filename( frame));
01566 check( header = cpl_propertylist_load( frame_name, 0));
01567 check( slit_min = xsh_pfits_get_rectify_space_min( header));
01568 xsh_free_propertylist( &header);
01569
01570 check( frame = cpl_frameset_get_frame( merge2d_frameset, 2));
01571 check( frame_name = cpl_frame_get_filename( frame));
01572 check( header = cpl_propertylist_load( frame_name, 0));
01573 check( slit_max = xsh_pfits_get_rectify_space_max( header));
01574 xsh_free_propertylist( &header);
01575
01576 slit_min += locifu_par->slitlow_edges_mask;
01577 slit_max -= locifu_par->slitup_edges_mask;
01578
01579 deg = locifu_par->bckg_deg;
01580
01581 check( result_frameset = cpl_frameset_new());
01582
01583 for( i = 0, slitlet = LOWER_IFU_SLITLET ; i < 3 ; i++, slitlet++ ) {
01584 cpl_frame * loc_frame = NULL;
01585 cpl_frame *merge2d_frame = NULL;
01586
01587 sprintf( fname ,"%s_LOCIFU_%s_%s.fits", prefix, SlitletName[slitlet],
01588 xsh_instrument_arm_tostring( instrument));
01589
01590 xsh_msg( "Localizing IFU in [%f,%f] slitlet %s, frame '%s'", slit_min, slit_max,
01591 SlitletName[slitlet], fname);
01592
01593 check( merge2d_frame = cpl_frameset_get_frame( merge2d_frameset, i));
01594
01595
01596 check( loc_frame = xsh_localize_ifu_slitlet( merge2d_frame, mask_frame,
01597 smooth_hsize, nscales,
01598 HF_skip, fname, cut_sigma_low, cut_sigma_up, cut_snr_low, cut_snr_up,
01599 slit_min, slit_max, deg, box_hsize, instrument));
01600
01601 check( cpl_frameset_insert( result_frameset, loc_frame));
01602 }
01603
01604 cleanup:
01605 if ( cpl_error_get_code() != CPL_ERROR_NONE ) {
01606 xsh_free_frameset( &result_frameset);
01607 xsh_free_propertylist( &header);
01608 }
01609 return result_frameset;
01610 }
01611