1.2.1. Implementation Classes

1.2.1.1. chase::Impl::ChASECPU

template<class T, typename MatrixType = chase::matrix::Matrix<T>>
class ChASECPU : public chase::ChaseBase<T>

CPU-based sequential Chase algorithm.

This class is responsible for solving generalized eigenvalue problems using CPU-based Lanczos iterations and other matrix operations.

Template Parameters:

T – The data type (e.g., float, double).

Public Functions

inline ChASECPU(std::size_t N, std::size_t nev, std::size_t nex, T *H, std::size_t ldh, T *V1, std::size_t ldv, chase::Base<T> *ritzv)

Constructs the ChASECPU object.

Initializes the matrices, vectors, and configuration necessary for the computation.

Parameters:
  • N – The size of the matrix (N x N).

  • nev – The number of eigenvalues to compute.

  • nex – The number of additional eigenvalues for extra space.

  • H – Pointer to the matrix of size N x N.

  • ldh – Leading dimension of H.

  • V1 – Pointer to the initial vector set of size N x (nev + nex).

  • ldv – Leading dimension of V1.

  • ritzv – Pointer to the Ritz values.

inline ChASECPU(std::size_t N, std::size_t nev, std::size_t nex, MatrixType *H, T *V1, std::size_t ldv, chase::Base<T> *ritzv)

Constructs the ChASECPU object.

Initializes the matrices, vectors, and configuration necessary for the computation.

Parameters:
  • N – The size of the matrix (N x N).

  • nev – The number of eigenvalues to compute.

  • nex – The number of additional eigenvalues for extra space.

  • H – Pointer to the matrix of size N x N.

  • ldh – Leading dimension of H.

  • V1 – Pointer to the initial vector set of size N x (nev + nex).

  • ldv – Leading dimension of V1.

  • ritzv – Pointer to the Ritz values.

ChASECPU(const ChASECPU&) = delete

Deleted copy constructor.

This class is not copyable, and the copy constructor is deleted to prevent object duplication.

inline ~ChASECPU()

Destructor for the ChASECPU class.

The destructor is defined to clean up the allocated memory (if any) and perform necessary cleanup tasks.

inline virtual std::size_t GetN() const override

Returns the size of the matrix.

This method returns the size of the matrix (number of rows or columns).

Returns:

The size of the matrix.

inline virtual std::size_t GetNev() override

Returns the number of eigenpairs to be computed.

This method returns the number of eigenpairs that need to be computed.

Returns:

The number of eigenpairs.

inline virtual std::size_t GetNex() override

Returns the value of nex, which is used for specifying the number of eigenpairs.

This method returns the size or number of eigenpairs required by the algorithm for the computation.

Returns:

The number of eigenpairs as a std::size_t value.

inline virtual chase::Base<T> *GetRitzv() override

Returns the computed Ritz values.

This method provides access to the Ritz values computed during the algorithm’s execution. These values are typically used for eigenvalue approximations and further analysis.

Returns:

A pointer to an array of Ritz values of type Base<T>.

inline virtual chase::Base<T> *GetResid() override

Returns the residuals of computed Ritz pairs.

This method returns the residuals associated with the Ritz values. Residuals are often used to measure the accuracy or convergence of the eigenvalue problem solution.

Returns:

A pointer to an array of residuals corresponding to the Ritz values of type Base<T>.

inline virtual ChaseConfig<T> &GetConfig() override

Returns the configuration parameters used in the algorithm.

This method provides access to the configuration object that contains various algorithm settings, parameters, and options.

Returns:

A reference to a ChaseConfig<T> class containing the configuration parameters.

inline virtual int get_nprocs() override

Returns the number of MPI processes used.

This method returns the number of MPI processes involved in the computation. If sequential ChASE is used, it returns 1, indicating a single process is running.

Returns:

The number of MPI processes used, represented as an integer.

inline void loadProblemFromFile(std::string filename)

Loads matrix data from a binary file.

This function reads the matrix data from a specified binary file and stores it in the internal matrix (Hmat_). The file is expected to contain the raw matrix data with a format that is compatible with the matrix’s internal representation.

Parameters:

filename – The path to the binary file containing the matrix data.

inline virtual void Output(std::string str) override

Print some intermediate infos during the solving procedure.

inline virtual bool checkSymmetryEasy() override

Checks for easy symmetry in the matrix.

This method performs a simple check for matrix symmetry.

Returns:

true if the matrix is symmetric, otherwise false.

inline virtual bool isSym() override

