/************************************************************************************* 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 << "Op: PVdag M "< 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 MGPreconditioner : 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 & _Aggregates; FineOperator & _FineOperator; FineSmoother & _PreSmoother; FineSmoother & _PostSmoother; CoarseOperator & _CoarseOperator; CoarseSolver & _CoarseSolve; int level; void Level(int lv) {level = lv; }; MGPreconditioner(Aggregates &Agg, FineOperator &Fine, FineSmoother &PreSmoother, FineSmoother &PostSmoother, CoarseOperator &CoarseOperator_, CoarseSolver &CoarseSolve_) : _Aggregates(Agg), _FineOperator(Fine), _PreSmoother(PreSmoother), _PostSmoother(PostSmoother), _CoarseOperator(CoarseOperator_), _CoarseSolve(CoarseSolve_), level(1) { } virtual void operator()(const FineField &in, FineField & out) { GridBase *CoarseGrid = _Aggregates.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 = 8; const int cb = 0 ; LatticeFermion prom(FGrid); typedef GeneralCoarsenedMatrix LittleDiracOperator; typedef LittleDiracOperator::CoarseVector CoarseVector; NextToNearestStencilGeometry5D geom(Coarse5d); std::cout< PVdagM_t; typedef ShiftedPVdagMLinearOperator ShiftedPVdagM_t; PVdagM_t PVdagM(Ddwf,Dpv); ShiftedPVdagM_t ShiftedPVdagM(2.0,Ddwf,Dpv); // Run power method on HOA?? // PowerMethod PM; PM(PVdagM,src); // Warning: This routine calls PVdagM.Op, not PVdagM.HermOp typedef Aggregation Subspace; Subspace AggregatesPD(Coarse5d,FGrid,cb); AggregatesPD.CreateSubspaceChebyshev(RNG5, PVdagM, nbasis, 4000.0, 2.0, 200, 200, 200, 0.0); LittleDiracOperator LittleDiracOpPV(geom,FGrid,Coarse5d); LittleDiracOpPV.CoarsenOperator(PVdagM,AggregatesPD); std::cout< subspace(nbasis,FGrid); subspace=AggregatesPD.subspace; Complex one(1.0); c_src = one; // 1 in every element for vector 1. blockPromote(c_src,err,subspace); prom=Zero(); for(int b=0;b simple; NonHermitianLinearOperator LinOpCoarse(LittleDiracOpPV); PrecGeneralisedConjugateResidualNonHermitian L2PGCR(1.0e-8, 100, LinOpCoarse,simple,10,10); L2PGCR.Level(2); c_res=Zero(); L2PGCR(c_src,c_res); //////////////////////////////////////// // Fine grid smoother //////////////////////////////////////// std::cout< simple_fine; // NonHermitianLinearOperator LinOpSmooth(PVdagM); PrecGeneralisedConjugateResidualNonHermitian SmootherGCR(0.01,10,ShiftedPVdagM,simple_fine,4,4); 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 MGPreconditioner TwoLevelMG; TwoLevelMG TwoLevelPrecon(AggregatesPD, PVdagM, SmootherGCR, SmootherGCR, LinOpCoarse, L2PGCR); PrecGeneralisedConjugateResidualNonHermitian L1PGCR(1.0e-8,1000,PVdagM,TwoLevelPrecon,8,8); L1PGCR.Level(1); f_res=Zero(); L1PGCR(f_src,f_res); std::cout<