00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #ifdef HAVE_CONFIG_H
00027 # include <config.h>
00028 #endif
00029
00030
00036
00039
00040
00041
00042
00043 #include <tests.h>
00044
00045 #include <xsh_error.h>
00046 #include <xsh_msg.h>
00047 #include <xsh_efficiency_response.h>
00048 #include <xsh_data_instrument.h>
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060 #include <cpl.h>
00061
00062
00063 #include <stdio.h>
00064 #include <stdlib.h>
00065 #include <math.h>
00066 #include <gsl/gsl_bspline.h>
00067 #include <gsl/gsl_multifit.h>
00068 #include <gsl/gsl_rng.h>
00069 #include <gsl/gsl_randist.h>
00070 #include <gsl/gsl_statistics.h>
00071
00072
00073
00074
00075
00076 #define MODULE_ID "XSH_BSPLINE"
00077
00078
00079 #define N 200
00080
00081
00082 #define NCOEFFS 12 // UVB 21 VIS 12 NIR 12
00083
00084 #define ORDER 6
00085
00086
00087 #define NBREAK (NCOEFFS +2 -ORDER)
00088
00096 int main( int argc, char **argv)
00097 {
00098
00099 int ret =0;
00100
00101
00102 size_t n = N;
00103 const size_t ncoeffs = NCOEFFS;
00104 const size_t nbreak = NBREAK;
00105 size_t i, j;
00106 gsl_bspline_workspace *bw;
00107 gsl_vector *B;
00108 double dy;
00109 gsl_rng *r;
00110 gsl_vector *c, *w;
00111 gsl_vector *x, *y;
00112 gsl_matrix *X, *cov;
00113 gsl_multifit_linear_workspace *mw;
00114 double chisq, Rsq, dof, tss;
00115 float* wave=NULL;
00116 float* response=NULL;
00117 double* pwav=NULL;
00118 double* pfit=NULL;
00119 int kh=0;
00120 double dump_factor=1.e10;
00121 cpl_table* tab_response=NULL;
00122 cpl_table* tab_resp_fit=NULL;
00123 HIGH_ABS_REGION * phigh;
00124 xsh_instrument* inst=NULL;
00125 const char* dir_uvb="/data3/xsh/spr/data_xsh_respon_slit_nod_uvb.sof/";
00126 const char* dir_vis="/data3/xsh/spr/data_xsh_respon_slit_nod_vis.sof/";
00127 const char* dir_nir="/data3/xsh/spr/data_xsh_respon_slit_nod_nir.sof/";
00128
00129 char fname[256];
00130
00131
00132
00133
00134
00135 TESTS_INIT(MODULE_ID);
00136 inst = xsh_instrument_new();
00137 xsh_instrument_set_arm(inst, XSH_ARM_NIR);
00138
00139 if( xsh_instrument_get_arm(inst) == XSH_ARM_UVB){
00140 sprintf(fname,"%sresponse_smooth.fits",dir_uvb);
00141 }
00142 else if( xsh_instrument_get_arm(inst) == XSH_ARM_VIS){
00143 sprintf(fname,"%sresponse_smooth.fits",dir_vis);
00144 }
00145 else if( xsh_instrument_get_arm(inst) == XSH_ARM_NIR){
00146 sprintf(fname,"%sresponse_smooth.fits",dir_nir);
00147 }
00148
00149 phigh=xsh_fill_high_abs_regions(inst,NULL);
00150
00151 tab_response=cpl_table_load(fname,1,0);
00152
00153 n=cpl_table_get_nrow(tab_response);
00154 wave=cpl_table_get_data_float(tab_response,"LAMBDA");
00155 response=cpl_table_get_data_float(tab_response,"FLUX");
00156
00157
00158
00159
00160
00161 gsl_rng_env_setup();
00162 r = gsl_rng_alloc(gsl_rng_default);
00163
00164
00165 bw = gsl_bspline_alloc(ORDER, nbreak);
00166 B = gsl_vector_alloc(ncoeffs);
00167
00168 x = gsl_vector_alloc(n);
00169 y = gsl_vector_alloc(n);
00170 X = gsl_matrix_alloc(n, ncoeffs);
00171 c = gsl_vector_alloc(ncoeffs);
00172 w = gsl_vector_alloc(n);
00173 cov = gsl_matrix_alloc(ncoeffs, ncoeffs);
00174 mw = gsl_multifit_linear_alloc(n, ncoeffs);
00175
00176 printf("#m=0,S=0\n");
00177
00178 for (i = 0; i < n; ++i)
00179 {
00180 double sigma;
00181 double xi = (double)wave[i];
00182 double yi = (double)response[i];
00183
00184 sigma = 0.1 * yi;
00185 dy = gsl_ran_gaussian(r, sigma);
00186 yi += dy;
00187
00188 gsl_vector_set(x, i, xi);
00189 gsl_vector_set(y, i, yi);
00190 gsl_vector_set(w, i, 1.0 / (sigma * sigma));
00191
00192
00193 }
00194
00195
00196 if ( phigh != NULL ) {
00197
00198 xsh_msg("Flag High Absorption Regions" ) ;
00199 for( kh = 0 ; phigh->lambda_min != 0. ; phigh++ ) {
00200 for( i = 0 ; i < n ; i++ ) {
00201 if( wave[i] >= phigh->lambda_min &&
00202 wave[i] <= phigh->lambda_max ) {
00203 gsl_vector_set(w, i, gsl_vector_get(w,i)/dump_factor );
00204 }
00205 }
00206 }
00207 }
00208
00209
00210
00211 gsl_bspline_knots_uniform(wave[0], wave[n-1], bw);
00212
00213
00214 for (i = 0; i < n; ++i)
00215 {
00216 double xi = gsl_vector_get(x, i);
00217
00218
00219 gsl_bspline_eval(xi, B, bw);
00220
00221
00222 for (j = 0; j < ncoeffs; ++j)
00223 {
00224 double Bj = gsl_vector_get(B, j);
00225 gsl_matrix_set(X, i, j, Bj);
00226 }
00227 }
00228
00229
00230 gsl_multifit_wlinear(X, w, y, c, cov, &chisq, mw);
00231
00232 dof = n - ncoeffs;
00233 tss = gsl_stats_wtss(w->data, 1, y->data, 1, y->size);
00234 Rsq = 1.0 - chisq / tss;
00235
00236 fprintf(stderr, "chisq/dof = %e, Rsq = %f\n",
00237 chisq / dof, Rsq);
00238
00239 tab_resp_fit=cpl_table_new(n);
00240 cpl_table_new_column(tab_resp_fit, "wave", CPL_TYPE_DOUBLE);
00241 cpl_table_new_column(tab_resp_fit, "fit", CPL_TYPE_DOUBLE);
00242 cpl_table_fill_column_window_double(tab_resp_fit,"wave",0,n,0);
00243 cpl_table_fill_column_window_double(tab_resp_fit,"fit",0,n,0);
00244 pwav=cpl_table_get_data_double(tab_resp_fit,"wave");
00245 pfit=cpl_table_get_data_double(tab_resp_fit,"fit");
00246
00247
00248 {
00249 double xi, yi, yerr;
00250 printf("#m=1,S=0\n");
00251 for (i = 0; i < n; i ++)
00252 {
00253 xi=(double)wave[i];
00254 gsl_bspline_eval(xi, B, bw);
00255 gsl_multifit_linear_est(B, c, cov, &yi, &yerr);
00256
00257 pwav[i]=xi;
00258 pfit[i]=yi;
00259 }
00260 }
00261 cpl_table_save(tab_resp_fit,NULL,NULL,"resp_fit.fits",CPL_IO_DEFAULT);
00262
00263 xsh_msg("ok3");
00264
00265
00266
00267 gsl_rng_free(r);
00268 gsl_bspline_free(bw);
00269 gsl_vector_free(B);
00270 gsl_vector_free(x);
00271 gsl_vector_free(y);
00272 gsl_matrix_free(X);
00273 gsl_vector_free(c);
00274 gsl_vector_free(w);
00275 gsl_matrix_free(cov);
00276 gsl_multifit_linear_free(mw);
00277 xsh_free_table(&tab_response);
00278 xsh_free_table(&tab_resp_fit);
00279 xsh_instrument_free(&inst);
00280 if (cpl_error_get_code() != CPL_ERROR_NONE) {
00281 xsh_error_dump(CPL_MSG_ERROR);
00282 ret= 1;
00283 }
00284 TEST_END();
00285 return ret ;
00286 }
00287