35 use,
intrinsic :: iso_fortran_env
74 integer,
public :: es_type
76 real(real64),
public :: tolerance
77 integer,
public :: es_maxiter
79 real(real64) :: imag_time
82 real(real64),
allocatable,
public :: diff(:, :)
83 integer,
public :: matvec
84 integer,
allocatable,
public :: converged(:)
87 type(preconditioner_t),
public :: pre
90 type(subspace_t) :: sdiag
92 integer :: rmmdiis_minimization_iter
94 logical,
public :: folded_spectrum
97 logical,
public :: orthogonalize_to_all
98 integer,
public :: conjugate_direction
99 logical,
public :: additional_terms
100 real(real64),
public :: energy_change_threshold
103 type(eigen_chebyshev_t),
public :: cheby_params
105 type(exponential_t) :: exponential_operator
111 integer,
public,
parameter :: &
121 subroutine eigensolver_init(eigens, namespace, gr, st, mc, space, deactivate_oracle)
122 type(eigensolver_t),
intent(out) :: eigens
123 type(namespace_t),
intent(in) :: namespace
124 type(grid_t),
intent(in) :: gr
125 type(states_elec_t),
intent(in) :: st
126 type(multicomm_t),
intent(in) :: mc
127 class(space_t),
intent(in) :: space
128 logical,
optional,
intent(in) :: deactivate_oracle
130 integer :: default_iter, default_es
131 real(real64) :: default_tol
170 if(st%parallel_in_states)
then
179 message(1) =
"The selected eigensolver is not parallel in states."
180 message(2) =
"Please use the rmmdiis or Chebyshev filtering eigensolvers."
189 select case(eigens%es_type)
204 call parse_variable(namespace,
'CGOrthogonalizeAll', .false., eigens%orthogonalize_to_all)
221 call parse_variable(namespace,
'CGDirection', option__cgdirection__fletcher, eigens%conjugate_direction)
234 call parse_variable(namespace,
'CGAdditionalTerms', .false., eigens%additional_terms)
235 if(eigens%additional_terms)
then
236 call messages_experimental(
"The additional terms for the CG eigensolver are not tested for all cases.")
254 call parse_variable(namespace,
'CGEnergyChangeThreshold', 0.1_real64, eigens%energy_change_threshold)
272 call parse_variable(namespace,
'EigensolverImaginaryTime', 0.1_real64, eigens%imag_time)
278 message(1) =
"Smearing of occupations is incompatible with imaginary time evolution."
296 call parse_variable(namespace,
'EigensolverMinimizationIter', 0, eigens%rmmdiis_minimization_iter)
313 eigens%cheby_params%n_lanczos)
353 eigens%cheby_params%optimize_degree)
367 eigens%cheby_params%bound_mixing)
379 eigens%cheby_params%n_iter)
398 eigens%folded_spectrum = .false.
407 default_tol = 1e-7_real64
409 call parse_variable(namespace,
'ConvRelDens', 1e-6_real64, default_tol)
410 default_tol = default_tol / 10.0_real64
412 call parse_variable(namespace,
'EigensolverTolerance', default_tol, eigens%tolerance)
428 call parse_variable(namespace,
'EigensolverMaxIter', default_iter, eigens%es_maxiter)
431 if(eigens%es_maxiter > default_iter)
then
432 call messages_write(
'You have specified a large number of eigensolver iterations (')
435 call messages_write(
'This is not a good idea as it might slow down convergence, even for', new_line = .
true.)
436 call messages_write(
'independent particles, as subspace diagonalization will not be used', new_line = .
true.)
445 safe_allocate(eigens%diff(1:st%nst, 1:st%nik))
446 eigens%diff(1:st%nst, 1:st%nik) = 0
448 safe_allocate(eigens%converged(1:st%nik))
449 eigens%converged(1:st%nik) = 0
452 call eigens%sdiag%init(namespace, st)
455 select case(eigens%es_type)
458 mem = (
m_two*eigens%es_maxiter -
m_one)*st%block_size*real(gr%np_part, real64)
462 mem = mem*16.0_real64
467 call messages_write(
' memory. This amount can be reduced by decreasing the value')
469 call messages_write(
' of the variable StatesBlockSize (currently set to ')
487 select case(eigens%es_type)
492 safe_deallocate_a(eigens%converged)
493 safe_deallocate_a(eigens%diff)
500 subroutine eigensolver_run(eigens, namespace, gr, st, hm, iter, conv, nstconv)
503 type(
grid_t),
intent(in) :: gr
506 integer,
intent(in) :: iter
507 logical,
optional,
intent(out) :: conv
508 integer,
optional,
intent(in) :: nstconv
511 integer :: ik, ist, nstconv_
513 logical :: conv_reduced
514 integer :: outcount, lmatvec
515 real(real64),
allocatable :: ldiff(:), leigenval(:)
516 real(real64),
allocatable :: ldiff_out(:), leigenval_out(:)
517 integer,
allocatable ::
lconv(:)
523 if(
present(conv)) conv = .false.
524 if(
present(nstconv))
then
536 ik_loop:
do ik = st%d%kpt%start, st%d%kpt%end
543 if (.not. eigens%folded_spectrum)
then
545 eigens%converged(ik) = 0
547 if(eigens%diff(ist, ik) < eigens%tolerance)
then
548 eigens%converged(ik) = ist
557 write(stdout,
'(1x)')
560 if(
present(conv)) conv = all(eigens%converged(st%d%kpt%start:st%d%kpt%end) >= nstconv_)
563 if (st%d%kpt%parallel)
then
564 if (
present(conv))
then
565 call st%d%kpt%mpi_grp%allreduce(conv, conv_reduced, 1, mpi_logical, mpi_land)
569 lmatvec = eigens%matvec
570 call st%d%kpt%mpi_grp%allreduce(lmatvec, eigens%matvec, 1, mpi_integer, mpi_sum)
572 safe_allocate(
lconv(1:st%d%kpt%nlocal))
573 lconv(1:st%d%kpt%nlocal) = eigens%converged(st%d%kpt%start:st%d%kpt%end)
575 assert(outcount == st%nik)
576 safe_deallocate_a(
lconv)
579 safe_allocate(ldiff(1:st%d%kpt%nlocal))
580 safe_allocate(leigenval(1:st%d%kpt%nlocal))
581 safe_allocate(ldiff_out(1:st%nik))
582 safe_allocate(leigenval_out(1:st%nik))
583 do ist = st%st_start, st%st_end
584 ldiff(1:st%d%kpt%nlocal) = eigens%diff(ist, st%d%kpt%start:st%d%kpt%end)
585 leigenval(1:st%d%kpt%nlocal) = st%eigenval(ist, st%d%kpt%start:st%d%kpt%end)
587 eigens%diff(ist, :) = ldiff_out
588 assert(outcount == st%nik)
589 call lmpi_gen_allgatherv(st%d%kpt%nlocal, leigenval, outcount, leigenval_out, st%d%kpt%mpi_grp)
590 st%eigenval(ist, :) = leigenval_out
591 assert(outcount == st%nik)
593 safe_deallocate_a(ldiff)
594 safe_deallocate_a(ldiff_out)
595 safe_deallocate_a(leigenval)
596 safe_deallocate_a(leigenval_out)
613 select case(this%es_type)
630 select case(this%es_type)
641#include "eigensolver_inc.F90"
642#include "eigen_plan_inc.F90"
643#include "eigen_evolution_inc.F90"
646#include "complex.F90"
647#include "eigensolver_inc.F90"
648#include "eigen_plan_inc.F90"
649#include "eigen_evolution_inc.F90"
This module implements batches of mesh functions.
This module implements common operations on batches of mesh functions.
type(debug_t), save, public debug
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
type(eigen_chebyshev_t), public default_chebyshev_params
Default Chebyshev input parameters Arguments 1 and 2 taken from 10.1016/j.jcp.2006....
subroutine eigensolver_run(eigens, namespace, gr, st, hm, iter, conv, nstconv)
subroutine zeigensolver_run(eigens, namespace, mesh, st, hm, iter, ik)
integer, parameter, public rs_evo
integer, parameter, public rs_cg
subroutine, public eigensolver_init(eigens, namespace, gr, st, mc, space, deactivate_oracle)
subroutine deigensolver_run(eigens, namespace, mesh, st, hm, iter, ik)
logical function eigensolver_parallel_in_states(this)
integer, parameter, public rs_chebyshev
logical function eigensolver_has_progress_bar(this)
subroutine, public eigensolver_end(eigens)
integer, parameter, public rs_rmmdiis
integer, parameter, public spinors
subroutine, public exponential_init(te, namespace, full_batch)
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
real(real64), parameter, public m_one
This module implements the underlying real-space grid.
A module to handle KS potential, without the external potential.
This module defines functions over batches of mesh functions.
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
character(len=512), private msg
subroutine, public messages_warning(no_lines, all_nodes, namespace)
subroutine, public messages_obsolete_variable(namespace, name, rep)
subroutine, public messages_new_line()
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_input_error(namespace, var, details, row, column)
subroutine, public messages_experimental(name, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
This module contains some common usage patterns of MPI routines.
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
type(mpi_grp_t), public mpi_world
This module handles the communicators for the various parallelization strategies.
logical function, public parse_is_defined(namespace, name)
subroutine, public preconditioner_end(this)
subroutine, public preconditioner_init(this, namespace, gr, mc, space)
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
integer, parameter, public smear_semiconductor
integer, parameter, public smear_fixed_occ
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
This module provides routines for communicating states when using states parallelization.
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
This module defines the unit system, used for input and output.
type(unit_t), public unit_megabytes
For large amounts of data (natural code units are bytes)
Description of the grid, containing information on derivatives, stencil, and symmetries.
The states_elec_t class contains all electronic wave functions.