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
00031
00036
00037
00038
00042
00043
00044
00045
00046 #include <math.h>
00047
00048 #include <xsh_data_dispersol.h>
00049 #include <xsh_utils.h>
00050 #include <xsh_error.h>
00051 #include <xsh_msg.h>
00052 #include <xsh_pfits.h>
00053 #include <xsh_dfs.h>
00054 #include <cpl.h>
00055 #include <xsh_utils_table.h>
00056 #include <xsh_data_spectralformat.h>
00057 #include <xsh_data_order.h>
00058
00059
00060
00061
00062
00063
00064
00076
00077 xsh_dispersol_list* xsh_dispersol_list_new( int size, int degx, int degy,
00078 xsh_instrument *instrument)
00079 {
00080 xsh_dispersol_list *result = NULL;
00081
00082
00083 XSH_ASSURE_NOT_ILLEGAL( size > 0);
00084 XSH_ASSURE_NOT_NULL( instrument);
00085
00086 XSH_CALLOC( result, xsh_dispersol_list, 1);
00087 result->size = size;
00088 result->degx = degx;
00089 result->degy = degy;
00090 check( result->binx = xsh_instrument_get_binx( instrument));
00091 check( result->biny = xsh_instrument_get_biny( instrument));
00092 XSH_CALLOC( result->list, xsh_dispersol, result->size);
00093 XSH_NEW_PROPERTYLIST( result->header);
00094
00095 cleanup:
00096 if ( cpl_error_get_code() != CPL_ERROR_NONE) {
00097 xsh_dispersol_list_free( &result);
00098 }
00099 return result;
00100 }
00101
00102
00103
00104
00113
00114 xsh_dispersol_list* xsh_dispersol_list_load( cpl_frame *frame,
00115 xsh_instrument *instrument)
00116 {
00117 xsh_dispersol_list *result = NULL;
00118 const char* tablename = NULL;
00119 cpl_table* table = NULL;
00120 int rows;
00121
00122 int degx, degy;
00123 char coefname[20];
00124 int i;
00125
00126
00127 XSH_ASSURE_NOT_ILLEGAL( frame);
00128
00129
00130 check( tablename = cpl_frame_get_filename(frame));
00131
00132 XSH_TABLE_LOAD( table, tablename);
00133 check( rows = cpl_table_get_nrow( table));
00134
00135
00136
00137 check( xsh_get_table_value( table, XSH_DISPERSOL_TABLE_COLNAME_DEGX,
00138 CPL_TYPE_INT, 0, °x));
00139
00140 check( xsh_get_table_value( table, XSH_DISPERSOL_TABLE_COLNAME_DEGY,
00141 CPL_TYPE_INT, 0, °y));
00142
00143 check( result = xsh_dispersol_list_new( rows/2, degx, degy, instrument));
00144
00145 for(i=0; i< rows/2; i++){
00146 cpl_polynomial *lambda_poly = NULL;
00147 cpl_polynomial *slit_poly = NULL;
00148 cpl_size pows[2];
00149 int k,l, absorder = 0;
00150 int j = 2*i;
00151
00152 check( xsh_get_table_value( table, XSH_DISPERSOL_TABLE_COLNAME_ORDER,
00153 CPL_TYPE_INT, j, &absorder));
00154 check( lambda_poly = cpl_polynomial_new( 2));
00155 check( slit_poly = cpl_polynomial_new( 2));
00156 for( k=0; k<= degx; k++){
00157 for( l=0; l<= degy; l++){
00158 double coef_lambda = 0.0, coef_slit = 0.0;
00159 sprintf(coefname,"C%d%d",k,l);
00160 check( xsh_get_table_value( table, coefname,CPL_TYPE_DOUBLE, j,
00161 &coef_lambda));
00162 check( xsh_get_table_value( table, coefname,CPL_TYPE_DOUBLE, j+1,
00163 &coef_slit));
00164 pows[0] =k;
00165 pows[1] = l;
00166 check( cpl_polynomial_set_coeff( lambda_poly, pows, coef_lambda));
00167 check( cpl_polynomial_set_coeff( slit_poly, pows, coef_slit));
00168 }
00169 }
00170 check( xsh_dispersol_list_add( result, i, absorder, lambda_poly, slit_poly));
00171 }
00172
00173 cleanup:
00174 xsh_free_table( &table);
00175 if ( cpl_error_get_code() != CPL_ERROR_NONE) {
00176 xsh_dispersol_list_free( &result);
00177 }
00178 return result;
00179 }
00180
00181
00182
00197
00198 void xsh_dispersol_list_add( xsh_dispersol_list *list, int idx,
00199 int absorder, cpl_polynomial *lambda_poly,
00200 cpl_polynomial *slit_poly)
00201 {
00202
00203 XSH_ASSURE_NOT_NULL( list);
00204 XSH_ASSURE_NOT_NULL( lambda_poly);
00205 XSH_ASSURE_NOT_NULL( slit_poly);
00206 XSH_ASSURE_NOT_ILLEGAL( idx >=0 && idx < list->size);
00207
00208 list->list[idx].absorder = absorder;
00209 list->list[idx].lambda_poly = lambda_poly;
00210 list->list[idx].slit_poly = slit_poly;
00211
00212 cleanup:
00213 return;
00214 }
00215
00216
00217
00224
00225 void xsh_dispersol_list_free( xsh_dispersol_list **list)
00226 {
00227 int i = 0;
00228
00229 if ( list && *list){
00230
00231 for( i = 0; i < (*list)->size; i++){
00232 xsh_free_polynomial( &(*list)->list[i].lambda_poly);
00233 xsh_free_polynomial( &(*list)->list[i].slit_poly);
00234 }
00235 if ((*list)->list){
00236 cpl_free( (*list)->list);
00237 }
00238 xsh_free_propertylist( &((*list)->header));
00239 cpl_free( *list);
00240 *list = NULL;
00241 }
00242 }
00243
00244
00245
00246
00261
00262 cpl_frame* xsh_dispersol_list_to_wavemap( xsh_dispersol_list *list,
00263 cpl_frame *order_frame,
00264 xsh_pre *pre,
00265 xsh_instrument *instr,
00266 const char* tag)
00267
00268 {
00269 cpl_frame *result = NULL;
00270 cpl_image *image = NULL;
00271 xsh_order_list *order_list = NULL;
00272 int nx, ny;
00273 cpl_vector* pos = NULL;
00274 double *data = NULL;
00275 int i, y;
00276 char filename[256];
00277 cpl_propertylist *wavemap_header = NULL;
00278
00279
00280 XSH_ASSURE_NOT_NULL( list);
00281 XSH_ASSURE_NOT_NULL( order_frame);
00282 XSH_ASSURE_NOT_NULL( pre);
00283 XSH_ASSURE_NOT_NULL( instr);
00284
00285
00286 check( order_list = xsh_order_list_load( order_frame, instr));
00287 nx = pre->nx;
00288 ny = pre->ny;
00289 pos = cpl_vector_new(2);
00290
00291
00292 check( image = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE));
00293 check( data = cpl_image_get_data_double( image));
00294 check( wavemap_header = cpl_propertylist_new());
00295
00296
00297 for( i=0; i < list->size; i++) {
00298 int starty, endy;
00299 double lambda_min, lambda_max;
00300 double lambda_starty_lo, lambda_starty_up, lambda_endy_lo, lambda_endy_up;
00301 double dxmin, dxmax ;
00302 int xmin, xmax;
00303 int absorder;
00304
00305
00306 absorder = order_list->list[i].absorder;
00307 check( starty = xsh_order_list_get_starty( order_list, i));
00308
00309 check( endy = xsh_order_list_get_endy( order_list, i));
00310
00311
00312 for ( y=starty ; y<=endy && y<ny; y++){
00313 int x;
00314
00315
00316 check( dxmin = xsh_order_list_eval( order_list,
00317 order_list->list[i].edglopoly,
00318 (double)y ) ) ;
00319 check( dxmax = xsh_order_list_eval( order_list,
00320 order_list->list[i].edguppoly,
00321 (double)y ) ) ;
00322 xmin = ceil( dxmin);
00323 xmax = floor( dxmax);
00324 for ( x=xmin ; x<=xmax && x<nx; x++){
00325 double lambda;
00326
00327 cpl_vector_set(pos, 0, (double)x);
00328 cpl_vector_set(pos, 1, (double)y);
00329
00330
00331 check( lambda = xsh_dispersol_list_eval(
00332 list, list->list[i].lambda_poly, pos));
00333
00334 data[x-1+(y-1)*nx] = (float)lambda;
00335
00336 }
00337
00338 }
00339
00340
00341 check( dxmin = xsh_order_list_eval( order_list,
00342 order_list->list[i].edglopoly,
00343 (double)starty));
00344 check( dxmax = xsh_order_list_eval( order_list,
00345 order_list->list[i].edguppoly,
00346 (double)starty));
00347 xmin = ceil( dxmin);
00348 xmax = floor( dxmax);
00349
00350 cpl_vector_set(pos, 0, (double)xmin);
00351 cpl_vector_set(pos, 1, (double)starty);
00352 check( lambda_starty_lo = xsh_dispersol_list_eval(
00353 list, list->list[i].lambda_poly, pos));
00354
00355 cpl_vector_set(pos, 0, (double)xmax);
00356 check( lambda_starty_up = xsh_dispersol_list_eval(
00357 list, list->list[i].lambda_poly, pos));
00358
00359 check( dxmin = xsh_order_list_eval( order_list,
00360 order_list->list[i].edglopoly,
00361 (double)endy));
00362 check( dxmax = xsh_order_list_eval( order_list,
00363 order_list->list[i].edguppoly,
00364 (double)endy));
00365 xmin = ceil( dxmin);
00366 xmax = floor( dxmax);
00367
00368 cpl_vector_set(pos, 0, (double)xmin);
00369 cpl_vector_set(pos, 1, (double)endy);
00370 check( lambda_endy_lo = xsh_dispersol_list_eval(
00371 list, list->list[i].lambda_poly, pos));
00372
00373 cpl_vector_set(pos, 0, (double)xmax);
00374 check( lambda_endy_up = xsh_dispersol_list_eval(
00375 list, list->list[i].lambda_poly, pos));
00376
00377 if ( lambda_starty_lo < lambda_endy_lo
00378 && lambda_starty_up < lambda_endy_up){
00379 if ( lambda_starty_lo > lambda_starty_up){
00380 lambda_min = lambda_starty_lo;
00381 }
00382 else{
00383 lambda_min = lambda_starty_up;
00384 }
00385 if ( lambda_endy_lo > lambda_endy_up){
00386 lambda_max = lambda_endy_up;
00387 }
00388 else{
00389 lambda_max = lambda_endy_lo;
00390 }
00391 }
00392 else if ( lambda_starty_lo > lambda_endy_lo
00393 && lambda_starty_up > lambda_endy_up){
00394 if ( lambda_starty_lo > lambda_starty_up){
00395 lambda_max = lambda_starty_up;
00396 }
00397 else{
00398 lambda_max = lambda_starty_lo;
00399 }
00400 if ( lambda_endy_lo > lambda_endy_up){
00401 lambda_min = lambda_endy_lo;
00402 }
00403 else{
00404 lambda_min = lambda_endy_up;
00405 }
00406 }
00407 else{
00408 xsh_msg_warning(
00409 "wavelength variation differs between upper(%f %f) and lower edges( %f %f)",
00410 lambda_starty_up, lambda_endy_up, lambda_endy_up, lambda_endy_lo);
00411 lambda_min = 0.0;
00412 lambda_max = 0.0;
00413 }
00414
00415
00416 check( xsh_pfits_set_wavemap_order_lambda_min(
00417 wavemap_header, absorder, lambda_min));
00418 check( xsh_pfits_set_wavemap_order_lambda_max(
00419 wavemap_header, absorder, lambda_max));
00420 }
00421
00422 sprintf(filename,"%s.fits",tag);
00423 check( xsh_pfits_set_pcatg( wavemap_header, tag));
00424
00425 check( cpl_image_save( image, filename, CPL_BPP_IEEE_FLOAT,
00426 wavemap_header, CPL_IO_DEFAULT));
00427
00428 check( result = xsh_frame_product(filename, tag,
00429 CPL_FRAME_TYPE_IMAGE, CPL_FRAME_GROUP_PRODUCT,
00430 CPL_FRAME_LEVEL_TEMPORARY));
00431
00432 cleanup:
00433 if (cpl_error_get_code() != CPL_ERROR_NONE){
00434 xsh_free_frame( &result);
00435 }
00436 xsh_free_image ( &image);
00437 xsh_free_propertylist( &wavemap_header);
00438 xsh_order_list_free( &order_list);
00439 xsh_free_vector( &pos);
00440
00441
00442 return result;
00443 }
00444
00445
00446
00447
00462
00463 cpl_frame* xsh_dispersol_list_to_slitmap( xsh_dispersol_list *list,
00464 cpl_frame *order_frame,
00465 xsh_pre *pre,
00466 xsh_instrument *instr,
00467 const char* tag)
00468 {
00469 cpl_frame *result = NULL;
00470 cpl_image *image = NULL;
00471 cpl_image *ifu_map = NULL;
00472 xsh_order_list *order_list = NULL;
00473 int nx, ny;
00474 cpl_vector* pos = NULL;
00475 double *data = NULL;
00476 int *imap = NULL;
00477 double crpix1=1.;
00478 double crpix2=1.;
00479 double crval1=1.;
00480 double crval2=1.;
00481 double cdelt1=1.;
00482 double cdelt2=1.;
00483
00484 int i, y;
00485 char filename[256];
00486 cpl_vector *med_low_vect = NULL, *med_up_vect = NULL, *med_cen_vect = NULL;
00487 cpl_vector *med_sliclo_vect = NULL, *med_slicup_vect = NULL;
00488 int nb_orders;
00489 double slit_up, slit_low, slit_cen;
00490 cpl_propertylist *slitmap_header = NULL;
00491
00492
00493 XSH_ASSURE_NOT_NULL( list);
00494 XSH_ASSURE_NOT_NULL( order_frame);
00495 XSH_ASSURE_NOT_NULL( pre);
00496 XSH_ASSURE_NOT_NULL( instr);
00497
00498
00499 check( order_list = xsh_order_list_load( order_frame, instr));
00500 nx = pre->nx;
00501 ny = pre->ny;
00502
00503
00504 XSH_ASSURE_NOT_ILLEGAL( order_list->bin_x == pre->binx);
00505 XSH_ASSURE_NOT_ILLEGAL( order_list->bin_x == list->binx);
00506
00507 pos = cpl_vector_new(2);
00508
00509 check( image = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE));
00510 check( data = cpl_image_get_data_double( image));
00511 check( ifu_map = cpl_image_new( nx, ny, CPL_TYPE_INT));
00512 check( imap = cpl_image_get_data_int( ifu_map));
00513 check( slitmap_header = cpl_propertylist_new());
00514
00515 nb_orders = list->size;
00516
00517 med_low_vect = cpl_vector_new( nb_orders);
00518 med_cen_vect = cpl_vector_new( nb_orders);
00519 med_up_vect = cpl_vector_new( nb_orders);
00520 if ( xsh_instrument_get_mode( instr) == XSH_MODE_IFU){
00521 med_sliclo_vect = cpl_vector_new( nb_orders);
00522 med_slicup_vect = cpl_vector_new( nb_orders);
00523 }
00524
00525
00526 for( i=0; i < nb_orders; i++) {
00527 int starty, endy, vect_size;
00528 int absorder;
00529 cpl_vector *low_vect = NULL, *up_vect = NULL, *cen_vect = NULL;
00530 cpl_vector *sliclo_vect = NULL, *slicup_vect = NULL;
00531
00532 check( starty = xsh_order_list_get_starty( order_list, i));
00533 check( endy = xsh_order_list_get_endy( order_list, i));
00534 absorder = order_list->list[i].absorder;
00535
00536
00537 int ymax = ( endy < ny ) ? endy:(ny-1);
00538
00539 vect_size = ymax-starty+1;
00540
00541
00542
00543
00544 if(vect_size<1) {
00545
00546 continue;
00547
00548 }
00549 low_vect = cpl_vector_new( vect_size);
00550 cen_vect = cpl_vector_new( vect_size);
00551 up_vect = cpl_vector_new( vect_size);
00552
00553 if ( xsh_instrument_get_mode( instr) == XSH_MODE_IFU){
00554 sliclo_vect = cpl_vector_new( vect_size);
00555 slicup_vect = cpl_vector_new( vect_size);
00556 }
00557
00558 for ( y=starty ; y<=ymax; y++){
00559
00560 double dxmin, dxmax, xcen;
00561 int xmin, xmax, x;
00562 double slit;
00563 double ypos;
00564
00565
00566 check( dxmin = xsh_order_list_eval( order_list,
00567 order_list->list[i].edglopoly,
00568 (double)y ) ) ;
00569 check( dxmax = xsh_order_list_eval( order_list,
00570 order_list->list[i].edguppoly,
00571 (double)y ) ) ;
00572 check( xcen = xsh_order_list_eval( order_list,
00573 order_list->list[i].cenpoly,
00574 (double)y ) ) ;
00575
00576 xmin = ceil( dxmin);
00577 xmax = floor( dxmax);
00578
00579 ypos = (double)y;
00580 cpl_vector_set(pos, 1, ypos);
00581
00582 for ( x=xmin ; x<=xmax && x<nx; x++){
00583 double xpos;
00584
00585 xpos = (double)x;
00586 cpl_vector_set(pos, 0, xpos);
00587 cpl_vector_set(pos, 1, ypos);
00588
00589 check( slit = xsh_dispersol_list_eval( list, list->list[i].slit_poly, pos));
00590
00591
00592 data[x-1+(y-1)*nx] = (float)slit;
00593 imap[x-1+(y-1)*nx] = 1;
00594 }
00595
00596 cpl_vector_set(pos, 0, xmin);
00597 cpl_vector_set(pos, 1, ypos);
00598 check( slit = xsh_dispersol_list_eval( list, list->list[i].slit_poly, pos));
00599 cpl_vector_set( low_vect, y-starty, slit);
00600
00601 cpl_vector_set(pos, 0, xmax);
00602 cpl_vector_set(pos, 1, ypos);
00603 check( slit = xsh_dispersol_list_eval( list, list->list[i].slit_poly, pos));
00604 cpl_vector_set( up_vect, y-starty, slit);
00605
00606 cpl_vector_set(pos, 0, xcen);
00607 cpl_vector_set(pos, 1, ypos);
00608 check( slit = xsh_dispersol_list_eval( list, list->list[i].slit_poly, pos));
00609 cpl_vector_set( cen_vect, y-starty, slit);
00610
00611 if ( xsh_instrument_get_mode( instr) == XSH_MODE_IFU){
00612 xmin = xsh_order_list_eval( order_list,
00613 order_list->list[i].sliclopoly,
00614 (double)y );
00615 xmax = xsh_order_list_eval( order_list,
00616 order_list->list[i].slicuppoly,
00617 (double)y );
00618 cpl_vector_set(pos, 0, xmin);
00619 cpl_vector_set(pos, 1, ypos);
00620 check( slit = xsh_dispersol_list_eval( list, list->list[i].slit_poly, pos));
00621 cpl_vector_set( sliclo_vect, y-starty, slit);
00622
00623 cpl_vector_set(pos, 0, xmax);
00624 cpl_vector_set(pos, 1, ypos);
00625 check( slit = xsh_dispersol_list_eval( list, list->list[i].slit_poly, pos));
00626 cpl_vector_set( slicup_vect, y-starty, slit);
00627 }
00628
00629 }
00630
00631 check( slit_up = cpl_vector_get_median( up_vect));
00632 check( slit_low = cpl_vector_get_median( low_vect));
00633 check( slit_cen = cpl_vector_get_median( cen_vect));
00634
00635
00636
00637 check( xsh_pfits_set_slitmap_order_cen( slitmap_header, absorder, slit_cen));
00638 check( cpl_vector_set( med_cen_vect, i, slit_cen));
00639
00640 if ( slit_up > slit_low){
00641 check( xsh_pfits_set_slitmap_order_edgup( slitmap_header, absorder, slit_up));
00642 check( xsh_pfits_set_slitmap_order_edglo( slitmap_header, absorder, slit_low));
00643 check( cpl_vector_set( med_low_vect, i, slit_low));
00644 check( cpl_vector_set( med_up_vect, i, slit_up));
00645 }
00646 else{
00647 check( xsh_pfits_set_slitmap_order_edgup( slitmap_header, absorder, slit_low));
00648 check( xsh_pfits_set_slitmap_order_edglo( slitmap_header, absorder, slit_up));
00649 check( cpl_vector_set( med_low_vect, i, slit_up));
00650 check( cpl_vector_set( med_up_vect, i, slit_low));
00651 }
00652
00653 if ( xsh_instrument_get_mode( instr) == XSH_MODE_IFU){
00654 check( slit_up = cpl_vector_get_median( slicup_vect));
00655 check( slit_low = cpl_vector_get_median( sliclo_vect));
00656 if ( slit_up > slit_low){
00657 check( xsh_pfits_set_slitmap_order_slicup( slitmap_header, absorder, slit_up));
00658 check( xsh_pfits_set_slitmap_order_sliclo( slitmap_header, absorder, slit_low));
00659 check( cpl_vector_set( med_sliclo_vect, i, slit_low));
00660 check( cpl_vector_set( med_slicup_vect, i, slit_up));
00661 }
00662 else{
00663 check( xsh_pfits_set_slitmap_order_slicup( slitmap_header, absorder, slit_low));
00664 check( xsh_pfits_set_slitmap_order_sliclo( slitmap_header, absorder, slit_up));
00665 check( cpl_vector_set( med_sliclo_vect, i, slit_up));
00666 check( cpl_vector_set( med_slicup_vect, i, slit_low));
00667 }
00668 xsh_free_vector( &slicup_vect);
00669 xsh_free_vector( &sliclo_vect);
00670 }
00671 xsh_free_vector( &up_vect);
00672 xsh_free_vector( &low_vect);
00673 xsh_free_vector( &cen_vect);
00674 }
00675
00676 check( slit_cen = cpl_vector_get_median( med_cen_vect));
00677 check( slit_up = cpl_vector_get_median( med_up_vect));
00678 check( slit_low = cpl_vector_get_median( med_low_vect));
00679
00680 check( xsh_pfits_set_slitmap_median_cen( slitmap_header, slit_cen));
00681 check( xsh_pfits_set_slitmap_median_edgup( slitmap_header, slit_up));
00682 check( xsh_pfits_set_slitmap_median_edglo( slitmap_header, slit_low));
00683
00684 if ( xsh_instrument_get_mode( instr) == XSH_MODE_IFU){
00685 check( slit_up = cpl_vector_get_median( med_slicup_vect));
00686 check( slit_low = cpl_vector_get_median( med_sliclo_vect));
00687 check( xsh_pfits_set_slitmap_median_slicup( slitmap_header, slit_up));
00688 check( xsh_pfits_set_slitmap_median_sliclo( slitmap_header, slit_low));
00689 }
00690
00691 sprintf(filename,"%s.fits",tag);
00692 xsh_pfits_set_pcatg( slitmap_header, tag);
00693 cdelt1=xsh_instrument_get_binx(instr);
00694 cdelt2=xsh_instrument_get_biny(instr);
00695 xsh_pfits_set_wcs(slitmap_header,crpix1,crval1,cdelt1,crpix2,crval2,cdelt2);
00696 check( cpl_image_save( image, filename, CPL_BPP_IEEE_FLOAT,
00697 slitmap_header, CPL_IO_DEFAULT));
00698
00699 if ( xsh_instrument_get_mode( instr) == XSH_MODE_IFU){
00700 check( cpl_image_save( ifu_map, filename, CPL_BPP_32_SIGNED,
00701 NULL, CPL_IO_EXTEND));
00702 }
00703 check( result = xsh_frame_product(filename,tag,
00704 CPL_FRAME_TYPE_IMAGE, CPL_FRAME_GROUP_PRODUCT,
00705 CPL_FRAME_LEVEL_TEMPORARY));
00706
00707 cleanup:
00708 if (cpl_error_get_code() != CPL_ERROR_NONE){
00709 xsh_free_frame( &result);
00710 }
00711
00712 xsh_free_image ( &image);
00713 xsh_free_image ( &ifu_map);
00714 xsh_free_propertylist( &slitmap_header);
00715 xsh_order_list_free( &order_list);
00716 xsh_free_vector( &pos);
00717 xsh_free_vector( &med_up_vect);
00718 xsh_free_vector( &med_low_vect);
00719 xsh_free_vector( &med_cen_vect);
00720 xsh_free_vector( &med_slicup_vect);
00721 xsh_free_vector( &med_sliclo_vect);
00722
00723 return result;
00724 }
00725
00726
00727
00740
00741 double
00742 xsh_dispersol_list_eval( xsh_dispersol_list *list,
00743 cpl_polynomial *poly,
00744 cpl_vector *pos)
00745 {
00746 double result=0;
00747 double y_bin, y_no_bin;
00748 double x_bin, x_no_bin;
00749
00750 XSH_ASSURE_NOT_NULL( list);
00751 XSH_ASSURE_NOT_NULL( pos);
00752 XSH_ASSURE_NOT_NULL( poly);
00753
00754 check( x_bin = cpl_vector_get(pos, 0));
00755 check( y_bin = cpl_vector_get(pos, 1));
00756
00757 x_no_bin = convert_bin_to_data( x_bin, list->binx);
00758 y_no_bin = convert_bin_to_data( y_bin, list->biny);
00759
00760 check( cpl_vector_set(pos, 0, x_no_bin));
00761 check( cpl_vector_set(pos, 1, y_no_bin));
00762
00763 check( result = cpl_polynomial_eval( poly, pos));
00764
00765 cleanup:
00766 return result;
00767 }
00768
00769
00770
00778
00779 cpl_frame*
00780 xsh_dispersol_list_save( xsh_dispersol_list *list, const char* tag)
00781 {
00782 cpl_frame *result = NULL;
00783 cpl_table *table = NULL;
00784 char coefname[20];
00785 int i,j,k;
00786 int nbcoefs =0;
00787 cpl_size pows[2];
00788 char filename[256];
00789
00790
00791 XSH_ASSURE_NOT_NULL( list);
00792
00793 nbcoefs = (list->degx+1)*(list->degy+1);
00794
00795 check( table = cpl_table_new( XSH_DISPERSOL_TABLE_NBCOL+nbcoefs));
00796 check(cpl_table_set_size(table, list->size*2));
00797
00798 check(
00799 cpl_table_new_column( table, XSH_DISPERSOL_TABLE_COLNAME_AXIS,
00800 CPL_TYPE_STRING));
00801 check(
00802 cpl_table_new_column( table, XSH_DISPERSOL_TABLE_COLNAME_ORDER,
00803 CPL_TYPE_INT));
00804 check(
00805 cpl_table_new_column( table, XSH_DISPERSOL_TABLE_COLNAME_DEGX,
00806 CPL_TYPE_INT));
00807 check(
00808 cpl_table_new_column( table, XSH_DISPERSOL_TABLE_COLNAME_DEGY,
00809 CPL_TYPE_INT));
00810
00811 for(i=0; i<= list->degx; i++){
00812 for(j = 0; j<= list->degy; j++){
00813 sprintf(coefname,"C%d%d",i,j);
00814 check( cpl_table_new_column( table, coefname, CPL_TYPE_DOUBLE));
00815 }
00816 }
00817
00818 for(i=0; i< list->size; i++){
00819 int i2 = i*2;
00820 check(cpl_table_set_string(table,XSH_DISPERSOL_TABLE_COLNAME_AXIS,
00821 i2, XSH_DISPERSOL_AXIS_LAMBDA));
00822 check(cpl_table_set_string(table,XSH_DISPERSOL_TABLE_COLNAME_AXIS,
00823 i2+1, XSH_DISPERSOL_AXIS_SLIT));
00824 check(cpl_table_set_int(table,XSH_DISPERSOL_TABLE_COLNAME_ORDER,
00825 i2, list->list[i].absorder));
00826 check(cpl_table_set_int(table,XSH_DISPERSOL_TABLE_COLNAME_ORDER,
00827 i2+1, list->list[i].absorder));
00828 check(cpl_table_set_int(table,XSH_DISPERSOL_TABLE_COLNAME_DEGX,
00829 i2, list->degx));
00830 check(cpl_table_set_int(table,XSH_DISPERSOL_TABLE_COLNAME_DEGX,
00831 i2+1, list->degx));
00832 check(cpl_table_set_int(table,XSH_DISPERSOL_TABLE_COLNAME_DEGY,
00833 i2, list->degy));
00834 check(cpl_table_set_int(table,XSH_DISPERSOL_TABLE_COLNAME_DEGY,
00835 i2+1, list->degy));
00836 for( j=0; j<= list->degx; j++){
00837 for( k=0; k<= list->degy; k++){
00838 double coef_lambda = 0.0, coef_slit = 0.0;
00839 sprintf(coefname,"C%d%d",j,k);
00840 pows[0] = j;
00841 pows[1] = k;
00842 check(coef_lambda =
00843 cpl_polynomial_get_coeff( list->list[i].lambda_poly, pows));
00844 check(cpl_table_set_double( table, coefname, i2, coef_lambda));
00845 check(coef_slit =
00846 cpl_polynomial_get_coeff( list->list[i].slit_poly, pows));
00847 check(cpl_table_set_double( table, coefname, i2+1, coef_slit));
00848 }
00849 }
00850 }
00851
00852 sprintf(filename, "%s.fits",tag);
00853 check( xsh_pfits_set_pcatg(list->header,tag));
00854 check(cpl_table_save(table, list->header, NULL,
00855 filename, CPL_IO_DEFAULT));
00856
00857 check( result= xsh_frame_product( filename, tag,
00858 CPL_FRAME_TYPE_TABLE,
00859 CPL_FRAME_GROUP_PRODUCT,
00860 CPL_FRAME_LEVEL_TEMPORARY));
00861
00862 cleanup:
00863 xsh_free_table( &table);
00864 if ( cpl_error_get_code() != CPL_ERROR_NONE) {
00865 }
00866 return result;
00867 }
00868
00869