X-shooter Pipeline Reference Manual 3.8.15
xsh_ngaussfdf.h
Go to the documentation of this file.
1/*
2 * File: xsh_ngaussdfd.h
3 * Author: sousasag
4 *
5 * Created on April 23, 2011, 12:28 PM
6 */
7/* expfit.c -- compute residual for exponential
8 +background model */
9
10struct data {
11 size_t n;
12 int para;
13 double * t;
14 double * y;
15 double * sigma;
16};
17
18int
19expb_f (const gsl_vector * x, void *data,
20 gsl_vector * f)
21{
22 size_t n = ((struct data *)data)->n;
23 double *t = ((struct data *)data)->t;
24 double *y = ((struct data *)data)->y;
25 double *sigma = ((struct data *) data)->sigma;
26 int para = ((struct data *)data)->para;
27
28 double ac[para];
29 int ijk;
30 for (ijk=0;ijk<para;ijk++)
31 ac[ijk]=gsl_vector_get (x, ijk);
32
33 size_t i;
34
35 for (i = 0; i < n; i++)
36 {
37 /* Model Yi = sum(i) : Ai * exp(-lambdai * (i-ci)*(i-ci))*/
38
39 double Yi=0.0;
40 int j;
41 for(j=0;j<para;j+=3)
42 {
43 Yi+= ac[j] * exp (-ac[j+1] * (t[i]-ac[j+2])*(t[i]-ac[j+2]) ) ;
44 }
45 gsl_vector_set (f, i, (Yi - y[i])/sigma[i]);
46
47 }
48
49 return GSL_SUCCESS;
50}
51
52int
53expb_df (const gsl_vector * x, void *data,
54 gsl_matrix * J)
55{
56 size_t n = ((struct data *)data)->n;
57 double *sigma = ((struct data *) data)->sigma;
58 double *t = ((struct data *)data)->t;
59 int para = ((struct data *)data)->para;
60
61 double ac[para];
62 int ijk;
63 for (ijk=0;ijk<para;ijk++)
64 ac[ijk]=gsl_vector_get (x, ijk);
65
66 size_t i;
67
68 for (i = 0; i < n; i++)
69 {
70 /* Jacobian matrix J(i,j) = dfi / dxj, */
71 /* where fi = (Yi - yi)/sigma[i], */
72 /* and the xj are the parameters (Ai,lambdai,ci,...) */
73
74 int j=0;
75 double ti = t[i];
76 double s = sigma[i];
77
78 for (j=0; j<para; j+=3)
79 {
80 double e1 = exp(-ac[j+1]*(ti-ac[j+2])*(ti-ac[j+2]) ) ;
81
82 gsl_matrix_set (J, i, j, e1/s);
83 gsl_matrix_set (J, i, j+1, - ac[j] * (ti-ac[j+2])*(ti-ac[j+2]) * e1/s);
84 gsl_matrix_set (J, i, j+2, ac[j] * ac[j+1] * 2. * (ti-ac[j+2]) * e1/s);
85 }
86
87 }
88 return GSL_SUCCESS;
89}
90
91int
92expb_fdf (const gsl_vector * x, void *data,
93 gsl_vector * f, gsl_matrix * J)
94{
95 expb_f (x, data, f);
96 expb_df (x, data, J);
97
98 return GSL_SUCCESS;
99}
static double sigma
int * y
int * x
static SimAnneal s
Definition: xsh_model_sa.c:99
double * y
Definition: xsh_ngaussfdf.h:14
size_t n
Definition: xsh_ngaussfdf.h:11
double * t
Definition: xsh_ngaussfdf.h:13
double * sigma
Definition: xsh_ngaussfdf.h:15
int para
Definition: xsh_ngaussfdf.h:12
int n
Definition: xsh_detmon_lg.c:92
int expb_fdf(const gsl_vector *x, void *data, gsl_vector *f, gsl_matrix *J)
Definition: xsh_ngaussfdf.h:92
int expb_f(const gsl_vector *x, void *data, gsl_vector *f)
Definition: xsh_ngaussfdf.h:19
int expb_df(const gsl_vector *x, void *data, gsl_matrix *J)
Definition: xsh_ngaussfdf.h:53