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
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): ",I0," of ",I0," 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 ! number of files to plot (2 for collinear magnetism, 1 otherwise)
57 if (ndmag == 1) then
58  nf=2
59 else
60  nf=1
61 end if
62 do f=1,nf
63  if (nf == 2) then
64  if (f == 1) then
65  open(50,file='FERMISURF_UP.bxsf',form='FORMATTED')
66  jst0=1; jst1=nstfv
67  else
68  open(50,file='FERMISURF_DN.bxsf',form='FORMATTED')
69  jst0=nstfv+1; jst1=2*nstfv
70  end if
71  else
72  open(50,file='FERMISURF.bxsf',form='FORMATTED')
73  jst0=1; jst1=nstsv
74  end if
75 ! find the range of eigenvalues which contribute to the Fermi surface (Lars)
76  ist0=jst1; ist1=jst0
77  do ist=jst0,jst1
78  e0=minval(evalsv(ist,:)); e1=maxval(evalsv(ist,:))
79 ! determine if the band crosses the Fermi energy
80  if ((e0 < efermi).and.(e1 > efermi)) then
81  ist0=min(ist0,ist); ist1=max(ist1,ist)
82  end if
83  end do
84  nst=ist1-ist0+1
85  write(50,'(" BEGIN_INFO")')
86  write(50,'(" # Band-XCRYSDEN-Structure-File for Fermi surface plotting")')
87  write(50,'(" # created by Elk version ",I0,".",I0,".",I0)') version
88  write(50,'(" # Launch as: xcrysden --bxsf FERMISURF(_UP/_DN).bxsf")')
89  write(50,'(" Fermi Energy: ",G18.10)') 0.d0
90  write(50,'(" END_INFO")')
91  write(50,'(" BEGIN_BLOCK_BANDGRID_3D")')
92  write(50, '(" band_energies")')
93  write(50,'(" BANDGRID_3D_BANDS")')
94  write(50,'(I4)') nst
95  write(50,'(3I6)') ngridk(:)+1
96  write(50,'(3G18.10)') 0.d0,0.d0,0.d0
97  do i=1,3
98  write(50,'(3G18.10)') bvec(:,i)
99  end do
100  do ist=ist0,ist1
101  write(50,'(" BAND: ",I4)') ist
102  do i1=0,ngridk(1)
103  j1=mod(i1,ngridk(1))
104  do i2=0,ngridk(2)
105  j2=mod(i2,ngridk(2))
106  do i3=0,ngridk(3)
107  j3=mod(i3,ngridk(3))
108  ik=ivkik(j1,j2,j3)
109  write(50,'(G18.10)') evalsv(ist,ik)-efermi
110  end do
111  end do
112  end do
113  end do
114  write(50,'(" END_BANDGRID_3D")')
115  write(50,'(" END_BLOCK_BANDGRID_3D")')
116  close(50)
117 end do
118 write(*,*)
119 write(*,'("Info(fermisurfbxsf):")')
120 if (ndmag == 1) then
121  write(*,'(" 3D Fermi surface data written to FERMISURF_UP.bxsf and &
122  &FERMISURF_DN.bxsf")')
123 else
124  write(*,'(" 3D Fermi surface data written to FERMISURF.bxsf")')
125 end if
126 write(*,'(" for plotting with XCrysDen (Fermi energy set to zero)")')
127 write(*,*)
128 write(*,'(" Launch as: xcrysden --bxsf FERMISURF(_UP/_DN).bxsf")')
129 ! restore original parameters
131 end subroutine
132 
integer nmatmax
Definition: modmain.f90:849
real(8) efermi
Definition: modmain.f90:901
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:913
integer nkpt
Definition: modmain.f90:461
integer ndmag
Definition: modmain.f90:238
integer reducek0
Definition: modmain.f90:455
Definition: modomp.f90:6
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
integer nstsv
Definition: modmain.f90:883
subroutine eveqn(ik, evalfv, evecfv, evecsv)
Definition: eveqn.f90:10
integer, dimension(3) ngridk
Definition: modmain.f90:448
real(8), dimension(3, 3) bvec
Definition: modmain.f90:16
subroutine init1
Definition: init1.f90:10
subroutine readstate
Definition: readstate.f90:10
integer, dimension(3), parameter version
Definition: modmain.f90:1281
subroutine freethd(nthd)
Definition: modomp.f90:112
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer reducek
Definition: modmain.f90:455
subroutine init0
Definition: init0.f90:10
integer nstfv
Definition: modmain.f90:881
integer nspnfv
Definition: modmain.f90:289
subroutine fermisurfbxsf