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