X-shooter Pipeline Reference Manual 3.8.15
xsh_calibrate_flux.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-05-03 12:05:47 $
23 * $Revision: 1.6 $
24*/
25
26#ifdef HAVE_CONFIG_H
27#include <config.h>
28#endif
29
30/*----------------------------------------------------------------------------*/
37/*----------------------------------------------------------------------------*/
40/*-----------------------------------------------------------------------------
41 Includes
42 -----------------------------------------------------------------------------*/
43
44#include <math.h>
45#include <xsh_drl.h>
46
47#include <xsh_data_rec.h>
49#include <xsh_data_pre.h>
50#include <xsh_data_order.h>
51#include <xsh_dfs.h>
52#include <xsh_pfits.h>
53#include <xsh_error.h>
54#include <xsh_msg.h>
55#include <xsh_fit.h>
56#include <xsh_ifu_defs.h>
58#include <xsh_data_atmos_ext.h>
59#include <cpl.h>
60#include <xsh_utils.h>
61#include <gsl/gsl_integration.h>
62#include <xsh_data_star_flux.h>
63#include <xsh_data_spectrum.h>
64
65#define USE_SPLINE
66#if defined(USE_SPLINE)
67#include <gsl/gsl_spline.h>
68#endif
69
70/*-----------------------------------------------------------------------------
71 Functions prototypes
72 -----------------------------------------------------------------------------*/
73
74
75
76/*-----------------------------------------------------------------------------
77 Implementation
78 -----------------------------------------------------------------------------*/
79#if defined(USE_SPLINE)
80static gsl_interp_accel * AcceleratorResp, * AcceleratorAtmos ;
81static gsl_spline * SplineResp, * SplineAtmos ;
82
83static void init_interpolate( double * x, double * yf, int nb,
84 gsl_spline ** spline,
85 gsl_interp_accel ** accel)
86{
87 int ok ;
88
89 //fout = fopen( "interpol.dat", "w" ) ;
90
91 *accel = gsl_interp_accel_alloc();
92 if ( *accel == NULL ) xsh_msg( "Accelerator = NULL" ) ;
93
94 *spline = gsl_spline_alloc( gsl_interp_cspline, nb);
95
96 if ( *spline == NULL ) xsh_msg( "Spline = NULL" ) ;
97
98 ok = gsl_spline_init( *spline, x, yf, nb ) ;
99
100 xsh_msg( "gsl_spline_init --> %d", ok ) ;
101
102 return ;
103}
104
105static double do_interpolation( double x, gsl_spline * spline,
106 gsl_interp_accel * accel )
107{
108 double y =0;
109
110 y = gsl_spline_eval( spline, x, accel );
111
112 return y ;
113}
114
115static void clear_interpolate( void )
116{
117 gsl_spline_free( SplineResp ) ;
118 gsl_spline_free( SplineAtmos ) ;
119 gsl_interp_accel_free( AcceleratorResp ) ;
120 gsl_interp_accel_free( AcceleratorAtmos ) ;
121
122}
123#else
124
125static int get_idx_by_lambda( double lambda, double * a_lambda, int from,
126 int to, double step )
127{
128 int idx =0;
129 double * plambda=NULL ;
130
131 for( idx = from, plambda = a_lambda+from ; idx<to ; idx++, plambda++ ) {
132 if ( lambda >= (*plambda - step/2) &&
133 lambda < (*plambda + step/2) ) {
134 return idx ;
135 }
136 }
137
138 return -1 ;
139}
140
141#endif
142
143static double myfunc (double x, void * params) {
144 double alpha = *(double *) params;
145 double thefunc = exp(alpha*x*x) ;
146
147 return thefunc ;
148}
149
150static double compute_Lx( double slit_width, double seeing )
151{
152 gsl_integration_workspace * w
153 = gsl_integration_workspace_alloc (10000);
154
155 double result, error;
156 double alpha = -1.0;
157
158 gsl_function F;
159 double limit ;
160
161 double Lx ;
162
163 F.function = &myfunc ;
164 F.params = &alpha;
165 limit = (slit_width/2.)/((seeing*sqrt(2.))/2.35482) ;
166
167 gsl_integration_qags (&F, 0, limit, 0, 1e-7, 10000,
168 w, &result, &error);
169
170
171 Lx = 1- (2/sqrt(M_PI))*result ;
172
173 return Lx ;
174}
175
177 xsh_star_flux_list * response_list,
178 xsh_atmos_ext_list * atmos_ext_list,
179 double airmass_ratio,
180 double Lx )
181{
182 xsh_spectrum * spectrum =NULL;
183 double * spectrum_flux = NULL ;
184 int * spectrum_qual = NULL ;
185 double spectrum_lambda_min=0;
186 //double spectrum_lambda_max=0 ;
187 int spectrum_nlambda=0, spectrum_nslit=0 ;
188 double lambda=0 ;
189 double lambda_step=0 ;
190
191 double * atmos_lambda = NULL ;
192 double * atmos_K = NULL ;
193 //double atmos_step=0 ;
194 int atmos_size=0 ;
195
196 double * response_lambda = NULL ;
197 double * response_flux = NULL ;
198 //double response_step=0 ;
199 int response_size=0 ;
200 int ipix=0 ;
201#if !defined(USE_SPLINE)
202 int atmidx = 0, respidx = 0 ;
203#endif
204
205 check( spectrum = xsh_spectrum_duplicate( spectrum_in ) ) ;
206
207 check( spectrum_nlambda = xsh_spectrum_get_size_lambda( spectrum ) ) ;
208 check( spectrum_nslit = xsh_spectrum_get_size_slit( spectrum ) ) ;
209 if ( spectrum_nslit > 1 ) {
210 xsh_msg( "2D Spectrum" ) ;
211 }
212 else {
213 xsh_msg( "1D Spectrum" ) ;
214 }
215
216 check( spectrum_flux = xsh_spectrum_get_flux( spectrum ) ) ;
217 check( spectrum_qual = xsh_spectrum_get_qual( spectrum ) ) ;
218 check( spectrum_lambda_min = xsh_spectrum_get_lambda_min( spectrum ) ) ;
219 //check( spectrum_lambda_max = xsh_spectrum_get_lambda_max( spectrum ) ) ;
221
222 check( atmos_lambda = xsh_atmos_ext_list_get_lambda( atmos_ext_list ) ) ;
223 check( atmos_K = xsh_atmos_ext_list_get_K( atmos_ext_list ) ) ;
224 if (atmos_ext_list != NULL) {
225 atmos_size = atmos_ext_list->size ;
226 }
227 //atmos_step = *(atmos_lambda+atmos_size-1) - *atmos_lambda ;
228#if defined(USE_SPLINE)
229 check( init_interpolate( atmos_lambda, atmos_K, atmos_size, &SplineAtmos,
230 &AcceleratorAtmos ) ) ;
231#endif
232
233 check( response_lambda = xsh_star_flux_list_get_lambda( response_list ) ) ;
234 check( response_flux = xsh_star_flux_list_get_flux( response_list ) ) ;
235 response_size = response_list->size ;
236 //response_step = *(response_lambda+response_size-1) - *response_lambda ;
237#if defined(USE_SPLINE)
238 check( init_interpolate( response_lambda, response_flux, response_size,
239 &SplineResp,
240 &AcceleratorResp ) ) ;
241#endif
242
243 lambda = spectrum_lambda_min ;
244 /* Loop over spectrum pixels */
245 for( ipix = 0 ; ipix < spectrum_nlambda ; ipix++, lambda += lambda_step ) {
246 double kvalue, vresp ;
247 int islit ;
248
249 /* Now correct with atmos extinction: how ? interpolate (how) ?
250 other ? or just use the value corresponding to the current lambda ? */
251 if ( atmos_ext_list != NULL ) {
252#if defined(USE_SPLINE)
253 double inter_k ;
254
255 inter_k = do_interpolation( lambda, SplineAtmos, AcceleratorAtmos ) ;
256 kvalue = pow( 10., -inter_k*0.4 ) ;
257#else
258 /* Trivial solution, use the table value in all the intervall */
259 atmidx = get_idx_by_lambda( lambda, atmos_lambda, atmidx, atmos_size,
260 atmos_step ) ;
261 kvalue = pow(10., -atmos_K[atmidx]*0.4 ) ;
262#endif
263 }
264 else kvalue = 1. ;
265 /* Now, more complicated, correct with response: how ? interpolate (how) ?
266 other ? Could do the same as for atmos ext ?
267 */
268#if defined(USE_SPLINE)
269 vresp = do_interpolation( lambda, SplineResp, AcceleratorResp ) ;
270#else
271 /* Trivial solution, use the table value in all the intervall */
272 respidx = get_idx_by_lambda( lambda, response_lambda, respidx,
273 response_size, response_step ) ;
274 vresp = response_flux[respidx] ;
275#endif
276 /* Loop over slit */
277 for( islit = 0 ; islit < spectrum_nslit ; islit++ ) {
278 int idx=0 ;
279
280 if ( spectrum_qual[idx] != QFLAG_GOOD_PIXEL )
281 continue ;
282 idx = ipix + islit*spectrum_nlambda ;
283 spectrum_flux[idx] *= (airmass_ratio*kvalue*vresp)/Lx ;
284 }
285 }
286
288
289 cleanup:
290 return spectrum ;
291}
292
293cpl_frame * xsh_calibrate_flux( cpl_frame * spectrum_frame,
294 cpl_frame * respon_frame,
295 cpl_frame * atmos_ext_frame,
296 const char * fname,
298{
299 cpl_frame * result = NULL ;
300 double s_airm_start, s_airm_end, r_airm_start, r_airm_end ;
301 double s_airmass = 1., r_airmass = 1., airmass_ratio = 1. ;
302 double slit_width, seeing_start, seeing_end, seeing ;
303 cpl_propertylist * s_header = NULL, * r_header = NULL ;
304 xsh_spectrum * spectrum = NULL;
305 xsh_star_flux_list * response_list = NULL ;
306 xsh_atmos_ext_list * atmos_ext_list = NULL ;
307 double Lx ;
308
309 XSH_ASSURE_NOT_NULL( spectrum_frame ) ;
310 XSH_ASSURE_NOT_NULL( respon_frame ) ;
312
313 /* First get AIRMASS's */
314 check( spectrum = xsh_spectrum_load( spectrum_frame)) ;
315 check( s_header = spectrum->flux_header ) ;
316 s_airm_start = xsh_pfits_get_airm_start( s_header ) ;
317 s_airm_end = xsh_pfits_get_airm_end( s_header ) ;
318 if ( cpl_error_get_code() != CPL_ERROR_NONE ) {
319 cpl_error_reset() ;
320 s_airmass = 1. ;
321 }
322 else s_airmass = (s_airm_start + s_airm_end)/2. ;
323
324 check( response_list = xsh_star_flux_list_load( respon_frame ) ) ;
325 check( r_header = response_list->header ) ;
326 r_airm_start = xsh_pfits_get_airm_start( r_header ) ;
327 r_airm_end = xsh_pfits_get_airm_end( r_header ) ;
328 if ( cpl_error_get_code() != CPL_ERROR_NONE ) {
329 cpl_error_reset() ;
330 r_airmass = 1. ;
331 }
332 else r_airmass = (r_airm_start + r_airm_end)/2. ;
333 airmass_ratio = s_airmass/r_airmass ;
334
335 /* Get slit width and seeing */
337 slit_width = WIDTH_SLIT_IFU ;
338 xsh_msg( "IFU Mode Slit width: %.2lf", slit_width ) ;
339 }
340 else {
341 check( slit_width = xsh_pfits_get_slit_width( s_header, instrument ) ) ;
342 xsh_msg( "SLIT Mode Slit width: %.2lf", slit_width ) ;
343 }
344 check( seeing_start = xsh_pfits_get_seeing_start( s_header ) ) ;
345 check( seeing_end = xsh_pfits_get_seeing_end( s_header ) ) ;
346 seeing = (seeing_start + seeing_end)/2. ;
347 xsh_msg_warning("SEEING=%g",seeing);
348
349 Lx = compute_Lx( slit_width, seeing ) ;
350
351 /* Load ATMOS EXT (if not NULL ) */
352 if ( atmos_ext_frame != NULL ) {
353 check( atmos_ext_list = xsh_atmos_ext_list_load( atmos_ext_frame ) ) ;
354 xsh_msg( "ATMOS EXT Loaded" ) ;
355 }
356
357 /* Same function for 1D/2D spectra */
358 check( do_calib_spectrum( spectrum, response_list, atmos_ext_list,
359 airmass_ratio, Lx ) ) ;
360 /* How about IFU ? */
361
362 check( result = xsh_phys_spectrum_save( spectrum, fname, instrument )) ;
363
364 cleanup:
365
366 return result ;
367}
static const double step
static xsh_instrument * instrument
static float lambda_step
static gsl_interp_accel * AcceleratorAtmos
static void init_interpolate(double *x, double *yf, int nb, gsl_spline **spline, gsl_interp_accel **accel)
static double do_interpolation(double x, gsl_spline *spline, gsl_interp_accel *accel)
cpl_frame * xsh_calibrate_flux(cpl_frame *spectrum_frame, cpl_frame *respon_frame, cpl_frame *atmos_ext_frame, const char *fname, xsh_instrument *instrument)
static xsh_spectrum * do_calib_spectrum(xsh_spectrum *spectrum_in, xsh_star_flux_list *response_list, xsh_atmos_ext_list *atmos_ext_list, double airmass_ratio, double Lx)
static gsl_spline * SplineAtmos
static gsl_spline * SplineResp
static gsl_interp_accel * AcceleratorResp
static void clear_interpolate(void)
static double compute_Lx(double slit_width, double seeing)
static double myfunc(double x, void *params)
int xsh_spectrum_get_size_lambda(xsh_spectrum *s)
Get lambda axis size of spectrum.
double xsh_spectrum_get_lambda_min(xsh_spectrum *s)
Get minimum lambda of spectrum.
xsh_spectrum * xsh_spectrum_load(cpl_frame *s1d_frame)
Load a 1D spectrum structure.
double xsh_spectrum_get_lambda_step(xsh_spectrum *s)
Get bin in lambda of spectrum.
cpl_frame * xsh_phys_spectrum_save(xsh_spectrum *s, const char *filename, xsh_instrument *instr)
save a spectrum
int * xsh_spectrum_get_qual(xsh_spectrum *s)
Get qual of spectrum.
xsh_spectrum * xsh_spectrum_duplicate(xsh_spectrum *org)
double * xsh_spectrum_get_flux(xsh_spectrum *s)
Get flux of spectrum.
int xsh_spectrum_get_size_slit(xsh_spectrum *s)
Get slit axis ize of spectrum.
#define check(COMMAND)
Definition: xsh_error.h:71
#define XSH_ASSURE_NOT_NULL(pointer)
Definition: xsh_error.h:99
XSH_MODE xsh_instrument_get_mode(xsh_instrument *i)
Get a mode on instrument structure.
int * y
int * x
#define xsh_msg_warning(...)
Print an warning message.
Definition: xsh_msg.h:88
#define xsh_msg(...)
Print a message on info level.
Definition: xsh_msg.h:121
double xsh_pfits_get_airm_end(const cpl_propertylist *plist)
find out the TEL AIRM END value
Definition: xsh_pfits.c:527
double xsh_pfits_get_seeing_start(const cpl_propertylist *plist)
find out the TEL AMBI START value (Seeing)
Definition: xsh_pfits.c:544
double xsh_pfits_get_slit_width(const cpl_propertylist *plist, xsh_instrument *instrument)
find out the INS OPTIx NAME value (the width of the slit)
Definition: xsh_pfits.c:580
double xsh_pfits_get_seeing_end(const cpl_propertylist *plist)
find out the TEL AMBI END value (Seeing)
Definition: xsh_pfits.c:561
double xsh_pfits_get_airm_start(const cpl_propertylist *plist)
find out the TEL AIRM START value
Definition: xsh_pfits.c:496
cpl_propertylist * flux_header
cpl_propertylist * header
#define QFLAG_GOOD_PIXEL
double * xsh_atmos_ext_list_get_K(xsh_atmos_ext_list *list)
xsh_atmos_ext_list * xsh_atmos_ext_list_load(cpl_frame *ext_frame)
double * xsh_atmos_ext_list_get_lambda(xsh_atmos_ext_list *list)
#define WIDTH_SLIT_IFU
@ XSH_MODE_IFU
double * xsh_star_flux_list_get_lambda(xsh_star_flux_list *list)
xsh_star_flux_list * xsh_star_flux_list_load(cpl_frame *star_frame)
double * xsh_star_flux_list_get_flux(xsh_star_flux_list *list)
#define M_PI
Definition: xsh_utils.h:43