/************************************************************************************* Grid physics library, www.github.com/paboyle/Grid Source file: ./lib/qcd/action/fermion/CayleyFermion5D.cc Copyright (C) 2015 Author: Peter Boyle Author: Peter Boyle Author: Peter Boyle Author: paboyle This program is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version. This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with this program; if not, write to the Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. See the full license in the file "LICENSE" in the top level distribution directory *************************************************************************************/ /* END LEGAL */ #include #include #include namespace Grid { namespace QCD { template CayleyFermion5D::CayleyFermion5D(GaugeField &_Umu, GridCartesian &FiveDimGrid, GridRedBlackCartesian &FiveDimRedBlackGrid, GridCartesian &FourDimGrid, GridRedBlackCartesian &FourDimRedBlackGrid, RealD _mass,RealD _M5,const ImplParams &p) : WilsonFermion5D(_Umu, FiveDimGrid, FiveDimRedBlackGrid, FourDimGrid, FourDimRedBlackGrid,_M5,p), mass(_mass) { } template void CayleyFermion5D::Dminus(const FermionField &psi, FermionField &chi) { int Ls=this->Ls; FermionField tmp_f(this->FermionGrid()); this->DW(psi,tmp_f,DaggerNo); for(int s=0;s void CayleyFermion5D::DminusDag(const FermionField &psi, FermionField &chi) { int Ls=this->Ls; FermionField tmp_f(this->FermionGrid()); this->DW(psi,tmp_f,DaggerYes); for(int s=0;s void CayleyFermion5D::CayleyReport(void) { this->Report(); std::vector latt = GridDefaultLatt(); RealD volume = this->Ls; for(int mu=0;mu_FourDimGrid->_Nprocessors; if ( M5Dcalls > 0 ) { std::cout << GridLogMessage << "#### M5D calls report " << std::endl; std::cout << GridLogMessage << "CayleyFermion5D Number of M5D Calls : " << M5Dcalls << std::endl; std::cout << GridLogMessage << "CayleyFermion5D ComputeTime/Calls : " << M5Dtime / M5Dcalls << " us" << std::endl; // Flops = 6.0*(Nc*Ns) *Ls*vol RealD mflops = 6.0*12*volume*M5Dcalls/M5Dtime/2; // 2 for red black counting std::cout << GridLogMessage << "Average mflops/s per call : " << mflops << std::endl; std::cout << GridLogMessage << "Average mflops/s per call per rank : " << mflops/NP << std::endl; } if ( MooeeInvCalls > 0 ) { std::cout << GridLogMessage << "#### MooeeInv calls report " << std::endl; std::cout << GridLogMessage << "CayleyFermion5D Number of MooeeInv Calls : " << MooeeInvCalls << std::endl; std::cout << GridLogMessage << "CayleyFermion5D ComputeTime/Calls : " << MooeeInvTime / MooeeInvCalls << " us" << std::endl; // Flops = MADD * Ls *Ls *4dvol * spin/colour/complex RealD mflops = 2.0*24*this->Ls*volume*MooeeInvCalls/MooeeInvTime/2; // 2 for red black counting std::cout << GridLogMessage << "Average mflops/s per call : " << mflops << std::endl; std::cout << GridLogMessage << "Average mflops/s per call per rank : " << mflops/NP << std::endl; } } template void CayleyFermion5D::CayleyZeroCounters(void) { this->ZeroCounters(); M5Dflops=0; M5Dcalls=0; M5Dtime=0; MooeeInvFlops=0; MooeeInvCalls=0; MooeeInvTime=0; } template void CayleyFermion5D::M5D (const FermionField &psi, FermionField &chi) { int Ls=this->Ls; std::vector diag (Ls,1.0); std::vector upper(Ls,-1.0); upper[Ls-1]=mass; std::vector lower(Ls,-1.0); lower[0] =mass; M5D(psi,chi,chi,lower,diag,upper); } template void CayleyFermion5D::Meooe5D (const FermionField &psi, FermionField &Din) { int Ls=this->Ls; std::vector diag = bs; std::vector upper= cs; std::vector lower= cs; upper[Ls-1]=-mass*upper[Ls-1]; lower[0] =-mass*lower[0]; M5D(psi,psi,Din,lower,diag,upper); } // FIXME Redunant with the above routine; check this and eliminate template void CayleyFermion5D::Meo5D (const FermionField &psi, FermionField &chi) { int Ls=this->Ls; std::vector diag = beo; std::vector upper(Ls); std::vector lower(Ls); for(int i=0;i void CayleyFermion5D::Mooee (const FermionField &psi, FermionField &chi) { int Ls=this->Ls; std::vector diag = bee; std::vector upper(Ls); std::vector lower(Ls); for(int i=0;i void CayleyFermion5D::MooeeDag (const FermionField &psi, FermionField &chi) { int Ls=this->Ls; std::vector diag = bee; std::vector upper(Ls); std::vector lower(Ls); for (int s=0;s void CayleyFermion5D::M5Ddag (const FermionField &psi, FermionField &chi) { int Ls=this->Ls; std::vector diag(Ls,1.0); std::vector upper(Ls,-1.0); std::vector lower(Ls,-1.0); upper[Ls-1]=-mass*upper[Ls-1]; lower[0] =-mass*lower[0]; M5Ddag(psi,chi,chi,lower,diag,upper); } template void CayleyFermion5D::MeooeDag5D (const FermionField &psi, FermionField &Din) { int Ls=this->Ls; std::vector diag =bs; std::vector upper=cs; std::vector lower=cs; for (int s=0;s RealD CayleyFermion5D::M (const FermionField &psi, FermionField &chi) { int Ls=this->Ls; FermionField Din(psi._grid); // Assemble Din Meooe5D(psi,Din); this->DW(Din,chi,DaggerNo); // ((b D_W + D_w hop terms +1) on s-diag axpby(chi,1.0,1.0,chi,psi); M5D(psi,chi); return(norm2(chi)); } template RealD CayleyFermion5D::Mdag (const FermionField &psi, FermionField &chi) { // Under adjoint //D1+ D1- P- -> D1+^dag P+ D2-^dag //D2- P+ D2+ P-D1-^dag D2+dag FermionField Din(psi._grid); // Apply Dw this->DW(psi,Din,DaggerYes); MeooeDag5D(Din,chi); M5Ddag(psi,chi); // ((b D_W + D_w hop terms +1) on s-diag axpby (chi,1.0,1.0,chi,psi); return norm2(chi); } // half checkerboard operations template void CayleyFermion5D::Meooe (const FermionField &psi, FermionField &chi) { int Ls=this->Ls; Meooe5D(psi,this->tmp()); if ( psi.checkerboard == Odd ) { this->DhopEO(this->tmp(),chi,DaggerNo); } else { this->DhopOE(this->tmp(),chi,DaggerNo); } } template void CayleyFermion5D::MeooeDag (const FermionField &psi, FermionField &chi) { // Apply 4d dslash if ( psi.checkerboard == Odd ) { this->DhopEO(psi,this->tmp(),DaggerYes); } else { this->DhopOE(psi,this->tmp(),DaggerYes); } MeooeDag5D(this->tmp(),chi); } template void CayleyFermion5D::Mdir (const FermionField &psi, FermionField &chi,int dir,int disp){ Meo5D(psi,this->tmp()); // Apply 4d dslash fragment this->DhopDir(this->tmp(),chi,dir,disp); } // force terms; five routines; default to Dhop on diagonal template void CayleyFermion5D::MDeriv (GaugeField &mat,const FermionField &U,const FermionField &V,int dag) { FermionField Din(V._grid); if ( dag == DaggerNo ) { // U d/du [D_w D5] V = U d/du DW D5 V Meooe5D(V,Din); this->DhopDeriv(mat,U,Din,dag); } else { // U d/du [D_w D5]^dag V = U D5^dag d/du DW^dag Y // implicit adj on U in call Meooe5D(U,Din); this->DhopDeriv(mat,Din,V,dag); } }; template void CayleyFermion5D::MoeDeriv(GaugeField &mat,const FermionField &U,const FermionField &V,int dag) { FermionField Din(V._grid); if ( dag == DaggerNo ) { // U d/du [D_w D5] V = U d/du DW D5 V Meooe5D(V,Din); this->DhopDerivOE(mat,U,Din,dag); } else { // U d/du [D_w D5]^dag V = U D5^dag d/du DW^dag Y // implicit adj on U in call Meooe5D(U,Din); this->DhopDerivOE(mat,Din,V,dag); } }; template void CayleyFermion5D::MeoDeriv(GaugeField &mat,const FermionField &U,const FermionField &V,int dag) { FermionField Din(V._grid); if ( dag == DaggerNo ) { // U d/du [D_w D5] V = U d/du DW D5 V Meooe5D(V,Din); this->DhopDerivEO(mat,U,Din,dag); } else { // U d/du [D_w D5]^dag V = U D5^dag d/du DW^dag Y // implicit adj on U in call Meooe5D(U,Din); this->DhopDerivEO(mat,Din,V,dag); } }; // Tanh template void CayleyFermion5D::SetCoefficientsTanh(Approx::zolotarev_data *zdata,RealD b,RealD c) { std::vector gamma(this->Ls); for(int s=0;sLs;s++) gamma[s] = zdata->gamma[s]; SetCoefficientsInternal(1.0,gamma,b,c); } //Zolo template void CayleyFermion5D::SetCoefficientsZolotarev(RealD zolo_hi,Approx::zolotarev_data *zdata,RealD b,RealD c) { std::vector gamma(this->Ls); for(int s=0;sLs;s++) gamma[s] = zdata->gamma[s]; SetCoefficientsInternal(zolo_hi,gamma,b,c); } //Zolo template void CayleyFermion5D::SetCoefficientsInternal(RealD zolo_hi,std::vector & gamma,RealD b,RealD c) { int Ls=this->Ls; /////////////////////////////////////////////////////////// // The Cayley coeffs (unprec) /////////////////////////////////////////////////////////// assert(gamma.size()==Ls); omega.resize(Ls); bs.resize(Ls); cs.resize(Ls); as.resize(Ls); // // Ts = ( [bs+cs]Dw )^-1 ( (bs+cs) Dw ) // -(g5 ------- -1 ) ( g5 --------- + 1 ) // ( {2+(bs-cs)Dw} ) ( 2+(bs-cs) Dw ) // // bs = 1/2( (1/omega_s + 1)*b + (1/omega - 1)*c ) = 1/2( 1/omega(b+c) + (b-c) ) // cs = 1/2( (1/omega_s - 1)*b + (1/omega + 1)*c ) = 1/2( 1/omega(b+c) - (b-c) ) // // bs+cs = 0.5*( 1/omega(b+c) + (b-c) + 1/omega(b+c) - (b-c) ) = 1/omega(b+c) // bs-cs = 0.5*( 1/omega(b+c) + (b-c) - 1/omega(b+c) + (b-c) ) = b-c // // So // // Ts = ( [b+c]Dw/omega_s )^-1 ( (b+c) Dw /omega_s ) // -(g5 ------- -1 ) ( g5 --------- + 1 ) // ( {2+(b-c)Dw} ) ( 2+(b-c) Dw ) // // Ts = ( [b+c]Dw )^-1 ( (b+c) Dw ) // -(g5 ------- -omega_s) ( g5 --------- + omega_s ) // ( {2+(b-c)Dw} ) ( 2+(b-c) Dw ) // double bpc = b+c; double bmc = b-c; for(int i=0; i < Ls; i++){ as[i] = 1.0; omega[i] = gamma[i]*zolo_hi; //NB reciprocal relative to Chroma NEF code assert(omega[i]!=Coeff_t(0.0)); bs[i] = 0.5*(bpc/omega[i] + bmc); cs[i] = 0.5*(bpc/omega[i] - bmc); } //////////////////////////////////////////////////////// // Constants for the preconditioned matrix Cayley form //////////////////////////////////////////////////////// bee.resize(Ls); cee.resize(Ls); beo.resize(Ls); ceo.resize(Ls); for(int i=0;iM5) +1.0); assert(bee[i]!=Coeff_t(0.0)); cee[i]=as[i]*(1.0-cs[i]*(4.0-this->M5)); beo[i]=as[i]*bs[i]; ceo[i]=-as[i]*cs[i]; } aee.resize(Ls); aeo.resize(Ls); for(int i=0;iMooeeInternalCompute(0,inv,MatpInv,MatmInv); this->MooeeInternalCompute(1,inv,MatpInvDag,MatmInvDag); } template void CayleyFermion5D::MooeeInternalCompute(int dag, int inv, Vector > & Matp, Vector > & Matm) { int Ls=this->Ls; GridBase *grid = this->FermionRedBlackGrid(); int LLs = grid->_rdimensions[0]; if ( LLs == Ls ) { return; // Not vectorised in 5th direction } Eigen::MatrixXcd Pplus = Eigen::MatrixXcd::Zero(Ls,Ls); Eigen::MatrixXcd Pminus = Eigen::MatrixXcd::Zero(Ls,Ls); for(int s=0;s::iscomplex() ) { sp[l] = PplusMat (l*istride+s1*ostride,s2); sm[l] = PminusMat(l*istride+s1*ostride,s2); } else { // if real scalar_type tmp; tmp = PplusMat (l*istride+s1*ostride,s2); sp[l] = scalar_type(tmp.real(),tmp.real()); tmp = PminusMat(l*istride+s1*ostride,s2); sm[l] = scalar_type(tmp.real(),tmp.real()); } } Matp[LLs*s2+s1] = Vp; Matm[LLs*s2+s1] = Vm; }} } FermOpTemplateInstantiate(CayleyFermion5D); GparityFermOpTemplateInstantiate(CayleyFermion5D); }}