/************************************************************************************* Grid physics library, www.github.com/paboyle/Grid Source file: ./tests/Test_padded_cell.cc Copyright (C) 2023 Author: Peter Boyle 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 #include #include #include using namespace std; using namespace Grid; template class PVdagMLinearOperator : public LinearOperatorBase { Matrix &_Mat; Matrix &_PV; public: PVdagMLinearOperator(Matrix &Mat,Matrix &PV): _Mat(Mat),_PV(PV){}; void OpDiag (const Field &in, Field &out) { assert(0); } void OpDir (const Field &in, Field &out,int dir,int disp) { assert(0); } void OpDirAll (const Field &in, std::vector &out){ assert(0); }; void Op (const Field &in, Field &out){ // std::cout << GridLogMessage<< "Op: PVdag M "< class MdagPVLinearOperator : public LinearOperatorBase { Matrix &_Mat; Matrix &_PV; public: MdagPVLinearOperator(Matrix &Mat,Matrix &PV): _Mat(Mat),_PV(PV){}; void OpDiag (const Field &in, Field &out) { assert(0); } void OpDir (const Field &in, Field &out,int dir,int disp) { assert(0); } void OpDirAll (const Field &in, std::vector &out){ assert(0); }; void Op (const Field &in, Field &out){ Field tmp(in.Grid()); // std::cout < class ShiftedPVdagMLinearOperator : public LinearOperatorBase { Matrix &_Mat; Matrix &_PV; RealD shift; public: ShiftedPVdagMLinearOperator(RealD _shift,Matrix &Mat,Matrix &PV): shift(_shift),_Mat(Mat),_PV(PV){}; void OpDiag (const Field &in, Field &out) { assert(0); } void OpDir (const Field &in, Field &out,int dir,int disp) { assert(0); } void OpDirAll (const Field &in, std::vector &out){ assert(0); }; void Op (const Field &in, Field &out){ // std::cout << "Op: PVdag M "< class MGPreconditionerSVD : public LinearFunction< Lattice > { public: using LinearFunction >::operator(); typedef Aggregation Aggregates; typedef typename Aggregation::FineField FineField; typedef typename Aggregation::CoarseVector CoarseVector; typedef typename Aggregation::CoarseMatrix CoarseMatrix; typedef LinearOperatorBase FineOperator; typedef LinearFunction FineSmoother; typedef LinearOperatorBase CoarseOperator; typedef LinearFunction CoarseSolver; Aggregates & _FineToCoarse; Aggregates & _CoarseToFine; FineOperator & _FineOperator; FineSmoother & _PreSmoother; FineSmoother & _PostSmoother; CoarseOperator & _CoarseOperator; CoarseSolver & _CoarseSolve; int level; void Level(int lv) {level = lv; }; MGPreconditionerSVD(Aggregates &FtoC, Aggregates &CtoF, FineOperator &Fine, FineSmoother &PreSmoother, FineSmoother &PostSmoother, CoarseOperator &CoarseOperator_, CoarseSolver &CoarseSolve_) : _FineToCoarse(FtoC), _CoarseToFine(CtoF), _FineOperator(Fine), _PreSmoother(PreSmoother), _PostSmoother(PostSmoother), _CoarseOperator(CoarseOperator_), _CoarseSolve(CoarseSolve_), level(1) { } virtual void operator()(const FineField &in, FineField & out) { GridBase *CoarseGrid = _FineToCoarse.CoarseGrid; // auto CoarseGrid = _CoarseOperator.Grid(); CoarseVector Csrc(CoarseGrid); CoarseVector Csol(CoarseGrid); FineField vec1(in.Grid()); FineField vec2(in.Grid()); std::cout< seeds4({1,2,3,4}); std::vector seeds5({5,6,7,8}); std::vector cseeds({5,6,7,8}); GridParallelRNG RNG5(FGrid); RNG5.SeedFixedIntegers(seeds5); GridParallelRNG RNG4(UGrid); RNG4.SeedFixedIntegers(seeds4); GridParallelRNG CRNG(Coarse5d);CRNG.SeedFixedIntegers(cseeds); LatticeFermion src(FGrid); random(RNG5,src); LatticeFermion result(FGrid); result=Zero(); LatticeFermion ref(FGrid); ref=Zero(); LatticeFermion tmp(FGrid); LatticeFermion err(FGrid); LatticeGaugeField Umu(UGrid); FieldMetaData header; std::string file("ckpoint_lat.4000"); NerscIO::readConfiguration(Umu,header,file); RealD mass=0.01; RealD M5=1.8; DomainWallFermionD Ddwf(Umu,*FGrid,*FrbGrid,*UGrid,*UrbGrid,mass,M5); DomainWallFermionD Dpv(Umu,*FGrid,*FrbGrid,*UGrid,*UrbGrid,1.0,M5); const int nbasis = 30; const int cb = 0 ; NextToNearestStencilGeometry5D geom(Coarse5d); std::cout< PVdagM_t; typedef MdagPVLinearOperator MdagPV_t; typedef ShiftedPVdagMLinearOperator ShiftedPVdagM_t; PVdagM_t PVdagM(Ddwf,Dpv); MdagPV_t MdagPV(Ddwf,Dpv); // ShiftedPVdagM_t ShiftedPVdagM(2.0,Ddwf,Dpv); // 355 // ShiftedPVdagM_t ShiftedPVdagM(1.0,Ddwf,Dpv); // 246 // ShiftedPVdagM_t ShiftedPVdagM(0.5,Ddwf,Dpv); // 183 // ShiftedPVdagM_t ShiftedPVdagM(0.25,Ddwf,Dpv); // 145 // ShiftedPVdagM_t ShiftedPVdagM(0.1,Ddwf,Dpv); // 134 // ShiftedPVdagM_t ShiftedPVdagM(0.1,Ddwf,Dpv); // 127 -- NULL space via inverse iteration // ShiftedPVdagM_t ShiftedPVdagM(0.1,Ddwf,Dpv); // 57 -- NULL space via inverse iteration; 3 iterations // ShiftedPVdagM_t ShiftedPVdagM(0.25,Ddwf,Dpv); // 57 , tighter inversion // ShiftedPVdagM_t ShiftedPVdagM(0.25,Ddwf,Dpv); // nbasis 20 -- 49 iters // ShiftedPVdagM_t ShiftedPVdagM(0.25,Ddwf,Dpv); // nbasis 20 -- 70 iters; asymmetric // ShiftedPVdagM_t ShiftedPVdagM(0.25,Ddwf,Dpv); // 58; Loosen coarse, tighten fine // ShiftedPVdagM_t ShiftedPVdagM(0.1,Ddwf,Dpv); // 56 ... // ShiftedPVdagM_t ShiftedPVdagM(0.1,Ddwf,Dpv); // 51 ... with 24 vecs // ShiftedPVdagM_t ShiftedPVdagM(0.1,Ddwf,Dpv); // 31 ... with 24 vecs and 2^4 blocking // ShiftedPVdagM_t ShiftedPVdagM(0.1,Ddwf,Dpv); // 43 ... with 16 vecs and 2^4 blocking, sloppier // ShiftedPVdagM_t ShiftedPVdagM(0.1,Ddwf,Dpv); // 35 ... with 20 vecs and 2^4 blocking // ShiftedPVdagM_t ShiftedPVdagM(0.1,Ddwf,Dpv); // 35 ... with 20 vecs and 2^4 blocking, looser coarse // ShiftedPVdagM_t ShiftedPVdagM(0.1,Ddwf,Dpv); // 64 ... with 20 vecs, Christoph setup, and 2^4 blocking, looser coarse ShiftedPVdagM_t ShiftedPVdagM(0.01,Ddwf,Dpv); // // Run power method on HOA?? PowerMethod PM; // PM(PVdagM,src); // PM(MdagPV,src); // Warning: This routine calls PVdagM.Op, not PVdagM.HermOp typedef Aggregation Subspace; Subspace V(Coarse5d,FGrid,cb); Subspace U(Coarse5d,FGrid,cb); // Breeds right singular vectors with call to HermOp (V) V.CreateSubspaceChebyshev(RNG5,PVdagM, nbasis, 4000.0,0.003, 500); // Breeds left singular vectors with call to HermOp (U) // U.CreateSubspaceChebyshev(RNG5,MdagPV, U.CreateSubspaceChebyshev(RNG5,PVdagM, nbasis, 4000.0,0.003, 500); typedef Aggregation CombinedSubspace; CombinedSubspace CombinedUV(Coarse5d,FGrid,cb); for(int b=0;b " < " < " < " < LittleDiracOperatorV; typedef LittleDiracOperatorV::CoarseVector CoarseVectorV; typedef GeneralCoarsenedMatrix LittleDiracOperator; typedef LittleDiracOperator::CoarseVector CoarseVector; V.Orthogonalise(); for(int b =0 ; b simple; NonHermitianLinearOperator LinOpCoarse(LittleDiracOpPV); // PrecGeneralisedConjugateResidualNonHermitian L2PGCR(1.0e-4, 100, LinOpCoarse,simple,10,10); PrecGeneralisedConjugateResidualNonHermitian L2PGCR(1.0e-2, 10, LinOpCoarse,simple,20,20); L2PGCR.Level(3); c_res=Zero(); L2PGCR(c_src,c_res); //////////////////////////////////////// // Fine grid smoother //////////////////////////////////////// std::cout< simple_fine; // NonHermitianLinearOperator LinOpSmooth(PVdagM); PrecGeneralisedConjugateResidualNonHermitian SmootherGCR(0.01,1,ShiftedPVdagM,simple_fine,16,16); SmootherGCR.Level(2); LatticeFermionD f_src(FGrid); LatticeFermionD f_res(FGrid); f_src = one; // 1 in every element for vector 1. f_res=Zero(); SmootherGCR(f_src,f_res); typedef MGPreconditionerSVD TwoLevelMG; TwoLevelMG TwoLevelPrecon(CombinedUV,CombinedUV, PVdagM, simple_fine, SmootherGCR, LinOpCoarse, L2PGCR); PrecGeneralisedConjugateResidualNonHermitian L1PGCR(1.0e-8,1000,PVdagM,TwoLevelPrecon,20,20); L1PGCR.Level(1); f_res=Zero(); L1PGCR(f_src,f_res); std::cout<