Octopus
eigen_chebyshev_oct_m Module Reference

Data Types

type  batch_pointer_t
 
type  eigen_chebyshev_t
 

Functions/Subroutines

subroutine, public dchebyshev_filter_solver (namespace, sdiag, mesh, st, hm, ik, subspace_tol, filter_params, scf_iter, prior_residuals)
 Driver for Chebyshev filter-based solver. More...
 
subroutine dfirstscf_iterative_chebyshev_filter (namespace, sdiag, mesh, st, hm, ik, tolerance, filter_params)
 Iterative application of Chebyshev filter, for use with the first SCF step. More...
 
subroutine dchebyshev_filter (namespace, mesh, st, hm, degree, bounds, ik)
 Chebyshev Filter. More...
 
subroutine, public zchebyshev_filter_solver (namespace, sdiag, mesh, st, hm, ik, subspace_tol, filter_params, scf_iter, prior_residuals)
 Driver for Chebyshev filter-based solver. More...
 
subroutine zfirstscf_iterative_chebyshev_filter (namespace, sdiag, mesh, st, hm, ik, tolerance, filter_params)
 Iterative application of Chebyshev filter, for use with the first SCF step. More...
 
subroutine zchebyshev_filter (namespace, mesh, st, hm, degree, bounds, ik)
 Chebyshev Filter. More...
 

Variables

type(eigen_chebyshev_t), public, protected default_chebyshev_params = eigen_chebyshev_t(5, 25, M_HALF, 5, .true.)
 Default Chebyshev input parameters Arguments 1 and 2 taken from 10.1016/j.jcp.2006.03.017 Argument 3. If this were a fixed value for the whole calculation, 10 would be reasonable (see 10.1016/j.jcp.2014.06.056) however as the default is to use this as a max value, which is optimised by Octopus, any large value is reasonable. Argument 4 set empirically. A value > 2 results in a more accurate density. Values > 6 show minimal to no improvement in the number of SCF steps. More...
 

Function/Subroutine Documentation

◆ dchebyshev_filter_solver()

subroutine, public eigen_chebyshev_oct_m::dchebyshev_filter_solver ( type(namespace_t), intent(in)  namespace,
type(subspace_t), intent(in)  sdiag,
class(mesh_t), intent(in)  mesh,
type(states_elec_t), intent(inout)  st,
type(hamiltonian_elec_t), intent(in)  hm,
integer, intent(in)  ik,
real(real64), intent(in)  subspace_tol,
type(eigen_chebyshev_t), intent(in)  filter_params,
integer, intent(in)  scf_iter,
real(real64), dimension(:), intent(in)  prior_residuals 
)

Driver for Chebyshev filter-based solver.

This routine is implemented according to Algorithm 4.1 in the paper: ""Chebyshev-filtered subspace iteration method free of sparse diagonalization for solving the Kohn–Sham equation"". http:

The scaled Chebyshev algorithm is always utilised, as we get an estimate of the lowest eigenvalue of H from the lowest Ritz value of the prior step. The reason for the scaling is to prevent a potential overflow, which may happen if the filter degree is large or if the smallest eigenvalue is mapped far away from −1.

For the first SCF step, a Chebyshev polynomial of filter_paramsdegree is applied filter_paramsn_iter times to a set of search vectors found from the initial guess for the density.

Parameters
[in]namespaceCalling namespace
[in]sdiagSubspace diagonalisation choice
[in]meshReal-space mesh
[in,out]stEigenstates
[in]hmHamiltonian
[in]ikk-point index
[in]subspace_tolSubspace iterative solver tolerance
[in]filter_paramsChebyshev filter parameters
[in]scf_iterSCF iteration
[in]prior_residualsEigenvalue residuals from prior SCF step

Definition at line 210 of file eigen_chebyshev.F90.

◆ dfirstscf_iterative_chebyshev_filter()

subroutine eigen_chebyshev_oct_m::dfirstscf_iterative_chebyshev_filter ( type(namespace_t), intent(in)  namespace,
type(subspace_t), intent(in)  sdiag,
class(mesh_t), intent(in)  mesh,
type(states_elec_t), intent(inout)  st,
type(hamiltonian_elec_t), intent(in)  hm,
integer, intent(in)  ik,
real(real64), intent(in)  tolerance,
type(eigen_chebyshev_t), intent(in)  filter_params 
)
private

Iterative application of Chebyshev filter, for use with the first SCF step.

The initial search vectors are found using the initial for the density (sum of atomic densities LCAO, etc), and returned in st: stgrouppsib.

Parameters
[in]namespaceCalling namespace
[in]sdiagSubspace diagonalisation choice
[in]meshReal-space mesh
[in,out]stInitial guess at search vectors, and much more
[in]hmHamiltonian
[in]ikk-point index
[in]toleranceSubspace iterative solver tolerance
[in]filter_paramsChebyshev filter parameters

Definition at line 266 of file eigen_chebyshev.F90.

◆ dchebyshev_filter()

subroutine eigen_chebyshev_oct_m::dchebyshev_filter ( type(namespace_t), intent(in)  namespace,
type(mesh_t), intent(in)  mesh,
type(states_elec_t), intent(inout), target  st,
class(hamiltonian_elec_t), intent(in)  hm,
integer, dimension(:), intent(in)  degree,
type(chebyshev_filter_bounds_t), intent(in)  bounds,
integer, intent(in)  ik 
)
private

Chebyshev Filter.

