]> git.treefish.org Git - phys/u1casc.git/blob - u1casc-ordinary/obs_diagcorr.hpp
seems to be working.
[phys/u1casc.git] / u1casc-ordinary / obs_diagcorr.hpp
1 #ifndef OBS_DIAGCORR_HPP
2 #define OBS_DIAGCORR_HPP
3
4 #include "latlib/o815/o815.h"
5
6 #include "latlib/writeout.h"
7
8 #include "latlib/obstat.hpp"
9
10 #include <iostream>
11 #include <complex>
12
13 #include <cmath>
14
15 #include <gsl/gsl_complex.h>
16 #include <gsl/gsl_complex_math.h>
17 #include <gsl/gsl_math.h>
18 #include <gsl/gsl_eigen.h>
19 #include <gsl/gsl_linalg.h>
20 #include <gsl/gsl_cblas.h>
21 #include <gsl/gsl_blas.h>
22
23 using namespace std;
24
25 class obs_diagcorr : public o815::obs {
26
27 public:
28   obs_diagcorr(o815 *_O815);
29   
30 private:
31   void _start();
32   void _meas(bool loadedobs, const int& nthmeas);
33   void _finish();
34
35   void corrCompute();
36   static complex<double> effMass(vector < complex<double> > *preCalculated, vector< complex<double> > *excludedMeas, int nmeas, void *para);
37   static void preEffMass(vector< vector < complex<double> > > *allVals, vector < complex<double> > *preCalculated, void *para);
38
39   sim *Sim;
40
41   int spatialV;
42
43   complex<double> *OM[16+1];
44
45   gsl_matrix_complex ***measurements;
46
47   gsl_vector_complex **measurements_disco;
48   
49   static void cdiag (gsl_vector *v, gsl_matrix_complex *evec, gsl_matrix_complex *m);
50 };
51
52 obs_diagcorr::obs_diagcorr(o815 *_O815) : o815::obs("diagcorr", 
53                                                     _O815->paraQ->getParaNames() + 
54                                                     "tsep"
55                                                     "...", 
56                                                     _O815, sizeof(complex<double>) * ( 16*_O815->comargs.lsize[1]/2 + 1*4 ) ) {
57   for (int ivar = 0; ivar<16+1; ivar++)
58     OM[ivar]  = (complex<double>*)( obsMem + ivar*sizeof(complex<double>)*( _O815->comargs.lsize[1]/2 ) );
59   
60   Sim = (sim*)O815->Sim;
61   spatialV = O815->comargs.lsize[0] * O815->comargs.lsize[0] * O815->comargs.lsize[0];
62
63   measurements = new gsl_matrix_complex**[O815->comargs.nmeas];
64   for (int imeas = 0; imeas<O815->comargs.nmeas; imeas++) {
65     measurements[imeas] = new gsl_matrix_complex*[_O815->comargs.lsize[1]/2];
66     for (int itsep=0; itsep<_O815->comargs.lsize[1]/2; itsep++)
67       measurements[imeas][itsep] = gsl_matrix_complex_alloc(4,4);
68   }
69
70   measurements_disco = new gsl_vector_complex*[O815->comargs.nmeas];
71   for (int imeas = 0; imeas<O815->comargs.nmeas; imeas++)
72     measurements_disco[imeas] = gsl_vector_complex_alloc(4);
73 }
74
75 void obs_diagcorr::_start() {
76   //  *out << O815->comargs.nmeas << endl;
77   
78   //*out << "OBS_test: start" << endl;
79 };
80
81 void obs_diagcorr::_meas(bool loadedobs, const int& nthmeas) {
82   if (!loadedobs)
83     corrCompute();
84
85   for (int icorr=0; icorr<4; icorr++)
86     for (int jcorr=0; jcorr<4; jcorr++)
87       for (int itsep=0; itsep<O815->comargs.lsize[1]/2; itsep++) {
88         gsl_matrix_complex_set(measurements[nthmeas][itsep], icorr, jcorr,
89                                gsl_complex_rect( OM[icorr*4+jcorr][itsep].real(), OM[icorr*4+jcorr][itsep].imag() ));
90       }
91   
92   for (int icorr=0; icorr<4; icorr++) {
93     gsl_vector_complex_set(measurements_disco[nthmeas], icorr,
94                            gsl_complex_rect( OM[16][icorr].real(), OM[16][icorr].imag() ) );
95   }
96 };
97
98 void obs_diagcorr::cdiag (gsl_vector *v, gsl_matrix_complex *evec, gsl_matrix_complex *m)
99 {
100   gsl_eigen_hermv_workspace *wspace = gsl_eigen_hermv_alloc(100);
101
102   gsl_eigen_hermv(m, v, evec, wspace);
103   gsl_eigen_hermv_sort(v, evec, GSL_EIGEN_SORT_VAL_ASC);
104   
105   gsl_eigen_hermv_free(wspace);
106 }
107
108 void obs_diagcorr::_finish() {
109   gsl_matrix_complex *totalval = gsl_matrix_complex_alloc(4,4);
110   gsl_matrix_complex *tmpmatrix = gsl_matrix_complex_alloc(4,4);
111   gsl_matrix_complex *tmpmatrix2 = gsl_matrix_complex_alloc(4,4);
112   gsl_vector *tmpvec = gsl_vector_alloc(4);
113   gsl_vector_complex *tmpvecc = gsl_vector_complex_alloc(4);
114   gsl_vector *tmpvec2 = gsl_vector_alloc(4);
115   gsl_vector *jackres = gsl_vector_alloc(4);
116   gsl_vector *jackerrornorm = gsl_vector_alloc(4);
117   gsl_vector_complex *totaldisco = gsl_vector_complex_alloc(4);
118   gsl_matrix_complex *evecres = gsl_matrix_complex_alloc(4,4);
119   gsl_matrix_complex *evecerrornorm = gsl_matrix_complex_alloc(4,4);
120   
121   gsl_vector_complex_set_zero(totaldisco);
122   for (int imeas=0; imeas<O815->comargs.nmeas; imeas++)
123     gsl_vector_complex_add( totaldisco, measurements_disco[imeas] );
124   
125   for (int itsep = 0; itsep < O815->comargs.lsize[1]/2; itsep++) {
126     gsl_matrix_complex_set_zero(totalval);
127     gsl_vector_set_zero(jackerrornorm);
128     gsl_matrix_complex_set_zero(evecerrornorm);
129
130     for (int imeas=0; imeas<O815->comargs.nmeas; imeas++)
131       gsl_matrix_complex_add( totalval, measurements[imeas][itsep] );
132     
133     gsl_matrix_complex_memcpy (tmpmatrix, totalval);
134     gsl_matrix_complex_scale (tmpmatrix, gsl_complex_rect(1.0/O815->comargs.nmeas, 0.0) );
135     for (int icorr=0; icorr<4; icorr++)
136       for (int jcorr=0; jcorr<4; jcorr++) {
137         gsl_complex discopart = gsl_complex_mul( gsl_vector_complex_get(totaldisco, icorr),
138                                                  gsl_complex_conjugate( gsl_vector_complex_get(totaldisco, jcorr) )
139                                                  );
140         discopart = gsl_complex_mul_real(discopart, pow(1.0/O815->comargs.nmeas,2));
141         gsl_matrix_complex_set( tmpmatrix, icorr, jcorr, 
142                                 gsl_complex_sub( gsl_matrix_complex_get( tmpmatrix, icorr, jcorr ),
143                                                  discopart
144                                                  )
145                                 );
146       }
147     cdiag(jackres, evecres, tmpmatrix);
148     
149     for (int imeas=0; imeas<O815->comargs.nmeas; imeas++) {
150       gsl_matrix_complex_memcpy (tmpmatrix, totalval);
151       gsl_matrix_complex_sub (tmpmatrix, measurements[imeas][itsep]);
152       gsl_matrix_complex_scale ( tmpmatrix, gsl_complex_rect(1.0/(O815->comargs.nmeas-1), 0.0) );
153       gsl_vector_complex_memcpy( tmpvecc, totaldisco );
154       gsl_vector_complex_sub( tmpvecc, measurements_disco[imeas] );
155       for (int icorr=0; icorr<4; icorr++)
156         for (int jcorr=0; jcorr<4; jcorr++) {
157           gsl_complex discopart = gsl_complex_mul( gsl_vector_complex_get(tmpvecc, icorr),
158                                                    gsl_complex_conjugate( gsl_vector_complex_get(tmpvecc, jcorr) )
159                                                    );
160           discopart = gsl_complex_mul_real(discopart, pow(1.0/(O815->comargs.nmeas-1),2));
161           gsl_matrix_complex_set( tmpmatrix, icorr, jcorr, 
162                                   gsl_complex_sub( gsl_matrix_complex_get( tmpmatrix, icorr, jcorr ),
163                                                    discopart
164                                                    )
165                                   );
166         }
167       cdiag(tmpvec, tmpmatrix2, tmpmatrix);
168       
169       gsl_vector_sub(tmpvec, jackres);
170       gsl_vector_memcpy(tmpvec2, tmpvec);
171       gsl_vector_mul(tmpvec, tmpvec2);
172       gsl_vector_add(jackerrornorm, tmpvec);
173
174       gsl_matrix_complex_sub(tmpmatrix2, evecres);
175       for (int i=0; i<4; i++)
176         for (int j=0; j<4; j++) {
177           gsl_complex tmperr =
178             gsl_complex_rect( pow(GSL_REAL( gsl_matrix_complex_get( tmpmatrix2, i, j ) ),2),
179                               pow(GSL_IMAG( gsl_matrix_complex_get( tmpmatrix2, i, j ) ),2));
180           gsl_matrix_complex_set(tmpmatrix, i, j, tmperr);
181         }
182       gsl_matrix_complex_add( evecerrornorm, tmpmatrix );
183     } 
184     gsl_vector_scale( jackerrornorm,
185                       (double)(O815->comargs.nmeas-1)/O815->comargs.nmeas );
186     gsl_matrix_complex_scale( evecerrornorm, gsl_complex_rect( (double)(O815->comargs.nmeas-1)/O815->comargs.nmeas, 0.0 ) );
187     
188     *out << O815->paraQ->getParaVals();
189     *out << "\t" << itsep;
190     for (int iev = 0; iev < 4; iev++) {
191       *out << "\t" << gsl_vector_get(jackres,iev)
192            << "\t" << sqrt(gsl_vector_get(jackerrornorm,iev));
193       for (int iinter=0; iinter<4; iinter++)
194         *out << "\t" << GSL_REAL(gsl_matrix_complex_get(evecres, iinter, iev))
195              << "\t" << sqrt( GSL_REAL(gsl_matrix_complex_get(evecerrornorm, iinter, iev)) );
196     }
197     *out << endl;
198   }
199
200   gsl_matrix_complex_free(totalval);
201   gsl_matrix_complex_free(tmpmatrix);
202   gsl_matrix_complex_free(evecres);
203   gsl_vector_free(tmpvec);
204   gsl_vector_complex_free(tmpvecc);
205   gsl_vector_free(tmpvec2);
206   gsl_vector_free(jackres);
207   gsl_vector_free(jackerrornorm);
208   gsl_vector_complex_free(totaldisco);
209 };
210
211 void obs_diagcorr::corrCompute()
212 {
213   complex<double> phislice[4][O815->comargs.lsize[1]];
214
215   for (int icorr=0; icorr<4; icorr++)
216     OM[16][icorr] = 0;
217   
218   for (int it = 0; it < O815->comargs.lsize[1]; it++) {
219     for (int icorr=0; icorr<4; icorr++)
220       phislice[icorr][it] = 0;
221     
222     for (int ix = 0; ix < spatialV; ix++) {
223       phislice[0][it] += conj(Sim->phi[ 0*Sim->lsize4 + it*spatialV + ix ]) * Sim->phi[ 0*Sim->lsize4 + it*spatialV + ix ];
224       phislice[1][it] += Sim->phi[ 0*Sim->lsize4 + it*spatialV + ix ] * Sim->phi[ 1*Sim->lsize4 + it*spatialV + ix ];
225       phislice[2][it] += conj(Sim->phi[ 0*Sim->lsize4 + it*spatialV + ix ]) * conj(Sim->phi[ 1*Sim->lsize4 + it*spatialV + ix ]);
226       phislice[3][it] += conj(Sim->phi[ 1*Sim->lsize4 + it*spatialV + ix ]) * Sim->phi[ 1*Sim->lsize4 + it*spatialV + ix ];
227     }
228
229     for (int icorr=0; icorr<4; icorr++) {
230       phislice[icorr][it] /= spatialV;
231       OM[16][icorr] += phislice[icorr][it];
232     }
233   }
234
235   for (int icorr = 0; icorr < 4; icorr++) {
236     for (int jcorr = 0; jcorr < 4; jcorr++) {
237       for (int itsep = 0; itsep < O815->comargs.lsize[1]/2; itsep++) {
238         OM[icorr*4+jcorr][itsep] = 0;
239         
240         for (int it = 0; it < O815->comargs.lsize[1]; it++)
241           OM[icorr*4+jcorr][itsep] += phislice[icorr][ (it+itsep)%O815->comargs.lsize[1] ] * conj( phislice[jcorr][it] );
242       
243         OM[icorr*4+jcorr][itsep] /= O815->comargs.lsize[1];
244       }
245     }
246     
247     OM[16][icorr] /= O815->comargs.lsize[1];
248   }
249 }
250
251 #endif