The Elk Code
hmlrad.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2016 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 
6 !BOP
7 ! !ROUTINE: hmlrad
8 ! !INTERFACE:
9 subroutine hmlrad
10 ! !USES:
11 use modmain
12 use modomp
13 ! !DESCRIPTION:
14 ! Calculates the radial Hamiltonian integrals of the APW and local-orbital
15 ! basis functions. In other words, for atom $\alpha$, it computes integrals of
16 ! the form
17 ! $$ h^{\alpha}_{qq';ll'l''m''}=\begin{cases}
18 ! \int_0^{R_i}u^{\alpha}_{q;l}(r)H u^{\alpha}_{q';l'}(r)r^2dr & l''=0 \\
19 ! \int_0^{R_i}u^{\alpha}_{q;l}(r)V^{\alpha}_{l''m''}(r)
20 ! u^{\alpha}_{q';l'}(r)r^2dr & l''>0 \end{cases} $$
21 ! where $u^{\alpha}_{q;l}$ is the $q$th APW radial function for angular
22 ! momentum $l$; $H$ is the Hamiltonian of the radial Schr\"{o}dinger equation;
23 ! and $V^{\alpha}_{l''m''}$ is the muffin-tin Kohn-Sham potential. Similar
24 ! integrals are calculated for APW-local-orbital and
25 ! local-orbital-local-orbital contributions.
26 !
27 ! !REVISION HISTORY:
28 ! Created December 2003 (JKD)
29 ! Updated for compressed muffin-tin functions, March 2016 (JKD)
30 !EOP
31 !BOC
32 implicit none
33 ! local variables
34 integer is,ias,nthd
35 integer nr,nri,iro,i0,i1
36 integer l1,l2,l3,lm2
37 integer io,jo,ilo,jlo
38 real(8) sm
39 ! begin loops over atoms and species
40 call holdthd(natmtot,nthd)
41 !$OMP PARALLEL DO DEFAULT(SHARED) &
42 !$OMP PRIVATE(is,nr,nri,iro) &
43 !$OMP PRIVATE(l1,l2,l3,io,jo,sm) &
44 !$OMP PRIVATE(lm2,i0,i1,ilo,jlo) &
45 !$OMP SCHEDULE(DYNAMIC) &
46 !$OMP NUM_THREADS(nthd)
47 do ias=1,natmtot
48  is=idxis(ias)
49  nr=nrmt(is)
50  nri=nrmti(is)
51  iro=nri+1
52 !---------------------------!
53 ! APW-APW integrals !
54 !---------------------------!
55  do l1=0,lmaxapw
56  do io=1,apword(l1,is)
57  do l3=0,lmaxapw
58  do jo=1,apword(l3,is)
59  if (l1 == l3) then
60  sm=sum((apwfr(1:nr,1,io,l1,ias)*apwfr(1:nr,2,jo,l1,ias) &
61  +apwfr(1:nr,2,io,l1,ias)*apwfr(1:nr,1,jo,l1,ias)) &
62  *wr2mt(1:nr,is))
63 ! add kinetic surface term
64  sm=sm+apwfr(nr,1,io,l1,ias)*apwdfr(jo,l1,ias) &
65  +apwdfr(io,l1,ias)*apwfr(nr,1,jo,l1,ias)
66  haa(1,jo,l1,io,l1,ias)=sm*y00i/2.d0
67  else
68  haa(1,jo,l3,io,l1,ias)=0.d0
69  end if
70  if (l1 >= l3) then
71  do l2=1,lmaxi
72  do lm2=l2**2+1,(l2+1)**2
73  i1=lmmaxi*(nri-1)+lm2
74  sm=sum(apwfr(1:nri,1,io,l1,ias)*apwfr(1:nri,1,jo,l3,ias) &
75  *wr2mt(1:nri,is)*vsmt(lm2:i1:lmmaxi,ias))
76  i0=i1+lmmaxi
77  i1=lmmaxo*(nr-iro)+i0
78  sm=sm+sum(apwfr(iro:nr,1,io,l1,ias)*apwfr(iro:nr,1,jo,l3,ias) &
79  *wr2mt(iro:nr,is)*vsmt(i0:i1:lmmaxo,ias))
80  haa(lm2,jo,l3,io,l1,ias)=sm
81  haa(lm2,io,l1,jo,l3,ias)=sm
82  end do
83  end do
84  do l2=lmaxi+1,lmaxo
85  do lm2=l2**2+1,(l2+1)**2
86  i0=lmmaxi*nri+lm2
87  i1=lmmaxo*(nr-iro)+i0
88  sm=sum(apwfr(iro:nr,1,io,l1,ias)*apwfr(iro:nr,1,jo,l3,ias) &
89  *wr2mt(iro:nr,is)*vsmt(i0:i1:lmmaxo,ias))
90  haa(lm2,jo,l3,io,l1,ias)=sm
91  haa(lm2,io,l1,jo,l3,ias)=sm
92  end do
93  end do
94  end if
95  end do
96  end do
97  end do
98  end do
99 !-------------------------------------!
100 ! local-orbital-APW integrals !
101 !-------------------------------------!
102  do ilo=1,nlorb(is)
103  l1=lorbl(ilo,is)
104  do l3=0,lmaxapw
105  do jo=1,apword(l3,is)
106  if (l1 == l3) then
107  sm=sum(lofr(1:nr,1,ilo,ias)*apwfr(1:nr,2,jo,l1,ias)*wr2mt(1:nr,is))
108  hloa(1,jo,l1,ilo,ias)=sm*y00i
109  else
110  hloa(1,jo,l3,ilo,ias)=0.d0
111  end if
112  do l2=1,lmaxi
113  do lm2=l2**2+1,(l2+1)**2
114  i1=lmmaxi*(nri-1)+lm2
115  sm=sum(lofr(1:nri,1,ilo,ias)*apwfr(1:nri,1,jo,l3,ias) &
116  *wr2mt(1:nri,is)*vsmt(lm2:i1:lmmaxi,ias))
117  i0=i1+lmmaxi
118  i1=lmmaxo*(nr-iro)+i0
119  sm=sm+sum(lofr(iro:nr,1,ilo,ias)*apwfr(iro:nr,1,jo,l3,ias) &
120  *wr2mt(iro:nr,is)*vsmt(i0:i1:lmmaxo,ias))
121  hloa(lm2,jo,l3,ilo,ias)=sm
122  end do
123  end do
124  do l2=lmaxi+1,lmaxo
125  do lm2=l2**2+1,(l2+1)**2
126  i0=lmmaxi*nri+lm2
127  i1=lmmaxo*(nr-iro)+i0
128  sm=sum(lofr(iro:nr,1,ilo,ias)*apwfr(iro:nr,1,jo,l3,ias) &
129  *wr2mt(iro:nr,is)*vsmt(i0:i1:lmmaxo,ias))
130  hloa(lm2,jo,l3,ilo,ias)=sm
131  end do
132  end do
133  end do
134  end do
135  end do
136 !-----------------------------------------------!
137 ! local-orbital-local-orbital integrals !
138 !-----------------------------------------------!
139  do ilo=1,nlorb(is)
140  l1=lorbl(ilo,is)
141  do jlo=1,nlorb(is)
142  l3=lorbl(jlo,is)
143  if (l1 == l3) then
144  sm=sum(lofr(1:nr,1,ilo,ias)*lofr(1:nr,2,jlo,ias)*wr2mt(1:nr,is))
145  hlolo(1,jlo,ilo,ias)=sm*y00i
146  else
147  hlolo(1,jlo,ilo,ias)=0.d0
148  end if
149  do l2=1,lmaxi
150  do lm2=l2**2+1,(l2+1)**2
151  i1=lmmaxi*(nri-1)+lm2
152  sm=sum(lofr(1:nri,1,ilo,ias)*lofr(1:nri,1,jlo,ias)*wr2mt(1:nri,is) &
153  *vsmt(lm2:i1:lmmaxi,ias))
154  i0=i1+lmmaxi
155  i1=lmmaxo*(nr-iro)+i0
156  sm=sm+sum(lofr(iro:nr,1,ilo,ias)*lofr(iro:nr,1,jlo,ias) &
157  *wr2mt(iro:nr,is)*vsmt(i0:i1:lmmaxo,ias))
158  hlolo(lm2,jlo,ilo,ias)=sm
159  end do
160  end do
161  do l2=lmaxi+1,lmaxo
162  do lm2=l2**2+1,(l2+1)**2
163  i0=lmmaxi*nri+lm2
164  i1=lmmaxo*(nr-iro)+i0
165  sm=sum(lofr(iro:nr,1,ilo,ias)*lofr(iro:nr,1,jlo,ias) &
166  *wr2mt(iro:nr,is)*vsmt(i0:i1:lmmaxo,ias))
167  hlolo(lm2,jlo,ilo,ias)=sm
168  end do
169  end do
170  end do
171  end do
172 ! end loops over atoms and species
173 end do
174 !$OMP END PARALLEL DO
175 call freethd(nthd)
176 end subroutine
177 !EOC
178 
subroutine hmlrad
Definition: hmlrad.f90:10
integer, dimension(maxspecies) nlorb
Definition: modmain.f90:786
real(8), dimension(:,:,:,:), allocatable lofr
Definition: modmain.f90:814
integer lmmaxo
Definition: modmain.f90:203
Definition: modomp.f90:6
integer lmaxo
Definition: modmain.f90:201
integer lmaxapw
Definition: modmain.f90:197
real(8), dimension(:,:,:,:,:,:), allocatable haa
Definition: modmain.f90:861
integer, dimension(0:maxlapw, maxspecies) apword
Definition: modmain.f90:758
real(8), dimension(:,:,:,:,:), allocatable hloa
Definition: modmain.f90:863
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer lmmaxi
Definition: modmain.f90:207
real(8), parameter y00i
Definition: modmain.f90:1237
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer natmtot
Definition: modmain.f90:40
real(8), dimension(:,:,:,:,:), allocatable apwfr
Definition: modmain.f90:774
integer, dimension(maxlorb, maxspecies) lorbl
Definition: modmain.f90:796
real(8), dimension(:,:,:,:), allocatable hlolo
Definition: modmain.f90:865
real(8), dimension(:,:), pointer, contiguous vsmt
Definition: modmain.f90:649
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
real(8), dimension(:,:,:), allocatable apwdfr
Definition: modmain.f90:778
real(8), dimension(:,:), allocatable wr2mt
Definition: modmain.f90:183
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150
integer lmaxi
Definition: modmain.f90:205