X-shooter Pipeline Reference Manual 3.8.15
xsh_data_dispersol.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-18 14:15:44 $
23 * $Revision: 1.38 $
24 * $Name: not supported by cvs2svn $
25 */
26
27#ifdef HAVE_CONFIG_H
28#include <config.h>
29#endif
30
31/*---------------------------------------------------------------------------*/
36/*---------------------------------------------------------------------------*/
37
38
42/*-----------------------------------------------------------------------------
43 Includes
44 ----------------------------------------------------------------------------*/
45
46#include <math.h>
47
48#include <xsh_data_dispersol.h>
49#include <xsh_utils.h>
50#include <xsh_error.h>
51#include <xsh_msg.h>
52#include <xsh_pfits.h>
53#include <xsh_dfs.h>
54#include <cpl.h>
55#include <xsh_utils_table.h>
57#include <xsh_data_order.h>
58
59
60/*----------------------------------------------------------------------------
61 Function implementation
62 ----------------------------------------------------------------------------*/
63
64/*--------------------------------------------------------------------------*/
76/*--------------------------------------------------------------------------*/
79{
80 xsh_dispersol_list *result = NULL;
81
82 /* check input parameters */
85
86 XSH_CALLOC( result, xsh_dispersol_list, 1);
87 result->size = size;
88 result->degx = degx;
89 result->degy = degy;
92 XSH_CALLOC( result->list, xsh_dispersol, result->size);
94
95 cleanup:
96 if ( cpl_error_get_code() != CPL_ERROR_NONE) {
98 }
99 return result;
100}
101/*--------------------------------------------------------------------------*/
102
103
104/*--------------------------------------------------------------------------*/
113/*--------------------------------------------------------------------------*/
116{
117 xsh_dispersol_list *result = NULL;
118 const char* tablename = NULL;
119 cpl_table* table = NULL;
120 int rows;
121 //int cols ;
122 int degx, degy;
123 char coefname[20];
124 int i;
125
126 /* check input parameters */
128
129 /* get table filename */
130 check( tablename = cpl_frame_get_filename(frame));
131
132 XSH_TABLE_LOAD( table, tablename);
133 check( rows = cpl_table_get_nrow( table));
134 //check( cols = cpl_table_get_ncol( table));
135
136
138 CPL_TYPE_INT, 0, &degx));
139
141 CPL_TYPE_INT, 0, &degy));
142
143 check( result = xsh_dispersol_list_new( rows/2, degx, degy, instrument));
144
145 for(i=0; i< rows/2; i++){
146 cpl_polynomial *lambda_poly = NULL;
147 cpl_polynomial *slit_poly = NULL;
148 cpl_size pows[2];
149 int k,l, absorder = 0;
150 int j = 2*i;
151
153 CPL_TYPE_INT, j, &absorder));
154 check( lambda_poly = cpl_polynomial_new( 2));
155 check( slit_poly = cpl_polynomial_new( 2));
156 for( k=0; k<= degx; k++){
157 for( l=0; l<= degy; l++){
158 double coef_lambda = 0.0, coef_slit = 0.0;
159 sprintf(coefname,"C%d%d",k,l);
160 check( xsh_get_table_value( table, coefname,CPL_TYPE_DOUBLE, j,
161 &coef_lambda));
162 check( xsh_get_table_value( table, coefname,CPL_TYPE_DOUBLE, j+1,
163 &coef_slit));
164 pows[0] =k;
165 pows[1] = l;
166 check( cpl_polynomial_set_coeff( lambda_poly, pows, coef_lambda));
167 check( cpl_polynomial_set_coeff( slit_poly, pows, coef_slit));
168 }
169 }
170 check( xsh_dispersol_list_add( result, i, absorder, lambda_poly, slit_poly));
171 }
172
173 cleanup:
174 xsh_free_table( &table);
175 if ( cpl_error_get_code() != CPL_ERROR_NONE) {
176 xsh_dispersol_list_free( &result);
177 }
178 return result;
179}
180/*--------------------------------------------------------------------------*/
181
182/*--------------------------------------------------------------------------*/
197/*--------------------------------------------------------------------------*/
199 int absorder, cpl_polynomial *lambda_poly,
200 cpl_polynomial *slit_poly)
201{
202 /* check input parameters */
203 XSH_ASSURE_NOT_NULL( list);
204 XSH_ASSURE_NOT_NULL( lambda_poly);
205 XSH_ASSURE_NOT_NULL( slit_poly);
206 XSH_ASSURE_NOT_ILLEGAL( idx >=0 && idx < list->size);
207
208 list->list[idx].absorder = absorder;
209 list->list[idx].lambda_poly = lambda_poly;
210 list->list[idx].slit_poly = slit_poly;
211
212 cleanup:
213 return;
214}
215/*--------------------------------------------------------------------------*/
216
217/*--------------------------------------------------------------------------*/
224/*--------------------------------------------------------------------------*/
226{
227 int i = 0;
228
229 if ( list && *list){
230 /* free the list */
231 for( i = 0; i < (*list)->size; i++){
232 xsh_free_polynomial( &(*list)->list[i].lambda_poly);
233 xsh_free_polynomial( &(*list)->list[i].slit_poly);
234 }
235 if ((*list)->list){
236 cpl_free( (*list)->list);
237 }
238 xsh_free_propertylist( &((*list)->header));
239 cpl_free( *list);
240 *list = NULL;
241 }
242}
243/*--------------------------------------------------------------------------*/
244
245
246/*--------------------------------------------------------------------------*/
261/*--------------------------------------------------------------------------*/
263 cpl_frame *order_frame,
264 xsh_pre *pre,
265 xsh_instrument *instr,
266 const char* tag)
267
268{
269 cpl_frame *result = NULL;
270 cpl_image *image = NULL;
271 xsh_order_list *order_list = NULL;
272 int nx, ny;
273 cpl_vector* pos = NULL;
274 double *data = NULL;
275 int i, y;
276 char filename[256];
277 cpl_propertylist *wavemap_header = NULL;
278
279 /* Check input */
280 XSH_ASSURE_NOT_NULL( list);
281 XSH_ASSURE_NOT_NULL( order_frame);
283 XSH_ASSURE_NOT_NULL( instr);
284
285 /* Load data */
286 check( order_list = xsh_order_list_load( order_frame, instr));
287 nx = pre->nx;
288 ny = pre->ny;
289 pos = cpl_vector_new(2);
290
291
292 check( image = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE));
293 check( data = cpl_image_get_data_double( image));
294 check( wavemap_header = cpl_propertylist_new());
295
296
297 for( i=0; i < list->size; i++) {
298 int starty, endy;
299 double lambda_min, lambda_max;
300 double lambda_starty_lo, lambda_starty_up, lambda_endy_lo, lambda_endy_up;
301 double dxmin, dxmax ;
302 int xmin, xmax;
303 int absorder;
304
305
306 absorder = order_list->list[i].absorder;
307 check( starty = xsh_order_list_get_starty( order_list, i));
308
309 check( endy = xsh_order_list_get_endy( order_list, i));
310
311
312 for ( y=starty ; y<=endy && y<ny; y++){
313 int x;
314
315 /* get xmin, xmax from order table*/
316 check( dxmin = xsh_order_list_eval( order_list,
317 order_list->list[i].edglopoly,
318 (double)y ) ) ;
319 check( dxmax = xsh_order_list_eval( order_list,
320 order_list->list[i].edguppoly,
321 (double)y ) ) ;
322 xmin = ceil( dxmin);
323 xmax = floor( dxmax);
324 for ( x=xmin ; x<=xmax && x<nx; x++){
325 double lambda;
326
327 cpl_vector_set(pos, 0, (double)x);
328 cpl_vector_set(pos, 1, (double)y);
329 /*calculate lambda(x,y)*/
330
332 list, list->list[i].lambda_poly, pos));
333 /*add lambda to image */
334 data[x-1+(y-1)*nx] = (float)lambda;
335
336 }
337
338 }
339
340 /* compute maxium and minimum of possible wavelength for this order */
341 check( dxmin = xsh_order_list_eval( order_list,
342 order_list->list[i].edglopoly,
343 (double)starty));
344 check( dxmax = xsh_order_list_eval( order_list,
345 order_list->list[i].edguppoly,
346 (double)starty));
347 xmin = ceil( dxmin);
348 xmax = floor( dxmax);
349
350 cpl_vector_set(pos, 0, (double)xmin);
351 cpl_vector_set(pos, 1, (double)starty);
352 check( lambda_starty_lo = xsh_dispersol_list_eval(
353 list, list->list[i].lambda_poly, pos));
354
355 cpl_vector_set(pos, 0, (double)xmax);
356 check( lambda_starty_up = xsh_dispersol_list_eval(
357 list, list->list[i].lambda_poly, pos));
358
359 check( dxmin = xsh_order_list_eval( order_list,
360 order_list->list[i].edglopoly,
361 (double)endy));
362 check( dxmax = xsh_order_list_eval( order_list,
363 order_list->list[i].edguppoly,
364 (double)endy));
365 xmin = ceil( dxmin);
366 xmax = floor( dxmax);
367
368 cpl_vector_set(pos, 0, (double)xmin);
369 cpl_vector_set(pos, 1, (double)endy);
370 check( lambda_endy_lo = xsh_dispersol_list_eval(
371 list, list->list[i].lambda_poly, pos));
372
373 cpl_vector_set(pos, 0, (double)xmax);
374 check( lambda_endy_up = xsh_dispersol_list_eval(
375 list, list->list[i].lambda_poly, pos));
376
377 if ( lambda_starty_lo < lambda_endy_lo
378 && lambda_starty_up < lambda_endy_up){
379 if ( lambda_starty_lo > lambda_starty_up){
380 lambda_min = lambda_starty_lo;
381 }
382 else{
383 lambda_min = lambda_starty_up;
384 }
385 if ( lambda_endy_lo > lambda_endy_up){
386 lambda_max = lambda_endy_up;
387 }
388 else{
389 lambda_max = lambda_endy_lo;
390 }
391 }
392 else if ( lambda_starty_lo > lambda_endy_lo
393 && lambda_starty_up > lambda_endy_up){
394 if ( lambda_starty_lo > lambda_starty_up){
395 lambda_max = lambda_starty_up;
396 }
397 else{
398 lambda_max = lambda_starty_lo;
399 }
400 if ( lambda_endy_lo > lambda_endy_up){
401 lambda_min = lambda_endy_lo;
402 }
403 else{
404 lambda_min = lambda_endy_up;
405 }
406 }
407 else{
409 "wavelength variation differs between upper(%f %f) and lower edges( %f %f)",
410 lambda_starty_up, lambda_endy_up, lambda_endy_up, lambda_endy_lo);
411 lambda_min = 0.0;
412 lambda_max = 0.0;
413 }
414
415 /* write KW */
417 wavemap_header, absorder, lambda_min));
419 wavemap_header, absorder, lambda_max));
420 }
421
422 sprintf(filename,"%s.fits",tag);
423 check( xsh_pfits_set_pcatg( wavemap_header, tag));
424
425 check( cpl_image_save( image, filename, CPL_BPP_IEEE_FLOAT,
426 wavemap_header, CPL_IO_DEFAULT));
427
428 check( result = xsh_frame_product(filename, tag,
429 CPL_FRAME_TYPE_IMAGE, CPL_FRAME_GROUP_PRODUCT,
430 CPL_FRAME_LEVEL_TEMPORARY));
431
432 cleanup:
433 if (cpl_error_get_code() != CPL_ERROR_NONE){
434 xsh_free_frame( &result);
435 }
436 xsh_free_image ( &image);
437 xsh_free_propertylist( &wavemap_header);
438 xsh_order_list_free( &order_list);
439 xsh_free_vector( &pos);
440
441
442 return result;
443}
444/*--------------------------------------------------------------------------*/
445
446
447/*--------------------------------------------------------------------------*/
462/*--------------------------------------------------------------------------*/
464 cpl_frame *order_frame,
465 xsh_pre *pre,
466 xsh_instrument *instr,
467 const char* tag)
468{
469 cpl_frame *result = NULL;
470 cpl_image *image = NULL;
471 cpl_image *ifu_map = NULL;
472 xsh_order_list *order_list = NULL;
473 int nx, ny;
474 cpl_vector* pos = NULL;
475 double *data = NULL;
476 int *imap = NULL;
477 double crpix1=1.;
478 double crpix2=1.;
479 double crval1=1.;
480 double crval2=1.;
481 double cdelt1=1.;
482 double cdelt2=1.;
483
484 int i, y;
485 char filename[256];
486 cpl_vector *med_low_vect = NULL, *med_up_vect = NULL, *med_cen_vect = NULL;
487 cpl_vector *med_sliclo_vect = NULL, *med_slicup_vect = NULL;
488 int nb_orders;
489 double slit_up, slit_low, slit_cen;
490 cpl_propertylist *slitmap_header = NULL;
491
492 /* Check input */
493 XSH_ASSURE_NOT_NULL( list);
494 XSH_ASSURE_NOT_NULL( order_frame);
496 XSH_ASSURE_NOT_NULL( instr);
497
498 /* Load data */
499 check( order_list = xsh_order_list_load( order_frame, instr));
500 nx = pre->nx;
501 ny = pre->ny;
502
503 /* Check that the binning are the same */
504 XSH_ASSURE_NOT_ILLEGAL( order_list->bin_x == pre->binx);
505 XSH_ASSURE_NOT_ILLEGAL( order_list->bin_x == list->binx);
506
507 pos = cpl_vector_new(2);
508
509 check( image = cpl_image_new( nx, ny, CPL_TYPE_DOUBLE));
510 check( data = cpl_image_get_data_double( image));
511 check( ifu_map = cpl_image_new( nx, ny, CPL_TYPE_INT));
512 check( imap = cpl_image_get_data_int( ifu_map));
513 check( slitmap_header = cpl_propertylist_new());
514
515 nb_orders = list->size;
516
517 med_low_vect = cpl_vector_new( nb_orders);
518 med_cen_vect = cpl_vector_new( nb_orders);
519 med_up_vect = cpl_vector_new( nb_orders);
520 if ( xsh_instrument_get_mode( instr) == XSH_MODE_IFU){
521 med_sliclo_vect = cpl_vector_new( nb_orders);
522 med_slicup_vect = cpl_vector_new( nb_orders);
523 }
524
525
526 for( i=0; i < nb_orders; i++) {
527 int starty, endy, vect_size;
528 int absorder;
529 cpl_vector *low_vect = NULL, *up_vect = NULL, *cen_vect = NULL;
530 cpl_vector *sliclo_vect = NULL, *slicup_vect = NULL;
531
532 check( starty = xsh_order_list_get_starty( order_list, i));
533 check( endy = xsh_order_list_get_endy( order_list, i));
534 absorder = order_list->list[i].absorder;
535
536 //vect_size = endy-starty+1;
537 int ymax = ( endy < ny ) ? endy:(ny-1);
538
539 vect_size = ymax-starty+1;
540/*
541 xsh_msg("starty=%d endy=%d ny=%d ymax=%d vect_size=%d",
542 starty,endy,ny,ymax,vect_size);
543*/
544 if(vect_size<1) {
545 /* exit in vcase of 0 or negative sixze of vector */
546 continue;
547
548 }
549 low_vect = cpl_vector_new( vect_size);
550 cen_vect = cpl_vector_new( vect_size);
551 up_vect = cpl_vector_new( vect_size);
552
553 if ( xsh_instrument_get_mode( instr) == XSH_MODE_IFU){
554 sliclo_vect = cpl_vector_new( vect_size);
555 slicup_vect = cpl_vector_new( vect_size);
556 }
557
558 for ( y=starty ; y<=ymax; y++){
559 //for ( y=starty ; y<=endy && y<ny; y++){
560 double dxmin, dxmax, xcen;
561 int xmin, xmax, x;
562 double slit;
563 double ypos;
564
565 /* get xmin, xmax */
566 check( dxmin = xsh_order_list_eval( order_list,
567 order_list->list[i].edglopoly,
568 (double)y ) ) ;
569 check( dxmax = xsh_order_list_eval( order_list,
570 order_list->list[i].edguppoly,
571 (double)y ) ) ;
572 check( xcen = xsh_order_list_eval( order_list,
573 order_list->list[i].cenpoly,
574 (double)y ) ) ;
575
576 xmin = ceil( dxmin);
577 xmax = floor( dxmax);
578
579 ypos = (double)y;
580 cpl_vector_set(pos, 1, ypos);
581
582 for ( x=xmin ; x<=xmax && x<nx; x++){
583 double xpos;
584
585 xpos = (double)x;
586 cpl_vector_set(pos, 0, xpos);
587 cpl_vector_set(pos, 1, ypos);
588 /*calculate lambda(x,y)*/
589 check( slit = xsh_dispersol_list_eval( list, list->list[i].slit_poly, pos));
590
591 /*add lambda to image */
592 data[x-1+(y-1)*nx] = (float)slit;
593 imap[x-1+(y-1)*nx] = 1;
594 }
595
596 cpl_vector_set(pos, 0, xmin);
597 cpl_vector_set(pos, 1, ypos);
598 check( slit = xsh_dispersol_list_eval( list, list->list[i].slit_poly, pos));
599 cpl_vector_set( low_vect, y-starty, slit);
600
601 cpl_vector_set(pos, 0, xmax);
602 cpl_vector_set(pos, 1, ypos);
603 check( slit = xsh_dispersol_list_eval( list, list->list[i].slit_poly, pos));
604 cpl_vector_set( up_vect, y-starty, slit);
605
606 cpl_vector_set(pos, 0, xcen);
607 cpl_vector_set(pos, 1, ypos);
608 check( slit = xsh_dispersol_list_eval( list, list->list[i].slit_poly, pos));
609 cpl_vector_set( cen_vect, y-starty, slit);
610
611 if ( xsh_instrument_get_mode( instr) == XSH_MODE_IFU){
612 xmin = xsh_order_list_eval( order_list,
613 order_list->list[i].sliclopoly,
614 (double)y );
615 xmax = xsh_order_list_eval( order_list,
616 order_list->list[i].slicuppoly,
617 (double)y );
618 cpl_vector_set(pos, 0, xmin);
619 cpl_vector_set(pos, 1, ypos);
620 check( slit = xsh_dispersol_list_eval( list, list->list[i].slit_poly, pos));
621 cpl_vector_set( sliclo_vect, y-starty, slit);
622
623 cpl_vector_set(pos, 0, xmax);
624 cpl_vector_set(pos, 1, ypos);
625 check( slit = xsh_dispersol_list_eval( list, list->list[i].slit_poly, pos));
626 cpl_vector_set( slicup_vect, y-starty, slit);
627 }
628
629 }
630
631 check( slit_up = cpl_vector_get_median( up_vect));
632 check( slit_low = cpl_vector_get_median( low_vect));
633 check( slit_cen = cpl_vector_get_median( cen_vect));
634 //xsh_msg("med up=%g low=%g cen=%g",slit_up,slit_low,slit_cen);
635
636 /* write some KW in header */
637 check( xsh_pfits_set_slitmap_order_cen( slitmap_header, absorder, slit_cen));
638 check( cpl_vector_set( med_cen_vect, i, slit_cen));
639
640 if ( slit_up > slit_low){
641 check( xsh_pfits_set_slitmap_order_edgup( slitmap_header, absorder, slit_up));
642 check( xsh_pfits_set_slitmap_order_edglo( slitmap_header, absorder, slit_low));
643 check( cpl_vector_set( med_low_vect, i, slit_low));
644 check( cpl_vector_set( med_up_vect, i, slit_up));
645 }
646 else{
647 check( xsh_pfits_set_slitmap_order_edgup( slitmap_header, absorder, slit_low));
648 check( xsh_pfits_set_slitmap_order_edglo( slitmap_header, absorder, slit_up));
649 check( cpl_vector_set( med_low_vect, i, slit_up));
650 check( cpl_vector_set( med_up_vect, i, slit_low));
651 }
652
653 if ( xsh_instrument_get_mode( instr) == XSH_MODE_IFU){
654 check( slit_up = cpl_vector_get_median( slicup_vect));
655 check( slit_low = cpl_vector_get_median( sliclo_vect));
656 if ( slit_up > slit_low){
657 check( xsh_pfits_set_slitmap_order_slicup( slitmap_header, absorder, slit_up));
658 check( xsh_pfits_set_slitmap_order_sliclo( slitmap_header, absorder, slit_low));
659 check( cpl_vector_set( med_sliclo_vect, i, slit_low));
660 check( cpl_vector_set( med_slicup_vect, i, slit_up));
661 }
662 else{
663 check( xsh_pfits_set_slitmap_order_slicup( slitmap_header, absorder, slit_low));
664 check( xsh_pfits_set_slitmap_order_sliclo( slitmap_header, absorder, slit_up));
665 check( cpl_vector_set( med_sliclo_vect, i, slit_up));
666 check( cpl_vector_set( med_slicup_vect, i, slit_low));
667 }
668 xsh_free_vector( &slicup_vect);
669 xsh_free_vector( &sliclo_vect);
670 }
671 xsh_free_vector( &up_vect);
672 xsh_free_vector( &low_vect);
673 xsh_free_vector( &cen_vect);
674 }
675
676 check( slit_cen = cpl_vector_get_median( med_cen_vect));
677 check( slit_up = cpl_vector_get_median( med_up_vect));
678 check( slit_low = cpl_vector_get_median( med_low_vect));
679
680 check( xsh_pfits_set_slitmap_median_cen( slitmap_header, slit_cen));
681 check( xsh_pfits_set_slitmap_median_edgup( slitmap_header, slit_up));
682 check( xsh_pfits_set_slitmap_median_edglo( slitmap_header, slit_low));
683
684 if ( xsh_instrument_get_mode( instr) == XSH_MODE_IFU){
685 check( slit_up = cpl_vector_get_median( med_slicup_vect));
686 check( slit_low = cpl_vector_get_median( med_sliclo_vect));
687 check( xsh_pfits_set_slitmap_median_slicup( slitmap_header, slit_up));
688 check( xsh_pfits_set_slitmap_median_sliclo( slitmap_header, slit_low));
689 }
690
691 sprintf(filename,"%s.fits",tag);
692 xsh_pfits_set_pcatg( slitmap_header, tag);
693 cdelt1=xsh_instrument_get_binx(instr);
694 cdelt2=xsh_instrument_get_biny(instr);
695 xsh_pfits_set_wcs(slitmap_header,crpix1,crval1,cdelt1,crpix2,crval2,cdelt2);
696 check( cpl_image_save( image, filename, CPL_BPP_IEEE_FLOAT,
697 slitmap_header, CPL_IO_DEFAULT));
698
699 if ( xsh_instrument_get_mode( instr) == XSH_MODE_IFU){
700 check( cpl_image_save( ifu_map, filename, CPL_BPP_32_SIGNED,
701 NULL, CPL_IO_EXTEND));
702 }
703 check( result = xsh_frame_product(filename,tag,
704 CPL_FRAME_TYPE_IMAGE, CPL_FRAME_GROUP_PRODUCT,
705 CPL_FRAME_LEVEL_TEMPORARY));
706
707 cleanup:
708 if (cpl_error_get_code() != CPL_ERROR_NONE){
709 xsh_free_frame( &result);
710 }
711
712 xsh_free_image ( &image);
713 xsh_free_image ( &ifu_map);
714 xsh_free_propertylist( &slitmap_header);
715 xsh_order_list_free( &order_list);
716 xsh_free_vector( &pos);
717 xsh_free_vector( &med_up_vect);
718 xsh_free_vector( &med_low_vect);
719 xsh_free_vector( &med_cen_vect);
720 xsh_free_vector( &med_slicup_vect);
721 xsh_free_vector( &med_sliclo_vect);
722
723 return result;
724}
725/*--------------------------------------------------------------------------*/
726
727/*--------------------------------------------------------------------------*/
740/*--------------------------------------------------------------------------*/
741double
743 cpl_polynomial *poly,
744 cpl_vector *pos)
745{
746 double result=0;
747 double y_bin, y_no_bin;
748 double x_bin, x_no_bin;
749
750 XSH_ASSURE_NOT_NULL( list);
752 XSH_ASSURE_NOT_NULL( poly);
753
754 check( x_bin = cpl_vector_get(pos, 0));
755 check( y_bin = cpl_vector_get(pos, 1));
756
757 x_no_bin = convert_bin_to_data( x_bin, list->binx);
758 y_no_bin = convert_bin_to_data( y_bin, list->biny);
759
760 check( cpl_vector_set(pos, 0, x_no_bin));
761 check( cpl_vector_set(pos, 1, y_no_bin));
762
763 check( result = cpl_polynomial_eval( poly, pos));
764
765 cleanup:
766 return result;
767}
768/*--------------------------------------------------------------------------*/
769
770/*--------------------------------------------------------------------------*/
778/*--------------------------------------------------------------------------*/
779cpl_frame*
781{
782 cpl_frame *result = NULL;
783 cpl_table *table = NULL;
784 char coefname[20];
785 int i,j,k;
786 int nbcoefs =0;
787 cpl_size pows[2];
788 char filename[256];
789
790 /* Check input parameters */
791 XSH_ASSURE_NOT_NULL( list);
792
793 nbcoefs = (list->degx+1)*(list->degy+1);
794 /* Create a table */
795 check( table = cpl_table_new( XSH_DISPERSOL_TABLE_NBCOL+nbcoefs));
796 check(cpl_table_set_size(table, list->size*2));
797
798 check(
799 cpl_table_new_column( table, XSH_DISPERSOL_TABLE_COLNAME_AXIS,
800 CPL_TYPE_STRING));
801 check(
802 cpl_table_new_column( table, XSH_DISPERSOL_TABLE_COLNAME_ORDER,
803 CPL_TYPE_INT));
804 check(
805 cpl_table_new_column( table, XSH_DISPERSOL_TABLE_COLNAME_DEGX,
806 CPL_TYPE_INT));
807 check(
808 cpl_table_new_column( table, XSH_DISPERSOL_TABLE_COLNAME_DEGY,
809 CPL_TYPE_INT));
810
811 for(i=0; i<= list->degx; i++){
812 for(j = 0; j<= list->degy; j++){
813 sprintf(coefname,"C%d%d",i,j);
814 check( cpl_table_new_column( table, coefname, CPL_TYPE_DOUBLE));
815 }
816 }
817
818 for(i=0; i< list->size; i++){
819 int i2 = i*2;
820 check(cpl_table_set_string(table,XSH_DISPERSOL_TABLE_COLNAME_AXIS,
822 check(cpl_table_set_string(table,XSH_DISPERSOL_TABLE_COLNAME_AXIS,
824 check(cpl_table_set_int(table,XSH_DISPERSOL_TABLE_COLNAME_ORDER,
825 i2, list->list[i].absorder));
826 check(cpl_table_set_int(table,XSH_DISPERSOL_TABLE_COLNAME_ORDER,
827 i2+1, list->list[i].absorder));
828 check(cpl_table_set_int(table,XSH_DISPERSOL_TABLE_COLNAME_DEGX,
829 i2, list->degx));
830 check(cpl_table_set_int(table,XSH_DISPERSOL_TABLE_COLNAME_DEGX,
831 i2+1, list->degx));
832 check(cpl_table_set_int(table,XSH_DISPERSOL_TABLE_COLNAME_DEGY,
833 i2, list->degy));
834 check(cpl_table_set_int(table,XSH_DISPERSOL_TABLE_COLNAME_DEGY,
835 i2+1, list->degy));
836 for( j=0; j<= list->degx; j++){
837 for( k=0; k<= list->degy; k++){
838 double coef_lambda = 0.0, coef_slit = 0.0;
839 sprintf(coefname,"C%d%d",j,k);
840 pows[0] = j;
841 pows[1] = k;
842 check(coef_lambda =
843 cpl_polynomial_get_coeff( list->list[i].lambda_poly, pows));
844 check(cpl_table_set_double( table, coefname, i2, coef_lambda));
845 check(coef_slit =
846 cpl_polynomial_get_coeff( list->list[i].slit_poly, pows));
847 check(cpl_table_set_double( table, coefname, i2+1, coef_slit));
848 }
849 }
850 }
851
852 sprintf(filename, "%s.fits",tag);
853 check( xsh_pfits_set_pcatg(list->header,tag));
854 check(cpl_table_save(table, list->header, NULL,
855 filename, CPL_IO_DEFAULT));
856 /* Create the frame */
857 check( result= xsh_frame_product( filename, tag,
858 CPL_FRAME_TYPE_TABLE,
859 CPL_FRAME_GROUP_PRODUCT,
860 CPL_FRAME_LEVEL_TEMPORARY));
861
862 cleanup:
863 xsh_free_table( &table);
864 if ( cpl_error_get_code() != CPL_ERROR_NONE) {
865 }
866 return result;
867}
868/*--------------------------------------------------------------------------*/
869
static int starty
static int endy
static xsh_instrument * instrument
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.
xsh_dispersol_list * xsh_dispersol_list_load(cpl_frame *frame, xsh_instrument *instrument)
Load a dispersion list from a frame.
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.
double xsh_dispersol_list_eval(xsh_dispersol_list *list, cpl_polynomial *poly, cpl_vector *pos)
Evaluate the polynomial according the binning.
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
double xsh_order_list_eval(xsh_order_list *list, cpl_polynomial *poly, double y)
Evaluate an order list poly.
#define XSH_ASSURE_NOT_ILLEGAL(cond)
Definition: xsh_error.h:107
#define check(COMMAND)
Definition: xsh_error.h:71
#define XSH_ASSURE_NOT_NULL(pointer)
Definition: xsh_error.h:99
int xsh_instrument_get_binx(xsh_instrument *instrument)
XSH_MODE xsh_instrument_get_mode(xsh_instrument *i)
Get a mode on instrument structure.
int xsh_instrument_get_biny(xsh_instrument *instrument)
int size
int * y
int * x
#define xsh_msg_warning(...)
Print an warning message.
Definition: xsh_msg.h:88
void xsh_pfits_set_pcatg(cpl_propertylist *plist, const char *value)
Write the PCATG value.
Definition: xsh_pfits.c:1008
void xsh_pfits_set_slitmap_median_cen(cpl_propertylist *plist, double value)
Definition: xsh_pfits.c:3881
void xsh_pfits_set_slitmap_median_slicup(cpl_propertylist *plist, double value)
Definition: xsh_pfits.c:3864
void xsh_pfits_set_slitmap_order_cen(cpl_propertylist *plist, int absorder, double value)
Definition: xsh_pfits.c:3794
cpl_error_code xsh_pfits_set_wcs(cpl_propertylist *header, const double crpix1, const double crval1, const double cdelt1, const double crpix2, const double crval2, const double cdelt2)
Definition: xsh_pfits.c:4364
void xsh_pfits_set_slitmap_median_edglo(cpl_propertylist *plist, double value)
Definition: xsh_pfits.c:3830
void xsh_pfits_set_slitmap_order_edgup(cpl_propertylist *plist, int absorder, double value)
Definition: xsh_pfits.c:3724
void xsh_pfits_set_slitmap_median_edgup(cpl_propertylist *plist, double value)
Definition: xsh_pfits.c:3813
void xsh_pfits_set_wavemap_order_lambda_max(cpl_propertylist *plist, int absorder, double value)
Definition: xsh_pfits.c:3997
void xsh_pfits_set_slitmap_order_slicup(cpl_propertylist *plist, int absorder, double value)
Definition: xsh_pfits.c:3759
void xsh_pfits_set_slitmap_median_sliclo(cpl_propertylist *plist, double value)
Definition: xsh_pfits.c:3847
void xsh_pfits_set_slitmap_order_edglo(cpl_propertylist *plist, int absorder, double value)
Definition: xsh_pfits.c:3742
void xsh_pfits_set_slitmap_order_sliclo(cpl_propertylist *plist, int absorder, double value)
Definition: xsh_pfits.c:3776
void xsh_pfits_set_wavemap_order_lambda_min(cpl_propertylist *plist, int absorder, double value)
Definition: xsh_pfits.c:3979
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_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 convert_bin_to_data(double bin_data, int binning)
Definition: xsh_utils.c:3228
void xsh_free_table(cpl_table **t)
Deallocate a table and set the pointer to NULL.
Definition: xsh_utils.c:2133
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.
cpl_propertylist * header
xsh_dispersol * list
cpl_polynomial * slit_poly
cpl_polynomial * lambda_poly
xsh_order * list
cpl_polynomial * edguppoly
cpl_polynomial * edglopoly
cpl_polynomial * slicuppoly
cpl_polynomial * sliclopoly
cpl_polynomial * cenpoly
int binx
Definition: xsh_data_pre.h:80
#define XSH_DISPERSOL_TABLE_COLNAME_DEGY
#define XSH_DISPERSOL_AXIS_SLIT
#define XSH_DISPERSOL_TABLE_NBCOL
#define XSH_DISPERSOL_TABLE_COLNAME_DEGX
#define XSH_DISPERSOL_TABLE_COLNAME_ORDER
#define XSH_DISPERSOL_AXIS_LAMBDA
#define XSH_DISPERSOL_TABLE_COLNAME_AXIS
@ XSH_MODE_IFU
int nx
int ny
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_NEW_PROPERTYLIST(POINTER)
Definition: xsh_utils.h:70
#define XSH_CALLOC(POINTER, TYPE, SIZE)
Definition: xsh_utils.h:56
#define XSH_TABLE_LOAD(TABLE, NAME)