X-shooter Pipeline Reference Manual 3.8.15
xsh_create_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-18 14:15:44 $
23 * $Revision: 1.56 $
24*/
25
26#ifdef HAVE_CONFIG_H
27#include <config.h>
28#endif
29
30/*---------------------------------------------------------------------------*/
39#include <xsh_drl.h>
40#include <xsh_pfits.h>
41#include <xsh_pfits_qc.h>
42#include <xsh_utils.h>
43#include <xsh_data_order.h>
44#include <xsh_error.h>
45#include <xsh_utils.h>
46#include <xsh_msg.h>
47#include <xsh_data_pre.h>
48#include <xsh_data_instrument.h>
49#include <xsh_data_order.h>
50#include <xsh_data_wavesol.h>
51#include <xsh_data_resid_tab.h>
52#include <xsh_data_wavemap.h>
54#include <xsh_model_io.h>
55#include <xsh_model_kernel.h>
56#include <cpl.h>
57
58
59
60/******************************************************************************/
61/*
62 @brief
63 Create maps from a dispersion solution
64 @param[in] dispsol_frame
65 The dispersion solution
66 @param(in] ordertab_frame
67 The order table
68 @param(in] pre_frame
69 The science image frame
70 @param[in] instrument
71 The data instrument structure
72 @param(out] wavemap_frame
73 The wave map frame
74 @param[out] slitmap_frame
75 The slit map frame
76*/
77/******************************************************************************/
78void xsh_create_map( cpl_frame *dispsol_frame, cpl_frame *ordertab_frame,
79 cpl_frame *pre_frame, xsh_instrument *instrument, cpl_frame **wavemap_frame,
80 cpl_frame **slitmap_frame,const char* rec_prefix)
81{
82 xsh_dispersol_list *dispersol_tab = NULL;
83 xsh_pre *pre = NULL;
84 char wavemap_tag[256];
85 char slitmap_tag[256];
86 /* Checking input parameters */
87 XSH_ASSURE_NOT_NULL( dispsol_frame);
88 XSH_ASSURE_NOT_NULL( ordertab_frame);
89 XSH_ASSURE_NOT_NULL( pre_frame);
91 XSH_ASSURE_NOT_NULL( wavemap_frame);
92 XSH_ASSURE_NOT_NULL( slitmap_frame);
93
94 /* Loading data */
95 check( pre = xsh_pre_load( pre_frame, instrument));
96 check( dispersol_tab = xsh_dispersol_list_load( dispsol_frame,
97 instrument));
98 sprintf(wavemap_tag,"%s_%s",
100 sprintf(slitmap_tag,"%s_%s",
102
103 check( *wavemap_frame = xsh_dispersol_list_to_wavemap( dispersol_tab,
104 ordertab_frame, pre, instrument,wavemap_tag));
105 check( *slitmap_frame = xsh_dispersol_list_to_slitmap( dispersol_tab,
106 ordertab_frame,
107 pre, instrument,
108 slitmap_tag));
109
110 cleanup:
111 xsh_dispersol_list_free( &dispersol_tab);
112 xsh_pre_free( &pre);
113 return;
114}
115/******************************************************************************/
116
117/******************************************************************************/
118/*
119 @brief
120 Create a wavemap using the physical model
121 @param[in] model_frame
122 The model configuration frame
123 @param[in] instrument
124 The instrument description
125*/
126/******************************************************************************/
127void xsh_create_model_map( cpl_frame* model_frame, xsh_instrument* instrument,
128 const char* wtag, const char* stag,
129 cpl_frame **wavemap_frame,
130 cpl_frame **slitmap_frame,const int save_tmp)
131{
132 xsh_xs_3 model_config;
133 //xsh_xs_3* p_xs_3;
134 //xsh_pre * pre_sci = NULL ;
135 XSH_ASSURE_NOT_NULL_MSG( model_frame,"If model-scenario is 0 make sure that the input model cfg has at least one parameter with Compute_Flag set to 1 and High_Limit>Low_limit");
137 XSH_ASSURE_NOT_NULL( wavemap_frame);
138 XSH_ASSURE_NOT_NULL( slitmap_frame);
139 XSH_ASSURE_NOT_NULL( wtag);
140 XSH_ASSURE_NOT_NULL( stag);
141
142 check( xsh_model_config_load_best( model_frame, &model_config));
143/*
144 check(xsh_model_temperature_update_frame(&config_model_frame,
145 ref_frame,
146 instr,&found_temp));
147*/
148
149 check( xsh_model_binxy( &model_config, instrument->binx,
150 instrument->biny));
151
152 check(xsh_model_maps_create(&model_config,instrument,wtag,stag,
153 wavemap_frame,slitmap_frame,save_tmp));
154
155 /*
156 xsh_msg("Generating dispersol frame...");
157 check( *dispsol_frame = xsh_model_dispsol_create(&model_config,
158 instrument,dtag));
159 */
160
161 cleanup:
162 // xsh_pre_free( &pre_sci);
163
164 return;
165}
166/******************************************************************************/
167
168/******************************************************************************/
169/*
170 @brief
171 Create a wavemap using the 2D solution
172 @param[in] pre_frame
173 The image frame
174 @param[in] wave_tab_2d_frame
175 The 2D solution
176 @param[in] order_tab_frame
177 The order tab frame
178 @param[in] spectral_format_frame
179 The spectral format frame
180 @param[in] dispsol_par
181 The dispersion solution parameters
182 @param[in] instrument
183 The instrument description
184 @param[in] vm_tag
185 The wavemap tag
186 @param[out] dispersol_frame
187 The dispersion solution frame
188 @param[out] slitmap_frame
189 The slitmap frame
190
191 @return
192 The wavemap
193*/
194/******************************************************************************/
195cpl_frame *
196xsh_create_poly_wavemap( cpl_frame *pre_frame,
197 cpl_frame *wave_tab_2d_frame,
198 cpl_frame *order_tab_frame,
199 cpl_frame *spectral_format_frame,
200 xsh_dispersol_param *dispsol_par,
202 const char * wm_tag,
203 cpl_frame **dispersol_frame,
204 cpl_frame** slitmap_frame)
205{
206 xsh_pre *pre = NULL;
207 xsh_spectralformat_list *spec_list = NULL;
208 xsh_wavesol *wave_tab_2d = NULL;
209 cpl_frame *result = NULL ;
210 float slit_step = 1.5;
211 float lambda_step = 0.1;
212 //float sol_min_lambda=0;
213 //float sol_max_lambda=0;
214 float sol_min_slit, sol_max_slit;
215 //int sol_min_order, sol_max_order;
216 int i, idx, size, slit_size;
217 double j, k;
218 double *vlambda = NULL, *vslit = NULL, *vorder = NULL;
219 double *pos_x = NULL, *pos_y = NULL;
220 /* take from create_wavemap*/
221 xsh_wavemap_list *wm_list = NULL ;
222 int binx;
223 //int biny ;
224 char wm_name[256];
225
226 /* Check input parameters */
227 XSH_ASSURE_NOT_NULL( wave_tab_2d_frame);
228 XSH_ASSURE_NOT_NULL( order_tab_frame);
229 XSH_ASSURE_NOT_NULL( spectral_format_frame);
230 XSH_ASSURE_NOT_NULL( dispsol_par);
232 XSH_ASSURE_NOT_NULL( dispersol_frame);
233
234 /* Load data */
235 if (pre_frame != NULL){
236 check( pre = xsh_pre_load( pre_frame, instrument));
237 }
238 /* Get binning (used by wavesol) */
240 //check( biny = xsh_instrument_get_biny( instrument )) ;
241
242 check( spec_list = xsh_spectralformat_list_load( spectral_format_frame,
243 instrument));
244 check( wave_tab_2d = xsh_wavesol_load( wave_tab_2d_frame,
245 instrument ));
246 /* Set binning, used in wavesol_eval_polx/y */
247 check( xsh_wavesol_set_bin_x( wave_tab_2d, binx ) ) ;
248 check( xsh_wavesol_set_bin_y( wave_tab_2d, binx ) ) ;
249
250 //sol_min_lambda = wave_tab_2d->min_lambda;
251 //sol_max_lambda = wave_tab_2d->max_lambda;
252 //sol_min_order = wave_tab_2d->min_order;
253 //sol_max_order = wave_tab_2d->max_order;
254 sol_min_slit = wave_tab_2d->min_slit;
255 sol_max_slit = wave_tab_2d->max_slit;
256 /*
257 xsh_msg("sol lambda %f,%f order %d,%d slit %f,%f",
258 sol_min_lambda,sol_max_lambda, sol_min_order, sol_max_order,
259 sol_min_slit, sol_max_slit);
260 */
261 /* evaluate size of grid points */
262 slit_size = (int)((sol_max_slit-sol_min_slit)/slit_step)+1;
263
264 size=0;
265 for( i=0; i< spec_list->size; i++){
266 float lambda_min, lambda_max;
267 int lambda_size;
268
269 lambda_min = spec_list->list[i].lambda_min_full;
270 lambda_max = spec_list->list[i].lambda_max_full;
271 XSH_ASSURE_NOT_ILLEGAL_MSG(lambda_max>= lambda_min,
272 "lambda_max< lambda_min!! Check input spectralformat table, columns WLMAXFUL and WLMINFUL");
273 lambda_size = (int)((lambda_max-lambda_min)/lambda_step)+1;
274 size += lambda_size*slit_size;
275 }
276 xsh_msg_dbg_medium( "size %d", size ) ;
277
278 idx = 0;
279 XSH_MALLOC( vorder, double, size);
280 XSH_MALLOC( vlambda, double, size);
281 XSH_MALLOC( vslit, double, size);
282 XSH_MALLOC( pos_x, double, size);
283 XSH_MALLOC( pos_y, double, size);
284
285 for( i=0; i< spec_list->size; i++){
286 double absorder;
287 float lambda_min, lambda_max;
288
289 absorder= (double)spec_list->list[i].absorder;
290 lambda_min = spec_list->list[i].lambda_min_full;
291 lambda_max = spec_list->list[i].lambda_max_full;
292 xsh_msg("order %f lambda %f-%f",absorder, lambda_min, lambda_max);
293
294 for(j=lambda_min; j <= lambda_max; j+=lambda_step){
295 for( k=sol_min_slit; k <= sol_max_slit; k+=slit_step){
296 double x, y;
297
298 check( x = xsh_wavesol_eval_polx( wave_tab_2d, j, absorder, k));
299 check( y = xsh_wavesol_eval_poly( wave_tab_2d, j, absorder, k));
300 vorder[idx] = absorder;
301 vlambda[idx] = j;
302 vslit[idx] = k;
303 pos_x[idx] = x;
304 pos_y[idx] = y;
305 idx++;
306 }
307 }
308 xsh_msg_dbg_medium( "i %d idx %d / %d", i, idx, size);
309 }
310
312 FILE *test = NULL;
313 int itest;
314
315 test = fopen( "wavemap_grid.log", "w");
316 for( itest=0; itest < size; itest++){
317 fprintf(test, "%f %f\n", pos_x[itest], pos_y[itest]);
318 }
319 fclose(test);
320 }
321
322 /* take from create_wavemap*/
324#if 0
325 check( xsh_wavemap_list_compute( vlambda, pos_x, pos_y, size,
326 vorder, wm_par, wm_list));
327#else
328 /* Had to replace 'size' by 'idx' because on MacOSX Intel, under certain
329 circumstances, the 'size' is too large by 8, though 'idx' is OK !!!
330 Actually on the Mac the following code gives a size of 766, instead of 765
331 on Linux !
332 Don't understand why !!!! Laurent
333 */
334 check( xsh_wavemap_list_compute_poly( vlambda, vslit, pos_x, pos_y, idx,
335 vorder, dispsol_par, wm_list));
336#endif
337 sprintf(wm_name,"%s.fits",wm_tag);
338 check( result = xsh_wavemap_list_save_poly( wm_list, order_tab_frame,
339 pre, instrument, wm_tag, dispersol_frame, slitmap_frame));
340 if ( pre != NULL){
341 //check( xsh_add_temporary_file( wm_name));
342 check( cpl_frame_set_tag( result,wm_tag));
343 check( cpl_frame_set_tag( *slitmap_frame,
345 }
346 cleanup:
347
348 XSH_FREE( vorder);
349 XSH_FREE( vlambda);
350 XSH_FREE( vslit);
351 XSH_FREE( pos_x);
352 XSH_FREE( pos_y);
353 xsh_pre_free( &pre);
354 xsh_spectralformat_list_free( &spec_list);
355 xsh_wavesol_free( &wave_tab_2d);
356 xsh_wavemap_list_free( &wm_list);
357 return result;
358}
359
360
361
362cpl_frame*
363xsh_create_dispersol_physmod(cpl_frame *pre_frame,
364 cpl_frame *order_tab_frame,
365 cpl_frame* mod_cfg_frame,
366 cpl_frame* wave_map_frame,
367 cpl_frame* slit_map_frame,
368 xsh_dispersol_param *dispsol_param,
369 cpl_frame* spectral_format_frame,
371 const int clean_tmp)
372{
373
374
375 xsh_pre *pre = NULL;
376 cpl_frame* disp_frame=NULL;
377 xsh_spectralformat_list *spec_list = NULL;
378 //int binx=0;
379 //int biny=0 ;
380 xsh_xs_3 model_config;
381 int slit_size=0;
382 int slit_step=1.5;
383
384
385 //float sol_min_lambda=0, sol_max_lambda=0;
386 float sol_min_slit=0, sol_max_slit=0;
387 int sol_min_order=0, sol_max_order=0;
388 cpl_image* wmap_ima=NULL;
389 cpl_image* smap_ima=NULL;
390 cpl_frame* loc_slit_map_frame=NULL;
391 const char* name=NULL;
392 int size=0;
393 int i=0;
394 int idx=0;
395 float lambda_step=0.1;
396
397 double *vlambda = NULL, *vslit = NULL, *vorder = NULL;
398 double *pos_x = NULL, *pos_y = NULL;
399 double j=0;
400 double k=0;
401
402 xsh_wavemap_list *wm_list = NULL ;
403 cpl_frame* result=NULL;
404 char wm_tag[256];
405 char wm_name[256];
406
407 /* Check input parameters */
408 XSH_ASSURE_NOT_NULL_MSG(mod_cfg_frame,"Null model cfg frame!");
409 XSH_ASSURE_NOT_NULL_MSG(spectral_format_frame,"Null spectral format frame!");
410 XSH_ASSURE_NOT_NULL_MSG( instrument,"Null instrument setting!");
411
412 if (pre_frame != NULL){
413 check( pre = xsh_pre_load( pre_frame, instrument));
414 }
415
416 /* Get binning (used by wavesol) */
417 //check( binx = xsh_instrument_get_binx( instrument )) ;
418 //check( biny = xsh_instrument_get_biny( instrument )) ;
419 check(spec_list=xsh_spectralformat_list_load(spectral_format_frame,
420 instrument));
421 check( xsh_model_config_load_best( mod_cfg_frame, &model_config));
422 check( xsh_model_binxy( &model_config, instrument->binx,
423 instrument->biny));
424
425 name=cpl_frame_get_filename(wave_map_frame);
426 wmap_ima=cpl_image_load(name,CPL_TYPE_FLOAT,0,0);
427
428 name=cpl_frame_get_filename(slit_map_frame);
429 smap_ima=cpl_image_load(name,CPL_TYPE_FLOAT,0,0);
430
431
432 //sol_min_lambda = cpl_image_get_min(wmap_ima);
433 //sol_max_lambda = cpl_image_get_max(wmap_ima);
434 xsh_free_image(&wmap_ima);
435
436 sol_min_slit = cpl_image_get_min(smap_ima);
437 sol_max_slit = cpl_image_get_max(smap_ima);
438 xsh_free_image(&smap_ima);
439
440 /* evaluate size of grid points */
441 slit_size = (int)((sol_max_slit-sol_min_slit)/slit_step)+1;
442
443 size=0;
444 for( i=0; i< spec_list->size; i++){
445 float lambda_min, lambda_max;
446 int lambda_size;
447
448 sol_min_order = (sol_min_order > spec_list->list[i].absorder) ? spec_list->list[i].absorder : sol_min_order;
449 sol_max_order = (sol_max_order < spec_list->list[i].absorder) ? spec_list->list[i].absorder : sol_max_order;
450
451 lambda_min = spec_list->list[i].lambda_min_full;
452 lambda_max = spec_list->list[i].lambda_max_full;
453 XSH_ASSURE_NOT_ILLEGAL_MSG(lambda_max>= lambda_min,
454 "lambda_max< lambda_min!! Check input spectralformat table, columns WLMAXFUL and WLMINFUL");
455 lambda_size = (int)((lambda_max-lambda_min)/lambda_step)+1;
456 size += lambda_size*slit_size;
457 }
458 xsh_msg_dbg_medium( "size %d", size ) ;
459
460
461 /* As in poly mode fill vector with positions */
462 /* Here we could directly use wave and slit maps, as we know results */
463 idx = 0;
464 XSH_MALLOC( vorder, double, size);
465 XSH_MALLOC( vlambda, double, size);
466 XSH_MALLOC( vslit, double, size);
467 XSH_MALLOC( pos_x, double, size);
468 XSH_MALLOC( pos_y, double, size);
469
470 for( i=0; i< spec_list->size; i++){
471 double absorder;
472 float lambda_min, lambda_max;
473 int morder=0;
474
475 absorder= (double)spec_list->list[i].absorder;
476 morder= spec_list->list[i].absorder;
477 lambda_min = spec_list->list[i].lambda_min_full;
478 lambda_max = spec_list->list[i].lambda_max_full;
479 xsh_msg("order %f lambda %f-%f",absorder, lambda_min, lambda_max);
480
481 for(j=lambda_min; j <= lambda_max; j+=lambda_step){
482 for( k=sol_min_slit; k <= sol_max_slit; k+=slit_step){
483 double x, y;
484 xsh_model_get_xy(&model_config,instrument,j,morder,k,&x,&y);
485 vorder[idx] = absorder;
486 vlambda[idx] = j;
487 vslit[idx] = k;
488 pos_x[idx] = x;
489 pos_y[idx] = y;
490 idx++;
491 }
492 }
493 xsh_msg_dbg_medium( "i %d idx %d / %d", i, idx, size);
494 }
495
497 FILE *test = NULL;
498 int itest;
499
500 test = fopen( "wavemap_grid.log", "w");
501 for( itest=0; itest < size; itest++){
502 fprintf(test, "%f %f\n", pos_x[itest], pos_y[itest]);
503 }
504 fclose(test);
505 }
506
507 /* generates wave list in order to have later the dispersion relation */
509
510 check( xsh_wavemap_list_compute_poly( vlambda, vslit, pos_x, pos_y, idx,
511 vorder, dispsol_param, wm_list));
512
513 sprintf(wm_tag,"WAVE_MAP_POLY_%s",xsh_instrument_arm_tostring(instrument));
514 sprintf(wm_name,"%s.fits",wm_tag);
515 check( result = xsh_wavemap_list_save_poly( wm_list, order_tab_frame,
516 pre, instrument, wm_tag,
517 &disp_frame,&loc_slit_map_frame));
518
519 if(clean_tmp) {
520 xsh_add_temporary_file(wm_name);
521 }
522 /* This generate recipe crash
523 if ( pre != NULL){
524 check( xsh_add_temporary_file( wm_name));
525 check( cpl_frame_set_tag( result,wm_tag));
526 check( cpl_frame_set_tag( &slit_map_frame,
527 XSH_GET_TAG_FROM_ARM( XSH_SLIT_MAP_POLY, instrument)));
528 }
529 */
530
531
532 cleanup:
533 xsh_free_image(&wmap_ima);
534 xsh_free_image(&smap_ima);
535 xsh_free_frame(&result);
536 xsh_free_frame(&loc_slit_map_frame);
537
538 XSH_FREE( vorder);
539 XSH_FREE( vlambda);
540 XSH_FREE( vslit);
541 XSH_FREE( pos_x);
542 XSH_FREE( pos_y);
543 xsh_pre_free( &pre);
544 xsh_spectralformat_list_free( &spec_list);
545 xsh_wavemap_list_free( &wm_list);
546 return disp_frame;
547}
548
549
550/*---------------------------------------------------------------------------*/
562/*---------------------------------------------------------------------------*/
563
564cpl_error_code
565xsh_wavemap_qc(cpl_frame* frm_map,const cpl_frame* frm_tab)
566{
567 int sx=0;
568 //int sy=0;
569 double* cx=NULL;
570 double* cy=NULL;
571 int wx=0;
572 int wy=0;
573 cpl_image* ima=NULL;
574 const char* name_tab=NULL;
575 const char* name_map=NULL;
576 char qc_wlen[40];
577 double* pima=NULL;
578 cpl_propertylist* plist=NULL;
579 cpl_table* tab=NULL;
580 cpl_table* ext=NULL;
581
582 double wlen=0;
583 int ord_min=0;
584 int ord_max=0;
585 int i=0;
586 int next=0;
587
588 XSH_ASSURE_NOT_NULL(frm_map);
589 XSH_ASSURE_NOT_NULL(frm_tab);
590 check(name_tab=cpl_frame_get_filename(frm_tab));
591 check(tab=cpl_table_load(name_tab,2,0));
592 check(ord_min=cpl_table_get_column_min(tab,"ABSORDER"));
593 check(ord_max=cpl_table_get_column_max(tab,"ABSORDER"));
594
595
596 name_map=cpl_frame_get_filename(frm_map);
597 ima=cpl_image_load(name_map,CPL_TYPE_DOUBLE,0,0);
598 pima=cpl_image_get_data_double(ima);
599 sx=cpl_image_get_size_x(ima);
600 //sy=cpl_image_get_size_y(ima);
601 plist=cpl_propertylist_load(name_map,0);
602 for(i=ord_min;i<=ord_max;i++) {
603 next=cpl_table_and_selected_int(tab,"ABSORDER",CPL_EQUAL_TO,i);;
604 ext=cpl_table_extract_selected(tab);
605 cx=cpl_table_get_data_double(ext,"CENTER_X");
606 cy=cpl_table_get_data_double(ext,"CENTER_Y");
607 wx=cx[next/2];
608 wy=cy[next/2];
609 wlen=pima[wx+sx*wy];
610 sprintf(qc_wlen,"%s%d",XSH_QC_WMAP_WAVEC,i);
611 cpl_propertylist_append_double(plist,qc_wlen,wlen);
612 xsh_free_table(&ext);
613 cpl_table_select_all(tab);
614
615 }
616
618
619 cleanup:
620 xsh_free_image(&ima);
621 xsh_free_table(&tab);
622 xsh_free_table(&ext);
623 xsh_free_propertylist(&plist);
624 return cpl_error_get_code();
625}
626
627
628/*---------------------------------------------------------------------------*/
640/*---------------------------------------------------------------------------*/
641
642cpl_error_code
643xsh_wavetab_qc(cpl_frame* frm_tab, const int is_poly)
644{
645
646 cpl_table* ext=NULL;
647 const char* name_tab=NULL;
648 cpl_table* tab=NULL;
649 cpl_table* tbl=NULL;
650 int nsel=0;
651 int ord_min=0;
652 int ord_max=0;
653 int i=0;
654 int j=0;
655 cpl_vector* loc_vec=NULL;
656 double* pvec=NULL;
657 double* ptab_all_ord=NULL;
658 double* py=NULL;
659 double ymin=0;
660 double ymax=0;
661 double ymin_all=FLT_MAX;
662 double ymax_all=FLT_MIN;
663
664 double ymed=0;
665 double yavg=0;
666 char qc_line[40];
667 cpl_propertylist* plist=NULL;
668 cpl_table* tab_all_ord=NULL;
669 int ymin_all_ord=0;
670 int ymax_all_ord=0;
671
672 XSH_ASSURE_NOT_NULL(frm_tab);
673
674 check(name_tab=cpl_frame_get_filename(frm_tab));
675 check(tab=cpl_table_load(name_tab,1,0));
676 check(ord_min=cpl_table_get_column_min(tab,"Order"));
677 check(ord_max=cpl_table_get_column_max(tab,"Order"));
678 check(plist=cpl_propertylist_load(name_tab,0));
679
680 check(nsel=cpl_table_and_selected_int(tab,"Slit_index",CPL_EQUAL_TO,4));
681 check(ext=cpl_table_extract_selected(tab));
682 tab_all_ord=cpl_table_new(ord_max-ord_min+1);
683 cpl_table_new_column(tab_all_ord,"data",CPL_TYPE_DOUBLE);
684 cpl_table_fill_column_window_double(tab_all_ord,"data",0,ord_max-ord_min+1,0);
685 check(ptab_all_ord=cpl_table_get_data_double(tab_all_ord,"data"));
686
687 for(i=ord_min;i<=ord_max;i++) {
688 check(nsel=cpl_table_and_selected_int(ext,"Order",CPL_EQUAL_TO,i));
689 check(tbl=cpl_table_extract_selected(ext));
690
691 if(nsel>1) {
692 check(loc_vec=cpl_vector_new(nsel-1));
693 check(pvec=cpl_vector_get_data(loc_vec));
694
695 if(is_poly) {
696 py=cpl_table_get_data_double(tbl,"Ypoly");
697 } else {
698 py=cpl_table_get_data_double(tbl,"Ythanneal");
699 }
700
701 for(j=0;j<nsel-1;j++) {
702 pvec[j]=fabs(py[j+1]-py[j]);
703
704 if(pvec[j]>ymax_all) {
705 ymax_all=pvec[j];
706 ymax_all_ord=i;
707 }
708 if(pvec[j]<ymin_all) {
709 ymin_all=pvec[j];
710 ymin_all_ord=i;
711 }
712
713
714 }
715
716 check(ymin=cpl_vector_get_min(loc_vec));
717 check(ymax=cpl_vector_get_max(loc_vec));
718 check(ymed=cpl_vector_get_median(loc_vec));
719 check(yavg=cpl_vector_get_mean(loc_vec));
720 check(ptab_all_ord[i-ord_min]=yavg);
721
722
723 //xsh_msg("ymin=%g ymax=%g ymed=%g yavg=%g",ymin,ymax,ymed,yavg);
724 sprintf(qc_line,"%s%d",XSH_QC_LINE_DIFMIN,i);
725 check(cpl_propertylist_append_double(plist,qc_line,ymin));
726 check(cpl_propertylist_set_comment(plist,qc_line,XSH_QC_LINE_DIFMIN_C));
727
728 sprintf(qc_line,"%s%d",XSH_QC_LINE_DIFMAX,i);
729 check(cpl_propertylist_append_double(plist,qc_line,ymax));
730 check(cpl_propertylist_set_comment(plist,qc_line,XSH_QC_LINE_DIFMAX_C));
731
732 sprintf(qc_line,"%s%d",XSH_QC_LINE_DIFMED,i);
733 check(cpl_propertylist_append_double(plist,qc_line,ymed));
734 check(cpl_propertylist_set_comment(plist,qc_line,XSH_QC_LINE_DIFMED_C));
735
736 sprintf(qc_line,"%s%d",XSH_QC_LINE_DIFAVG,i);
737 check(cpl_propertylist_append_double(plist,qc_line,yavg));
738 check(cpl_propertylist_set_comment(plist,qc_line,XSH_QC_LINE_DIFAVG_C));
739
740 xsh_free_vector(&loc_vec);
741 } else {
742 xsh_msg_warning("Too few values of 'Slit_index=4' for Order=%d",i);
743 xsh_msg_warning("No %s (and similar) QC parameters can be generated for order %d",XSH_QC_LINE_DIFMIN,i);
744 check(cpl_table_set_invalid(tab_all_ord,"data",i-ord_min));
745 }
746 check(cpl_table_select_all(ext));
747 xsh_free_table(&tbl);
748 }
749 cpl_table_erase_invalid(tab_all_ord);
750
751
752 check(yavg=cpl_table_get_column_mean(tab_all_ord,"data"));
753 check(ymax=cpl_table_get_column_max(tab_all_ord,"data"));
754 check(ymed=cpl_table_get_column_median(tab_all_ord,"data"));
755 check(ymin=cpl_table_get_column_min(tab_all_ord,"data"));
756 /*
757 xsh_msg("All: ymin_all=%g ymin_all_ord=%g ymax_all=%g ymax_all_ord=%g ymed=%g yavg=%g",
758 ymin_all,ymin_all_ord, ymax_all,ymax_all_ord, ymed,yavg);
759 */
760 sprintf(qc_line,"%s",XSH_QC_LINE_DIFMIN);
761 check(cpl_propertylist_append_double(plist,qc_line,ymin_all));
762 check(cpl_propertylist_set_comment(plist,qc_line,XSH_QC_LINE_DIFMIN_C));
763
764 sprintf(qc_line,"%s",XSH_QC_LINE_DIFMIN_ORD);
765 check(cpl_propertylist_append_int(plist,qc_line,ymin_all_ord));
766 check(cpl_propertylist_set_comment(plist,qc_line,XSH_QC_LINE_DIFMIN_C));
767
768 sprintf(qc_line,"%s",XSH_QC_LINE_DIFMAX);
769 check(cpl_propertylist_append_double(plist,qc_line,ymax_all));
770 check(cpl_propertylist_set_comment(plist,qc_line,XSH_QC_LINE_DIFMAX_C));
771
772 sprintf(qc_line,"%s",XSH_QC_LINE_DIFMAX_ORD);
773 check(cpl_propertylist_append_int(plist,qc_line,ymax_all_ord));
774 check(cpl_propertylist_set_comment(plist,qc_line,XSH_QC_LINE_DIFMAX_C));
775
776 sprintf(qc_line,"%s",XSH_QC_LINE_DIFMED);
777 check(cpl_propertylist_append_double(plist,qc_line,ymed));
778 check(cpl_propertylist_set_comment(plist,qc_line,XSH_QC_LINE_DIFMED_C));
779
780 sprintf(qc_line,"%s",XSH_QC_LINE_DIFAVG);
781 check(cpl_propertylist_append_double(plist,qc_line,yavg));
782 check(cpl_propertylist_set_comment(plist,qc_line,XSH_QC_LINE_DIFAVG_C));
783
784 check(cpl_table_select_all(tab));
785 check(cpl_table_save(tab, plist, NULL,name_tab, CPL_IO_DEFAULT));
786
787 cleanup:
788 xsh_free_table(&tab);
789 xsh_free_table(&ext);
790 xsh_free_table(&tbl);
791 xsh_free_vector(&loc_vec);
792 xsh_free_table(&tab_all_ord);
793 xsh_free_propertylist(&plist);
794
795
796 return cpl_error_get_code();
797
798}
799
static xsh_instrument * instrument
int binx
static float slit_step
static float lambda_step
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.
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_pre * xsh_pre_load(cpl_frame *frame, xsh_instrument *instr)
Load a xsh_pre structure from a frame.
Definition: xsh_data_pre.c:849
void xsh_pre_free(xsh_pre **pre)
Free a xsh_pre structure.
Definition: xsh_data_pre.c:823
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_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_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.
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
xsh_wavemap_list * xsh_wavemap_list_create(xsh_instrument *instr)
create an empty order list
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.
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_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
#define XSH_ASSURE_NOT_NULL_MSG(pointer, msg)
Definition: xsh_error.h:103
#define check(COMMAND)
Definition: xsh_error.h:71
#define XSH_ASSURE_NOT_ILLEGAL_MSG(cond, msg)
Definition: xsh_error.h:111
#define XSH_ASSURE_NOT_NULL(pointer)
Definition: xsh_error.h:99
const char * xsh_instrument_arm_tostring(xsh_instrument *i)
Get the string associated with an arm.
int xsh_instrument_get_binx(xsh_instrument *instrument)
cpl_error_code xsh_model_config_load_best(cpl_frame *config_frame, xsh_xs_3 *p_xs_3)
Load the config model table and fill the struct.
Definition: xsh_model_io.c:174
void xsh_model_get_xy(xsh_xs_3 *p_xs_3, xsh_instrument *instr, double lambda_nm, int morder, double ent_slit_pos, double *x, double *y)
Compute the detector location (floating point pixels) of a given wavelength/entrance slit position.
void xsh_model_binxy(xsh_xs_3 *p_xs_3, int bin_X, int bin_Y)
corrects model for detector's binning
cpl_error_code xsh_model_maps_create(xsh_xs_3 *p_xs_3, xsh_instrument *instr, const char *wtag, const char *stag, cpl_frame **wmap_frame, cpl_frame **smap_frame, const int save_tmp)
Compute the wavelength and slit maps.
int size
int * y
int * x
#define xsh_msg_warning(...)
Print an warning message.
Definition: xsh_msg.h:88
#define xsh_msg_dbg_medium(...)
Definition: xsh_msg.h:44
#define xsh_msg(...)
Print a message on info level.
Definition: xsh_msg.h:121
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
int xsh_debug_level_get(void)
get debug level
Definition: xsh_utils.c:3142
cpl_error_code xsh_update_pheader_in_image_multi(cpl_frame *frame, const cpl_propertylist *pheader)
Update FITS header.
Definition: xsh_utils.c:4258
void xsh_free_frame(cpl_frame **f)
Deallocate a frame and set the pointer to NULL.
Definition: xsh_utils.c:2269
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
void xsh_add_temporary_file(const char *name)
Add temporary file to temprary files list.
Definition: xsh_utils.c:1432
cpl_error_code xsh_wavetab_qc(cpl_frame *frm_tab, const int is_poly)
Monitor min/max/med/avg distance between detected lines on each ordee.
cpl_frame * xsh_create_dispersol_physmod(cpl_frame *pre_frame, cpl_frame *order_tab_frame, cpl_frame *mod_cfg_frame, cpl_frame *wave_map_frame, cpl_frame *slit_map_frame, xsh_dispersol_param *dispsol_param, cpl_frame *spectral_format_frame, xsh_instrument *instrument, const int clean_tmp)
cpl_error_code xsh_wavemap_qc(cpl_frame *frm_map, const cpl_frame *frm_tab)
Monitor Flux level along the orders traces given by an input table.
void xsh_create_map(cpl_frame *dispsol_frame, cpl_frame *ordertab_frame, cpl_frame *pre_frame, xsh_instrument *instrument, cpl_frame **wavemap_frame, cpl_frame **slitmap_frame, const char *rec_prefix)
void xsh_create_model_map(cpl_frame *model_frame, xsh_instrument *instrument, const char *wtag, const char *stag, cpl_frame **wavemap_frame, cpl_frame **slitmap_frame, const int save_tmp)
cpl_frame * xsh_create_poly_wavemap(cpl_frame *pre_frame, cpl_frame *wave_tab_2d_frame, cpl_frame *order_tab_frame, cpl_frame *spectral_format_frame, xsh_dispersol_param *dispsol_par, xsh_instrument *instrument, const char *wm_tag, cpl_frame **dispersol_frame, cpl_frame **slitmap_frame)
#define XSH_SLIT_MAP_POLY
Definition: xsh_dfs.h:151
#define XSH_WAVE_MAP_POLY
Definition: xsh_dfs.h:918
#define XSH_GET_TAG_FROM_ARM(TAG, instr)
Definition: xsh_dfs.h:1548
#define XSH_QC_LINE_DIFMIN_ORD
Definition: xsh_pfits_qc.h:94
#define XSH_QC_WMAP_WAVEC
Definition: xsh_pfits_qc.h:92
#define XSH_QC_LINE_DIFMAX_C
Definition: xsh_pfits_qc.h:100
#define XSH_QC_LINE_DIFMED_C
Definition: xsh_pfits_qc.h:102
#define XSH_QC_LINE_DIFAVG
Definition: xsh_pfits_qc.h:103
#define XSH_QC_LINE_DIFMIN
Definition: xsh_pfits_qc.h:97
#define XSH_QC_LINE_DIFMAX
Definition: xsh_pfits_qc.h:99
#define XSH_QC_LINE_DIFMAX_ORD
Definition: xsh_pfits_qc.h:95
#define XSH_QC_LINE_DIFAVG_C
Definition: xsh_pfits_qc.h:104
#define XSH_QC_LINE_DIFMED
Definition: xsh_pfits_qc.h:101
#define XSH_QC_LINE_DIFMIN_C
Definition: xsh_pfits_qc.h:98
#define XSH_FREE(POINTER)
Definition: xsh_utils.h:92
@ XSH_DEBUG_LEVEL_MEDIUM
Definition: xsh_utils.h:138
#define XSH_MALLOC(POINTER, TYPE, SIZE)
Definition: xsh_utils.h:49