Returns if the matrix is Hermitian or not.

This method simply returns if the matrix is Hermitian or not

Returns:

true if the matrix is Hermitian, otherwise false.

inline virtual bool checkPseudoHermicityEasy() override

Checks for easy pseudo-hermicity in the matrix.

This method performs a simple check for matrix pseudo-hermicity

Returns:

true if the matrix is pseudo-hermicity, otherwise false.

inline virtual bool isPseudoHerm() override

Returns if the matrix is Pseudo-Hermitian or not.

This method simply returns if the matrix is Pseudo-Hermitian or not

Returns:

true if the matrix is Pseudo-Hermitian, otherwise false.

inline virtual void symOrHermMatrix(char uplo) override

Sets matrix type for symmetric or Hermitian matrix.

This method defines whether the matrix is symmetric or Hermitian based on the specified input.

Parameters:

uplo – Specifies whether the upper or lower triangular part of the matrix is stored.

inline virtual void Start() override

Indicates the start of an eigenproblem solution.

This method marks the beginning of the eigenproblem-solving process.

inline virtual void initVecs(bool random) override

Initializes vectors randomly, if necessary.

This method initializes the vectors either randomly or from previously computed Ritz values, depending on the input flag.

Parameters:

random – Flag to determine whether to initialize vectors randomly.

inline virtual void Lanczos(std::size_t m, chase::Base<T> *upperb) override

Estimates the upper bound of the user-interested spectrum using Lanczos eigensolver.

This method estimates the upper bound of the spectrum of interest using the Lanczos eigensolver.

Parameters:
  • m – Number of iterative steps for the Lanczos method.

  • upperb – Pointer to store the estimated upper bound of the spectrum.

inline virtual void Lanczos(std::size_t M, std::size_t numvec, chase::Base<T> *upperb, chase::Base<T> *ritzv, chase::Base<T> *Tau, chase::Base<T> *ritzV) override

Performs Lanczos eigensolver with additional parameters for eigenvector and Tau storage.

This method performs the Lanczos eigensolver with more advanced parameters, storing eigenvectors and Tau values for further computation.

Parameters:
  • M – Number of iterations for the Lanczos algorithm.

  • numvec – Number of vectors to be used in the Lanczos procedure.

  • upperb – Pointer to store the upper bound of the spectrum.

  • ritzv – Pointer to store the Ritz values.

  • Tau – Pointer to store the Tau values.

  • ritzV – Pointer to store the Ritz vectors.

inline virtual void LanczosDos(std::size_t idx, std::size_t m, T *ritzVc) override

Lanczos method for Density of States (DOS) computation.

This method performs a Lanczos procedure for computing the density of states.

Parameters:
  • idx – Index for the specific Lanczos vector.

  • m – Number of iterations for the Lanczos method.

  • ritzVc – Pointer to store the Ritz vectors.

inline virtual void Shift(T c, bool isunshift = false) override

Shifts the diagonal of matrix A with a specified value c.

This method applies a shift to the diagonal of the matrix A by the specified amount c. Optionally, the shift can be undone if the isunshift flag is set to true.

Parameters:
  • c – The shift value applied to the diagonal of matrix A.

  • isunshift – If true, undoes the shift applied to the diagonal.

inline virtual void HEMM(std::size_t block, T alpha, T beta, std::size_t offset) override

Performs a matrix operation of the form (V1 = \alpha \cdot H \cdot V2 + \beta \cdot V1).

This method computes the matrix operation for the input vectors (V1) and (V2) using scalars (\alpha) and (\beta), where (H) is a matrix and (V1) and (V2) are input vectors.

Parameters:
  • nev – Number of eigenpairs involved in the operation.

  • alpha – Scalar that scales the result of (H \cdot V2).

  • beta – Scalar that scales the result of (V1).

  • offset – Column offset for the operation.

inline virtual void QR(std::size_t fixednev, chase::Base<T> cond) override

Performs QR factorization with optional Householder transformation.

This method performs a QR factorization on the matrix, optionally using Householder transformations for the computation.

Parameters:
  • fixednev – Total number of converged eigenpairs before the QR factorization.

  • cond – Condition number for controlling the factorization process.

inline virtual void RR(chase::Base<T> *ritzv, std::size_t block) override

Performs Rayleigh-Ritz projection for eigenproblem reduction.

This method projects the eigenproblem into a smaller subspace, solves the smaller problem, and reconstructs the eigenvectors.

