ChaseMpi

template<template<typename> class MF, class T>
class chase::mpi::ChaseMpi : public chase::Chase<T>

A derived class of Chase to implement ChASE based on MPI and Dense Linear Algebra (DLA) routines.

This is a calls derived from the Chase class which plays the role of interface for the kernels used by the library.

  • All members of the Chase class are virtual functions. These functions are re-implemented in the ChaseMpi class.

  • All the members functions of ChaseMpi, which are the implementation of the virtual functions in class Chase, are implemented using the DLA routines provided by the class ChaseMpiDLAInterface.

  • The DLA functions in ChaseMpiDLAInterface are also vritual functions, which are differently implemented targeting different computing architectures (sequential/parallel, CPU/GPU, shared-memory/distributed-memory, etc). In the class ChaseMpi, the calling of DLA functions are indeed calling their implementations from different derived classes. Thus this ChaseMpi class is able to have customized implementation for various architectures.

  • For the implementation of the class ChaseMpi targeting distributed-memory platforms based on MPI, the setup of MPI environment and communication scheme, and the distribution of data (matrix, vectors) across MPI nodes are following the ChaseMpiProperties class, the distribution of matrix can be either Block or Block-Cyclic scheme.

    Template Parameters
    • MF: A class derived from ChaseMpiDLAInterface, which indicates the selected implementation of DLA that to be used by ChaseMpi Object.

    • T: the scalar type used for the application. ChASE is templated for real and complex numbers with both Single Precision and Double Precision, thus T can be one of float, double, std::complex<float> and std::complex<double>.

Public Functions

ChaseMpi(std::size_t N, std::size_t nev, std::size_t nex, T *V1 = nullptr, Base<T> *ritzv = nullptr, T *H = nullptr, T *V2 = nullptr, Base<T> *resid = nullptr)

A constructor of the ChaseMpi class which gives an implenentation of ChASE for shared-memory architecture, without MPI.

The private members of this classes are initialized by the parameters of this constructor.

  • For the variable N_, it is initialized by the first parameter of this constructor N.

  • The variables nev_ and and nex_ are initialized by the parameters of this constructor nev and nex, respectively.

  • The variable rank_ is set to be 0 since non MPI is supported. The variable locked_ is initially set to be 0.

  • The variable config_ is setup by the constructor of ChaseConfig which takes the parameters N, nev and nex.

  • The variable matrices_ is setup directly by the constructor of ChaseMpiMatrices which takes the paramter H, N, nev, nex, V1, ritzv, V2 and resid.

  • The variable dla_ is initialized by MF which is a derived class of ChaseMpiDLAInterface. In ChASE, the candidates of MF for non-MPI case are the classes ChaseMpiDLABlaslapackSeq, ChaseMpiDLABlaslapackSeqInplace and ChaseMpiDLACudaSeq.

Parameters
  • N: size of the square matrix defining the eigenproblem.

  • nev: Number of desired extremal eigenvalues.

  • nex: Number of eigenvalues augmenting the search space. Usually a relatively small fraction of nev.

  • V1: a pointer to a rectangular matrix of size N * (nev+nex). If V1 is not provided by the user, it will be internally allocated in ChaseMpiMatrices class. After the solving step, the first nev columns of V1 are overwritten by the desired Ritz vectors.

  • ritzv: a pointer to an array to store the computed Ritz values. If ritzv is not provided by the user, it will be internally allocated in ChaseMpiMatrices class. Its minimal size should be nev+nex.

  • H: a pointer to the memory which stores local part of matrix on each MPI rank. If H is not provided by the user, it will be internally allocated in ChaseMpiMatrices class.

  • V2: a pointer to anther rectangular matrix of size N * (nev+nex). f V2 is not provided by the user, it will be internally allocated in ChaseMpiMatrices class.

  • resid: a pointer to an array to store the residual of computed eigenpairs. If resid is not provided by the user, it will be internally allocated in ChaseMpiMatrices class. Its minimal size should be nev+nex.

ChaseMpi(ChaseMpiProperties<T> *properties, T *V1 = nullptr, Base<T> *ritzv = nullptr, T *V2 = nullptr, Base<T> *resid = nullptr)

A constructor of the ChaseMpi class which gives an implenentation of ChASE for distributed-memory architecture, with the support of MPI. In this case, the buffer for matrix H is allocated internally.

