The Elk Code
elnes.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 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 elnes
7 use modmain
8 use modomp
9 use modtest
10 implicit none
11 ! local variables
12 integer ik,jk,ikq,isym,nsk(3)
13 integer ist,jst,iw,n,nthd
14 real(8) vgqc(3),gqc
15 real(8) vkql(3),v(3)
16 real(8) q,wd,dw,w,t1
17 ! allocatable arrays
18 real(8), allocatable :: jlgqr(:,:),ddcs(:),e(:,:,:),f(:,:,:)
19 complex(8), allocatable :: ylmgq(:),sfacgq(:),expmt(:,:),emat(:,:)
20 ! initialise universal variables
21 call init0
22 call init1
23 call init2
24 ! check q-vector is commensurate with k-point grid
25 v(:)=dble(ngridk(:))*vecql(:)
26 v(:)=abs(v(:)-nint(v(:)))
27 if ((v(1) > epslat).or.(v(2) > epslat).or.(v(3) > epslat)) then
28  write(*,*)
29  write(*,'("Error(elnes): q-vector incommensurate with k-point grid")')
30  write(*,'(" ngridk :",3(X,I0))') ngridk
31  write(*,'(" vecql : ",3G18.10)') vecql
32  write(*,*)
33  stop
34 end if
35 ! read in the density and potentials from file
36 call readstate
37 ! find the new linearisation energies
38 call linengy
39 ! generate the APW radial functions
40 call genapwfr
41 ! generate the local-orbital radial functions
42 call genlofr
43 ! get the second-variational eigenvalues and occupation numbers from file
44 do ik=1,nkpt
45  call getevalsv(filext,ik,vkl(:,ik),evalsv(:,ik))
46  call getoccsv(filext,ik,vkl(:,ik),occsv(:,ik))
47 end do
48 ! generate the phase factor function exp(iq.r) in the muffin-tins
49 allocate(jlgqr(njcmax,nspecies))
50 allocate(ylmgq(lmmaxo),sfacgq(natmtot))
51 allocate(expmt(npcmtmax,natmtot))
52 call gengqf(1,vecqc,vgqc,gqc,jlgqr,ylmgq,sfacgq)
53 call genexpmt(1,jlgqr,ylmgq,1,sfacgq,expmt)
54 deallocate(jlgqr,ylmgq,sfacgq)
55 allocate(e(nstsv,nstsv,nkptnr),f(nstsv,nstsv,nkptnr))
56 e(:,:,:)=0.d0
57 f(:,:,:)=0.d0
58 ! begin parallel loop over non-reduced k-points
59 call holdthd(nkptnr,nthd)
60 !$OMP PARALLEL DEFAULT(SHARED) &
61 !$OMP PRIVATE(emat,jk,vkql,isym) &
62 !$OMP PRIVATE(ikq,ist,jst,t1) &
63 !$OMP NUM_THREADS(nthd)
64 allocate(emat(nstsv,nstsv))
65 !$OMP DO SCHEDULE(DYNAMIC)
66 do ik=1,nkptnr
67 !$OMP CRITICAL(elnes_)
68  write(*,'("Info(elnes): ",I0," of ",I0," k-points")') ik,nkptnr
69 !$OMP END CRITICAL(elnes_)
70 ! equivalent reduced k-point
71  jk=ivkik(ivk(1,ik),ivk(2,ik),ivk(3,ik))
72 ! k+q-vector in lattice coordinates
73  vkql(:)=vkl(:,ik)+vecql(:)
74 ! index to k+q-vector
75  call findkpt(vkql,isym,ikq)
76 ! compute < i,k+q | exp(iq.r) | j,k > matrix elements
77  call genexpmat(vkl(:,ik),expmt,emat)
78 ! add to the double differential scattering cross-section
79  do jst=1,nstsv
80  if (evalsv(jst,jk) < emaxelnes) then
81  do ist=1,nstsv
82  e(ist,jst,ik)=evalsv(ist,ikq)-evalsv(jst,jk)
83  t1=dble(emat(ist,jst))**2+aimag(emat(ist,jst))**2
84  f(ist,jst,ik)=t1*occsv(jst,jk)*(occmax-occsv(ist,ikq))
85  end do
86  end if
87  end do
88 end do
89 !$OMP END DO
90 deallocate(emat)
91 !$OMP END PARALLEL
92 call freethd(nthd)
93 ! number of subdivisions used for interpolation
94 nsk(:)=max(ngrkf/ngridk(:),1)
95 n=nstsv*nstsv
96 ! integrate over the Brillouin zone
97 allocate(ddcs(nwplot))
98 call brzint(nswplot,ngridk,nsk,ivkiknr,nwplot,wplot,n,n,e,f,ddcs)
99 q=sqrt(vecqc(1)**2+vecqc(2)**2+vecqc(3)**2)
100 t1=2.d0/(omega*occmax)
101 if (q > epslat) t1=t1/q**4
102 ddcs(:)=t1*ddcs(:)
103 open(50,file='ELNES.OUT',form='FORMATTED')
104 wd=wplot(2)-wplot(1)
105 dw=wd/dble(nwplot)
106 do iw=1,nwplot
107  w=dw*dble(iw-1)+wplot(1)
108  write(50,'(2G18.10)') w,ddcs(iw)
109 end do
110 close(50)
111 write(*,*)
112 write(*,'("Info(elnes):")')
113 write(*,'(" ELNES double differential cross-section written to ELNES.OUT")')
114 ! write ELNES distribution to test file
115 call writetest(140,'ELNES cross-section',nv=nwplot,tol=1.d-2,rva=ddcs)
116 deallocate(e,f,ddcs,expmt)
117 end subroutine
118 
integer ngrkf
Definition: modmain.f90:1073
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition: modtest.f90:16
subroutine genexpmat(vpl, expmt, emat)
Definition: genexpmat.f90:7
integer npcmtmax
Definition: modmain.f90:216
character(256) filext
Definition: modmain.f90:1299
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:919
integer lmmaxo
Definition: modmain.f90:203
integer nkpt
Definition: modmain.f90:461
subroutine genlofr
Definition: genlofr.f90:10
subroutine gengqf(ng, vqpc, vgqc, gqc, jlgqr, ylmgq, sfacgq)
Definition: gengqf.f90:7
real(8) omega
Definition: modmain.f90:20
subroutine getevalsv(fext, ikp, vpl, evalsv_)
Definition: getevalsv.f90:7
Definition: modomp.f90:6
integer, dimension(:,:,:), allocatable ivkik
Definition: modmain.f90:467
subroutine brzint(nsm, ngridk, nsk, ivkik, nw, wint, n, ld, e, f, g)
Definition: brzint.f90:10
integer nkptnr
Definition: modmain.f90:463
real(8), dimension(3) vecql
Definition: modmain.f90:1102
subroutine linengy
Definition: linengy.f90:10
integer nstsv
Definition: modmain.f90:887
real(8), dimension(2) wplot
Definition: modmain.f90:1077
real(8) occmax
Definition: modmain.f90:899
subroutine elnes
Definition: elnes.f90:7
subroutine init2
Definition: init2.f90:7
integer, dimension(3) ngridk
Definition: modmain.f90:448
real(8), dimension(:,:), allocatable occsv
Definition: modmain.f90:903
subroutine init1
Definition: init1.f90:10
subroutine genexpmt(ngp, jlgpr, ylmgp, ld, sfacgp, expmt)
Definition: genexpmt.f90:7
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
subroutine genapwfr
Definition: genapwfr.f90:10
real(8) epslat
Definition: modmain.f90:24
subroutine getoccsv(fext, ikp, vpl, occsvp)
Definition: getoccsv.f90:7
subroutine readstate
Definition: readstate.f90:10
integer nswplot
Definition: modmain.f90:1075
integer nspecies
Definition: modmain.f90:34
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer njcmax
Definition: modmain.f90:1171
real(8) emaxelnes
Definition: modmain.f90:1104
subroutine init0
Definition: init0.f90:10
integer natmtot
Definition: modmain.f90:40
subroutine findkpt(vpl, isym, ik)
Definition: findkpt.f90:7
integer nwplot
Definition: modmain.f90:1071
integer, dimension(:,:), allocatable ivk
Definition: modmain.f90:465
real(8), dimension(3) vecqc
Definition: modmain.f90:1102
integer, dimension(:,:,:), allocatable ivkiknr
Definition: modmain.f90:469