CHROMA
Classes | Typedefs | Functions

Classes

class  Chroma::LinAlg::Vectors< T >
 Hold vectors. More...
 
class  Chroma::LinAlg::RitzPairs< T >
 Holds eigenvalues and eigenvectors. More...
 
class  Chroma::LinAlg::Matrix< T >
 This is a square matrix. More...
 
class  Chroma::LinAlg::VectorArrays< T >
 Hold vectors. More...
 
class  Chroma::LinAlg::RitzPairsArray< T >
 Holds eigenvalues and arrays of eigenvectors for use in 5D work. More...
 
class  Chroma::LinOpSysSolverMDWFArray
 AVP's DWF Solver interface. More...
 
struct  Chroma::MultiSysSolverCGParams
 Params for CG inverter. More...
 
class  Chroma::LinOpMultiSysSolverMR< T >
 Solve a MR system. Here, the operator is NOT assumed to be hermitian. More...
 
class  Chroma::MdagMMultiSysSolverCG< T >
 Solve a CG2 system. Here, the operator is NOT assumed to be hermitian. More...
 
class  Chroma::MdagMMultiSysSolverCGAccumulate< T >
 Solve a CG2 system. Here, the operator is NOT assumed to be hermitian. More...
 
class  Chroma::MdagMMultiSysSolverCGAccumulateArray< T >
 Solve a CG system. Here, the operator is NOT assumed to be hermitian. More...
 
class  Chroma::MdagMMultiSysSolverCGArray< T >
 Solve a CG system. Here, the operator is NOT assumed to be hermitian. More...
 
class  Chroma::MdagMMultiSysSolverCGChronoClover
 Solve a CG2 system. Here, the operator is NOT assumed to be hermitian. More...
 
struct  Chroma::MultiSysSolverMRParams
 Params for MR inverter. More...
 
class  Chroma::LinOpSysSolverQOPMG
 Solve a M*psi=chi linear system using the external QDP multigrid inverter. More...
 
class  Chroma::MdagMSysSolverQOPMG
 Solve a M*psi=chi linear system using the external QDP multigrid inverter. More...
 
struct  Chroma::SysSolverQOPMGParams
 Parameters for the external QDP multigrid inverter. More...
 
struct  Chroma::SysSolverBiCGStabParams
 Params for BiCGStab inverter. More...
 
struct  Chroma::SysSolverCGParams
 Params for CG inverter. More...
 
struct  Chroma::SysSolverEigCGParams
 Params for EigCG inverter. More...
 
struct  Chroma::SysSolverFGMRESDRParams
 Params for FGMRESDR inverter. More...
 
class  Chroma::LinOpSysSolverBiCGStab< T >
 Solve a M*psi=chi linear system by BICGSTAB. More...
 
class  Chroma::LinOpSysSolverBiCRStab< T >
 Solve a M*psi=chi linear system by BICGSTAB. More...
 
class  Chroma::LinOpSysSolverCG< T >
 Solve a M*psi=chi linear system by CG2. More...
 
class  Chroma::LinOpSysSolverCGArray< T >
 Solve a M*psi=chi linear system by CG2. More...
 
class  Chroma::LinOpSysSolverCGTiming< T >
 Solve a M*psi=chi linear system by CG2. More...
 
class  Chroma::LinOpSysSolverEigCG< T >
 Solve a M*psi=chi linear system by EigCG with eigenvectors. More...
 
class  Chroma::LinOpSysSolverEigCGArray< T >
 Solve a M*psi=chi linear system by CG2 with eigenvectors. More...
 
class  Chroma::LinOpSysSolverFGMRESDR
 Solve a M*psi=chi linear system by FGMRESDR. More...
 
class  Chroma::LinOpSysSolverIBiCGStab< T >
 Solve a M*psi=chi linear system by IBICGSTAB. More...
 
class  Chroma::LinOpSysSolverMR< T >
 Solve a M*psi=chi linear system by MR. More...
 
class  Chroma::LinOpSysSolverOptEigBiCG< T >
 Solve a M*psi=chi linear system by biCG with eigenvectors. More...
 
class  Chroma::LinOpSysSolverReliableBiCGStabClover
 Solve a system using Richardson iteration. More...
 
class  Chroma::LinOpSysSolverReliableCGClover
 Solve a system using Richardson iteration. More...
 
class  Chroma::LinOpSysSolverReliableIBiCGStabClover
 Solve a system using Richardson iteration. More...
 
class  Chroma::LinOpSysSolverRichardsonClover
 Solve a system using Richardson iteration. More...
 
class  Chroma::MdagMSysSolverBiCGStab< T >
 Solve a BiCGStab system. Here, the operator is NOT assumed to be hermitian. More...
 
class  Chroma::MdagMSysSolverCG< T >
 Solve a CG2 system. Here, the operator is NOT assumed to be hermitian. More...
 
class  Chroma::MdagMSysSolverCGArray< T >
 Solve a CG2 system. Here, the operator is NOT assumed to be hermitian. More...
 
class  Chroma::MdagMSysSolverCGLFClover
 Solve a system using CG iteration. More...
 
class  Chroma::MdagMSysSolverCGTimings< T >
 Solve a CG2 system. Here, the operator is NOT assumed to be hermitian. More...
 
class  Chroma::MdagMSysSolverQDPEigCG< T >
 Solve a M*psi=chi linear system by CG2 with eigenvectors. More...
 
class  Chroma::MdagMSysSolverIBiCGStab< T >
 Solve a IBiCGStab system. Here, the operator is NOT assumed to be hermitian. More...
 
class  Chroma::MdagMSysSolverMR< T >
 Solve a MR system. Here, the operator is NOT assumed to be hermitian. More...
 
class  Chroma::MdagMSysSolverOptEigCG< T >
 Solve a M*psi=chi linear system by CG2 with eigenvectors. More...
 
class  Chroma::MdagMSysSolverReliableBiCGStabClover
 Solve a system using Richardson iteration. More...
 
class  Chroma::MdagMSysSolverReliableCGClover
 Solve a system using Richardson iteration. More...
 
class  Chroma::MdagMSysSolverReliableIBiCGStabClover
 Solve a system using Richardson iteration. More...
 
class  Chroma::MdagMSysSolverRichardsonClover
 Solve a system using Richardson iteration. More...
 
struct  Chroma::SysSolverMRParams
 Params for MR inverter. More...
 
struct  Chroma::SysSolverOptEigBiCGParams
 Params for EigBiCG inverter. More...
 
struct  Chroma::SysSolverOptEigCGParams
 Params for EigCG inverter. More...
 
class  Chroma::PolyPrecSysSolverCG< T >
 Solve a PolyPrec*psi=chi linear system by CG1. More...
 

Typedefs

typedef SingletonHolder< ObjectFactory< LinOpMultiSystemSolver< LatticeFermion >, std::string, TYPELIST_3(XMLReader &, const std::string &, Handle< LinearOperator< LatticeFermion > >), LinOpMultiSystemSolver< LatticeFermion > *(*)(XMLReader &, const std::string &, Handle< LinearOperator< LatticeFermion > >), StringFactoryError > > Chroma::TheLinOpFermMultiSystemSolverFactory
 LinOp system solver factory (foundry) More...
 
typedef SingletonHolder< ObjectFactory< LinOpMultiSystemSolver< LatticeStaggeredFermion >, std::string, TYPELIST_3(XMLReader &, const std::string &, Handle< LinearOperator< LatticeStaggeredFermion > >), LinOpMultiSystemSolver< LatticeStaggeredFermion > *(*)(XMLReader &, const std::string &, Handle< LinearOperator< LatticeStaggeredFermion > >), StringFactoryError > > Chroma::TheLinOpStagFermMultiSystemSolverFactory
 LinOp system solver factory (foundry) More...
 
typedef SingletonHolder< ObjectFactory< MdagMMultiSystemSolverAccumulate< LatticeFermion >, std::string, TYPELIST_3(XMLReader &, const std::string &, Handle< LinearOperator< LatticeFermion > >), MdagMMultiSystemSolverAccumulate< LatticeFermion > *(*)(XMLReader &, const std::string &, Handle< LinearOperator< LatticeFermion > >), StringFactoryError > > Chroma::TheMdagMFermMultiSystemSolverAccumulateFactory
 MdagM system solver factory (foundry) More...
 
typedef SingletonHolder< ObjectFactory< MdagMMultiSystemSolverAccumulateArray< LatticeFermion >, std::string, TYPELIST_3(XMLReader &, const std::string &, Handle< LinearOperatorArray< LatticeFermion > >), MdagMMultiSystemSolverAccumulateArray< LatticeFermion > *(*)(XMLReader &, const std::string &, Handle< LinearOperatorArray< LatticeFermion > >), StringFactoryError > > Chroma::TheMdagMFermMultiSystemSolverAccumulateArrayFactory
 MdagM system solver factory (foundry) More...
 
typedef SingletonHolder< ObjectFactory< MdagMMultiSystemSolverAccumulate< LatticeStaggeredFermion >, std::string, TYPELIST_3(XMLReader &, const std::string &, Handle< LinearOperator< LatticeStaggeredFermion > >), MdagMMultiSystemSolverAccumulate< LatticeStaggeredFermion > *(*)(XMLReader &, const std::string &, Handle< LinearOperator< LatticeStaggeredFermion > >), StringFactoryError > > Chroma::TheMdagMStagFermMultiSystemSolverAccumulateFactory
 MdagM system solver factory (foundry) More...
 
typedef SingletonHolder< ObjectFactory< MdagMMultiSystemSolver< LatticeFermion >, std::string, TYPELIST_4(XMLReader &, const std::string &, FSHandle, Handle< LinearOperator< LatticeFermion > >), MdagMMultiSystemSolver< LatticeFermion > *(*)(XMLReader &, const std::string &, FSHandle, Handle< LinearOperator< LatticeFermion > >), StringFactoryError > > Chroma::TheMdagMFermMultiSystemSolverFactory
 MdagM system solver factory (foundry) More...
 
typedef SingletonHolder< ObjectFactory< MdagMMultiSystemSolverArray< LatticeFermion >, std::string, TYPELIST_4(XMLReader &, const std::string &, FSHandle, Handle< LinearOperatorArray< LatticeFermion > >), MdagMMultiSystemSolverArray< LatticeFermion > *(*)(XMLReader &, const std::string &, FSHandle, Handle< LinearOperatorArray< LatticeFermion > >), StringFactoryError > > Chroma::TheMdagMFermMultiSystemSolverArrayFactory
 MdagM system solver factory (foundry) More...
 
typedef SingletonHolder< ObjectFactory< MdagMMultiSystemSolver< LatticeStaggeredFermion >, std::string, TYPELIST_3(XMLReader &, const std::string &, Handle< LinearOperator< LatticeStaggeredFermion > >), MdagMMultiSystemSolver< LatticeStaggeredFermion > *(*)(XMLReader &, const std::string &, Handle< LinearOperator< LatticeStaggeredFermion > >), StringFactoryError > > Chroma::TheMdagMStagFermMultiSystemSolverFactory
 MdagM system solver factory (foundry) More...
 
typedef SingletonHolder< ObjectFactory< LinOpSystemSolver< LatticeFermion >, std::string, TYPELIST_4(XMLReader &, const std::string &, FSHandle, Handle< LinearOperator< LatticeFermion > >), LinOpSystemSolver< LatticeFermion > *(*)(XMLReader &, const std::string &, FSHandle, Handle< LinearOperator< LatticeFermion > >), StringFactoryError > > Chroma::TheLinOpFermSystemSolverFactory
 LinOp system solver factory (foundry) More...
 
typedef SingletonHolder< ObjectFactory< LinOpSystemSolverArray< LatticeFermion >, std::string, TYPELIST_4(XMLReader &, const std::string &, FSHandle, Handle< LinearOperatorArray< LatticeFermion > >), LinOpSystemSolverArray< LatticeFermion > *(*)(XMLReader &, const std::string &, FSHandle, Handle< LinearOperatorArray< LatticeFermion > >), StringFactoryError > > Chroma::TheLinOpFermSystemSolverArrayFactory
 LinOp system solver factory (foundry) More...
 
typedef SingletonHolder< ObjectFactory< LinOpSystemSolver< LatticeStaggeredFermion >, std::string, TYPELIST_3(XMLReader &, const std::string &, Handle< LinearOperator< LatticeStaggeredFermion > >), LinOpSystemSolver< LatticeStaggeredFermion > *(*)(XMLReader &, const std::string &, Handle< LinearOperator< LatticeStaggeredFermion > >), StringFactoryError > > Chroma::TheLinOpStagFermSystemSolverFactory
 LinOp system solver factory (foundry) More...
 
typedef SingletonHolder< ObjectFactory< MdagMSystemSolver< LatticeFermion >, std::string, TYPELIST_4(XMLReader &, const std::string &, FactoryEnv::FSHandle, Handle< LinearOperator< LatticeFermion > >), MdagMSystemSolver< LatticeFermion > *(*)(XMLReader &, const std::string &, FactoryEnv::FSHandle, Handle< LinearOperator< LatticeFermion > >), StringFactoryError > > Chroma::TheMdagMFermSystemSolverFactory
 MdagM system solver factory (foundry) More...
 
typedef SingletonHolder< ObjectFactory< MdagMSystemSolverArray< LatticeFermion >, std::string, TYPELIST_4(XMLReader &, const std::string &, FactoryEnv::FSHandle, Handle< LinearOperatorArray< LatticeFermion > >), MdagMSystemSolverArray< LatticeFermion > *(*)(XMLReader &, const std::string &, FactoryEnv::FSHandle, Handle< LinearOperatorArray< LatticeFermion > >), StringFactoryError > > Chroma::TheMdagMFermSystemSolverArrayFactory
 MdagM system solver factory (foundry) More...
 
typedef SingletonHolder< ObjectFactory< MdagMSystemSolver< LatticeStaggeredFermion >, std::string, TYPELIST_3(XMLReader &, const std::string &, Handle< LinearOperator< LatticeStaggeredFermion > >), MdagMSystemSolver< LatticeStaggeredFermion > *(*)(XMLReader &, const std::string &, Handle< LinearOperator< LatticeStaggeredFermion > >), StringFactoryError > > Chroma::TheMdagMStagFermSystemSolverFactory
 MdagM system solver factory (foundry) More...
 
