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