The Elk Code
fermisurfbxsf.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2009 F. Cricchio, F. Bultmark and L. Nordstrom.
3
! This file is distributed under the terms of the GNU General Public License.
4
! See the file COPYING for license details.
5
6
subroutine
fermisurfbxsf
7
use
modmain
8
use
modomp
9
implicit none
10
! local variables
11
integer
ik,nst,ist
12
integer
ist0,ist1,jst0,jst1
13
integer
i1,i2,i3,j1,j2,j3
14
integer
nf,f,i,nthd
15
real(8)
e0,e1
16
! allocatable arrays
17
real(8)
,
allocatable
:: evalfv(:,:)
18
complex(8)
,
allocatable
:: evecfv(:,:,:),evecsv(:,:)
19
! initialise universal variables
20
call
init0
21
! no k-point reduction for the collinear case
22
reducek0
=
reducek
23
if
(
ndmag
== 1)
reducek
=0
24
call
init1
25
! read density and potentials from file
26
call
readstate
27
! Fourier transform Kohn-Sham potential to G-space
28
call
genvsig
29
! find the new linearisation energies
30
call
linengy
31
! generate the APW and local-orbital radial functions and integrals
32
call
genapwlofr
33
! generate the spin-orbit coupling radial functions
34
call
gensocfr
35
! begin parallel loop over reduced k-points set
36
call
holdthd
(
nkpt
,nthd)
37
!$OMP PARALLEL DEFAULT(SHARED) &
38
!$OMP PRIVATE(evalfv,evecfv,evecsv) &
39
!$OMP NUM_THREADS(nthd)
40
allocate
(evalfv(
nstfv
,
nspnfv
))
41
allocate
(evecfv(
nmatmax
,
nstfv
,
nspnfv
))
42
allocate
(evecsv(
nstsv
,
nstsv
))
43
!$OMP DO SCHEDULE(DYNAMIC)
44
do
ik=1,
nkpt
45
!$OMP CRITICAL(fermisurfbxsf_)
46
write
(*,
'("Info(fermisurfbxsf): ",I6," of ",I6," k-points")'
) ik,
nkpt
47
!$OMP END CRITICAL(fermisurfbxsf_)
48
! solve the first- and second-variational eigenvalue equations
49
call
eveqn
(ik,evalfv,evecfv,evecsv)
50
! end loop over reduced k-points set
51
end do
52
!$OMP END DO
53
deallocate
(evalfv,evecfv,evecsv)
54
!$OMP END PARALLEL
55
call
freethd
(nthd)
56
! if iterative diagonalisation is used the eigenvalues must be reordered
57
if
(
tefvit
.and.(.not.
spinpol
))
then
58
! allocate(idx(nstsv),e(nstsv))
59
do
ik=1,
nkpt
60
call
sort
(
nstsv
,
evalsv
(:,ik))
61
end do
62
end if
63
! number of files to plot (2 for collinear magnetism, 1 otherwise)
64
if
(
ndmag
== 1)
then
65
nf=2
66
else
67
nf=1
68
end if
69
do
f=1,nf
70
if
(nf == 2)
then
71
if
(f == 1)
then
72
open
(50,file=
'FERMISURF_UP.bxsf'
,form=
'FORMATTED'
)
73
jst0=1; jst1=
nstfv
74
else
75
open
(50,file=
'FERMISURF_DN.bxsf'
,form=
'FORMATTED'
)
76
jst0=
nstfv
+1; jst1=2*
nstfv
77
end if
78
else
79
open
(50,file=
'FERMISURF.bxsf'
,form=
'FORMATTED'
)
80
jst0=1; jst1=
nstsv
81
end if
82
! find the range of eigenvalues which contribute to the Fermi surface (Lars)
83
ist0=jst1; ist1=jst0
84
do
ist=jst0,jst1
85
e0=minval(
evalsv
(ist,:)); e1=maxval(
evalsv
(ist,:))
86
! determine if the band crosses the Fermi energy
87
if
((e0 <
efermi
).and.(e1 >
efermi
))
then
88
ist0=min(ist0,ist); ist1=max(ist1,ist)
89
end if
90
end do
91
nst=ist1-ist0+1
92
write
(50,
'(" BEGIN_INFO")'
)
93
write
(50,
'(" # Band-XCRYSDEN-Structure-File for Fermi surface plotting")'
)
94
write
(50,
'(" # created by Elk version ",I0,".",I0,".",I0)'
)
version
95
write
(50,
'(" # Launch as: xcrysden --bxsf FERMISURF(_UP/_DN).bxsf")'
)
96
write
(50,
'(" Fermi Energy: ",G18.10)'
) 0.d0
97
write
(50,
'(" END_INFO")'
)
98
write
(50,
'(" BEGIN_BLOCK_BANDGRID_3D")'
)
99
write
(50,
'(" band_energies")'
)
100
write
(50,
'(" BANDGRID_3D_BANDS")'
)
101
write
(50,
'(I4)'
) nst
102
write
(50,
'(3I6)'
)
ngridk
(:)+1
103
write
(50,
'(3G18.10)'
) 0.d0,0.d0,0.d0
104
do
i=1,3
105
write
(50,
'(3G18.10)'
)
bvec
(:,i)
106
end do
107
do
ist=ist0,ist1
108
write
(50,
'(" BAND: ",I4)'
) ist
109
do
i1=0,
ngridk
(1)
110
j1=mod(i1,
ngridk
(1))
111
do
i2=0,
ngridk
(2)
112
j2=mod(i2,
ngridk
(2))
113
do
i3=0,
ngridk
(3)
114
j3=mod(i3,
ngridk
(3))
115
ik=
ivkik
(j1,j2,j3)
116
write
(50,
'(G18.10)'
)
evalsv
(ist,ik)-
efermi
117
end do
118
end do
119
end do
120
end do
121
write
(50,
'(" END_BANDGRID_3D")'
)
122
write
(50,
'(" END_BLOCK_BANDGRID_3D")'
)
123
close
(50)
124
end do
125
write
(*,*)
126
write
(*,
'("Info(fermisurfbxsf):")'
)
127
if
(
ndmag
== 1)
then
128
write
(*,
'(" 3D Fermi surface data written to FERMISURF_UP.bxsf and &
129
&FERMISURF_DN.bxsf")'
)
130
else
131
write
(*,
'(" 3D Fermi surface data written to FERMISURF.bxsf")'
)
132
end if
133
write
(*,
'(" for plotting with XCrysDen (Fermi energy set to zero)")'
)
134
write
(*,*)
135
write
(*,
'(" Launch as: xcrysden --bxsf FERMISURF(_UP/_DN).bxsf")'
)
136
! restore original parameters
137
reducek
=
reducek0
138
end subroutine
139
modmain::nmatmax
integer nmatmax
Definition:
modmain.f90:853
modmain::efermi
real(8) efermi
Definition:
modmain.f90:907
modmain::evalsv
real(8), dimension(:,:), allocatable evalsv
Definition:
modmain.f90:921
modmain::spinpol
logical spinpol
Definition:
modmain.f90:228
modmain::nkpt
integer nkpt
Definition:
modmain.f90:461
modmain::tefvit
logical tefvit
Definition:
modmain.f90:873
modmain::ndmag
integer ndmag
Definition:
modmain.f90:238
modmain::reducek0
integer reducek0
Definition:
modmain.f90:455
modomp
Definition:
modomp.f90:6
modmain::ivkik
integer, dimension(:,:,:), allocatable ivkik
Definition:
modmain.f90:467
genvsig
subroutine genvsig
Definition:
genvsig.f90:10
genapwlofr
subroutine genapwlofr
Definition:
genapwlofr.f90:7
gensocfr
subroutine gensocfr
Definition:
gensocfr.f90:7
linengy
subroutine linengy
Definition:
linengy.f90:10
modmain
Definition:
modmain.f90:6
modmain::nstsv
integer nstsv
Definition:
modmain.f90:889
eveqn
subroutine eveqn(ik, evalfv, evecfv, evecsv)
Definition:
eveqn.f90:10
modmain::ngridk
integer, dimension(3) ngridk
Definition:
modmain.f90:448
modmain::bvec
real(8), dimension(3, 3) bvec
Definition:
modmain.f90:16
init1
subroutine init1
Definition:
init1.f90:10
readstate
subroutine readstate
Definition:
readstate.f90:10
modmain::version
integer, dimension(3), parameter version
Definition:
modmain.f90:1289
modomp::freethd
subroutine freethd(nthd)
Definition:
modomp.f90:106
modomp::holdthd
subroutine holdthd(nloop, nthd)
Definition:
modomp.f90:78
sort
subroutine sort(n, x)
Definition:
sort.f90:10
modmain::reducek
integer reducek
Definition:
modmain.f90:455
init0
subroutine init0
Definition:
init0.f90:10
modmain::nstfv
integer nstfv
Definition:
modmain.f90:887
modmain::nspnfv
integer nspnfv
Definition:
modmain.f90:289
fermisurfbxsf
subroutine fermisurfbxsf
Definition:
fermisurfbxsf.f90:7
fermisurfbxsf.f90
Generated by
1.8.14