typedef SingletonHolder< ObjectFactory< PolyPrecSystemSolver< LatticeFermion >, std::string, TYPELIST_3(XMLReader &, const std::string &, Handle< LinearOperator< LatticeFermion > >), PolyPrecSystemSolver< LatticeFermion > *(*)(XMLReader &, const std::string &, Handle< LinearOperator< LatticeFermion > >), StringFactoryError > > Chroma::ThePolyPrecFermSystemSolverFactory
 PolyPrec system solver factory (foundry) More...
 

Functions

template<typename T >
void Chroma::InvRelCG1_a (const LinearOperator< T > &A, const T &chi, T &psi, const Real &RsdCG, int MaxCG, int &n_count)
 Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator. More...
 
template<>
void Chroma::InvRelCG1 (const LinearOperator< T > &A, const T &chi, T &psi, const Real &RsdCG, int MaxCG, int &n_count)
 Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator. More...
 
template<typename T >
void Chroma::InvRelCG2_a (const LinearOperator< T > &M, const T &chi, T &psi, const Real &RsdCG, int MaxCG, int &n_count)
 Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator. More...
 
template<>
void Chroma::InvRelCG2 (const LinearOperator< T > &M, const T &chi, T &psi, const Real &RsdCG, int MaxCG, int &n_count)
 Relaxed Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator. More...
 
template<typename T >
SystemSolverResults_t Chroma::InvCG1_a (const LinearOperator< T > &A, const T &chi, T &psi, const Real &RsdCG, int MaxCG, int MinCG=0)
 Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator. More...
 
template<>
SystemSolverResults_t Chroma::InvCG1 (const LinearOperator< T > &A, const T &chi, T &psi, const Real &RsdCG, int MaxCG, int MinCG=0)
 Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator. More...
 
template<typename T >
void Chroma::InvCG1_a (const LinearOperatorArray< T > &A, const multi1d< T > &chi, multi1d< T > &psi, const Real &RsdCG, int MaxCG, int &n_count)
 Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator. More...
 
template<>
void Chroma::InvCG1 (const LinearOperatorArray< T > &A, const multi1d< T > &chi, multi1d< T > &psi, const Real &RsdCG, int MaxCG, int &n_count)
 Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator. More...
 
template<typename T , typename RT >
SystemSolverResults_t Chroma::InvCG2_a (const LinearOperator< T > &M, const T &chi, T &psi, const RT &RsdCG, int MaxCG)
 Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator. More...
 
template<typename T , typename C >
SystemSolverResults_t Chroma::InvCG2_a (const C &M, const multi1d< T > &chi, multi1d< T > &psi, const Real &RsdCG, int MaxCG)
 Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator. More...
 
template<typename T >
SystemSolverResults_t Chroma::InvCG2_timings_a (const LinearOperator< T > &M, const T &chi, T &psi, int MaxCG)
 Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator. More...
 
template<>
SystemSolverResults_t Chroma::InvCG2_timings (const LinearOperator< T > &M, const T &chi, T &psi, int MaxCG)
 Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator. More...
 
template<typename T >
void Chroma::MInvRelCG_a (const LinearOperator< T > &A, const T &chi, multi1d< T > &psi, const multi1d< Real > &shifts, const multi1d< Real > &RsdCG, int MaxCG, int &n_count)
 Multishift Conjugate-Gradient (CG1) algorithm for a Linear Operator. More...
 
template<typename T >
void Chroma::MInvCG_a (const LinearOperator< T > &A, const T &chi, multi1d< T > &psi, const multi1d< Real > &shifts, const multi1d< Real > &RsdCG, int MaxCG, int &n_count)
 Multishift Conjugate-Gradient (CG1) algorithm for a Linear Operator. More...
 
template<>
void Chroma::MInvCG (const LinearOperator< LatticeFermion > &M, const LatticeFermion &chi, multi1d< LatticeFermion > &psi, const multi1d< Real > &shifts, const multi1d< Real > &RsdCG, int MaxCG, int &n_count)
 
template<>
void Chroma::MInvCG (const DiffLinearOperator< LatticeFermion, multi1d< LatticeColorMatrix >, multi1d< LatticeColorMatrix > > &M, const LatticeFermion &chi, multi1d< LatticeFermion > &psi, const multi1d< Real > &shifts, const multi1d< Real > &RsdCG, int MaxCG, int &n_count)
 
template<typename T , typename R >
void Chroma::MInvCG2_a (const LinearOperator< T > &M, const T &chi, multi1d< T > &psi, const multi1d< R > &shifts, const multi1d< R > &RsdCG, int MaxCG, int &n_count)
 Multishift Conjugate-Gradient (CG1) algorithm for a Linear Operator. More...
 
void Chroma::MInvCG2 (const LinearOperator< LatticeFermionF > &M, const LatticeFermionF &chi, multi1d< LatticeFermionF > &psi, const multi1d< RealF > &shifts, const multi1d< RealF > &RsdCG, int MaxCG, int &n_count)
 
void Chroma::MInvCG2 (const LinearOperator< LatticeFermionD > &M, const LatticeFermionD &chi, multi1d< LatticeFermionD > &psi, const multi1d< RealD > &shifts, const multi1d< RealD > &RsdCG, int MaxCG, int &n_count)
 
void Chroma::MInvCG2 (const DiffLinearOperator< LatticeFermionF, multi1d< LatticeColorMatrixF >, multi1d< LatticeColorMatrixF > > &M, const LatticeFermionF &chi, multi1d< LatticeFermionF > &psi, const multi1d< RealF > &shifts, const multi1d< RealF > &RsdCG, int MaxCG, int &n_count)
 
template<typename T >
void Chroma::MInvCG2 (const LinearOperator< T > &M, const T &chi, multi1d< T > &psi, const multi1d< Real > &shifts, const multi1d< Real > &RsdCG, int MaxCG, int &n_count)
 
template<typename T , typename P , typename Q >
void Chroma::MInvCG2 (const DiffLinearOperator< T, P, Q > &M, const T &chi, multi1d< T > &psi, const multi1d< Real > &shifts, const multi1d< Real > &RsdCG, int MaxCG, int &n_count)
 
template<typename T >
void Chroma::MInvCG2Accum_a (const LinearOperator< T > &M, const T &chi, T &psi, const Real &norm, const multi1d< Real > &residues, const multi1d< Real > &poles, const Real &RsdCG, int MaxCG, int &n_count)
 Multishift Conjugate-Gradient (CG1) algorithm for a Linear Operator. More...
 
template<>
void Chroma::MInvCG2Accum (const LinearOperator< LatticeFermion > &M, const LatticeFermion &chi, LatticeFermion &psi, const Real &norm, const multi1d< Real > &residues, const multi1d< Real > &poles, const Real &RsdCG, int MaxCG, int &n_count)
 
template<>
void Chroma::MInvCG2Accum (const DiffLinearOperator< LatticeFermion, multi1d< LatticeColorMatrix >, multi1d< LatticeColorMatrix > > &M, const LatticeFermion &chi, LatticeFermion &psi, const Real &norm, const multi1d< Real > &residues, const multi1d< Real > &poles, const Real &RsdCG, int MaxCG, int &n_count)
 
template<typename T , typename C >
void Chroma::MInvCGAccum_a (const C &A, const multi1d< T > &chi, multi1d< T > &psi, const Real &norm, const multi1d< Real > &residues, const multi1d< Real > &poles, const Real &RsdCG, const int MaxCG, int &n_count)
 Multishift Conjugate-Gradient (CG1) algorithm for a Linear Operator. More...
 
template<>
void Chroma::MInvCGAccum (const LinearOperatorArray< LatticeFermion > &M, const multi1d< LatticeFermion > &chi, multi1d< LatticeFermion > &psi, const Real &norm, const multi1d< Real > &residues, const multi1d< Real > &poles, const Real &RsdCG, const int MaxCG, int &n_count)
 
template<>
void Chroma::MInvCGAccum (const DiffLinearOperatorArray< LatticeFermion, multi1d< LatticeColorMatrix >, multi1d< LatticeColorMatrix > > &M, const multi1d< LatticeFermion > &chi, multi1d< LatticeFermion > &psi, const Real &norm, const multi1d< Real > &residues, const multi1d< Real > &poles, const Real &RsdCG, const int MaxCG, int &n_count)
 
template<typename T , typename C >
void Chroma::MInvCG_a (const C &A, const multi1d< T > &chi, multi1d< multi1d< T > > &psi, const multi1d< Real > &shifts, const multi1d< Real > &RsdCG, int MaxCG, int &n_count)
 Multishift Conjugate-Gradient (CG1) algorithm for a Linear Operator. More...
 
template<>
void Chroma::MInvCG (const LinearOperatorArray< LatticeFermion > &M, const multi1d< LatticeFermion > &chi, multi1d< multi1d< LatticeFermion > > &psi, const multi1d< Real > &shifts, const multi1d< Real > &RsdCG, int MaxCG, int &n_count)
 
template<>
void Chroma::MInvCG (const DiffLinearOperatorArray< LatticeFermion, multi1d< LatticeColorMatrix >, multi1d< LatticeColorMatrix > > &M, const multi1d< LatticeFermion > &chi, multi1d< multi1d< LatticeFermion > > &psi, const multi1d< Real > &shifts, const multi1d< Real > &RsdCG, int MaxCG, int &n_count)
 
void Chroma::read (XMLReader &xml, const std::string &path, MultiSysSolverCGParams &param)
 
void Chroma::write (XMLWriter &xml, const std::string &path, const MultiSysSolverCGParams &param)
 
void Chroma::read (XMLReader &xml, const std::string &path, MultiSysSolverMRParams &param)
 
void Chroma::write (XMLWriter &xml, const std::string &path, const MultiSysSolverMRParams &param)
 
template<typename T >
void Chroma::normGramSchmidt_T (multi1d< T > &vec, int f, int t, const Subset &sub)
 Gram-Schmidt with normalization. More...
 
template<typename T >
void Chroma::normGramSchmidtArray_T (multi2d< T > &vec, int f, int t, const Subset &sub)
 Gram-Schmidt with normalization. More...
 
void Chroma::read (XMLReader &xml, const std::string &path, SysSolverQOPMGParams &param)
 
void Chroma::write (XMLWriter &xml, const std::string &path, const SysSolverQOPMGParams &param)
 
void Chroma::read (XMLReader &xml, const std::string &path, SysSolverBiCGStabParams &param)
 
void Chroma::write (XMLWriter &xml, const std::string &path, const SysSolverBiCGStabParams &param)
 
void Chroma::read (XMLReader &xml, const std::string &path, SysSolverCGParams &param)
 
void Chroma::write (XMLWriter &xml, const std::string &path, const SysSolverCGParams &param)
 
void Chroma::read (XMLReader &xml, const std::string &path, SysSolverEigCGParams &param)
 
void Chroma::write (XMLWriter &xml, const std::string &path, const SysSolverEigCGParams &param)
 
void Chroma::read (XMLReader &xml, const std::string &path, SysSolverFGMRESDRParams &p)
 
void Chroma::write (XMLWriter &xml, const std::string &path, const SysSolverFGMRESDRParams &p)
 
void Chroma::read (XMLReader &xml, const std::string &path, SysSolverMRParams &param)
 
void Chroma::write (XMLWriter &xml, const std::string &path, const SysSolverMRParams &param)
 
void Chroma::read (XMLReader &xml, const std::string &path, SysSolverOptEigBiCGParams &param)
 
void Chroma::write (XMLWriter &xml, const std::string &path, const SysSolverOptEigBiCGParams &param)
 
void Chroma::read (XMLReader &xml, const std::string &path, SysSolverOptEigCGParams &param)
 
void Chroma::write (XMLWriter &xml, const std::string &path, const SysSolverOptEigCGParams &param)
 
void InvCG2EvenOddPrecWilsLinOp (const WilsonDslash &D, const LFerm &chi, LFerm &psi, const LScal &mass, const LScal &RsdCG, int MaxCG, int &n_count)
 Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator. More...
 
void InvCG2EvenOddPrecWilsLinOpTHack (const WilsonDslash &D, const LFerm &chi, LFerm &psi, const LScal &mass, const LScal &RsdCG, int MaxCG, int &n_count)
 Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator. More...
 
void WlInvCG2 (const LinearOperator &M, const LatticeFermion &chi, LatticeFermion &psi, const Real &RsdCG, int MaxCG, int &n_count)
 Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator. More...
 
template<typename T >
SystemSolverResults_t Chroma::InvBiCGStab (const LinearOperator< T > &A, const T &chi, T &psi, const Real &RsdBiCGStab, int MaxBiCGStab, enum PlusMinus isign)
 Bi-CG stabilized. More...
 
template<typename T >
SystemSolverResults_t Chroma::InvBiCRStab (const LinearOperator< T > &A, const T &chi, T &psi, const Real &RsdBiCGStab, int MaxBiCGStab, enum PlusMinus isign)
 Bi-CG stabilized. More...
 
SystemSolverResults_t Chroma::InvCG2 (const LinearOperator< LatticeFermionF > &M, const LatticeFermionF &chi, LatticeFermionF &psi, const Real &RsdCG, int MaxCG)
 Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator. More...
 
SystemSolverResults_t Chroma::InvCG2 (const LinearOperator< LatticeFermionD > &M, const LatticeFermionD &chi, LatticeFermionD &psi, const Real &RsdCG, int MaxCG)
 
SystemSolverResults_t Chroma::InvCG2 (const LinearOperator< LatticeStaggeredFermionF > &M, const LatticeStaggeredFermionF &chi, LatticeStaggeredFermionF &psi, const Real &RsdCG, int MaxCG)
 
SystemSolverResults_t Chroma::InvCG2 (const LinearOperator< LatticeStaggeredFermionD > &M, const LatticeStaggeredFermionD &chi, LatticeStaggeredFermionD &psi, const Real &RsdCG, int MaxCG)
 
SystemSolverResults_t Chroma::InvCG2 (const LinearOperatorArray< LatticeFermionF > &M, const multi1d< LatticeFermionF > &chi, multi1d< LatticeFermionF > &psi, const Real &RsdCG, int MaxCG)
 Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator. More...
 
SystemSolverResults_t Chroma::InvCG2 (const LinearOperatorArray< LatticeFermionD > &M, const multi1d< LatticeFermionD > &chi, multi1d< LatticeFermionD > &psi, const Real &RsdCG, int MaxCG)
 
template<typename T >
SystemSolverResults_t Chroma::InvIBiCGStab (const LinearOperator< T > &A, const T &chi, T &psi, const Real &RsdBiCGStab, int MaxBiCGStab, enum PlusMinus isign)
 Bi-CG stabilized. More...
 
