The Elk Code
olpistl.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: olpistl
8 ! !INTERFACE:
9 pure subroutine olpistl(ngp,igpig,ld,o)
10 ! !USES:
11 use modmain
12 ! !INPUT/OUTPUT PARAMETERS:
13 ! ngp : number of G+p-vectors (in,integer)
14 ! igpig : index from G+p-vectors to G-vectors (in,integer(ngkmax))
15 ! ld : leading dimension of o (in,integer)
16 ! o : overlap matrix (inout,complex(*))
17 ! !DESCRIPTION:
18 ! Computes the interstitial contribution to the overlap matrix for the APW
19 ! basis functions. The overlap is given by
20 ! $$ O^{\rm I}({\bf G+k,G'+k})=\tilde{\Theta}({\bf G-G'}), $$
21 ! where $\tilde{\Theta}$ is the characteristic function. See routine
22 ! {\tt gencfun}.
23 !
24 ! !REVISION HISTORY:
25 ! Created April 2003 (JKD)
26 !EOP
27 !BOC
28 implicit none
29 ! arguments
30 integer, intent(in) :: ngp,igpig(ngkmax),ld
31 complex(8), intent(out) :: o(ld,*)
32 ! local variables
33 integer ig,j1,j2,j3,i,j
34 do j=1,ngp
35  ig=igpig(j)
36  j1=ivg(1,ig); j2=ivg(2,ig); j3=ivg(3,ig)
37  do i=1,j
38  ig=igpig(i)
39  ig=ivgig(ivg(1,ig)-j1,ivg(2,ig)-j2,ivg(3,ig)-j3)
40  o(i,j)=cfunig(ig)
41  end do
42 end do
43 end subroutine
44 !EOC
45 
integer ngkmax
Definition: modmain.f90:499
integer, dimension(:,:,:), allocatable ivgig
Definition: modmain.f90:402
integer, dimension(:,:), allocatable ivg
Definition: modmain.f90:400
pure subroutine olpistl(ngp, igpig, ld, o)
Definition: olpistl.f90:10
complex(8), dimension(:), allocatable cfunig
Definition: modmain.f90:434