X-shooter Pipeline Reference Manual 3.8.15
test-xsh_response_fit_bspline.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: 2013-05-14 09:55:30 $
23 * $Revision: 1.3 $
24 */
25
26#ifdef HAVE_CONFIG_H
27# include <config.h>
28#endif
29
30/*-------------------------------------------------------------------------*/
36/*-------------------------------------------------------------------------*/
39/*--------------------------------------------------------------------------
40 Includes
41 --------------------------------------------------------------------------*/
42
43#include <tests.h>
44
45#include <xsh_error.h>
46#include <xsh_msg.h>
48#include <xsh_data_instrument.h>
49
50/*
51#include <xsh_data_pre.h>
52#include <xsh_data_rec.h>
53#include <xsh_data_localization.h>
54#include <xsh_drl.h>
55#include <xsh_pfits.h>
56#include <xsh_utils.h>
57
58#include <xsh_badpixelmap.h>
59#include <getopt.h>
60*/
61#include <cpl.h>
62
63
64#include <stdio.h>
65#include <stdlib.h>
66#include <math.h>
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>
72
73/*--------------------------------------------------------------------------
74 Defines
75 --------------------------------------------------------------------------*/
76
77#define MODULE_ID "XSH_RESPONSE_FIT_BSPLINE"
78
79/* number of data points to fit */
80#define N 200
81
82/* number of fit coefficients */
83#define NCOEFFS 12 // UVB 21 VIS 12 NIR 12
84
85#define ORDER 6
86
87/* nbreak = ncoeffs +2 -ORDER = ncoeffs -2 since ORDER = 4 */
88#define NBREAK (NCOEFFS +2 -ORDER)
89
97int main( int argc, char **argv)
98{
99 /* Declarations */
100 int ret =0;
101
102
103 size_t n = N;
104 const size_t ncoeffs = NCOEFFS;
105 const size_t nbreak = NBREAK;
106 size_t i, j;
107 gsl_bspline_workspace *bw;
108 gsl_vector *B;
109 double dy;
110 gsl_rng *r;
111 gsl_vector *c, *w;
112 gsl_vector *x, *y;
113 gsl_matrix *X, *cov;
114 gsl_multifit_linear_workspace *mw;
115 double chisq, Rsq, dof, tss;
116 float* wave=NULL;
117 float* response=NULL;
118 double* pwav=NULL;
119 double* pfit=NULL;
120 int kh=0;
121 double dump_factor=1.e10;
122 cpl_table* tab_response=NULL;
123 cpl_table* tab_resp_fit=NULL;
124 HIGH_ABS_REGION * phigh;
125 xsh_instrument* inst=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/";
129
130 char fname[256];
131 const int activate=0;
134
135 cpl_ensure(activate > 0, CPL_ERROR_ILLEGAL_INPUT,NULL);
136 XSH_ASSURE_NOT_NULL_MSG(fname,"Null input response");
137 /* Initialize libraries */
138
139 inst = xsh_instrument_new();
141
143 sprintf(fname,"%sresponse_smooth.fits",dir_uvb);
144 }
145 else if( xsh_instrument_get_arm(inst) == XSH_ARM_VIS){
146 sprintf(fname,"%sresponse_smooth.fits",dir_vis);
147 }
148 else if( xsh_instrument_get_arm(inst) == XSH_ARM_NIR){
149 sprintf(fname,"%sresponse_smooth.fits",dir_nir);
150 }
151
153 phigh=xsh_fill_high_abs_regions(inst,NULL);
154
155 tab_response=cpl_table_load(fname,1,0);
156
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");
160
161 /*
162 for (i = 0; i < n; ++i) {
163 xsh_msg("wave=%g resp=%g",wave[i],response[i]);
164 }
165 */
166 gsl_rng_env_setup();
167 r = gsl_rng_alloc(gsl_rng_default);
168
169 /* allocate a cubic bspline workspace (ORDER = 4) */
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);
179
180 printf("#m=0,S=0\n");
181 /* this is the data to be fitted */
182 for (i = 0; i < n; ++i)
183 {
184 double sigma;
185 double xi = (double)wave[i];
186 double yi = (double)response[i];
187
188 sigma = 0.1 * yi;
189 dy = gsl_ran_gaussian(r, sigma);
190 yi += dy;
191
192 gsl_vector_set(x, i, xi);
193 gsl_vector_set(y, i, yi);
194 gsl_vector_set(w, i, 1.0 / (sigma * sigma));
195
196 //printf("%f %g\n", xi, yi);
197 }
198
199 if ( phigh != NULL ) {
200
201 xsh_msg("Flag High Absorption Regions" ) ;
202 for( kh = 0 ; phigh->lambda_min != 0. ; phigh++ ) {
203 for( i = 0 ; i < n ; i++ ) {
204 if( wave[i] >= phigh->lambda_min &&
205 wave[i] <= phigh->lambda_max ) {
206 gsl_vector_set(w, i, gsl_vector_get(w,i)/dump_factor );
207 }
208 }
209 }
210 }
211
212 /* use uniform breakpoints on [0, 15] */
213 gsl_bspline_knots_uniform(wave[0], wave[n-1], bw);
214
215 /* construct the fit matrix X */
216 for (i = 0; i < n; ++i)
217 {
218 double xi = gsl_vector_get(x, i);
219
220 /* compute B_j(xi) for all j */
221 gsl_bspline_eval(xi, B, bw);
222
223 /* fill in row i of X */
224 for (j = 0; j < ncoeffs; ++j)
225 {
226 double Bj = gsl_vector_get(B, j);
227 gsl_matrix_set(X, i, j, Bj);
228 }
229 }
230
231 /* do the fit */
232 gsl_multifit_wlinear(X, w, y, c, cov, &chisq, mw);
233
234 dof = n - ncoeffs;
235 tss = gsl_stats_wtss(w->data, 1, y->data, 1, y->size);
236 Rsq = 1.0 - chisq / tss;
237
238 fprintf(stderr, "chisq/dof = %e, Rsq = %f\n",
239 chisq / dof, Rsq);
240
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");
248
249 /* output the smoothed curve */
250 {
251 double xi, yi, yerr;
252 printf("#m=1,S=0\n");
253 for (i = 0; i < n; i ++)
254 {
255 xi=(double)wave[i];
256 gsl_bspline_eval(xi, B, bw);
257 gsl_multifit_linear_est(B, c, cov, &yi, &yerr);
258 //printf("%f %g\n", xi, yi);
259 pwav[i]=xi;
260 pfit[i]=yi;
261 }
262 }
263 cpl_table_save(tab_resp_fit,NULL,NULL,"resp_fit.fits",CPL_IO_DEFAULT);
264
265 gsl_rng_free(r);
266 gsl_bspline_free(bw);
267 gsl_vector_free(B);
268 gsl_vector_free(x);
269 gsl_vector_free(y);
270 gsl_matrix_free(X);
271 gsl_vector_free(c);
272 gsl_vector_free(w);
273 gsl_matrix_free(cov);
274 gsl_multifit_linear_free(mw);
275 xsh_free_table(&tab_response);
276 xsh_free_table(&tab_resp_fit);
277 xsh_instrument_free(&inst);
278
279 cleanup:
280
281 if (cpl_error_get_code() != CPL_ERROR_NONE) {
282 xsh_error_dump(CPL_MSG_ERROR);
283 ret= 1;
284 }
286 TEST_END();
287 return ret ;
288}
289
int main()
Unit test of xsh_bspline_interpol.
static double sigma
#define MODULE_ID
#define XSH_ASSURE_NOT_NULL_MSG(pointer, msg)
Definition: xsh_error.h:103
#define xsh_error_dump(level)
Definition: xsh_error.h:92
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
int * y
int * x
#define xsh_msg(...)
Print a message on info level.
Definition: xsh_msg.h:121
void xsh_free_table(cpl_table **t)
Deallocate a table and set the pointer to NULL.
Definition: xsh_utils.c:2133
XSH_INSTRCONFIG * config
#define TESTS_CLEAN_WORKSPACE(DRL_ID)
Definition: tests.h:139
#define TESTS_INIT_WORKSPACE(DRL_ID)
Definition: tests.h:133
#define TEST_END()
Definition: tests.h:111
#define TESTS_INIT(DRL_ID)
Definition: tests.h:105
@ XSH_ARM_UVB
@ XSH_ARM_NIR
@ XSH_ARM_VIS
int n
Definition: xsh_detmon_lg.c:92
HIGH_ABS_REGION * xsh_fill_high_abs_regions(xsh_instrument *instrument, cpl_frame *high_abs_frame)
#define XSH_MALLOC(POINTER, TYPE, SIZE)
Definition: xsh_utils.h:49