template<typename T , typename C >
SystemSolverResults_t Chroma::InvMR_a (const C &M, const T &chi, T &psi, const Real &MRovpar, const Real &RsdMR, int MaxMR, enum PlusMinus isign)
 Minimal-residual (MR) algorithm for a generic Linear Operator. More...
 
template<>
SystemSolverResults_t Chroma::InvMR (const LinearOperator< T > &M, const T &chi, T &psi, const Real &MRovpar, const Real &RsdMR, int MaxMR, enum PlusMinus isign)
 Minimal-residual (MR) algorithm for a generic Linear Operator. More...
 
template<>
SystemSolverResults_t Chroma::InvMR (const LinearOperator< LatticeStaggeredFermion > &M, const LatticeStaggeredFermion &chi, LatticeStaggeredFermion &psi, const Real &MRovpar, const Real &RsdMR, int MaxMR, enum PlusMinus isign)
 
template<typename T >
void Chroma::MInvMR_a (const LinearOperator< T > &A, const T &chi, multi1d< T > &psi, const multi1d< Real > &shifts, const multi1d< Real > &RsdMR, int MaxMR, int &n_count)
 Multishift Conjugate-Gradient (CG1) algorithm for a Linear Operator. More...
 
template<>
void Chroma::MInvMR (const LinearOperator< T > &A, const T &chi, multi1d< T > &psi, const multi1d< Real > &shifts, const multi1d< Real > &RsdCG, int MaxCG, int &n_count)
 Multishift Conjugate-Gradient (CG1) algorithm for a Linear Operator. More...
 
template<>
void Chroma::MInvMR (const DiffLinearOperator< LatticeFermion, multi1d< LatticeColorMatrix >, multi1d< LatticeColorMatrix > > &M, const LatticeFermion &chi, multi1d< LatticeFermion > &psi, const multi1d< Real > &shifts, const multi1d< Real > &RsdMR, int MaxMR, int &n_count)
 
void Chroma::normGramSchmidt (multi1d< LatticeFermionF > &vec, int f, int t, const Subset &sub)
 Gram-Schmidt with normalization. More...
 
void Chroma::normGramSchmidt (multi1d< LatticeFermionD > &vec, int f, int t, const Subset &sub)
 
void Chroma::normGramSchmidt (multi1d< LatticeStaggeredFermionF > &vec, int f, int t, const Subset &sub)
 
void Chroma::normGramSchmidt (multi1d< LatticeStaggeredFermionD > &vec, int f, int t, const Subset &sub)
 
void Chroma::normGramSchmidt (multi2d< LatticeFermionF > &vec, int f, int t, const Subset &sub)
 
void Chroma::normGramSchmidt (multi2d< LatticeFermionD > &vec, int f, int t, const Subset &sub)
 
SystemSolverResults_t Chroma::InvBiCGStabReliable (const LinearOperator< LatticeFermionF > &A, const LatticeFermionF &chi, LatticeFermionF &psi, const Real &RsdBiCGStab, const Real &Delta, int MaxBiCGStab, enum PlusMinus isign)
 Bi-CG stabilized. More...
 
SystemSolverResults_t Chroma::InvBiCGStabReliable (const LinearOperator< LatticeFermionD > &A, const LatticeFermionD &chi, LatticeFermionD &psi, const Real &RsdBiCGStab, const Real &Delta, int MaxBiCGStab, enum PlusMinus isign)
 
SystemSolverResults_t Chroma::InvBiCGStabReliable (const LinearOperator< LatticeFermionD > &A, const LinearOperator< LatticeFermionF > &AF, const LatticeFermionD &chi, LatticeFermionD &psi, const Real &RsdBiCGStab, const Real &Delta, int MaxBiCGStab, enum PlusMinus isign)
 
SystemSolverResults_t Chroma::InvCGReliable (const LinearOperator< LatticeFermionF > &A, const LatticeFermionF &chi, LatticeFermionF &psi, const Real &RsdCG, const Real &Delta, int MaxCG)
 Bi-CG stabilized. More...
 
SystemSolverResults_t Chroma::InvCGReliable (const LinearOperator< LatticeFermionD > &A, const LatticeFermionD &chi, LatticeFermionD &psi, const Real &RsdCG, const Real &Delta, int MaxCG)
 
SystemSolverResults_t Chroma::InvCGReliable (const LinearOperator< LatticeFermionD > &A, const LinearOperator< LatticeFermionF > &AF, const LatticeFermionD &chi, LatticeFermionD &psi, const Real &RsdCG, const Real &Delta, int MaxCG)
 
SystemSolverResults_t Chroma::InvIBiCGStabReliable (const LinearOperator< LatticeFermionF > &A, const LatticeFermionF &chi, LatticeFermionF &psi, const Real &RsdBiCGStab, const Real &Delta, int MaxBiCGStab, enum PlusMinus isign)
 Bi-CG stabilized. More...
 
SystemSolverResults_t Chroma::InvIBiCGStabReliable (const LinearOperator< LatticeFermionD > &A, const LatticeFermionD &chi, LatticeFermionD &psi, const Real &RsdBiCGStab, const Real &Delta, int MaxBiCGStab, enum PlusMinus isign)
 
SystemSolverResults_t Chroma::InvIBiCGStabReliable (const LinearOperator< LatticeFermionD > &A, const LinearOperator< LatticeFermionF > &AF, const LatticeFermionD &chi, LatticeFermionD &psi, const Real &RsdBiCGStab, const Real &Delta, int MaxBiCGStab, enum PlusMinus isign)
 

Detailed Description

Various inverters for fermion linear operators

Typedef Documentation

◆ TheLinOpFermMultiSystemSolverFactory

typedef SingletonHolder< ObjectFactory<LinOpMultiSystemSolver<LatticeFermion>, std::string, TYPELIST_3(XMLReader&, const std::string&, Handle< LinearOperator<LatticeFermion> >), LinOpMultiSystemSolver<LatticeFermion>* (*)(XMLReader&, const std::string&, Handle< LinearOperator<LatticeFermion> >), StringFactoryError> > Chroma::TheLinOpFermMultiSystemSolverFactory

LinOp system solver factory (foundry)

Definition at line 27 of file multi_syssolver_linop_factory.h.

◆ TheLinOpFermSystemSolverArrayFactory

typedef SingletonHolder< ObjectFactory<LinOpSystemSolverArray<LatticeFermion>, std::string, TYPELIST_4(XMLReader&, const std::string&, FSHandle, Handle< LinearOperatorArray<LatticeFermion> >), LinOpSystemSolverArray<LatticeFermion>* (*)(XMLReader&, const std::string&, FSHandle, Handle< LinearOperatorArray<LatticeFermion> >), StringFactoryError> > Chroma::TheLinOpFermSystemSolverArrayFactory

LinOp system solver factory (foundry)

Definition at line 77 of file syssolver_linop_factory.h.

◆ TheLinOpFermSystemSolverFactory

typedef SingletonHolder< ObjectFactory<LinOpSystemSolver<LatticeFermion>, std::string, TYPELIST_4(XMLReader&, const std::string&, FSHandle, Handle< LinearOperator<LatticeFermion> >), LinOpSystemSolver<LatticeFermion>* (*)(XMLReader&, const std::string&, FSHandle, Handle< LinearOperator<LatticeFermion> >), StringFactoryError> > Chroma::TheLinOpFermSystemSolverFactory

LinOp system solver factory (foundry)

Definition at line 39 of file syssolver_linop_factory.h.

◆ TheLinOpStagFermMultiSystemSolverFactory

typedef SingletonHolder< ObjectFactory<LinOpMultiSystemSolver<LatticeStaggeredFermion>, std::string, TYPELIST_3(XMLReader&, const std::string&, Handle< LinearOperator<LatticeStaggeredFermion> >), LinOpMultiSystemSolver<LatticeStaggeredFermion>* (*)(XMLReader&, const std::string&, Handle< LinearOperator<LatticeStaggeredFermion> >), StringFactoryError> > Chroma::TheLinOpStagFermMultiSystemSolverFactory

LinOp system solver factory (foundry)

Definition at line 55 of file multi_syssolver_linop_factory.h.

◆ TheLinOpStagFermSystemSolverFactory

typedef SingletonHolder< ObjectFactory<LinOpSystemSolver<LatticeStaggeredFermion>, std::string, TYPELIST_3(XMLReader&, const std::string&, Handle< LinearOperator<LatticeStaggeredFermion> >), LinOpSystemSolver<LatticeStaggeredFermion>* (*)(XMLReader&, const std::string&, Handle< LinearOperator<LatticeStaggeredFermion> >), StringFactoryError> > Chroma::TheLinOpStagFermSystemSolverFactory

LinOp system solver factory (foundry)

Definition at line 91 of file syssolver_linop_factory.h.

◆ TheMdagMFermMultiSystemSolverAccumulateArrayFactory

typedef SingletonHolder< ObjectFactory<MdagMMultiSystemSolverAccumulateArray<LatticeFermion>, std::string, TYPELIST_3(XMLReader&, const std::string&, Handle< LinearOperatorArray<LatticeFermion> >), MdagMMultiSystemSolverAccumulateArray<LatticeFermion>* (*)(XMLReader&, const std::string&, Handle< LinearOperatorArray<LatticeFermion> >), StringFactoryError> > Chroma::TheMdagMFermMultiSystemSolverAccumulateArrayFactory

MdagM system solver factory (foundry)

Definition at line 40 of file multi_syssolver_mdagm_accumulate_factory.h.

◆ TheMdagMFermMultiSystemSolverAccumulateFactory

typedef SingletonHolder< ObjectFactory<MdagMMultiSystemSolverAccumulate<LatticeFermion>, std::string, TYPELIST_3(XMLReader&, const std::string&, Handle< LinearOperator<LatticeFermion> >), MdagMMultiSystemSolverAccumulate<LatticeFermion>* (*)(XMLReader&, const std::string&, Handle< LinearOperator<LatticeFermion> >), StringFactoryError> > Chroma::TheMdagMFermMultiSystemSolverAccumulateFactory

MdagM system solver factory (foundry)

Definition at line 27 of file multi_syssolver_mdagm_accumulate_factory.h.

◆ TheMdagMFermMultiSystemSolverArrayFactory

typedef SingletonHolder< ObjectFactory<MdagMMultiSystemSolverArray<LatticeFermion>, std::string, TYPELIST_4(XMLReader&, const std::string&, FSHandle, Handle< LinearOperatorArray<LatticeFermion> >), MdagMMultiSystemSolverArray<LatticeFermion>* (*)(XMLReader&, const std::string&, FSHandle, Handle< LinearOperatorArray<LatticeFermion> >), StringFactoryError> > Chroma::TheMdagMFermMultiSystemSolverArrayFactory

MdagM system solver factory (foundry)

Definition at line 47 of file multi_syssolver_mdagm_factory.h.

◆ TheMdagMFermMultiSystemSolverFactory

typedef SingletonHolder< ObjectFactory<MdagMMultiSystemSolver<LatticeFermion>, std::string, TYPELIST_4(XMLReader&, const std::string&, FSHandle, Handle< LinearOperator<LatticeFermion> >), MdagMMultiSystemSolver<LatticeFermion>* (*)(XMLReader&, const std::string&, FSHandle, Handle< LinearOperator<LatticeFermion> >), StringFactoryError> > Chroma::TheMdagMFermMultiSystemSolverFactory

MdagM system solver factory (foundry)

Definition at line 33 of file multi_syssolver_mdagm_factory.h.

◆ TheMdagMFermSystemSolverArrayFactory

typedef SingletonHolder< ObjectFactory<MdagMSystemSolverArray<LatticeFermion>, std::string, TYPELIST_4(XMLReader&, const std::string&, FactoryEnv::FSHandle, Handle< LinearOperatorArray<LatticeFermion> >), MdagMSystemSolverArray<LatticeFermion>* (*)(XMLReader&, const std::string&, FactoryEnv::FSHandle, Handle< LinearOperatorArray<LatticeFermion> >), StringFactoryError> > Chroma::TheMdagMFermSystemSolverArrayFactory

MdagM system solver factory (foundry)

Definition at line 75 of file syssolver_mdagm_factory.h.

◆ TheMdagMFermSystemSolverFactory

typedef SingletonHolder< ObjectFactory<MdagMSystemSolver<LatticeFermion>, std::string, TYPELIST_4(XMLReader&, const std::string&, FactoryEnv::FSHandle, Handle< LinearOperator<LatticeFermion> >), MdagMSystemSolver<LatticeFermion>* (*)(XMLReader&, const std::string&, FactoryEnv::FSHandle, Handle< LinearOperator<LatticeFermion> >), StringFactoryError> > Chroma::TheMdagMFermSystemSolverFactory

MdagM system solver factory (foundry)

Definition at line 39 of file syssolver_mdagm_factory.h.

◆ TheMdagMStagFermMultiSystemSolverAccumulateFactory

typedef SingletonHolder< ObjectFactory<MdagMMultiSystemSolverAccumulate<LatticeStaggeredFermion>, std::string, TYPELIST_3(XMLReader&, const std::string&, Handle< LinearOperator<LatticeStaggeredFermion> >), MdagMMultiSystemSolverAccumulate<LatticeStaggeredFermion>* (*)(XMLReader&, const std::string&, Handle< LinearOperator<LatticeStaggeredFermion> >), StringFactoryError> > Chroma::TheMdagMStagFermMultiSystemSolverAccumulateFactory

MdagM system solver factory (foundry)

Definition at line 53 of file multi_syssolver_mdagm_accumulate_factory.h.

◆ TheMdagMStagFermMultiSystemSolverFactory

typedef SingletonHolder< ObjectFactory<MdagMMultiSystemSolver<LatticeStaggeredFermion>, std::string, TYPELIST_3(XMLReader&, const std::string&, Handle< LinearOperator<LatticeStaggeredFermion> >), MdagMMultiSystemSolver<LatticeStaggeredFermion>* (*)(XMLReader&, const std::string&, Handle< LinearOperator<LatticeStaggeredFermion> >), StringFactoryError> > Chroma::TheMdagMStagFermMultiSystemSolverFactory

MdagM system solver factory (foundry)

Definition at line 60 of file multi_syssolver_mdagm_factory.h.

◆ TheMdagMStagFermSystemSolverFactory

typedef SingletonHolder< ObjectFactory<MdagMSystemSolver<LatticeStaggeredFermion>, std::string, TYPELIST_3(XMLReader&, const std::string&, Handle< LinearOperator<LatticeStaggeredFermion> >), MdagMSystemSolver<LatticeStaggeredFermion>* (*)(XMLReader&, const std::string&, Handle< LinearOperator<LatticeStaggeredFermion> >), StringFactoryError> > Chroma::TheMdagMStagFermSystemSolverFactory

