X-shooter Pipeline Reference Manual 3.8.15
test-xsh_data_spectrum_merge_3d.c
Go to the documentation of this file.
1/* *
2 * This file is part of the ESO X-shooter Pipeline *
3 * Copyright (C) 2006 European Southern Observatory *
4 * *
5 * This library is free software; you can redistribute it and/or modify *
6 * it under the terms of the GNU General Public License as published by *
7 * the Free Software Foundation; either version 2 of the License, or *
8 * (at your option) any later version. *
9 * *
10 * This program is distributed in the hope that it will be useful, *
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
13 * GNU General Public License for more details. *
14 * *
15 * You should have received a copy of the GNU General Public License *
16 * along with this program; if not, write to the Free Software *
17 * Foundation, 51 Franklin St, Fifth Floor, Boston, MA 02111-1307 USA *
18 * */
19
20/*
21 * $Author: amodigli $
22 * $Date: 2011-12-02 14:13:14 $
23 * $Revision: 1.9 $
24 */
25
26#ifdef HAVE_CONFIG_H
27# include <config.h>
28#endif
29
30/*-------------------------------------------------------------------------*/
35/*-------------------------------------------------------------------------*/
38/*--------------------------------------------------------------------------
39 Includes
40 --------------------------------------------------------------------------*/
41
42#include <tests.h>
43#include <xsh_data_pre.h>
44#include <xsh_error.h>
45#include <xsh_msg.h>
46#include <xsh_data_instrument.h>
47#include <xsh_data_rec.h>
49#include <xsh_drl.h>
50#include <xsh_pfits.h>
51#include <xsh_data_pre_3d.h>
52#include <xsh_badpixelmap.h>
53
54#include <cpl.h>
55#include <math.h>
56
57#include <getopt.h>
58
59/*--------------------------------------------------------------------------
60 Defines
61 --------------------------------------------------------------------------*/
62
63#define MODULE_ID "XSH_DATA_SPECTRUM_MERGE_3D"
64
65static const char * Options = "" ;
66
67enum {
69};
70
71static struct option long_options[] = {
72 {"lambda-chunk-hsize", required_argument, 0, ZCHUNK_HSIZE_OPT},
73 {0, 0, 0, 0}
74};
75
76static void Help( void )
77{
78 puts( "Creates for each slitlet:");
79 puts( " - a median profile along the slit");
80 puts( " - a localization as a function of wavelength");
81
82 puts( "Usage: test_xsh_data_spectrum_merge_3d [options] CUBE SKY_LINE");
83
84 puts( "Options" ) ;
85 puts( " --lambda-chunk-hsize=<n> : Chunk half size in pixel [50]");
86
87 puts( "\nInput Files" ) ;
88 puts( "CUBE : data cube (tag = *_IFU_MERGE3D_IFU_UVB)" );
89 puts( "SKY_LINE [OPTIONAL] : sky lines (tag = SKY_LINE_LIST_arm)" );
90
91 puts( "\nOutput Files");
92 puts( "1) MEDIAN_CUBE.fits : median profile in fits image");
93 puts( "2) median_profile_slitbin_*.dat : ASCII profile along the slit direction for every slitlet by collapsing in wavelength");
94 puts( "3) profile_lambda_*.dat : ASCII file containing for each chunk around lref gaussian central fitting positions");
95 puts( "4) gsl_locfit_*.dat : ASCII central slit position for each slitlet");
96 TEST_END();
97}
98
99
100static void HandleOptions( int argc, char **argv,
101 int *zchunk_hsize)
102{
103 int opt ;
104 int option_index = 0;
105
106 while (( opt = getopt_long (argc, argv, Options,
107 long_options, &option_index)) != EOF )
108 switch ( opt ) {
109 case ZCHUNK_HSIZE_OPT:
110 *zchunk_hsize = atoi( optarg);
111 break;
112 default: Help() ; exit( 0);
113 }
114}
115/*--------------------------------------------------------------------------
116 Implementation
117 --------------------------------------------------------------------------*/
118
126int main( int argc, char **argv)
127{
128 int ret;
129 /* Declarations */
130 const char* sci_name = NULL;
131 cpl_frame *sci_frame = NULL;
132
133 const char* sky_name = NULL;
134 cpl_frame *sky_frame = NULL;
135
136 xsh_pre_3d *cube = NULL;
137 xsh_image_3d *cube_data_img = NULL;
138 int cube_x, cube_y, cube_z;
139 float *cube_data = NULL;
140 double *median = NULL;
141 int median_size =0;
142 double *coadd =NULL;
143 int zchunk_hsize = 100;
144 cpl_vector *med_vect = NULL;
145 cpl_vector *coadd_vect = NULL;
146 cpl_vector *positions = NULL;
147 int x,y,z;
148 double med_z;
149 cpl_image *med_img = NULL;
150 double *med_img_data = NULL;
151 FILE *profil_file = NULL;
152 char filename[256];
153 int zmed;
154 cpl_propertylist *header = NULL;
155 double crval1, crval2, crval3;
156 double crpix1, crpix2, crpix3;
157 double cdelt1, cdelt2, cdelt3;
158 double x0= 0.0,sigma =0.0, area=0.0, offset=0.0;
159 int bad_fit;
161
162 /* Initialize libraries */
164 cpl_msg_set_level( CPL_MSG_DEBUG);
166
167 HandleOptions( argc, argv, &zchunk_hsize);
168
169 /* Analyse parameters */
170 if (argc-optind >= 1 ) {
171 sci_name = argv[optind];
172 if (argc-optind >= 2 ) {
173 sky_name = argv[optind+1];
174 }
175 }
176 else{
177 Help();
178 exit(0);
179 }
180
181 TESTS_XSH_FRAME_CREATE( sci_frame, "CUBE_arm", sci_name);
182 if ( sky_name != NULL){
183 TESTS_XSH_FRAME_CREATE( sky_frame, "SKY_LINE_LIST_arm", sky_name);
184 }
185
186 /* USE load function */
187 xsh_msg("CUBE : %s",
188 cpl_frame_get_filename( sci_frame));
189 if ( sky_frame != NULL){
190 xsh_msg("SKY LINE LIST : %s",
191 cpl_frame_get_filename( sky_frame));
192 }
193 else{
194 xsh_msg("Don't use sky line list!!");
195 }
196
197 check( instrument = create_instrument( sci_name));
199
200 check( cube = xsh_pre_3d_load( sci_frame));
201 header = cube->data_header;
202
203 check( cube_data_img = xsh_pre_3d_get_data( cube));
204 check( cube_x = xsh_image_3d_get_size_x( cube_data_img));
205 check( cube_y = xsh_image_3d_get_size_y( cube_data_img));
206 check( cube_z = xsh_image_3d_get_size_z( cube_data_img));
207
208 xsh_msg(" Cube SLITLET*SLIT*LAMBDAS (%dx%dx%d)", cube_x, cube_y, cube_z);
209
210 check( crpix1 = xsh_pfits_get_crpix1( header));
211 check( crval1 = xsh_pfits_get_crval1( header));
212 check( cdelt1 = xsh_pfits_get_cdelt1( header));
213
214 check( crpix2 = xsh_pfits_get_crpix2( header));
215 check( crval2 = xsh_pfits_get_crval2( header));
216 check( cdelt2 = xsh_pfits_get_cdelt2( header));
217
218 check( crpix3 = xsh_pfits_get_crpix3( header));
219 check( crval3 = xsh_pfits_get_crval3( header));
220 check( cdelt3 = xsh_pfits_get_cdelt3( header));
221
222 xsh_msg(" SLITLET pix %f ref %f with step %f", crpix1, crval1, cdelt1);
223 xsh_msg(" SLIT pix %f ref %f with step %f", crpix2, crval2, cdelt2);
224 xsh_msg(" LAMBDAS pix %f ref %f with step %f", crpix3, crval3, cdelt3);
225
226 check( cube_data = (float*)xsh_image_3d_get_data( cube_data_img));
227
228 /* test1 */
229 XSH_MALLOC( median, double, cube_z);
230 check( med_vect = cpl_vector_wrap( cube_z, median));
231 check( med_img = cpl_image_new( cube_x, cube_y, CPL_TYPE_DOUBLE));
232 check( med_img_data = cpl_image_get_data_double( med_img));
233
234 for( y=0; y<cube_y; y++){
235 for( x=0; x<cube_x; x++){
236 for(z=0; z< cube_z; z++){
237 median[z] = cube_data[cube_x*cube_y*z+cube_x*y+x];
238 }
239 check (med_z = cpl_vector_get_median( med_vect));
240 med_img_data[y*cube_x+x] = med_z;
241 }
242 }
243 check( cpl_image_save( med_img, "MEDIAN_CUBE.fits", CPL_BPP_IEEE_DOUBLE,
244 NULL, CPL_IO_CREATE));
245
246 sprintf( filename, "median_profile_slitbin_%.3f.dat", cdelt2);
247 profil_file = fopen( filename, "w");
248 fprintf( profil_file, "#slit UP CEN LO\n");
249
250 for( y=0; y<cube_y; y++){
251 fprintf( profil_file, "%f ", crval2+y*cdelt2);
252 for( x=0; x<cube_x; x++){
253 fprintf( profil_file, "%f ", med_img_data[y*cube_x+x]);
254 }
255 fprintf( profil_file, "\n");
256 }
257 fclose( profil_file);
258
259 XSH_MALLOC( coadd, double, cube_y);
260 check( coadd_vect = cpl_vector_wrap( cube_y, coadd));
261 check( positions = cpl_vector_new( cube_y));
262
263 for( y=0; y<cube_y; y++){
264 cpl_vector_set( positions, y, y);
265 cpl_vector_set( coadd_vect, y, med_img_data[y*cube_x+1]);
266 }
267 cpl_vector_fit_gaussian( positions, NULL, coadd_vect, NULL,CPL_FIT_ALL,&x0,&sigma,&area,&offset,NULL,NULL,NULL);
268
269 if (cpl_error_get_code() != CPL_ERROR_NONE){
271 }
272 else{
273 sprintf( filename, "gaussian_fit.dat");
274 profil_file = fopen( filename, "w");
275 fprintf( profil_file, "#slit cen gauss_cen\n");
276
277 for( y=0; y<cube_y; y++){
278 double sval = 0.0;
279 double cen_val = 0.0;
280 double gauss_val = 0.0;
281 double z;
282
283 z = (y-x0)/(sigma*XSH_MATH_SQRT_2);
284 gauss_val = area / sqrt(2 *M_PI*sigma*sigma)*exp(-(z*z))+offset;
285 sval = crval2+y*cdelt2;
286 cen_val = cpl_vector_get( coadd_vect, y);
287
288 fprintf( profil_file, "%f %f %f\n", sval, cen_val, gauss_val);
289 }
290
291 fclose( profil_file);
292 xsh_msg(" Median cube central slitlet gaussian fit %f arcsec\n", crval2+x0*cdelt2);
293 }
294
295 /* test 2 */
296 sprintf( filename, "profile_lambda_%d.dat", zchunk_hsize);
297
298 profil_file = fopen( filename, "w");
299
300 median_size = 0;
301 bad_fit = 0;
302 xsh_unwrap_vector( &med_vect);
303
304 fprintf( profil_file, "#lref up cen low\n");
305 for(z=zchunk_hsize; z< cube_z; z+=zchunk_hsize*2){
306 double cen[3];
307
308 for( x=0; x<cube_x; x++){
309
310 for( y=0; y<cube_y; y++){
311 coadd[y]=0;
312 }
313 for(zmed=z-zchunk_hsize; zmed< z+zchunk_hsize; zmed++){
314 for( y=0; y<cube_y; y++){
315 coadd[y] += cube_data[cube_x*cube_y*zmed+cube_x*y+x];
316 }
317 }
318 cpl_vector_fit_gaussian( positions, NULL, coadd_vect, NULL,CPL_FIT_ALL,&x0,&sigma,&area,&offset,NULL,NULL,NULL);
319 if (cpl_error_get_code() != CPL_ERROR_NONE){
321 cen[x]=-1;
322 if (x == 1){
323 bad_fit++;
324 }
325 }
326 else{
327 cen[x] = crval2+x0*cdelt2;
328 if (x == 1){
329 median[median_size] = cen[1];
330 median_size++;
331 }
332 }
333 }
334 fprintf( profil_file, "%f %f %f %f\n", crval3+z*cdelt3, cen[0], cen[1], cen[2]);
335 }
336 fclose( profil_file);
337 check( med_vect = cpl_vector_wrap( median_size, median));
338 med_z = cpl_vector_get_median( med_vect);
339 xsh_msg(" Statistics of gaussian fit bad %d lines, good %d lines",
340 bad_fit, median_size);
341 xsh_msg(" Median from gaussian fit %f arcsec\n", med_z);
342
343 /* GSL loc fit */
344 for( x=0; x< cube_x; x++){
345 double init_par[6];
346 double fit_errs[6];
347 int status;
348 float slitcen, fit_err;
349
350 sprintf( filename, "gsl_localize_pos%d.dat", x);
351 profil_file = fopen( filename, "w");
352 fprintf( profil_file, "#wavelength slit_fit fit_err\n");
353
354 for( z=0; z < cube_z; z++){
355 for( y=0; y < cube_y; y++){
356 cpl_vector_set( coadd_vect, y, cube_data[cube_x*cube_y*z+cube_x*y+x]);
357 }
358 xsh_gsl_init_gaussian_fit( positions, coadd_vect, init_par);
359 xsh_gsl_fit_gaussian( positions, coadd_vect, 0, init_par, fit_errs, &status);
360 slitcen = crval2+init_par[4]*cdelt2;
361 fit_err = fit_errs[4]*cdelt2;
362 fprintf( profil_file, "%f %f %f\n", crval3+z*cdelt3, slitcen, fit_err);
363 }
364 fclose( profil_file);
365 }
366
367 check( xsh_center_cube( sci_frame, sky_frame, zchunk_hsize*2,
368 instrument));
369
370 cleanup:
371 if (cpl_error_get_code() != CPL_ERROR_NONE) {
372 xsh_error_dump(CPL_MSG_ERROR);
373 ret = 1;
374 }
376 xsh_free_frame( &sci_frame);
377 xsh_pre_3d_free( &cube);
378 XSH_FREE( median);
379 XSH_FREE( coadd);
380 xsh_unwrap_vector( &coadd_vect);
381 xsh_unwrap_vector( &med_vect);
382 xsh_free_image( &med_img);
383 TEST_END();
384 return ret;
385}
386
int main()
Unit test of xsh_bspline_interpol.
static double sigma
static xsh_instrument * instrument
static const char * Options
static void Help(void)
static void HandleOptions(int argc, char **argv, int *zchunk_hsize)
static struct option long_options[]
xsh_instrument * create_instrument(const char *filename)
Definition: tests.c:622
void * xsh_image_3d_get_data(xsh_image_3d *img_3d)
int xsh_image_3d_get_size_x(xsh_image_3d *img_3d)
int xsh_image_3d_get_size_y(xsh_image_3d *img_3d)
int xsh_image_3d_get_size_z(xsh_image_3d *img_3d)
void xsh_pre_3d_free(xsh_pre_3d **pre_3d)
xsh_pre_3d * xsh_pre_3d_load(cpl_frame *frame)
Load a xsh_pre_3d structure from a frame.
xsh_image_3d * xsh_pre_3d_get_data(xsh_pre_3d *pre_3d)
Get data.
#define check(COMMAND)
Definition: xsh_error.h:71
#define xsh_error_dump(level)
Definition: xsh_error.h:92
#define xsh_error_reset()
Definition: xsh_error.h:87
void xsh_center_cube(cpl_frame *cube_frame, cpl_frame *sky_line_frame, int chunk_size, xsh_instrument *instrument)
Shift a cube to center object at 0 arcsec.
Definition: xsh_format.c:343
void xsh_instrument_set_mode(xsh_instrument *i, XSH_MODE mode)
Set a mode on instrument structure.
void xsh_instrument_free(xsh_instrument **instrument)
free an instrument structure
int * y
int * x
#define xsh_msg(...)
Print a message on info level.
Definition: xsh_msg.h:121
double xsh_pfits_get_crval2(const cpl_propertylist *plist)
find out the crval2
Definition: xsh_pfits.c:1926
double xsh_pfits_get_cdelt2(const cpl_propertylist *plist)
find out the cdelt2
Definition: xsh_pfits.c:2216
double xsh_pfits_get_cdelt3(const cpl_propertylist *plist)
find out the cdelt3
Definition: xsh_pfits.c:2235
double xsh_pfits_get_crpix1(const cpl_propertylist *plist)
find out the CRPIX1 value
Definition: xsh_pfits.c:1965
double xsh_pfits_get_crpix3(const cpl_propertylist *plist)
find out the CRPIX3 value
Definition: xsh_pfits.c:1999
double xsh_pfits_get_cdelt1(const cpl_propertylist *plist)
find out the cdelt1
Definition: xsh_pfits.c:2196
double xsh_pfits_get_crpix2(const cpl_propertylist *plist)
find out the CRPIX2 value
Definition: xsh_pfits.c:1982
double xsh_pfits_get_crval1(const cpl_propertylist *plist)
find out the crval1
Definition: xsh_pfits.c:1907
double xsh_pfits_get_crval3(const cpl_propertylist *plist)
find out the crval3
Definition: xsh_pfits.c:1946
void xsh_unwrap_vector(cpl_vector **v)
Unwrap a vector and set the pointer to NULL.
Definition: xsh_utils.c:2345
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
int xsh_debug_level_set(int level)
set debug level
Definition: xsh_utils.c:3125
void xsh_gsl_init_gaussian_fit(cpl_vector *xpos_vect, cpl_vector *ypos_vect, double *init_par)
Definition: xsh_utils.c:6862
void xsh_gsl_fit_gaussian(cpl_vector *xpos_vect, cpl_vector *ypos_vect, int deg, double *params, double *errs, int *status)
Definition: xsh_utils.c:6936
cpl_propertylist * data_header
#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
@ XSH_MODE_IFU
#define XSH_MATH_SQRT_2
Definition: xsh_drl.h:567
#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
#define M_PI
Definition: xsh_utils.h:43