pwdft  0.1
PW-DFT code for education
plot Module Reference

Functions/Subroutines

subroutine band_plot ()
 Output band dispersion. More...
 
subroutine fermi_plot ()
 Output FermiSurfer file. More...
 

Function/Subroutine Documentation

◆ band_plot()

subroutine plot::band_plot

Output band dispersion.

Definition at line 32 of file plot.F90.

32  !
33  use constant, only : htr2ev
34  use k_point, only : kvec, nk
35  use kohn_sham, only : nbnd, eval, ef
36  !
37  integer :: fo = 20, ibnd, ik
38  real(8) :: xk(nk), dkvec(3), dxk
39  !
40  write(*,*) " Output band.dat"
41  open(fo, file="band.dat")
42  !
43  xk(1) = 0.0d0
44  do ik = 2, nk
45  dkvec(1:3) = kvec(1:3,ik) - kvec(1:3,ik-1)
46  dxk = sqrt(dot_product(dkvec(1:3), dkvec(1:3)))
47  xk(ik) = xk(ik-1) + dxk
48  end do
49  !
50  write(fo,*) "# k energy[eV]"
51  do ibnd = 1, nbnd
52  do ik = 1, nk
53  write(fo,*) xk(ik), (eval(ibnd, ik)-ef)*htr2ev
54  end do
55  write(fo,*) ""
56  end do
57  !
58  close(fo)
59  !

References kohn_sham::ef, kohn_sham::eval, constant::htr2ev, k_point::kvec, kohn_sham::nbnd, and k_point::nk.

Referenced by pwdft().

Here is the caller graph for this function:

◆ fermi_plot()

subroutine plot::fermi_plot

Output FermiSurfer file.

Definition at line 65 of file plot.F90.

65  !
66  use atm_spec, only : bvec
67  use k_point, only : kgrd, nk, kvec
68  use kohn_sham, only : eval, nbnd, evec, ef
69  use gvec, only : g_wf
70  use constant, only : htr2ev
71  !
72  integer :: ik, ibnd, fo = 20, ii, ipw
73  real(8) :: kgv(g_wf%npw,3), avf(nbnd,nk), vf(3)
74  !
75  do ik = 1, nk
76  do ipw = 1, g_wf%npw
77  kgv(ipw,1:3) = kvec(1:3,ik) + matmul(bvec(1:3,1:3), dble(g_wf%mill(1:3,g_wf%map(ipw))))
78  end do
79  do ibnd = 1, nbnd
80  do ii = 1, 3
81  vf(ii) = dble(sum(conjg(evec(1:g_wf%npw,ibnd,ik)) &
82  & * evec(1:g_wf%npw,ibnd,ik) &
83  & * kgv(1:g_wf%npw,ii ) ))
84  end do
85  avf(ibnd,ik) = sqrt(dot_product(vf(1:3),vf(1:3)))
86  end do
87  end do
88  !
89  write(*,*) " Output vf.frmsf"
90  open(fo, file = "vf.frmsf")
91  write(fo,*) kgrd(3), kgrd(2), kgrd(1)
92  write(fo,*) 1
93  write(fo,*) nbnd
94  write(fo,*) real(bvec(1:3,3))
95  write(fo,*) real(bvec(1:3,2))
96  write(fo,*) real(bvec(1:3,1))
97  do ibnd = 1, nbnd
98  do ik = 1, nk
99  write(fo,*) real((eval(ibnd,ik) - ef)*htr2ev)
100  end do
101  end do
102  do ibnd = 1, nbnd
103  do ik = 1, nk
104  write(fo,*) real(avf(ibnd,ik))
105  end do
106  end do
107  close(fo)
108  !

References atm_spec::bvec, kohn_sham::ef, kohn_sham::eval, kohn_sham::evec, gvec::g_wf, constant::htr2ev, k_point::kgrd, k_point::kvec, kohn_sham::nbnd, and k_point::nk.

Referenced by pwdft().

Here is the caller graph for this function:
constant
Definition: constant.F90:23
k_point
Definition: k_point.F90:23
gvec
Definition: gvec.F90:23
atm_spec
Definition: atm_spec.F90:23
kohn_sham
Definition: kohn_sham.F90:23