The Elk Code
writew90mmn.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2015 Manh Duc Le, 2017-18 Arsenii Gerasimov, Yaroslav Kvashnin
3 ! and Lars Nordstrom. This file is distributed under the terms of the GNU
4 ! General Public License. See the file COPYING for license details.
5 
6 subroutine writew90mmn
7 use modmain
8 use modw90
9 implicit none
10 ! local variables
11 integer ik,jk,ist,jst,i
12 integer ngp(nspnfv),ngpq(nspnfv)
13 real(8) vl(3),vc(3),q,vkql(3)
14 character(256) fname
15 ! automatic arrays
16 complex(8) ylmq(lmmaxo),sfacq(natmtot)
17 ! allocatable arrays
18 integer, allocatable :: igpig(:,:),igpqig(:,:)
19 real(8), allocatable :: jlqr(:,:)
20 complex(8), allocatable :: wfmt(:,:,:,:),wfir(:,:,:)
21 complex(8), allocatable :: wfmtq(:,:,:,:),wfgpq(:,:,:)
22 complex(8), allocatable :: expqmt(:,:),oq(:,:)
23 ! allocate local arrays
24 allocate(igpig(ngkmax,nspnfv),igpqig(ngkmax,nspnfv))
25 allocate(jlqr(njcmax,nspecies))
26 allocate(wfmt(npcmtmax,natmtot,nspinor,num_bands))
27 allocate(wfir(ngtc,nspinor,num_bands))
28 allocate(wfmtq(npcmtmax,natmtot,nspinor,num_bands))
29 allocate(wfgpq(ngkmax,nspinor,num_bands))
30 allocate(expqmt(npcmtmax,natmtot))
31 allocate(oq(num_bands,num_bands))
32 fname=trim(seedname)//'.mmn'
33 open(50,file=trim(fname),action='WRITE',form='FORMATTED')
34 write(50,'("Generated by Elk version ",I0,".",I0,".",I0)') version
35 write(50,'(3I8)') num_bands,nkptnr,nntot
36 do ik=1,nkptnr
37  call genwfsvp(.false.,.false.,num_bands,idxw90,ngdgc,igfc,vkl(:,ik),ngp, &
38  igpig,wfmt,ngtc,wfir)
39  do i=1,nntot
40  jk=nnlist(ik,i)
41 ! q-vector in lattice coordinates
42  vl(:)=dble(nncell(:,ik,i))+vkl(:,jk)-vkl(:,ik)
43 ! q-vector in Cartesian coordinates
44  call r3mv(bvec,vl,vc)
45 ! q-vector length
46  q=sqrt(vc(1)**2+vc(2)**2+vc(3)**2)
47 ! generate phase factor function exp(iq.r) in the muffin-tins
48  call genjlgpr(1,q,jlqr)
49  call genylmv(.true.,lmaxo,vc,ylmq)
50  call gensfacgp(1,vc,1,sfacq)
51  call genexpmt(1,jlqr,ylmq,1,sfacq,expqmt)
52 ! k+q-vector in lattice coordinates
53  vkql(:)=vkl(:,ik)+vl(:)
54 ! generate the wavefunctions at k+q
55  call genwfsvp(.false.,.true.,num_bands,idxw90,ngdgc,igfc,vkql,ngpq,igpqig, &
56  wfmtq,ngkmax,wfgpq)
57 ! determine the overlap matrix
58  call genolpq(num_bands,expqmt,ngpq,igpqig,wfmt,wfir,wfmtq,wfgpq,oq)
59 ! write overlap matrix to file
60  write(50,'(5I8)') ik,jk,nncell(:,ik,i)
61  do jst=1,num_bands
62  do ist=1,num_bands
63  write(50,'(2G18.10)') dble(oq(jst,ist)),-aimag(oq(jst,ist))
64  end do
65  end do
66  end do
67 end do
68 close(50)
69 write(*,*)
70 write(*,'("Info(writew90mmn): created file ",A)') trim(fname)
71 deallocate(igpig,igpqig,jlqr)
72 deallocate(wfmt,wfir,wfmtq,wfgpq)
73 deallocate(expqmt,oq)
74 end subroutine
75 
integer num_bands
Definition: modw90.f90:20
integer npcmtmax
Definition: modmain.f90:216
pure subroutine gensfacgp(ngp, vgpc, ld, sfacgp)
Definition: gensfacgp.f90:10
subroutine genolpq(nst, expqmt, ngpq, igpqig, wfmt, wfir, wfmtq, wfgpq, oq)
Definition: genolpq.f90:7
integer ngtc
Definition: modmain.f90:392
subroutine writew90mmn
Definition: writew90mmn.f90:7
integer ngkmax
Definition: modmain.f90:499
integer lmaxo
Definition: modmain.f90:201
integer nkptnr
Definition: modmain.f90:463
pure subroutine genylmv(t4pil, lmax, v, ylm)
Definition: genylmv.f90:10
character(256) seedname
Definition: modw90.f90:12
integer, dimension(:,:,:), allocatable nncell
Definition: modw90.f90:36
Definition: modw90.f90:6
integer, dimension(:), allocatable igfc
Definition: modmain.f90:410
integer, dimension(:), allocatable idxw90
Definition: modw90.f90:22
integer nspinor
Definition: modmain.f90:267
integer nntot
Definition: modw90.f90:32
real(8), dimension(3, 3) bvec
Definition: modmain.f90:16
subroutine genexpmt(ngp, jlgpr, ylmgp, ld, sfacgp, expmt)
Definition: genexpmt.f90:7
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
integer, dimension(3) ngdgc
Definition: modmain.f90:388
subroutine genjlgpr(ngp, gpc, jlgpr)
Definition: genjlgpr.f90:7
integer, dimension(3), parameter version
Definition: modmain.f90:1289
integer nspecies
Definition: modmain.f90:34
integer njcmax
Definition: modmain.f90:1173
pure subroutine r3mv(a, x, y)
Definition: r3mv.f90:10
subroutine genwfsvp(tsh, tgp, nst, idx, ngridg_, igfft_, vpl, ngp, igpig, wfmt, ld, wfir)
Definition: genwfsvp.f90:7
integer, dimension(:,:), allocatable nnlist
Definition: modw90.f90:34