w90_wan_ham Module

This module contain operations on the Hamiltonian in the WF basis


Uses

  • module~~w90_wan_ham~~UsesGraph module~w90_wan_ham w90_wan_ham module~w90_constants w90_constants module~w90_wan_ham->module~w90_constants

Used by

  • module~~w90_wan_ham~~UsedByGraph module~w90_wan_ham w90_wan_ham proc~calctdfanddos calcTDFandDOS proc~calctdfanddos->module~w90_wan_ham proc~berry_get_sc_klist berry_get_sc_klist proc~berry_get_sc_klist->module~w90_wan_ham proc~berry_get_shc_klist berry_get_shc_klist proc~berry_get_shc_klist->module~w90_wan_ham proc~berry_get_imfgh_klist berry_get_imfgh_klist proc~berry_get_imfgh_klist->module~w90_wan_ham proc~k_slice k_slice proc~k_slice->module~w90_wan_ham proc~dos_main dos_main proc~dos_main->module~w90_wan_ham proc~gyrotropic_get_k_list gyrotropic_get_k_list proc~gyrotropic_get_k_list->module~w90_wan_ham module~w90_geninterp w90_geninterp module~w90_geninterp->module~w90_wan_ham proc~berry_get_kubo_k berry_get_kubo_k proc~berry_get_kubo_k->module~w90_wan_ham program~postw90 postw90 program~postw90->module~w90_geninterp

Contents


Subroutines

private subroutine wham_get_D_h_a(delHH_a, UU, eig, ef, D_h_a)

Compute D^H_a=UU^dag.del_a UU (a=alpha,beta), using Eq.(24) of WYSV06

Arguments

Type IntentOptional AttributesName
complex(kind=dp), intent(in), dimension(:, :):: delHH_a
complex(kind=dp), intent(in), dimension(:, :):: UU
real(kind=dp), intent(in), dimension(:):: eig
real(kind=dp), intent(in) :: ef
complex(kind=dp), intent(out), dimension(:, :):: D_h_a

public subroutine wham_get_D_h(delHH, UU, eig, D_h)

Compute D^H_a=UU^dag.del_a UU (a=x,y,z) using Eq.(24) of WYSV06

Arguments

Type IntentOptional AttributesName
complex(kind=dp), intent(in), dimension(:, :, :):: delHH
complex(kind=dp), intent(in), dimension(:, :):: UU
real(kind=dp), intent(in), dimension(:):: eig
complex(kind=dp), intent(out), dimension(:, :, :):: D_h

public subroutine wham_get_D_h_P_value(delHH, UU, eig, D_h)

Compute D^H_a=UU^dag.del_a UU (a=x,y,z) using Eq.(24) of WYSV06

Arguments

Type IntentOptional AttributesName
complex(kind=dp), intent(in), dimension(:, :, :):: delHH
complex(kind=dp), intent(in), dimension(:, :):: UU
real(kind=dp), intent(in), dimension(:):: eig
complex(kind=dp), intent(out), dimension(:, :, :):: D_h

private subroutine wham_get_JJp_JJm_list(delHH, UU, eig, JJp_list, JJm_list, occ)

Arguments

Type IntentOptional AttributesName
complex(kind=dp), intent(inout), dimension(:, :):: delHH
complex(kind=dp), intent(in), dimension(:, :):: UU
real(kind=dp), intent(in), dimension(:):: eig
complex(kind=dp), intent(out), dimension(:, :, :):: JJp_list
complex(kind=dp), intent(out), dimension(:, :, :):: JJm_list
real(kind=dp), intent(in), optional dimension(:):: occ

public subroutine wham_get_occ_mat_list(UU, f_list, g_list, eig, occ)

Occupation matrix f, and g=1-f for a list of Fermi energies

Arguments

Type IntentOptional AttributesName
complex(kind=dp), intent(in), dimension(:, :):: UU
complex(kind=dp), intent(out), dimension(:, :, :):: f_list
complex(kind=dp), intent(out), dimension(:, :, :):: g_list
real(kind=dp), intent(in), optional dimension(:):: eig
real(kind=dp), intent(in), optional dimension(:):: occ