The private members of this classes are initialized by the parameters of this constructor.

  • For the variable N_, nev_ and nex_ are initialized by the first parameter of this constructor properties_.

  • The variable rank_ is initialized by MPI_Comm_rank. The variable locked_ is initially set to be 0.

  • The variable config_ is setup by the constructor of ChaseConfig which takes the parameters N, nev and nex.

  • The variable matrices_ is constructed by the create_matrices function defined in ChaseMpiProperties.

  • The variable dla_ for the distributed-memory is initialized by the constructor of ChaseMpiDLA which takes both properties_ and MF. In MPI case, the implementation of are split into two classses:

Parameters
  • properties: an object of ChaseMpiProperties which setups the MPI environment and data distribution scheme for ChaseMpi targeting distributed-memory systems.

  • V1: a pointer to a rectangular matrix of size m * (nev+nex). m can be obtained through ChaseMpiProperties::get_m(). V1 is partially distributed within each column communicator and is redundant among different column communicator. If V1 is not provided by the user, it will be internally allocated in ChaseMpiMatrices class. After the solving step, the first nev columns of V1 are overwritten by the desired Ritz vectors.

  • ritzv: a pointer to an array to store the computed Ritz values. If ritzv is not provided by the user, it will be internally allocated in ChaseMpiMatrices class. Its minimal size should be nev+nex.

  • V2: a pointer to anther rectangular matrix of size n * (nev+nex). If V2 is not provided by the user, it will be internally allocated in ChaseMpiMatrices class. n can be obtained through ChaseMpiProperties::get_n(). V2 is partially distributed within each row communicator and is redundant among different row communicator.

  • resid: a pointer to an array to store the residual of computed eigenpairs. If resid is not provided by the user, it will be internally allocated in ChaseMpiMatrices class. Its minimal size should be nev+nex.

ChaseMpi(ChaseMpiProperties<T> *properties, T *H, std::size_t ldh, T *V1, Base<T> *ritzv, T *V2 = nullptr, Base<T> *resid = nullptr)

A constructor of the ChaseMpi class which gives an implenentation of ChASE for distributed-memory architecture, with the support of MPI. In this case, the buffer for matrix H should be allocated externally and provided by users.

The private members of this classes are initialized by the parameters of this constructor.

  • For the variable N_, nev_ and nex_ are initialized by the first parameter of this constructor properties_.

  • The variable rank_ is initialized by MPI_Comm_rank. The variable locked_ is initially set to be 0.

  • The variable config_ is setup by the constructor of ChaseConfig which takes the parameters N, nev and nex.

  • The variable matrices_ is constructed by the create_matrices function defined in ChaseMpiProperties.

  • The variable dla_ for the distributed-memory is initialized by the constructor of ChaseMpiDLA which takes both properties_ and MF. In MPI case, the implementation of are split into two classses:

Parameters
  • properties: an object of ChaseMpiProperties which setups the MPI environment and data distribution scheme for ChaseMpi targeting distributed-memory systems.

  • H: a pointer to a local matrix of the distributed Symmetric/Hermtitian matrix to be diagonalised. It should be provided and allocated by users. The minimal dimension is mxn, in which m and n can be obtained through ChaseMpiProperties::get_m() and ChaseMpiProperties::get_n().

  • ldh: leading dimension of local matrix H, it should be >= ChaseMpiProperties::get_m()

  • V1: a pointer to a rectangular matrix of size m * (nev+nex). m can be obtained through ChaseMpiProperties::get_m(). V1 is partially distributed within each column communicator and is redundant among different column communicator. If V1 is not provided by the user, it will be internally allocated in ChaseMpiMatrices class. After the solving step, the first nev columns of V1 are overwritten by the desired Ritz vectors.

  • ritzv: a pointer to an array to store the computed Ritz values. If ritzv is not provided by the user, it will be internally allocated in ChaseMpiMatrices class. Its minimal size should be nev+nex.

  • V2: a pointer to anther rectangular matrix of size n * (nev+nex). V2 is partially distributed within each row communicator and is redundant among different row communicator. If V2 is not provided by the user, it will be internally allocated in ChaseMpiMatrices class. m can be obtained through ChaseMpiProperties::get_n().

  • resid: a pointer to an array to store the residual of computed eigenpairs. If resid is not provided by the user, it will be internally allocated in ChaseMpiMatrices class. Its minimal size should be nev+nex.