Parameters:
  • ritzv – Array to store the computed eigenvalues.

  • block – The number of non-converged eigenpairs used for the small eigenproblem.

inline virtual void Sort(chase::Base<T> *ritzv, chase::Base<T> *residLast, chase::Base<T> *resid) override

Sort the ritz values, the residuals, and the eigenvectors in ascending order.

This method sorts the ritz values, the residuals, and the ritz vectors in ascending order of the ritz values.

Parameters:
  • ritzv – Array to store the computed eigenvalues.

  • residLast – Array to store the residuals of the previous iteration.

  • resid – Array to store the residuals of the current iteration.

inline virtual void Resd(chase::Base<T> *ritzv, chase::Base<T> *resd, std::size_t fixednev) override

Computes residuals for unconverged eigenpairs.

This method computes the residuals of the eigenpairs that have not yet converged.

Parameters:
  • ritzv – Array containing the eigenvalues.

  • resd – Array to store the computed residuals for each eigenpair.

  • fixednev – Number of converged eigenpairs, which helps determine the number of non-converged eigenpairs.

inline virtual void Swap(std::size_t i, std::size_t j) override

Swaps two columns in a matrix for Chebyshev filtering.

This method swaps the specified columns of a matrix to facilitate Chebyshev filtering.

Parameters:
  • i – One of the column indices to be swapped.

  • j – Another column index to be swapped.

inline virtual void Lock(std::size_t new_converged) override

Locks the newly converged eigenvectors.

This method locks the newly converged eigenvectors and increments the locked count.

Parameters:

new_converged – The number of newly converged eigenvectors.

inline virtual void End() override

Indicates the end of an eigenproblem solution.

This method marks the conclusion of the eigenproblem-solving process.

1.2.1.2. chase::Impl::ChASEGPU

template<class T, typename MatrixType = chase::matrix::Matrix<T, chase::platform::GPU>>
class ChASEGPU : public chase::ChaseBase<T>

GPU-based sequential Chase algorithm.

This class is responsible for solving generalized eigenvalue problems using GPU-based Lanczos iterations and other matrix operations. It uses CUDA, CUBLAS, and CUSOLVER libraries for efficient matrix operations and eigenvalue solvers.

Template Parameters:

T – The data type (e.g., float, double).

Public Functions

inline ChASEGPU(std::size_t N, std::size_t nev, std::size_t nex, T *H, std::size_t ldh, T *V1, std::size_t ldv, chase::Base<T> *ritzv)

Constructor for the ChASEGPU class.

Initializes matrices and GPU resources, including CUDA streams, CUBLAS, and CUSOLVER handles. Allocates memory for the GPU matrices and computes necessary buffer sizes for matrix operations.

Parameters:
  • N – The size of the matrix.

  • nev – The number of eigenvalues to compute.

  • nex – The number of extra vectors.

  • H – Pointer to the matrix ( H ).

  • ldh – The leading dimension of matrix ( H ).

  • V1 – Pointer to the matrix ( V_1 ).

  • ldv – The leading dimension of matrix ( V_1 ).

  • ritzv – Pointer to the Ritz values vector.

inline ChASEGPU(std::size_t N, std::size_t nev, std::size_t nex, MatrixType *H, T *V1, std::size_t ldv, chase::Base<T> *ritzv)
inline void CUBLAS_INIT()
ChASEGPU(const ChASEGPU&) = delete
inline ~ChASEGPU()

Destructor for the ChASEGPU class.

Frees all dynamically allocated memory and destroys CUDA resources (CUBLAS, CUSOLVER, CUDA buffers).

inline virtual std::size_t GetN() const override

Returns the size of the matrix.

This method returns the size of the matrix (number of rows or columns).

Returns:

The size of the matrix.

inline virtual std::size_t GetNev() override

Returns the number of eigenpairs to be computed.

This method returns the number of eigenpairs that need to be computed.

Returns:

The number of eigenpairs.

inline virtual std::size_t GetNex() override

Returns the value of nex, which is used for specifying the number of eigenpairs.

This method returns the size or number of eigenpairs required by the algorithm for the computation.

Returns:

The number of eigenpairs as a std::size_t value.

inline virtual chase::Base<T> *GetRitzv() override

Returns the computed Ritz values.

This method provides access to the Ritz values computed during the algorithm’s execution. These values are typically used for eigenvalue approximations and further analysis.

Returns:

A pointer to an array of Ritz values of type Base<T>.

inline virtual chase::Base<T> *GetResid() override

Returns the residuals of computed Ritz pairs.

