X-shooter Pipeline Reference Manual 3.8.15
test-xsh_gaussian_fit.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-04-27 17:35:03 $
23 * $Revision: 1.5 $
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#include <stdlib.h>
45#include <stdio.h>
46
47
48#include <xsh_data_pre.h>
49#include <xsh_error.h>
50#include <xsh_msg.h>
51#include <xsh_data_instrument.h>
52#include <xsh_data_spectrum.h>
54#include <xsh_drl.h>
55#include <xsh_pfits.h>
56#include <xsh_parameters.h>
57#include <xsh_badpixelmap.h>
58#include <xsh_utils_ifu.h>
59#include <xsh_utils.h>
60#include <cpl.h>
61#include <math.h>
62
63#include <string.h>
64#include <getopt.h>
65
66/*--------------------------------------------------------------------------
67 Defines
68 --------------------------------------------------------------------------*/
69
70#define MODULE_ID "XSH_GAUSSIAN_FIT"
71
72enum {
74} ;
75
76static struct option long_options[] = {
77 {"help", 0, 0, HELP_OPT},
78 {0, 0, 0, 0}
79};
80
81static void Help( void )
82{
83 puts( "Unitary test of xsh_gaussian_fit");
84 puts( "Usage: test_xsh_gaussian_fit [options] DATA_FILE");
85
86 puts( "Options" ) ;
87 puts( " --help : What you see" ) ;
88 puts( "\nInput Files" ) ;
89 puts( "DATA_FILE :ASCII x,y file y column is use for decomposition");
90 TEST_END();
91}
92
93
94static void HandleOptions( int argc, char **argv)
95{
96 int opt ;
97 int option_index = 0;
98
99 while (( opt = getopt_long (argc, argv,
100 "help",
101 long_options, &option_index)) != EOF ){
102
103 switch ( opt ) {
104 default:
105 Help(); exit(-1);
106 }
107 }
108 return;
109}
110
111
112int main( int argc, char **argv)
113{
114 /* Declarations */
115 int ret = 0 ;
116 const char *file_name = NULL;
117 FILE* file = NULL;
118 char line[200];
119 double xpos[10000];
120 double ypos[10000];
121 int i, size=0;
122 int max = 10000;
123 char res_name[256];
124
125 /* gsl */
126 double cpl_area=0, cpl_sigma=0, cpl_x0=0, cpl_offset=0;
127 int status;
128
129 cpl_vector *pos_vect = NULL;
130 cpl_vector *data_vect = NULL;
131
132
133 double init_par[6];
134 double errs_par[6];
135
136 /* Initialize libraries */
138
139 cpl_msg_set_level(CPL_MSG_DEBUG);
141
142 /* Analyse parameters */
143 HandleOptions( argc, argv);
144
145 if ( (argc - optind) > 0 ) {
146 file_name = argv[optind];
147 }
148 else {
149 Help();
150 exit( 0);
151 }
152
153 xsh_msg("---Input Files");
154 xsh_msg("File : %s ", file_name);
155
156 /* READ data */
157 file = fopen( file_name, "r");
158 if (file != NULL){
159 size=0;
160 while ( size < max && fgets( line, 200, file)){
161 char col1[20];
162 char col2[20];
163
164 if ( line[0] != '#'){
165 /* Note: update the string format (the %__s) if col1(2) changes size.
166 * Remember: field width must equal 1 less than sizeof(col1(2)). */
167 sscanf(line,"%19s %19s", col1, col2);
168 xpos[size] = atof(col1);
169 ypos[size] = atof(col2);
170 size++;
171 }
172 }
173 }
174 if(file != NULL) {
175 fclose(file);
176 }
177
178 check( pos_vect = cpl_vector_wrap( size, xpos));
179 check( data_vect = cpl_vector_wrap( size, ypos));
180
181 cpl_vector_fit_gaussian( pos_vect, NULL, data_vect,
182 NULL,CPL_FIT_ALL,&cpl_x0,&cpl_sigma,&cpl_area,&cpl_offset,NULL,NULL,NULL);
183
184 printf( "CPL FIT area %f x0 %f sigma %f offset %f\n", cpl_area, cpl_x0, cpl_sigma, cpl_offset);
185
186 xsh_gsl_init_gaussian_fit( pos_vect, data_vect, init_par);
187 check( xsh_gsl_fit_gaussian( pos_vect, data_vect, 0, init_par, errs_par, &status));
188
189
190 printf ("area = %.5f +/- %.5f\n", init_par[0], errs_par[0]);
191 printf ("a = %.5f +/- %.5f\n", init_par[1], errs_par[1]);
192 printf ("b = %.5f +/- %.5f\n", init_par[2], errs_par[2]);
193 printf ("c = %.5f +/- %.5f\n", init_par[3], errs_par[3]);
194 printf ("x0 = %.5f +/- %.5f\n", init_par[4], errs_par[4]);
195 printf ("sigma = %.5f +/- %.5f\n", init_par[5], errs_par[5]);
196
197 sprintf( res_name, "GSL_FIT_%s.dat", file_name);
198
199 file = fopen( res_name, "w+");
200 fprintf( file, "#x y\n");
201
202 for( i=0; i< size; i++){
203 double step = xpos[1]-xpos[0];
204 double j;
205
206 for (j=0; j< step; j+=step/10.){
207 double t = xpos[i]+j;
208
209 double area = init_par[0];
210 double a = init_par[1];
211 double b = init_par[2];
212 double c = init_par[3];
213 double x0 = init_par[4];
214 double sigma =init_par[5];
215 double height = area/sqrt(2*M_PI*sigma*sigma);
216
217 double W = ((t-x0)*(t-x0))/(2*sigma*sigma);
218
219 double Yi = height*exp(-W)+a+b*t+c*t*t;
220
221
222 fprintf( file, "%f %f\n", t, Yi);
223 }
224 }
225 fclose(file);
226
227 sprintf( res_name, "CPL_FIT_%s.dat", file_name);
228
229 file = fopen( res_name, "w+");
230 fprintf( file, "#x y\n");
231
232 for( i=0; i< size; i++){
233 double step = xpos[1]-xpos[0];
234 double j;
235
236 for (j=0; j< step; j+=step/10.){
237 double t = xpos[i]+j;
238
239 double area = cpl_area;
240 double a = cpl_offset;
241 double b = 0;
242 double c = 0;
243 double x0 = cpl_x0;
244 double sigma = cpl_sigma;
245 double height = area/sqrt(2*M_PI*sigma*sigma);
246
247 double W = ((t-x0)*(t-x0))/(2*sigma*sigma);
248
249 double Yi = height*exp(-W)+a+b*t+c*t*t;
250
251
252 fprintf( file, "%f %f\n", t, Yi);
253 }
254 }
255 fclose(file);
256
257
258
259
260 cleanup:
261 if (cpl_error_get_code() != CPL_ERROR_NONE) {
262 xsh_error_dump(CPL_MSG_ERROR);
263 ret=1;
264 }
265 xsh_unwrap_vector( &pos_vect);
266 xsh_unwrap_vector( &data_vect);
267 TEST_END();
268 return ret ;
269}
270
int main()
Unit test of xsh_bspline_interpol.
static const double step
static void HandleOptions(int argc, char **argv)
#define MODULE_ID
static struct option long_options[]
static double sigma
#define check(COMMAND)
Definition: xsh_error.h:71
#define xsh_error_dump(level)
Definition: xsh_error.h:92
int size
#define xsh_msg(...)
Print a message on info level.
Definition: xsh_msg.h:121
void xsh_unwrap_vector(cpl_vector **v)
Unwrap a vector and set the pointer to NULL.
Definition: xsh_utils.c:2345
int xsh_debug_level_set(int level)
set debug level
Definition: xsh_utils.c:3125
void xsh_gsl_init_gaussian_fit(cpl_vector *xpos_vect, cpl_vector *ypos_vect, double *init_par)
Definition: xsh_utils.c:6862
void xsh_gsl_fit_gaussian(cpl_vector *xpos_vect, cpl_vector *ypos_vect, int deg, double *params, double *errs, int *status)
Definition: xsh_utils.c:6936
static void Help(void)
#define TEST_END()
Definition: tests.h:111
#define TESTS_INIT(DRL_ID)
Definition: tests.h:105
#define max(a, b)
@ XSH_DEBUG_LEVEL_MEDIUM
Definition: xsh_utils.h:138
#define M_PI
Definition: xsh_utils.h:43