(git:d18deda)
Loading...
Searching...
No Matches
cp_fm_diag Module Reference

used for collecting some of the diagonalization schemes available for cp_fm_type. cp_fm_power also moved here as it is very related More...

Functions/Subroutines

subroutine, public diag_init (diag_lib, fallback_applied, elpa_kernel, elpa_neigvec_min_input, elpa_qr, elpa_print, elpa_qr_unsafe, dlaf_neigvec_min_input, eps_check_diag_input)
 Setup the diagonalization library to be used.
 
subroutine, public diag_finalize ()
 Finalize the diagonalization library.
 
subroutine, public choose_eigv_solver (matrix, eigenvectors, eigenvalues, info)
 Choose the Eigensolver depending on which library is available ELPA seems to be unstable for small systems.
 
subroutine, public cp_fm_syevd (matrix, eigenvectors, eigenvalues, info)
 Computes all eigenvalues and vectors of a real symmetric matrix significantly faster than syevx, scales also much better. Needs workspace to allocate all the eigenvectors.
 
subroutine, public cp_fm_syevx (matrix, eigenvectors, eigenvalues, neig, work_syevx)
 compute eigenvalues and optionally eigenvectors of a real symmetric matrix using scalapack. If eigenvectors are required this routine will replicate a full matrix on each CPU... if more than a handful of vectors are needed, use cp_fm_syevd instead
 
subroutine, public cp_fm_svd (matrix_a, matrix_eigvl, matrix_eigvr_t, eigval, info)
 decomposes a quadratic matrix into its singular value decomposition
 
subroutine, public cp_fm_power (matrix, work, exponent, threshold, n_dependent, verbose, eigvals)
 ...
 
subroutine, public cp_fm_block_jacobi (matrix, eigenvectors, eigval, thresh, start_sec_block)
 ...
 
subroutine, public cp_fm_geeig (amatrix, bmatrix, eigenvectors, eigenvalues, work)
 General Eigenvalue Problem AX = BXE Single option version: Cholesky decomposition of B.
 
subroutine, public cp_fm_geeig_canon (amatrix, bmatrix, eigenvectors, eigenvalues, work, epseig)
 General Eigenvalue Problem AX = BXE Use canonical diagonalization : U*s**(-1/2)
 

Variables

real(kind=dp), parameter, public eps_check_diag_default = 5.0E-14_dp
 
integer, save, public diag_type = 0
 
integer, parameter, public fm_diag_type_scalapack = 101
 
integer, parameter, public fm_diag_type_elpa = 102
 
integer, parameter, public fm_diag_type_cusolver = 103
 
integer, parameter, public fm_diag_type_dlaf = 104
 
integer, parameter, public fm_diag_type_default = FM_DIAG_TYPE_SCALAPACK
 

Detailed Description

used for collecting some of the diagonalization schemes available for cp_fm_type. cp_fm_power also moved here as it is very related

Note
first version : most routines imported
History
  • unused Jacobi routines removed, cosmetics (05.04.06,MK)
Author
Joost VandeVondele (2003-08)

Function/Subroutine Documentation

◆ diag_init()

subroutine, public cp_fm_diag::diag_init ( character(len=*), intent(in)  diag_lib,
logical, intent(out)  fallback_applied,
integer, intent(in)  elpa_kernel,
integer, intent(in)  elpa_neigvec_min_input,
logical, intent(in)  elpa_qr,
logical, intent(in)  elpa_print,
logical, intent(in)  elpa_qr_unsafe,
integer, intent(in)  dlaf_neigvec_min_input,
real(kind=dp), intent(in)  eps_check_diag_input 
)

Setup the diagonalization library to be used.

Parameters
diag_libdiag_library flag from GLOBAL section in input
fallback_applied.TRUE. if support for the requested library was not compiled-in and fallback to ScaLAPACK was applied, .FALSE. otherwise.
elpa_kernelinteger that determines which ELPA kernel to use for diagonalization
elpa_neigvec_min_input...
elpa_qrlogical that determines if ELPA should try to use QR to accelerate the diagonalization procedure of suitably sized matrices
elpa_printlogical that determines if information about the ELPA diagonalization should be printed
elpa_qr_unsafelogical that enables potentially unsafe ELPA options
dlaf_neigvec_min_input...
eps_check_diag_input...
History
  • Add support for DLA-Future (05.09.2023, RMeli)
Author
MI 11.2013

Definition at line 150 of file cp_fm_diag.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ diag_finalize()

subroutine, public cp_fm_diag::diag_finalize

Finalize the diagonalization library.

Definition at line 208 of file cp_fm_diag.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ choose_eigv_solver()

subroutine, public cp_fm_diag::choose_eigv_solver ( type(cp_fm_type), intent(in)  matrix,
type(cp_fm_type), intent(in)  eigenvectors,
real(kind=dp), dimension(:), intent(out)  eigenvalues,
integer, intent(out), optional  info 
)

Choose the Eigensolver depending on which library is available ELPA seems to be unstable for small systems.

Parameters
matrix...
eigenvectors...
eigenvalues...
info...
info If present returns error code and prevents program stops.
Works currently only for cp_fm_syevd with ScaLAPACK. Other solvers will end the program regardless of PRESENT(info).
History
  • Do not use ELPA for small matrices and use instead ScaLAPACK as fallback (10.05.2021, MK)

Definition at line 228 of file cp_fm_diag.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cp_fm_syevd()

subroutine, public cp_fm_diag::cp_fm_syevd ( type(cp_fm_type), intent(in)  matrix,
type(cp_fm_type), intent(in)  eigenvectors,
real(kind=dp), dimension(:), intent(out)  eigenvalues,
integer, intent(out), optional  info 
)

