void mul_one_pm_imu_inv(spinor * const l, const double _sign, const int N){ #ifdef OMP #pragma omp parallel { #endif _Complex double ALIGN z,w; int ix; double sign=-1.; spinor *r; su3_vector ALIGN phi1; double ALIGN nrm = 1./(1.+g_mu*g_mu); if(_sign < 0.){ sign = 1.; } z = nrm + (sign * nrm * g_mu) * I; w = conj(z); /************ loop over all lattice sites ************/ #ifdef OMP #pragma omp for #endif for(ix = 0; ix < N; ix++){ r=l + ix; /* Multiply the spinorfield with the inverse of 1+imu\gamma_5 */ #if ( defined SSE2 || defined SSE3 ) _prefetch_spinor((r+predist)); _sse_load_up(r->s0); _sse_vector_cmplx_mul(z); _sse_store_nt_up(r->s0); _sse_load_up(r->s1); _sse_vector_cmplx_mul_two(); _sse_store_nt_up(r->s1); _sse_load_up(r->s2); _sse_vector_cmplx_mul(w); _sse_store_nt_up(r->s2); _sse_load_up(r->s3); _sse_vector_cmplx_mul_two(); _sse_store_nt_up(r->s3); #else _complex_times_vector(phi1, z, r->s0); _vector_assign(r->s0, phi1); _complex_times_vector(phi1, z, r->s1); _vector_assign(r->s1, phi1); _complex_times_vector(phi1, w, r->s2); _vector_assign(r->s2, phi1); _complex_times_vector(phi1, w, r->s3); _vector_assign(r->s3, phi1); #endif } #ifdef OMP } /* OpenMP closing brace */ #endif }
void assign_mul_one_pm_imu(spinor * const l, spinor * const k, const double _sign, const int N){ #ifdef OMP #pragma omp parallel { #endif _Complex double z,w; int ix; double sign = 1.; spinor *r, *s; 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++){ s=l+ix; r=k+ix; /* Multiply the spinorfield with of 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_store_nt_up(s->s0); _sse_load_up(r->s1); _sse_vector_cmplx_mul_two(); _sse_store_nt_up(s->s1); _sse_load_up(r->s2); _sse_vector_cmplx_mul(w); _sse_store_nt_up(s->s2); _sse_load_up(r->s3); _sse_vector_cmplx_mul_two(); _sse_store_nt_up(s->s3); #else _complex_times_vector(s->s0, z, r->s0); _complex_times_vector(s->s1, z, r->s1); _complex_times_vector(s->s2, w, r->s2); _complex_times_vector(s->s3, w, r->s3); #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 }
void Hopping_Matrix(const int ieo, spinor * const l, spinor * const k){ int ix, i; su3 * restrict U ALIGN; static spinor rs; spinor * restrict s ALIGN; halfspinor ** phi ALIGN; #if defined OPTERON const int predist=2; #else const int predist=1; #endif #ifdef _KOJAK_INST #pragma pomp inst begin(hoppingmatrix) #endif #ifdef _GAUGE_COPY if(g_update_gauge_copy) { update_backward_gauge(); } #endif /* We will run through the source vector now */ /* instead of the solution vector */ s = k; _prefetch_spinor(s); if(ieo == 0) { U = g_gauge_field_copy[0][0]; } else { U = g_gauge_field_copy[1][0]; } phi = NBPointer[ieo]; _prefetch_su3(U); /**************** loop over all lattice sites ******************/ ix=0; for(i = 0; i < (VOLUME)/2; i++){ /*********************** direction +0 ************************/ _prefetch_su3(U+predist); _sse_load((*s).s0); _sse_load_up((*s).s2); _sse_vector_add(); _sse_su3_multiply((*U)); _sse_vector_cmplx_mul(ka0); _sse_store_nt_up((*phi[ix]).s0); _sse_load((*s).s1); _sse_load_up((*s).s3); _sse_vector_add(); _sse_su3_multiply((*U)); _sse_vector_cmplx_mul(ka0); _sse_store_nt_up((*phi[ix]).s1); U++; ix++; /*********************** direction -0 ************************/ _sse_load((*s).s0); _sse_load_up((*s).s2); _sse_vector_sub(); _sse_store_nt((*phi[ix]).s0); _sse_load((*s).s1); _sse_load_up((*s).s3); _sse_vector_sub(); _sse_store_nt((*phi[ix]).s1); ix++; /*********************** direction +1 ************************/ _prefetch_su3(U+predist); _sse_load((*s).s0); /*next not needed?*/ _sse_load_up((*s).s3); _sse_vector_i_mul(); _sse_vector_add(); _sse_su3_multiply((*U)); _sse_vector_cmplx_mul(ka1); _sse_store_nt_up((*phi[ix]).s0); _sse_load((*s).s1); _sse_load_up((*s).s2); _sse_vector_i_mul(); _sse_vector_add(); _sse_su3_multiply((*U)); _sse_vector_cmplx_mul(ka1); _sse_store_nt_up((*phi[ix]).s1); ix++; U++; /*********************** direction -1 ************************/ _sse_load((*s).s0); _sse_load_up((*s).s3); _sse_vector_i_mul(); _sse_vector_sub(); _sse_store_nt((*phi[ix]).s0); _sse_load((*s).s1); _sse_load_up((*s).s2); _sse_vector_i_mul(); _sse_vector_sub(); _sse_store_nt((*phi[ix]).s1); ix++; /*********************** direction +2 ************************/ _prefetch_su3(U+predist); _sse_load((*s).s0); _sse_load_up((*s).s3); _sse_vector_add(); _sse_su3_multiply((*U)); _sse_vector_cmplx_mul(ka2); _sse_store_nt_up((*phi[ix]).s0); _sse_load((*s).s1); _sse_load_up((*s).s2); _sse_vector_sub(); _sse_su3_multiply((*U)); _sse_vector_cmplx_mul(ka2); _sse_store_nt_up((*phi[ix]).s1); ix++; U++; /*********************** direction -2 ************************/ _sse_load((*s).s0); _sse_load_up((*s).s3); _sse_vector_sub(); _sse_store_nt((*phi[ix]).s0); _sse_load((*s).s1); _sse_load_up((*s).s2); _sse_vector_add(); _sse_store_nt((*phi[ix]).s1); ix++; /*********************** direction +3 ************************/ _prefetch_su3(U+predist); _prefetch_spinor(s+1); _sse_load((*s).s0); _sse_load_up((*s).s2); _sse_vector_i_mul(); _sse_vector_add(); _sse_su3_multiply((*U)); _sse_vector_cmplx_mul(ka3); _sse_store_nt_up((*phi[ix]).s0); _sse_load((*s).s1); _sse_load_up((*s).s3); _sse_vector_i_mul(); _sse_vector_sub(); _sse_su3_multiply((*U)); _sse_vector_cmplx_mul(ka3); _sse_store_nt_up((*phi[ix]).s1); ix++; U++; /*********************** direction -3 ************************/ _sse_load((*s).s0); _sse_load_up((*s).s2); _sse_vector_i_mul(); _sse_vector_sub(); _sse_store_nt((*phi[ix]).s0); _sse_load((*s).s1); _sse_load_up((*s).s3); _sse_vector_i_mul(); _sse_vector_add(); _sse_store_nt((*phi[ix]).s1); ix++; s++; } # 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]; } _prefetch_su3(U); /* Now we sum up and expand to a full spinor */ ix = 0; 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 ************************/ _prefetch_su3(U+predist); _sse_load((*phi[ix]).s0); _sse_su3_inverse_multiply((*U)); _sse_vector_cmplxcg_mul(ka0); _sse_load(rs.s0); _sse_vector_add(); _sse_store(rs.s0); _sse_load(rs.s2); _sse_vector_sub(); _sse_store(rs.s2); _sse_load((*phi[ix]).s1); _sse_su3_inverse_multiply((*U)); _sse_vector_cmplxcg_mul(ka0); _sse_load(rs.s1); _sse_vector_add(); _sse_store(rs.s1); _sse_load(rs.s3); _sse_vector_sub(); _sse_store(rs.s3); ix++; U++; /*********************** direction +1 ************************/ _sse_load_up((*phi[ix]).s0); _sse_load(rs.s0); _sse_vector_add(); _sse_store(rs.s0); _sse_load(rs.s3); _sse_vector_i_mul(); _sse_vector_sub(); _sse_store(rs.s3); _sse_load_up((*phi[ix]).s1); _sse_load(rs.s1); _sse_vector_add(); _sse_store(rs.s1); _sse_load(rs.s2); _sse_vector_i_mul(); _sse_vector_sub(); _sse_store(rs.s2); ix++; /*********************** direction -1 ************************/ _prefetch_su3(U+predist); _sse_load((*phi[ix]).s0); _sse_su3_inverse_multiply((*U)); _sse_vector_cmplxcg_mul(ka1); _sse_load(rs.s0); _sse_vector_add(); _sse_store(rs.s0); _sse_load(rs.s3); _sse_vector_i_mul(); _sse_vector_add(); _sse_store(rs.s3); _sse_load((*phi[ix]).s1); _sse_su3_inverse_multiply((*U)); _sse_vector_cmplxcg_mul(ka1); _sse_load(rs.s1); _sse_vector_add(); _sse_store(rs.s1); _sse_load(rs.s2); _sse_vector_i_mul(); _sse_vector_add(); _sse_store(rs.s2); ix++; U++; /*********************** direction +2 ************************/ _sse_load_up((*phi[ix]).s0); _sse_load(rs.s0); _sse_vector_add(); _sse_store(rs.s0); _sse_load(rs.s3); _sse_vector_add(); _sse_store(rs.s3); _sse_load_up((*phi[ix]).s1); _sse_load(rs.s1); _sse_vector_add(); _sse_store(rs.s1); _sse_load(rs.s2); _sse_vector_sub(); _sse_store(rs.s2); ix++; /*********************** direction -2 ************************/ _prefetch_su3(U+predist); _sse_load((*phi[ix]).s0); _sse_su3_inverse_multiply((*U)); _sse_vector_cmplxcg_mul(ka2); _sse_load(rs.s0); _sse_vector_add(); _sse_store(rs.s0); _sse_load(rs.s3); _sse_vector_sub(); _sse_store(rs.s3); _sse_load((*phi[ix]).s1); _sse_su3_inverse_multiply((*U)); _sse_vector_cmplxcg_mul(ka2); _sse_load(rs.s1); _sse_vector_add(); _sse_store(rs.s1); _sse_load(rs.s2); _sse_vector_add(); _sse_store(rs.s2); ix++; U++; /*********************** direction +3 ************************/ _sse_load_up((*phi[ix]).s0); _sse_load(rs.s0); _sse_vector_add(); _sse_store(rs.s0); _sse_load(rs.s2); _sse_vector_i_mul(); _sse_vector_sub(); _sse_store(rs.s2); _sse_load_up((*phi[ix]).s1); _sse_load(rs.s1); _sse_vector_add(); _sse_store(rs.s1); _sse_load(rs.s3); _sse_vector_i_mul(); _sse_vector_add(); _sse_store(rs.s3); ix++; /*********************** direction -3 ************************/ _prefetch_su3(U+predist); _prefetch_spinor(s+1); _sse_load((*phi[ix]).s0); _sse_su3_inverse_multiply((*U)); _sse_vector_cmplxcg_mul(ka3); _sse_load(rs.s0); _sse_vector_add(); _sse_store_nt((*s).s0); _sse_load(rs.s2); _sse_vector_i_mul(); _sse_vector_add(); _sse_store_nt((*s).s2); _sse_load((*phi[ix]).s1); _sse_su3_inverse_multiply((*U)); _sse_vector_cmplxcg_mul(ka3); _sse_load(rs.s1); _sse_vector_add(); _sse_store_nt((*s).s1); _sse_load(rs.s3); _sse_vector_i_mul(); _sse_vector_sub(); _sse_store_nt((*s).s3); ix++; U++; s++; } #ifdef _KOJAK_INST #pragma pomp inst end(hoppingmatrix) #endif }