The Elk Code
 
Loading...
Searching...
No Matches
dhmlrad.f90
Go to the documentation of this file.
1
2! Copyright (C) 2013 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 dhmlrad
7use modmain
8use modphonon
9implicit none
10! local variables
11integer is,ias
12integer nr,nri,iro,i0,i1
13integer l1,l2,l3,lm2
14integer io,jo,ilo,jlo
15complex(8) zsm
16! begin loops over atoms and species
17do ias=1,natmtot
18 is=idxis(ias)
19 nr=nrmt(is)
20 nri=nrmti(is)
21 iro=nri+1
22!---------------------------!
23! APW-APW integrals !
24!---------------------------!
25 do l1=0,lmaxapw
26 do io=1,apword(l1,is)
27 do l3=0,lmaxapw
28 do jo=1,apword(l3,is)
29 do l2=0,lmaxi
30 do lm2=l2**2+1,(l2+1)**2
31 i1=lmmaxi*(nri-1)+lm2
32 zsm=sum(apwfr(1:nri,1,io,l1,ias)*apwfr(1:nri,1,jo,l3,ias) &
33 *wr2mt(1:nri,is)*dvsmt(lm2:i1:lmmaxi,ias))
34 i0=i1+lmmaxi
35 i1=lmmaxo*(nr-iro)+i0
36 zsm=zsm+sum(apwfr(iro:nr,1,io,l1,ias)*apwfr(iro:nr,1,jo,l3,ias) &
37 *wr2mt(iro:nr,is)*dvsmt(i0:i1:lmmaxo,ias))
38 dhaa(lm2,jo,l3,io,l1,ias)=zsm
39 end do
40 end do
41 do l2=lmaxi+1,lmaxo
42 do lm2=l2**2+1,(l2+1)**2
43 i0=lmmaxi*nri+lm2
44 i1=lmmaxo*(nr-iro)+i0
45 zsm=sum(apwfr(iro:nr,1,io,l1,ias)*apwfr(iro:nr,1,jo,l3,ias) &
46 *wr2mt(iro:nr,is)*dvsmt(i0:i1:lmmaxo,ias))
47 dhaa(lm2,jo,l3,io,l1,ias)=zsm
48 end do
49 end do
50 end do
51 end do
52 end do
53 end do
54!-------------------------------------!
55! local-orbital-APW integrals !
56!-------------------------------------!
57 do ilo=1,nlorb(is)
58 l1=lorbl(ilo,is)
59 do l3=0,lmaxapw
60 do io=1,apword(l3,is)
61 do l2=0,lmaxi
62 do lm2=l2**2+1,(l2+1)**2
63 i1=lmmaxi*(nri-1)+lm2
64 zsm=sum(lofr(1:nri,1,ilo,ias)*apwfr(1:nri,1,io,l3,ias) &
65 *wr2mt(1:nri,is)*dvsmt(lm2:i1:lmmaxi,ias))
66 i0=i1+lmmaxi
67 i1=lmmaxo*(nr-iro)+i0
68 zsm=zsm+sum(lofr(iro:nr,1,ilo,ias)*apwfr(iro:nr,1,io,l3,ias) &
69 *wr2mt(iro:nr,is)*dvsmt(i0:i1:lmmaxo,ias))
70 dhloa(lm2,io,l3,ilo,ias)=zsm
71 end do
72 end do
73 do l2=lmaxi+1,lmaxo
74 do lm2=l2**2+1,(l2+1)**2
75 i0=lmmaxi*nri+lm2
76 i1=lmmaxo*(nr-iro)+i0
77 zsm=sum(lofr(iro:nr,1,ilo,ias)*apwfr(iro:nr,1,io,l3,ias) &
78 *wr2mt(iro:nr,is)*dvsmt(i0:i1:lmmaxo,ias))
79 dhloa(lm2,io,l3,ilo,ias)=zsm
80 end do
81 end do
82 end do
83 end do
84 end do
85!-----------------------------------------------!
86! local-orbital-local-orbital integrals !
87!-----------------------------------------------!
88 do ilo=1,nlorb(is)
89 l1=lorbl(ilo,is)
90 do jlo=1,nlorb(is)
91 l3=lorbl(jlo,is)
92 do l2=0,lmaxi
93 do lm2=l2**2+1,(l2+1)**2
94 i1=lmmaxi*(nri-1)+lm2
95 zsm=sum(lofr(1:nri,1,ilo,ias)*lofr(1:nri,1,jlo,ias)*wr2mt(1:nri,is) &
96 *dvsmt(lm2:i1:lmmaxi,ias))
97 i0=i1+lmmaxi
98 i1=lmmaxo*(nr-iro)+i0
99 zsm=zsm+sum(lofr(iro:nr,1,ilo,ias)*lofr(iro:nr,1,jlo,ias) &
100 *wr2mt(iro:nr,is)*dvsmt(i0:i1:lmmaxo,ias))
101 dhlolo(lm2,jlo,ilo,ias)=zsm
102 end do
103 end do
104 do l2=lmaxi+1,lmaxo
105 do lm2=l2**2+1,(l2+1)**2
106 i0=lmmaxi*nri+lm2
107 i1=lmmaxo*(nr-iro)+i0
108 zsm=sum(lofr(iro:nr,1,ilo,ias)*lofr(iro:nr,1,jlo,ias) &
109 *wr2mt(iro:nr,is)*dvsmt(i0:i1:lmmaxo,ias))
110 dhlolo(lm2,jlo,ilo,ias)=zsm
111 end do
112 end do
113 end do
114 end do
115! end loops over atoms and species
116end do
117end subroutine
118
subroutine dhmlrad
Definition dhmlrad.f90:7
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
integer lmmaxi
Definition modmain.f90:207
integer, dimension(0:maxlapw, maxspecies) apword
Definition modmain.f90:758
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer natmtot
Definition modmain.f90:40
integer lmaxapw
Definition modmain.f90:197
real(8), dimension(:,:,:,:,:), allocatable apwfr
Definition modmain.f90:774
integer lmaxo
Definition modmain.f90:201
real(8), dimension(:,:), allocatable wr2mt
Definition modmain.f90:183
integer lmaxi
Definition modmain.f90:205
integer, dimension(maxspecies) nlorb
Definition modmain.f90:786
integer lmmaxo
Definition modmain.f90:203
real(8), dimension(:,:,:,:), allocatable lofr
Definition modmain.f90:814
integer, dimension(maxlorb, maxspecies) lorbl
Definition modmain.f90:796
complex(8), dimension(:,:,:,:), allocatable dhlolo
complex(8), dimension(:,:), pointer, contiguous dvsmt
complex(8), dimension(:,:,:,:,:), allocatable dhloa
complex(8), dimension(:,:,:,:,:,:), allocatable dhaa