The Elk Code
putkmat.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007-2010 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 putkmat(tfv,tvclcr,ik,vmt,vir,bmt,bir)
7 use modmain
8 use modmpi
9 use modramdisk
10 implicit none
11 ! arguments
12 logical, intent(in) :: tfv,tvclcr
13 integer, intent(in) :: ik
14 real(8), intent(in) :: vmt(npcmtmax,natmtot),vir(ngtc)
15 real(8), intent(in) :: bmt(npcmtmax,natmtot,ndmag),bir(ngtc,ndmag)
16 ! local variables
17 integer jst,ispn,recl
18 ! automatic arrays
19 complex(8) kmat(nstsv,nstsv),a(nstsv,nstsv)
20 ! allocatable arrays
21 complex(8), allocatable :: apwalm(:,:,:,:,:),evecfv(:,:,:),evecsv(:,:)
22 complex(4), allocatable :: wfmt(:,:,:,:),wfgk(:,:,:)
23 ! get the eigenvalues/vectors from file for input reduced k-point
24 allocate(evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv))
25 call getevalsv(filext,ik,vkl(:,ik),evalsv(:,ik))
26 call getevecfv(filext,ik,vkl(:,ik),vgkl(:,:,:,ik),evecfv)
27 call getevecsv(filext,ik,vkl(:,ik),evecsv)
28 ! find the matching coefficients
29 allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv))
30 do ispn=1,nspnfv
31  call match(ngk(ispn,ik),vgkc(:,:,ispn,ik),gkc(:,ispn,ik), &
32  sfacgk(:,:,ispn,ik),apwalm(:,:,:,:,ispn))
33 end do
34 ! calculate the wavefunctions for all states of the input k-point
35 allocate(wfmt(npcmtmax,natmtot,nspinor,nstsv),wfgk(ngkmax,nspinor,nstsv))
36 call genwfsv_sp(.false.,.true.,nstsv,[0],ngridg,igfft,ngk(:,ik),igkig(:,:,ik), &
37  apwalm,evecfv,evecsv,wfmt,ngkmax,wfgk)
38 deallocate(apwalm,evecfv)
39 ! compute Kohn-Sham potential matrix elements
40 if (spinpol) then
41  call genvbmatk(vmt,vir,bmt,bir,ngk(:,ik),igkig(:,:,ik),wfmt,ngkmax,wfgk,kmat)
42 else
43  call genvmatk(vmt,vir,ngk(:,ik),igkig(:,:,ik),wfmt,ngkmax,wfgk,kmat)
44 end if
45 deallocate(wfgk)
46 ! negate the potential matrix elements because we have to subtract them and add
47 ! the second-variational eigenvalues along the diagonal
48 do jst=1,nstsv
49  kmat(1:nstsv,jst)=-kmat(1:nstsv,jst)
50  kmat(jst,jst)=kmat(jst,jst)+evalsv(jst,ik)
51 end do
52 ! add the Coulomb core matrix elements if required
53 if (tvclcr) call vclcore(wfmt,kmat)
54 ! rotate kinetic matrix elements to first-variational basis if required
55 if (tfv) then
56  call zgemm('N','C',nstsv,nstsv,nstsv,zone,kmat,nstsv,evecsv,nstsv,zzero,a, &
57  nstsv)
58  call zgemm('N','N',nstsv,nstsv,nstsv,zone,evecsv,nstsv,a,nstsv,zzero,kmat, &
59  nstsv)
60 end if
61 deallocate(evecsv,wfmt)
62 !$OMP CRITICAL(u220)
63 ! write to RAM disk if required
64 if (ramdisk) then
65  call putrd('KMAT.OUT',ik,v1=vkl(1:3,ik),n1=nstsv,nzv=nstsv*nstsv,zva=kmat)
66 end if
67 ! write to disk if required
68 if (wrtdsk) then
69 ! determine the record length
70  inquire(iolength=recl) vkl(1:3,1),nstsv,kmat
71  open(220,file='KMAT.OUT',form='UNFORMATTED',access='DIRECT',recl=recl)
72  write(220,rec=ik) vkl(1:3,ik),nstsv,kmat
73  close(220)
74 end if
75 !$OMP END CRITICAL(u220)
76 end subroutine
77 
integer nmatmax
Definition: modmain.f90:853
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
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:921
logical spinpol
Definition: modmain.f90:228
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 getevalsv(fext, ikp, vpl, evalsv_)
Definition: getevalsv.f90:7
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), parameter zone
Definition: modmain.f90:1240
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition: modmain.f90:509
subroutine vclcore(wfmt, vmat)
Definition: vclcore.f90:7
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
logical wrtdsk
Definition: modramdisk.f90:15
integer nspinor
Definition: modmain.f90:267
real(8), dimension(:,:,:,:), allocatable vgkc
Definition: modmain.f90:505
complex(8), parameter zzero
Definition: modmain.f90:1240
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
subroutine genvbmatk(vmt, vir, bmt, bir, ngp, igpig, wfmt, ld, wfgp, vbmat)
Definition: genvbmatk.f90:7
integer, dimension(:,:,:), allocatable igkig
Definition: modmain.f90:501
integer nstfv
Definition: modmain.f90:887
subroutine putkmat(tfv, tvclcr, ik, vmt, vir, bmt, bir)
Definition: putkmat.f90:7
integer nspnfv
Definition: modmain.f90:289