The Elk Code
 
Loading...
Searching...
No Matches
vclcore.f90
Go to the documentation of this file.
1
2! Copyright (C) 2012 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 vclcore(wfmt,vmat)
7use modmain
8use modomp
9implicit none
10! arguments
11complex(4), intent(in) :: wfmt(npcmtmax,natmtot,nspinor,nstsv)
12complex(8), intent(inout) :: vmat(nstsv,nstsv)
13! local variables
14integer ist1,ist2,ist3
15integer is,ia,ias,m,nthd
16integer nrc,nrci,npc
17! automatic arrays
18complex(4) wfcr(npcmtmax,2),cfmt(npcmtmax)
19complex(8) v(nstsv,nstsv)
20! allocatable arrays
21complex(4), allocatable :: crhomt(:,:)
22! external functions
23complex(8), external :: zcfmtinp
24allocate(crhomt(npcmtmax,nstsv))
25! zero the upper triangular matrix elements
26do ist2=1,nstsv
27 v(1:ist2,ist2)=0.d0
28end do
29call holdthd(nstsv,nthd)
30!$OMP PARALLEL DEFAULT(SHARED) &
31!$OMP PRIVATE(cfmt,is,ia,ias) &
32!$OMP PRIVATE(nrc,nrci,npc) &
33!$OMP PRIVATE(ist1,ist2,ist3,m) &
34!$OMP NUM_THREADS(nthd)
35do is=1,nspecies
36 nrc=nrcmt(is)
37 nrci=nrcmti(is)
38 npc=npcmt(is)
39 do ia=1,natoms(is)
40 ias=idxas(ia,is)
41 do ist3=1,nstsp(is)
42 if (spcore(ist3,is)) then
43 do m=-ksp(ist3,is),ksp(ist3,is)-1
44! generate the core wavefunction in spherical coordinates (pass in m-1/2)
45!$OMP SINGLE
46 call wavefcr(.false.,lradstp,is,ia,ist3,m,npcmtmax,wfcr)
47!$OMP END SINGLE
48!$OMP DO SCHEDULE(DYNAMIC)
49 do ist1=1,nstsv
50! calculate the complex overlap density in spherical harmonics
51 if (spinpol) then
52 call crho2(npc,wfcr,wfcr(:,2),wfmt(:,ias,1,ist1), &
53 wfmt(:,ias,2,ist1),cfmt)
54 else
55 call crho1(npc,wfcr,wfmt(:,ias,1,ist1),cfmt)
56 end if
57 call cfsht(nrc,nrci,cfmt,crhomt(:,ist1))
58 end do
59!$OMP END DO
60!$OMP DO SCHEDULE(DYNAMIC)
61 do ist2=1,nstsv
62 call cpotclmt(nrc,nrci,nrcmtmax,rlcmt(:,:,is),wprcmt(:,:,is), &
63 crhomt(:,ist2),cfmt)
64 do ist1=1,ist2
65 v(ist1,ist2)=v(ist1,ist2)-zcfmtinp(nrc,nrci,wr2cmt(:,is), &
66 crhomt(:,ist1),cfmt)
67 end do
68 end do
69!$OMP END DO
70 end do
71 end if
72 end do
73 end do
74end do
75!$OMP END PARALLEL
76call freethd(nthd)
77do ist1=1,nstsv
78! set the lower triangular part of the matrix
79 do ist2=1,ist1-1
80 v(ist1,ist2)=conjg(v(ist2,ist1))
81 end do
82! make the diagonal elements real
83 v(ist1,ist1)=dble(v(ist1,ist1))
84end do
85! scale the Coulomb matrix elements in the case of a hybrid functional
86if (hybrid) v(1:nstsv,1:nstsv)=hybridc*v(1:nstsv,1:nstsv)
87! add to input matrix
88vmat(1:nstsv,1:nstsv)=vmat(1:nstsv,1:nstsv)+v(1:nstsv,1:nstsv)
89deallocate(crhomt)
90return
91
92contains
93
94pure subroutine crho1(n,wf1,wf2,crho)
95implicit none
96integer, intent(in) :: n
97complex(4), intent(in) :: wf1(n),wf2(n)
98complex(4), intent(out) :: crho(n)
99crho(1:n)=conjg(wf1(1:n))*wf2(1:n)
100end subroutine
101
102pure subroutine crho2(n,wf11,wf12,wf21,wf22,crho)
103implicit none
104integer, intent(in) :: n
105complex(4), intent(in) :: wf11(n),wf12(n),wf21(n),wf22(n)
106complex(4), intent(out) :: crho(n)
107crho(1:n)=conjg(wf11(1:n))*wf21(1:n)+conjg(wf12(1:n))*wf22(1:n)
108end subroutine
109
110end subroutine
111
subroutine cfsht(nr, nri, cfmt1, cfmt2)
Definition cfsht.f90:7
pure subroutine cpotclmt(nr, nri, ld, rl, wpr, crhomt, cvclmt)
Definition cpotclmt.f90:7
pure subroutine crho2(n, wf11, wf12, wf21, wf22, crho)
Definition exxengy.f90:86
pure subroutine crho1(n, wf1, wf2, crho)
Definition exxengyk.f90:174
real(8) hybridc
Definition modmain.f90:1151
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
logical, dimension(maxstsp, maxspecies) spcore
Definition modmain.f90:127
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
logical spinpol
Definition modmain.f90:228
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
logical hybrid
Definition modmain.f90:1149
integer lradstp
Definition modmain.f90:171
integer, dimension(maxspecies) nstsp
Definition modmain.f90:113
real(8), dimension(:,:,:), allocatable rlcmt
Definition modmain.f90:181
real(8), dimension(:,:), allocatable wr2cmt
Definition modmain.f90:189
real(8), dimension(:,:,:), allocatable wprcmt
Definition modmain.f90:191
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
integer, dimension(maxstsp, maxspecies) ksp
Definition modmain.f90:125
integer nrcmtmax
Definition modmain.f90:175
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
integer nspecies
Definition modmain.f90:34
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine vclcore(wfmt, vmat)
Definition vclcore.f90:7
pure subroutine wavefcr(tsh, lrstp, is, ia, ist, m, ld, wfcr)
Definition wavefcr.f90:7
pure complex(8) function zcfmtinp(nr, nri, wr, cfmt1, cfmt2)
Definition zcfmtinp.f90:7