/************************************************************************************* Runs the Krylov-Schur algorithm on a (pre-conditioned) domain-wall fermion operator to determine part of its spectrum. Usage : $ ./Example_spec_kryschur Nm = Maximum size of approximation subspace. Nk = Size of truncation subspace maxiter = Maximum number of iterations. Nstop = Stop when Nstop eigenvalues have converged. inFile = Gauge configuration to read in. outDir = Directory to write output to. rf = (Optional) RitzFilter to sort with. Takes in any string in {EvalNormSmall, EvalNormLarge, EvalReSmall, EvalReLarge, EvalImSmall, EvalImLarge} Output: ${outDir}/evals.txt = Contains all eigenvalues. Each line is formatted as `$idx $eval $ritz`, where: - $idx is the index of the eigenvalue. - $eval is the eigenvalue, formated as "(re,im)". - $ritz is the Ritz estimate of the eigenvalue (deviation from being a true eigenvalue) ${outDir}/evec${idx} = Eigenvector $idx written out in SCIDAC format (if LIME is enabled). Grid physics library, www.github.com/paboyle/Grid Source file: ./tests/Test_padded_cell.cc Copyright (C) 2023 Author: Peter Boyle Author: Patrick Oare 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 #include #include #include using namespace std; using namespace Grid; template void writeFile(T& in, std::string const fname){ #ifdef HAVE_LIME // Ref: https://github.com/paboyle/Grid/blob/feature/scidac-wp1/tests/debug/Test_general_coarse_hdcg_phys48.cc#L111 std::cout << Grid::GridLogMessage << "Writes to: " << fname << std::endl; Grid::emptyUserRecord record; Grid::ScidacWriter WR(in.Grid()->IsBoss()); WR.open(fname); WR.writeScidacFieldRecord(in,record,0); // Lexico WR.close(); #endif } /** * Writes the eigensystem of a Krylov Schur object to a directory. * * Parameters * ---------- * std::string path * Directory to write to. */ template void writeEigensystem(KrylovSchur KS, std::string outDir) { int Nk = KS.getNk(); std::cout << GridLogMessage << "Writing output to directory: " << outDir << std::endl; // Write evals std::string evalPath = outDir + "/evals.txt"; std::ofstream fEval; fEval.open(evalPath); Eigen::VectorXcd evals = KS.getEvals(); std::vector ritz = KS.getRitzEstimates(); for (int i = 0; i < Nk; i++) { // write eigenvalues and Ritz estimates fEval << i << " " << evals(i) << " " << ritz[i]; if (i < Nk - 1) { fEval << "\n"; } } fEval.close(); // Write evecs std::vector evecs = KS.getEvecs(); for (int i = 0; i < Nk; i++) { std::string fName = outDir + "/evec" + std::to_string(i); writeFile(evecs[i], fName); // using method from Grid/HMC/ComputeWilsonFlow.cc } } // Hermitize a DWF operator by squaring it template class SquaredLinearOperator : public LinearOperatorBase { public: Matrix &_Mat; public: SquaredLinearOperator(Matrix &Mat): _Mat(Mat) {}; 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 is overloaded as HermOp" << std::endl; HermOp(in, out); } void AdjOp (const Field &in, Field &out){ HermOp(in, out); } void _Op (const Field &in, Field &out){ // std::cout << "Op: M "< 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 ShiftedComplexPVdagMLinearOperator : public LinearOperatorBase { Matrix &_Mat; Matrix &_PV; ComplexD shift; public: ShiftedComplexPVdagMLinearOperator(ComplexD _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 "< std::string NmStr = argv[1]; std::string NkStr = argv[2]; std::string maxIterStr = argv[3]; std::string NstopStr = argv[4]; std::string file = argv[5]; std::string outDir = argv[6]; RitzFilter RF; if (argc == 8) { std::string rf = argv[7]; RF = selectRitzFilter(rf); } else { RF = EvalReSmall; } std::cout << "Sorting eigenvalues using " << rfToString(RF) << std::endl; //const int Ls=16; const int Ls = 8; // GridCartesian * UGrid = SpaceTimeGrid::makeFourDimGrid(GridDefaultLatt(), GridDefaultSimd(Nd,vComplex::Nsimd()),GridDefaultMpi()); //std::vector lat_size {16, 16, 16, 32}; std::vector lat_size {8, 8, 8, 8}; std::cout << "Lattice size: " << lat_size << std::endl; GridCartesian * UGrid = SpaceTimeGrid::makeFourDimGrid(lat_size, GridDefaultSimd(Nd,vComplex::Nsimd()), GridDefaultMpi()); GridRedBlackCartesian * UrbGrid = SpaceTimeGrid::makeFourDimRedBlackGrid(UGrid); GridCartesian * FGrid = SpaceTimeGrid::makeFiveDimGrid(Ls,UGrid); GridRedBlackCartesian * FrbGrid = SpaceTimeGrid::makeFiveDimRedBlackGrid(Ls,UGrid); std::vector seeds4({1,2,3,4}); std::vector seeds5({5,6,7,8}); GridParallelRNG RNG5(FGrid); RNG5.SeedFixedIntegers(seeds5); GridParallelRNG RNG4(UGrid); RNG4.SeedFixedIntegers(seeds4); LatticeFermion src(FGrid); random(RNG5,src); LatticeGaugeField Umu(UGrid); FieldMetaData header; NerscIO::readConfiguration(Umu,header,file); // RealD mass=0.01; RealD mass=0.001; RealD M5=1.8; DomainWallFermionD Ddwf(Umu,*FGrid,*FrbGrid,*UGrid,*UrbGrid,mass,M5); DomainWallFermionD Dpv(Umu,*FGrid,*FrbGrid,*UGrid,*UrbGrid,1.0,M5); std::cout< PVdagM_t; typedef ShiftedPVdagMLinearOperator ShiftedPVdagM_t; typedef ShiftedComplexPVdagMLinearOperator ShiftedComplexPVdagM_t; PVdagM_t PVdagM(Ddwf, Dpv); ShiftedPVdagM_t ShiftedPVdagM(0.1,Ddwf,Dpv); SquaredLinearOperator Dsq (Ddwf); NonHermitianLinearOperator DLinOp (Ddwf); int Nm = std::stoi(NmStr); int Nk = std::stoi(NkStr); int maxIter = std::stoi(maxIterStr); int Nstop = std::stoi(NstopStr); std::cout << GridLogMessage << "Runnning Krylov Schur. Nm = " << Nm << ", Nk = " << Nk << ", maxIter = " << maxIter << ", Nstop = " << Nstop << std::endl; KrylovSchur KrySchur (PVdagM, FGrid, 1e-8, RF); // use preconditioned PV^\dag D_{dwf} // KrylovSchur KrySchur (DLinOp, FGrid, 1e-8, RF); // use D_{dwf} KrySchur(src, maxIter, Nm, Nk, Nstop); std::cout<