1 #ifndef OBS_DIAGCORR_HPP
2 #define OBS_DIAGCORR_HPP
4 #include "latlib/o815/o815.h"
6 #include "latlib/writeout.h"
8 #include "latlib/obstat.hpp"
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>
22 class obs_diagcorr : public o815::obs {
25 obs_diagcorr(o815 *_O815);
29 void _meas(bool loadedobs, const int& nthmeas);
33 static complex<double> effMass(vector < complex<double> > *preCalculated, vector< complex<double> > *excludedMeas, int nmeas, void *para);
34 static void preEffMass(vector< vector < complex<double> > > *allVals, vector < complex<double> > *preCalculated, void *para);
40 complex<double> *OM[16+1];
42 gsl_matrix_complex ***measurements;
44 gsl_vector_complex **measurements_disco;
46 static void cdiag (gsl_vector *v, gsl_matrix_complex *m);
49 obs_diagcorr::obs_diagcorr(o815 *_O815) : o815::obs("diagcorr",
50 _O815->paraQ->getParaNames() +
53 _O815, sizeof(complex<double>) * ( 16*_O815->comargs.lsize[1]/2 + 1*4 ) ) {
54 for (int ivar = 0; ivar<16+1; ivar++)
55 OM[ivar] = (complex<double>*)( obsMem + ivar*sizeof(complex<double>)*( _O815->comargs.lsize[1]/2 ) );
57 Sim = (sim*)O815->Sim;
58 spatialV = O815->comargs.lsize[0] * O815->comargs.lsize[0] * O815->comargs.lsize[0];
60 measurements = new gsl_matrix_complex**[O815->comargs.nmeas];
61 for (int imeas = 0; imeas<O815->comargs.nmeas; imeas++) {
62 measurements[imeas] = new gsl_matrix_complex*[_O815->comargs.lsize[1]/2];
63 for (int itsep=0; itsep<_O815->comargs.lsize[1]/2; itsep++)
64 measurements[imeas][itsep] = gsl_matrix_complex_alloc(4,4);
67 measurements_disco = new gsl_vector_complex*[O815->comargs.nmeas];
68 for (int imeas = 0; imeas<O815->comargs.nmeas; imeas++)
69 measurements_disco[imeas] = gsl_vector_complex_alloc(4);
72 void obs_diagcorr::_start() {
73 // *out << O815->comargs.nmeas << endl;
75 //*out << "OBS_test: start" << endl;
78 void obs_diagcorr::_meas(bool loadedobs, const int& nthmeas) {
82 for (int icorr=0; icorr<4; icorr++)
83 for (int jcorr=0; jcorr<4; jcorr++)
84 for (int itsep=0; itsep<O815->comargs.lsize[1]/2; itsep++) {
85 gsl_matrix_complex_set(measurements[nthmeas][itsep], icorr, jcorr,
86 gsl_complex_rect( OM[icorr*4+jcorr][itsep].real(), OM[icorr*4+jcorr][itsep].imag() ));
89 for (int icorr=0; icorr<4; icorr++) {
90 gsl_vector_complex_set(measurements_disco[nthmeas], icorr,
91 gsl_complex_rect( OM[16][icorr].real(), OM[16][icorr].imag() ) );
95 void obs_diagcorr::cdiag (gsl_vector *v, gsl_matrix_complex *m)
97 gsl_eigen_herm_workspace *wspace = gsl_eigen_herm_alloc(100);
98 gsl_eigen_herm(m, v, wspace);
99 gsl_eigen_herm_free(wspace);
102 void obs_diagcorr::_finish() {
103 gsl_matrix_complex *totalval = gsl_matrix_complex_alloc(4,4);
104 gsl_matrix_complex *tmpmatrix = gsl_matrix_complex_alloc(4,4);
105 gsl_vector *tmpvec = gsl_vector_alloc(4);
106 gsl_vector_complex *tmpvecc = gsl_vector_complex_alloc(4);
107 gsl_vector *tmpvec2 = gsl_vector_alloc(4);
108 gsl_vector *jackres = gsl_vector_alloc(4);
109 gsl_vector *jackerrornorm = gsl_vector_alloc(4);
110 gsl_vector_complex *totaldisco = gsl_vector_complex_alloc(4);
112 gsl_vector_complex_set_zero(totaldisco);
113 for (int imeas=0; imeas<O815->comargs.nmeas; imeas++)
114 gsl_vector_complex_add( totaldisco, measurements_disco[imeas] );
116 for (int itsep = 0; itsep < O815->comargs.lsize[1]/2; itsep++) {
117 gsl_matrix_complex_set_zero(totalval);
118 gsl_vector_set_zero(jackerrornorm);
120 for (int imeas=0; imeas<O815->comargs.nmeas; imeas++)
121 gsl_matrix_complex_add( totalval, measurements[imeas][itsep] );
123 gsl_matrix_complex_memcpy (tmpmatrix, totalval);
124 gsl_matrix_complex_scale (tmpmatrix, gsl_complex_rect(1.0/O815->comargs.nmeas, 0.0) );
125 for (int icorr=0; icorr<4; icorr++)
126 for (int jcorr=0; jcorr<4; jcorr++) {
127 gsl_complex discopart = gsl_complex_mul( gsl_vector_complex_get(totaldisco, icorr),
128 gsl_complex_conjugate( gsl_vector_complex_get(totaldisco, jcorr) )
130 discopart = gsl_complex_mul_real(discopart, pow(1.0/O815->comargs.nmeas,2));
131 gsl_matrix_complex_set( tmpmatrix, icorr, jcorr,
132 gsl_complex_sub( gsl_matrix_complex_get( tmpmatrix, icorr, jcorr ),
137 cdiag(jackres, tmpmatrix);
139 for (int imeas=0; imeas<O815->comargs.nmeas; imeas++) {
140 gsl_matrix_complex_memcpy (tmpmatrix, totalval);
141 gsl_matrix_complex_sub (tmpmatrix, measurements[imeas][itsep]);
142 gsl_matrix_complex_scale ( tmpmatrix, gsl_complex_rect(1.0/(O815->comargs.nmeas-1), 0.0) );
143 gsl_vector_complex_memcpy( tmpvecc, totaldisco );
144 gsl_vector_complex_sub( tmpvecc, measurements_disco[imeas] );
145 for (int icorr=0; icorr<4; icorr++)
146 for (int jcorr=0; jcorr<4; jcorr++) {
147 gsl_complex discopart = gsl_complex_mul( gsl_vector_complex_get(tmpvecc, icorr),
148 gsl_complex_conjugate( gsl_vector_complex_get(tmpvecc, jcorr) )
150 discopart = gsl_complex_mul_real(discopart, pow(1.0/(O815->comargs.nmeas-1),2));
151 gsl_matrix_complex_set( tmpmatrix, icorr, jcorr,
152 gsl_complex_sub( gsl_matrix_complex_get( tmpmatrix, icorr, jcorr ),
157 cdiag(tmpvec, tmpmatrix);
158 gsl_vector_sub(tmpvec, jackres);
159 gsl_vector_memcpy(tmpvec2, tmpvec);
160 gsl_vector_mul(tmpvec, tmpvec2);
161 gsl_vector_add(jackerrornorm, tmpvec);
163 gsl_vector_scale( jackerrornorm,
164 (double)(O815->comargs.nmeas-1)/O815->comargs.nmeas );
166 *out << O815->paraQ->getParaVals();
167 *out << "\t" << itsep;
168 for (int iev = 0; iev < 4; iev++) {
169 *out << "\t" << gsl_vector_get(jackres,iev)
170 << "\t" << sqrt(gsl_vector_get(jackerrornorm,iev));
175 gsl_matrix_complex_free(totalval);
176 gsl_matrix_complex_free(tmpmatrix);
177 gsl_vector_free(tmpvec);
178 gsl_vector_complex_free(tmpvecc);
179 gsl_vector_free(tmpvec2);
180 gsl_vector_free(jackres);
181 gsl_vector_free(jackerrornorm);
182 gsl_vector_complex_free(totaldisco);
185 void obs_diagcorr::corrCompute()
187 complex<double> phislice[4][O815->comargs.lsize[1]];
189 for (int icorr=0; icorr<4; icorr++)
192 for (int it = 0; it < O815->comargs.lsize[1]; it++) {
193 for (int icorr=0; icorr<4; icorr++)
194 phislice[icorr][it] = 0;
196 for (int ix = 0; ix < spatialV; ix++) {
197 phislice[0][it] += conj(Sim->phi[ 0*Sim->lsize4 + it*spatialV + ix ]) * Sim->phi[ 0*Sim->lsize4 + it*spatialV + ix ];
198 phislice[1][it] += Sim->phi[ 0*Sim->lsize4 + it*spatialV + ix ] * Sim->phi[ 1*Sim->lsize4 + it*spatialV + ix ];
199 phislice[2][it] += conj(Sim->phi[ 0*Sim->lsize4 + it*spatialV + ix ]) * conj(Sim->phi[ 1*Sim->lsize4 + it*spatialV + ix ]);
200 phislice[3][it] += conj(Sim->phi[ 1*Sim->lsize4 + it*spatialV + ix ]) * Sim->phi[ 1*Sim->lsize4 + it*spatialV + ix ];
203 for (int icorr=0; icorr<4; icorr++) {
204 phislice[icorr][it] /= spatialV;
205 OM[16][icorr] += phislice[icorr][it];
209 for (int icorr = 0; icorr < 4; icorr++) {
210 for (int jcorr = 0; jcorr < 4; jcorr++) {
211 for (int itsep = 0; itsep < O815->comargs.lsize[1]/2; itsep++) {
212 OM[icorr*4+jcorr][itsep] = 0;
214 for (int it = 0; it < O815->comargs.lsize[1]; it++)
215 OM[icorr*4+jcorr][itsep] += phislice[icorr][ (it+itsep)%O815->comargs.lsize[1] ] * conj( phislice[jcorr][it] );
217 OM[icorr*4+jcorr][itsep] /= O815->comargs.lsize[1];
221 OM[16][icorr] /= O815->comargs.lsize[1];