void sw_spinor(const int ieo, spinor * const kk, spinor * const ll) { int ioff; int icx; int x; spinor *r,*s; static su3 v0,v1,v2,v3; static su3 u0,u1,u2,u3; static su3 lswp[4],lswm[4]; if(ieo == 0) { ioff=0; } else { ioff=(VOLUME+RAND)/2; } /************************ loop over half of the lattice sites ***********/ for(icx = ioff; icx < (VOLUME/2+ioff); icx++) { x = g_eo2lexic[icx]; r = kk + icx - ioff; s = ll + icx - ioff; _vector_tensor_vector(v0,(*r).s0,(*s).s0); _vector_tensor_vector(v1,(*r).s0,(*s).s1); _vector_tensor_vector(v2,(*r).s1,(*s).s1); _vector_tensor_vector(v3,(*r).s1,(*s).s0); _vector_tensor_vector(u0,(*r).s2,(*s).s2); _vector_tensor_vector(u1,(*r).s2,(*s).s3); _vector_tensor_vector(u2,(*r).s3,(*s).s3); _vector_tensor_vector(u3,(*r).s3,(*s).s2); /* compute the insertion matrix */ _su3_plus_su3(lswp[0],u0,v0); _su3_plus_su3(lswp[1],u1,v1); _su3_plus_su3(lswp[2],u2,v2); _su3_plus_su3(lswp[3],u3,v3); _su3_minus_su3(lswm[0],u0,v0); _su3_minus_su3(lswm[1],u1,v1); _su3_minus_su3(lswm[2],u2,v2); _su3_minus_su3(lswm[3],u3,v3); /* add up the swm[0] and swp[0] */ _su3_acc(swm[x][0], lswm[0]); _su3_acc(swm[x][1], lswm[1]); _su3_acc(swm[x][2], lswm[2]); _su3_acc(swm[x][3], lswm[3]); _su3_acc(swp[x][0], lswp[0]); _su3_acc(swp[x][1], lswp[1]); _su3_acc(swp[x][2], lswp[2]); _su3_acc(swp[x][3], lswp[3]); } return; }
void sw_term(su3 ** const gf, const double kappa, const double c_sw) { int k,l; int x,xpk,xpl,xmk,xml,xpkml,xplmk,xmkml; su3 *w1,*w2,*w3,*w4; double ka_csw_8 = kappa*c_sw/8.; static su3 v1,v2,plaq; static su3 fkl[4][4]; static su3 magnetic[4],electric[4]; static su3 aux; /* compute the clover-leave */ /* l __ __ | | | | |__| |__| __ __ | | | | |__| |__| k */ for(x = 0; x < VOLUME; x++) { for(k = 0; k < 4; k++) { for(l = k+1; l < 4; l++) { xpk=g_iup[x][k]; xpl=g_iup[x][l]; xmk=g_idn[x][k]; xml=g_idn[x][l]; xpkml=g_idn[xpk][l]; xplmk=g_idn[xpl][k]; xmkml=g_idn[xml][k]; w1=&gf[x][k]; w2=&gf[xpk][l]; w3=&gf[xpl][k]; w4=&gf[x][l]; _su3_times_su3(v1,*w1,*w2); _su3_times_su3(v2,*w4,*w3); _su3_times_su3d(plaq,v1,v2); w1=&gf[x][l]; w2=&gf[xplmk][k]; w3=&gf[xmk][l]; w4=&gf[xmk][k]; _su3_times_su3d(v1,*w1,*w2); _su3d_times_su3(v2,*w3,*w4); _su3_times_su3_acc(plaq,v1,v2); w1=&gf[xmk][k]; w2=&gf[xmkml][l]; w3=&gf[xmkml][k]; w4=&gf[xml][l]; _su3_times_su3(v1,*w2,*w1); _su3_times_su3(v2,*w3,*w4); _su3d_times_su3_acc(plaq,v1,v2); w1=&gf[xml][l]; w2=&gf[xml][k]; w3=&gf[xpkml][l]; w4=&gf[x][k]; _su3d_times_su3(v1,*w1,*w2); _su3_times_su3d(v2,*w3,*w4); _su3_times_su3_acc(plaq,v1,v2); _su3_dagger(v2,plaq); _su3_minus_su3(fkl[k][l],plaq,v2); } } // this is the one in flavour and colour space // twisted mass term is treated in clover, sw_inv and // clover_gamma5 _su3_one(sw[x][0][0]); _su3_one(sw[x][2][0]); _su3_one(sw[x][0][1]); _su3_one(sw[x][2][1]); for(k = 1; k < 4; k++) { _su3_assign(electric[k], fkl[0][k]); } _su3_assign(magnetic[1], fkl[2][3]); _su3_minus_assign(magnetic[2], fkl[1][3]); _su3_assign(magnetic[3], fkl[1][2]); /* upper left block 6x6 matrix */ _itimes_su3_minus_su3(aux,electric[3],magnetic[3]); _su3_refac_acc(sw[x][0][0],ka_csw_8,aux); _itimes_su3_minus_su3(aux,electric[1],magnetic[1]); _su3_minus_su3(v2,electric[2],magnetic[2]); _su3_acc(aux,v2); _real_times_su3(sw[x][1][0],ka_csw_8,aux); _itimes_su3_minus_su3(aux,magnetic[3],electric[3]); _su3_refac_acc(sw[x][2][0],ka_csw_8,aux); /* lower right block 6x6 matrix */ _itimes_su3_plus_su3(aux,electric[3],magnetic[3]); _su3_refac_acc(sw[x][0][1],(-ka_csw_8),aux); _itimes_su3_plus_su3(aux,electric[1],magnetic[1]); _su3_plus_su3(v2,electric[2],magnetic[2]); _su3_acc(aux,v2); _real_times_su3(sw[x][1][1],(-ka_csw_8),aux); _itimes_su3_plus_su3(aux,magnetic[3],electric[3]); _su3_refac_acc(sw[x][2][1],ka_csw_8,aux); } return; }
void sw_all(hamiltonian_field_t * const hf, const double kappa, const double c_sw) { int k,l; int x,xpk,xpl,xmk,xml,xpkml,xplmk,xmkml; su3 *w1,*w2,*w3,*w4; double ka_csw_8 = kappa*c_sw/8.; static su3 v1,v2,vv1,vv2,plaq; static su3 vis[4][4]; for(x = 0; x < VOLUME; x++) { _minus_itimes_su3_plus_su3(vis[0][1],swm[x][1],swm[x][3]); _su3_minus_su3(vis[0][2],swm[x][1],swm[x][3]); _itimes_su3_minus_su3(vis[0][3],swm[x][2],swm[x][0]); _minus_itimes_su3_plus_su3(vis[2][3],swp[x][1],swp[x][3]); _su3_minus_su3(vis[1][3],swp[x][3],swp[x][1]); _itimes_su3_minus_su3(vis[1][2],swp[x][2],swp[x][0]); // project to the traceless anti-hermitian part _su3_dagger(v1,vis[0][1]); _su3_minus_su3(vis[0][1],vis[0][1],v1); _su3_dagger(v1,vis[0][2]); _su3_minus_su3(vis[0][2],vis[0][2],v1); _su3_dagger(v1,vis[0][3]); _su3_minus_su3(vis[0][3],vis[0][3],v1); _su3_dagger(v1,vis[2][3]); _su3_minus_su3(vis[2][3],vis[2][3],v1); _su3_dagger(v1,vis[1][3]); _su3_minus_su3(vis[1][3],vis[1][3],v1); _su3_dagger(v1,vis[1][2]); _su3_minus_su3(vis[1][2],vis[1][2],v1); for(k = 0; k < 4; k++) { for(l = k+1; l < 4; l++) { xpk=g_iup[x][k]; xpl=g_iup[x][l]; xmk=g_idn[x][k]; xml=g_idn[x][l]; xpkml=g_idn[xpk][l]; xplmk=g_idn[xpl][k]; xmkml=g_idn[xml][k]; w1=&hf->gaugefield[x][k]; w2=&hf->gaugefield[xpk][l]; w3=&hf->gaugefield[xpl][k]; /*dag*/ w4=&hf->gaugefield[x][l]; /*dag*/ _su3_times_su3(v1,*w1,*w2); _su3_times_su3(v2,*w4,*w3); _su3_times_su3d(plaq,v1,v2); _su3_times_su3(vv1,plaq,vis[k][l]); _trace_lambda_mul_add_assign(hf->derivative[x][k], -2.*ka_csw_8, vv1); _su3d_times_su3(vv2,*w1,vv1); _su3_times_su3(vv1,vv2,*w1); _trace_lambda_mul_add_assign(hf->derivative[xpk][l], -2.*ka_csw_8, vv1); _su3_times_su3(vv2,vis[k][l],plaq); _su3_dagger(vv1,vv2); _trace_lambda_mul_add_assign(hf->derivative[x][l], -2.*ka_csw_8, vv1); _su3d_times_su3(vv2,*w4,vv1); _su3_times_su3(vv1,vv2,*w4); _trace_lambda_mul_add_assign(hf->derivative[xpl][k], -2.*ka_csw_8, vv1); w1=&hf->gaugefield[x][l]; w2=&hf->gaugefield[xplmk][k]; /*dag*/ w3=&hf->gaugefield[xmk][l]; /*dag*/ w4=&hf->gaugefield[xmk][k]; _su3_times_su3d(v1,*w1,*w2); _su3d_times_su3(v2,*w3,*w4); _su3_times_su3(plaq,v1,v2); _su3_times_su3(vv1,plaq,vis[k][l]); _trace_lambda_mul_add_assign(hf->derivative[x][l], -2.*ka_csw_8, vv1); _su3_dagger(vv1,v1); _su3_times_su3d(vv2,vv1,vis[k][l]); _su3_times_su3d(vv1,vv2,v2); _trace_lambda_mul_add_assign(hf->derivative[xplmk][k], -2.*ka_csw_8, vv1); _su3_times_su3(vv2,*w3,vv1); _su3_times_su3d(vv1,vv2,*w3); _trace_lambda_mul_add_assign(hf->derivative[xmk][l], -2.*ka_csw_8, vv1); _su3_dagger(vv2,vv1); _trace_lambda_mul_add_assign(hf->derivative[xmk][k], -2.*ka_csw_8, vv2); w1=&hf->gaugefield[xmk][k]; /*dag*/ w2=&hf->gaugefield[xmkml][l]; /*dag*/ w3=&hf->gaugefield[xmkml][k]; w4=&hf->gaugefield[xml][l]; _su3_times_su3(v1,*w2,*w1); _su3_times_su3(v2,*w3,*w4); _su3_times_su3d(vv1,*w1,vis[k][l]); _su3_times_su3d(vv2,vv1,v2); _su3_times_su3(vv1,vv2,*w2); _trace_lambda_mul_add_assign(hf->derivative[xmk][k], -2.*ka_csw_8, vv1); _su3_times_su3(vv2,*w2,vv1); _su3_times_su3d(vv1,vv2,*w2); _trace_lambda_mul_add_assign(hf->derivative[xmkml][l], -2.*ka_csw_8, vv1); _su3_dagger(vv2,vv1); _trace_lambda_mul_add_assign(hf->derivative[xmkml][k], -2.*ka_csw_8, vv2); _su3d_times_su3(vv1,*w3,vv2); _su3_times_su3(vv2,vv1,*w3); _trace_lambda_mul_add_assign(hf->derivative[xml][l], -2.*ka_csw_8, vv2); w1=&hf->gaugefield[xml][l]; /*dag*/ w2=&hf->gaugefield[xml][k]; w3=&hf->gaugefield[xpkml][l]; w4=&hf->gaugefield[x][k]; /*dag*/ _su3d_times_su3(v1,*w1,*w2); _su3_times_su3d(v2,*w3,*w4); _su3_times_su3d(vv1,*w1,vis[k][l]); _su3_times_su3d(vv2,vv1,v2); _su3_times_su3d(vv1,vv2,*w2); _trace_lambda_mul_add_assign(hf->derivative[xml][l], -2.*ka_csw_8, vv1); _su3_dagger(vv2,vv1); _trace_lambda_mul_add_assign(hf->derivative[xml][k], -2.*ka_csw_8, vv2); _su3d_times_su3(vv1,*w2,vv2); _su3_times_su3(vv2,vv1,*w2); _trace_lambda_mul_add_assign(hf->derivative[xpkml][l], -2.*ka_csw_8, vv2); _su3_dagger(vv2,v2); _su3_times_su3d(vv1,vv2,v1); _su3_times_su3d(vv2,vv1,vis[k][l]); _trace_lambda_mul_add_assign(hf->derivative[x][k], -2.*ka_csw_8, vv2); } } } return; }
void sw_deriv(const int ieo, const double mu) { int ioff; int x; double fac = 1.0000; static su3 lswp[4],lswm[4]; /* convention: Tr clover-leaf times insertion */ if(ieo == 0) { ioff=0; } else { ioff = (VOLUME+RAND)/2; } if(fabs(mu) > 0.) fac = 0.5; for(int icx = ioff, icy=0; icx < (VOLUME/2+ioff); icx++, icy++) { x = g_eo2lexic[icx]; /* compute the insertion matrix */ _su3_plus_su3(lswp[0],sw_inv[icy][0][1],sw_inv[icy][0][0]); _su3_plus_su3(lswp[1],sw_inv[icy][1][1],sw_inv[icy][1][0]); _su3_plus_su3(lswp[2],sw_inv[icy][2][1],sw_inv[icy][2][0]); _su3_plus_su3(lswp[3],sw_inv[icy][3][1],sw_inv[icy][3][0]); _su3_minus_su3(lswm[0],sw_inv[icy][0][1],sw_inv[icy][0][0]); _su3_minus_su3(lswm[1],sw_inv[icy][1][1],sw_inv[icy][1][0]); _su3_minus_su3(lswm[2],sw_inv[icy][2][1],sw_inv[icy][2][0]); _su3_minus_su3(lswm[3],sw_inv[icy][3][1],sw_inv[icy][3][0]); /* add up to swm[] and swp[] */ _su3_refac_acc(swm[x][0], fac, lswm[0]); _su3_refac_acc(swm[x][1], fac, lswm[1]); _su3_refac_acc(swm[x][2], fac, lswm[2]); _su3_refac_acc(swm[x][3], fac, lswm[3]); _su3_refac_acc(swp[x][0], fac, lswp[0]); _su3_refac_acc(swp[x][1], fac, lswp[1]); _su3_refac_acc(swp[x][2], fac, lswp[2]); _su3_refac_acc(swp[x][3], fac, lswp[3]); if(fabs(mu) > 0.) { /* compute the insertion matrix */ _su3_plus_su3(lswp[0],sw_inv[icy+VOLUME/2][0][1],sw_inv[icy+VOLUME/2][0][0]); _su3_plus_su3(lswp[1],sw_inv[icy+VOLUME/2][1][1],sw_inv[icy+VOLUME/2][1][0]); _su3_plus_su3(lswp[2],sw_inv[icy+VOLUME/2][2][1],sw_inv[icy+VOLUME/2][2][0]); _su3_plus_su3(lswp[3],sw_inv[icy+VOLUME/2][3][1],sw_inv[icy+VOLUME/2][3][0]); _su3_minus_su3(lswm[0],sw_inv[icy+VOLUME/2][0][1],sw_inv[icy+VOLUME/2][0][0]); _su3_minus_su3(lswm[1],sw_inv[icy+VOLUME/2][1][1],sw_inv[icy+VOLUME/2][1][0]); _su3_minus_su3(lswm[2],sw_inv[icy+VOLUME/2][2][1],sw_inv[icy+VOLUME/2][2][0]); _su3_minus_su3(lswm[3],sw_inv[icy+VOLUME/2][3][1],sw_inv[icy+VOLUME/2][3][0]); /* add up to swm[] and swp[] */ _su3_refac_acc(swm[x][0], fac, lswm[0]); _su3_refac_acc(swm[x][1], fac, lswm[1]); _su3_refac_acc(swm[x][2], fac, lswm[2]); _su3_refac_acc(swm[x][3], fac, lswm[3]); _su3_refac_acc(swp[x][0], fac, lswp[0]); _su3_refac_acc(swp[x][1], fac, lswp[1]); _su3_refac_acc(swp[x][2], fac, lswp[2]); _su3_refac_acc(swp[x][3], fac, lswp[3]); } } return; }
int update_tm(double *plaquette_energy, double *rectangle_energy, char * filename, const int return_check, const int acctest, const int traj_counter) { su3 *v, *w; static int ini_g_tmp = 0; int accept, i=0, j=0, iostatus=0; double yy[1]; double dh, expmdh, ret_dh=0., ret_gauge_diff=0., tmp; double atime=0., etime=0.; double ks = 0., kc = 0., ds, tr, ts, tt; char tmp_filename[50]; /* Energy corresponding to the Gauge part */ double new_plaquette_energy=0., new_rectangle_energy = 0.; /* Energy corresponding to the Momenta part */ double enep=0., enepx=0., ret_enep = 0.; /* Energy corresponding to the pseudo fermion part(s) */ FILE * datafile=NULL, * ret_check_file=NULL; hamiltonian_field_t hf; paramsXlfInfo *xlfInfo; hf.gaugefield = g_gauge_field; hf.momenta = moment; hf.derivative = df0; hf.update_gauge_copy = g_update_gauge_copy; hf.update_gauge_energy = g_update_gauge_energy; hf.update_rectangle_energy = g_update_rectangle_energy; hf.traj_counter = traj_counter; integrator_set_fields(&hf); strcpy(tmp_filename, ".conf.tmp"); if(ini_g_tmp == 0) { ini_g_tmp = init_gauge_tmp(VOLUME); if(ini_g_tmp != 0) { exit(-1); } ini_g_tmp = 1; } atime = gettime(); /* * here the momentum and spinor fields are initialized * and their respective actions are calculated */ /* * copy the gauge field to gauge_tmp */ #ifdef OMP #pragma omp parallel for private(w,v) #endif for(int ix=0;ix<VOLUME;ix++) { for(int mu=0;mu<4;mu++) { v=&hf.gaugefield[ix][mu]; w=&gauge_tmp[ix][mu]; _su3_assign(*w,*v); } } /* heatbath for all monomials */ for(i = 0; i < Integrator.no_timescales; i++) { for(j = 0; j < Integrator.no_mnls_per_ts[i]; j++) { monomial_list[ Integrator.mnls_per_ts[i][j] ].hbfunction(Integrator.mnls_per_ts[i][j], &hf); } } if(Integrator.monitor_forces) monitor_forces(&hf); /* initialize the momenta */ enep = random_su3adj_field(reproduce_randomnumber_flag, hf.momenta); g_sloppy_precision = 1; /* run the trajectory */ if(Integrator.n_int[Integrator.no_timescales-1] > 0) { Integrator.integrate[Integrator.no_timescales-1](Integrator.tau, Integrator.no_timescales-1, 1); } g_sloppy_precision = 0; /* compute the final energy contributions for all monomials */ dh = 0.; for(i = 0; i < Integrator.no_timescales; i++) { for(j = 0; j < Integrator.no_mnls_per_ts[i]; j++) { dh += monomial_list[ Integrator.mnls_per_ts[i][j] ].accfunction(Integrator.mnls_per_ts[i][j], &hf); } } enepx = moment_energy(hf.momenta); if (!bc_flag) { /* if PBC */ new_plaquette_energy = measure_gauge_action( (const su3**) hf.gaugefield); if(g_rgi_C1 > 0. || g_rgi_C1 < 0.) { new_rectangle_energy = measure_rectangles( (const su3**) hf.gaugefield); } } if(g_proc_id == 0 && g_debug_level > 3) printf("called moment_energy: dh = %1.10e\n", (enepx - enep)); /* Compute the energy difference */ dh = dh + (enepx - enep); if(g_proc_id == 0 && g_debug_level > 3) { printf("called momenta_acc dH = %e\n", (enepx - enep)); } expmdh = exp(-dh); /* the random number is only taken at node zero and then distributed to the other sites */ ranlxd(yy,1); if(g_proc_id==0) { #ifdef MPI for(i = 1; i < g_nproc; i++) { MPI_Send(&yy[0], 1, MPI_DOUBLE, i, 31, MPI_COMM_WORLD); } #endif } #ifdef MPI else{ MPI_Recv(&yy[0], 1, MPI_DOUBLE, 0, 31, MPI_COMM_WORLD, &status); } #endif accept = (!acctest | (expmdh > yy[0])); if(g_proc_id == 0) { fprintf(stdout, "# Trajectory is %saccepted.\n", (accept ? "" : "not ")); } /* Here a reversibility test is performed */ /* The trajectory is integrated back */ if(return_check) { if(g_proc_id == 0) { fprintf(stdout, "# Performing reversibility check.\n"); } if(accept) { /* save gauge file to disk before performing reversibility check */ xlfInfo = construct_paramsXlfInfo((*plaquette_energy)/(6.*VOLUME*g_nproc), -1); // Should write this to temporary file first, and then check if(g_proc_id == 0 && g_debug_level > 0) { fprintf(stdout, "# Writing gauge field to file %s.\n", tmp_filename); } if((iostatus = write_gauge_field( tmp_filename, 64, xlfInfo) != 0 )) { /* Writing failed directly */ fprintf(stderr, "Error %d while writing gauge field to %s\nAborting...\n", iostatus, tmp_filename); exit(-2); } /* There is double writing of the gauge field, also in hmc_tm.c in this case */ /* No reading back check needed here, as reading back is done further down */ if(g_proc_id == 0 && g_debug_level > 0) { fprintf(stdout, "# Writing done.\n"); } free(xlfInfo); } g_sloppy_precision = 1; /* run the trajectory back */ Integrator.integrate[Integrator.no_timescales-1](-Integrator.tau, Integrator.no_timescales-1, 1); g_sloppy_precision = 0; /* compute the energy contributions from the pseudo-fermions */ ret_dh = 0.; for(i = 0; i < Integrator.no_timescales; i++) { for(j = 0; j < Integrator.no_mnls_per_ts[i]; j++) { ret_dh += monomial_list[ Integrator.mnls_per_ts[i][j] ].accfunction(Integrator.mnls_per_ts[i][j], &hf); } } ret_enep = moment_energy(hf.momenta); /* Compute the energy difference */ ret_dh += ret_enep - enep ; /* Compute Differences in the fields */ ks = 0.; kc = 0.; #ifdef OMP #pragma omp parallel private(w,v,tt,tr,ts,ds,ks,kc) { int thread_num = omp_get_thread_num(); #endif su3 ALIGN v0; #ifdef OMP #pragma omp for #endif for(int ix = 0; ix < VOLUME; ++ix) { for(int mu = 0; mu < 4; ++mu) { v=&hf.gaugefield[ix][mu]; w=&gauge_tmp[ix][mu]; _su3_minus_su3(v0, *v, *w); _su3_square_norm(ds, v0); tr = sqrt(ds) + kc; ts = tr + ks; tt = ts-ks; ks = ts; kc = tr-tt; } } kc=ks+kc; #ifdef OMP g_omp_acc_re[thread_num] = kc; } /* OpenMP parallel section closing brace */ /* sum up contributions from thread-local kahan summations */ for(int k = 0; k < omp_num_threads; ++k) ret_gauge_diff += g_omp_acc_re[k]; #else ret_gauge_diff = kc; #endif #ifdef MPI tmp = ret_gauge_diff; MPI_Reduce(&tmp, &ret_gauge_diff, 1, MPI_DOUBLE, MPI_SUM, 0, MPI_COMM_WORLD); #endif /* compute the total H */ tmp = enep; for(i = 0; i < Integrator.no_timescales; i++) { for(j = 0; j < Integrator.no_mnls_per_ts[i]; j++) { tmp += monomial_list[ Integrator.mnls_per_ts[i][j] ].energy0; } } /* Output */ if(g_proc_id == 0) { ret_check_file = fopen("return_check.data","a"); fprintf(ret_check_file,"ddh = %1.4e ddU= %1.4e ddh/H = %1.4e\n", ret_dh, ret_gauge_diff/4./((double)(VOLUME*g_nproc))/3., ret_dh/tmp); fclose(ret_check_file); } if(accept) { /* Read back gauge field */ if(g_proc_id == 0 && g_debug_level > 0) { fprintf(stdout, "# Trying to read gauge field from file %s.\n", tmp_filename); } if((iostatus = read_gauge_field(tmp_filename) != 0)) { fprintf(stderr, "Error %d while reading gauge field from %s\nAborting...\n", iostatus, tmp_filename); exit(-2); } if(g_proc_id == 0 && g_debug_level > 0) { fprintf(stdout, "# Reading done.\n"); } } if(g_proc_id == 0) { fprintf(stdout, "# Reversibility check done.\n"); } } /* end of reversibility check */ if(accept) { *plaquette_energy = new_plaquette_energy; *rectangle_energy = new_rectangle_energy; /* put the links back to SU(3) group */ if (!bc_flag) { /* periodic boundary conditions */ #ifdef OMP #pragma omp parallel for private(v) #endif for(int ix=0;ix<VOLUME;ix++) { for(int mu=0;mu<4;mu++) { v=&hf.gaugefield[ix][mu]; restoresu3_in_place(v); } } } } else { /* reject: copy gauge_tmp to hf.gaugefield */ #ifdef OMP #pragma omp parallel for private(w) private(v) #endif for(int ix=0;ix<VOLUME;ix++) { for(int mu=0;mu<4;mu++){ v=&hf.gaugefield[ix][mu]; w=&gauge_tmp[ix][mu]; _su3_assign(*v,*w); } } } hf.update_gauge_copy = 1; g_update_gauge_copy = 1; hf.update_gauge_energy = 1; g_update_gauge_energy = 1; hf.update_rectangle_energy = 1; g_update_rectangle_energy = 1; #ifdef MPI xchange_gauge(hf.gaugefield); #endif etime=gettime(); /* printing data in the .data file */ if(g_proc_id==0) { datafile = fopen(filename, "a"); if (!bc_flag) { /* if Periodic Boundary Conditions */ fprintf(datafile, "%.8d %14.12f %14.12f %e ", traj_counter, (*plaquette_energy)/(6.*VOLUME*g_nproc), dh, expmdh); } for(i = 0; i < Integrator.no_timescales; i++) { for(j = 0; j < Integrator.no_mnls_per_ts[i]; j++) { if(monomial_list[ Integrator.mnls_per_ts[i][j] ].type != GAUGE && monomial_list[ Integrator.mnls_per_ts[i][j] ].type != SFGAUGE && monomial_list[ Integrator.mnls_per_ts[i][j] ].type != NDPOLY && monomial_list[ Integrator.mnls_per_ts[i][j] ].type != NDCLOVER && monomial_list[ Integrator.mnls_per_ts[i][j] ].type != CLOVERNDTRLOG && monomial_list[ Integrator.mnls_per_ts[i][j] ].type != CLOVERTRLOG ) { fprintf(datafile,"%d %d ", monomial_list[ Integrator.mnls_per_ts[i][j] ].iter0, monomial_list[ Integrator.mnls_per_ts[i][j] ].iter1); } } } fprintf(datafile, "%d %e", accept, etime-atime); if(g_rgi_C1 > 0. || g_rgi_C1 < 0) { fprintf(datafile, " %e", (*rectangle_energy)/(12*VOLUME*g_nproc)); } fprintf(datafile, "\n"); fflush(datafile); fclose(datafile); } return(accept); }