void mul_one_sub_mul_gamma5(spinor * const l, spinor * const k, spinor * const j){ #ifdef OMP #pragma omp parallel { #endif spinor *r, *s, *t; /************ loop over all lattice sites ************/ #ifdef OMP #pragma omp for #endif for(int ix = 0; ix < (VOLUME/2); ++ix) { r = k+ix; s = j+ix; t = l+ix; /* Subtract s and store the result in t */ /* multiply with gamma5 included by */ /* reversed order of s and r (2&3) */ _vector_sub(t->s0, r->s0, s->s0); _vector_sub(t->s1, r->s1, s->s1); _vector_sub(t->s2, s->s2, r->s2); _vector_sub(t->s3, s->s3, r->s3); } #ifdef OMP } /* OpenMP closing brace */ #endif }
void mul_one_pm_imu_sub_mul_gamma5(spinor * const l, spinor * const k, spinor * const j, const double _sign){ #ifdef OMP #pragma omp parallel { #endif _Complex double z,w; int ix; double sign=1.; spinor *r, *s, *t; su3_vector ALIGN phi1, phi2, phi3, phi4; if(_sign < 0.){ sign = -1.; } z = 1. + (sign * g_mu) * I; w = conj(z); /************ loop over all lattice sites ************/ #ifdef OMP #pragma omp for #endif for(ix = 0; ix < (VOLUME/2); ix++){ r = k+ix; s = j+ix; t = l+ix; /* Multiply the spinorfield with 1+imu\gamma_5 */ _complex_times_vector(phi1, z, r->s0); _complex_times_vector(phi2, z, r->s1); _complex_times_vector(phi3, w, r->s2); _complex_times_vector(phi4, w, r->s3); /* Subtract s and store the result in t */ /* multiply with gamma5 included by */ /* reversed order of s and phi3|4 */ _vector_sub(t->s0, phi1, s->s0); _vector_sub(t->s1, phi2, s->s1); _vector_sub(t->s2, s->s2, phi3); _vector_sub(t->s3, s->s3, phi4); } #ifdef OMP } /* OpenMP closing brace */ #endif }
void deriv_Sb_D_psi(spinor * const l, spinor * const k, hamiltonian_field_t * const hf, const double factor) { #ifdef BGL __alignx(16, l); __alignx(16, k); #endif /* for parallelization */ #ifdef MPI xchange_lexicfield(k); xchange_lexicfield(l); #endif #ifdef OMP #define static #pragma omp parallel { #endif int ix,iy; su3 * restrict up ALIGN; su3 * restrict um ALIGN; static su3 v1,v2; static su3_vector psia,psib,phia,phib; static spinor rr; /* spinor * restrict r ALIGN; */ spinor * restrict sp ALIGN; spinor * restrict sm ALIGN; #ifdef OMP #undef static #endif #ifdef _KOJAK_INST #pragma pomp inst begin(derivSb) #endif #ifdef XLC #pragma disjoint(*sp, *sm, *up, *um) #endif /************** loop over all lattice sites ****************/ #ifdef OMP #pragma omp for #endif for(ix = 0; ix < (VOLUME); ix++){ rr = (*(l + ix)); /* rr=g_spinor_field[l][icx-ioff]; */ /*multiply the left vector with gamma5*/ _vector_minus_assign(rr.s2, rr.s2); _vector_minus_assign(rr.s3, rr.s3); /*********************** direction +0 ********************/ iy=g_iup[ix][0]; sp = k + iy; up=&hf->gaugefield[ix][0]; _vector_add(psia,(*sp).s0,(*sp).s2); _vector_add(psib,(*sp).s1,(*sp).s3); _vector_add(phia,rr.s0,rr.s2); _vector_add(phib,rr.s1,rr.s3); _vector_tensor_vector_add(v1, phia, psia, phib, psib); _su3_times_su3d(v2,*up,v1); _complex_times_su3(v1,ka0,v2); _trace_lambda_mul_add_assign_nonlocal(hf->derivative[ix][0], 2.*factor, v1); /************** direction -0 ****************************/ iy=g_idn[ix][0]; sm = k + iy; um=&hf->gaugefield[iy][0]; _vector_sub(psia,(*sm).s0,(*sm).s2); _vector_sub(psib,(*sm).s1,(*sm).s3); _vector_sub(phia,rr.s0,rr.s2); _vector_sub(phib,rr.s1,rr.s3); _vector_tensor_vector_add(v1, psia, phia, psib, phib); _su3_times_su3d(v2,*um,v1); _complex_times_su3(v1,ka0,v2); _trace_lambda_mul_add_assign_nonlocal(hf->derivative[iy][0], 2.*factor, v1); /*************** direction +1 **************************/ iy=g_iup[ix][1]; sp = k + iy; up=&hf->gaugefield[ix][1]; _vector_i_add(psia,(*sp).s0,(*sp).s3); _vector_i_add(psib,(*sp).s1,(*sp).s2); _vector_i_add(phia,rr.s0,rr.s3); _vector_i_add(phib,rr.s1,rr.s2); _vector_tensor_vector_add(v1, phia, psia, phib, psib); _su3_times_su3d(v2,*up,v1); _complex_times_su3(v1,ka1,v2); _trace_lambda_mul_add_assign_nonlocal(hf->derivative[ix][1], 2.*factor, v1); /**************** direction -1 *************************/ iy=g_idn[ix][1]; sm = k + iy; um=&hf->gaugefield[iy][1]; _vector_i_sub(psia,(*sm).s0,(*sm).s3); _vector_i_sub(psib,(*sm).s1,(*sm).s2); _vector_i_sub(phia,rr.s0,rr.s3); _vector_i_sub(phib,rr.s1,rr.s2); _vector_tensor_vector_add(v1, psia, phia, psib, phib); _su3_times_su3d(v2,*um,v1); _complex_times_su3(v1,ka1,v2); _trace_lambda_mul_add_assign_nonlocal(hf->derivative[iy][1], 2.*factor, v1); /*************** direction +2 **************************/ iy=g_iup[ix][2]; sp = k + iy; up=&hf->gaugefield[ix][2]; _vector_add(psia,(*sp).s0,(*sp).s3); _vector_sub(psib,(*sp).s1,(*sp).s2); _vector_add(phia,rr.s0,rr.s3); _vector_sub(phib,rr.s1,rr.s2); _vector_tensor_vector_add(v1, phia, psia, phib, psib); _su3_times_su3d(v2,*up,v1); _complex_times_su3(v1,ka2,v2); _trace_lambda_mul_add_assign_nonlocal(hf->derivative[ix][2], 2.*factor, v1); /***************** direction -2 ************************/ iy=g_idn[ix][2]; sm = k + iy; um=&hf->gaugefield[iy][2]; _vector_sub(psia,(*sm).s0,(*sm).s3); _vector_add(psib,(*sm).s1,(*sm).s2); _vector_sub(phia,rr.s0,rr.s3); _vector_add(phib,rr.s1,rr.s2); _vector_tensor_vector_add(v1, psia, phia, psib, phib); _su3_times_su3d(v2,*um,v1); _complex_times_su3(v1,ka2,v2); _trace_lambda_mul_add_assign_nonlocal(hf->derivative[iy][2], 2.*factor, v1); /****************** direction +3 ***********************/ iy=g_iup[ix][3]; sp = k + iy; up=&hf->gaugefield[ix][3]; _vector_i_add(psia,(*sp).s0,(*sp).s2); _vector_i_sub(psib,(*sp).s1,(*sp).s3); _vector_i_add(phia,rr.s0,rr.s2); _vector_i_sub(phib,rr.s1,rr.s3); _vector_tensor_vector_add(v1, phia, psia, phib, psib); _su3_times_su3d(v2,*up,v1); _complex_times_su3(v1,ka3,v2); _trace_lambda_mul_add_assign_nonlocal(hf->derivative[ix][3], 2.*factor, v1); /***************** direction -3 ************************/ iy=g_idn[ix][3]; sm = k + iy; um=&hf->gaugefield[iy][3]; _vector_i_sub(psia,(*sm).s0,(*sm).s2); _vector_i_add(psib,(*sm).s1,(*sm).s3); _vector_i_sub(phia,rr.s0,rr.s2); _vector_i_add(phib,rr.s1,rr.s3); _vector_tensor_vector_add(v1, psia, phia, psib, phib); _su3_times_su3d(v2,*um,v1); _complex_times_su3(v1,ka3,v2); _trace_lambda_mul_add_assign_nonlocal(hf->derivative[iy][3], 2.*factor, v1); /****************** end of loop ************************/ } #ifdef _KOJAK_INST #pragma pomp inst end(derivSb) #endif #ifdef OMP } /*OpenMP closing brace */ #endif }
void mul_one_pm_imu_sub_mul(spinor * const l, spinor * const k, spinor * const j, const double _sign, const int N){ #ifdef OMP #pragma omp parallel { #endif _Complex double z,w; int ix; double sign=1.; spinor *r, *s, *t; #if (!defined SSE2 && !defined SSE3) su3_vector ALIGN phi1, phi2, phi3, phi4; #endif if(_sign < 0.){ sign = -1.; } z = 1. + (sign * g_mu) * I; w = conj(z); /************ loop over all lattice sites ************/ #ifdef OMP #pragma omp for #endif for(ix = 0; ix < N; ix++){ r = k+ix; s = j+ix; t = l+ix; /* Multiply the spinorfield with 1+imu\gamma_5 */ #if (defined SSE2 || defined SSE3) _prefetch_spinor((r+predist)); _prefetch_spinor((s+predist)); _sse_load_up(r->s0); _sse_vector_cmplx_mul(z); _sse_load(s->s0); _sse_vector_sub_up(); _sse_store_nt_up(t->s0); _sse_load_up(r->s1); _sse_vector_cmplx_mul_two(); _sse_load(s->s1); _sse_vector_sub_up(); _sse_store_nt_up(t->s1); _sse_load_up(r->s2); _sse_vector_cmplx_mul(w); _sse_load(s->s2); _sse_vector_sub_up(); _sse_store_nt_up(t->s2); _sse_load_up(r->s3); _sse_vector_cmplx_mul_two(); _sse_load(s->s3); _sse_vector_sub_up(); _sse_store_nt_up(t->s3); #else _complex_times_vector(phi1, z, r->s0); _complex_times_vector(phi2, z, r->s1); _complex_times_vector(phi3, w, r->s2); _complex_times_vector(phi4, w, r->s3); /* Subtract s and store the result in t */ _vector_sub(t->s0, phi1, s->s0); _vector_sub(t->s1, phi2, s->s1); _vector_sub(t->s2, phi3, s->s2); _vector_sub(t->s3, phi4, s->s3); #endif } #ifdef OMP } /* OpenMP closing brace */ #endif }
/* for ieo=0, k resides on odd sites and l on even sites */ void Hopping_Matrix(int ieo, spinor * const l, spinor * const k){ int ix,iy; int ioff,ioff2,icx,icy; su3 * restrict up, * restrict um; spinor * restrict r, * restrict sp, * restrict sm; spinor temp; #ifdef _GAUGE_COPY if(g_update_gauge_copy) { update_backward_gauge(); } #endif /* for parallelization */ # if (defined MPI && !(defined _NO_COMM)) xchange_field(k, ieo); # endif if(k == l){ printf("Error in H_psi (simple.c):\n"); printf("Arguments k and l must be different\n"); printf("Program aborted\n"); exit(1); } if(ieo == 0){ ioff = 0; } else{ ioff = (VOLUME+RAND)/2; } ioff2 = (VOLUME+RAND)/2-ioff; /**************** loop over all lattice sites ****************/ for (icx = ioff; icx < (VOLUME/2 + ioff); icx++){ ix=g_eo2lexic[icx]; r=l+(icx-ioff); /*********************** direction +0 ************************/ iy=g_iup[ix][0]; icy=g_lexic2eosub[iy]; sp=k+icy; # if ((defined _GAUGE_COPY)) up=&g_gauge_field_copy[icx][0]; # else up=&g_gauge_field[ix][0]; # endif _vector_add(psi,(*sp).s0,(*sp).s2); _su3_multiply(chi,(*up),psi); _complex_times_vector(psi,ka0,chi); _vector_assign(temp.s0,psi); _vector_assign(temp.s2,psi); _vector_add(psi,(*sp).s1,(*sp).s3); _su3_multiply(chi,(*up),psi); _complex_times_vector(psi,ka0,chi); _vector_assign(temp.s1,psi); _vector_assign(temp.s3,psi); /*********************** direction -0 ************************/ iy=g_idn[ix][0]; icy=g_lexic2eosub[iy]; sm=k+icy; # if ((defined _GAUGE_COPY)) um = up+1; # else um=&g_gauge_field[iy][0]; # endif _vector_sub(psi,(*sm).s0,(*sm).s2); _su3_inverse_multiply(chi,(*um),psi); _complexcjg_times_vector(psi,ka0,chi); _vector_add_assign(temp.s0,psi); _vector_sub_assign(temp.s2,psi); _vector_sub(psi,(*sm).s1,(*sm).s3); _su3_inverse_multiply(chi,(*um),psi); _complexcjg_times_vector(psi,ka0,chi); _vector_add_assign(temp.s1,psi); _vector_sub_assign(temp.s3,psi); /*********************** direction +1 ************************/ iy=g_iup[ix][1]; icy=g_lexic2eosub[iy]; sp=k+icy; # if ((defined _GAUGE_COPY)) up=um+1; # else up+=1; # endif _vector_i_add(psi,(*sp).s0,(*sp).s3); _su3_multiply(chi,(*up),psi); _complex_times_vector(psi,ka1,chi); _vector_add_assign(temp.s0,psi); _vector_i_sub_assign(temp.s3,psi); _vector_i_add(psi,(*sp).s1,(*sp).s2); _su3_multiply(chi,(*up),psi); _complex_times_vector(psi,ka1,chi); _vector_add_assign(temp.s1,psi); _vector_i_sub_assign(temp.s2,psi); /*********************** direction -1 ************************/ iy=g_idn[ix][1]; icy=g_lexic2eosub[iy]; sm=k+icy; # ifndef _GAUGE_COPY um=&g_gauge_field[iy][1]; # else um=up+1; # endif _vector_i_sub(psi,(*sm).s0,(*sm).s3); _su3_inverse_multiply(chi,(*um),psi); _complexcjg_times_vector(psi,ka1,chi); _vector_add_assign(temp.s0,psi); _vector_i_add_assign(temp.s3,psi); _vector_i_sub(psi,(*sm).s1,(*sm).s2); _su3_inverse_multiply(chi,(*um),psi); _complexcjg_times_vector(psi,ka1,chi); _vector_add_assign(temp.s1,psi); _vector_i_add_assign(temp.s2,psi); /*********************** direction +2 ************************/ iy=g_iup[ix][2]; icy=g_lexic2eosub[iy]; sp=k+icy; # if ((defined _GAUGE_COPY)) up=um+1; # else up+=1; # endif _vector_add(psi,(*sp).s0,(*sp).s3); _su3_multiply(chi,(*up),psi); _complex_times_vector(psi,ka2,chi); _vector_add_assign(temp.s0,psi); _vector_add_assign(temp.s3,psi); _vector_sub(psi,(*sp).s1,(*sp).s2); _su3_multiply(chi,(*up),psi); _complex_times_vector(psi,ka2,chi); _vector_add_assign(temp.s1,psi); _vector_sub_assign(temp.s2,psi); /*********************** direction -2 ************************/ iy=g_idn[ix][2]; icy=g_lexic2eosub[iy]; sm=k+icy; # ifndef _GAUGE_COPY um = &g_gauge_field[iy][2]; # else um = up +1; # endif _vector_sub(psi,(*sm).s0,(*sm).s3); _su3_inverse_multiply(chi,(*um),psi); _complexcjg_times_vector(psi,ka2,chi); _vector_add_assign(temp.s0,psi); _vector_sub_assign(temp.s3,psi); _vector_add(psi,(*sm).s1,(*sm).s2); _su3_inverse_multiply(chi,(*um),psi); _complexcjg_times_vector(psi,ka2,chi); _vector_add_assign(temp.s1,psi); _vector_add_assign(temp.s2,psi); /*********************** direction +3 ************************/ iy=g_iup[ix][3]; icy=g_lexic2eosub[iy]; sp=k+icy; # if ((defined _GAUGE_COPY)) up=um+1; # else up+=1; # endif _vector_i_add(psi,(*sp).s0,(*sp).s2); _su3_multiply(chi,(*up),psi); _complex_times_vector(psi,ka3,chi); _vector_add_assign(temp.s0,psi); _vector_i_sub_assign(temp.s2,psi); _vector_i_sub(psi,(*sp).s1,(*sp).s3); _su3_multiply(chi,(*up),psi); _complex_times_vector(psi,ka3,chi); _vector_add_assign(temp.s1,psi); _vector_i_add_assign(temp.s3,psi); /*********************** direction -3 ************************/ iy=g_idn[ix][3]; icy=g_lexic2eosub[iy]; sm=k+icy; # ifndef _GAUGE_COPY um = &g_gauge_field[iy][3]; # else um = up+1; # endif _vector_i_sub(psi,(*sm).s0,(*sm).s2); _su3_inverse_multiply(chi,(*um),psi); _complexcjg_times_vector(psi,ka3,chi); _vector_add((*r).s0, temp.s0, psi); _vector_i_add((*r).s2, temp.s2, psi); _vector_i_add(psi,(*sm).s1,(*sm).s3); _su3_inverse_multiply(chi,(*um),psi); _complexcjg_times_vector(psi,ka3,chi); _vector_add((*r).s1, temp.s1, psi); _vector_i_sub((*r).s3, temp.s3, psi); /************************ end of loop ************************/ } }
void deriv_Sb_D_psi(spinor * const l, spinor * const k) { /* const int l, const int k){ */ int ix,iy; su3 * restrict up ALIGN; su3 * restrict um ALIGN; /* su3adj * restrict ddd; */ /* static su3adj der; */ static su3 v1,v2; static su3_vector psia,psib,phia,phib; static spinor rr; /* spinor * restrict r ALIGN; */ spinor * restrict sp ALIGN; spinor * restrict sm ALIGN; #ifdef _KOJAK_INST #pragma pomp inst begin(derivSb) #endif #ifdef XLC #pragma disjoint(*sp, *sm, *up, *um) #endif #ifdef BGL __alignx(16, l); __alignx(16, k); #endif /* for parallelization */ #ifdef MPI xchange_lexicfield(k); xchange_lexicfield(l); #endif /************** loop over all lattice sites ****************/ for(ix = 0; ix < (VOLUME); ix++){ rr = (*(l + ix)); /* rr=g_spinor_field[l][icx-ioff]; */ /*multiply the left vector with gamma5*/ _vector_minus_assign(rr.s2, rr.s2); _vector_minus_assign(rr.s3, rr.s3); /*********************** direction +0 ********************/ iy=g_iup[ix][0]; sp = k + iy; /* sp=&g_spinor_field[k][icy]; */ up=&g_gauge_field[ix][0]; _vector_add(psia,(*sp).s0,(*sp).s2); _vector_add(psib,(*sp).s1,(*sp).s3); _vector_add(phia,rr.s0,rr.s2); _vector_add(phib,rr.s1,rr.s3); _vector_tensor_vector_add(v1, phia, psia, phib, psib); /* _vector_tensor_vector(v1,phia,psia); */ /* _vector_tensor_vector(v2,phib,psib); */ /* _su3_plus_su3(v1,v1,v2); */ _su3_times_su3d(v2,*up,v1); _complex_times_su3(v1,ka0,v2); _trace_lambda_add_assign(df0[ix][0], v1); /* _trace_lambda(der,v1); */ /* ddd=&df0[ix][0]; */ /* _add_su3adj(*ddd,der); */ /************** direction -0 ****************************/ iy=g_idn[ix][0]; sm = k + iy; /* sm=&g_spinor_field[k][icy]; */ um=&g_gauge_field[iy][0]; _vector_sub(psia,(*sm).s0,(*sm).s2); _vector_sub(psib,(*sm).s1,(*sm).s3); _vector_sub(phia,rr.s0,rr.s2); _vector_sub(phib,rr.s1,rr.s3); _vector_tensor_vector_add(v1, psia, phia, psib, phib); /* _vector_tensor_vector(v1,psia,phia); */ /* _vector_tensor_vector(v2,psib,phib); */ /* _su3_plus_su3(v1,v1,v2); */ _su3_times_su3d(v2,*um,v1); _complex_times_su3(v1,ka0,v2); _trace_lambda_add_assign(df0[iy][0], v1); /* _trace_lambda(der,v1); */ /* ddd=&df0[iy][0]; */ /* _add_su3adj(*ddd,der); */ /*************** direction +1 **************************/ iy=g_iup[ix][1]; sp = k + iy; /* sp=&g_spinor_field[k][icy]; */ up=&g_gauge_field[ix][1]; _vector_i_add(psia,(*sp).s0,(*sp).s3); _vector_i_add(psib,(*sp).s1,(*sp).s2); _vector_i_add(phia,rr.s0,rr.s3); _vector_i_add(phib,rr.s1,rr.s2); _vector_tensor_vector_add(v1, phia, psia, phib, psib); /* _vector_tensor_vector(v1,phia,psia); */ /* _vector_tensor_vector(v2,phib,psib); */ /* _su3_plus_su3(v1,v1,v2); */ _su3_times_su3d(v2,*up,v1); _complex_times_su3(v1,ka1,v2); _trace_lambda_add_assign(df0[ix][1], v1); /* _trace_lambda(der,v1); */ /* ddd=&df0[ix][1]; */ /* _add_su3adj(*ddd,der); */ /**************** direction -1 *************************/ iy=g_idn[ix][1]; sm = k + iy; /* sm=&g_spinor_field[k][icy]; */ um=&g_gauge_field[iy][1]; _vector_i_sub(psia,(*sm).s0,(*sm).s3); _vector_i_sub(psib,(*sm).s1,(*sm).s2); _vector_i_sub(phia,rr.s0,rr.s3); _vector_i_sub(phib,rr.s1,rr.s2); _vector_tensor_vector_add(v1, psia, phia, psib, phib); /* _vector_tensor_vector(v1,psia,phia); */ /* _vector_tensor_vector(v2,psib,phib); */ /* _su3_plus_su3(v1,v1,v2); */ _su3_times_su3d(v2,*um,v1); _complex_times_su3(v1,ka1,v2); _trace_lambda_add_assign(df0[iy][1], v1); /* _trace_lambda(der,v1); */ /* ddd=&df0[iy][1]; */ /* _add_su3adj(*ddd,der); */ /*************** direction +2 **************************/ iy=g_iup[ix][2]; sp = k + iy; /* sp=&g_spinor_field[k][icy]; */ up=&g_gauge_field[ix][2]; _vector_add(psia,(*sp).s0,(*sp).s3); _vector_sub(psib,(*sp).s1,(*sp).s2); _vector_add(phia,rr.s0,rr.s3); _vector_sub(phib,rr.s1,rr.s2); _vector_tensor_vector_add(v1, phia, psia, phib, psib); /* _vector_tensor_vector(v1,phia,psia); */ /* _vector_tensor_vector(v2,phib,psib); */ /* _su3_plus_su3(v1,v1,v2); */ _su3_times_su3d(v2,*up,v1); _complex_times_su3(v1,ka2,v2); _trace_lambda_add_assign(df0[ix][2], v1); /* _trace_lambda(der,v1); */ /* ddd=&df0[ix][2]; */ /* _add_su3adj(*ddd,der); */ /***************** direction -2 ************************/ iy=g_idn[ix][2]; sm = k + iy; /* sm=&g_spinor_field[k][icy]; */ um=&g_gauge_field[iy][2]; _vector_sub(psia,(*sm).s0,(*sm).s3); _vector_add(psib,(*sm).s1,(*sm).s2); _vector_sub(phia,rr.s0,rr.s3); _vector_add(phib,rr.s1,rr.s2); _vector_tensor_vector_add(v1, psia, phia, psib, phib); /* _vector_tensor_vector(v1,psia,phia); */ /* _vector_tensor_vector(v2,psib,phib); */ /* _su3_plus_su3(v1,v1,v2); */ _su3_times_su3d(v2,*um,v1); _complex_times_su3(v1,ka2,v2); _trace_lambda_add_assign(df0[iy][2], v1); /* _trace_lambda(der,v1); */ /* ddd=&df0[iy][2]; */ /* _add_su3adj(*ddd,der); */ /****************** direction +3 ***********************/ iy=g_iup[ix][3]; sp = k + iy; /* sp=&g_spinor_field[k][icy]; */ up=&g_gauge_field[ix][3]; _vector_i_add(psia,(*sp).s0,(*sp).s2); _vector_i_sub(psib,(*sp).s1,(*sp).s3); _vector_i_add(phia,rr.s0,rr.s2); _vector_i_sub(phib,rr.s1,rr.s3); _vector_tensor_vector_add(v1, phia, psia, phib, psib); /* _vector_tensor_vector(v1,phia,psia); */ /* _vector_tensor_vector(v2,phib,psib); */ /* _su3_plus_su3(v1,v1,v2); */ _su3_times_su3d(v2,*up,v1); _complex_times_su3(v1,ka3,v2); _trace_lambda_add_assign(df0[ix][3], v1); /* _trace_lambda(der,v1); */ /* ddd=&df0[ix][3]; */ /* _add_su3adj(*ddd,der); */ /***************** direction -3 ************************/ iy=g_idn[ix][3]; sm = k + iy; /* sm=&g_spinor_field[k][icy]; */ um=&g_gauge_field[iy][3]; _vector_i_sub(psia,(*sm).s0,(*sm).s2); _vector_i_add(psib,(*sm).s1,(*sm).s3); _vector_i_sub(phia,rr.s0,rr.s2); _vector_i_add(phib,rr.s1,rr.s3); _vector_tensor_vector_add(v1, psia, phia, psib, phib); /* _vector_tensor_vector(v1,psia,phia); */ /* _vector_tensor_vector(v2,psib,phib); */ /* _su3_plus_su3(v1,v1,v2); */ _su3_times_su3d(v2,*um,v1); _complex_times_su3(v1,ka3,v2); _trace_lambda_add_assign(df0[iy][3], v1); /* _trace_lambda(der,v1); */ /* ddd=&df0[iy][3]; */ /* _add_su3adj(*ddd,der); */ /****************** end of loop ************************/ } #ifdef _KOJAK_INST #pragma pomp inst end(derivSb) #endif }
/* for ieo=0, k resides on odd sites and l on even sites */ void Hopping_Matrix(const int ieo, spinor * const l, spinor * const k){ int i,ix; su3 * restrict U ALIGN; spinor * restrict s ALIGN; spinor rs; static su3_vector psi, chi, psi2, chi2; halfspinor * restrict * phi ALIGN; halfspinor32 * restrict * phi32 ALIGN; #ifdef _KOJAK_INST #pragma pomp inst begin(hoppingmatrix) #endif #ifdef XLC #pragma disjoint(*l, *k, *U, *s) #endif #ifdef _GAUGE_COPY if(g_update_gauge_copy) { update_backward_gauge(); } #endif if(k == l){ printf("Error in H_psi (simple.c):\n"); printf("Arguments k and l must be different\n"); printf("Program aborted\n"); exit(1); } s = k; if(ieo == 0) { U = g_gauge_field_copy[0][0]; } else { U = g_gauge_field_copy[1][0]; } if(g_sloppy_precision == 1 && g_sloppy_precision_flag == 1) { phi32 = NBPointer32[ieo]; /**************** loop over all lattice sites ****************/ ix=0; for(i = 0; i < (VOLUME)/2; i++){ _vector_assign(rs.s0, (*s).s0); _vector_assign(rs.s1, (*s).s1); _vector_assign(rs.s2, (*s).s2); _vector_assign(rs.s3, (*s).s3); s++; /*********************** direction +0 ************************/ _vector_add(psi, rs.s0, rs.s2); _su3_multiply(chi,(*U),psi); _complex_times_vector((*phi32[ix]).s0, ka0, chi); _vector_add(psi, rs.s1, rs.s3); _su3_multiply(chi,(*U),psi); _complex_times_vector((*phi32[ix]).s1, ka0, chi); U++; ix++; /*********************** direction -0 ************************/ _vector_sub((*phi32[ix]).s0, rs.s0, rs.s2); _vector_sub((*phi32[ix]).s1, rs.s1, rs.s3); ix++; /*********************** direction +1 ************************/ _vector_i_add(psi, rs.s0, rs.s3); _su3_multiply(chi, (*U), psi); _complex_times_vector((*phi32[ix]).s0, ka1, chi); _vector_i_add(psi, rs.s1, rs.s2); _su3_multiply(chi, (*U), psi); _complex_times_vector((*phi32[ix]).s1, ka1, chi); U++; ix++; /*********************** direction -1 ************************/ _vector_i_sub((*phi32[ix]).s0, rs.s0, rs.s3); _vector_i_sub((*phi32[ix]).s1, rs.s1, rs.s2); ix++; /*********************** direction +2 ************************/ _vector_add(psi, rs.s0, rs.s3); _su3_multiply(chi,(*U),psi); _complex_times_vector((*phi32[ix]).s0, ka2, chi); _vector_sub(psi, rs.s1, rs.s2); _su3_multiply(chi,(*U),psi); _complex_times_vector((*phi32[ix]).s1, ka2, chi); U++; ix++; /*********************** direction -2 ************************/ _vector_sub((*phi32[ix]).s0, rs.s0, rs.s3); _vector_add((*phi32[ix]).s1, rs.s1, rs.s2); ix++; /*********************** direction +3 ************************/ _vector_i_add(psi, rs.s0, rs.s2); _su3_multiply(chi, (*U), psi); _complex_times_vector((*phi32[ix]).s0, ka3, chi); _vector_i_sub(psi, rs.s1, rs.s3); _su3_multiply(chi,(*U),psi); _complex_times_vector((*phi32[ix]).s1, ka3, chi); U++; ix++; /*********************** direction -3 ************************/ _vector_i_sub((*phi32[ix]).s0, rs.s0, rs.s2); _vector_i_add((*phi32[ix]).s1, rs.s1, rs.s3); ix++; /************************ end of loop ************************/ } # if (defined MPI && !defined _NO_COMM) xchange_halffield32(); # endif s = l; phi32 = NBPointer32[2 + ieo]; if(ieo == 0) { U = g_gauge_field_copy[1][0]; } else { U = g_gauge_field_copy[0][0]; } ix = 0; for(i = 0; i < (VOLUME)/2; i++){ /*********************** direction +0 ************************/ _vector_assign(rs.s0, (*phi32[ix]).s0); _vector_assign(rs.s2, (*phi32[ix]).s0); _vector_assign(rs.s1, (*phi32[ix]).s1); _vector_assign(rs.s3, (*phi32[ix]).s1); ix++; /*********************** direction -0 ************************/ _vector_assign(psi, (*phi32[ix]).s0); _su3_inverse_multiply(chi,(*U), psi); _complexcjg_times_vector(psi,ka0,chi); _vector_add_assign(rs.s0, psi); _vector_sub_assign(rs.s2, psi); _vector_assign(psi, (*phi32[ix]).s1); _su3_inverse_multiply(chi,(*U), psi); _complexcjg_times_vector(psi,ka0,chi); _vector_add_assign(rs.s1, psi); _vector_sub_assign(rs.s3, psi); ix++; U++; /*********************** direction +1 ************************/ _vector_add_assign(rs.s0, (*phi32[ix]).s0); _vector_i_sub_assign(rs.s3, (*phi32[ix]).s0); _vector_add_assign(rs.s1, (*phi32[ix]).s1); _vector_i_sub_assign(rs.s2, (*phi32[ix]).s1); ix++; /*********************** direction -1 ************************/ _vector_assign(psi, (*phi32[ix]).s0); _su3_inverse_multiply(chi,(*U), psi); _complexcjg_times_vector(psi,ka1,chi); _vector_add_assign(rs.s0, psi); _vector_i_add_assign(rs.s3, psi); _vector_assign(psi, (*phi32[ix]).s1); _su3_inverse_multiply(chi,(*U), psi); _complexcjg_times_vector(psi,ka1,chi); _vector_add_assign(rs.s1, psi); _vector_i_add_assign(rs.s2, psi); U++; ix++; /*********************** direction +2 ************************/ _vector_add_assign(rs.s0, (*phi32[ix]).s0); _vector_add_assign(rs.s3, (*phi32[ix]).s0); _vector_add_assign(rs.s1, (*phi32[ix]).s1); _vector_sub_assign(rs.s2, (*phi32[ix]).s1); ix++; /*********************** direction -2 ************************/ _vector_assign(psi, (*phi32[ix]).s0); _su3_inverse_multiply(chi,(*U), psi); _complexcjg_times_vector(psi,ka2,chi); _vector_add_assign(rs.s0, psi); _vector_sub_assign(rs.s3, psi); _vector_assign(psi, (*phi32[ix]).s1); _su3_inverse_multiply(chi, (*U), psi); _complexcjg_times_vector(psi,ka2,chi); _vector_add_assign(rs.s1, psi); _vector_add_assign(rs.s2, psi); U++; ix++; /*********************** direction +3 ************************/ _vector_add_assign(rs.s0, (*phi32[ix]).s0); _vector_i_sub_assign(rs.s2, (*phi32[ix]).s0); _vector_add_assign(rs.s1, (*phi32[ix]).s1); _vector_i_add_assign(rs.s3, (*phi32[ix]).s1); ix++; /*********************** direction -3 ************************/ _vector_assign(psi, (*phi32[ix]).s0); _su3_inverse_multiply(chi,(*U), psi); _complexcjg_times_vector(psi,ka3,chi); _vector_add((*s).s0, rs.s0, psi); _vector_i_add((*s).s2, rs.s2, psi); _vector_assign(psi, (*phi32[ix]).s1); _su3_inverse_multiply(chi,(*U), psi); _complexcjg_times_vector(psi,ka3,chi); _vector_add((*s).s1, rs.s1, psi); _vector_i_sub((*s).s3, rs.s3, psi); U++; ix++; s++; } } else { phi = NBPointer[ieo]; /**************** loop over all lattice sites ****************/ ix=0; /* #pragma ivdep*/ for(i = 0; i < (VOLUME)/2; i++){ _vector_assign(rs.s0, (*s).s0); _vector_assign(rs.s1, (*s).s1); _vector_assign(rs.s2, (*s).s2); _vector_assign(rs.s3, (*s).s3); s++; /*********************** direction +0 ************************/ _vector_add(psi, rs.s0, rs.s2); _vector_add(psi2, rs.s1, rs.s3); _su3_multiply(chi,(*U),psi); _su3_multiply(chi2,(*U),psi2); _complex_times_vector((*phi[ix]).s0, ka0, chi); _complex_times_vector((*phi[ix]).s1, ka0, chi2); U++; ix++; /*********************** direction -0 ************************/ _vector_sub((*phi[ix]).s0, rs.s0, rs.s2); _vector_sub((*phi[ix]).s1, rs.s1, rs.s3); ix++; /*********************** direction +1 ************************/ _vector_i_add(psi, rs.s0, rs.s3); _vector_i_add(psi2, rs.s1, rs.s2); _su3_multiply(chi, (*U), psi); _su3_multiply(chi2, (*U), psi2); _complex_times_vector((*phi[ix]).s0, ka1, chi); _complex_times_vector((*phi[ix]).s1, ka1, chi2); U++; ix++; /*********************** direction -1 ************************/ _vector_i_sub((*phi[ix]).s0, rs.s0, rs.s3); _vector_i_sub((*phi[ix]).s1, rs.s1, rs.s2); ix++; /*********************** direction +2 ************************/ _vector_add(psi, rs.s0, rs.s3); _vector_sub(psi2, rs.s1, rs.s2); _su3_multiply(chi,(*U),psi); _su3_multiply(chi2,(*U),psi2); _complex_times_vector((*phi[ix]).s0, ka2, chi); _complex_times_vector((*phi[ix]).s1, ka2, chi2); U++; ix++; /*********************** direction -2 ************************/ _vector_sub((*phi[ix]).s0, rs.s0, rs.s3); _vector_add((*phi[ix]).s1, rs.s1, rs.s2); ix++; /*********************** direction +3 ************************/ _vector_i_add(psi, rs.s0, rs.s2); _vector_i_sub(psi2, rs.s1, rs.s3); _su3_multiply(chi, (*U), psi); _su3_multiply(chi2,(*U),psi2); _complex_times_vector((*phi[ix]).s0, ka3, chi); _complex_times_vector((*phi[ix]).s1, ka3, chi2); U++; ix++; /*********************** direction -3 ************************/ _vector_i_sub((*phi[ix]).s0, rs.s0, rs.s2); _vector_i_add((*phi[ix]).s1, rs.s1, rs.s3); ix++; /************************ end of loop ************************/ } # if (defined MPI && !defined _NO_COMM) xchange_halffield(); # endif s = l; phi = NBPointer[2 + ieo]; if(ieo == 0) { U = g_gauge_field_copy[1][0]; } else { U = g_gauge_field_copy[0][0]; } ix = 0; /* #pragma ivdep */ for(i = 0; i < (VOLUME)/2; i++){ /*********************** direction +0 ************************/ _vector_assign(rs.s0, (*phi[ix]).s0); _vector_assign(rs.s2, (*phi[ix]).s0); _vector_assign(rs.s1, (*phi[ix]).s1); _vector_assign(rs.s3, (*phi[ix]).s1); ix++; /*********************** direction -0 ************************/ _su3_inverse_multiply(chi,(*U),(*phi[ix]).s0); _su3_inverse_multiply(chi2,(*U),(*phi[ix]).s1); _complexcjg_times_vector(psi,ka0,chi); _complexcjg_times_vector(psi2,ka0,chi2); _vector_add_assign(rs.s0, psi); _vector_sub_assign(rs.s2, psi); _vector_add_assign(rs.s1, psi2); _vector_sub_assign(rs.s3, psi2); ix++; U++; /*********************** direction +1 ************************/ _vector_add_assign(rs.s0, (*phi[ix]).s0); _vector_i_sub_assign(rs.s3, (*phi[ix]).s0); _vector_add_assign(rs.s1, (*phi[ix]).s1); _vector_i_sub_assign(rs.s2, (*phi[ix]).s1); ix++; /*********************** direction -1 ************************/ _su3_inverse_multiply(chi,(*U), (*phi[ix]).s0); _su3_inverse_multiply(chi2, (*U), (*phi[ix]).s1); _complexcjg_times_vector(psi,ka1,chi); _complexcjg_times_vector(psi2,ka1,chi2); _vector_add_assign(rs.s0, psi); _vector_i_add_assign(rs.s3, psi); _vector_add_assign(rs.s1, psi2); _vector_i_add_assign(rs.s2, psi2); U++; ix++; /*********************** direction +2 ************************/ _vector_add_assign(rs.s0, (*phi[ix]).s0); _vector_add_assign(rs.s3, (*phi[ix]).s0); _vector_add_assign(rs.s1, (*phi[ix]).s1); _vector_sub_assign(rs.s2, (*phi[ix]).s1); ix++; /*********************** direction -2 ************************/ _su3_inverse_multiply(chi,(*U), (*phi[ix]).s0); _su3_inverse_multiply(chi2, (*U), (*phi[ix]).s1); _complexcjg_times_vector(psi,ka2,chi); _complexcjg_times_vector(psi2,ka2,chi2); _vector_add_assign(rs.s0, psi); _vector_sub_assign(rs.s3, psi); _vector_add_assign(rs.s1, psi2); _vector_add_assign(rs.s2, psi2); U++; ix++; /*********************** direction +3 ************************/ _vector_add_assign(rs.s0, (*phi[ix]).s0); _vector_i_sub_assign(rs.s2, (*phi[ix]).s0); _vector_add_assign(rs.s1, (*phi[ix]).s1); _vector_i_add_assign(rs.s3, (*phi[ix]).s1); ix++; /*********************** direction -3 ************************/ _su3_inverse_multiply(chi,(*U), (*phi[ix]).s0); _su3_inverse_multiply(chi2, (*U), (*phi[ix]).s1); _complexcjg_times_vector(psi,ka3,chi); _complexcjg_times_vector(psi2,ka3,chi2); _vector_add((*s).s0, rs.s0, psi); _vector_i_add((*s).s2, rs.s2, psi); _vector_add((*s).s1, rs.s1, psi2); _vector_i_sub((*s).s3, rs.s3, psi2); U++; ix++; s++; } } #ifdef _KOJAK_INST #pragma pomp inst end(hoppingmatrix) #endif }