MdagM system solver factory (foundry)

Definition at line 88 of file syssolver_mdagm_factory.h.

◆ ThePolyPrecFermSystemSolverFactory

typedef SingletonHolder< ObjectFactory<PolyPrecSystemSolver<LatticeFermion>, std::string, TYPELIST_3(XMLReader&, const std::string&, Handle< LinearOperator<LatticeFermion> >), PolyPrecSystemSolver<LatticeFermion>* (*)(XMLReader&, const std::string&, Handle< LinearOperator<LatticeFermion> >), StringFactoryError> > Chroma::ThePolyPrecFermSystemSolverFactory

PolyPrec system solver factory (foundry)

Definition at line 26 of file syssolver_polyprec_factory.h.

Function Documentation

◆ InvBiCGStab()

template<typename T >
SystemSolverResults_t Chroma::InvBiCGStab ( const LinearOperator< T > &  A,
const T chi,
T psi,
const Real &  RsdBiCGStab,
int  MaxBiCGStab,
enum PlusMinus  isign 
)

Bi-CG stabilized.

◆ InvBiCGStabReliable() [1/3]

SystemSolverResults_t Chroma::InvBiCGStabReliable ( const LinearOperator< LatticeFermionD > &  A,
const LatticeFermionD &  chi,
LatticeFermionD &  psi,
const Real &  RsdBiCGStab,
const Real &  Delta,
int  MaxBiCGStab,
enum PlusMinus  isign 
)

Definition at line 309 of file reliable_bicgstab.cc.

References Chroma::A(), Chroma::chi(), Chroma::isign, and Chroma::psi.

◆ InvBiCGStabReliable() [2/3]

SystemSolverResults_t Chroma::InvBiCGStabReliable ( const LinearOperator< LatticeFermionD > &  A,
const LinearOperator< LatticeFermionF > &  AF,
const LatticeFermionD &  chi,
LatticeFermionD &  psi,
const Real &  RsdBiCGStab,
const Real &  Delta,
int  MaxBiCGStab,
enum PlusMinus  isign 
)

Definition at line 323 of file reliable_bicgstab.cc.

References Chroma::A(), Chroma::chi(), Chroma::isign, and Chroma::psi.

◆ InvBiCGStabReliable() [3/3]

SystemSolverResults_t Chroma::InvBiCGStabReliable ( const LinearOperator< LatticeFermionF > &  A,
const LatticeFermionF &  chi,
LatticeFermionF &  psi,
const Real &  RsdBiCGStab,
const Real &  Delta,
int  MaxBiCGStab,
enum PlusMinus  isign 
)

◆ InvBiCRStab()

template<typename T >
SystemSolverResults_t Chroma::InvBiCRStab ( const LinearOperator< T > &  A,
const T chi,
T psi,
const Real &  RsdBiCGStab,
int  MaxBiCGStab,
enum PlusMinus  isign 
)

Bi-CG stabilized.

◆ InvCG1() [1/2]

template<>
SystemSolverResults_t Chroma::InvCG1 ( const LinearOperator< T > &  A,
const T chi,
T psi,
const Real &  RsdCG,
int  MaxCG,
int  MinCG = 0 
)

Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations

        Chi  =  A . Psi

where A is hermitian

Algorithm:

Psi[0] := initial guess; Linear interpolation (argument) r[0] := Chi - A. Psi[0] ; Initial residual p[1] := r[0] ; Initial direction IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k-1]|**2 / <p[k],A p[k]> ; Psi[k] += a[k] p[k] ; New solution std::vector r[k] -= a[k] A p[k] ; New residual IF |r[k]| <= RsdCG |Chi| THEN RETURN; Converged? b[k+1] := |r[k]|**2 / |r[k-1]|**2 ; p[k+1] := r[k] + b[k+1] p[k]; New direction

Arguments:

Parameters
ALinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGCG residual accuracy (Read)
MaxCGMaximum CG iterations (Read)
Returns
res System solver results

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > ap Temporary for M.p

Subroutines:


  • A Apply matrix M or M to std::vector

Operations:

2 A + 2 Nc Ns + N_Count ( 2 A + 10 Nc Ns )

Definition at line 215 of file invcg1.cc.

References Chroma::A(), chi, Chroma::InvCG1_a(), MaxCG, psi, and Chroma::RsdCG.

Referenced by Chroma::InvGMRESR_CG_a(), main(), Chroma::Ovlap4DQprop::operator()(), Chroma::PolyPrecSysSolverCG< T >::operator()(), and Chroma::EvenOddFermActQprop< T, P, Q >::operator()().

◆ InvCG1() [2/2]

template<>
void Chroma::InvCG1 ( const LinearOperatorArray< T > &  A,
const multi1d< T > &  chi,
multi1d< T > &  psi,
const Real &  RsdCG,
int  MaxCG,
int &  n_count 
)

Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations

        Chi  =  A . Psi

where A is hermitian

Algorithm:

Psi[0] := initial guess; Linear interpolation (argument) r[0] := Chi - A . Psi[0] ; Initial residual p[1] := r[0] ; Initial direction IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k-1]|**2 / <A p[k],p[k]> ; Psi[k] += a[k] p[k] ; New solution std::vector r[k] -= a[k] A. p[k] ; New residual IF |r[k]| <= RsdCG |Chi| THEN RETURN; Converged? b[k+1] := |r[k]|**2 / |r[k-1]|**2 ; p[k+1] := r[k] + b[k+1] p[k]; New direction

Arguments:

Parameters
MLinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGCG residual accuracy (Read)
MaxCGMaximum CG iterations (Read)
n_countNumber of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Mp Temporary for M.p

Subroutines:


Operations:

2 A + 2 Nc Ns + N_Count ( 2 A + 10 Nc Ns )

Definition at line 216 of file invcg1_array.cc.

References Chroma::A(), chi, Chroma::InvCG1_a(), MaxCG, n_count, psi, and Chroma::RsdCG.

◆ InvCG1_a() [1/2]

template<typename T >
SystemSolverResults_t Chroma::InvCG1_a ( const LinearOperator< T > &  A,
const T chi,
T psi,
const Real &  RsdCG,
int  MaxCG,
int  MinCG = 0 
)

Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations

        Chi  =  A . Psi

where A is Hermitian Positive Definite

Algorithm:

Psi[0] := initial guess; Linear interpolation (argument) r[0] := Chi - A. Psi[0] ; Initial residual p[1] := r[0] ; Initial direction IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k-1]|**2 / <p[k],Ap[k]> ; Psi[k] += a[k] p[k] ; New solution std::vector r[k] -= a[k] A. p[k] ; New residual IF |r[k]| <= RsdCG |Chi| THEN RETURN; Converged? b[k+1] := |r[k]|**2 / |r[k-1]|**2 ; p[k+1] := r[k] + b[k+1] p[k]; New direction

Arguments:

Parameters
MLinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGCG residual accuracy (Rea/Write)
MaxCGMaximum CG iterations (Read)
Returns
res System solver results

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > ap Temporary for A.p

Subroutines:


  • A Apply matrix M or M to std::vector

Operations:

2 A + 2 Nc Ns + N_Count ( 2 A + 10 Nc Ns )

Definition at line 69 of file invcg1.cc.

References Chroma::a, Chroma::A(), Chroma::b, c, chi, Chroma::cp, Chroma::d, END_CODE, Chroma::k, MaxCG, Chroma::SystemSolverResults_t::n_count, Chroma::p, Chroma::PLUS, psi, Chroma::QDP_error_exit(), r(), Chroma::SystemSolverResults_t::resid, Chroma::rsd_sq, Chroma::RsdCG, s, and START_CODE.

Referenced by Chroma::InvCG1().

◆ InvCG1_a() [2/2]

template<typename T >
void Chroma::InvCG1_a ( const LinearOperatorArray< T > &  A,
const multi1d< T > &  chi,
multi1d< T > &  psi,
const Real &  RsdCG,
int  MaxCG,
int &  n_count 
)

Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations

        Chi  =  A . Psi

where A is hermitian

Algorithm:

Psi[0] := initial guess; Linear interpolation (argument) r[0] := Chi - M^dag . M . Psi[0] ; Initial residual p[1] := r[0] ; Initial direction IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k-1]|**2 / <Mp[k],Mp[k]> ; Psi[k] += a[k] p[k] ; New solution std::vector r[k] -= a[k] M^dag . M . p[k] ; New residual IF |r[k]| <= RsdCG |Chi| THEN RETURN; Converged? b[k+1] := |r[k]|**2 / |r[k-1]|**2 ; p[k+1] := r[k] + b[k+1] p[k]; New direction

Arguments:

Parameters
MLinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGCG residual accuracy (Read)
MaxCGMaximum CG iterations (Read)
n_countNumber of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Mp Temporary for M.p

Subroutines:


  • A Apply matrix M or M to std::vector

Operations:

2 A + 2 Nc Ns + N_Count ( 2 A + 10 Nc Ns )

Definition at line 67 of file invcg1_array.cc.

References Chroma::a, Chroma::A(), Chroma::b, c, chi, Chroma::cp, Chroma::d, END_CODE, i, Chroma::k, MaxCG, n, n_count, Chroma::p, Chroma::PLUS, psi, r(), Chroma::rsd_sq, Chroma::RsdCG, s, and START_CODE.

Referenced by Chroma::InvCG1().

◆ InvCG2() [1/6]

SystemSolverResults_t Chroma::InvCG2 ( const LinearOperator< LatticeFermionD > &  M,
const LatticeFermionD &  chi,
LatticeFermionD &  psi,
const Real &  RsdCG,
int  MaxCG 
)

Definition at line 258 of file invcg2.cc.

References chi, MaxCG, psi, and Chroma::RsdCG.

◆ InvCG2() [2/6]

SystemSolverResults_t Chroma::InvCG2 ( const LinearOperator< LatticeFermionF > &  M,
const LatticeFermionF &  chi,
LatticeFermionF &  psi,
const Real &  RsdCG,
int  MaxCG 
)

Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations

        Chi  =  A . Psi

where A = M^dag . M

Algorithm:

Psi[0] := initial guess; Linear interpolation (argument) r[0] := Chi - M^dag . M . Psi[0] ; Initial residual p[1] := r[0] ; Initial direction IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k-1]|**2 / <Mp[k],Mp[k]> ; Psi[k] += a[k] p[k] ; New solution std::vector r[k] -= a[k] M^dag . M . p[k] ; New residual IF |r[k]| <= RsdCG |Chi| THEN RETURN; Converged? b[k+1] := |r[k]|**2 / |r[k-1]|**2 ; p[k+1] := r[k] + b[k+1] p[k]; New direction

Arguments:

Parameters
MLinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGCG residual accuracy (Read)
MaxCGMaximum CG iterations (Read)
n_countNumber of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Mp Temporary for M.p

Subroutines:


  • A Apply matrix M or M to std::vector

Operations:

2 A + 2 Nc Ns + N_Count ( 2 A + 10 Nc Ns )

Definition at line 240 of file invcg2.cc.

References chi, MaxCG, psi, and Chroma::RsdCG.

Referenced by main(), Chroma::UnprecDWF4DLinOp< T >::operator()(), Chroma::UnprecDWFTransfLinOp::operator()(), Chroma::UnprecDWFTransfMdagMLinOp::operator()(), Chroma::UnprecPDWF4DLinOp< T, P, Q >::operator()(), Chroma::UnprecPPDWF4DLinOp< T, P, Q >::operator()(), Chroma::Ovlap4DQprop::operator()(), Chroma::LinOpSysSolverCGArray< T >::operator()(), Chroma::MdagMSysSolverCGArray< T >::operator()(), Chroma::HtContFrac5DQprop< T, P, Q >::operator()(), Chroma::PrecOvExt5DQprop< T, P, Q >::operator()(), Chroma::ContFrac5DQprop< T, P, Q >::operator()(), Chroma::OvHTCFZ5DQprop< T >::operator()(), Chroma::OvExt5DQprop< T >::operator()(), Chroma::OvUnprecCF5DQprop< T >::operator()(), Chroma::LinOpSysSolverCG< T >::operator()(), Chroma::MdagMSysSolverBiCGStab< T >::operator()(), Chroma::MdagMSysSolverCG< T >::operator()(), Chroma::MdagMSysSolverCGLFClover::operator()(), and Chroma::MdagMSysSolverIBiCGStab< T >::operator()().

◆ InvCG2() [3/6]

SystemSolverResults_t Chroma::InvCG2 ( const LinearOperator< LatticeStaggeredFermionD > &  M,
const LatticeStaggeredFermionD &  chi,
LatticeStaggeredFermionD &  psi,
const Real &  RsdCG,
int  MaxCG 
)

Definition at line 286 of file invcg2.cc.

References chi, MaxCG, psi, and Chroma::RsdCG.

◆ InvCG2() [4/6]

SystemSolverResults_t Chroma::InvCG2 ( const LinearOperator< LatticeStaggeredFermionF > &  M,
const LatticeStaggeredFermionF &  chi,
LatticeStaggeredFermionF &  psi,
const Real &  RsdCG,
int  MaxCG 
)

Definition at line 275 of file invcg2.cc.

References chi, MaxCG, psi, and Chroma::RsdCG.

◆ InvCG2() [5/6]

SystemSolverResults_t Chroma::InvCG2 ( const LinearOperatorArray< LatticeFermionD > &  M,
const multi1d< LatticeFermionD > &  chi,
multi1d< LatticeFermionD > &  psi,
const Real &  RsdCG,
int  MaxCG 
)

◆ InvCG2() [6/6]

SystemSolverResults_t Chroma::InvCG2 ( const LinearOperatorArray< LatticeFermionF > &  M,
const multi1d< LatticeFermionF > &  chi,
multi1d< LatticeFermionF > &  psi,
const Real &  RsdCG,
int  MaxCG 
)

Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations

        Chi  =  A . Psi

where A = M^dag . M

Algorithm:

Psi[0] := initial guess; Linear interpolation (argument) r[0] := Chi - M^dag . M . Psi[0] ; Initial residual p[1] := r[0] ; Initial direction IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k-1]|**2 / <Mp[k],Mp[k]> ; Psi[k] += a[k] p[k] ; New solution std::vector r[k] -= a[k] M^dag . M . p[k] ; New residual IF |r[k]| <= RsdCG |Chi| THEN RETURN; Converged? b[k+1] := |r[k]|**2 / |r[k-1]|**2 ; p[k+1] := r[k] + b[k+1] p[k]; New direction