ChaseMpi(const ChaseMpi&) = delete

It prevents the copy operation of the constructor of ChaseMpi.

ChaseConfig<T> &GetConfig() override

This member function implements the virtual one declared in Chase class.

Return

config_: an object of ChaseConfig class which defines the eigensolver.

std::size_t GetN() const override

This member function implements the virtual one declared in Chase class. Returns the rank of matrix A which is distributed within 2D MPI grid.

Return

N_: the rank of matrix A.

std::size_t GetNev() override

This member function implements the virtual one declared in Chase class. Returns Number of desired extremal eigenpairs, which was set by users.

Return

nev_: Number of desired extremal eigenpairs.

std::size_t GetNex() override

This member function implements the virtual one declared in Chase class. Returns the Increment of the search subspace so that its total size, which was set by users.

Return

nex_: Increment of the search subspace.

Base<T> *GetRitzv() override

This member function implements the virtual one declared in Chase class. Returns the array which stores the computed Ritz values.

Return

ritzv_: an array stores the computed Ritz values.

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

This member function implements the virtual one declared in Chase class. This member function shifts the diagonal of matrix A with a shift value c. It is implemented by calling the different implementation in the derived class of ChaseMpiDLAInterface. For the construction of ChaseMpi with MPI, this function is naturally in parallel across the MPI nodes, since the shift operation only takes places on selected local matrices which stores in a distributed manner.

Parameters
  • c: the shift value on the diagonal of matrix A.

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

