67#include <gsl/gsl_bspline.h>
68#include <gsl/gsl_multifit.h>
69#include <gsl/gsl_rng.h>
70#include <gsl/gsl_randist.h>
71#include <gsl/gsl_statistics.h>
77#define MODULE_ID "XSH_RESPONSE_FIT_BSPLINE"
88#define NBREAK (NCOEFFS +2 -ORDER)
97int main(
int argc,
char **argv)
104 const size_t ncoeffs =
NCOEFFS;
105 const size_t nbreak =
NBREAK;
107 gsl_bspline_workspace *bw;
114 gsl_multifit_linear_workspace *mw;
115 double chisq, Rsq, dof, tss;
117 float* response=NULL;
121 double dump_factor=1.e10;
122 cpl_table* tab_response=NULL;
123 cpl_table* tab_resp_fit=NULL;
126 const char* dir_uvb=
"/data3/xsh/spr/data_xsh_respon_slit_nod_uvb.sof/";
127 const char* dir_vis=
"/data3/xsh/spr/data_xsh_respon_slit_nod_vis.sof/";
128 const char* dir_nir=
"/data3/xsh/spr/data_xsh_respon_slit_nod_nir.sof/";
131 const int activate=0;
135 cpl_ensure(activate > 0, CPL_ERROR_ILLEGAL_INPUT,NULL);
143 sprintf(fname,
"%sresponse_smooth.fits",dir_uvb);
146 sprintf(fname,
"%sresponse_smooth.fits",dir_vis);
149 sprintf(fname,
"%sresponse_smooth.fits",dir_nir);
155 tab_response=cpl_table_load(fname,1,0);
157 n=cpl_table_get_nrow(tab_response);
158 wave=cpl_table_get_data_float(tab_response,
"LAMBDA");
159 response=cpl_table_get_data_float(tab_response,
"FLUX");
167 r = gsl_rng_alloc(gsl_rng_default);
170 bw = gsl_bspline_alloc(
ORDER, nbreak);
171 B = gsl_vector_alloc(ncoeffs);
172 x = gsl_vector_alloc(
n);
173 y = gsl_vector_alloc(
n);
174 X = gsl_matrix_alloc(
n, ncoeffs);
175 c = gsl_vector_alloc(ncoeffs);
176 w = gsl_vector_alloc(
n);
177 cov = gsl_matrix_alloc(ncoeffs, ncoeffs);
178 mw = gsl_multifit_linear_alloc(
n, ncoeffs);
180 printf(
"#m=0,S=0\n");
182 for (i = 0; i <
n; ++i)
185 double xi = (double)wave[i];
186 double yi = (double)response[i];
189 dy = gsl_ran_gaussian(r,
sigma);
192 gsl_vector_set(
x, i, xi);
193 gsl_vector_set(
y, i, yi);
199 if ( phigh != NULL ) {
201 xsh_msg(
"Flag High Absorption Regions" ) ;
202 for( kh = 0 ; phigh->
lambda_min != 0. ; phigh++ ) {
203 for( i = 0 ; i <
n ; i++ ) {
206 gsl_vector_set(w, i, gsl_vector_get(w,i)/dump_factor );
213 gsl_bspline_knots_uniform(wave[0], wave[
n-1], bw);
216 for (i = 0; i <
n; ++i)
218 double xi = gsl_vector_get(
x, i);
221 gsl_bspline_eval(xi, B, bw);
224 for (j = 0; j < ncoeffs; ++j)
226 double Bj = gsl_vector_get(B, j);
227 gsl_matrix_set(X, i, j, Bj);
232 gsl_multifit_wlinear(X, w,
y, c, cov, &chisq, mw);
235 tss = gsl_stats_wtss(w->data, 1,
y->data, 1,
y->size);
236 Rsq = 1.0 - chisq / tss;
238 fprintf(stderr,
"chisq/dof = %e, Rsq = %f\n",
241 tab_resp_fit=cpl_table_new(
n);
242 cpl_table_new_column(tab_resp_fit,
"wave", CPL_TYPE_DOUBLE);
243 cpl_table_new_column(tab_resp_fit,
"fit", CPL_TYPE_DOUBLE);
244 cpl_table_fill_column_window_double(tab_resp_fit,
"wave",0,
n,0);
245 cpl_table_fill_column_window_double(tab_resp_fit,
"fit",0,
n,0);
246 pwav=cpl_table_get_data_double(tab_resp_fit,
"wave");
247 pfit=cpl_table_get_data_double(tab_resp_fit,
"fit");
252 printf(
"#m=1,S=0\n");
253 for (i = 0; i <
n; i ++)
256 gsl_bspline_eval(xi, B, bw);
257 gsl_multifit_linear_est(B, c, cov, &yi, &yerr);
263 cpl_table_save(tab_resp_fit,NULL,NULL,
"resp_fit.fits",CPL_IO_DEFAULT);
266 gsl_bspline_free(bw);
273 gsl_matrix_free(cov);
274 gsl_multifit_linear_free(mw);
281 if (cpl_error_get_code() != CPL_ERROR_NONE) {
int main()
Unit test of xsh_bspline_interpol.
#define XSH_ASSURE_NOT_NULL_MSG(pointer, msg)
#define xsh_error_dump(level)
void xsh_instrument_set_arm(xsh_instrument *i, XSH_ARM arm)
Set an arm on instrument structure.
XSH_ARM xsh_instrument_get_arm(xsh_instrument *i)
Get an arm on instrument structure.
void xsh_instrument_free(xsh_instrument **instrument)
free an instrument structure
xsh_instrument * xsh_instrument_new(void)
create new instrument structure
#define xsh_msg(...)
Print a message on info level.
void xsh_free_table(cpl_table **t)
Deallocate a table and set the pointer to NULL.
#define TESTS_CLEAN_WORKSPACE(DRL_ID)
#define TESTS_INIT_WORKSPACE(DRL_ID)
#define TESTS_INIT(DRL_ID)
HIGH_ABS_REGION * xsh_fill_high_abs_regions(xsh_instrument *instrument, cpl_frame *high_abs_frame)
#define XSH_MALLOC(POINTER, TYPE, SIZE)