pwdft  0.1
PW-DFT code for education
k_point Module Reference

Functions/Subroutines

subroutine ksum_rho ()
 
subroutine ksum_dos ()
 DOS. More...
 

Variables

integer, save nk
 
integer, dimension(3), save kgrd
 
integer, dimension(:,:), allocatable, save tetra
 
real(8), dimension(:,:), allocatable, save kvec
 

Function/Subroutine Documentation

◆ ksum_rho()

subroutine k_point::ksum_rho

Definition at line 35 of file k_point.F90.

35  !
36  use atm_spec, only : bvec, nelec
37  use kohn_sham, only : nbnd, ef, eval, evec
38  use fftw_wrapper, only : fft_g2r_w
39  use rho_v, only : rho
40  use gvec, only : g_rh, g_wf
41  use libtetrabz, only : libtetrabz_fermieng
42  !
43  integer :: ik, ibnd
44  real(8) :: occ(nbnd,nk)
45  complex(8) :: psir(g_rh%nr)
46  !
47  call libtetrabz_fermieng(2,bvec,nbnd,kgrd,eval,kgrd,occ,ef,nelec*0.5d0)
48  !
49  rho(1:g_rh%nr) = 0.0d0
50  do ik = 1, nk
51  do ibnd = 1, nbnd
52  call fft_g2r_w(evec(1:g_wf%npw,ibnd,ik), psir)
53  rho(1:g_rh%nr) = rho(1:g_rh%nr) &
54  & + dble(conjg(psir(1:g_rh%nr))*psir(1:g_rh%nr))*occ(ibnd,ik)*2.0d0
55  end do
56  end do
57  !

References atm_spec::bvec, kohn_sham::ef, kohn_sham::eval, kohn_sham::evec, fftw_wrapper::fft_g2r_w(), gvec::g_rh, gvec::g_wf, kgrd, kohn_sham::nbnd, atm_spec::nelec, and rho_v::rho.

Referenced by scf::kohn_sham_eq().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ksum_dos()

subroutine k_point::ksum_dos

DOS.

Definition at line 63 of file k_point.F90.

63  !
64  use constant, only : htr2ev
65  use atm_spec, only : bvec
66  use kohn_sham, only : nbnd, ef, eval
67  use libtetrabz, only : libtetrabz_dos, libtetrabz_intdos
68  !
69  integer,parameter :: ne = 500
70  integer :: ie, fo = 20
71  real(8) :: e0(ne), dos(ne), intdos(ne), wght(ne,nbnd,nk), &
72  & emax, emin, de
73  !
74  emax = maxval(eval(1:nbnd,1:nk))
75  emin = minval(eval(1:nbnd,1:nk))
76  de = (emax - emin) / dble(ne)
77  do ie = 1, ne
78  e0(ie) = emin + de * (ie -1)
79  end do
80  !
81  call libtetrabz_dos(2,bvec,nbnd,kgrd,eval,kgrd,wght,ne,e0)
82  !
83  do ie = 1, ne
84  dos(ie) = sum(wght(ie,1:nbnd,1:nk))
85  end do
86  !
87  call libtetrabz_intdos(2,bvec,nbnd,kgrd,eval,kgrd,wght,ne,e0)
88  !
89  do ie = 1, ne
90  intdos(ie) = sum(wght(ie,1:nbnd,1:nk))
91  end do
92  !
93  write(*,*) " Output dos.dat"
94  open(fo, file = "dos.dat")
95  !
96  write(fo,*) "# Energy[eV] DOS[eV^-1] IntDOS"
97  !
98  do ie = 1, ne
99  write(fo,*) (e0(ie)-ef)*htr2ev, dos(ie)*2.0d0/htr2ev, intdos(ie)*2.0d0
100  end do
101  !
102  close(fo)
103  !

References atm_spec::bvec, kohn_sham::ef, kohn_sham::eval, constant::htr2ev, kgrd, and kohn_sham::nbnd.

Referenced by pwdft().

Here is the caller graph for this function:

Variable Documentation

◆ nk

integer, save k_point::nk

Definition at line 27 of file k_point.F90.

27  integer,save :: nk

Referenced by plot::band_plot(), plot::fermi_plot(), energy::kinetic(), scf::kohn_sham_eq(), pwdft(), and stdin::read_stdin().

◆ kgrd

integer, dimension(3), save k_point::kgrd

Definition at line 28 of file k_point.F90.

28  integer,save :: kgrd(3)

Referenced by plot::fermi_plot(), energy::kinetic(), ksum_dos(), ksum_rho(), and stdin::read_stdin().

◆ tetra

integer, dimension(:,:), allocatable, save k_point::tetra

Definition at line 29 of file k_point.F90.

29  integer,allocatable,save :: tetra(:,:)

◆ kvec

real(8), dimension(:,:), allocatable, save k_point::kvec

Definition at line 30 of file k_point.F90.

30  real(8),allocatable,save :: kvec(:,:)

Referenced by plot::band_plot(), plot::fermi_plot(), energy::kinetic(), scf::kohn_sham_eq(), and stdin::read_stdin().

constant
Definition: constant.F90:23
libtetrabz
Definition: libtetrabz.F90:23
rho_v::rho
real(8), dimension(:), allocatable rho
(g_rhnr) Charge density
Definition: rho_v.F90:27
kohn_sham::ef
real(8), save ef
Fermi energy [Htr].
Definition: kohn_sham.F90:31
rho_v
Definition: rho_v.F90:23
fftw_wrapper
Definition: fftw_wrapper.F90:23
gvec
Definition: gvec.F90:23
k_point::kgrd
integer, dimension(3), save kgrd
Definition: k_point.F90:28
atm_spec::bvec
real(8), dimension(3, 3), save bvec
Unit reciplocal lattice vector [Bohr^-1].
Definition: atm_spec.F90:50
atm_spec::nelec
real(8), save nelec
Number of electrons per u.c.
Definition: atm_spec.F90:50
fftw_wrapper::fft_g2r_w
subroutine fft_g2r_w(wfG, wfR)
Uncompress -> W(G) e^{iGr} -> w(r)
Definition: fftw_wrapper.F90:121
atm_spec
Definition: atm_spec.F90:23
kohn_sham
Definition: kohn_sham.F90:23