1.1.3. chase::ChaseConfig
-
template<class T>
class ChaseConfig A class to set up all the parameters of the eigensolver.
Besides setting up the standard parameters such as size of the matrix
N_defining the eigenproblem, number of wanted eigenvaluesnev_, the public functions of this class initialize all internal parameters and allow the experienced user to set up the values of parameters of core functionalities (e.g. lanczos DoS). The aim is to influence the behavior of the library in special cases when the default values of the parameters return a suboptimal efficiency in terms of performance and/or accuracy.Public Functions
-
inline ChaseConfig(std::size_t _N, std::size_t _nev, std::size_t _nex)
Constructor for the ChaseConfig class.
Requires the explicit values for the initalization of the size
N_of the matrix A, the number of sought after extremal eigenvaluesnev_, and the number of extra eigenvaluenex_which defines, together withnev_, the search space. All the other private members of the class are initialized using default values either specified directly (e.g.max_iter_) or specified using a function that is part of the ChaseConfig namespace (e.g.initMaxDeg).- 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.
-
inline bool UseApprox() const
Returns the value of the
approx_flag.The value of the
approx_flag indicates if ChASE has been used with the engagement of approximate solutions as it is typical when solving for sequences of eigenvalue problems occurring in Density Functional Theory.- Returns:
The value of the
approx_flag.
-
inline bool DoOptimization() const
Sets the
approx_flag to eithertrueorfalse.This function is used to change the value of
approx_so that the eigensolver can use approximate solutions inputed through a matrix of initial vectors.- Parameters:
flag – A boolean parameter which admits either a
trueorfalsevalue.
-
inline void SetApprox(bool flag)
Sets the
optimization_flag to eithertrueorfalse.This function is used to change the value of
optimization_so that the eigensolver minimizes the number of FLOPs needed to reach convergence for the entire sought after subspace of the spectrum.- Parameters:
flag – A boolean parameter which admits either a
trueorfalsevalue.
-
inline void SetOpt(bool flag)
Returns the value of the
optimization_flag.The value of the
optimization_flag indicates when ChASE computes a polynomial degree optimized for each single desired eigenpairs. The optimization minimizes the number of operations required for the eigenpairs to have a residual which is just below the specified tolerance threshold.- Returns:
The value of the
optimization_flag
-
inline std::size_t GetMaxDeg() const
Returns the integer value of the maximum degree used by the polynomial filter.
The value of
max_deg_indicates the upper bound for the degree of the polynomial for any of the vectors filtered. Such bound is important to avoid potential numerical instabilities that may occur and impede the convergence of the eigenpairs, expecially the one close to the upper end of the desired subspace of the spectrum.- Returns:
The value of the maximum degree of the Chebyshev filter.
-
inline void SetMaxDeg(std::size_t _maxDeg)
Sets the maximum value of the degree of the Chebyshev filter.
When
optimization_ == 'true', the Chebyshev filter degree is computed automatically. Because the computed values could be quite large for eigenvectors at the end of the sought after spectrum, a maximum value is set to avoid numerical instabilities that may trigger eigenpairs divergence.- Parameters:
_maxDeg – This value should be set by the expert user. It is set to 36 by default. It can be lowered in case of the onset of early instabilities but it should not be lower than 20-25 to avoid the filter becomes ineffective. It can be increased whenever it is known there is a spectral gap between the value of
nev_and the value ofnev_ + nex_. It is strongly suggested to never exceed the value of 70.
-
inline std::size_t GetDegExtra() const
Returns the extra degree added to the polynomial filter.
When
optimization_ == 'true', each vector is filtered with a polynomial of a given calculated degree. Because the degree is predicted based on an heuristic fomula, such degree is augmented by a small value to ensure that the residual of the corresponding vector will be almost always lower than the required threshold tolerance.- Returns:
The extra value added to a computed vector of polynomial degrees.
-
inline void SetDegExtra(std::size_t degExtra)
Sets the value of the extra degree added to the polynomial filter.
The value of
degExtrashould be a single digit number usually not exceeding 5 or 6. The expert user can modify the default value (which is 2) in those cases where the heuristic to automatically compute the vector of optimal degrees seems to underestimate the value of the degree necessary for the eigenpairs to be declared converged.- Parameters:
degExtra – Value of the extra degree.
-
inline std::size_t GetMaxIter() const
Returns the value of the maximum number of subspace iterations allowed within ChASE.
In order to avoid that the eigensolver would runoff unchecked, ChASE is given a upper bound on the number of subspace iteration it can execute. This is a normal safety mechanism when the algorithm fails to converge some eigenpairs whose residuals would continue to oscillate without above the tolerance threshold.
- Returns:
The value of the maximum number of subspace iterations allowed.
-
inline void SetMaxIter(std::size_t maxIter)
Sets the value of the maximum number of subspace iterations within ChASE.
Typically ChASE requires a number of single digit iterations to converge. In extreme cases such number can grow up to 10 or 12. An increasing number of iterations usually implies that the eigensolver is not used for the intended purpose or that the spectrum of eigenproblem tackled is particularly problematic. The default value of
maxIteris set to 25.- Parameters:
maxIter – Value of the maximum number of subspace iterations.
-
inline std::size_t GetDeg() const
Returns the degree of the Chebyshev filter used by ChASE.
The value returned is the degree used by the filter when it is called (when
optimization_ == 'true'this value is used only the first time the filter is called)- Returns:
The value used by the Chebyshev filter
-
inline void SetDeg(std::size_t _deg)
Set the value of the initial degree of the Chebyshev filter.
Depending if the
optimization_parameter is set tofalseortrue, the value of_degis used by the Chebyshev filter respectively every time or just the first time it is called.- Parameters:
_deg – Value set by the expert user and should in general be between 10 and 25. The default value is 20. If a odd value is inserted, the function makes it even. This is necessary due to the swapping of the matrix of vectors within the filter. It is strongly suggested to avoid values above the higher between 40 and the value returned by GetMaxDeg().
-
inline double GetTol() const
Returns the threshold value of the eigenpair’s residual tolerance.
The value of the tolerance is used as threshold for all the residuals of the desired eigenpaits. Whenever an eigenpair’s residual decreases below such a value it is declared as converged, and is consequently deflated and locked. For an eigenpair \((\lambda_i,V_i)\), the residual in ChASE is defined as the Euclidean norm: \(||AV_i-\lambda_iV_i||_2\).
- Returns:
The value of the
tol_parameter.
-
inline void SetTol(double _tol)
Sets the value of the threshold of the eigenpair’s residual tolerance.
The value of the tolerance should be set carefully keeping in mind that the residual of the eigenpairs is limited by the accuracy of the dense eigensolver used within the Rayleigh-Ritz procedure. As such it should hardly be set below \(1e-14\) in double precision. As a rule of thumb a minimum value of 1e-04 and 1e-08 should be used respectively in single and double precision.
- Parameters:
_tol – A type double number usually specified in scientific notation (e.g. 1e-10).
-
inline std::size_t GetLanczosIter() const
Returns the number of Lanczos iterations executed by ChASE.
In order to estimate the spectral bounds, ChASE executes a limited number of Lanczos steps. These steps are then used to compute a spectral estimate based on the Density of State (DoS) algorithm.
- Returns:
The total number of the Lanczos iterations used by the DoS algorithm.
-
inline void SetLanczosIter(std::size_t lanczosIter)
Sets the number of Lanczos iterations executed by ChASE.
For the DoS algorithm to work effectively without overburdening the eigensolver, the number of Lanczos iteration should be not less than 10 but also no more than 100 . ChASE does not need very precise spectral estimates because at each iteration such estimates are automatically improved by the approximate spectrum computed. This is the reason why the default value of
lanczos_iter_is 25.- Parameters:
lanczosIter – Value of the total number of Lanczos iterations executed by ChASE.
-
inline std::size_t GetNumLanczos() const
Returns the number of stochastic vectors used for the spectral estimates.
After having executed a number of Lanczos steps, ChASE uses a cheap and efficient estimator to calculate the value of the upper extremum of the search space. Such an estimator uses a small number of stochastic vectors.
- Returns:
Number of stochastic vectors used by ChASE for the spectral estimate.
-
inline void SetNumLanczos(std::size_t lanczosIter)
Sets the number of stochastic vectors used for the spectral estimates.
The stochastic estimator used by ChASE is based on a cheap and efficient DoS algorithm. Because ChASE does not need precise estimates of the upper extremum of the search space, the number of vectors used is quite small. The default value used is 4. The expert user can change the value to a larger number (it is not suggested to use a smaller value) and pay a slightly higher computing cost. It is not suggested to exceed a value for
num_lanczos_higher than 20.- Parameters:
lanczosIter – Number of stochastic vectors used by ChASE.
-
inline std::size_t GetN() const
Returns the size of the eigenproblem.
This function returns the size of the matrix defining the standard or the generalized eigenvalue problem.
- Returns:
Rank of the matrix.
-
inline std::size_t GetNev() const
Returns the number of desired eigenpairs.
The number of sought after eigenpairs as also specified in the constructor ChASEConfig.
- Returns:
Number of desired eigenpairs.
-
inline std::size_t GetNex() const
Returns the number of extra eigenpairs that are used to augment the search subspace.
ChASE effectively uses an enlarged subspace to improve the convergence of the algorithm. With a very small value of
nex_, the eigenpairs at the end of the desired interval may have hard time to converge. By including a small but substantial number of extra values (in most case no more than 20% ofnev_), ensures that ChASE converges smoothly without stagnating.- Returns:
Number of extra eigenpairs augmenting the search space.
-
inline void SetCholQR(bool flag)
Sets the
cholqr_flag to eithertrueorfalse.This function is used to change the value of
cholqr_so either flexible CholQR (true) or Househoulder QR (false) will be used.- Parameters:
flag – A boolean parameter which admits either a
trueorfalsevalue.
-
inline bool DoCholQR()
Return the value of
cholqr_
-
inline void EnableSymCheck(bool flag)
-
inline bool DoSymCheck()
-
inline float GetDecayingRate() const
Returns the decaying rate for the polynomial lower bound.
The lower bound of the chebyshev lower bound is set based on an approximation of the eigenvalues by few iterations of lanczos. It might be better to use under estimation of the lower bound in certain cases, for instance if the target eigenvalues are packed.
-
inline void SetDecayingRate(float decayingRate)
Sets the decaying rate for the polynomial lower bound.
The lower bound of the chebyshev lower bound is set based on an approximation of the eigenvalues by few iterations of lanczos. It might be better to use under estimation of the lower bound in certain cases, for instance if the target eigenvalues are packed.
-
inline bool UseClusterAwareDegrees() const
Returns whether cluster-aware degree optimization is enabled.
When enabled, the algorithm detects clusters of eigenvalues and adjusts polynomial degrees accordingly to improve convergence for clustered eigenvalues.
-
inline void SetClusterAwareDegrees(bool flag)
Sets whether to use cluster-aware degree optimization.
This enables/disables the cluster detection algorithm that adjusts polynomial degrees based on eigenvalue clustering patterns.
-
inline float GetUpperbScaleRate() const
Returns the scale rate for upperb based on its sign.
This variable controls how upperb is scaled based on its sign. For positive upperb, it’s multiplied by this rate. For negative upperb, it’s multiplied by (2 - rate). Default value is 1.2.
-
inline void SetUpperbScaleRate(float upperbScaleRate)
Sets the scale rate for upperb based on its sign.
This variable controls how upperb is scaled based on its sign. For positive upperb, it’s multiplied by this rate. For negative upperb, it’s multiplied by (2 - rate). Default value is 1.2.
-
inline ChaseConfig(std::size_t _N, std::size_t _nev, std::size_t _nex)