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 use modomp
10 implicit none
11 ! local variables
12 integer ik,jk,ist,jst,i,nthd
13 integer ngp(nspnfv),ngpq(nspnfv)
14 real(8) vl(3),vc(3),q,vkql(3)
15 character(256) fname
16 ! automatic arrays
17 integer igpig(ngkmax,nspnfv),igpqig(ngkmax,nspnfv)
18 real(8) jlqr(njcmax,nspecies)
19 complex(8) ylmq(lmmaxo),sfacq(natmtot)
20 ! allocatable arrays
21 complex(8), allocatable :: wfmt(:,:,:,:),wfir(:,:,:)
22 complex(8), allocatable :: wfmtq(:,:,:,:),wfgpq(:,:,:)
23 complex(8), allocatable :: expqmt(:,:),oq(:,:,:,:)
24 allocate(oq(num_bands,num_bands,nntot,nkptnr))
25 ! begin parallel loop over k-points
26 call holdthd(nkptnr,nthd)
27 !$OMP PARALLEL DEFAULT(SHARED) &
28 !$OMP PRIVATE(igpig,igpqig,jlqr,ylmq,sfacq) &
29 !$OMP PRIVATE(wfmt,wfir,wfmtq,wfgpq,expqmt) &
30 !$OMP PRIVATE(ngp,ngpq,i,jk,vl,vc,q,vkql) &
31 !$OMP NUM_THREADS(nthd)
32 allocate(wfmt(npcmtmax,natmtot,nspinor,num_bands))
33 allocate(wfir(ngtc,nspinor,num_bands))
34 allocate(wfmtq(npcmtmax,natmtot,nspinor,num_bands))
35 allocate(wfgpq(ngkmax,nspinor,num_bands))
36 allocate(expqmt(npcmtmax,natmtot))
37 !$OMP DO SCHEDULE(DYNAMIC)
38 do ik=1,nkptnr
39 ! generate the wavefunctions at k
40  call genwfsvp(.false.,.false.,num_bands,idxw90,ngdgc,igfc,vkl(:,ik),ngp, &
41  igpig,wfmt,ngtc,wfir)
42  do i=1,nntot
43  jk=nnlist(ik,i)
44 ! q-vector in lattice coordinates
45  vl(1:3)=dble(nncell(1:3,ik,i))+vkl(1:3,jk)-vkl(1:3,ik)
46 ! q-vector in Cartesian coordinates
47  call r3mv(bvec,vl,vc)
48 ! q-vector length
49  q=sqrt(vc(1)**2+vc(2)**2+vc(3)**2)
50 ! generate phase factor function exp(iq⋅r) in the muffin-tins
51  call genjlgpr(1,q,jlqr)
52  call genylmv(.true.,lmaxo,vc,ylmq)
53  call gensfacgp(1,vc,1,sfacq)
54  call genexpmt(1,jlqr,ylmq,1,sfacq,expqmt)
55 ! k+q-vector in lattice coordinates
56  vkql(1:3)=vkl(1:3,ik)+vl(1:3)
57 ! generate the wavefunctions at k+q
58  call genwfsvp(.false.,.true.,num_bands,idxw90,ngdgc,igfc,vkql,ngpq,igpqig, &
59  wfmtq,ngkmax,wfgpq)
60 ! determine the overlap matrix
61  call genolpq(num_bands,expqmt,ngpq,igpqig,wfmt,wfir,wfmtq,wfgpq, &
62  oq(:,:,i,ik))
63  end do
64 end do
65 !$OMP END DO
66 deallocate(wfmt,wfir,wfmtq,wfgpq,expqmt)
67 !$OMP END PARALLEL
68 call freethd(nthd)
69 ! write overlap matrix to file
70 fname=trim(seedname)//'.mmn'
71 open(50,file=trim(fname),action='WRITE',form='FORMATTED')
72 write(50,'("Generated by Elk version ",I0,".",I0,".",I0)') version
73 write(50,'(3I8)') num_bands,nkptnr,nntot
74 do ik=1,nkptnr
75  do i=1,nntot
76  jk=nnlist(ik,i)
77  write(50,'(5I8)') ik,jk,nncell(1:3,ik,i)
78  do jst=1,num_bands
79  do ist=1,num_bands
80 ! lower half of matrix required
81  write(50,'(2G18.10)') conjg(oq(jst,ist,i,ik))
82  end do
83  end do
84  end do
85 end do
86 close(50)
87 write(*,*)
88 write(*,'("Info(writew90mmn): created file ",A)') trim(fname)
89 deallocate(oq)
90 end subroutine
91 
integer num_bands
Definition: modw90.f90:22
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
Definition: modomp.f90:6
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:14
integer, dimension(:,:,:), allocatable nncell
Definition: modw90.f90:46
Definition: modw90.f90:6
integer, dimension(:), allocatable igfc
Definition: modmain.f90:410
integer, dimension(:), allocatable idxw90
Definition: modw90.f90:24
integer nspinor
Definition: modmain.f90:267
integer nntot
Definition: modw90.f90:42
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:1276
subroutine freethd(nthd)
Definition: modomp.f90:112
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
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:44