The Elk Code
Loading...
Searching...
No Matches
writew90eig.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2015 Manh Duc Le, 2017-18 Arsenii Gerasimov, Yaroslav Kvashnin
3
! and Lars Nordstrom. This file is distributed under the terms of the GNU
4
! General Public License. See the file COPYING for license details.
5
6
subroutine
writew90eig
7
use
modmain
8
use
modw90
9
implicit none
10
! local variables
11
integer
ik,jk,ist,i
12
real
(8) t1
13
character(256)
fname
14
fname=trim(
seedname
)//
'.eig'
15
open
(50,file=trim(fname),action=
'WRITE'
,form=
'FORMATTED'
)
16
! loop over non-reduced k-points
17
do
ik=1,
nkptnr
18
! equivalent reduced k-point
19
jk=
ivkik
(
ivk
(1,ik),
ivk
(2,ik),
ivk
(3,ik))
20
do
i=1,
num_bands
21
ist=
idxw90
(i)
22
t1=
evalsv
(ist,jk)-
efermi
23
write
(50,
'(2I6,G18.10)'
) i,ik,t1*
ha_ev
24
end do
25
end do
26
close
(50)
27
write
(*,*)
28
write
(*,
'("Info(writew90eig): created file ",A)'
) trim(fname)
29
end subroutine
30
modmain
Definition
modmain.f90:6
modmain::efermi
real(8) efermi
Definition
modmain.f90:904
modmain::nkptnr
integer nkptnr
Definition
modmain.f90:463
modmain::ivk
integer, dimension(:,:), allocatable ivk
Definition
modmain.f90:465
modmain::ha_ev
real(8), parameter ha_ev
Definition
modmain.f90:1256
modmain::ivkik
integer, dimension(:,:,:), allocatable ivkik
Definition
modmain.f90:467
modmain::evalsv
real(8), dimension(:,:), allocatable evalsv
Definition
modmain.f90:918
modw90
Definition
modw90.f90:6
modw90::idxw90
integer, dimension(:), allocatable idxw90
Definition
modw90.f90:22
modw90::num_bands
integer num_bands
Definition
modw90.f90:20
modw90::seedname
character(256) seedname
Definition
modw90.f90:12
writew90eig
subroutine writew90eig
Definition
writew90eig.f90:7
writew90eig.f90
Generated by
1.9.8