Filter an eigenspectrum by an m-degree Chebyshev polynomial that dampens on the interval [a,b]. Based on algorithm 3.2 of ""Chebyshev-filtered subspace iteration method free of sparse diagonalization for solving the Kohn–Sham equation"". http:

Application of the simple or scaled filter depends upon the choice of bounds. In the simple case, sigma will reduce to 1.

Parameters
[in]namespaceCalling namespace
[in]meshReal-space mesh
[in]hmHamiltonian
[in]degreeChebyshev polynomial degree, per block
[in]boundsPolynomial filter bounds of the subspace to dampen
[in]ikk-point index
[in,out]stKS containing input eigenvectors and

Definition at line 372 of file eigen_chebyshev.F90.

◆ zchebyshev_filter_solver()

subroutine, public eigen_chebyshev_oct_m::zchebyshev_filter_solver ( type(namespace_t), intent(in)  namespace,
type(subspace_t), intent(in)  sdiag,
class(mesh_t), intent(in)  mesh,
type(states_elec_t), intent(inout)  st,
type(hamiltonian_elec_t), intent(in)  hm,
integer, intent(in)  ik,
real(real64), intent(in)  subspace_tol,
type(eigen_chebyshev_t), intent(in)  filter_params,
integer, intent(in)  scf_iter,
real(real64), dimension(:), intent(in)  prior_residuals 
)

Driver for Chebyshev filter-based solver.

This routine is implemented according to Algorithm 4.1 in the paper: ""Chebyshev-filtered subspace iteration method free of sparse diagonalization for solving the Kohn–Sham equation"". http:

The scaled Chebyshev algorithm is always utilised, as we get an estimate of the lowest eigenvalue of H from the lowest Ritz value of the prior step. The reason for the scaling is to prevent a potential overflow, which may happen if the filter degree is large or if the smallest eigenvalue is mapped far away from −1.

For the first SCF step, a Chebyshev polynomial of filter_paramsdegree is applied filter_paramsn_iter times to a set of search vectors found from the initial guess for the density.

Parameters
[in]namespaceCalling namespace
[in]sdiagSubspace diagonalisation choice
[in]meshReal-space mesh
[in,out]stEigenstates
[in]hmHamiltonian
[in]ikk-point index
[in]subspace_tolSubspace iterative solver tolerance
[in]filter_paramsChebyshev filter parameters
[in]scf_iterSCF iteration
[in]prior_residualsEigenvalue residuals from prior SCF step

Definition at line 547 of file eigen_chebyshev.F90.

◆ zfirstscf_iterative_chebyshev_filter()

subroutine eigen_chebyshev_oct_m::zfirstscf_iterative_chebyshev_filter ( type(namespace_t), intent(in)  namespace,
type(subspace_t), intent(in)  sdiag,
class(mesh_t), intent(in)  mesh,
type(states_elec_t), intent(inout)  st,
type(hamiltonian_elec_t), intent(in)  hm,
integer, intent(in)  ik,
real(real64), intent(in)  tolerance,
type(eigen_chebyshev_t), intent(in)  filter_params 
)
private

Iterative application of Chebyshev filter, for use with the first SCF step.

The initial search vectors are found using the initial for the density (sum of atomic densities LCAO, etc), and returned in st: stgrouppsib.

Parameters
[in]namespaceCalling namespace
[in]sdiagSubspace diagonalisation choice
[in]meshReal-space mesh
[in,out]stInitial guess at search vectors, and much more
[in]hmHamiltonian
[in]ikk-point index
[in]toleranceSubspace iterative solver tolerance
[in]filter_paramsChebyshev filter parameters

Definition at line 603 of file eigen_chebyshev.F90.

◆ zchebyshev_filter()

subroutine eigen_chebyshev_oct_m::zchebyshev_filter ( type(namespace_t), intent(in)  namespace,
type(mesh_t), intent(in)  mesh,
type(states_elec_t), intent(inout), target  st,
class(hamiltonian_elec_t), intent(in)  hm,
integer, dimension(:), intent(in)  degree,
type(chebyshev_filter_bounds_t), intent(in)  bounds,
integer, intent(in)  ik 
)
private

Chebyshev Filter.

Filter an eigenspectrum by an m-degree Chebyshev polynomial that dampens on the interval [a,b]. Based on algorithm 3.2 of ""Chebyshev-filtered subspace iteration method free of sparse diagonalization for solving the Kohn–Sham equation"". http:

Application of the simple or scaled filter depends upon the choice of bounds. In the simple case, sigma will reduce to 1.

Parameters
[in]namespaceCalling namespace
[in]meshReal-space mesh
[in]hmHamiltonian
[in]degreeChebyshev polynomial degree, per block
[in]boundsPolynomial filter bounds of the subspace to dampen
[in]ikk-point index
[in,out]stKS containing input eigenvectors and

Definition at line 709 of file eigen_chebyshev.F90.

Variable Documentation

◆ default_chebyshev_params

type(eigen_chebyshev_t), public, protected eigen_chebyshev_oct_m::default_chebyshev_params = eigen_chebyshev_t(5, 25, M_HALF, 5, .true.)

Default Chebyshev input parameters Arguments 1 and 2 taken from 10.1016/j.jcp.2006.03.017 Argument 3. If this were a fixed value for the whole calculation, 10 would be reasonable (see 10.1016/j.jcp.2014.06.056) however as the default is to use this as a max value, which is optimised by Octopus, any large value is reasonable. Argument 4 set empirically. A value > 2 results in a more accurate density. Values > 6 show minimal to no improvement in the number of SCF steps.

Definition at line 149 of file eigen_chebyshev.F90.