pwdft  0.1
PW-DFT code for education
k_point.F90
Go to the documentation of this file.
1 !
2 ! Copyright (c) 2018 Mitsuaki Kawamura
3 !
4 ! Permission is hereby granted, free of charge, to any person obtaining a
5 ! copy of this software and associated documentation files (the
6 ! "Software"), to deal in the Software without restriction, including
7 ! without limitation the rights to use, copy, modify, merge, publish,
8 ! distribute, sublicense, and/or sell copies of the Software, and to
9 ! permit persons to whom the Software is furnished to do so, subject to
10 ! the following conditions:
11 !
12 ! The above copyright notice and this permission notice shall be included
13 ! in all copies or substantial portions of the Software.
14 !
15 ! THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
16 ! OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF
17 ! MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT.
18 ! IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY
19 ! CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT,
20 ! TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE
21 ! SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
22 !
23 module k_point
24  !
25  implicit none
26  !
27  integer,save :: nk
28  integer,save :: kgrd(3)
29  integer,allocatable,save :: tetra(:,:)
30  real(8),allocatable,save :: kvec(:,:)
31  !
32 contains
33  !
34  subroutine ksum_rho()
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  !
58  end subroutine ksum_rho
62  subroutine ksum_dos()
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  !
104  end subroutine ksum_dos
105  !
106 end module k_point
k_point::tetra
integer, dimension(:,:), allocatable, save tetra
Definition: k_point.F90:29
constant
Definition: constant.F90:23
libtetrabz
Definition: libtetrabz.F90:23
gvec::g_wf
type(fft), save g_wf
Definition: gvec.F90:38
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
kohn_sham::evec
complex(8), dimension(:,:,:), allocatable, save evec
(g_wfnpw,nbnd,nk) Kohn-Sham eigenvector (orbital)
Definition: kohn_sham.F90:35
fftw_wrapper
Definition: fftw_wrapper.F90:23
k_point
Definition: k_point.F90:23
k_point::ksum_rho
subroutine ksum_rho()
Definition: k_point.F90:35
gvec
Definition: gvec.F90:23
kohn_sham::nbnd
integer, save nbnd
Number of bands.
Definition: kohn_sham.F90:29
k_point::nk
integer, save nk
Definition: k_point.F90:27
k_point::kgrd
integer, dimension(3), save kgrd
Definition: k_point.F90:28
kohn_sham::eval
real(8), dimension(:,:), allocatable, save eval
(nbnd,nk) Kohn-Sham eigenvalue (energy)
Definition: kohn_sham.F90:33
gvec::g_rh
type(fft), save g_rh
Definition: gvec.F90:37
constant::htr2ev
real(8), parameter htr2ev
Definition: constant.F90:27
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
k_point::kvec
real(8), dimension(:,:), allocatable, save kvec
Definition: k_point.F90:30
atm_spec
Definition: atm_spec.F90:23
k_point::ksum_dos
subroutine ksum_dos()
DOS.
Definition: k_point.F90:63
kohn_sham
Definition: kohn_sham.F90:23