This method returns the residuals associated with the Ritz values. Residuals are often used to measure the accuracy or convergence of the eigenvalue problem solution.

Returns:

A pointer to an array of residuals corresponding to the Ritz values of type Base<T>.

inline virtual ChaseConfig<T> &GetConfig() override

Returns the configuration parameters used in the algorithm.

This method provides access to the configuration object that contains various algorithm settings, parameters, and options.

Returns:

A reference to a ChaseConfig<T> class containing the configuration parameters.

inline virtual int get_nprocs() override

Returns the number of MPI processes used.

This method returns the number of MPI processes involved in the computation. If sequential ChASE is used, it returns 1, indicating a single process is running.

Returns:

The number of MPI processes used, represented as an integer.

inline void loadProblemFromFile(std::string filename)
inline virtual void Output(std::string str) override

Print some intermediate infos during the solving procedure.

inline virtual bool checkSymmetryEasy() override

Checks for easy symmetry in the matrix.

This method performs a simple check for matrix symmetry.

Returns:

true if the matrix is symmetric, otherwise false.

inline virtual bool isSym() override

Returns if the matrix is Hermitian or not.

This method simply returns if the matrix is Hermitian or not

Returns:

true if the matrix is Hermitian, otherwise false.

inline virtual bool checkPseudoHermicityEasy() override

Checks for easy pseudo-hermicity in the matrix.

This method performs a simple check for matrix pseudo-hermicity

Returns:

true if the matrix is pseudo-hermicity, otherwise false.

inline virtual bool isPseudoHerm() override

Returns if the matrix is Pseudo-Hermitian or not.

This method simply returns if the matrix is Pseudo-Hermitian or not

Returns:

true if the matrix is Pseudo-Hermitian, otherwise false.

inline virtual void symOrHermMatrix(char uplo) override

Sets matrix type for symmetric or Hermitian matrix.

This method defines whether the matrix is symmetric or Hermitian based on the specified input.

Parameters:

uplo – Specifies whether the upper or lower triangular part of the matrix is stored.

inline virtual void Start() override

Indicates the start of an eigenproblem solution.

This method marks the beginning of the eigenproblem-solving process.

inline virtual void initVecs(bool random) override

Initializes vectors randomly, if necessary.

This method initializes the vectors either randomly or from previously computed Ritz values, depending on the input flag.

Parameters:

random – Flag to determine whether to initialize vectors randomly.

inline virtual void Lanczos(std::size_t M, chase::Base<T> *upperb) override

Estimates the upper bound of the user-interested spectrum using Lanczos eigensolver.

This method estimates the upper bound of the spectrum of interest using the Lanczos eigensolver.

Parameters:
  • m – Number of iterative steps for the Lanczos method.

  • upperb – Pointer to store the estimated upper bound of the spectrum.

inline virtual void Lanczos(std::size_t M, std::size_t numvec, chase::Base<T> *upperb, chase::Base<T> *ritzv, chase::Base<T> *Tau, chase::Base<T> *ritzV) override

Performs Lanczos eigensolver with additional parameters for eigenvector and Tau storage.

This method performs the Lanczos eigensolver with more advanced parameters, storing eigenvectors and Tau values for further computation.

Parameters:
  • M – Number of iterations for the Lanczos algorithm.

  • numvec – Number of vectors to be used in the Lanczos procedure.

  • upperb – Pointer to store the upper bound of the spectrum.

  • ritzv – Pointer to store the Ritz values.

  • Tau – Pointer to store the Tau values.

  • ritzV – Pointer to store the Ritz vectors.

inline virtual void LanczosDos(std::size_t idx, std::size_t m, T *ritzVc) override

Lanczos method for Density of States (DOS) computation.

This method performs a Lanczos procedure for computing the density of states.

Parameters:
  • idx – Index for the specific Lanczos vector.

  • m – Number of iterations for the Lanczos method.

  • ritzVc – Pointer to store the Ritz vectors.

inline virtual void Shift(T c, bool isunshift = false) override

Shifts the diagonal of matrix A with a specified value c.

This method applies a shift to the diagonal of the matrix A by the specified amount c. Optionally, the shift can be undone if the isunshift flag is set to true.

Parameters:
  • c – The shift value applied to the diagonal of matrix A.

  • isunshift – If true, undoes the shift applied to the diagonal.

inline virtual void HEMM(std::size_t block, T alpha, T beta, std::size_t offset) override

Performs a matrix operation of the form (V1 = \alpha \cdot H \cdot V2 + \beta \cdot V1).

