The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine writew90mmn
7use modmain
8use modw90
9implicit none
10! local variables
11integer ik,jk,ist,jst,i
12integer ngp(nspnfv),ngpq(nspnfv)
13real(8) vl(3),vc(3),q,vkql(3)
14character(256) fname
15! automatic arrays
16complex(8) ylmq(lmmaxo),sfacq(natmtot)
17! allocatable arrays
18integer, allocatable :: igpig(:,:),igpqig(:,:)
19real(8), allocatable :: jlqr(:,:)
20complex(8), allocatable :: wfmt(:,:,:,:),wfir(:,:,:)
21complex(8), allocatable :: wfmtq(:,:,:,:),wfgpq(:,:,:)
22complex(8), allocatable :: expqmt(:,:),oq(:,:)
23! allocate local arrays
24allocate(igpig(ngkmax,nspnfv),igpqig(ngkmax,nspnfv))
25allocate(jlqr(njcmax,nspecies))
26allocate(wfmt(npcmtmax,natmtot,nspinor,num_bands))
27allocate(wfir(ngtc,nspinor,num_bands))
28allocate(wfmtq(npcmtmax,natmtot,nspinor,num_bands))
29allocate(wfgpq(ngkmax,nspinor,num_bands))
30allocate(expqmt(npcmtmax,natmtot))
31allocate(oq(num_bands,num_bands))
32fname=trim(seedname)//'.mmn'
33open(50,file=trim(fname),action='WRITE',form='FORMATTED')
34write(50,'("Generated by Elk version ",I0,".",I0,".",I0)') version
35write(50,'(3I8)') num_bands,nkptnr,nntot
36do 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
67end do
68close(50)
69write(*,*)
70write(*,'("Info(writew90mmn): created file ",A)') trim(fname)
71deallocate(igpig,igpqig,jlqr)
72deallocate(wfmt,wfir,wfmtq,wfgpq)
73deallocate(expqmt,oq)
74end subroutine
75
subroutine genexpmt(ngp, jlgpr, ylmgp, ld, sfacgp, expmt)
Definition genexpmt.f90:7
subroutine genjlgpr(ngp, gpc, jlgpr)
Definition genjlgpr.f90:7
subroutine genolpq(nst, expqmt, ngpq, igpqig, wfmt, wfir, wfmtq, wfgpq, oq)
Definition genolpq.f90:7
pure subroutine gensfacgp(ngp, vgpc, ld, sfacgp)
Definition gensfacgp.f90:10
subroutine genwfsvp(tsh, tgp, nst, idx, ngridg_, igfft_, vpl, ngp, igpig, wfmt, ld, wfir)
Definition genwfsvp.f90:7
pure subroutine genylmv(t4pil, lmax, v, ylm)
Definition genylmv.f90:10
integer njcmax
Definition modmain.f90:1170
integer nkptnr
Definition modmain.f90:463
integer nspinor
Definition modmain.f90:267
integer, dimension(3) ngdgc
Definition modmain.f90:388
real(8), dimension(3, 3) bvec
Definition modmain.f90:16
integer, dimension(:), allocatable igfc
Definition modmain.f90:410
integer ngtc
Definition modmain.f90:392
integer npcmtmax
Definition modmain.f90:216
integer, dimension(3), parameter version
Definition modmain.f90:1288
integer lmaxo
Definition modmain.f90:201
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
integer ngkmax
Definition modmain.f90:499
integer nspecies
Definition modmain.f90:34
integer nntot
Definition modw90.f90:32
integer, dimension(:,:), allocatable nnlist
Definition modw90.f90:34
integer, dimension(:), allocatable idxw90
Definition modw90.f90:22
integer num_bands
Definition modw90.f90:20
integer, dimension(:,:,:), allocatable nncell
Definition modw90.f90:36
character(256) seedname
Definition modw90.f90:12
pure subroutine r3mv(a, x, y)
Definition r3mv.f90:10
subroutine writew90mmn