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