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