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