This method computes the matrix operation for the input vectors (V1) and (V2) using scalars (\alpha) and (\beta), where (H) is a matrix and (V1) and (V2) are input vectors.

Parameters:
  • nev – Number of eigenpairs involved in the operation.

  • alpha – Scalar that scales the result of (H \cdot V2).

  • beta – Scalar that scales the result of (V1).

  • offset – Column offset for the operation.

inline virtual void QR(std::size_t fixednev, chase::Base<T> cond) override

Performs QR factorization with optional Householder transformation.

This method performs a QR factorization on the matrix, optionally using Householder transformations for the computation.

Parameters:
  • fixednev – Total number of converged eigenpairs before the QR factorization.

  • cond – Condition number for controlling the factorization process.

inline virtual void RR(chase::Base<T> *ritzv, std::size_t block) override

Performs Rayleigh-Ritz projection for eigenproblem reduction.

This method projects the eigenproblem into a smaller subspace, solves the smaller problem, and reconstructs the eigenvectors.

Parameters:
  • ritzv – Array to store the computed eigenvalues.

  • block – The number of non-converged eigenpairs used for the small eigenproblem.

inline virtual void Sort(chase::Base<T> *ritzv, chase::Base<T> *residLast, chase::Base<T> *resid) override

Sort the ritz values, the residuals, and the eigenvectors in ascending order.

This method sorts the ritz values, the residuals, and the ritz vectors in ascending order of the ritz values.

Parameters:
  • ritzv – Array to store the computed eigenvalues.

  • residLast – Array to store the residuals of the previous iteration.

  • resid – Array to store the residuals of the current iteration.

inline virtual void Resd(chase::Base<T> *ritzv, chase::Base<T> *resd, std::size_t fixednev) override

Computes residuals for unconverged eigenpairs.

This method computes the residuals of the eigenpairs that have not yet converged.

Parameters:
  • ritzv – Array containing the eigenvalues.

  • resd – Array to store the computed residuals for each eigenpair.

  • fixednev – Number of converged eigenpairs, which helps determine the number of non-converged eigenpairs.

inline virtual void Swap(std::size_t i, std::size_t j) override

Swaps two columns in a matrix for Chebyshev filtering.

This method swaps the specified columns of a matrix to facilitate Chebyshev filtering.

Parameters:
  • i – One of the column indices to be swapped.

  • j – Another column index to be swapped.

inline virtual void Lock(std::size_t new_converged) override

Locks the newly converged eigenvectors.

This method locks the newly converged eigenvectors and increments the locked count.

Parameters:

new_converged – The number of newly converged eigenvectors.

inline virtual void End() override

Indicates the end of an eigenproblem solution.

This method marks the conclusion of the eigenproblem-solving process.

1.2.1.3. chase::Impl::pChASECPU

template<typename MatrixType, typename InputMultiVectorType, typename BackendType = chase::grid::backend::MPI>
class pChASECPU : public chase::ChaseBase<MatrixType::value_type>

CPU-based parallel Chase algorithm with MPI.

This class solves generalized eigenvalue problems using a parallel implementation of the Chase algorithm. It leverages MPI for parallel processing and works with matrix and multi-vector data types.

Template Parameters:

Public Functions

inline pChASECPU(std::size_t nev, std::size_t nex, MatrixType *H, InputMultiVectorType *V, chase::Base<T> *ritzv)

Constructor for the pChASECPU class.

Initializes the Chase algorithm with the given matrix, eigenvector, and Ritz values. The matrix and eigenvectors must be mapped to the same MPI grid. The constructor also verifies that the matrix is square and sets up the necessary data structures for the algorithm.

Parameters:
  • nev – The number of eigenvalues to compute.

  • nex – The number of additional eigenvalues.

  • H – Pointer to the matrix to solve.

  • V – Pointer to the input multi-vector of eigenvectors.

  • ritzv – Pointer to the Ritz values used in the algorithm.

pChASECPU(const pChASECPU&) = delete

Deleted copy constructor.

The copy constructor is deleted to prevent copying of this class, as it manages unique resources such as MPI communication.

inline ~pChASECPU()

Destructor for the pChASECPU class.

Cleans up any resources allocated by the constructor, including matrix data and MPI-specific information.

inline virtual std::size_t GetN() const override

Returns the size of the matrix.

This method returns the size of the matrix (number of rows or columns).

Returns:

The size of the matrix.

inline virtual std::size_t GetNev() override

Returns the number of eigenpairs to be computed.