Arguments:

Parameters
MLinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGCG residual accuracy (Read)
MaxCGMaximum CG iterations (Read)
n_countNumber of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Mp Temporary for M.p

Subroutines:


  • A Apply matrix M or M to std::vector

Operations:

2 A + 2 Nc Ns + N_Count ( 2 A + 10 Nc Ns )

Definition at line 281 of file invcg2_array.cc.

References Chroma::chi(), Chroma::InvCG2_a(), Chroma::MaxCG, Chroma::psi, and Chroma::RsdCG.

◆ InvCG2_a() [1/2]

template<typename T , typename C >
SystemSolverResults_t Chroma::InvCG2_a ( const C &  M,
const multi1d< T > &  chi,
multi1d< T > &  psi,
const Real &  RsdCG,
int  MaxCG 
)

Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations

        Chi  =  A . Psi

where A = M^dag . M

Algorithm:

Psi[0] := initial guess; Linear interpolation (argument) r[0] := Chi - M^dag . M . Psi[0] ; Initial residual p[1] := r[0] ; Initial direction IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k-1]|**2 / <Mp[k],Mp[k]> ; Psi[k] += a[k] p[k] ; New solution std::vector r[k] -= a[k] M^dag . M . p[k] ; New residual IF |r[k]| <= RsdCG |Chi| THEN RETURN; Converged? b[k+1] := |r[k]|**2 / |r[k-1]|**2 ; p[k+1] := r[k] + b[k+1] p[k]; New direction

Arguments:

Parameters
MLinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGCG residual accuracy (Read)
MaxCGMaximum CG iterations (Read)
Returns
res System solver results

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Mp Temporary for M.p

Subroutines:


  • A Apply matrix M or M to std::vector

Operations:

2 A + 2 Nc Ns + N_Count ( 2 A + 10 Nc Ns )

Definition at line 68 of file invcg2_array.cc.

References Chroma::a, Chroma::b, Chroma::c, Chroma::chi(), Chroma::cp, Chroma::d, Chroma::END_CODE(), Chroma::i, Chroma::k, Chroma::MaxCG, Chroma::MINUS, Chroma::mp(), n, Chroma::SystemSolverResults_t::n_count, Chroma::p, Chroma::PLUS, Chroma::psi, Chroma::r, Chroma::SystemSolverResults_t::resid, Chroma::rsd_sq, Chroma::RsdCG, Chroma::START_CODE(), and Chroma::zero.

◆ InvCG2_a() [2/2]

template<typename T , typename RT >
SystemSolverResults_t Chroma::InvCG2_a ( const LinearOperator< T > &  M,
const T chi,
T psi,
const RT &  RsdCG,
int  MaxCG 
)

Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations

        Chi  =  A . Psi

where A = M^dag . M

Algorithm:

Psi[0] := initial guess; Linear interpolation (argument) r[0] := Chi - M^dag . M . Psi[0] ; Initial residual p[1] := r[0] ; Initial direction IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k-1]|**2 / <Mp[k],Mp[k]> ; Psi[k] += a[k] p[k] ; New solution std::vector r[k] -= a[k] M^dag . M . p[k] ; New residual IF |r[k]| <= RsdCG |Chi| THEN RETURN; Converged? b[k+1] := |r[k]|**2 / |r[k-1]|**2 ; p[k+1] := r[k] + b[k+1] p[k]; New direction

Arguments:

Parameters
MLinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGCG residual accuracy (Read)
MaxCGMaximum CG iterations (Read)
Returns
res System solver results

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Mp Temporary for M.p

Subroutines:


  • A Apply matrix M or M to std::vector

Operations:

2 A + 2 Nc Ns + N_Count ( 2 A + 10 Nc Ns )

Definition at line 72 of file invcg2.cc.

References Chroma::a, Chroma::b, c, chi, Chroma::cp, Chroma::d, END_CODE, Chroma::k, MaxCG, Chroma::MINUS, Chroma::mp(), Chroma::SystemSolverResults_t::n_count, Chroma::LinearOperator< T >::nFlops(), Chroma::p, Chroma::PLUS, psi, r(), Chroma::SystemSolverResults_t::resid, Chroma::rsd_sq, Chroma::RsdCG, s, START_CODE, and Chroma::LinearOperator< T >::subset().

Referenced by Chroma::InvCG2().

◆ InvCG2_timings()

template<>
SystemSolverResults_t Chroma::InvCG2_timings ( const LinearOperator< T > &  M,
const T chi,
T psi,
int  MaxCG 
)

Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations

        Chi  =  A . Psi

where A = M^dag . M

Algorithm:

Psi[0] := initial guess; Linear interpolation (argument) r[0] := Chi - M^dag . M . Psi[0] ; Initial residual p[1] := r[0] ; Initial direction IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k-1]|**2 / <Mp[k],Mp[k]> ; Psi[k] += a[k] p[k] ; New solution std::vector r[k] -= a[k] M^dag . M . p[k] ; New residual IF |r[k]| <= RsdCG |Chi| THEN RETURN; Converged? b[k+1] := |r[k]|**2 / |r[k-1]|**2 ; p[k+1] := r[k] + b[k+1] p[k]; New direction

Arguments:

Parameters
MLinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGCG residual accuracy (Read)
MaxCGMaximum CG iterations (Read)
n_countNumber of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Mp Temporary for M.p

Subroutines:


  • A Apply matrix M or M to std::vector

Operations:

2 A + 2 Nc Ns + N_Count ( 2 A + 10 Nc Ns )

Definition at line 184 of file invcg2_timing_hacks.cc.

References Chroma::chi(), Chroma::InvCG2_timings_a(), Chroma::MaxCG, and Chroma::psi.

Referenced by Chroma::LinOpSysSolverCGTiming< T >::operator()(), and Chroma::MdagMSysSolverCGTimings< T >::operator()().

◆ InvCG2_timings_a()

template<typename T >
SystemSolverResults_t Chroma::InvCG2_timings_a ( const LinearOperator< T > &  M,
const T chi,
T psi,
int  MaxCG 
)

Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations

        Chi  =  A . Psi

where A = M^dag . M

Algorithm:

Psi[0] := initial guess; Linear interpolation (argument) r[0] := Chi - M^dag . M . Psi[0] ; Initial residual p[1] := r[0] ; Initial direction IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k-1]|**2 / <Mp[k],Mp[k]> ; Psi[k] += a[k] p[k] ; New solution std::vector r[k] -= a[k] M^dag . M . p[k] ; New residual IF |r[k]| <= RsdCG |Chi| THEN RETURN; Converged? b[k+1] := |r[k]|**2 / |r[k-1]|**2 ; p[k+1] := r[k] + b[k+1] p[k]; New direction

Arguments:

Parameters
MLinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGCG residual accuracy (Read)
MaxCGMaximum CG iterations (Read)
n_countNumber of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Mp Temporary for M.p

Subroutines:


  • A Apply matrix M or M to std::vector

Operations:

2 A + 2 Nc Ns + N_Count ( 2 A + 10 Nc Ns )

Definition at line 70 of file invcg2_timing_hacks.cc.

References Chroma::a, Chroma::b, Chroma::c, Chroma::chi(), Chroma::cp, Chroma::d, Chroma::END_CODE(), Chroma::k, Chroma::MaxCG, Chroma::MINUS, Chroma::mp(), Chroma::SystemSolverResults_t::n_count, Chroma::LinearOperator< T >::nFlops(), Chroma::p, Chroma::PLUS, Chroma::psi, Chroma::r, Chroma::SystemSolverResults_t::resid, Chroma::s(), Chroma::START_CODE(), and Chroma::LinearOperator< T >::subset().

Referenced by Chroma::InvCG2_timings().

◆ InvCG2EvenOddPrecWilsLinOp()

void InvCG2EvenOddPrecWilsLinOp ( const WilsonDslash &  D,
const LFerm chi,
LFerm psi,
const LScal mass,
const LScal RsdCG,
int  MaxCG,
int &  n_count 
)

Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations

        Chi  =  A . Psi

where A = M^dag . M

Algorithm:

Psi[0] := initial guess; Linear interpolation (argument) r[0] := Chi - M^dag . M . Psi[0] ; Initial residual p[1] := r[0] ; Initial direction IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k-1]|**2 / <Mp[k],Mp[k]> ; Psi[k] += a[k] p[k] ; New solution std::vector r[k] -= a[k] M^dag . M . p[k] ; New residual IF |r[k]| <= RsdCG |Chi| THEN RETURN; Converged? b[k+1] := |r[k]|**2 / |r[k-1]|**2 ; p[k+1] := r[k] + b[k+1] p[k]; New direction

Arguments:

Parameters
MLinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGCG residual accuracy (Read)
MaxCGMaximum CG iterations (Read)
n_countNumber of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Mp Temporary for M.p

Subroutines:


  • A Apply matrix M or M to std::vector

Operations:

2 A + 2 Nc Ns + N_Count ( 2 A + 10 Nc Ns )

Definition at line 66 of file invcg2_prec_wilson.cc.

References Chroma::a, AT_REAL, Chroma::b, c, chi, Chroma::cp, Chroma::d, Chroma::DOUBLE, FIRST_ELEM, Chroma::k, Chroma::mass, MaxCG, Chroma::MINUS, Chroma::mp(), n_count, Nd, Chroma::p, Chroma::PLUS, psi, Chroma::QDP_error_exit(), r(), Chroma::rsd_sq, Chroma::RsdCG, s, sum, tmp2, tmp3, and vaxpy3_norm().

◆ InvCG2EvenOddPrecWilsLinOpTHack()

void InvCG2EvenOddPrecWilsLinOpTHack ( const WilsonDslash &  D,
const LFerm chi,
LFerm psi,
const LScal mass,
const LScal RsdCG,
int  MaxCG,
int &  n_count 
)

Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations

        Chi  =  A . Psi

where A = M^dag . M

Algorithm:

Psi[0] := initial guess; Linear interpolation (argument) r[0] := Chi - M^dag . M . Psi[0] ; Initial residual p[1] := r[0] ; Initial direction IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k-1]|**2 / <Mp[k],Mp[k]> ; Psi[k] += a[k] p[k] ; New solution std::vector r[k] -= a[k] M^dag . M . p[k] ; New residual IF |r[k]| <= RsdCG |Chi| THEN RETURN; Converged? b[k+1] := |r[k]|**2 / |r[k-1]|**2 ; p[k+1] := r[k] + b[k+1] p[k]; New direction

Arguments:

Parameters
MLinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGCG residual accuracy (Read)
MaxCGMaximum CG iterations (Read)
n_countNumber of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Mp Temporary for M.p

Subroutines:


  • A Apply matrix M or M to std::vector

Operations:

2 A + 2 Nc Ns + N_Count ( 2 A + 10 Nc Ns )

Definition at line 66 of file invcg2_timing_hacks_2.cc.

References Chroma::a, AT_REAL, Chroma::b, c, chi, Chroma::cp, Chroma::d, Chroma::DOUBLE, FIRST_ELEM, Chroma::k, Chroma::mass, Chroma::MINUS, Chroma::mp(), n_count, Nd, Chroma::p, Chroma::PLUS, psi, r(), Chroma::rsd_sq, Chroma::RsdCG, s, sum, tmp2, tmp3, and vaxpy3_norm().

◆ InvCGReliable() [1/3]

SystemSolverResults_t Chroma::InvCGReliable ( const LinearOperator< LatticeFermionD > &  A,
const LatticeFermionD &  chi,
LatticeFermionD &  psi,
const Real &  RsdCG,
const Real &  Delta,
int  MaxCG 
)

Definition at line 206 of file reliable_cg.cc.

References Chroma::A(), Chroma::chi(), Chroma::MaxCG, Chroma::psi, and Chroma::RsdCG.

◆ InvCGReliable() [2/3]

SystemSolverResults_t Chroma::InvCGReliable ( const LinearOperator< LatticeFermionD > &  A,
const LinearOperator< LatticeFermionF > &  AF,
const LatticeFermionD &  chi,
LatticeFermionD &  psi,
const Real &  RsdCG,
const Real &  Delta,
int  MaxCG 
)

Definition at line 218 of file reliable_cg.cc.

References Chroma::A(), Chroma::chi(), Chroma::MaxCG, Chroma::psi, and Chroma::RsdCG.

◆ InvCGReliable() [3/3]

SystemSolverResults_t Chroma::InvCGReliable ( const LinearOperator< LatticeFermionF > &  A,
const LatticeFermionF &  chi,
LatticeFermionF &  psi,
const Real &  RsdCG,
const Real &  Delta,
int  MaxCG 
)

◆ InvIBiCGStab()

template<typename T >
SystemSolverResults_t Chroma::InvIBiCGStab ( const LinearOperator< T > &  A,
const T chi,
T psi,
const Real &  RsdBiCGStab,
int  MaxBiCGStab,
enum PlusMinus  isign 
)

Bi-CG stabilized.

◆ InvIBiCGStabReliable() [1/3]

SystemSolverResults_t Chroma::InvIBiCGStabReliable ( const LinearOperator< LatticeFermionD > &  A,
const LatticeFermionD &  chi,
LatticeFermionD &  psi,
const Real &  RsdBiCGStab,
const Real &  Delta,
int  MaxBiCGStab,
enum PlusMinus  isign 
)

Definition at line 451 of file reliable_ibicgstab.cc.

References Chroma::A(), Chroma::chi(), Chroma::isign, and Chroma::psi.

◆ InvIBiCGStabReliable() [2/3]

SystemSolverResults_t Chroma::InvIBiCGStabReliable ( const LinearOperator< LatticeFermionD > &  A,
const LinearOperator< LatticeFermionF > &  AF,
const LatticeFermionD &  chi,
LatticeFermionD &  psi,
const Real &  RsdBiCGStab,
const Real &  Delta,
int  MaxBiCGStab,
enum PlusMinus  isign 
)

Definition at line 465 of file reliable_ibicgstab.cc.

References Chroma::A(), Chroma::chi(), Chroma::isign, and Chroma::psi.

◆ InvIBiCGStabReliable() [3/3]

SystemSolverResults_t Chroma::InvIBiCGStabReliable ( const LinearOperator< LatticeFermionF > &  A,
const LatticeFermionF &  chi,
LatticeFermionF &  psi,
const Real &  RsdBiCGStab,
const Real &  Delta,
int  MaxBiCGStab,
enum PlusMinus  isign 
)

◆ InvMR() [1/2]

