|
QUDA
0.9.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 | solve (ColorSpinorField &out, ColorSpinorField &in) |
| bool | convergence (const double &r2, const double &hq2, const double &r2_tol, const double &hq_tol) |
| bool | convergenceHQ (const double &r2, const double &hq2, const double &r2_tol, const double &hq_tol) |
| bool | convergenceL2 (const double &r2, const double &hq2, const double &r2_tol, const double &hq_tol) |
| void | PrintStats (const char *, int k, const double &r2, const double &b2, const double &hq2) |
| void | PrintSummary (const char *name, int k, const double &r2, const double &b2) |
| virtual double | flops () const |
Static Public Member Functions | |
| static Solver * | create (SolverParam ¶m, DiracMatrix &mat, DiracMatrix &matSloppy, DiracMatrix &matPrecon, TimeProfile &profile) |
| static double | stopping (const double &tol, const double &b2, QudaResidualType residual_type) |
Protected Attributes | |
| SolverParam & | param |
| TimeProfile & | profile |
Definition at line 325 of file invert_quda.h.
|
inline |
Definition at line 332 of file invert_quda.h.
|
inlinevirtual |
Definition at line 333 of file invert_quda.h.
| bool quda::Solver::convergence | ( | const double & | r2, |
| const double & | hq2, | ||
| const double & | r2_tol, | ||
| const double & | hq_tol | ||
| ) |
Test for solver convergence
| r2 | L2 norm squared of the residual |
| hq2 | Heavy quark residual |
| r2_tol | Solver L2 tolerance |
| hq_tol | Solver heavy-quark tolerance |
Definition at line 139 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::IncEigCG::eigCGsolve(), quda::CG::operator()(), quda::MPCG::operator()(), quda::PreconCG::operator()(), quda::BiCGstab::operator()(), quda::SimpleBiCGstab::operator()(), quda::MPBiCGstab::operator()(), quda::BiCGstabL::operator()(), quda::GCR::operator()(), quda::GMResDR::operator()(), and quda::CG::solve().

| bool quda::Solver::convergenceHQ | ( | const double & | r2, |
| const double & | hq2, | ||
| const double & | r2_tol, | ||
| const double & | hq_tol | ||
| ) |
Test for HQ solver convergence – ignore L2 residual
| r2 | L2 norm squared of the residual |
| hq2 | Heavy quark residual |
| r2_tol | Solver L2 tolerance |
| hq_tol | Solver heavy-quark tolerance |
Definition at line 156 of file solver.cpp.
References param, QUDA_HEAVY_QUARK_RESIDUAL, and quda::SolverParam::residual_type.
Referenced by quda::CG::operator()().

| bool quda::Solver::convergenceL2 | ( | const double & | r2, |
| const double & | hq2, | ||
| const double & | r2_tol, | ||
| const double & | hq_tol | ||
| ) |
Test for L2 solver convergence – ignore HQ residual
| r2 | L2 norm squared of the residual |
| hq2 | Heavy quark residual |
| r2_tol | Solver L2 tolerance |
| hq_tol | Solver heavy-quark tolerance |
Definition at line 167 of file solver.cpp.
References param, QUDA_L2_ABSOLUTE_RESIDUAL, QUDA_L2_RELATIVE_RESIDUAL, and quda::SolverParam::residual_type.
Referenced by quda::CG::operator()().

|
static |
Solver factory
Definition at line 13 of file solver.cpp.
References errorQuda, quda::SolverParam::inv_type, quda::SolverParam::inv_type_precondition, mat(), quda::SolverParam::maxiter, quda::multigrid_solver::mg, param, quda::SolverParam::precision_precondition, quda::SolverParam::precision_sloppy, quda::SolverParam::preconditioner, profile, QUDA_BICGSTAB_INVERTER, QUDA_BICGSTABL_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::createSmoother(), quda::MG::generateNullVectors(), invertMultiSrcQuda(), invertQuda(), and quda::MG::MG().


|
inlinevirtual |
Return flops
Reimplemented in quda::MG.
Definition at line 399 of file invert_quda.h.
Referenced by quda::GCR::operator()().

