File: putkmat.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (100 lines) | stat: -rw-r--r-- 3,171 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100

! Copyright (C) 2007-2010 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

subroutine putkmat(tfv,tvclcr,ik,vmt,vir)
use modmain
use modmpi
implicit none
! arguments
logical, intent(in) :: tfv,tvclcr
integer, intent(in) :: ik
real(8), intent(in) :: vmt(npcmtmax,natmtot),vir(ngtot)
! local variables
integer ist,ispn,recl,i
! automatic arrays
integer idx(nstsv)
! allocatable arrays
complex(8), allocatable :: apwalm(:,:,:,:,:),evecfv(:,:,:),evecsv(:,:)
complex(8), allocatable :: wfmt(:,:,:,:),wfgk(:,:,:)
complex(8), allocatable :: kmat(:,:),a(:,:)
! get the eigenvalues/vectors from file for input reduced k-point
allocate(evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv))
call getevalsv(filext,ik,vkl(:,ik),evalsv(:,ik))
call getevecfv(filext,ik,vkl(:,ik),vgkl(:,:,:,ik),evecfv)
call getevecsv(filext,ik,vkl(:,ik),evecsv)
! find the matching coefficients
allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv))
do ispn=1,nspnfv
  call match(ngk(ispn,ik),gkc(:,ispn,ik),tpgkc(:,:,ispn,ik), &
   sfacgk(:,:,ispn,ik),apwalm(:,:,:,:,ispn))
end do
! index to all states
do ist=1,nstsv
  idx(ist)=ist
end do
! calculate the wavefunctions for all states of the input k-point
allocate(wfmt(npcmtmax,natmtot,nspinor,nstsv),wfgk(ngkmax,nspinor,nstsv))
call genwfsv(.false.,.true.,nstsv,idx,ngridg,igfft,ngk(:,ik),igkig(:,:,ik), &
 apwalm,evecfv,evecsv,wfmt,ngkmax,wfgk)
deallocate(apwalm,evecfv)
! compute Kohn-Sham potential matrix elements
allocate(kmat(nstsv,nstsv))
if (spinpol) then
  call genvbmatk(vmt,vir,bsmt,bsir,ngk(:,ik),igkig(:,:,ik),wfmt,ngkmax,wfgk, &
   kmat)
else
  call genvmatk(vmt,vir,ngk(:,ik),igkig(:,:,ik),wfmt,ngkmax,wfgk,kmat)
end if
deallocate(wfgk)
! add the DFT+U matrix elements if required
call genvmmtsv(wfmt,kmat)
! negate the potential matrix elements because we have to subtract them
kmat(:,:)=-kmat(:,:)
! add second-variational eigenvalues along the diagonal
do ist=1,nstsv
  kmat(ist,ist)=kmat(ist,ist)+evalsv(ist,ik)
end do
! add scissor correction if required
if (scissor.ne.0.d0) then
  do ist=1,nstsv
    if (evalsv(ist,ik).gt.efermi) then
      kmat(ist,ist)=kmat(ist,ist)+scissor
    end if
  end do
end if
! add the Coulomb core matrix elements if required
if (tvclcr) call vclcore(wfmt,kmat)
! rotate kinetic matrix elements to first-variational basis if required
if (tfv) then
  allocate(a(nstsv,nstsv))
  call zgemm('N','C',nstsv,nstsv,nstsv,zone,kmat,nstsv,evecsv,nstsv,zzero,a, &
   nstsv)
  call zgemm('N','N',nstsv,nstsv,nstsv,zone,evecsv,nstsv,a,nstsv,zzero,kmat, &
   nstsv)
  deallocate(a)
end if
deallocate(evecsv,wfmt)
! determine the record length
inquire(iolength=recl) vkl(:,1),nstsv,kmat
!$OMP CRITICAL(u140)
do i=1,2
  open(140,file='KMAT.OUT',form='UNFORMATTED',access='DIRECT',recl=recl,err=10)
  write(140,rec=ik,err=10) vkl(:,ik),nstsv,kmat
  close(140)
  exit
10 continue
  if (i.eq.2) then
    write(*,*)
    write(*,'("Error(putkmat): unable to write to KMAT.OUT")')
    write(*,*)
    stop
  end if
  close(140)
end do
!$OMP END CRITICAL(u140)
deallocate(kmat)
return
end subroutine