template<>
SystemSolverResults_t Chroma::InvMR ( const LinearOperator< LatticeStaggeredFermion > &  M,
const LatticeStaggeredFermion &  chi,
LatticeStaggeredFermion &  psi,
const Real &  MRovpar,
const Real &  RsdMR,
int  MaxMR,
enum PlusMinus  isign 
)

Definition at line 200 of file invmr.cc.

References chi, Chroma::InvMR_a(), Chroma::isign, psi, and RsdMR.

◆ InvMR() [2/2]

template<>
SystemSolverResults_t Chroma::InvMR ( const LinearOperator< T > &  M,
const T chi,
T psi,
const Real &  MRovpar,
const Real &  RsdMR,
int  MaxMR,
enum PlusMinus  isign 
)

Minimal-residual (MR) algorithm for a generic Linear Operator.

This subroutine uses the Minimal Residual (MR) algorithm to determine the solution of the set of linear equations. Here we allow M to be nonhermitian.

        Chi  =  M . Psi 

Algorithm:

Psi[0] Argument r[0] := Chi - M . Psi[0] ; Initial residual IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO MR iterations a[k-1] := <M.r[k-1],r[k-1]> / <M.r[k-1],M.r[k-1]> ; ap[k-1] := MRovpar * a[k] ; Overrelaxtion step Psi[k] += ap[k-1] r[k-1] ; New solution std::vector r[k] -= ap[k-1] A . r[k-1] ; New residual IF |r[k]| <= RsdCG |Chi| THEN RETURN; Converged?

Arguments:

Parameters
MLinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGMR residual accuracy (Read)
MRovparOverrelaxation parameter (Read)
MaxCGMaximum MR iterations (Read)

Local Variables:

r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k MR iteration counter a a[k] d < M.r[k], M.r[k] > R_Aux Temporary for M.Psi Mr Temporary for M.r

Global Variables:

MaxCG Maximum number of MR iterations allowed RsdCG Maximum acceptable MR residual (relative to source)

Subroutines:

M Apply matrix to std::vector

Definition at line 185 of file invmr.cc.

References chi, Chroma::InvMR_a(), Chroma::isign, psi, and RsdMR.

Referenced by Chroma::Ovlap4DQprop::operator()(), Chroma::LinOpSysSolverMR< T >::operator()(), and Chroma::MdagMSysSolverMR< T >::operator()().

◆ InvMR_a()

template<typename T , typename C >
SystemSolverResults_t Chroma::InvMR_a ( const C &  M,
const T chi,
T psi,
const Real &  MRovpar,
const Real &  RsdMR,
int  MaxMR,
enum PlusMinus  isign 
)

Minimal-residual (MR) algorithm for a generic Linear Operator.

This subroutine uses the Minimal Residual (MR) algorithm to determine the solution of the set of linear equations. Here we allow M to be nonhermitian.

        Chi  =  M . Psi 

Algorithm:

Psi[0] Argument r[0] := Chi - M . Psi[0] ; Initial residual IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO MR iterations a[k-1] := <M.r[k-1],r[k-1]> / <M.r[k-1],M.r[k-1]> ; ap[k-1] := MRovpar * a[k] ; Overrelaxtion step Psi[k] += ap[k-1] r[k-1] ; New solution std::vector r[k] -= ap[k-1] A . r[k-1] ; New residual IF |r[k]| <= RsdCG |Chi| THEN RETURN; Converged?

Arguments:

Parameters
MLinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGMR residual accuracy (Read)
MRovparOverrelaxation parameter (Read)
MaxMRMaximum MR iterations (Read)

Local Variables:

r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k MR iteration counter a a[k] d < M.r[k], M.r[k] > R_Aux Temporary for M.Psi Mr Temporary for M.r

Global Variables:

MaxMR Maximum number of MR iterations allowed RsdCG Maximum acceptable MR residual (relative to source)

Subroutines:

M Apply matrix to std::vector

Definition at line 66 of file invmr.cc.

References Chroma::a, c, chi, Chroma::cp, Chroma::d, END_CODE, Chroma::InlinePropAndMatElemDistillation2Env::local::innerProduct(), Chroma::isign, Chroma::k, Chroma::SystemSolverResults_t::n_count, psi, r(), Chroma::SystemSolverResults_t::resid, Chroma::rsd_sq, RsdMR, s, and START_CODE.

Referenced by Chroma::InvMR().

◆ InvRelCG1()

template<>
void Chroma::InvRelCG1 ( const LinearOperator< T > &  A,
const T chi,
T psi,
const Real &  RsdCG,
int  MaxCG,
int &  n_count 
)

Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations

        Chi  =  A . Psi

where A is hermitian

Algorithm:

Psi[0] := initial guess; Linear interpolation (argument) r[0] := Chi - A. Psi[0] ; Initial residual p[1] := r[0] ; Initial direction c := cp := || r[0] ||^2 zeta := 1/c;

IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations inner_tol := RsdCG*||chi||*||p||*sqrt(zeta)

q := A(inner_tol) p

a[k] := |r[k-1]|**2 / <q,p[k]> ; Psi[k] += a[k] p[k] ; New solution std::vector r[k] -= a[k] q ; New residual c := <r , r> zeta := zeta + 1/c

b[k+1] := |r[k]|**2 / |r[k-1]|**2 = c/cp; p[k+1] := r[k] + b[k+1] p[k]; New direction

cp := c IF ( c < RsdCG^2 || chi || ) RETURN

Arguments:

Parameters
ALinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGCG residual accuracy (Read)
MaxCGMaximum CG iterations (Read)
n_countNumber of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > ap Temporary for M.p

Subroutines:


  • A Apply matrix M or M to std::vector

Operations:

2 A + 2 Nc Ns + N_Count ( 2 A + 10 Nc Ns )

Definition at line 144 of file inv_rel_cg1.cc.

References Chroma::A(), Chroma::chi(), Chroma::InvRelCG1_a(), Chroma::MaxCG, Chroma::n_count, Chroma::psi, and Chroma::RsdCG.

Referenced by Chroma::InvRelGMRESR_CG_a(), and Chroma::Ovlap4DQprop::operator()().

◆ InvRelCG1_a()

template<typename T >
void Chroma::InvRelCG1_a ( const LinearOperator< T > &  A,
const T chi,
T psi,
const Real &  RsdCG,
int  MaxCG,
int &  n_count 
)

Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations

        Chi  =  A . Psi

Definition at line 20 of file inv_rel_cg1.cc.

References Chroma::a, Chroma::A(), Chroma::b, Chroma::c, Chroma::chi(), Chroma::cp, Chroma::d, Chroma::END_CODE(), Chroma::k, Chroma::MaxCG, Chroma::n_count, Chroma::p, Chroma::PLUS, Chroma::psi, q, Chroma::r, Chroma::rsd_sq, Chroma::RsdCG, Chroma::s(), Chroma::START_CODE(), Chroma::tmp, tmp2, and Chroma::zero.

Referenced by Chroma::InvRelCG1().

◆ InvRelCG2()

template<>
void Chroma::InvRelCG2 ( const LinearOperator< T > &  M,
const T chi,
T psi,
const Real &  RsdCG,
int  MaxCG,
int &  n_count 
)

Relaxed Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations

        Chi  =  M^{dag}M . Psi

Algorithm:

Psi[0] := initial guess; Linear interpolation (argument) r[0] := Chi - M^dag . M . Psi[0] ; Initial residual p[1] := r[0] ; Initial direction c = cp := || r[0] ||^2
zeta := 1/c;

IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged?

FOR k FROM 1 TO MaxCG DO CG iterations

inner_tol := RsdCG*||chi||*||p||*sqrt(zeta)/2;
q := M^{dag}(tol) M(tol) p;
a[k] := c / <q , p>
Psi[k] += a[k] p[k] ;             New solution std::vector
r[k] -= a[k] q;                        New residual
c := || r[k]^2 ||
zeta = zeta + 1/c;
b[k+1] := |r[k]|**2 / |r[k-1]|**2 = c/cp;
p[k+1] := r[k] + b[k+1] p[k];          New direction
cp := c;
IF |r[k]| <= RsdCG |Chi| THEN RETURN;  Converged?

Arguments:

Parameters
MApproxLinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGCG residual accuracy (Read)
MaxCGMaximum CG iterations (Read)
n_countNumber of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Mp Temporary for M.p

Subroutines:


  • A Apply matrix M or M to std::vector

Operations:

2 A + 2 Nc Ns + N_Count ( 2 A + 10 Nc Ns )

Definition at line 175 of file inv_rel_cg2.cc.

References Chroma::chi(), Chroma::InvRelCG2_a(), Chroma::MaxCG, Chroma::n_count, Chroma::psi, and Chroma::RsdCG.

Referenced by Chroma::Ovlap4DQprop::operator()().

◆ InvRelCG2_a()

template<typename T >
void Chroma::InvRelCG2_a ( const LinearOperator< T > &  M,
const T chi,
T psi,
const Real &  RsdCG,
int  MaxCG,
int &  n_count 
)

Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations

        Chi  =  A . Psi

where A = M^dag . M

Algorithm:

Psi[0] := initial guess; Linear interpolation (argument) r[0] := Chi - M^dag . M . Psi[0] ; Initial residual p[1] := r[0] ; Initial direction IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k-1]|**2 / <Mp[k],Mp[k]> ; Psi[k] += a[k] p[k] ; New solution std::vector r[k] -= a[k] M^dag . M . p[k] ; New residual IF |r[k]| <= RsdCG |Chi| THEN RETURN; Converged? b[k+1] := |r[k]|**2 / |r[k-1]|**2 ; p[k+1] := r[k] + b[k+1] p[k]; New direction

Arguments:

Parameters
MLinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGCG residual accuracy (Read)
MaxCGMaximum CG iterations (Read)
n_countNumber of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Mp Temporary for M.p

Subroutines:


  • A Apply matrix M or M to std::vector

Operations:

2 A + 2 Nc Ns + N_Count ( 2 A + 10 Nc Ns )

Definition at line 64 of file inv_rel_cg2.cc.

References Chroma::a, Chroma::b, Chroma::c, Chroma::chi(), Chroma::cp, Chroma::d, Chroma::k, Chroma::MaxCG, Chroma::MINUS, Chroma::mp(), Chroma::n_count, Chroma::p, Chroma::PLUS, Chroma::psi, Chroma::r, Chroma::rsd_sq, Chroma::RsdCG, Chroma::s(), Chroma::LinearOperator< T >::subset(), and Chroma::zero.

Referenced by Chroma::InvRelCG2().

◆ MInvCG() [1/4]

template<>
void Chroma::MInvCG ( const DiffLinearOperator< LatticeFermion, multi1d< LatticeColorMatrix >, multi1d< LatticeColorMatrix > > &  M,
const LatticeFermion &  chi,
multi1d< LatticeFermion > &  psi,
const multi1d< Real > &  shifts,
const multi1d< Real > &  RsdCG,
int  MaxCG,
int &  n_count 
)

◆ MInvCG() [2/4]

template<>
void Chroma::MInvCG ( const DiffLinearOperatorArray< LatticeFermion, multi1d< LatticeColorMatrix >, multi1d< LatticeColorMatrix > > &  M,
const multi1d< LatticeFermion > &  chi,
multi1d< multi1d< LatticeFermion > > &  psi,
const multi1d< Real > &  shifts,
const multi1d< Real > &  RsdCG,
int  MaxCG,
int &  n_count 
)

Definition at line 453 of file minvcg_array.cc.

References chi, MaxCG, Chroma::MInvCG_a(), n_count, psi, and Chroma::RsdCG.

◆ MInvCG() [3/4]

template<>
void Chroma::MInvCG ( const LinearOperator< LatticeFermion > &  M,
const LatticeFermion &  chi,
multi1d< LatticeFermion > &  psi,
const multi1d< Real > &  shifts,
const multi1d< Real > &  RsdCG,
int  MaxCG,
int &  n_count 
)

◆ MInvCG() [4/4]

template<>
void Chroma::MInvCG ( const LinearOperatorArray< LatticeFermion > &  M,
const multi1d< LatticeFermion > &  chi,
multi1d< multi1d< LatticeFermion > > &  psi,
const multi1d< Real > &  shifts,
const multi1d< Real > &  RsdCG,
int  MaxCG,
int &  n_count 
)

Definition at line 439 of file minvcg_array.cc.

References chi, MaxCG, Chroma::MInvCG_a(), n_count, psi, and Chroma::RsdCG.

◆ MInvCG2() [1/5]

void Chroma::MInvCG2 ( const DiffLinearOperator< LatticeFermionF, multi1d< LatticeColorMatrixF >, multi1d< LatticeColorMatrixF > > &  M,
const LatticeFermionF &  chi,
multi1d< LatticeFermionF > &  psi,
const multi1d< RealF > &  shifts,
const multi1d< RealF > &  RsdCG,
int  MaxCG,
int &  n_count 
)

◆ MInvCG2() [2/5]

template<typename T , typename P , typename Q >
void Chroma::MInvCG2 ( const DiffLinearOperator< T, P, Q > &  M,
const T chi,
multi1d< T > &  psi,
const multi1d< Real > &  shifts,
const multi1d< Real > &  RsdCG,
int  MaxCG,
int &  n_count 
)

◆ MInvCG2() [3/5]

void Chroma::MInvCG2 ( const LinearOperator< LatticeFermionD > &  M,
const LatticeFermionD &  chi,
multi1d< LatticeFermionD > &  psi,
const multi1d< RealD > &  shifts,
const multi1d< RealD > &  RsdCG,
int  MaxCG,
int &  n_count 
)

◆ MInvCG2() [4/5]

void Chroma::MInvCG2 ( const LinearOperator< LatticeFermionF > &  M,
const LatticeFermionF &  chi,
multi1d< LatticeFermionF > &  psi,
const multi1d< RealF > &  shifts,
const multi1d< RealF > &  RsdCG,
int  MaxCG,
int &  n_count 
)

◆ MInvCG2() [5/5]

template<typename T >
void Chroma::MInvCG2 ( const LinearOperator< T > &  M,
const T chi,
multi1d< T > &  psi,
const multi1d< Real > &  shifts,
const multi1d< Real > &  RsdCG,
int  MaxCG,
int &  n_count 
)

◆ MInvCG2_a()

template<typename T , typename R >
void Chroma::MInvCG2_a ( const LinearOperator< T > &  M,
const T chi,
multi1d< T > &  psi,
const multi1d< R > &  shifts,
const multi1d< R > &  RsdCG,
int  MaxCG,
int &  n_count 
)

