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
00035
00036
00040
00041
00042
00043 #include <string.h>
00044 #include <math.h>
00045 #include <xsh_utils_table.h>
00046 #include <xsh_parameters.h>
00047 #include <xsh_data_pre.h>
00048 #include <xsh_data_order.h>
00049 #include <xsh_data_dispersol.h>
00050 #include <xsh_utils_wrappers.h>
00051 #include <xsh_data_wavemap.h>
00052 #include <xsh_error.h>
00053 #include <xsh_msg.h>
00054 #include <xsh_pfits.h>
00055 #include <xsh_dfs.h>
00056 #include <math.h>
00057 #include <gsl/gsl_multifit.h>
00058 #include <cpl.h>
00059
00060
00061
00062
00069 void
00070 xsh_wavemap_list_dump( xsh_wavemap_list * list,
00071 const char * fname )
00072 {
00073 int i ;
00074 FILE * fout = NULL;
00075
00076 XSH_ASSURE_NOT_NULL( list ) ;
00077 if ( fname == NULL ) fout = stdout ;
00078 else fout = fopen( fname, "w" ) ;
00079 XSH_ASSURE_NOT_NULL( fout ) ;
00080
00081 fprintf( fout, "Wavemap List. Nb of orders: %d\n", list->size ) ;
00082 for( i = 0 ; i<list->size ; i++ ) {
00083 fprintf( fout, " Entry %2d: Order %d, Ndata: %d\n",
00084 i, list->list[i].order, list->list[i].sky_size);
00085 }
00086
00087 cleanup:
00088 if ( fname != NULL && fout != NULL ) fclose( fout ) ;
00089 return ;
00090 }
00091
00092
00093
00100
00101 xsh_wavemap_list*
00102 xsh_wavemap_list_create( xsh_instrument* instr)
00103 {
00104 xsh_wavemap_list* result = NULL;
00105 XSH_INSTRCONFIG * config = NULL;
00106
00107 int i, size = 0;
00108
00109
00110 XSH_ASSURE_NOT_NULL( instr);
00111
00112
00113
00114 check( config = xsh_instrument_get_config( instr));
00115 size = config->orders;
00116
00117 XSH_CALLOC(result, xsh_wavemap_list, 1 );
00118
00119 result->size = size;
00120
00121 XSH_ASSURE_NOT_ILLEGAL(result->size > 0);
00122 result->instrument = instr;
00123 XSH_CALLOC(result->list, xsh_wavemap, result->size);
00124 XSH_NEW_PROPERTYLIST(result->header);
00125
00126
00127 for( i = 0 ; i<result->size ; i++ ) {
00128 result->list[i].tcheb_pol_lambda = NULL ;
00129 }
00130
00131 cleanup:
00132 if (cpl_error_get_code() != CPL_ERROR_NONE) {
00133 xsh_wavemap_list_free(&result);
00134 }
00135 return result;
00136 }
00137
00147
00148 void
00149 xsh_wavemap_list_set_max_size( xsh_wavemap_list * list, int idx,
00150 int absorder, int max_size)
00151 {
00152 xsh_wavemap * pwavemap = NULL ;
00153
00154
00155 XSH_ASSURE_NOT_NULL(list) ;
00156 XSH_ASSURE_NOT_ILLEGAL( idx >= 0 && idx < list->size && max_size > 0);
00157 pwavemap = &list->list[idx] ;
00158 XSH_ASSURE_NOT_NULL( pwavemap);
00159
00160 pwavemap->order = absorder;
00161 pwavemap->sky_size = max_size;
00162 pwavemap->object_size = max_size;
00163 XSH_CALLOC( pwavemap->sky, wavemap_item, max_size);
00164 XSH_CALLOC( pwavemap->object, wavemap_item, max_size);
00165
00166 cleanup:
00167 return ;
00168 }
00169
00170
00175
00176 void
00177 xsh_wavemap_list_free(xsh_wavemap_list** list)
00178 {
00179 int i = 0;
00180 xsh_wavemap_list *plist = NULL;
00181
00182 if (list != NULL && *list != NULL ) {
00183 plist = *list ;
00184
00185 for (i = 0; i < plist->size; i++) {
00186 xsh_wavemap *pr = &(plist->list[i]) ;
00187
00188 xsh_msg_dbg_high( "Freeing order index %d", i ) ;
00189 if ( pr == NULL ) continue ;
00190 xsh_msg_dbg_high( " Abs Order: %d", pr->order ) ;
00191 cpl_free( pr->sky);
00192 cpl_free( pr->object);
00193 if (pr->pol_lambda != NULL){
00194 xsh_free_polynomial( &(pr->pol_lambda));
00195 }
00196 if (pr->pol_slit != NULL){
00197 xsh_free_polynomial( &(pr->pol_slit));
00198 }
00199 xsh_free_polynomial( &(pr->tcheb_pol_lambda));
00200 }
00201 if ((*list)->list){
00202 cpl_free((*list)->list);
00203 }
00204 xsh_free_propertylist(&((*list)->header));
00205 cpl_free(*list);
00206 *list = NULL;
00207 }
00208 }
00209
00210 static void
00211 lambda_fit( double * lambda, double * xpos, double * ypos,
00212 int size, xsh_wavemap_list * wmap,
00213 xsh_dispersol_param * dispsol_param, int order)
00214 {
00215 double xmin, xmax, ymin, ymax, lmin, lmax;
00216 double *tcheb_lambda = NULL, *tcheb_xpos = NULL, *tcheb_ypos = NULL;
00217 int nbcoefs, xdeg, ydeg;
00218 cpl_vector *tcheb_coef_x = NULL;
00219 cpl_vector *tcheb_coef_y = NULL;
00220 int i, j, k ;
00221 double chisq;
00222 gsl_matrix *X = NULL, *cov = NULL;
00223 gsl_vector *y = NULL, *c = NULL;
00224 gsl_multifit_linear_workspace *work = NULL;
00225 cpl_polynomial* result = NULL;
00226 cpl_size pows[3];
00227
00228 XSH_ASSURE_NOT_NULL( lambda);
00229 XSH_ASSURE_NOT_NULL( xpos);
00230 XSH_ASSURE_NOT_NULL( ypos);
00231 XSH_ASSURE_NOT_NULL( wmap);
00232 XSH_ASSURE_NOT_NULL( dispsol_param);
00233 XSH_ASSURE_NOT_ILLEGAL( size > 0);
00234 XSH_ASSURE_NOT_ILLEGAL( order >= 0 && order < wmap->size);
00235
00236 xsh_msg_dbg_low( " Entering lambda_fit. Order %d, size: %d",
00237 wmap->list[order].order, size);
00238 xdeg = dispsol_param->deg_x;
00239 ydeg = dispsol_param->deg_y;
00240 nbcoefs = (xdeg+1)*(ydeg+1) ;
00241
00242 XSH_ASSURE_NOT_ILLEGAL( size >= nbcoefs);
00243
00244 wmap->degx = xdeg ;
00245 wmap->degy = ydeg ;
00246
00247
00248 check(xsh_tools_min_max(size, lambda, &lmin, &lmax));
00249 check(xsh_tools_min_max(size, xpos, &xmin, &xmax));
00250 check(xsh_tools_min_max(size, ypos, &ymin, &ymax));
00251
00252 wmap->list[order].xmin = xmin ;
00253 wmap->list[order].xmax = xmax ;
00254 wmap->list[order].ymin = ymin ;
00255 wmap->list[order].ymax = ymax ;
00256 wmap->list[order].lambda_min = lmin ;
00257 wmap->list[order].lambda_max = lmax ;
00258
00259 xsh_msg_dbg_medium(
00260 " Min,Max. Lambda: %.6lf,%.6lf - X: %.6lf,%.6lf - Y: %.6lf,%.6lf",
00261 lmin, lmax, xmin, xmax, ymin, ymax);
00262
00263
00264 XSH_CALLOC(tcheb_lambda, double, size);
00265 XSH_CALLOC(tcheb_xpos, double, size);
00266 XSH_CALLOC(tcheb_ypos, double, size);
00267
00268 check(xsh_tools_tchebitchev_transform_tab(size, lambda, lmin, lmax,
00269 tcheb_lambda));
00270 check(xsh_tools_tchebitchev_transform_tab(size, xpos, xmin, xmax,
00271 tcheb_xpos));
00272 check(xsh_tools_tchebitchev_transform_tab(size, ypos, ymin, ymax,
00273 tcheb_ypos));
00274 X = gsl_matrix_alloc( size, nbcoefs);
00275 y = gsl_vector_alloc( size);
00276 c = gsl_vector_alloc( nbcoefs);
00277 cov = gsl_matrix_alloc( nbcoefs, nbcoefs);
00278
00279 for (i = 0; i < size; i++) {
00280 check( tcheb_coef_x = xsh_tools_tchebitchev_poly_eval(
00281 xdeg, tcheb_xpos[i]));
00282 check( tcheb_coef_y = xsh_tools_tchebitchev_poly_eval(
00283 ydeg, tcheb_ypos[i]));
00284
00285 for(j = 0; j < dispsol_param->deg_x +1 ; j++) {
00286 for(k = 0; k < dispsol_param->deg_y +1 ; k++){
00287 double vx, vy;
00288
00289 check( vx = cpl_vector_get( tcheb_coef_x, j));
00290 check( vy = cpl_vector_get( tcheb_coef_y, k));
00291 gsl_matrix_set (X, i, k+j*(xdeg+1), vx*vy );
00292 }
00293 }
00294 gsl_vector_set (y, i, tcheb_lambda[i]);
00295 xsh_free_vector(&tcheb_coef_x);
00296 xsh_free_vector(&tcheb_coef_y);
00297
00298 }
00299 xsh_msg_dbg_medium( " GSL Initialized" ) ;
00300
00301
00302 XSH_ASSURE_NOT_ILLEGAL_MSG( size >nbcoefs,
00303 "Not enough Points vs Number of Coeffs" ) ;
00304 work = gsl_multifit_linear_alloc(size, nbcoefs);
00305 gsl_multifit_linear(X, y, c, cov, &chisq, work);
00306
00307 xsh_msg_dbg_medium( " GSL Done" ) ;
00308
00309 wmap->list[order].tcheb_pol_lambda = cpl_polynomial_new( 2);
00310
00311 result = wmap->list[order].tcheb_pol_lambda ;
00312
00313 for(i = 0; i < (xdeg +1) ; i++){
00314 for(j = 0; j < (ydeg +1) ; j++){
00315 pows[0] = j;
00316 pows[1] = i;
00317 cpl_polynomial_set_coeff(result, pows,
00318 gsl_vector_get(c, j+i*(ydeg+1)));
00319 }
00320 }
00321 xsh_msg_dbg_medium( " GSL Finished" ) ;
00322
00323 gsl_multifit_linear_free (work);
00324 gsl_matrix_free (X);
00325 gsl_vector_free (y);
00326 gsl_vector_free (c);
00327 gsl_matrix_free (cov);
00328
00329 cleanup:
00330 XSH_FREE( tcheb_lambda);
00331 XSH_FREE( tcheb_xpos);
00332 XSH_FREE( tcheb_ypos);
00333 return ;
00334 }
00335
00336
00352
00353
00354
00355
00356 void
00357 xsh_wavemap_list_compute_poly( double * vlambda,
00358 double * vslit,
00359 double * xpos,
00360 double * ypos,
00361 int nitems,
00362 double * orders,
00363 xsh_dispersol_param *dispsol_param,
00364 xsh_wavemap_list *wmap)
00365 {
00366
00367
00368
00369
00370 int i, ordersize, nborder ;
00371 double *vx = NULL, *vy = NULL, *vl = NULL, *vs = NULL;
00372
00373
00374 xsh_msg( "Entering xsh_wavemap_compute" ) ;
00375
00376 XSH_ASSURE_NOT_NULL( vlambda);
00377 XSH_ASSURE_NOT_NULL( vslit);
00378 XSH_ASSURE_NOT_NULL( xpos);
00379 XSH_ASSURE_NOT_NULL( ypos);
00380 XSH_ASSURE_NOT_NULL( orders);
00381 XSH_ASSURE_NOT_NULL( wmap);
00382 XSH_ASSURE_NOT_NULL( dispsol_param);
00383 XSH_ASSURE_NOT_ILLEGAL( nitems > 0);
00384 XSH_ASSURE_NOT_ILLEGAL( wmap->size);
00385
00386 xsh_msg( " X degree = %d, Y degree = %d", dispsol_param->deg_x,
00387 dispsol_param->deg_y) ;
00388
00389
00390 ordersize = 0;
00391 nborder = 0;
00392
00393 wmap->degx = dispsol_param->deg_x;
00394 wmap->degy = dispsol_param->deg_y;
00395
00396 xsh_msg("Compute POLY for lambda");
00397 XSH_REGDEBUG( "nitems %d ", nitems);
00398
00399 for( i = 1 ; i<=nitems ; i++ ) {
00400
00401 if ( i < nitems && orders[i-1] == orders[i] ) {
00402 ordersize++;
00403 }
00404 else {
00405 cpl_vector* posx = NULL;
00406 cpl_vector* posy = NULL;
00407 cpl_vector* lvalues = NULL;
00408 cpl_vector* svalues = NULL;
00409 cpl_bivector *positions = NULL;
00410 double mse =0.0;
00411
00412 ordersize++;
00413 int i_minus_ordersize=i-ordersize;
00414 XSH_MALLOC( vx, double, ordersize);
00415 memcpy( vx, &(xpos[i_minus_ordersize]), ordersize*sizeof(double));
00416
00417 XSH_MALLOC( vy, double, ordersize);
00418 memcpy( vy, &(ypos[i_minus_ordersize]), ordersize*sizeof(double));
00419
00420 XSH_MALLOC( vl, double, ordersize);
00421 memcpy( vl, &(vlambda[i_minus_ordersize]), ordersize*sizeof(double) ) ;
00422
00423 XSH_MALLOC( vs, double, ordersize);
00424 memcpy( vs, &( vslit[i_minus_ordersize]), ordersize*sizeof(double) ) ;
00425 wmap->list[nborder].sky_size = ordersize ;
00426
00427
00428 wmap->list[nborder].order = orders[i-1] ;
00429 xsh_msg_dbg_high( "Order: %d, Size: %d", wmap->list[nborder].order,
00430 ordersize );
00431
00432 posx = cpl_vector_wrap( ordersize, vx);
00433 posy = cpl_vector_wrap( ordersize, vy);
00434
00435 positions = cpl_bivector_wrap_vectors( posx, posy);
00436 lvalues = cpl_vector_wrap( ordersize, vl);
00437 svalues = cpl_vector_wrap( ordersize, vs);
00438 cpl_size loc_degx=(cpl_size)dispsol_param->deg_x;
00439 wmap->list[nborder].pol_lambda =
00440 xsh_polynomial_fit_2d_create( positions, lvalues, &loc_degx,
00441 &mse);
00442
00443 wmap->list[nborder].pol_slit =
00444 xsh_polynomial_fit_2d_create( positions, svalues, &loc_degx,
00445 &mse);
00446
00447 cpl_bivector_unwrap_vectors( positions);
00448
00449 cpl_vector_unwrap( posx);
00450 cpl_vector_unwrap( posy);
00451 cpl_vector_unwrap( lvalues);
00452 cpl_vector_unwrap( svalues);
00453
00454 XSH_FREE( vx);
00455 XSH_FREE( vy);
00456 XSH_FREE( vl);
00457 XSH_FREE( vs);
00458 nborder++ ;
00459 ordersize = 0 ;
00460
00461 }
00462 }
00463
00464 cleanup:
00465 XSH_FREE( vx);
00466 XSH_FREE( vy);
00467 XSH_FREE( vl);
00468 XSH_FREE( vs);
00469 return ;
00470 }
00471
00472
00486
00487 void
00488 xsh_wavemap_list_compute( double * vlambda,
00489 double * xpos,
00490 double * ypos,
00491 int nitems,
00492 double * orders,
00493 xsh_dispersol_param * dispsol_param,
00494 xsh_wavemap_list * wmap )
00495 {
00496 int i, ordersize, nborder ;
00497 double *vx = NULL, *vy = NULL, *vl = NULL;
00498
00499 xsh_msg( "Entering xsh_wavemap_compute" ) ;
00500
00501 XSH_ASSURE_NOT_NULL( vlambda ) ;
00502 XSH_ASSURE_NOT_NULL( xpos ) ;
00503 XSH_ASSURE_NOT_NULL( ypos ) ;
00504 XSH_ASSURE_NOT_NULL( orders ) ;
00505 XSH_ASSURE_NOT_NULL( wmap ) ;
00506 XSH_ASSURE_NOT_NULL( dispsol_param ) ;
00507 XSH_ASSURE_NOT_ILLEGAL( nitems > 0 ) ;
00508 XSH_ASSURE_NOT_ILLEGAL( wmap->size ) ;
00509
00510 xsh_msg( " X degree = %d, Y degree = %d", dispsol_param->deg_x,
00511 dispsol_param->deg_y) ;
00512
00513
00514
00515
00516
00517 ordersize = 0;
00518 nborder = 0;
00519
00520 XSH_REGDEBUG( "nitems %d ", nitems);
00521
00522 for( i = 1 ; i<=nitems ; i++ ) {
00523
00524 if ( i < nitems && orders[i-1] == orders[i] ) {
00525 ordersize++;
00526 }
00527 else {
00528 ordersize++;
00529 XSH_MALLOC( vx, double, ordersize);
00530 memcpy( vx, &(xpos[i-ordersize]), ordersize*sizeof(double));
00531
00532 XSH_MALLOC( vy, double, ordersize);
00533 memcpy( vy, &(ypos[i-ordersize]), ordersize*sizeof(double));
00534
00535 XSH_MALLOC( vl, double, ordersize);
00536 memcpy( vl, &(vlambda[i-ordersize]), ordersize*sizeof(double) ) ;
00537
00538 wmap->list[nborder].sky_size = ordersize;
00539
00540 wmap->list[nborder].order = orders[i-1] ;
00541 xsh_msg_dbg_high( "Order: %d, Size: %d", wmap->list[nborder].order,
00542 ordersize );
00543
00544 check( lambda_fit( vl, vx, vy, ordersize, wmap, dispsol_param, nborder ) ) ;
00545 XSH_FREE( vx ) ;
00546 XSH_FREE( vy ) ;
00547 XSH_FREE( vl ) ;
00548 nborder++ ;
00549 ordersize = 0 ;
00550 }
00551 }
00552
00553 cleanup:
00554 XSH_FREE( vx ) ;
00555 XSH_FREE( vy ) ;
00556 XSH_FREE( vl ) ;
00557 return ;
00558 }
00559
00560
00570
00571
00572 static double
00573 xsh_wavemap_list_eval_lambda( xsh_wavemap_list * wmap,
00574 double x, double y, int ord )
00575 {
00576 double tcheb_x = 0.0, tcheb_y = 0.0;
00577 double lambda = -1.;
00578 double res=0.0;
00579 int i, j ;
00580 cpl_size pows[2] ;
00581 int degx, degy ;
00582 cpl_vector *tcheb_coef_x = NULL;
00583 cpl_vector *tcheb_coef_y = NULL;
00584
00585 XSH_ASSURE_NOT_NULL( wmap);
00586 degx = wmap->degx;
00587 degy = wmap->degy;
00588 XSH_ASSURE_NOT_ILLEGAL( degx >0 && degy > 0);
00589
00590 xsh_msg_dbg_high( " eval_lambda: degx: %d, min: %lf, max: %lf",
00591 degx, wmap->list[ord].xmin,
00592 wmap->list[ord].xmax ) ;
00593
00594 if ( wmap->list[ord].xmin >= wmap->list[ord].xmax ||
00595 wmap->list[ord].lambda_min >= wmap->list[ord].lambda_max ||
00596 wmap->list[ord].ymin >= wmap->list[ord].ymax ) {
00597 return 0. ;
00598 }
00599
00600 check( tcheb_x = xsh_tools_tchebitchev_transform( x, wmap->list[ord].xmin,
00601 wmap->list[ord].xmax));
00602 check( tcheb_y = xsh_tools_tchebitchev_transform( y,
00603 wmap->list[ord].ymin, wmap->list[ord].ymax));
00604 check( tcheb_coef_x = xsh_tools_tchebitchev_poly_eval(
00605 degx, tcheb_x));
00606 check( tcheb_coef_y = xsh_tools_tchebitchev_poly_eval(
00607 degy, tcheb_y));
00608
00609 for ( i=0; i < degx+1; i++) {
00610 for( j=0; j < degy+1; j++) {
00611 double vx, vy;
00612 double coef;
00613
00614 check( vx = cpl_vector_get( tcheb_coef_x, i));
00615 check( vy = cpl_vector_get( tcheb_coef_y, j));
00616 pows[0] = j;
00617 pows[1] = i;
00618 check( coef = cpl_polynomial_get_coeff( wmap->list[ord].tcheb_pol_lambda,
00619 pows));
00620 res += coef*vx*vy ;
00621 }
00622 }
00623 check( lambda = xsh_tools_tchebitchev_reverse_transform( res,
00624 wmap->list[ord].lambda_min, wmap->list[ord].lambda_max));
00625
00626 cleanup:
00627 xsh_free_vector( &tcheb_coef_x);
00628 xsh_free_vector( &tcheb_coef_y);
00629 return lambda ;
00630 }
00631
00644 cpl_frame *
00645 xsh_wavemap_list_save( xsh_wavemap_list * wmap,
00646 cpl_frame * order_frame,
00647 xsh_pre * pre,
00648 xsh_instrument * instr,
00649 const char * prefix)
00650 {
00651 cpl_frame * result = NULL ;
00652 int order, x, y ;
00653 int start_y, end_y ;
00654 int xmin, xmax ;
00655 int nx, ny ;
00656 double lambda ;
00657 cpl_image * wimg = NULL ;
00658 double * dimg = NULL ;
00659 cpl_propertylist * wheader = NULL ;
00660 char * final_name = NULL;
00661 xsh_order_list * order_list = NULL ;
00662
00663 XSH_ASSURE_NOT_NULL( wmap ) ;
00664 XSH_ASSURE_NOT_NULL( order_frame ) ;
00665 XSH_ASSURE_NOT_NULL( pre ) ;
00666 XSH_ASSURE_NOT_NULL( prefix ) ;
00667 XSH_ASSURE_NOT_NULL( instr ) ;
00668
00669 final_name = xsh_stringcat_any( prefix, ".fits", NULL ) ;
00670 xsh_msg( "Entering xsh_wavemap_save, file \"%s\"", final_name ) ;
00671
00672
00673
00674
00675
00676
00677
00678
00679 check(order_list=xsh_order_list_load(order_frame, instr )) ;
00680
00681
00682
00683
00684
00685 nx = pre->nx ;
00686 ny = pre->ny ;
00687 xsh_msg( "Image size:%d,%d", nx, ny ) ;
00688
00689 check( wimg = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE ) ) ;
00690 check( dimg = cpl_image_get_data_double( wimg ) ) ;
00691
00692 for( order = 0 ; order < wmap->size ; order++ ) {
00693 int count = 0 ;
00694
00695
00696 if ( wmap->list[order].tcheb_pol_lambda == NULL ) {
00697 xsh_msg( "Order %d: NULL Polynome", order ) ;
00698 continue ;
00699 }
00700 if ( wmap->list[order].sky_size <= 0 ) {
00701 xsh_msg( "NOT ENOUGH DATA FOR ORDER %d",
00702 order_list->list[order].absorder ) ;
00703 continue ;
00704 }
00705
00706 start_y = xsh_order_list_get_starty( order_list, order ) ;
00707 end_y = xsh_order_list_get_endy( order_list, order ) ;
00708 xsh_msg_dbg_low( " Order %d, Ymin: %d, Ymax: %d",
00709 order, start_y, end_y ) ;
00710
00711 for ( y = start_y ; y < end_y ; y++ ) {
00712 double dxmin, dxmax ;
00713
00714
00715 check( dxmin = cpl_polynomial_eval_1d(order_list->list[order].edglopoly,
00716 (double)y, NULL ) ) ;
00717 check( dxmax = cpl_polynomial_eval_1d(order_list->list[order].edguppoly,
00718 (double)y, NULL ) ) ;
00719 xmin = floor(dxmin) ;
00720 xmax = ceil(dxmax) ;
00721 for (x = xmin ; x<xmax ; x++ ) {
00722 count++ ;
00723
00724 lambda = xsh_wavemap_list_eval_lambda( wmap, (double)x,
00725 (double)y, order);
00726 if ( cpl_error_get_code() != CPL_ERROR_NONE){
00727 lambda =0.0;
00728 xsh_error_reset();
00729 }
00730
00731 *(dimg + x + (y*nx)) = (float)lambda ;
00732 }
00733 }
00734 xsh_msg_dbg_high( " %d points to calculate", count ) ;
00735 }
00736
00737
00738 check( wheader = cpl_propertylist_duplicate( pre->data_header ) ) ;
00739
00740 check( cpl_image_save( wimg, final_name, CPL_BPP_IEEE_FLOAT,
00741 wheader, CPL_IO_DEFAULT ) ) ;
00742 check(result=xsh_frame_product(final_name,
00743 "TAG",
00744 CPL_FRAME_TYPE_IMAGE,
00745 CPL_FRAME_GROUP_PRODUCT,
00746 CPL_FRAME_LEVEL_TEMPORARY));
00747
00748 cleanup:
00749 xsh_order_list_free( &order_list);
00750 XSH_FREE( final_name ) ;
00751 xsh_free_image ( &wimg ) ;
00752 xsh_free_propertylist( &wheader);
00753
00754 return result ;
00755 }
00756
00769 cpl_frame *
00770 xsh_wavemap_list_save2( xsh_wavemap_list * wmap,
00771 xsh_order_list * order_list,
00772 xsh_pre * pre,
00773 xsh_instrument * instr,
00774 const char * prefix)
00775 {
00776 cpl_frame * result = NULL ;
00777 int order, x, y ;
00778 int start_y, end_y ;
00779 int xmin, xmax ;
00780 int nx, ny ;
00781 double lambda ;
00782 cpl_image * wimg = NULL ;
00783 double * dimg = NULL ;
00784 cpl_propertylist * wheader = NULL ;
00785 char * final_name = NULL;
00786
00787 XSH_ASSURE_NOT_NULL( wmap ) ;
00788 XSH_ASSURE_NOT_NULL( order_list ) ;
00789 XSH_ASSURE_NOT_NULL( pre ) ;
00790 XSH_ASSURE_NOT_NULL( prefix ) ;
00791 XSH_ASSURE_NOT_NULL( instr ) ;
00792
00793 final_name = xsh_stringcat_any( prefix, ".fits", NULL ) ;
00794 xsh_msg( "Entering xsh_wavemap_save, file \"%s\"", final_name ) ;
00795
00796
00797
00798
00799
00800
00801
00802
00803
00804
00805
00806
00807
00808
00809 nx = pre->nx ;
00810 ny = pre->ny ;
00811 xsh_msg( "Image size:%d,%d", nx, ny ) ;
00812
00813 check( wimg = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE ) ) ;
00814 check( dimg = cpl_image_get_data_double( wimg ) ) ;
00815
00816 for( order = 0 ; order < wmap->size ; order++ ) {
00817 int count = 0 ;
00818
00819
00820 if ( wmap->list[order].tcheb_pol_lambda == NULL ) {
00821 xsh_msg( "Order %d: NULL Polynome", order ) ;
00822 continue ;
00823 }
00824 if ( wmap->list[order].sky_size <= 0 ) {
00825 xsh_msg( "NOT ENOUGH DATA FOR ORDER %d",
00826 order_list->list[order].absorder ) ;
00827 continue ;
00828 }
00829
00830 start_y = xsh_order_list_get_starty( order_list, order ) ;
00831 end_y = xsh_order_list_get_endy( order_list, order ) ;
00832 xsh_msg_dbg_low( " Order %d, Ymin: %d, Ymax: %d",
00833 order, start_y, end_y ) ;
00834
00835 for ( y = start_y ; y < end_y ; y++ ) {
00836 double dxmin, dxmax ;
00837
00838
00839 check( dxmin = cpl_polynomial_eval_1d(order_list->list[order].edglopoly,
00840 (double)y, NULL ) ) ;
00841 check( dxmax = cpl_polynomial_eval_1d(order_list->list[order].edguppoly,
00842 (double)y, NULL ) ) ;
00843 xmin = floor(dxmin) ;
00844 xmax = ceil(dxmax) ;
00845 for (x = xmin ; x<xmax ; x++ ) {
00846 count++ ;
00847
00848 lambda = xsh_wavemap_list_eval_lambda( wmap, (double)x,
00849 (double)y, order);
00850 if ( cpl_error_get_code() != CPL_ERROR_NONE){
00851 lambda =0.0;
00852 xsh_error_reset();
00853 }
00854
00855 *(dimg + x + (y*nx)) = (float)lambda ;
00856 }
00857 }
00858 xsh_msg_dbg_high( " %d points to calculate", count ) ;
00859 }
00860
00861
00862 check( wheader = cpl_propertylist_duplicate( pre->data_header ) ) ;
00863
00864 check( cpl_image_save( wimg, final_name, CPL_BPP_IEEE_FLOAT,
00865 wheader, CPL_IO_DEFAULT ) ) ;
00866 check(result=xsh_frame_product(final_name,
00867 "TAG",
00868 CPL_FRAME_TYPE_IMAGE,
00869 CPL_FRAME_GROUP_PRODUCT,
00870 CPL_FRAME_LEVEL_TEMPORARY));
00871
00872 cleanup:
00873 xsh_order_list_free( &order_list);
00874 XSH_FREE( final_name ) ;
00875 xsh_free_image ( &wimg ) ;
00876 xsh_free_propertylist( &wheader);
00877
00878 return result ;
00879 }
00880
00881
00882
00896
00897 cpl_frame*
00898 xsh_wavemap_list_save_poly( xsh_wavemap_list* wmap,
00899 cpl_frame *order_frame,
00900 xsh_pre *pre,
00901 xsh_instrument *instr,
00902 const char *prefix,
00903 cpl_frame **dispersol_frame,
00904 cpl_frame **slitmap_frame)
00905 {
00906 cpl_frame *result = NULL;
00907 const char *slit_tag=NULL;
00908 int order;
00909 xsh_dispersol_list *dispersol_list = NULL;
00910 const char* disp_tag=NULL;
00911
00912 XSH_ASSURE_NOT_NULL( wmap);
00913 XSH_ASSURE_NOT_NULL( order_frame);
00914 XSH_ASSURE_NOT_NULL( prefix);
00915 XSH_ASSURE_NOT_NULL( dispersol_frame);
00916 XSH_ASSURE_NOT_NULL( instr);
00917
00918
00919 check( dispersol_list = xsh_dispersol_list_new( wmap->size,
00920 wmap->degx, wmap->degy, instr));
00921
00922 for( order = 0 ; order < wmap->size ; order++ ) {
00923
00924 check( xsh_dispersol_list_add( dispersol_list, order,
00925 wmap->list[order].order, wmap->list[order].pol_lambda,
00926 wmap->list[order].pol_slit));
00927 wmap->list[order].pol_lambda = NULL;
00928 wmap->list[order].pol_slit = NULL;
00929 }
00930
00931 if (pre != NULL){
00932
00933 check( result = xsh_dispersol_list_to_wavemap( dispersol_list,
00934 order_frame, pre,
00935 instr,prefix));
00936
00937 slit_tag = XSH_GET_TAG_FROM_ARM( XSH_SLIT_MAP_POLY, instr);
00938
00939
00940 check (*slitmap_frame = xsh_dispersol_list_to_slitmap( dispersol_list,
00941 order_frame, pre,
00942 instr, slit_tag));
00943 }
00944
00945 if(strstr(cpl_frame_get_tag(order_frame),"AFC") != NULL) {
00946 disp_tag = XSH_GET_TAG_FROM_ARM( XSH_DISP_TAB_AFC,instr);
00947 } else {
00948 disp_tag = XSH_GET_TAG_FROM_ARM( XSH_DISP_TAB,instr);
00949 }
00950 check(*dispersol_frame = xsh_dispersol_list_save( dispersol_list,disp_tag));
00951
00952 cleanup:
00953 if (cpl_error_get_code() != CPL_ERROR_NONE){
00954 xsh_free_frame( &result);
00955 }
00956 xsh_dispersol_list_free( &dispersol_list);
00957
00958
00959 return result ;
00960 }
00961