Computes all eigenvalues and vectors of a real symmetric matrix significantly faster than syevx, scales also much better. Needs workspace to allocate all the eigenvectors.

Parameters
matrix...
eigenvectors...
eigenvalues...
info...
matrix is supposed to be in upper triangular form, and overwritten by this routine
info If present returns error code and prevents program stops.
Works currently only for scalapack. Other solvers will end the program regardless of PRESENT(info).

Definition at line 433 of file cp_fm_diag.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cp_fm_syevx()

subroutine, public cp_fm_diag::cp_fm_syevx ( type(cp_fm_type), intent(in)  matrix,
type(cp_fm_type), intent(in), optional  eigenvectors,
real(kind=dp), dimension(:), intent(out)  eigenvalues,
integer, intent(in), optional  neig,
real(kind=dp), intent(in), optional  work_syevx 
)

compute eigenvalues and optionally eigenvectors of a real symmetric matrix using scalapack. If eigenvectors are required this routine will replicate a full matrix on each CPU... if more than a handful of vectors are needed, use cp_fm_syevd instead

Parameters
matrix...
eigenvectors...
eigenvalues...
neig...
work_syevx...
matrix is supposed to be in upper triangular form, and overwritten by this routine
neig is the number of vectors needed (default all) work_syevx evec calculation only, is the fraction of the working buffer allowed (1.0 use full buffer) reducing this saves time, but might cause the routine to fail

Definition at line 677 of file cp_fm_diag.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cp_fm_svd()

subroutine, public cp_fm_diag::cp_fm_svd ( type(cp_fm_type), intent(in)  matrix_a,
type(cp_fm_type), intent(inout)  matrix_eigvl,
type(cp_fm_type), intent(inout)  matrix_eigvr_t,
real(kind=dp), dimension(:), intent(inout), pointer  eigval,
integer, intent(out), optional  info 
)

decomposes a quadratic matrix into its singular value decomposition

Parameters
matrix_a...
matrix_eigvl...
matrix_eigvr_t...
eigval...
info...
Author
Maximilian Graml

Definition at line 915 of file cp_fm_diag.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cp_fm_power()

subroutine, public cp_fm_diag::cp_fm_power ( type(cp_fm_type), intent(in)  matrix,
type(cp_fm_type), intent(in)  work,
real(kind=dp), intent(in)  exponent,
real(kind=dp), intent(in)  threshold,
integer, intent(out)  n_dependent,
logical, intent(in), optional  verbose,
real(kind=dp), dimension(2), intent(out), optional  eigvals 
)

...

Parameters
matrix...
work...
exponent...
threshold...
n_dependent...
verbose...
eigvals...

Definition at line 1037 of file cp_fm_diag.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cp_fm_block_jacobi()

subroutine, public cp_fm_diag::cp_fm_block_jacobi ( type(cp_fm_type), intent(in)  matrix,
type(cp_fm_type), intent(in)  eigenvectors,
real(kind=dp), dimension(:), intent(in)  eigval,
real(kind=dp), intent(in)  thresh,
integer, intent(in)  start_sec_block 
)

...

Parameters
matrix...
eigenvectors...
eigval...
thresh...
start_sec_block...

Definition at line 1203 of file cp_fm_diag.F.

Here is the caller graph for this function:

◆ cp_fm_geeig()

subroutine, public cp_fm_diag::cp_fm_geeig ( type(cp_fm_type), intent(in)  amatrix,
type(cp_fm_type), intent(in)  bmatrix,
type(cp_fm_type), intent(in)  eigenvectors,
real(kind=dp), dimension(:)  eigenvalues,
type(cp_fm_type), intent(in)  work 
)

General Eigenvalue Problem AX = BXE Single option version: Cholesky decomposition of B.

Parameters
amatrix...
bmatrix...
eigenvectors...
eigenvalues...
work...

Definition at line 1403 of file cp_fm_diag.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cp_fm_geeig_canon()

subroutine, public cp_fm_diag::cp_fm_geeig_canon ( type(cp_fm_type), intent(in)  amatrix,
type(cp_fm_type), intent(in)  bmatrix,
type(cp_fm_type), intent(in)  eigenvectors,
real(kind=dp), dimension(:), intent(out)  eigenvalues,
type(cp_fm_type), intent(in)  work,
real(kind=dp), intent(in)  epseig 
)

General Eigenvalue Problem AX = BXE Use canonical diagonalization : U*s**(-1/2)

Parameters
amatrix...
bmatrix...
eigenvectors...
eigenvalues...
work...
epseig...

Definition at line 1456 of file cp_fm_diag.F.

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ eps_check_diag_default

real(kind=dp), parameter, public cp_fm_diag::eps_check_diag_default = 5.0E-14_dp

Definition at line 86 of file cp_fm_diag.F.

◆ diag_type

integer, save, public cp_fm_diag::diag_type = 0

Definition at line 90 of file cp_fm_diag.F.

◆ fm_diag_type_scalapack

integer, parameter, public cp_fm_diag::fm_diag_type_scalapack = 101

Definition at line 104 of file cp_fm_diag.F.

◆ fm_diag_type_elpa

integer, parameter, public cp_fm_diag::fm_diag_type_elpa = 102

Definition at line 104 of file cp_fm_diag.F.

◆ fm_diag_type_cusolver

integer, parameter, public cp_fm_diag::fm_diag_type_cusolver = 103

Definition at line 104 of file cp_fm_diag.F.

◆ fm_diag_type_dlaf

integer, parameter, public cp_fm_diag::fm_diag_type_dlaf = 104

Definition at line 104 of file cp_fm_diag.F.

◆ fm_diag_type_default

integer, parameter, public cp_fm_diag::fm_diag_type_default = FM_DIAG_TYPE_SCALAPACK

Definition at line 115 of file cp_fm_diag.F.