Multishift Conjugate-Gradient (CG1) algorithm for a Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations Method used is described in Jegerlehner, hep-lat/9708029

We are searching in a subspace orthogonal to the eigenvectors EigVec of A. The source chi is assumed to already be orthogonal!

        Chi  =  M^\dag M . Psi

Algorithm:

Psi[0] := 0; Zeroed r[0] := Chi; Initial residual p[1] := Chi ; Initial direction b[0] := |r[0]|**2 / <M p[0], M p[0]> ; z[0] := 1 / (1 - (shift - shift(0))*b) bs[0] := b[0] * z[0]
r[1] += b[k] M^\dag M . p[0] ; New residual Psi[1] = - b[k] p[k] ; Starting solution std::vector IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k]|**2 / |r[k-1]|**2 ; p[k] := r[k] + a[k] p[k-1]; New direction b[k+1] := |r[k]|**2 / <M p[k], Mp[k]> ; r[k+1] += b[k+1] M^\dag M . p[k] ; New residual Psi[k+1] -= b[k+1] p[k] ; New solution std::vector IF |[k+1]| <= RsdCG |Chi| THEN RETURN; Converged?

Arguments:

A Hermitian linear operator (Read) Chi Source (Read) Psi array of solutions (Write) shifts shifts of form A + mass (Read) RsdCG residual accuracy (Read/Write) n_count Number of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Ap Temporary for M.p

MaxCG Maximum number of CG iterations allowed

Subroutines: A Apply matrix hermitian A to std::vector

Definition at line 75 of file minvcg2.cc.

References Chroma::a, Chroma::b, Chroma::c, Chroma::chi(), Chroma::chi_norm, Chroma::cp, Chroma::d, Chroma::END_CODE(), Chroma::i, isz, iz, Chroma::k, Chroma::MaxCG, Chroma::MINUS, Chroma::n_count, Chroma::LinearOperator< T >::nFlops(), Chroma::p, Chroma::PLUS, Chroma::psi, Chroma::QDP_error_exit(), Chroma::r, Chroma::rsd_sq, Chroma::RsdCG, Chroma::s(), Chroma::START_CODE(), Chroma::LinearOperator< T >::subset(), z, and Chroma::zero.

Referenced by Chroma::MInvCG2().

◆ MInvCG2Accum() [1/2]

template<>
void Chroma::MInvCG2Accum ( const DiffLinearOperator< LatticeFermion, multi1d< LatticeColorMatrix >, multi1d< LatticeColorMatrix > > &  M,
const LatticeFermion &  chi,
LatticeFermion &  psi,
const Real &  norm,
const multi1d< Real > &  residues,
const multi1d< Real > &  poles,
const Real &  RsdCG,
int  MaxCG,
int &  n_count 
)

◆ MInvCG2Accum() [2/2]

template<>
void Chroma::MInvCG2Accum ( const LinearOperator< LatticeFermion > &  M,
const LatticeFermion &  chi,
LatticeFermion &  psi,
const Real &  norm,
const multi1d< Real > &  residues,
const multi1d< Real > &  poles,
const Real &  RsdCG,
int  MaxCG,
int &  n_count 
)

◆ MInvCG2Accum_a()

template<typename T >
void Chroma::MInvCG2Accum_a ( const LinearOperator< T > &  M,
const T chi,
T psi,
const Real &  norm,
const multi1d< Real > &  residues,
const multi1d< Real > &  poles,
const Real &  RsdCG,
int  MaxCG,
int &  n_count 
)

Multishift Conjugate-Gradient (CG1) algorithm for a Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations Method used is described in Jegerlehner, hep-lat/9708029

We are searching in a subspace orthogonal to the eigenvectors EigVec of A. The source chi is assumed to already be orthogonal!

        Chi  =  M^\dag M . Psi

Algorithm:

Psi[0] := 0; Zeroed r[0] := Chi; Initial residual p[1] := Chi ; Initial direction b[0] := |r[0]|**2 / <M p[0], M p[0]> ; z[0] := 1 / (1 - (shift - shift(0))*b) bs[0] := b[0] * z[0]
r[1] += b[k] M^\dag M . p[0] ; New residual Psi[1] = - b[k] p[k] ; Starting solution std::vector IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k]|**2 / |r[k-1]|**2 ; p[k] := r[k] + a[k] p[k-1]; New direction b[k+1] := |r[k]|**2 / <M p[k], Mp[k]> ; r[k+1] += b[k+1] M^\dag M . p[k] ; New residual Psi[k+1] -= b[k+1] p[k] ; New solution std::vector IF |[k+1]| <= RsdCG |Chi| THEN RETURN; Converged?

Arguments:

A Hermitian linear operator (Read) Chi Source (Read) Psi array of solutions (Write) shifts shifts of form A + mass (Read) RsdCG residual accuracy (Read/Write) n_count Number of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Ap Temporary for M.p

MaxCG Maximum number of CG iterations allowed

Subroutines: A Apply matrix hermitian A to std::vector

Definition at line 72 of file minvcg2_accum.cc.

References Chroma::a, Chroma::b, Chroma::c, Chroma::chi(), Chroma::chi_norm, Chroma::cp, Chroma::d, Chroma::END_CODE(), Chroma::i, isz, iz, Chroma::k, Chroma::MaxCG, Chroma::MINUS, Chroma::n_count, Chroma::LinearOperator< T >::nFlops(), norm, Chroma::p, Chroma::PLUS, Chroma::psi, Chroma::QDP_error_exit(), Chroma::r, Chroma::rsd_sq, Chroma::RsdCG, Chroma::s(), Chroma::START_CODE(), Chroma::LinearOperator< T >::subset(), Chroma::tmp, z, and Chroma::zero.

Referenced by Chroma::MInvCG2Accum().

◆ MInvCG_a() [1/2]

template<typename T , typename C >
void Chroma::MInvCG_a ( const C &  A,
const multi1d< T > &  chi,
multi1d< multi1d< T > > &  psi,
const multi1d< Real > &  shifts,
const multi1d< Real > &  RsdCG,
int  MaxCG,
int &  n_count 
)

Multishift Conjugate-Gradient (CG1) algorithm for a Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations Method used is described in Jegerlehner, hep-lat/9708029

We are searching in a subspace orthogonal to the eigenvectors EigVec of A. The source chi is assumed to already be orthogonal!

        Chi  =  A . Psi

Algorithm:

Psi[0] := 0; Zeroed r[0] := Chi; Initial residual p[1] := Chi ; Initial direction b[0] := |r[0]|**2 / <p[0],Ap[0]> ; z[0] := 1 / (1 - (shift - shift(0))*b) bs[0] := b[0] * z[0]
r[1] += b[k] A . p[0] ; New residual Psi[1] = - b[k] p[k] ; Starting solution std::vector IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k]|**2 / |r[k-1]|**2 ; p[k] := r[k] + a[k] p[k-1]; New direction b[k+1] := |r[k]|**2 / <p[k],Ap[k]> ; r[k+1] += b[k+1] A . p[k] ; New residual Psi[k+1] -= b[k+1] p[k] ; New solution std::vector IF |[k+1]| <= RsdCG |Chi| THEN RETURN; Converged?

Arguments:

A Hermitian linear operator (Read) Chi Source (Read) Psi array of solutions (Write) shifts shifts of form A + mass (Read) RsdCG residual accuracy (Read/Write) n_count Number of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Ap Temporary for M.p

MaxCG Maximum number of CG iterations allowed

Subroutines: A Apply matrix hermitian A to std::vector

Definition at line 74 of file minvcg_array.cc.

References Chroma::a, Chroma::A(), Chroma::b, c, chi, chi_norm, Chroma::cp, Chroma::d, END_CODE, i, isz, iz, Chroma::k, MaxCG, n, n_count, Chroma::p, Chroma::PLUS, psi, r(), Chroma::rsd_sq, Chroma::RsdCG, s, START_CODE, z, and Chroma::zero.

Referenced by Chroma::MInvCG().

◆ MInvCG_a() [2/2]

template<typename T >
void Chroma::MInvCG_a ( const LinearOperator< T > &  A,
const T chi,
multi1d< T > &  psi,
const multi1d< Real > &  shifts,
const multi1d< Real > &  RsdCG,
int  MaxCG,
int &  n_count 
)

Multishift Conjugate-Gradient (CG1) algorithm for a Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations Method used is described in Jegerlehner, hep-lat/9708029

We are searching in a subspace orthogonal to the eigenvectors EigVec of A. The source chi is assumed to already be orthogonal!

        Chi  =  A . Psi

Algorithm:

Psi[0] := 0; Zeroed r[0] := Chi; Initial residual p[1] := Chi ; Initial direction b[0] := |r[0]|**2 / <p[0],Ap[0]> ; z[0] := 1 / (1 - (shift - shift(0))*b) bs[0] := b[0] * z[0]
r[1] += b[k] A . p[0] ; New residual Psi[1] = - b[k] p[k] ; Starting solution std::vector IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k]|**2 / |r[k-1]|**2 ; p[k] := r[k] + a[k] p[k-1]; New direction b[k+1] := |r[k]|**2 / <p[k],Ap[k]> ; r[k+1] += b[k+1] A . p[k] ; New residual Psi[k+1] -= b[k+1] p[k] ; New solution std::vector IF |[k+1]| <= RsdCG |Chi| THEN RETURN; Converged?

Arguments:

A Hermitian linear operator (Read) Chi Source (Read) Psi array of solutions (Write) shifts shifts of form A + mass (Read) RsdCG residual accuracy (Read/Write) n_count Number of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Ap Temporary for M.p

MaxCG Maximum number of CG iterations allowed

Subroutines: A Apply matrix hermitian A to std::vector

Definition at line 72 of file minvcg.cc.

References Chroma::a, Chroma::A(), Chroma::b, Chroma::c, Chroma::chi(), Chroma::chi_norm, Chroma::cp, Chroma::d, Chroma::END_CODE(), Chroma::i, isz, iz, Chroma::k, Chroma::MaxCG, Chroma::n_count, Chroma::p, Chroma::PLUS, Chroma::psi, Chroma::QDP_error_exit(), Chroma::r, Chroma::rsd_sq, Chroma::RsdCG, Chroma::s(), Chroma::START_CODE(), z, and Chroma::zero.

Referenced by Chroma::MInvCG().

◆ MInvCGAccum() [1/2]

template<>
void Chroma::MInvCGAccum ( const DiffLinearOperatorArray< LatticeFermion, multi1d< LatticeColorMatrix >, multi1d< LatticeColorMatrix > > &  M,
const multi1d< LatticeFermion > &  chi,
multi1d< LatticeFermion > &  psi,
const Real &  norm,
const multi1d< Real > &  residues,
const multi1d< Real > &  poles,
const Real &  RsdCG,
const int  MaxCG,
int &  n_count 
)

Definition at line 477 of file minvcg_accumulate_array.cc.

References chi, MaxCG, Chroma::MInvCGAccum_a(), n_count, norm, psi, and Chroma::RsdCG.

◆ MInvCGAccum() [2/2]

template<>
void Chroma::MInvCGAccum ( const LinearOperatorArray< LatticeFermion > &  M,
const multi1d< LatticeFermion > &  chi,
multi1d< LatticeFermion > &  psi,
const Real &  norm,
const multi1d< Real > &  residues,
const multi1d< Real > &  poles,
const Real &  RsdCG,
const int  MaxCG,
int &  n_count 
)

◆ MInvCGAccum_a()

template<typename T , typename C >
void Chroma::MInvCGAccum_a ( const C &  A,
const multi1d< T > &  chi,
multi1d< T > &  psi,
const Real &  norm,
const multi1d< Real > &  residues,
const multi1d< Real > &  poles,
const Real &  RsdCG,
const int  MaxCG,
int &  n_count 
)

Multishift Conjugate-Gradient (CG1) algorithm for a Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations Method used is described in Jegerlehner, hep-lat/9708029

We are searching in a subspace orthogonal to the eigenvectors EigVec of A. The source chi is assumed to already be orthogonal!

        Chi  =  A . Psi

Algorithm:

Psi[0] := 0; Zeroed r[0] := Chi; Initial residual p[1] := Chi ; Initial direction b[0] := |r[0]|**2 / <p[0],Ap[0]> ; z[0] := 1 / (1 - (shift - shift(0))*b) bs[0] := b[0] * z[0]
r[1] += b[k] A . p[0] ; New residual Psi[1] = - b[k] p[k] ; Starting solution std::vector IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k]|**2 / |r[k-1]|**2 ; p[k] := r[k] + a[k] p[k-1]; New direction b[k+1] := |r[k]|**2 / <p[k],Ap[k]> ; r[k+1] += b[k+1] A . p[k] ; New residual Psi[k+1] -= b[k+1] p[k] ; New solution std::vector IF |[k+1]| <= RsdCG |Chi| THEN RETURN; Converged?

Arguments:

A Hermitian linear operator (Read) Chi Source (Read) Psi array of solutions (Write) shifts shifts of form A + mass (Read) RsdCG residual accuracy (Read/Write) n_count Number of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Ap Temporary for M.p

MaxCG Maximum number of CG iterations allowed

Subroutines: A Apply matrix hermitian A to std::vector

Definition at line 74 of file minvcg_accumulate_array.cc.

References Chroma::a, Chroma::A(), Chroma::b, c, chi, chi_norm, Chroma::cp, Chroma::d, END_CODE, i, isz, iz, Chroma::k, MaxCG, n, n_count, norm, Chroma::p, Chroma::PLUS, psi, r(), Chroma::rsd_sq, Chroma::RsdCG, s, START_CODE, tmp, z, and Chroma::zero.

Referenced by Chroma::MInvCGAccum().

◆ MInvMR() [1/2]

template<>
void Chroma::MInvMR ( const DiffLinearOperator< LatticeFermion, multi1d< LatticeColorMatrix >, multi1d< LatticeColorMatrix > > &  M,
const LatticeFermion &  chi,
multi1d< LatticeFermion > &  psi,
const multi1d< Real > &  shifts,
const multi1d< Real > &  RsdMR,
int  MaxMR,
int &  n_count 
)

Definition at line 272 of file minvmr.cc.

References Chroma::chi(), Chroma::MInvMR_a(), Chroma::n_count, Chroma::psi, and RsdMR.

◆ MInvMR() [2/2]

template<>
void Chroma::MInvMR ( const LinearOperator< LatticeFermion > &  M,
const LatticeFermion &  chi,
multi1d< LatticeFermion > &  psi,
const multi1d< Real > &  shifts,
const multi1d< Real > &  RsdCG,
int  MaxCG,
int &  n_count 
)

