The Elk Code
olplolo.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 pure subroutine olplolo(is,ias,ngp,ld,o)
7 use modmain
8 implicit none
9 ! arguments
10 integer, intent(in) :: is,ias,ngp,ld
11 complex(8), intent(out) :: o(ld,*)
12 ! local variables
13 integer ilo,jlo,l,lm,i,j
14 do jlo=1,nlorb(is)
15  l=lorbl(jlo,is)
16  do ilo=1,jlo
17  if (lorbl(ilo,is) /= l) cycle
18  do lm=l**2+1,(l+1)**2
19  i=ngp+idxlo(lm,ilo,ias)
20  j=ngp+idxlo(lm,jlo,ias)
21  o(i,j)=ololo(ilo,jlo,ias)
22  end do
23  end do
24 end do
25 end subroutine
26 
integer, dimension(maxspecies) nlorb
Definition: modmain.f90:786
integer, dimension(:,:,:), allocatable idxlo
Definition: modmain.f90:855
real(8), dimension(:,:,:), allocatable ololo
Definition: modmain.f90:859
pure subroutine olplolo(is, ias, ngp, ld, o)
Definition: olplolo.f90:7
integer, dimension(maxlorb, maxspecies) lorbl
Definition: modmain.f90:796