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
00031
00038
00041
00042
00043
00044
00045 #include <math.h>
00046 #include <xsh_drl.h>
00047 #include <xsh_drl_check.h>
00048
00049 #include <xsh_badpixelmap.h>
00050 #include <xsh_data_pre.h>
00051 #include <xsh_data_order.h>
00052 #include <xsh_data_wavemap.h>
00053 #include <xsh_data_localization.h>
00054 #include <xsh_data_rec.h>
00055 #include <xsh_dfs.h>
00056 #include <xsh_pfits.h>
00057 #include <xsh_error.h>
00058 #include <xsh_utils_image.h>
00059 #include <xsh_msg.h>
00060 #include <xsh_fit.h>
00061
00062 #include <cpl.h>
00063
00064 #include <math.h>
00065 #include <gsl/gsl_bspline.h>
00066 #include <gsl/gsl_linalg.h>
00067 #ifdef DARWIN
00068 #include <vecLib/clapack.h>
00069 #else
00070 #include <f2c.h>
00071 #include <clapack.h>
00072 #endif
00073
00074 #define REGDEBUG_MEDIAN_SPLINE 0
00075
00076
00077
00078 static void fit_spline( xsh_wavemap_list * wave_list, int idx,
00079 int nbkpts, int order, int niter,
00080 float kappa, float ron, float gain,
00081 const char * sarm );
00082 static xsh_wavemap_list* xsh_wavemap_list_new( cpl_image *wavemap,
00083 cpl_image *slitmap, xsh_localization *list,
00084 xsh_order_list *order_table, xsh_pre *pre_sci, int nbkpts,cpl_frame* usr_defined_break_points_frame,
00085 xsh_subtract_sky_single_param* sky_par, xsh_instrument *instrument);
00086 static xsh_rec_list* xsh_wavemap_list_build_sky( xsh_wavemap_list* list,
00087 xsh_pre* pre_mflat,
00088 xsh_instrument *instrument);
00089
00090 static int wave_compare( const void * un, const void * deux ) ;
00091
00092
00093
00094
00095
00096
00097
00098 static cpl_frame* xsh_wavelist_subtract_sky( xsh_pre * pre_sci,
00099 xsh_pre* pre_mflat,
00100 xsh_rec_list* sky_list,
00101 xsh_wavemap_list * wave_list,
00102 xsh_instrument* instrument,
00103 const char* rec_prefix);
00104
00105 static int xsh_bspline_eval_all(const double x, gsl_vector *B, size_t *idx,
00106 gsl_bspline_workspace *w, size_t *start_i);
00107
00108 static size_t xsh_bspline_find_interval(const double x, int *flag,
00109 gsl_bspline_workspace *w, size_t *start_i);
00110
00111
00112
00113
00114
00115
00116
00117
00118 void
00119 fit_spline( xsh_wavemap_list * wave_list, int idx,
00120 int nbkpts, int bs_order, int niter,
00121 float kappa, float ron, float gain,
00122 const char * sarm )
00123 {
00124 float sqrtgain= sqrtf(gain);
00125
00126 gsl_bspline_workspace *bw;
00127 gsl_matrix *Bn, *Bn_full;
00128 gsl_vector *Bkpts;
00129 int *Bidx,*Bidx_full;
00130 float *y,*yf;
00131 double *lambda_fit;
00132
00133 double* pw=NULL;
00134 double* pf=NULL;
00135 double* ps=NULL;
00136 double* px=NULL;
00137 double* py=NULL;
00138 double* pfit=NULL;
00139 double* perr=NULL;
00140
00141 double start, end ;
00142 int ncoeffs ;
00143 int i, j ;
00144 int ii,jj,kk;
00145 int iii;
00146 int order, nitem ;
00147 wavemap_item * pentry ;
00148
00149 float *Xtp,*c;
00150
00151
00152 size_t idxb=0;
00153 size_t start_i;
00154
00155 int ord=bs_order;
00156 double *bwB_ptr;
00157 double *Bkpts_ptr;
00158
00159 double *Bn_ptr;
00160 int Bn_strd;
00161
00162 double *Bn_full_ptr;
00163 int Bn_full_strd;
00164 int level;
00165 FILE* fout = NULL;
00166 cpl_table* debug_fit=NULL;
00167
00168 char fname[128], dirname[128], cmd[128];
00169
00170
00171
00172 pentry = wave_list->list[idx].sky;
00173 XSH_ASSURE_NOT_NULL( pentry ) ;
00174 order = wave_list->list[idx].order ;
00175 nitem = wave_list->list[idx].sky_size;
00176 start = wave_list->list[idx].lambda_min ;
00177 end = wave_list->list[idx].lambda_max ;
00178
00179 xsh_msg_dbg_medium( "Fit Spline(%d), Order %d, bs_order: %d, ndata: %d, nbkpts: %d",
00180 idx, order, bs_order, nitem, nbkpts ) ;
00181 xsh_msg_dbg_medium( " Lambda Min: %lf, Max: %lf", start, end ) ;
00182
00183 bw = gsl_bspline_alloc(ord, nbkpts);
00184 ncoeffs = gsl_bspline_ncoeffs(bw) ;
00185
00186 Bkpts = gsl_vector_alloc(nbkpts);
00187 Bn = gsl_matrix_alloc(ord,nitem);
00188 Bidx = calloc(nitem,sizeof(int));
00189
00190 Bn_full = gsl_matrix_alloc(ord,nitem);
00191 Bidx_full = calloc(nitem,sizeof(int));
00192
00193 Xtp=calloc(ncoeffs*ord,sizeof(float));
00194 c=calloc(ncoeffs,sizeof(float));
00195
00196 y=calloc(nitem,sizeof(float));
00197 yf=calloc(nitem,sizeof(float));
00198
00199 lambda_fit=calloc(nitem,sizeof(double));
00200
00201 bwB_ptr=gsl_vector_ptr(bw->B,0);
00202 Bkpts_ptr=gsl_vector_ptr(Bkpts,0);
00203
00204 Bn_ptr=gsl_matrix_ptr(Bn,0,0);
00205 Bn_strd=Bn->tda;
00206
00207 Bn_full_ptr=gsl_matrix_ptr(Bn_full,0,0);
00208 Bn_full_strd=Bn_full->tda;
00209
00210 level = xsh_debug_level_get() ;
00211
00212
00213 debug_fit=cpl_table_new(nitem);
00214
00215 cpl_table_new_column(debug_fit,"WAVE",CPL_TYPE_DOUBLE);
00216 cpl_table_new_column(debug_fit,"FLUX",CPL_TYPE_DOUBLE);
00217 cpl_table_new_column(debug_fit,"SIGMA",CPL_TYPE_DOUBLE);
00218 cpl_table_new_column(debug_fit,"FIT",CPL_TYPE_DOUBLE);
00219 cpl_table_new_column(debug_fit,"ERR",CPL_TYPE_DOUBLE);
00220 cpl_table_new_column(debug_fit,"X",CPL_TYPE_DOUBLE);
00221 cpl_table_new_column(debug_fit,"Y",CPL_TYPE_DOUBLE);
00222
00223 cpl_table_fill_column_window_double(debug_fit,"WAVE",0,nitem,0.);
00224 cpl_table_fill_column_window_double(debug_fit,"FLUX",0,nitem,0.);
00225 cpl_table_fill_column_window_double(debug_fit,"SIGMA",0,nitem,0.);
00226 cpl_table_fill_column_window_double(debug_fit,"FIT",0,nitem,0.);
00227 cpl_table_fill_column_window_double(debug_fit,"ERR",0,nitem,0.);
00228 cpl_table_fill_column_window_double(debug_fit,"X",0,nitem,0.);
00229 cpl_table_fill_column_window_double(debug_fit,"Y",0,nitem,0.);
00230
00231 pw=cpl_table_get_data_double(debug_fit,"WAVE");
00232 pf=cpl_table_get_data_double(debug_fit,"FLUX");
00233 ps=cpl_table_get_data_double(debug_fit,"SIGMA");
00234 pfit=cpl_table_get_data_double(debug_fit,"FIT");
00235 perr=cpl_table_get_data_double(debug_fit,"ERR");
00236 px=cpl_table_get_data_double(debug_fit,"X");
00237 py=cpl_table_get_data_double(debug_fit,"Y");
00238
00239
00240
00241 for (iii=0;iii<niter;iii++){
00242 int nfit=0, err;
00243 char uplo='U';
00244
00245
00246
00247 #ifdef DARWIN
00248 __CLPK_integer clpk_n=ncoeffs;
00249 __CLPK_integer kd=ord-1;
00250 __CLPK_integer nrhs=1;
00251 __CLPK_integer ldab=ord;
00252 __CLPK_integer ldb=ncoeffs;
00253 __CLPK_integer info=0;
00254 #else
00255 integer clpk_n=ncoeffs;
00256 integer kd=ord-1;
00257 integer nrhs=1;
00258 integer ldab=ord;
00259 integer ldb=ncoeffs;
00260 integer info=0;
00261 #endif
00262
00263
00264
00265 pentry = wave_list->list[idx].sky;
00266
00267 for(i=0;i<ncoeffs*ord;i++)
00268 Xtp[i]=0;
00269 for(i=0;i<ncoeffs;i++)
00270 c[i]=0;
00271 if(iii==0)
00272 gsl_bspline_knots_uniform(start, end, bw);
00273 else{
00274 for ( i = 0 ; i<nitem ; i++, pentry++ ) {
00275 if(fabs(yf[i]-pentry->flux) < kappa*(ron*ron+sqrtgain*sqrt(fabs(yf[i])))){
00276 lambda_fit[nfit]=pentry->lambda;
00277 nfit++;
00278 }
00279 }
00280 for(i=0;i<nbkpts-1;i++)
00281 Bkpts_ptr[i]=lambda_fit[i*(nfit/nbkpts)];
00282 Bkpts_ptr[0]=start;
00283 Bkpts_ptr[nbkpts-1]=end;
00284 gsl_bspline_knots(Bkpts, bw);
00285 }
00286
00287 pentry = wave_list->list[idx].sky;
00288 nfit=0;
00289 start_i=bw->k - 1;
00290
00291
00292 for ( i = 0 ; i<nitem ; i++, pentry++ ) {
00293 if ( isnan( (double)pentry->lambda ) || isnan( (double)pentry->flux ) ||
00294 isnan( (double)pentry->sigma ) ) {
00295 xsh_msg( "Nan in lambda, flux or sigma" ) ;
00296 }
00297
00298 xsh_bspline_eval_all( (double)pentry->lambda, bw->B, &idxb,bw, &start_i);
00299 Bidx_full[i]=idxb;
00300 for (j = 0; j < ord; ++j)
00301 Bn_full_ptr[j*Bn_full_strd+i] = bwB_ptr[j];
00302
00303 if(iii==0 || fabs(yf[i]-pentry->flux) < kappa*(ron*ron+sqrtgain*sqrt(fabs(yf[i])))){
00304 y[nfit]=(float)pentry->flux;
00305 Bidx[nfit]=idxb;
00306 for (j = 0; j < ord; ++j)
00307 Bn_ptr[j*Bn_strd+nfit] = bwB_ptr[j];
00308 nfit++;
00309 }
00310 }
00311 if (nfit > 0) {
00312 xsh_msg_debug("niter %d nfit: %d %f %f", iii, nfit,y[0],y[nfit-1]);
00313 }
00314 else{
00315 xsh_msg_debug("niter %d nfit: %d ", iii, nfit);
00316 }
00317
00318
00319
00320 for(ii=0;ii<nfit;ii++){
00321 for(jj=0;jj<ord;jj++){
00322 for(kk=jj;kk<ord;kk++){
00323 idxb=Bidx[ii]-(ord-1);
00324 Xtp[(idxb+kk)*ord+(ord-1)+(idxb+jj)-(idxb+kk)]+=Bn_ptr[jj*Bn_strd+ii]*Bn_ptr[kk*Bn_strd+ii];
00325 }
00326 }
00327 }
00328
00329 for(ii=0;ii<nfit;ii++){
00330 for(jj=0;jj<ord;jj++){
00331 idxb=Bidx[ii]-(ord-1);
00332 c[idxb+jj]+=Bn_ptr[jj*Bn_strd+ii]*y[ii];
00333
00334 }
00335 }
00336
00337
00338 err=spbtrf_(&uplo,&clpk_n,&kd,Xtp,&ldab,&info);
00339 err=spbtrs_(&uplo,&clpk_n,&kd,&nrhs,Xtp,&ldab,c,&ldb,&info);
00340
00341
00342
00343 for(ii=0;ii<nitem;ii++)
00344 yf[ii]=0;
00345
00346 for(ii=0;ii<nitem;ii++){
00347 for(jj=0;jj<ord;jj++){
00348 idxb=Bidx_full[ii]-(ord-1);
00349
00350
00351
00352 if(!isnan(c[idxb+jj])) {
00353 yf[ii]+=Bn_full_ptr[jj*Bn_full_strd+ii]*c[idxb+jj];
00354 }
00355
00356
00357
00358
00359
00360
00361
00362
00363
00364 }
00365 }
00366 xsh_msg_debug("end of iter");
00367 }
00368
00369
00370
00371 if ( level >= XSH_DEBUG_LEVEL_LOW ) {
00372
00373 xsh_msg_dbg_medium( "Output fitted data" ) ;
00374 sprintf( dirname, "Wave_%s", sarm ) ;
00375 if ( access( dirname, 0 ) != 0 ) {
00376 sprintf( cmd, "mkdir %s", dirname);
00377 system( cmd);
00378 }
00379 sprintf( fname, "%s/wave_%02d.dat", dirname, order ) ;
00380 fout = fopen( fname, "w" ) ;
00381 }
00382
00383 pentry = wave_list->list[idx].sky;
00384
00385 for (i = 0; i < nitem ; i++, pentry++ ) {
00386 double xi, yi=0, yerr=0;
00387
00388 xi = pentry->lambda ;
00389 yi=yf[i];
00390
00391 pentry->fitted = yi ;
00392 pentry->fit_err = yerr ;
00393
00394 if ( level >= XSH_DEBUG_LEVEL_LOW ) {
00395 fprintf( fout, "%f %f %f %d %d ", pentry->lambda, pentry->flux,
00396 pentry->sigma, pentry->ix, pentry->iy );
00397 fprintf( fout, "%lf %lf\n", yi, yerr);
00398 }
00399
00400
00401 pw[i]=pentry->lambda;
00402 pf[i]=pentry->flux;
00403 ps[i]=pentry->sigma;
00404
00405 pfit[i]=pentry->fitted;
00406 perr[i]=pentry->fit_err;
00407
00408 px[i]=pentry->ix;
00409 py[i]=pentry->iy;
00410
00411
00412 }
00413 if (level >= XSH_DEBUG_LEVEL_HIGH) {
00414 sprintf(fname,"debug_fit%02d.fits",order);
00415 check(cpl_table_save(debug_fit,NULL,NULL,fname,CPL_IO_DEFAULT));
00416 }
00417 xsh_free_table(&debug_fit);
00418
00419
00420 if ( fout ) fclose( fout ) ;
00421
00422 free(Xtp);
00423 free(c);
00424 free(y);
00425 free(yf);
00426 free(Bidx);
00427 free(Bidx_full);
00428 free(lambda_fit);
00429
00430 gsl_bspline_free(bw);
00431 gsl_vector_free( Bkpts);
00432 gsl_matrix_free(Bn);
00433 gsl_matrix_free(Bn_full);
00434 cleanup:
00435 xsh_free_table( &debug_fit);
00436 return ;
00437 }
00438
00439
00440
00441
00442 static int wave_compare( const void * un, const void * deux )
00443 {
00444 wavemap_item * one = (wavemap_item *)un ;
00445 wavemap_item * two = (wavemap_item *)deux ;
00446
00447 if ( one->lambda < two->lambda ) return -1 ;
00448 else if ( one->lambda == two->lambda ) return 0 ;
00449 return 1 ;
00450 }
00451
00452
00453
00458
00459 double xsh_nbkpts_uvb[XSH_ORDERS_UVB]={2.0,2.0,2.0,2.0,
00460 2.0,2.0,2.0,2.0,
00461 2.0,2.0,2.0,2.0};
00462
00463
00464
00465
00466 double xsh_nbkpts_vis[XSH_ORDERS_VIS]={1.5,0.8,0.8,0.7,
00467 0.8,0.8,0.7,0.7,
00468 0.7,0.75,0.7,0.6,
00469 0.6,0.6,0.6};
00470
00471
00472
00473
00474
00475
00476
00477
00478
00479
00480
00481
00482 double xsh_nbkpts_nir[XSH_ORDERS_NIR]={0.30,0.28,0.28,0.7,
00483 0.9,0.8,0.9,0.9,
00484 0.5,0.7,1.0,0.5,
00485 0.7,0.7,0.7,0.4};
00486
00487
00488
00489
00490
00491
00492 static xsh_wavemap_list* xsh_wavemap_list_new( cpl_image *wavemap,
00493 cpl_image *slitmap, xsh_localization *list,
00494 xsh_order_list *order_table, xsh_pre *pre_sci, int nbkpts,cpl_frame* usr_defined_break_points_frame,
00495 xsh_subtract_sky_single_param *sky_par, xsh_instrument *instrument)
00496 {
00497 cpl_image *sci = NULL;
00498 xsh_wavemap_list *wave_list = NULL;
00499 float *pflux = NULL ;
00500 float *perrs = NULL ;
00501 int *pqual = NULL ;
00502 double *plambda = NULL ;
00503 double *pslit = NULL;
00504 int i, ix, iy, nx, ny, max_size, sky_size, object_size;
00505 int iorder ;
00506 int level =0;
00507 int nrow=0;
00508 const double* pfactor=NULL;
00509 float slit_edges_mask = 0.5;
00510
00511 int median_hsize = 7;
00512 double nbkpts_arm=0;
00513 int nbkpts_ord=0;
00514 const char * name=NULL;
00515 cpl_table* usr_defined_break_points_tab=NULL;
00516 int wmap_xsize_diff=0;
00517 double obj_slit_min=0, obj_slit_max=0;
00518 double pos1=0, hheight1=0;
00519 double pos2=0, hheight2=0;
00520 double sky_slit_min=0, sky_slit_max=0;
00521
00522
00523 XSH_ASSURE_NOT_NULL( wavemap);
00524 XSH_ASSURE_NOT_NULL( slitmap);
00525 XSH_ASSURE_NOT_NULL( order_table);
00526 XSH_ASSURE_NOT_NULL( pre_sci);
00527 XSH_ASSURE_NOT_NULL( sky_par);
00528 XSH_ASSURE_NOT_NULL( instrument);
00529
00530
00531 median_hsize = sky_par->median_hsize;
00532 slit_edges_mask = sky_par->slit_edges_mask;
00533 pos1 = sky_par->pos1;
00534 hheight1 = sky_par->hheight1;
00535
00536 pos2 = sky_par->pos2;
00537 hheight2 = sky_par->hheight2;
00538
00539
00540 check( sci = xsh_pre_get_data( pre_sci));
00541 nx = xsh_pre_get_nx( pre_sci);
00542 ny = xsh_pre_get_ny( pre_sci);
00543
00544
00545 check( pflux = cpl_image_get_data_float( xsh_pre_get_data( pre_sci)));
00546 check( perrs = cpl_image_get_data_float( xsh_pre_get_errs( pre_sci)));
00547 check( pqual = cpl_image_get_data_int( xsh_pre_get_qual( pre_sci)));
00548 check( plambda = cpl_image_get_data_double( wavemap));
00549 check( pslit = cpl_image_get_data_double( slitmap));
00550
00551
00552 if( xsh_instrument_get_arm(instrument) == XSH_ARM_NIR){
00553 wmap_xsize_diff=pre_sci->nx-cpl_image_get_size_x(slitmap);
00554 }
00555
00556
00557
00558 check( wave_list = xsh_wavemap_list_create( instrument));
00559
00560
00561
00562 if ( (hheight1 == 0) && (hheight2 == 0) && (list != NULL) ){
00563 obj_slit_min = cpl_polynomial_eval_1d( list->edglopoly,
00564 600, NULL);
00565 obj_slit_max = cpl_polynomial_eval_1d( list->edguppoly,
00566 600, NULL);
00567 xsh_msg_dbg_medium("Limit in slit given by localization %f => %f",
00568 obj_slit_min,obj_slit_max);
00569 }
00570 else {
00571 if (hheight1 > 0){
00572 sky_slit_min = pos1 - hheight1;
00573 sky_slit_max = pos1 + hheight1;
00574
00575 if (hheight2 > 0){
00576 double s2_min, s2_max, s1_min, s1_max;
00577
00578 s1_min = pos1 - hheight1;
00579 s1_max = pos1 + hheight1;
00580
00581 s2_min = pos2 - hheight2;
00582 s2_max = pos2 + hheight2;
00583
00584 if ( s2_min < s1_min){
00585 sky_slit_min = s2_min;
00586 }
00587 else{
00588 sky_slit_min = s1_min;
00589 }
00590
00591 if ( s2_max > s1_max){
00592 sky_slit_max = s2_max;
00593 }
00594 else{
00595 sky_slit_max = s1_max;
00596 }
00597
00598 if ( s2_min > s1_max){
00599 obj_slit_min = s1_max;
00600 obj_slit_max = s2_min;
00601 }
00602
00603 if ( s1_min > s2_max){
00604 obj_slit_min = s2_max;
00605 obj_slit_max = s1_min;
00606 }
00607 }
00608 }
00609 else{
00610 sky_slit_min = pos2 - hheight2;
00611 sky_slit_max = pos2 + hheight2;
00612 }
00613 }
00614
00615
00616 for ( iorder = 0; iorder < order_table->size; iorder++){
00617 int miny, maxy, minx, maxx , diff_minmax_x;
00618 double lambda ;
00619 wavemap_item *psky = NULL;
00620 wavemap_item *pobject = NULL;
00621 int order;
00622 double lambda_min = 10000., lambda_max = 0. ;
00623 double sky_slit_min_edge=0, sky_slit_max_edge=0;
00624
00625
00626 miny = xsh_order_list_get_starty( order_table, iorder);
00627 maxy = xsh_order_list_get_endy( order_table, iorder);
00628 order = order_table->list[iorder].absorder ;
00629 xsh_msg_dbg_medium( "Order %d, miny %d, maxy %d", order, miny, maxy ) ;
00630
00631
00632
00633 max_size = 0;
00634
00635 for( iy = miny ; iy < maxy ; iy++){
00636 float slit_val;
00637
00638 check( minx = xsh_order_list_eval_int( order_table,
00639 order_table->list[iorder].edglopoly, iy)+1);
00640 check( maxx = xsh_order_list_eval_int( order_table,
00641 order_table->list[iorder].edguppoly, iy)-1);
00642
00643 slit_val = pslit[minx+iy*(nx-wmap_xsize_diff)];
00644
00645 if (slit_val > sky_slit_max_edge){
00646 sky_slit_max_edge = slit_val;
00647 }
00648 if (slit_val < sky_slit_min_edge){
00649 sky_slit_min_edge = slit_val;
00650 }
00651
00652 slit_val = pslit[maxx+iy*(nx-wmap_xsize_diff)];
00653
00654 if (slit_val > sky_slit_max_edge){
00655 sky_slit_max_edge = slit_val;
00656 }
00657 if (slit_val < sky_slit_min_edge){
00658 sky_slit_min_edge = slit_val;
00659 }
00660
00661 diff_minmax_x = maxx - minx;
00662 max_size += diff_minmax_x+1;
00663 }
00664
00665 XSH_ASSURE_NOT_ILLEGAL( sky_slit_min_edge <= sky_slit_max_edge);
00666 XSH_CMP_INT( max_size, >, 0, "Not enough points for the sky order %d",
00667 ,order);
00668 sky_slit_min_edge = sky_slit_min_edge+slit_edges_mask;
00669 sky_slit_max_edge = sky_slit_max_edge-slit_edges_mask;
00670
00671 if ( (hheight1 == 0) && (hheight1 == 0)){
00672 sky_slit_min = sky_slit_min_edge;
00673 sky_slit_max = sky_slit_max_edge;
00674 }
00675 else{
00676 if (sky_slit_min_edge > sky_slit_min){
00677 sky_slit_min = sky_slit_min_edge;
00678 }
00679 if (sky_slit_max_edge < sky_slit_max){
00680 sky_slit_max = sky_slit_max_edge;
00681 }
00682 }
00683 xsh_msg_dbg_medium( "SKY data slit range : [%f,%f],[%f,%f]",
00684 sky_slit_min, obj_slit_min, obj_slit_max, sky_slit_max);
00685
00686
00687 check( xsh_wavemap_list_set_max_size( wave_list, iorder, order,
00688 max_size));
00689
00690
00691 psky = wave_list->list[iorder].sky;
00692 pobject = wave_list->list[iorder].object;
00693 sky_size = 0;
00694 object_size = 0;
00695
00696 for( iy = miny ; iy < maxy ; iy++ ) {
00697 check( minx = xsh_order_list_eval_int(order_table,
00698 order_table->list[iorder].edglopoly, iy)+1);
00699 check( maxx = xsh_order_list_eval_int(order_table,
00700 order_table->list[iorder].edguppoly, iy)-1);
00701
00702 for( ix = minx ; ix <= maxx; ix++){
00703 double derrs;
00704 int idx;
00705 int sdx;
00706 float slit_val;
00707
00708 idx = ix+iy*nx;
00709 sdx = ix+iy*(nx-wmap_xsize_diff);
00710
00711
00712 if ( pqual[idx] <= XSH_GOOD_PIXEL_LEVEL){
00713 lambda = plambda[sdx];
00714 if ( lambda == 0 ) {
00715 continue;
00716 }
00717 if ( lambda < lambda_min ) {
00718 lambda_min = lambda ;
00719 }
00720 if ( lambda > lambda_max ) {
00721 lambda_max = lambda ;
00722 }
00723 slit_val = pslit[sdx];
00724
00725 if ( ((slit_val >= obj_slit_min) && (slit_val <= obj_slit_max))||
00726 (slit_val < sky_slit_min) || (slit_val > sky_slit_max)){
00727
00728 pobject->lambda = lambda;
00729 pobject->slit = slit_val;
00730 pobject->flux = pflux[idx];
00731 derrs = perrs[idx];
00732 pobject->sigma = (double)1./(derrs*derrs);
00733 pobject->ix = ix;
00734 pobject->iy = iy;
00735 object_size++;
00736 pobject++;
00737 }
00738 else{
00739 psky->lambda = lambda ;
00740 psky->slit = slit_val;
00741 psky->flux = *(pflux + ix + (iy*nx)) ;
00742
00743
00744 derrs = *(perrs + ix + (iy*nx)) ;
00745 psky->sigma = (double)1./(derrs*derrs) ;
00746 psky->ix=ix;
00747 psky->iy=iy;
00748 sky_size++ ;
00749 psky++ ;
00750 }
00751 }
00752 }
00753 }
00754 XSH_CMP_INT( sky_size, >, 0, "ERROR! On order %d sky_size 0. Make sure the order edge table is proper on that order. Probably now it over-estimates corresponding order size",,order);
00755
00756 wave_list->list[iorder].sky_size = sky_size;
00757 wave_list->list[iorder].object_size = object_size;
00758 wave_list->list[iorder].lambda_min = lambda_min ;
00759 wave_list->list[iorder].lambda_max = lambda_max ;
00760
00761
00762 qsort( wave_list->list[iorder].sky, sky_size, sizeof( wavemap_item),
00763 wave_compare);
00764 qsort( wave_list->list[iorder].object, object_size, sizeof( wavemap_item),
00765 wave_compare);
00766
00767 level = xsh_debug_level_get();
00768
00769 if ( level >= XSH_DEBUG_LEVEL_MEDIUM ) {
00770 FILE* debug = NULL;
00771 char debug_name[256];
00772 wavemap_item *pentry = NULL;
00773
00774 pentry = wave_list->list[iorder].sky;
00775 sprintf( debug_name, "Wave_sky_data_%d.log", order);
00776 debug = fopen( debug_name, "w");
00777 fprintf( debug, "# lambda flux slit x y\n");
00778 for( i=0; i< sky_size; i++){
00779 fprintf(debug, "%f %f %f %d %d\n", pentry->lambda, pentry->flux,
00780 pentry->slit, pentry->ix, pentry->iy);
00781 pentry++;
00782 }
00783 fclose( debug);
00784
00785 pentry = wave_list->list[iorder].object;
00786 sprintf( debug_name, "Wave_object_data_%d.log", order);
00787 debug = fopen( debug_name, "w");
00788 fprintf( debug, "# lambda flux slit\n");
00789 for( i=0; i< object_size; i++){
00790 fprintf(debug, "%f %f %f\n", pentry->lambda, pentry->flux, pentry->slit);
00791 pentry++;
00792 }
00793 fclose( debug);
00794 }
00795
00796 if(usr_defined_break_points_frame!=NULL){
00797
00798 name=cpl_frame_get_filename(usr_defined_break_points_frame);
00799 usr_defined_break_points_tab=cpl_table_load(name,1,0);
00800 nrow=cpl_table_get_nrow(usr_defined_break_points_tab);
00801 check(pfactor=cpl_table_get_data_double_const(usr_defined_break_points_tab,"FACTOR"));
00802 switch(xsh_instrument_get_arm(instrument)) {
00803 case XSH_ARM_UVB:XSH_ASSURE_NOT_ILLEGAL_MSG(nrow==XSH_ORDERS_UVB,
00804 "Wrong number of factors for single frame sky subtraction table");
00805 xsh_nbkpts_uvb[iorder]=pfactor[iorder];
00806 break;
00807 case XSH_ARM_VIS:XSH_ASSURE_NOT_ILLEGAL_MSG(nrow==XSH_ORDERS_VIS,
00808 "Wrong number of factors for single frame sky subtraction table");
00809 xsh_msg("iorder=%d factor=%f",iorder,pfactor[iorder]);
00810 xsh_nbkpts_vis[iorder]=pfactor[iorder];
00811 break;
00812 case XSH_ARM_NIR:XSH_ASSURE_NOT_ILLEGAL_MSG(nrow==XSH_ORDERS_NIR,
00813 "Wrong number of factors for single frame sky subtraction table");
00814 xsh_nbkpts_nir[iorder]=pfactor[iorder];
00815 break;
00816 default:xsh_msg("Arm not supported"); break;
00817 }
00818
00819 }
00820
00821
00822 if ( sky_par->method == BSPLINE_METHOD){
00823 switch(xsh_instrument_get_arm(instrument)) {
00824 case XSH_ARM_UVB:nbkpts_arm=xsh_nbkpts_uvb[iorder];break;
00825 case XSH_ARM_VIS:nbkpts_arm=xsh_nbkpts_vis[iorder];break;
00826 case XSH_ARM_NIR:nbkpts_arm=xsh_nbkpts_nir[iorder];break;
00827 default:xsh_msg("Arm not supported"); break;
00828 }
00829 if (sky_par->bspline_sampling == FINE) {
00830 nbkpts_ord=(int)(nbkpts*nbkpts_arm+0.5);
00831 } else {
00832 nbkpts_ord=nbkpts;
00833 }
00834 check( fit_spline(wave_list,iorder,nbkpts_ord,
00835 sky_par->bezier_spline_order,
00836 sky_par->niter, sky_par->kappa,
00837 sky_par->ron, sky_par->gain,
00838 xsh_instrument_arm_tostring(instrument)));
00839 }
00840 else{
00841 cpl_vector* median_vector = NULL;
00842 int mediani;
00843 double median;
00844 int nb_points = 2*median_hsize+1;
00845
00846 median_vector = cpl_vector_new( nb_points);
00847 psky = wave_list->list[iorder].sky;
00848
00849 for( i=median_hsize; i< (sky_size-median_hsize); i++){
00850 double median_err=0.0;
00851
00852 for(mediani=0; mediani< nb_points; mediani++){
00853 cpl_vector_set( median_vector, mediani, psky[i-median_hsize+mediani].flux);
00854 median_err += 1.0/psky[i-median_hsize+mediani].sigma;
00855 }
00856 median = cpl_vector_get_median(median_vector);
00857 psky[i].fitted = median;
00858 psky[i].fit_err = sqrt(
00859 M_PI/2.0*((double) nb_points/((double) nb_points- 1.))) *
00860 (1./(double)nb_points) * sqrt (median_err);
00861 }
00862 xsh_free_vector( &median_vector);
00863
00864 #if REGDEBUG_MEDIAN_SPLINE
00865 {
00866 for( i=0; i< sky_size; i++){
00867 psky[i].flux = psky[i].fitted;
00868 }
00869 XSH_REGDEBUG("fit median data");
00870
00871 switch(xsh_instrument_get_arm(instrument)) {
00872 case XSH_ARM_UVB:nbkpts_arm=xsh_nbkpts_uvb[iorder];break;
00873 case XSH_ARM_VIS:nbkpts_arm=xsh_nbkpts_vis[iorder];break;
00874 case XSH_ARM_NIR:nbkpts_arm=xsh_nbkpts_nir[iorder];break;
00875 default:xsh_msg("Arm not supported"); break;
00876 }
00877 if (sky_par->bspline_sampling == FINE) {
00878 nbkpts_ord=(int)(nbkpts*nbkpts_arm+0.5);
00879 } else {
00880 nbkpts_ord=nbkpts;
00881 }
00882 check( fit_spline(wave_list,iorder,nbkpts_ord,
00883 sky_par->bezier_spline_order,
00884 sky_par->niter, sky_par->kappa,
00885 sky_par->ron, sky_par->gain,
00886 xsh_instrument_arm_tostring(instrument)));
00887 }
00888 #endif
00889 }
00890 }
00891 cleanup :
00892 return wave_list ;
00893 }
00894
00895
00896
00902
00903 xsh_rec_list*
00904 xsh_wavemap_list_build_sky( xsh_wavemap_list* wave_list,
00905 xsh_pre* pre_mflat,
00906 xsh_instrument* instr)
00907 {
00908 xsh_rec_list* sky_list = NULL;
00909 int i, j;
00910 int level;
00911
00912
00913 XSH_ASSURE_NOT_NULL( wave_list);
00914 XSH_ASSURE_NOT_NULL( instr);
00915
00916 xsh_msg("Build sky model");
00917 level = xsh_debug_level_get();
00918
00919 check( sky_list = xsh_rec_list_create_with_size( wave_list->size, instr));
00920
00921
00922 for ( i=0; i < wave_list->size ; i++){
00923 wavemap_item *psky = NULL;
00924 int order;
00925 int sky_size;
00926 double *lambdas = NULL;
00927 float *sky_flux = NULL;
00928 float *sky_err = NULL;
00929 float* pflat=NULL;
00930 int sx=0;
00931 psky = wave_list->list[i].sky;
00932 sky_size = wave_list->list[i].sky_size;
00933 order = wave_list->list[i].order;
00934
00935 check( xsh_rec_list_set_data_size( sky_list, i, order, sky_size, 1));
00936 check( lambdas = xsh_rec_list_get_lambda( sky_list, i));
00937 check( sky_flux = xsh_rec_list_get_data1( sky_list, i));
00938 check( sky_err = xsh_rec_list_get_errs1( sky_list, i));
00939
00940 if ( pre_mflat!=NULL) {
00941 check( pflat = cpl_image_get_data_float( xsh_pre_get_data( pre_mflat)));
00942 sx=cpl_image_get_size_x(xsh_pre_get_data( pre_mflat));
00943
00944 for( j=0; j< sky_size; j++){
00945 lambdas[j] = psky->lambda;
00946 sky_flux[j] = psky->fitted*pflat[psky->ix+(psky->iy)*sx];
00947 sky_err[j] = psky->fit_err;
00948 psky++;
00949 }
00950 }
00951 else {
00952 for( j=0; j< sky_size; j++){
00953 lambdas[j] = psky->lambda;
00954 sky_flux[j] = psky->fitted;
00955 sky_err[j] = psky->fit_err;
00956 psky++;
00957 }
00958 }
00959
00960 if ( level >= XSH_DEBUG_LEVEL_MEDIUM ) {
00961 FILE* debug = NULL;
00962 char debug_name[256];
00963
00964 sprintf( debug_name, "fitted_data_sky_%d.log", order);
00965 debug = fopen( debug_name, "w");
00966 fprintf( debug, "# lambda flux err x y slit\n");
00967 psky = wave_list->list[i].sky;
00968
00969 for( j=0; j <sky_size; j++){
00970 fprintf(debug, "%f %f %f %d %d %f\n", psky->lambda, psky->fitted,
00971 psky->fit_err, psky->ix, psky->iy, psky->slit);
00972 psky++;
00973 }
00974 fclose( debug);
00975 }
00976 }
00977
00978 cleanup:
00979 if ( cpl_error_get_code() != CPL_ERROR_NONE){
00980 xsh_rec_list_free( &sky_list);
00981 }
00982 return sky_list;
00983 }
00984
00985
00986
00987
00988
00989
00990
00991
00992
00993
00994
00995
00996
00997
00998
00999
01000
01001
01002
01003
01004
01005
01006
01007
01008
01009
01010
01011
01012
01013
01014
01015
01016
01017
01018
01019
01024
01025 static cpl_frame* xsh_wavelist_subtract_sky( xsh_pre * pre_sci,
01026 xsh_pre* pre_mflat,
01027 xsh_rec_list* sky_list,
01028 xsh_wavemap_list * wave_list,
01029 xsh_instrument* instrument,
01030 const char* rec_prefix)
01031 {
01032 int level;
01033 float* perrs = NULL ;
01034 float* pflux = NULL ;
01035 float* pflat = NULL ;
01036 int* pqual = NULL ;
01037 int i, nx, ny ;
01038 xsh_pre* pre_res = NULL ;
01039 cpl_image* res_image = NULL ;
01040 cpl_frame* res_frame = NULL ;
01041 int iorder;
01042 char tag[80];
01043 char fname[80];
01044
01045
01046
01047 XSH_ASSURE_NOT_NULL( pre_sci);
01048 XSH_ASSURE_NOT_NULL( sky_list);
01049 XSH_ASSURE_NOT_NULL( wave_list);
01050
01051 check( pre_res = xsh_pre_duplicate( pre_sci));
01052 check( res_image = xsh_pre_get_data( pre_res));
01053 check( pflux = cpl_image_get_data_float( res_image));
01054 check( perrs = cpl_image_get_data_float( xsh_pre_get_errs( pre_res)));
01055 check( pqual = cpl_image_get_data_int( xsh_pre_get_qual( pre_res)));
01056
01057 if(pre_mflat != NULL) {
01058 check( pflat = cpl_image_get_data_float( xsh_pre_get_data( pre_mflat)));
01059 }
01060
01061 nx = xsh_pre_get_nx( pre_res);
01062 ny = xsh_pre_get_ny( pre_res);
01063
01064 for( iorder=0; iorder < wave_list->size; iorder++){
01065 wavemap_item* psky = NULL;
01066 wavemap_item* pobject = NULL;
01067 int order, sky_size, object_size;
01068 double* sky_lambdas = NULL;
01069 int sky_lambdas_idx = 0;
01070 int sky_lambdas_size;
01071 float* sky_flux = NULL;
01072 float* sky_err = NULL;
01073 float mflat=0;
01074 order = wave_list->list[iorder].order ;
01075 psky = wave_list->list[iorder].sky;
01076 sky_size = wave_list->list[iorder].sky_size;
01077
01078 check( sky_lambdas_size = xsh_rec_list_get_nlambda( sky_list, iorder));
01079 check( sky_lambdas = xsh_rec_list_get_lambda( sky_list, iorder));
01080 check( sky_flux = xsh_rec_list_get_data1( sky_list, iorder));
01081 check( sky_err = xsh_rec_list_get_errs1( sky_list, iorder));
01082
01083 xsh_msg_dbg_medium( "Subtract Sky - Order %d", order);
01084
01085 for( i=0; i < sky_size; i++) {
01086 int x,y;
01087 float flux, fitted, err, fitted_err;
01088
01089 x = psky->ix;
01090 y = psky->iy;
01091 flux = pflux[x+y*nx];
01092 fitted = psky->fitted;
01093 fitted_err = psky->fit_err;
01094 if(pre_mflat != NULL) {
01095 mflat=pflat[x+nx*y];
01096 pflux[x+y*nx] = flux-fitted*mflat;
01097 } else {
01098 pflux[x+y*nx] = flux-fitted;
01099 }
01100 err = perrs[x+y*nx];
01101 perrs[x+y*nx] = sqrt(fitted_err*fitted_err + err*err);
01102 psky++;
01103 }
01104 object_size = wave_list->list[iorder].object_size;
01105 pobject = wave_list->list[iorder].object;
01106
01107 xsh_msg_dbg_medium( "Subtract Sky on object - Order %d size %d", order, object_size);
01108
01109 for( i=0; i< object_size; i++){
01110 int x,y;
01111 float flux, err, fitted, fitted_err;
01112 float lambda, lambda_min, lambda_max;
01113 float flux_min, flux_max, err_min, err_max;
01114 float cte_min, cte_max;
01115
01116 x = pobject->ix;
01117 y = pobject->iy;
01118 lambda = pobject->lambda;
01119 flux = pflux[x+y*nx];
01120
01121 while ( (sky_lambdas_idx < sky_lambdas_size) &&
01122 ( sky_lambdas[sky_lambdas_idx]<= lambda)){
01123 sky_lambdas_idx++;
01124 }
01125
01126 if ( sky_lambdas_idx >= (sky_lambdas_size-1)){
01127 break;
01128 }
01129 lambda_max = sky_lambdas[sky_lambdas_idx];
01130
01131 if ( sky_lambdas_idx == 0){
01132 pobject++;
01133 continue;
01134 }
01135 lambda_min = sky_lambdas[sky_lambdas_idx-1];
01136
01137 flux_min = sky_flux[sky_lambdas_idx-1];
01138 flux_max = sky_flux[sky_lambdas_idx];
01139 err_min = sky_err[sky_lambdas_idx-1];
01140 err_max = sky_err[sky_lambdas_idx];
01141 cte_min = 1-(lambda-lambda_min)/(lambda_max-lambda_min);
01142 cte_max = (lambda-lambda_min)/(lambda_max-lambda_min);
01143 fitted = flux_min*cte_min+flux_max*cte_max;
01144 fitted_err = sqrt( err_min*err_min*cte_min + err_max*err_max*cte_max);
01145 pobject->fitted = fitted;
01146 pobject->fit_err = fitted_err;
01147 pflux[x+y*nx] = flux-fitted;
01148 err = perrs[x+y*nx];
01149 perrs[x+y*nx] = sqrt(err*err + fitted_err*fitted_err);
01150 pobject++;
01151 }
01152 }
01153
01154 level = xsh_debug_level_get();
01155
01156 if ( level >= XSH_DEBUG_LEVEL_MEDIUM ) {
01157 for( iorder=0; iorder < wave_list->size; iorder++){
01158 wavemap_item* pobject = NULL;
01159 int j, order, object_size;
01160 FILE* debug = NULL;
01161 char debug_name[256];
01162
01163
01164 order = wave_list->list[iorder].order ;
01165 object_size = wave_list->list[iorder].object_size;
01166 pobject = wave_list->list[iorder].object;
01167
01168
01169 sprintf( debug_name, "fitted_data_obj_%d.log", order);
01170 debug = fopen( debug_name, "w");
01171 fprintf( debug, "# lambda flux err x y slit\n");
01172
01173 for( j=0; j <object_size; j++){
01174 fprintf(debug, "%f %f %f %d %d %f\n", pobject->lambda, pobject->fitted,
01175 pobject->fit_err, pobject->ix, pobject->iy, pobject->slit);
01176 pobject++;
01177 }
01178 fclose( debug);
01179 }
01180 }
01181
01182 sprintf(tag,"%s_SUB_SKY_%s",rec_prefix,
01183 xsh_instrument_arm_tostring(instrument)) ;
01184 sprintf( fname, "%s.fits", tag);
01185
01186 check( res_frame = xsh_pre_save( pre_res, fname, tag, 0 ) );
01187 check(xsh_frame_config(fname,tag,CPL_FRAME_TYPE_IMAGE,CPL_FRAME_GROUP_CALIB,
01188 CPL_FRAME_LEVEL_FINAL,&res_frame));
01189
01190 cleanup:
01191 xsh_pre_free( &pre_res);
01192 return res_frame ;
01193 }
01194
01195
01196
01208 static cpl_error_code
01209 xsh_mflat_normalize(xsh_order_list * order_table, xsh_pre* pre_mflat)
01210 {
01211
01212 float* pflux=NULL;
01213 float* psmo=NULL;
01214 int* pqual=NULL;
01215 int iorder=0;
01216 int abs_order=0;
01217 float flux_max=FLT_MIN;
01218 int minx=0;
01219 int miny=0;
01220 int maxx=0;
01221 int maxy=0;
01222 int sx=0;
01223 int sy=0;
01224
01225
01226 int iy=0;
01227 int ix=0;
01228 int rad=2;
01229 float flux=0;
01230 cpl_image* data=NULL;
01231 cpl_image* smo=NULL;
01232
01233 XSH_ASSURE_NOT_NULL( order_table);
01234 XSH_ASSURE_NOT_NULL( pre_mflat);
01235 data=xsh_pre_get_data(pre_mflat);
01236
01237 smo=xsh_image_smooth_median_x(data,rad);
01238
01239 check( pflux = cpl_image_get_data_float( data));
01240 check( psmo = cpl_image_get_data_float( smo));
01241 check( pqual = cpl_image_get_data_int ( xsh_pre_get_qual( pre_mflat )));
01242
01243 sx = xsh_pre_get_nx( pre_mflat);
01244 sy = xsh_pre_get_ny( pre_mflat);
01245
01246 for ( iorder = 0; iorder < order_table->size; iorder++){
01247 miny = xsh_order_list_get_starty( order_table, iorder);
01248 maxy = xsh_order_list_get_endy( order_table, iorder);
01249 abs_order = order_table->list[iorder].absorder ;
01250
01251 for( iy = miny ; iy < maxy ; iy++ ) {
01252 check( minx = xsh_order_list_eval_int(order_table,
01253 order_table->list[iorder].edglopoly, iy)-3);
01254 check( maxx = xsh_order_list_eval_int(order_table,
01255 order_table->list[iorder].edguppoly, iy)+2);
01256
01257
01258 flux_max=FLT_MIN;
01259 for( ix = minx ; ix <= maxx; ix++){
01260 if(pqual[iy*sx+ix]==0) {
01261 flux=psmo[iy*sx+ix];
01262 flux_max = (flux>flux_max)? flux : flux_max;
01263 }
01264 }
01265
01266 for( ix = minx ; ix <= maxx; ix++){
01267 pflux[iy*sx+ix]/=flux_max;
01268 }
01269 }
01270 }
01271
01272 cleanup:
01273 xsh_free_image(&smo);
01274 return cpl_error_get_code();
01275
01276 }
01277
01312
01313 cpl_frame *
01314 xsh_subtract_sky_single (cpl_frame *sci_frame,
01315 cpl_frame *order_table_frame,
01316 cpl_frame *slitmap_frame,
01317 cpl_frame *wavemap_frame,
01318 cpl_frame *loc_table_frame,
01319 cpl_frame *mflat_frame,
01320 cpl_frame* usr_defined_break_points_frame,
01321 xsh_instrument *instrument,
01322 int nbkpts,
01323 xsh_subtract_sky_single_param* sky_par,
01324 cpl_frame **sky_spectrum,
01325 cpl_frame** sky_spectrum_eso,
01326 const char* rec_prefix)
01327 {
01328 cpl_frame *res_frame = NULL ;
01329 const char* wavemap_name = NULL;
01330 cpl_image *img_wavemap = NULL;
01331 const char* slitmap_name = NULL;
01332 cpl_image *img_slitmap = NULL;
01333 xsh_pre *pre_sci = NULL;
01334 xsh_order_list *order_table = NULL;
01335 xsh_wavemap_list *wave_list = NULL;
01336 xsh_localization *local_list = NULL;
01337 xsh_rec_list *sky_list = NULL;
01338 cpl_frame *sky_frame = NULL;
01339 cpl_frame *sky_frame2 = NULL;
01340 char fname[80];
01341 char tag[80];
01342 const char* mflat_name=NULL;
01343 xsh_pre* pre_mflat=NULL;
01344
01345
01346 XSH_ASSURE_NOT_NULL_MSG( sci_frame,"Required science frame is missing");
01347 XSH_ASSURE_NOT_NULL_MSG( order_table_frame,"Required order table frame is missing");
01348 XSH_ASSURE_NOT_NULL_MSG( slitmap_frame, "Required slitmap frame is missing, provide it or set compute-map to TRUE");
01349 XSH_ASSURE_NOT_NULL_MSG( wavemap_frame,"Required wavemap frame is missing");
01350 XSH_ASSURE_NOT_NULL_MSG( instrument,"Instrument setting undefined");
01351 XSH_ASSURE_NOT_ILLEGAL( nbkpts > 1);
01352 XSH_ASSURE_NOT_NULL_MSG( sky_par,"Undefined input sky parameters");
01353
01354
01355 xsh_msg_dbg_low( "method %s edges slit mask %f",
01356 SKY_METHOD_PRINT( sky_par->method), sky_par->slit_edges_mask);
01357 if (sky_par->method == BSPLINE_METHOD){
01358 xsh_msg_dbg_medium("start niter=%d kappa=%f ron=%f gain=%f",
01359 sky_par->niter, sky_par->kappa, sky_par->ron, sky_par->gain);
01360 }
01361 else{
01362 xsh_msg_dbg_low("median_hsize %d", sky_par->median_hsize);
01363 }
01364
01365 if ( loc_table_frame == NULL ) {
01366 xsh_msg( "Subtract sky single no localization");
01367 }
01368 else {
01369 xsh_msg( "Subtract sky single using localization");
01370 check( local_list = xsh_localization_load( loc_table_frame));
01371 }
01372
01373 check( slitmap_name = cpl_frame_get_filename( slitmap_frame));
01374 check( img_slitmap = cpl_image_load( slitmap_name, CPL_TYPE_DOUBLE,
01375 0, 0));
01376
01377 check( pre_sci = xsh_pre_load( sci_frame, instrument));
01378
01379 if( sky_par->gain < 0.) {
01380 if( xsh_instrument_get_arm(instrument) == XSH_ARM_NIR){
01381 sky_par->gain=2.12;
01382 } else {
01383 sky_par->gain=xsh_pfits_get_gain(pre_sci->data_header);
01384 }
01385 }
01386
01387
01388 if( sky_par->ron < 0.) {
01389 if( xsh_instrument_get_arm(instrument) == XSH_ARM_NIR){
01390 sky_par->ron=10;
01391 } else {
01392 sky_par->ron=xsh_pfits_get_ron(pre_sci->data_header);
01393 }
01394 }
01395
01396 check( order_table = xsh_order_list_load (order_table_frame, instrument));
01397
01398 check( wavemap_name = cpl_frame_get_filename( wavemap_frame));
01399 check( img_wavemap = cpl_image_load( wavemap_name, CPL_TYPE_DOUBLE,
01400 0, 0));
01401
01402
01403 if(mflat_frame!=NULL) {
01404 xsh_msg("--subtract_sky_single : normalize master flat");
01405 mflat_name=cpl_frame_get_filename( mflat_frame);
01406 pre_mflat=xsh_pre_load( mflat_frame, instrument);
01407 check( xsh_mflat_normalize( order_table, pre_mflat));
01408
01409 check(cpl_image_save( xsh_pre_get_data( pre_mflat),"mflat_norm.fits",
01410 CPL_BPP_IEEE_FLOAT,NULL,CPL_IO_DEFAULT));
01411 }
01412
01413 check( wave_list = xsh_wavemap_list_new( img_wavemap, img_slitmap,
01414 local_list, order_table, pre_sci, nbkpts, usr_defined_break_points_frame,sky_par, instrument));
01415
01416
01417 check (sky_list = xsh_wavemap_list_build_sky(wave_list,pre_mflat,instrument));
01418
01419 sprintf(tag,"%s_DRL_SKY_ORD1D_%s", rec_prefix,
01420 xsh_instrument_arm_tostring(instrument));
01421 sprintf(fname,"%s.fits",tag);
01422
01423 check( sky_frame = xsh_rec_list_save( sky_list,fname,tag, CPL_TRUE));
01424
01425 if ( sky_spectrum != NULL){
01426 *sky_spectrum = cpl_frame_duplicate( sky_frame);
01427 }
01428
01429 sprintf( tag,"%s_SKY_ORD1D_%s", rec_prefix,
01430 xsh_instrument_arm_tostring(instrument));
01431 sprintf(fname,"%s.fits",tag);
01432
01433 check( sky_frame2 = xsh_rec_list_save2( sky_list,fname,tag));
01434
01435 if ( sky_spectrum != NULL){
01436 *sky_spectrum_eso = cpl_frame_duplicate( sky_frame2);
01437 }
01438
01439
01440 check( res_frame = xsh_wavelist_subtract_sky( pre_sci, pre_mflat,sky_list,
01441 wave_list,instrument,
01442 rec_prefix));
01443 cleanup:
01444 xsh_free_frame( &sky_frame);
01445 xsh_free_frame( &sky_frame2);
01446 xsh_free_image( &img_slitmap);
01447 xsh_free_image( &img_wavemap);
01448 xsh_pre_free( &pre_mflat);
01449 xsh_pre_free( &pre_sci);
01450 xsh_localization_free( &local_list);
01451 xsh_order_list_free( &order_table);
01452 xsh_wavemap_list_free( &wave_list);
01453 xsh_rec_list_free( &sky_list);
01454 return res_frame;
01455 }
01456
01457
01458
01459
01460 static int
01461 xsh_bspline_eval_all(const double x, gsl_vector *B, size_t *idx,
01462 gsl_bspline_workspace *w, size_t *start_i)
01463 {
01464 if (B->size != w->k)
01465 {
01466 GSL_ERROR("B vector not of length k", GSL_EBADLEN);
01467 }
01468 else
01469 {
01470 size_t i;
01471 size_t j;
01472 size_t ii;
01473 int flg = 0;
01474 double saved;
01475 double term;
01476
01477 double *knot_ptr;
01478 double *deltar_ptr;
01479 double *deltal_ptr;
01480 double *B_ptr;
01481
01482 knot_ptr=gsl_vector_ptr(w->knots,0);
01483 deltar_ptr=gsl_vector_ptr(w->deltar,0);
01484 deltal_ptr=gsl_vector_ptr(w->deltal,0);
01485 B_ptr=gsl_vector_ptr(B,0);
01486
01487 i = xsh_bspline_find_interval(x, &flg, w, start_i);
01488
01489 if (flg == -1)
01490 {
01491 GSL_ERROR("x outside of knot interval", GSL_EINVAL);
01492 }
01493 else if (flg == 1)
01494 {
01495 if (x <= knot_ptr[i] + GSL_DBL_EPSILON)
01496 {
01497 --i;
01498 }
01499 else
01500 {
01501 GSL_ERROR("x outside of knot interval", GSL_EINVAL);
01502 }
01503 }
01504
01505 *idx = i;
01506
01507 B_ptr[0]=1.0;
01508
01509 for (j = 0; j < w->k - 1; ++j)
01510 {
01511 deltar_ptr[j]=knot_ptr[i + j + 1] - x;
01512 deltal_ptr[j]=x - knot_ptr[i - j];
01513
01514 saved = 0.0;
01515
01516 for (ii = 0; ii <= j; ++ii)
01517 {
01518 term = B_ptr[ii] /
01519 (deltar_ptr[ii] + deltal_ptr[j - ii]);
01520
01521 B_ptr[ii]=saved + deltar_ptr[ii] * term;
01522
01523 saved = deltal_ptr[j - ii] * term;
01524 }
01525
01526 B_ptr[j + 1]=saved;
01527 }
01528
01529 return GSL_SUCCESS;
01530 }
01531 }
01532
01533
01534
01535
01536 static size_t
01537 xsh_bspline_find_interval(const double x, int *flg,
01538 gsl_bspline_workspace *w, size_t *start_i)
01539 {
01540 size_t i;
01541
01542 double *knot_ptr;
01543
01544 knot_ptr=gsl_vector_ptr(w->knots,0);
01545
01546 if (x < knot_ptr[0])
01547 {
01548 *flg = -1;
01549 return 0;
01550 }
01551
01552
01553 for (i = *start_i; i < w->k + w->l - 1; ++i)
01554 {
01555 const double ti = knot_ptr[i];
01556 const double tip1 = knot_ptr[i + 1];
01557
01558
01559
01560
01561
01562
01563 if (ti <= x && x < tip1)
01564 break;
01565
01566 if (ti < x && x == tip1 &&
01567 tip1 == knot_ptr[ w->k + w->l - 1])
01568 break;
01569
01570 }
01571
01572 if (i == w->k + w->l - 1)
01573 *flg = 1;
01574 else
01575 *flg = 0;
01576
01577 *start_i=i;
01578 return i;
01579 }
01580
01581
01582
01583
01584
01585
01586
01587 cpl_frame* xsh_save_sky_model( cpl_frame* obj_frame, cpl_frame* sub_sky_frame,
01588 const char* sky_tag,xsh_instrument* instrument)
01589 {
01590 char sky_name[80];
01591 cpl_frame* result=NULL;
01592
01593 sprintf(sky_name,"%s.fits",sky_tag);
01594 check(result=cpl_frame_duplicate(sub_sky_frame));
01595 result=xsh_pre_frame_subtract( obj_frame, sub_sky_frame, sky_name, instrument);
01596
01597 cpl_frame_set_filename(result,sky_name);
01598 cpl_frame_set_tag(result,sky_tag);
01599
01600 cleanup:
01601
01602 return result;
01603 }
01604
01605
01606
01607
01608 cpl_frame * xsh_add_sky_model( cpl_frame *subsky_frame, cpl_frame *sky_frame,
01609 xsh_instrument * instrument, const char* prefix)
01610 {
01611 char result_tag[256];
01612 char result_name[256];
01613 cpl_frame *result = NULL;
01614 xsh_pre *pre_sci = NULL;
01615 const char *sky_name = NULL;
01616 cpl_image *sky_img = NULL;
01617
01618 sprintf( result_tag, "%s_OBJ_AND_SKY_NOCRH_%s", prefix,
01619 xsh_instrument_arm_tostring(instrument));
01620 sprintf( result_name, "%s.fits", result_tag);
01621
01622
01623 check( pre_sci = xsh_pre_load( subsky_frame, instrument));
01624
01625 check( sky_name = cpl_frame_get_filename(sky_frame));
01626 check( sky_img = cpl_image_load( sky_name, CPL_TYPE_FLOAT, 0, 0));
01627 check( cpl_image_add( pre_sci->data, sky_img));
01628 check( result = xsh_pre_save( pre_sci, result_name, result_tag, 1));
01629
01630 cleanup:
01631 xsh_free_image( &sky_img);
01632 xsh_pre_free( &pre_sci);
01633
01634 return result;
01635 }
01636
01637