The Elk Code
 
Loading...
Searching...
No Matches
writeexpmat.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
6subroutine writeexpmat
7use modmain
8implicit none
9! local variables
10integer nk,ik,jk,i,j
11real(8) vgqc(3),gqc
12real(8) a,b
13! allocatable arrays
14real(8), allocatable :: jlgqr(:,:)
15complex(8), allocatable :: ylmgq(:),sfacgq(:)
16complex(8), allocatable :: expmt(:,:),emat(:,:)
17! initialise universal variables
18call init0
19call init1
20call init2
21! read in the density and potentials from file
22call readstate
23! find the new linearisation energies
24call linengy
25! generate the APW radial functions
26call genapwfr
27! generate the local-orbital radial functions
28call genlofr
29! generate the phase factor function exp(iq.r) in the muffin-tins
30allocate(jlgqr(njcmax,nspecies))
31allocate(ylmgq(lmmaxo),sfacgq(natmtot))
32allocate(expmt(npcmtmax,natmtot))
33call gengqf(1,vecqc,vgqc,gqc,jlgqr,ylmgq,sfacgq)
34call genexpmt(1,jlgqr,ylmgq,1,sfacgq,expmt)
35expmt(:,:)=omega*expmt(:,:)
36deallocate(jlgqr,ylmgq,sfacgq)
37! number of k-points to write out
38if (kstlist(1,1) < 1) then
39 nk=nkpt
40else
41 nk=nkstlist
42end if
43open(50,file='EXPIQR.OUT',form='FORMATTED')
44write(50,*)
45write(50,'("q-vector (lattice coordinates) :")')
46write(50,'(3G18.10)') vecql
47write(50,'("q-vector (Cartesian coordinates) :")')
48write(50,'(3G18.10)') vecqc
49write(50,*)
50write(50,'(I8," : number of k-points")') nk
51write(50,'(I6," : number of states per k-point")') nstsv
52allocate(emat(nstsv,nstsv))
53do jk=1,nk
54 if (kstlist(1,1) < 1) then
55 ik=jk
56 else
57 ik=kstlist(1,jk)
58 end if
59 if ((ik < 1).or.(ik > nkpt)) then
60 write(*,*)
61 write(*,'("Error(writeexpiqr): k-point out of range : ",I8)') ik
62 write(*,*)
63 stop
64 end if
65 write(50,*)
66 write(50,'(" k-point (lattice coordinates) :")')
67 write(50,'(3G18.10)') vkl(:,ik)
68 write(50,*)
69 write(50,'(" k-point (Cartesian coordinates) :")')
70 write(50,'(3G18.10)') vkc(:,ik)
71 call genexpmat(vkl(:,ik),expmt,emat)
72 do i=1,nstsv
73 write(50,*)
74 write(50,²'(I6," : state i; state j, <...>, |<...>| below")') i
75 do j=1,nstsv
76 a=dble(emat(i,j))
77 b=aimag(emat(i,j))
78 write(50,'(I6,3G18.10)') j,a,b,a**2+b**2
79 end do
80 end do
81! end loop over k-points
82end do
83close(50)
84write(*,*)
85write(*,'("Info(writeexpmat)")')
86write(*,'(" < i,k+q | exp(iq.r) | j,k > matrix elements written to &
87 &EXPIQR.OUT")')
88write(*,'(" for the q-vector in vecql and all k-points in kstlist")')
89deallocate(expmt,emat)
90end subroutine
91
subroutine genapwfr
Definition genapwfr.f90:10
subroutine genexpmat(vpl, expmt, emat)
Definition genexpmat.f90:7
subroutine genexpmt(ngp, jlgpr, ylmgp, ld, sfacgp, expmt)
Definition genexpmt.f90:7
subroutine gengqf(ng, vqpc, vgqc, gqc, jlgqr, ylmgq, sfacgq)
Definition gengqf.f90:7
subroutine genlofr
Definition genlofr.f90:10
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine init2
Definition init2.f90:7
subroutine linengy
Definition linengy.f90:10
integer njcmax
Definition modmain.f90:1170
real(8), dimension(3) vecql
Definition modmain.f90:1101
real(8) omega
Definition modmain.f90:20
real(8), dimension(3) vecqc
Definition modmain.f90:1101
integer nkpt
Definition modmain.f90:461
integer natmtot
Definition modmain.f90:40
integer nkstlist
Definition modmain.f90:924
integer npcmtmax
Definition modmain.f90:216
integer nstsv
Definition modmain.f90:886
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8), dimension(:,:), allocatable vkc
Definition modmain.f90:473
integer lmmaxo
Definition modmain.f90:203
integer nspecies
Definition modmain.f90:34
integer, dimension(2, maxkst) kstlist
Definition modmain.f90:926
subroutine readstate
Definition readstate.f90:10
subroutine writeexpmat