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
00027 #ifdef HAVE_CONFIG_H
00028 #include <config.h>
00029 #endif
00030
00035
00038
00039
00040
00041
00042 #include <math.h>
00043 #include <xsh_drl.h>
00044
00045 #include <xsh_badpixelmap.h>
00046 #include <xsh_data_rec.h>
00047 #include <xsh_data_pre.h>
00048 #include <xsh_utils_wrappers.h>
00049 #include <xsh_data_order.h>
00050 #include <xsh_data_wavesol.h>
00051 #include <xsh_data_spectralformat.h>
00052 #include <xsh_data_localization.h>
00053 #include <xsh_dfs.h>
00054 #include <xsh_pfits.h>
00055 #include <xsh_error.h>
00056 #include <xsh_msg.h>
00057 #include <xsh_fit.h>
00058 #include <xsh_model_io.h>
00059 #include <xsh_model_kernel.h>
00060 #include <xsh_rectify.h>
00061 #include <cpl.h>
00062 #include <gsl/gsl_sf_erf.h>
00063 #include <xsh_blaze.h>
00064
00065
00066
00067 #define SLIT_USEFUL_WINDOW_FACTOR 0.80
00068 #define FIT_FWHM_LIMIT 30
00069 #define OPT_EXTRACT_SLIT_SIZE 80
00070
00071 #define FLUX_MODE 0
00072 #define WAVEMAP_MODE 1
00073
00074 #define REGDEBUG_PIXELSIZE 0
00075 #define REGDEBUG_INTEGRATE 0
00076
00077
00078
00079
00080
00081 static void xsh_image_gaussian_fit_y( cpl_image* img, int chunk_size,
00082 int deg_poly, int oversample,
00083 cpl_polynomial** center, cpl_polynomial** height,
00084 cpl_polynomial** width, cpl_polynomial** offset);
00085
00086 static cpl_image* xsh_image_create_gaussian_image( cpl_image *src,
00087 cpl_polynomial* centerp, cpl_polynomial* heightp, cpl_polynomial* widthp,
00088 cpl_polynomial* offsetp);
00089
00090 static cpl_image* xsh_image_create_model_image( cpl_image *src_img,
00091 double *x_data, double *y_data, double kappa, int niter, double frac_min);
00092
00093 static cpl_vector* xsh_vector_integrate( int biny, int absorder, int oversample,
00094 cpl_vector *ref_pos, double step,
00095 cpl_vector *ref_values, cpl_vector *err_values, cpl_vector *qual_values,
00096 cpl_vector *new_pos,
00097 cpl_vector **spectrum_err, cpl_vector **spectrum_qual);
00098
00099
00100
00101
00102
00103
00104 static void xsh_interpolate_spectrum( int biny, int oversample,
00105 int absorder, double lambda_step,
00106 cpl_vector *init_pos,
00107 cpl_vector *std_flux, cpl_vector *std_err, cpl_vector *std_qual,
00108 cpl_vector *opt_flux, cpl_vector *opt_err, cpl_vector *opt_qual,
00109 cpl_vector **res_pos,
00110 cpl_vector **res_std_flux, cpl_vector **res_std_err, cpl_vector **res_std_qual,
00111 cpl_vector **res_opt_flux, cpl_vector **res_opt_err, cpl_vector **res_opt_qual)
00112 {
00113
00114 int init_size;
00115 double lambda_min, lambda_max;
00116 double binlambda_min, binlambda_max;
00117 int mult_min, mult_max;
00118 int res_size, i;
00119
00120 XSH_ASSURE_NOT_NULL( init_pos);
00121 XSH_ASSURE_NOT_NULL( std_flux);
00122 XSH_ASSURE_NOT_NULL( std_err);
00123 XSH_ASSURE_NOT_NULL( std_qual);
00124 XSH_ASSURE_NOT_NULL( opt_flux);
00125 XSH_ASSURE_NOT_NULL( opt_err);
00126 XSH_ASSURE_NOT_NULL( opt_qual);
00127 XSH_ASSURE_NOT_NULL( res_pos);
00128 XSH_ASSURE_NOT_NULL( res_std_flux);
00129 XSH_ASSURE_NOT_NULL( res_std_err);
00130 XSH_ASSURE_NOT_NULL( res_std_qual);
00131 XSH_ASSURE_NOT_NULL( res_opt_flux);
00132 XSH_ASSURE_NOT_NULL( res_opt_err);
00133 XSH_ASSURE_NOT_NULL( res_opt_qual);
00134
00135 check( init_size = cpl_vector_get_size( init_pos));
00136
00137 check( lambda_min = cpl_vector_get( init_pos, 0));
00138 check( lambda_max = cpl_vector_get( init_pos, init_size-1));
00139
00140
00141 mult_min = (int) ceil(lambda_min/lambda_step)+1;
00142 binlambda_min = mult_min*lambda_step;
00143
00144 mult_max = (int) floor(lambda_max/lambda_step)-1;
00145 binlambda_max = mult_max*lambda_step;
00146
00147 xsh_msg("lambda_min %f lambda_max %f TEST %f %f : step %f", lambda_min, lambda_max, binlambda_min, binlambda_max, lambda_step);
00148
00149 res_size = mult_max-mult_min+1;
00150
00151 check( *res_pos = cpl_vector_new( res_size));
00152
00153 for( i=0; i< res_size; i++){
00154 check( cpl_vector_set( *res_pos, i,binlambda_min+i*lambda_step));
00155 }
00156
00157 check( *res_std_flux = xsh_vector_integrate(
00158 biny, absorder, oversample, init_pos, lambda_step,
00159 std_flux, std_err, std_qual, *res_pos, res_std_err, res_std_qual));
00160 check( *res_opt_flux = xsh_vector_integrate(
00161 biny, absorder, oversample, init_pos, lambda_step,
00162 opt_flux, opt_err, opt_qual, *res_pos, res_opt_err, res_opt_qual));
00163
00164 cleanup:
00165 if (cpl_error_get_code() != CPL_ERROR_NONE){
00166 xsh_free_vector( res_pos);
00167 xsh_free_vector( res_std_flux);
00168 xsh_free_vector( res_std_err);
00169 xsh_free_vector( res_std_qual);
00170 xsh_free_vector( res_opt_flux);
00171 xsh_free_vector( res_opt_err);
00172 xsh_free_vector( res_opt_qual);
00173 }
00174 return;
00175 }
00176
00177
00178 cpl_image* xsh_optextract_produce_model( cpl_image* s2Dby1D_img,
00179 int method, int chunk_ovsamp_size, int deg_poly, int oversample,
00180 double *extract_x_data, double *extract_y_data, int abs_order,
00181 double kappa, int niter, double frac_min){
00182
00183 cpl_polynomial *center_gauss_poly = NULL;
00184 cpl_polynomial *height_gauss_poly = NULL;
00185 cpl_polynomial *offset_gauss_poly = NULL;
00186 cpl_polynomial *width_gauss_poly = NULL;
00187 cpl_image *model_img = NULL;
00188 int model_x, model_y, ix, iy;
00189 double *model_data = NULL;
00190
00191 XSH_ASSURE_NOT_NULL( s2Dby1D_img);
00192 XSH_ASSURE_NOT_NULL( extract_x_data);
00193 XSH_ASSURE_NOT_NULL( extract_y_data);
00194
00195 xsh_msg_dbg_high( "Produce model for order %d", abs_order);
00196
00197 if ( method == GAUSS_METHOD){
00198 check( xsh_image_gaussian_fit_y( s2Dby1D_img, chunk_ovsamp_size,
00199 deg_poly, oversample,
00200 ¢er_gauss_poly, &height_gauss_poly, &width_gauss_poly,
00201 &offset_gauss_poly));
00202
00203 check( model_img = xsh_image_create_gaussian_image( s2Dby1D_img,
00204 center_gauss_poly,
00205 height_gauss_poly, width_gauss_poly, offset_gauss_poly));
00206
00207 #if REGDEBUG_GAUSSIAN
00208 if ( xsh_debug_level_get() >= XSH_DEBUG_LEVEL_MEDIUM) {
00209 char poly_name[256];
00210 FILE* poly_file = NULL;
00211 int nx = cpl_image_get_size_x( model_img);
00212 int i;
00213
00214 sprintf( poly_name, "gaussian_center_poly_%d.dat", abs_order);
00215 poly_file = fopen( poly_name, "w");
00216
00217 for(i=0; i< nx; i++){
00218 fprintf( poly_file, "%d %d\n", i+1,
00219 (int) cpl_polynomial_eval_1d( center_gauss_poly, i, NULL)+1);
00220 }
00221 fclose( poly_file);
00222
00223 sprintf( poly_name, "gaussian_sigma_poly_%d.dat", abs_order);
00224 poly_file = fopen( poly_name, "w");
00225
00226 for(i=0; i< nx; i++){
00227 fprintf( poly_file, "%d %d\n", i+1,
00228 (int) cpl_polynomial_eval_1d( width_gauss_poly, i, NULL)+1);
00229 }
00230 fclose( poly_file);
00231 }
00232 #endif
00233 }
00234 else{
00235 check( model_img = xsh_image_create_model_image( s2Dby1D_img,
00236 extract_x_data, extract_y_data, kappa, niter, frac_min));
00237 }
00238
00239
00240 model_x = cpl_image_get_size_x( model_img);
00241 model_y = cpl_image_get_size_y( model_img);
00242 model_data = cpl_image_get_data_double( model_img);
00243
00244 for( ix=0; ix < model_x; ix++){
00245 double psum =0.0;
00246
00247 for( iy=0; iy < model_y; iy++){
00248 int ipos = iy*model_x+ix;
00249
00250 if ( model_data[ipos] < 0){
00251 model_data[ipos] = 0;
00252 }
00253 psum += model_data[ipos];
00254 }
00255
00256 for( iy=0; iy < model_y; iy++){
00257 int ipos = iy*model_x+ix;
00258
00259 model_data[ipos] /= psum;
00260 }
00261 }
00262
00263 cleanup:
00264 if (cpl_error_get_code() != CPL_ERROR_NONE){
00265 xsh_free_image( &model_img);
00266 }
00267 xsh_free_polynomial( ¢er_gauss_poly);
00268 xsh_free_polynomial( &width_gauss_poly);
00269 xsh_free_polynomial( &height_gauss_poly);
00270 xsh_free_polynomial( &offset_gauss_poly);
00271 return model_img;
00272 }
00273
00274
00275
00276
00277
00278
00279 static void xsh_object_localize( cpl_frame *slitmap_frame,
00280 cpl_frame *loc_frame, int oversample, int box_hsize,
00281 int nlambdas, int ny_extract, double *extract_x_data,
00282 double *extract_y_data, int *ymin, int *ymax)
00283 {
00284 int u;
00285 double scen=0.0;
00286 cpl_image *slitmap = NULL;
00287 float *slit_data = NULL;
00288 int slit_nx;
00289 int ny_extract_min =-1, ny_extract_max=-1, ny_extract_cen=-1;
00290 int mid_lambda;
00291 double diff_s_min = 50.0;
00292 xsh_localization *loc = NULL;
00293 int ref_ov=0;
00294
00295 XSH_ASSURE_NOT_NULL( slitmap_frame);
00296 XSH_ASSURE_NOT_NULL( loc_frame);
00297 XSH_ASSURE_NOT_NULL( extract_x_data);
00298 XSH_ASSURE_NOT_NULL( extract_y_data);
00299 XSH_ASSURE_NOT_NULL( ymin);
00300 XSH_ASSURE_NOT_NULL( ymax);
00301
00302 check( slitmap = cpl_image_load( cpl_frame_get_filename( slitmap_frame),
00303 CPL_TYPE_FLOAT, 0, 0));
00304 check( loc = xsh_localization_load( loc_frame));
00305
00306 check( slit_nx = cpl_image_get_size_x( slitmap));
00307 check( slit_data = cpl_image_get_data_float( slitmap));
00308
00309 check( scen = cpl_polynomial_eval_1d( loc->cenpoly, 0, NULL));
00310
00311 mid_lambda = 0.5*nlambdas;
00312
00313 for( u=1; u < ny_extract; u++){
00314 int x,y;
00315 double s;
00316 double diff;
00317
00318 x = (int) (extract_x_data[u*nlambdas+mid_lambda]+0.5);
00319 y = (int) (extract_y_data[u*nlambdas+mid_lambda]+0.5);
00320 s = slit_data[x+y*slit_nx];
00321
00322 diff = fabs(s-scen);
00323
00324 if (diff < diff_s_min){
00325 ny_extract_cen= u;
00326 diff_s_min = diff;
00327 }
00328 }
00329
00330
00331 ref_ov += floor(0.5*(oversample-1));
00332
00333 ny_extract_min= ny_extract_cen -box_hsize-ref_ov;
00334 ny_extract_max= ny_extract_cen +box_hsize+ref_ov;
00335
00336 if ( ny_extract_min < 0){
00337 ny_extract_min = 0;
00338 }
00339 if ( ny_extract_max >= ny_extract){
00340 ny_extract_max = ny_extract-1;
00341 }
00342 *ymin = ny_extract_min;
00343 *ymax = ny_extract_max;
00344
00345 cleanup:
00346 xsh_free_image( &slitmap);
00347 xsh_localization_free( &loc);
00348 return;
00349 }
00350
00351
00352 static int xsh_interpolate_linear(float *fluxtab, float *errtab, int *qualtab,
00353 int nx, int ny, float pos_x,
00354 float pos_y, double *flux, double *err, int* qual, int mode){
00355
00356 float f00=0.0, f10=0.0, f01=0.0, f11=0.0;
00357 float e00=0.0, e10=0.0, e01=0.0, e11=0.0;
00358 int intx, inty;
00359 float fx=0.0, fy=0.0;
00360
00361 int qual_comb = -1;
00362 double A1, A2, A3, A4;
00363 int ret=0;
00364
00365 intx = (int) pos_x;
00366 inty = (int) pos_y;
00367
00368 XSH_ASSURE_NOT_ILLEGAL( intx >= 0 && intx <nx);
00369 XSH_ASSURE_NOT_ILLEGAL( inty >= 0 && inty <ny);
00370 XSH_ASSURE_NOT_NULL( flux);
00371 XSH_ASSURE_NOT_NULL( err);
00372 int pix=intx+inty*nx;
00373 int pix_plus_1=pix+1;
00374
00375 f00 = fluxtab[pix];
00376 e00 = errtab[pix];
00377 qual_comb = qualtab[pix];
00378
00379 if ( (intx +1) < nx){
00380 f10 = fluxtab[pix_plus_1];
00381 e10 = errtab[pix_plus_1];
00382 fx = pos_x-intx;
00383 qual_comb |= qualtab[pix_plus_1];
00384 }
00385 if ( (inty +1) < ny){
00386 f01 = fluxtab[pix+nx];
00387 e01 = errtab[pix+nx];
00388 fy = pos_y-inty;
00389 qual_comb |= qualtab[pix+nx];
00390 }
00391 if ( ((intx +1) < nx) && ((inty +1) < ny)){
00392 f11 = fluxtab[pix_plus_1+nx];
00393 e11 = errtab[pix_plus_1+nx];
00394 qual_comb |= qualtab[pix_plus_1+nx];
00395 }
00396
00397 if ( mode == WAVEMAP_MODE &&
00398 ( f00 == 0.0 || f10 == 0.0 || f01 == 0 || f11 == 0)){
00399 xsh_msg_dbg_medium("pixel %f, %f at zero, interpolate with "\
00400 "(%d,%d)%f, (%d,%d)%f (%d,%d)%f, (%d,%d)%f", pos_x, pos_y,
00401 intx, inty, f00, intx+1, inty, f10, intx, inty+1, f01,
00402 intx+1, inty+1, f11);
00403 ret = 1;
00404 }
00405 A1 = (1-fx)*(1-fy);
00406 A2 = fx*(1-fy);
00407 A3 = (1-fx)*fy;
00408 A4 = fx*fy;
00409
00410 *flux = f00*A1 + f10*A2 + f01*A3 + f11*A4;
00411
00412 *err = sqrt( A1*A1*e00*e00+A2*A2*e10*e10+A3*A3*e01*e01+A4*A4*e11*e11);
00413 *qual = qual_comb;
00414
00415 cleanup:
00416 return ret;
00417 }
00418
00419
00420
00455 static void xsh_wavemap_lambda_range( cpl_frame *wavemap_frame,
00456 cpl_frame *slitmap_frame,
00457 int starty,
00458 int endy, int oversample, xsh_order_list *order_list, int iorder,double *xtab,
00459 double *ytab, int order, xsh_spectralformat_list *spectralformat,
00460 double *lambdastab, double *slitstab, int* sizetab, xsh_instrument *instr)
00461 {
00462 cpl_image* wavemap = NULL;
00463 cpl_image *slitmap = NULL;
00464 int* qual = NULL;
00465 float *wavemap_data = NULL;
00466 float *slitmap_data = NULL;
00467 int i, nx=0, ny=0, qual_val=0;
00468 double lambda=0.0, lambda_next=0.0, lambda_old=0, err=0;
00469 double lambda_sf_max=0;
00470 double slit= -1;
00471 float y_ovs=0, x=0, y=0;
00472 int size=0;
00473 int ycen=0;
00474 double xcen=0;
00475 double lambda_min=0, lambda_max=0;
00476 const char* wave_filename = NULL;
00477 const char* slitmap_name = NULL;
00478
00479 XSH_ASSURE_NOT_NULL( wavemap_frame);
00480 XSH_ASSURE_NOT_NULL( sizetab);
00481 XSH_ASSURE_NOT_NULL( order_list);
00482 XSH_ASSURE_NOT_NULL( xtab);
00483 XSH_ASSURE_NOT_NULL( ytab);
00484 XSH_ASSURE_NOT_NULL( spectralformat);
00485 XSH_ASSURE_NOT_NULL( lambdastab);
00486
00487 check( wave_filename = cpl_frame_get_filename( wavemap_frame));
00488 check( slitmap_name = cpl_frame_get_filename( slitmap_frame));
00489 xsh_msg_dbg_medium("Wave map name %s", wave_filename);
00490 xsh_msg_dbg_medium("Slit map name %s", slitmap_name);
00491
00492 check( wavemap = cpl_image_load( wave_filename,
00493 CPL_TYPE_FLOAT, 0, 0));
00494 check( slitmap = cpl_image_load( slitmap_name,
00495 CPL_TYPE_FLOAT, 0, 0));
00496
00497 check( nx = cpl_image_get_size_x( wavemap));
00498 check( ny = cpl_image_get_size_y( wavemap));
00499
00500 XSH_CALLOC(qual, int, nx*ny);
00501 check( wavemap_data = cpl_image_get_data_float( wavemap));
00502 check( slitmap_data = cpl_image_get_data_float( slitmap));
00503
00504 check( lambda_min = xsh_spectralformat_list_get_lambda_min(
00505 spectralformat, order));
00506 check( lambda_sf_max = xsh_spectralformat_list_get_lambda_max(
00507 spectralformat, order));
00508
00509 xsh_msg_dbg_high( "SPECTRAL FORMAT lambda min %f max %f", lambda_min, lambda_sf_max);
00510
00511 lambda_max = lambda_sf_max;
00512
00513
00514 if ( xsh_instrument_get_arm( instr) !=XSH_ARM_NIR){
00515
00516 ycen = starty-1;
00517
00518 xsh_msg_dbg_high("DBG ycen %d", ycen);
00519
00520 while( !((lambda >= lambda_min) && (lambda_next > lambda)) ){
00521 double start_x_next;
00522
00523 ycen++;
00524 check( xcen = xsh_order_list_eval( order_list, order_list->list[iorder].cenpoly,
00525 ycen));
00526 xsh_msg_dbg_high("DBG xcen %f ycen %d", xcen, ycen);
00527 check( xsh_interpolate_linear( wavemap_data, wavemap_data, qual,
00528 nx, ny, xcen-1, ycen-1, &lambda, &err, &qual_val, WAVEMAP_MODE));
00529 xsh_msg_dbg_high("interLambda %f", lambda);
00530
00531 check( start_x_next= xsh_order_list_eval( order_list, order_list->list[iorder].cenpoly,
00532 ycen+1));
00533 check( xsh_interpolate_linear( wavemap_data, wavemap_data, qual,
00534 nx, ny, start_x_next-1, ycen, &lambda_next, &err, &qual_val, WAVEMAP_MODE));
00535 }
00536 y_ovs=ycen*oversample+1;
00537 }
00538 else{
00539
00540 ycen = endy+1;
00541
00542
00543
00544 while( !((lambda >= lambda_min) && (lambda_next > lambda)) ){
00545 double start_x_next;
00546
00547 ycen--;
00548 check( xcen = xsh_order_list_eval( order_list, order_list->list[iorder].cenpoly,
00549 ycen));
00550
00551 check( xsh_interpolate_linear( wavemap_data, wavemap_data, qual,
00552 nx, ny, xcen-1, ycen-1, &lambda, &err, &qual_val, WAVEMAP_MODE));
00553
00554
00555 check( start_x_next= xsh_order_list_eval( order_list, order_list->list[iorder].cenpoly,
00556 ycen-1));
00557 check( xsh_interpolate_linear( wavemap_data, wavemap_data, qual,
00558 nx, ny, start_x_next-1, ycen-2, &lambda_next, &err, &qual_val, WAVEMAP_MODE));
00559 }
00560 y_ovs=ycen*oversample-1;
00561 }
00562
00563
00564 size=0;
00565 x= (float) xcen;
00566 y =(float) ycen;
00567 lambda_old = -1;
00568
00569 xsh_msg_dbg_high("first X,Y %f,%f lambda %f in [%f,%f]", x,y, lambda, lambda_min, lambda_max);
00570
00571 if ( xsh_instrument_get_arm( instr) !=XSH_ARM_NIR){
00572
00573 while( (lambda >= lambda_min) && (lambda > lambda_old) &&
00574 (lambda <= lambda_max) &&
00575 ( y_ovs <= (endy * oversample)) ){
00576
00577 xsh_msg_dbg_high("size %d lambda %f ,x %f y %f", size, lambda, x,y);
00578
00579 lambdastab[size] = lambda;
00580 xtab[size] = x;
00581 ytab[size] = y;
00582 size++;
00583 lambda_old = lambda;
00584 y = (float) (y_ovs) / (float) oversample;
00585
00586 check( x = xsh_order_list_eval( order_list, order_list->list[iorder].cenpoly,
00587 y));
00588
00589 check( xsh_interpolate_linear( wavemap_data, wavemap_data, qual,
00590 nx, ny, x-1, y-1, &lambda, &err, &qual_val, WAVEMAP_MODE));
00591 y_ovs++;
00592 }
00593 }
00594 else{
00595 while( (lambda >= lambda_min) &&
00596 (lambda <= lambda_max) &&
00597 ( y_ovs >= (starty * oversample)) ){
00598
00599 xsh_msg_dbg_high("size %d lambda %f ,x %f y %f", size, lambda, x,y);
00600
00601 lambdastab[size] = lambda;
00602 xtab[size] = x;
00603 ytab[size] = y;
00604 size++;
00605 lambda_old = lambda;
00606 y = (float) (y_ovs) / (float) oversample;
00607
00608 check( x = xsh_order_list_eval( order_list, order_list->list[iorder].cenpoly,
00609 y));
00610
00611 check( xsh_interpolate_linear( wavemap_data, wavemap_data, qual,
00612 nx, ny, x-1, y-1, &lambda, &err, &qual_val, WAVEMAP_MODE));
00613 y_ovs--;
00614
00615
00616 }
00617 }
00618
00619 *sizetab = size;
00620
00621 for (i=0; i< size; i++){
00622 x = xtab[i];
00623 y = ytab[i];
00624 check( xsh_interpolate_linear( slitmap_data, slitmap_data, qual,
00625 nx, ny, x-1, y-1, &slit, &err, &qual_val, WAVEMAP_MODE));
00626 slitstab[i] = slit;
00627 }
00628
00629
00630
00631
00632
00633
00634
00635
00636
00637
00638
00639
00640
00641 cleanup:
00642 XSH_FREE( qual);
00643 xsh_free_image( &wavemap);
00644 return;
00645 }
00646
00647
00648
00649
00664
00665 static void xsh_vector_divide_poly( cpl_vector *vector, double oversample,
00666 cpl_polynomial *poly, int shift, xsh_instrument *instr)
00667 {
00668 int i, size;
00669
00670 XSH_ASSURE_NOT_NULL( vector);
00671 XSH_ASSURE_NOT_NULL( poly);
00672
00673 check( size = cpl_vector_get_size( vector));
00674
00675 for( i=0; i< size; i++){
00676 double flux, val;
00677
00678 check( val = cpl_polynomial_eval_1d( poly,
00679 (double)i/oversample+shift, NULL));
00680 if ( xsh_instrument_get_arm( instr) == XSH_ARM_NIR){
00681 check( val = cpl_polynomial_eval_1d( poly,
00682 (double)(size-1-i)/oversample+shift, NULL));
00683 }
00684 else{
00685 check( val = cpl_polynomial_eval_1d( poly,
00686 (double)i/oversample+shift, NULL));
00687 }
00688 check( flux = cpl_vector_get( vector, i));
00689 check( cpl_vector_set( vector, i, flux/val));
00690 }
00691
00692 cleanup:
00693 return;
00694 }
00695
00696
00697
00719
00720 static void xsh_image_gaussian_fit_y( cpl_image* img, int chunk_size,
00721 int deg_poly, int oversample,
00722 cpl_polynomial** center, cpl_polynomial** height,
00723 cpl_polynomial** width, cpl_polynomial** offset)
00724 {
00725
00726 int nx, ny, nb_chunk, i;
00727 cpl_vector *center_gauss = NULL;
00728 cpl_vector *height_gauss = NULL;
00729 cpl_vector *width_gauss = NULL;
00730 cpl_vector *offset_gauss = NULL;
00731 cpl_vector *xpos_gauss = NULL;
00732 cpl_vector *chunk_vector = NULL;
00733 cpl_vector *chunk_pos_vector = NULL;
00734 cpl_vector *median_vect = NULL;
00735 double *data = NULL;
00736 double *center_gauss_data = NULL;
00737 double *height_gauss_data = NULL;
00738 double *width_gauss_data = NULL;
00739 double *offset_gauss_data = NULL;
00740 double *xpos_gauss_data = NULL;
00741 double *median_data = NULL;
00742
00743 int nbfit = 0;
00744
00745 XSH_ASSURE_NOT_NULL( img);
00746 XSH_ASSURE_NOT_ILLEGAL( chunk_size >=0);
00747 XSH_ASSURE_NOT_ILLEGAL( deg_poly >=0);
00748 XSH_ASSURE_NOT_NULL( center);
00749 XSH_ASSURE_NOT_NULL( height);
00750 XSH_ASSURE_NOT_NULL( width);
00751 XSH_ASSURE_NOT_NULL( offset);
00752
00753 check( nx = cpl_image_get_size_x( img));
00754 check( ny = cpl_image_get_size_y( img));
00755 check( data = cpl_image_get_data_double( img));
00756
00757 if ( (chunk_size >= nx) || (chunk_size == 0)){
00758 chunk_size = nx;
00759 }
00760
00761 nb_chunk = nx / chunk_size;
00762
00763 XSH_MALLOC( center_gauss_data, double, nb_chunk*2);
00764 XSH_MALLOC( height_gauss_data, double, nb_chunk*2);
00765 XSH_MALLOC( width_gauss_data, double, nb_chunk*2);
00766 XSH_MALLOC( offset_gauss_data, double, nb_chunk*2);
00767 XSH_MALLOC( xpos_gauss_data, double, nb_chunk*2);
00768 XSH_CALLOC( median_data, double, chunk_size);
00769
00770 check( chunk_vector = cpl_vector_new( ny));
00771 check( chunk_pos_vector = cpl_vector_new( ny));
00772
00773 xsh_msg_dbg_medium( "nx %d nb_chunks %d of size %d", nx, nb_chunk,
00774 chunk_size);
00775
00776 for( i=0; i< nb_chunk; i++){
00777 double x0, sigma, area, offset_val, height_val, fwhm;
00778 double med;
00779 int j;
00780 int ideb, ilast;
00781
00782 ideb = i*chunk_size;
00783 ilast = ideb+chunk_size;
00784
00785 for( j=0; j< ny; j++){
00786 int l, val_nb=0;
00787
00788 for( l=ideb; l< ilast; l++){
00789 val_nb++;
00790 median_data[l-ideb] = data[j*nx+l];
00791 }
00792 check( median_vect = cpl_vector_wrap( val_nb, median_data));
00793 check( med = cpl_vector_get_median( median_vect));
00794 check( cpl_vector_set( chunk_vector, j, med));
00795 check( cpl_vector_set( chunk_pos_vector, j, j));
00796 }
00797 #if REGDEBUG_OPTEXTRACT
00798 #if REGDEBUG_GAUSSIAN
00799
00800 {
00801 FILE* regdebug = NULL;
00802 char filename[256];
00803
00804 sprintf(filename, "gaussian_points_%d.dat",i);
00805 regdebug = fopen(filename,"w");
00806 for( j=0; j< ny; j++){
00807 double pos, val;
00808
00809 pos = cpl_vector_get( chunk_pos_vector,j);
00810 val = cpl_vector_get( chunk_vector,j);
00811
00812 fprintf(regdebug, "%f %f\n", pos, val);
00813 }
00814 fclose(regdebug);
00815 }
00816 #endif
00817 #endif
00818 cpl_vector_fit_gaussian( chunk_pos_vector, NULL, chunk_vector, NULL,
00819 CPL_FIT_ALL, &x0, &sigma, &area, &offset_val, NULL, NULL, NULL);
00820
00821 if ( cpl_error_get_code() == CPL_ERROR_NONE){
00822 height_val = area / sqrt(2*CPL_MATH_PI*sigma*sigma);
00823 fwhm = CPL_MATH_FWHM_SIG*sigma;
00824 if ( fwhm < (oversample*FIT_FWHM_LIMIT)){
00825 center_gauss_data[nbfit] = x0;
00826 height_gauss_data[nbfit] = height_val;
00827 width_gauss_data[nbfit] = sigma;
00828 xpos_gauss_data[nbfit] = ideb;
00829 offset_gauss_data[nbfit] = offset_val;
00830 nbfit++;
00831 center_gauss_data[nbfit] = x0;
00832 height_gauss_data[nbfit] = height_val;
00833 width_gauss_data[nbfit] = sigma;
00834 xpos_gauss_data[nbfit] = ilast;
00835 offset_gauss_data[nbfit] = offset_val;
00836 nbfit++;
00837 }
00838 else{
00839 xsh_msg("WARNING chunk %d-%d so big FWHM %f expected < %d",
00840 ideb, ilast, fwhm, oversample*FIT_FWHM_LIMIT);
00841 }
00842 }
00843 else {
00844 xsh_msg("WARNING chunk %d-%d don't converge with gaussian fit",
00845 ideb, ilast);
00846 xsh_error_reset();
00847 }
00848 }
00849
00850
00851 if ( nbfit == 0){
00852 xsh_msg("WARNING nbfit == 0 force poly degree to 0");
00853 nbfit = 1;
00854 deg_poly = 0;
00855 xpos_gauss_data[0] = 0;
00856 center_gauss_data[0] = ny/2.0;
00857 height_gauss_data[0] = 1;
00858 width_gauss_data[0] = ny;
00859 offset_gauss_data[0] = 0.0;
00860 }
00861 else if ( nbfit <= deg_poly){
00862 xsh_msg("WARNING %d (nbfit) < %d (poly degree) force poly degree to %d",
00863 nbfit, deg_poly, nbfit-1);
00864 deg_poly = nbfit-1;
00865 }
00866
00867 check(xpos_gauss = cpl_vector_wrap( nbfit, xpos_gauss_data));
00868 check( center_gauss = cpl_vector_wrap( nbfit, center_gauss_data));
00869 check( height_gauss = cpl_vector_wrap( nbfit, height_gauss_data));
00870 check( width_gauss = cpl_vector_wrap( nbfit, width_gauss_data));
00871 check( offset_gauss = cpl_vector_wrap( nbfit, offset_gauss_data));
00872
00873 check( *center = xsh_polynomial_fit_1d_create( xpos_gauss,
00874 center_gauss, deg_poly, NULL));
00875 check( *height = xsh_polynomial_fit_1d_create( xpos_gauss,
00876 height_gauss, deg_poly, NULL));
00877 check( *width = xsh_polynomial_fit_1d_create( xpos_gauss,
00878 width_gauss, deg_poly, NULL));
00879 check( *offset = xsh_polynomial_fit_1d_create( xpos_gauss,
00880 offset_gauss, deg_poly, NULL));
00881
00882 cleanup:
00883 if ( cpl_error_get_code() != CPL_ERROR_NONE){
00884 xsh_free_polynomial( center);
00885 xsh_free_polynomial( height);
00886 xsh_free_polynomial( width);
00887 xsh_free_polynomial( offset);
00888 }
00889 xsh_unwrap_vector( ¢er_gauss);
00890 xsh_unwrap_vector( &width_gauss);
00891 xsh_unwrap_vector( &height_gauss);
00892 xsh_unwrap_vector( &offset_gauss);
00893 xsh_unwrap_vector( &xpos_gauss);
00894 xsh_unwrap_vector( &median_vect);
00895 XSH_FREE( center_gauss_data);
00896 XSH_FREE( width_gauss_data);
00897 XSH_FREE( height_gauss_data);
00898 XSH_FREE( offset_gauss_data);
00899 XSH_FREE( xpos_gauss_data);
00900 XSH_FREE( median_data);
00901 xsh_free_vector( &chunk_vector);
00902 xsh_free_vector( &chunk_pos_vector);
00903 return;
00904 }
00905
00906
00907 static cpl_image* xsh_image_create_model_image( cpl_image *src_img,
00908 double *x_data, double *y_data, double kappa, int niter, double frac_min)
00909 {
00910 cpl_image *result = NULL;
00911 int nx,ny;
00912 int i,j;
00913 double *data = NULL;
00914 double *src_data = NULL;
00915 double *pos = NULL;
00916 double *line = NULL;
00917 cpl_vector *pos_vect = NULL;
00918 cpl_vector *line_vect = NULL;
00919 int size;
00920 cpl_polynomial *fit = NULL;
00921
00922
00923 XSH_ASSURE_NOT_NULL( src_img);
00924 check( nx = cpl_image_get_size_x( src_img));
00925 check( ny = cpl_image_get_size_y( src_img));
00926 check( result = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE));
00927 check( data = cpl_image_get_data_double( result));
00928 check( src_data = cpl_image_get_data_double( src_img));
00929 XSH_CALLOC( pos, double, nx);
00930 XSH_CALLOC( line, double, nx);
00931 double before=cpl_test_get_walltime();
00932 xsh_msg_debug("setting timing");
00933
00934 for(j=0; j<ny; j++){
00935 double val;
00936 double chisq;
00937 int *flags;
00938 int j_nx=j*nx;
00939 size = 0;
00940 for(i=0; i<nx; i++){
00941 double diffx, diffy;
00942 int pix=i+j_nx;
00943 diffx = x_data[pix]-(int)x_data[pix];
00944 if(diffx<0.2) {
00945 diffy = y_data[pix]-(int)y_data[pix];
00946
00947 if ( diffy < 0.2){
00948 pos[size] = i;
00949 line[size] = src_data[pix];
00950 size++;
00951 }
00952 }
00953 }
00954 check( pos_vect = cpl_vector_wrap( size, pos));
00955 check( line_vect = cpl_vector_wrap( size, line));
00956
00957 check( xsh_array_clip_poly1d( pos_vect, line_vect, kappa, niter,
00958 frac_min, 3, &fit, &chisq, &flags));
00959
00960 for(i=0; i<nx; i++){
00961 check( val = cpl_polynomial_eval_1d( fit, i, NULL));
00962 data[i+j_nx] = fabs(val);
00963 }
00964 xsh_unwrap_vector( &pos_vect);
00965 xsh_unwrap_vector( &line_vect);
00966 xsh_free_polynomial( &fit);
00967 XSH_FREE( flags);
00968 }
00969 double after=cpl_test_get_walltime();
00970 xsh_msg_debug("setting timing");
00971 xsh_msg("Time spent %g %g %g",before,after,after-before);
00972
00973 cleanup:
00974 if ( cpl_error_get_code() != CPL_ERROR_NONE){
00975 xsh_free_image( &result);
00976 xsh_unwrap_vector( &pos_vect);
00977 xsh_unwrap_vector( &line_vect);
00978 xsh_free_polynomial( &fit);
00979 }
00980 XSH_FREE( pos);
00981 XSH_FREE( line);
00982 return result;
00983 }
00984
00985 static cpl_image* xsh_image_create_gaussian_image( cpl_image *src,
00986 cpl_polynomial* centerp, cpl_polynomial* heightp, cpl_polynomial* widthp,
00987 cpl_polynomial* offsetp)
00988 {
00989 cpl_image *result = NULL;
00990 int nx,ny;
00991 int i,j;
00992 double *data = NULL;
00993
00994
00995 XSH_ASSURE_NOT_NULL( src);
00996 XSH_ASSURE_NOT_NULL( centerp);
00997 XSH_ASSURE_NOT_NULL( heightp);
00998 XSH_ASSURE_NOT_NULL( widthp);
00999 XSH_ASSURE_NOT_NULL( offsetp);
01000
01001 check( nx = cpl_image_get_size_x( src));
01002 check( ny = cpl_image_get_size_y( src));
01003 check( result = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE));
01004 check( data = cpl_image_get_data_double( result));
01005
01006 for(i=0; i<nx; i++){
01007 double x0, height, width, offset, z;
01008
01009
01010
01011 check( x0 = cpl_polynomial_eval_1d( centerp, i, NULL));
01012 check( height = cpl_polynomial_eval_1d( heightp, i, NULL));
01013 check( width = cpl_polynomial_eval_1d( widthp, i, NULL));
01014 check( offset = cpl_polynomial_eval_1d( offsetp, i, NULL));
01015
01016 double fct=width*XSH_MATH_SQRT_2;
01017 double inv_fct=1.0/fct;
01018
01019
01020
01021 for(j=0; j<ny; j++){
01022 z = (j-x0)*inv_fct;
01023
01024
01025 data[i+j*nx] = height*exp(-(z*z))+offset;
01026
01027
01028
01029
01030 }
01031
01032 }
01033
01034
01035 cleanup:
01036 if ( cpl_error_get_code() != CPL_ERROR_NONE){
01037 xsh_free_image( &result);
01038 }
01039 return result;
01040 }
01041
01058
01059
01060 static cpl_vector* xsh_image_extract_standard( cpl_image* img,
01061 cpl_image *err_img, cpl_image *qual_img,
01062 cpl_vector **err_v, cpl_vector **qual_v)
01063 {
01064 cpl_vector *result = NULL;
01065 int nx, ny, i;
01066 double *data = NULL;
01067 double *err = NULL;
01068 int *qual = NULL;
01069 double* pres=NULL;
01070 double* perr=NULL;
01071 double* pqua=NULL;
01072
01073
01074 XSH_ASSURE_NOT_NULL( img);
01075 XSH_ASSURE_NOT_NULL( err_img);
01076 XSH_ASSURE_NOT_NULL( qual_img);
01077 XSH_ASSURE_NOT_NULL( err_v);
01078 XSH_ASSURE_NOT_NULL( qual_v);
01079
01080
01081 check (nx = cpl_image_get_size_x( img));
01082 check (ny = cpl_image_get_size_y( img));
01083 check( data = cpl_image_get_data_double( img));
01084 check( err = cpl_image_get_data_double( err_img));
01085 check( qual = cpl_image_get_data_int( qual_img));
01086
01087
01088 check( result = cpl_vector_new( nx));
01089 check( *err_v = cpl_vector_new( nx));
01090 check( *qual_v = cpl_vector_new( nx));
01091 pres=cpl_vector_get_data(result);
01092 perr=cpl_vector_get_data(*err_v);
01093 pqua=cpl_vector_get_data(*qual_v);
01094
01095 for( i=0; i< nx; i++){
01096 int j;
01097 double pix_val = 0.0, err_val = 0.0;
01098 int qual_val = 0;
01099
01100 for( j=0; j< ny; j++){
01101 int pix=j*nx+i;
01102 pix_val += data[pix];
01103 err_val += err[pix]*err[pix];
01104 qual_val |= qual[pix];
01105 }
01106
01107 pres[i]=pix_val;
01108 perr[i]=sqrt(err_val);
01109 pqua[i]=(double)qual_val;
01110
01111 }
01112
01113 cleanup:
01114 if ( cpl_error_get_code() != CPL_ERROR_NONE){
01115 xsh_free_vector( &result);
01116 xsh_free_vector( err_v);
01117 xsh_free_vector( qual_v);
01118 }
01119 return result;
01120 }
01121
01122
01123
01140
01141 static cpl_image* xsh_image_divide_1D( cpl_image *image, cpl_image *err_img,
01142 cpl_vector *s1D, cpl_vector *err_s1D, cpl_image **err_2dby1D_img)
01143 {
01144 cpl_image *result = NULL;
01145 int i,j, nx, ny, size;
01146 double *data = NULL;
01147 double *err = NULL;
01148 double *data_res= NULL;
01149 double *err_res = NULL;
01150
01151
01152 XSH_ASSURE_NOT_NULL( image);
01153 XSH_ASSURE_NOT_NULL( err_img);
01154 XSH_ASSURE_NOT_NULL( s1D);
01155 XSH_ASSURE_NOT_NULL( err_s1D);
01156 XSH_ASSURE_NOT_NULL( err_2dby1D_img);
01157
01158 check( size = cpl_vector_get_size( s1D));
01159 check( nx = cpl_image_get_size_x( image));
01160 check( ny = cpl_image_get_size_y( image));
01161 XSH_ASSURE_NOT_ILLEGAL( size == nx);
01162
01163 check( data = cpl_image_get_data_double( image));
01164 check( err = cpl_image_get_data_double( err_img));
01165 check( result = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE));
01166 check( *err_2dby1D_img = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE));
01167 check( data_res = cpl_image_get_data_double( result));
01168 check( err_res = cpl_image_get_data_double( *err_2dby1D_img));
01169
01170 for(i=0; i< nx; i++){
01171 double d, d2, e2, d1, e1;
01172
01173 check( d2 = cpl_vector_get( s1D, i));
01174 check( e2 = cpl_vector_get( err_s1D, i));
01175
01176 if ( d2 == 0.0){
01177 if (i > 0 && i < (nx-1)){
01178 d2 = 0.5*(cpl_vector_get( s1D, i-1)+cpl_vector_get( s1D, i+1));
01179 e2 = 0.5*sqrt(cpl_vector_get( err_s1D, i-1)+cpl_vector_get( err_s1D, i+1));
01180 }
01181 }
01182 for( j=0; j <ny; j++){
01183 int pix=i+j*nx;
01184 d1 = data[pix];
01185 e1 = err[pix];
01186 d = d1/d2;
01187
01188 data_res[pix] = d;
01189 err_res[pix] = d * sqrt( (e1*e1)/(d1*d1) + (e2*e2)/(d2*d2) );
01190 }
01191 }
01192
01193 cleanup:
01194 if ( cpl_error_get_code() != CPL_ERROR_NONE){
01195 xsh_free_image( &result);
01196 xsh_free_image( err_2dby1D_img);
01197 }
01198 return result;
01199 }
01200
01201
01202
01203
01228
01229
01230 static cpl_vector* xsh_image_extract_optimal( cpl_image* img,
01231 cpl_image *errs_img, cpl_image *qual_img, xsh_opt_extract_param* opt_par,
01232 cpl_image *model_img,const int decode_bp,
01233 cpl_image** corr_img, cpl_vector **err_v, cpl_vector **qual_v)
01234 {
01235 int nb_rejected =-1, iter=1, clip_niter = 1, nbtotal_rejected=0;
01236 double clip_kappa = 3, clip_frac=0, clip_frac_rej =0.4;
01237 int *rejected = NULL;
01238 double *ab = NULL;
01239 cpl_vector *median_vect = NULL;
01240
01241 cpl_vector *result = NULL;
01242 int nx, ny, i;
01243 double *data = NULL;
01244 double *errs = NULL;
01245 int *qual_data = NULL;
01246 cpl_image *corr_image = NULL;
01247 double *corr_data = NULL;
01248 double *model_data = NULL;
01249
01250
01251 XSH_ASSURE_NOT_NULL( img);
01252 XSH_ASSURE_NOT_NULL( errs_img);
01253 XSH_ASSURE_NOT_NULL( qual_img);
01254 XSH_ASSURE_NOT_NULL( model_img);
01255 XSH_ASSURE_NOT_NULL( corr_img);
01256 XSH_ASSURE_NOT_NULL( err_v);
01257 XSH_ASSURE_NOT_NULL( qual_v);
01258 XSH_ASSURE_NOT_NULL( opt_par);
01259
01260
01261 clip_kappa = opt_par->clip_kappa;
01262 clip_niter = opt_par->clip_niter;
01263 clip_frac = opt_par->clip_frac;
01264
01265
01266 check (nx = cpl_image_get_size_x( img));
01267 check (ny = cpl_image_get_size_y( img));
01268 check( data = cpl_image_get_data_double( img));
01269 check( errs = cpl_image_get_data_double( errs_img));
01270 check( qual_data = cpl_image_get_data_int( qual_img));
01271
01272 check( corr_image = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE));
01273 check( corr_data = cpl_image_get_data_double( corr_image));
01274 check( model_data = cpl_image_get_data_double( model_img));
01275
01276
01277 check( result = cpl_vector_new( nx));
01278 check( *err_v = cpl_vector_new( nx));
01279 check( *qual_v = cpl_vector_new( nx));
01280
01281 XSH_CALLOC( rejected, int, ny);
01282 XSH_CALLOC( ab, double, nx);
01283
01284
01285 for( i=0; i< nx; i++){
01286 int j;
01287 double pix_val = 0.0;
01288 double f_err = 0.0;
01289 double f_num = 0.0;
01290 double f_denom = 0.0;
01291 double sum_MP = 0.0;
01292
01293 double err_val=0.0;
01294 int qual_val = 0;
01295
01296 int count = 0;
01297 double med_flux;
01298 double sum, sum_err;
01299 int sum_qual = 0;
01300
01301 for( j=0; j< ny; j++){
01302 rejected[j] = 1;
01303 }
01304
01305
01306
01307 count = 0;
01308 nbtotal_rejected=0;
01309 sum = 0;
01310 sum_err = 0;
01311
01312 for( j=0; j< ny; j++){
01313 int qual;
01314 double flux = 0.0, err_val = 0.0;
01315 double model_flux = 0.0;
01316
01317 flux = data[i+j*nx];
01318 err_val = errs[i+j*nx];
01319 qual = qual_data[i+j*nx];
01320 model_flux = model_data[i+j*nx];
01321
01322 if ( ((qual & decode_bp ) == 0) && model_flux != 0){
01323
01324 ab[count] = flux/model_flux;
01325 count++;
01326 qual_val |= qual;
01327 }
01328 else{
01329 rejected[j] =0;
01330 nbtotal_rejected++;
01331 sum += flux;
01332 sum_err += err_val*err_val;
01333 sum_qual |= qual;
01334 }
01335 }
01336 clip_frac = nbtotal_rejected/(double)ny;
01337 xsh_msg_dbg_medium("Bad pixel fraction %f for %d", clip_frac, i);
01338
01339 if (clip_frac < 1){
01340 check( median_vect = cpl_vector_wrap( count, ab));
01341 check( med_flux = cpl_vector_get_median( median_vect));
01342 xsh_unwrap_vector( &median_vect);
01343 }
01344 else{
01345 med_flux = sum;
01346 sum_err = sqrt( sum_err);
01347 }
01348
01349 iter=1;
01350 nb_rejected =-1;
01351
01352
01353 while( iter <= clip_niter && nb_rejected != 0 && count != 0
01354 && clip_frac <= clip_frac_rej){
01355 nb_rejected = 0;
01356 int nbdiff =0;
01357 double sigma;
01358
01359 xsh_msg_dbg_medium("Sigma clipping : Iteration %d/%d", iter, clip_niter);
01360
01361
01362 for( j=0; j< ny; j++){
01363 double flux = 0.0;
01364 double model_flux = 0.0;
01365 double diff = 0.0;
01366
01367 model_flux = model_data[i+j*nx];
01368 flux = data[i+j*nx];
01369
01370
01371 if ( rejected[j] == 1){
01372 diff = flux-med_flux*model_flux;
01373 ab[nbdiff] = diff;
01374 nbdiff++;
01375 }
01376 }
01377 check( median_vect = cpl_vector_wrap( nbdiff, ab));
01378 check( sigma = cpl_vector_get_stdev( median_vect));
01379 xsh_unwrap_vector( &median_vect);
01380 xsh_msg_dbg_medium( "Sigma %f", sigma);
01381
01382
01383 for( j=0; j< ny; j++){
01384 double flux = 0.0;
01385 double model_flux = 0.0;
01386 double diff = 0.0;
01387
01388 model_flux = model_data[i+j*nx];
01389 flux = data[i+j*nx];
01390 diff = flux-med_flux*model_flux;
01391
01392 if ( rejected[j] == 1 && fabs( diff) > clip_kappa*sigma){
01393 nb_rejected++;
01394 rejected[j] = 0;
01395 qual_val = QFLAG_COSMIC_RAY_REMOVED;
01396 }
01397 }
01398 nbtotal_rejected += nb_rejected;
01399
01400 count = 0;
01401
01402 for( j=0; j< ny; j++){
01403 double flux = 0.0;
01404 double model_flux = 0.0;
01405
01406 flux = data[i+j*nx];
01407 model_flux = model_data[i+j*nx];
01408
01409 if ( rejected[j] == 1 && model_flux != 0){
01410 ab[count] = flux/model_flux;
01411 count++;
01412 }
01413 }
01414 if ( count != 0){
01415 check( median_vect = cpl_vector_wrap( count, ab));
01416 check( med_flux = cpl_vector_get_median( median_vect));
01417 xsh_unwrap_vector( &median_vect);
01418 }
01419 clip_frac = nbtotal_rejected/(double)ny;
01420 iter++;
01421 }
01422 xsh_msg_dbg_medium("Fraction rejected %f", clip_frac);
01423
01424
01425 for( j=0; j< ny; j++){
01426 int qual;
01427
01428 qual = qual_data[i+j*nx];
01429 if ( ( (qual & decode_bp) == 0 ) && rejected[j] == 1 ){
01430 corr_data[i+j*nx] = data[i+j*nx];
01431 }
01432 else{
01433 corr_data[i+j*nx] = 0;
01434 }
01435 }
01436
01437
01438 if ( clip_frac < 1){
01439 for( j=0; j< ny; j++){
01440 double p, p2;
01441 double variance;
01442
01443 p = model_data[i+j*nx];
01444 p2 = p*p;
01445 variance = errs[i+j*nx]*errs[i+j*nx];
01446
01447 sum_MP += rejected[j]*p;
01448
01449 f_err += rejected[j]*p/variance;
01450
01451 f_num += (double)rejected[j]*p*corr_data[i+j*nx]/variance;
01452 f_denom += (double)rejected[j]*p2/variance;
01453 }
01454 pix_val = f_num / f_denom;
01455 err_val = sqrt( sum_MP / f_denom);
01456 }
01457 else{
01458 pix_val = sum;
01459 err_val = sum_err;
01460 }
01461 check( cpl_vector_set( result, i, pix_val));
01462 check( cpl_vector_set( *err_v, i, err_val));
01463 check( cpl_vector_set( *qual_v, i, (double)qual_val));
01464 }
01465
01466
01467 #if REGDEBUG_N
01468 if ( xsh_debug_level_get() >= XSH_DEBUG_LEVEL_MEDIUM) {
01469 char std_spectrum_name[256];
01470 FILE* std_spectrum_file = NULL;
01471
01472 sprintf( std_spectrum_name, "n.dat");
01473 std_spectrum_file = fopen( std_spectrum_name, "w");
01474
01475 fprintf( std_spectrum_file, "#index n\n");
01476
01477 for(i=0; i< nx; i++){
01478 fprintf(std_spectrum_file, "%d %f\n", i, n[i]);
01479 }
01480 fclose( std_spectrum_file);
01481
01482 sprintf( std_spectrum_name, "rejfrac.dat");
01483 std_spectrum_file = fopen( std_spectrum_name, "w");
01484
01485 fprintf( std_spectrum_file, "#index rejfrac\n");
01486
01487 for(i=0; i< nx; i++){
01488 fprintf(std_spectrum_file, "%d %f\n", i, rej[i]);
01489 }
01490 fclose( std_spectrum_file);
01491 }
01492 #endif
01493 *corr_img = corr_image;
01494
01495 cleanup:
01496 if ( cpl_error_get_code() != CPL_ERROR_NONE){
01497 xsh_free_vector( &result);
01498 xsh_free_vector( err_v);
01499 xsh_free_vector( qual_v);
01500 }
01501 XSH_FREE( rejected);
01502 XSH_FREE( ab);
01503 return result;
01504 }
01505
01506
01533
01534
01535 static cpl_vector* xsh_vector_integrate( int biny,
01536 int oversample, int absorder,
01537 cpl_vector *ref_pos, double step,
01538 cpl_vector *ref_values, cpl_vector *err_values, cpl_vector *qual_values,
01539 cpl_vector *new_pos,
01540 cpl_vector **spectrum_err, cpl_vector **spectrum_qual)
01541 {
01542
01543 int new_pos_size, ref_pos_size, pixel_pos_size;
01544 int i, j;
01545 cpl_vector *result=NULL;
01546 cpl_vector *pixel_pos_vect = NULL, *pixel_size_vect = NULL;
01547 cpl_polynomial *pixel_size_poly = NULL;
01548 double hstep;
01549 double lambda_bin_factor =0.0;
01550
01551 XSH_ASSURE_NOT_NULL( ref_pos);
01552 XSH_ASSURE_NOT_NULL( ref_values);
01553 XSH_ASSURE_NOT_NULL( err_values);
01554 XSH_ASSURE_NOT_NULL( qual_values);
01555 XSH_ASSURE_NOT_NULL( new_pos);
01556 XSH_ASSURE_NOT_NULL( spectrum_err);
01557 XSH_ASSURE_NOT_NULL( spectrum_qual);
01558
01559 check( new_pos_size = cpl_vector_get_size( new_pos));
01560 check( ref_pos_size = cpl_vector_get_size( ref_pos));
01561 check( result = cpl_vector_new( new_pos_size));
01562
01563 check( *spectrum_err = cpl_vector_new( new_pos_size));
01564 check( *spectrum_qual = cpl_vector_new( new_pos_size));
01565
01566
01567 pixel_pos_size = ref_pos_size-2*oversample;
01568 check( pixel_pos_vect = cpl_vector_new( pixel_pos_size));
01569 check( pixel_size_vect = cpl_vector_new( pixel_pos_size));
01570
01571 for( j=oversample; j< ref_pos_size-oversample; j++){
01572 double xA, xB, xC;
01573 double sizeAB, sizeBC, size;
01574
01575 check( xA = cpl_vector_get( ref_pos, j-oversample));
01576 check( xB = cpl_vector_get( ref_pos, j));
01577 check( xC = cpl_vector_get( ref_pos, j+oversample));
01578 sizeAB = xB-xA;
01579 sizeBC = xC-xB;
01580 size = (sizeAB+sizeBC)/2.0;
01581 check( cpl_vector_set( pixel_pos_vect, j-oversample, xB));
01582 check( cpl_vector_set( pixel_size_vect, j-oversample, size));
01583 }
01584
01585 check( pixel_size_poly = xsh_polynomial_fit_1d_create( pixel_pos_vect,
01586 pixel_size_vect, 1, NULL));
01587
01588 j =1;
01589 hstep = step/2.0;
01590
01591 for( i=0; i< new_pos_size; i++){
01592 double xA=0, xB=0, xC=0;
01593 double deltaA1=0, deltaA2=0;
01594
01595
01596 double size_min=0, pos_min=0, frac_min=0, sigy_min=0, flux_min=0;
01597 double size_max=0, pos_max=0, frac_max=0, sigy_max=0, flux_max=0;
01598
01599 double sig_min=0, sig_max=0;
01600 int qual_min=0, qual_max=0, qual=0;
01601
01602
01603 double y=0, ymin=0, ymax=0;
01604
01605 double xmin=0, xmax=0;
01606
01607 double x=0, sig_y = 26;
01608 int j_start=0, j_end=0, k=0;
01609
01610
01611 double flux=0, err=0, pixel_size_comp=0;
01612 double nb_x=0;
01613
01614 check( x= cpl_vector_get( new_pos, i));
01615 xmin = x-hstep;
01616 xmax = x+hstep;
01617 check( xA = cpl_vector_get( ref_pos, j-1));
01618 check( xB = cpl_vector_get( ref_pos, j));
01619 deltaA2 = (xB-xA)/2.0;
01620
01621 if ( j>1){
01622 check( xC = cpl_vector_get( ref_pos, j-2));
01623 deltaA1 = (xA-xC)/2.0;
01624 }
01625 else{
01626 deltaA1 = deltaA2;
01627 }
01628
01629 while( !(xmin >= (xA-deltaA1) && xmin <= (xA+deltaA2)) ){
01630 j++;
01631 check( xA = cpl_vector_get( ref_pos, j-1));
01632 check( xB = cpl_vector_get( ref_pos, j));
01633 check( xC = cpl_vector_get( ref_pos, j-2));
01634 deltaA1 = (xA-xC)/2.0;
01635 deltaA2 = (xB-xA)/2.0;
01636 }
01637
01638 j_start = j-1;
01639 size_min = deltaA1+deltaA2;
01640 pos_min = xA+deltaA2;
01641 frac_min = (pos_min-xmin)/size_min;
01642 check( ymin = cpl_vector_get( ref_values, j_start));
01643 check( sigy_min = cpl_vector_get( err_values, j_start));
01644 check( qual_min = cpl_vector_get( qual_values, j_start));
01645
01646 flux_min = frac_min*ymin;
01647 sig_min = frac_min*frac_min*sigy_min*sigy_min;
01648
01649 while( !(xmax >= (xA-deltaA1) && xmax <= (xA+deltaA2)) ){
01650 j++;
01651 check( xA = cpl_vector_get( ref_pos, j-1));
01652 check( xC = cpl_vector_get( ref_pos, j-2));
01653
01654 deltaA1 = (xA-xC)/2.0;
01655
01656 if ( j < ref_pos_size){
01657 check( xB = cpl_vector_get( ref_pos, j));
01658 deltaA2 = (xB-xA)/2.0;
01659 }
01660 else{
01661 deltaA2 = deltaA1;
01662 }
01663 }
01664
01665 j_end = j-1;
01666 j = j_end;
01667
01668 if ( j_start == j_end){
01669 flux_min =0;
01670 frac_min = 0;
01671 frac_max = (xmax-xmin)/size_min;
01672 }
01673 else{
01674 size_max = deltaA1+deltaA2;
01675 pos_max = xA-deltaA1;
01676 frac_max = (xmax-pos_max)/size_max;
01677 }
01678 check( ymax = cpl_vector_get( ref_values, j_end));
01679 check( sigy_max = cpl_vector_get( err_values, j_end));
01680 check( qual_max = cpl_vector_get( qual_values, j_end));
01681
01682 flux_max = frac_max*ymax;
01683 sig_max = frac_max*frac_max*sigy_max*sigy_max;
01684
01685 xsh_msg_dbg_high( "xmin %f pos_min %f diff %f size %f frac*size %f", xmin , pos_min, pos_min-xmin, size_min, frac_min*size_min);
01686 xsh_msg_dbg_high( "xmax %f pos_max %f diff %f size %f frac*size %f", xmax , pos_max, xmax-pos_max, size_max, frac_max*size_max);
01687
01688 nb_x = j_end-j_start-1;
01689 if ( nb_x > 0){
01690 nb_x += frac_min+frac_max;
01691 }
01692 else{
01693 nb_x = frac_min+frac_max;;
01694 }
01695
01696 qual = qual_min;
01697 qual |= qual_max;
01698 flux = flux_min+flux_max;
01699 err = sig_min+sig_max;
01700
01701 xsh_msg_dbg_high("nb_x %f size %f %f cont_min %f max_cont %f", nb_x, size_min, size_max, frac_min, frac_max);
01702
01703 for( k=j_start+1; k <= j_end-1; k++){
01704 check( y = cpl_vector_get( ref_values, k));
01705 check( sig_y = cpl_vector_get( err_values, k));
01706 flux +=y;
01707 err += sig_y*sig_y;
01708 qual |= (int)xsh_round_double(cpl_vector_get( qual_values, k));
01709 }
01710 err = sqrt( err);
01711
01712
01713 check( pixel_size_comp = cpl_polynomial_eval_1d( pixel_size_poly, x, NULL));
01714
01715 lambda_bin_factor = 1.0/(step*biny);
01716
01717 flux *= pixel_size_comp*lambda_bin_factor;
01718 err *= pixel_size_comp*lambda_bin_factor;
01719
01720 check( cpl_vector_set( result, i, flux));
01721 check( cpl_vector_set( *spectrum_err, i, err));
01722 check( cpl_vector_set( *spectrum_qual, i, qual));
01723 }
01724
01725 #if REGDEBUG_PIXELSIZE
01726 if ( xsh_debug_level_get() >= XSH_DEBUG_LEVEL_MEDIUM) {
01727 FILE* debug_file = NULL;
01728 char debug_name[256];
01729
01730 sprintf( debug_name, "pixel_size_%d.dat", absorder);
01731 xsh_msg("Create %s", debug_name);
01732
01733 debug_file = fopen( debug_name, "w");
01734
01735 fprintf( debug_file,"# wavelength size fit\n");
01736
01737 for( j=0; j< pixel_pos_size; j++){
01738 double pos, size, pol;
01739
01740 check( pos = cpl_vector_get( pixel_pos_vect, j));
01741 check( size = cpl_vector_get( pixel_size_vect, j));
01742 check( pol = cpl_polynomial_eval_1d( pixel_size_poly, pos, NULL));
01743 fprintf( debug_file ,"%f %f %f\n", pos, size, pol);
01744 }
01745 fclose( debug_file);
01746 }
01747 #endif
01748
01749 #if REGDEBUG_INTEGRATE
01750 if ( xsh_debug_level_get() >= XSH_DEBUG_LEVEL_MEDIUM) {
01751 FILE* debug_file = NULL;
01752 char debug_name[256];
01753
01754 sprintf( debug_name, "integrate_flux_%d.dat", absorder);
01755 xsh_msg("Create %s", debug_name);
01756
01757 debug_file = fopen( debug_name, "w");
01758
01759 fprintf( debug_file,"# wavelength flux err qual\n");
01760
01761 for( j=0; j< new_pos_size; j++){
01762 double xA, yA, sig_yA, qual_yA;
01763
01764 check( xA = cpl_vector_get( new_pos, j));
01765 check( yA = cpl_vector_get( result, j));
01766 check( sig_yA = cpl_vector_get( *spectrum_err, j));
01767 check( qual_yA = cpl_vector_get( *spectrum_qual, j));
01768 fprintf( debug_file ,"%f %f %f %f\n", xA, yA, sig_yA, qual_yA);
01769 }
01770
01771 fclose( debug_file);
01772 }
01773 #endif
01774
01775
01776 cleanup:
01777 xsh_free_vector( &pixel_pos_vect);
01778 xsh_free_vector( &pixel_size_vect);
01779 xsh_free_polynomial( &pixel_size_poly);
01780 if( cpl_error_get_code() != CPL_ERROR_NONE){
01781 xsh_free_vector( &result);
01782 xsh_free_vector( spectrum_err);
01783 xsh_free_vector( spectrum_qual);
01784 }
01785 return result;
01786 }
01787
01788
01789
01790
01831
01832 void xsh_opt_extract(cpl_frame *sci_frame, cpl_frame *orderlist_frame,
01833 cpl_frame *wavesol_frame, cpl_frame *model_frame, cpl_frame *wavemap_frame,
01834 cpl_frame *slitmap_frame, cpl_frame *loc_frame,
01835 cpl_frame *spectralformat_frame, cpl_frame *masterflat_frame,
01836 xsh_instrument *instrument, xsh_opt_extract_param* opt_extract_par,
01837 const char* rec_prefix, cpl_frame** orderext1d_frame,
01838 cpl_frame** orderoxt1d_frame, cpl_frame** orderoxt1d_eso_frame,
01839 cpl_frame** qc_subextract_frame, cpl_frame** qc_s2ddiv1d_frame,
01840 cpl_frame** qc_model_frame, cpl_frame** qc_weight_frame)
01841
01842 {
01843
01844 check(
01845 xsh_opt_extract_orders( sci_frame, orderlist_frame, wavesol_frame,
01846 model_frame, wavemap_frame, slitmap_frame, loc_frame,spectralformat_frame,
01847 masterflat_frame, instrument, opt_extract_par, 0, 100, rec_prefix,
01848 orderext1d_frame, orderoxt1d_frame, orderoxt1d_eso_frame, qc_subextract_frame,
01849 qc_s2ddiv1d_frame, qc_model_frame, qc_weight_frame));
01850
01851 cleanup: return;
01852 }
01853
01854
01897
01898 void
01899 xsh_opt_extract_orders( cpl_frame *sci_frame,
01900 cpl_frame *orderlist_frame,
01901 cpl_frame *wavesol_frame,
01902 cpl_frame *model_frame,
01903 cpl_frame *wavemap_frame,
01904 cpl_frame *slitmap_frame,
01905 cpl_frame *loc_frame,
01906 cpl_frame *spectralformat_frame,
01907 cpl_frame *masterflat_frame,
01908 xsh_instrument *instrument,
01909 xsh_opt_extract_param* opt_extract_par,
01910 int min_index, int max_index,
01911 const char* rec_prefix,
01912 cpl_frame** orderext1d_frame,
01913 cpl_frame** orderoxt1d_frame,
01914 cpl_frame** res_frame_ext,
01915 cpl_frame** qc_subextract_frame,
01916 cpl_frame** qc_s2ddiv1d_frame,
01917 cpl_frame** qc_model_frame,
01918 cpl_frame** qc_weight_frame)
01919 {
01920 xsh_pre *sci_pre = NULL;
01921 xsh_order_list *order_list = NULL;
01922 cpl_image *blaze_img = NULL;
01923 float *sci_pre_data = NULL;
01924 float *sci_pre_errs = NULL;
01925 int *sci_pre_qual = NULL;
01926 xsh_wavesol *wavesol = NULL;
01927 xsh_spectralformat_list *spectralformat_list= NULL;
01928
01929 int iorder, i;
01930 double *slits_cen = NULL;
01931 double *lambdas = NULL;
01932 double *pos_cen_x = NULL;
01933 double *pos_cen_y = NULL;
01934 int nlambdas = 0;
01935 cpl_vector *extract_lambda_vect = NULL;
01936 cpl_vector *result_lambda_vect = NULL;
01937 int result_lambda_vect_size;
01938
01939 cpl_vector *div_flux_vect = NULL;
01940 cpl_vector *div_err_vect = NULL;
01941
01942 cpl_vector *std_flux_vect = NULL;
01943 cpl_vector *std_err_vect = NULL;
01944 cpl_vector *std_qual_vect = NULL;
01945
01946 cpl_vector *opt_flux_vect = NULL;
01947 cpl_vector *opt_err_vect = NULL;
01948 cpl_vector *opt_qual_vect = NULL;
01949
01950 cpl_image *extract_img = NULL;
01951 cpl_image *extract_errs_img = NULL;
01952 cpl_image *extract_qual_img = NULL;
01953
01954 cpl_image *x_img = NULL;
01955 cpl_image *y_img = NULL;
01956
01957 cpl_image *sub_extract_img = NULL;
01958 cpl_image *sub_extract_errs_img = NULL;
01959 cpl_image *sub_extract_qual_img = NULL;
01960
01961 cpl_image *sub_x_img = NULL;
01962 cpl_image *sub_y_img = NULL;
01963
01964 cpl_image *s2Dby1D_img = NULL;
01965 cpl_image *s2Dby1D_err_img = NULL;
01966 cpl_image *weight_img = NULL;
01967 double *extract_data = NULL;
01968 double *extract_errs = NULL;
01969 int *extract_qual = NULL;
01970 double *extract_x_data = NULL;
01971 double *extract_y_data = NULL;
01972
01973 int nx_extract, ny_extract;
01974 int box_hsize;
01975 int niter=1, iiter;
01976 int chunk_ovsamp_size;
01977 int extract_slit_hsize=0;
01978 int deg_poly = 3;
01979 cpl_vector* std_extract_vect = NULL;
01980 cpl_vector* std_err_extract_vect = NULL;
01981 cpl_vector* std_qual_extract_vect = NULL;
01982
01983 cpl_vector *opt_extract_vect = NULL;
01984 cpl_vector *opt_err_extract_vect = NULL;
01985 cpl_vector *opt_qual_extract_vect = NULL;
01986
01987 double lambda_step=0.0;
01988 xsh_rec_list *orderext1d = NULL;
01989 xsh_rec_list *orderoxt1d = NULL;
01990 char tag_drl[256];
01991 char tag[256];
01992
01993 char filename[256];
01994 char filename_drl[256];
01995
01996 cpl_image *model_img = NULL;
01997 int rec_min_index, rec_max_index;
01998 double slit_down, slit_up;
01999 double lambda_min, lambda_max;
02000 int blaze_shift;
02001 xsh_xs_3 config_model;
02002 double conad = 0.0;
02003 double square_oversample;
02004
02005 int biny;
02006
02007
02008 char qc_extname[256];
02009 cpl_propertylist *qc_header = NULL;
02010 char qc_subextract_name[256];
02011 char qc_s2ddiv1d_name[256];
02012 char qc_model_name[256];
02013 char qc_weight_name[256];
02014 int mode = CPL_IO_DEFAULT;
02015
02016
02017 int decode_bp=instrument->decode_bp;
02018
02019 XSH_ASSURE_NOT_NULL( sci_frame);
02020 XSH_ASSURE_NOT_NULL( orderlist_frame);
02021 XSH_ASSURE_NOT_NULL( wavemap_frame);
02022 XSH_ASSURE_NOT_NULL( slitmap_frame);
02023 XSH_ASSURE_NOT_NULL( loc_frame);
02024 XSH_ASSURE_NOT_NULL( spectralformat_frame);
02025 XSH_ASSURE_NOT_NULL( masterflat_frame);
02026 XSH_ASSURE_NOT_NULL( instrument);
02027 XSH_ASSURE_NOT_NULL( opt_extract_par);
02028 XSH_ASSURE_NOT_NULL( orderext1d_frame);
02029 XSH_ASSURE_NOT_NULL( orderoxt1d_frame);
02030
02031
02032 check( sci_pre = xsh_pre_load( sci_frame, instrument));
02033
02034
02035 square_oversample = opt_extract_par->oversample*opt_extract_par->oversample;
02036
02037 conad = sci_pre->conad;
02038
02039
02040 biny = instrument->biny;
02041 if (model_frame == NULL){
02042 XSH_ASSURE_NOT_NULL( wavesol_frame);
02043 check( wavesol = xsh_wavesol_load( wavesol_frame, instrument));
02044 }
02045 else{
02046 check( xsh_model_config_load_best( model_frame, &config_model));
02047 check( xsh_model_binxy( &config_model, instrument->binx, instrument->biny));
02048 }
02049
02050 lambda_step = opt_extract_par->lambda_step;
02051 niter = opt_extract_par->niter;
02052
02053 check( sci_pre_data = cpl_image_get_data_float( sci_pre->data));
02054 check( sci_pre_errs = cpl_image_get_data_float( sci_pre->errs));
02055 check( sci_pre_qual = cpl_image_get_data_int( sci_pre->qual));
02056
02057 check( order_list = xsh_order_list_load( orderlist_frame, instrument));
02058
02059 check( spectralformat_list = xsh_spectralformat_list_load(
02060 spectralformat_frame, instrument));
02061
02062 nx_extract = sci_pre->ny*opt_extract_par->oversample;
02063 ny_extract = (OPT_EXTRACT_SLIT_SIZE * opt_extract_par->oversample)+1;
02064
02065 XSH_MALLOC( lambdas, double, nx_extract);
02066 XSH_MALLOC( slits_cen, double, nx_extract);
02067 XSH_MALLOC( pos_cen_x, double, nx_extract);
02068 XSH_MALLOC( pos_cen_y, double, nx_extract);
02069 XSH_MALLOC( extract_data, double, nx_extract*ny_extract);
02070 XSH_MALLOC( extract_errs, double, nx_extract*ny_extract);
02071 XSH_MALLOC( extract_qual, int, nx_extract*ny_extract);
02072
02073 XSH_MALLOC( extract_x_data, double, nx_extract*ny_extract);
02074 XSH_MALLOC( extract_y_data, double, nx_extract*ny_extract);
02075
02076 chunk_ovsamp_size = opt_extract_par->chunk_size*opt_extract_par->oversample;
02077 box_hsize = opt_extract_par->box_hsize*opt_extract_par->oversample;
02078
02079 xsh_msg( "Create and multiply by blaze function");
02080
02081 xsh_msg_dbg_low("---opt_extract : create blaze");
02082 check( blaze_img = xsh_create_blaze( masterflat_frame,
02083 order_list, instrument));
02084
02085 xsh_msg_dbg_low("---opt_extract : multiply by blaze");
02086 check( xsh_pre_multiply_image( sci_pre, blaze_img));
02087
02088 #if REGDEBUG_BLAZE
02089 if ( xsh_debug_level_get() >= XSH_DEBUG_LEVEL_MEDIUM) {
02090 cpl_frame *norm_frame = NULL;
02091
02092 check( norm_frame = xsh_pre_save( sci_pre,"OPTEXT_MULT_BLAZE.fits",
02093 "OPTEXT_MULT_BLAZE",1));
02094 xsh_free_frame( &norm_frame);
02095 }
02096 #endif
02097
02098
02099 check( orderext1d = xsh_rec_list_create( instrument));
02100 check( orderoxt1d = xsh_rec_list_create( instrument));
02101
02102 if (min_index >= 0){
02103 rec_min_index = min_index;
02104 }
02105 else{
02106 rec_min_index = 0;
02107 }
02108 if ( max_index < orderext1d->size){
02109 rec_max_index = max_index;
02110 }
02111 else{
02112 rec_max_index = orderext1d->size-1;
02113 }
02114
02115 check( xsh_get_slit_edges( slitmap_frame, &slit_down, &slit_up,
02116 NULL, NULL, instrument));
02117
02118
02119 XSH_NEW_PROPERTYLIST( qc_header);
02120
02121 for( iorder= rec_min_index; iorder <= rec_max_index; iorder++) {
02122 int abs_order,start_y, end_y;
02123 double cos_tilt=0, sin_tilt=0;
02124 int ny_extract_min=-1, ny_extract_max = -1;
02125
02126
02127
02128 abs_order = order_list->list[iorder].absorder;
02129 check( start_y = xsh_order_list_get_starty( order_list, iorder));
02130 check( end_y = xsh_order_list_get_endy( order_list, iorder));
02131 xsh_msg("Extracting order %d", abs_order);
02132 xsh_msg_dbg_low("---opt_extract : order %d [%d --> %d]", abs_order, start_y, end_y);
02133
02134 check( xsh_wavemap_lambda_range( wavemap_frame, slitmap_frame, start_y, end_y,
02135 opt_extract_par->oversample, order_list, iorder,
02136 pos_cen_x, pos_cen_y, abs_order, spectralformat_list, lambdas, slits_cen,
02137 &nlambdas, instrument));
02138
02139 #if REGDEBUG_OPTEXTRACT
02140
02141 {
02142 FILE *test_file;
02143 char test_name[256];
02144 float x, y, lambda, slit;
02145 double x_slit_min=0, y_slit_min=0;
02146 double x_slit_max=0, y_slit_max=0;
02147 double x_slit_cen=0, y_slit_cen=0;
02148 double diffx, diffy;
02149
02150 if (model_frame == NULL){
02151 sprintf( test_name, "poly%dx%d_lambda_range_%d.dat", binx, biny, abs_order);
02152 }
02153 else{
02154 sprintf( test_name, "model%dx%d_lambda_range_%d.dat", binx, biny, abs_order);
02155 }
02156 xsh_msg("Write %s", test_name);
02157 test_file = fopen( test_name, "w");
02158 fprintf( test_file, "# x y lambda slit slit_min slit_max x_slit_min y_slit_min x_slit_cen y_slit_cen x_slit_max y_slit_max diffx diffy\n");
02159
02160 for(i=0; i< nlambdas; i++){
02161 x = pos_cen_x[i];
02162 y = pos_cen_y[i];
02163 lambda = lambdas[i];
02164 slit = slits_cen[i];
02165
02166 if (model_frame == NULL){
02167 check( x_slit_min = xsh_wavesol_eval_polx( wavesol, lambda,(double) abs_order,
02168 slit_down));
02169 check( x_slit_max = xsh_wavesol_eval_polx( wavesol, lambda, (double)abs_order,
02170 slit_up));
02171 check( y_slit_min = xsh_wavesol_eval_poly( wavesol, lambda, (double)abs_order,
02172 slit_down));
02173 check( y_slit_max = xsh_wavesol_eval_poly( wavesol, lambda, (double)abs_order,
02174 slit_up));
02175 check( x_slit_cen = xsh_wavesol_eval_polx( wavesol, lambda,(double) abs_order,
02176 slit));
02177 check( y_slit_cen = xsh_wavesol_eval_poly( wavesol, lambda, (double)abs_order,
02178 slit));
02179 }
02180 else{
02181 check( xsh_model_get_xy( &config_model, instrument, lambda, (double) abs_order,
02182 slit_down, &x_slit_min, &y_slit_min));
02183 check( xsh_model_get_xy( &config_model, instrument, lambda, (double) abs_order,
02184 slit_up, &x_slit_max, &y_slit_max));
02185 check( xsh_model_get_xy( &config_model, instrument, lambda, (double) abs_order,
02186 slit, &x_slit_cen, &y_slit_cen));
02187 }
02188 diffx = x-x_slit_cen;
02189 diffy = y-y_slit_cen;
02190
02191 fprintf( test_file, "%f %f %f %f %f %f %f %f %f %f %f %f %f %f\n", x, y, lambda, slit, slit_down, slit_up,
02192 x_slit_min, y_slit_min, x_slit_cen, y_slit_cen, x_slit_max, y_slit_max, diffx, diffy);
02193 }
02194 fclose( test_file);
02195 }
02196 #endif
02197
02198
02199
02200 xsh_msg_dbg_low("---opt_extract : (%f --> %f) in %d", lambdas[0], lambdas[nlambdas-1],
02201 nlambdas);
02202
02203 xsh_msg_dbg_low("---opt_extract : s (%f --> %f)", slit_down, slit_up);
02204
02205 for(i=0; i< nlambdas; i++){
02206 float x, y, lambda, slit;
02207 double x_slit_min=0, y_slit_min=0;
02208 double x_slit_max=0, y_slit_max=0;
02209 float x_min, x_max;
02210 float y_min, y_max;
02211 int l;
02212
02213 x = pos_cen_x[i];
02214 y = pos_cen_y[i];
02215 lambda = lambdas[i];
02216 slit = slits_cen[i];
02217 xsh_msg_dbg_high(" x y lambda slit => %f %f %f %f", x, y, lambda, slit);
02218
02219 if ( ( i%opt_extract_par->oversample) == 0){
02220 xsh_msg_dbg_high("evaluate tilt lambda order %f %f", lambda, (double) abs_order);
02221
02222 if (model_frame == NULL){
02223 check( x_slit_min = xsh_wavesol_eval_polx( wavesol, lambda,(double) abs_order,
02224 slit_down));
02225 check( x_slit_max = xsh_wavesol_eval_polx( wavesol, lambda, (double)abs_order,
02226 slit_up));
02227 check( y_slit_min = xsh_wavesol_eval_poly( wavesol, lambda, (double)abs_order,
02228 slit_down));
02229 check( y_slit_max = xsh_wavesol_eval_poly( wavesol, lambda, (double)abs_order,
02230 slit_up));
02231 }
02232 else{
02233
02234 check( xsh_model_get_xy( &config_model, instrument, lambda, (double) abs_order,
02235 slit_down, &x_slit_min, &y_slit_min));
02236 check( xsh_model_get_xy( &config_model, instrument, lambda, (double) abs_order,
02237 slit_up, &x_slit_max, &y_slit_max));
02238 }
02239 if ( x_slit_min < x_slit_max){
02240 x_min = x_slit_min;
02241 x_max = x_slit_max;
02242 }
02243 else{
02244 x_min = x_slit_max;
02245 x_max = x_slit_min;
02246 }
02247 if ( y_slit_min < y_slit_max){
02248 y_min = y_slit_min;
02249 y_max = y_slit_max;
02250 }
02251 else{
02252 y_min = y_slit_max;
02253 y_max = y_slit_min;
02254 }
02255 xsh_msg_dbg_high("min(%f,%f) max(%f,%f)", x_min, y_min, x_max, y_max);
02256
02257 if ( x_min >= 1 && x_max <= sci_pre->nx && (x_min <= x) &&
02258 ( x <= x_max) &&
02259 y_min >= 1 && y_max <= sci_pre->ny){
02260 double tilt_x, tilt_y;
02261 double hypotenuse_tilt;
02262
02263 if ((y_min >= y) || (y >= y_max)){
02264 xsh_msg_warning("y out of bounds : %f [%f,%f]", y, y_min,y_max);
02265 }
02266 tilt_x = x_slit_max-x_slit_min;
02267 tilt_y = y_slit_max-y_slit_min;
02268 xsh_msg_dbg_high("tilt_x %f tilt_y %f", tilt_x, tilt_y);
02269
02270 if ( extract_slit_hsize <= 0){
02271 extract_slit_hsize = (int)(fabs(tilt_x)/2.0*SLIT_USEFUL_WINDOW_FACTOR)+1;
02272 ny_extract = extract_slit_hsize*opt_extract_par->oversample*2+1;
02273 xsh_msg_dbg_high("SLIT hsize %d global %d", extract_slit_hsize,ny_extract);
02274 }
02275 hypotenuse_tilt = sqrt( tilt_y*tilt_y+tilt_x*tilt_x);
02276 cos_tilt = tilt_x/hypotenuse_tilt;
02277 sin_tilt = tilt_y/hypotenuse_tilt;
02278 }
02279 }
02280
02281 for (l=0; l< ny_extract; l++){
02282 float xs, fy, fx;
02283 double flux=-1, err=-1;
02284 int qual_val=0;
02285
02286 xs = (float)l/(float)opt_extract_par->oversample-
02287 extract_slit_hsize;
02288
02289
02290 fy = xs*sin_tilt+y-1;
02291
02292 if ( (fy >= 0) && (fy < sci_pre->ny) ){
02293 int pixel=i+l*nlambdas;
02294
02295 fx = xs*cos_tilt+x-1;
02296
02297 check(xsh_interpolate_linear( sci_pre_data, sci_pre_errs, sci_pre_qual, sci_pre->nx,
02298 sci_pre->ny, fx, fy, &flux, &err, &qual_val, FLUX_MODE));
02299
02300 extract_x_data[pixel] = fx;
02301 extract_y_data[pixel] = fy;
02302
02303
02304 flux *= conad/square_oversample;
02305 extract_data[pixel] = flux;
02306 extract_errs[pixel] = err * conad /
02307 opt_extract_par->oversample;
02308 extract_qual[pixel] = qual_val;
02309 }
02310 else{
02311 break;
02312 }
02313 }
02314 }
02315
02316
02317
02318 check( xsh_object_localize( slitmap_frame, loc_frame, opt_extract_par->oversample,
02319 box_hsize, nlambdas,
02320 ny_extract, extract_x_data,
02321 extract_y_data, &ny_extract_min, &ny_extract_max));
02322
02323 check( extract_img = cpl_image_wrap_double( nlambdas,
02324 ny_extract, extract_data));
02325 check( extract_errs_img = cpl_image_wrap_double( nlambdas,
02326 ny_extract, extract_errs));
02327 check( extract_qual_img = cpl_image_wrap_int( nlambdas,
02328 ny_extract, extract_qual));
02329
02330 #if REGDEBUG_EXTRACT
02331 { char name[256];
02332 sprintf( name, "extract_%d.fits", abs_order);
02333 cpl_image_save( extract_img, name, CPL_BPP_IEEE_DOUBLE,
02334 NULL,CPL_IO_DEFAULT);
02335 xsh_msg("EXTRACT save %s",name);
02336 sprintf( name, "extract_errs_%d.fits", abs_order);
02337 cpl_image_save( extract_errs_img, name, CPL_BPP_IEEE_DOUBLE,
02338 NULL,CPL_IO_DEFAULT);
02339 sprintf( name, "extract_qual_%d.fits", abs_order);
02340 cpl_image_save( extract_qual_img, name, XSH_PRE_QUAL_BPP,
02341 NULL,CPL_IO_DEFAULT);
02342 }
02343 #endif
02344
02345
02346 check( x_img = cpl_image_wrap_double( nlambdas,
02347 ny_extract, extract_x_data));
02348 check( y_img = cpl_image_wrap_double( nlambdas,
02349 ny_extract, extract_y_data));
02350
02351 #if REGDEBUG_EXTRACT_XY
02352 { char name[256];
02353
02354 sprintf( name, "extract_x_%d.fits", abs_order);
02355 xsh_msg("EXTRACT_XY save %s",name);
02356 cpl_image_save( x_img, name, CPL_BPP_IEEE_DOUBLE, NULL,CPL_IO_DEFAULT);
02357 sprintf( name, "extract_y_%d.fits", abs_order);
02358 xsh_msg("EXTRACT_XY save %s",name);
02359 cpl_image_save( y_img, name, CPL_BPP_IEEE_DOUBLE, NULL,CPL_IO_DEFAULT);
02360 }
02361 #endif
02362
02363
02364 check( sub_extract_img = cpl_image_extract( extract_img, 1,
02365 ny_extract_min+1, nlambdas, ny_extract_max+1));
02366 check( sub_extract_errs_img = cpl_image_extract( extract_errs_img, 1,
02367 ny_extract_min+1, nlambdas, ny_extract_max+1));
02368 check( sub_extract_qual_img = cpl_image_extract( extract_qual_img, 1,
02369 ny_extract_min+1, nlambdas, ny_extract_max+1));
02370
02371
02372
02373 sprintf( qc_subextract_name, "sub_extract.fits");
02374
02375 sprintf( qc_extname, "ORD%d_FLUX", abs_order);
02376 check( xsh_pfits_set_extname ( qc_header, qc_extname));
02377 cpl_image_save( sub_extract_img, qc_subextract_name, CPL_BPP_IEEE_DOUBLE,
02378 qc_header,mode);
02379
02380 sprintf( qc_extname, "ORD%d_ERRS", abs_order);
02381 check( xsh_pfits_set_extname ( qc_header, qc_extname));
02382 cpl_image_save( sub_extract_errs_img, qc_subextract_name, CPL_BPP_IEEE_DOUBLE,
02383 qc_header, CPL_IO_EXTEND);
02384
02385 sprintf( qc_extname, "ORD%d_QUAL", abs_order);
02386 check( xsh_pfits_set_extname ( qc_header, qc_extname));
02387 cpl_image_save( sub_extract_qual_img, qc_subextract_name, XSH_PRE_QUAL_BPP,
02388 qc_header, CPL_IO_EXTEND);
02389
02390
02391
02392 check( sub_x_img = cpl_image_extract( x_img, 1,
02393 ny_extract_min+1, nlambdas, ny_extract_max+1));
02394 check( sub_y_img = cpl_image_extract( y_img, 1,
02395 ny_extract_min+1, nlambdas, ny_extract_max+1));
02396
02397
02398 #if REGDEBUG_SUBEXTRACT_XY
02399 { char name[256];
02400
02401 sprintf( name, "sub_x_%d.fits", abs_order);
02402 xsh_msg("EXTRACT_XY save %s",name);
02403 cpl_image_save( sub_x_img, name, CPL_BPP_IEEE_DOUBLE, NULL,CPL_IO_DEFAULT);
02404 sprintf( name, "sub_y_%d.fits", abs_order);
02405 xsh_msg("EXTRACT_XY save %s",name);
02406 cpl_image_save( sub_y_img, name, CPL_BPP_IEEE_DOUBLE, NULL,CPL_IO_DEFAULT);
02407 }
02408 #endif
02409
02410 check( extract_lambda_vect = cpl_vector_wrap( nlambdas, lambdas));
02411
02412
02413 check( std_extract_vect = xsh_image_extract_standard( sub_extract_img,
02414 sub_extract_errs_img, sub_extract_qual_img, &std_err_extract_vect,
02415 &std_qual_extract_vect));
02416
02417 #if REGDEBUG_EXTRACT_STD
02418 {
02419 char std_spectrum_name[256];
02420 FILE* std_spectrum_file = NULL;
02421 int nx = cpl_vector_get_size( std_extract_vect);
02422
02423 sprintf( std_spectrum_name, "extract_std_%d.dat", abs_order);
02424 std_spectrum_file = fopen( std_spectrum_name, "w");
02425
02426 fprintf(std_spectrum_file, "#index lambda flux err qual\n");
02427
02428 for(i=0; i< nx; i++){
02429 double lambdav;
02430 double fluxv;
02431 double lambda_err;
02432 double qual_v;
02433
02434 lambdav = cpl_vector_get( extract_lambda_vect, i);
02435 fluxv = cpl_vector_get( std_extract_vect, i);
02436 lambda_err = cpl_vector_get( std_err_extract_vect, i);
02437 qual_v = cpl_vector_get( std_qual_extract_vect, i);
02438 fprintf(std_spectrum_file, "%d %f %f %f %f\n", i, lambdav, fluxv,
02439 lambda_err, qual_v);
02440 }
02441 fclose( std_spectrum_file);
02442 }
02443 #endif
02444
02445
02446 div_flux_vect = std_extract_vect;
02447 div_err_vect = std_err_extract_vect;
02448
02449 for( iiter=1; iiter <= niter; iiter++){
02450 double kappa = 3.0, frac_min=0.7;
02451 int niter = 2;
02452
02453 xsh_free_image( &s2Dby1D_img);
02454 xsh_free_image( &s2Dby1D_err_img);
02455 xsh_free_image( &model_img);
02456 xsh_free_image( &weight_img);
02457
02458 xsh_msg_dbg_low("---opt_extract : Do optimal extraction %d / %d", iiter, niter);
02459
02460 check( s2Dby1D_img = xsh_image_divide_1D( sub_extract_img,
02461 sub_extract_errs_img, div_flux_vect, div_err_vect,
02462 &s2Dby1D_err_img));
02463
02464 sprintf( qc_s2ddiv1d_name, "s2Dby1D.fits");
02465 sprintf( qc_extname, "ORD%d_FLUX", abs_order);
02466 check( xsh_pfits_set_extname ( qc_header, qc_extname));
02467 cpl_image_save( s2Dby1D_img, qc_s2ddiv1d_name, CPL_BPP_IEEE_DOUBLE, qc_header,
02468 mode);
02469 sprintf( qc_extname, "ORD%d_ERRS", abs_order);
02470 cpl_image_save( s2Dby1D_err_img, qc_s2ddiv1d_name, CPL_BPP_IEEE_DOUBLE, qc_header,
02471 CPL_IO_EXTEND);
02472
02473 check( model_img = xsh_optextract_produce_model( s2Dby1D_img,
02474 opt_extract_par->method, chunk_ovsamp_size, deg_poly,
02475 opt_extract_par->oversample, extract_x_data, extract_y_data,
02476 abs_order, kappa, niter, frac_min));
02477
02478 sprintf( qc_model_name, "model.fits");
02479 sprintf( qc_extname, "ORD%d_FLUX", abs_order);
02480 check( xsh_pfits_set_extname ( qc_header, qc_extname));
02481 cpl_image_save( model_img, qc_model_name, CPL_BPP_IEEE_DOUBLE, qc_header, mode);
02482
02483
02484 xsh_free_vector( &opt_extract_vect);
02485 xsh_free_vector( &opt_err_extract_vect);
02486 xsh_free_vector( &opt_qual_extract_vect);
02487 check( opt_extract_vect = xsh_image_extract_optimal( sub_extract_img,
02488 sub_extract_errs_img, sub_extract_qual_img, opt_extract_par,
02489 model_img, decode_bp,&weight_img, &opt_err_extract_vect, &opt_qual_extract_vect));
02490 div_flux_vect = opt_extract_vect;
02491 div_err_vect = opt_err_extract_vect;
02492 }
02493
02494 sprintf( qc_weight_name, "weight.fits");
02495 sprintf( qc_extname, "ORD%d_FLUX", abs_order);
02496 check( xsh_pfits_set_extname ( qc_header, qc_extname));
02497 cpl_image_save( weight_img, qc_weight_name, CPL_BPP_IEEE_DOUBLE, qc_header,
02498 mode);
02499
02500
02501 if ( xsh_instrument_get_arm(instrument) != XSH_ARM_NIR){
02502 blaze_shift = (int)pos_cen_y[0];
02503 }
02504 else{
02505 int nx = cpl_vector_get_size( std_extract_vect);
02506
02507 blaze_shift = (int)pos_cen_y[nx-1];
02508 }
02509
02510 check( xsh_vector_divide_poly( std_extract_vect,
02511 opt_extract_par->oversample,
02512 order_list->list[iorder].blazepoly, blaze_shift, instrument));
02513
02514 check( xsh_vector_divide_poly( std_err_extract_vect,
02515 opt_extract_par->oversample,
02516 order_list->list[iorder].blazepoly, blaze_shift, instrument));
02517
02518 #if REGDEBUG_BLAZE
02519 if ( xsh_debug_level_get() >= XSH_DEBUG_LEVEL_MEDIUM) {
02520 char std_spectrum_name[256];
02521 FILE* std_spectrum_file = NULL;
02522 int nx = cpl_vector_get_size( std_extract_vect);
02523
02524 sprintf( std_spectrum_name, "spectrum_divbyblaze_std_%d.dat", abs_order);
02525 std_spectrum_file = fopen( std_spectrum_name, "w");
02526
02527 fprintf( std_spectrum_file, "#index lambda flux err qual\n");
02528 for(i=0; i< nx; i++){
02529 double lambdav;
02530 double fluxv, flux_err, qual_v;
02531
02532 lambdav = cpl_vector_get( extract_lambda_vect, i);
02533 fluxv = cpl_vector_get( std_extract_vect, i);
02534 flux_err = cpl_vector_get( std_err_extract_vect, i);
02535 qual_v = cpl_vector_get( std_qual_extract_vect, i);
02536 fprintf(std_spectrum_file, "%d %f %f %f %f\n", i, lambdav, fluxv,
02537 flux_err, qual_v);
02538 }
02539 fclose( std_spectrum_file);
02540 }
02541 #endif
02542
02543 check( xsh_vector_divide_poly( opt_extract_vect,
02544 opt_extract_par->oversample,
02545 order_list->list[iorder].blazepoly, blaze_shift, instrument));
02546
02547 check( xsh_vector_divide_poly( opt_err_extract_vect,
02548 opt_extract_par->oversample,
02549 order_list->list[iorder].blazepoly, blaze_shift, instrument));
02550
02551
02552 check( xsh_interpolate_spectrum( biny, abs_order, opt_extract_par->oversample,
02553 opt_extract_par->lambda_step,
02554 extract_lambda_vect,
02555 std_extract_vect, std_err_extract_vect, std_qual_extract_vect,
02556 opt_extract_vect, opt_err_extract_vect, opt_qual_extract_vect,
02557 &result_lambda_vect,
02558 &std_flux_vect, &std_err_vect, &std_qual_vect,
02559 &opt_flux_vect, &opt_err_vect, &opt_qual_vect));
02560
02561 result_lambda_vect_size = cpl_vector_get_size( result_lambda_vect);
02562
02563
02564 if ( xsh_debug_level_get() >= XSH_DEBUG_LEVEL_MEDIUM ) {
02565 char std_spectrum_name[256];
02566 FILE* std_spectrum_file = NULL;
02567
02568 sprintf( std_spectrum_name, "spectrum_std_%d.dat", abs_order);
02569 std_spectrum_file = fopen( std_spectrum_name, "w");
02570
02571 fprintf(std_spectrum_file, "#lambda flux err qual sn %d\n",
02572 opt_extract_par->oversample);
02573 for(i=0; i< result_lambda_vect_size; i++){
02574 double lambdav, flux, err;
02575 double qual;
02576
02577 lambdav = cpl_vector_get( result_lambda_vect, i);
02578 flux = cpl_vector_get( std_flux_vect, i);
02579 err = cpl_vector_get( std_err_vect, i);
02580 qual = cpl_vector_get( std_qual_vect, i);
02581 fprintf(std_spectrum_file, "%f %f %f %f %f\n", lambdav, flux / conad,
02582 err / conad, qual, flux/err);
02583 }
02584 fclose( std_spectrum_file);
02585
02586 sprintf( std_spectrum_name, "spectrum_opt_%d.dat", abs_order);
02587 std_spectrum_file = fopen( std_spectrum_name, "w");
02588
02589 fprintf(std_spectrum_file, "#lambda flux err qual sn %d\n",
02590 opt_extract_par->oversample);
02591 for(i=0; i< result_lambda_vect_size; i++){
02592 double lambdav, flux, err;
02593 double qual;
02594
02595 lambdav = cpl_vector_get( result_lambda_vect, i);
02596 flux = cpl_vector_get( opt_flux_vect, i);
02597 err = cpl_vector_get( opt_err_vect, i);
02598 qual = cpl_vector_get( opt_qual_vect, i);
02599
02600 fprintf(std_spectrum_file, "%f %f %f %f %f\n", lambdav, flux / conad,
02601 err / conad, qual, flux/err);
02602 }
02603 fclose( std_spectrum_file);
02604 }
02605
02606 xsh_rec_list_set_data_size( orderext1d, iorder, abs_order,
02607 result_lambda_vect_size, 1);
02608 xsh_rec_list_set_data_size( orderoxt1d, iorder, abs_order,
02609 result_lambda_vect_size, 1);
02610 for(i=0; i< result_lambda_vect_size; i++){
02611 double lambdav, extflux,oxtflux;
02612 double exterr, oxterr;
02613 int extqual, oxtqual;
02614
02615 lambdav = cpl_vector_get( result_lambda_vect, i);
02616 extflux = cpl_vector_get( std_flux_vect, i);
02617 oxtflux = cpl_vector_get( opt_flux_vect, i);
02618
02619 exterr = cpl_vector_get( std_err_vect, i);
02620 oxterr = cpl_vector_get( opt_err_vect, i);
02621
02622 extqual = (int) cpl_vector_get( std_qual_vect, i);
02623 oxtqual = (int) cpl_vector_get( opt_qual_vect, i);
02624
02625 orderext1d->list[iorder].lambda[i] = lambdav;
02626 orderoxt1d->list[iorder].lambda[i] = lambdav;
02627
02628
02629 orderext1d->list[iorder].data1[i] = extflux / conad;
02630 orderoxt1d->list[iorder].data1[i] = oxtflux / conad;
02631
02632 orderext1d->list[iorder].errs1[i] = exterr / conad;
02633 orderoxt1d->list[iorder].errs1[i] = oxterr / conad;
02634
02635 orderext1d->list[iorder].qual1[i] = extqual;
02636 orderoxt1d->list[iorder].qual1[i] = oxtqual;
02637 mode = CPL_IO_EXTEND;
02638
02639 }
02640
02641 xsh_free_image( &s2Dby1D_img);
02642 xsh_free_image( &s2Dby1D_err_img);
02643 xsh_free_image( &model_img);
02644 xsh_free_image( &weight_img);
02645 xsh_unwrap_image( &extract_img);
02646 xsh_unwrap_image( &extract_errs_img);
02647 xsh_unwrap_image( &extract_qual_img);
02648 xsh_unwrap_image( &x_img);
02649 xsh_unwrap_image( &y_img);
02650 xsh_free_image( &sub_extract_img);
02651 xsh_free_image( &sub_extract_errs_img);
02652 xsh_free_image( &sub_extract_qual_img);
02653 xsh_free_image( &sub_x_img);
02654 xsh_free_image( &sub_y_img);
02655
02656 xsh_free_vector( &std_extract_vect);
02657 xsh_free_vector( &std_err_extract_vect);
02658 xsh_free_vector( &std_qual_extract_vect);
02659
02660 xsh_free_vector( &opt_extract_vect);
02661 xsh_free_vector( &opt_err_extract_vect);
02662 xsh_free_vector( &opt_qual_extract_vect);
02663
02664 xsh_unwrap_vector( &extract_lambda_vect);
02665 xsh_free_vector( &result_lambda_vect);
02666 xsh_free_vector( &std_flux_vect);
02667 xsh_free_vector( &std_err_vect);
02668 xsh_free_vector( &std_qual_vect);
02669 xsh_free_vector( &opt_flux_vect);
02670 xsh_free_vector( &opt_err_vect);
02671 xsh_free_vector( &opt_qual_vect);
02672 }
02673
02674
02675
02676 check( cpl_propertylist_append( orderext1d->header, sci_pre->data_header));
02677
02678 check( xsh_pfits_set_rectify_bin_lambda( orderext1d->header, lambda_step));
02679 check( xsh_pfits_set_rectify_bin_space( orderext1d->header, 0.0));
02680
02681 check( lambda_min = xsh_rec_list_get_lambda_min( orderext1d));
02682 check( xsh_pfits_set_rectify_lambda_min( orderext1d->header, lambda_min));
02683 check( lambda_max = xsh_rec_list_get_lambda_max( orderext1d));
02684 check( xsh_pfits_set_rectify_lambda_max( orderext1d->header, lambda_max));
02685 sprintf(tag,"%s_%s",rec_prefix,XSH_GET_TAG_FROM_ARM(XSH_ORDER_EXT1D,instrument));
02686
02687 check( xsh_pfits_set_pcatg( orderext1d->header, tag));
02688 check( xsh_pfits_set_rectify_space_min( orderext1d->header, 0.0)) ;
02689 check( xsh_pfits_set_rectify_space_max( orderext1d->header, 0.0));
02690 sprintf(filename,"%s.fits",tag);
02691 check( *orderext1d_frame = xsh_rec_list_save( orderext1d,
02692 filename, tag, CPL_TRUE));
02693
02694
02695
02696
02697 check( cpl_propertylist_append( orderoxt1d->header, sci_pre->data_header));
02698
02699 check( xsh_pfits_set_rectify_bin_lambda( orderoxt1d->header, lambda_step));
02700 check( xsh_pfits_set_rectify_bin_space( orderoxt1d->header, 0.0));
02701 check( lambda_min = xsh_rec_list_get_lambda_min( orderoxt1d));
02702 check( xsh_pfits_set_rectify_lambda_min(orderoxt1d->header, lambda_min));
02703 check( lambda_max = xsh_rec_list_get_lambda_max( orderoxt1d));
02704 check( xsh_pfits_set_rectify_lambda_max( orderoxt1d->header, lambda_max));
02705 sprintf(tag,"%s_%s",rec_prefix,XSH_GET_TAG_FROM_ARM(XSH_ORDER_OXT1D,instrument));
02706 sprintf(tag_drl,"%s_DRL",tag);
02707
02708 check( xsh_pfits_set_pcatg( orderoxt1d->header, tag_drl));
02709 check( xsh_pfits_set_rectify_space_min( orderoxt1d->header, 0.0)) ;
02710 check( xsh_pfits_set_rectify_space_max( orderoxt1d->header, 0.0));
02711 sprintf(filename,"%s.fits",tag);
02712 sprintf(filename_drl,"DRL_%s.fits",filename);
02713 xsh_add_temporary_file(filename);
02714 xsh_add_temporary_file(filename_drl);
02715
02716 check( *orderoxt1d_frame = xsh_rec_list_save( orderoxt1d, filename_drl,
02717 tag_drl,CPL_TRUE));
02718
02719 check( *res_frame_ext = xsh_rec_list_save2( orderoxt1d,filename,tag));
02720
02721
02722 XSH_ASSURE_NOT_NULL( qc_subextract_frame);
02723 XSH_ASSURE_NOT_NULL( qc_s2ddiv1d_frame);
02724 XSH_ASSURE_NOT_NULL( qc_model_frame);
02725 XSH_ASSURE_NOT_NULL( qc_weight_frame);
02726
02727 *qc_subextract_frame = cpl_frame_new();
02728
02729 sprintf(tag,"%s_OXT_SUBEXTRACT",rec_prefix);
02730 check ((cpl_frame_set_filename ( *qc_subextract_frame, qc_subextract_name),
02731 cpl_frame_set_tag( *qc_subextract_frame, tag),
02732 cpl_frame_set_type ( *qc_subextract_frame, CPL_FRAME_TYPE_IMAGE) )) ;
02733
02734 *qc_s2ddiv1d_frame = cpl_frame_new();
02735 sprintf(tag,"%s_OXT_S2DDIV1D",rec_prefix);
02736 check ((cpl_frame_set_filename ( *qc_s2ddiv1d_frame, qc_s2ddiv1d_name),
02737 cpl_frame_set_tag( *qc_s2ddiv1d_frame, tag),
02738 cpl_frame_set_type ( *qc_s2ddiv1d_frame, CPL_FRAME_TYPE_IMAGE) )) ;
02739
02740 *qc_model_frame = cpl_frame_new();
02741 sprintf(tag,"%s_OXT_MODEL",rec_prefix);
02742 check ((cpl_frame_set_filename ( *qc_model_frame, qc_model_name),
02743 cpl_frame_set_tag( *qc_model_frame, tag),
02744 cpl_frame_set_type ( *qc_model_frame, CPL_FRAME_TYPE_IMAGE) ));
02745
02746 *qc_weight_frame = cpl_frame_new();
02747 sprintf(tag,"%s_OXT_WEIGHT",rec_prefix);
02748 check ((cpl_frame_set_filename ( *qc_weight_frame, qc_weight_name),
02749 cpl_frame_set_tag( *qc_weight_frame, tag),
02750 cpl_frame_set_type ( *qc_weight_frame, CPL_FRAME_TYPE_IMAGE) ));
02751
02752 check( xsh_add_temporary_file( qc_subextract_name));
02753 check( xsh_add_temporary_file( qc_s2ddiv1d_name));
02754 check( xsh_add_temporary_file( qc_model_name));
02755 check( xsh_add_temporary_file( qc_weight_name));
02756
02757 cleanup:
02758 if ( cpl_error_get_code() != CPL_ERROR_NONE) {
02759 xsh_unwrap_image( &sub_extract_img);
02760 xsh_unwrap_image( &sub_extract_errs_img);
02761 xsh_unwrap_image( &sub_extract_qual_img);
02762
02763 xsh_free_vector( &std_extract_vect);
02764 xsh_free_vector( &std_err_extract_vect);
02765 xsh_free_vector( &std_qual_extract_vect);
02766
02767 xsh_free_vector( &opt_extract_vect);
02768 xsh_free_vector( &opt_err_extract_vect);
02769 xsh_free_vector( &opt_qual_extract_vect);
02770
02771 xsh_unwrap_vector( &extract_lambda_vect);
02772
02773 xsh_free_vector( &result_lambda_vect);
02774 xsh_free_vector( &std_flux_vect);
02775 xsh_free_vector( &std_err_vect);
02776 xsh_free_vector( &std_qual_vect);
02777 xsh_free_vector( &opt_flux_vect);
02778 xsh_free_vector( &opt_err_vect);
02779 xsh_free_vector( &opt_qual_vect);
02780 }
02781 xsh_pre_free( &sci_pre);
02782 xsh_order_list_free( &order_list);
02783 xsh_spectralformat_list_free( &spectralformat_list);
02784 xsh_wavesol_free( &wavesol);
02785 XSH_FREE( lambdas);
02786 XSH_FREE( pos_cen_x);
02787 XSH_FREE( pos_cen_y);
02788 XSH_FREE( extract_data);
02789 XSH_FREE( extract_x_data);
02790 XSH_FREE( extract_y_data);
02791 XSH_FREE( extract_errs);
02792 XSH_FREE( extract_qual);
02793 xsh_free_image( &blaze_img);
02794 xsh_rec_list_free( &orderext1d);
02795 xsh_rec_list_free( &orderoxt1d);
02796 xsh_free_propertylist( &qc_header);
02797 return;
02798 }
02799
02800
02801
02802