The Elk Code
 
Loading...
Searching...
No Matches
exxengy.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2006 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine exxengy
7use modmain
8use modmpi
9use modomp
10implicit none
11! local variables
12integer ik,ist,jst,is,ia
13integer nrc,nrci,npc
14integer m1,m2,nthd
15complex(8) z1
16! automatic arrays
17complex(4) wfcr1(npcmtmax,2),wfcr2(npcmtmax,2)
18complex(4) crhomt(npcmtmax),cvclmt(npcmtmax)
19! external functions
20complex(8), external :: zcfmtinp
21! zero the exchange energy
22engyx=0.d0
23!--------------------------------------------------!
24! val-val-val and val-cr-val contributions !
25!--------------------------------------------------!
26call holdthd(nkpt/np_mpi,nthd)
27!$OMP PARALLEL DEFAULT(SHARED) &
28!$OMP NUM_THREADS(nthd)
29!$OMP DO SCHEDULE(DYNAMIC)
30do ik=1,nkpt
31! distribute among MPI processes
32 if (mod(ik-1,np_mpi) /= lp_mpi) cycle
33!$OMP CRITICAL(exxengy_)
34 write(*,'("Info(exxengy): ",I6," of ",I6," k-points")') ik,nkpt
35!$OMP END CRITICAL(exxengy_)
36 call exxengyk(ik)
37end do
38!$OMP END DO
39!$OMP END PARALLEL
40call freethd(nthd)
41! add energies from each process and redistribute
42call mpi_allreduce(mpi_in_place,engyx,1,mpi_double_precision,mpi_sum,mpicom, &
43 ierror)
44!-----------------------------------!
45! core-core-core contribution !
46!-----------------------------------!
47! begin loops over atoms and species
48do is=1,nspecies
49 nrc=nrcmt(is)
50 nrci=nrcmti(is)
51 npc=npcmt(is)
52 do ia=1,natoms(is)
53 do jst=1,nstsp(is)
54 if (spcore(jst,is)) then
55 do m2=-ksp(jst,is),ksp(jst,is)-1
56! generate the core wavefunction in spherical coordinates (pass in m-1/2)
57 call wavefcr(.false.,lradstp,is,ia,jst,m2,npcmtmax,wfcr2)
58 do ist=1,nstsp(is)
59 if (spcore(ist,is)) then
60 do m1=-ksp(ist,is),ksp(ist,is)-1
61 call wavefcr(.false.,lradstp,is,ia,ist,m1,npcmtmax,wfcr1)
62! calculate the complex overlap density
63 call crho2(npc,wfcr1,wfcr1(:,2),wfcr2,wfcr2(:,2),crhomt)
64 call cfshtip(nrc,nrci,crhomt)
65! calculate the Coulomb potential
66 call cpotclmt(nrc,nrci,nrcmtmax,rlcmt(:,:,is),wprcmt(:,:,is), &
67 crhomt,cvclmt)
68 z1=zcfmtinp(nrc,nrci,wr2cmt(:,is),crhomt,cvclmt)
69 engyx=engyx-0.5d0*z1%re
70 end do
71! end loop over ist
72 end if
73 end do
74 end do
75! end loop over jst
76 end if
77 end do
78! end loops over atoms and species
79 end do
80end do
81return
82
83contains
84
85pure subroutine crho2(n,wf11,wf12,wf21,wf22,crho)
86implicit none
87integer, intent(in) :: n
88complex(4), intent(in) :: wf11(n),wf12(n),wf21(n),wf22(n)
89complex(4), intent(out) :: crho(n)
90crho(1:n)=conjg(wf11(1:n))*wf21(1:n)+conjg(wf12(1:n))*wf22(1:n)
91end subroutine
92
93end subroutine
94
subroutine cfshtip(nr, nri, cfmt)
Definition cfshtip.f90:7
pure subroutine cpotclmt(nr, nri, ld, rl, wpr, crhomt, cvclmt)
Definition cpotclmt.f90:7
subroutine exxengy
Definition exxengy.f90:7
pure subroutine crho2(n, wf11, wf12, wf21, wf22, crho)
Definition exxengy.f90:86
subroutine exxengyk(ikp)
Definition exxengyk.f90:7
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
logical, dimension(maxstsp, maxspecies) spcore
Definition modmain.f90:127
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
integer lradstp
Definition modmain.f90:171
integer, dimension(maxspecies) nstsp
Definition modmain.f90:113
real(8), dimension(:,:,:), allocatable rlcmt
Definition modmain.f90:181
integer nkpt
Definition modmain.f90:461
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
real(8) engyx
Definition modmain.f90:972
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
integer lp_mpi
Definition modmpi.f90:15
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
integer np_mpi
Definition modmpi.f90:13
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
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