61#include <gsl/gsl_integration.h>
66#if defined(USE_SPLINE)
67#include <gsl/gsl_spline.h>
79#if defined(USE_SPLINE)
85 gsl_interp_accel ** accel)
91 *accel = gsl_interp_accel_alloc();
92 if ( *accel == NULL )
xsh_msg(
"Accelerator = NULL" ) ;
94 *spline = gsl_spline_alloc( gsl_interp_cspline, nb);
96 if ( *spline == NULL )
xsh_msg(
"Spline = NULL" ) ;
98 ok = gsl_spline_init( *spline,
x, yf, nb ) ;
100 xsh_msg(
"gsl_spline_init --> %d", ok ) ;
106 gsl_interp_accel * accel )
110 y = gsl_spline_eval( spline,
x, accel );
125static int get_idx_by_lambda(
double lambda,
double * a_lambda,
int from,
126 int to,
double step )
129 double * plambda=NULL ;
131 for( idx = from, plambda = a_lambda+from ; idx<to ; idx++, plambda++ ) {
132 if ( lambda >= (*plambda -
step/2) &&
133 lambda < (*plambda +
step/2) ) {
143static double myfunc (
double x,
void * params) {
144 double alpha = *(
double *) params;
145 double thefunc = exp(alpha*
x*
x) ;
152 gsl_integration_workspace * w
153 = gsl_integration_workspace_alloc (10000);
155 double result, error;
165 limit = (slit_width/2.)/((seeing*sqrt(2.))/2.35482) ;
167 gsl_integration_qags (&F, 0, limit, 0, 1e-7, 10000,
171 Lx = 1- (2/sqrt(
M_PI))*result ;
179 double airmass_ratio,
183 double * spectrum_flux = NULL ;
184 int * spectrum_qual = NULL ;
185 double spectrum_lambda_min=0;
187 int spectrum_nlambda=0, spectrum_nslit=0 ;
191 double * atmos_lambda = NULL ;
192 double * atmos_K = NULL ;
196 double * response_lambda = NULL ;
197 double * response_flux = NULL ;
199 int response_size=0 ;
201#if !defined(USE_SPLINE)
202 int atmidx = 0, respidx = 0 ;
209 if ( spectrum_nslit > 1 ) {
224 if (atmos_ext_list != NULL) {
225 atmos_size = atmos_ext_list->
size ;
228#if defined(USE_SPLINE)
235 response_size = response_list->
size ;
237#if defined(USE_SPLINE)
243 lambda = spectrum_lambda_min ;
245 for( ipix = 0 ; ipix < spectrum_nlambda ; ipix++, lambda +=
lambda_step ) {
246 double kvalue, vresp ;
251 if ( atmos_ext_list != NULL ) {
252#if defined(USE_SPLINE)
256 kvalue = pow( 10., -inter_k*0.4 ) ;
259 atmidx = get_idx_by_lambda( lambda, atmos_lambda, atmidx, atmos_size,
261 kvalue = pow(10., -atmos_K[atmidx]*0.4 ) ;
268#if defined(USE_SPLINE)
272 respidx = get_idx_by_lambda( lambda, response_lambda, respidx,
273 response_size, response_step ) ;
274 vresp = response_flux[respidx] ;
277 for( islit = 0 ; islit < spectrum_nslit ; islit++ ) {
282 idx = ipix + islit*spectrum_nlambda ;
283 spectrum_flux[idx] *= (airmass_ratio*kvalue*vresp)/Lx ;
294 cpl_frame * respon_frame,
295 cpl_frame * atmos_ext_frame,
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 ;
318 if ( cpl_error_get_code() != CPL_ERROR_NONE ) {
322 else s_airmass = (s_airm_start + s_airm_end)/2. ;
328 if ( cpl_error_get_code() != CPL_ERROR_NONE ) {
332 else r_airmass = (r_airm_start + r_airm_end)/2. ;
333 airmass_ratio = s_airmass/r_airmass ;
338 xsh_msg(
"IFU Mode Slit width: %.2lf", slit_width ) ;
342 xsh_msg(
"SLIT Mode Slit width: %.2lf", slit_width ) ;
346 seeing = (seeing_start + seeing_end)/2. ;
352 if ( atmos_ext_frame != NULL ) {
354 xsh_msg(
"ATMOS EXT Loaded" ) ;
359 airmass_ratio, Lx ) ) ;
static xsh_instrument * instrument
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 XSH_ASSURE_NOT_NULL(pointer)
XSH_MODE xsh_instrument_get_mode(xsh_instrument *i)
Get a mode on instrument structure.
#define xsh_msg_warning(...)
Print an warning message.
#define xsh_msg(...)
Print a message on info level.
double xsh_pfits_get_airm_end(const cpl_propertylist *plist)
find out the TEL AIRM END value
double xsh_pfits_get_seeing_start(const cpl_propertylist *plist)
find out the TEL AMBI START value (Seeing)
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)
double xsh_pfits_get_seeing_end(const cpl_propertylist *plist)
find out the TEL AMBI END value (Seeing)
double xsh_pfits_get_airm_start(const cpl_propertylist *plist)
find out the TEL AIRM START value
cpl_propertylist * flux_header
cpl_propertylist * header
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)
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)