(git:6a2e663)
qs_scf_diagonalization Module Reference

Different diagonalization schemes that can be used for the iterative solution of the eigenvalue problem. More...

Functions/Subroutines

subroutine, public general_eigenproblem (scf_env, mos, matrix_ks, matrix_s, scf_control, scf_section, diis_step)
 the inner loop of scf, specific to diagonalization with S matrix basically, in goes the ks matrix out goes a new p matrix More...
 
subroutine, public do_general_diag (scf_env, mos, matrix_ks, matrix_s, scf_control, scf_section, diis_step)
 ... More...
 
subroutine, public do_general_diag_kp (matrix_ks, matrix_s, kpoints, scf_env, scf_control, update_p, diis_step, diis_error, qs_env)
 Kpoint diagonalization routine Transforms matrices to kpoint, distributes kpoint groups, performs general diagonalization (no storgae of overlap decomposition), stores MOs, calculates occupation numbers, calculates density matrices in kpoint representation, transforms density matrices to real space. More...
 
subroutine, public do_scf_diag_subspace (qs_env, scf_env, subspace_env, mos, rho, ks_env, scf_section, scf_control)
 inner loop within MOS subspace, to refine occupation and density, before next diagonalization of the Hamiltonian More...
 
subroutine, public diag_subspace_allocate (subspace_env, qs_env, mos)
 ... More...
 
subroutine, public do_special_diag (scf_env, mos, matrix_ks, scf_control, scf_section, diis_step)
 the inner loop of scf, specific to diagonalization without S matrix basically, in goes the ks matrix out goes a new p matrix More...
 
subroutine, public do_ot_diag (scf_env, mos, matrix_ks, matrix_s, scf_control, scf_section, diis_step)
 the inner loop of scf, specific to iterative diagonalization using OT with S matrix; basically, in goes the ks matrix out goes a new p matrix More...
 
subroutine, public do_roks_diag (scf_env, mos, matrix_ks, matrix_s, scf_control, scf_section, diis_step, orthogonal_basis)
 Solve a set restricted open Kohn-Sham (ROKS) equations based on the alpha and beta Kohn-Sham matrices from unrestricted Kohn-Sham. More...
 
subroutine, public do_block_krylov_diag (scf_env, mos, matrix_ks, scf_control, scf_section, check_moconv_only)
 iterative diagonalization using the block Krylov-space approach More...
 
subroutine, public do_block_davidson_diag (qs_env, scf_env, mos, matrix_ks, matrix_s, scf_control, scf_section, check_moconv_only)
 iterative diagonalization using the block davidson space approach More...
 

Detailed Description

Different diagonalization schemes that can be used for the iterative solution of the eigenvalue problem.

History
started from routines previously located in the qs_scf module 05.2009

Function/Subroutine Documentation

◆ general_eigenproblem()

subroutine, public qs_scf_diagonalization::general_eigenproblem ( type(qs_scf_env_type), pointer  scf_env,
type(mo_set_type), dimension(:), intent(in)  mos,
type(dbcsr_p_type), dimension(:), pointer  matrix_ks,
type(dbcsr_p_type), dimension(:), pointer  matrix_s,
type(scf_control_type), pointer  scf_control,
type(section_vals_type), pointer  scf_section,
logical, intent(inout)  diis_step 
)

the inner loop of scf, specific to diagonalization with S matrix basically, in goes the ks matrix out goes a new p matrix

Parameters
scf_env...
mos...
matrix_ks...
matrix_s...
scf_control...
scf_section...
diis_step...
History
03.2006 created [Joost VandeVondele]

Definition at line 148 of file qs_scf_diagonalization.F.

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

◆ do_general_diag()

subroutine, public qs_scf_diagonalization::do_general_diag ( type(qs_scf_env_type), pointer  scf_env,
type(mo_set_type), dimension(:), intent(inout)  mos,
type(dbcsr_p_type), dimension(:), pointer  matrix_ks,
type(dbcsr_p_type), dimension(:), pointer  matrix_s,
type(scf_control_type), pointer  scf_control,
type(section_vals_type), pointer  scf_section,
logical, intent(inout)  diis_step 
)

...

Parameters
scf_env...
mos...
matrix_ks...
matrix_s...
scf_control...
scf_section...
diis_step...

Definition at line 318 of file qs_scf_diagonalization.F.

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

◆ do_general_diag_kp()

subroutine, public qs_scf_diagonalization::do_general_diag_kp ( type(dbcsr_p_type), dimension(:, :), pointer  matrix_ks,
type(dbcsr_p_type), dimension(:, :), pointer  matrix_s,
type(kpoint_type), pointer  kpoints,
type(qs_scf_env_type), pointer  scf_env,
type(scf_control_type), pointer  scf_control,
logical, intent(in)  update_p,
logical, intent(inout)  diis_step,
real(dp), intent(inout), optional  diis_error,
type(qs_environment_type), optional, pointer  qs_env 
)

Kpoint diagonalization routine Transforms matrices to kpoint, distributes kpoint groups, performs general diagonalization (no storgae of overlap decomposition), stores MOs, calculates occupation numbers, calculates density matrices in kpoint representation, transforms density matrices to real space.

Parameters
matrix_ksKohn-sham matrices (RS indices, global)
matrix_sOverlap matrices (RS indices, global)
kpointsKpoint environment
scf_envSCF environment
scf_controlSCF control variables
update_p...
diis_step...
diis_error...
qs_env...
History
08.2014 created [JGH]

Definition at line 373 of file qs_scf_diagonalization.F.

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