This method returns the number of eigenpairs that need to be computed.

Returns:

The number of eigenpairs.

inline virtual std::size_t GetNex() override

Returns the value of nex, which is used for specifying the number of eigenpairs.

This method returns the size or number of eigenpairs required by the algorithm for the computation.

Returns:

The number of eigenpairs as a std::size_t value.

inline virtual chase::Base<T> *GetRitzv() override

Returns the computed Ritz values.

This method provides access to the Ritz values computed during the algorithm’s execution. These values are typically used for eigenvalue approximations and further analysis.

Returns:

A pointer to an array of Ritz values of type Base<T>.

inline virtual chase::Base<T> *GetResid() override

Returns the residuals of computed Ritz pairs.

This method returns the residuals associated with the Ritz values. Residuals are often used to measure the accuracy or convergence of the eigenvalue problem solution.

Returns:

A pointer to an array of residuals corresponding to the Ritz values of type Base<T>.

inline virtual ChaseConfig<T> &GetConfig() override

Returns the configuration parameters used in the algorithm.

This method provides access to the configuration object that contains various algorithm settings, parameters, and options.

Returns:

A reference to a ChaseConfig<T> class containing the configuration parameters.

inline virtual int get_nprocs() override

Returns the number of MPI processes used.

This method returns the number of MPI processes involved in the computation. If sequential ChASE is used, it returns 1, indicating a single process is running.

Returns:

The number of MPI processes used, represented as an integer.

inline int get_rank()
inline void loadProblemFromFile(std::string filename)
inline void saveProblemToFile(std::string filename)
inline virtual void Output(std::string str) override

Print some intermediate infos during the solving procedure.

inline virtual bool checkSymmetryEasy() override

Checks for easy symmetry in the matrix.

This method performs a simple check for matrix symmetry.

Returns:

true if the matrix is symmetric, otherwise false.

inline virtual bool isSym() override

Returns if the matrix is Hermitian or not.

This method simply returns if the matrix is Hermitian or not

Returns:

true if the matrix is Hermitian, otherwise false.

inline virtual bool checkPseudoHermicityEasy() override

Checks for easy pseudo-hermicity in the matrix.

This method performs a simple check for matrix pseudo-hermicity

Returns:

true if the matrix is pseudo-hermicity, otherwise false.

inline virtual bool isPseudoHerm() override

Returns if the matrix is Pseudo-Hermitian or not.

This method simply returns if the matrix is Pseudo-Hermitian or not

Returns:

true if the matrix is Pseudo-Hermitian, otherwise false.

inline virtual void symOrHermMatrix(char uplo) override

Sets matrix type for symmetric or Hermitian matrix.

This method defines whether the matrix is symmetric or Hermitian based on the specified input.

Parameters:

uplo – Specifies whether the upper or lower triangular part of the matrix is stored.

inline virtual void Start() override

Indicates the start of an eigenproblem solution.

This method marks the beginning of the eigenproblem-solving process.

inline virtual void initVecs(bool random) override

Initializes vectors randomly, if necessary.

This method initializes the vectors either randomly or from previously computed Ritz values, depending on the input flag.

Parameters:

random – Flag to determine whether to initialize vectors randomly.

inline void Lanczos(std::size_t m, chase::Base<T> *upperb) override
inline void Lanczos(std::size_t M, std::size_t numvec, chase::Base<T> *upperb, chase::Base<T> *ritzv, chase::Base<T> *Tau, chase::Base<T> *ritzV) override
inline virtual void LanczosDos(std::size_t idx, std::size_t m, T *ritzVc) override

Lanczos method for Density of States (DOS) computation.

This method performs a Lanczos procedure for computing the density of states.

Parameters:
  • idx – Index for the specific Lanczos vector.

  • m – Number of iterations for the Lanczos method.

  • ritzVc – Pointer to store the Ritz vectors.

inline virtual void Shift(T c, bool isunshift = false) override

Shifts the diagonal of matrix A with a specified value c.

This method applies a shift to the diagonal of the matrix A by the specified amount c. Optionally, the shift can be undone if the isunshift flag is set to true.

Parameters:
  • c – The shift value applied to the diagonal of matrix A.

  • isunshift – If true, undoes the shift applied to the diagonal.

inline virtual void HEMM(std::size_t block, T alpha, T beta, std::size_t offset) override

Performs a matrix operation of the form (V1 = \alpha \cdot H \cdot V2 + \beta \cdot V1).

