The Elk Code
fermisurf.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2008 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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 fermisurf
7 use modmain
8 use modomp
9 implicit none
10 ! local variables
11 integer ik,ist,np,nfp,i
12 integer ist0,ist1,jst0,jst1
13 integer l,i1,i2,i3,nthd
14 real(8) e0,e1,v(3),x,fn
15 ! allocatable arrays
16 real(8), allocatable :: evalfv(:,:),vpc(:,:)
17 complex(8), allocatable :: evecfv(:,:,:),evecsv(:,:)
18 ! external functions
19 real(8), external :: sdelta
20 ! initialise universal variables
21 call init0
22 ! no k-point reduction for the collinear case
24 if (ndmag == 1) reducek=0
25 call init1
26 ! read density and potentials from file
27 call readstate
28 ! Fourier transform Kohn-Sham potential to G-space
29 call genvsig
30 ! find the new linearisation energies
31 call linengy
32 ! generate the APW and local-orbital radial functions and integrals
33 call genapwlofr
34 ! generate the spin-orbit coupling radial functions
35 call gensocfr
36 ! begin parallel loop over reduced k-points set
37 call holdthd(nkpt,nthd)
38 !$OMP PARALLEL DEFAULT(SHARED) &
39 !$OMP PRIVATE(evalfv,evecfv,evecsv) &
40 !$OMP NUM_THREADS(nthd)
41 allocate(evalfv(nstfv,nspnfv))
42 allocate(evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv))
43 !$OMP DO SCHEDULE(DYNAMIC)
44 do ik=1,nkpt
45 !$OMP CRITICAL(fermisurf_)
46  write(*,'("Info(fermisurf): ",I0," of ",I0," k-points")') ik,nkpt
47 !$OMP END CRITICAL(fermisurf_)
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 ! generate the plotting point grid in Cartesian coordinates (this has the same
57 ! arrangement as the k-point grid)
58 np=np3d(1)*np3d(2)*np3d(3)
59 allocate(vpc(3,np))
60 call plotpt3d(vpc)
61 do i=1,np
62  v(:)=vpc(:,i)
63  call r3mv(bvec,v,vpc(:,i))
64 end do
65 ! number of files to plot (2 for collinear magnetism, 1 otherwise)
66 if (ndmag == 1) then
67  nfp=2
68 else
69  nfp=1
70 end if
71 do l=1,nfp
72  if (nfp == 2) then
73  if (l == 1) then
74  open(50,file='FERMISURF_UP.OUT',form='FORMATTED',action='WRITE')
75  jst0=1; jst1=nstfv
76  else
77  open(50,file='FERMISURF_DN.OUT',form='FORMATTED',action='WRITE')
78  jst0=nstfv+1; jst1=2*nstfv
79  end if
80  else
81  open(50,file='FERMISURF.OUT',form='FORMATTED',action='WRITE')
82  jst0=1; jst1=nstsv
83  end if
84  if ((task == 100).or.(task == 101)) then
85 ! find the range of eigenvalues which contribute to the Fermi surface (Lars)
86  ist0=jst1; ist1=jst0
87  do ist=jst0,jst1
88  e0=minval(evalsv(ist,:)); e1=maxval(evalsv(ist,:))
89 ! determine if the band crosses the Fermi energy
90  if ((e0 < efermi).and.(e1 > efermi)) then
91  ist0=min(ist0,ist); ist1=max(ist1,ist)
92  end if
93  end do
94  else
95 ! use all eigenvalues
96  ist0=jst0; ist1=jst1
97  end if
98  if ((task == 100).or.(task == 103)) then
99  write(50,'(3I6," : grid size")') np3d(:)
100  i=0
101  do i3=0,ngridk(3)-1
102  do i2=0,ngridk(2)-1
103  do i1=0,ngridk(1)-1
104  i=i+1
105  ik=ivkik(i1,i2,i3)
106  if (task == 100) then
107 ! write product of eigenstates minus the Fermi energy
108  fn=product(evalsv(ist0:ist1,ik)-efermi)
109  else
110 ! write single smooth delta function at the Fermi energy
111  fn=0.d0
112  do ist=ist0,ist1
113  x=(evalsv(ist,ik)-efermi)/swidth
114  fn=fn+sdelta(stype,x)/swidth
115  end do
116  end if
117  write(50,'(4G18.10)') vpc(:,i),fn
118  end do
119  end do
120  end do
121  else
122  write(50,'(4I6," : grid size, number of states")') np3d(:),ist1-ist0+1
123  i=0
124  do i3=0,ngridk(3)-1
125  do i2=0,ngridk(2)-1
126  do i1=0,ngridk(1)-1
127  i=i+1
128  ik=ivkik(i1,i2,i3)
129  write(50,'(3G18.10)',advance='NO') vpc(:,i)
130  do ist=ist0,ist1
131  if (task == 101) then
132 ! write the eigenvalues minus the Fermi energy separately
133  write(50,'(F14.8)',advance='NO') evalsv(ist,ik)-efermi
134  else
135 ! write separate smooth delta functions at the Fermi energy
136  x=(evalsv(ist,ik)-efermi)/swidth
137  write(50,'(F14.8)',advance='NO') sdelta(stype,x)/swidth
138  end if
139  end do
140  write(50,*)
141  end do
142  end do
143  end do
144  end if
145  close(50)
146 end do
147 write(*,*)
148 write(*,'("Info(fermisurf):")')
149 if (ndmag == 1) then
150  write(*,'(" 3D Fermi surface data written to FERMISURF_UP.OUT and &
151  &FERMISURF_DN.OUT")')
152 else
153  write(*,'(" 3D Fermi surface data written to FERMISURF.OUT")')
154 end if
155 if (task == 100) then
156  write(*,'(" in terms of the product of eigenvalues minus the Fermi energy")')
157 else if (task == 101) then
158  write(*,'(" in terms of separate eigenvalues minus the Fermi energy")')
159 else if (task == 103) then
160  write(*,'(" in terms of a smooth delta function at the Fermi energy")')
161 else
162  write(*,'(" in terms of separate delta functions at the Fermi energy")')
163 end if
164 deallocate(vpc)
165 ! restore original parameters
167 end subroutine
168 
integer nmatmax
Definition: modmain.f90:849
real(8) efermi
Definition: modmain.f90:901
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:913
integer task
Definition: modmain.f90:1291
integer nkpt
Definition: modmain.f90:461
integer ndmag
Definition: modmain.f90:238
integer reducek0
Definition: modmain.f90:455
Definition: modomp.f90:6
real(8) swidth
Definition: modmain.f90:889
integer, dimension(:,:,:), allocatable ivkik
Definition: modmain.f90:467
subroutine genvsig
Definition: genvsig.f90:10
subroutine genapwlofr
Definition: genapwlofr.f90:7
subroutine gensocfr
Definition: gensocfr.f90:10
subroutine linengy
Definition: linengy.f90:10
pure subroutine plotpt3d(vpl)
Definition: plotpt3d.f90:7
integer nstsv
Definition: modmain.f90:883
subroutine eveqn(ik, evalfv, evecfv, evecsv)
Definition: eveqn.f90:10
subroutine fermisurf
Definition: fermisurf.f90:7
integer, dimension(3) ngridk
Definition: modmain.f90:448
real(8), dimension(3, 3) bvec
Definition: modmain.f90:16
subroutine init1
Definition: init1.f90:10
integer stype
Definition: modmain.f90:885
subroutine readstate
Definition: readstate.f90:10
subroutine freethd(nthd)
Definition: modomp.f90:112
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer, dimension(3) np3d
Definition: modmain.f90:1126
integer reducek
Definition: modmain.f90:455
subroutine init0
Definition: init0.f90:10
pure subroutine r3mv(a, x, y)
Definition: r3mv.f90:10
integer nstfv
Definition: modmain.f90:881
integer nspnfv
Definition: modmain.f90:289