Octopus
exponential_oct_m Module Reference

Data Types

type  exponential_t
 

Functions/Subroutines

subroutine, public exponential_init (te, namespace, full_batch)
 
subroutine, public exponential_copy (teo, tei)
 
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. More...
 
subroutine exponential_apply_batch (te, namespace, mesh, hm, psib, deltat, psib2, deltat2, vmagnus, imag_time, inh_psib)
 This routine performs the operation: More...
 
subroutine exponential_taylor_series_batch (te, namespace, mesh, hm, psib, deltat, psib2, deltat2, vmagnus, inh_psib, phik_shift)
 
subroutine exponential_lanczos_batch (te, namespace, mesh, hm, psib, deltat, vmagnus, inh_psib)
 
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. More...
 
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. More...
 
subroutine operate_batch (hm, namespace, mesh, psib, hpsib, vmagnus)
 
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 inhomogeneous term in the Hamiltonian hm. More...
 
subroutine exponential_apply_phi_batch (te, namespace, mesh, hm, psib, deltat, k, vmagnus)
 

Variables

integer, parameter, public exp_lanczos = 2
 
integer, parameter, public exp_taylor = 3
 
integer, parameter, public exp_chebyshev = 4
 

Function/Subroutine Documentation

◆ exponential_init()

subroutine, public exponential_oct_m::exponential_init ( type(exponential_t), intent(out)  te,
type(namespace_t), intent(in)  namespace,
logical, intent(in), optional  full_batch 
)

Definition at line 176 of file exponential.F90.

◆ exponential_copy()

subroutine, public exponential_oct_m::exponential_copy ( type(exponential_t), intent(inout)  teo,
type(exponential_t), intent(in)  tei 
)

Definition at line 314 of file exponential.F90.

◆ exponential_apply_single()

subroutine exponential_oct_m::exponential_apply_single ( class(exponential_t), intent(inout)  te,
type(namespace_t), intent(in)  namespace,
class(mesh_t), intent(in)  mesh,
type(hamiltonian_elec_t), intent(inout)  hm,
complex(real64), dimension(:, :), intent(inout), contiguous  zpsi,
integer, intent(in)  ist,
integer, intent(in)  ik,
real(real64), intent(in)  deltat,
real(real64), dimension(mesh%np, hm%d%nspin, 2), intent(in), optional  vmagnus,
logical, intent(in), optional  imag_time 
)
private

Wrapper to batchified routine for applying exponential to an array.

Definition at line 330 of file exponential.F90.

◆ exponential_apply_batch()

subroutine exponential_oct_m::exponential_apply_batch ( class(exponential_t), intent(inout)  te,
type(namespace_t), intent(in)  namespace,
class(mesh_t), intent(in)  mesh,
class(hamiltonian_abst_t), intent(inout)  hm,
class(batch_t), intent(inout)  psib,
real(real64), intent(in)  deltat,
class(batch_t), intent(inout), optional  psib2,
real(real64), intent(in), optional  deltat2,
real(real64), dimension(:,:,:), intent(in), optional  vmagnus,
logical, intent(in), optional  imag_time,
class(batch_t), intent(inout), optional  inh_psib 
)
private

This routine performs the operation:

\[ \exp{-i*\Delta t*hm(t)}|\psi> <-- |\psi> \]

If imag_time is present and is set to true, it performs instead:

\[ \exp{ \Delta t*hm(t)}|\psi> <-- |\psi> \]

If an inhomogeneous term is present, the operation is:

\[ \exp{-i*\Delta t*hm(t)}|\psi> + \Delta t*\phi_1{-i*\Delta t*hm(t)}|inh\psi> <-- |\psi> \]

where:

\[ \phi_1(x) = (e^x - 1)/x \]

Parameters
[in,out]inh_psibinhomogeneous term

Definition at line 395 of file exponential.F90.

◆ exponential_taylor_series_batch()

subroutine exponential_oct_m::exponential_taylor_series_batch ( type(exponential_t), intent(inout)  te,
type(namespace_t), intent(in)  namespace,
class(mesh_t), intent(in)  mesh,
class(hamiltonian_abst_t), intent(inout)  hm,
class(batch_t), intent(inout)  psib,
complex(real64), intent(in)  deltat,
class(batch_t), intent(inout), optional  psib2,
complex(real64), intent(in), optional  deltat2,
real(real64), dimension(:,:,:), intent(in), optional  vmagnus,
class(batch_t), intent(inout), optional  inh_psib,
integer, intent(in), optional  phik_shift 
)
private
Parameters
[in,out]inh_psibinhomogeneous term
[in]phik_shiftshift in the Taylor expansion coefficients for phi_k

Definition at line 511 of file exponential.F90.

◆ exponential_lanczos_batch()

