QUDA  1.0.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
Public Member Functions | Private Member Functions | Private Attributes | List of all members
quda::BiCGstabL Class Reference

#include <invert_quda.h>

Inheritance diagram for quda::BiCGstabL:
Inheritance graph
[legend]
Collaboration diagram for quda::BiCGstabL:
Collaboration graph
[legend]

Public Member Functions

 BiCGstabL (DiracMatrix &mat, DiracMatrix &matSloppy, SolverParam &param, TimeProfile &profile)
 
virtual ~BiCGstabL ()
 
void operator() (ColorSpinorField &out, ColorSpinorField &in)
 
- Public Member Functions inherited from quda::Solver
 Solver (SolverParam &param, TimeProfile &profile)
 
virtual ~Solver ()
 
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
 

Private Member Functions

int reliable (double &rNorm, double &maxrx, double &maxrr, const double &r2, const double &delta)
 Current residual, in BiCG language. More...
 
void computeTau (Complex **tau, double *sigma, std::vector< ColorSpinorField *> r, int begin, int size, int j)
 
void updateR (Complex **tau, std::vector< ColorSpinorField *> r, int begin, int size, int j)
 
void orthoDir (Complex **tau, double *sigma, std::vector< ColorSpinorField *> r, int j, int pipeline)
 
void updateUend (Complex *gamma, std::vector< ColorSpinorField *> u, int nKrylov)
 
void updateXRend (Complex *gamma, Complex *gamma_prime, Complex *gamma_prime_prime, std::vector< ColorSpinorField *> r, ColorSpinorField &x, int nKrylov)
 

Private Attributes

DiracMatrixmat
 
const DiracMatrixmatSloppy
 
int nKrylov
 
Complex rho0
 
Complex rho1
 
Complex alpha
 
Complex omega
 
Complex beta
 
Complexgamma
 
Complexgamma_prime
 
Complexgamma_prime_prime
 
Complex ** tau
 
double * sigma
 
ColorSpinorFieldr_fullp
 
ColorSpinorFieldyp
 Full precision residual. More...
 
ColorSpinorFieldtempp
 Full precision temporary. More...
 
std::vector< ColorSpinorField * > r
 Sloppy temporary vector. More...
 
std::vector< ColorSpinorField * > u
 
ColorSpinorFieldx_sloppy_saved_p
 
ColorSpinorFieldr0_saved_p
 Sloppy solution vector. More...
 
ColorSpinorFieldr_sloppy_saved_p
 Shadow residual, in BiCG language. More...
 
bool init
 
std::string solver_name
 

Additional Inherited Members

