The Elk Code
genvmat.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
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 genvmat(vmt,vir,vmat)
7 ! generates potential matrix elements for all states and k-points
8 use modmain
9 use modmpi
10 use modomp
11 implicit none
12 ! arguments
13 real(8), intent(in) :: vmt(npmtmax,natmtot),vir(ngtot)
14 complex(8), intent(out) :: vmat(nstsv,nstsv,nkpt)
15 ! local variables
16 integer ik,ispn,is,ias
17 integer nrc,nrci,n,lp,nthd
18 ! automatic arrays
19 real(8) vmtc(npcmtmax,natmtot),virc(ngtc)
20 ! allocatable arrays
21 complex(8), allocatable :: apwalm(:,:,:,:,:),evecfv(:,:),evecsv(:,:)
22 complex(4), allocatable :: wfmt(:,:,:,:),wfgk(:,:,:)
23 call holdthd(natmtot,nthd)
24 !$OMP PARALLEL DO DEFAULT(SHARED) &
25 !$OMP PRIVATE(is,nrc,nrci) &
26 !$OMP SCHEDULE(DYNAMIC) &
27 !$OMP NUM_THREADS(nthd)
28 do ias=1,natmtot
29  is=idxis(ias)
30  nrc=nrcmt(is)
31  nrci=nrcmti(is)
32 ! convert muffin-tin potential to a coarse radial mesh
33  call rfmtftoc(nrc,nrci,vmt(:,ias),vmtc(:,ias))
34 ! convert potential to spherical coordinates
35  call rbshtip(nrc,nrci,vmtc(:,ias))
36 ! multiply by the radial integration weights
37  call rfcmtwr(nrc,nrci,wr2cmt(:,is),vmtc(:,ias))
38 end do
39 !$OMP END PARALLEL DO
40 call freethd(nthd)
41 ! multiply by characteristic function and convert to coarse grid
42 call rfirftoc(vir,virc)
43 ! loop over k-points
44 call holdthd(nkpt/np_mpi,nthd)
45 !$OMP PARALLEL DEFAULT(SHARED) &
46 !$OMP PRIVATE(apwalm,evecfv,evecsv) &
47 !$OMP PRIVATE(wfmt,wfgk,ispn) &
48 !$OMP NUM_THREADS(nthd)
49 allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv))
50 allocate(evecfv(nmatmax,nstfv),evecsv(nstsv,nstsv))
51 allocate(wfmt(npcmtmax,natmtot,nspinor,nstsv),wfgk(ngkmax,nspinor,nstsv))
52 !$OMP DO SCHEDULE(DYNAMIC)
53 do ik=1,nkpt
54 ! distribute among MPI processes
55  if (mod(ik-1,np_mpi) /= lp_mpi) cycle
56 ! find the matching coefficients
57  do ispn=1,nspnfv
58  call match(ngk(ispn,ik),vgkc(:,:,ispn,ik),gkc(:,ispn,ik), &
59  sfacgk(:,:,ispn,ik),apwalm(:,:,:,:,ispn))
60  end do
61 ! get the eigenvectors from file
62  call getevecfv(filext,ik,vkl(:,ik),vgkl(:,:,:,ik),evecfv)
63  call getevecsv(filext,ik,vkl(:,ik),evecsv)
64 ! calculate the wavefunctions for all states of the input k-point
65  call genwfsv_sp(.false.,.true.,nstsv,[0],ngridg,igfft,ngk(:,ik), &
66  igkig(:,:,ik),apwalm,evecfv,evecsv,wfmt,ngkmax,wfgk)
67  call genvmatk(vmtc,virc,ngk(:,ik),igkig(:,:,ik),wfmt,ngkmax,wfgk,vmat(:,:,ik))
68 end do
69 !$OMP END DO
70 deallocate(apwalm,evecfv,evecsv,wfmt,wfgk)
71 !$OMP END PARALLEL
72 call freethd(nthd)
73 ! broadcast matrix elements to every process
74 if (np_mpi > 1) then
75  n=nstsv*nstsv
76  do ik=1,nkpt
77  lp=mod(ik-1,np_mpi)
78  call mpi_bcast(vmat(:,:,ik),n,mpi_double_complex,lp,mpicom,ierror)
79  end do
80 end if
81 end subroutine
82 
integer nmatmax
Definition: modmain.f90:853
subroutine rbshtip(nr, nri, rfmt)
Definition: rbshtip.f90:7
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
integer, dimension(3) ngridg
Definition: modmain.f90:386
character(256) filext
Definition: modmain.f90:1301
integer lmmaxapw
Definition: modmain.f90:199
integer ngkmax
Definition: modmain.f90:499
subroutine getevecfv(fext, ikp, vpl, vgpl, evecfv)
Definition: getevecfv.f90:10
subroutine match(ngp, vgpc, gpc, sfacgp, apwalm)
Definition: match.f90:10
Definition: modomp.f90:6
subroutine genvmat(vmt, vir, vmat)
Definition: genvmat.f90:7
pure subroutine rfcmtwr(nr, nri, wr, rfmt)
Definition: rfcmtwr.f90:7
integer np_mpi
Definition: modmpi.f90:13
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition: modmain.f90:509
integer, dimension(:,:), allocatable ngk
Definition: modmain.f90:497
subroutine genvmatk(vmt, vir, ngp, igpig, wfmt, ld, wfgp, vmat)
Definition: genvmatk.f90:7
integer, dimension(:), allocatable igfft
Definition: modmain.f90:406
subroutine genwfsv_sp(tsh, tgp, nst, idx, ngridg_, igfft_, ngp, igpig, apwalm, evecfv, evecsv, wfmt, ld, wfir)
Definition: genwfsv_sp.f90:8
real(8), dimension(:,:,:,:), allocatable vgkl
Definition: modmain.f90:503
integer nspinor
Definition: modmain.f90:267
real(8), dimension(:,:,:,:), allocatable vgkc
Definition: modmain.f90:505
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
integer apwordmax
Definition: modmain.f90:760
real(8), dimension(:,:), allocatable wr2cmt
Definition: modmain.f90:189
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
Definition: modmpi.f90:6
real(8), dimension(:,:,:), allocatable gkc
Definition: modmain.f90:507
integer lp_mpi
Definition: modmpi.f90:15
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
pure subroutine rfmtftoc(nrc, nrci, rfmt, rfcmt)
Definition: rfmtftoc.f90:7
integer, dimension(:,:,:), allocatable igkig
Definition: modmain.f90:501
subroutine rfirftoc(rfir, rfirc)
Definition: rfirftoc.f90:7
integer nstfv
Definition: modmain.f90:887
integer mpicom
Definition: modmpi.f90:11
integer nspnfv
Definition: modmain.f90:289
integer ierror
Definition: modmpi.f90:19