00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #ifdef HAVE_CONFIG_H
00027 #include <config.h>
00028 #endif
00029
00030
00037
00040
00041
00042
00043
00044 #include <math.h>
00045 #include <xsh_drl.h>
00046
00047 #include <xsh_baryvel.h>
00048 #include <xsh_badpixelmap.h>
00049 #include <xsh_data_rec.h>
00050 #include <xsh_data_localization.h>
00051 #include <xsh_data_pre.h>
00052 #include <xsh_data_order.h>
00053 #include <xsh_data_wavesol.h>
00054 #include <xsh_data_spectralformat.h>
00055 #include <xsh_dfs.h>
00056 #include <xsh_pfits.h>
00057 #include <xsh_error.h>
00058 #include <xsh_msg.h>
00059 #include <xsh_fit.h>
00060 #include <xsh_model_io.h>
00061 #include <xsh_model_kernel.h>
00062 #include <xsh_ifu_defs.h>
00063 #include <xsh_data_dispersol.h>
00064 #include <cpl.h>
00065 #include <xsh_utils.h>
00066 #include <xsh_utils_ifu.h>
00067 #include <xsh_utils_table.h>
00068 #include <xsh_rectify.h>
00069 #include <xsh_model_utils.h>
00070
00071
00072
00073
00074
00075
00076
00077 #define SLIT_FRAC_STEP 50
00078
00079 static cpl_vector * rec_profile = NULL ;
00080 static cpl_vector * err_profile = NULL ;
00081
00082 static double compute_shift_with_localization( cpl_frame *loc_frame,
00083 cpl_frame *loc0_frame);
00084
00085 static double compute_shift_with_kw( cpl_propertylist * header,
00086 xsh_rectify_param *rectify_par,
00087 double **ref_ra, double **ref_dec, int flag);
00088
00089 static cpl_frame*
00090 xsh_shift( cpl_frame *rec_frame, xsh_instrument *instrument,
00091 const char *fname, double slit_shift, cpl_frame** res_frame_ext);
00092
00093 static void xsh_frame_set_shiftifu_ref( cpl_propertylist *header,
00094 cpl_frame *shift_frame);
00095
00096
00097
00098
00099
00100 static void get_errors( float *err, int ix, int iy, int nx, int ny,
00101 double radius, float *err_val, int* qual_val)
00102 {
00103 int xmin, xmax, ymin, ymax ;
00104 int i, j, x, y, n ;
00105 double sum2 ;
00106 double *err_vect = NULL;
00107 int err_vect_size = 0;
00108
00109 XSH_ASSURE_NOT_NULL( err);
00110 XSH_ASSURE_NOT_NULL( err_val);
00111 XSH_ASSURE_NOT_NULL( qual_val);
00112
00113
00114
00115
00116
00117 n = 0 ;
00118 sum2 = 0 ;
00119 xmin = ix - radius ;
00120 xmax = ix + radius ;
00121 ymin = iy - radius ;
00122 ymax = iy + radius ;
00123
00124
00125 err_vect_size = (2*radius+1)*(2*radius+1);
00126
00127 XSH_CALLOC( err_vect, double, err_vect_size);
00128
00129 for( y = ymin ; y < ymax ; y++ ){
00130 for( x = xmin ; x < xmax ; x++ ) {
00131 if ( y < 0 || y >= ny || x < 0 || x >= nx ) continue;
00132 err_vect[n] = err[x+y*nx];
00133 n++;
00134 }
00135 }
00136
00137 for(i=0; i< n; i++){
00138 for(j=i; j< n; j++){
00139 if (i==j){
00140 sum2 += err_vect[i]*err_vect[j];
00141 }
00142 else{
00143 sum2 += 2*err_vect[i]*err_vect[j];
00144 }
00145 }
00146 }
00147
00148 if (n != 0){
00149 *err_val = sqrt( sum2 )/(double)n ;
00150 }
00151 else{
00152 *err_val = 1;
00153 *qual_val = QFLAG_OUTSIDE_DATA_RANGE;
00154 }
00155 cleanup:
00156 XSH_FREE( err_vect);
00157 }
00158
00159
00160
00161 static void
00162 xsh_rec_list_rectify( xsh_rec_list *rec_list, int iorder, int irec,
00163 xsh_pre *sci_pre,cpl_image* var_img,
00164 double fx, double fy, double radius, cpl_vector *profile,
00165 cpl_vector* profile2,
00166 double mult, int *good, int *bad, int *out)
00167 {
00168 int nx,ny;
00169 float *rec_flux = NULL;
00170 float *rec_errs = NULL;
00171 int *rec_qual = NULL;
00172 int ix,iy;
00173
00174 XSH_ASSURE_NOT_NULL( rec_list);
00175 XSH_ASSURE_NOT_NULL( sci_pre);
00176 XSH_ASSURE_NOT_NULL( good);
00177 XSH_ASSURE_NOT_NULL( bad);
00178 XSH_ASSURE_NOT_NULL( out);
00179 XSH_ASSURE_NOT_NULL( profile);
00180
00181 check( rec_flux = xsh_rec_list_get_data1( rec_list, iorder));
00182 check( rec_errs = xsh_rec_list_get_errs1( rec_list, iorder));
00183 check( rec_qual = xsh_rec_list_get_qual1( rec_list, iorder));
00184
00185 check( nx = xsh_pre_get_nx( sci_pre));
00186 check( ny = xsh_pre_get_ny( sci_pre));
00187
00188 ix = floor( fx - 0.5) ;
00189 iy = floor( fy - 0.5) ;
00190
00191 if ( ix <0 || ix >= nx || iy < 0 || iy >= ny ) {
00192 xsh_msg_dbg_high( "OUT_OF_RANGE at %d,%d", ix, iy);
00193 *out +=1;
00194 rec_flux[irec] = 0;
00195 rec_errs[irec] = 1;
00196 rec_qual[irec] = QFLAG_OUTSIDE_DATA_RANGE;
00197 }
00198 else{
00199 cpl_image *sci_img = NULL;
00200 cpl_image *err_img = NULL;
00201
00202 float *data = NULL;
00203 float *errs = NULL;
00204 int *qual = NULL;
00205 int imqual_val;
00206 float imflux_val, imerr_val;
00207
00208 check( sci_img = xsh_pre_get_data( sci_pre));
00209 check( err_img = xsh_pre_get_errs( sci_pre));
00210 check( data = cpl_image_get_data_float( sci_img));
00211 check( qual = cpl_image_get_data_int( xsh_pre_get_qual( sci_pre)));
00212 check( errs = cpl_image_get_data_float( err_img));
00213
00214 int ix_iynx=ix+iy*nx;
00215 imqual_val = qual[ix_iynx];
00216 imflux_val = data[ix_iynx];
00217 imerr_val = errs[ix_iynx];
00218
00219 if ( radius == 0.){
00220 rec_flux[irec] = imflux_val;
00221 rec_errs[irec] = imerr_val;
00222 rec_qual[irec] = imqual_val;
00223 }
00224 else {
00225 double flux_val=0, confidence=0 ;
00226 float err_val=0;
00227 int qual_val = QFLAG_GOOD_PIXEL;
00228 check( flux_val = cpl_image_get_interpolated( sci_img, fx, fy,
00229 profile, radius, profile, radius, &confidence));
00230
00231
00232
00233
00234
00235 check( err_val = cpl_image_get_interpolated( var_img, fx, fy,
00236 profile2, radius, profile2, radius, &confidence));
00237
00238 if (confidence >= 0.5){
00239 rec_flux[irec] = flux_val;
00240 rec_errs[irec] = sqrt(err_val);
00241 rec_qual[irec] = qual_val;
00242 *good +=1;
00243 }
00244 else if (confidence > 0){
00245 if (flux_val == 0){
00246 XSH_REGDEBUG("fx %f fy %f confidence %f", fx,fy,confidence);
00247 }
00248 rec_flux[irec] = flux_val;
00249 rec_errs[irec] = sqrt(err_val);
00250 rec_qual[irec] = QFLAG_MISSING_DATA;
00251 *bad +=1;
00252 }
00253 else{
00254 rec_flux[irec] = 0;
00255 rec_errs[irec] = 1;
00256 rec_qual[irec] = QFLAG_MISSING_DATA;
00257 *bad += 1;
00258 }
00259
00260 }
00261 }
00262 rec_flux[irec] *= mult;
00263 rec_errs[irec] *= mult;
00264
00265 cleanup:
00266 return;
00267 }
00268
00269
00270
00291
00292
00293 static inline void
00294 fill_rectified( xsh_pre *pre_sci,
00295 xsh_rec_list *rec_list,
00296 int idx,
00297 xsh_wavesol *wavesol,
00298 xsh_xs_3 *model_config,
00299 xsh_instrument *instrument,
00300 xsh_dispersol_list *disp_list,
00301 float slit_min,
00302 float slit_max,
00303 double lambda_min,
00304 int skip_low,
00305 int skip_up,
00306 xsh_rectify_param *rectify_par,
00307 double slit_shift,
00308 cpl_frame *slit_shiftab_frame)
00309 {
00310
00311 int ns, nl;
00312 int nslit, nlambda ;
00313 float * slit = NULL ;
00314 double * lambda = NULL ;
00315 float * flux1 = NULL ;
00316 float * errs1 = NULL ;
00317 int * qual1 = NULL ;
00318 int order, last_slit, last_lambda;
00319 double *x_data = NULL;
00320 double *y_data = NULL;
00321
00322 int good_pixel = 0, out_of_bound = 0, bad_pixel = 0 ;
00323 double *slit_shift_data = NULL;
00324 const char *slit_shifttab_name = NULL;
00325 cpl_table *shift_tab = NULL;
00326 int shift_size;
00327 double *slit_shifttab_data = NULL;
00328 double *wave_shifttab_data = NULL;
00329 cpl_image* var_img=NULL;
00330
00331
00332 XSH_ASSURE_NOT_NULL( pre_sci);
00333 XSH_ASSURE_NOT_NULL( rec_list);
00334 check( slit = xsh_rec_list_get_slit( rec_list, idx));
00335 XSH_ASSURE_NOT_NULL( slit);
00336 check( lambda = xsh_rec_list_get_lambda( rec_list, idx));
00337 XSH_ASSURE_NOT_NULL( lambda);
00338 check( flux1 = xsh_rec_list_get_data1( rec_list, idx));
00339 check( errs1 = xsh_rec_list_get_errs1( rec_list, idx));
00340 check( qual1 = xsh_rec_list_get_qual1( rec_list, idx));
00341
00342 var_img=cpl_image_duplicate(xsh_pre_get_errs( pre_sci));
00343 cpl_image_multiply(var_img,var_img);
00344
00345
00346 if ( rec_profile == NULL && rectify_par->rectif_radius != 0 ) {
00347 check( rec_profile = cpl_vector_new( CPL_KERNEL_DEF_SAMPLES));
00348
00349 assure(rectify_par->rectif_radius > 0,CPL_ERROR_ILLEGAL_INPUT,
00350 "rectify-radius must be positive");
00351 check( cpl_vector_fill_kernel_profile( rec_profile,
00352 rectify_par->kernel_type, rectify_par->rectif_radius));
00353 err_profile=cpl_vector_duplicate(rec_profile);
00354 cpl_vector_multiply(err_profile,err_profile);
00355 }
00356
00357
00358 check( order = xsh_rec_list_get_order( rec_list, idx));
00359 check( nslit = xsh_rec_list_get_nslit( rec_list, idx));
00360 check( nlambda = xsh_rec_list_get_nlambda( rec_list, idx));
00361
00362 last_slit = nslit-1;
00363 last_lambda = nlambda-1;
00364
00365 XSH_CALLOC( x_data, double, nslit*nlambda);
00366 XSH_CALLOC( y_data, double, nslit*nlambda);
00367
00368
00369 for( ns = 0 ; ns < nslit ; ns++ ) {
00370
00371
00372
00373
00374 rec_list->list[idx].slit[ns] = slit_min +
00375 (double)ns*rectify_par->rectif_bin_space;
00376 }
00377
00378 for( nl = 0 ; nl < nlambda ; nl++ ) {
00379
00380 rec_list->list[idx].lambda[nl] = lambda_min +
00381 (double)nl*rectify_par->rectif_bin_lambda;
00382 }
00383
00384 xsh_msg_dbg_low( "Order %d slit (%f,%f):%d lambda (%f,%f):%d",
00385 order, rec_list->list[idx].slit[0],
00386 rec_list->list[idx].slit[last_slit], nslit,
00387 rec_list->list[idx].lambda[0],
00388 rec_list->list[idx].lambda[last_lambda], nlambda);
00389
00390
00391 XSH_CALLOC( slit_shift_data, double, nlambda);
00392 if( slit_shiftab_frame != NULL){
00393 check( slit_shifttab_name = cpl_frame_get_filename( slit_shiftab_frame));
00394 xsh_msg_dbg_low("Load slit shift tab %s", slit_shifttab_name);
00395 XSH_TABLE_LOAD( shift_tab, slit_shifttab_name);
00396 check( wave_shifttab_data = cpl_table_get_data_double( shift_tab,
00397 XSH_SHIFTIFU_COLNAME_WAVELENGTH));
00398 check( slit_shifttab_data = cpl_table_get_data_double( shift_tab,
00399 XSH_SHIFTIFU_COLNAME_SHIFTSLIT));
00400 shift_size = cpl_table_get_nrow( shift_tab);
00401
00402 for( nl = 0 ; nl < nlambda ; nl++ ) {
00403 double wave = lambda[nl];
00404 double cshift;
00405 cshift = xsh_data_interpolate( wave, shift_size, wave_shifttab_data,
00406 slit_shifttab_data);
00407 slit_shift_data[nl] = cshift;
00408 }
00409
00410 if ( xsh_debug_level_get() >= XSH_DEBUG_LEVEL_MEDIUM) {
00411 FILE *shift_file = NULL;
00412 char shift_name[256];
00413
00414 sprintf( shift_name, "shift_%d_%.1f_%.1f.dat", order, slit_min, slit_max);
00415 shift_file = fopen( shift_name, "w+");
00416
00417 fprintf( shift_file, "# wave delta_s\n");
00418
00419 for( nl = 0 ; nl < nlambda ; nl++ ) {
00420 fprintf ( shift_file, "%f %f\n", lambda[nl], slit_shift_data[nl]);
00421 }
00422
00423 fclose( shift_file);
00424 }
00425 }
00426
00427
00428 if ( wavesol != NULL){
00429 for( ns = skip_low ; ns < (nslit - skip_up) ; ns++ ) {
00430 double xslit = slit[ns]+slit_shift;
00431
00432 for( nl = 0 ; nl < nlambda ; nl++ ) {
00433 int irec;
00434 double xlambda = lambda[nl];
00435 double xslitrec;
00436
00437 irec= nl+ns*nlambda;
00438 xslitrec = xslit+slit_shift_data[nl];
00439
00440 check( x_data[irec] = xsh_wavesol_eval_polx( wavesol, xlambda, order,
00441 xslitrec));
00442 check( y_data[irec] = xsh_wavesol_eval_poly( wavesol, xlambda, order,
00443 xslitrec));
00444 }
00445 }
00446 }
00447 else{
00448 XSH_ASSURE_NOT_NULL( model_config);
00449 for( ns = skip_low ; ns < (nslit - skip_up) ; ns++ ) {
00450 double xslit = slit[ns]+slit_shift;
00451
00452 for( nl = 0 ; nl < nlambda ; nl++ ) {
00453 int irec;
00454 double pm_x=0, pm_y=0;
00455 double xlambda = lambda[nl];
00456 double xslitrec;
00457
00458 irec= nl+ns*nlambda;
00459 xslitrec = xslit+slit_shift_data[nl];
00460 check( xsh_model_get_xy( model_config, instrument, xlambda, order,
00461 xslitrec,
00462 &pm_x, &pm_y));
00463 x_data[irec] = convert_data_to_bin( pm_x, instrument->binx );
00464 y_data[irec] = convert_data_to_bin( pm_y, instrument->biny );
00465 }
00466 }
00467 }
00468
00469 if ( xsh_debug_level_get() > XSH_DEBUG_LEVEL_MEDIUM) {
00470 FILE *rec_file = NULL;
00471 char recfile_name[256];
00472
00473 sprintf( recfile_name, "rec_%d_%.1f_%.1f.dat", order, slit_min ,slit_max);
00474 rec_file = fopen( recfile_name, "w+");
00475
00476 fprintf( rec_file, "#irec wave slit x y\n");
00477
00478 for( ns = skip_low ; ns < (nslit - skip_up) ; ns++ ) {
00479 double xslit = slit[ns]+slit_shift;
00480
00481 for( nl = 0 ; nl < nlambda ; nl++ ) {
00482 int irec;
00483 double xlambda = lambda[nl];
00484 double xslitrec;
00485
00486 irec= nl+ns*nlambda;
00487 xslitrec = xslit+slit_shift_data[nl];
00488 fprintf ( rec_file, "%d %f %f %f %f\n", irec, xlambda, xslitrec, x_data[irec],
00489 y_data[irec]);
00490 }
00491 }
00492 fclose( rec_file);
00493 }
00494
00495 if (disp_list == NULL) {
00496
00497 for( ns = skip_low ; ns < (nslit - skip_up) ; ns++ ) {
00498 for( nl = 0 ; nl < nlambda ; nl++ ) {
00499 double fx=0, fy=0;
00500 int irec=0;
00501
00502 irec= nl+ns*nlambda;
00503 fx = x_data[irec];
00504 fy = y_data[irec];
00505
00506 check( xsh_rec_list_rectify( rec_list, idx, irec, pre_sci,var_img,
00507 fx, fy,
00508 rectify_par->rectif_radius,
00509 rec_profile,err_profile,
00510 1, &good_pixel, &bad_pixel,
00511 &out_of_bound));
00512 }
00513
00514 }
00515 }
00516 else{
00517 cpl_polynomial *plambdadx = NULL, *plambdady = NULL;
00518 cpl_polynomial *pslitdx = NULL, *pslitdy = NULL;
00519 cpl_vector *val = NULL;
00520 int abs_order=0;
00521
00522 check( plambdadx = cpl_polynomial_duplicate(
00523 disp_list->list[idx].lambda_poly));
00524 check( plambdady = cpl_polynomial_duplicate(
00525 disp_list->list[idx].lambda_poly));
00526 check( pslitdx = cpl_polynomial_duplicate(
00527 disp_list->list[idx].slit_poly));
00528 check( pslitdy = cpl_polynomial_duplicate(
00529 disp_list->list[idx].slit_poly));
00530
00531 abs_order = disp_list->list[idx].absorder;
00532
00533 check( cpl_polynomial_derivative( plambdadx, 0));
00534 check( cpl_polynomial_derivative( plambdady, 1));
00535 check( cpl_polynomial_derivative( pslitdx, 0));
00536 check( cpl_polynomial_derivative( pslitdy, 1));
00537
00538 check( val = cpl_vector_new(2));
00539 double fx, fy;
00540 double ldx, ldy, sdx, sdy, abs_determinant;
00541 int irec;
00542 int ns_nlambda=0;
00543 double bin_size=rectify_par->rectif_bin_space*rectify_par->rectif_bin_lambda;
00544 for( ns = skip_low ; ns < (nslit - skip_up) ; ns++ ) {
00545 ns_nlambda=ns*nlambda;
00546 for( nl = 0 ; nl < nlambda ; nl++ ) {
00547
00548 irec= nl+ns_nlambda;
00549 fx = x_data[irec];
00550 fy = y_data[irec];
00551
00552 check( cpl_vector_set( val, 0, fx));
00553 check( cpl_vector_set( val, 1, fy));
00554
00555 check( ldx = cpl_polynomial_eval( plambdadx, val));
00556 check( ldy = cpl_polynomial_eval( plambdady, val));
00557 check( sdx = cpl_polynomial_eval( pslitdx, val));
00558 check( sdy = cpl_polynomial_eval( pslitdy, val));
00559 abs_determinant = bin_size/(fabs(ldx*sdy-sdx*ldy));
00560
00561 check( xsh_rec_list_rectify( rec_list, idx, irec, pre_sci,var_img,
00562 fx, fy,
00563 rectify_par->rectif_radius,
00564 rec_profile, err_profile,abs_determinant,
00565 &good_pixel, &bad_pixel, &out_of_bound));
00566 }
00567 }
00568 xsh_free_polynomial( &plambdadx);
00569 xsh_free_polynomial( &plambdady);
00570 xsh_free_polynomial( &pslitdx);
00571 xsh_free_polynomial( &pslitdy);
00572 xsh_free_vector( &val);
00573 }
00574
00575 xsh_msg_dbg_low( " Good pixels: %d, Out of bound: %d, Bad pixels: %d",
00576 good_pixel, out_of_bound, bad_pixel ) ;
00577 cleanup:
00578 xsh_free_image(&var_img);
00579 xsh_free_vector( &rec_profile);
00580 xsh_free_vector( &err_profile);
00581 XSH_FREE( x_data);
00582 XSH_FREE( y_data);
00583 XSH_FREE( slit_shift_data);
00584 XSH_TABLE_FREE( shift_tab);
00585 return ;
00586 }
00587
00588
00589
00590
00600 void
00601 xsh_get_slit_edges( cpl_frame *slitmap_frame, double *sdown, double *sup,
00602 double *sldown, double *slup, xsh_instrument *instrument)
00603 {
00604 const char *slitmap_name = NULL;
00605 cpl_propertylist *slitmap_header = NULL;
00606
00607 if ( slitmap_frame != NULL) {
00608 XSH_ASSURE_NOT_NULL( sdown);
00609 XSH_ASSURE_NOT_NULL( sup);
00610
00611 check( slitmap_name = cpl_frame_get_filename( slitmap_frame));
00612 check( slitmap_header = cpl_propertylist_load( slitmap_name, 0));
00613
00614 *sdown = xsh_pfits_get_slitmap_median_edglo( slitmap_header);
00615 if ( cpl_error_get_code() != CPL_ERROR_NONE){
00616 xsh_msg_warning(
00617 "Keyword 'MEDIAN EDGLO SLIT' not found in SLIT_MAP %s. Using default value %f",
00618 slitmap_name, MIN_SLIT);
00619 xsh_error_reset();
00620 }
00621
00622 *sup = xsh_pfits_get_slitmap_median_edgup( slitmap_header);
00623 if ( cpl_error_get_code() != CPL_ERROR_NONE){
00624 xsh_msg_warning(
00625 "Keyword 'MEDIAN EDGUP SLIT' not found in SLIT_MAP %s. Using default value %f",
00626 slitmap_name, MAX_SLIT);
00627 xsh_error_reset();
00628 }
00629 if ( xsh_instrument_get_mode(instrument) == XSH_MODE_IFU){
00630 XSH_ASSURE_NOT_NULL( sldown);
00631 XSH_ASSURE_NOT_NULL( slup);
00632 *sldown = xsh_pfits_get_slitmap_median_sliclo( slitmap_header);
00633 if ( cpl_error_get_code() != CPL_ERROR_NONE){
00634 xsh_msg_warning(
00635 "Keyword 'MEDIAN SLICLO SLIT' not found in SLIT_MAP %s. Using default value %f",
00636 slitmap_name, SLITLET_CEN_CENTER-2);
00637 xsh_error_reset();
00638 }
00639 *slup = xsh_pfits_get_slitmap_median_slicup( slitmap_header);
00640 if ( cpl_error_get_code() != CPL_ERROR_NONE){
00641 xsh_msg_warning(
00642 "Keyword 'MEDIAN SLICUP SLIT' not found in SLIT_MAP %s. Using default value %f",
00643 slitmap_name, SLITLET_CEN_CENTER+2);
00644 xsh_error_reset();
00645 }
00646 }
00647 }
00648 else {
00649 xsh_msg_warning( "No provided SLIT_MAP. Using default values: 'MEDIAN EDGLO SLIT' %f - 'MEDIAN EDGUP SLIT' %f",
00650 MIN_SLIT, MAX_SLIT);
00651 *sdown = MIN_SLIT;
00652 *sup = MAX_SLIT;
00653 if ( xsh_instrument_get_mode(instrument) == XSH_MODE_IFU){
00654 xsh_msg_warning( "Using default values: 'MEDIAN SLICLO SLIT' %f - 'MEDIAN SLICUP SLIT' %f",
00655 SLITLET_CEN_CENTER-2, SLITLET_CEN_CENTER+2);
00656 xsh_error_reset();
00657 XSH_ASSURE_NOT_NULL( sldown);
00658 XSH_ASSURE_NOT_NULL( slup);
00659 *sldown = SLITLET_CEN_CENTER-2;
00660 *slup = SLITLET_CEN_CENTER+2;
00661 }
00662 }
00663
00664 if ( xsh_instrument_get_mode(instrument) == XSH_MODE_IFU){
00665 xsh_msg( "IFU limits: slitlet DOWN [%f %f], size: %f arcsec",
00666 *sdown, *sldown, *sldown-*sdown);
00667 xsh_msg( "IFU limits: slitlet CEN [%f %f], size: %f arcsec",
00668 *sldown, *slup, *slup-*sldown);
00669 xsh_msg( "IFU limits: slitlet UP [%f %f], size: %f arcsec",
00670 *slup, *sup, *sup-*slup);
00671 }
00672 else{
00673 xsh_msg( "SLIT limits: [%f %f], total size: %f arcsec", *sdown, *sup, *sup-*sdown);
00674 }
00675 cleanup:
00676 xsh_free_propertylist( &slitmap_header);
00677 return ;
00678 }
00679
00680
00688 void xsh_rec_slit_size( xsh_rectify_param *rectify_par,
00689 double *slit_min, int *nslit, XSH_MODE mode)
00690 {
00691 double smin= MIN_SLIT, smax=MAX_SLIT;
00692 double slit_step;
00693
00694 XSH_ASSURE_NOT_NULL( rectify_par);
00695 XSH_ASSURE_NOT_NULL( slit_min);
00696 XSH_ASSURE_NOT_NULL( nslit);
00697
00698 slit_step = rectify_par->rectif_bin_space;
00699
00700 if ( mode == XSH_MODE_SLIT){
00701 if ( rectify_par->rectify_full_slit != 1 ) {
00702 xsh_msg_warning(" Option not READY go to FULL_SLIT");
00703 }
00704 *nslit = (smax - smin)/slit_step;
00705 *slit_min = smin;
00706 smax = smin+(*nslit-1)*slit_step;
00707
00708 xsh_msg( "SLIT : (%.3f,%.3f) used only (%.3f,%.3f) in %d elts",
00709 MIN_SLIT, MAX_SLIT, smin, smax, *nslit);
00710 }
00711
00712 cleanup:
00713 return;
00714 }
00715
00716
00717 static void
00718 adjust_lambdas( xsh_spectralformat_list * spec_list,
00719 xsh_rectify_param * rectify_par )
00720 {
00721 int i ;
00722 double step = rectify_par->rectif_bin_lambda;
00723
00724 xsh_msg_dbg_high( "Adjust Lambdas: %d orders", spec_list->size ) ;
00725
00726
00727 for ( i = 0 ; i<spec_list->size ; i++ ) {
00728 double min, max ;
00729 int imin, imax ;
00730
00731 min = spec_list->list[i].lambda_min ;
00732 max = spec_list->list[i].lambda_max ;
00733 imin = ceil( min/step) ;
00734 min = (double)imin * step;
00735 imax = floor( max / step);
00736 max = (double)imax*step ;
00737 xsh_msg_dbg_high( "--- test lambdamin - %d : %.4lf -> %.4lf, %lf -> %lf", i,
00738 spec_list->list[i].lambda_min, min,
00739 spec_list->list[i].lambda_max, max ) ;
00740 spec_list->list[i].lambda_min = min;
00741 spec_list->list[i].lambda_max = max;
00742 }
00743 }
00744
00745 static float
00746 get_step_slit( float * slit, int nslit )
00747 {
00748 float result = 0. ;
00749 float min, max ;
00750
00751 max = *(slit+nslit-1) ;
00752 min = *slit ;
00753 result = (max-min)/(float)nslit ;
00754
00755 xsh_msg( " Step Slit = %f", result ) ;
00756 return result ;
00757
00758 }
00759
00760
00761
00794
00795 cpl_frame *
00796 xsh_rectify( cpl_frame *sci_frame,
00797 cpl_frame * orderlist_frame,
00798 cpl_frame *wavesol_frame,
00799 cpl_frame * model_frame,
00800 xsh_instrument *instrument,
00801 xsh_rectify_param *rectify_par,
00802 cpl_frame *spectralformat_frame,
00803 cpl_frame *disp_tab_frame,
00804 const char * res_name,
00805 cpl_frame** res_frame_ext,
00806 cpl_frame** res_frame_tab,
00807 const char* rec_prefix)
00808 {
00809 cpl_frame *result = NULL;
00810 xsh_order_list *orderlist = NULL ;
00811 int nslit;
00812 double slit_min;
00813 char tag[80];
00814
00815 XSH_ASSURE_NOT_NULL( orderlist_frame);
00816
00817 check( orderlist = xsh_order_list_load ( orderlist_frame, instrument));
00818
00819 sprintf(tag,"%s_%s",rec_prefix,
00820 XSH_GET_TAG_FROM_ARM( XSH_ORDER2D, instrument));
00821
00822 xsh_rec_slit_size( rectify_par, &slit_min, &nslit, XSH_MODE_SLIT);
00823
00824 check( result = xsh_rectify_orders( sci_frame, orderlist,
00825 wavesol_frame, model_frame, instrument, rectify_par, spectralformat_frame,
00826 disp_tab_frame, res_name, tag,
00827 res_frame_ext, res_frame_tab,
00828 0, 100, slit_min, nslit,0, NULL));
00829
00830 cleanup:
00831 xsh_order_list_free( &orderlist);
00832 return result;
00833 }
00834
00835
00857 cpl_frame * xsh_rectify_and_shift( cpl_frame *sci_frame,
00858 cpl_frame * orderlist_frame,
00859 cpl_frame *wavesol_frame,
00860 cpl_frame * model_frame,
00861 xsh_instrument *instrument,
00862 xsh_rectify_param *rectify_par,
00863 cpl_frame *spectralformat_frame,
00864 cpl_frame * loc_frame,
00865 cpl_frame * loc0_frame,
00866 double *throw_shift,
00867 cpl_frame *disp_tab_frame,
00868 const char * res_name,
00869 cpl_frame** res_frame_ext,
00870 cpl_frame** res_frame_tab)
00871 {
00872 cpl_frame *result = NULL, *shift_rec = NULL;
00873 cpl_frame * shift_rec_ext = NULL;
00874 xsh_order_list *orderlist = NULL ;
00875 int nslit;
00876 double slit_min;
00877
00878 const char *tag = NULL;
00879 cpl_propertylist *header = NULL;
00880 double bin_space;
00881 double shift_arcsec=0, shift_pix=0, shift_bin_arcsec=0;
00882 double rectify_shift;
00883 #if 0
00884
00885 cpl_frame * result0 = NULL;
00886 char * res0_name = NULL ;
00887 #endif
00888 xsh_msg( "===> xsh_rectify_and_shift" ) ;
00889
00890 XSH_ASSURE_NOT_NULL( orderlist_frame);
00891 XSH_ASSURE_NOT_NULL( rectify_par);
00892
00893 check( orderlist = xsh_order_list_load ( orderlist_frame, instrument));
00894
00895 tag = XSH_GET_TAG_FROM_ARM( XSH_ORDER2D, instrument);
00896 xsh_rec_slit_size( rectify_par, &slit_min, &nslit, XSH_MODE_SLIT);
00897 bin_space = rectify_par->rectif_bin_space;
00898
00899 if ( throw_shift == NULL) {
00900 if ( loc0_frame != NULL){
00901 XSH_ASSURE_NOT_NULL( loc_frame ) ;
00902 xsh_msg("Compute shift with localization");
00903 shift_arcsec = compute_shift_with_localization( loc0_frame, loc_frame);
00904 }
00905 else{
00906 xsh_msg("Fixed shift with localization to 0");
00907 shift_arcsec = 0;
00908 }
00909 }
00910 else {
00911 shift_arcsec = *throw_shift;
00912 xsh_msg("Use throw shift %f", shift_arcsec);
00913 }
00914
00915 shift_pix = xsh_round_double(shift_arcsec/bin_space);
00916 shift_bin_arcsec = shift_pix*bin_space;
00917 rectify_shift = shift_bin_arcsec-shift_arcsec;
00918
00919 xsh_msg(" Mesured Shift for rectify : %f", rectify_shift);
00920
00921 check( shift_rec = xsh_rectify_orders( sci_frame, orderlist,
00922 wavesol_frame, model_frame, instrument,
00923 rectify_par, spectralformat_frame,
00924 disp_tab_frame, res_name, tag,
00925 &shift_rec_ext, res_frame_tab,
00926 0, 100, slit_min, nslit, rectify_shift,
00927 NULL));
00928
00929 check( result = xsh_shift( shift_rec, instrument, res_name, shift_bin_arcsec,
00930 res_frame_ext));
00931
00932 cleanup:
00933 xsh_free_frame( &shift_rec);
00934 xsh_free_frame( &shift_rec_ext);
00935 xsh_order_list_free( &orderlist);
00936 xsh_free_propertylist( &header);
00937 return result;
00938 }
00939
00940
00941
00942
00943 cpl_frameset* xsh_rectify_ifu(cpl_frame *sci_frame,
00944 cpl_frame *orderlist_frame, cpl_frameset *wavesol_frameset,
00945 cpl_frameset *shiftifu_frameset,
00946 cpl_frame *model_frame, xsh_instrument *instrument,
00947 xsh_rectify_param *rectify_par, cpl_frame *spectralformat_frame,
00948 cpl_frame * slitmap_frame,
00949 cpl_frameset** rec_frameset_ext,
00950 cpl_frameset** rec_frameset_tab,
00951 const char* rec_prefix )
00952 {
00953 cpl_frameset *result = NULL;
00954 xsh_order_list *orderlist = NULL ;
00955
00956 XSH_ASSURE_NOT_NULL( orderlist_frame);
00957 check( orderlist = xsh_order_list_load ( orderlist_frame, instrument));
00958
00959
00960 XSH_REGDEBUG("TODO : ADD disp_tab frameset, res_frame_ext frameset");
00961
00962 check( result = xsh_rectify_orders_ifu( sci_frame, orderlist,
00963 wavesol_frameset, shiftifu_frameset, model_frame, instrument, rectify_par,
00964 spectralformat_frame,
00965 slitmap_frame,rec_frameset_ext,rec_frameset_tab, 0, 100, rec_prefix));
00966
00967 cleanup:
00968 xsh_order_list_free( &orderlist);
00969 return result;
00970 }
00971
00972
00973
01005
01006 cpl_frame*
01007 xsh_rectify_orders( cpl_frame *sci_frame,
01008 xsh_order_list *orderlist,
01009 cpl_frame *wavesol_frame,
01010 cpl_frame *model_frame,
01011 xsh_instrument *instrument,
01012 xsh_rectify_param *rectify_par,
01013 cpl_frame *spectralformat_frame,
01014 cpl_frame *disp_tab_frame,
01015 const char *res_name,
01016 const char* tag,
01017 cpl_frame** res_frame_ext,
01018 cpl_frame** res_frame_tab,
01019 int min_index,
01020 int max_index,
01021 double slit_min,
01022 int nslit,
01023 double slit_shift,
01024 cpl_frame *slitshift_tab)
01025 {
01026 xsh_rec_list * rec_list = NULL ;
01027 xsh_xs_3 config_model ;
01028
01029 int solution_type = 0;
01030 xsh_wavesol * wavesol = NULL ;
01031 cpl_frame * res_frame = NULL ;
01032 xsh_pre * pre_sci = NULL ;
01033 int i, nlambda, order ;
01034 float slit_max ;
01035 xsh_spectralformat_list *spec_list = NULL ;
01036 char tag_drl[256];
01037 char res_name_drl[256];
01038 cpl_mask* qual_mask = NULL;
01039 xsh_dispersol_list* disp_list = NULL;
01040 int rec_min_index, rec_max_index;
01041 const char* solution_type_name[2] = { "POLY", "MODEL"};
01042 int found_line=true;
01043 double barycor=0;
01044 double helicor=0;
01045
01046 XSH_ASSURE_NOT_NULL( sci_frame);
01047 XSH_ASSURE_NOT_NULL( orderlist);
01048 XSH_ASSURE_NOT_NULL( spectralformat_frame);
01049 XSH_ASSURE_NOT_NULL( rectify_par);
01050 XSH_ASSURE_NOT_NULL( instrument);
01051 XSH_ASSURE_NOT_ILLEGAL( min_index <= max_index);
01052
01053 if (rectify_par->conserve_flux){
01054 assure( disp_tab_frame != NULL,CPL_ERROR_ILLEGAL_INPUT,
01055 "If -rectify-conserve-flux=TRUE "
01056 "you must provide an input dispersion table "
01057 "tagged as DISP_TAB_%s or DISP_TAB_AFC_%s .",
01058 xsh_instrument_arm_tostring(instrument),
01059 xsh_instrument_arm_tostring(instrument));
01060 check( disp_list = xsh_dispersol_list_load( disp_tab_frame, instrument));
01061 }
01062
01063 check( pre_sci = xsh_pre_load( sci_frame, instrument));
01064 check( rec_list = xsh_rec_list_create( instrument));
01065
01066 check( xsh_baryvel(pre_sci->data_header, &barycor,&helicor));
01067 cpl_propertylist_append_double(pre_sci->data_header,XSH_QC_VRAD_BARYCOR,barycor);
01068 cpl_propertylist_set_comment(pre_sci->data_header,XSH_QC_VRAD_BARYCOR,
01069 XSH_QC_VRAD_BARYCOR_C);
01070 cpl_propertylist_append_double(pre_sci->data_header,XSH_QC_VRAD_HELICOR,helicor);
01071 cpl_propertylist_set_comment(pre_sci->data_header,XSH_QC_VRAD_HELICOR,
01072 XSH_QC_VRAD_HELICOR_C);
01073
01074 if ( model_frame != NULL) {
01075 solution_type = XSH_RECTIFY_TYPE_MODEL;
01076 check(xsh_model_temperature_update_frame(&model_frame,sci_frame,
01077 instrument,&found_line));
01078 check( xsh_model_config_load_best( model_frame, &config_model));
01079 }
01080 else if ( wavesol_frame != NULL ) {
01081 solution_type = XSH_RECTIFY_TYPE_POLY;
01082 check( wavesol = xsh_wavesol_load( wavesol_frame, instrument));
01083 }
01084 else {
01085 xsh_error_msg(
01086 "Undefined solution type (POLY or MODEL). See your input sof");
01087 }
01088
01089 xsh_msg("===== Solution type %s", solution_type_name[solution_type]);
01090
01091 check( spec_list = xsh_spectralformat_list_load( spectralformat_frame,
01092 instrument));
01093
01094
01095 xsh_msg_dbg_high( "Kernel: %s, Radius: %lf", rectify_par->rectif_kernel,
01096 rectify_par->rectif_radius);
01097 xsh_msg_dbg_high( "Slit Binning: %f, Lambda Binning: %f",
01098 rectify_par->rectif_bin_space,
01099 rectify_par->rectif_bin_lambda);
01100
01101
01102 adjust_lambdas( spec_list, rectify_par);
01103
01104
01105 rec_list->nslit = nslit;
01106 rec_list->slit_min = slit_min;
01107 slit_max = slit_min+rectify_par->rectif_bin_space*nslit;
01108 rec_list->slit_max = slit_max;
01109
01110 if (min_index >= 0){
01111 rec_min_index = min_index;
01112 }
01113 else{
01114 rec_min_index = 0;
01115 }
01116 if ( max_index < rec_list->size){
01117 rec_max_index = max_index;
01118 }
01119 else{
01120 rec_max_index = rec_list->size-1;
01121 }
01122
01123 for( i = rec_min_index; i <= rec_max_index; i++ ) {
01124 double lambda_min, lambda_max ;
01125 double n;
01126
01127 order = spec_list->list[i].absorder ;
01128 rec_list->list[i].nslit = nslit ;
01129 lambda_min = spec_list->list[i].lambda_min ;
01130 lambda_max = spec_list->list[i].lambda_max ;
01131
01132
01133 n = (lambda_max-lambda_min)/rectify_par->rectif_bin_lambda+1.0;
01134 nlambda = (int)xsh_round_double(n);
01135 assure(nlambda>0,CPL_ERROR_ILLEGAL_INPUT,
01136 "Number of wavelength sampling points is %d. Must be > 0. "
01137 "You may have to decrease rectify-bin-lambda param value.",
01138 nlambda);
01139
01140 assure(nslit>0,CPL_ERROR_ILLEGAL_INPUT,
01141 "Number of spatial sampling points is %d. Must be > 0. "
01142 "You may have to decrease rectify-bin-slit param value.",
01143 nslit);
01144
01145 rec_list->list[i].nlambda = nlambda ;
01146
01147
01148 check( xsh_rec_list_set_data_size( rec_list, i, order, nlambda, nslit));
01149
01150
01151
01152
01153
01154 xsh_msg_dbg_medium( "Order %d, Nslit: %d, Nlambda: %d", order, nslit, nlambda ) ;
01155 xsh_msg_dbg_medium( " Lambda (%f,%f) Slit(%f,%f)",
01156 lambda_min, lambda_max , slit_min, slit_max);
01157
01158 check( qual_mask = xsh_pre_get_bpmap( pre_sci));
01159 check( cpl_image_reject_from_mask( pre_sci->data, qual_mask));
01160
01161 check( fill_rectified( pre_sci, rec_list, i, wavesol,
01162 &config_model, instrument, disp_list,
01163 slit_min, slit_max, lambda_min, 0, 0, rectify_par,
01164 slit_shift, slitshift_tab));
01165 }
01166
01167 rec_list->slit_min = rec_list->list[rec_min_index].slit[0];
01168 rec_list->slit_max = rec_list->list[rec_min_index].slit[nslit-1];
01169
01170 xsh_msg_dbg_medium( "Saving Rectified Frame '%s'", res_name);
01171
01172 sprintf( tag_drl ,"%s_DRL", tag);
01173 sprintf( res_name_drl ,"DRL_%s", res_name);
01174
01175
01176 check( xsh_rec_list_update_header( rec_list, pre_sci, rectify_par, tag_drl));
01177
01178 if ( slitshift_tab != NULL){
01179 xsh_frame_set_shiftifu_ref( rec_list->header, slitshift_tab);
01180 }
01181
01182 check( res_frame = xsh_rec_list_save( rec_list, res_name_drl, tag_drl, CPL_TRUE));
01183 xsh_add_temporary_file(res_name_drl);
01184
01185 if ( res_frame_ext != NULL){
01186
01187 check( *res_frame_ext = xsh_rec_list_save2(rec_list,res_name,tag));
01188 sprintf( tag_drl ,"%s_TAB", tag);
01189 sprintf( res_name_drl ,"TAB_%s", res_name);
01190 check( *res_frame_tab = xsh_rec_list_save_table(rec_list,res_name_drl,tag_drl,
01191 CPL_TRUE));
01192
01193 xsh_add_temporary_file(res_name_drl);
01194 }
01195
01196 cleanup :
01197 xsh_dispersol_list_free( &disp_list);
01198 xsh_rec_list_free( &rec_list);
01199 xsh_wavesol_free( &wavesol);
01200 xsh_spectralformat_list_free( &spec_list);
01201 xsh_pre_free( &pre_sci);
01202 return res_frame;
01203 }
01204
01205
01206 static void xsh_get_shift_ref( cpl_frameset *set,
01207 double *down, double *cen, double *up)
01208 {
01209 cpl_frame *offsetdown_frame = NULL;
01210 cpl_frame *offsetcen_frame = NULL;
01211 cpl_frame *offsetup_frame = NULL;
01212 const char *offsetdown_name = NULL;
01213 const char *offsetcen_name = NULL;
01214 const char *offsetup_name = NULL;
01215 cpl_propertylist *offsetdown_list = NULL;
01216 cpl_propertylist *offsetcen_list = NULL;
01217 cpl_propertylist *offsetup_list = NULL;
01218
01219 XSH_ASSURE_NOT_NULL( down);
01220 XSH_ASSURE_NOT_NULL( cen);
01221 XSH_ASSURE_NOT_NULL( up);
01222
01223 check( offsetdown_frame = cpl_frameset_get_frame( set, 0));
01224 check( offsetcen_frame = cpl_frameset_get_frame( set, 1));
01225 check( offsetup_frame = cpl_frameset_get_frame( set, 2));
01226 check( offsetdown_name = cpl_frame_get_filename( offsetdown_frame));
01227 check( offsetcen_name = cpl_frame_get_filename( offsetcen_frame));
01228 check( offsetup_name = cpl_frame_get_filename( offsetup_frame));
01229 check( offsetdown_list = cpl_propertylist_load( offsetdown_name, 0));
01230 check( offsetcen_list = cpl_propertylist_load( offsetcen_name, 0));
01231 check( offsetup_list = cpl_propertylist_load( offsetup_name, 0));
01232
01233 check( *down = xsh_pfits_get_shiftifu_slitref( offsetdown_list));
01234 check( *cen = xsh_pfits_get_shiftifu_slitref( offsetcen_list));
01235 check( *up = xsh_pfits_get_shiftifu_slitref( offsetup_list));
01236
01237 cleanup:
01238 xsh_free_propertylist( &offsetdown_list);
01239 xsh_free_propertylist( &offsetcen_list);
01240 xsh_free_propertylist( &offsetup_list);
01241 return;
01242 }
01243
01244
01245
01246
01247 static void xsh_frame_set_shiftifu_ref( cpl_propertylist *header,
01248 cpl_frame *shift_frame)
01249 {
01250 cpl_propertylist *shift_header = NULL;
01251 const char *shift_name = NULL;
01252 double waveref, slitref;
01253
01254 XSH_ASSURE_NOT_NULL( header);
01255 XSH_ASSURE_NOT_NULL( shift_frame);
01256
01257 check( shift_name = cpl_frame_get_filename( shift_frame));
01258 check( shift_header = cpl_propertylist_load( shift_name, 0));
01259
01260 check( waveref = xsh_pfits_get_shiftifu_lambdaref( shift_header));
01261 check( slitref = xsh_pfits_get_shiftifu_slitref( shift_header));
01262
01263 check( xsh_pfits_set_shiftifu_lambdaref( header, waveref));
01264 check( xsh_pfits_set_shiftifu_slitref( header, slitref));
01265
01266 cleanup:
01267 xsh_free_propertylist( &shift_header);
01268 return;
01269 }
01270
01271
01272
01287 void xsh_compute_slitlet_limits( cpl_frameset *shift_set, double sdown,
01288 double sldown, double slup, double sup, double slit_bin,
01289 double *slitmin_tab, int *nslit_tab, double *slitcen_tab)
01290 {
01291 double offset_down, offset_cen, offset_up;
01292 double cen_slitlet_down_size;
01293 double cen_slitlet_up_size;
01294 double down_slitlet_down_size;
01295 double down_slitlet_up_size;
01296 double up_slitlet_down_size;
01297 double up_slitlet_up_size;
01298 double slitlet_up_size;
01299 double slitlet_down_size;
01300 double cen_ref_down, cen_ref_up;
01301 double cen_pixpos = 0.0;
01302 int ref_slit_size;
01303 int down_nslit;
01304 int up_nslit;
01305
01306 XSH_ASSURE_NOT_NULL( shift_set);
01307 XSH_ASSURE_NOT_NULL( slitmin_tab);
01308 XSH_ASSURE_NOT_NULL( nslit_tab);
01309 XSH_ASSURE_NOT_NULL( slitcen_tab);
01310
01311 check( xsh_get_shift_ref( shift_set, &offset_down, &offset_cen,
01312 &offset_up));
01313
01314
01315 xsh_msg( "Shift reference: down %f arcsec, center %f arcsec, up %f arcsec",
01316 offset_down, offset_cen, offset_up);
01317
01318 slitcen_tab[0] = offset_down;
01319 slitcen_tab[1] = offset_cen;
01320 slitcen_tab[2] =offset_up;
01321
01322 down_slitlet_down_size = offset_down-sdown;
01323 down_slitlet_up_size = sldown-offset_down;
01324
01325 xsh_msg_dbg_medium("In down slitlet [%f,%f] size lo %f up %f", sdown, sldown,
01326 down_slitlet_down_size, down_slitlet_up_size);
01327
01328 cen_slitlet_down_size = offset_cen-sldown;
01329 cen_slitlet_up_size = slup-offset_cen;
01330
01331 xsh_msg_dbg_medium("In cen slitlet [%f,%f] size lo %f up %f", sldown, slup,
01332 cen_slitlet_down_size, cen_slitlet_up_size);
01333
01334 up_slitlet_down_size = offset_up-slup;
01335 up_slitlet_up_size = sup-offset_up;
01336
01337 xsh_msg_dbg_medium("In up slitlet [%f,%f] size lo %f up %f", slup, sup,
01338 up_slitlet_down_size, up_slitlet_up_size);
01339
01340
01341 if ( cen_slitlet_down_size < up_slitlet_up_size){
01342 slitlet_down_size = cen_slitlet_down_size;
01343 }
01344 else{
01345 slitlet_down_size = up_slitlet_up_size;
01346 }
01347 if ( slitlet_down_size > down_slitlet_up_size){
01348 slitlet_down_size = down_slitlet_up_size;
01349 }
01350
01351 if ( cen_slitlet_up_size < up_slitlet_down_size){
01352 slitlet_up_size = cen_slitlet_up_size;
01353 }
01354 else{
01355 slitlet_up_size = up_slitlet_down_size;
01356 }
01357 if ( slitlet_up_size > down_slitlet_down_size){
01358 slitlet_up_size = down_slitlet_down_size;
01359 }
01360
01361 xsh_msg_dbg_medium( "Slitlet size DOWN %f UP %f", slitlet_down_size, slitlet_up_size);
01362
01363
01364 cen_ref_down = offset_cen-slitlet_down_size;
01365 cen_ref_up = offset_cen+slitlet_up_size;
01366 if (cen_ref_down > 0){
01367 down_nslit = ceil( cen_ref_down/slit_bin);
01368 }
01369 else{
01370 down_nslit = floor( cen_ref_down/slit_bin);
01371 }
01372 if ( cen_ref_up > 0){
01373 up_nslit = ceil( cen_ref_up/slit_bin);
01374 }
01375 else{
01376 up_nslit = floor( cen_ref_up/slit_bin);
01377 }
01378
01379 xsh_msg( "Adjust central reference slitlet [%f %f] with bin %f arcsec: [%f %f]",
01380 cen_ref_down, cen_ref_up, slit_bin, down_nslit*slit_bin, up_nslit*slit_bin);
01381
01382 cen_pixpos = (offset_cen-down_nslit*slit_bin)/slit_bin;
01383 ref_slit_size = up_nslit-down_nslit+1;
01384 xsh_msg( "Center position in pixel %f (Total %d pix)", cen_pixpos, ref_slit_size);
01385 slitmin_tab[1] = down_nslit*slit_bin;
01386 nslit_tab[1] = ref_slit_size;
01387
01388 slitmin_tab[0] = offset_down-(ref_slit_size-1-cen_pixpos)*slit_bin;
01389 nslit_tab[0] = ref_slit_size;
01390
01391 slitmin_tab[2] = offset_up-(ref_slit_size-1-cen_pixpos)*slit_bin;
01392 nslit_tab[2] = ref_slit_size;
01393
01394 xsh_msg("Prepare the cube bin %f arcsec", slit_bin);
01395 xsh_msg("DOWN [%f, %f]", slitmin_tab[0], slitmin_tab[0]+nslit_tab[0]*slit_bin);
01396 xsh_msg("CEN [%f, %f]", slitmin_tab[1], slitmin_tab[1]+nslit_tab[1]*slit_bin);
01397 xsh_msg("UP [%f, %f]", slitmin_tab[2], slitmin_tab[2]+nslit_tab[2]*slit_bin);
01398
01399 cleanup:
01400 return;
01401 }
01402
01403
01404
01427
01428 cpl_frameset *
01429 xsh_rectify_orders_ifu(cpl_frame *sci_frame,
01430 xsh_order_list *orderlist,
01431 cpl_frameset *wavesol_frameset,
01432 cpl_frameset *shift_frameset,
01433 cpl_frame *model_frame,
01434 xsh_instrument *instrument,
01435 xsh_rectify_param *rectify_par,
01436 cpl_frame *spectralformat_frame,
01437 cpl_frame * slitmap_frame,
01438 cpl_frameset ** res_frameset_ext,
01439 cpl_frameset ** res_frameset_tab,
01440 int min_index,
01441 int max_index,
01442 const char* rec_prefix)
01443 {
01444 cpl_frameset *res_frameset = NULL ;
01445
01446 cpl_frame *wavesolup_frame = NULL;
01447 cpl_frame *wavesolcen_frame = NULL;
01448 cpl_frame *wavesoldown_frame = NULL;
01449 cpl_frame *shiftup_frame = NULL;
01450 cpl_frame *shiftcen_frame = NULL;
01451 cpl_frame *shiftdown_frame = NULL;
01452
01453 int slitlet=0;
01454 double sdown=0, sldown=0, slup=0, sup=0;
01455 double slitlet_down_l=0, slitlet_down_u=0;
01456
01457 double slitlet_cen_l=0, slitlet_cen_u=0;
01458 double slitlet_cen_u_bin=0;
01459 double slitlet_up_u=0;
01460 double down_offset=0.0;
01461 int nslit=0;
01462 double bin_space;
01463 XSH_ARM arm ;
01464 int nslit_tab[3];
01465 double slitmin_tab[3];
01466 double slitcen_tab[3];
01467 int cut_edges = 0;
01468
01469 XSH_ASSURE_NOT_NULL( sci_frame);
01470 XSH_ASSURE_NOT_NULL( orderlist);
01471 XSH_ASSURE_NOT_NULL( rectify_par);
01472 XSH_ASSURE_NOT_NULL( instrument);
01473 XSH_ASSURE_NOT_NULL( spectralformat_frame);
01474
01475 bin_space = rectify_par->rectif_bin_space;
01476
01477
01478 arm = xsh_instrument_get_arm( instrument ) ;
01479
01480 if ( wavesol_frameset != NULL){
01481 check( wavesoldown_frame = cpl_frameset_get_frame( wavesol_frameset, 0));
01482 check( wavesolcen_frame = cpl_frameset_get_frame( wavesol_frameset, 1));
01483 check( wavesolup_frame = cpl_frameset_get_frame( wavesol_frameset, 2));
01484 }
01485
01486 if ( shift_frameset != NULL){
01487 check( shiftdown_frame = cpl_frameset_get_frame( shift_frameset, 0));
01488 check( shiftcen_frame = cpl_frameset_get_frame( shift_frameset, 1));
01489 check( shiftup_frame = cpl_frameset_get_frame( shift_frameset, 2));
01490 }
01491
01492
01493 check( res_frameset = cpl_frameset_new());
01494
01495 check( xsh_get_slit_edges( slitmap_frame, &sdown,
01496 &sup, &sldown, &slup, instrument));
01497
01498 if ( shift_frameset != NULL){
01499 check( xsh_compute_slitlet_limits( shift_frameset, sdown,
01500 sldown, slup, sup, bin_space,
01501 slitmin_tab, nslit_tab, slitcen_tab));
01502 }
01503 else{
01504 slitlet_cen_l = sldown;
01505 slitlet_cen_u = slup;
01506
01507 slitlet_down_l = sdown;
01508 slitlet_down_u = slitlet_cen_l;
01509 slitlet_up_u = sup;
01510
01511
01512 slitmin_tab[1] = slitlet_cen_l;
01513 nslit_tab[1] = ceil((slitlet_cen_u-slitlet_cen_l)/rectify_par->rectif_bin_space);
01514 slitlet_cen_u_bin = slitmin_tab[1]+nslit_tab[1]*rectify_par->rectif_bin_space;
01515
01516
01517 nslit_tab[0] = nslit_tab[1];
01518 slitmin_tab[0] = slitlet_down_u-down_offset-
01519 nslit_tab[0]*rectify_par->rectif_bin_space;
01520
01521
01522 slitmin_tab[2] = slitlet_cen_u_bin;
01523 nslit_tab[2] = nslit_tab[1];
01524 }
01525
01526 for ( slitlet = LOWER_IFU_SLITLET ; slitlet <= UPPER_IFU_SLITLET ;
01527 slitlet++ ) {
01528
01529 cpl_frame *wavesol_frame = NULL;
01530 cpl_frame *shift_frame = NULL;
01531 cpl_frame *res_frame = NULL;
01532 cpl_frame *resext_frame = NULL;
01533 cpl_frame *restab_frame = NULL;
01534 double slitlet_slit_min = 0.0;
01535 char tag[80];
01536 char res_name[80];
01537 double rec_shift=0.0;
01538 double rec_min=0.0;
01539 const char *slitlet_name = NULL;
01540
01541 switch( slitlet ) {
01542 case LOWER_IFU_SLITLET:
01543
01544 slitlet_name = "DOWN";
01545 nslit = nslit_tab[0];
01546 slitlet_slit_min = slitmin_tab[0];
01547 wavesol_frame = wavesoldown_frame;
01548 shift_frame = shiftdown_frame;
01549 sprintf(tag,"%s_%s",rec_prefix,
01550 XSH_GET_TAG_FROM_ARM(XSH_ORDER2D_DOWN_IFU,instrument));
01551 sprintf(res_name,"%s.fits",tag);
01552 break ;
01553
01554 case CENTER_IFU_SLITLET:
01555 slitlet_name = "CEN";
01556 slitlet_slit_min = slitmin_tab[1];
01557 nslit = nslit_tab[1];
01558 wavesol_frame = wavesolcen_frame;
01559 shift_frame = shiftcen_frame;
01560 sprintf(tag,"%s_%s",rec_prefix,
01561 XSH_GET_TAG_FROM_ARM(XSH_ORDER2D_CEN_IFU,instrument));
01562 sprintf(res_name,"%s.fits",tag);
01563 break ;
01564
01565 case UPPER_IFU_SLITLET:
01566 slitlet_name = "UP";
01567 slitlet_slit_min = slitmin_tab[2];
01568 nslit = nslit_tab[2];
01569
01570 wavesol_frame = wavesolup_frame;
01571 shift_frame = shiftup_frame;
01572 sprintf(tag,"%s_%s",rec_prefix,
01573 XSH_GET_TAG_FROM_ARM(XSH_ORDER2D_UP_IFU,instrument));
01574 sprintf(res_name,"%s.fits",tag);
01575 break ;
01576 }
01577
01578 rec_min = slitlet_slit_min;
01579 rec_shift = 0;
01580
01581 xsh_msg("Cut edges with %d pix", cut_edges);
01582 rec_min += cut_edges*rectify_par->rectif_bin_space;
01583 nslit -= cut_edges*2;
01584
01585 XSH_CMP_INT( nslit, >, 0, "Check size in slit",);
01586
01587
01588 xsh_msg( "%s [%f %f] in %d pix of %f arcsec" ,slitlet_name,
01589 slitlet_slit_min,
01590 slitlet_slit_min+nslit*rectify_par->rectif_bin_space,
01591 nslit, rectify_par->rectif_bin_space);
01592
01593 check( res_frame = xsh_rectify_orders( sci_frame, orderlist,
01594 wavesol_frame, model_frame, instrument, rectify_par,
01595 spectralformat_frame, NULL, res_name, tag, &resext_frame, &restab_frame,
01596 min_index, max_index, rec_min, nslit, rec_shift,shift_frame));
01597
01598 check( cpl_frameset_insert( res_frameset, res_frame));
01599
01600 if(res_frameset_ext !=NULL && res_frameset_tab !=NULL) {
01601 check( cpl_frameset_insert(*res_frameset_ext, resext_frame));
01602 check( cpl_frameset_insert(*res_frameset_tab, restab_frame));
01603 }
01604 }
01605
01606 cleanup:
01607 return res_frameset ;
01608 }
01609
01610
01611
01612
01620
01621 static double
01622 compute_shift_with_localization( cpl_frame *loc_frame,
01623 cpl_frame *loc0_frame)
01624 {
01625 double shift_a=0, shift_b=0, shift;
01626 xsh_localization *loc = NULL, *loc0 = NULL;
01627
01628 XSH_ASSURE_NOT_NULL( loc_frame);
01629 XSH_ASSURE_NOT_NULL( loc0_frame);
01630
01631 check( loc = xsh_localization_load( loc_frame));
01632 check( loc0 = xsh_localization_load( loc0_frame));
01633
01634
01635 check( shift_a = cpl_polynomial_eval_1d( loc0->cenpoly, 0., NULL));
01636 check( shift_b = cpl_polynomial_eval_1d( loc->cenpoly, 0., NULL));
01637
01638 xsh_msg_dbg_medium("Shift A %f B %f", shift_a, shift_b)
01639 shift = shift_b-shift_a;
01640
01641 xsh_msg( "Shift (localization) = %lf", shift);
01642
01643 cleanup:
01644 xsh_localization_free( &loc);
01645 xsh_localization_free( &loc0);
01646 return shift;
01647 }
01648
01649
01650
01651
01652
01653
01654
01655
01656
01657
01658
01659
01660
01661
01662
01663
01664
01665
01666
01667
01668
01669
01670
01671
01672
01673
01674
01675
01676
01677
01678
01679
01680
01681
01682
01683
01684
01685
01686
01687
01688
01689
01690
01691
01692
01693
01694
01695
01696
01697
01698
01699
01700
01701
01702
01703
01704
01705
01706
01707
01708
01709
01710
01711
01712
01713
01714
01715
01716
01717
01718
01719
01720
01721
01722
01723
01724
01725
01726
01727
01728
01729
01730
01731
01732
01733
01734
01735
01736
01737
01738
01739
01740
01741
01742
01743
01744
01745
01746
01747
01748
01749
01750
01751
01752
01753
01754
01755
01756
01757
01758
01759
01760
01761
01762
01763
01764
01765
01766
01767
01768
01769
01770
01771
01772
01773
01774
01775
01776
01777
01778
01779
01780
01781
01782
01783
01784
01785
01786
01787
01788
01789
01790
01791
01792
01793
01794
01795
01796
01797
01798
01799 static double
01800 compute_shift_with_kw( cpl_propertylist *header,
01801 xsh_rectify_param *rectify_par,
01802 double **ref_ra,
01803 double **ref_dec,
01804 int flag)
01805 {
01806 double ref_ra_cumoff, ref_dec_cumoff;
01807 double ra_cumoff, dec_cumoff;
01808 double ra_reloff, dec_reloff, ra_off, dec_off;
01809 double posang=0;
01810
01811 double shift_arcsec=0;
01812
01813
01814
01815 xsh_msg_dbg_high( "==> compute_shift_with_kw" ) ;
01816 XSH_ASSURE_NOT_NULL( header);
01817 XSH_ASSURE_NOT_NULL( rectify_par);
01818
01819 check( posang = xsh_pfits_get_posang( header));
01820 posang = posang/180.0* M_PI;
01821
01822 if ( *ref_ra != NULL){
01823 ref_ra_cumoff = **ref_ra;
01824 }
01825 else{
01826 check( ref_ra_cumoff = xsh_pfits_get_ra_cumoffset( header));
01827 XSH_MALLOC( *ref_ra, double, 1);
01828 **ref_ra = ref_ra_cumoff;
01829 }
01830 if ( *ref_dec != NULL){
01831 ref_dec_cumoff = **ref_dec;
01832 }
01833 else{
01834 check( ref_dec_cumoff = xsh_pfits_get_dec_cumoffset( header));
01835 XSH_MALLOC( *ref_dec, double, 1);
01836 **ref_dec = ref_dec_cumoff;
01837 }
01838
01839 if ( flag == 0){
01840 check( ra_cumoff = xsh_pfits_get_ra_cumoffset( header));
01841 check( dec_cumoff = xsh_pfits_get_dec_cumoffset( header));
01842 }
01843 else{
01844 check( ra_cumoff = xsh_pfits_get_b_ra_cumoffset( header));
01845 check( dec_cumoff = xsh_pfits_get_b_dec_cumoffset( header));
01846 }
01847
01848 check( ra_reloff = xsh_pfits_get_b_ra_reloffset( header));
01849 check( dec_reloff = xsh_pfits_get_b_dec_reloffset( header));
01850
01851
01852 xsh_msg( "POSANG %f (rad) REF CUM_(RA DEC) : %f %f ", posang, ref_ra_cumoff,
01853 ref_dec_cumoff);
01854 xsh_msg( "OBJ CUM_(RA DEC) : %f %f ", ra_cumoff, dec_cumoff);
01855 xsh_msg( "REL_(RA DEC) : %f %f", ra_reloff, dec_reloff);
01856
01857 ra_off = ra_cumoff-ref_ra_cumoff;
01858 dec_off = dec_cumoff-ref_dec_cumoff;
01859
01860 xsh_msg( "COMPUTE OFF_(RA DEC) : %f %f", ra_off, dec_off);
01861
01862 shift_arcsec = cos(-posang)*dec_off+sin(-posang)*ra_off;
01863
01864 xsh_msg( "COMPUTE shift in arcsec %f :", shift_arcsec);
01865
01866 cleanup:
01867 return shift_arcsec;
01868 }
01869
01870
01871
01886 cpl_frame*
01887 shift_with_kw( cpl_frame *rec_frame,
01888 xsh_instrument *instrument,
01889 xsh_rectify_param *rectify_par,
01890 const char *fname,
01891 cpl_frame** res_frame_ext,
01892 double **ref_ra,
01893 double **ref_dec,
01894 int flag)
01895 {
01896 cpl_frame *result = NULL ;
01897 int shift_pix;
01898 double bin_space, shift_arcsec;
01899 const char *filename = NULL;
01900 cpl_propertylist* header = NULL;
01901
01902
01903 XSH_ASSURE_NOT_NULL( rec_frame);
01904 XSH_ASSURE_NOT_NULL( instrument);
01905 XSH_ASSURE_NOT_NULL( fname);
01906 XSH_ASSURE_NOT_NULL( res_frame_ext);
01907 XSH_ASSURE_NOT_NULL( ref_ra);
01908 XSH_ASSURE_NOT_NULL( ref_dec);
01909
01910 check( filename = cpl_frame_get_filename( rec_frame));
01911 check( header = cpl_propertylist_load( filename, 0));
01912
01913 check( bin_space = xsh_pfits_get_rectify_bin_space( header));
01914
01915 check( shift_arcsec = compute_shift_with_kw( header,
01916 rectify_par, ref_ra, ref_dec, flag));
01917
01918
01919 shift_pix = xsh_round_double(shift_arcsec/bin_space);
01920
01921 shift_arcsec = shift_pix*bin_space;
01922
01923 xsh_msg( "SHIFT A-->B : %f in arcsec", shift_arcsec);
01924
01925 check( result = xsh_shift( rec_frame, instrument, fname, shift_arcsec,
01926 res_frame_ext));
01927
01928 cleanup:
01929 if ( cpl_error_get_code() !=CPL_ERROR_NONE){
01930 xsh_free_frame( &result);
01931 }
01932 xsh_free_propertylist( &header);
01933 return result ;
01934 }
01935
01936
01937
01938
01939 static cpl_frame*
01940 xsh_shift( cpl_frame *rec_frame,
01941 xsh_instrument *instrument,
01942 const char *fname,
01943 double slit_shift,
01944 cpl_frame** res_frame_ext)
01945 {
01946 xsh_rec_list *rec_list = NULL ;
01947 cpl_frame *result = NULL ;
01948 int nb_orders, order = 0 ;
01949 float *slit = NULL ;
01950
01951 char* fname_drl=NULL;
01952 char* tag_drl=NULL;
01953 char* tag=NULL;
01954 int nslit=0;
01955
01956 XSH_ASSURE_NOT_NULL( rec_frame);
01957 XSH_ASSURE_NOT_NULL( instrument);
01958 XSH_ASSURE_NOT_NULL( fname);
01959 XSH_ASSURE_NOT_NULL( res_frame_ext);
01960
01961 check( rec_list = xsh_rec_list_load( rec_frame, instrument));
01962
01963 nb_orders = rec_list->size ;
01964 check( nslit = xsh_rec_list_get_nslit( rec_list, 0));
01965
01966
01967 for( order = 0 ; order < nb_orders ; order++ ) {
01968 int ns;
01969
01970 check( slit = xsh_rec_list_get_slit( rec_list, order));
01971 for( ns = 0 ; ns < nslit ; ns++ ) {
01972 slit[ns]+=slit_shift;
01973 }
01974 }
01975 rec_list->slit_min = slit[0];
01976 rec_list->slit_max = slit[nslit-1];
01977
01978 check( xsh_pfits_set_rectify_space_min( rec_list->header,
01979 rec_list->slit_min));
01980 check( xsh_pfits_set_rectify_space_max( rec_list->header,
01981 rec_list->slit_max));
01982
01983
01984 tag= xsh_stringcat_any( cpl_frame_get_tag( rec_frame ), "", NULL ) ;
01985 tag_drl = xsh_stringcat_any( cpl_frame_get_tag( rec_frame ), "_DRL", NULL ) ;
01986 fname_drl = xsh_stringcat_any( "DRL_", fname, NULL);
01987
01988 check( *res_frame_ext=xsh_rec_list_save2( rec_list, fname,tag)) ;
01989 check( result = xsh_rec_list_save( rec_list, fname_drl,
01990 tag_drl, CPL_TRUE));
01991 xsh_add_temporary_file(fname_drl);
01992
01993 cleanup:
01994 if ( cpl_error_get_code() !=CPL_ERROR_NONE){
01995 xsh_free_frame( &result);
01996 }
01997 XSH_FREE( fname_drl);
01998 XSH_FREE( tag_drl);
01999 XSH_FREE( tag);
02000 xsh_rec_list_free( &rec_list);
02001 return result ;
02002 }
02003
02020 cpl_frame * xsh_shift_rectified( cpl_frame * rec_frame,
02021 cpl_frame * loc_frame,
02022 cpl_frame * loc0_frame,
02023 const char * file_name,
02024 xsh_combine_nod_param * combine_nod_param,
02025 xsh_rectify_param * rectif_par,
02026 xsh_instrument * instrument,
02027 cpl_frame** res_frame_ext )
02028 {
02029 cpl_frame * result = NULL ;
02030 double *ra = NULL, *dec = NULL;
02031
02032 xsh_msg( "===> xsh_shift_rectified" ) ;
02033
02034 XSH_ASSURE_NOT_NULL( combine_nod_param ) ;
02035 XSH_ASSURE_NOT_NULL( rec_frame ) ;
02036
02037 check( result = shift_with_kw( rec_frame, instrument,
02038 rectif_par, file_name, res_frame_ext, &ra, &dec, 1));
02039
02040 cleanup:
02041 return result ;
02042 }
02043
02044