The Elk Code
putpmat.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2014 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 putpmat(ik)
7 use modmain
8 use modmpi
9 use modramdisk
10 implicit none
11 ! arguments
12 integer, intent(in) :: ik
13 ! local variables
14 integer ispn,recl
15 ! allocatable arrays
16 complex(8), allocatable :: apwalm(:,:,:,:,:),evecfv(:,:,:),evecsv(:,:)
17 complex(8), allocatable :: wfmt(:,:,:,:),wfgk(:,:,:),pmat(:,:,:)
18 ! get the eigenvectors from file
19 allocate(evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv))
20 call getevecfv(filext,ik,vkl(:,ik),vgkl(:,:,:,ik),evecfv)
21 call getevecsv(filext,ik,vkl(:,ik),evecsv)
22 ! find the matching coefficients
23 allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv))
24 do ispn=1,nspnfv
25  call match(ngk(ispn,ik),vgkc(:,:,ispn,ik),gkc(:,ispn,ik), &
26  sfacgk(:,:,ispn,ik),apwalm(:,:,:,:,ispn))
27 end do
28 ! calculate the wavefunctions for all states
29 allocate(wfmt(npcmtmax,natmtot,nspinor,nstsv),wfgk(ngkmax,nspinor,nstsv))
30 call genwfsv(.true.,.true.,nstsv,[0],ngridg,igfft,ngk(:,ik),igkig(:,:,ik), &
31  apwalm,evecfv,evecsv,wfmt,ngkmax,wfgk)
32 deallocate(evecfv,evecsv,apwalm)
33 ! calculate the momentum matrix elements
34 allocate(pmat(nstsv,nstsv,3))
35 call genpmatk(ngk(:,ik),igkig(:,:,ik),vgkc(:,:,:,ik),wfmt,wfgk,pmat)
36 deallocate(wfmt,wfgk)
37 ! write the matrix elements in the second-variational basis
38 !$OMP CRITICAL(u230)
39 ! write to RAM disk if required
40 if (ramdisk) then
41  call putrd('PMAT.OUT',ik,v1=vkl(1:3,ik),n1=nstsv,nzv=nstsv*nstsv*3,zva=pmat)
42 end if
43 ! write to disk if required
44 if (wrtdsk) then
45 ! determine the record length
46  inquire(iolength=recl) vkl(1:3,1),nstsv,pmat
47  open(230,file='PMAT.OUT',form='UNFORMATTED',access='DIRECT',recl=recl)
48  write(230,rec=ik) vkl(1:3,ik),nstsv,pmat
49  close(230)
50 end if
51 !$OMP END CRITICAL(u230)
52 deallocate(pmat)
53 end subroutine
54 
integer nmatmax
Definition: modmain.f90:853
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
integer npcmtmax
Definition: modmain.f90:216
integer, dimension(3) ngridg
Definition: modmain.f90:386
character(256) filext
Definition: modmain.f90:1301
subroutine genwfsv(tsh, tgp, nst, idx, ngridg_, igfft_, ngp, igpig, apwalm, evecfv, evecsv, wfmt, ld, wfir)
Definition: genwfsv.f90:11
integer lmmaxapw
Definition: modmain.f90:199
logical ramdisk
Definition: modramdisk.f90:9
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
subroutine putrd(fname, irec, n1, n2, n3, v1, v2, nrv, rva, nzv, zva)
Definition: modramdisk.f90:120
type(file_t), dimension(:), allocatable, private file
Definition: modramdisk.f90:29
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition: modmain.f90:509
subroutine putpmat(ik)
Definition: putpmat.f90:7
subroutine genpmatk(ngp, igpig, vgpc, wfmt, wfgp, pmat)
Definition: genpmatk.f90:10
integer nstsv
Definition: modmain.f90:889
integer, dimension(:,:), allocatable ngk
Definition: modmain.f90:497
integer, dimension(:), allocatable igfft
Definition: modmain.f90:406
real(8), dimension(:,:,:,:), allocatable vgkl
Definition: modmain.f90:503
logical wrtdsk
Definition: modramdisk.f90:15
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
Definition: modmpi.f90:6
real(8), dimension(:,:,:), allocatable gkc
Definition: modmain.f90:507
integer natmtot
Definition: modmain.f90:40
integer, dimension(:,:,:), allocatable igkig
Definition: modmain.f90:501
integer nstfv
Definition: modmain.f90:887
integer nspnfv
Definition: modmain.f90:289