2024-07-25 10:27:17 +02:00

1477 lines
56 KiB
Fortran
Executable File

************************************************************************
subroutine ecp(iang,a,ax,ay,az,jang,b,bx,by,bz,cx,cy,cz,work,ints,
$kangmax,ecpprim,ecpcoef,ecpn,ecpexp,nprimmax,nangmax,ptol,nicart,
$njcart,harmtab,ylmrst,dfact,icts,jcts,nispher,njspher,t4,dero,dx,
$dy,dz,intd,lid,ljd,lkd,bnm)
************************************************************************
* Calculate ECP integarls
* L. E. McMurchie and E. R. Davidson, J. Comp. Phys. 44, 289 (1981)
* R. M. Pitzer and N. W. Winter, Int. J. Quantum Chem. 40, 773 (1991)
* M. Kallay, unpublished (2015)
************************************************************************
implicit none
integer iang,jang,kangmax,nprimmax,nangmax,eang,l,nprim,dero
integer ecpn(nprimmax,0:nangmax),ecpprim(0:nangmax),nicart,njcart
integer harmtab(2,-6*nangmax:6*nangmax,0:6*nangmax),dx,dy,dz
integer njmax,i,nmax,lmax,mmax,lammaxi,lammaxj,nispher,njspher
real*8 a,ax,ay,az,b,bx,by,bz,work(*),cx,cy,cz,ints(njcart,nicart)
real*8 ecpexp(nprimmax,0:nangmax),ptol,ylmrst(*),icts(*),jcts(*)
real*8 ecpcoef(nprimmax,0:nangmax),dfact(-1:9*nangmax+1)
real*8 intd(njspher,nispher),t4(njcart,nispher)
real*8 bnm(0:2*nangmax+1,0:2*nangmax+1)
logical lid,ljd,lkd
C Calculate maximimum values
njmax=0
do eang=0,kangmax
do nprim=1,ecpprim(eang)
njmax=max(njmax,ecpn(nprim,eang))
enddo
enddo
nmax=iang+jang+njmax+dero
lmax=kangmax+10+dero
mmax=2*lmax+1
lammaxi=iang+lmax+1
lammaxj=jang+lmax+1
C Initialize derivative integrals
c if(dero.gt.0) ints=0.d0
C Type 1 integrals
c ints=0.d0 !szemet!!!szemet!!!szemet!!!szemet
if(dero.eq.0) then
call type1(work,ints,nmax,work((nmax+1)**2+1),iang,jang,
$a,b,ax,ay,az,bx,by,bz,cx,cy,cz,ecpprim,ecpcoef,ecpexp,ecpn,
$harmtab,ylmrst,nangmax,dfact,0,0,0,0,0,0,0,iang,jang,0.d0,1.d0,
$ptol,bnm)
else
if(lid) then
call type1(work,ints,nmax,work((nmax+1)**2+1),iang+dero,jang,
$a,b,ax,ay,az,bx,by,bz,cx,cy,cz,ecpprim,ecpcoef,ecpexp,ecpn,
$harmtab,ylmrst,nangmax,dfact,dero,dx,dy,dz,0,0,0,iang,jang,a,1.d0,
$ptol,bnm)
endif
if(ljd) then
call type1(work,ints,nmax,work((nmax+1)**2+1),iang,jang+dero,
$a,b,ax,ay,az,bx,by,bz,cx,cy,cz,ecpprim,ecpcoef,ecpexp,ecpn,
$harmtab,ylmrst,nangmax,dfact,dero,0,0,0,dx,dy,dz,iang,jang,b,1.d0,
$ptol,bnm)
endif
if(lkd) then
call type1(work,ints,nmax,work((nmax+1)**2+1),iang+dero,jang,
$a,b,ax,ay,az,bx,by,bz,cx,cy,cz,ecpprim,ecpcoef,ecpexp,ecpn,
$harmtab,ylmrst,nangmax,dfact,dero,dx,dy,dz,0,0,0,iang,jang,a,-1d0,
$ptol,bnm)
call type1(work,ints,nmax,work((nmax+1)**2+1),iang,jang+dero,
$a,b,ax,ay,az,bx,by,bz,cx,cy,cz,ecpprim,ecpcoef,ecpexp,ecpn,
$harmtab,ylmrst,nangmax,dfact,dero,0,0,0,dx,dy,dz,iang,jang,b,-1d0,
$ptol,bnm)
endif
endif
c
C Type 2 integrals
if(dero.eq.0) then
call type2(work,work((iang+1)*mmax*lammaxi+1),ints,lammaxi,
$lammaxj,iang+jang,
$work((iang+1)*mmax*lammaxi+(jang+1)*mmax*lammaxj+1),iang,jang,a,b,
$ax,ay,az,bx,by,bz,cx,cy,cz,kangmax,ecpprim,ecpcoef,ecpexp,ecpn,
$nprimmax,lmax,harmtab,ylmrst,nangmax,dfact,0,0,0,0,0,0,iang,jang,
$1.d0,ptol,bnm)
else
if(lid) then
call type2(work,work((iang+dero+1)*mmax*lammaxi+1),ints,
$lammaxi,lammaxj,iang+jang+dero,
$work((iang+dero+1)*mmax*lammaxi+(jang+1)*mmax*lammaxj+1),
$iang+dero,jang,a,b,ax,ay,az,bx,by,bz,cx,cy,cz,kangmax,ecpprim,
$ecpcoef,ecpexp,ecpn,nprimmax,lmax,harmtab,ylmrst,nangmax,dfact,
$dero,dx,dy,dz,1,0,iang,jang,1.d0,ptol,bnm)
endif
if(ljd) then
call type2(work,work((iang+1)*mmax*lammaxi+1),ints,
$lammaxi,lammaxj,iang+jang+dero,
$work((iang+1)*mmax*lammaxi+(jang+dero+1)*mmax*lammaxj+1),
$iang,jang+dero,a,b,ax,ay,az,bx,by,bz,cx,cy,cz,kangmax,ecpprim,
$ecpcoef,ecpexp,ecpn,nprimmax,lmax,harmtab,ylmrst,nangmax,dfact,
$dero,dx,dy,dz,0,1,iang,jang,1.d0,ptol,bnm)
endif
if(lkd) then
call type2(work,work((iang+dero+1)*mmax*lammaxi+1),ints,
$lammaxi,lammaxj,iang+jang+dero,
$work((iang+dero+1)*mmax*lammaxi+(jang+1)*mmax*lammaxj+1),
$iang+dero,jang,a,b,ax,ay,az,bx,by,bz,cx,cy,cz,kangmax,ecpprim,
$ecpcoef,ecpexp,ecpn,nprimmax,lmax,harmtab,ylmrst,nangmax,dfact,
$dero,dx,dy,dz,1,0,iang,jang,-1.d0,ptol,bnm)
call type2(work,work((iang+1)*mmax*lammaxi+1),ints,
$lammaxi,lammaxj,iang+jang+dero,
$work((iang+1)*mmax*lammaxi+(jang+dero+1)*mmax*lammaxj+1),
$iang,jang+dero,a,b,ax,ay,az,bx,by,bz,cx,cy,cz,kangmax,ecpprim,
$ecpcoef,ecpexp,ecpn,nprimmax,lmax,harmtab,ylmrst,nangmax,dfact,
$dero,dx,dy,dz,0,1,iang,jang,-1.d0,ptol,bnm)
endif
C Spherical contraction in the case of derivatives
C ints(jcart,icart) -> t4(jcart,ishper)
c call dgemm('n','t',njcart,nispher,nicart,1.d0,ints,njcart,icts,
c $nispher,0.d0,t4,njcart)
C t4(jcart,ishper) -> intd(jshper,ishper)
c call dgemm('n','n',njspher,nispher,njcart,-1.d0,jcts,njspher,t4,
c $njcart,1.d0,intd,njspher)
endif
c write(6,"(a5,5i3)") 'ints2'
c write(6,"(7f14.8)") dabs(ints)
C
return
end
C
************************************************************************
subroutine type1(ang,ints,nmax,qnl,iang,jang,aa,bb,ax,ay,az,bx,by,
$bz,cx,cy,cz,ecpprim,ecpcoef,ecpexp,ecpn,harmtab,ylmrst,nangmax,
$dfact,dero,dxa,dya,dza,dxb,dyb,dzb,diang,djang,ab,kfc,ptol,bnm)
************************************************************************
C Calculate type 1 ECP integrals
************************************************************************
implicit none
integer ecpprim,a,b,c,d,e,f,lam,lamup,iang,jang,ii,dero
integer ecpn(ecpprim),lammax,nmax,n,nprim,iprim,na,la,ma,nb,lb,mb
integer nangmax,harmtab(2,-6*nangmax:6*nangmax,0:6*nangmax)
integer diang,djang,dxa,dya,dza,dxb,dyb,dzb
real*8 ints(*),fac,cax,cay,caz,cbx,cby,cbz,kx,ky,kz,aa,bb,cx,cy,cz
real*8 rc,ambx,amby,ambz,rr,rr2,abp,sum,fac1,ylmrst(*),ab
real*8 ang(0:nmax,0:nmax),dfact(-1:*),ax,ay,az,kfc,ptol
real*8 qnl(0:nmax,0:nmax),t,temp,prr,pprr2,alpha,p,ar,prr2
real*8 ecpcoef(ecpprim),ecpexp(ecpprim),qnlam,bx,by,bz
real*8 bnm(0:2*nangmax+1,0:2*nangmax+1)
logical lid,ljd,lkd
C
p=aa+bb
cax=cx-ax
cay=cy-ay
caz=cz-az
cbx=cx-bx
cby=cy-by
cbz=cz-bz
ambx=ax-bx
amby=ay-by
ambz=az-bz
abp=(aa-bb)/p
ar=aa*bb*(ambx**2+amby**2+ambz**2)/p
kx=-2.d0*(aa*cax+bb*cbx)
ky=-2.d0*(aa*cay+bb*cby)
kz=-2.d0*(aa*caz+bb*cbz)
rc=kx*kx+ky*ky+kz*kz
if(rc.eq.0.d0) then
rr=0.d0
prr2=0.d0
pprr2=0.d0
prr=0.d0
lammax=0
lamup=0
else
rc=dsqrt(rc)
kx=kx/rc
ky=ky/rc
kz=kz/rc
rr2=(0.5d0*(ax+bx+abp*ambx)-cx)**2
$ +(0.5d0*(ay+by+abp*amby)-cy)**2
$ +(0.5d0*(az+bz+abp*ambz)-cz)**2
rr=dsqrt(rr2)
prr2=p*rr2
pprr2=p*prr2
prr=2.d0*p*rr
lammax=iang+jang
lamup=3*(iang+jang)
endif
C Calculate radial integrals
qnl(0:nmax,0:lammax)=0.d0
do iprim=1,ecpprim
nprim=ecpn(iprim)
alpha=p+ecpexp(iprim)
temp=ar+ecpexp(iprim)*prr2/alpha
if(temp.le.50.d0) then
fac=22.27331198732683138088d0*dexp(-temp)*ecpcoef(iprim)
if(prr.eq.0.d0) then
do n=0,nmax-mod(nmax,2),2
qnl(n,0)=qnl(n,0)+fac*qnlam(nprim+n,0,prr,alpha)
enddo
else
do lam=0,lammax
qnl(lam,lam)=qnl(lam,lam)+fac*
$qnlam(nprim+lam,lam,prr,alpha)
enddo
endif
endif
enddo
C Recursion, PW p. 777
if(prr.ne.0.d0) then
temp=dble(2*lammax-1)
do lam=lammax-2,0,-1
fac=temp/prr
do n=lam+2,lammax-mod(lammax-lam+1,2)+1,2
qnl(n,lam)=qnl(n,lam+2)+fac*qnl(n-1,lam+1)
enddo
temp=temp-2.d0
enddo
endif
C Prescreening
temp=0.d0
do n=0,nmax
do lam=0,lammax
temp=max(temp,dabs(qnl(n,lam)))
enddo
enddo
if(1d5*dble((nmax+1)*(lammax+1))*temp.lt.ptol) return
C Loop over Cartesian components
ii=0
do na=0,diang
do la=0,diang-na
ma=diang-na-la
do nb=0,djang
do lb=0,djang-nb
mb=djang-nb-lb
ii=ii+1
C Calculate angular integrals
ang(0:nmax,0:lammax)=0.d0
if(dero.eq.0) then
call omega_ijk_lam(na,la,ma,nb,lb,mb,kx,ky,kz,cax,cay,
$caz,cbx,cby,cbz,lamup,nmax,lammax,ang,harmtab,ylmrst,ylmrst,
$nangmax,dfact,1.d0,bnm)
else
call omega_ijk_lam(na-dxa,la-dya,ma-dza,nb-dxb,lb-dyb,
$mb-dzb,kx,ky,kz,cax,cay,caz,cbx,cby,cbz,lamup,nmax,lammax,ang,
$harmtab,ylmrst,ylmrst,nangmax,dfact,
$kfc*dble(na*dxa+la*dya+ma*dza+nb*dxb+lb*dyb+mb*dzb),bnm)
call omega_ijk_lam(na+dxa,la+dya,ma+dza,nb+dxb,lb+dyb,
$mb+dzb,kx,ky,kz,cax,cay,caz,cbx,cby,cbz,lamup,nmax,lammax,ang,
$harmtab,ylmrst,ylmrst,nangmax,dfact,-2.d0*kfc*ab,bnm)
endif
C Multiply angular and radial integrals
sum=0.d0
do lam=0,lammax
do n=lam,min(nmax,nmax-mod(nmax+1-lam,2)+1),2
sum=sum+ang(n,lam)*qnl(n,lam)
enddo
enddo
ints(ii)=ints(ii)+sum
enddo
enddo
enddo
enddo
C
return
end
C
************************************************************************
subroutine omega_ijk_lam(na,la,ma,nb,lb,mb,kx,ky,kz,cax,cay,caz,
$cbx,cby,cbz,lamup,nmax,lammax,ang,harmtab,dylmrst,iylmrst,nangmax,
$dfact,dfc,bnm)
************************************************************************
* Calculate Omega^IJK_lam (MD Eq. 28)
************************************************************************
implicit none
integer i,j,k,lam,r,s,t,mu,a,b,c,d,e,f,lamup,nmax,lammax,na,la,ma
integer nangmax,nb,lb,mb,ii,nn,n
integer harmtab(2,-6*nangmax:6*nangmax,0:6*nangmax)
integer*8 iylmrst(*)
real*8 fac,kx,ky,kz,om,om1,ylmrst,ang(0:nmax,0:*),dfc
real*8 cax,cay,caz,cbx,cby,cbz,fa,fb,fc,fd,fe,ff,dylmrst(*)
real*8 dfact(-1:9*nangmax+1),bnm(0:2*nangmax+1,0:2*nangmax+1)
C
do a=0,na
fa=bnm(na,a)*cax**(na-a)*dfc
do b=0,la
fb=bnm(la,b)*cay**(la-b)*fa
do c=0,ma
fc=bnm(ma,c)*caz**(ma-c)*fb
do d=0,nb
fd=bnm(nb,d)*cbx**(nb-d)*fc
i=a+d
do e=0,lb
fe=bnm(lb,e)*cby**(lb-e)*fd
j=b+e
do f=0,mb
ff=bnm(mb,f)*cbz**(mb-f)*fe
k=c+f
do lam=0,min(i+j+k,lamup)
if(mod(i+j+k-lam,2).eq.0) then
om=0.d0
do mu=-lam,lam
fac=0.d0
om1=0.d0
ii=harmtab(1,mu,lam)
nn=harmtab(2,mu,lam)
do n=1,nn
ylmrst=dylmrst(ii)
ii=ii+1
r=iylmrst(ii)
ii=ii+1
s=iylmrst(ii)
ii=ii+1
t=lam-r-s
fac=fac+ylmrst*kx**r*ky**s*kz**t
if(mod(i+r,2)+mod(j+s,2)+mod(k+t,2).eq.0)
$om1=om1+
$ylmrst*dfact(i+r-1)*dfact(j+s-1)*dfact(k+t-1)/dfact(i+r+j+s+k+t+1)
enddo
om=om+fac*om1
enddo
ang(i+j+k,lam)=ang(i+j+k,lam)+om*ff
endif
enddo
enddo
enddo
enddo
enddo
enddo
enddo
C
return
end
C
************************************************************************
subroutine omega_abc_lamlm(na,la,ma,kx,ky,kz,cx,cy,cz,iang,l,
$lamiup,lamupi,ang,lmax,harmtab,dylmrst,iylmrst,nangmax,dfact,dfc,
$bnm)
************************************************************************
* Calculate Omega^abc_lamlm (MD Eq. 29)
************************************************************************
implicit none
integer a,b,c,lam,l,m,r,s,t,mu,u,v,w,na,la,ma,lamiup,lamupi,lmax
integer iang,nangmax,ii,jj,nn,mm,iylmrst(*),n,jjj,mmm
integer harmtab(2,-6*nangmax:6*nangmax,0:6*nangmax)
real*8 fac,kx,ky,kz,om,om1,ylmrst,ylmuvw,cx,cy,cz,dfc
real*8 ang(0:iang,-lmax:lmax,0:*),fa,fb,fc,dylmrst(*)
real*8 dfact(-1:9*nangmax+1),bnm(0:2*nangmax+1,0:2*nangmax+1)
C
do a=0,na
fa=bnm(na,a)*cx**(na-a)*dfc
do b=0,la
fb=bnm(la,b)*cy**(la-b)*fa
do c=0,ma
fc=bnm(ma,c)*cz**(ma-c)*fb
do lam=max(l-a-b-c,0),min(l+a+b+c,lamupi)
if(mod(l+a+b+c-lam,2).eq.0) then
do m=-l,l
om=0.d0
jjj=harmtab(1,m,l)
mmm=harmtab(2,m,l)
do mu=-lam,lam
fac=0.d0
om1=0.d0
ii=harmtab(1,mu,lam)
nn=harmtab(2,mu,lam)
do n=1,nn
ylmrst=dylmrst(ii)
ii=ii+1
r=iylmrst(ii)
ii=ii+1
s=iylmrst(ii)
ii=ii+1
t=lam-r-s
fac=fac+ylmrst*kx**r*ky**s*kz**t
jj=jjj
do mm=1,mmm
ylmuvw=dylmrst(jj)
jj=jj+1
u=iylmrst(jj)
jj=jj+1
v=iylmrst(jj)
jj=jj+1
w=l-u-v
if(mod(a+r+u,2)+mod(b+s+v,2)+mod(c+t+w,2).eq.0)
$om1=om1+ylmrst*ylmuvw*dfact(a+r+u-1)*dfact(b+s+v-1)*dfact(c+t+w-1)
$/dfact(a+r+u+b+s+v+c+t+w+1)
enddo
enddo
om=om+fac*om1
enddo
ang(a+b+c,m,lam)=ang(a+b+c,m,lam)+om*fc
enddo
endif
enddo
enddo
enddo
enddo
C
return
end
C
************************************************************************
real*8 function qnlam(n,l,k,alpha)
************************************************************************
C Calculate Q^N_lambda(k,alpha) integrals (MD Eq. 19)
************************************************************************
implicit none
integer n,l,a,b
real*8 ddfact
real*8 k,zmin(0:8),z,alpha,fac,sum,fct,phi
data zmin/31.d0,28.d0,25.d0,23.d0,22.d0,20.d0,19.d0,18.d0,15.d0/
C
z=k**2/(4.d0*alpha)
if(mod(n+l,2).eq.0.and.n.gt.l) then
C Confluent hypergeometric series (MD Eq. 33)
fac=(0.5d0*k/alpha)**l*ddfact(l+n-1)/
$(2.d0*(2.d0*alpha)**((n-l)/2)*dsqrt(alpha)*ddfact(2*l+1))
b=2*l+3
sum=1.d0
fct=-1.d0
phi=sum
do a=l-n+2,-1,2
sum=sum*dble(a)*z/(dble(b)*fct)
fct=fct-1.d0
phi=phi+sum
b=b+2
enddo
else if(z.ge.zmin(min(n,8))) then
C Asymptotic series (MD Eq. 34)
fac=0.25d0*(k/(2.d0*alpha))**(n-2)/alpha**1.5d0
sum=1.d0
fct=1.d0
phi=1.d0
a=1-l-n
b=l-n+2
do while(dabs(sum/phi).gt.1d-12)
sum=sum*0.25d0*dble(a*b)/(fct*z)
phi=phi+sum
a=a+2
b=b+2
fct=fct+1.d0
enddo
else
C Confluent hypergeometric series (MD Eq. 36)
fac=dexp(-z)*(0.5d0*k/alpha)**l/(2.d0*alpha)**((n-l+1)/2)
if(mod(n+l,2).eq.0) then
fac=0.5d0*fac/dsqrt(alpha)
else
fac=0.56418958354775628694d0*fac
endif
a=l+n-1
b=2*l+1
sum=ddfact(a)/ddfact(b)
phi=sum
fct=0.d0
do while(dabs(sum/phi).gt.1d-14)
a=a+2
b=b+2
fct=fct+1.d0
sum=sum*z*dble(a)/(fct*dble(b))
phi=phi+sum
enddo
endif
qnlam=fac*phi
C
return
end
C
************************************************************************
real*8 function qnnlam(n,la,lb,ka,kb,alpha,shift)
************************************************************************
C Calculate Q^N_lambdalambda'(k_a,k_b,alpha) integrals (MD Eq. 25)
************************************************************************
implicit none
integer n,la,lb,npt,i
real*8 ka,kb,alpha,x,r,wt,sum,klz,shift,sumo
C
npt=1
r=0.66765696312261313405d0
wt=4.28571428571428558740d0
sum=0.5d0*
$wt*r**n*dexp((ka+kb)*r-alpha*r*r+shift)*klz(la,ka*r)*klz(lb,kb*r)
sumo=sum+1.d0
do while((dabs(sumo-sum).gt.1d-16.or.npt.lt.15).and.npt.lt.1024)
sumo=sum
sum=0.5d0*sum
npt=2*npt+1
do i=1,npt,2
x=dble(i)/dble(npt+1)
r=-5.d0*dlog(1.d0-x**3)
wt=15.d0*x**2/((1.d0-x**3)*dble(npt+1))
sum=sum+
$wt*r**n*dexp((ka+kb)*r-alpha*r*r+shift)*klz(la,ka*r)*klz(lb,kb*r)
enddo
enddo
qnnlam=sum/dsqrt(3.1415926535897932384626433832795028841d0)
C
return
end
C
************************************************************************
subroutine type2(angi,angj,ints,lammaxi,lammaxj,nmax,qnl,iang,
$jang,aa,bb,ax,ay,az,bx,by,bz,cx,cy,cz,kangmax,ecpprim,ecpcoef,
$ecpexp,ecpn,nprimmax,lmax,harmtab,ylmrst,nangmax,dfact,dero,
$dx,dy,dz,ddi,ddj,diang,djang,kfc,ptol,bnm)
************************************************************************
C Calculate type 2 ECP integrals
************************************************************************
implicit none
integer lamupi,lamupj,ii,lmax,kangmax,lin,l,lamilo,lamiup,dero
integer nprimmax,ecpprim(*),ecpn(nprimmax,*),lamjlo,lamjup
integer na,la,ma,nb,lb,mb,llo,lup,lami,lli,nli,lamj,llj,nlj,n,m
integer nlamilo,nlamiup,nlm,nmax,iang,jang,lammaxj,iprim,nprim
integer nangmax,harmtab(2,-6*nangmax:6*nangmax,0:6*nangmax),lamip
integer lamjp,lammaxi,dx,dy,dz,ddi,ddj,diang,djang
real*8 qnlam,qnnlam,rc,arg,fac1,fac2,sht,alpha,fac12,fac22,kfc
real*8 aa,bb,ax,ay,az,bx,by,bz,cx,cy,cz,aarr2,ylmrst(*),ptol
real*8 kax,kay,kaz,kbx,kby,kbz,rca,rcb,cax,cay,caz,cbx,cby,cbz,sum
real*8 ecpcoef(nprimmax,*),ecpexp(nprimmax,*),ints(*),fac,temp,p
real*8 qnl(0:nmax,0:lammaxi,0:lammaxj),angi(0:iang,-lmax:lmax,0:*)
real*8 angj(0:jang,-lmax:lmax,0:*),ka,kb
real*8 dfact(-1:9*nangmax+1),bnm(0:2*nangmax+1,0:2*nangmax+1)
C
p=aa+bb
cax=cx-ax
cay=cy-ay
caz=cz-az
cbx=cx-bx
cby=cy-by
cbz=cz-bz
rca=dsqrt(cax*cax+cay*cay+caz*caz)
rcb=dsqrt(cbx*cbx+cby*cby+cbz*cbz)
if(rca+rcb.eq.0.d0) then
llo=mod(iang,2)
lup=min(kangmax-1,iang,jang)
if(llo.ne.mod(jang,2)) return
lin=2
else if(rca.eq.0.d0) then
llo=mod(iang,2)
lup=min(kangmax-1,iang)
lin=2
else if(rcb.eq.0.d0) then
llo=mod(jang,2)
lup=min(kangmax-1,jang)
lin=2
else
llo=0
lup=kangmax-1
lin=1
endif
kax=-cax/rca
kay=-cay/rca
kaz=-caz/rca
kbx=-cbx/rcb
kby=-cby/rcb
kbz=-cbz/rcb
aarr2=(rca-rcb)**2*aa*bb/p
if(rca.eq.0.d0) then
lamupi=0
ka=0.d0
else
lamupi=kangmax+iang-1
ka=2.d0*aa*rca
endif
if(rcb.eq.0.d0) then
lamupj=0
kb=0.d0
else
lamupj=kangmax+jang-1
kb=2.d0*bb*rcb
endif
C Loop over l
do l=llo,lup,lin
lamilo=max(l-iang,0)
lamiup=min(lamupi,l+iang)
lamjlo=max(l-jang,0)
lamjup=min(lamupj,l+jang)
C Calculate radial integrals
qnl(0:nmax,0:lamiup,0:lamjup)=0.d0
C Loop over terms in ECP
do iprim=1,ecpprim(l+2)
nprim=ecpn(iprim,l+2)
alpha=p+ecpexp(iprim,l+2)
rc=0.5d0*(ka+kb)/alpha
arg=alpha*rc**2
sht=aarr2+ecpexp(iprim,l+2)*arg/p
if(sht.lt.50.d0) then
fac2=22.27331198732683138088d0*ecpcoef(iprim,l+2)
fac1=fac2*dexp(-sht)
C Calculate a subset of integrals, MD type 1 radial algorithm + grid
if(ka.eq.0.d0.and.kb.eq.0.d0) then
qnl(nmax,0,0)=qnl(nmax,0,0)+
$fac1*qnlam(nprim+nmax,0,0.d0,alpha)
else if(ka.eq.0.d0) then
do lamj=l,lamjup
qnl(lamj-l+iang,0,lamj)=qnl(lamj-l+iang,0,lamj)+
$fac1*qnlam(nprim+lamj-l+iang,lamj,kb,alpha)
enddo
else if(kb.eq.0.d0) then
do lami=l,lamiup
qnl(lami-l+jang,lami,0)=qnl(lami-l+jang,lami,0)+
$fac1*qnlam(nprim+lami-l+jang,lami,ka,alpha)
enddo
else if(nprim.eq.2) then !PW p. 778
call qnlbess(ka,kb,alpha,lammaxi,lamiup,lamjup,nmax,l,
$fac1,qnl,dfact,bnm,nangmax)
else !grid
do lamip=0,lamiup-l
do lamjp=0,lamjup-l
n=lamip+lamjp
qnl(n,lamip+l,lamjp+l)=qnl(n,lamip+l,lamjp+l)+
$fac2*qnnlam(n+nprim,lamip+l,lamjp+l,ka,kb,alpha,-sht-arg)
enddo
enddo
endif
endif
enddo
C Calculate the rest by recursion, PW p. 778
if(ka.eq.0.d0.and.kb.ne.0.d0) then
fac2=dble(2*lamjup-1)
do lamj=lamjup-2,lamjlo,-1
fac22=fac2/kb
do n=iang+iabs(lamj-l+1)+1,
$iang+jang-mod(jang-iabs(lamj-l),2),2
qnl(n,0,lamj)=qnl(n,0,lamj+2)+fac22*qnl(n-1,0,lamj+1)
enddo
fac2=fac2-2.d0
enddo
else if(ka.ne.0.d0.and.kb.eq.0.d0) then
fac1=dble(2*lamiup-1)
do lami=lamiup-2,lamilo,-1
fac12=fac1/ka
do n=jang+iabs(lami-l+1)+1,
$iang+jang-mod(iang-iabs(lami-l),2),2
qnl(n,lami,0)=qnl(n,lami+2,0)+fac12*qnl(n-1,lami+1,0)
enddo
fac1=fac1-2.d0
enddo
else if(ka.ne.0.d0.and.kb.ne.0.d0) then
fac1=dble(2*lamiup+3)
do lami=lamiup,lamilo,-1
lli=iabs(l-lami)+1
fac2=dble(2*lamjup+3)
fac12=fac1/ka
do lamj=lamjup,lamjlo,-1
llj=iabs(l-lamj)
fac22=fac2/kb
do n=lli+llj-1,(nmax-mod(iang-lli+1,2))-mod(jang-llj,2),2
if(lami+lamj-n.eq.2*l) cycle
if(lami.gt.lamiup-2.or.n+1.le.iabs(l-lami-2)+llj) then
qnl(n,lami,lamj)=qnl(n ,lami,lamj+2)+
$ fac22*qnl(n-1,lami,lamj+1)
else
qnl(n,lami,lamj)=qnl(n ,lami+2,lamj)+
$ fac12*qnl(n-1,lami+1,lamj)
endif
enddo
fac2=fac2-2.d0
enddo
fac1=fac1-2.d0
enddo
endif
C Prescreening
temp=0.d0
do n=0,nmax
do lami=0,lamiup
do lamj=0,lamjup
temp=max(temp,dabs(qnl(n,lami,lamj)))
enddo
enddo
enddo
if(1d6*dble((nmax+1)*(lamiup+1)*(lamjup+1))*temp.lt.ptol)
$goto 1234
C Loop over Cartesian components for bra
ii=0
do na=0,diang
do la=0,diang-na
ma=diang-na-la
C Calculate angular integrals for bra
angi(0:iang,-l:l,0:lamiup)=0.d0
if(ddi.eq.0) then
call omega_abc_lamlm(na,la,ma,kax,kay,kaz,cax,cay,caz,
$iang,l,lamiup,lamupi,angi,lmax,harmtab,ylmrst,ylmrst,nangmax,
$dfact,1.d0,bnm)
else
call omega_abc_lamlm(na-dx,la-dy,ma-dz,kax,kay,kaz,cax,
$cay,caz,iang,l,lamiup,lamupi,angi,lmax,harmtab,ylmrst,ylmrst,
$nangmax,dfact,kfc*dble(na*dx+la*dy+ma*dz),bnm)
call omega_abc_lamlm(na+dx,la+dy,ma+dz,kax,kay,kaz,cax,
$cay,caz,iang,l,lamiup,lamupi,angi,lmax,harmtab,ylmrst,ylmrst,
$nangmax,dfact,-2.d0*kfc*aa,bnm)
endif
C Loop over Cartesian components for ket
do nb=0,djang
do lb=0,djang-nb
mb=djang-nb-lb
C Calculate angular integrals for ket
angj(0:jang,-l:l,0:lamjup)=0.d0
if(ddj.eq.0) then
call omega_abc_lamlm(nb,lb,mb,kbx,kby,kbz,cbx,cby,cbz,
$jang,l,lamjup,lamupj,angj,lmax,harmtab,ylmrst,ylmrst,nangmax,
$dfact,1.d0,bnm)
else
call omega_abc_lamlm(nb-dx,lb-dy,mb-dz,kbx,kby,kbz,
$cbx,cby,cbz,jang,l,lamjup,lamupj,angj,lmax,harmtab,ylmrst,ylmrst,
$nangmax,dfact,kfc*dble(nb*dx+lb*dy+mb*dz),bnm)
call omega_abc_lamlm(nb+dx,lb+dy,mb+dz,kbx,kby,kbz,
$cbx,cby,cbz,jang,l,lamjup,lamupj,angj,lmax,harmtab,ylmrst,ylmrst,
$nangmax,dfact,-2.d0*kfc*bb,bnm)
endif
C Multiply angular and radial integrals
ii=ii+1
sum=0.d0
do lami=lamilo,lamiup
lli=iabs(l-lami)
nli=iang-mod(iang-lli,2)
do lamj=lamjlo,lamjup
llj=iabs(l-lamj)
nlj=jang-mod(jang-llj,2)
do n=lli+llj,nli+nlj,2
nlamilo=max(lli,n-nlj)
nlamiup=min(nli,n-llj)
temp=0.d0
do m=-l,l
do nlm=nlamilo,nlamiup,2
temp=temp+angi(nlm,m,lami)*angj(n-nlm,m,lamj)
enddo
enddo
sum=sum+temp*qnl(n,lami,lamj)
enddo
enddo
enddo
ints(ii)=ints(ii)+sum
enddo
enddo
enddo
enddo
1234 continue
enddo
C
return
end
C
************************************************************************
subroutine genharmtab(harmtab,dylmrst,iylmrst,nangmax,ii)
************************************************************************
* Generate coefficients for harmonic polynomials
************************************************************************
implicit none
integer nangmax,lam,mu,ii,nn,r,s,t
integer harmtab(2,-6*nangmax:6*nangmax,0:6*nangmax)
real*8 dylmrst(*),y,normcoeff,getcf,nrm
integer*8 iylmrst(*)
C
ii=1
do lam=0,6*nangmax
do mu=-lam,lam
nrm=dsqrt(dble(2*lam+1))*normcoeff(lam,mu)
harmtab(1,mu,lam)=ii
nn=0
do r=0,lam
do s=0,lam-r
t=lam-r-s
y=nrm*getcf(r,s,t,lam,mu)
if(dabs(y).gt.1d-20) then
nn=nn+1
dylmrst(ii)=y
ii=ii+1
iylmrst(ii)=r
ii=ii+1
iylmrst(ii)=s
ii=ii+1
endif
enddo
enddo
harmtab(2,mu,lam)=nn
enddo
enddo
C
return
end
C
************************************************************************
real*8 function klz(l,z)
************************************************************************
C Calculate modified spherical Bessel function scaled by e^-z, K_l(z)
************************************************************************
implicit none
integer l,k,j
real*8 z,rlz,fac,ddfact
C
if(z.lt.1d-7) then
klz=(1.d0-z)*z**l/ddfact(2*l+1)
else if(z.gt.16.d0) then
rlz=1.d0
fac=1.d0
do k=1,l
fac=fac*dble((l+k)*(l-k+1))/dble(k)
rlz=rlz+fac/(-2.d0*z)**k
enddo
klz=0.5d0*rlz/z
else
fac=dexp(-z)/ddfact(2*l+1)
rlz=fac
j=0
do while(fac.gt.1d-16)
j=j+1
fac=0.5d0*fac*z*z/dble(j*(2*j+2*l+1))
rlz=rlz+fac
enddo
klz=z**l*rlz
endif
C
return
end
C
************************************************************************
subroutine ecpint(natoms,nangmax,ncontrmax,nprimmax,ncartmax,
$nsphermax,nang,nangmin,ncontr,nprim,gexp,gcoef,coord,ctostr,
$nbasis,smat,work,gcn,nshrange,cartg,itol,necpatoms,ecpn,ecpnang,
$ecpprim,ecpexp,ecpcoef,ecpatoms,inti,intj,intf,ints,intc,dero,xyz,
$diatoms,bnm)
************************************************************************
* Calculate ECP integrals
************************************************************************
implicit none
integer natoms,nangmax,ncontrmax,nprimmax,ncartmax,iatoms,jatoms,j
integer nang(natoms),ncontr(0:nangmax,natoms),iang,jang,nsphermax
integer nprim(0:nangmax,natoms),iprim,jprim,i,l,iintc,kkatoms,dero
integer nicontr,niprim,nicart,njcontr,njprim,njcart,nispher
integer njspher,jntcc,katoms,ii,nmax,nmboys,jj,jangmax,icart,jcart
integer icontr,jcontr,ispher,jspher,jjntc,nbasis,necpatoms,idfact
integer nangmin(natoms),niprimmax,njprimmax,itn,jtn,iitn,jjtn
integer ecpn(nprimmax,0:nangmax,necpatoms),ecpnang(necpatoms),iecp
integer ecpprim(0:nangmax,necpatoms),ecpatoms(necpatoms),iylmrst
integer gcn(2,ncontrmax,0:nangmax,natoms),dx,dy,dz,xyz
integer nshrange(2,0:nangmax,natoms),diatoms
integer*8 binom
real*8 inti(*),intj(*),intf(*),ints(*),intc(*),itol,ptol
real*8 ecpexp(nprimmax,0:nangmax,necpatoms)
real*8 ecpcoef(nprimmax,0:nangmax,necpatoms)
real*8 gexp(nprimmax,0:nangmax,natoms),coord(3,natoms)
real*8 gcoef(nprimmax,ncontrmax,0:nangmax,natoms)
real*8 cx,cy,cz,work(*),ax,bx,ay,by,az,bz,a,b
real*8 ctostr(ncartmax**2,0:nangmax),smat(nbasis,nbasis),ddfact
real*8 bnm(0:2*(nangmax+dero)+1,0:2*(nangmax+dero)+1)
logical cartg
C Generate constants for ECP calculation
iylmrst=2*(12*(nangmax+dero)+1)*(6*(nangmax+dero)+1)+1
call genharmtab(work,work(iylmrst),work(iylmrst),nangmax+dero,
$iecp)
iecp=iecp+iylmrst
idfact=iecp
do i=-1,9*(nangmax+dero)+1
work(iecp)=ddfact(i)
iecp=iecp+1
enddo
do i=0,2*(nangmax+dero)+1
do j=0,i
bnm(i,j)=dble(binom(i,j))
enddo
enddo
dx=0
dy=0
dz=0
if(dero.gt.0) then
if(xyz.eq.1) then
dx=1
else if(xyz.eq.2) then
dy=1
else
dz=1
endif
endif
C Loop over atoms
do iatoms=1,natoms
ax=coord(1,iatoms)
ay=coord(2,iatoms)
az=coord(3,iatoms)
do jatoms=1,iatoms
bx=coord(1,jatoms)
by=coord(2,jatoms)
bz=coord(3,jatoms)
C Prescreening nuclear attraction integrals
do iang=nangmin(iatoms),nang(iatoms)
iitn=nshrange(1,iang,iatoms)
nicontr=ncontr(iang,iatoms)
niprim=nprim(iang,iatoms)
nicart=(iang+1)*(iang+2)/2
nispher=2*iang+1
if(cartg) nispher=nicart
jangmax=nang(jatoms)
if(iatoms.eq.jatoms) jangmax=iang
do jang=nangmin(jatoms),jangmax
jjtn=nshrange(1,jang,jatoms)
njcontr=ncontr(jang,jatoms)
njprim=nprim(jang,jatoms)
njcart=(jang+1)*(jang+2)/2
njspher=2*jang+1
if(cartg) njspher=njcart
ptol=itol/dble(max(1,niprim*njprim*natoms))
iintc=1
do iprim=1,niprim
a=gexp(iprim,iang,iatoms)
jjntc=1
do jprim=1,njprim
b=gexp(jprim,jang,jatoms)
ints(1:nicart*njcart)=0.d0
if(dero.eq.0)then
C Calculate core potential integrals
do kkatoms=1,necpatoms
katoms=ecpatoms(kkatoms)
cx=coord(1,katoms)
cy=coord(2,katoms)
cz=coord(3,katoms)
c write(6,"(a4,5i3)") 'ints',iang,jang,iprim,jprim,katoms
call ecp(iang,a,ax,ay,az,jang,b,bx,by,bz,cx,cy,cz,
$work(iecp),ints,ecpnang(kkatoms),ecpprim(0,kkatoms),
$ecpcoef(1,0,kkatoms),ecpn(1,0,kkatoms),ecpexp(1,0,kkatoms),
$nprimmax,nangmax,ptol,nicart,njcart,work,work(iylmrst),
$work(idfact),ctostr(1,iang),ctostr(1,jang),nispher,njspher,work,
$0,0,0,0,ints,.false.,.false.,.false.,bnm)
enddo
C Calculate core potential integral derivatives
else
do kkatoms=1,necpatoms
katoms=ecpatoms(kkatoms)
if(iatoms.eq.diatoms.or.jatoms.eq.diatoms.or.
$ katoms.eq.diatoms) then
cx=coord(1,katoms)
cy=coord(2,katoms)
cz=coord(3,katoms)
call ecp(iang,a,ax,ay,az,jang,b,bx,by,bz,cx,cy,
$cz,work(iecp),ints,ecpnang(kkatoms),ecpprim(0,kkatoms),
$ecpcoef(1,0,kkatoms),ecpn(1,0,kkatoms),ecpexp(1,0,kkatoms),
$nprimmax,nangmax+dero,ptol,nicart,njcart,work,work(iylmrst),
$work(idfact),ctostr(1,iang),ctostr(1,jang),nispher,njspher,work,
$dero,dx,dy,dz,ints,iatoms.eq.diatoms,
$jatoms.eq.diatoms,katoms.eq.diatoms,bnm)
endif
enddo
ints(1:nicart*njcart)=-ints(1:nicart*njcart)
endif
C Convert Cartesian GTOs to spherical ones for i
C intss(jcart,icart) -> intcs(jcart,ispher)
call sphconi(intc,ints,ctostr(1,iang),nicart,njcart,
$nispher,iang,cartg)
C Convert Cartesian GTOs to spherical ones for j
C intcs(jcart,ispher) -> intjs(jspher,ispher,jprim)
call sphconj(intj(jjntc),intc,ctostr(1,jang),njcart,
$nispher,njspher,jang,cartg)
jjntc=jjntc+nispher*njspher
enddo ! jprim
C Multiply with contrcation coefficients for j
C intj(jspher,ispher,jprim) -> inti(jspher,ispher,jcontr,iprim)
call pricon(inti(iintc),intj,gcoef(1,1,jang,jatoms),
$nispher*njspher,njcontr,nprimmax,gcn(1,1,jang,jatoms))
iintc=iintc+nispher*njspher*njcontr
enddo ! iprim
C Multiply with contrcation coefficients for i
C inti(jspher,ispher,jcontr,iprim) -> intf(jspher,ispher,jcontr,icontr)
call pricon(intf,inti,gcoef(1,1,iang,iatoms),
$nispher*njspher*njcontr,nicontr,nprimmax,gcn(1,1,iang,iatoms))
C Save integrals to memory
l=0
itn=iitn
do icontr=1,nicontr
jtn=jjtn
do jcontr=1,njcontr
ii=itn
do ispher=1,nispher
ii=ii+1
jj=jtn
do jspher=1,njspher
jj=jj+1
l=l+1
smat(ii,jj)=smat(ii,jj)+intf(l)
c write(6,"(2i3,f14.5)") ii,jj,intf(l)
enddo
enddo
jtn=jtn+njspher
enddo
itn=itn+nispher
enddo
C
enddo !jang
enddo !iang
enddo !jatoms
enddo !iatoms
C
return
end
C
************************************************************************
real*8 function normcoeff(l,m)
************************************************************************
* Calculate norm for transformation matrix for a given l,m
************************************************************************
implicit none
integer l,m
integer*8 fact
C
if(m.eq.0) then
normcoeff=
$dsqrt(2.d0*dble(fact(l+iabs(m)))*dble(fact(l-iabs(m)))/2.d0)/
$(2.d0**iabs(m)*dble(fact(l)))
else
normcoeff=
$dsqrt(2.d0*dble(fact(l+iabs(m)))*dble(fact(l-iabs(m))) )/
$(2.d0**iabs(m)*dble(fact(l)))
endif
c write(6,*)normcoeff,l,m
C
return
end
C
************************************************************************
subroutine sphconi(intc,ints,ctostr,nicart,njcart,nispher,iang,
$cartg)
************************************************************************
* Spherical contraction for i
************************************************************************
implicit none
integer nicart,njcart,nispher,iang
real*8 intc(njcart,nispher),ints(njcart,nicart)
real*8 ctostr(nispher,nicart)
logical cartg
C
if(iang.le.1.or.cartg) then !s, p
intc=ints
else if(iang.eq.2) then!d
intc(1:njcart, 1)= ints(1:njcart, 5)*1.73205080756887697113d0
intc(1:njcart, 2)= ints(1:njcart, 2)*1.73205080756887697113d0
intc(1:njcart, 3)= ints(1:njcart, 1)-
$ (ints(1:njcart, 3)+
$ ints(1:njcart, 6))*0.5d0
intc(1:njcart, 4)= ints(1:njcart, 4)*1.73205080756887697113d0
intc(1:njcart, 5)=(ints(1:njcart, 6)-
$ ints(1:njcart, 3))*0.86602540378443848557d0
else if(iang.eq.3) then!f
intc(1:njcart, 1)= ints(1:njcart, 9)*2.37170824512628453107d0-
$ ints(1:njcart, 4)*0.79056941504209476967d0
intc(1:njcart, 2)= ints(1:njcart, 6)*3.87298334620741657730d0
intc(1:njcart, 3)= ints(1:njcart, 2)*2.44948974278317788134d0-
$ (ints(1:njcart, 4)+
$ ints(1:njcart, 9))*0.61237243569579447033d0
intc(1:njcart, 4)= ints(1:njcart, 1)-
$ (ints(1:njcart, 3)+
$ ints(1:njcart, 8))*1.5d0
intc(1:njcart, 5)= ints(1:njcart, 5)*2.44948974278317788134d0-
$ (ints(1:njcart, 7)+
$ ints(1:njcart,10))*0.61237243569579447033d0
intc(1:njcart, 6)=(ints(1:njcart, 8)-
$ ints(1:njcart, 3))*1.93649167310370828865d0
intc(1:njcart, 7)= ints(1:njcart,10)*0.79056941504209476967d0-
$ ints(1:njcart, 7)*2.37170824512628453107d0
else if(iang.eq.4) then!g
intc(1:njcart, 1)=(ints(1:njcart,14)-
$ ints(1:njcart, 9))*2.95803989154980806475d0
intc(1:njcart, 2)= ints(1:njcart,11)*6.27495019900556627590d0-
$ ints(1:njcart, 4)*2.09165006633518890666d0
intc(1:njcart, 3)= ints(1:njcart, 7)*6.70820393249936941515d0-
$ (ints(1:njcart, 9)+
$ ints(1:njcart,14))*1.11803398874989490253d0
intc(1:njcart, 4)= ints(1:njcart, 2)*3.16227766016837907870d0-
$ (ints(1:njcart, 4)+
$ ints(1:njcart,11))*2.37170824512628453107d0
intc(1:njcart, 5)= ints(1:njcart, 1)-
$ (ints(1:njcart, 3)+
$ ints(1:njcart,10))*3.d0+
$ (ints(1:njcart, 5)+
$ ints(1:njcart,15))*0.375d0+
$ ints(1:njcart,12)*0.75d0
intc(1:njcart, 6)= ints(1:njcart, 6)*3.16227766016837907870d0-
$ (ints(1:njcart, 8)+
$ ints(1:njcart,13))*2.37170824512628453107d0
intc(1:njcart, 7)=(ints(1:njcart,10)-
$ ints(1:njcart, 3))*3.35410196624968470758d0+
$ (ints(1:njcart, 5)-
$ ints(1:njcart,15))*0.55901699437494745126d0
intc(1:njcart, 8)= ints(1:njcart,13)*2.09165006633518890666d0-
$ ints(1:njcart, 8)*6.27495019900556627590d0
intc(1:njcart, 9)=(ints(1:njcart, 5)+
$ ints(1:njcart,15))*0.73950997288745201619d0-
$ ints(1:njcart,12)*4.43705983732471231917d0
else if(iang.eq.5) then!h
intc(1:njcart, 1)= ints(1:njcart, 6)*0.70156076002011391601d0-
$ ints(1:njcart,15)*7.01560760020113871605d0+
$ ints(1:njcart,20)*3.50780380010056935802d0
intc(1:njcart, 2)=(ints(1:njcart,17)-
$ ints(1:njcart,10))*8.87411967464942463835d0
intc(1:njcart, 3)= ints(1:njcart, 6)*0.52291251658379711564d0-
$ ints(1:njcart, 4)*4.18330013267037692515d0+
$ ints(1:njcart,13)*12.54990039801113255180d0-
$ ints(1:njcart,15)*1.04582503316759423129d0-
$ ints(1:njcart,20)*1.56873754975139156898d0
intc(1:njcart, 4)= ints(1:njcart, 8)*10.24695076595959974952d0-
$ (ints(1:njcart,10)+
$ ints(1:njcart,17))*5.12347538297979987476d0
intc(1:njcart, 5)= ints(1:njcart, 2)*3.87298334620741613321d0-
$ (ints(1:njcart, 4)+
$ ints(1:njcart,13))*5.80947501931112419982d0+
$ ints(1:njcart,15)*0.96824583655185403330d0+
$ (ints(1:njcart, 6)+
$ ints(1:njcart,20))*0.48412291827592701665d0
intc(1:njcart, 6)= ints(1:njcart, 1)-
$ (ints(1:njcart, 3)+
$ ints(1:njcart,12))*5.d0+
$ ints(1:njcart,14)*3.75d0+
$ (ints(1:njcart, 5)+
$ ints(1:njcart,19))*1.875d0
intc(1:njcart, 7)= ints(1:njcart, 7)*3.87298334620741613321d0-
$ (ints(1:njcart, 9)+
$ ints(1:njcart,16))*5.80947501931112419982d0+
$ ints(1:njcart,18)*0.96824583655185403330d0+
$ (ints(1:njcart,11)+
$ ints(1:njcart,21))*0.48412291827592701665d0
intc(1:njcart, 8)=(ints(1:njcart,12)-
$ ints(1:njcart, 3))*5.12347538297979987476d0+
$ (ints(1:njcart, 5)-
$ ints(1:njcart,19))*2.56173769148989993738d0
intc(1:njcart, 9)= ints(1:njcart,11)*1.56873754975139156898d0-
$ ints(1:njcart, 9)*12.54990039801113255180d0+
$ ints(1:njcart,16)*4.18330013267037692515d0+
$ ints(1:njcart,18)*1.04582503316759423129d0-
$ ints(1:njcart,21)*0.52291251658379711564d0
intc(1:njcart,10)=(ints(1:njcart, 5)+
$ ints(1:njcart,19))*2.21852991866235615959d0-
$ ints(1:njcart,14)*13.31117951197413695752d0
intc(1:njcart,11)= ints(1:njcart,11)*3.50780380010056935802d0-
$ ints(1:njcart,18)*7.01560760020113871605d0+
$ ints(1:njcart,21)*0.70156076002011391601d0
else
intc=matmul(ints,transpose(ctostr))
endif
C
return
end
C
************************************************************************
subroutine sphconj(intj,intc,ctostr,njcart,nispher,njspher,jang,
$cartg)
************************************************************************
* Spherical contraction for j
************************************************************************
implicit none
integer njcart,nispher,njspher,jang
real*8 intj(njspher,nispher),intc(njcart,nispher)
real*8 ctostr(njspher,njcart)
logical cartg
C
if(jang.le.1.or.cartg) then !s, p
intj=intc
else if(jang.eq.2) then!d
intj(1,1:nispher)= intc( 5,1:nispher)*1.73205080756887697113d0
intj(2,1:nispher)= intc( 2,1:nispher)*1.73205080756887697113d0
intj(3,1:nispher)= intc( 1,1:nispher)-
$ (intc( 3,1:nispher)+
$ intc( 6,1:nispher))*0.5d0
intj(4,1:nispher)= intc( 4,1:nispher)*1.73205080756887697113d0
intj(5,1:nispher)=(intc( 6,1:nispher)-
$ intc( 3,1:nispher))*0.86602540378443848557d0
else if(jang.eq.3) then!f
intj(1,1:nispher)= intc( 9,1:nispher)*2.37170824512628453107d0-
$ intc( 4,1:nispher)*0.79056941504209476967d0
intj(2,1:nispher)= intc( 6,1:nispher)*3.87298334620741657730d0
intj(3,1:nispher)= intc( 2,1:nispher)*2.44948974278317788134d0-
$ (intc( 4,1:nispher)+
$ intc( 9,1:nispher))*0.61237243569579447033d0
intj(4,1:nispher)= intc( 1,1:nispher)-
$ (intc( 3,1:nispher)+
$ intc( 8,1:nispher))*1.5d0
intj(5,1:nispher)= intc( 5,1:nispher)*2.44948974278317788134d0-
$ (intc( 7,1:nispher)+
$ intc(10,1:nispher))*0.61237243569579447033d0
intj(6,1:nispher)=(intc( 8,1:nispher)-
$ intc( 3,1:nispher))*1.93649167310370828865d0
intj(7,1:nispher)= intc(10,1:nispher)*0.79056941504209476967d0-
$ intc( 7,1:nispher)*2.37170824512628453107d0
else if(jang.eq.4) then!g
intj(1,1:nispher)=(intc(14,1:nispher)-
$ intc( 9,1:nispher))*2.95803989154980806475d0
intj(2,1:nispher)= intc(11,1:nispher)*6.27495019900556627590d0-
$ intc( 4,1:nispher)*2.09165006633518890666d0
intj(3,1:nispher)= intc( 7,1:nispher)*6.70820393249936941515d0-
$ (intc( 9,1:nispher)+
$ intc(14,1:nispher))*1.11803398874989490253d0
intj(4,1:nispher)= intc( 2,1:nispher)*3.16227766016837907870d0-
$ (intc( 4,1:nispher)+
$ intc(11,1:nispher))*2.37170824512628453107d0
intj(5,1:nispher)= intc( 1,1:nispher)-
$ (intc( 3,1:nispher)+
$ intc(10,1:nispher))*3.d0+
$ (intc( 5,1:nispher)+
$ intc(15,1:nispher))*0.375d0+
$ intc(12,1:nispher)*0.75d0
intj(6,1:nispher)= intc( 6,1:nispher)*3.16227766016837907870d0-
$ (intc( 8,1:nispher)+
$ intc(13,1:nispher))*2.37170824512628453107d0
intj(7,1:nispher)=(intc(10,1:nispher)-
$ intc( 3,1:nispher))*3.35410196624968470758d0+
$ (intc( 5,1:nispher)-
$ intc(15,1:nispher))*0.55901699437494745126d0
intj(8,1:nispher)= intc(13,1:nispher)*2.09165006633518890666d0-
$ intc( 8,1:nispher)*6.27495019900556627590d0
intj(9,1:nispher)=(intc( 5,1:nispher)+
$ intc(15,1:nispher))*0.73950997288745201619d0-
$ intc(12,1:nispher)*4.43705983732471231917d0
else if(jang.eq.5) then!h
intj(1,1:nispher)= intc( 6,1:nispher)*0.70156076002011391601d0-
$ intc(15,1:nispher)*7.01560760020113871605d0+
$ intc(20,1:nispher)*3.50780380010056935802d0
intj(2,1:nispher)=(intc(17,1:nispher)-
$ intc(10,1:nispher))*8.87411967464942463835d0
intj(3,1:nispher)= intc( 6,1:nispher)*0.52291251658379711564d0-
$ intc( 4,1:nispher)*4.18330013267037692515d0+
$ intc(13,1:nispher)*12.54990039801113255180d0-
$ intc(15,1:nispher)*1.04582503316759423129d0-
$ intc(20,1:nispher)*1.56873754975139156898d0
intj(4,1:nispher)= intc( 8,1:nispher)*10.24695076595959974952d0-
$ (intc(10,1:nispher)+
$ intc(17,1:nispher))*5.12347538297979987476d0
intj(5,1:nispher)= intc( 2,1:nispher)*3.87298334620741613321d0-
$ (intc( 4,1:nispher)+
$ intc(13,1:nispher))*5.80947501931112419982d0+
$ intc(15,1:nispher)*0.96824583655185403330d0+
$ (intc( 6,1:nispher)+
$ intc(20,1:nispher))*0.48412291827592701665d0
intj(6,1:nispher)= intc( 1,1:nispher)-
$ (intc( 3,1:nispher)+
$ intc(12,1:nispher))*5.d0+
$ intc(14,1:nispher)*3.75d0+
$ (intc( 5,1:nispher)+
$ intc(19,1:nispher))*1.875d0
intj(7,1:nispher)= intc( 7,1:nispher)*3.87298334620741613321d0-
$ (intc( 9,1:nispher)+
$ intc(16,1:nispher))*5.80947501931112419982d0+
$ intc(18,1:nispher)*0.96824583655185403330d0+
$ (intc(11,1:nispher)+
$ intc(21,1:nispher))*0.48412291827592701665d0
intj(8,1:nispher)=(intc(12,1:nispher)-
$ intc( 3,1:nispher))*5.12347538297979987476d0+
$ (intc( 5,1:nispher)-
$ intc(19,1:nispher))*2.56173769148989993738d0
intj(9,1:nispher)= intc(11,1:nispher)*1.56873754975139156898d0-
$ intc( 9,1:nispher)*12.54990039801113255180d0+
$ intc(16,1:nispher)*4.18330013267037692515d0+
$ intc(18,1:nispher)*1.04582503316759423129d0-
$ intc(21,1:nispher)*0.52291251658379711564d0
intj(10,1:nispher)=(intc( 5,1:nispher)+
$ intc(19,1:nispher))*2.21852991866235615959d0-
$ intc(14,1:nispher)*13.31117951197413695752d0
intj(11,1:nispher)=intc(11,1:nispher)*3.50780380010056935802d0-
$ intc(18,1:nispher)*7.01560760020113871605d0+
$ intc(21,1:nispher)*0.70156076002011391601d0
else
intj=matmul(ctostr,intc)
endif
C
return
end
C
************************************************************************
subroutine pricon(intf,inti,gcoefi,n,nicontr,nprimmax,gcni)
************************************************************************
* Primitive contraction
************************************************************************
implicit none
integer n,nicontr,nprimmax,i,icontr,iprim,gcni(2,nicontr)
real*8 intf(n,nicontr),inti(n,*),gcoefi(nprimmax,nicontr)
C
intf=0.d0
do icontr=1,nicontr
do iprim=gcni(1,icontr),gcni(2,icontr)
intf(1:n,icontr)=
$ intf(1:n,icontr)+inti(1:n,iprim)*gcoefi(iprim,icontr)
enddo
enddo
C
return
end
C
************************************************************************
real*8 function getcf(i,j,k,l,m)
************************************************************************
implicit none
integer i,j,k,l,m,t,u,ii,jj,kk,v,binom,tmax,vmax
real*8 gcf,isum
real*8 vv,vm,normcoeff
C
if(m.lt.0) then
vm=0.5d0
if(mod(m,2).eq.0) then
vmax=iabs(m)/2-1
else
vmax=(iabs(m)-1)/2
endif
else
vm=0.d0
if(mod(m,2).eq.0) then
vmax=iabs(m)/2
else
vmax=(iabs(m)-1)/2
endif
endif
isum=0
if(mod(l-iabs(m),2).eq.0) then
tmax=(l-iabs(m))/2
else
tmax=(l-iabs(m)-1)/2
endif
do t=0,tmax
do u=0,t
do v=0,vmax
vv=dble(v)+vm
ii=2*t+iabs(m)-idnint(2.d0*(dble(u)+vv))
jj=idnint(2.d0*(dble(u)+vv))
kk=l-2*t-iabs(m)
if(i.eq.ii.and.j.eq.jj.and.k.eq.kk) then
isum=isum+(-1.d0)**(dble(t)+vv-vm)*
$dble(binom(l,t))*dble(binom(l-t,iabs(m)+t))*dble(binom(t,u))*
$dble(binom(iabs(m),idnint(2.d0*vv)))/4.d0**t
C write(6,*)'ijk m: ',i,j,k,m
C write(6,*)'tuv: ',t,u,v,
C $ 'Coeff:',gcf(l,m,t,u,v),
C $ 'Norm: ', normcoeff(l,m)
C write(6,*)
endif
enddo
enddo
enddo
getcf=isum
c getcf=isum*normcoeff(l,m)
c write(6,*)normcoeff(l,m)
C
return
end
C
************************************************************************
subroutine qnlbess(ka,kb,alpha,lammaxi,laup,muup,nmax,nu,fac1,qnl,
$dfact,bnm,nangmax)
************************************************************************
* Calculate type 2 radial integrals using Bessel function expansion
* Pitzer and Winter, p. 778
************************************************************************
implicit none
integer lammaxi,laup,muup,nmax,nu,u,t,la,mu,nangmax,binom
real*8 ka,kb,kata,kbta,tmp,faca,faca1,faca2,facu,facu1,talpha
real*8 alpha,fac1,facm,facm1,facm2,fact,fact1,fact2,facu2
real*8 qnl(0:nmax,0:lammaxi,0:*),blx(0:laup+muup-nu)
real*8 afac(0:muup+1),bfac(0:laup+1),afac1(0:muup+1),dfact(-1:*)
real*8 bfac1(0:laup+1),afac2(0:muup+1),bbfac(0:muup+1)
real*8 bnm(0:2*nangmax+1,0:2*nangmax+1)
C
talpha=2.d0*alpha
kata=ka/talpha
kbta=kb/talpha
tmp=ka*kata
afac(0)=1.d0
do mu=1,muup-nu
afac(mu)=tmp*afac(mu-1)
enddo
tmp=kb*kbta
bfac(0)=1.d0
do la=1,laup-nu
bfac(la)=tmp*bfac(la-1)
enddo
bbfac(0)=kbta**nu
do mu=1,muup-nu
bbfac(mu)=kbta*bbfac(mu-1)
enddo
call ilx(ka*kbta,blx,laup+muup-nu)
do la=0,laup+muup-nu
blx(la)=blx(la)/dfact(2*la-1)
enddo
faca=ka**nu*fac1*talpha**(-1.5d0)/dsqrt(2.d0)/dble(2*nu+1)
faca1=dble(2*nu+3)
faca2=faca1
do la=0,laup-nu
facu=1.d0
facu1=dble(2*nu+1)
facu2=dble(2*(la+nu)+3)
do u=0,muup-nu
afac1(u)=facu*afac(u)
facu=facu*facu1/facu2
facu1=facu1+2.d0
facu2=facu2+2.d0
enddo
do t=0,la
bfac1(t)=bnm(la,t)*bfac(t)
enddo
facm=faca
facm1=dble(2*(la+nu)+3)
facm2=dble(2*nu+3)
do mu=0,muup-nu
do u=0,mu
afac2(u)=bnm(mu,u)*afac1(u)
enddo
tmp=0.d0
fact=1.d0
fact1=dble(2*nu+1)
fact2=dble(2*(mu+nu)+3)
do t=0,la
do u=0,mu
tmp=tmp+fact*afac2(u)*bfac1(t)*blx(t+u+nu)
enddo
fact=fact*fact1/fact2
fact1=fact1+2.d0
fact2=fact2+2.d0
enddo
qnl(la+mu,la+nu,mu+nu)=
$ qnl(la+mu,la+nu,mu+nu)+facm*bbfac(mu)*tmp
facm=facm*facm1/facm2
facm1=facm1+2.d0
facm2=facm2+2.d0
enddo
faca=kata*faca*faca1/faca2
faca1=faca1+2.d0
faca2=faca2+2.d0
enddo
C
return
end
C
************************************************************************
subroutine ilx(x,blx,lmax)
************************************************************************
* Compute scaled spherical modified Bessel functions of the first kind
* Wahl, Cade, and Roothaan, JCP 41, 2578 (1964), p. 2S93
************************************************************************
implicit none
integer lmax,l
real*8 x,blx(0:lmax),x2,em2x,fac
real*8 pj,pjm1,qj,qjm1,a,aj,ajpjm1,ajqjm1
C
x2=x*x
C Calculate b_l(x) by simple recurrence relation for large x (VI.4)
if(x.gt.dble(abs(3*lmax-1))) then
em2x=dexp(-2.d0*x)
blx(0)=0.5d0*(1.d0-em2x)/x
if(lmax.eq.0) return
blx(1)=1.5d0*(1.d0+em2x+(em2x-1.d0)/x)/x2
if(lmax.eq.1) return
do l=2,lmax
blx(l)=dble((2*l+1)*(2*l-1))*(blx(l-2)-blx(l-1))/x2
enddo
else
C Caclulate rho_l(x) by continued fraction (VI.8), double recursion alg.
l=lmax
pj=1.d0
pjm1=0.d0
qj=1.d0
qjm1=1.d0
a=1.d0
l=lmax-1
do while(dabs(a).gt.1d-16)
l=l+1
aj=x2/dble((2*l+3)*(2*l+1))
ajpjm1=aj*pjm1
pjm1=pj
pj=pj+ajpjm1
ajqjm1=aj*qjm1
qjm1=qj
qj=qj+ajqjm1
a=aj*a
enddo
blx(lmax)=pj/qj
C Caclulate rho_l(x) by downward recursion
l=lmax
do while(l.ne.0)
l=l-1
fac=dble((2*l+3)*(2*l+1))
blx(l)=fac/(fac+x2*blx(l+1))
enddo
blx(0)=blx(0)/(1.d0+x*blx(0))
C Caclulate b_l(x) by VI.6
do l=1,lmax
blx(l)=blx(l)*blx(l-1)
enddo
endif
C
return
end
C