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
00041
00042
00043
00044 #include <xsh_data_wavesol.h>
00045 #include <xsh_data_spectralformat.h>
00046 #include <xsh_utils.h>
00047 #include <xsh_error.h>
00048 #include <xsh_msg.h>
00049 #include <xsh_pfits.h>
00050 #include <xsh_dfs.h>
00051 #include <math.h>
00052 #include <gsl/gsl_multifit.h>
00053 #include <cpl.h>
00054 #include <xsh_utils_table.h>
00055
00056
00057
00058
00059
00072
00073 xsh_wavesol* xsh_wavesol_create( cpl_frame* spectral_format_frame,
00074 xsh_detect_arclines_param* params, xsh_instrument *instrument)
00075 {
00076 xsh_wavesol* res = NULL;
00077 xsh_spectralformat_list *spec_list = NULL;
00078 int i;
00079 float sol_lambda_min=0, sol_lambda_max=0;
00080 int sol_absorder_min=0, sol_absorder_max=0;
00081
00082
00083 XSH_ASSURE_NOT_NULL( params);
00084 XSH_ASSURE_NOT_NULL( spectral_format_frame);
00085 XSH_ASSURE_NOT_NULL( instrument ) ;
00086
00087
00088 XSH_MALLOC(res, xsh_wavesol, 1);
00089 res->deg_slit = params->wavesol_deg_slit;
00090 res->deg_order = params->wavesol_deg_order;
00091 res->deg_lambda = params->wavesol_deg_lambda;
00092 res->nbcoefs = (res->deg_order+1)*(res->deg_slit+1)*(res->deg_lambda+1);
00093 xsh_msg_dbg_high("nbcoef %d deg_lambda %d deg_n %d deg_s %d",res->nbcoefs,
00094 res->deg_lambda, res->deg_order, res->deg_slit);
00095 res->polx = cpl_polynomial_new(3);
00096 res->poly = cpl_polynomial_new(3);
00097 res->dim = cpl_vector_new(3);
00098 res->header = cpl_propertylist_new();
00099
00100
00101 res->bin_x = xsh_instrument_get_binx( instrument ) ;
00102 res->bin_y = xsh_instrument_get_biny( instrument ) ;
00103
00104
00105 check( spec_list = xsh_spectralformat_list_load( spectral_format_frame,
00106 instrument));
00107 sol_lambda_min = spec_list->list[0].lambda_min_full;
00108 sol_lambda_max = spec_list->list[0].lambda_max_full;
00109 sol_absorder_min = spec_list->list[0].absorder;
00110 sol_absorder_max = spec_list->list[0].absorder;
00111
00112 for( i=1; i< spec_list->size; i++){
00113 int absorder;
00114 float lambda_min, lambda_max;
00115 absorder = spec_list->list[i].absorder;
00116 lambda_min = spec_list->list[i].lambda_min_full;
00117 lambda_max = spec_list->list[i].lambda_max_full;
00118 if ( absorder > sol_absorder_max){
00119 sol_absorder_max = absorder;
00120 }
00121 if ( absorder < sol_absorder_min){
00122 sol_absorder_min = absorder;
00123 }
00124 if ( lambda_min < sol_lambda_min){
00125 sol_lambda_min = lambda_min;
00126 }
00127 if ( lambda_max > sol_lambda_max){
00128 sol_lambda_max = lambda_max;
00129 }
00130 }
00131 xsh_msg_dbg_high("Order range %d-%d",sol_absorder_min, sol_absorder_max);
00132 xsh_msg_dbg_high("Lambda range %f-%f", sol_lambda_min, sol_lambda_max);
00133 res->min_lambda = sol_lambda_min;
00134 res->max_lambda = sol_lambda_max;
00135 res->min_order = sol_absorder_min;
00136 res->max_order = sol_absorder_max;
00137
00138 cleanup:
00139 xsh_spectralformat_list_free( &spec_list);
00140 return res;
00141 }
00142
00149
00150 xsh_wavesol * xsh_wavesol_duplicate( xsh_wavesol * org )
00151 {
00152 xsh_wavesol * res = NULL ;
00153
00154 XSH_MALLOC( res, xsh_wavesol, 1 ) ;
00155 res->bin_x = org->bin_x ;
00156 res->bin_y = org->bin_y ;
00157 res->deg_slit = org->deg_slit ;
00158 res->deg_order = org->deg_order ;
00159 res->deg_lambda = org->deg_lambda ;
00160 res->min_lambda = org->min_lambda ;
00161 res->max_lambda = org->max_lambda ;
00162 res->min_order = org->min_order ;
00163 res->max_order = org->max_order ;
00164 res->min_slit = org->min_slit ;
00165 res->max_slit = org->max_slit ;
00166 res->min_x = org->min_x ;
00167 res->max_x = org->max_x ;
00168 res->min_y = org->min_y ;
00169 res->max_y = org->max_y ;
00170 res->nbcoefs = org->nbcoefs ;
00171 res->polx = cpl_polynomial_duplicate( org->polx ) ;
00172 res->poly = cpl_polynomial_duplicate( org->poly ) ;
00173 res->dim = cpl_vector_duplicate( org->dim ) ;
00174 res->header = cpl_propertylist_duplicate( org->header ) ;
00175
00176 cleanup:
00177 return res ;
00178 }
00179
00186 void xsh_wavesol_add_poly( xsh_wavesol * to, xsh_wavesol * from )
00187 {
00188 int i, j, k ;
00189 int pows[3];
00190
00191 pows[0] = 0;
00192 i = 0 ;
00193
00194 for(j = 0; j<= from->deg_order; j++)
00195 for(k = 0; k <= from->deg_lambda; k++) {
00196 double vto = 0.0, vfrom = 0.0;
00197
00198 pows[0] = i ;
00199 pows[1] = j;
00200 pows[2] = k;
00201 xsh_msg_dbg_high( "Add_poly: %d %d %d", i, j, k ) ;
00202 check(vfrom = cpl_polynomial_get_coeff(from->poly, pows));
00203 check(vto = cpl_polynomial_get_coeff(to->poly, pows));
00204 vto += vfrom;
00205 check( cpl_polynomial_set_coeff( to->poly, pows, vto));
00206 }
00207
00208 cleanup:
00209 return ;
00210 }
00211
00212
00218
00219 void xsh_wavesol_free(xsh_wavesol** w)
00220 {
00221 if (w && *w) {
00222 xsh_free_polynomial(&(*w)->polx);
00223 xsh_free_polynomial(&(*w)->poly);
00224 xsh_free_vector(&(*w)->dim);
00225 xsh_free_propertylist(&(*w)->header);
00226 cpl_free(*w);
00227 *w = NULL;
00228 }
00229 }
00230
00231
00232
00239
00240 cpl_polynomial* xsh_wavesol_get_poly(xsh_wavesol* sol)
00241 {
00242 cpl_polynomial* res = NULL;
00243
00244 XSH_ASSURE_NOT_NULL(sol);
00245 res = sol->poly;
00246
00247 cleanup:
00248 return res;
00249 }
00250
00251
00252
00259
00260
00261 void xsh_wavesol_set_type(xsh_wavesol * wsol, enum wavesol_type type)
00262 {
00263 XSH_ASSURE_NOT_NULL( wsol);
00264 wsol->type = type;
00265
00266 cleanup:
00267 return;
00268 }
00269
00270
00271
00272
00279 enum wavesol_type xsh_wavesol_get_type(xsh_wavesol *wsol)
00280 {
00281 enum wavesol_type type = XSH_WAVESOL_GUESS;
00282 XSH_ASSURE_NOT_NULL( wsol);
00283 type = wsol->type;
00284
00285 cleanup:
00286 return type;
00287 }
00288
00289
00296
00297 cpl_polynomial* xsh_wavesol_get_polx(xsh_wavesol* sol)
00298 {
00299 cpl_polynomial* res = NULL;
00300
00301 XSH_ASSURE_NOT_NULL(sol);
00302 res = sol->polx;
00303
00304 cleanup:
00305 return res;
00306 }
00307
00308
00314
00315 cpl_propertylist* xsh_wavesol_get_header(xsh_wavesol* sol)
00316 {
00317 cpl_propertylist * res = NULL;
00318
00319 XSH_ASSURE_NOT_NULL(sol);
00320 res = sol->header;
00321 cleanup:
00322 return res;
00323 }
00324
00325
00335
00336 double xsh_wavesol_eval_polx(xsh_wavesol* sol, double lambda, double order,
00337 double slit)
00338 {
00339 double tcheb_slit = 0.0, tcheb_order, tcheb_lambda;
00340 double pos = 0.0, res = 0.0;
00341 int i, j, k;
00342 int pows[3];
00343 cpl_vector *tcheb_coef_lambda = NULL;
00344 cpl_vector *tcheb_coef_order = NULL;
00345 cpl_vector *tcheb_coef_slit = NULL;
00346
00347
00348 XSH_ASSURE_NOT_NULL( sol);
00349
00350 if (sol->deg_slit > 0){
00351 check(tcheb_slit = xsh_tools_tchebitchev_transform(slit, sol->min_slit,
00352 sol->max_slit));
00353 }
00354 check(tcheb_order = xsh_tools_tchebitchev_transform(order, sol->min_order,
00355 sol->max_order));
00356 check(tcheb_lambda = xsh_tools_tchebitchev_transform(lambda,
00357 sol->min_lambda, sol->max_lambda));
00358 if ( tcheb_lambda < -1. ) tcheb_lambda = -1. ;
00359 if ( tcheb_lambda > 1. ) tcheb_lambda = 1. ;
00360
00361 check( tcheb_coef_lambda = xsh_tools_tchebitchev_poly_eval(
00362 sol->deg_lambda, tcheb_lambda));
00363 check( tcheb_coef_order = xsh_tools_tchebitchev_poly_eval(
00364 sol->deg_order, tcheb_order));
00365 check( tcheb_coef_slit = xsh_tools_tchebitchev_poly_eval(
00366 sol->deg_slit, tcheb_slit));
00367
00368 for(i = 0; i < sol->deg_lambda +1 ; i++) {
00369 for(j = 0; j < sol->deg_order +1 ; j++){
00370 for(k = 0; k < sol->deg_slit +1 ; k++){
00371 double coef;
00372 double vlambda ;
00373 double vslit ;
00374 double vorder ;
00375
00376 check( vlambda = cpl_vector_get( tcheb_coef_lambda, i));
00377 check( vorder = cpl_vector_get( tcheb_coef_order, j));
00378 check( vslit = cpl_vector_get( tcheb_coef_slit, k));
00379 pows[0] = k;
00380 pows[1] = j;
00381 pows[2] = i;
00382 check( coef = cpl_polynomial_get_coeff(sol->polx, pows));
00383 res += coef*vslit*vorder*vlambda;
00384 }
00385 }
00386 }
00387 check(pos = xsh_tools_tchebitchev_reverse_transform(res, sol->min_x,
00388 sol->max_x));
00389
00390
00391 pos = convert_data_to_bin( pos, sol->bin_x ) ;
00392
00393 cleanup:
00394 xsh_free_vector( &tcheb_coef_lambda);
00395 xsh_free_vector( &tcheb_coef_order);
00396 xsh_free_vector( &tcheb_coef_slit);
00397 return pos;
00398 }
00399
00400
00401
00411
00412 double xsh_wavesol_eval_poly(xsh_wavesol* sol, double lambda, double order,
00413 double slit)
00414 {
00415 double tcheb_slit = 0.0, tcheb_order, tcheb_lambda;
00416 double pos = 0.0, res = 0.0;
00417 int i, j, k;
00418 int pows[3];
00419 cpl_vector *tcheb_coef_lambda = NULL;
00420 cpl_vector *tcheb_coef_order = NULL;
00421 cpl_vector *tcheb_coef_slit = NULL;
00422
00423
00424 XSH_ASSURE_NOT_NULL(sol);
00425
00426 if (sol->deg_slit > 0){
00427 check(tcheb_slit = xsh_tools_tchebitchev_transform(slit, sol->min_slit,
00428 sol->max_slit));
00429 }
00430 check(tcheb_order = xsh_tools_tchebitchev_transform(order, sol->min_order,
00431 sol->max_order));
00432 check(tcheb_lambda = xsh_tools_tchebitchev_transform(lambda,
00433 sol->min_lambda,
00434 sol->max_lambda));
00435 if ( tcheb_lambda < -1. ) tcheb_lambda = -1. ;
00436 if ( tcheb_lambda > 1. ) tcheb_lambda = 1. ;
00437
00438 check( tcheb_coef_lambda = xsh_tools_tchebitchev_poly_eval(
00439 sol->deg_lambda, tcheb_lambda));
00440 check( tcheb_coef_order = xsh_tools_tchebitchev_poly_eval(
00441 sol->deg_order, tcheb_order));
00442 check( tcheb_coef_slit = xsh_tools_tchebitchev_poly_eval(
00443 sol->deg_slit, tcheb_slit));
00444
00445 for(i = 0; i < sol->deg_lambda +1 ; i++) {
00446 for(j = 0; j < sol->deg_order +1 ; j++){
00447 for(k = 0; k < sol->deg_slit +1 ; k++){
00448 double coef;
00449 double vlambda ;
00450 double vslit ;
00451 double vorder ;
00452
00453 check( vlambda = cpl_vector_get( tcheb_coef_lambda, i));
00454 check( vorder = cpl_vector_get( tcheb_coef_order, j));
00455 check( vslit = cpl_vector_get( tcheb_coef_slit, k));
00456 pows[0] = k;
00457 pows[1] = j;
00458 pows[2] = i;
00459 check( coef = cpl_polynomial_get_coeff(sol->poly, pows));
00460 res += coef*vslit*vorder*vlambda;
00461 }
00462 }
00463 }
00464 check(pos = xsh_tools_tchebitchev_reverse_transform(res, sol->min_y,
00465 sol->max_y));
00466
00467
00468 pos = convert_data_to_bin( pos, sol->bin_y ) ;
00469
00470 cleanup:
00471 xsh_free_vector( &tcheb_coef_lambda);
00472 xsh_free_vector( &tcheb_coef_order);
00473 xsh_free_vector( &tcheb_coef_slit);
00474 return pos;
00475 }
00476
00477
00478
00479
00493
00494
00495 void xsh_wavesol_compute(xsh_wavesol* sol, int size,
00496 double* pos,double *posmin, double *posmax, double* lambda,
00497 double* order, double* slit, cpl_polynomial* result)
00498 {
00499 int i, j, k, l;
00500 double chisq;
00501 gsl_matrix *X = NULL, *cov = NULL;
00502 gsl_vector *y = NULL, *c = NULL;
00503 gsl_multifit_linear_workspace *work = NULL;
00504 int pows[3];
00505 double *tcheb_pos = NULL, *tcheb_lambda = NULL;
00506 double *tcheb_order = NULL, *tcheb_slit = NULL;
00507
00508
00509 XSH_ASSURE_NOT_NULL(sol);
00510 XSH_ASSURE_NOT_ILLEGAL( sol->nbcoefs < size);
00511
00512
00513 check(xsh_tools_min_max(size, slit, &(sol->min_slit), &(sol->max_slit)));
00514
00515
00516 sol->min_slit *= XSH_SLIT_RANGE;
00517 sol->max_slit *= XSH_SLIT_RANGE;
00518
00519 check(xsh_tools_min_max(size, pos, posmin, posmax));
00520
00521
00522 XSH_CALLOC(tcheb_lambda, double, size);
00523 XSH_CALLOC(tcheb_order, double, size);
00524 XSH_CALLOC(tcheb_slit, double, size);
00525 XSH_CALLOC(tcheb_pos, double, size);
00526
00527 check(xsh_tools_tchebitchev_transform_tab( size, lambda, sol->min_lambda,
00528 sol->max_lambda, tcheb_lambda));
00529 check(xsh_tools_tchebitchev_transform_tab( size, order, sol->min_order,
00530 sol->max_order, tcheb_order));
00531 if (sol->deg_slit > 0){
00532 check( xsh_tools_tchebitchev_transform_tab( size, slit, sol->min_slit,
00533 sol->max_slit, tcheb_slit));
00534 }
00535 check(xsh_tools_tchebitchev_transform_tab( size, pos, *posmin, *posmax,
00536 tcheb_pos));
00537
00538 X = gsl_matrix_alloc (size, sol->nbcoefs);
00539 y = gsl_vector_alloc (size);
00540 c = gsl_vector_alloc (sol->nbcoefs);
00541 cov = gsl_matrix_alloc (sol->nbcoefs, sol->nbcoefs);
00542
00543 for (i = 0; i < size; i++) {
00544 for(j = 0; j < sol->deg_lambda +1 ; j++) {
00545 for(k = 0; k < sol->deg_order +1 ; k++){
00546 for(l = 0; l < sol->deg_slit +1 ; l++){
00547 double vslit = cos(l*acos(tcheb_slit[i]));
00548 double vorder = cos(k*acos(tcheb_order[i]));
00549 double vlambda = cos(j*acos(tcheb_lambda[i]));
00550
00551 gsl_matrix_set (X, i, l+k*(sol->deg_slit+1)+j*(sol->deg_order+1)*
00552 (sol->deg_slit+1), vslit*vorder*vlambda);
00553 }
00554 }
00555 }
00556 gsl_vector_set (y, i, tcheb_pos[i]);
00557 }
00558
00559 work = gsl_multifit_linear_alloc(size,sol->nbcoefs);
00560 gsl_multifit_linear(X, y, c, cov, &chisq, work);
00561
00562 for(j=0; j< sol->deg_lambda+1; j++){
00563 for(k = 0; k < sol->deg_order +1 ; k++){
00564 for(l = 0; l < sol->deg_slit +1 ; l++){
00565 pows[0] = l;
00566 pows[1] = k;
00567 pows[2] = j;
00568 cpl_polynomial_set_coeff(result, pows, gsl_vector_get(c,
00569 l+k*(sol->deg_slit+1)+j*(sol->deg_order+1)*(sol->deg_slit+1)));
00570 }
00571 }
00572 }
00573 gsl_multifit_linear_free (work);
00574 gsl_matrix_free (X);
00575 gsl_vector_free (y);
00576 gsl_vector_free (c);
00577 gsl_matrix_free (cov);
00578
00579 cleanup:
00580 XSH_FREE(tcheb_pos);
00581 XSH_FREE(tcheb_lambda);
00582 XSH_FREE(tcheb_order);
00583 XSH_FREE(tcheb_slit);
00584 }
00585
00601 void xsh_wavesol_residual(xsh_wavesol* sol, xsh_wavesol * adj, int size,
00602 double* new_pos, double* lambda,
00603 double* order, double* slit,
00604 cpl_polynomial* result, char axis)
00605 {
00606 int i, j, k, l;
00607 double chisq;
00608 gsl_matrix *X = NULL, *cov = NULL;
00609 gsl_vector *y = NULL, *c = NULL;
00610 gsl_multifit_linear_workspace *work = NULL;
00611 int pows[3];
00612 double *tcheb_pos = NULL, *tcheb_lambda = NULL;
00613 double *tcheb_order = NULL, *tcheb_slit = NULL;
00614 double * tcheb_new = NULL ;
00615 double posmin, posmax ;
00616
00617 XSH_ASSURE_NOT_NULL(sol);
00618 XSH_ASSURE_NOT_NULL(adj);
00619
00620
00621 sol->min_lambda = adj->min_lambda ;
00622 sol->max_lambda = adj->max_lambda ;
00623 sol->min_order = adj->min_order ;
00624 sol->max_order = adj->max_order ;
00625 sol->min_slit = adj->min_slit ;
00626 sol->max_slit = adj->max_slit ;
00627 sol->min_y = adj->min_y ;
00628 sol->max_y = adj->max_y ;
00629 sol->min_x = adj->min_x ;
00630 sol->max_x = adj->max_x ;
00631
00632 if ( axis == 'y' ){
00633 posmin = adj->min_y ;
00634 posmax = adj->max_y ;
00635 }
00636 else{
00637 posmin = adj->min_x ;
00638 posmax = adj->max_x ;
00639 }
00640
00641
00642 XSH_CALLOC(tcheb_lambda, double, size);
00643 XSH_CALLOC(tcheb_order, double, size);
00644 XSH_CALLOC(tcheb_slit, double, size);
00645 XSH_CALLOC(tcheb_pos, double, size);
00646 XSH_CALLOC(tcheb_new, double, size);
00647
00648 check(xsh_tools_tchebitchev_transform_tab(size, lambda, sol->min_lambda,
00649 sol->max_lambda, tcheb_lambda));
00650 check(xsh_tools_tchebitchev_transform_tab(size, order, sol->min_order,
00651 sol->max_order,tcheb_order));
00652 if (sol->deg_slit > 0){
00653 check(xsh_tools_tchebitchev_transform_tab(size, slit, sol->min_slit,
00654 sol->max_slit,
00655 tcheb_slit));
00656 }
00657 check(xsh_tools_tchebitchev_transform_tab(size, new_pos, posmin, posmax,
00658 tcheb_new));
00659
00660
00661
00662
00663
00664 {
00665 int indx ;
00666 double temp0, temp1, a, b ;
00667
00668 a = 2/(posmax-posmin);
00669 b = 1-2*posmax/(posmax-posmin);
00670
00671 for( indx = 0 ; indx < size ; indx++ ) {
00672 if ( axis == 'y' ) {
00673 temp0 = xsh_wavesol_eval_poly( adj, lambda[indx], order[indx],
00674 slit[indx] ) ;
00675 }
00676 else
00677 temp0 = xsh_wavesol_eval_polx( adj, lambda[indx], order[indx],
00678 slit[indx] ) ;
00679
00680 temp1 = (temp0*a) + b ;
00681 tcheb_pos[indx] = tcheb_new[indx] - temp1 ;
00682 }
00683 }
00684
00685 X = gsl_matrix_alloc (size, sol->nbcoefs);
00686 y = gsl_vector_alloc (size);
00687 c = gsl_vector_alloc (sol->nbcoefs);
00688 cov = gsl_matrix_alloc (sol->nbcoefs, sol->nbcoefs);
00689
00690 for (i = 0; i < size; i++) {
00691 for(j = 0; j < sol->deg_lambda +1 ; j++) {
00692 for(k = 0; k < sol->deg_order +1 ; k++) {
00693 for(l = 0; l < sol->deg_slit +1 ; l++) {
00694 double vslit = cos(l*acos(tcheb_slit[i]));
00695 double vorder = cos(k*acos(tcheb_order[i]));
00696 double vlambda = cos(j*acos(tcheb_lambda[i]));
00697
00698 gsl_matrix_set (X, i, l+k*(sol->deg_slit+1)+j*(sol->deg_order+1)*
00699 (sol->deg_slit+1), vslit*vorder*vlambda);
00700 }
00701 }
00702 }
00703 gsl_vector_set (y, i, tcheb_pos[i]);
00704 }
00705
00706 XSH_ASSURE_NOT_ILLEGAL_MSG( size > sol->nbcoefs,
00707 "Not enough Points vs Number of Coeffs" ) ;
00708 work = gsl_multifit_linear_alloc(size, sol->nbcoefs);
00709 gsl_multifit_linear(X, y, c, cov, &chisq, work);
00710
00711 for(j=0; j< sol->deg_lambda+1; j++){
00712 for(k = 0; k < sol->deg_order +1 ; k++){
00713 for(l = 0; l < sol->deg_slit +1 ; l++){
00714 pows[0] = l;
00715 pows[1] = k;
00716 pows[2] = j;
00717 cpl_polynomial_set_coeff(result, pows, gsl_vector_get(c,
00718 l+k*(sol->deg_slit+1)+j*(sol->deg_order+1)*(sol->deg_slit+1)));
00719 }
00720 }
00721 }
00722 gsl_multifit_linear_free (work);
00723 gsl_matrix_free (X);
00724 gsl_vector_free (y);
00725 gsl_vector_free (c);
00726 gsl_matrix_free (cov);
00727
00728 cleanup:
00729 XSH_FREE(tcheb_pos);
00730 XSH_FREE(tcheb_lambda);
00731 XSH_FREE(tcheb_order);
00732 XSH_FREE(tcheb_slit);
00733 XSH_FREE(tcheb_new);
00734 }
00735
00736
00746
00747
00748 cpl_frame*
00749 xsh_wavesol_save(xsh_wavesol *w, cpl_table* trace, const char* filename,const char* tag)
00750 {
00751 cpl_table* table = NULL;
00752 cpl_frame * result = NULL ;
00753 int i, j, k;
00754 char coefname[20];
00755 int pows[3];
00756
00757 XSH_ASSURE_NOT_NULL(w);
00758 XSH_ASSURE_NOT_NULL(trace);
00759 XSH_ASSURE_NOT_NULL(filename);
00760
00761
00762 check(table = cpl_table_new(XSH_WAVESOL_TABLE_NB_COL+w->nbcoefs));
00763 check(cpl_table_set_size(table,XSH_WAVESOL_TABLE_NB_ROWS));
00764
00765 check(
00766 cpl_table_new_column(table,XSH_WAVESOL_TABLE_COLNAME_AXIS,
00767 CPL_TYPE_STRING));
00768 check(
00769 cpl_table_new_column(table,XSH_WAVESOL_TABLE_COLNAME_DEGSLIT,
00770 CPL_TYPE_INT));
00771 check(
00772 cpl_table_new_column(table,XSH_WAVESOL_TABLE_COLNAME_DEGORDER,
00773 CPL_TYPE_INT));
00774 check(
00775 cpl_table_new_column(table,XSH_WAVESOL_TABLE_COLNAME_DEGLAMBDA,
00776 CPL_TYPE_INT));
00777 check(cpl_table_set_string(table,XSH_WAVESOL_TABLE_COLNAME_AXIS,
00778 0,"X"));
00779 check(cpl_table_set_string(table,XSH_WAVESOL_TABLE_COLNAME_AXIS,
00780 1,"Y"));
00781 check(cpl_table_set_int(table,XSH_WAVESOL_TABLE_COLNAME_DEGSLIT,
00782 0,w->deg_slit));
00783 check(cpl_table_set_int(table,XSH_WAVESOL_TABLE_COLNAME_DEGSLIT,
00784 1, w->deg_slit));
00785 check(cpl_table_set_int(table,XSH_WAVESOL_TABLE_COLNAME_DEGORDER,
00786 0,w->deg_order));
00787 check(cpl_table_set_int(table,XSH_WAVESOL_TABLE_COLNAME_DEGORDER,
00788 1, w->deg_order));
00789 check(cpl_table_set_int(table,XSH_WAVESOL_TABLE_COLNAME_DEGLAMBDA,
00790 0,w->deg_lambda));
00791 check(cpl_table_set_int(table,XSH_WAVESOL_TABLE_COLNAME_DEGLAMBDA,
00792 1, w->deg_lambda));
00793
00794 for(i=0; i<= w->deg_slit; i++){
00795 for(j = 0; j<= w->deg_order; j++){
00796 for(k = 0; k <= w->deg_lambda; k++){
00797 double vx = 0.0, vy = 0.0;
00798 pows[0] = i;
00799 pows[1] = j;
00800 pows[2] = k;
00801 sprintf(coefname,"C%d%d%d",i,j,k);
00802 check(cpl_table_new_column(table,coefname, CPL_TYPE_DOUBLE));
00803 check(vx = cpl_polynomial_get_coeff(w->polx, pows));
00804 check(vy = cpl_polynomial_get_coeff(w->poly, pows));
00805 check(cpl_table_set_double( table, coefname, 0, vx));
00806 check(cpl_table_set_double( table, coefname, 1, vy));
00807 }
00808 }
00809 }
00810
00811
00812 check(xsh_pfits_set_wavesol_lambda_min(w->header, w->min_lambda));
00813 check(xsh_pfits_set_wavesol_lambda_max(w->header, w->max_lambda));
00814 check(xsh_pfits_set_wavesol_order_min(w->header, w->min_order));
00815 check(xsh_pfits_set_wavesol_order_max(w->header, w->max_order));
00816 check(xsh_pfits_set_wavesol_slit_min(w->header, w->min_slit));
00817 check(xsh_pfits_set_wavesol_slit_max(w->header, w->max_slit));
00818 check(xsh_pfits_set_wavesol_x_min(w->header, w->min_x));
00819 check(xsh_pfits_set_wavesol_x_max(w->header, w->max_x));
00820 check(xsh_pfits_set_wavesol_y_min(w->header, w->min_y));
00821 check(xsh_pfits_set_wavesol_y_max(w->header, w->max_y));
00822 check( xsh_pfits_set_pcatg( w->header, tag));
00823 check(cpl_table_save(table, w->header, NULL, filename, CPL_IO_DEFAULT));
00824 check(cpl_table_save(trace, NULL, NULL, filename, CPL_IO_EXTEND));
00825
00826
00827 check(result=xsh_frame_product(filename,
00828 tag,
00829 CPL_FRAME_TYPE_TABLE,
00830 CPL_FRAME_GROUP_PRODUCT,
00831 CPL_FRAME_LEVEL_TEMPORARY));
00832
00833 cleanup:
00834 XSH_TABLE_FREE( table);
00835 return result;
00836 }
00837
00838
00839
00847
00848
00849 xsh_wavesol * xsh_wavesol_load( cpl_frame * frame,
00850 xsh_instrument * instrument )
00851 {
00852 cpl_table* table = NULL;
00853 const char* tablename = NULL;
00854 xsh_wavesol * result = NULL ;
00855 int rows, cols ;
00856
00857 XSH_ASSURE_NOT_NULL(frame);
00858
00859
00860 check(tablename = cpl_frame_get_filename(frame));
00861
00862 XSH_TABLE_LOAD( table, tablename);
00863
00864 XSH_CALLOC( result, xsh_wavesol, 1 ) ;
00865
00866 check(result->header = cpl_propertylist_load(tablename,0));
00867
00868 check( rows = cpl_table_get_nrow( table ) ) ;
00869 check( cols = cpl_table_get_ncol( table ) ) ;
00870 xsh_msg_dbg_high( " === Wavesol Rows = %d, Columns = %d", rows, cols ) ;
00871
00872
00873 result->bin_x = 1;
00874 result->bin_y = 1;
00875
00876 check( result->min_lambda = xsh_pfits_get_wavesol_lambda_min(
00877 result->header));
00878 check( result->max_lambda = xsh_pfits_get_wavesol_lambda_max(
00879 result->header));
00880 check( result->min_order = xsh_pfits_get_wavesol_order_min(
00881 result->header));
00882 check( result->max_order = xsh_pfits_get_wavesol_order_max(
00883 result->header));
00884 check( result->min_slit = xsh_pfits_get_wavesol_slit_min(
00885 result->header));
00886 check( result->max_slit = xsh_pfits_get_wavesol_slit_max(
00887 result->header));
00888 check( result->min_x = xsh_pfits_get_wavesol_x_min(
00889 result->header));
00890 check( result->max_x = xsh_pfits_get_wavesol_x_max(
00891 result->header));
00892 check( result->min_y = xsh_pfits_get_wavesol_y_min(
00893 result->header));
00894 check( result->max_y = xsh_pfits_get_wavesol_y_max(
00895 result->header));
00896
00897
00898 {
00899 int i, j, k ;
00900
00901
00902 check( xsh_get_table_value( table, XSH_WAVESOL_TABLE_COLNAME_DEGSLIT,
00903 CPL_TYPE_INT, 0, &result->deg_slit ) ) ;
00904 check( xsh_get_table_value( table, XSH_WAVESOL_TABLE_COLNAME_DEGORDER,
00905 CPL_TYPE_INT, 0, &result->deg_order ) ) ;
00906 check( xsh_get_table_value( table, XSH_WAVESOL_TABLE_COLNAME_DEGLAMBDA,
00907 CPL_TYPE_INT, 0, &result->deg_lambda ) ) ;
00908
00909 result->nbcoefs = (result->deg_order+1)*(result->deg_slit+1)*
00910 (result->deg_lambda+1);
00911 xsh_msg_dbg_high( "nbcoef %d deg_lambda %d deg_n %d deg_s %d", result->nbcoefs,
00912 result->deg_lambda, result->deg_order, result->deg_slit);
00913 result->polx = cpl_polynomial_new(3);
00914 result->poly = cpl_polynomial_new(3);
00915 result->dim = cpl_vector_new(3);
00916
00917
00918 for( i = 0 ; i<=result->deg_slit ; i++ ) {
00919 for( j = 0 ; j<=result->deg_order ; j++ ) {
00920 for( k = 0 ; k<=result->deg_lambda ; k++ ) {
00921 char coefname[16] ;
00922 double coef ;
00923 int pows[3] ;
00924
00925 pows[0] = i ;
00926 pows[1] = j ;
00927 pows[2] = k ;
00928 sprintf( coefname, "C%d%d%d", i, j, k ) ;
00929
00930 check( xsh_get_table_value( table, coefname, CPL_TYPE_DOUBLE,
00931 0, &coef ) ) ;
00932
00933 check( cpl_polynomial_set_coeff( result->polx, pows, coef ) ) ;
00934 check( xsh_get_table_value( table, coefname, CPL_TYPE_DOUBLE,
00935 1, &coef ) ) ;
00936
00937 check( cpl_polynomial_set_coeff( result->poly, pows, coef ) ) ;
00938 }
00939 }
00940 }
00941 }
00942
00943
00944 result->bin_x = xsh_instrument_get_binx( instrument ) ;
00945 result->bin_y = xsh_instrument_get_biny( instrument ) ;
00946
00947 cleanup:
00948 if (cpl_error_get_code () != CPL_ERROR_NONE) {
00949 xsh_error_msg("can't load Wavesol frame %s", tablename );
00950 xsh_wavesol_free(&result);
00951 }
00952 XSH_TABLE_FREE( table) ;
00953
00954 return result ;
00955 }
00956
00957 void xsh_wavesol_dump( xsh_wavesol * wsol, const char * fname, int nb )
00958 {
00959 FILE * fout = NULL ;
00960 int i, j, k, l ;
00961 int pows[3] ;
00962
00963 if ( fname != NULL ) fout = fopen( fname, "w" ) ;
00964 l = 0 ;
00965
00966 for( i = 0 ; i <= wsol->deg_slit ; i++ )
00967 for( j = 0 ; j <= wsol->deg_lambda ; j++ )
00968 for( k = 0 ; k <= wsol->deg_order ; k++ ) {
00969 double coeff ;
00970 pows[0] = i ;
00971 pows[1] = j;
00972 pows[2] = k;
00973 check(coeff = cpl_polynomial_get_coeff( wsol->poly, pows));
00974 if ( fout == NULL )
00975 xsh_msg( " %d%d%d; %lf", i, j, k, coeff ) ;
00976 else
00977 fprintf( fout, "%d%d%d: %lf\n", i, j, k, coeff ) ;
00978 l++ ;
00979 if ( nb != 0 && l >= nb ) goto cleanup ;
00980 }
00981 cleanup:
00982 if ( fout != NULL ) fclose( fout ) ;
00983 }
00984
00985 cpl_table*
00986 xsh_wavesol_trace( xsh_wavesol * wsol, double* lambda,
00987 double* order, double* slit,int size)
00988 {
00989 cpl_table* table = NULL ;
00990 int i;
00991 double* pw=NULL;
00992 double* po=NULL;
00993 double* px=NULL;
00994 double* py=NULL;
00995 double* ps=NULL;
00996
00997
00998
00999 XSH_ASSURE_NOT_NULL( wsol);
01000 XSH_ASSURE_NOT_NULL( lambda);
01001 XSH_ASSURE_NOT_NULL( order);
01002 XSH_ASSURE_NOT_NULL( slit);
01003
01004
01005 table=cpl_table_new(size);
01006 cpl_table_new_column(table,"WAVELENGTH",CPL_TYPE_DOUBLE);
01007 cpl_table_new_column(table,"ORDER",CPL_TYPE_DOUBLE);
01008 cpl_table_new_column(table,"X",CPL_TYPE_DOUBLE);
01009 cpl_table_new_column(table,"Y",CPL_TYPE_DOUBLE);
01010 cpl_table_new_column(table,"S",CPL_TYPE_DOUBLE);
01011
01012 cpl_table_fill_column_window(table,"WAVELENGTH",0,size,0);
01013 cpl_table_fill_column_window(table,"ORDER",0,size,0);
01014 cpl_table_fill_column_window(table,"X",0,size,0);
01015 cpl_table_fill_column_window(table,"Y",0,size,0);
01016 cpl_table_fill_column_window(table,"S",0,size,0);
01017
01018 po=cpl_table_get_data_double(table,"ORDER");
01019 pw=cpl_table_get_data_double(table,"WAVELENGTH");
01020 px=cpl_table_get_data_double(table,"X");
01021 py=cpl_table_get_data_double(table,"Y");
01022 ps=cpl_table_get_data_double(table,"S");
01023
01024 for(i = 0 ; i<size; i++) {
01025 pw[i] = lambda[i];
01026 po[i] = order[i];
01027 ps[i] = slit[i];
01028 check( px[i] = xsh_wavesol_eval_polx(wsol,pw[i],po[i],ps[i]));
01029 check( py[i] = xsh_wavesol_eval_poly(wsol,pw[i],po[i],ps[i]));
01030 }
01031
01032 cleanup:
01033 return table;
01034 }
01035
01044 void xsh_wavesol_set_bin_x( xsh_wavesol * wsol, int bin )
01045 {
01046 XSH_ASSURE_NOT_NULL( wsol);
01047 wsol->bin_x = bin;
01048
01049 cleanup:
01050 return;
01051 }
01052
01061 void xsh_wavesol_set_bin_y( xsh_wavesol * wsol, int bin )
01062 {
01063 XSH_ASSURE_NOT_NULL( wsol);
01064 wsol->bin_y = bin;
01065
01066 cleanup:
01067 return;
01068 }
01069
01080 void xsh_wavesol_apply_shift( xsh_wavesol *wsol, float shift_x, float shift_y)
01081 {
01082
01083 XSH_ASSURE_NOT_NULL( wsol);
01084
01085 wsol->min_x = wsol->min_x+shift_x;
01086 wsol->max_x = wsol->max_x+shift_x;
01087 wsol->min_y = wsol->min_y+shift_y;
01088 wsol->max_y = wsol->max_y+shift_y;
01089
01090 cleanup:
01091 return;
01092 }