X-shooter Pipeline Reference Manual 3.8.15
xsh_data_wavemap.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: 2012-12-01 18:59:24 $
23 * $Revision: 1.57 $
24 */
25
26#ifdef HAVE_CONFIG_H
27#include <config.h>
28#endif
29#define REGDEBUG_BSPLINE 0
30/*---------------------------------------------------------------------------*/
35/*---------------------------------------------------------------------------*/
36
40/*-----------------------------------------------------------------------------
41 Includes
42 ----------------------------------------------------------------------------*/
43#include <string.h>
44#include <math.h>
45#include <xsh_utils_table.h>
46#include <xsh_parameters.h>
47#include <xsh_data_pre.h>
48#include <xsh_data_order.h>
49#include <xsh_data_dispersol.h>
50#include <xsh_utils_wrappers.h>
51#include <xsh_data_wavemap.h>
52#include <xsh_error.h>
53#include <xsh_msg.h>
54#include <xsh_pfits.h>
55#include <xsh_dfs.h>
56#include <math.h>
57#include <gsl/gsl_multifit.h>
58#include <cpl.h>
59
60/*----------------------------------------------------------------------------
61 Function implementation
62 ----------------------------------------------------------------------------*/
69void
71 const char * fname )
72{
73 int i ;
74 FILE * fout = NULL;
75
76 XSH_ASSURE_NOT_NULL( list ) ;
77 if ( fname == NULL ) fout = stdout ;
78 else fout = fopen( fname, "w" ) ;
79 XSH_ASSURE_NOT_NULL( fout ) ;
80
81 fprintf( fout, "Wavemap List. Nb of orders: %d\n", list->size ) ;
82 for( i = 0 ; i<list->size ; i++ ) {
83 fprintf( fout, " Entry %2d: Order %d, Ndata: %d\n",
84 i, list->list[i].order, list->list[i].sky_size);
85 }
86
87 cleanup:
88 if ( fname != NULL && fout != NULL ) fclose( fout ) ;
89 return ;
90}
91
92
93/*---------------------------------------------------------------------------*/
100/*---------------------------------------------------------------------------*/
103{
104 xsh_wavemap_list* result = NULL;
105 XSH_INSTRCONFIG * config = NULL;
106 //XSH_ARM xsh_arm;
107 int i, size = 0;
108
109 /* check input parameters */
110 XSH_ASSURE_NOT_NULL( instr);
111
112 /* Get arm */
113 //xsh_arm = xsh_instrument_get_arm( instr) ;
114 check( config = xsh_instrument_get_config( instr));
115 size = config->orders;
116
117 XSH_CALLOC(result, xsh_wavemap_list, 1 );
118
119 result->size = size;
120
121 XSH_ASSURE_NOT_ILLEGAL(result->size > 0);
122 result->instrument = instr;
123 XSH_CALLOC(result->list, xsh_wavemap, result->size);
125
126 /* Create polynomials for each order */
127 for( i = 0 ; i<result->size ; i++ ) {
128 result->list[i].tcheb_pol_lambda = NULL ;
129 }
130
131 cleanup:
132 if (cpl_error_get_code() != CPL_ERROR_NONE) {
133 xsh_wavemap_list_free(&result);
134 }
135 return result;
136}
137/*---------------------------------------------------------------------------*/
147/*---------------------------------------------------------------------------*/
148void
150 int absorder, int max_size)
151{
152 xsh_wavemap * pwavemap = NULL ;
153
154 /* Check input parameters */
155 XSH_ASSURE_NOT_NULL(list) ;
156 XSH_ASSURE_NOT_ILLEGAL( idx >= 0 && idx < list->size && max_size > 0);
157 pwavemap = &list->list[idx] ;
158 XSH_ASSURE_NOT_NULL( pwavemap);
159
160 pwavemap->order = absorder;
161 pwavemap->sky_size = max_size;
162 pwavemap->object_size = max_size;
163 pwavemap->all_size = max_size;
164 XSH_CALLOC( pwavemap->sky, wavemap_item, max_size);
165 XSH_CALLOC( pwavemap->object, wavemap_item, max_size);
166 XSH_CALLOC( pwavemap->all, wavemap_item, max_size);
167
168 cleanup:
169 return ;
170}
171
172/*---------------------------------------------------------------------------*/
177/*---------------------------------------------------------------------------*/
178void
180{
181 int i = 0;
182 xsh_wavemap_list *plist = NULL;
183
184 if (list != NULL && *list != NULL ) {
185 plist = *list ;
186 /* free the list */
187 for (i = 0; i < plist->size; i++) {
188 xsh_wavemap *pr = &(plist->list[i]) ;
189
190 xsh_msg_dbg_high( "Freeing order index %d", i ) ;
191 if ( pr == NULL ) continue ;
192 xsh_msg_dbg_high( " Abs Order: %d", pr->order ) ;
193 cpl_free( pr->sky);
194 cpl_free( pr->object);
195 cpl_free( pr->all);
196 if (pr->pol_lambda != NULL){
198 }
199 if (pr->pol_slit != NULL){
201 }
203 }
204 if ((*list)->list){
205 cpl_free((*list)->list);
206 }
207 xsh_free_propertylist(&((*list)->header));
208 cpl_free(*list);
209 *list = NULL;
210 }
211}
212
213static void
214lambda_fit( double * lambda, double * xpos, double * ypos,
215 int size, xsh_wavemap_list * wmap,
216 xsh_dispersol_param * dispsol_param, int order)
217{
218 double xmin, xmax, ymin, ymax, lmin, lmax;
219 double *tcheb_lambda = NULL, *tcheb_xpos = NULL, *tcheb_ypos = NULL;
220 int nbcoefs, xdeg, ydeg;
221 cpl_vector *tcheb_coef_x = NULL;
222 cpl_vector *tcheb_coef_y = NULL;
223 int i, j, k ;
224 double chisq;
225 gsl_matrix *X = NULL, *cov = NULL;
226 gsl_vector *y = NULL, *c = NULL;
227 gsl_multifit_linear_workspace *work = NULL;
228 cpl_polynomial* result = NULL;
229 cpl_size pows[3];
230
231 XSH_ASSURE_NOT_NULL( lambda);
232 XSH_ASSURE_NOT_NULL( xpos);
233 XSH_ASSURE_NOT_NULL( ypos);
234 XSH_ASSURE_NOT_NULL( wmap);
235 XSH_ASSURE_NOT_NULL( dispsol_param);
237 XSH_ASSURE_NOT_ILLEGAL( order >= 0 && order < wmap->size);
238
239 xsh_msg_dbg_low( " Entering lambda_fit. Order %d, size: %d",
240 wmap->list[order].order, size);
241 xdeg = dispsol_param->deg_x;
242 ydeg = dispsol_param->deg_y;
243 nbcoefs = (xdeg+1)*(ydeg+1) ;
244
245 XSH_ASSURE_NOT_ILLEGAL( size >= nbcoefs);
246
247 wmap->degx = xdeg ;
248 wmap->degy = ydeg ;
249
250 /* Calculate min and max (necessary to transform into Tchebitchev space) */
251 check(xsh_tools_min_max(size, lambda, &lmin, &lmax));
252 check(xsh_tools_min_max(size, xpos, &xmin, &xmax));
253 check(xsh_tools_min_max(size, ypos, &ymin, &ymax));
254
255 wmap->list[order].xmin = xmin ;
256 wmap->list[order].xmax = xmax ;
257 wmap->list[order].ymin = ymin ;
258 wmap->list[order].ymax = ymax ;
259 wmap->list[order].lambda_min = lmin ;
260 wmap->list[order].lambda_max = lmax ;
261
263 " Min,Max. Lambda: %.6lf,%.6lf - X: %.6lf,%.6lf - Y: %.6lf,%.6lf",
264 lmin, lmax, xmin, xmax, ymin, ymax);
265
266 /* Transform to tchebitchev space */
267 XSH_CALLOC(tcheb_lambda, double, size);
268 XSH_CALLOC(tcheb_xpos, double, size);
269 XSH_CALLOC(tcheb_ypos, double, size);
270
272 tcheb_lambda));
274 tcheb_xpos));
276 tcheb_ypos));
277 X = gsl_matrix_alloc( size, nbcoefs);
278 y = gsl_vector_alloc( size);
279 c = gsl_vector_alloc( nbcoefs);
280 cov = gsl_matrix_alloc( nbcoefs, nbcoefs);
281
282 for (i = 0; i < size; i++) {
284 xdeg, tcheb_xpos[i]));
286 ydeg, tcheb_ypos[i]));
287
288 for(j = 0; j < dispsol_param->deg_x +1 ; j++) {
289 for(k = 0; k < dispsol_param->deg_y +1 ; k++){
290 double vx, vy;
291
292 check( vx = cpl_vector_get( tcheb_coef_x, j));
293 check( vy = cpl_vector_get( tcheb_coef_y, k));
294 gsl_matrix_set (X, i, k+j*(xdeg+1), vx*vy );
295 }
296 }
297 gsl_vector_set (y, i, tcheb_lambda[i]);
298 xsh_free_vector(&tcheb_coef_x);
299 xsh_free_vector(&tcheb_coef_y);
300
301 }
302 xsh_msg_dbg_medium( " GSL Initialized" ) ;
303
304
306 "Not enough Points vs Number of Coeffs" ) ;
307 work = gsl_multifit_linear_alloc(size, nbcoefs);
308 gsl_multifit_linear(X, y, c, cov, &chisq, work);
309
310 xsh_msg_dbg_medium( " GSL Done" ) ;
311
312 wmap->list[order].tcheb_pol_lambda = cpl_polynomial_new( 2);
313
314 result = wmap->list[order].tcheb_pol_lambda ;
315
316 for(i = 0; i < (xdeg +1) ; i++){
317 for(j = 0; j < (ydeg +1) ; j++){
318 pows[0] = j;
319 pows[1] = i;
320 cpl_polynomial_set_coeff(result, pows,
321 gsl_vector_get(c, j+i*(ydeg+1)));
322 }
323 }
324 xsh_msg_dbg_medium( " GSL Finished" ) ;
325
326 gsl_multifit_linear_free (work);
327 gsl_matrix_free (X);
328 gsl_vector_free (y);
329 gsl_vector_free (c);
330 gsl_matrix_free (cov);
331
332 cleanup:
333 XSH_FREE( tcheb_lambda);
334 XSH_FREE( tcheb_xpos);
335 XSH_FREE( tcheb_ypos);
336 return ;
337}
338
339/*---------------------------------------------------------------------------*/
355/*---------------------------------------------------------------------------*/
356
357
358
359void
361 double * vslit,
362 double * xpos,
363 double * ypos,
364 int nitems,
365 double * orders,
366 xsh_dispersol_param *dispsol_param,
367 xsh_wavemap_list *wmap)
368{
369 /*
370 * calculer par ordre ==> separer en fonction de l'ordre (vorderdata)
371 * pour chaque ordre, fitter lambda en fonction de X et Y (tchebitchev)
372 * */
373 int i, ordersize, nborder ;
374 double *vx = NULL, *vy = NULL, *vl = NULL, *vs = NULL;
375
376
377 xsh_msg( "Entering xsh_wavemap_compute" ) ;
378
379 XSH_ASSURE_NOT_NULL( vlambda);
380 XSH_ASSURE_NOT_NULL( vslit);
381 XSH_ASSURE_NOT_NULL( xpos);
382 XSH_ASSURE_NOT_NULL( ypos);
383 XSH_ASSURE_NOT_NULL( orders);
384 XSH_ASSURE_NOT_NULL( wmap);
385 XSH_ASSURE_NOT_NULL( dispsol_param);
386 XSH_ASSURE_NOT_ILLEGAL( nitems > 0);
388
389 xsh_msg( " X degree = %d, Y degree = %d", dispsol_param->deg_x,
390 dispsol_param->deg_y) ;
391
392 /* split by order */
393 ordersize = 0;
394 nborder = 0;
395
396 wmap->degx = dispsol_param->deg_x;
397 wmap->degy = dispsol_param->deg_y;
398
399 xsh_msg("Compute POLY for lambda");
400 XSH_REGDEBUG( "nitems %d ", nitems);
401
402 for( i = 1 ; i<=nitems ; i++ ) {
403
404 if ( i < nitems && orders[i-1] == orders[i] ) {
405 ordersize++;
406 }
407 else {
408 cpl_vector* posx = NULL;
409 cpl_vector* posy = NULL;
410 cpl_vector* lvalues = NULL;
411 cpl_vector* svalues = NULL;
412 cpl_bivector *positions = NULL;
413 double mse =0.0;
414
415 ordersize++;
416 int i_minus_ordersize=i-ordersize;
417 XSH_MALLOC( vx, double, ordersize);
418 memcpy( vx, &(xpos[i_minus_ordersize]), ordersize*sizeof(double));
419
420 XSH_MALLOC( vy, double, ordersize);
421 memcpy( vy, &(ypos[i_minus_ordersize]), ordersize*sizeof(double));
422
423 XSH_MALLOC( vl, double, ordersize);
424 memcpy( vl, &(vlambda[i_minus_ordersize]), ordersize*sizeof(double) ) ;
425
426 XSH_MALLOC( vs, double, ordersize);
427 memcpy( vs, &( vslit[i_minus_ordersize]), ordersize*sizeof(double) ) ;
428 wmap->list[nborder].sky_size = ordersize ;
429
430 /* Now fit lambda */
431 wmap->list[nborder].order = orders[i-1] ;
432 xsh_msg_dbg_high( "Order: %d, Size: %d", wmap->list[nborder].order,
433 ordersize );
434
435 posx = cpl_vector_wrap( ordersize, vx);
436 posy = cpl_vector_wrap( ordersize, vy);
437
438 positions = cpl_bivector_wrap_vectors( posx, posy);
439 lvalues = cpl_vector_wrap( ordersize, vl);
440 svalues = cpl_vector_wrap( ordersize, vs);
441 cpl_size loc_degx=(cpl_size)dispsol_param->deg_x;
442 wmap->list[nborder].pol_lambda =
443 xsh_polynomial_fit_2d_create( positions, lvalues, &loc_degx,
444 &mse);
445
446 wmap->list[nborder].pol_slit =
447 xsh_polynomial_fit_2d_create( positions, svalues, &loc_degx,
448 &mse);
449
450 cpl_bivector_unwrap_vectors( positions);
451
452 cpl_vector_unwrap( posx);
453 cpl_vector_unwrap( posy);
454 cpl_vector_unwrap( lvalues);
455 cpl_vector_unwrap( svalues);
456
457 XSH_FREE( vx);
458 XSH_FREE( vy);
459 XSH_FREE( vl);
460 XSH_FREE( vs);
461 nborder++ ;
462 ordersize = 0 ;
463
464 }
465 }
466
467 cleanup:
468 XSH_FREE( vx);
469 XSH_FREE( vy);
470 XSH_FREE( vl);
471 XSH_FREE( vs);
472 return ;
473}
474
475/*---------------------------------------------------------------------------*/
489/*---------------------------------------------------------------------------*/
490void
491xsh_wavemap_list_compute( double * vlambda,
492 double * xpos,
493 double * ypos,
494 int nitems,
495 double * orders,
496 xsh_dispersol_param * dispsol_param,
497 xsh_wavemap_list * wmap )
498{
499 int i, ordersize, nborder ;
500 double *vx = NULL, *vy = NULL, *vl = NULL;
501
502 xsh_msg( "Entering xsh_wavemap_compute" ) ;
503
504 XSH_ASSURE_NOT_NULL( vlambda ) ;
505 XSH_ASSURE_NOT_NULL( xpos ) ;
506 XSH_ASSURE_NOT_NULL( ypos ) ;
507 XSH_ASSURE_NOT_NULL( orders ) ;
508 XSH_ASSURE_NOT_NULL( wmap ) ;
509 XSH_ASSURE_NOT_NULL( dispsol_param ) ;
510 XSH_ASSURE_NOT_ILLEGAL( nitems > 0 ) ;
512
513 xsh_msg( " X degree = %d, Y degree = %d", dispsol_param->deg_x,
514 dispsol_param->deg_y) ;
515
516 /*
517 split data by order (according to orders array)
518 then fit a tchebitchev polynomial Lambda = T(X,Y) per order
519 */
520 ordersize = 0;
521 nborder = 0;
522
523 XSH_REGDEBUG( "nitems %d ", nitems);
524
525 for( i = 1 ; i<=nitems ; i++ ) {
526
527 if ( i < nitems && orders[i-1] == orders[i] ) {
528 ordersize++;
529 }
530 else {
531 ordersize++;
532 XSH_MALLOC( vx, double, ordersize);
533 memcpy( vx, &(xpos[i-ordersize]), ordersize*sizeof(double));
534
535 XSH_MALLOC( vy, double, ordersize);
536 memcpy( vy, &(ypos[i-ordersize]), ordersize*sizeof(double));
537
538 XSH_MALLOC( vl, double, ordersize);
539 memcpy( vl, &(vlambda[i-ordersize]), ordersize*sizeof(double) ) ;
540
541 wmap->list[nborder].sky_size = ordersize;
542 /* Now fit lambda */
543 wmap->list[nborder].order = orders[i-1] ;
544 xsh_msg_dbg_high( "Order: %d, Size: %d", wmap->list[nborder].order,
545 ordersize );
546
547 check( lambda_fit( vl, vx, vy, ordersize, wmap, dispsol_param, nborder ) ) ;
548 XSH_FREE( vx ) ;
549 XSH_FREE( vy ) ;
550 XSH_FREE( vl ) ;
551 nborder++ ;
552 ordersize = 0 ;
553 }
554 }
555
556 cleanup:
557 XSH_FREE( vx ) ;
558 XSH_FREE( vy ) ;
559 XSH_FREE( vl ) ;
560 return ;
561}
562
563/*---------------------------------------------------------------------------*/
573/*---------------------------------------------------------------------------*/
574
575static double
577 double x, double y, int ord )
578{
579 double tcheb_x = 0.0, tcheb_y = 0.0;
580 double lambda = -1.;
581 double res=0.0;
582 int i, j ;
583 cpl_size pows[2] ;
584 int degx, degy ;
585 cpl_vector *tcheb_coef_x = NULL;
586 cpl_vector *tcheb_coef_y = NULL;
587
588 XSH_ASSURE_NOT_NULL( wmap);
589 degx = wmap->degx;
590 degy = wmap->degy;
591 XSH_ASSURE_NOT_ILLEGAL( degx >0 && degy > 0);
592
593 xsh_msg_dbg_high( " eval_lambda: degx: %d, min: %lf, max: %lf",
594 degx, wmap->list[ord].xmin,
595 wmap->list[ord].xmax ) ;
596 // Bug fixing
597 if ( wmap->list[ord].xmin >= wmap->list[ord].xmax ||
598 wmap->list[ord].lambda_min >= wmap->list[ord].lambda_max ||
599 wmap->list[ord].ymin >= wmap->list[ord].ymax ) {
600 return 0. ;
601 }
602
603 check( tcheb_x = xsh_tools_tchebitchev_transform( x, wmap->list[ord].xmin,
604 wmap->list[ord].xmax));
606 wmap->list[ord].ymin, wmap->list[ord].ymax));
608 degx, tcheb_x));
610 degy, tcheb_y));
611
612 for ( i=0; i < degx+1; i++) {
613 for( j=0; j < degy+1; j++) {
614 double vx, vy;
615 double coef;
616
617 check( vx = cpl_vector_get( tcheb_coef_x, i));
618 check( vy = cpl_vector_get( tcheb_coef_y, j));
619 pows[0] = j;
620 pows[1] = i;
621 check( coef = cpl_polynomial_get_coeff( wmap->list[ord].tcheb_pol_lambda,
622 pows));
623 res += coef*vx*vy ;
624 }
625 }
627 wmap->list[ord].lambda_min, wmap->list[ord].lambda_max));
628
629 cleanup:
630 xsh_free_vector( &tcheb_coef_x);
631 xsh_free_vector( &tcheb_coef_y);
632 return lambda ;
633}
634
647cpl_frame *
649 cpl_frame * order_frame,
650 xsh_pre * pre,
651 xsh_instrument * instr,
652 const char * prefix)
653{
654
655 cpl_frame * result = NULL ;
656 int order, x, y ;
657 int start_y, end_y ;
658 int xmin, xmax ;
659 int nx, ny ;
660 double lambda ;
661 cpl_image * wimg = NULL ;
662 double * dimg = NULL ;
663 cpl_propertylist * wheader = NULL ;
664 char * final_name = NULL;
665 xsh_order_list * order_list = NULL ;
666
667 XSH_ASSURE_NOT_NULL( wmap ) ;
668 XSH_ASSURE_NOT_NULL( order_frame ) ;
669 XSH_ASSURE_NOT_NULL( pre ) ;
670 XSH_ASSURE_NOT_NULL( prefix ) ;
671 XSH_ASSURE_NOT_NULL( instr ) ;
672
673 final_name = xsh_stringcat_any( prefix, ".fits", (void*)NULL ) ;
674 xsh_msg( "Entering xsh_wavemap_save, file \"%s\"", final_name ) ;
675 /*
676 create the wave map image filled with zeroes.
677 Fill the orders with lambda calculated from wave map polynomial
678 */
679 /*
680 Load the order table into the order_list
681 */
682 check(order_list=xsh_order_list_load(order_frame, instr )) ;
683 /*
684 Create a new image according to parameters included in instrument
685 structure
686 */
687
688 nx = pre->nx ;
689 ny = pre->ny ;
690
691 xsh_msg/*_dbg_medium*/( "Image size:%d,%d", nx, ny ) ;
692
693 check( wimg = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE ) ) ;
694 check( dimg = cpl_image_get_data_double( wimg ) ) ;
695
696 for( order = 0 ; order < wmap->size ; order++ ) {
697 int count = 0 ;
698
699 // Fix
700 if ( wmap->list[order].tcheb_pol_lambda == NULL ) {
701 xsh_msg( "Order %d: NULL Polynome", order ) ;
702 continue ;
703 }
704 if ( wmap->list[order].sky_size <= 0 ) {
705 xsh_msg( "NOT ENOUGH DATA FOR ORDER %d",
706 order_list->list[order].absorder ) ;
707 continue ;
708 }
709 /* get start_y, end_y from order list */
710 start_y = xsh_order_list_get_starty( order_list, order ) ;
711 end_y = xsh_order_list_get_endy( order_list, order ) ;
712 xsh_msg_dbg_low( " Order %d, Ymin: %d, Ymax: %d",
713 order, start_y, end_y ) ;
714
715 for ( y = start_y ; y < end_y ; y++ ) {
716 double dxmin, dxmax ;
717
718 /* get xmin, xmax */
719 check( dxmin = cpl_polynomial_eval_1d(order_list->list[order].edglopoly,
720 (double)y, NULL ) ) ;
721 check( dxmax = cpl_polynomial_eval_1d(order_list->list[order].edguppoly,
722 (double)y, NULL ) ) ;
723 xmin = floor(dxmin) ;
724 xmax = ceil(dxmax) ;
725 for (x = xmin ; x<xmax ; x++ ) {
726 count++ ;
727 /*calculate lambda(x,y)*/
728 lambda = xsh_wavemap_list_eval_lambda( wmap, (double)x,
729 (double)y, order);
730 if ( cpl_error_get_code() != CPL_ERROR_NONE){
731 lambda =0.0;
733 }
734 /*add lambda to image */
735 *(dimg + x + (y*nx)) = (float)lambda ;
736 }
737 }
738 xsh_msg_dbg_high( " %d points to calculate", count ) ;
739 }
740
741 /* Now save image */
742 check( wheader = cpl_propertylist_duplicate( pre->data_header ) ) ;
743 check( cpl_image_save( wimg, final_name, CPL_BPP_IEEE_FLOAT,
744 wheader, CPL_IO_DEFAULT ) ) ;
745 check(result=xsh_frame_product(final_name,
746 "TAG",
747 CPL_FRAME_TYPE_IMAGE,
748 CPL_FRAME_GROUP_PRODUCT,
749 CPL_FRAME_LEVEL_TEMPORARY));
750
751 cleanup:
752 xsh_order_list_free( &order_list);
753 XSH_FREE( final_name ) ;
754 xsh_free_image ( &wimg ) ;
755 xsh_free_propertylist( &wheader);
756
757 return result ;
758}
759
760
761
762cpl_error_code
764 xsh_instrument * instr,
765 const int abs_ord, const int sid,const int iter)
766{
767
768 xsh_pre* sky_model = NULL ;
769 xsh_pre* fit_model = NULL ;
770 cpl_image* sky_wmap = NULL;
771 cpl_image* sky_smap = NULL;
772
773 float* pdata_sky =NULL;
774 float* perrs_sky =NULL;
775 int* pqual_sky =NULL;
776 float* pdata_fit =NULL;
777 float* perrs_fit =NULL;
778 int* pqual_fit =NULL;
779
780 float* psky_wmap =NULL;
781 float* psky_smap =NULL;
782
783 int order;
784 int sx, sy ;
785 //char fname[80];
786 char * final_name = NULL;
787
788 XSH_ASSURE_NOT_NULL( smap ) ;
789 XSH_ASSURE_NOT_NULL( instr ) ;
790
791 sx = instr->config->nx /instr->binx;
792 sy = instr->config->ny /instr->biny;
793 xsh_msg( "Image size:%d,%d", sx, sy ) ;
794
795 /* create a model of the object and the sky as a PRE format frame */
796 sky_model = xsh_pre_new(sx,sy);
797 fit_model = xsh_pre_new(sx,sy);
798
799 pdata_sky = cpl_image_get_data_float(xsh_pre_get_data(sky_model));
800 perrs_sky = cpl_image_get_data_float(xsh_pre_get_errs(sky_model));
801 pqual_sky = cpl_image_get_data_int(xsh_pre_get_qual(sky_model));
802
803 pdata_fit = cpl_image_get_data_float(xsh_pre_get_data(fit_model));
804 perrs_fit = cpl_image_get_data_float(xsh_pre_get_errs(fit_model));
805 pqual_fit = cpl_image_get_data_int(xsh_pre_get_qual(fit_model));
806
807
808 /* associate to object and sky the corresponding wave and slit maps */
809 sky_wmap = cpl_image_new(sx,sy,XSH_PRE_DATA_TYPE);
810 sky_smap = cpl_image_new(sx,sy,XSH_PRE_DATA_TYPE);
811
812 psky_wmap = cpl_image_get_data_float(sky_wmap);
813 psky_smap = cpl_image_get_data_float(sky_smap);
814
815 int ix, iy,pix;
816 for( order = 0 ; order < smap->size ; order++ ) {
817
818 wavemap_item * psky = smap->list[order].sky;
819 int sky_size = smap->list[order].sky_size;
820
821 for (int k = 0; k < sky_size; k++) {
822 ix = psky->ix;
823 iy = psky->iy;
824 pix = iy*sx+ix;
825 pdata_sky[pix]=psky->flux;
826 perrs_sky[pix]=psky->sigma;
827 pqual_sky[pix]=psky->qual;
828 pdata_fit[pix]=psky->fitted;
829 perrs_fit[pix]=psky->fit_err;
830 pqual_fit[pix]=psky->qual;
831
832 psky_wmap[pix]=psky->lambda;
833 psky_smap[pix]=psky->slit;
834 psky++;
835 }
836
837 }
838
839#if REGDEBUG_BSPLINE
840 cpl_frame* frm = NULL;
841 sprintf(fname,"sky_model_ord_%2.2d_slice_%2.2d_iter_%2.2d.fits",abs_ord,sid,iter);
842 frm = xsh_pre_save(sky_model,fname, "SKY_MODEL",0);
843 xsh_free_frame(&frm);
844 sprintf(fname,"sky_fit_ord_%2.2d_slice_%2.2d_iter_%2.2d.fits",abs_ord,sid,iter);
845 frm = xsh_pre_save(fit_model,fname, "SKY_FIT",0);
846 xsh_free_frame(&frm);
847 sprintf(fname,"sky_wmap_%2.2d.fits",iter);
848 cpl_image_save(sky_wmap,"sky_wmap.fits", XSH_PRE_DATA_BPP, NULL, CPL_IO_DEFAULT);
849 sprintf(fname,"sky_smap_%2.2d.fits",iter);
850 cpl_image_save(sky_smap,"sky_smap.fits", XSH_PRE_DATA_BPP, NULL, CPL_IO_DEFAULT);
851#endif
852
853 cleanup:
854 cpl_free(final_name);
855 xsh_pre_free(&fit_model);
856 xsh_pre_free(&sky_model);
857 xsh_free_image(&sky_wmap);
858 xsh_free_image(&sky_smap);
859
860 return cpl_error_get_code() ;
861}
871cpl_error_code
873 xsh_instrument * instr,
874 const int iter)
875{
876
877 xsh_pre* obj_model = NULL ;
878 cpl_image* obj_wmap = NULL;
879 cpl_image* obj_smap = NULL;
880
881 float* pdata_obj =NULL;
882 float* perrs_obj =NULL;
883 int* pqual_obj =NULL;
884
885 float* pobj_wmap =NULL;
886 float* pobj_smap =NULL;
887
888 int order;
889 int sx, sy ;
890 //char fname[80];
891
892 XSH_ASSURE_NOT_NULL( omap ) ;
893 XSH_ASSURE_NOT_NULL( instr ) ;
894
895 sx = instr->config->nx / instr->binx;
896 sy = instr->config->ny / instr->biny;
897
898 /* create a model of the object and the sky as a PRE format frame */
899 obj_model = xsh_pre_new(sx,sy);
900
901 pdata_obj = cpl_image_get_data_float(xsh_pre_get_data(obj_model));
902 perrs_obj = cpl_image_get_data_float(xsh_pre_get_errs(obj_model));
903 pqual_obj = cpl_image_get_data_int(xsh_pre_get_qual(obj_model));
904
905 /* associate to object and sky the corresponding wave and slit maps */
906 obj_wmap = cpl_image_new(sx,sy,XSH_PRE_DATA_TYPE);
907 obj_smap = cpl_image_new(sx,sy,XSH_PRE_DATA_TYPE);
908
909 pobj_wmap = cpl_image_get_data_float(obj_wmap);
910 pobj_smap = cpl_image_get_data_float(obj_smap);
911
912 int ix, iy,pix;
913 for( order = 0 ; order < omap->size ; order++ ) {
914
915 wavemap_item * pobj = omap->list[order].object;
916 int obj_size = omap->list[order].object_size;
917
918 for (int k = 0; k < obj_size; k++) {
919 ix = pobj->ix;
920 iy = pobj->iy;
921 pix = iy*sx+ix;
922
923 pdata_obj[pix]=pobj->flux;
924 perrs_obj[pix]=pobj->sigma;
925 pqual_obj[pix]=pobj->qual;
926 pobj_wmap[pix]=pobj->lambda;
927 pobj_smap[pix]=pobj->slit;
928
929 pobj++;
930 }
931
932
933 }
934#if REGDEBUG_BSPLINE
935 cpl_frame* frm = NULL;
936 sprintf(fname,"object_model_%2.2d.fits",iter);
937 frm = xsh_pre_save(obj_model,fname, "OBJECT_MODEL",0);
938 xsh_free_frame(&frm);
939
940 sprintf(fname,fname,iter);
941 cpl_image_save(obj_wmap,"object_wmap.fits", XSH_PRE_DATA_BPP, NULL, CPL_IO_DEFAULT);
942 sprintf(fname,"object_smap_%2.2d.fits",iter);
943 cpl_image_save(obj_smap,fname, XSH_PRE_DATA_BPP, NULL, CPL_IO_DEFAULT);
944#endif
945
946 cleanup:
947 xsh_pre_free(&obj_model);
948 xsh_free_image(&obj_wmap);
949 xsh_free_image(&obj_smap);
950
951 return cpl_error_get_code() ;
952}
953
954
964cpl_error_code
966 xsh_instrument * instr,
967 const char * prefix)
968{
969
970 xsh_pre* obj_model = NULL ;
971 xsh_pre* sky_model = NULL ;
972 cpl_image* obj_wmap = NULL;
973 cpl_image* obj_smap = NULL;
974 cpl_image* sky_wmap = NULL;
975 cpl_image* sky_smap = NULL;
976
977 float* pdata_obj =NULL;
978 float* perrs_obj =NULL;
979 int* pqual_obj =NULL;
980
981 float* pdata_sky =NULL;
982 float* perrs_sky =NULL;
983 int* pqual_sky =NULL;
984
985 float* pobj_wmap =NULL;
986 float* pobj_smap =NULL;
987
988 float* psky_wmap =NULL;
989 float* psky_smap =NULL;
990
991 int order;
992 int sx, sy ;
993 char * final_name = NULL;
994
995 XSH_ASSURE_NOT_NULL( wmap ) ;
996 XSH_ASSURE_NOT_NULL( prefix ) ;
997 XSH_ASSURE_NOT_NULL( instr ) ;
998
999 final_name = xsh_stringcat_any( prefix, ".fits", (void*)NULL ) ;
1000 xsh_msg( "Entering xsh_wavemap_save, file \"%s\"", final_name ) ;
1001
1002 sx = instr->config->nx /instr->binx;
1003 sy = instr->config->ny /instr->biny;
1004 xsh_msg( "Image size:%d,%d", sx, sy ) ;
1005
1006 /* create a model of the object and the sky as a PRE format frame */
1007 obj_model = xsh_pre_new(sx,sy);
1008 sky_model = xsh_pre_new(sx,sy);
1009
1010 pdata_obj = cpl_image_get_data_float(xsh_pre_get_data(obj_model));
1011 perrs_obj = cpl_image_get_data_float(xsh_pre_get_errs(obj_model));
1012 pqual_obj = cpl_image_get_data_int(xsh_pre_get_qual(obj_model));
1013
1014 pdata_sky = cpl_image_get_data_float(xsh_pre_get_data(sky_model));
1015 perrs_sky = cpl_image_get_data_float(xsh_pre_get_errs(sky_model));
1016 pqual_sky = cpl_image_get_data_int(xsh_pre_get_qual(sky_model));
1017
1018 /* associate to object and sky the corresponding wave and slit maps */
1019 obj_wmap = cpl_image_new(sx,sy,XSH_PRE_DATA_TYPE);
1020 obj_smap = cpl_image_new(sx,sy,XSH_PRE_DATA_TYPE);
1021 sky_wmap = cpl_image_new(sx,sy,XSH_PRE_DATA_TYPE);
1022 sky_smap = cpl_image_new(sx,sy,XSH_PRE_DATA_TYPE);
1023
1024 pobj_wmap = cpl_image_get_data_float(obj_wmap);
1025 pobj_smap = cpl_image_get_data_float(obj_smap);
1026
1027 psky_wmap = cpl_image_get_data_float(sky_wmap);
1028 psky_smap = cpl_image_get_data_float(sky_smap);
1029
1030 int ix, iy,pix;
1031 for( order = 0 ; order < wmap->size ; order++ ) {
1032
1033 wavemap_item * psky = wmap->list[order].sky;
1034 wavemap_item * pobj = wmap->list[order].object;
1035 int sky_size = wmap->list[order].sky_size;
1036 int obj_size = wmap->list[order].object_size;
1037
1038 for (int k = 0; k < obj_size; k++) {
1039 ix = pobj->ix;
1040 iy = pobj->iy;
1041 pix = iy*sx+ix;
1042 pdata_obj[pix]=pobj->flux;
1043 perrs_obj[pix]=pobj->sigma;
1044 pqual_obj[pix]=pobj->qual;
1045 pobj_wmap[pix]=pobj->lambda;
1046 pobj_smap[pix]=pobj->slit;
1047 pobj++;
1048 }
1049
1050 for (int k = 0; k < sky_size; k++) {
1051 ix = psky->ix;
1052 iy = psky->iy;
1053 pix = iy*sx+ix;
1054 pdata_sky[pix]=psky->flux;
1055 perrs_sky[pix]=psky->sigma;
1056 pqual_sky[pix]=psky->qual;
1057 psky_wmap[pix]=psky->lambda;
1058 psky_smap[pix]=psky->slit;
1059 psky++;
1060 }
1061
1062 }
1063#if REGDEBUG_BSPLINE
1064 cpl_frame* frm = NULL;
1065 frm = xsh_pre_save(obj_model,"object_model.fits", "OBJECT_MODEL",0);
1066 xsh_free_frame(&frm);
1067 frm = xsh_pre_save(sky_model,"sky_model.fits", "SKY_MODEL",0);
1068 xsh_free_frame(&frm);
1069
1070 cpl_image_save(obj_wmap,"object_wmap.fits", XSH_PRE_DATA_BPP, NULL, CPL_IO_DEFAULT);
1071 cpl_image_save(obj_smap,"object_smap.fits", XSH_PRE_DATA_BPP, NULL, CPL_IO_DEFAULT);
1072 cpl_image_save(sky_wmap,"sky_wmap.fits", XSH_PRE_DATA_BPP, NULL, CPL_IO_DEFAULT);
1073 cpl_image_save(sky_smap,"sky_smap.fits", XSH_PRE_DATA_BPP, NULL, CPL_IO_DEFAULT);
1074#endif
1075
1076 cleanup:
1077 xsh_pre_free(&obj_model);
1078 xsh_pre_free(&sky_model);
1079 xsh_free_image(&obj_wmap);
1080 xsh_free_image(&obj_smap);
1081 xsh_free_image(&sky_wmap);
1082 xsh_free_image(&sky_smap);
1083
1084 return cpl_error_get_code() ;
1085}
1086
1099cpl_frame *
1101 xsh_order_list * order_list,
1102 xsh_pre * pre,
1103 xsh_instrument * instr,
1104 const char * prefix)
1105{
1106 cpl_frame * result = NULL ;
1107 int order, x, y ;
1108 int start_y, end_y ;
1109 int xmin, xmax ;
1110 int nx, ny ;
1111 double lambda ;
1112 cpl_image * wimg = NULL ;
1113 double * dimg = NULL ;
1114 cpl_propertylist * wheader = NULL ;
1115 char * final_name = NULL;
1116
1117 XSH_ASSURE_NOT_NULL( wmap ) ;
1118 XSH_ASSURE_NOT_NULL( order_list ) ;
1119 XSH_ASSURE_NOT_NULL( pre ) ;
1120 XSH_ASSURE_NOT_NULL( prefix ) ;
1121 XSH_ASSURE_NOT_NULL( instr ) ;
1122
1123 final_name = xsh_stringcat_any( prefix, ".fits", (void*)NULL ) ;
1124 xsh_msg( "Entering xsh_wavemap_save, file \"%s\"", final_name ) ;
1125
1126 /*
1127 create the wave map image filled with zeroes.
1128 Fill the orders with lambda calculated from wave map polynomial
1129 */
1130 /*
1131 Load the order table into the order_list
1132 */
1133
1134 /*
1135 Create a new image according to parameters included in instrument
1136 structure
1137 */
1138
1139 nx = pre->nx ;
1140 ny = pre->ny ;
1141 xsh_msg/*_dbg_medium*/( "Image size:%d,%d", nx, ny ) ;
1142
1143 check( wimg = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE ) ) ;
1144 check( dimg = cpl_image_get_data_double( wimg ) ) ;
1145
1146 for( order = 0 ; order < wmap->size ; order++ ) {
1147 int count = 0 ;
1148
1149 /*
1150 // Fix
1151 if ( wmap->list[order].tcheb_pol_lambda == NULL ) {
1152 xsh_msg( "Order %d: NULL Polynome", order ) ;
1153 continue ;
1154 }
1155 if ( wmap->list[order].sky_size <= 0 ) {
1156 xsh_msg( "NOT ENOUGH DATA FOR ORDER %d",
1157 order_list->list[order].absorder ) ;
1158 continue ;
1159 }
1160 */
1161 /* get start_y, end_y from order list */
1162 start_y = xsh_order_list_get_starty( order_list, order ) ;
1163 end_y = xsh_order_list_get_endy( order_list, order ) ;
1164 xsh_msg_dbg_low( " Order %d, Ymin: %d, Ymax: %d",
1165 order, start_y, end_y ) ;
1166 xsh_msg( " Order %d, Ymin: %d, Ymax: %d",
1167 order, start_y, end_y ) ;
1168 for ( y = start_y ; y < end_y ; y++ ) {
1169 double dxmin, dxmax ;
1170
1171 /* get xmin, xmax */
1172 check( dxmin = cpl_polynomial_eval_1d(order_list->list[order].edglopoly,
1173 (double)y, NULL ) ) ;
1174 check( dxmax = cpl_polynomial_eval_1d(order_list->list[order].edguppoly,
1175 (double)y, NULL ) ) ;
1176 xmin = floor(dxmin) ;
1177 xmax = ceil(dxmax) ;
1178 xsh_msg( " Order %d, Xmin: %d, Xmax: %d",
1179 order, xmin, xmax ) ;
1180 for (x = xmin ; x<xmax ; x++ ) {
1181 count++ ;
1182 /*calculate lambda(x,y)*/
1183 lambda = xsh_wavemap_list_eval_lambda( wmap, (double)x,
1184 (double)y, order);
1185 if ( cpl_error_get_code() != CPL_ERROR_NONE){
1186 lambda =0.0;
1187 xsh_msg("degx=%d degy=%d", wmap->degx, wmap->degy);
1189 //xsh_error_reset();
1190 }
1191 /*add lambda to image */
1192 *(dimg + x + (y*nx)) = (float)lambda ;
1193 }
1194 }
1195 xsh_msg_dbg_high( " %d points to calculate", count ) ;
1196 }
1197
1198 /* Now save image */
1199 check( wheader = cpl_propertylist_duplicate( pre->data_header ) ) ;
1200
1201 check( cpl_image_save( wimg, final_name, CPL_BPP_IEEE_FLOAT,
1202 wheader, CPL_IO_DEFAULT ) ) ;
1203 check(result=xsh_frame_product(final_name,
1204 "TAG",
1205 CPL_FRAME_TYPE_IMAGE,
1206 CPL_FRAME_GROUP_PRODUCT,
1207 CPL_FRAME_LEVEL_TEMPORARY));
1208
1209 cleanup:
1210 xsh_order_list_free( &order_list);
1211 XSH_FREE( final_name ) ;
1212 xsh_free_image ( &wimg ) ;
1213 xsh_free_propertylist( &wheader);
1214
1215 return result ;
1216}
1217
1218
1219/*****************************************************************************/
1233/*****************************************************************************/
1234cpl_frame*
1236 cpl_frame *order_frame,
1237 xsh_pre *pre,
1238 xsh_instrument *instr,
1239 const char *prefix,
1240 cpl_frame **dispersol_frame,
1241 cpl_frame **slitmap_frame)
1242{
1243 cpl_frame *result = NULL;
1244 const char *slit_tag=NULL;
1245 int order;
1246 xsh_dispersol_list *dispersol_list = NULL;
1247 const char* disp_tag=NULL;
1248
1249 XSH_ASSURE_NOT_NULL( wmap);
1250 XSH_ASSURE_NOT_NULL( order_frame);
1251 XSH_ASSURE_NOT_NULL( prefix);
1252 XSH_ASSURE_NOT_NULL( dispersol_frame);
1253 XSH_ASSURE_NOT_NULL( instr);
1254
1255
1256 check( dispersol_list = xsh_dispersol_list_new( wmap->size,
1257 wmap->degx, wmap->degy, instr));
1258
1259 for( order = 0 ; order < wmap->size ; order++ ) {
1260 /* get start_y, end_y from order list */
1261 check( xsh_dispersol_list_add( dispersol_list, order,
1262 wmap->list[order].order, wmap->list[order].pol_lambda,
1263 wmap->list[order].pol_slit));
1264 wmap->list[order].pol_lambda = NULL;
1265 wmap->list[order].pol_slit = NULL;
1266 }
1267
1268 if (pre != NULL){
1269
1270 check( result = xsh_dispersol_list_to_wavemap( dispersol_list,
1271 order_frame, pre,
1272 instr,prefix));
1273
1274 slit_tag = XSH_GET_TAG_FROM_ARM( XSH_SLIT_MAP_POLY, instr);
1275
1276
1277 check (*slitmap_frame = xsh_dispersol_list_to_slitmap( dispersol_list,
1278 order_frame, pre,
1279 instr, slit_tag));
1280 }
1281
1282 if(strstr(cpl_frame_get_tag(order_frame),"AFC") != NULL) {
1283 disp_tag = XSH_GET_TAG_FROM_ARM( XSH_DISP_TAB_AFC,instr);
1284 } else {
1285 disp_tag = XSH_GET_TAG_FROM_ARM( XSH_DISP_TAB,instr);
1286 }
1287 check(*dispersol_frame = xsh_dispersol_list_save( dispersol_list,disp_tag));
1288
1289 cleanup:
1290 if (cpl_error_get_code() != CPL_ERROR_NONE){
1291 xsh_free_frame( &result);
1292 }
1293 xsh_dispersol_list_free( &dispersol_list);
1294
1295
1296 return result ;
1297}
1298/*****************************************************************************/
cpl_frame * xsh_dispersol_list_save(xsh_dispersol_list *list, const char *tag)
Save a dispersion list on the disk.
void xsh_dispersol_list_add(xsh_dispersol_list *list, int idx, int absorder, cpl_polynomial *lambda_poly, cpl_polynomial *slit_poly)
Add a dispersion solution in the list.
void xsh_dispersol_list_free(xsh_dispersol_list **list)
Free the dispersion list.
cpl_frame * xsh_dispersol_list_to_slitmap(xsh_dispersol_list *list, cpl_frame *order_frame, xsh_pre *pre, xsh_instrument *instr, const char *tag)
Save a SLIT MAP image.
xsh_dispersol_list * xsh_dispersol_list_new(int size, int degx, int degy, xsh_instrument *instrument)
Create a new dispersion solution list.
cpl_frame * xsh_dispersol_list_to_wavemap(xsh_dispersol_list *list, cpl_frame *order_frame, xsh_pre *pre, xsh_instrument *instr, const char *tag)
Save a WAVE MAP image.
xsh_order_list * xsh_order_list_load(cpl_frame *frame, xsh_instrument *instr)
load an order list from a frame
int xsh_order_list_get_starty(xsh_order_list *list, int i)
get position on Y axis of first pixel detected on order
int xsh_order_list_get_endy(xsh_order_list *list, int i)
get position on Y axis of last pixel detected on order
void xsh_order_list_free(xsh_order_list **list)
free memory associated to an order_list
cpl_image * xsh_pre_get_data(xsh_pre *pre)
Get data.
cpl_image * xsh_pre_get_qual(xsh_pre *pre)
Get qual.
xsh_pre * xsh_pre_new(int nx, int ny)
Create new PRE image.
void xsh_pre_free(xsh_pre **pre)
Free a xsh_pre structure.
Definition: xsh_data_pre.c:823
cpl_image * xsh_pre_get_errs(xsh_pre *pre)
Get errs.
cpl_frame * xsh_pre_save(const xsh_pre *pre, const char *filename, const char *tag, int temp)
Save PRE on disk.
void xsh_wavemap_list_compute_poly(double *vlambda, double *vslit, double *xpos, double *ypos, int nitems, double *orders, xsh_dispersol_param *dispsol_param, xsh_wavemap_list *wmap)
compute a wave-map-list
void xsh_wavemap_list_dump(xsh_wavemap_list *list, const char *fname)
Dump main info about a rec table (for each order of the list)
void xsh_wavemap_list_free(xsh_wavemap_list **list)
free memory associated to a wavemap_list
cpl_frame * xsh_wavemap_list_save_poly(xsh_wavemap_list *wmap, cpl_frame *order_frame, xsh_pre *pre, xsh_instrument *instr, const char *prefix, cpl_frame **dispersol_frame, cpl_frame **slitmap_frame)
Save the wave_map slit_map and disp_tab.
cpl_error_code xsh_wavemap_list_sky_image_save(xsh_wavemap_list *smap, xsh_instrument *instr, const int abs_ord, const int sid, const int iter)
cpl_error_code xsh_wavemap_list_object_image_save(xsh_wavemap_list *omap, xsh_instrument *instr, const int iter)
static double xsh_wavemap_list_eval_lambda(xsh_wavemap_list *wmap, double x, double y, int ord)
eval the polynomial wave solution in x,y,ord
cpl_frame * xsh_wavemap_list_save2(xsh_wavemap_list *wmap, xsh_order_list *order_list, xsh_pre *pre, xsh_instrument *instr, const char *prefix)
void xsh_wavemap_list_compute(double *vlambda, double *xpos, double *ypos, int nitems, double *orders, xsh_dispersol_param *dispsol_param, xsh_wavemap_list *wmap)
compute a wave-map-list
cpl_error_code xsh_wavemap_list_save4debug(xsh_wavemap_list *wmap, xsh_instrument *instr, const char *prefix)
static void lambda_fit(double *lambda, double *xpos, double *ypos, int size, xsh_wavemap_list *wmap, xsh_dispersol_param *dispsol_param, int order)
void xsh_wavemap_list_set_max_size(xsh_wavemap_list *list, int idx, int absorder, int max_size)
set max size of wavemap
xsh_wavemap_list * xsh_wavemap_list_create(xsh_instrument *instr)
create an empty order list
cpl_frame * xsh_wavemap_list_save(xsh_wavemap_list *wmap, cpl_frame *order_frame, xsh_pre *pre, xsh_instrument *instr, const char *prefix)
#define XSH_REGDEBUG(...)
Definition: xsh_error.h:132
#define XSH_ASSURE_NOT_ILLEGAL(cond)
Definition: xsh_error.h:107
#define check(COMMAND)
Definition: xsh_error.h:71
#define xsh_error_reset()
Definition: xsh_error.h:87
#define XSH_ASSURE_NOT_ILLEGAL_MSG(cond, msg)
Definition: xsh_error.h:111
#define XSH_ASSURE_NOT_NULL(pointer)
Definition: xsh_error.h:99
XSH_INSTRCONFIG * xsh_instrument_get_config(xsh_instrument *i)
Get the instrument default set of keywords.
int size
int * y
int * x
#define xsh_msg_dbg_medium(...)
Definition: xsh_msg.h:44
#define xsh_msg(...)
Print a message on info level.
Definition: xsh_msg.h:121
#define xsh_msg_dbg_low(...)
Definition: xsh_msg.h:48
#define xsh_msg_dbg_high(...)
Definition: xsh_msg.h:40
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
void xsh_free_image(cpl_image **i)
Deallocate an image and set the pointer to NULL.
Definition: xsh_utils.c:2116
void xsh_free_frame(cpl_frame **f)
Deallocate a frame and set the pointer to NULL.
Definition: xsh_utils.c:2269
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
char * xsh_stringcat_any(const char *s,...)
Concatenate an arbitrary number of strings.
Definition: xsh_utils.c:1925
void xsh_free_propertylist(cpl_propertylist **p)
Deallocate a property list and set the pointer to NULL.
Definition: xsh_utils.c:2179
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
XSH_INSTRCONFIG * config
xsh_order * list
cpl_polynomial * edguppoly
cpl_polynomial * edglopoly
cpl_propertylist * data_header
Definition: xsh_data_pre.h:66
xsh_instrument * instrument
xsh_wavemap * list
cpl_propertylist * header
wavemap_item * sky
cpl_polynomial * pol_lambda
cpl_polynomial * pol_slit
wavemap_item * all
cpl_polynomial * tcheb_pol_lambda
wavemap_item * object
#define XSH_PRE_DATA_TYPE
Definition: xsh_data_pre.h:42
#define XSH_PRE_DATA_BPP
Definition: xsh_data_pre.h:43
int nx
int ny
int order
Definition: xsh_detmon_lg.c:80
int xsh_print_rec_status(const int val)
Check if an error has happened and returns error kind and location.
Definition: xsh_dfs.c:877
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_SLIT_MAP_POLY
Definition: xsh_dfs.h:151
#define XSH_DISP_TAB
Definition: xsh_dfs.h:897
#define XSH_GET_TAG_FROM_ARM(TAG, instr)
Definition: xsh_dfs.h:1548
#define XSH_DISP_TAB_AFC
Definition: xsh_dfs.h:902
#define XSH_NEW_PROPERTYLIST(POINTER)
Definition: xsh_utils.h:70
#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
cpl_polynomial * xsh_polynomial_fit_2d_create(cpl_bivector *xy_pos, cpl_vector *values, cpl_size *degree, double *mse)