This package contains a Modern Fortran implementation of the Davidson diagonalization algorithms. Different schemas are available to compute the correction.
Available correction methods are:
- DPR: Diagonal-Preconditioned-Residue
- GJD: Generalized Jacobi Davidson
The Davidson method is suitable for diagonal-dominant hermitian matrices, that are quite common in certain scientific problems like electronic structure. The Davidson method could be not practical for other kinds of hermitian matrices.
The following program call the eigensolver
subroutine from the davidson
module and computes
the lowest 3 eigenvalues and corresponding eigenvectors, using the GJD method with a tolerance
of 1e-8
and 100
maximum iteration.
program main
use numeric_kinds, only: dp
use davidson, only: generalized_eigensolver, generate_diagonal_dominant
implicit none
integer, parameter :: dim = 50
integer, parameter :: lowest = 3
complex(dp), dimension(dim, dim) :: matrix, second_matrix
real(dp), dimension(lowest) :: eigenvalues
complex(dp), dimension(dim, lowest) :: eigenvectors
real(dp) :: tolerance
integer:: max_dim_subspace, max_iterations, lowest
matrix = generate_diagonal_dominant(dim, 1d-4)
second_matrix = generate_diagonal_dominant(dim, 1d-4, 1.0_dp)
max_iterations = 1000
max_dim_subspace = 20
tolerance = 1d-8
call generalized_eigensolver(matrix, eigenvalues, eigenvectors, lowest, "GJD", max_iterations, &
tolerance, final_iterations, max_dim_subspace, second_matrix)
print *, eigenvalues
print *, eigenvectors
end program main
The helper generate_diagonal_dominant
function generates a diagonal dominant
matrix with entries to the diagonal close to row number (i=1, number_of_rows)
and random number of the order 1e-4
on the off-diagonal entries.
Variables:
matrix
(in) matrix to diagonalizeeigenvalues
(out) resulting eigenvalueseigenvectors
(out) resulting eigenvectorslowest
(in) number of eigenvalues to computemethod
(in) Either "DPR" or "GJD"max_iterations
(in) maximum number of iterationstolerance
(in) Numerical tolerance for convergencefinal_iterations
(output) returns the number of iterations that were needed to convergemax_dim_subspace
(in, optional) Dimension of the subspace of searchsecond_matrix
(in, optional) Optional matrix to compute the generalized eigenvalue problem
- Davidson diagonalization method and its applications to electronic structure calculations
- Numerical Methods for Large Eigenvalue Problem
To compile execute:
cmake -H. -Bbuild && cmake --build build
To use another compiler (e.g. ifort):
cmake -H. -Bbuild -DCMAKE_Fortran_COMPILER=ifort && cmake --build build
To run the test:
cmake -H. -Bbuild -DENABLE_TEST=ON && cmake --build build
cd build && ctest -V
To Debug compile as:
cmake -H. -Bbuild -DCMAKE_BUILD_TYPE=Debug && cmake --build build
This packages assumes that you have installed the following packages:
Optionally, If an MKL library is available the package will try to find it.