Multishift Conjugate-Gradient (CG1) algorithm for a Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations Method used is described in Jegerlehner, hep-lat/9708029

We are searching in a subspace orthogonal to the eigenvectors EigVec of A. The source chi is assumed to already be orthogonal!

        Chi  =  A . Psi

Algorithm:

Psi[0] := 0; Zeroed r[0] := Chi; Initial residual p[1] := Chi ; Initial direction b[0] := |r[0]|**2 / <p[0],Ap[0]> ; z[0] := 1 / (1 - (shift - shift(0))*b) bs[0] := b[0] * z[0]
r[1] += b[k] A . p[0] ; New residual Psi[1] = - b[k] p[k] ; Starting solution std::vector IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k]|**2 / |r[k-1]|**2 ; p[k] := r[k] + a[k] p[k-1]; New direction b[k+1] := |r[k]|**2 / <p[k],Ap[k]> ; r[k+1] += b[k+1] A . p[k] ; New residual Psi[k+1] -= b[k+1] p[k] ; New solution std::vector IF |[k+1]| <= RsdCG |Chi| THEN RETURN; Converged?

Arguments:

A Hermitian linear operator (Read) Chi Source (Read) Psi array of solutions (Write) shifts shifts of form A + mass (Read) RsdCG residual accuracy (Read/Write) n_count Number of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Ap Temporary for M.p

MaxCG Maximum number of CG iterations allowed

Subroutines: A Apply matrix hermitian A to std::vector

Definition at line 258 of file minvmr.cc.

References Chroma::chi(), Chroma::MInvMR_a(), Chroma::n_count, Chroma::psi, and RsdMR.

Referenced by Chroma::OverlapFermActBase::multiQprop(), and Chroma::LinOpMultiSysSolverMR< T >::operator()().

◆ MInvMR_a()

template<typename T >
void Chroma::MInvMR_a ( const LinearOperator< T > &  A,
const T chi,
multi1d< T > &  psi,
const multi1d< Real > &  shifts,
const multi1d< Real > &  RsdMR,
int  MaxMR,
int &  n_count 
)

Multishift Conjugate-Gradient (CG1) algorithm for a Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations Method used is described in Jegerlehner, hep-lat/9708029

We are searching in a subspace orthogonal to the eigenvectors EigVec of A. The source chi is assumed to already be orthogonal!

        Chi  =  A . Psi

Algorithm:

Psi[0] := 0; Zeroed r[0] := Chi; Initial residual p[1] := Chi ; Initial direction b[0] := |r[0]|**2 / <p[0],Ap[0]> ; z[0] := 1 / (1 - (shift - shift(0))*b) bs[0] := b[0] * z[0]
r[1] += b[k] A . p[0] ; New residual Psi[1] = - b[k] p[k] ; Starting solution std::vector IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k]|**2 / |r[k-1]|**2 ; p[k] := r[k] + a[k] p[k-1]; New direction b[k+1] := |r[k]|**2 / <p[k],Ap[k]> ; r[k+1] += b[k+1] A . p[k] ; New residual Psi[k+1] -= b[k+1] p[k] ; New solution std::vector IF |[k+1]| <= RsdCG |Chi| THEN RETURN; Converged?

Arguments:

A Hermitian linear operator (Read) Chi Source (Read) Psi array of solutions (Write) shifts shifts of form A + mass (Read) RsdCG residual accuracy (Read/Write) n_count Number of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Ap Temporary for M.p

MaxCG Maximum number of CG iterations allowed

Subroutines: A Apply matrix hermitian A to std::vector

Definition at line 74 of file minvmr.cc.

References Chroma::a, Chroma::A(), Chroma::c, Chroma::cb, Chroma::chi(), Chroma::chi_norm, Chroma::cp, Chroma::d, Chroma::END_CODE(), FILL(), Chroma::GramSchm(), Chroma::i, Chroma::InlinePropAndMatElemDistillation2Env::local::innerProduct(), isz, Chroma::k, Chroma::mass, Chroma::n_count, Chroma::Ncb, NOperEig, Chroma::PLUS, PRINTF(), Chroma::psi, Chroma::QDP_error_exit(), Chroma::r, RsdMR, Chroma::s(), Chroma::START_CODE(), and Chroma::zero.

Referenced by Chroma::MInvMR().

◆ MInvRelCG_a()

template<typename T >
void Chroma::MInvRelCG_a ( const LinearOperator< T > &  A,
const T chi,
multi1d< T > &  psi,
const multi1d< Real > &  shifts,
const multi1d< Real > &  RsdCG,
int  MaxCG,
int &  n_count 
)

Multishift Conjugate-Gradient (CG1) algorithm for a Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations Method used is described in Jegerlehner, hep-lat/9708029

We are searching in a subspace orthogonal to the eigenvectors EigVec of A. The source chi is assumed to already be orthogonal!

        Chi  =  A . Psi

Algorithm:

Psi[0] := 0; Zeroed r[0] := Chi; Initial residual p[1] := Chi ; Initial direction b[0] := |r[0]|**2 / <p[0],Ap[0]> ; z[0] := 1 / (1 - (shift - shift(0))*b) bs[0] := b[0] * z[0]
r[1] += b[k] A . p[0] ; New residual Psi[1] = - b[k] p[k] ; Starting solution std::vector IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k]|**2 / |r[k-1]|**2 ; p[k] := r[k] + a[k] p[k-1]; New direction b[k+1] := |r[k]|**2 / <p[k],Ap[k]> ; r[k+1] += b[k+1] A . p[k] ; New residual Psi[k+1] -= b[k+1] p[k] ; New solution std::vector IF |[k+1]| <= RsdCG |Chi| THEN RETURN; Converged?

Arguments:

A Hermitian linear operator (Read) Chi Source (Read) Psi array of solutions (Write) shifts shifts of form A + mass (Read) RsdCG residual accuracy (Read/Write) n_count Number of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Ap Temporary for M.p

MaxCG Maximum number of CG iterations allowed

Subroutines: A Apply matrix hermitian A to std::vector

Definition at line 70 of file minv_rel_cg.cc.

References Chroma::a, Chroma::A(), Chroma::b, Chroma::c, Chroma::chi(), Chroma::chi_norm, Chroma::cp, Chroma::d, Chroma::END_CODE(), Chroma::i, Chroma::InlinePropAndMatElemDistillation2Env::local::innerProduct(), isz, iz, Chroma::k, Chroma::MaxCG, Chroma::n_count, Chroma::p, Chroma::PLUS, Chroma::psi, Chroma::QDP_error_exit(), Chroma::r, Chroma::rsd_sq, Chroma::RsdCG, Chroma::s(), Chroma::START_CODE(), z, and Chroma::zero.

Referenced by Chroma::MInvRelCG().

◆ normGramSchmidt() [1/6]

void Chroma::normGramSchmidt ( multi1d< LatticeFermionD > &  vec,
int  f,
int  t,
const Subset &  sub 
)

Definition at line 88 of file norm_gram_schm.cc.

References t.

◆ normGramSchmidt() [2/6]

void Chroma::normGramSchmidt ( multi1d< LatticeFermionF > &  vec,
int  f,
int  t,
const Subset &  sub 
)

Gram-Schmidt with normalization.

Definition at line 80 of file norm_gram_schm.cc.

References t.

Referenced by Chroma::InvEigCG2Env::InvEigCG2_T(), and Chroma::InvEigCG2ArrayEnv::old_InvEigCG2_T().

◆ normGramSchmidt() [3/6]

void Chroma::normGramSchmidt ( multi1d< LatticeStaggeredFermionD > &  vec,
int  f,
int  t,
const Subset &  sub 
)

Definition at line 104 of file norm_gram_schm.cc.

References t.

◆ normGramSchmidt() [4/6]

void Chroma::normGramSchmidt ( multi1d< LatticeStaggeredFermionF > &  vec,
int  f,
int  t,
const Subset &  sub 
)

Definition at line 96 of file norm_gram_schm.cc.

References t.

◆ normGramSchmidt() [5/6]

void Chroma::normGramSchmidt ( multi2d< LatticeFermionD > &  vec,
int  f,
int  t,
const Subset &  sub 
)

Definition at line 122 of file norm_gram_schm.cc.

References Chroma::normGramSchmidtArray_T(), and t.

◆ normGramSchmidt() [6/6]

void Chroma::normGramSchmidt ( multi2d< LatticeFermionF > &  vec,
int  f,
int  t,
const Subset &  sub 
)

Definition at line 114 of file norm_gram_schm.cc.

References Chroma::normGramSchmidtArray_T(), and t.

◆ normGramSchmidt_T()

template<typename T >
void Chroma::normGramSchmidt_T ( multi1d< T > &  vec,
int  f,
int  t,
const Subset &  sub 
)

Gram-Schmidt with normalization.

Definition at line 14 of file norm_gram_schm.cc.

References Chroma::i, Chroma::in, Chroma::InlinePropAndMatElemDistillation2Env::local::innerProduct(), Chroma::k, and t.

◆ normGramSchmidtArray_T()

template<typename T >
void Chroma::normGramSchmidtArray_T ( multi2d< T > &  vec,
int  f,
int  t,
const Subset &  sub 
)

Gram-Schmidt with normalization.

Definition at line 45 of file norm_gram_schm.cc.

References Chroma::i, Chroma::in, Chroma::InlinePropAndMatElemDistillation2Env::local::innerProduct(), Chroma::k, Chroma::s(), and t.

Referenced by Chroma::normGramSchmidt().

◆ read() [1/10]

void Chroma::read ( XMLReader &  xml,
const std::string &  path,
MultiSysSolverCGParams param 
)

◆ read() [2/10]

void Chroma::read ( XMLReader &  xml,
const std::string &  path,
MultiSysSolverMRParams param 
)

◆ read() [3/10]

void Chroma::read ( XMLReader &  xml,
const std::string &  path,
SysSolverBiCGStabParams param 
)

◆ read() [4/10]

void Chroma::read ( XMLReader &  xml,
const std::string &  path,
SysSolverCGParams param 
)

◆ read() [5/10]

void Chroma::read ( XMLReader &  xml,
const std::string &  path,
SysSolverEigCGParams param 
)

◆ read() [6/10]

void Chroma::read ( XMLReader &  xml,
const std::string &  path,
SysSolverFGMRESDRParams p 
)

◆ read() [7/10]

void Chroma::read ( XMLReader &  xml,
const std::string &  path,
SysSolverMRParams param 
)

◆ read() [8/10]

void Chroma::read ( XMLReader &  xml,
const std::string &  path,
SysSolverOptEigBiCGParams param 
)

◆ read() [9/10]

void Chroma::read ( XMLReader &  xml,
const std::string &  path,
SysSolverOptEigCGParams param 
)

◆ read() [10/10]

void Chroma::read ( XMLReader &  xml,
const std::string &  path,
SysSolverQOPMGParams param 
)

◆ WlInvCG2()

void WlInvCG2 ( const LinearOperator M,
const LatticeFermion &  chi,
LatticeFermion &  psi,
const Real &  RsdCG,
int  MaxCG,
int &  n_count 
)

Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.

This subroutine uses the Conjugate Gradient (CG) algorithm to find the solution of the set of linear equations

        Chi  =  A . Psi

where A = M^dag . M

Algorithm:

Psi[0] := initial guess; Linear interpolation (argument) r[0] := Chi - M^dag . M . Psi[0] ; Initial residual p[1] := r[0] ; Initial direction IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged? FOR k FROM 1 TO MaxCG DO CG iterations a[k] := |r[k-1]|**2 / <Mp[k],Mp[k]> ; Psi[k] += a[k] p[k] ; New solution std::vector r[k] -= a[k] M^dag . M . p[k] ; New residual IF |r[k]| <= RsdCG |Chi| THEN RETURN; Converged? b[k+1] := |r[k]|**2 / |r[k-1]|**2 ; p[k+1] := r[k] + b[k+1] p[k]; New direction

Arguments:

Parameters
MLinear Operator (Read)
chiSource (Read)
psiSolution (Modify)
RsdCGCG residual accuracy (Read)
MaxCGMaximum CG iterations (Read)
n_countNumber of CG iteration (Write)

Local Variables:

p Direction std::vector r Residual std::vector cp | r[k] |**2 c | r[k-1] |**2 k CG iteration counter a a[k] b b[k+1] d < p[k], A.p[k] > Mp Temporary for M.p

Subroutines:


  • A Apply matrix M or M to std::vector

Operations:

2 A + 2 Nc Ns + N_Count ( 2 A + 10 Nc Ns )

Definition at line 69 of file t_wlinvcg.cc.

References Chroma::a, Chroma::b, Chroma::c, Chroma::chi(), Chroma::cp, Chroma::d, Chroma::k, Chroma::MaxCG, Chroma::MINUS, Chroma::mp(), Chroma::n_count, Chroma::p, Chroma::PLUS, Chroma::psi, Chroma::QDP_error_exit(), Chroma::r, Chroma::rsd_sq, Chroma::RsdCG, Chroma::s(), and Chroma::LinearOperator< T >::subset().

◆ write() [1/10]

void Chroma::write ( XMLWriter &  xml,
const std::string &  path,
const MultiSysSolverCGParams param 
)

◆ write() [2/10]

void Chroma::write ( XMLWriter &  xml,
const std::string &  path,
const MultiSysSolverMRParams param 
)

◆ write() [3/10]

void Chroma::write ( XMLWriter &  xml,
const std::string &  path,
const SysSolverBiCGStabParams param 
)

◆ write() [4/10]

void Chroma::write ( XMLWriter &  xml,
const std::string &  path,
const SysSolverCGParams param 
)

◆ write() [5/10]

void Chroma::write ( XMLWriter &  xml,
const std::string &  path,
const SysSolverEigCGParams param 
)

◆ write() [6/10]

void Chroma::write ( XMLWriter &  xml,
const std::string &  path,
const SysSolverFGMRESDRParams p 
)

◆ write() [7/10]

void Chroma::write ( XMLWriter &  xml,
const std::string &  path,
const SysSolverMRParams param 
)

◆ write() [8/10]

void Chroma::write ( XMLWriter &  xml,
const std::string &  path,
const SysSolverOptEigBiCGParams param 
)

◆ write() [9/10]

void Chroma::write ( XMLWriter &  xml,
const std::string &  path,
const SysSolverOptEigCGParams param 
)

◆ write() [10/10]

void Chroma::write ( XMLWriter &  xml,
const std::string &  path,
const SysSolverQOPMGParams param 
)