X-shooter Pipeline Reference Manual 3.8.15
test-xsh_geom_ifu.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-01-16 21:09:42 $
23 * $Revision: 1.11 $
24 */
25
26#ifdef HAVE_CONFIG_H
27# include <config.h>
28#endif
29
30/*-------------------------------------------------------------------------*/
36/*-------------------------------------------------------------------------*/
39/*--------------------------------------------------------------------------
40 Includes
41 --------------------------------------------------------------------------*/
42
43#include <tests.h>
44
45#include <xsh_data_pre.h>
46#include <xsh_error.h>
47#include <xsh_msg.h>
48#include <xsh_data_instrument.h>
49#include <xsh_data_the_map.h>
52#include <xsh_utils_table.h>
53
54#include <xsh_drl.h>
55#include <xsh_pfits.h>
56
57#include <xsh_badpixelmap.h>
58
59#include <cpl.h>
60#include <math.h>
61
62#include <getopt.h>
63
64/*--------------------------------------------------------------------------
65 Defines
66 --------------------------------------------------------------------------*/
67
68#define MODULE_ID "XSH_GEOM_IFU"
69
70/*--------------------------------------------------------------------------
71 Implementation
72 --------------------------------------------------------------------------*/
73enum {
77} ;
78
79static const char * Options = "" ;
80
81static struct option long_options[] = {
82 {"kernel", required_argument, 0, KERNEL_OPT},
83 {"radius", required_argument, 0, RADIUS_OPT},
84 {"bin-lambda", required_argument, 0, BIN_LAMBDA_OPT},
85 {"bin-space", required_argument, 0, BIN_SPACE_OPT},
86 {"order-min", required_argument, 0, MIN_ORDER_OPT},
87 {"order-max", required_argument, 0, MAX_ORDER_OPT},
88 {"slit-min",required_argument, 0, SLIT_MIN_OPT},
89 {"slit-n",required_argument, 0, NSLIT_OPT},
90 {"merge-method", required_argument, 0, MERGE_METHOD_OPT},
91 {"lambda-ref", required_argument, 0, LAMBDA_REF_OPT},
92 {"help", 0, 0, HELP_OPT},
93 {0, 0, 0, 0}
94};
95
96static void Help( void )
97{
98 puts( "Unitary test of xsh_geom_ifu" ) ;
99 puts( "Usage: test_xsh_geom_ifu [options] <input_files>" ) ;
100 puts( "Options" ) ;
101 puts( " --kernel=<name> : TANH, SINC, SINC2, LANCZOS, HAMMING, HANN [TANH]");
102 puts( " --radius=<nn> : Radius [4]" ) ;
103 puts( " --bin-lambda=<n> : Bin in Lambda [0.1]" ) ;
104 puts( " --bin-space=<n> : Bin in Slit [0.1]" ) ;
105 puts( " --order-min=<n> : Minimum abs order" );
106 puts( " --order-max=<n> : Maximum abs order" );
107 puts( " --slit-min=<n> : Minimum slit to rectify" );
108 puts( " --slit-n=<n> : Number of pixels in slit rectified frame" );
109 puts( " --lambda-ref=<n> : Lambda reference for creating SLICE OFFSET");
110 puts( " --help : What you see" ) ;
111 puts( "\nInput Files" ) ;
112 puts( "The input files argument MUST be in this order:" ) ;
113 puts( " 1. Science frame in PRE format" ) ;
114 puts( " 2. SOF\n" ) ;
115 TEST_END();
116}
117
118static void HandleOptions( int argc, char **argv,
119 xsh_rectify_param *rectify_par, xsh_localize_obj_param *locobj_par,
120 int *order_min, int *order_max,
121 double *slit_min, int *nslit, double *lambda_ref)
122{
123 int opt ;
124 int option_index = 0;
125
126 while (( opt = getopt_long (argc, argv, Options,
127 long_options, &option_index)) != EOF )
128 switch ( opt ) {
129 case KERNEL_OPT:
130 strcpy( rectify_par->rectif_kernel, optarg);
131 if ( strcmp( optarg, RECTIFY_KERNEL_PRINT( CPL_KERNEL_TANH)) == 0){
132 rectify_par->kernel_type = CPL_KERNEL_TANH;
133 }
134 else if ( strcmp( optarg, RECTIFY_KERNEL_PRINT( CPL_KERNEL_TANH)) == 0){
135 rectify_par->kernel_type = CPL_KERNEL_TANH;
136 }
137 else if ( strcmp( optarg, RECTIFY_KERNEL_PRINT( CPL_KERNEL_SINC)) == 0){
138 rectify_par->kernel_type = CPL_KERNEL_SINC;
139 }
140 else if ( strcmp( optarg, RECTIFY_KERNEL_PRINT( CPL_KERNEL_SINC2)) == 0){
141 rectify_par->kernel_type = CPL_KERNEL_SINC2;
142 }
143 else if ( strcmp( optarg, RECTIFY_KERNEL_PRINT( CPL_KERNEL_LANCZOS)) == 0){
144 rectify_par->kernel_type = CPL_KERNEL_LANCZOS;
145 }
146 else if ( strcmp( optarg, RECTIFY_KERNEL_PRINT( CPL_KERNEL_HAMMING)) == 0){
147 rectify_par->kernel_type = CPL_KERNEL_HAMMING;
148 }
149 else if ( strcmp( optarg, RECTIFY_KERNEL_PRINT( CPL_KERNEL_HANN)) == 0){
150 rectify_par->kernel_type = CPL_KERNEL_HANN;
151 }
152 else{
153 xsh_msg("Invalid kernel option %s", optarg);
154 Help();
155 exit(0);
156 }
157 break;
158 case RADIUS_OPT:
159 rectify_par->rectif_radius = atof( optarg);
160 break ;
161 case BIN_LAMBDA_OPT:
162 sscanf( optarg, "%64lf", &rectify_par->rectif_bin_lambda ) ;
163 break ;
164 case BIN_SPACE_OPT:
165 sscanf( optarg, "%64lf", &rectify_par->rectif_bin_space ) ;
166 break ;
167 case MIN_ORDER_OPT:
168 sscanf( optarg, "%64d", order_min);
169 break;
170 case MAX_ORDER_OPT:
171 sscanf( optarg, "%64d", order_max);
172 break;
173 case SLIT_MIN_OPT:
174 sscanf( optarg, "%64lf", slit_min) ;
175 break ;
176 case LAMBDA_REF_OPT:
177 sscanf( optarg, "%64lf", lambda_ref);
178 break ;
179 case NSLIT_OPT:
180 sscanf( optarg, "%64d", nslit);
181 break ;
182 default: Help() ; exit( 0 ) ;
183 }
184}
185
194int main( int argc, char **argv)
195{
196 /* Declarations */
197 int ret = 0 ;
198
200
201 const char *sof_name = NULL;
202 cpl_frameset *set = NULL;
203 const char * sci_name = NULL ;
204 xsh_localize_obj_param loc_obj_par = {10, 0.1, 1, 0, LOC_GAUSSIAN_METHOD, 0, 0.5, 2.0, 3, FALSE};;
205
206 cpl_frameset *wavesol_frameset = NULL;
207 cpl_frameset *rect_frameset = NULL;
208 cpl_frameset *loc_frameset = NULL;
209 cpl_frame *sci_frame = NULL;
210 cpl_frame *orderlist_frame = NULL;
211 cpl_frame *slitmap_frame = NULL;
212 cpl_frame *model_frame = NULL;
213 cpl_frame *spectralformat_frame = NULL;
214 cpl_frame *recorder_frame = NULL ;
215 cpl_frame *recordereso_frame = NULL ;
216 cpl_frame *recordertab_frame = NULL ;
217 cpl_frame *cube = NULL;
218 int order_min = -1, order_max = -1;
219 int rec_min_index=-1, rec_max_index=-1;
220 double lambda_ref = -1;
221 xsh_order_list* order_list = NULL;
222 xsh_rectify_param rectify_par;
223 int nslit_c, nslit=-100;
224 double slit_min_c, slit_min=-100;
226 const int decode_bp=2147483647;
227 /* Initialize libraries */
229 cpl_msg_set_level(CPL_MSG_DEBUG);
231
232 /* Set rectify default params */
233 strcpy( rectify_par.rectif_kernel, "default");
234 rectify_par.kernel_type = CPL_KERNEL_TANH;
235 rectify_par.rectif_radius = 4 ;
236 rectify_par.rectif_bin_lambda = 0.1 ;
237 rectify_par.rectif_bin_space = 0.1 ;
238 rectify_par.conserve_flux = FALSE;
239
240 loc_obj_par.loc_deg_poly = 0;
241
242 HandleOptions( argc, argv, &rectify_par, &loc_obj_par, &order_min, &order_max,
243 &slit_min, &nslit, &lambda_ref);
244
245 if ( (argc - optind) >=2 ) {
246 sci_name = argv[optind];
247 sof_name = argv[optind+1];
248 }
249 else {
250 Help();
251 exit(0);
252 }
253 /* Create frameset from sof */
254 check( set = sof_to_frameset( sof_name));
255 /* Validate frame set */
258
259 check( spectralformat_frame = xsh_find_spectral_format( set, instrument));
260 check( orderlist_frame = xsh_find_order_tab_edges( set, instrument));
261 check( slitmap_frame = xsh_find_slitmap( set, instrument));
262
263 TESTS_XSH_FRAME_CREATE( sci_frame, "SLIT_STARE_DIVIDED_arm", sci_name);
265
266 xsh_msg("SCI : %s",
267 cpl_frame_get_filename( sci_frame));
268 xsh_msg("ORDERLIST : %s",
269 cpl_frame_get_filename( orderlist_frame));
270 xsh_msg("SPECTRALFORMAT : %s",
271 cpl_frame_get_filename( spectralformat_frame));
272 if ( slitmap_frame != NULL){
273 xsh_msg("SLITMAP : %s",
274 cpl_frame_get_filename( slitmap_frame));
275 }
276
277 wavesol_frameset = xsh_find_wave_tab_ifu( set, instrument);
278 cpl_error_reset();
279 if ( wavesol_frameset == NULL){
280 model_frame = xsh_find_model_config_tab( set, instrument);
281 xsh_msg("MODEL : %s",
282 cpl_frame_get_filename( model_frame));
283 }
284 else{
285 int iframe;
286 int size;
287
288 size = cpl_frameset_get_size( wavesol_frameset);
289 for(iframe=0; iframe < size; iframe++){
290 cpl_frame *frame = cpl_frameset_get_frame(wavesol_frameset, iframe);
291 xsh_msg("WAVESOL : %s",
292 cpl_frame_get_filename( frame));
293 }
294 }
295
296 xsh_msg(" Parameters ");
297 xsh_msg(" kernel %s", RECTIFY_KERNEL_PRINT( rectify_par.kernel_type));
298 xsh_msg(" radius %f", rectify_par.rectif_radius);
299 xsh_msg(" bin-space %f", rectify_par.rectif_bin_space);
300 xsh_msg(" bin-lambda %f", rectify_par.rectif_bin_lambda);
301
302
303 check( order_list = xsh_order_list_load ( orderlist_frame, instrument));
304
305 if ( order_min != -1) {
306 check( rec_min_index = xsh_order_list_get_index_by_absorder( order_list,
307 order_min));
308 xsh_msg("Order min %d => index %d", order_min, rec_min_index);
309 }
310 else{
311 rec_min_index = 0;
312 }
313
314 if ( order_max != -1) {
315 check( rec_max_index = xsh_order_list_get_index_by_absorder( order_list,
316 order_max));
317 xsh_msg("Order max %d => index %d", order_max, rec_max_index);
318 }
319 else{
320 rec_max_index = order_list->size;
321 }
322
323
324 /* Now rectify the frame */
325
326 check( xsh_rec_slit_size( &rectify_par,
327 &slit_min_c, &nslit_c, mode));
328
329 if (slit_min == -100){
330 slit_min = slit_min_c;
331 }
332 if ( nslit == -100){
333 nslit = nslit_c;
334 }
335
336 xsh_msg("SLIT min = %f and has size %d", slit_min, nslit);
337 xsh_msg("Lambda ref %f", lambda_ref);
338 if (lambda_ref < 0.0){
339 check( rect_frameset = xsh_rectify_orders_ifu( sci_frame, order_list,
340 wavesol_frameset, NULL, model_frame, instrument, &rectify_par,
341 spectralformat_frame, slitmap_frame,NULL,NULL,rec_min_index,
342 rec_max_index, "test"));
343 }
344 else{
345 double lambda_ref_min;
346 double lambda_ref_max;
347 const char * tablename = NULL;
348 cpl_table *table = NULL;
349 cpl_frame *spectralformat2_frame = NULL;
350 xsh_spectralformat_list* list = NULL;
351 cpl_vector* orders = NULL;
352 int size, absorder, iorder;
353
354 lambda_ref_min = lambda_ref-rectify_par.rectif_bin_lambda*10;
355 lambda_ref_max = lambda_ref+rectify_par.rectif_bin_lambda*10;
356
357 xsh_msg( "Do cube for reference %f ==> %f", lambda_ref_min, lambda_ref_max);
358 check( list = xsh_spectralformat_list_load( spectralformat_frame,
359 instrument));
360 check( orders = xsh_spectralformat_list_get_orders( list, lambda_ref));
361 size = cpl_vector_get_size( orders);
362
363 tablename = cpl_frame_get_filename( spectralformat_frame);
364 XSH_TABLE_LOAD( table, tablename);
365
366 for( iorder=0; iorder < size; iorder++){
367 absorder = cpl_vector_get( orders, iorder);
368 xsh_msg("In Order %d", absorder);
369 cpl_table_set_int(table, XSH_SPECTRALFORMAT_TABLE_COLNAME_ORDER, iorder, absorder);
370 cpl_table_set_float(table, XSH_SPECTRALFORMAT_TABLE_COLNAME_WLMIN, iorder, lambda_ref_min);
371 cpl_table_set_float(table, XSH_SPECTRALFORMAT_TABLE_COLNAME_WLMAX, iorder, lambda_ref_max);
372 }
373 cpl_table_save( table, NULL,NULL,"LAMBDA_REF_FORMAT.fits",CPL_IO_DEFAULT);
374 TESTS_XSH_FRAME_CREATE( spectralformat2_frame, "SPEC_arm", "LAMBDA_REF_FORMAT.fits");
375
376 check( rect_frameset = xsh_rectify_orders_ifu( sci_frame, order_list,
377 wavesol_frameset, NULL, model_frame, instrument, &rectify_par,
378 spectralformat2_frame, slitmap_frame,NULL,NULL, 0,
379 size-1, "test"));
380 }
381 /* Localize object */
382 xsh_msg( "---Localize the objects for IFU");
383 check( loc_frameset = xsh_localize_obj_ifu( rect_frameset,
384 NULL, instrument, &loc_obj_par, NULL));
385
386 cleanup:
387 if (cpl_error_get_code() != CPL_ERROR_NONE) {
388 xsh_error_dump(CPL_MSG_ERROR);
389 ret = 1;
390 }
391 xsh_order_list_free( &order_list);
392 xsh_free_frame( &sci_frame);
393 xsh_free_frameset( &rect_frameset);
394 xsh_free_frameset( &loc_frameset);
395 xsh_free_frameset( &wavesol_frameset);
396 xsh_free_frameset( &set);
398 xsh_free_frame( &recorder_frame);
399 xsh_free_frame( &recordereso_frame);
400 xsh_free_frame( &recordertab_frame);
401 xsh_free_frame( &cube);
402 TEST_END();
403 return ret;
404}
int main()
Unit test of xsh_bspline_interpol.
static const char * Options
static void HandleOptions(int argc, char **argv, xsh_rectify_param *rectify_par, xsh_localize_obj_param *locobj_par, int *order_min, int *order_max, double *slit_min, int *nslit, double *lambda_ref)
#define MODULE_ID
static struct option long_options[]
@ KERNEL_OPT
@ HELP_OPT
@ MAX_ORDER_OPT
@ RADIUS_OPT
@ BIN_SPACE_OPT
@ MERGE_METHOD_OPT
@ NSLIT_OPT
@ MIN_ORDER_OPT
@ SLIT_MIN_OPT
@ BIN_LAMBDA_OPT
@ LAMBDA_REF_OPT
static char mode[32]
static xsh_instrument * instrument
cpl_frameset * sof_to_frameset(const char *sof_name)
Definition: tests.c:576
int xsh_order_list_get_index_by_absorder(xsh_order_list *list, double absorder)
xsh_order_list * xsh_order_list_load(cpl_frame *frame, xsh_instrument *instr)
load an order list from a frame
void xsh_order_list_free(xsh_order_list **list)
free memory associated to an order_list
cpl_vector * xsh_spectralformat_list_get_orders(xsh_spectralformat_list *list, float lambda)
Returns list of absolute orders containing lambda.
xsh_spectralformat_list * xsh_spectralformat_list_load(cpl_frame *frame, xsh_instrument *instr)
Load a spectralformat list from a frame.
#define check(COMMAND)
Definition: xsh_error.h:71
#define xsh_error_dump(level)
Definition: xsh_error.h:92
XSH_MODE xsh_instrument_get_mode(xsh_instrument *i)
Get a mode on instrument structure.
void xsh_instrument_free(xsh_instrument **instrument)
free an instrument structure
void xsh_instrument_set_decode_bp(xsh_instrument *i, const int decode_bp)
Set bad pixel code.
cpl_frameset * xsh_localize_obj_ifu(cpl_frameset *sci_frame, cpl_frame *skymask_frame, xsh_instrument *instrument, xsh_localize_obj_param *loc_obj_par, xsh_slit_limit_param *slit_limit_param)
int size
#define xsh_msg(...)
Print a message on info level.
Definition: xsh_msg.h:121
cpl_frameset * xsh_rectify_orders_ifu(cpl_frame *sci_frame, xsh_order_list *orderlist, cpl_frameset *wavesol_frameset, cpl_frameset *shift_frameset, cpl_frame *model_frame, xsh_instrument *instrument, xsh_rectify_param *rectify_par, cpl_frame *spectralformat_frame, cpl_frame *slitmap_frame, cpl_frameset **rec_frameset_ext, cpl_frameset **rec_frameset_tab, int min_index, int max_index, const char *rec_prefix)
Definition: xsh_rectify.c:1514
void xsh_rec_slit_size(xsh_rectify_param *rectify_par, double *slit_min, int *nslit, XSH_MODE mode)
rectify frame
Definition: xsh_rectify.c:758
void xsh_free_frame(cpl_frame **f)
Deallocate a frame and set the pointer to NULL.
Definition: xsh_utils.c:2269
void xsh_free_frameset(cpl_frameset **f)
Deallocate a frame set and set the pointer to NULL.
Definition: xsh_utils.c:2254
int xsh_debug_level_set(int level)
set debug level
Definition: xsh_utils.c:3125
cpl_kernel kernel_type
static void Help(void)
#define TESTS_XSH_FRAME_CREATE(frame, tag, name)
Definition: tests.h:123
#define TEST_END()
Definition: tests.h:111
#define TESTS_INIT(DRL_ID)
Definition: tests.h:105
#define XSH_SPECTRALFORMAT_TABLE_COLNAME_ORDER
#define XSH_SPECTRALFORMAT_TABLE_COLNAME_WLMIN
#define XSH_SPECTRALFORMAT_TABLE_COLNAME_WLMAX
cpl_frame * xsh_find_spectral_format(cpl_frameset *frames, xsh_instrument *instr)
Find spectral format frame.
Definition: xsh_dfs.c:4318
cpl_frame * xsh_find_order_tab_edges(cpl_frameset *frames, xsh_instrument *instr)
Find an order tab EDGES.
Definition: xsh_dfs.c:3595
cpl_frameset * xsh_find_wave_tab_ifu(cpl_frameset *frames, xsh_instrument *instrument)
Find wave tab ARC (for IFU 3 frames)
Definition: xsh_dfs.c:3756
cpl_frame * xsh_find_slitmap(cpl_frameset *frames, xsh_instrument *instr)
Find a slit map.
Definition: xsh_dfs.c:3673
xsh_instrument * xsh_dfs_set_groups(cpl_frameset *set)
Set the group as RAW or CALIB in a frameset and return the instrument detected.
Definition: xsh_dfs.c:1046
cpl_frame * xsh_find_model_config_tab(cpl_frameset *frames, xsh_instrument *instr)
Find a model configuration table frame.
Definition: xsh_dfs.c:3957
#define RECTIFY_KERNEL_PRINT(method)
@ LOC_GAUSSIAN_METHOD
@ XSH_DEBUG_LEVEL_MEDIUM
Definition: xsh_utils.h:138
#define XSH_TABLE_LOAD(TABLE, NAME)