private subroutine wham_get_deleig_a(deleig_a, eig, delHH_a, UU)

Band derivatives dE/dk_a

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(out) :: deleig_a(num_wann)
real(kind=dp), intent(in) :: eig(num_wann)
complex(kind=dp), intent(in), dimension(:, :):: delHH_a
complex(kind=dp), intent(in), dimension(:, :):: UU

public subroutine wham_get_eig_deleig(kpt, eig, del_eig, HH, delHH, UU)

Given a k point, this function returns eigenvalues E and derivatives of the eigenvalues dE/dk_a, using wham_get_deleig_a

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(in), dimension(3):: kpt

the three coordinates of the k point vector (in relative coordinates)

real(kind=dp), intent(out) :: eig(num_wann)

the calculated eigenvalues at kpt

real(kind=dp), intent(out) :: del_eig(num_wann,3)

the calculated derivatives of the eigenvalues at kpt [first component: band; second component: 1,2,3 for the derivatives along the three k directions]

complex(kind=dp), intent(out), dimension(:, :):: HH

the Hamiltonian matrix at kpt

complex(kind=dp), intent(out), dimension(:, :, :):: delHH

the delHH matrix (derivative of H) at kpt

complex(kind=dp), intent(out), dimension(:, :):: UU

the rotation matrix that gives the eigenvectors of HH

public subroutine wham_get_eig_deleig_TB_conv(kpt, eig, del_eig, delHH, UU)

Given a k point, this function returns eigenvalues E and derivatives of the eigenvalues dE/dk_a, using wham_get_deleig_a

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(in), dimension(3):: kpt

the three coordinates of the k point vector (in relative coordinates)

real(kind=dp), intent(in) :: eig(num_wann)
real(kind=dp), intent(out) :: del_eig(num_wann,3)
complex(kind=dp), intent(in), dimension(:, :, :):: delHH

the delHH matrix (derivative of H) at kpt

complex(kind=dp), intent(in), dimension(:, :):: UU

the rotation matrix that gives the eigenvectors of HH

public subroutine wham_get_eig_UU_HH_JJlist(kpt, eig, UU, HH, JJp_list, JJm_list, occ)

Wrapper routine used to reduce number of Fourier calls

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(in), dimension(3):: kpt
real(kind=dp), intent(out) :: eig(num_wann)
complex(kind=dp), intent(out), dimension(:, :):: UU
complex(kind=dp), intent(out), dimension(:, :):: HH
complex(kind=dp), intent(out), dimension(:, :, :, :):: JJp_list
complex(kind=dp), intent(out), dimension(:, :, :, :):: JJm_list
real(kind=dp), intent(in), optional dimension(:):: occ

public subroutine wham_get_eig_UU_HH_AA_sc_TB_conv(kpt, eig, UU, HH, HH_da, HH_dadb)

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(in), dimension(3):: kpt
real(kind=dp), intent(out) :: eig(num_wann)
complex(kind=dp), intent(out), dimension(:, :):: UU
complex(kind=dp), intent(out), dimension(:, :):: HH
complex(kind=dp), intent(out), dimension(:, :, :):: HH_da
complex(kind=dp), intent(out), dimension(:, :, :, :):: HH_dadb

public subroutine wham_get_eig_UU_HH_AA_sc(kpt, eig, UU, HH, HH_da, HH_dadb)

Wrapper routine used to reduce number of Fourier calls

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(in), dimension(3):: kpt
real(kind=dp), intent(out) :: eig(num_wann)
complex(kind=dp), intent(out), dimension(:, :):: UU
complex(kind=dp), intent(out), dimension(:, :):: HH
complex(kind=dp), intent(out), dimension(:, :, :):: HH_da
complex(kind=dp), intent(out), dimension(:, :, :, :):: HH_dadb