X-shooter Pipeline Reference Manual 3.8.15
xsh_data_wavesol.c
Go to the documentation of this file.
1/* *
2 * This file is part of the ESO X-shooter Pipeline *
3 * Copyright (C) 2006 European Southern Observatory *
4 * *
5 * This library is free software; you can redistribute it and/or modify *
6 * it under the terms of the GNU General Public License as published by *
7 * the Free Software Foundation; either version 2 of the License, or *
8 * (at your option) any later version. *
9 * *
10 * This program is distributed in the hope that it will be useful, *
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
13 * GNU General Public License for more details. *
14 * *
15 * You should have received a copy of the GNU General Public License *
16 * along with this program; if not, write to the Free Software *
17 * Foundation, 51 Franklin St, Fifth Floor, Boston, MA 02111-1307 USA *
18 * */
19
20/*
21 * $Author: amodigli $
22 * $Date: 2011-12-09 10:10:12 $
23 * $Revision: 1.59 $
24 * $Name: not supported by cvs2svn $
25 */
26
27#ifdef HAVE_CONFIG_H
28#include <config.h>
29#endif
30
31/*---------------------------------------------------------------------------*/
36/*---------------------------------------------------------------------------*/
37
41/*-----------------------------------------------------------------------------
42 Includes
43 ----------------------------------------------------------------------------*/
44#include <xsh_data_wavesol.h>
46#include <xsh_utils.h>
47#include <xsh_error.h>
48#include <xsh_msg.h>
49#include <xsh_pfits.h>
50#include <xsh_dfs.h>
51#include <math.h>
52#include <gsl/gsl_multifit.h>
53#include <cpl.h>
54#include <xsh_utils_table.h>
55/*----------------------------------------------------------------------------
56 Function implementation
57 ----------------------------------------------------------------------------*/
58
59/*---------------------------------------------------------------------------*/
72/*---------------------------------------------------------------------------*/
73xsh_wavesol* xsh_wavesol_create( cpl_frame* spectral_format_frame,
75{
76 xsh_wavesol* res = NULL;
77 xsh_spectralformat_list *spec_list = NULL;
78 int i;
79 float sol_lambda_min=0, sol_lambda_max=0;
80 int sol_absorder_min=0, sol_absorder_max=0;
81
82 /* Check parameters */
83 XSH_ASSURE_NOT_NULL( params);
84 XSH_ASSURE_NOT_NULL( spectral_format_frame);
86
87 /* Allocate memory */
88 XSH_MALLOC(res, xsh_wavesol, 1);
89 res->deg_slit = params->wavesol_deg_slit;
90 res->deg_order = params->wavesol_deg_order;
91 res->deg_lambda = params->wavesol_deg_lambda;
92 res->nbcoefs = (res->deg_order+1)*(res->deg_slit+1)*(res->deg_lambda+1);
93 xsh_msg_dbg_high("nbcoef %d deg_lambda %d deg_n %d deg_s %d",res->nbcoefs,
94 res->deg_lambda, res->deg_order, res->deg_slit);
95 res->polx = cpl_polynomial_new(3);
96 res->poly = cpl_polynomial_new(3);
97 res->dim = cpl_vector_new(3);
98 res->header = cpl_propertylist_new();
99
100 /* Add binx,y */
103
104 /* compute range of lambda and orders */
105 check( spec_list = xsh_spectralformat_list_load( spectral_format_frame,
106 instrument));
107 sol_lambda_min = spec_list->list[0].lambda_min_full;
108 sol_lambda_max = spec_list->list[0].lambda_max_full;
109 sol_absorder_min = spec_list->list[0].absorder;
110 sol_absorder_max = spec_list->list[0].absorder;
111
112 for( i=1; i< spec_list->size; i++){
113 int absorder;
114 float lambda_min, lambda_max;
115 absorder = spec_list->list[i].absorder;
116 lambda_min = spec_list->list[i].lambda_min_full;
117 lambda_max = spec_list->list[i].lambda_max_full;
118 if ( absorder > sol_absorder_max){
119 sol_absorder_max = absorder;
120 }
121 if ( absorder < sol_absorder_min){
122 sol_absorder_min = absorder;
123 }
124 if ( lambda_min < sol_lambda_min){
125 sol_lambda_min = lambda_min;
126 }
127 if ( lambda_max > sol_lambda_max){
128 sol_lambda_max = lambda_max;
129 }
130 }
131 xsh_msg_dbg_high("Order range %d-%d",sol_absorder_min, sol_absorder_max);
132 xsh_msg_dbg_high("Lambda range %f-%f", sol_lambda_min, sol_lambda_max);
133 res->min_lambda = sol_lambda_min;
134 res->max_lambda = sol_lambda_max;
135 res->min_order = sol_absorder_min;
136 res->max_order = sol_absorder_max;
137
138 cleanup:
139 xsh_spectralformat_list_free( &spec_list);
140 return res;
141}
142
149/*---------------------------------------------------------------------------*/
151{
152 xsh_wavesol * res = NULL ;
153
154 XSH_MALLOC( res, xsh_wavesol, 1 ) ;
155 res->bin_x = org->bin_x ;
156 res->bin_y = org->bin_y ;
157 res->deg_slit = org->deg_slit ;
158 res->deg_order = org->deg_order ;
159 res->deg_lambda = org->deg_lambda ;
160 res->min_lambda = org->min_lambda ;
161 res->max_lambda = org->max_lambda ;
162 res->min_order = org->min_order ;
163 res->max_order = org->max_order ;
164 res->min_slit = org->min_slit ;
165 res->max_slit = org->max_slit ;
166 res->min_x = org->min_x ;
167 res->max_x = org->max_x ;
168 res->min_y = org->min_y ;
169 res->max_y = org->max_y ;
170 res->nbcoefs = org->nbcoefs ;
171 res->polx = cpl_polynomial_duplicate( org->polx ) ;
172 res->poly = cpl_polynomial_duplicate( org->poly ) ;
173 res->dim = cpl_vector_duplicate( org->dim ) ;
174 res->header = cpl_propertylist_duplicate( org->header ) ;
175
176 cleanup:
177 return res ;
178}
179
187{
188 int i, j, k ;
189 cpl_size pows[3];
190
191 pows[0] = 0;
192 i = 0 ;
193
194 for(j = 0; j<= from->deg_order; j++)
195 for(k = 0; k <= from->deg_lambda; k++) {
196 double vto = 0.0, vfrom = 0.0;
197
198 pows[0] = i ;
199 pows[1] = j;
200 pows[2] = k;
201 xsh_msg_dbg_high( "Add_poly: %d %d %d", i, j, k ) ;
202 check(vfrom = cpl_polynomial_get_coeff(from->poly, pows));
203 check(vto = cpl_polynomial_get_coeff(to->poly, pows));
204 vto += vfrom;
205 check( cpl_polynomial_set_coeff( to->poly, pows, vto));
206 }
207
208 cleanup:
209 return ;
210}
211
212/*---------------------------------------------------------------------------*/
218/*---------------------------------------------------------------------------*/
220{
221 if (w && *w) {
222 xsh_free_polynomial(&(*w)->polx);
223 xsh_free_polynomial(&(*w)->poly);
224 xsh_free_vector(&(*w)->dim);
225 xsh_free_propertylist(&(*w)->header);
226 cpl_free(*w);
227 *w = NULL;
228 }
229}
230
231
232/*---------------------------------------------------------------------------*/
239/*---------------------------------------------------------------------------*/
241{
242 cpl_polynomial* res = NULL;
243
245 res = sol->poly;
246
247 cleanup:
248 return res;
249}
250
251
252/*---------------------------------------------------------------------------*/
259/*---------------------------------------------------------------------------*/
260
262{
263 XSH_ASSURE_NOT_NULL( wsol);
264 wsol->type = type;
265
266 cleanup:
267 return;
268}
269
270
271
272/*---------------------------------------------------------------------------*/
280{
282 XSH_ASSURE_NOT_NULL( wsol);
283 type = wsol->type;
284
285 cleanup:
286 return type;
287}
288
289/*---------------------------------------------------------------------------*/
296/*---------------------------------------------------------------------------*/
298{
299 cpl_polynomial* res = NULL;
300
302 res = sol->polx;
303
304 cleanup:
305 return res;
306}
307
308/*---------------------------------------------------------------------------*/
314/*---------------------------------------------------------------------------*/
315cpl_propertylist* xsh_wavesol_get_header(xsh_wavesol* sol)
316{
317 cpl_propertylist * res = NULL;
318
320 res = sol->header;
321 cleanup:
322 return res;
323}
324
325/*---------------------------------------------------------------------------*/
335/*---------------------------------------------------------------------------*/
336double xsh_wavesol_eval_polx(xsh_wavesol* sol, double lambda, double order,
337 double slit)
338{
339 double tcheb_slit = 0.0, tcheb_order, tcheb_lambda;
340 double pos = 0.0, res = 0.0;
341 int i, j, k;
342 cpl_size pows[3];
343 cpl_vector *tcheb_coef_lambda = NULL;
344 cpl_vector *tcheb_coef_order = NULL;
345 cpl_vector *tcheb_coef_slit = NULL;
346
347 /* Check input parameters */
349
350 if (sol->deg_slit > 0){
351 check(tcheb_slit = xsh_tools_tchebitchev_transform(slit, sol->min_slit,
352 sol->max_slit));
353 }
355 sol->max_order));
356 check(tcheb_lambda = xsh_tools_tchebitchev_transform(lambda,
357 sol->min_lambda, sol->max_lambda));
358 if ( tcheb_lambda < -1. ) tcheb_lambda = -1. ;
359 if ( tcheb_lambda > 1. ) tcheb_lambda = 1. ;
360
361 check( tcheb_coef_lambda = xsh_tools_tchebitchev_poly_eval(
362 sol->deg_lambda, tcheb_lambda));
363 check( tcheb_coef_order = xsh_tools_tchebitchev_poly_eval(
364 sol->deg_order, tcheb_order));
365 check( tcheb_coef_slit = xsh_tools_tchebitchev_poly_eval(
366 sol->deg_slit, tcheb_slit));
367
368 for(i = 0; i < sol->deg_lambda +1 ; i++) {
369 for(j = 0; j < sol->deg_order +1 ; j++){
370 for(k = 0; k < sol->deg_slit +1 ; k++){
371 double coef;
372 double vlambda ;
373 double vslit ;
374 double vorder ;
375
376 check( vlambda = cpl_vector_get( tcheb_coef_lambda, i));
377 check( vorder = cpl_vector_get( tcheb_coef_order, j));
378 check( vslit = cpl_vector_get( tcheb_coef_slit, k));
379 pows[0] = k;
380 pows[1] = j;
381 pows[2] = i;
382 check( coef = cpl_polynomial_get_coeff(sol->polx, pows));
383 res += coef*vslit*vorder*vlambda;
384 }
385 }
386 }
388 sol->max_x));
389
390 /* Convert with binning */
391 pos = convert_data_to_bin( pos, sol->bin_x ) ;
392
393 cleanup:
394 xsh_free_vector( &tcheb_coef_lambda);
395 xsh_free_vector( &tcheb_coef_order);
396 xsh_free_vector( &tcheb_coef_slit);
397 return pos;
398}
399
400
401/*---------------------------------------------------------------------------*/
411/*---------------------------------------------------------------------------*/
412double xsh_wavesol_eval_poly(xsh_wavesol* sol, double lambda, double order,
413 double slit)
414{
415 double tcheb_slit = 0.0, tcheb_order, tcheb_lambda;
416 double pos = 0.0, res = 0.0;
417 int i, j, k;
418 cpl_size pows[3];
419 cpl_vector *tcheb_coef_lambda = NULL;
420 cpl_vector *tcheb_coef_order = NULL;
421 cpl_vector *tcheb_coef_slit = NULL;
422
423
425
426 if (sol->deg_slit > 0){
427 check(tcheb_slit = xsh_tools_tchebitchev_transform(slit, sol->min_slit,
428 sol->max_slit));
429 }
431 sol->max_order));
432 check(tcheb_lambda = xsh_tools_tchebitchev_transform(lambda,
433 sol->min_lambda,
434 sol->max_lambda));
435 if ( tcheb_lambda < -1. ) tcheb_lambda = -1. ;
436 if ( tcheb_lambda > 1. ) tcheb_lambda = 1. ;
437
438 check( tcheb_coef_lambda = xsh_tools_tchebitchev_poly_eval(
439 sol->deg_lambda, tcheb_lambda));
440 check( tcheb_coef_order = xsh_tools_tchebitchev_poly_eval(
441 sol->deg_order, tcheb_order));
442 check( tcheb_coef_slit = xsh_tools_tchebitchev_poly_eval(
443 sol->deg_slit, tcheb_slit));
444
445 for(i = 0; i < sol->deg_lambda +1 ; i++) {
446 for(j = 0; j < sol->deg_order +1 ; j++){
447 for(k = 0; k < sol->deg_slit +1 ; k++){
448 double coef;
449 double vlambda ;
450 double vslit ;
451 double vorder ;
452
453 check( vlambda = cpl_vector_get( tcheb_coef_lambda, i));
454 check( vorder = cpl_vector_get( tcheb_coef_order, j));
455 check( vslit = cpl_vector_get( tcheb_coef_slit, k));
456 pows[0] = k;
457 pows[1] = j;
458 pows[2] = i;
459 check( coef = cpl_polynomial_get_coeff(sol->poly, pows));
460 res += coef*vslit*vorder*vlambda;
461 }
462 }
463 }
465 sol->max_y));
466
467 /* Convert with binning */
468 pos = convert_data_to_bin( pos, sol->bin_y ) ;
469
470 cleanup:
471 xsh_free_vector( &tcheb_coef_lambda);
472 xsh_free_vector( &tcheb_coef_order);
473 xsh_free_vector( &tcheb_coef_slit);
474 return pos;
475}
476
477
478
479/*---------------------------------------------------------------------------*/
493/*---------------------------------------------------------------------------*/
494
496 double* pos,double *posmin, double *posmax, double* lambda,
497 double* order, double* slit, cpl_polynomial* result)
498{
499 int i, j, k, l;
500 double chisq;
501 gsl_matrix *X = NULL, *cov = NULL;
502 gsl_vector *y = NULL, *c = NULL;
503 gsl_multifit_linear_workspace *work = NULL;
504 cpl_size pows[3];
505 double *tcheb_pos = NULL, *tcheb_lambda = NULL;
506 double *tcheb_order = NULL, *tcheb_slit = NULL;
507
508
511
512 /* search min and max */
513 check(xsh_tools_min_max(size, slit, &(sol->min_slit), &(sol->max_slit)));
514
515 /* extrapolate limit in slit */
516 sol->min_slit *= XSH_SLIT_RANGE;
517 sol->max_slit *= XSH_SLIT_RANGE;
518
519 check(xsh_tools_min_max(size, pos, posmin, posmax));
520
521 /* transform in tchebitchev plan */
522 XSH_CALLOC(tcheb_lambda, double, size);
523 XSH_CALLOC(tcheb_order, double, size);
524 XSH_CALLOC(tcheb_slit, double, size);
525 XSH_CALLOC(tcheb_pos, double, size);
526
528 sol->max_lambda, tcheb_lambda));
530 sol->max_order, tcheb_order));
531 if (sol->deg_slit > 0){
533 sol->max_slit, tcheb_slit));
534 }
535 check(xsh_tools_tchebitchev_transform_tab( size, pos, *posmin, *posmax,
536 tcheb_pos));
537
538 X = gsl_matrix_alloc (size, sol->nbcoefs);
539 y = gsl_vector_alloc (size);
540 c = gsl_vector_alloc (sol->nbcoefs);
541 cov = gsl_matrix_alloc (sol->nbcoefs, sol->nbcoefs);
542
543 for (i = 0; i < size; i++) {
544 for(j = 0; j < sol->deg_lambda +1 ; j++) {
545 for(k = 0; k < sol->deg_order +1 ; k++){
546 for(l = 0; l < sol->deg_slit +1 ; l++){
547 double vslit = cos(l*acos(tcheb_slit[i]));
548 double vorder = cos(k*acos(tcheb_order[i]));
549 double vlambda = cos(j*acos(tcheb_lambda[i]));
550
551 gsl_matrix_set (X, i, l+k*(sol->deg_slit+1)+j*(sol->deg_order+1)*
552 (sol->deg_slit+1), vslit*vorder*vlambda);
553 }
554 }
555 }
556 gsl_vector_set (y, i, tcheb_pos[i]);
557 }
558
559 work = gsl_multifit_linear_alloc(size,sol->nbcoefs);
560 gsl_multifit_linear(X, y, c, cov, &chisq, work);
561
562 for(j=0; j< sol->deg_lambda+1; j++){
563 for(k = 0; k < sol->deg_order +1 ; k++){
564 for(l = 0; l < sol->deg_slit +1 ; l++){
565 pows[0] = l;
566 pows[1] = k;
567 pows[2] = j;
568 cpl_polynomial_set_coeff(result, pows, gsl_vector_get(c,
569 l+k*(sol->deg_slit+1)+j*(sol->deg_order+1)*(sol->deg_slit+1)));
570 }
571 }
572 }
573 gsl_multifit_linear_free (work);
574 gsl_matrix_free (X);
575 gsl_vector_free (y);
576 gsl_vector_free (c);
577 gsl_matrix_free (cov);
578
579 cleanup:
580 XSH_FREE(tcheb_pos);
581 XSH_FREE(tcheb_lambda);
582 XSH_FREE(tcheb_order);
583 XSH_FREE(tcheb_slit);
584}
585
602 double* new_pos, double* lambda,
603 double* order, double* slit,
604 cpl_polynomial* result, char axis)
605{
606 int i, j, k, l;
607 double chisq;
608 gsl_matrix *X = NULL, *cov = NULL;
609 gsl_vector *y = NULL, *c = NULL;
610 gsl_multifit_linear_workspace *work = NULL;
611 cpl_size pows[3];
612 double *tcheb_pos = NULL, *tcheb_lambda = NULL;
613 double *tcheb_order = NULL, *tcheb_slit = NULL;
614 double * tcheb_new = NULL ;
615 double posmin, posmax ;
616
619
620 /* set min and max */
621 sol->min_lambda = adj->min_lambda ;
622 sol->max_lambda = adj->max_lambda ;
623 sol->min_order = adj->min_order ;
624 sol->max_order = adj->max_order ;
625 sol->min_slit = adj->min_slit ;
626 sol->max_slit = adj->max_slit ;
627 sol->min_y = adj->min_y ;
628 sol->max_y = adj->max_y ;
629 sol->min_x = adj->min_x ;
630 sol->max_x = adj->max_x ;
631
632 if ( axis == 'y' ){
633 posmin = adj->min_y ;
634 posmax = adj->max_y ;
635 }
636 else{
637 posmin = adj->min_x ;
638 posmax = adj->max_x ;
639 }
640
641 /* transform in tchebitchev space */
642 XSH_CALLOC(tcheb_lambda, double, size);
643 XSH_CALLOC(tcheb_order, double, size);
644 XSH_CALLOC(tcheb_slit, double, size);
645 XSH_CALLOC(tcheb_pos, double, size);
646 XSH_CALLOC(tcheb_new, double, size);
647
649 sol->max_lambda, tcheb_lambda));
651 sol->max_order,tcheb_order));
652 if (sol->deg_slit > 0){
654 sol->max_slit,
655 tcheb_slit));
656 }
657 check(xsh_tools_tchebitchev_transform_tab(size, new_pos, posmin, posmax,
658 tcheb_new));
659
660 /* Calculate residuals in Tcheb Space
661 tcheb_new = fitted (new) positions in tcheb space
662 tcheb_pos = residual in tcheb space
663 */
664 {
665 int indx ;
666 double temp0, temp1, a, b ;
667
668 a = 2/(posmax-posmin);
669 b = 1-2*posmax/(posmax-posmin);
670
671 for( indx = 0 ; indx < size ; indx++ ) {
672 if ( axis == 'y' ) {
673 temp0 = xsh_wavesol_eval_poly( adj, lambda[indx], order[indx],
674 slit[indx] ) ;
675 }
676 else
677 temp0 = xsh_wavesol_eval_polx( adj, lambda[indx], order[indx],
678 slit[indx] ) ;
679
680 temp1 = (temp0*a) + b ;
681 tcheb_pos[indx] = tcheb_new[indx] - temp1 ;
682 }
683 }
684
685 X = gsl_matrix_alloc (size, sol->nbcoefs);
686 y = gsl_vector_alloc (size);
687 c = gsl_vector_alloc (sol->nbcoefs);
688 cov = gsl_matrix_alloc (sol->nbcoefs, sol->nbcoefs);
689
690 for (i = 0; i < size; i++) {
691 for(j = 0; j < sol->deg_lambda +1 ; j++) {
692 for(k = 0; k < sol->deg_order +1 ; k++) {
693 for(l = 0; l < sol->deg_slit +1 ; l++) {
694 double vslit = cos(l*acos(tcheb_slit[i]));
695 double vorder = cos(k*acos(tcheb_order[i]));
696 double vlambda = cos(j*acos(tcheb_lambda[i]));
697
698 gsl_matrix_set (X, i, l+k*(sol->deg_slit+1)+j*(sol->deg_order+1)*
699 (sol->deg_slit+1), vslit*vorder*vlambda);
700 }
701 }
702 }
703 gsl_vector_set (y, i, tcheb_pos[i]);
704 }
705
707 "Not enough Points vs Number of Coeffs" ) ;
708 work = gsl_multifit_linear_alloc(size, sol->nbcoefs);
709 gsl_multifit_linear(X, y, c, cov, &chisq, work);
710
711 for(j=0; j< sol->deg_lambda+1; j++){
712 for(k = 0; k < sol->deg_order +1 ; k++){
713 for(l = 0; l < sol->deg_slit +1 ; l++){
714 pows[0] = l;
715 pows[1] = k;
716 pows[2] = j;
717 cpl_polynomial_set_coeff(result, pows, gsl_vector_get(c,
718 l+k*(sol->deg_slit+1)+j*(sol->deg_order+1)*(sol->deg_slit+1)));
719 }
720 }
721 }
722 gsl_multifit_linear_free (work);
723 gsl_matrix_free (X);
724 gsl_vector_free (y);
725 gsl_vector_free (c);
726 gsl_matrix_free (cov);
727
728 cleanup:
729 XSH_FREE(tcheb_pos);
730 XSH_FREE(tcheb_lambda);
731 XSH_FREE(tcheb_order);
732 XSH_FREE(tcheb_slit);
733 XSH_FREE(tcheb_new);
734}
735
736/*---------------------------------------------------------------------------*/
746/*---------------------------------------------------------------------------*/
747
748cpl_frame*
749xsh_wavesol_save(xsh_wavesol *w, cpl_table* trace, const char* filename,const char* tag)
750{
751 cpl_table* table = NULL;
752 cpl_frame * result = NULL ;
753 int i, j, k;
754 char coefname[20];
755 cpl_size pows[3];
756
758 XSH_ASSURE_NOT_NULL(trace);
759 XSH_ASSURE_NOT_NULL(filename);
760
761 /* create a table */
762 check(table = cpl_table_new(XSH_WAVESOL_TABLE_NB_COL+w->nbcoefs));
763 check(cpl_table_set_size(table,XSH_WAVESOL_TABLE_NB_ROWS));
764 /* create column names */
765 check(
766 cpl_table_new_column(table,XSH_WAVESOL_TABLE_COLNAME_AXIS,
767 CPL_TYPE_STRING));
768 check(
769 cpl_table_new_column(table,XSH_WAVESOL_TABLE_COLNAME_DEGSLIT,
770 CPL_TYPE_INT));
771 check(
772 cpl_table_new_column(table,XSH_WAVESOL_TABLE_COLNAME_DEGORDER,
773 CPL_TYPE_INT));
774 check(
775 cpl_table_new_column(table,XSH_WAVESOL_TABLE_COLNAME_DEGLAMBDA,
776 CPL_TYPE_INT));
777 check(cpl_table_set_string(table,XSH_WAVESOL_TABLE_COLNAME_AXIS,
778 0,"X"));
779 check(cpl_table_set_string(table,XSH_WAVESOL_TABLE_COLNAME_AXIS,
780 1,"Y"));
781 check(cpl_table_set_int(table,XSH_WAVESOL_TABLE_COLNAME_DEGSLIT,
782 0,w->deg_slit));
783 check(cpl_table_set_int(table,XSH_WAVESOL_TABLE_COLNAME_DEGSLIT,
784 1, w->deg_slit));
785 check(cpl_table_set_int(table,XSH_WAVESOL_TABLE_COLNAME_DEGORDER,
786 0,w->deg_order));
787 check(cpl_table_set_int(table,XSH_WAVESOL_TABLE_COLNAME_DEGORDER,
788 1, w->deg_order));
789 check(cpl_table_set_int(table,XSH_WAVESOL_TABLE_COLNAME_DEGLAMBDA,
790 0,w->deg_lambda));
791 check(cpl_table_set_int(table,XSH_WAVESOL_TABLE_COLNAME_DEGLAMBDA,
792 1, w->deg_lambda));
793
794 for(i=0; i<= w->deg_slit; i++){
795 for(j = 0; j<= w->deg_order; j++){
796 for(k = 0; k <= w->deg_lambda; k++){
797 double vx = 0.0, vy = 0.0;
798 pows[0] = i;
799 pows[1] = j;
800 pows[2] = k;
801 sprintf(coefname,"C%d%d%d",i,j,k);
802 check(cpl_table_new_column(table,coefname, CPL_TYPE_DOUBLE));
803 check(vx = cpl_polynomial_get_coeff(w->polx, pows));
804 check(vy = cpl_polynomial_get_coeff(w->poly, pows));
805 check(cpl_table_set_double( table, coefname, 0, vx));
806 check(cpl_table_set_double( table, coefname, 1, vy));
807 }
808 }
809 }
810
811 /* write min max in header */
822 check( xsh_pfits_set_pcatg( w->header, tag));
823 check(cpl_table_save(table, w->header, NULL, filename, CPL_IO_DEFAULT));
824 check(cpl_table_save(trace, NULL, NULL, filename, CPL_IO_EXTEND));
825
826 /* Create the frame */
827 check(result=xsh_frame_product(filename,
828 tag,
829 CPL_FRAME_TYPE_TABLE,
830 CPL_FRAME_GROUP_PRODUCT,
831 CPL_FRAME_LEVEL_TEMPORARY));
832
833 cleanup:
834 XSH_TABLE_FREE( table);
835 return result;
836}
837
838
839/*---------------------------------------------------------------------------*/
847/*---------------------------------------------------------------------------*/
848
849xsh_wavesol * xsh_wavesol_load( cpl_frame * frame,
851{
852 cpl_table* table = NULL;
853 const char* tablename = NULL;
854 xsh_wavesol * result = NULL ;
855 int rows, cols ;
856
857 XSH_ASSURE_NOT_NULL(frame);
858
859 /* get table filename */
860 check(tablename = cpl_frame_get_filename(frame));
861
862 XSH_TABLE_LOAD( table, tablename);
863
864 XSH_CALLOC( result, xsh_wavesol, 1 ) ;
865
866 check(result->header = cpl_propertylist_load(tablename,0));
867
868 check( rows = cpl_table_get_nrow( table ) ) ;
869 check( cols = cpl_table_get_ncol( table ) ) ;
870 xsh_msg_dbg_high( " === Wavesol Rows = %d, Columns = %d", rows, cols ) ;
871
872 /* For now fix the binning to 1 */
873 result->bin_x = 1;
874 result->bin_y = 1;
875
877 result->header));
879 result->header));
881 result->header));
883 result->header));
885 result->header));
887 result->header));
889 result->header));
891 result->header));
893 result->header));
895 result->header));
896
897 /* Now populate the structure from the table */
898 {
899 int i, j, k ;
900
901 /* Get degrees */
903 CPL_TYPE_INT, 0, &result->deg_slit ) ) ;
905 CPL_TYPE_INT, 0, &result->deg_order ) ) ;
907 CPL_TYPE_INT, 0, &result->deg_lambda ) ) ;
908 /* Create polynomials */
909 result->nbcoefs = (result->deg_order+1)*(result->deg_slit+1)*
910 (result->deg_lambda+1);
911 xsh_msg_dbg_high( "nbcoef %d deg_lambda %d deg_n %d deg_s %d", result->nbcoefs,
912 result->deg_lambda, result->deg_order, result->deg_slit);
913 result->polx = cpl_polynomial_new(3);
914 result->poly = cpl_polynomial_new(3);
915 result->dim = cpl_vector_new(3);
916
917 /* Get and set coeffs */
918 for( i = 0 ; i<=result->deg_slit ; i++ ) {
919 for( j = 0 ; j<=result->deg_order ; j++ ) {
920 for( k = 0 ; k<=result->deg_lambda ; k++ ) {
921 char coefname[16] ;
922 double coef ;
923 cpl_size pows[3] ;
924
925 pows[0] = i ;
926 pows[1] = j ;
927 pows[2] = k ;
928 sprintf( coefname, "C%d%d%d", i, j, k ) ;
929 /* Get coeffs */
930 check( xsh_get_table_value( table, coefname, CPL_TYPE_DOUBLE,
931 0, &coef ) ) ;
932 /* Set at the right place */
933 check( cpl_polynomial_set_coeff( result->polx, pows, coef ) ) ;
934 check( xsh_get_table_value( table, coefname, CPL_TYPE_DOUBLE,
935 1, &coef ) ) ;
936 /* Set at the right place */
937 check( cpl_polynomial_set_coeff( result->poly, pows, coef ) ) ;
938 }
939 }
940 }
941 }
942
943 /* Set binning */
946
947 cleanup:
948 if (cpl_error_get_code () != CPL_ERROR_NONE) {
949 xsh_error_msg("can't load Wavesol frame %s", tablename );
950 xsh_wavesol_free(&result);
951 }
952 XSH_TABLE_FREE( table) ;
953
954 return result ;
955}
956
957void xsh_wavesol_dump( xsh_wavesol * wsol, const char * fname, int nb )
958{
959 FILE * fout = NULL ;
960 int i, j, k, l ;
961 cpl_size pows[3] ;
962
963 if ( fname != NULL ) fout = fopen( fname, "w" ) ;
964 l = 0 ;
965
966 for( i = 0 ; i <= wsol->deg_slit ; i++ )
967 for( j = 0 ; j <= wsol->deg_lambda ; j++ )
968 for( k = 0 ; k <= wsol->deg_order ; k++ ) {
969 double coeff ;
970 pows[0] = i ;
971 pows[1] = j;
972 pows[2] = k;
973 check(coeff = cpl_polynomial_get_coeff( wsol->poly, pows));
974 if ( fout == NULL )
975 xsh_msg( " %d%d%d; %lf", i, j, k, coeff ) ;
976 else
977 fprintf( fout, "%d%d%d: %lf\n", i, j, k, coeff ) ;
978 l++ ;
979 if ( nb != 0 && l >= nb ) goto cleanup ;
980 }
981 cleanup:
982 if ( fout != NULL ) fclose( fout ) ;
983}
984
985cpl_table*
986xsh_wavesol_trace( xsh_wavesol * wsol, double* lambda,
987 double* order, double* slit,int size)
988{
989 cpl_table* table = NULL ;
990 int i;
991 double* pw=NULL;
992 double* po=NULL;
993 double* px=NULL;
994 double* py=NULL;
995 double* ps=NULL;
996
997
998 /* Check parameters */
999 XSH_ASSURE_NOT_NULL( wsol);
1000 XSH_ASSURE_NOT_NULL( lambda);
1002 XSH_ASSURE_NOT_NULL( slit);
1003
1004
1005 table=cpl_table_new(size);
1006 cpl_table_new_column(table,"WAVELENGTH",CPL_TYPE_DOUBLE);
1007 cpl_table_new_column(table,"ORDER",CPL_TYPE_DOUBLE);
1008 cpl_table_new_column(table,"X",CPL_TYPE_DOUBLE);
1009 cpl_table_new_column(table,"Y",CPL_TYPE_DOUBLE);
1010 cpl_table_new_column(table,"S",CPL_TYPE_DOUBLE);
1011
1012 cpl_table_fill_column_window(table,"WAVELENGTH",0,size,0);
1013 cpl_table_fill_column_window(table,"ORDER",0,size,0);
1014 cpl_table_fill_column_window(table,"X",0,size,0);
1015 cpl_table_fill_column_window(table,"Y",0,size,0);
1016 cpl_table_fill_column_window(table,"S",0,size,0);
1017
1018 po=cpl_table_get_data_double(table,"ORDER");
1019 pw=cpl_table_get_data_double(table,"WAVELENGTH");
1020 px=cpl_table_get_data_double(table,"X");
1021 py=cpl_table_get_data_double(table,"Y");
1022 ps=cpl_table_get_data_double(table,"S");
1023
1024 for(i = 0 ; i<size; i++) {
1025 pw[i] = lambda[i];
1026 po[i] = order[i];
1027 ps[i] = slit[i];
1028 check( px[i] = xsh_wavesol_eval_polx(wsol,pw[i],po[i],ps[i]));
1029 check( py[i] = xsh_wavesol_eval_poly(wsol,pw[i],po[i],ps[i]));
1030 }
1031
1032 cleanup:
1033 return table;
1034}
1035
1044void xsh_wavesol_set_bin_x( xsh_wavesol * wsol, int bin )
1045{
1046 XSH_ASSURE_NOT_NULL( wsol);
1047 wsol->bin_x = bin;
1048
1049 cleanup:
1050 return;
1051}
1052
1061void xsh_wavesol_set_bin_y( xsh_wavesol * wsol, int bin )
1062{
1063 XSH_ASSURE_NOT_NULL( wsol);
1064 wsol->bin_y = bin;
1065
1066 cleanup:
1067 return;
1068}
1069
1080void xsh_wavesol_apply_shift( xsh_wavesol *wsol, float shift_x, float shift_y)
1081{
1082
1083 XSH_ASSURE_NOT_NULL( wsol);
1084
1085 wsol->min_x = wsol->min_x+shift_x;
1086 wsol->max_x = wsol->max_x+shift_x;
1087 wsol->min_y = wsol->min_y+shift_y;
1088 wsol->max_y = wsol->max_y+shift_y;
1089
1090 cleanup:
1091 return;
1092}
static xsh_instrument * instrument
xsh_spectralformat_list * xsh_spectralformat_list_load(cpl_frame *frame, xsh_instrument *instr)
Load a spectralformat list from a frame.
void xsh_spectralformat_list_free(xsh_spectralformat_list **list)
Free memory associated to an spactralformat_list.
void xsh_wavesol_apply_shift(xsh_wavesol *wsol, float shift_x, float shift_y)
Apply a shift on X and Y to wave solution.
void xsh_wavesol_set_bin_y(xsh_wavesol *wsol, int bin)
Set the bin of wave table in y.
void xsh_wavesol_set_bin_x(xsh_wavesol *wsol, int bin)
Set the bin of wave table in x.
cpl_table * xsh_wavesol_trace(xsh_wavesol *wsol, double *lambda, double *order, double *slit, int size)
void xsh_wavesol_compute(xsh_wavesol *sol, int size, double *pos, double *posmin, double *posmax, double *lambda, double *order, double *slit, cpl_polynomial *result)
compute a wavelength solution
cpl_polynomial * xsh_wavesol_get_polx(xsh_wavesol *sol)
get the solution in X
void xsh_wavesol_dump(xsh_wavesol *wsol, const char *fname, int nb)
xsh_wavesol * xsh_wavesol_load(cpl_frame *frame, xsh_instrument *instrument)
load a wavelength solution
double xsh_wavesol_eval_poly(xsh_wavesol *sol, double lambda, double order, double slit)
eval the polynomial solution in Y
void xsh_wavesol_residual(xsh_wavesol *sol, xsh_wavesol *adj, int size, double *new_pos, double *lambda, double *order, double *slit, cpl_polynomial *result, char axis)
xsh_wavesol * xsh_wavesol_create(cpl_frame *spectral_format_frame, xsh_detect_arclines_param *params, xsh_instrument *instrument)
Create a new wavelength solution structure.
cpl_polynomial * xsh_wavesol_get_poly(xsh_wavesol *sol)
get the solution in Y
void xsh_wavesol_free(xsh_wavesol **w)
free wavelength solution structure
double xsh_wavesol_eval_polx(xsh_wavesol *sol, double lambda, double order, double slit)
eval the polynomial solution in X
enum wavesol_type xsh_wavesol_get_type(xsh_wavesol *wsol)
get the type of the wave table
void xsh_wavesol_set_type(xsh_wavesol *wsol, enum wavesol_type type)
set the type of the wave table
cpl_propertylist * xsh_wavesol_get_header(xsh_wavesol *sol)
get header of the table
xsh_wavesol * xsh_wavesol_duplicate(xsh_wavesol *org)
duplicate a wavelength solution structure
void xsh_wavesol_add_poly(xsh_wavesol *to, xsh_wavesol *from)
cpl_frame * xsh_wavesol_save(xsh_wavesol *w, cpl_table *trace, const char *filename, const char *tag)
save a wavelength solution
#define XSH_ASSURE_NOT_ILLEGAL(cond)
Definition: xsh_error.h:107
#define check(COMMAND)
Definition: xsh_error.h:71
#define XSH_ASSURE_NOT_ILLEGAL_MSG(cond, msg)
Definition: xsh_error.h:111
#define xsh_error_msg(...)
Definition: xsh_error.h:94
#define XSH_ASSURE_NOT_NULL(pointer)
Definition: xsh_error.h:99
int xsh_instrument_get_binx(xsh_instrument *instrument)
int xsh_instrument_get_biny(xsh_instrument *instrument)
int size
int * y
#define xsh_msg(...)
Print a message on info level.
Definition: xsh_msg.h:121
#define xsh_msg_dbg_high(...)
Definition: xsh_msg.h:40
double xsh_pfits_get_wavesol_lambda_min(cpl_propertylist *plist)
find out the min lambda
Definition: xsh_pfits.c:2919
double xsh_pfits_get_wavesol_lambda_max(cpl_propertylist *plist)
find out the wavesol max lambda
Definition: xsh_pfits.c:2938
void xsh_pfits_set_pcatg(cpl_propertylist *plist, const char *value)
Write the PCATG value.
Definition: xsh_pfits.c:1008
void xsh_pfits_set_wavesol_x_min(cpl_propertylist *plist, double value)
WRITE the min x.
Definition: xsh_pfits.c:2850
double xsh_pfits_get_wavesol_order_min(cpl_propertylist *plist)
find out the min order
Definition: xsh_pfits.c:2957
double xsh_pfits_get_wavesol_y_max(cpl_propertylist *plist)
find out the wavesol max y position
Definition: xsh_pfits.c:3128
double xsh_pfits_get_wavesol_slit_min(cpl_propertylist *plist)
find out the min slit
Definition: xsh_pfits.c:2995
void xsh_pfits_set_wavesol_y_min(cpl_propertylist *plist, double value)
WRITE the min y.
Definition: xsh_pfits.c:2885
double xsh_pfits_get_wavesol_x_max(cpl_propertylist *plist)
find out the wavesol max x position
Definition: xsh_pfits.c:3090
void xsh_pfits_set_wavesol_y_max(cpl_propertylist *plist, double value)
WRITE the max y position.
Definition: xsh_pfits.c:2902
double xsh_pfits_get_wavesol_y_min(cpl_propertylist *plist)
find out the min y position
Definition: xsh_pfits.c:3109
void xsh_pfits_set_wavesol_slit_min(cpl_propertylist *plist, double value)
WRITE the min slit.
Definition: xsh_pfits.c:2781
void xsh_pfits_set_wavesol_lambda_max(cpl_propertylist *plist, double value)
WRITE the max lambda.
Definition: xsh_pfits.c:2730
void xsh_pfits_set_wavesol_slit_max(cpl_propertylist *plist, double value)
WRITE the max slit.
Definition: xsh_pfits.c:2798
void xsh_pfits_set_wavesol_order_max(cpl_propertylist *plist, double value)
WRITE the max order.
Definition: xsh_pfits.c:2764
double xsh_pfits_get_wavesol_x_min(cpl_propertylist *plist)
find out the min x position
Definition: xsh_pfits.c:3071
double xsh_pfits_get_wavesol_slit_max(cpl_propertylist *plist)
find out the wavesol max slit
Definition: xsh_pfits.c:3014
void xsh_pfits_set_wavesol_lambda_min(cpl_propertylist *plist, double value)
WRITE the min lambda.
Definition: xsh_pfits.c:2713
void xsh_pfits_set_wavesol_x_max(cpl_propertylist *plist, double value)
WRITE the max x position.
Definition: xsh_pfits.c:2867
double xsh_pfits_get_wavesol_order_max(cpl_propertylist *plist)
find out the wavesol max order
Definition: xsh_pfits.c:2976
void xsh_pfits_set_wavesol_order_min(cpl_propertylist *plist, double value)
WRITE the min order.
Definition: xsh_pfits.c:2747
void xsh_free_polynomial(cpl_polynomial **p)
Deallocate a polynomial and set the pointer to NULL.
Definition: xsh_utils.c:2194
void xsh_free_vector(cpl_vector **v)
Deallocate a vector and set the pointer to NULL.
Definition: xsh_utils.c:2284
void xsh_tools_min_max(int size, double *tab, double *min, double *max)
computes min & max in ab array
Definition: xsh_utils.c:2799
double convert_data_to_bin(double data, int binning)
Definition: xsh_utils.c:3246
double xsh_tools_tchebitchev_transform(double pos, double min, double max)
computes Tchebitchev transformation
Definition: xsh_utils.c:2917
cpl_vector * xsh_tools_tchebitchev_poly_eval(int n, double X)
Compute tchebitchev Tn(X) first coefficient for tchebitchev polynomial.
Definition: xsh_utils.c:2836
double xsh_tools_tchebitchev_reverse_transform(double pos, double min, double max)
computes reverse Tchebitchev transformation
Definition: xsh_utils.c:2954
void xsh_free_propertylist(cpl_propertylist **p)
Deallocate a property list and set the pointer to NULL.
Definition: xsh_utils.c:2179
cpl_error_code xsh_get_table_value(const cpl_table *table, const char *colname, cpl_type coltype, int i, void *result)
Read a table value from a fits table.
void xsh_tools_tchebitchev_transform_tab(int size, double *pos, double min, double max, double *tcheb_pos)
computes Tchebitchev transformation
Definition: xsh_utils.c:2881
cpl_polynomial * polx
cpl_propertylist * header
cpl_vector * dim
enum wavesol_type type
cpl_polynomial * poly
#define XSH_WAVESOL_TABLE_COLNAME_DEGSLIT
#define XSH_WAVESOL_TABLE_COLNAME_DEGLAMBDA
#define XSH_WAVESOL_TABLE_NB_COL
#define XSH_WAVESOL_TABLE_COLNAME_AXIS
wavesol_type
@ XSH_WAVESOL_GUESS
#define XSH_SLIT_RANGE
#define XSH_WAVESOL_TABLE_NB_ROWS
#define XSH_WAVESOL_TABLE_COLNAME_DEGORDER
int order
Definition: xsh_detmon_lg.c:80
cpl_frame * xsh_frame_product(const char *fname, const char *tag, cpl_frame_type type, cpl_frame_group group, cpl_frame_level level)
Creates a frame with given characteristics.
Definition: xsh_dfs.c:930
#define XSH_FREE(POINTER)
Definition: xsh_utils.h:92
#define XSH_MALLOC(POINTER, TYPE, SIZE)
Definition: xsh_utils.h:49
#define XSH_CALLOC(POINTER, TYPE, SIZE)
Definition: xsh_utils.h:56
#define XSH_TABLE_LOAD(TABLE, NAME)
#define XSH_TABLE_FREE(TABLE)