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