1.1.4. Performance Classes
1.1.4.1. chase::ChasePerfData
-
template<class T>
class ChasePerfData ChASE class for collecting data relative to FLOPs, timings, etc.
The ChasePerfData class collects and handles information relative to the execution of the eigensolver. It collects information about
Number of subspace iterations
Number of filtered vectors
Timings of each main algorithmic procedure (Lanczos, Filter, etc.)
Number of FLOPs executed
The number of iterations and filtered vectors can be used to monitor the behavior of the algorithm as it attempts to converge all the desired eigenpairs. The timings and number of FLOPs are use to measure performance, especially parallel performance. The timings are stored in a vector of objects derived by the class template
std::chrono::duration.Public Types
Public Functions
-
inline ChasePerfData()
-
inline void Reset()
-
inline std::size_t get_iter_count()
Returns the number of total subspace iterations executed by ChASE.
The S in ChASE stands for Subspace iteration. The main engine under the hood of ChASE is a loop enveloping all the main routines executed by the code. Because of this structure, ChASE is a truly iterative algorithm based on subspace filtering. Counting the number of times such a loop is repeated gives a measure of the effectiveness of the algorithm and it is usually a non-linear function of the spectral distribution. For example, when using the flag
approximate_ = 'true'to solve a sequence of eigenproblems, one can observe that the number of subspace iteration decreases as a function of sequences index.- Returns:
The total number of subspace iterations.
-
inline std::size_t get_filtered_vecs()
Returns the cumulative number of times each column vector is filtered by one degree.
The most computationally expensive routine of ChASE is the Chebyshev filter. Within the filter a matrix of vectors V is filtered with a varying degree each time a subspace iteration is executed. This counter return the total number of times each vector in V goes through a filtering step. For instance, when the flag
optim_ = false, such a number roughly corresponds to rank(V) x degree x iter_count. When theoptim_is set totruesuch a calculation is quite more complicated. Roughly speaking, this counter is useful to monitor the convergence ration of the filtered vectors and together withget_iter_countconvey the effectiveness of the algorithm.- Returns:
Cumulative number of filtered vectors.
-
inline std::vector<std::chrono::duration<double>> get_timings()
-
inline std::size_t get_flops(std::size_t N)
Returns the total number of FLOPs executed by ChASE.
When measuring performance, it is fundamental to understand how many operations a routine executes against the total time to solutions. This counter returns the total amount of operations executed by ChASE and can be used to extract the performance of ChASE and compare it with theoretical peak performance of the platform where the code is executed.
- Parameters:
N – Size of the eigenproblem matrix
- Returns:
The total number of operations executed by ChASE.
-
inline std::size_t get_filter_flops(std::size_t N)
Returns the total number of FLOPs of the Chebyshev filter.
Similar to
get_flops, this counter return the total number of operations executed by the Chebyshev filter alone. Since the filter is the routine that executes, on average, 80% of the total FLOPs of ChASE, this counter is a good indicator of the performance of the entire algorithm. Because the filter executes almost exclusively BLAS-3 operations, this counter is quite useful to monitor how well the filter is close to the peak performance of the platform where ChASE is executed. This can be quite useful to fine tune the use of the computational resources used.- Parameters:
N – Size of the eigenproblem matrix
- Returns:
The total number of operations executed by the polynomial filter.
-
inline void set_nprocs(int nProcs)
-
inline void add_iter_count(std::size_t add)
-
inline void add_iter_blocksize(std::size_t nevex)
-
inline void add_filtered_vecs(std::size_t add)
-
inline auto getTimePoint() -> TimePointType
-
inline void print(std::size_t N = 0)
Print function outputting counters and timings for all routines.
It prints by default ( for N = 0) in the order,
size of the eigenproblem
total number of subspace iterations executed
total number of filtered vectors
time-to-solution of the following 6 main sections of the ChASE algorithm:
Total time-to-solution
Estimates of the spectral bounds based on Lanczos,
Chebyshev filter,
QR decomposition,
Raleygh-Ritz procedure including the solution of the reduced dense problem,
Computation of the eigenpairs residuals
When the parameter
Nis set to be a number else than zero, the function returns total FLOPs and filter FLOPs, respectively.- Parameters:
N – Control parameter. By default equal to 0.
1.1.4.2. chase::PerformanceDecoratorChase
-
template<class T>
class PerformanceDecoratorChase : public chase::ChaseBase<T> A derived class used to extract performance and configuration data.
This is a class 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 PerformanceDecoratorChase class. All derived members that provide an interface to computational kernels are reimplemented by decorating the original function with time pointers which are members of the ChasePerfData class. All derived members that provide an interface to input or output data are called without any specific decoration. In addition to the virtual member of the Chase class, the PerformanceDecoratorChase class has also among its public members a reference to an object of type ChasePerfData. When using Chase to solve an eigenvalue problem, the members of the PerformanceDecoratorChase are called instead of the virtual functions members of the Chase class. In this way, all parameters and counters are automatically invoked and returned in the correct order.
See also
Chase
Public Functions
-
inline virtual void initVecs(bool random)
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 Shift(T c, bool isunshift = false)
Shifts the diagonal of matrix
Awith a specified valuec.This method applies a shift to the diagonal of the matrix
Aby the specified amountc. Optionally, the shift can be undone if theisunshiftflag 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 nev, T alpha, T beta, std::size_t offset)
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, Base<T> cond)
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(Base<T> *ritzv, std::size_t block)
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(Base<T> *ritzv, Base<T> *residLast, Base<T> *resid)
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(Base<T> *ritzv, Base<T> *resd, std::size_t fixednev)
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 Lanczos(std::size_t m, Base<T> *upperb)
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 idx, Base<T> *upperb, Base<T> *ritzv, Base<T> *Tau, Base<T> *ritzV)
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)
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 Swap(std::size_t i, std::size_t j)
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)
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 bool checkSymmetryEasy()
Checks for easy symmetry in the matrix.
This method performs a simple check for matrix symmetry.
- Returns:
trueif the matrix is symmetric, otherwisefalse.
-
inline virtual bool isSym()
Returns if the matrix is Hermitian or not.
This method simply returns if the matrix is Hermitian or not
- Returns:
trueif the matrix is Hermitian, otherwisefalse.
-
inline virtual bool checkPseudoHermicityEasy()
Checks for easy pseudo-hermicity in the matrix.
This method performs a simple check for matrix pseudo-hermicity
- Returns:
trueif the matrix is pseudo-hermicity, otherwisefalse.
-
inline virtual bool isPseudoHerm()
Returns if the matrix is Pseudo-Hermitian or not.
This method simply returns if the matrix is Pseudo-Hermitian or not
- Returns:
trueif the matrix is Pseudo-Hermitian, otherwisefalse.
-
inline virtual void symOrHermMatrix(char uplo)
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()
Indicates the start of an eigenproblem solution.
This method marks the beginning of the eigenproblem-solving process.
-
inline virtual int get_nprocs()
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 virtual void End()
Indicates the end of an eigenproblem solution.
This method marks the conclusion of the eigenproblem-solving process.
-
inline virtual std::size_t GetN() const
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()
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()
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_tvalue.
-
inline virtual Base<T> *GetRitzv()
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 Base<T> *GetResid()
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()
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 ChasePerfData<T> &GetPerfData()
-
inline virtual void Output(std::string str)
Prints intermediate information during the solving procedure.
This method is used to output intermediate information or debug data during the execution of the algorithm. It is typically enabled when debugging or logging is required.
- Parameters:
str – A string containing the information to be printed.
-
inline virtual void initVecs(bool random)