subroutine exponential_oct_m::exponential_lanczos_batch ( type(exponential_t), intent(inout)  te,
type(namespace_t), intent(in)  namespace,
class(mesh_t), intent(in)  mesh,
class(hamiltonian_abst_t), intent(inout)  hm,
class(batch_t), intent(inout)  psib,
complex(real64), intent(in)  deltat,
real(real64), dimension(:,:,:), intent(in), optional  vmagnus,
class(batch_t), intent(in), optional  inh_psib 
)
private
Parameters
[in]inh_psibinhomogeneous term

Definition at line 596 of file exponential.F90.

◆ exponential_lanczos_function_batch()

subroutine, public exponential_oct_m::exponential_lanczos_function_batch ( type(exponential_t), intent(inout)  te,
type(namespace_t), intent(in)  namespace,
class(mesh_t), intent(in)  mesh,
class(hamiltonian_abst_t), intent(inout)  hm,
class(batch_t), intent(inout)  psib,
complex(real64), intent(in)  deltat,
  fun,
real(real64), dimension(:,:,:), intent(in), optional  vmagnus 
)

Compute fun(H) psib, i.e. the application of a function of the Hamiltonian to a batch.

Usually, the function is an exponential, or a related function.

Some details of the implementation can be understood from Saad, Y. (1992). Analysis of some Krylov subspace approximations to the matrix exponential operator. SIAM Journal on Numerical Analysis, 29(1), 209-228.

A pdf can be accessed [here](https: or [via the DOI](https: Equation numbers below refer to this paper.

Definition at line 635 of file exponential.F90.

◆ exponential_cheby_batch()

subroutine exponential_oct_m::exponential_cheby_batch ( type(exponential_t), intent(inout)  te,
type(namespace_t), intent(in)  namespace,
class(mesh_t), intent(in)  mesh,
class(hamiltonian_abst_t), intent(inout)  hm,
class(batch_t), intent(inout)  psib,
real(real64), intent(in)  deltat,
class(chebyshev_function_t), intent(in)  chebyshev_function,
real(real64), dimension(:,:,:), intent(in), optional  vmagnus 
)
private

Calculates the exponential of the Hamiltonian through an expansion in Chebyshev polynomials.

It uses the algorithm as outlined in [Hochbruck, M. & Ostermann, A.: Exponential Runge–Kutta methods for parabolic problems. Applied Numerical Mathematics 53, 323–339 (2005)](https: Section 4.1. See also [Kosloff, J. Phys. Chem. (1988), 92, 2087-2100](http: Section III.1 and [Kosloff, Annu. Rev. Phys. Chem. (1994), 45, 145–178](https: equations 7.1ff and especially figure 2 for the convergence properties of the coefficients.

Definition at line 792 of file exponential.F90.

◆ operate_batch()

subroutine exponential_oct_m::operate_batch ( class(hamiltonian_abst_t), intent(inout)  hm,
type(namespace_t), intent(in)  namespace,
class(mesh_t), intent(in)  mesh,
class(batch_t), intent(inout)  psib,
class(batch_t), intent(inout)  hpsib,
real(real64), dimension(:, :, :), intent(in), optional  vmagnus 
)
private

Definition at line 896 of file exponential.F90.

◆ exponential_apply_all()

subroutine, public exponential_oct_m::exponential_apply_all ( type(exponential_t), intent(inout)  te,
type(namespace_t), intent(in)  namespace,
class(mesh_t), intent(inout)  mesh,
type(hamiltonian_elec_t), intent(inout)  hm,
type(states_elec_t), intent(inout)  st,
real(real64), intent(in)  deltat,
integer, intent(inout), optional  order 
)

Note that this routine not only computes the exponential, but also an extra term if there is a inhomogeneous term in the Hamiltonian hm.

Definition at line 919 of file exponential.F90.

◆ exponential_apply_phi_batch()

subroutine exponential_oct_m::exponential_apply_phi_batch ( class(exponential_t), intent(inout)  te,
type(namespace_t), intent(in)  namespace,
class(mesh_t), intent(in)  mesh,
class(hamiltonian_abst_t), intent(inout)  hm,
class(batch_t), intent(inout)  psib,
real(real64), intent(in)  deltat,
integer, intent(in)  k,
real(real64), dimension(:,:,:), intent(in), optional  vmagnus 
)
private

Definition at line 1008 of file exponential.F90.

Variable Documentation

◆ exp_lanczos

integer, parameter, public exponential_oct_m::exp_lanczos = 2

Definition at line 154 of file exponential.F90.

◆ exp_taylor

integer, parameter, public exponential_oct_m::exp_taylor = 3

Definition at line 154 of file exponential.F90.

◆ exp_chebyshev

integer, parameter, public exponential_oct_m::exp_chebyshev = 4

Definition at line 154 of file exponential.F90.