QUDA
1.0.0
|
#include <invert_quda.h>
Public Member Functions | |
Solver (SolverParam ¶m, TimeProfile &profile) | |
virtual | ~Solver () |
virtual void | operator() (ColorSpinorField &out, ColorSpinorField &in)=0 |
virtual void | blocksolve (ColorSpinorField &out, ColorSpinorField &in) |
bool | convergence (double r2, double hq2, double r2_tol, double hq_tol) |
bool | convergenceHQ (double r2, double hq2, double r2_tol, double hq_tol) |
Test for HQ solver convergence – ignore L2 residual. More... | |
bool | convergenceL2 (double r2, double hq2, double r2_tol, double hq_tol) |
Test for L2 solver convergence – ignore HQ residual. More... | |
void | PrintStats (const char *name, int k, double r2, double b2, double hq2) |
Prints out the running statistics of the solver (requires a verbosity of QUDA_VERBOSE) More... | |
void | PrintSummary (const char *name, int k, double r2, double b2, double r2_tol, double hq_tol) |
Prints out the summary of the solver convergence (requires a verbosity of QUDA_SUMMARIZE). Assumes SolverParam.true_res and SolverParam.true_res_hq has been set. More... | |
void | constructDeflationSpace (const ColorSpinorField &meta, const DiracMatrix &mat, bool svd) |
Constructs the deflation space. More... | |
virtual double | flops () const |
Static Public Member Functions | |
static Solver * | create (SolverParam ¶m, DiracMatrix &mat, DiracMatrix &matSloppy, DiracMatrix &matPrecon, TimeProfile &profile) |
static double | stopping (double tol, double b2, QudaResidualType residual_type) |
Set the solver L2 stopping condition. More... | |
Public Attributes | |
EigenSolver * | eig_solve |
bool | deflate_init = false |
std::vector< ColorSpinorField * > | defl_tmp1 |
std::vector< ColorSpinorField * > | defl_tmp2 |
Protected Attributes | |
SolverParam & | param |
TimeProfile & | profile |
int | node_parity |
Definition at line 460 of file invert_quda.h.
quda::Solver::Solver | ( | SolverParam & | param, |
TimeProfile & | profile | ||
) |
Definition at line 13 of file solver.cpp.
References commCoords(), and node_parity.
|
virtual |
Definition at line 24 of file solver.cpp.
References eig_solve.
|
virtual |
Reimplemented in quda::CG.
Definition at line 198 of file solver.cpp.
References quda::ColorSpinorField::Component(), quda::SolverParam::num_src, param, quda::SolverParam::true_res, quda::SolverParam::true_res_hq, quda::SolverParam::true_res_hq_offset, and quda::SolverParam::true_res_offset.
Referenced by invertMultiSrcQuda().
void quda::Solver::constructDeflationSpace | ( | const ColorSpinorField & | meta, |
const DiracMatrix & | mat, | ||
bool | svd | ||
) |
Constructs the deflation space.
[in] | meta | A sample ColorSpinorField with which to instantiate the eigensolver |
[in] | mat | The operator to eigensolve |
[in] | Whether | to compute the SVD |
Definition at line 159 of file solver.cpp.
References quda::EigenSolver::computeSVD(), quda::EigenSolver::create(), quda::ColorSpinorParam::create, quda::ColorSpinorField::Create(), csParam, defl_tmp1, defl_tmp2, deflate_init, quda::SolverParam::eig_param, eig_solve, quda::SolverParam::evals, quda::SolverParam::evecs, QudaEigParam_s::nConv, param, quda::SolverParam::precision_sloppy, profile, QUDA_INVALID_PRECISION, quda::QUDA_PROFILE_INIT, QUDA_ZERO_FIELD_CREATE, and quda::ColorSpinorParam::setPrecision().
Referenced by quda::CACG::create(), quda::CAGCR::create(), quda::CG::operator()(), and quda::GCR::operator()().
bool quda::Solver::convergence | ( | double | r2, |
double | hq2, | ||
double | r2_tol, | ||
double | hq_tol | ||
) |
for solver convergence
[in] | r2 | L2 norm squared of the residual |
[in] | hq2 | Heavy quark residual |
[in] | r2_tol | Solver L2 tolerance |
[in] | hq_tol | Solver heavy-quark tolerance |
Definition at line 223 of file solver.cpp.
References param, QUDA_HEAVY_QUARK_RESIDUAL, QUDA_L2_ABSOLUTE_RESIDUAL, QUDA_L2_RELATIVE_RESIDUAL, and quda::SolverParam::residual_type.
Referenced by quda::CG::blocksolve(), quda::IncEigCG::eigCGsolve(), quda::CG::operator()(), quda::CG3::operator()(), quda::CG3NE::operator()(), quda::MPCG::operator()(), quda::PreconCG::operator()(), quda::BiCGstab::operator()(), quda::SimpleBiCGstab::operator()(), quda::MPBiCGstab::operator()(), quda::BiCGstabL::operator()(), quda::GCR::operator()(), quda::CACG::operator()(), quda::CAGCR::operator()(), and quda::GMResDR::operator()().
bool quda::Solver::convergenceHQ | ( | double | r2, |
double | hq2, | ||
double | r2_tol, | ||
double | hq_tol | ||
) |
Test for HQ solver convergence – ignore L2 residual.
[in] | r2 | L2 norm squared of the residual |
[in] | hq2 | Heavy quark residual |
[in] | r2_tol | Solver L2 tolerance |
Definition at line 237 of file solver.cpp.
References param, QUDA_HEAVY_QUARK_RESIDUAL, and quda::SolverParam::residual_type.
Referenced by quda::CG::blocksolve(), quda::CG::operator()(), quda::CG3::operator()(), and quda::CG3NE::operator()().
bool quda::Solver::convergenceL2 | ( | double | r2, |
double | hq2, | ||
double | r2_tol, | ||
double | hq_tol | ||
) |
Test for L2 solver convergence – ignore HQ residual.
[in] | r2 | L2 norm squared of the residual |
[in] | hq2 | Heavy quark residual |
[in] | r2_tol | Solver L2 tolerance |
[in] | hq_tol | Solver heavy-quark tolerance |
Definition at line 246 of file solver.cpp.
References param, QUDA_L2_ABSOLUTE_RESIDUAL, QUDA_L2_RELATIVE_RESIDUAL, and quda::SolverParam::residual_type.
Referenced by quda::CG::blocksolve(), and quda::CG::operator()().
|
static |
Solver factory
Definition at line 33 of file solver.cpp.
References errorQuda, quda::SolverParam::inv_type, quda::SolverParam::inv_type_precondition, quda::multigrid_solver::mg, quda::SolverParam::mg_instance, quda::SolverParam::precision_precondition, quda::SolverParam::precision_sloppy, quda::SolverParam::preconditioner, QUDA_BICGSTAB_INVERTER, QUDA_BICGSTABL_INVERTER, QUDA_CA_CG_INVERTER, QUDA_CA_CGNE_INVERTER, QUDA_CA_CGNR_INVERTER, QUDA_CA_GCR_INVERTER, QUDA_CG3_INVERTER, QUDA_CG3NE_INVERTER, QUDA_CG3NR_INVERTER, QUDA_CG_INVERTER, QUDA_CGNE_INVERTER, QUDA_CGNR_INVERTER, QUDA_EIGCG_INVERTER, QUDA_GCR_INVERTER, QUDA_GMRESDR_INVERTER, QUDA_INC_EIGCG_INVERTER, QUDA_MG_INVERTER, QUDA_MPBICGSTAB_INVERTER, QUDA_MPCG_INVERTER, QUDA_MR_INVERTER, QUDA_PCG_INVERTER, QUDA_SD_INVERTER, QUDA_XSD_INVERTER, and quda::report().
Referenced by quda::MG::createCoarseSolver(), quda::MG::createSmoother(), quda::MG::generateNullVectors(), invertMultiSrcQuda(), and invertQuda().
|
inlinevirtual |
Return flops
Reimplemented in quda::MG.
Definition at line 563 of file invert_quda.h.
Referenced by quda::GCR::operator()().
|
pure virtual |
Implemented in quda::GMResDR, quda::IncEigCG, quda::PreconditionedSolver, quda::XSD, quda::SD, quda::CAGCR, quda::CACGNR, quda::CACGNE, quda::CACG, quda::MR, quda::GCR, quda::BiCGstabL, quda::MPBiCGstab, quda::SimpleBiCGstab, quda::BiCGstab, quda::PreconCG, quda::MPCG, quda::CGNR, quda::CGNE, quda::CG3NE, quda::CG3, quda::CG, and quda::MG.
void quda::Solver::PrintStats | ( | const char * | name, |
int | k, | ||
double | r2, | ||
double | b2, | ||
double | hq2 | ||
) |
Prints out the running statistics of the solver (requires a verbosity of QUDA_VERBOSE)
[in] | name | Name of solver that called this |
[in] | k | iteration count |
[in] | r2 | L2 norm squared of the residual |
[in] | hq2 | Heavy quark residual |
Definition at line 256 of file solver.cpp.
References errorQuda, getVerbosity(), param, printfQuda, QUDA_HEAVY_QUARK_RESIDUAL, QUDA_VERBOSE, quda::SolverParam::residual_type, and quda::sqrt().
Referenced by quda::CG::blocksolve(), quda::IncEigCG::eigCGsolve(), quda::CG::operator()(), quda::CG3::operator()(), quda::CG3NE::operator()(), quda::MPCG::operator()(), quda::PreconCG::operator()(), quda::BiCGstab::operator()(), quda::SimpleBiCGstab::operator()(), quda::MPBiCGstab::operator()(), quda::BiCGstabL::operator()(), quda::GCR::operator()(), quda::CACG::operator()(), quda::CAGCR::operator()(), and quda::GMResDR::operator()().
void quda::Solver::PrintSummary | ( | const char * | name, |
int | k, | ||
double | r2, | ||
double | b2, | ||
double | r2_tol, | ||
double | hq_tol | ||
) |
Prints out the summary of the solver convergence (requires a verbosity of QUDA_SUMMARIZE). Assumes SolverParam.true_res and SolverParam.true_res_hq has been set.
[in] | name | Name of solver that called this |
[in] | k | iteration count |
[in] | r2 | L2 norm squared of the residual |
[in] | hq2 | Heavy quark residual |
[in] | r2_tol | Solver L2 tolerance |
[in] | hq_tol | Solver heavy-quark tolerance |
Definition at line 270 of file solver.cpp.
References quda::SolverParam::compute_true_res, getVerbosity(), param, printfQuda, QUDA_HEAVY_QUARK_RESIDUAL, QUDA_SUMMARIZE, quda::SolverParam::residual_type, quda::sqrt(), quda::SolverParam::true_res, and quda::SolverParam::true_res_hq.
Referenced by quda::CG::blocksolve(), quda::IncEigCG::eigCGsolve(), quda::CG::operator()(), quda::CG3::operator()(), quda::CG3NE::operator()(), quda::CGNE::operator()(), quda::CGNR::operator()(), quda::MPCG::operator()(), quda::BiCGstab::operator()(), quda::SimpleBiCGstab::operator()(), quda::MPBiCGstab::operator()(), quda::BiCGstabL::operator()(), quda::GCR::operator()(), quda::CACG::operator()(), quda::CACGNE::operator()(), quda::CACGNR::operator()(), quda::CAGCR::operator()(), quda::IncEigCG::operator()(), and quda::GMResDR::operator()().
|
static |
Set the solver L2 stopping condition.
[in] | Desired | solver tolerance |
[in] | b2 | L2 norm squared of the source vector |
[in] | residual_type | The type of residual we want to solve for |
Definition at line 206 of file solver.cpp.
References QUDA_L2_ABSOLUTE_RESIDUAL, QUDA_L2_RELATIVE_RESIDUAL, and tol.
Referenced by quda::CG::blocksolve(), quda::IncEigCG::eigCGsolve(), quda::CG::operator()(), quda::CG3::operator()(), quda::CG3NE::operator()(), quda::CGNE::operator()(), quda::CGNR::operator()(), quda::MPCG::operator()(), quda::PreconCG::operator()(), quda::BiCGstab::operator()(), quda::SimpleBiCGstab::operator()(), quda::MPBiCGstab::operator()(), quda::BiCGstabL::operator()(), quda::GCR::operator()(), quda::CACG::operator()(), quda::CACGNE::operator()(), quda::CACGNR::operator()(), quda::CAGCR::operator()(), and quda::MultiShiftCG::operator()().
std::vector<ColorSpinorField *> quda::Solver::defl_tmp1 |
Definition at line 547 of file invert_quda.h.
Referenced by constructDeflationSpace(), quda::CG::operator()(), quda::GCR::operator()(), quda::CACG::operator()(), quda::CAGCR::operator()(), quda::CACG::~CACG(), quda::CAGCR::~CAGCR(), quda::CG::~CG(), and quda::GCR::~GCR().
std::vector<ColorSpinorField *> quda::Solver::defl_tmp2 |
Definition at line 548 of file invert_quda.h.
Referenced by constructDeflationSpace(), quda::GCR::operator()(), quda::CACG::operator()(), quda::CAGCR::operator()(), quda::CACG::~CACG(), quda::CAGCR::~CAGCR(), quda::CG::~CG(), and quda::GCR::~GCR().
bool quda::Solver::deflate_init = false |
Definition at line 546 of file invert_quda.h.
Referenced by constructDeflationSpace(), quda::CACG::create(), quda::CAGCR::create(), quda::CG::operator()(), quda::GCR::operator()(), quda::CACG::~CACG(), quda::CAGCR::~CAGCR(), quda::CG::~CG(), and quda::GCR::~GCR().
EigenSolver* quda::Solver::eig_solve |
Deflation objects
Definition at line 545 of file invert_quda.h.
Referenced by constructDeflationSpace(), quda::MG::generateEigenVectors(), quda::CG::operator()(), quda::GCR::operator()(), quda::CACG::operator()(), quda::CAGCR::operator()(), and ~Solver().
|
protected |
Definition at line 465 of file invert_quda.h.
Referenced by quda::MR::operator()(), and Solver().
|
protected |
Definition at line 463 of file invert_quda.h.
Referenced by blocksolve(), quda::CG::blocksolve(), quda::CACG::compute_alpha(), quda::CACG::compute_beta(), constructDeflationSpace(), convergence(), quda::MultiShiftSolver::convergence(), convergenceHQ(), convergenceL2(), quda::CACG::create(), quda::CAGCR::create(), quda::IncEigCG::eigCGsolve(), quda::GMResDR::FlexArnoldiProcedure(), quda::IncEigCG::initCGsolve(), quda::CG::operator()(), quda::CG3::operator()(), quda::CG3NE::operator()(), quda::CGNE::operator()(), quda::CGNR::operator()(), quda::MPCG::operator()(), quda::PreconCG::operator()(), quda::BiCGstab::operator()(), quda::SimpleBiCGstab::operator()(), quda::MPBiCGstab::operator()(), quda::BiCGstabL::operator()(), quda::GCR::operator()(), quda::MR::operator()(), quda::CACG::operator()(), quda::CACGNE::operator()(), quda::CACGNR::operator()(), quda::CAGCR::operator()(), quda::SD::operator()(), quda::IncEigCG::operator()(), quda::GMResDR::operator()(), PrintStats(), PrintSummary(), quda::IncEigCG::RestartVT(), quda::GMResDR::RestartVZH(), quda::CAGCR::solve(), quda::GMResDR::UpdateSolution(), quda::IncEigCG::UpdateVm(), quda::CACG::~CACG(), quda::CAGCR::~CAGCR(), quda::CG::~CG(), quda::CG3::~CG3(), quda::CG3NE::~CG3NE(), quda::GCR::~GCR(), quda::GMResDR::~GMResDR(), quda::MR::~MR(), quda::SD::~SD(), and quda::XSD::~XSD().
|
protected |
Definition at line 464 of file invert_quda.h.
Referenced by quda::CG::blocksolve(), quda::CACG::compute_alpha(), quda::CACG::compute_beta(), constructDeflationSpace(), quda::CACG::create(), quda::CAGCR::create(), quda::CG::operator()(), quda::CG3::operator()(), quda::CG3NE::operator()(), quda::MPCG::operator()(), quda::PreconCG::operator()(), quda::BiCGstab::operator()(), quda::SimpleBiCGstab::operator()(), quda::MPBiCGstab::operator()(), quda::BiCGstabL::operator()(), quda::GCR::operator()(), quda::MR::operator()(), quda::CACG::operator()(), quda::CAGCR::operator()(), quda::CAGCR::solve(), quda::BiCGstab::~BiCGstab(), quda::BiCGstabL::~BiCGstabL(), quda::CACG::~CACG(), quda::CAGCR::~CAGCR(), quda::CG::~CG(), quda::GCR::~GCR(), quda::MR::~MR(), quda::PreconCG::~PreconCG(), quda::SD::~SD(), and quda::XSD::~XSD().