- Static Public Member Functions inherited from quda::Solver
static Solvercreate (SolverParam &param, 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 inherited from quda::Solver
EigenSolvereig_solve
 
bool deflate_init = false
 
std::vector< ColorSpinorField * > defl_tmp1
 
std::vector< ColorSpinorField * > defl_tmp2
 
- Protected Attributes inherited from quda::Solver
SolverParamparam
 
TimeProfileprofile
 
int node_parity
 

Detailed Description

Definition at line 755 of file invert_quda.h.

Constructor & Destructor Documentation

◆ BiCGstabL()

quda::BiCGstabL::BiCGstabL ( DiracMatrix mat,
DiracMatrix matSloppy,
SolverParam param,
TimeProfile profile 
)

Definition at line 255 of file inv_bicgstabl_quda.cpp.

References gamma, gamma_prime, gamma_prime_prime, nKrylov, r, sigma, solver_name, tau, and u.

◆ ~BiCGstabL()

quda::BiCGstabL::~BiCGstabL ( )
virtual

Member Function Documentation

◆ computeTau()

void quda::BiCGstabL::computeTau ( Complex **  tau,
double *  sigma,
std::vector< ColorSpinorField *>  r,
int  begin,
int  size,
int  j 
)
private

Internal routines for pipelined Gram-Schmidt. Made to not conflict with GCR's implementation.

Definition at line 21 of file inv_bicgstabl_quda.cpp.

References quda::blas::cDotProduct(), and quda::size.

Referenced by orthoDir().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator()()

void quda::BiCGstabL::operator() ( ColorSpinorField out,
ColorSpinorField in 
)
virtual

Implements quda::Solver.

Definition at line 322 of file inv_bicgstabl_quda.cpp.

References alpha, quda::dslash::aux_worker, beta, quda::BICGSTABL_UPDATE_R, quda::BICGSTABL_UPDATE_U, quda::blas::caxpby(), quda::blas::caxpy(), quda::blas::cDotProduct(), quda::blas::cDotProductNormA(), quda::SolverParam::compute_null_vector, quda::SolverParam::compute_true_res, quda::Solver::convergence(), quda::blas::copy(), quda::ColorSpinorParam::create, quda::ColorSpinorField::Create(), csParam, quda::SolverParam::delta, errorQuda, quda::blas::flops, quda::DiracMatrix::flops(), gamma, gamma_prime, gamma_prime_prime, quda::DiracMatrix::getStencilSteps(), getVerbosity(), quda::SolverParam::gflops, quda::SolverParam::heavy_quark_check, quda::blas::HeavyQuarkResidualNorm(), init, quda::SolverParam::is_preconditioner, quda::SolverParam::iter, quda::TimeProfile::Last(), mat, matSloppy, quda::SolverParam::maxiter, nKrylov, quda::blas::norm2(), omega, orthoDir(), quda::Solver::param, pipeline, quda::SolverParam::pipeline, quda::LatticeField::Precision(), quda::SolverParam::precision_sloppy, quda::SolverParam::preserve_source, printfQuda, quda::Solver::PrintStats(), quda::Solver::PrintSummary(), quda::Solver::profile, QUDA_COMPUTE_NULL_VECTOR_NO, QUDA_HEAVY_QUARK_RESIDUAL, QUDA_PRESERVE_SOURCE_NO, quda::QUDA_PROFILE_COMPUTE, quda::QUDA_PROFILE_EPILOGUE, quda::QUDA_PROFILE_FREE, quda::QUDA_PROFILE_PREAMBLE, QUDA_USE_INIT_GUESS_YES, QUDA_VERBOSE, QUDA_ZERO_FIELD_CREATE, r, r0_saved_p, r_fullp, r_sloppy_saved_p, reliable(), quda::SolverParam::residual_type, rho0, rho1, quda::SolverParam::secs, quda::ColorSpinorParam::setPrecision(), sigma, solver_name, quda::sqrt(), quda::Solver::stopping(), tau, tempp, quda::SolverParam::tol, quda::SolverParam::tol_hq, quda::SolverParam::true_res, quda::SolverParam::true_res_hq, u, updateUend(), updateXRend(), quda::SolverParam::use_init_guess, quda::SolverParam::use_sloppy_partial_accumulator, warningQuda, x_sloppy_saved_p, quda::blas::xmyNorm(), quda::blas::xpy(), quda::blas::xpyHeavyQuarkResidualNorm(), yp, and quda::blas::zero().

Here is the call graph for this function:

◆ orthoDir()

void quda::BiCGstabL::orthoDir ( Complex **  tau,
double *  sigma,
std::vector< ColorSpinorField *>  r,
int  j,
int  pipeline 
)
private

Definition at line 57 of file inv_bicgstabl_quda.cpp.

References quda::blas::caxpy(), quda::blas::caxpyDotzy(), quda::blas::cDotProduct(), computeTau(), pipeline, and updateR().

Referenced by operator()().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ reliable()

int quda::BiCGstabL::reliable ( double &  rNorm,
double &  maxrx,
double &  maxrr,
const double &  r2,
const double &  delta 
)
private

Current residual, in BiCG language.

Internal routine for reliable updates. Made to not conflict with BiCGstab's implementation.

Definition at line 308 of file inv_bicgstabl_quda.cpp.

References quda::sqrt(), and updateR().

Referenced by operator()().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ updateR()

void quda::BiCGstabL::updateR ( Complex **  tau,
std::vector< ColorSpinorField *>  r,
int  begin,
int  size,
int  j 
)
private

Definition at line 40 of file inv_bicgstabl_quda.cpp.

References quda::blas::caxpy(), and quda::size.

Referenced by orthoDir(), and reliable().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ updateUend()

void quda::BiCGstabL::updateUend ( Complex gamma,
std::vector< ColorSpinorField *>  u,
int  nKrylov 
)
private

Definition at line 108 of file inv_bicgstabl_quda.cpp.

References quda::blas::caxpy(), and nKrylov.

Referenced by operator()().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ updateXRend()

void quda::BiCGstabL::updateXRend ( Complex gamma,
Complex gamma_prime,
Complex gamma_prime_prime,
std::vector< ColorSpinorField *>  r,
ColorSpinorField x,
int  nKrylov 
)
private

Definition at line 125 of file inv_bicgstabl_quda.cpp.

References quda::blas::caxpyBxpz(), and nKrylov.

Referenced by operator()().

Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ alpha

Complex quda::BiCGstabL::alpha
private

Definition at line 767 of file invert_quda.h.

Referenced by operator()().

◆ beta

Complex quda::BiCGstabL::beta
private

Definition at line 767 of file invert_quda.h.

Referenced by operator()().

◆ gamma

Complex* quda::BiCGstabL::gamma
private

Definition at line 768 of file invert_quda.h.

Referenced by BiCGstabL(), operator()(), and ~BiCGstabL().

◆ gamma_prime

Complex * quda::BiCGstabL::gamma_prime
private

Definition at line 768 of file invert_quda.h.

Referenced by BiCGstabL(), operator()(), and ~BiCGstabL().

◆ gamma_prime_prime

Complex * quda::BiCGstabL::gamma_prime_prime
private

Definition at line 768 of file invert_quda.h.

Referenced by BiCGstabL(), operator()(), and ~BiCGstabL().

◆ init

bool quda::BiCGstabL::init
private

Solver uses lazy allocation: this flag determines whether we have allocated or not.

Definition at line 805 of file invert_quda.h.

Referenced by operator()(), and ~BiCGstabL().

◆ mat

DiracMatrix& quda::BiCGstabL::mat
private

Definition at line 758 of file invert_quda.h.

Referenced by operator()().

◆ matSloppy

const DiracMatrix& quda::BiCGstabL::matSloppy
private

Definition at line 759 of file invert_quda.h.

Referenced by operator()().

◆ nKrylov

int quda::BiCGstabL::nKrylov
private

The size of the Krylov space that BiCGstabL uses.

Definition at line 764 of file invert_quda.h.

Referenced by BiCGstabL(), operator()(), updateUend(), updateXRend(), and ~BiCGstabL().

◆ omega

Complex quda::BiCGstabL::omega
private

Definition at line 767 of file invert_quda.h.

Referenced by operator()().

◆ r

std::vector<ColorSpinorField*> quda::BiCGstabL::r
private

Sloppy temporary vector.

Definition at line 778 of file invert_quda.h.

Referenced by BiCGstabL(), operator()(), and ~BiCGstabL().

◆ r0_saved_p

ColorSpinorField* quda::BiCGstabL::r0_saved_p
private

Sloppy solution vector.

Definition at line 783 of file invert_quda.h.

Referenced by operator()(), and ~BiCGstabL().

◆ r_fullp

ColorSpinorField* quda::BiCGstabL::r_fullp
private

Definition at line 774 of file invert_quda.h.

Referenced by operator()(), and ~BiCGstabL().

◆ r_sloppy_saved_p

ColorSpinorField* quda::BiCGstabL::r_sloppy_saved_p
private

Shadow residual, in BiCG language.

Definition at line 784 of file invert_quda.h.

Referenced by operator()(), and ~BiCGstabL().

◆ rho0

Complex quda::BiCGstabL::rho0
private

Definition at line 767 of file invert_quda.h.

Referenced by operator()().

◆ rho1

Complex quda::BiCGstabL::rho1
private

Definition at line 767 of file invert_quda.h.

Referenced by operator()().

◆ sigma

double* quda::BiCGstabL::sigma
private

Definition at line 770 of file invert_quda.h.

Referenced by BiCGstabL(), operator()(), and ~BiCGstabL().

◆ solver_name

std::string quda::BiCGstabL::solver_name
private

Definition at line 807 of file invert_quda.h.

Referenced by BiCGstabL(), and operator()().

◆ tau

Complex** quda::BiCGstabL::tau
private

Definition at line 769 of file invert_quda.h.

Referenced by BiCGstabL(), operator()(), and ~BiCGstabL().

◆ tempp

ColorSpinorField* quda::BiCGstabL::tempp
private

Full precision temporary.

Definition at line 777 of file invert_quda.h.

Referenced by operator()(), and ~BiCGstabL().

◆ u

std::vector<ColorSpinorField*> quda::BiCGstabL::u
private

Definition at line 779 of file invert_quda.h.

Referenced by BiCGstabL(), operator()(), and ~BiCGstabL().

◆ x_sloppy_saved_p

ColorSpinorField* quda::BiCGstabL::x_sloppy_saved_p
private

Definition at line 782 of file invert_quda.h.

Referenced by operator()(), and ~BiCGstabL().

◆ yp

ColorSpinorField* quda::BiCGstabL::yp
private

Full precision residual.

Definition at line 775 of file invert_quda.h.

Referenced by operator()(), and ~BiCGstabL().


The documentation for this class was generated from the following files: