The Elk Code
olprad.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 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: olprad
8 ! !INTERFACE:
9 subroutine olprad
10 ! !USES:
11 use modmain
12 use modomp
13 ! !DESCRIPTION:
14 ! Calculates the radial overlap integrals of the APW and local-orbital basis
15 ! functions. In other words, for atom $\alpha$, it computes integrals of the
16 ! form
17 ! $$ o^{\alpha}_{qp}=\int_0^{R_i}u^{\alpha}_{q;l_p}(r)v^{\alpha}_p(r)r^2dr $$
18 ! and
19 ! $$ o^{\alpha}_{pp'}=\int_0^{R_i}v^{\alpha}_p(r)v^{\alpha}_{p'}(r)r^2dr,
20 ! \quad l_p=l_{p'} $$
21 ! where $u^{\alpha}_{q;l}$ is the $q$th APW radial function for angular
22 ! momentum $l$; and $v^{\alpha}_p$ is the $p$th local-orbital radial function
23 ! and has angular momentum $l_p$.
24 !
25 ! !REVISION HISTORY:
26 ! Created November 2003 (JKD)
27 !EOP
28 !BOC
29 implicit none
30 ! local variables
31 integer is,ias,nr,nthd
32 integer ilo,jlo,l,io
33 ! loop over atoms
34 call holdthd(natmtot,nthd)
35 !$OMP PARALLEL DO DEFAULT(SHARED) &
36 !$OMP PRIVATE(is,nr,ilo,jlo,l,io) &
37 !$OMP SCHEDULE(DYNAMIC) &
38 !$OMP NUM_THREADS(nthd)
39 do ias=1,natmtot
40  is=idxis(ias)
41  nr=nrmt(is)
42 ! loop over local-orbitals
43  do ilo=1,nlorb(is)
44  l=lorbl(ilo,is)
45 !-------------------------------------!
46 ! APW-local-orbital integrals !
47 !-------------------------------------!
48  do io=1,apword(l,is)
49  oalo(io,ilo,ias)=sum(apwfr(1:nr,1,io,l,ias)*lofr(1:nr,1,ilo,ias) &
50  *wr2mt(1:nr,is))
51  end do
52 !-----------------------------------------------!
53 ! local-orbital-local-orbital integrals !
54 !-----------------------------------------------!
55  do jlo=1,nlorb(is)
56  if (lorbl(jlo,is) == l) then
57  ololo(ilo,jlo,ias)=sum(lofr(1:nr,1,ilo,ias)*lofr(1:nr,1,jlo,ias) &
58  *wr2mt(1:nr,is))
59  end if
60  end do
61  end do
62 end do
63 !$OMP END PARALLEL DO
64 call freethd(nthd)
65 end subroutine
66 !EOC
67 
integer, dimension(maxspecies) nlorb
Definition: modmain.f90:786
real(8), dimension(:,:,:,:), allocatable lofr
Definition: modmain.f90:814
real(8), dimension(:,:,:), allocatable oalo
Definition: modmain.f90:857
Definition: modomp.f90:6
real(8), dimension(:,:,:), allocatable ololo
Definition: modmain.f90:859
integer, dimension(0:maxlapw, maxspecies) apword
Definition: modmain.f90:758
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
subroutine olprad
Definition: olprad.f90:10
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 wr2mt
Definition: modmain.f90:183
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150