This method computes the matrix operation for the input vectors (V1) and (V2) using scalars (\alpha) and (\beta), where (H) is a matrix and (V1) and (V2) are input vectors.

Parameters:
  • nev – Number of eigenpairs involved in the operation.

  • alpha – Scalar that scales the result of (H \cdot V2).

  • beta – Scalar that scales the result of (V1).

  • offset – Column offset for the operation.

inline void QR(std::size_t fixednev, chase::Base<T> cond) override
inline void RR(chase::Base<T> *ritzv, std::size_t block) override
inline void Sort(chase::Base<T> *ritzv, chase::Base<T> *residLast, chase::Base<T> *resid) override
inline void Resd(chase::Base<T> *ritzv, chase::Base<T> *resd, std::size_t fixednev) override
inline virtual void Swap(std::size_t i, std::size_t j) override

Swaps two columns in a matrix for Chebyshev filtering.

This method swaps the specified columns of a matrix to facilitate Chebyshev filtering.

Parameters:
  • i – One of the column indices to be swapped.

  • j – Another column index to be swapped.

inline virtual void Lock(std::size_t new_converged) override

Locks the newly converged eigenvectors.

This method locks the newly converged eigenvectors and increments the locked count.

Parameters:

new_converged – The number of newly converged eigenvectors.

inline virtual void End() override

Indicates the end of an eigenproblem solution.

This method marks the conclusion of the eigenproblem-solving process.

1.2.1.4. chase::Impl::pChASEGPU

template<typename MatrixType, typename InputMultiVectorType, typename BackendType = chase::grid::backend::NCCL>
class pChASEGPU : public chase::ChaseBase<MatrixType::value_type>

GPU-based parallel Chase algorithm with NCCL.

This class provides an implementation of the Chase algorithm for solving generalized eigenvalue problems on GPU-based platforms. It leverages NVIDIA’s NCCL for efficient GPU-to-GPU communication and cuSolver and cuBLAS for GPU-based linear algebra operations. Designed for distributed GPU environments, it operates on matrix and multi-vector data types with MPI parallelism.

Template Parameters:

Public Functions

inline pChASEGPU(std::size_t nev, std::size_t nex, MatrixType *H, InputMultiVectorType *V, chase::Base<T> *ritzv)

Constructs the pChASEGPU object.

Initializes the Chase algorithm parameters, matrix data, multi-vector data, and GPU-specific resources. Sets up memory allocations for intermediate data and configures cuBLAS and cuSolver for GPU computation.

Parameters:
  • nev – Number of eigenvalues to compute.

  • nex – Number of extra vectors for iterative refinement.

  • H – Pointer to the Hamiltonian matrix.

  • V – Pointer to the input multi-vector (typically for eigenvectors).

  • ritzv – Pointer to the Ritz values, used for storing eigenvalues.

pChASEGPU(const pChASEGPU&) = delete
inline ~pChASEGPU()
inline virtual std::size_t GetN() const override

Returns the size of the matrix.

This method returns the size of the matrix (number of rows or columns).

Returns:

The size of the matrix.

inline virtual std::size_t GetNev() override

Returns the number of eigenpairs to be computed.

This method returns the number of eigenpairs that need to be computed.

Returns:

The number of eigenpairs.

inline virtual std::size_t GetNex() override

Returns the value of nex, which is used for specifying the number of eigenpairs.

This method returns the size or number of eigenpairs required by the algorithm for the computation.

Returns:

The number of eigenpairs as a std::size_t value.

inline virtual chase::Base<T> *GetRitzv() override

Returns the computed Ritz values.

This method provides access to the Ritz values computed during the algorithm’s execution. These values are typically used for eigenvalue approximations and further analysis.

Returns:

A pointer to an array of Ritz values of type Base<T>.

inline virtual chase::Base<T> *GetResid() override

Returns the residuals of computed Ritz pairs.

This method returns the residuals associated with the Ritz values. Residuals are often used to measure the accuracy or convergence of the eigenvalue problem solution.

Returns:

A pointer to an array of residuals corresponding to the Ritz values of type Base<T>.

inline virtual ChaseConfig<T> &GetConfig() override

Returns the configuration parameters used in the algorithm.

This method provides access to the configuration object that contains various algorithm settings, parameters, and options.

Returns:

A reference to a ChaseConfig<T> class containing the configuration parameters.

inline virtual int get_nprocs() override

Returns the number of MPI processes used.

This method returns the number of MPI processes involved in the computation. If sequential ChASE is used, it returns 1, indicating a single process is running.