This member function implements the virtual one declared in Chase class. This member function computes \(V1 = alpha * H*V2 + beta * V1\) or \(V2 = alpha * H'*V1 + beta * V2\). This operation starts with an offset offset of column,

Parameters
  • block: the number of vectors in V1 and V2 for this HEMM operation.

  • alpha: a scalar of type T which scales H*V1 or H*V2`.

  • beta: a scalar of type T which scales V1 or V2, respectively.

  • offset: the offset of column which the HEMM starts from.

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

This member function performs a QR factorization with an explicit construction of the unitary matrix Q. After the explicit construction of Q, its first locked_ number of vectors are overwritten by the converged eigenvectors stored previously.

CholQR is used in default, and it can be disabled by adding a environment variable CHASE_DISABLE_CHOLQR=1. When CholQR is disabled, ScaLAPACK Householder QR is used whenever possible,

Parameters
  • fixednev: total number of converged eigenpairs before this time QR factorization.

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

This member function implements the virtual one declared in Chase class.

This function performs the Rayleigh-Ritz projection, which projects the eigenproblem to be a small one, then solves the small problem and reconstructs the eigenvectors.

Parameters
  • ritzv: a pointer to the array which stores the computed eigenvalues.

  • block: the number of non-converged eigenpairs, which determines the size of small eigenproblem.

void Resd(Base<T> *ritzv, Base<T> *resid, std::size_t fixednev) override

This member function implements the virtual one declared in Chase class. This member function computes the residuals of unconverged eigenpairs.

Parameters
  • ritzv: an array stores the eigenvalues.

  • resid: a pointer to an array which stores the residual for each eigenpairs.

  • fixednev: number of converged eigenpairs. Thus the number of non-converged one which perform this operation of computing residual is nev_ + nex_ + fixednev_.

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

This member function implements the virtual one declared in Chase class. It swaps the two matrices of vectors used in the Chebyschev filter

Parameters
  • i: one of the column index to be swapped

  • j: another of the column index to be swapped rectangular matrices approxV_ and workspace_.

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

This member function implements the virtual one declared in Chase class. It estimates the upper bound of user-interested spectrum by Lanczos eigensolver

Parameters
  • m: the iterative steps for Lanczos eigensolver.

  • upperb: a pointer to the upper bound estimated by Lanczos eigensolver.

void Lanczos(std::size_t M, std::size_t idx, Base<T> *upperb, Base<T> *ritzv, Base<T> *Tau, Base<T> *ritzV) override

This member function implements the virtual one declared in Chase class. It estimates the upper bound of user-interested spectrum by Lanczos eigensolver

void Lock(std::size_t new_converged) override

This member function implements the virtual one declared in Chase class. It locks the new_converged eigenvectors, which makes locked_ += new_converged.

Parameters
  • new_converged: number of newly converged eigenpairs in the present iterative step.

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

This member function implements the virtual one declared in Chase class. It estimates the spectral distribution of eigenvalues.

T *GetMatrixPtr()

Return

H_: A pointer to the memory allocated to store (local part if applicable) of matrix A.

Base<T> *GetResid() override

This member function implements the virtual one declared in Chase class.

Return

resid_: a pointer to the memory allocated to store the residual of each computed eigenpair.

Private Members

const std::size_t N_

Global size of the matrix A defining the eigenproblem.

  • For the constructor of class ChaseMpi without MPI, this variable is initialized by the constructor using the value of the first of its input parameters N.

  • For the constructor of class ChaseMpi with MPI, this variable is initialized by properties_, a pointer to an object of ChaseMpiProperties. This variable is private, it can be access by the member function GetN().

const std::size_t nev_

Number of desired extremal eigenpairs.

  • For the constructor of class ChaseMpi without MPI, this variable is initialized by the constructor using the value of the second of its input parameters nev.

  • For the constructor of class ChaseMpi with MPI, this variable is initialized by properties, a pointer to an object of ChaseMpiProperties. This variable is private, it can be access by the member function GetNev().

const std::size_t nex_

Increment of the search subspace so that its total size is nev+nex.

  • For the constructor of class ChaseMpi without MPI, this variable is initialized by the constructor using the value of the third of its input parameters nex.

  • For the constructor of class ChaseMpi with MPI, this variable is initialized by properties, a pointer to an object of ChaseMpiProperties. This variable is private, it can be access by the member function GetNex().

int rank_

The rank of each MPI node with the working MPI communicator.

  • For the constructor of class ChaseMpi without MPI, this variable is 0.

  • For the constructor of class ChaseMpi with MPI, this variable is gotten through MPI function MPI_Comm_rank.

std::size_t locked_

The number of eigenvectors to be locked, which indicates the number of eigenpairs converged into acceptable tolerance.

  • For the constructor of class ChaseMpi without MPI, this variable is initialized 0.

  • For the constructor of class ChaseMpi with MPI, this variable is initialized also 0.

  • During the process of ChASE solving eigenproblem, the value of this variable will increase with more and more eigenpairs computed.

Base<T> *ritzv_

A pointer to the memory allocated to store the computed eigenvalues. The values inside are always real since the matrix of eigenproblem is Hermitian (symmetric).

  • For the constructor of class ChaseMpi without MPI, the memory is allocated directly by ChaseMpiMatrices of size nex_ + nex_.

  • For the constructor of class ChaseMpi with MPI, this variable is a pointer to a vector of size nex_ + nex_, which are identical on each MPI node. It is initalized within the construction of ChaseMpiMatrices.

  • The final converged eigenvalues will be stored in it.

  • This variable is private, and it can be assessed through the member function GetRitzv().

Base<T> *resid_

A pointer to the memory allocated to store the residual of each computed eigenpair. The values inside are always real.

  • For the constructor of class ChaseMpi without MPI, the memory is allocated directly by ChaseMpiMatrices of size nex_ + nex_.

  • For the constructor of class ChaseMpi with MPI, this variable is a pointer to a vector of size nex_ + nex_, which are identical on each MPI node. It is initalized within the construction of ChaseMpiMatrices.

  • The residuals of final converged eigenvalues will be stored in it.

  • This variable is private, and it can be assessed through the member function GetResid().

std::unique_ptr<ChaseMpiProperties<T>> properties_

A smart pointer to an object of ChaseMpiProperties class which is used by ChaseMpi class. This variable only matters for the constructor of class ChaseMpi with MPI.

ChaseMpiMatrices<T> matrices_

An object of ChaseMpiMatrices to setup the matrices and vectors used by ChaseMpi class.

std::unique_ptr<ChaseMpiDLAInterface<T>> dla_

A smart pointer to an object of ChaseMpiDLAInterface class which is used by ChaseMpi class.

  • For the constructor of class ChaseMpi without MPI, this variable is initialized directly by the template parameter MF.

  • For the constructor of class ChaseMpi with MPI, this variable is implemented by considering properties_ and MF. It is initalized within the construction of ChaseMpiMatrices.

ChaseConfig<T> config_

An object of ChaseConfig class which setup all the parameters of ChASE, these parameters are either provided by users, or using the default values. This variable is initialized by the constructor of ChaseConfig class by using N, nev and nex, which are the parameters of the constructors of ChaseMpi.