|
pure virtual |
| void quda::Solver::PrintStats | ( | const char * | name, |
| int | k, | ||
| const double & | r2, | ||
| const double & | b2, | ||
| const double & | hq2 | ||
| ) |
Prints out the running statistics of the solver (requires a verbosity of QUDA_VERBOSE)
Definition at line 179 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::IncEigCG::eigCGsolve(), quda::CG::operator()(), quda::MPCG::operator()(), quda::PreconCG::operator()(), quda::BiCGstab::operator()(), quda::SimpleBiCGstab::operator()(), quda::MPBiCGstab::operator()(), quda::BiCGstabL::operator()(), quda::GCR::operator()(), quda::GMResDR::operator()(), and quda::CG::solve().


Prints out the summary of the solver convergence (requires a versbosity of QUDA_SUMMARIZE). Assumes SolverParam.true_res and SolverParam.true_res_hq has been set
Definition at line 194 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::IncEigCG::eigCGsolve(), quda::CG::operator()(), quda::CGNR::operator()(), quda::MPCG::operator()(), quda::BiCGstab::operator()(), quda::SimpleBiCGstab::operator()(), quda::MPBiCGstab::operator()(), quda::BiCGstabL::operator()(), quda::GCR::operator()(), quda::IncEigCG::operator()(), quda::GMResDR::operator()(), and quda::CG::solve().


|
virtual |
Reimplemented in quda::CG.
Definition at line 114 of file solver.cpp.
References fused_exterior_ndeg_tm_dslash_cuda_gen::i, in, quda::SolverParam::num_src, out, param, quda::SolverParam::true_res, quda::SolverParam::true_res_hq, quda::SolverParam::true_res_hq_offset, and quda::SolverParam::true_res_offset.
Referenced by quda::MG::generateNullVectors().

|
static |
Set the solver stopping condition
| b2 | L2 norm squared of the source vector |
Definition at line 122 of file solver.cpp.
References QUDA_L2_ABSOLUTE_RESIDUAL, QUDA_L2_RELATIVE_RESIDUAL, and tol.
Referenced by quda::IncEigCG::eigCGsolve(), quda::CG::operator()(), quda::MPCG::operator()(), quda::PreconCG::operator()(), quda::BiCGstab::operator()(), quda::SimpleBiCGstab::operator()(), quda::MPBiCGstab::operator()(), quda::BiCGstabL::operator()(), quda::GCR::operator()(), quda::MultiShiftCG::operator()(), and quda::CG::solve().

|
protected |
Definition at line 328 of file invert_quda.h.
Referenced by convergence(), convergenceHQ(), convergenceL2(), create(), quda::IncEigCG::eigCGsolve(), quda::GMResDR::FlexArnoldiProcedure(), quda::GCR::GCR(), quda::GMResDR::GMResDR(), quda::IncEigCG::IncEigCG(), quda::IncEigCG::initCGsolve(), quda::CG::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::SD::operator()(), quda::IncEigCG::operator()(), quda::GMResDR::operator()(), quda::PreconCG::PreconCG(), PrintStats(), PrintSummary(), quda::IncEigCG::RestartVT(), quda::GMResDR::RestartVZH(), solve(), quda::CG::solve(), quda::GMResDR::UpdateSolution(), quda::IncEigCG::UpdateVm(), quda::XSD::XSD(), quda::GCR::~GCR(), quda::GMResDR::~GMResDR(), quda::MR::~MR(), quda::SD::~SD(), and quda::XSD::~XSD().
|
protected |
Definition at line 329 of file invert_quda.h.
Referenced by create(), quda::GCR::GCR(), quda::CG::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::PreconCG::PreconCG(), quda::CG::solve(), quda::XSD::XSD(), quda::BiCGstab::~BiCGstab(), quda::BiCGstabL::~BiCGstabL(), quda::GCR::~GCR(), quda::MR::~MR(), quda::PreconCG::~PreconCG(), quda::SD::~SD(), and quda::XSD::~XSD().
1.8.14