void cloverndpoly_heatbath(const int id, hamiltonian_field_t * const hf) { int j; monomial * mnl = &monomial_list[id]; spinor *up0, *dn0, *up1, *dn1, *dummy; double atime, etime; atime = gettime(); ndpoly_set_global_parameter(mnl, 0); g_mu3 = 0.; init_sw_fields(); sw_term((const su3**)hf->gaugefield, mnl->kappa, mnl->c_sw); sw_invert_nd(mnl->mubar*mnl->mubar - mnl->epsbar*mnl->epsbar); // we measure before trajectory! if((mnl->rec_ev != 0) && (hf->traj_counter%mnl->rec_ev == 0)) { phmc_compute_ev(hf->traj_counter-1, id, &Qsw_pm_ndbipsi); } mnl->energy0 = 0.; random_spinor_field_eo(g_chi_up_spinor_field[0], mnl->rngrepro, RN_GAUSS); mnl->energy0 = square_norm(g_chi_up_spinor_field[0], VOLUME/2, 1); random_spinor_field_eo(g_chi_dn_spinor_field[0], mnl->rngrepro, RN_GAUSS); mnl->energy0 += square_norm(g_chi_dn_spinor_field[0], VOLUME/2, 1); Qsw_ndpsi(g_chi_up_spinor_field[1], g_chi_dn_spinor_field[1], g_chi_up_spinor_field[0], g_chi_dn_spinor_field[0]); up0 = g_chi_up_spinor_field[0]; up1 = g_chi_up_spinor_field[1]; dn0 = g_chi_dn_spinor_field[0]; dn1 = g_chi_dn_spinor_field[1]; for(j = 1; j < (mnl->MDPolyDegree); j++){ Qsw_tau1_sub_const_ndpsi(up0, dn0, up1, dn1, mnl->MDPolyRoots[mnl->MDPolyDegree-2+j]); dummy = up1; up1 = up0; up0 = dummy; dummy = dn1; dn1 = dn0; dn0 = dummy; } Ptilde_ndpsi(up0, dn0, mnl->PtildeCoefs, mnl->PtildeDegree, up1, dn1, &Qsw_pm_ndpsi); assign(mnl->pf, up0, VOLUME/2); assign(mnl->pf2, dn0, VOLUME/2); etime = gettime(); if(g_proc_id == 0) { if(g_debug_level > 1) { printf("# Time for %s monomial heatbath: %e s\n", mnl->name, etime-atime); } if(g_debug_level > 3) { printf("called cloverndpoly_heatbath for id %d energy %f\n", id, mnl->energy0); } } return; }
double cloverndpoly_acc(const int id, hamiltonian_field_t * const hf) { int j; monomial * mnl = &monomial_list[id]; spinor *up0, *dn0, *up1, *dn1, *dummy; double atime, etime; atime = gettime(); ndpoly_set_global_parameter(mnl, 0); g_mu3 = 0.; sw_term((const su3**) hf->gaugefield, mnl->kappa, mnl->c_sw); sw_invert_nd(mnl->mubar*mnl->mubar - mnl->epsbar*mnl->epsbar); mnl->energy1 = 0.; up0 = g_chi_up_spinor_field[0]; up1 = g_chi_up_spinor_field[1]; dn0 = g_chi_dn_spinor_field[0]; dn1 = g_chi_dn_spinor_field[1]; /* This is needed if we consider only "1" in eq. 9 */ assign(up0, mnl->pf , VOLUME/2); assign(dn0, mnl->pf2, VOLUME/2); for(j = 1; j <= (mnl->MDPolyDegree-1); j++) { Qsw_tau1_sub_const_ndpsi(up1, dn1, up0, dn0, mnl->MDPolyRoots[j-1]); dummy = up1; up1 = up0; up0 = dummy; dummy = dn1; dn1 = dn0; dn0 = dummy; /* result always in up0 and dn0 */ } mnl->energy1 = square_norm(up0, VOLUME/2, 1); mnl->energy1 += square_norm(dn0, VOLUME/2, 1); etime = gettime(); if(g_proc_id == 0) { if(g_debug_level > 1) { printf("# Time for %s monomial acc step: %e s\n", mnl->name, etime-atime); } if(g_debug_level > 3) { printf("called cloverndpoly_acc for id %d dH = %1.10e\n", id, mnl->energy1 - mnl->energy0); } } return(mnl->energy1 - mnl->energy0); }
void op_invert(const int op_id, const int index_start, const int write_prop) { operator * optr = &operator_list[op_id]; double atime = 0., etime = 0., nrm1 = 0., nrm2 = 0.; int i; optr->iterations = 0; optr->reached_prec = -1.; g_kappa = optr->kappa; boundary(g_kappa); atime = gettime(); if(optr->type == TMWILSON || optr->type == WILSON || optr->type == CLOVER) { g_mu = optr->mu; g_c_sw = optr->c_sw; if(optr->type == CLOVER) { if (g_cart_id == 0 && g_debug_level > 1) { printf("#\n# csw = %e, computing clover leafs\n", g_c_sw); } init_sw_fields(VOLUME); sw_term( (const su3**) g_gauge_field, optr->kappa, optr->c_sw); /* this must be EE here! */ /* to match clover_inv in Qsw_psi */ sw_invert(EE, optr->mu); } for(i = 0; i < 2; i++) { if (g_cart_id == 0) { printf("#\n# 2 kappa mu = %e, kappa = %e, c_sw = %e\n", g_mu, g_kappa, g_c_sw); } if(optr->type != CLOVER) { if(use_preconditioning){ g_precWS=(void*)optr->precWS; } else { g_precWS=NULL; } optr->iterations = invert_eo( optr->prop0, optr->prop1, optr->sr0, optr->sr1, optr->eps_sq, optr->maxiter, optr->solver, optr->rel_prec, 0, optr->even_odd_flag,optr->no_extra_masses, optr->extra_masses, optr->id ); /* check result */ M_full(g_spinor_field[4], g_spinor_field[5], optr->prop0, optr->prop1); } else { optr->iterations = invert_clover_eo(optr->prop0, optr->prop1, optr->sr0, optr->sr1, optr->eps_sq, optr->maxiter, optr->solver, optr->rel_prec, &g_gauge_field, &Qsw_pm_psi, &Qsw_minus_psi); /* check result */ Msw_full(g_spinor_field[4], g_spinor_field[5], optr->prop0, optr->prop1); } diff(g_spinor_field[4], g_spinor_field[4], optr->sr0, VOLUME / 2); diff(g_spinor_field[5], g_spinor_field[5], optr->sr1, VOLUME / 2); nrm1 = square_norm(g_spinor_field[4], VOLUME / 2, 1); nrm2 = square_norm(g_spinor_field[5], VOLUME / 2, 1); optr->reached_prec = nrm1 + nrm2; /* convert to standard normalisation */ /* we have to mult. by 2*kappa */ if (optr->kappa != 0.) { mul_r(optr->prop0, (2*optr->kappa), optr->prop0, VOLUME / 2); mul_r(optr->prop1, (2*optr->kappa), optr->prop1, VOLUME / 2); } if (optr->solver != CGMMS && write_prop) /* CGMMS handles its own I/O */ optr->write_prop(op_id, index_start, i); if(optr->DownProp) { optr->mu = -optr->mu; } else break; } } else if(optr->type == DBTMWILSON || optr->type == DBCLOVER) { g_mubar = optr->mubar; g_epsbar = optr->epsbar; g_c_sw = 0.; if(optr->type == DBCLOVER) { g_c_sw = optr->c_sw; if (g_cart_id == 0 && g_debug_level > 1) { printf("#\n# csw = %e, computing clover leafs\n", g_c_sw); } init_sw_fields(VOLUME); sw_term( (const su3**) g_gauge_field, optr->kappa, optr->c_sw); sw_invert_nd(optr->mubar*optr->mubar-optr->epsbar*optr->epsbar); } for(i = 0; i < SourceInfo.no_flavours; i++) { if(optr->type != DBCLOVER) { optr->iterations = invert_doublet_eo( optr->prop0, optr->prop1, optr->prop2, optr->prop3, optr->sr0, optr->sr1, optr->sr2, optr->sr3, optr->eps_sq, optr->maxiter, optr->solver, optr->rel_prec); } else { optr->iterations = invert_cloverdoublet_eo( optr->prop0, optr->prop1, optr->prop2, optr->prop3, optr->sr0, optr->sr1, optr->sr2, optr->sr3, optr->eps_sq, optr->maxiter, optr->solver, optr->rel_prec); } g_mu = optr->mubar; if(optr->type != DBCLOVER) { M_full(g_spinor_field[DUM_DERI+1], g_spinor_field[DUM_DERI+2], optr->prop0, optr->prop1); } else { Msw_full(g_spinor_field[DUM_DERI+1], g_spinor_field[DUM_DERI+2], optr->prop0, optr->prop1); } assign_add_mul_r(g_spinor_field[DUM_DERI+1], optr->prop2, -optr->epsbar, VOLUME/2); assign_add_mul_r(g_spinor_field[DUM_DERI+2], optr->prop3, -optr->epsbar, VOLUME/2); g_mu = -g_mu; if(optr->type != DBCLOVER) { M_full(g_spinor_field[DUM_DERI+3], g_spinor_field[DUM_DERI+4], optr->prop2, optr->prop3); } else { Msw_full(g_spinor_field[DUM_DERI+3], g_spinor_field[DUM_DERI+4], optr->prop2, optr->prop3); } assign_add_mul_r(g_spinor_field[DUM_DERI+3], optr->prop0, -optr->epsbar, VOLUME/2); assign_add_mul_r(g_spinor_field[DUM_DERI+4], optr->prop1, -optr->epsbar, VOLUME/2); diff(g_spinor_field[DUM_DERI+1], g_spinor_field[DUM_DERI+1], optr->sr0, VOLUME/2); diff(g_spinor_field[DUM_DERI+2], g_spinor_field[DUM_DERI+2], optr->sr1, VOLUME/2); diff(g_spinor_field[DUM_DERI+3], g_spinor_field[DUM_DERI+3], optr->sr2, VOLUME/2); diff(g_spinor_field[DUM_DERI+4], g_spinor_field[DUM_DERI+4], optr->sr3, VOLUME/2); nrm1 = square_norm(g_spinor_field[DUM_DERI+1], VOLUME/2, 1); nrm1 += square_norm(g_spinor_field[DUM_DERI+2], VOLUME/2, 1); nrm1 += square_norm(g_spinor_field[DUM_DERI+3], VOLUME/2, 1); nrm1 += square_norm(g_spinor_field[DUM_DERI+4], VOLUME/2, 1); optr->reached_prec = nrm1; g_mu = g_mu1; /* For standard normalisation */ /* we have to mult. by 2*kappa */ mul_r(g_spinor_field[DUM_DERI], (2*optr->kappa), optr->prop0, VOLUME/2); mul_r(g_spinor_field[DUM_DERI+1], (2*optr->kappa), optr->prop1, VOLUME/2); mul_r(g_spinor_field[DUM_DERI+2], (2*optr->kappa), optr->prop2, VOLUME/2); mul_r(g_spinor_field[DUM_DERI+3], (2*optr->kappa), optr->prop3, VOLUME/2); /* the final result should be stored in the convention used in */ /* hep-lat/0606011 */ /* this requires multiplication of source with */ /* (1+itau_2)/sqrt(2) and the result with (1-itau_2)/sqrt(2) */ mul_one_pm_itau2(optr->prop0, optr->prop2, g_spinor_field[DUM_DERI], g_spinor_field[DUM_DERI+2], -1., VOLUME/2); mul_one_pm_itau2(optr->prop1, optr->prop3, g_spinor_field[DUM_DERI+1], g_spinor_field[DUM_DERI+3], -1., VOLUME/2); /* write propagator */ if(write_prop) optr->write_prop(op_id, index_start, i); mul_r(optr->prop0, 1./(2*optr->kappa), g_spinor_field[DUM_DERI], VOLUME/2); mul_r(optr->prop1, 1./(2*optr->kappa), g_spinor_field[DUM_DERI+1], VOLUME/2); mul_r(optr->prop2, 1./(2*optr->kappa), g_spinor_field[DUM_DERI+2], VOLUME/2); mul_r(optr->prop3, 1./(2*optr->kappa), g_spinor_field[DUM_DERI+3], VOLUME/2); /* mirror source, but not for volume sources */ if(i == 0 && SourceInfo.no_flavours == 2 && SourceInfo.type != 1) { if (g_cart_id == 0) { fprintf(stdout, "# Inversion done in %d iterations, squared residue = %e!\n", optr->iterations, optr->reached_prec); } mul_one_pm_itau2(g_spinor_field[DUM_DERI], g_spinor_field[DUM_DERI+2], optr->sr0, optr->sr2, -1., VOLUME/2); mul_one_pm_itau2(g_spinor_field[DUM_DERI+1], g_spinor_field[DUM_DERI+3], optr->sr1, optr->sr3, -1., VOLUME/2); mul_one_pm_itau2(optr->sr0, optr->sr2, g_spinor_field[DUM_DERI+2], g_spinor_field[DUM_DERI], +1., VOLUME/2); mul_one_pm_itau2(optr->sr1, optr->sr3, g_spinor_field[DUM_DERI+3], g_spinor_field[DUM_DERI+1], +1., VOLUME/2); } /* volume sources need only one inversion */ else if(SourceInfo.type == 1) i++; } } else if(optr->type == OVERLAP) { g_mu = 0.; m_ov=optr->m; eigenvalues(&optr->no_ev, 5000, optr->ev_prec, 0, optr->ev_readwrite, nstore, optr->even_odd_flag); /* ov_check_locality(); */ /* index_jd(&optr->no_ev_index, 5000, 1.e-12, optr->conf_input, nstore, 4); */ ov_n_cheby=optr->deg_poly; if(use_preconditioning==1) g_precWS=(void*)optr->precWS; else g_precWS=NULL; if(g_debug_level > 3) ov_check_ginsparg_wilson_relation_strong(); invert_overlap(op_id, index_start); if(write_prop) optr->write_prop(op_id, index_start, 0); } etime = gettime(); if (g_cart_id == 0 && g_debug_level > 0) { fprintf(stdout, "# Inversion done in %d iterations, squared residue = %e!\n", optr->iterations, optr->reached_prec); fprintf(stdout, "# Inversion done in %1.2e sec. \n", etime - atime); } return; }
void ndratcor_heatbath(const int id, hamiltonian_field_t * const hf) { monomial * mnl = &monomial_list[id]; double atime, etime, delta; spinor * up0, * dn0, * up1, * dn1, * tup, * tdn, * Zup, * Zdn; double coefs[6] = {1./4., -3./32., 7./128., -77./2048., 231./8192., -1463./65536.}; // series of (1+x)^(1/4) double coefs_check[6] = {1./2., -1./8., 1./16., -5./128., 7./256., -21./1024.}; // series of (1+x)^(1/2) atime = gettime(); nd_set_global_parameter(mnl); g_mu3 = 0.; mnl->iter0 = 0; if(mnl->type == NDCLOVERRATCOR) { init_sw_fields(); sw_term((const su3**)hf->gaugefield, mnl->kappa, mnl->c_sw); sw_invert_nd(mnl->mubar*mnl->mubar - mnl->epsbar*mnl->epsbar); copy_32_sw_fields(); } // we measure before the trajectory! if((mnl->rec_ev != 0) && (hf->traj_counter%mnl->rec_ev == 0)) { if(mnl->type != NDCLOVERRAT) phmc_compute_ev(hf->traj_counter-1, id, &Qtm_pm_ndbipsi); else phmc_compute_ev(hf->traj_counter-1, id, &Qsw_pm_ndbipsi); } // the Gaussian distributed random fields mnl->energy0 = 0.; random_spinor_field_eo(mnl->pf, mnl->rngrepro, RN_GAUSS); mnl->energy0 = square_norm(mnl->pf, VOLUME/2, 1); random_spinor_field_eo(mnl->pf2, mnl->rngrepro, RN_GAUSS); mnl->energy0 += square_norm(mnl->pf2, VOLUME/2, 1); mnl->solver_params.max_iter = mnl->maxiter; mnl->solver_params.squared_solver_prec = mnl->accprec; mnl->solver_params.no_shifts = mnl->rat.np; mnl->solver_params.shifts = mnl->rat.mu; mnl->solver_params.type = mnl->solver; mnl->solver_params.M_ndpsi = &Qtm_pm_ndpsi; mnl->solver_params.M_ndpsi32 = &Qtm_pm_ndpsi_32; if(mnl->type == NDCLOVERRATCOR) { mnl->solver_params.M_ndpsi = &Qsw_pm_ndpsi; mnl->solver_params.M_ndpsi32 = &Qsw_pm_ndpsi_32; } mnl->solver_params.sdim = VOLUME/2; mnl->solver_params.rel_prec = g_relative_precision_flag; // apply B to the random field to generate pseudo-fermion fields up0 = mnl->w_fields[0]; dn0 = mnl->w_fields[1]; up1 = mnl->w_fields[2]; dn1 = mnl->w_fields[3]; Zup = mnl->w_fields[4]; Zdn = mnl->w_fields[5]; apply_Z_ndpsi(up0, dn0, mnl->pf, mnl->pf2, id, hf, &(mnl->solver_params)); // computing correction to energy1 delta = coefs_check[0]*(scalar_prod_r(mnl->pf, up0, VOLUME/2, 1) + scalar_prod_r(mnl->pf2, dn0, VOLUME/2, 1)); if(g_debug_level > 2 && g_proc_id == 0) printf("# NDRATCOR heatbath: c_%d*(R * Z^%d * R) = %e\n", 1, 1, delta); // debug for showing that the old check was giving a smaller delta if(g_debug_level > 3) { double delta_old = square_norm(up0, VOLUME/2, 1) + square_norm(dn0, VOLUME/2, 1); if(g_proc_id == 0) { printf("# NDRATCOR old check: || Z^%d * R ||^2 = %e\n", 1, delta_old); printf("# NDRATCOR new check: (c_%d*(R * Z^%d * R))^2 = %e\n", 1, 1, delta*delta); } } if(delta*delta > mnl->accprec) { assign_add_mul_r(mnl->pf, up0, coefs[0], VOLUME/2); assign_add_mul_r(mnl->pf2, dn0, coefs[0], VOLUME/2); // saving first application assign(Zup, up0, VOLUME/2); assign(Zdn, dn0, VOLUME/2); for(int i = 2; i < 8; i++) { // computing next order correction to energy1 delta = coefs_check[i-1]*(scalar_prod_r(Zup, up0, VOLUME/2, 1) + scalar_prod_r(Zup, dn0, VOLUME/2, 1)); if(g_debug_level > 2 && g_proc_id == 0) printf("# NDRATCOR heatbath: c_%d*(R * Z^%d * R) = %e\n", i, i, delta); // debug for showing that the old check was giving a smaller delta if(g_debug_level > 3) { double delta_old = square_norm(up0, VOLUME/2, 1) + square_norm(dn0, VOLUME/2, 1); if(g_proc_id == 0) { printf("# NDRATCOR old check: || Z^%d * R ||^2 = %e\n", 1, delta_old); printf("# NDRATCOR new check: (c_%d*(R * Z^%d * R))^2 = %e\n", 1, 1, delta*delta); } } if(delta*delta < mnl->accprec) break; apply_Z_ndpsi(up1, dn1, up0, dn0, id, hf, &(mnl->solver_params)); assign_add_mul_r(mnl->pf, up1, coefs[i-1], VOLUME/2); assign_add_mul_r(mnl->pf2, dn1, coefs[i-1], VOLUME/2); tup = up0; tdn = dn0; up0 = up1; dn0 = dn1; up1 = tup; dn1 = tdn; } } etime = gettime(); if(g_proc_id == 0) { if(g_debug_level > 1) { printf("# Time for %s monomial heatbath: %e s\n", mnl->name, etime-atime); } if(g_debug_level > 3) { printf("called ndratcor_heatbath for id %d energy %f\n", id, mnl->energy0); } } return; }
double ndratcor_acc(const int id, hamiltonian_field_t * const hf) { monomial * mnl = &monomial_list[id]; double atime, etime, delta; spinor * up0, * dn0, * up1, * dn1, * tup, * tdn; double coefs[6] = {-1./2., 3./8., -5./16., 35./128., -63./256., 231./1024.}; atime = gettime(); nd_set_global_parameter(mnl); g_mu3 = 0.; if(mnl->type == NDCLOVERRATCOR) { sw_term((const su3**) hf->gaugefield, mnl->kappa, mnl->c_sw); sw_invert_nd(mnl->mubar*mnl->mubar - mnl->epsbar*mnl->epsbar); copy_32_sw_fields(); } mnl->energy1 = square_norm(mnl->pf, VOLUME/2, 1) + square_norm(mnl->pf2, VOLUME/2, 1); mnl->solver_params.max_iter = mnl->maxiter; mnl->solver_params.squared_solver_prec = mnl->accprec; mnl->solver_params.no_shifts = mnl->rat.np; mnl->solver_params.shifts = mnl->rat.mu; mnl->solver_params.type = mnl->solver; mnl->solver_params.M_ndpsi = &Qtm_pm_ndpsi; mnl->solver_params.M_ndpsi32 = &Qtm_pm_ndpsi_32; if(mnl->type == NDCLOVERRATCOR) { mnl->solver_params.M_ndpsi = &Qsw_pm_ndpsi; mnl->solver_params.M_ndpsi32 = &Qsw_pm_ndpsi_32; } mnl->solver_params.sdim = VOLUME/2; mnl->solver_params.rel_prec = g_relative_precision_flag; // apply (Q R)^(-1) to pseudo-fermion fields up0 = mnl->w_fields[0]; dn0 = mnl->w_fields[1]; up1 = mnl->w_fields[2]; dn1 = mnl->w_fields[3]; apply_Z_ndpsi(up0, dn0, mnl->pf, mnl->pf2, id, hf, &(mnl->solver_params)); delta = coefs[0]*(scalar_prod_r(mnl->pf, up0, VOLUME/2, 1) + scalar_prod_r(mnl->pf2, dn0, VOLUME/2, 1)); mnl->energy1 += delta; if(g_debug_level > 2 && g_proc_id == 0) printf("# NDRATCOR acc step: c_%d*(phi * Z^%d * phi) = %e\n", 1, 1, delta); for(int i = 2; i < 8; i++) { if(delta*delta < mnl->accprec) break; delta = coefs[i-1]*(square_norm(up0, VOLUME/2, 1) + square_norm(dn0, VOLUME/2, 1)); mnl->energy1 += delta; if(g_debug_level > 2 && g_proc_id == 0) printf("# NDRATCOR acc step: c_%d*(phi * Z^%d * phi) = %e\n", i, i, delta); i++; //incrementing i if(delta*delta < mnl->accprec) break; apply_Z_ndpsi(up1, dn1, up0, dn0, id, hf, &(mnl->solver_params)); delta = coefs[i-1]*(scalar_prod_r(up0, up1, VOLUME/2, 1) + scalar_prod_r(dn0, dn1, VOLUME/2, 1)); mnl->energy1 += delta; if(g_debug_level > 2 && g_proc_id == 0) printf("# NDRATCOR acc step: c_%d*(phi * Z^%d * phi) = %e\n", i, i, delta); tup = up0; tdn = dn0; up0 = up1; dn0 = dn1; up1 = tup; dn1 = tdn; } etime = gettime(); if(g_proc_id == 0) { if(g_debug_level > 1) { printf("# Time for %s monomial acc step: %e s\n", mnl->name, etime-atime); } if(g_debug_level > 3) { // shoud be 3 printf("called ndratcor_acc for id %d dH = %1.10e\n", id, mnl->energy1 - mnl->energy0); } } return(mnl->energy1 - mnl->energy0); }
void cloverndpoly_derivative(const int id, hamiltonian_field_t * const hf) { int j, k; monomial * mnl = &monomial_list[id]; double atime, etime; atime = gettime(); for(int i = 0; i < VOLUME; i++) { for(int mu = 0; mu < 4; mu++) { _su3_zero(swm[i][mu]); _su3_zero(swp[i][mu]); } } ndpoly_set_global_parameter(mnl, 0); // we compute the clover term (1 + T_ee(oo)) for all sites x sw_term( (const su3**) hf->gaugefield, mnl->kappa, mnl->c_sw); // we invert it for the even sites only sw_invert_nd(mnl->mubar*mnl->mubar - mnl->epsbar*mnl->epsbar); mnl->forcefactor = -phmc_Cpol*mnl->EVMaxInv; /* Recall: The GAMMA_5 left of delta M_eo is done in deriv_Sb !!! */ /* Here comes the definitions for the chi_j fields */ /* from j=0 (chi_0 = phi) ..... to j = n-1 */ /* in g_chi_up_spinor_field[0] (g_chi_dn_spinor_field[0] we expect */ /* to find the phi field, the pseudo fermion field */ /* i.e. must be equal to mnl->pf (mnl->pf2) */ assign(g_chi_up_spinor_field[0], mnl->pf, VOLUME/2); assign(g_chi_dn_spinor_field[0], mnl->pf2, VOLUME/2); for(k = 1; k < (mnl->MDPolyDegree-1); k++) { Qsw_tau1_sub_const_ndpsi(g_chi_up_spinor_field[k], g_chi_dn_spinor_field[k], g_chi_up_spinor_field[k-1], g_chi_dn_spinor_field[k-1], mnl->MDPolyRoots[k-1]); } /* Here comes the remaining fields chi_k ; k=n,...,2n-1 */ /*They are evaluated step-by-step overwriting the same field (mnl->MDPolyDegree)*/ assign(g_chi_up_spinor_field[mnl->MDPolyDegree], g_chi_up_spinor_field[mnl->MDPolyDegree-2], VOLUME/2); assign(g_chi_dn_spinor_field[mnl->MDPolyDegree], g_chi_dn_spinor_field[mnl->MDPolyDegree-2], VOLUME/2); for(j = (mnl->MDPolyDegree-1); j > 0; j--) { assign(g_chi_up_spinor_field[mnl->MDPolyDegree-1], g_chi_up_spinor_field[mnl->MDPolyDegree], VOLUME/2); assign(g_chi_dn_spinor_field[mnl->MDPolyDegree-1], g_chi_dn_spinor_field[mnl->MDPolyDegree], VOLUME/2); Qsw_tau1_sub_const_ndpsi(g_chi_up_spinor_field[mnl->MDPolyDegree], g_chi_dn_spinor_field[mnl->MDPolyDegree], g_chi_up_spinor_field[mnl->MDPolyDegree-1], g_chi_dn_spinor_field[mnl->MDPolyDegree-1], mnl->MDPolyRoots[2*mnl->MDPolyDegree-3-j]); /* Get the even parts of the (j-1)th chi_spinors */ H_eo_sw_ndpsi(mnl->w_fields[0], mnl->w_fields[1], g_chi_up_spinor_field[j-1], g_chi_dn_spinor_field[j-1]); /* \delta M_eo sandwitched by chi[j-1]_e^\dagger and chi[2N-j]_o */ deriv_Sb(EO, mnl->w_fields[0], g_chi_up_spinor_field[mnl->MDPolyDegree], hf, mnl->forcefactor);/* UP */ deriv_Sb(EO, mnl->w_fields[1], g_chi_dn_spinor_field[mnl->MDPolyDegree], hf, mnl->forcefactor);/* DN */ /* Get the even parts of the (2N-j)-th chi_spinors */ H_eo_sw_ndpsi(mnl->w_fields[2], mnl->w_fields[3], g_chi_up_spinor_field[mnl->MDPolyDegree], g_chi_dn_spinor_field[mnl->MDPolyDegree]); /* \delta M_oe sandwitched by chi[j-1]_o^\dagger and chi[2N-j]_e */ deriv_Sb(OE, g_chi_up_spinor_field[j-1], mnl->w_fields[2], hf, mnl->forcefactor); deriv_Sb(OE, g_chi_dn_spinor_field[j-1], mnl->w_fields[3], hf, mnl->forcefactor); // even/even sites sandwiched by gamma_5 Y_e and gamma_5 X_e sw_spinor(EE, mnl->w_fields[3], mnl->w_fields[0], mnl->forcefactor); // odd/odd sites sandwiched by gamma_5 Y_o and gamma_5 X_o sw_spinor(OO, g_chi_up_spinor_field[j-1], g_chi_dn_spinor_field[mnl->MDPolyDegree], mnl->forcefactor); // even/even sites sandwiched by gamma_5 Y_e and gamma_5 X_e sw_spinor(EE, mnl->w_fields[2], mnl->w_fields[1], mnl->forcefactor); // odd/odd sites sandwiched by gamma_5 Y_o and gamma_5 X_o sw_spinor(OO, g_chi_dn_spinor_field[j-1], g_chi_up_spinor_field[mnl->MDPolyDegree], mnl->forcefactor); } // trlog part does not depend on the normalisation of the polynomial sw_deriv_nd(EE); sw_all(hf, mnl->kappa, mnl->c_sw); etime = gettime(); if(g_debug_level > 1 && g_proc_id == 0) { printf("# Time for %s monomial derivative: %e s\n", mnl->name, etime-atime); } return; }