4 #define BESSELCACHE 10000
7 #include <gsl/gsl_rng.h>
8 #include <gsl/gsl_errno.h>
9 #include <gsl/gsl_integration.h>
12 #include <boost/math/special_functions/bessel.hpp>
14 #include "latlib/neigh.h"
16 class sim : public o815::sim {
22 double kappa[2], lambda[2], Mu[2], beta;
23 int wArg(int x, int i);
24 double WF(int nom, int denom, int iflav);
31 int lineUpdate(int x, int vardir);
32 int cubeUpdate(int x, int orient);
33 void getCube(int cube[3], int orient);
34 int plaqUpdate(int x, int mu, int nu, int i);
35 int iPlaq(int i, int j);
36 int lpUpdate(int x, int mu, int i);
38 static double wkern(double t, void *params);
39 double wGsl(int n, int iflav);
44 bcache besselCache[BESSELCACHE];
55 gsl_integration_cquad_workspace *w;
56 double rhoLine(int xstart, int vardir, int modifier);
57 double rhoPlaq_Fac(int x, int mu, int i, int modifier);
58 double rhoPlaq(int x, int mu, int nu, int i, int modifier);
59 double rhoLp(int x, int mu, int i, int modifier);
60 bool neverComputedB, neverComputedW[2];
61 double rhoCube(int x, int *cube, int modifier);
62 double wLastLambda[2], wLastKappa[2];
63 double lastBesselBeta;
64 int ipow(int base, int expo);
65 wcache wCache[2][WCACHE*WCACHE];
66 double WFl0(int nom, int denom, int iflav);
69 sim::sim(o815 *_O815) : o815::sim( _O815,
71 _O815->comargs.lsize[0]*_O815->comargs.lsize[0]*_O815->comargs.lsize[0]*_O815->comargs.lsize[1]*(2*4+2*4+6) ) {
75 lsize4 = _O815->comargs.lsize[0]*_O815->comargs.lsize[0]*_O815->comargs.lsize[0]*_O815->comargs.lsize[1];
77 nb = new neigh(4, _O815->comargs.lsize[0], _O815->comargs.lsize[0], _O815->comargs.lsize[0], _O815->comargs.lsize[1]);
80 lp = (int*)(confMem + sizeof(int)*2*lsize4*4);
81 p = (int*)(confMem + sizeof(int)*2*lsize4*4*2);
83 neverComputedB = true;
84 neverComputedW[0] = true;
85 neverComputedW[1] = true;
87 w = gsl_integration_cquad_workspace_alloc(10000);
89 gettimeofday(&tv, NULL);
90 rangsl = gsl_rng_alloc(gsl_rng_ranlxs0);
91 gsl_rng_set(rangsl, 1000000 * tv.tv_sec + tv.tv_usec);
94 void sim::_makeSweep()
98 for( int x=0; x<lsize4; x++) {
99 for(int orient=0; orient<4; orient++)
100 cubeUpdate(x, orient);
102 for(int mu=0; mu<3; mu++)
103 for(int nu=mu+1; nu<4; nu++)
104 for(int i=0; i<2; i++)
105 plaqUpdate(x, mu, nu, i);
107 for(int mu=0; mu<4; mu++)
108 for(int i=0; i<2; i++)
113 void sim::_newParas() {
114 kappa[0] = (*O815->paraQ)["kappaone"];
115 kappa[1] = (*O815->paraQ)["kappatwo"];
116 lambda[0] = (*O815->paraQ)["lambdaone"];
117 lambda[1] = (*O815->paraQ)["lambdatwo"];
118 beta = (*O815->paraQ)["beta"];
119 Mu[0] = (*O815->paraQ)["muone"];
120 Mu[1] = (*O815->paraQ)["mutwo"];
122 memset(confMem, 0, sizeof(int)*lsize4*(2*4+2*4+6));
127 for (int vardir = 0; vardir < 4; vardir++)
128 for (int x0 = 0; x0 < O815->comargs.lsize[0]; x0++) {
129 for (int x1 = 0; x1 < O815->comargs.lsize[0]; x1++) {
130 for (int x2 = 0; x2 < O815->comargs.lsize[0]; x2++) {
131 for (int x3 = 0; x3 < O815->comargs.lsize[1]; x3++) {
133 + x1 * O815->comargs.lsize[0]
134 + x2 * O815->comargs.lsize[0] * O815->comargs.lsize[0]
135 + x3 * O815->comargs.lsize[0] * O815->comargs.lsize[0] * O815->comargs.lsize[0], vardir);
146 int sim::lineUpdate(int x, int vardir)
148 int modifier = 1 - 2 * gsl_rng_uniform_int(rangsl, 2);
150 if ( gsl_rng_uniform(rangsl) < rhoLine(x, vardir, modifier) )
152 const int stepsize = ipow(O815->comargs.lsize[0], vardir);
154 for (int xvar = 0; xvar < O815->comargs.lsize[vardir/3]; xvar++)
155 for (int i = 0; i < 2; i++)
156 l[ i*lsize4*4 + (x+xvar*stepsize)*4 + vardir ] += modifier;
164 int sim::cubeUpdate(int x, int orient)
166 int modifier = 1 - 2 * gsl_rng_uniform_int(rangsl, 2);
168 getCube(cube, orient);
170 if ( gsl_rng_uniform(rangsl) < rhoCube(x, cube, modifier) )
172 p[ x*6 + iPlaq( cube[0],cube[1] ) ] += modifier;
173 p[ x*6 + iPlaq( cube[0],cube[2] ) ] -= modifier;
174 p[ x*6 + iPlaq( cube[1],cube[2] ) ] += modifier;
175 p[ (*nb)[x*8+cube[0]]*6 + iPlaq( cube[1],cube[2] ) ] -= modifier;
176 p[ (*nb)[x*8+cube[1]]*6 + iPlaq( cube[0],cube[2] ) ] += modifier;
177 p[ (*nb)[x*8+cube[2]]*6 + iPlaq( cube[0],cube[1] ) ] -= modifier;
185 void sim::getCube(int cube[3], int orient)
188 for (int i=0; i<4; i++) {
197 int sim::plaqUpdate(int x, int mu, int nu, int i)
199 int modifier = 1 - 2 * gsl_rng_uniform_int(rangsl, 2);
200 int lmodsign = 1 - 2*i;
202 if ( gsl_rng_uniform(rangsl) < rhoPlaq(x, mu, nu, i, modifier) )
204 p[ x*6 + iPlaq(mu,nu) ] = p[ x*6 + iPlaq(mu,nu) ] + modifier;
205 l[ i*lsize4*4 + x*4 + mu ] = l[ i*lsize4*4 + x*4 + mu ] - lmodsign*modifier;
206 l[ i*lsize4*4 + x*4 + nu ] = l[ i*lsize4*4 + x*4 + nu ] + lmodsign*modifier;
207 l[ i*lsize4*4 + (*nb)[x*8+nu]*4 + mu ] = l[ i*lsize4*4 + (*nb)[x*8+nu]*4 + mu ] + lmodsign*modifier;
208 l[ i*lsize4*4 + (*nb)[x*8+mu]*4 + nu ] = l[ i*lsize4*4 + (*nb)[x*8+mu]*4 + nu ] - lmodsign*modifier;
215 int sim::iPlaq(int i, int j)
217 if ( i == 0 ) return j-1;
218 else if ( i == 1 ) return j+1;
222 int sim::lpUpdate(int x, int mu, int i)
224 int modifier = 1 - 2 * gsl_rng_uniform_int(rangsl, 2);
226 if ( (lp[ i*lsize4*4 + x*4 + mu ] + modifier) < 0 )
229 if ( gsl_rng_uniform(rangsl) < rhoLp(x, mu, i, modifier) ) {
230 lp[ i*lsize4*4 + x*4 + mu ] += modifier;
239 /* http://www.boost.org/doc/libs/1_35_0/libs/math/doc/sf_and_dist/html/math_toolkit/special/bessel/mbessel.html */
243 if ( beta != lastBesselBeta || neverComputedB )
245 memset(besselCache, 0, sizeof(bcache)*BESSELCACHE);
246 lastBesselBeta = beta;
247 neverComputedB = false;
250 if ( p < BESSELCACHE )
252 if ( besselCache[p].set ) return besselCache[p].val;
255 besselCache[p].val = boost::math::cyl_bessel_i(p, beta);
256 besselCache[p].set = true;
257 return besselCache[p].val;
262 return boost::math::cyl_bessel_i(p, beta);
266 int sim::wArg(int x, int i)
270 for (int nu=0; nu<4; nu++)
271 warg += abs( l[ i*lsize4*4 + x*4 + nu ] )
272 + 2 * ( lp[ i*lsize4*4 + x*4 + nu ] + lp[ i*lsize4*4 + (*nb)[x*8+nu+4]*4 + nu ] )
273 + abs( l[ i*lsize4*4 + (*nb)[x*8+nu+4]*4 + nu ] );
278 double sim::WF(int nom, int denom, int iflav)
280 if ( nom == denom ) return 1;
282 if ( wLastLambda[iflav] != lambda[iflav] || wLastKappa[iflav] != kappa[iflav] || neverComputedW[iflav] )
284 memset(wCache[iflav], 0, sizeof(wcache)*WCACHE*WCACHE);
285 wLastLambda[iflav] = lambda[iflav];
286 wLastKappa[iflav] = kappa[iflav];
287 neverComputedW[iflav] = false;
290 if ( nom/2 < WCACHE && denom/2 < WCACHE )
292 if ( wCache[iflav][ (nom/2)*WCACHE + denom/2 ].set )
293 return wCache[iflav][ (nom/2)*WCACHE + denom/2 ].val;
296 /* integral can be solved analytically for lambda==0 */
298 wCache[iflav][ (nom/2)*WCACHE + denom/2 ].val = WFl0(nom, denom, iflav);
300 wCache[iflav][ (nom/2)*WCACHE + denom/2 ].val = wGsl(nom, iflav) / wGsl(denom, iflav);
302 wCache[iflav][ (nom/2)*WCACHE + denom/2 ].set = true;
304 return wCache[iflav][ (nom/2)*WCACHE + denom/2 ].val;
310 return WFl0(nom, denom, iflav);
312 return wGsl(nom, iflav) / wGsl(denom, iflav);
316 double sim::wkern(double t, void *params)
318 kernpara *para = (kernpara*)params;
319 double wkern = pow((1-t)/t, para->n + 1) * exp( -para->kappa*pow((1-t)/t, 2) -para->lambda*pow((1-t)/t, 4) ) / pow(t,2);
323 double sim::wGsl(int n, int iflav)
331 para.kappa = kappa[iflav];
332 para.lambda = lambda[iflav];
334 gslF.function = &wkern;
337 status = gsl_integration_cquad(&gslF, 0, 1, 0, 1e-7, w, &result, NULL, NULL);
342 double sim::rhoLine(int xstart, int vardir, int modifier)
347 const int stepsize = ipow(O815->comargs.lsize[0],vardir);
348 for(int xvar=0; xvar<O815->comargs.lsize[vardir/3]; xvar++)
349 for(int i=0; i<2; i++)
351 const int x = xstart + xvar*stepsize;
353 if( abs( l[ i*lsize4*4 + x*4 + vardir ] + modifier ) > abs( l[ i*lsize4*4 + x*4 + vardir ] ) )
354 rho /= abs( l[ i*lsize4*4 + x*4 + vardir ] + modifier ) + lp[ i*lsize4*4 + x*4 + vardir ];
356 rho *= abs( l[ i*lsize4*4 + x*4 + vardir ] ) + lp[ i*lsize4*4 + x*4 + vardir ];
360 - abs( l[ i*lsize4*4 + x*4 + vardir ] )
361 + abs( l[ i*lsize4*4 + x*4 + vardir ] + modifier )
362 - abs( l[ i*lsize4*4 + (*nb)[x*8+vardir+4]*4 + vardir ] )
363 + abs( l[ i*lsize4*4 + (*nb)[x*8+vardir+4]*4 + vardir ] + modifier ),
367 rho *= exp(modifier*Mu[i]);
373 double sim::rhoPlaq(int x, int mu, int nu, int i, int modifier)
376 int lmodsign = 1 - 2*i;
379 rho *= I( p[ x*6 + iPlaq(mu,nu) ] + modifier ) / I( p[ x*6 + iPlaq(mu,nu) ] );
383 - abs( l[ i*lsize4*4 + x*4 + mu] ) - abs( l[ i*lsize4*4 + x*4 + nu ] )
384 + abs( l[ i*lsize4*4 + x*4 + mu] - lmodsign*modifier ) + abs( l[ i*lsize4*4 + x*4 + nu ] + lmodsign*modifier ),
387 tmpwarg = wArg((*nb)[x*8+mu],i);
389 - abs( l[ i*lsize4*4 + x*4 + mu] ) - abs( l[ i*lsize4*4 + (*nb)[x*8+mu]*4 + nu ] )
390 + abs( l[ i*lsize4*4 + x*4 + mu] - lmodsign*modifier ) + abs( l[ i*lsize4*4 + (*nb)[x*8+mu]*4 + nu ] - lmodsign*modifier ),
393 tmpwarg = wArg((*nb)[(*nb)[x*8+mu]*8+nu],i);
395 - abs( l[ i*lsize4*4 + (*nb)[x*8+mu]*4 + nu ] ) - abs( l[ i*lsize4*4 + (*nb)[x*8+nu]*4 + mu ] )
396 + abs( l[ i*lsize4*4 + (*nb)[x*8+mu]*4 + nu ] - lmodsign*modifier ) + abs( l[ i*lsize4*4 + (*nb)[x*8+nu]*4 + mu ] + lmodsign*modifier ),
399 tmpwarg = wArg((*nb)[x*8+nu],i);
401 - abs( l[ i*lsize4*4 + (*nb)[x*8+nu]*4 + mu ] ) - abs( l[ i*lsize4*4 + x*4 + nu ] )
402 + abs( l[ i*lsize4*4 + (*nb)[x*8+nu]*4 + mu ] + lmodsign*modifier ) + abs( l[ i*lsize4*4 + x*4 + nu ] + lmodsign*modifier ),
405 rho *= rhoPlaq_Fac(x, mu, i, -lmodsign*modifier);
406 rho *= rhoPlaq_Fac(x, nu, i, +lmodsign*modifier);
407 rho *= rhoPlaq_Fac((*nb)[x*8+mu], nu, i, -lmodsign*modifier);
408 rho *= rhoPlaq_Fac((*nb)[x*8+nu], mu, i, +lmodsign*modifier);
413 double sim::rhoPlaq_Fac(int x, int mu, int i, int modifier)
415 if( abs( l[ i*lsize4*4 + x*4 + mu ] + modifier ) > abs( l[ i*lsize4*4 + x*4 + mu] ) )
416 return 1.0 / ( abs( l[ i*lsize4*4 + x*4 + mu ] + modifier ) + lp[ i*lsize4*4 + x*4 + mu ] );
418 return ( abs( l[ i*lsize4*4 + x*4 + mu ] ) + lp[ i*lsize4*4 + x*4 + mu ] );
421 double sim::rhoLp(int x, int mu, int i, int modifier)
427 rho /= ( abs(l[ i*lsize4*4 + x*4 + mu ]) + lp[ i*lsize4*4 + x*4 + mu ] + modifier ) * ( lp[ i*lsize4*4 + x*4 + mu ] + modifier );
429 rho *= ( abs(l[ i*lsize4*4 + x*4 + mu ]) + lp[ i*lsize4*4 + x*4 + mu ] ) * lp[ i*lsize4*4 + x*4 + mu ];
432 rho *= WF( tmpwarg + 2*modifier, tmpwarg, i );
433 tmpwarg = wArg((*nb)[x*8+mu],i);
434 rho *= WF( tmpwarg + 2*modifier, tmpwarg, i );
439 double sim::rhoCube(int x, int *cube, int modifier)
442 rho *= I( p[ x*6 + iPlaq( cube[0],cube[1] ) ] + modifier ) / I( p[ x*6 + iPlaq( cube[0],cube[1] ) ] );
443 rho *= I( p[ x*6 + iPlaq( cube[0],cube[2] ) ] - modifier ) / I( p[ x*6 + iPlaq( cube[0],cube[2] ) ] );
444 rho *= I( p[ x*6 + iPlaq( cube[1],cube[2] ) ] + modifier ) / I( p[ x*6 + iPlaq( cube[1],cube[2] ) ] );
445 rho *= I( p[ (*nb)[x*8+cube[0]]*6 + iPlaq( cube[1],cube[2] ) ] - modifier ) / I( p[ (*nb)[x*8+cube[0]]*6 + iPlaq( cube[1],cube[2] ) ] );
446 rho *= I( p[ (*nb)[x*8+cube[1]]*6 + iPlaq( cube[0],cube[2] ) ] + modifier ) / I( p[ (*nb)[x*8+cube[1]]*6 + iPlaq( cube[0],cube[2] ) ] );
447 rho *= I( p[ (*nb)[x*8+cube[2]]*6 + iPlaq( cube[0],cube[1] ) ] - modifier ) / I( p[ (*nb)[x*8+cube[2]]*6 + iPlaq( cube[0],cube[1] ) ] );
451 int sim::ipow(int base, int expo)
462 double sim::WFl0(int nom, int denom, int iflav)
467 for (int i=nom/2; i>denom/2; i--)
469 return wfl0 * pow(kappa[iflav], 0.5*(denom-nom));
472 for (int i=denom/2; i>nom/2; i--)
474 return (1.0 / wfl0) * pow(kappa[iflav], 0.5*(denom-nom));