Returns:

The number of MPI processes used, represented as an integer.

inline int get_rank()
inline void loadProblemFromFile(std::string filename)
inline void saveProblemToFile(std::string filename)
inline virtual void Output(std::string str) override

Print some intermediate infos during the solving procedure.

inline virtual bool checkSymmetryEasy() override

Checks for easy symmetry in the matrix.

This method performs a simple check for matrix symmetry.

Returns:

true if the matrix is symmetric, otherwise false.

inline virtual bool isSym() override

Returns if the matrix is Hermitian or not.

This method simply returns if the matrix is Hermitian or not

Returns:

true if the matrix is Hermitian, otherwise false.

inline virtual bool checkPseudoHermicityEasy() override

Checks for easy pseudo-hermicity in the matrix.

This method performs a simple check for matrix pseudo-hermicity

Returns:

true if the matrix is pseudo-hermicity, otherwise false.

inline virtual bool isPseudoHerm() override

Returns if the matrix is Pseudo-Hermitian or not.

This method simply returns if the matrix is Pseudo-Hermitian or not

Returns:

true if the matrix is Pseudo-Hermitian, otherwise false.

inline virtual void symOrHermMatrix(char uplo) override

Sets matrix type for symmetric or Hermitian matrix.

This method defines whether the matrix is symmetric or Hermitian based on the specified input.

Parameters:

uplo – Specifies whether the upper or lower triangular part of the matrix is stored.

inline virtual void Start() override

Indicates the start of an eigenproblem solution.

This method marks the beginning of the eigenproblem-solving process.

inline virtual void initVecs(bool random) override

Initializes vectors randomly, if necessary.

This method initializes the vectors either randomly or from previously computed Ritz values, depending on the input flag.

Parameters:

random – Flag to determine whether to initialize vectors randomly.

inline void Lanczos(std::size_t m, chase::Base<T> *upperb) override
inline void Lanczos(std::size_t M, std::size_t numvec, chase::Base<T> *upperb, chase::Base<T> *ritzv, chase::Base<T> *Tau, chase::Base<T> *ritzV) override
inline virtual void LanczosDos(std::size_t idx, std::size_t m, T *ritzVc) override

Lanczos method for Density of States (DOS) computation.

This method performs a Lanczos procedure for computing the density of states.

Parameters:
  • idx – Index for the specific Lanczos vector.

  • m – Number of iterations for the Lanczos method.

  • ritzVc – Pointer to store the Ritz vectors.

inline virtual void Shift(T c, bool isunshift = false) override

Shifts the diagonal of matrix A with a specified value c.

This method applies a shift to the diagonal of the matrix A by the specified amount c. Optionally, the shift can be undone if the isunshift flag is set to true.

Parameters:
  • c – The shift value applied to the diagonal of matrix A.

  • isunshift – If true, undoes the shift applied to the diagonal.

inline virtual void HEMM(std::size_t block, T alpha, T beta, std::size_t offset) override

Performs a matrix operation of the form (V1 = \alpha \cdot H \cdot V2 + \beta \cdot V1).

This method computes the matrix operation for the input vectors (V1) and (V2) using scalars (\alpha) and (\beta), where (H) is a matrix and (V1) and (V2) are input vectors.

Parameters:
  • nev – Number of eigenpairs involved in the operation.

  • alpha – Scalar that scales the result of (H \cdot V2).

  • beta – Scalar that scales the result of (V1).

  • offset – Column offset for the operation.

inline void QR(std::size_t fixednev, chase::Base<T> cond) override
inline void RR(chase::Base<T> *ritzv, std::size_t block) override
inline void Sort(chase::Base<T> *ritzv, chase::Base<T> *residLast, chase::Base<T> *resid) override
inline void Resd(chase::Base<T> *ritzv, chase::Base<T> *resd, std::size_t fixednev) override
inline virtual void Swap(std::size_t i, std::size_t j) override

Swaps two columns in a matrix for Chebyshev filtering.

This method swaps the specified columns of a matrix to facilitate Chebyshev filtering.

Parameters:
  • i – One of the column indices to be swapped.

  • j – Another column index to be swapped.

inline virtual void Lock(std::size_t new_converged) override

Locks the newly converged eigenvectors.

This method locks the newly converged eigenvectors and increments the locked count.

Parameters:

new_converged – The number of newly converged eigenvectors.

inline virtual void End() override

Indicates the end of an eigenproblem solution.

This method marks the conclusion of the eigenproblem-solving process.