33 use,
intrinsic :: iso_fortran_env
61 integer,
public,
parameter :: &
68 integer,
public :: exp_method
69 real(real64) :: lanczos_tol
70 real(real64) :: chebyshev_tol
71 integer,
public :: exp_order
73 logical,
public :: full_batch = .false.
84 type(exponential_t),
intent(out) :: te
85 type(namespace_t),
intent(in) :: namespace
86 logical,
optional,
intent(in) :: full_batch
145 select case (te%exp_method)
158 call parse_variable(namespace,
'TDChebyshevTol', 1e-10_real64, te%chebyshev_tol)
172 call parse_variable(namespace,
'TDLanczosTol', 1e-5_real64, te%lanczos_tol)
189 call parse_variable(namespace,
'TDExpOrder', default__tdexporder, te%exp_order)
194 te%arnoldi_gs = option__arnoldiorthogonalization__cgs
195 if (te%exp_method == exp_lanczos)
then
210 call parse_variable(namespace,
'ArnoldiOrthogonalization', option__arnoldiorthogonalization__cgs, &
227 teo%exp_method = tei%exp_method
228 teo%lanczos_tol = tei%lanczos_tol
229 teo%exp_order = tei%exp_order
230 teo%arnoldi_gs = tei%arnoldi_gs
240 class(
mesh_t),
intent(in) :: mesh
242 integer,
intent(in) :: ist
243 integer,
intent(in) :: ik
244 complex(real64),
contiguous,
intent(inout) :: zpsi(:, :)
245 real(real64),
intent(in) :: deltat
246 real(real64),
optional,
intent(in) :: vmagnus(mesh%np, hm%d%nspin, 2)
247 logical,
optional,
intent(in) :: imag_time
250 complex(real64),
allocatable :: zpsi_inh(:, :)
256 if (hm%phase%is_allocated())
then
257 call hm%phase%apply_to_single(zpsi, mesh%np, hm%d%dim, ik, .false.)
263 safe_allocate(zpsi_inh(1:mesh%np_part, 1:hm%d%dim))
265 call wfs_elec_init(inh_psib, hm%d%dim, ist, ist, zpsi_inh, ik)
266 call te%apply_batch(namespace, mesh, hm, psib, deltat, &
267 vmagnus=vmagnus, imag_time=imag_time, inh_psib=inh_psib)
269 safe_deallocate_a(zpsi_inh)
271 call te%apply_batch(namespace, mesh, hm, psib, deltat, &
272 vmagnus=vmagnus, imag_time=imag_time)
277 if (hm%phase%is_allocated())
then
278 call hm%phase%apply_to_single(zpsi, mesh%np, hm%d%dim, ik, .
true.)
302 subroutine exponential_apply_batch(te, namespace, mesh, hm, psib, deltat, psib2, deltat2, vmagnus, imag_time, inh_psib)
305 class(
mesh_t),
intent(in) :: mesh
307 class(
batch_t),
intent(inout) :: psib
308 real(real64),
intent(in) :: deltat
309 class(
batch_t),
optional,
intent(inout) :: psib2
310 real(real64),
optional,
intent(in) :: deltat2
311 real(real64),
optional,
intent(in) :: vmagnus(:,:,:)
312 logical,
optional,
intent(in) :: imag_time
313 class(
batch_t),
optional,
intent(inout) :: inh_psib
315 complex(real64) :: deltat_, deltat2_
317 logical :: imag_time_
324 assert(
present(psib2) .eqv.
present(deltat2))
325 if (
present(inh_psib))
then
326 assert(inh_psib%nst == psib%nst)
329 if (
present(vmagnus))
then
332 assert(
size(vmagnus, 1) >= mesh%np)
333 assert(
size(vmagnus, 2) == hm%d%nspin)
334 assert(
size(vmagnus, 3) == 2)
336 write(
message(1),
'(a)')
'Magnus operators only implemented for electrons at the moment'
345 deltat_ = -
m_zi*deltat
346 if (
present(deltat2)) deltat2_ =
m_zi*deltat2
348 deltat_ = cmplx(deltat,
m_zero, real64)
349 if (
present(deltat2)) deltat2_ = cmplx(deltat2,
m_zero, real64)
352 if (.not. hm%is_hermitian() .and. te%exp_method ==
exp_chebyshev)
then
353 write(
message(1),
'(a)')
'The Chebyshev expansion cannot be used for non-Hermitian operators.'
354 write(
message(2),
'(a)')
'Please use the Lanczos exponentiation scheme ("TDExponentialMethod = lanczos")'
355 write(
message(3),
'(a)')
'or the Taylor expansion ("TDExponentialMethod = taylor") method.'
359 select case (te%exp_method)
362 if (
present(deltat2))
then
364 psib2, deltat2_, vmagnus)
369 if (
present(inh_psib))
then
370 if (
present(deltat2))
then
372 psib2, deltat2_, vmagnus, inh_psib)
375 vmagnus=vmagnus, inh_psib=inh_psib)
380 if (
present(psib2))
call psib%copy_data_to(mesh%np, psib2)
382 if (
present(inh_psib))
then
385 if (
present(psib2))
then
387 if (
present(inh_psib))
then
393 if (
present(inh_psib))
then
394 write(
message(1),
'(a)')
'Chebyshev exponential ("TDExponentialMethod = chebyshev")'
395 write(
message(2),
'(a)')
'with inhomogeneous term is not implemented'
402 chebyshev_function =>
chebyshev_exp_t(hm%spectral_half_span, hm%spectral_middle_point, deltat)
404 if (
present(psib2))
call psib%copy_data_to(mesh%np, psib2)
406 if (
present(psib2))
then
409 deallocate(chebyshev_function)
418 subroutine exponential_taylor_series_batch(te, namespace, mesh, hm, psib, deltat, psib2, deltat2, vmagnus, inh_psib, phik_shift)
421 class(
mesh_t),
intent(in) :: mesh
423 class(
batch_t),
intent(inout) :: psib
424 complex(real64),
intent(in) :: deltat
425 class(
batch_t),
optional,
intent(inout) :: psib2
426 complex(real64),
optional,
intent(in) :: deltat2
427 real(real64),
optional,
intent(in) :: vmagnus(:,:,:)
428 class(
batch_t),
optional,
intent(inout) :: inh_psib
429 integer,
optional,
intent(in) :: phik_shift
431 complex(real64) :: zfact, zfact2
432 integer :: iter, denom, phik_shift_
433 logical :: zfact_is_real
434 class(
batch_t),
allocatable :: psi1b, hpsi1b
439 call psib%clone_to(psi1b)
440 call psib%clone_to(hpsi1b)
444 zfact_is_real = abs(deltat-real(deltat, real64)) <
m_epsilon
446 if (
present(psib2))
call psib%copy_data_to(mesh%np, psib2)
448 if (
present(inh_psib))
then
450 call batch_axpy(mesh%np, real(zfact, real64), inh_psib, psib)
452 if (
present(psib2))
then
453 zfact2 = zfact2*deltat2
454 call batch_axpy(mesh%np, real(zfact2, real64), inh_psib, psib2)
461 do iter = 1, te%exp_order
462 denom = iter+phik_shift_
463 if (
present(inh_psib)) denom = denom + 1
464 zfact = zfact*(-
m_zi*deltat)/denom
465 if (
present(deltat2)) zfact2 = zfact2*(-
m_zi*deltat2)/denom
466 zfact_is_real = .not. zfact_is_real
473 call operate_batch(hm, namespace, mesh, psi1b, hpsi1b, vmagnus)
475 if (
present(inh_psib))
then
476 call operate_batch(hm, namespace, mesh, inh_psib, hpsi1b, vmagnus)
478 call operate_batch(hm, namespace, mesh, psib, hpsi1b, vmagnus)
482 if (zfact_is_real)
then
483 call batch_axpy(mesh%np, real(zfact, real64), hpsi1b, psib)
484 if (
present(psib2))
call batch_axpy(mesh%np, real(zfact2, real64), hpsi1b, psib2)
487 if (
present(psib2))
call batch_axpy(mesh%np, zfact2, hpsi1b, psib2)
490 if (iter /= te%exp_order)
call hpsi1b%copy_data_to(mesh%np, psi1b)
496 safe_deallocate_a(psi1b)
497 safe_deallocate_a(hpsi1b)
506 class(
mesh_t),
intent(in) :: mesh
508 class(
batch_t),
intent(inout) :: psib
509 complex(real64),
intent(in) :: deltat
510 real(real64),
optional,
intent(in) :: vmagnus(:,:,:)
511 class(
batch_t),
optional,
intent(in) :: inh_psib
513 class(
batch_t),
allocatable :: tmpb
517 if (
present(inh_psib))
then
518 call inh_psib%clone_to(tmpb, copy_data=.
true.)
521 call batch_axpy(mesh%np, real(deltat, real64), tmpb, psib)
523 safe_deallocate_a(tmpb)
546 class(
mesh_t),
intent(in) :: mesh
548 class(
batch_t),
intent(inout) :: psib
549 complex(real64),
intent(in) :: deltat
551 complex(real64) function fun(z)
553 complex(real64),
intent(in) :: z
556 real(real64),
optional,
intent(in) :: vmagnus(:,:,:)
558 integer :: iter, l, ii, ist, max_initialized
559 complex(real64),
allocatable :: hamilt(:,:,:), expo(:,:,:)
560 real(real64),
allocatable :: beta(:), res(:), norm(:)
561 integer,
parameter :: max_order = 200
567 if (te%exp_method /= exp_lanczos)
then
568 message(1) =
"The exponential method needs to be set to Lanzcos (TDExponentialMethod=lanczos)."
572 safe_allocate(beta(1:psib%nst))
573 safe_allocate(res(1:psib%nst))
574 safe_allocate(norm(1:psib%nst))
575 safe_allocate(vb(1:max_order))
576 call psib%clone_to(vb(1)%p)
579 call psib%copy_data_to(mesh%np, vb(1)%p, async=.
true.)
582 if (te%full_batch) beta = norm2(beta)
585 if (all(abs(beta) <= 1.0e-12_real64))
then
586 safe_deallocate_a(beta)
587 safe_deallocate_a(res)
588 safe_deallocate_a(norm)
590 safe_deallocate_a(vb)
598 safe_allocate(hamilt(1:max_order+1, 1:max_order+1, 1:psib%nst))
599 safe_allocate( expo(1:max_order+1, 1:max_order+1, 1:psib%nst))
602 do iter = 1, max_order-1
603 call psib%clone_to(vb(iter + 1)%p)
604 max_initialized = iter + 1
607 call operate_batch(hm, namespace, mesh, vb(iter)%p, vb(iter+1)%p, vmagnus)
610 if (hm%is_hermitian())
then
612 hamilt(1:max(l-1, 1), iter, 1:psib%nst) =
m_zero
616 hamilt(iter, 1:iter-2, 1:psib%nst) =
m_zero
622 normalize = .false., overlap = hamilt(l:iter, iter, 1:psib%nst), norm = hamilt(iter + 1, iter, 1:psib%nst), &
623 gs_scheme = te%arnoldi_gs, full_batch = te%full_batch)
630 res(ii) = abs(hamilt(iter + 1, iter, ii) * abs(expo(iter, 1, ii)))
634 if (all(abs(hamilt(iter + 1, iter, :)) < 1.0e4_real64*
m_epsilon))
exit
639 if (abs(hamilt(iter + 1, iter, ist)) >= 1.0e4_real64 *
m_epsilon)
then
640 norm(ist) =
m_one / abs(hamilt(iter + 1, iter, ist))
644 call batch_scal(mesh%np, norm, vb(iter+1)%p, a_full = .false.)
646 if (iter > 3 .and. all(res < te%lanczos_tol))
exit
650 if (iter == max_order)
then
651 write(
message(1),
'(a,i5,a,es9.2)')
'Lanczos exponential expansion did not converge after ', iter, &
652 ' iterations. Residual: ', maxval(res)
655 write(
message(1),
'(a,i5)')
'Debug: Lanczos exponential iterations: ', iter
661 call batch_scal(mesh%np, expo(1,1,1:psib%nst), psib, a_full = .false.)
664 call batch_axpy(mesh%np, beta(1:psib%nst)*expo(ii,1,1:psib%nst), vb(ii)%p, psib, a_full = .false.)
669 do ii = 1, max_initialized
673 safe_deallocate_a(vb)
674 safe_deallocate_a(hamilt)
675 safe_deallocate_a(expo)
676 safe_deallocate_a(beta)
677 safe_deallocate_a(res)
678 safe_deallocate_a(norm)
703 class(
mesh_t),
intent(in) :: mesh
705 class(
batch_t),
intent(inout) :: psib
706 real(real64),
intent(in) :: deltat
708 real(real64),
optional,
intent(in) :: vmagnus(:,:,:)
710 integer :: j, order_needed
711 complex(real64) :: coefficient
712 complex(real64),
allocatable :: coefficients(:)
713 real(real64) :: error
714 class(
batch_t),
allocatable,
target :: psi0, psi1, psi2
715 class(
batch_t),
pointer :: psi_n, psi_n1, psi_n2
716 integer,
parameter :: max_order = 200
721 call psib%clone_to(psi0)
722 call psib%clone_to(psi1)
723 call psib%clone_to(psi2)
724 call psib%copy_data_to(mesh%np, psi0)
726 order_needed = max_order
728 error = chebyshev_function%get_error(j)
729 if (error >
m_zero .and. error < te%chebyshev_tol)
then
735 call chebyshev_function%get_coefficients(j, coefficients)
738 call batch_scal(mesh%np, coefficients(0), psib)
742 call batch_axpy(mesh%np, -hm%spectral_middle_point, psi0, psi1)
745 call batch_axpy(mesh%np, coefficients(1), psi1, psib)
752 do j = 2, order_needed
754 call operate_batch(hm, namespace, mesh, psi_n1, psi_n, vmagnus)
755 call batch_axpy(mesh%np, -hm%spectral_middle_point, psi_n1, psi_n)
760 call batch_axpy(mesh%np, coefficients(j), psi_n, psib)
763 if (mod(j, 3) == 2)
then
767 else if (mod(j, 3) == 1)
then
778 if (order_needed == max_order)
then
779 write(
message(1),
'(a,i5,a,es9.2)')
'Chebyshev exponential expansion did not converge after ', j, &
780 ' iterations. Coefficient: ', coefficient
783 write(
message(1),
'(a,i5)')
'Debug: Chebyshev exponential iterations: ', j
790 safe_deallocate_a(psi0)
791 safe_deallocate_a(psi1)
792 safe_deallocate_a(psi2)
794 safe_deallocate_a(coefficients)
802 subroutine operate_batch(hm, namespace, mesh, psib, hpsib, vmagnus)
805 class(
mesh_t),
intent(in) :: mesh
806 class(
batch_t),
intent(inout) :: psib
807 class(
batch_t),
intent(inout) :: hpsib
808 real(real64),
optional,
intent(in) :: vmagnus(:, :, :)
812 if (
present(vmagnus))
then
813 call hm%zmagnus_apply(namespace, mesh, psib, hpsib, vmagnus)
815 call hm%zapply(namespace, mesh, psib, hpsib)
828 class(
mesh_t),
intent(inout) :: mesh
831 real(real64),
intent(in) :: deltat
832 integer,
optional,
intent(inout) :: order
835 real(real64) :: zfact
847 do i = 1, te%exp_order
848 zfact = zfact * deltat / i
856 do ik = st%d%kpt%start, st%d%kpt%end
857 do ib = st%group%block_start, st%group%block_end
859 call batch_axpy(mesh%np, -
m_zi, hst1%group%psib(ib, ik), st1%group%psib(ib, ik))
860 call batch_axpy(mesh%np, zfact, st1%group%psib(ib, ik), st%group%psib(ib, ik))
878 do ik = st%d%kpt%start, st%d%kpt%end
879 do ib = st%group%block_start, st%group%block_end
880 call batch_axpy(mesh%np, deltat, st1%group%psib(ib, ik), st%group%psib(ib, ik))
885 do i = 1, te%exp_order
886 zfact = zfact * deltat / (i+1)
894 do ik = st%d%kpt%start, st%d%kpt%end
895 do ib = st%group%block_start, st%group%block_end
897 call batch_axpy(mesh%np, -
m_zi, hst1%group%psib(ib, ik), st1%group%psib(ib, ik))
898 call batch_axpy(mesh%np, deltat * zfact, st1%group%psib(ib, ik), st%group%psib(ib, ik))
909 if (
present(order)) order = te%exp_order*st%nik*st%nst
917 class(
mesh_t),
intent(in) :: mesh
919 class(
batch_t),
intent(inout) :: psib
920 real(real64),
intent(in) :: deltat
921 integer,
intent(in) :: k
922 real(real64),
optional,
intent(in) :: vmagnus(:,:,:)
925 complex(real64) :: deltat_
930 if (
present(vmagnus))
then
933 assert(
size(vmagnus, 1) >= mesh%np)
934 assert(
size(vmagnus, 2) == hm%d%nspin)
935 assert(
size(vmagnus, 3) == 2)
937 write(
message(1),
'(a)')
'Magnus operators only implemented for electrons at the moment'
942 if (.not. hm%is_hermitian() .and. te%exp_method ==
exp_chebyshev)
then
943 write(
message(1),
'(a)')
'The Chebyshev expansion for the exponential will only converge if the imaginary'
944 write(
message(2),
'(a)')
'eigenvalues are small enough compared to the span of the real eigenvalues,'
945 write(
message(3),
'(a)')
'i.e., for ratios smaller than about 1e-3.'
946 write(
message(4),
'(a)')
'The Lanczos method ("TDExponentialMethod = lanczos") is guaranteed to'
947 write(
message(5),
'(a)')
'always converge in this case.'
951 deltat_ = cmplx(deltat,
m_zero, real64)
953 select case (te%exp_method)
960 else if (k == 2)
then
963 write(
message(1),
'(a)')
'Lanczos expansion not implemented for phi_k, k > 2'
970 else if (k == 2)
then
973 write(
message(1),
'(a)')
'Chebyshev expansion not implemented for phi_k, k > 2'
977 deallocate(chebyshev_function)
batchified version of the BLAS axpy routine:
scale a batch by a constant or vector
subroutine, public accel_finish()
This module implements batches of mesh functions.
This module implements common operations on batches of mesh functions.
subroutine, public batch_set_zero(this, np, async)
fill all mesh functions of the batch with zero
This module contains interfaces for BLAS routines You should not use these routines directly....
subroutine, public exponential_copy(teo, tei)
subroutine, public exponential_apply_all(te, namespace, mesh, hm, st, deltat, order)
Note that this routine not only computes the exponential, but also an extra term if there is a inhomo...
subroutine operate_batch(hm, namespace, mesh, psib, hpsib, vmagnus)
subroutine exponential_apply_phi_batch(te, namespace, mesh, hm, psib, deltat, k, vmagnus)
subroutine, public exponential_init(te, namespace, full_batch)
integer, parameter, public exp_taylor
subroutine exponential_lanczos_batch(te, namespace, mesh, hm, psib, deltat, vmagnus, inh_psib)
subroutine exponential_apply_single(te, namespace, mesh, hm, zpsi, ist, ik, deltat, vmagnus, imag_time)
Wrapper to batchified routine for applying exponential to an array.
subroutine exponential_taylor_series_batch(te, namespace, mesh, hm, psib, deltat, psib2, deltat2, vmagnus, inh_psib, phik_shift)
subroutine exponential_apply_batch(te, namespace, mesh, hm, psib, deltat, psib2, deltat2, vmagnus, imag_time, inh_psib)
This routine performs the operation:
subroutine, public exponential_lanczos_function_batch(te, namespace, mesh, hm, psib, deltat, fun, vmagnus)
Compute fun(H) psib, i.e. the application of a function of the Hamiltonian to a batch.
integer, parameter, public exp_chebyshev
subroutine exponential_cheby_batch(te, namespace, mesh, hm, psib, deltat, chebyshev_function, vmagnus)
Calculates the exponential of the Hamiltonian through an expansion in Chebyshev polynomials.
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
complex(real64), parameter, public m_zi
real(real64), parameter, public m_epsilon
complex(real64), parameter, public m_z1
real(real64), parameter, public m_one
This module defines an abstract class for Hamiltonians.
subroutine, public zhamiltonian_elec_apply_all(hm, namespace, mesh, st, hst)
pure logical function, public hamiltonian_elec_inh_term(hm)
subroutine, public zlalg_matrix_function(n, factor, a, fun_a, fun, hermitian)
This routine calculates a function of a matrix by using an eigenvalue decomposition.
This module is intended to contain "only mathematical" functions and procedures.
complex(real64) pure function, public phi2(z)
Compute phi2(z) = (phi1(z)-1)/z = (exp(z) - z - 1)/z^2.
complex(real64) pure function, public exponential(z)
Wrapper for exponential.
complex(real64) pure function, public phi1(z)
Compute phi1(z) = (exp(z)-1)/z.
This module defines functions over batches of mesh functions.
subroutine, public mesh_batch_nrm2(mesh, aa, nrm2, reduce)
Calculate the norms (norm2, not the square!) of a batch of mesh functions.
subroutine, public zmesh_batch_orthogonalization(mesh, nst, psib, phib, normalize, overlap, norm, gs_scheme, full_batch)
Orthonormalizes states of phib to the orbitals of nst batches of psi.
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
subroutine, public messages_warning(no_lines, all_nodes, namespace)
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_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
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.
subroutine, public states_elec_end(st)
finalize the states_elec_t object
subroutine, public states_elec_copy(stout, stin, exclude_wfns, exclude_eigenval, special)
make a (selective) copy of a states_elec_t object
type(type_t), parameter, public type_cmplx
Class defining batches of mesh functions.
The abstract Hamiltonian class defines a skeleton for specific implementations.
Describes mesh distribution to nodes.
The states_elec_t class contains all electronic wave functions.
batches of electronic states