Octopus
kdotp_calc_oct_m Module Reference

Functions/Subroutines

character(len=100) function, public kdotp_wfs_tag (dir, dir2)
 
subroutine, public calc_band_velocity (namespace, space, gr, st, hm, pert, velocity)
 Computes the k-point and band-resolved velocity. More...
 
subroutine, public dcalc_eff_mass_inv (namespace, space, gr, st, hm, lr, pert, eff_mass_inv, degen_thres)
 Computes the effective mass tensor. More...
 
subroutine, public dkdotp_add_occ (namespace, space, gr, st, hm, pert, kdotp_lr, degen_thres)
 add projection onto occupied states, by sum over states More...
 
subroutine, public zcalc_eff_mass_inv (namespace, space, gr, st, hm, lr, pert, eff_mass_inv, degen_thres)
 Computes the effective mass tensor. More...
 
subroutine, public zkdotp_add_occ (namespace, space, gr, st, hm, pert, kdotp_lr, degen_thres)
 add projection onto occupied states, by sum over states More...
 

Function/Subroutine Documentation

◆ kdotp_wfs_tag()

character(len=100) function, public kdotp_calc_oct_m::kdotp_wfs_tag ( integer, intent(in)  dir,
integer, intent(in), optional  dir2 
)

Definition at line 151 of file kdotp_calc.F90.

◆ calc_band_velocity()

subroutine, public kdotp_calc_oct_m::calc_band_velocity ( type(namespace_t), intent(in)  namespace,
class(space_t), intent(in)  space,
type(grid_t), intent(in)  gr,
type(states_elec_t), intent(in)  st,
type(hamiltonian_elec_t), intent(inout)  hm,
class(perturbation_t), intent(inout)  pert,
real(real64), dimension(:,:,:), intent(out)  velocity 
)

Computes the k-point and band-resolved velocity.

This is done using the following formula

\[ v = (dE_nk/dk)/hbar = -Im < u_nk | -i grad | u_nk > \]

This is identically zero for real wavefunctions.

Definition at line 173 of file kdotp_calc.F90.

◆ dcalc_eff_mass_inv()

subroutine, public kdotp_calc_oct_m::dcalc_eff_mass_inv ( type(namespace_t), intent(in)  namespace,
class(space_t), intent(in)  space,
type(grid_t), intent(in)  gr,
type(states_elec_t), intent(in)  st,
type(hamiltonian_elec_t), intent(inout)  hm,
type(lr_t), dimension(:,:), intent(in)  lr,
class(perturbation_t), intent(inout)  pert,
real(real64), dimension(:,:,:,:), intent(out)  eff_mass_inv,
real(real64), intent(in)  degen_thres 
)

Computes the effective mass tensor.

This is obtained using the formula

\[ m^-1[ij] = <psi0|H2ij|psi0> + 2*Re<psi0|H'i|psi'j> \]

for each state, spin, and k-point The off-diagonal elements are not correct in a degenerate subspace

Parameters
[in]lr(1, pdim)
[out]eff_mass_inv(pdim, pdim, nik, nst)

Definition at line 301 of file kdotp_calc.F90.

◆ dkdotp_add_occ()

subroutine, public kdotp_calc_oct_m::dkdotp_add_occ ( type(namespace_t), intent(in)  namespace,
class(space_t), intent(in)  space,
type(grid_t), intent(in)  gr,
type(states_elec_t), intent(in)  st,
type(hamiltonian_elec_t), intent(inout)  hm,
class(perturbation_t), intent(in)  pert,
type(lr_t), intent(inout)  kdotp_lr,
real(real64), intent(in)  degen_thres 
)

add projection onto occupied states, by sum over states

Definition at line 399 of file kdotp_calc.F90.

◆ zcalc_eff_mass_inv()

subroutine, public kdotp_calc_oct_m::zcalc_eff_mass_inv ( type(namespace_t), intent(in)  namespace,
class(space_t), intent(in)  space,
type(grid_t), intent(in)  gr,
type(states_elec_t), intent(in)  st,
type(hamiltonian_elec_t), intent(inout)  hm,
type(lr_t), dimension(:,:), intent(in)  lr,
class(perturbation_t), intent(inout)  pert,
real(real64), dimension(:,:,:,:), intent(out)  eff_mass_inv,
real(real64), intent(in)  degen_thres 
)

Computes the effective mass tensor.

This is obtained using the formula

\[ m^-1[ij] = <psi0|H2ij|psi0> + 2*Re<psi0|H'i|psi'j> \]

for each state, spin, and k-point The off-diagonal elements are not correct in a degenerate subspace

Parameters
[in]lr(1, pdim)
[out]eff_mass_inv(pdim, pdim, nik, nst)

Definition at line 540 of file kdotp_calc.F90.

◆ zkdotp_add_occ()

subroutine, public kdotp_calc_oct_m::zkdotp_add_occ ( type(namespace_t), intent(in)  namespace,
class(space_t), intent(in)  space,
type(grid_t), intent(in)  gr,
type(states_elec_t), intent(in)  st,
type(hamiltonian_elec_t), intent(inout)  hm,
class(perturbation_t), intent(in)  pert,
type(lr_t), intent(inout)  kdotp_lr,
real(real64), intent(in)  degen_thres 
)

add projection onto occupied states, by sum over states

Definition at line 638 of file kdotp_calc.F90.