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
00037
00040
00041
00042
00043
00044 #include <math.h>
00045 #include <xsh_drl.h>
00046
00047 #include <xsh_utils_table.h>
00048 #include <xsh_badpixelmap.h>
00049 #include <xsh_utils_wrappers.h>
00050 #include <xsh_data_pre.h>
00051 #include <xsh_dfs.h>
00052 #include <xsh_pfits.h>
00053 #include <xsh_error.h>
00054 #include <xsh_msg.h>
00055 #include <xsh_fit.h>
00056 #include <xsh_data_instrument.h>
00057 #include <xsh_data_localization.h>
00058 #include <xsh_data_spectrum.h>
00059 #include <xsh_data_rec.h>
00060 #include <xsh_ifu_defs.h>
00061 #include <xsh_parameters.h>
00062 #include <cpl.h>
00063
00064
00065
00066
00067
00068 static void chunk_coadd( double* data, double* flux, int* qual, int nx, int ny,
00069 int ifirst, int ilast, int *skymask);
00070
00071 static void tab_minmax_get_index( double* data, int size,
00072 xsh_slit_limit_param * slit_limit_par,
00073 int* min, int* max);
00074
00075 static xsh_localization* xsh_localize_obj_auto( cpl_frame * rec_frame,
00076 cpl_frame *skymask_frame,
00077 xsh_instrument* instrument, xsh_localize_obj_param * loc_obj_par,
00078 xsh_slit_limit_param * slit_limit_par);
00079
00080 static xsh_localization* xsh_localize_obj_manual(
00081 xsh_instrument *instrument, xsh_localize_obj_param *loc_obj_par);
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105 static void chunk_coadd( double* data, double *flux, int* qual, int nx, int ny,
00106 int ifirst, int ilast, int *skymask)
00107 {
00108 int i,j;
00109 int* nbgoodpixels = NULL;
00110
00111 XSH_ASSURE_NOT_NULL(data);
00112 XSH_ASSURE_NOT_NULL(flux);
00113 XSH_ASSURE_NOT_NULL(qual);
00114 XSH_ASSURE_NOT_NULL(skymask);
00115 XSH_ASSURE_NOT_ILLEGAL(ifirst >= 0 && ifirst <= ilast && ilast < nx);
00116
00117 XSH_CALLOC( nbgoodpixels, int, ny);
00118
00119 for(j=0; j<ny; j++){
00120 data[j] = 0;
00121 }
00122
00123 for(i=ifirst; i<=ilast; i++){
00124 if (skymask[i] == 0){
00125 for( j=0 ; j<ny ; j++) {
00126 if (qual[i+j*nx] <= XSH_GOOD_PIXEL_LEVEL){
00127 nbgoodpixels[j]++;
00128 data[j] +=flux[i+j*nx];
00129 }
00130 else{
00131 xsh_msg_dbg_high("bad pixel at %d %d",i,j);
00132 }
00133 }
00134 }
00135 else{
00136 xsh_msg_dbg_medium("Mask by sky lines : column %d", i);
00137 }
00138 }
00139
00140 for(j=0; j<ny; j++){
00141 if ( nbgoodpixels[j] != 0 ) data[j] = data[j] / (float)(nbgoodpixels[j]);
00142 else data[j] = 0. ;
00143 }
00144
00145 cleanup:
00146 XSH_FREE(nbgoodpixels);
00147 return;
00148 }
00149
00150 static void tab_minmax_get_index( double* data, int size,
00151 xsh_slit_limit_param * slit_limit_par,
00152 int* min, int* max)
00153 {
00154 int i, maxi, mini ;
00155 int slit0, slit1 ;
00156
00157 XSH_ASSURE_NOT_NULL(data);
00158 XSH_ASSURE_NOT_NULL(min);
00159 XSH_ASSURE_NOT_NULL(max);
00160
00161 if ( slit_limit_par == NULL ) {
00162 slit0 = 0 ;
00163 slit1 = size ;
00164 }
00165 else {
00166 slit0 = slit_limit_par->min_slit_idx + 1 ;
00167 slit1 = slit_limit_par->max_slit_idx ;
00168 }
00169 maxi = 0 ;
00170 mini = 0 ;
00171
00172 for(i = slit0; i<slit1; i++) {
00173 if (data[i]>data[maxi]){
00174 maxi = i;
00175 }
00176 else if (data[i]<data[mini]){
00177 mini = i;
00178 }
00179 }
00180 *min = mini;
00181 *max = maxi;
00182
00183 cleanup:
00184 return;
00185 }
00186
00187
00188
00189
00190
00191
00209
00210 cpl_frame* xsh_localize_obj (cpl_frame *rec_frame,
00211 cpl_frame *skymask_frame,
00212 xsh_instrument* instrument,
00213 xsh_localize_obj_param * loc_obj_par,
00214 xsh_slit_limit_param * slitlimit_par,
00215 const char * fname)
00216 {
00217 cpl_frame *merge_frame = NULL;
00218 cpl_frame *res_frame = NULL;
00219 xsh_localization *loc_list = NULL;
00220 xsh_merge_param merge_par = {MEAN_MERGE_METHOD};
00221 const char* rec_prefix = "LOCALIZE";
00222
00223
00224
00225 XSH_ASSURE_NOT_NULL( loc_obj_par);
00226 XSH_ASSURE_NOT_NULL( instrument);
00227
00228 xsh_msg_dbg_medium( "Entering xsh_localize_obj") ;
00229
00230 if ( rec_frame != NULL){
00231 check( merge_frame = xsh_merge_ord( rec_frame, instrument, &merge_par,
00232 rec_prefix));
00233 }
00234
00235
00236 xsh_msg_dbg_medium("method %s",LOCALIZE_METHOD_PRINT( loc_obj_par->method));
00237
00238 if ( loc_obj_par->method == LOC_MANUAL_METHOD){
00239 xsh_msg_dbg_medium("slit position %f slit half height %f",
00240 loc_obj_par->slit_position, loc_obj_par->slit_hheight);
00241 check( loc_list = xsh_localize_obj_manual( instrument, loc_obj_par));
00242 }
00243 else{
00244 check( loc_list = xsh_localize_obj_auto( merge_frame, skymask_frame,
00245 instrument,
00246 loc_obj_par, slitlimit_par));
00247 }
00248
00249 xsh_msg_dbg_low( "Saving Localization Table" ) ;
00250 if ( fname == NULL ) {
00251 fname = "LOCALIZATION_TABLE.fits";
00252 }
00253 check( res_frame = xsh_localization_save( loc_list,
00254 fname, instrument));
00255
00256 cleanup:
00257 xsh_localization_free( &loc_list);
00258 xsh_free_frame( &merge_frame);
00259 return res_frame ;
00260 }
00261
00262
00263 static xsh_localization* xsh_localize_obj_manual(
00264 xsh_instrument *instrument, xsh_localize_obj_param *loc_obj_par)
00265 {
00266 xsh_localization *loc_list = NULL;
00267 int deg_poly = 0;
00268 double slit_cen, slit_up, slit_low;
00269 int pows = 0;
00270
00271
00272 XSH_ASSURE_NOT_NULL( loc_obj_par);
00273 XSH_ASSURE_NOT_NULL( instrument);
00274
00275
00276 xsh_msg_dbg_medium("slit position %f slit half height %f",
00277 loc_obj_par->slit_position, loc_obj_par->slit_hheight);
00278
00279 check( loc_list = xsh_localization_create());
00280
00281 slit_cen = loc_obj_par->slit_position;
00282 slit_up = slit_cen + loc_obj_par->slit_hheight;
00283 slit_low = slit_cen - loc_obj_par->slit_hheight;
00284
00285 loc_list->pol_degree = deg_poly;
00286 check( loc_list->cenpoly = cpl_polynomial_new( 1));
00287 check( loc_list->edglopoly = cpl_polynomial_new( 1));
00288 check( loc_list->edguppoly = cpl_polynomial_new( 1));
00289 check( cpl_polynomial_set_coeff( loc_list->cenpoly,
00290 &pows, slit_cen));
00291 check( cpl_polynomial_set_coeff( loc_list->edglopoly,
00292 &pows, slit_low));
00293 check( cpl_polynomial_set_coeff( loc_list->edguppoly,
00294 &pows, slit_up));
00295
00296 cleanup:
00297 return loc_list;
00298 }
00299
00300
00316
00317 static xsh_localization*
00318 xsh_localize_obj_auto( cpl_frame *merge_frame,
00319 cpl_frame *skymask_frame,
00320 xsh_instrument *instrument,
00321 xsh_localize_obj_param *loc_obj_par,
00322 xsh_slit_limit_param *slitlimit_par)
00323 {
00324 xsh_localization *loc_list = NULL;
00325 xsh_spectrum *spectrum = NULL;
00326 double* slit_center = NULL;
00327 double* slit_upper = NULL;
00328 double* slit_lower = NULL;
00329 double* chunk_center = NULL;
00330 cpl_vector* slit_center_v = NULL;
00331 cpl_vector* slit_upper_v = NULL;
00332 cpl_vector* slit_lower_v = NULL;
00333 cpl_vector* chunk_center_v = NULL;
00334
00335 #if 0
00336 int i = 0;
00337 int skip_order_nb = 0;
00338 int first_order_index = 0;
00339 int slit_nod, slit_lim ;
00340 #endif
00341
00342 int level;
00343
00344 int nslit, nlambda, chunk_size;
00345 double* slit_tab = NULL;
00346 double *flux = NULL;
00347 double *errs = NULL;
00348 int *qual = NULL;
00349 int ichunk, skip_chunk=0, nb_chunk;
00350 cpl_polynomial *cen_poly = NULL;
00351 int loc_degree;
00352 double lambda_min, lambda_step;
00353 double slit_min, slit_step, kappa=2.0;
00354 int islit, iter, niter=5, nbrej=-1;
00355 cpl_vector *dispslit = NULL;
00356 cpl_vector *gausspos = NULL;
00357 cpl_vector *gaussval = NULL;
00358 int *sky_mask = NULL;
00359 cpl_table *skymask_table = NULL;
00360 const char* skymask_name = NULL;
00361
00362
00363 XSH_ASSURE_NOT_NULL( loc_obj_par);
00364 XSH_ASSURE_NOT_NULL( instrument);
00365 XSH_ASSURE_NOT_NULL( merge_frame);
00366
00367
00368 kappa = loc_obj_par->kappa;
00369 niter = loc_obj_par->niter;
00370 loc_degree = loc_obj_par->loc_deg_poly;
00371 xsh_msg_dbg_medium( "Entering xsh_localize_obj_auto") ;
00372 xsh_msg_dbg_medium("Localize deg_poly %d chunk %d Thresh %f kappa %f niter %d",
00373 loc_degree, loc_obj_par->loc_chunk_nb, loc_obj_par->loc_thresh,
00374 kappa, niter);
00375 if ( loc_obj_par->use_skymask){
00376 XSH_ASSURE_NOT_NULL( skymask_frame);
00377 check( skymask_name = cpl_frame_get_filename( skymask_frame));
00378 xsh_msg_dbg_medium("Sky mask %s", skymask_name);
00379 }
00380
00381 level = xsh_debug_level_get();
00382
00383
00384 check( spectrum = xsh_spectrum_load( merge_frame, instrument));
00385 check( nslit = xsh_spectrum_get_size_slit( spectrum));
00386 check( nlambda = xsh_spectrum_get_size_lambda( spectrum));
00387 lambda_min = spectrum->lambda_min;
00388 lambda_step = spectrum->lambda_step;
00389 slit_min = spectrum->slit_min;
00390 slit_step = spectrum->slit_step;
00391
00392 check( flux = cpl_image_get_data_double( spectrum->flux));
00393 check( errs = cpl_image_get_data_double( spectrum->errs));
00394 check( qual = cpl_image_get_data_int( spectrum->qual));
00395
00396
00397 XSH_CALLOC( chunk_center, double, loc_obj_par->loc_chunk_nb);
00398 XSH_CALLOC( slit_center, double, loc_obj_par->loc_chunk_nb);
00399 XSH_CALLOC( slit_upper, double, loc_obj_par->loc_chunk_nb);
00400 XSH_CALLOC( slit_lower, double, loc_obj_par->loc_chunk_nb);
00401
00402
00403 XSH_CALLOC( sky_mask, int, nlambda);
00404
00405 if ( loc_obj_par->use_skymask){
00406 float *skymask_data = NULL;
00407 int irow, nrow;
00408 double fwhm =0.0, sky_min, sky_max;
00409 int isky_min, isky_max, imask;
00410 double width, resolution;
00411
00412 XSH_TABLE_LOAD( skymask_table, skymask_name);
00413
00414 check( xsh_sort_table_1( skymask_table, "WAVELENGTH", CPL_FALSE));
00415 check( skymask_data = cpl_table_get_data_float( skymask_table,
00416 "WAVELENGTH"));
00417 check( nrow = cpl_table_get_nrow( skymask_table));
00418
00419 xsh_msg_dbg_low("lambda min %f, step %f", lambda_min, lambda_step);
00420
00421 for( irow=0; irow < nrow; irow++){
00422 check( width = xsh_pfits_get_slit_width( spectrum->flux_header, instrument));
00423 resolution = xsh_resolution_get( instrument, width);
00424 fwhm = skymask_data[irow]/resolution;
00425 sky_min = skymask_data[irow]-fwhm;
00426 sky_max = skymask_data[irow]+fwhm;
00427 isky_min = (int)xsh_round_double((sky_min-lambda_min)/lambda_step);
00428 isky_max = (int)xsh_round_double((sky_max-lambda_min)/lambda_step);
00429 for( imask=isky_min; imask <=isky_max; imask++){
00430 sky_mask[imask] = 1;
00431 }
00432 }
00433 if (level >= XSH_DEBUG_LEVEL_MEDIUM){
00434 FILE *mask_file = NULL;
00435 char mask_name[256];
00436 int idbg=0;
00437
00438 sprintf( mask_name, "skymask.reg");
00439 mask_file = fopen( mask_name, "w");
00440
00441 fprintf( mask_file,"# Region file format: DS9 version 4.1\n");
00442 fprintf( mask_file,"global color=green dashlist=8 3 width=1 font=\"helvetica 10 normal\"\
00443 select=1 highlite=1 dash=0 fixed=0 edit=1 move=1 delete=1 include=1 source=1\n");
00444 fprintf(mask_file,"physical\n");
00445
00446 for(idbg=0; idbg< nlambda; idbg++){
00447
00448 if (sky_mask[idbg] == 1){
00449 fprintf( mask_file, "line(%d,%d,%d,%d)\n",idbg+1, nslit, idbg+1, 1);
00450 }
00451 }
00452 fclose( mask_file);
00453 }
00454 }
00455
00456 check( loc_list = xsh_localization_create());
00457
00458 #if 0
00459
00460 {
00461 float arcsec ;
00462 switch( xsh_instrument_get_arm( instrument ) ) {
00463 case XSH_ARM_NIR:
00464 arcsec = XSH_ARCSEC_NIR ;
00465 break ;
00466 case XSH_ARM_UVB:
00467 arcsec = XSH_ARCSEC_UVB ;
00468 break ;
00469 case XSH_ARM_VIS:
00470 arcsec = XSH_ARCSEC_VIS ;
00471 break ;
00472 default:
00473 arcsec = 0.15 ;
00474 break ;
00475 }
00476 slit_nod = loc_obj_par->nod_step/arcsec ;
00477 xsh_msg_dbg_medium( "nod_step: %lf, slit_nod: %d",
00478 loc_obj_par->nod_step, slit_nod ) ;
00479 }
00480 #endif
00481 XSH_CALLOC( slit_tab, double, nslit);
00482 check( gaussval = cpl_vector_wrap( nslit, slit_tab));
00483 check( gausspos = cpl_vector_new( nslit));
00484 for( islit=0; islit < nslit; islit++){
00485 cpl_vector_set( gausspos, islit, islit);
00486 }
00487
00488 chunk_size = (nlambda-1) / (double)(loc_obj_par->loc_chunk_nb);
00489
00490
00491
00492
00493 for( ichunk=0; ichunk< loc_obj_par->loc_chunk_nb; ichunk++){
00494
00495
00496
00497
00498
00499 int ifirst=0,ilast=0,icenter=0;
00500 double slit_cen=0, slit_low=0, slit_up=0;
00501 double threshold = 0.0;
00502 int is_valid_chunk = 1;
00503
00504 ifirst = ichunk*chunk_size;
00505 ilast = (ichunk+1)*chunk_size;
00506 icenter = (int)((ifirst+ilast)/2);
00507
00508 xsh_msg_dbg_medium("%d-%d", ifirst,ilast);
00509
00510
00511 chunk_coadd( slit_tab, flux, qual, nlambda, nslit, ifirst, ilast, sky_mask);
00512
00513 if (level >= XSH_DEBUG_LEVEL_MEDIUM){
00514 FILE *coadd_file = NULL;
00515 char coadd_name[256];
00516 int icoadd=0;
00517 const char* filename = cpl_frame_get_filename( merge_frame);
00518
00519 sprintf( coadd_name, "%s_%d_%d.dat",filename,ifirst,ilast);
00520 coadd_file = fopen( coadd_name, "w");
00521
00522 for(icoadd=0; icoadd<nslit; icoadd++){
00523 fprintf(coadd_file, "%d %f\n",icoadd,slit_tab[icoadd]);
00524 }
00525 fclose( coadd_file);
00526 }
00527
00528 if ( loc_obj_par->method == LOC_GAUSSIAN_METHOD){
00529 double cenpos=0, sigma=0,area=0, offset=0;
00530 double kappa=3.0;
00531
00532 xsh_msg_dbg_medium("Using GAUSSIAN method to fit chunk %d_%d", ifirst,ilast);
00533 cpl_vector_fit_gaussian( gausspos, NULL, gaussval, NULL, CPL_FIT_ALL, &cenpos,
00534 &sigma,
00535 &area, &offset, NULL,
00536 NULL, NULL);
00537 if( cpl_error_get_code() == CPL_ERROR_CONTINUE ){
00538 xsh_msg_dbg_low("CONTINUE to fit %d_%d x0 %f sigma %f", ifirst,ilast,cenpos,sigma);
00539 }
00540 if( cpl_error_get_code() == CPL_ERROR_NONE ){
00541 slit_cen = cenpos;
00542 slit_low = cenpos-kappa*sigma;
00543 slit_up = cenpos+kappa*sigma;
00544 if (slit_low < 0) slit_low = 0;
00545 if (slit_up >= nslit) slit_up = nslit-1;
00546 }
00547 else{
00548 xsh_msg_dbg_low("FAILED to fit %d_%d", ifirst,ilast);
00549 is_valid_chunk=0;
00550 xsh_error_reset();
00551 }
00552 }
00553 else
00554 {
00555 int ymax=0, ymin=0;
00556 int ylow=0, yup=0;
00557 tab_minmax_get_index( slit_tab, nslit, slitlimit_par, &ymin, &ymax);
00558 if (ymin == ymax){
00559 xsh_msg_warning( "No maximum found in chunk (skip it)");
00560 is_valid_chunk=0;
00561 }
00562 else{
00563
00564
00565 xsh_msg_dbg_medium("ymin [%d] = %f ymax [%d] = %f",ymin,
00566 slit_tab[ymin], ymax, slit_tab[ymax]);
00567
00568 threshold = slit_tab[ymin]+(slit_tab[ymax]-slit_tab[ymin])*
00569 loc_obj_par->loc_thresh;
00570 xsh_msg_dbg_medium("Threshold at %f",threshold);
00571
00572
00573 yup = ymax+1;
00574 while( yup < nslit && (slit_tab[yup] >= threshold)) {
00575 yup++;
00576 }
00577
00578
00579 ylow = ymax-1;
00580 while( ylow > 0 && (slit_tab[ylow] >= threshold)) {
00581 ylow--;
00582 }
00583
00584 slit_cen = ymax;
00585 slit_low = ylow;
00586 slit_up = yup;
00587 }
00588 }
00589
00590 if ( is_valid_chunk == 1){
00591 chunk_center[ichunk-skip_chunk] = lambda_min+lambda_step*icenter;
00592 slit_center[ichunk-skip_chunk] = slit_min+slit_step*(0.5+slit_cen);
00593 slit_lower[ichunk-skip_chunk] = slit_min+slit_step*(0.5+slit_low);
00594 slit_upper[ichunk-skip_chunk] = slit_min+slit_step*(0.5+slit_up);
00595 }
00596 else{
00597 skip_chunk++;
00598 }
00599 }
00600
00601 nb_chunk = loc_obj_par->loc_chunk_nb-skip_chunk;
00602
00603 XSH_ASSURE_NOT_ILLEGAL(loc_obj_par->loc_deg_poly < nb_chunk);
00604
00605 check( chunk_center_v = cpl_vector_wrap(
00606 nb_chunk, chunk_center));
00607 check( slit_center_v = cpl_vector_wrap(
00608 nb_chunk, slit_center));
00609
00610 check( cen_poly =
00611 xsh_polynomial_fit_1d_create( chunk_center_v, slit_center_v,
00612 loc_degree, NULL));
00613
00614
00615 if ( niter > 0){
00616 xsh_msg_dbg_low("Doing sigma clipping");
00617 iter = 0;
00618
00619 while ( (loc_degree < (nb_chunk-nbrej)) && iter < niter
00620 && nbrej != 0){
00621 double sigma_med;
00622
00623 nbrej = 0;
00624 xsh_msg_dbg_medium( " *** NBITER = %d / %d ***", iter+1, niter);
00625 dispslit = cpl_vector_new( nb_chunk);
00626
00627 for(ichunk = 0; ichunk < nb_chunk; ichunk++){
00628 double slit_fit;
00629 double lambda;
00630 double slit_diff;
00631
00632 lambda = chunk_center[ichunk];
00633
00634 check( slit_fit = cpl_polynomial_eval_1d( cen_poly,
00635 lambda, NULL));
00636
00637 slit_diff = slit_center[ichunk]-slit_fit;
00638 xsh_msg_dbg_low("slit_center %f FIT %f, DIFF %d %f", slit_center[ichunk], slit_fit, ichunk, slit_diff);
00639 check( cpl_vector_set( dispslit, ichunk, slit_diff));
00640 }
00641
00642 check( sigma_med = cpl_vector_get_stdev( dispslit));
00643 xsh_msg_dbg_medium(" kappa %f SIGMA MEDIAN = %f", kappa, sigma_med);
00644
00645 for(ichunk = 0; ichunk < nb_chunk; ichunk++){
00646 if ( fabs(cpl_vector_get( dispslit, ichunk)) > (kappa * sigma_med) ){
00647 nbrej++;
00648 }
00649 else{
00650 chunk_center[ichunk-nbrej] = chunk_center[ichunk];
00651 slit_center[ichunk-nbrej] = slit_center[ichunk];
00652 slit_lower[ichunk-nbrej] = slit_lower[ichunk];
00653 slit_upper[ichunk-nbrej] = slit_upper[ichunk];
00654 }
00655 }
00656
00657 xsh_msg_dbg_medium(" Removed %d points", nbrej);
00658
00659 nb_chunk -= nbrej;
00660
00661
00662 xsh_unwrap_vector( &chunk_center_v);
00663 xsh_unwrap_vector( &slit_center_v);
00664 xsh_free_polynomial( &cen_poly);
00665
00666 check( chunk_center_v = cpl_vector_wrap(
00667 nb_chunk, chunk_center));
00668 check( slit_center_v = cpl_vector_wrap(
00669 nb_chunk, slit_center));
00670
00671 check( cen_poly =
00672 xsh_polynomial_fit_1d_create( chunk_center_v, slit_center_v,
00673 loc_degree, NULL));
00674 xsh_free_vector( &dispslit);
00675 iter++;
00676 }
00677 }
00678
00679 if (level >= XSH_DEBUG_LEVEL_MEDIUM){
00680 FILE *debug_file = NULL;
00681 char debug_name[256];
00682 int idebug=0;
00683 const char* filename = cpl_frame_get_filename( merge_frame);
00684
00685 sprintf( debug_name, "%s_points.dat",filename);
00686 debug_file = fopen( debug_name, "w");
00687
00688 fprintf( debug_file,"#chunk_pos slit_low slit_cen slit_up\n");
00689
00690 for(idebug=0; idebug<nb_chunk; idebug++){
00691 fprintf( debug_file, "%f %f %f %f\n",chunk_center[idebug],
00692 slit_lower[idebug], slit_center[idebug], slit_upper[idebug]);
00693 }
00694 fclose( debug_file);
00695 }
00696
00697
00698 check( slit_lower_v = cpl_vector_wrap(
00699 nb_chunk, slit_lower));
00700 check( slit_upper_v = cpl_vector_wrap(
00701 nb_chunk, slit_upper));
00702 loc_list->pol_degree = loc_degree;
00703
00704 check(loc_list->cenpoly = cpl_polynomial_duplicate( cen_poly));
00705 check(loc_list->edglopoly =
00706 xsh_polynomial_fit_1d_create( chunk_center_v, slit_lower_v,
00707 loc_degree, NULL));
00708 check(loc_list->edguppoly =
00709 xsh_polynomial_fit_1d_create( chunk_center_v, slit_upper_v,
00710 loc_degree, NULL));
00711
00712 cleanup:
00713 xsh_unwrap_vector( &chunk_center_v);
00714 xsh_unwrap_vector( &slit_center_v);
00715 xsh_unwrap_vector( &slit_upper_v);
00716 xsh_unwrap_vector( &slit_lower_v);
00717 xsh_unwrap_vector( &gaussval);
00718 xsh_free_vector( &gausspos);
00719 XSH_FREE( chunk_center);
00720 XSH_FREE( slit_center);
00721 XSH_FREE( slit_upper);
00722 XSH_FREE( slit_lower);
00723 XSH_FREE( slit_tab);
00724 XSH_FREE( sky_mask);
00725 xsh_free_table( &skymask_table);
00726 xsh_free_polynomial( &cen_poly);
00727 xsh_spectrum_free( &spectrum);
00728 return loc_list;
00729
00730 }
00731
00732 static int comp_lambda( const void * one, const void * two )
00733 {
00734 float * un = (float *)one, * deux = (float *)two ;
00735
00736 if ( *un < *deux ) return -1 ;
00737 else if ( *un == *deux ) return 0 ;
00738 else return 1 ;
00739 }
00740
00741 static cpl_frame * concat_rec( cpl_frame * rect_frame,
00742 xsh_instrument * instrument,
00743 int slitlet, int * nb_orders )
00744 {
00745 cpl_frame * result = NULL ;
00746 xsh_rec_list * res_list = NULL, * rect_list = NULL ;
00747 int norders, i, ns, nlambda, depth ;
00748 float * res_slit = NULL, * rect_slit = NULL ;
00749 double * res_lambda = NULL, * rect_lambda = NULL, * plambda ;
00750 const char * tag = "ORDER2D_DOWN_IFU_VIS" ;
00751 char * fname = NULL ;
00752 float * res_data = NULL, * res_errs = NULL ;
00753 int * res_qual = NULL ;
00754 float * pdata, * perrs, * pdata1, * perrs1 ;
00755 int * pqual, * pqual1 ;
00756 float * rect_data = NULL, * rect_errs = NULL ;
00757 int * rect_qual = NULL ;
00758
00759 XSH_ASSURE_NOT_NULL( rect_frame ) ;
00760
00761 check( rect_list = xsh_rec_list_load( rect_frame, instrument ) ) ;
00762
00763 check( res_list = xsh_rec_list_create_with_size( 1, instrument ) ) ;
00764
00765 ns = xsh_rec_list_get_nslit( rect_list, 0);
00766 norders = rect_list->size ;
00767 *nb_orders = norders ;
00768 for( i = 0, nlambda = 0 ; i<norders ; i++ ) {
00769 nlambda += xsh_rec_list_get_nlambda( rect_list, i);
00770 xsh_msg_dbg_medium( "Nlambda[%d] = %d, sum = %d", ns,
00771 rect_list->list[i].nlambda, nlambda);
00772 }
00773 xsh_msg_dbg_medium( "Single Rectify list: ns = %d, Nlambda = %d",
00774 ns, nlambda ) ;
00775
00776
00777 check( xsh_rec_list_set_data_size( res_list, 0, 16, nlambda, ns ) ) ;
00778
00779 res_slit = xsh_rec_list_get_slit( res_list, 0 ) ;
00780 rect_slit = xsh_rec_list_get_slit( rect_list, 0 ) ;
00781 memcpy( res_slit, rect_slit, ns*sizeof( float ) ) ;
00782
00783
00784 res_lambda = xsh_rec_list_get_lambda( res_list, 0 ) ;
00785 plambda = res_lambda ;
00786 for( i = 0 ; i<norders ; i++ ) {
00787 int nl ;
00788
00789 rect_lambda = xsh_rec_list_get_lambda( rect_list, i ) ;
00790 nl = xsh_rec_list_get_nlambda( rect_list, i) ;
00791 memcpy( plambda, rect_lambda, nl*sizeof(double) ) ;;
00792 plambda += nl ;
00793 }
00794
00795
00796 depth = ns*nlambda ;
00797 res_data = xsh_rec_list_get_data1( res_list, 0 ) ;
00798 pdata = res_data ;
00799 pdata1 = res_data ;
00800 res_errs = xsh_rec_list_get_errs1( res_list, 0 ) ;
00801 perrs = res_errs ;
00802 perrs1 = res_errs ;
00803 res_qual = xsh_rec_list_get_qual1( res_list, 0 ) ;
00804 pqual = res_qual ;
00805 pqual1 = pqual ;
00806 for( i = 0 ; i<norders ; i++ ) {
00807 int nl ;
00808 float * rdata, * pdata2 ;
00809 float * rerrs, * perrs2 ;
00810 int * rqual, * pqual2 ;
00811 int ks, kl ;
00812
00813 check( rect_data = xsh_rec_list_get_data1( rect_list, i ) ) ;
00814 check( rect_errs = xsh_rec_list_get_errs1( rect_list, i ) ) ;
00815 check( rect_qual = xsh_rec_list_get_qual1( rect_list, i ) ) ;
00816 nl = xsh_rec_list_get_nlambda( rect_list, i) ;
00817 xsh_msg_dbg_medium( "Order %d, nlambdas: %d [%d]", i, nl, nlambda ) ;
00818 rdata = rect_data ;
00819 pdata2 = pdata1 ;
00820 rerrs = rect_errs ;
00821 perrs2 = perrs1 ;
00822 rqual = rect_qual ;
00823 pqual2 = pqual1 ;
00824 for ( ks = 0 ; ks < ns ; ks++ ) {
00825 for ( kl = 0 ; kl < nl ; kl++ ) {
00826 *pdata++ = *rdata++ ;
00827 *perrs++ = *rerrs++ ;
00828 *pqual++ = *rqual++ ;
00829 }
00830 pdata2 += nlambda ;
00831 pdata = pdata2 ;
00832 perrs2 += nlambda ;
00833 perrs = perrs2 ;
00834 pqual2 += nlambda ;
00835 pqual = pqual2 ;
00836 }
00837 pdata1 += nl ;
00838 pdata = pdata1 ;
00839 perrs1 += nl ;
00840 perrs = perrs1 ;
00841 pqual1 += nl ;
00842 pqual = pqual1 ;
00843 }
00844
00845
00846 {
00847 int * idx ;
00848 int ks ;
00849
00850 idx = xsh_sort( res_lambda, nlambda, sizeof( double ), comp_lambda ) ;
00851 pdata = res_data ;
00852 perrs = res_errs ;
00853 pqual = res_qual ;
00854 for( ks = 0 ; ks < ns ; ks++ ) {
00855 xsh_reindex_float( pdata, idx, nlambda ) ;
00856 pdata += nlambda ;
00857 xsh_reindex_float( perrs, idx, nlambda ) ;
00858 perrs += nlambda ;
00859 xsh_reindex_int( pqual, idx, nlambda ) ;
00860 pqual += nlambda ;
00861 }
00862 }
00863
00864 fname = xsh_stringcat_any( "SINGLE_RECTIFY_", SlitletName[slitlet],
00865 ".fits", NULL ) ;
00866 check( result = xsh_rec_list_save( res_list, fname, tag, 1 ) ) ;
00867 XSH_FREE( fname ) ;
00868 fname = xsh_stringcat_any( "MULTIPLE_RECTIFY_", SlitletName[slitlet],
00869 ".fits", NULL ) ;
00870 check( rect_frame = xsh_rec_list_save( rect_list,fname, tag, 1 ) ) ;
00871 cleanup:
00872 xsh_rec_list_free( &res_list ) ;
00873 xsh_rec_list_free( &rect_list ) ;
00874 return result ;
00875 }
00876
00877
00878
00879 cpl_frameset * xsh_localize_obj_ifu( cpl_frameset *rec_frameset,
00880 cpl_frame *skymask_frame,
00881 xsh_instrument * instrument,
00882 xsh_localize_obj_param * locobj_par,
00883 xsh_slit_limit_param *slitlimit_par)
00884 {
00885 int i, slitlet;
00886 cpl_frameset *result_frameset = NULL;
00887 char fname[256];
00888
00889 XSH_ASSURE_NOT_NULL( rec_frameset);
00890 XSH_ASSURE_NOT_NULL( instrument);
00891 XSH_ASSURE_NOT_NULL( locobj_par);
00892
00893 check( result_frameset = cpl_frameset_new());
00894
00895 for( i = 0, slitlet = LOWER_IFU_SLITLET ; i < 3 ; i++, slitlet++ ) {
00896 cpl_frame * loc_frame = NULL;
00897 cpl_frame *rec_frame = NULL;
00898
00899 sprintf( fname ,"LOCALIZE_TABLE_%s_IFU_%s.fits", SlitletName[slitlet],
00900 xsh_instrument_arm_tostring( instrument));
00901
00902 xsh_msg( "Localizing slitlet %s, frame '%s'", SlitletName[slitlet], fname);
00903
00904 check( rec_frame = cpl_frameset_get_frame( rec_frameset, i));
00905
00906 check( loc_frame = xsh_localize_obj( rec_frame, skymask_frame, instrument,
00907 locobj_par,
00908 slitlimit_par, fname));
00909 check( cpl_frameset_insert( result_frameset, loc_frame));
00910 }
00911
00912 cleanup:
00913 if ( cpl_error_get_code() != CPL_ERROR_NONE ) {
00914 xsh_free_frameset( &result_frameset);
00915 }
00916 return result_frameset;
00917 }
00918
00919