◆ do_scf_diag_subspace()

subroutine, public qs_scf_diagonalization::do_scf_diag_subspace ( type(qs_environment_type), pointer  qs_env,
type(qs_scf_env_type), pointer  scf_env,
type(subspace_env_type), pointer  subspace_env,
type(mo_set_type), dimension(:), intent(inout)  mos,
type(qs_rho_type), pointer  rho,
type(qs_ks_env_type), pointer  ks_env,
type(section_vals_type), pointer  scf_section,
type(scf_control_type), pointer  scf_control 
)

inner loop within MOS subspace, to refine occupation and density, before next diagonalization of the Hamiltonian

Parameters
qs_env...
scf_env...
subspace_env...
mos...
rho...
ks_env...
scf_section...
scf_control...
History
09.2009 created [MI]
Note
it is assumed that when diagonalization is used, also some mixing procedure is active

Definition at line 725 of file qs_scf_diagonalization.F.

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

◆ diag_subspace_allocate()

subroutine, public qs_scf_diagonalization::diag_subspace_allocate ( type(subspace_env_type), pointer  subspace_env,
type(qs_environment_type), pointer  qs_env,
type(mo_set_type), dimension(:), intent(in)  mos 
)

...

Parameters
subspace_env...
qs_env...
mos...

Definition at line 934 of file qs_scf_diagonalization.F.

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

◆ do_special_diag()

subroutine, public qs_scf_diagonalization::do_special_diag ( type(qs_scf_env_type), pointer  scf_env,
type(mo_set_type), dimension(:), intent(inout)  mos,
type(dbcsr_p_type), dimension(:), pointer  matrix_ks,
type(scf_control_type), pointer  scf_control,
type(section_vals_type), pointer  scf_section,
logical, intent(inout)  diis_step 
)

the inner loop of scf, specific to diagonalization without S matrix basically, in goes the ks matrix out goes a new p matrix

Parameters
scf_env...
mos...
matrix_ks...
scf_control...
scf_section...
diis_step...
History
03.2006 created [Joost VandeVondele]

Definition at line 1003 of file qs_scf_diagonalization.F.

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

◆ do_ot_diag()

subroutine, public qs_scf_diagonalization::do_ot_diag ( type(qs_scf_env_type), pointer  scf_env,
type(mo_set_type), dimension(:), intent(inout)  mos,
type(dbcsr_p_type), dimension(:), pointer  matrix_ks,
type(dbcsr_p_type), dimension(:), pointer  matrix_s,
type(scf_control_type), pointer  scf_control,
type(section_vals_type), pointer  scf_section,
logical, intent(inout)  diis_step 
)

the inner loop of scf, specific to iterative diagonalization using OT with S matrix; basically, in goes the ks matrix out goes a new p matrix

Parameters
scf_env...
mos...
matrix_ks...
matrix_s...
scf_control...
scf_section...
diis_step...
History
10.2008 created [JGH]

Definition at line 1099 of file qs_scf_diagonalization.F.

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

◆ do_roks_diag()

subroutine, public qs_scf_diagonalization::do_roks_diag ( type(qs_scf_env_type), pointer  scf_env,
type(mo_set_type), dimension(:), intent(in)  mos,
type(dbcsr_p_type), dimension(:), pointer  matrix_ks,
type(dbcsr_p_type), dimension(:), pointer  matrix_s,
type(scf_control_type), pointer  scf_control,
type(section_vals_type), pointer  scf_section,
logical, intent(inout)  diis_step,
logical, intent(in)  orthogonal_basis 
)

Solve a set restricted open Kohn-Sham (ROKS) equations based on the alpha and beta Kohn-Sham matrices from unrestricted Kohn-Sham.

Parameters
scf_env...
mos...
matrix_ks...
matrix_s...
scf_control...
scf_section...
diis_step...
orthogonal_basis...
History
04.2006 created [MK] Revised (01.05.06,MK)
Note
this is only a high-spin ROKS.

Definition at line 1205 of file qs_scf_diagonalization.F.

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

◆ do_block_krylov_diag()

subroutine, public qs_scf_diagonalization::do_block_krylov_diag ( type(qs_scf_env_type), pointer  scf_env,
type(mo_set_type), dimension(:), intent(inout)  mos,
type(dbcsr_p_type), dimension(:), pointer  matrix_ks,
type(scf_control_type), pointer  scf_control,
type(section_vals_type), pointer  scf_section,
logical, intent(in), optional  check_moconv_only 
)

iterative diagonalization using the block Krylov-space approach

Parameters
scf_env...
mos...
matrix_ks...
scf_control...
scf_section...
check_moconv_only...

Definition at line 1479 of file qs_scf_diagonalization.F.

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

◆ do_block_davidson_diag()

subroutine, public qs_scf_diagonalization::do_block_davidson_diag ( type(qs_environment_type), pointer  qs_env,
type(qs_scf_env_type), pointer  scf_env,
type(mo_set_type), dimension(:), intent(inout)  mos,
type(dbcsr_p_type), dimension(:), pointer  matrix_ks,
type(dbcsr_p_type), dimension(:), pointer  matrix_s,
type(scf_control_type), pointer  scf_control,
type(section_vals_type), pointer  scf_section,
logical, intent(in), optional  check_moconv_only 
)

iterative diagonalization using the block davidson space approach

Parameters
qs_env...
scf_env...
mos...
matrix_ks...
matrix_s...
scf_control...
scf_section...
check_moconv_only...

Definition at line 1679 of file qs_scf_diagonalization.F.

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