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

3796 lines
127 KiB
Fortran

C***********************************************************************
subroutine propcore(r8heap,i4heap,dscfa,dscfb,dscfaao,dscfbao,
$dscftao,ca,cb,s,xdip,focka,fockb,
$da,da2,db,db2,fockmoa,fockmob,iprimea,
$iprimeb,xinta,xintb,xintrohf,xintrohf2,xintaoa,xintaob,xintfa,
$xintfb,xintfrohf,gmtx,gmtxuhf,gmtxrohf,daoa,
$daob,ha,hb,ipiv,ipivuhf,ipivrohf,ipivrohf2,
$zvect,xvect,scftype,core,calctype,nos,natoms,
$coord,atnum,sder,hder,tdera,tderb,wr,wi,vl,vr,work,gmtxrohf_2,
$xintrohf_2,xintrohf2_2,g2aa,g2bb,g2ab,vint2aa,vint2bb,
$vint2ab,maxw,g2int,g2int2,gaoa,gaob,gaoab,tdera2,tderb2,lcatoms,
$symat,pggen,efield,atsymbol,grdens,dft,intpos,natrange,bfat,atchg,
$popul,ncent,qmmm) !HB
C***********************************************************************
C Driver for property calculations
C***********************************************************************
#include "MRCCCOMMON"
integer maxw,nblock,noccblock,nosblock,nosshift,nblockal,nblockbe
integer nblockcore,ndoublyoccblock,nblockvirt,nvirtshift,nblockao
integer lcatoms(natoms),symat(natoms,nir),natrange(2,natoms)
integer ii,jj,nn,i,j,k,l,m,n,p,q,r,z,t,a,e,f,ff,g,h,hh
integer info,ipiv(nvirtal*nal),info2,info3
integer mm,ipivuhf(nvirtal*nal+nvirtbe*nbe),nos,lambda,sigma
integer ipivrohf((nvirt+nos)*(nocc+nos))
integer ipivrohf2((nvirt+nos)*(nocc+nos)-nos*nos)
integer ncent
real*8 r8heap(*),fa,fi,fe,fm,fj,fp,fq,chfx,pggen(3,3,nir)
real*8 dscfa(nbasis,nbasis),dscfaao(nbasis,nbasis)
real*8 dscfb(nbasis,nbasis),dscfbao(nbasis,nbasis)
real*8 dscftao(nbasis,nbasis),efield(3,ncent),deco(natoms)
real*8 ca(nbasis,nbasis),s(nbasis,nbasis),sum,fd,fg,sum2
real*8 cb(nbasis,nbasis),focka(nbasis,nbasis)
real*8 sumscf(nbasis,nbasis)
real*8 xdip(nbasis,nbasis)
real*8 fockb(nbasis,nbasis)
real*8 da(nbasis,nbasis),da2(nbasis,nbasis)
real*8 db(nbasis,nbasis),db2(nbasis,nbasis)
real*8 fockmoa(nbasis,nbasis),fockmob(nbasis,nbasis)
real*8 iprimea(nbasis,nbasis),iprimeb(nbasis,nbasis)
real*8 xinta(nvirtal,nal),xintb(nvirtbe,nbe)
real*8 xintrohf(nvirt+nos,nocc+nos)
real*8 gmtx(nvirtal,nal,nvirtal,nal)
real*8 xintrohf2(nvirt+nos,nocc+nos)
real*8 xintrohf_2((nvirt+nos)*(nocc+nos))
real*8 xintrohf2_2((nvirt+nos)*(nocc+nos)),xintaoa(nbasis,nbasis)
real*8 xintaob(nbasis,nbasis),xintfa(nal-ncore,ncore)
real*8 xintfb(nbe-ncore,ncore),xintfrohf(nocc-ncore,ncore)
real*8 gmtxuhf(nvirtal*nal+nvirtbe*nbe,nvirtal*nal+nvirtbe*nbe)
real*8 gmtxrohf(nvirt+nos,nocc+nos,nvirt+nos,nocc+nos)
real*8 gmtxrohf_2((nvirt+nos)*(nocc+nos),(nvirt+nos)*(nocc+nos))
real*8 daoa(nbasis,nbasis),daob(nbasis,nbasis),dipx,dipy,dipz
real*8 nucdip(3),nucsec(3,3),ha(nbasis,nbasis)
real*8 hb(nbasis,nbasis),qmtx(3,3),zvect(nvirtal*nal+nvirtbe*nbe)
real*8 xvect(nvirtal*nal+nvirtbe*nbe)
real*8 wr((nocc+nos)*(nvirt+nos)-(nos*nos))
real*8 wi((nocc+nos)*(nvirt+nos)-(nos*nos))
real*8 vl(1,(nocc+nos)*(nvirt+nos)-(nos*nos))
real*8 vr(1,(nocc+nos)*(nvirt+nos)-(nos*nos))
real*8 work(4*((nocc+nos)*(nvirt+nos)-(nos*nos)))
real*8 g2aa(nbasis,nbasis,nbasis,maxw)
real*8 g2bb(nbasis,nbasis,nbasis,maxw)
real*8 g2ab(nbasis,nbasis,nbasis,maxw)
real*8 gaoa(nbasis,nbasis,nbasis,maxw)
real*8 gaob(nbasis,nbasis,nbasis,maxw)
real*8 gaoab(nbasis,nbasis,nbasis,maxw)
real*8 g2int(nbasis,nbasis,nbasis,maxw)
real*8 g2int2(nbasis,nbasis,nbasis,maxw)
real*8 vint2aa(nbasis,nbasis,nbasis,maxw)
real*8 vint2bb(nbasis,nbasis,nbasis,maxw)
real*8 vint2ab(nbasis,nbasis,nbasis,maxw)
integer*4 i4heap(*)
character*2 atsymbol(natoms)
character*5 scftype
character*6 core
character*8 grdens,popul,qmmm
character*15 cscr15
character*16 calctype
character*32 dft
C SCF gradient
integer natoms,datoms,xyz,atnum(ncent),mu,nu,iatoms
real*8 sder(nbasis,nbasis,3),hder(nbasis,nbasis,3),coord(3,ncent)
real*8 tdera(nbasis,nbasis,3),tderb(nbasis,nbasis,3),vnn,rab3,sum1
real*8 tdera2(nbasis,nbasis,3),tderb2(nbasis,nbasis,3)
real*8 sum3,grads(3,natoms,5),tegrad(3),tegrad2(3),sum4,ddot
C QMMM stuff
real*8 finalener,selfener
C Deco
integer teintf(10),ipos,maxm,nninteg,iposlo,npos,ip,imem1,ninteg
integer intpos(3,(nbasis+1)*nbasis/2) !Equivalent to tdera
integer bfat(nbasis) !Equivalent to tdera2
integer*2 ii2,kk2
real*8 nucrepat,atchg(ncent),sum5,sum6
common/memcom/ imem1
C core issues...
call getkey('core',4,core,6)
if(core.eq.'corr ') ncore=0
nal=nal-ncore
nbe=nbe-ncore
nbasis=nbasis-ncore
nvirtal=nbasis-nal
nvirtbe=nbasis-nbe
c nocc=min(nal,nbe)
nvirt=min(nvirtal,nvirtbe)
call getvar('coord ',coord)
call getvar('atnum ',atnum)
call getvar('chfx ',chfx)
if(qmmm.eq.'amber') then
call getenergy(finalener,cscr15)
call getvar('selfenergy',selfener)
endif
C Set factors for UHF/RHF
if(scftype.eq.'uhf ') then
fd=0.5d0
fg=0.25d0
else
fd=1.d0
fg=1.d0
endif
C Get nuclear contributions
call getvar('nucdip ',nucdip)
call getvar('nucsec ',nucsec)
C Read SCF densities
call dfillzero(dscfaao,(nbasis+ncore)**2)
open(scrfile2,file='SCFDENSITIES',form='UNFORMATTED')
rewind(scrfile2)
call roeint(r8heap,i4heap,dscfaao,scrfile2,nbasis+ncore)
c Cases UHF and ROHF
if(scftype.ne.'rhf ') then
call dfillzero(dscfbao,(nbasis+ncore)**2)
call roeint(r8heap,i4heap,dscfbao,scrfile2,nbasis+ncore)
endif
c Create total AO density
if(scftype.eq.'uhf ') then
call dcopy((nbasis+ncore)**2,dscfaao,1,dscftao,1)
call daxpy((nbasis+ncore)**2,1.d0,dscfbao,1,dscftao,1)
endif
C Read MO coefficients from mocoeffile
call dfillzero(ca,(nbasis+ncore)**2)
open(mocoeffile,file='MOCOEF',form='UNFORMATTED')
rewind(mocoeffile)
read(mocoeffile) nn,(r8heap(jj),jj=1,2*nn)
do ii=1,nn
i=i4heap((ii-1)*4+3)
j=i4heap((ii-1)*4+4)
ca(i,j)=r8heap(ii*2-1)
enddo
if(scftype.eq.'uhf ') then
call dfillzero(cb,(nbasis+ncore)**2)
read(mocoeffile) nn,(r8heap(jj),jj=1,2*nn)
do ii=1,nn
i=i4heap((ii-1)*4+3)
j=i4heap((ii-1)*4+4)
cb(i,j)=r8heap(ii*2-1)
enddo
endif
close(mocoeffile)
C Read overlap integrals
open(oeintfile,file='OEINT',form='UNFORMATTED')
call roeint(r8heap,i4heap,s,oeintfile,nbasis+ncore)
close(oeintfile)
C Read AO Fock-matrix
open(scrfile1,file='FOCK',form='UNFORMATTED')
rewind(scrfile1)
call roeint(r8heap,i4heap,focka,scrfile1,nbasis+ncore)
if(scftype.eq.'uhf ')
$call roeint(r8heap,i4heap,fockb,scrfile1,nbasis+ncore)
close(scrfile1)
C Read one- and two-electron integrals !szemet
c Only in correlated calculations
if(calctype.ne.'scf ') then
call dfillzero(ha,(nbasis+ncore)**2)
C frozen core
if(core.ne.'corr '.and.ncore.ne.0) then
open(unit=scrfile1,file='PLUSCORE.55')
rewind(scrfile1)
C RHF and ROHF cases
if(scftype.ne.'uhf ') then
1016 read(scrfile1,*) sum,i,j,k,l
if(l.eq.0) then
ha(i,j)=sum
ha(j,i)=sum
goto 1016
endif
1017 read(scrfile1,*) sum,i,j,k,l
if(l.ne.0) then
goto 1017
endif
else
C UHF case
C aaaa list
1018 read(scrfile1,*) sum,i,j,k,l
if(i.ne.0) then
ha(i,j)=sum
ha(j,i)=sum
goto 1018
endif
1019 read(scrfile1,*) sum,i,j,k,l
if(l.ne.0) then
goto 1019
endif
call dfillzero(hb,(nbasis+ncore)**2)
C bbbb list
1020 read(scrfile1,*) sum,i,j,k,l
if(i.ne.0) then
hb(i,j)=sum
hb(j,i)=sum
goto 1020
endif
1021 read(scrfile1,*) sum,i,j,k,l
if(l.ne.0) then
goto 1021
endif
c abab list
1022 read(scrfile1,*) sum,k,l,i,j
if(l.ne.0) then
goto 1022
endif
endif
nuc=sum
close(scrfile1)
else
C correlated core
read(inp,*) i,j
read(inp,*) (i,j=1,nbasis)
read(inp,*) i
C aaaa list
1004 read(inp,*) sum,i,j,k,l
if(l.ne.0) then
goto 1004
endif
if(scftype.eq.'uhf ') then
call dfillzero(hb,nbasis**2)
C bbbb list
1008 read(inp,*) sum,i,j,k,l
if(l.ne.0) then
goto 1008
endif
c abab list
1009 read(inp,*) sum,i,j,k,l
if(l.ne.0) then
goto 1009
endif
endif
C aa list
if(scftype.ne.'uhf ') then
ha(i,j)=sum
ha(j,i)=sum
endif
1007 read(inp,*) sum,i,j,k,l
if(i.ne.0) then
ha(i,j)=sum
ha(j,i)=sum
goto 1007
endif
if(scftype.eq.'uhf ') then
1010 read(inp,*) sum,i,j,k,l
if(i.ne.0) then
hb(i,j)=sum
hb(j,i)=sum
goto 1010
endif
read(inp,*) sum,i,j,k,l
endif
nuc=sum
close(inp)
endif
C Read density matrices !szemet
call dfillzero(da,(nbasis+ncore)**2)
open(unit=densfile,file='CCDENSITIES')
rewind(densfile)
c aaaa list
1005 read(densfile,*) sum,i,k,j,l
if(l.ne.0) then
goto 1005
endif
if(scftype.eq.'uhf ') then
c bbbb list
1011 read(densfile,*) sum,i,k,j,l
if(l.ne.0) then
goto 1011
endif
C abab list
1012 read(densfile,*) sum,i,k,j,l
if(l.ne.0) then
goto 1012
endif
C Ujabb verzio rrr
C aa list
call dfillzero(db,(nbasis+ncore)**2)
1013 read(densfile,*) sum,i,k,j,l
da(i,k)=sum
da(k,i)=sum
if(i.gt.0) then
goto 1013
else
db(i,k)=sum
db(k,i)=sum
endif
C bb list
do
read(densfile,*,end=1006) sum,i,k,j,l
db(i,k)=sum
db(k,i)=sum
enddo
1006 continue
close(densfile)
else
da(i,k)=sum
da(k,i)=sum
do
read(densfile,*,end=1014) sum,i,k,j,l
da(i,k)=sum
da(k,i)=sum
enddo
1014 continue
close(densfile)
endif
C **********************************************************************
C frozen core - shift the density (add the core)
C **********************************************************************
if(core.ne.'corr '.and.ncore.ne.0) then
C one electron
call dfillzero(da2,(nbasis+ncore)**2)
do p=1,nbasis
do q=1,nbasis
da2(p+ncore,q+ncore)=da(p,q)
enddo
enddo
do f=1,ncore
da2(f,f)=2.d0*fd
enddo
call dcopy((nbasis+ncore)**2,da2,1,da,1)
C UHF case
if(scftype.eq.'uhf ') then
call dfillzero(db2,(nbasis+ncore)**2)
do p=1,nbasis
do q=1,nbasis
db2(p+ncore,q+ncore)=db(p,q)
enddo
enddo
do f=1,ncore
db2(f,f)=1.d0
enddo
call dcopy((nbasis+ncore)**2,db2,1,db,1)
endif
endif
c end of correlated section
endif
C Transform aa-Fock-matrix into MO-basis
call dgemm('N','N',nbasis+ncore,nbasis+ncore,nbasis+ncore,1.d0,
$focka,nbasis+ncore,ca,nbasis+ncore,0.d0,da2,nbasis+ncore)
call dgemm('T','N',nbasis+ncore,nbasis+ncore,nbasis+ncore,1.d0,
$ca,nbasis+ncore,da2,nbasis+ncore,0.d0,fockmoa,nbasis+ncore)
C Transform bb-Fock-matrix into MO-basis
if(scftype.eq.'uhf ') then
call dgemm('N','N',nbasis+ncore,nbasis+ncore,nbasis+ncore,1.d0,
$fockb,nbasis+ncore,cb,nbasis+ncore,0.d0,db2,nbasis+ncore)
call dgemm('T','N',nbasis+ncore,nbasis+ncore,nbasis+ncore,1.d0,
$cb,nbasis+ncore,db2,nbasis+ncore,0.d0,fockmob,nbasis+ncore)
endif
C correlated section
if(calctype.ne.'scf ') then
C **********************************************************************
C Build the "I prime" intermediate
C **********************************************************************
c reseting the matrices
call dfillzero(iprimea,(nbasis+ncore)**2)
if(scftype.eq.'rhf ')
$call dfillzero(gmtx,(nvirtal*(nal+ncore))**2)
if(scftype.eq.'uhf ') then
call dfillzero(iprimeb,(nbasis+ncore)**2)
call dfillzero(gmtxuhf,(nvirtal*(nal+ncore)+nvirtbe*(nbe+ncore))
$**2)
endif
if(scftype.eq.'rohf ')
$call dfillzero(gmtxrohf,(nvirt+nos)**2*(nocc+nos)**2)
nosshift=0
c one electron part*****************************************************
c RHF case
do p=1,nbasis+ncore
do q=1,nbasis+ncore
sum=0.d0
do t=1,nbasis+ncore
sum=sum-0.5d0*
$(ha(p,t)*da(q,t)+ha(t,p)*da(t,q))
enddo
iprimea(p,q)=sum
enddo
enddo
c UHF case
if(scftype.eq.'uhf ') then
do p=1,nbasis+ncore
do q=1,nbasis+ncore
sum=0.d0
do t=1,nbasis+ncore
sum=sum-0.5d0*
$(hb(p,t)*db(q,t)+hb(t,p)*db(t,q))
enddo
iprimeb(p,q)=sum
enddo
enddo
endif
c two electron part*****************************************************
c open the density files and integral lists
if(core.ne.'corr '.and.ncore.ne.0)
$ open(unit=scrfile1,file='PLUSCORE.55')
open(unit=densfile,file='CCDENSITIES')
c loop on the blocks
do ii=1,(nbasis+ncore)/maxw+1
c block sizes
c RHF case
if(ii*maxw.le.nbasis+ncore) then
nblock=maxw
else
nblock=mod(nbasis+ncore,maxw)
endif
c UHF case
c alpha
if(ii*maxw.le.nal+ncore) then
nblockal=maxw
else
if((ii-1)*maxw.gt.nal+ncore) then
nblockal=0
else
nblockal=mod(nal+ncore,maxw)
endif
endif
c beta
if(ii*maxw.le.nbe+ncore) then
nblockbe=maxw
else
if((ii-1)*maxw.gt.nbe+ncore) then
nblockbe=0
else
nblockbe=mod(nbe+ncore,maxw)
endif
endif
c ***read the integrals and densities needed***
C INTEGRALS
call dfillzero(vint2aa,maxw*(nbasis+ncore)**3)
C frozen core
if(core.ne.'corr '.and.ncore.ne.0) then
rewind(scrfile1)
C RHF and ROHF cases
if(scftype.ne.'uhf ') then
1032 read(scrfile1,*) sum2,i,j,k,l
if(l.eq.0) goto 1032
if((i.gt.(ii-1)*maxw).and.(i.le.(ii-1)*maxw+nblock)) then
vint2aa(l,j,k,i-(ii-1)*maxw)=sum2
vint2aa(k,j,l,i-(ii-1)*maxw)=sum2
endif
if((j.gt.(ii-1)*maxw).and.(j.le.(ii-1)*maxw+nblock)) then
vint2aa(k,i,l,j-(ii-1)*maxw)=sum2
vint2aa(l,i,k,j-(ii-1)*maxw)=sum2
endif
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblock)) then
vint2aa(j,l,i,k-(ii-1)*maxw)=sum2
vint2aa(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblock)) then
vint2aa(i,k,j,l-(ii-1)*maxw)=sum2
vint2aa(j,k,i,l-(ii-1)*maxw)=sum2
endif
1033 read(scrfile1,*) sum2,i,j,k,l
if((i.gt.(ii-1)*maxw).and.(i.le.(ii-1)*maxw+nblock)) then
vint2aa(l,j,k,i-(ii-1)*maxw)=sum2
vint2aa(k,j,l,i-(ii-1)*maxw)=sum2
endif
if((j.gt.(ii-1)*maxw).and.(j.le.(ii-1)*maxw+nblock)) then
vint2aa(k,i,l,j-(ii-1)*maxw)=sum2
vint2aa(l,i,k,j-(ii-1)*maxw)=sum2
endif
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblock)) then
vint2aa(j,l,i,k-(ii-1)*maxw)=sum2
vint2aa(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblock)) then
vint2aa(i,k,j,l-(ii-1)*maxw)=sum2
vint2aa(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1033
else
C UHF case
C aaaa list
1034 read(scrfile1,*) sum2,i,j,k,l
if(i.ne.0) goto 1034
1035 read(scrfile1,*) sum2,i,j,k,l
if((i.gt.(ii-1)*maxw).and.(i.le.(ii-1)*maxw+nblock)) then
vint2aa(l,j,k,i-(ii-1)*maxw)=sum2
vint2aa(k,j,l,i-(ii-1)*maxw)=sum2
endif
if((j.gt.(ii-1)*maxw).and.(j.le.(ii-1)*maxw+nblock)) then
vint2aa(k,i,l,j-(ii-1)*maxw)=sum2
vint2aa(l,i,k,j-(ii-1)*maxw)=sum2
endif
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblock)) then
vint2aa(j,l,i,k-(ii-1)*maxw)=sum2
vint2aa(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblock)) then
vint2aa(i,k,j,l-(ii-1)*maxw)=sum2
vint2aa(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1035
call dfillzero(vint2bb,maxw*(nbasis+ncore)**3)
call dfillzero(vint2ab,maxw*(nbasis+ncore)**3)
C bbbb list
1036 read(scrfile1,*) sum2,i,j,k,l
if(i.ne.0) goto 1036
1037 read(scrfile1,*) sum2,i,j,k,l
if((i.gt.(ii-1)*maxw).and.(i.le.(ii-1)*maxw+nblock)) then
vint2bb(l,j,k,i-(ii-1)*maxw)=sum2
vint2bb(k,j,l,i-(ii-1)*maxw)=sum2
endif
if((j.gt.(ii-1)*maxw).and.(j.le.(ii-1)*maxw+nblock)) then
vint2bb(k,i,l,j-(ii-1)*maxw)=sum2
vint2bb(l,i,k,j-(ii-1)*maxw)=sum2
endif
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblock)) then
vint2bb(j,l,i,k-(ii-1)*maxw)=sum2
vint2bb(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblock)) then
vint2bb(i,k,j,l-(ii-1)*maxw)=sum2
vint2bb(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1037
c abab list
1038 read(scrfile1,*) sum2,k,l,i,j
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblock)) then
vint2ab(j,l,i,k-(ii-1)*maxw)=sum2
vint2ab(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblock)) then
vint2ab(i,k,j,l-(ii-1)*maxw)=sum2
vint2ab(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1038
endif
nuc=sum2
else
C correlated core
rewind(inp)
read(inp,*) i,j
read(inp,*) (i,j=1,nbasis+ncore)
read(inp,*) i
C aaaa list
1023 read(inp,*) sum2,i,j,k,l
if((i.gt.(ii-1)*maxw).and.(i.le.(ii-1)*maxw+nblock)) then
vint2aa(l,j,k,i-(ii-1)*maxw)=sum2
vint2aa(k,j,l,i-(ii-1)*maxw)=sum2
endif
if((j.gt.(ii-1)*maxw).and.(j.le.(ii-1)*maxw+nblock)) then
vint2aa(k,i,l,j-(ii-1)*maxw)=sum2
vint2aa(l,i,k,j-(ii-1)*maxw)=sum2
endif
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblock)) then
vint2aa(j,l,i,k-(ii-1)*maxw)=sum2
vint2aa(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblock)) then
vint2aa(i,k,j,l-(ii-1)*maxw)=sum2
vint2aa(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1023
c UHF case
if(scftype.eq.'uhf ') then
call dfillzero(vint2bb,maxw*(nbasis+ncore)**3)
call dfillzero(vint2ab,maxw*(nbasis+ncore)**3)
C bbbb list
1025 read(inp,*) sum2,i,j,k,l
if((i.gt.(ii-1)*maxw).and.(i.le.(ii-1)*maxw+nblock)) then
vint2bb(l,j,k,i-(ii-1)*maxw)=sum2
vint2bb(k,j,l,i-(ii-1)*maxw)=sum2
endif
if((j.gt.(ii-1)*maxw).and.(j.le.(ii-1)*maxw+nblock)) then
vint2bb(k,i,l,j-(ii-1)*maxw)=sum2
vint2bb(l,i,k,j-(ii-1)*maxw)=sum2
endif
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblock)) then
vint2bb(j,l,i,k-(ii-1)*maxw)=sum2
vint2bb(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblock)) then
vint2bb(i,k,j,l-(ii-1)*maxw)=sum2
vint2bb(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1025
c abab list
1026 read(inp,*) sum2,i,j,k,l
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblock)) then
vint2ab(j,l,i,k-(ii-1)*maxw)=sum2
vint2ab(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblock)) then
vint2ab(i,k,j,l-(ii-1)*maxw)=sum2
vint2ab(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1026
endif
endif
c DENSITIES
call dfillzero(g2aa,maxw*(nbasis+ncore)**3)
rewind(densfile)
c aaaa list
1024 read(densfile,*) sum2,i,k,j,l
if((i+ncore.gt.(ii-1)*maxw).and.
$(i+ncore.le.(ii-1)*maxw+nblock)) then
g2aa(l+ncore,j+ncore,k+ncore,i-(ii-1)*maxw+ncore)=sum2
g2aa(k+ncore,j+ncore,l+ncore,i-(ii-1)*maxw+ncore)=sum2
endif
if((j+ncore.gt.(ii-1)*maxw).and.
$(j+ncore.le.(ii-1)*maxw+nblock)) then
g2aa(k+ncore,i+ncore,l+ncore,j-(ii-1)*maxw+ncore)=sum2
g2aa(l+ncore,i+ncore,k+ncore,j-(ii-1)*maxw+ncore)=sum2
endif
if((k+ncore.gt.(ii-1)*maxw).and.
$(k+ncore.le.(ii-1)*maxw+nblock)) then
g2aa(j+ncore,l+ncore,i+ncore,k-(ii-1)*maxw+ncore)=sum2
g2aa(i+ncore,l+ncore,j+ncore,k-(ii-1)*maxw+ncore)=sum2
endif
if((l+ncore.gt.(ii-1)*maxw).and.
$(l+ncore.le.(ii-1)*maxw+nblock)) then
g2aa(i+ncore,k+ncore,j+ncore,l-(ii-1)*maxw+ncore)=sum2
g2aa(j+ncore,k+ncore,i+ncore,l-(ii-1)*maxw+ncore)=sum2
endif
if(l.ne.0) goto 1024
c build the frozen parts
if(core.ne.'corr '.and.ncore.ne.0) then
do f=1,ncore
if((scftype.ne.'uhf ').and.(f.gt.(ii-1)*maxw).and.
$(f.le.(ii-1)*maxw+nblock)) g2aa(f,f,f,f-(ii-1)*maxw)=2.d0
enddo
do f=1,ncore
do g=1,ncore
if(g.ne.f) then
if((f.gt.(ii-1)*maxw).and.(f.le.(ii-1)*maxw+nblock))
$g2aa(g,f,g,f-(ii-1)*maxw)=4.d0*fg
if((g.gt.(ii-1)*maxw).and.(g.le.(ii-1)*maxw+nblock)) then
g2aa(g,f,f,g-(ii-1)*maxw)=-1.d0*fd
g2aa(f,f,g,g-(ii-1)*maxw)=-1.d0*fd
endif
endif
enddo
enddo
do f=1,ncore
do p=ncore+1,ncore+nbasis
do q=ncore+1,ncore+nbasis
if((f.gt.(ii-1)*maxw).and.(f.le.(ii-1)*maxw+nblock)) then
g2aa(p,f,q,f-(ii-1)*maxw)=2.d0*da(p,q)*fd
g2aa(p,q,f,f-(ii-1)*maxw)=-0.5d0*da(p,q)
g2aa(f,q,p,f-(ii-1)*maxw)=-0.5d0*da(p,q)
endif
if((q.gt.(ii-1)*maxw).and.(q.le.(ii-1)*maxw+nblock)) then
g2aa(f,p,f,q-(ii-1)*maxw)=2.d0*da(p,q)*fd
g2aa(p,f,f,q-(ii-1)*maxw)=-0.5d0*da(p,q)
g2aa(f,f,p,q-(ii-1)*maxw)=-0.5d0*da(p,q)
endif
enddo
enddo
enddo
endif
c UHF case
if(scftype.eq.'uhf ') then
call dfillzero(g2bb,maxw*(nbasis+ncore)**3)
call dfillzero(g2ab,maxw*(nbasis+ncore)**3)
c bbbb list
1027 read(densfile,*) sum2,i,k,j,l
if((i+ncore.gt.(ii-1)*maxw).and.
$(i+ncore.le.(ii-1)*maxw+nblock)) then
g2bb(l+ncore,j+ncore,k+ncore,i-(ii-1)*maxw+ncore)=sum2
g2bb(k+ncore,j+ncore,l+ncore,i-(ii-1)*maxw+ncore)=sum2
endif
if((j+ncore.gt.(ii-1)*maxw).and.
$(j+ncore.le.(ii-1)*maxw+nblock)) then
g2bb(k+ncore,i+ncore,l+ncore,j-(ii-1)*maxw+ncore)=sum2
g2bb(l+ncore,i+ncore,k+ncore,j-(ii-1)*maxw+ncore)=sum2
endif
if((k+ncore.gt.(ii-1)*maxw).and.
$(k+ncore.le.(ii-1)*maxw+nblock)) then
g2bb(j+ncore,l+ncore,i+ncore,k-(ii-1)*maxw+ncore)=sum2
g2bb(i+ncore,l+ncore,j+ncore,k-(ii-1)*maxw+ncore)=sum2
endif
if((l+ncore.gt.(ii-1)*maxw).and.
$(l+ncore.le.(ii-1)*maxw+nblock)) then
g2bb(i+ncore,k+ncore,j+ncore,l-(ii-1)*maxw+ncore)=sum2
g2bb(j+ncore,k+ncore,i+ncore,l-(ii-1)*maxw+ncore)=sum2
endif
if(l.ne.0) goto 1027
c build the frozen parts
if(core.ne.'corr '.and.ncore.ne.0) then
do f=1,ncore
do g=1,ncore
if(g.ne.f) then
if((f.gt.(ii-1)*maxw).and.(f.le.(ii-1)*maxw+nblock))
$g2bb(g,f,g,f-(ii-1)*maxw)=1.d0
if((g.gt.(ii-1)*maxw).and.(g.le.(ii-1)*maxw+nblock)) then
g2bb(g,f,f,g-(ii-1)*maxw)=-0.5d0
g2bb(f,f,g,g-(ii-1)*maxw)=-0.5d0
endif
endif
enddo
enddo
do f=1,ncore
do p=ncore+1,ncore+nbasis
do q=ncore+1,ncore+nbasis
if((f.gt.(ii-1)*maxw).and.(f.le.(ii-1)*maxw+nblock)) then
g2bb(p,f,q,f-(ii-1)*maxw)=db(p,q)
g2bb(p,q,f,f-(ii-1)*maxw)=-0.5d0*db(p,q)
g2bb(f,q,p,f-(ii-1)*maxw)=-0.5d0*db(p,q)
endif
if((q.gt.(ii-1)*maxw).and.(q.le.(ii-1)*maxw+nblock)) then
g2bb(f,p,f,q-(ii-1)*maxw)=db(p,q)
g2bb(p,f,f,q-(ii-1)*maxw)=-0.5d0*db(p,q)
g2bb(f,f,p,q-(ii-1)*maxw)=-0.5d0*db(p,q)
endif
enddo
enddo
enddo
endif
C abab list
1028 read(densfile,*) sum2,i,k,j,l
if((k+ncore.gt.(ii-1)*maxw).and.
$(k+ncore.le.(ii-1)*maxw+nblock)) then
g2ab(j+ncore,l+ncore,i+ncore,k-(ii-1)*maxw+ncore)=sum2
g2ab(i+ncore,l+ncore,j+ncore,k-(ii-1)*maxw+ncore)=sum2
endif
if((l+ncore.gt.(ii-1)*maxw).and.
$(l+ncore.le.(ii-1)*maxw+nblock)) then
g2ab(i+ncore,k+ncore,j+ncore,l-(ii-1)*maxw+ncore)=sum2
g2ab(j+ncore,k+ncore,i+ncore,l-(ii-1)*maxw+ncore)=sum2
endif
if(l.ne.0) goto 1028
call dscal(maxw*(nbasis+ncore)**3,fd,g2ab,1)
c build the frozen parts
if(core.ne.'corr '.and.ncore.ne.0) then
do f=1,ncore
if((f.gt.(ii-1)*maxw).and.(f.le.(ii-1)*maxw+nblock))
$g2ab(f,f,f,f-(ii-1)*maxw)=1.d0
enddo
do f=1,ncore
do g=1,ncore
if((g.ne.f).and.(f.gt.(ii-1)*maxw).and.
$(f.le.(ii-1)*maxw+nblock)) g2ab(g,f,g,f-(ii-1)*maxw)=1.d0
enddo
enddo
do f=1,ncore
do p=ncore+1,ncore+nbasis
do q=ncore+1,ncore+nbasis
if((f.gt.(ii-1)*maxw).and.(f.le.(ii-1)*maxw+nblock))
$g2ab(p,f,q,f-(ii-1)*maxw)=da(p,q)
if((q.gt.(ii-1)*maxw).and.(q.le.(ii-1)*maxw+nblock))
$g2ab(f,p,f,q-(ii-1)*maxw)=db(p,q)
enddo
enddo
enddo
endif
endif
c BUILD THE MATRIX
C UHF case
if(scftype.eq.'uhf ') then
C I prime-aa
call dgemm('n','t',nbasis+ncore,nbasis+ncore,
$nblock*(nbasis+ncore)**2,-1.d0,vint2aa,nbasis+ncore,g2aa,
$nbasis+ncore,1.d0,iprimea,nbasis+ncore)
call dgemm('n','t',nbasis+ncore,nbasis+ncore,
$nblock*(nbasis+ncore)**2,-1.d0,vint2ab,nbasis+ncore,g2ab,
$nbasis+ncore,1.d0,iprimea,nbasis+ncore)
C I prime-bb
call dgemm('n','t',nbasis+ncore,nbasis+ncore,
$nblock*(nbasis+ncore)**2,-1.d0,vint2bb,nbasis+ncore,g2bb,
$nbasis+ncore,1.d0,iprimeb,nbasis+ncore)
do t=1,nblock
do z=1,nbasis+ncore
call dgemm('t','n',nbasis+ncore,nbasis+ncore,
$nbasis+ncore,-1.d0,vint2ab(1,1,z,t),nbasis+ncore,g2ab(1,1,z,t),
$nbasis+ncore,1.d0,iprimeb,nbasis+ncore)
enddo
enddo
else
C RHF case
C I prime-aa
call dgemm('n','t',nbasis+ncore,nbasis+ncore,
$nblock*(nbasis+ncore)**2,-1.d0,vint2aa,nbasis+ncore,g2aa,
$nbasis+ncore,1.d0,iprimea,nbasis+ncore)
endif
C **********************************************************************
C Build the "G" hypermatrix of the Z-vector method
C **********************************************************************
C RHF case
if(scftype.eq.'rhf ') then
do a=1,nvirtal
do i=1,nal+ncore
c the last index is virtual
if((ii-1)*maxw.ge.nal+ncore) then
do m=1,nal+ncore
do e=1,nblock
gmtx(a,i,e+(ii-1)*maxw-nal-ncore,m)
$=gmtx(a,i,e+(ii-1)*maxw-nal-ncore,m)+4.d0*
$vint2aa(a+nal+ncore,m,i,e)-
$vint2aa(a+nal+ncore,i,m,e)-
$vint2aa(i,a+nal+ncore,m,e)
enddo
enddo
endif
c the last index is either virtual or not
if(((ii-1)*maxw.lt.nal+ncore).and.((ii-1)*maxw+nblock.gt.
$nal+ncore)) then
noccblock=nal+ncore-(ii-1)*maxw
c virtual
do m=1,nal+ncore
do e=noccblock+1,nblock
gmtx(a,i,e+(ii-1)*maxw-nal-ncore,m)
$=gmtx(a,i,e+(ii-1)*maxw-nal-ncore,m)+4.d0*
$vint2aa(a+nal+ncore,m,i,e)-
$vint2aa(a+nal+ncore,i,m,e)-
$vint2aa(i,a+nal+ncore,m,e)
enddo
enddo
endif
enddo
enddo
endif
C UHF case
if(scftype.eq.'uhf ') then
c aaaa-block
c the last index is not virtual
if((ii-1)*maxw+nblockal.le.nal+ncore) then
q=1+(ii-1)*maxw*nvirtal
do m=1,nblockal
do e=1,nvirtal
p=1
do i=1,nal+ncore
do a=1,nvirtal
gmtxuhf(p,q)=gmtxuhf(p,q)+2.d0*
$vint2aa(a+nal+ncore,e+nal+ncore,i,m)-
$vint2aa(i,a+nal+ncore,e+nal+ncore,m)-
$vint2aa(a+nal+ncore,i,e+nal+ncore,m)
p=p+1
enddo
enddo
q=q+1
enddo
enddo
endif
c the last index is either virtual or not
if(((ii-1)*maxw.lt.nal+ncore).and.((ii-1)*maxw+nblockal.gt.
$nal+ncore)) then
noccblock=nal+ncore-(ii-1)*maxw
c not virtual
q=1+(ii-1)*maxw*nvirtal
do m=1,noccblock
do e=1,nvirtal
p=1
do i=1,nal+ncore
do a=1,nvirtal
gmtxuhf(p,q)=gmtxuhf(p,q)+2.d0*
$vint2aa(a+nal+ncore,e+nal+ncore,i,m)-
$vint2aa(i,a+nal+ncore,e+nal+ncore,m)-
$vint2aa(a+nal+ncore,i,e+nal+ncore,m)
p=p+1
enddo
enddo
q=q+1
enddo
enddo
endif
C aabb-block
c the last index is not virtual
if((ii-1)*maxw+nblockbe.le.nbe+ncore) then
q=1+nvirtal*(nal+ncore)+(ii-1)*maxw*nvirtbe
do m=1,nblockbe
do e=1,nvirtbe
p=1
do i=1,nal+ncore
do a=1,nvirtal
gmtxuhf(p,q)=gmtxuhf(p,q)+2.d0*
$vint2ab(a+nal+ncore,e+nbe+ncore,i,m)
p=p+1
enddo
enddo
q=q+1
enddo
enddo
endif
c the last index is either virtual or not
if(((ii-1)*maxw.lt.nbe+ncore).and.((ii-1)*maxw+nblockbe.gt.
$nbe+ncore)) then
noccblock=nbe+ncore-(ii-1)*maxw
c not virtual
q=1+nvirtal*(nal+ncore)+(ii-1)*maxw*nvirtbe
do m=1,noccblock
do e=1,nvirtbe
p=1
do i=1,nal+ncore
do a=1,nvirtal
gmtxuhf(p,q)=gmtxuhf(p,q)+2.d0*
$vint2ab(a+nal+ncore,e+nbe+ncore,i,m)
p=p+1
enddo
enddo
q=q+1
enddo
enddo
endif
c bbaa-block
c the last index is not virtual
if((ii-1)*maxw+nblockbe.le.nbe+ncore) then
p=1+nvirtal*(nal+ncore)+(ii-1)*maxw*nvirtbe
do i=1,nblockbe
do a=1,nvirtbe
q=1
do m=1,nal+ncore
do e=1,nvirtal
gmtxuhf(p,q)=gmtxuhf(p,q)+2.d0*
$vint2ab(e+nal+ncore,a+nbe+ncore,m,i)
q=q+1
enddo
enddo
p=p+1
enddo
enddo
endif
c the last index is either virtual or not
if(((ii-1)*maxw.lt.nbe+ncore).and.((ii-1)*maxw+nblockbe.gt.
$nbe+ncore)) then
noccblock=nbe+ncore-(ii-1)*maxw
c not virtual
p=1+nvirtal*(nal+ncore)+(ii-1)*maxw*nvirtbe
do i=1,noccblock
do a=1,nvirtbe
q=1
do m=1,nal+ncore
do e=1,nvirtal
gmtxuhf(p,q)=gmtxuhf(p,q)+2.d0*
$vint2ab(e+nal+ncore,a+nbe+ncore,m,i)
q=q+1
enddo
enddo
p=p+1
enddo
enddo
endif
c bbbb-block
c the last index is not virtual
if((ii-1)*maxw+nblockbe.le.nbe+ncore) then
q=1+nvirtal*(nal+ncore)+(ii-1)*maxw*nvirtbe
do m=1,nblockbe
do e=1,nvirtbe
p=1+nvirtal*(nal+ncore)
do i=1,nbe+ncore
do a=1,nvirtbe
gmtxuhf(p,q)=gmtxuhf(p,q)+2.d0*
$vint2bb(a+nbe+ncore,e+nbe+ncore,i,m)-
$vint2bb(i,a+nbe+ncore,e+nbe+ncore,m)-
$vint2bb(a+nbe+ncore,i,e+nbe+ncore,m)
p=p+1
enddo
enddo
q=q+1
enddo
enddo
endif
c the last index is either virtual or not
if(((ii-1)*maxw.lt.nbe+ncore).and.((ii-1)*maxw+nblockbe.gt.
$nbe+ncore)) then
noccblock=nbe+ncore-(ii-1)*maxw
c not virtual
q=1+nvirtal*(nal+ncore)+(ii-1)*maxw*nvirtbe
do m=1,noccblock
do e=1,nvirtbe
p=1+nvirtal*(nal+ncore)
do i=1,nbe+ncore
do a=1,nvirtbe
gmtxuhf(p,q)=gmtxuhf(p,q)+2.d0*
$vint2bb(a+nbe+ncore,e+nbe+ncore,i,m)-
$vint2bb(i,a+nbe+ncore,e+nbe+ncore,m)-
$vint2bb(a+nbe+ncore,i,e+nbe+ncore,m)
p=p+1
enddo
enddo
q=q+1
enddo
enddo
endif
endif
C ROHF case
if(scftype.eq.'rohf ') then
c the last index is 'i'********
info=0
c the last index is not virtual
if((ii-1)*maxw+nblock.le.nocc+nos) then
noccblock=maxw
info=1
else
c the last index is either virtual or not
if(((ii-1)*maxw.lt.nocc+nos).and.((ii-1)*maxw+nblock.gt.
$nocc+nos)) then
noccblock=nocc+nos-(ii-1)*maxw
info=1
endif
endif
c BUILD THE MATRIX
if(info.eq.1) then
do a=1,nvirt+nos
do i=1,noccblock
do e=1,nvirt+nos
do m=1,nocc+nos
c os-os
if((e.le.nos.and.m.gt.nocc).or.(a.le.nos.and.i+(ii-1)*maxw
$.gt.nocc)) then
else
c virt-occ
if((e.gt.nos.and.m.le.nocc)) then
gmtxrohf(e,m,a,i+(ii-1)*maxw)=
$gmtxrohf(e,m,a,i+(ii-1)*maxw)+4.d0*vint2aa(e+nocc,a+nocc,m,i)-
$vint2aa(a+nocc,m,e+nocc,i)-vint2aa(a+nocc,e+nocc,m,i)
endif
c Occupation number
if(i+(ii-1)*maxw.le.nocc) then
fi=2.d0
else
fi=1.d0
endif
if(a.le.nos) then
fa=1.d0
else
fa=0.d0
endif
c virt-os
if((e.gt.nos.and.m.gt.nocc)) then
gmtxrohf(e,m,a,i+(ii-1)*maxw)=
$gmtxrohf(e,m,a,i+(ii-1)*maxw)+2.d0*vint2aa(e+nocc,a+nocc,m,i)
$-(3.d0-fa-fi)*0.5d0*(vint2aa(a+nocc,m,e+nocc,i)
$+vint2aa(a+nocc,e+nocc,m,i))
endif
c os-occ
if((e.le.nos.and.m.le.nocc)) then
gmtxrohf(e,m,a,i+(ii-1)*maxw)=
$gmtxrohf(e,m,a,i+(ii-1)*maxw)+2.d0*vint2aa(e+nocc,a+nocc,m,i)
$+(1.d0-fa-fi)*0.5d0*(vint2aa(a+nocc,m,e+nocc,i)
$+vint2aa(a+nocc,e+nocc,m,i))
endif
endif
enddo
enddo
enddo
enddo
endif
c the last index is 'z'********
info=0
c the last index is either doubly or half occupied (or virtual)
if(((ii-1)*maxw.lt.nocc).and.((ii-1)*maxw+nblock.gt.
$nocc)) then
nosblock=min(nos,(ii-1)*maxw+nblock-nocc)
nosshift=nocc-(ii-1)*maxw
info=1
else
c the last index is half occupied
if(((ii-1)*maxw.ge.nocc).and.((ii-1)*maxw+nblock.le.
$nocc+nos)) then
nosblock=maxw
nosshift=0
info=1
else
c the last index is either half occupied or virtual
if(((ii-1)*maxw.lt.nocc+nos).and.((ii-1)*maxw+nblock.gt.
$nocc+nos)) then
nosblock=nocc+nos-(ii-1)*maxw
nosshift=0
info=1
endif
endif
endif
if(info.eq.1) then
do a=1,nvirt+nos
do i=1,nocc+nos
do e=1,nvirt+nos
do m=1,nocc+nos
c os-os
if((e.le.nos.and.m.gt.nocc).or.(a.le.nos.and.i
$.gt.nocc)) then
else
c virt-occ
if((e.gt.nos.and.m.le.nocc)) then
if(a.eq.e) then
sum=0.d0
do z=nosshift+1,nosshift+nosblock
sum=sum+vint2aa(m,i,z+(ii-1)*maxw,z)
enddo
gmtxrohf(e,m,a,i)=gmtxrohf(e,m,a,i)+sum
endif
if(i.eq.e+nocc) then
sum=0.d0
do z=nosshift+1,nosshift+nosblock
sum=sum+vint2aa(a+nocc,m,z+(ii-1)*maxw,z)
enddo
gmtxrohf(e,m,a,i)=gmtxrohf(e,m,a,i)+sum
endif
if(a+nocc.eq.m) then
sum=0.d0
do z=nosshift+1,nosshift+nosblock
sum=sum+vint2aa(e+nocc,i,z+(ii-1)*maxw,z)
enddo
gmtxrohf(e,m,a,i)=gmtxrohf(e,m,a,i)+sum
endif
if(i.eq.m) then
sum=0.d0
do z=nosshift+1,nosshift+nosblock
sum=sum+vint2aa(a+nocc,e+nocc,z+(ii-1)*maxw,z)
enddo
gmtxrohf(e,m,a,i)=gmtxrohf(e,m,a,i)+sum
endif
endif
c virt-os
if((e.gt.nos.and.m.gt.nocc)) then
if(a.eq.e) then
sum=0.d0
do z=nosshift+1,nosshift+nosblock
sum=sum+0.5d0*vint2aa(m,i,z+(ii-1)*maxw,z)
enddo
gmtxrohf(e,m,a,i)=gmtxrohf(e,m,a,i)+sum
endif
if(i.eq.e+nocc) then
sum=0.d0
do z=nosshift+1,nosshift+nosblock
sum=sum+0.5d0*vint2aa(a+nocc,m,z+(ii-1)*maxw,z)
enddo
gmtxrohf(e,m,a,i)=gmtxrohf(e,m,a,i)+sum
endif
if(a+nocc.eq.m) then
sum=0.d0
do z=nosshift+1,nosshift+nosblock
sum=sum+0.5d0*vint2aa(e+nocc,i,z+(ii-1)*maxw,z)
enddo
gmtxrohf(e,m,a,i)=gmtxrohf(e,m,a,i)+sum
endif
if(i.eq.m) then
sum=0.d0
do z=nosshift+1,nosshift+nosblock
sum=sum+0.5d0*vint2aa(a+nocc,e+nocc,z+(ii-1)*maxw,z)
enddo
gmtxrohf(e,m,a,i)=gmtxrohf(e,m,a,i)+sum
endif
endif
c os-occ
if((e.le.nos.and.m.le.nocc)) then
if(a.eq.e) then
sum=0.d0
do z=nosshift+1,nosshift+nosblock
sum=sum+0.5d0*vint2aa(m,i,z+(ii-1)*maxw,z)
enddo
gmtxrohf(e,m,a,i)=gmtxrohf(e,m,a,i)+sum
endif
if(i.eq.e+nocc) then
sum=0.d0
do z=nosshift+1,nosshift+nosblock
sum=sum+0.5d0*vint2aa(a+nocc,m,z+(ii-1)*maxw,z)
enddo
gmtxrohf(e,m,a,i)=gmtxrohf(e,m,a,i)+sum
endif
if(a+nocc.eq.m) then
sum=0.d0
do z=nosshift+1,nosshift+nosblock
sum=sum+0.5d0*vint2aa(e+nocc,i,z+(ii-1)*maxw,z)
enddo
gmtxrohf(e,m,a,i)=gmtxrohf(e,m,a,i)+sum
endif
if(i.eq.m) then
sum=0.d0
do z=nosshift+1,nosshift+nosblock
sum=sum+0.5d0*vint2aa(a+nocc,e+nocc,z+(ii-1)*maxw,z)
enddo
gmtxrohf(e,m,a,i)=gmtxrohf(e,m,a,i)+sum
endif
endif
endif
enddo
enddo
enddo
enddo
endif
endif
c end of the loop on the blocks
enddo
c close the density and integral files
if(core.ne.'corr '.and.ncore.ne.0) close(scrfile1)
close(densfile)
c diagonal part of the G matrices
c RHF case
if(scftype.eq.'rhf ') then
do a=1,nvirtal
do i=1,nal+ncore
do e=1,nvirtal
do m=1,nal+ncore
if(a.eq.e.and.m.eq.i)
$gmtx(a,i,e,m)=gmtx(a,i,e,m)+fockmoa(a+nal+ncore,a+nal+ncore)-
$fockmoa(i,i)
enddo
enddo
enddo
enddo
endif
c UHF case
if(scftype.eq.'uhf ') then
c aaaa-block
q=1
do m=1,nal+ncore
do e=1,nvirtal
p=1
do i=1,nal+ncore
do a=1,nvirtal
if(a.eq.e.and.m.eq.i) gmtxuhf(p,q)=gmtxuhf(p,q)+
$fockmoa(a+nal+ncore,a+nal+ncore)-fockmoa(i,i)
p=p+1
enddo
enddo
q=q+1
enddo
enddo
c bbbb-block
q=1+nvirtal*(nal+ncore)
do m=1,nbe+ncore
do e=1,nvirtbe
p=1+nvirtal*(nal+ncore)
do i=1,nbe+ncore
do a=1,nvirtbe
if(a.eq.e.and.m.eq.i) gmtxuhf(p,q)=gmtxuhf(p,q)+
$fockmob(a+nbe+ncore,a+nbe+ncore)-fockmob(i,i)
p=p+1
enddo
enddo
q=q+1
enddo
enddo
endif
c ROHF case
if(scftype.eq.'rohf ') then
do a=1,nvirt+nos
do i=1,nocc+nos
do e=1,nvirt+nos
do m=1,nocc+nos
c os-os
if((e.le.nos.and.m.gt.nocc).or.(a.le.nos.and.i
$.gt.nocc)) then
else
if(a.eq.e.and.m.eq.i)
$gmtxrohf(a,i,e,m)=gmtxrohf(a,i,e,m)+fockmoa(a+nocc,a+nocc)-
$fockmoa(i,i)
endif
enddo
enddo
enddo
enddo
endif
C **********************************************************************
C Build the "X" intermediate
C **********************************************************************
C X-aa
if(scftype.ne.'rohf ') then
call dfillzero(xinta,nvirtal*(nal+ncore))
do a=1,nvirtal
do i=1,nal+ncore
xinta(a,i)=iprimea(a+nal+ncore,i)-iprimea(i,a+nal+ncore)
enddo
enddo
C X-bb (UHF case)
if(scftype.eq.'uhf ') then
call dfillzero(xintb,nvirtbe*(nbe+ncore))
do a=1,nvirtbe
do i=1,nbe+ncore
xintb(a,i)=iprimeb(a+nbe+ncore,i)-iprimeb(i,a+nbe+ncore)
enddo
enddo
endif
C ROHF case
else
call dfillzero(xintrohf,(nvirt+nos)*(nocc+nos))
do a=1,nvirt+nos
do i=1,nocc+nos
if ((a.le.nos).and.(i.gt.nocc)) then
else
xintrohf(a,i)=iprimea(a+nocc,i)-iprimea(i,a+nocc)
endif
enddo
enddo
endif
C FROZEN CORE
c f-i blocks
if(core.ne.'corr '.and.ncore.ne.0) then
call dfillzero(xintfa,nal*ncore)
do i=1,nal
do f=1,ncore
xintfa(i,f)=iprimea(i+ncore,f)-iprimea(f,i+ncore)
enddo
enddo
c UHF case
if(scftype.eq.'uhf ') then
call dfillzero(xintfb,nbe*ncore)
do i=1,nbe
do f=1,ncore
xintfb(i,f)=iprimeb(i+ncore,f)-iprimeb(f,i+ncore)
enddo
enddo
endif
c ROHF case
if(scftype.eq.'rohf ') then
call dfillzero(xintfrohf,(nocc-ncore)*ncore)
do i=1,nocc-ncore
do f=1,ncore
xintfrohf(i,f)=iprimea(i+ncore,f)-iprimea(f,i+ncore)
enddo
enddo
endif
c extra terms in the gradient expression
c core Hamiltonian part - added to the one-electron density
if(scftype.ne.'rohf ') then
do i=1,nal
do f=1,ncore
da(i+ncore,f)=da(i+ncore,f)-
$1.d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)
da(f,i+ncore)=da(f,i+ncore)-
$1.d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)
enddo
enddo
c UHF case
if(scftype.eq.'uhf ') then
do i=1,nbe
do f=1,ncore
db(i+ncore,f)=db(i+ncore,f)-
$1.d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)
db(f,i+ncore)=db(f,i+ncore)-
$1.d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)
enddo
enddo
endif
c ROHF case
else
do i=1,nocc-ncore
do f=1,ncore
da(i+ncore,f)=da(i+ncore,f)-
$1.d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)
da(f,i+ncore)=da(f,i+ncore)-
$1.d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)
enddo
enddo
endif
C **********************************************************************
C X_al -> Y_al
C **********************************************************************
open(unit=scrfile1,file='PLUSCORE.55')
nvirtshift=0
nosshift=0
c loop on the blocks
do ii=1,max(nal+ncore,nbe+ncore)/maxw+1
c block sizes
c alpha
if(ii*maxw.le.nal+ncore) then
nblockal=maxw
else
if((ii-1)*maxw.gt.nal+ncore) then
nblockal=0
else
nblockal=mod(nal+ncore,maxw)
endif
endif
c beta
if(ii*maxw.le.nbe+ncore) then
nblockbe=maxw
else
if((ii-1)*maxw.gt.nbe+ncore) then
nblockbe=0
else
nblockbe=mod(nbe+ncore,maxw)
endif
endif
c core
if(ii*maxw.le.ncore) then
nblockcore=maxw
else
if((ii-1)*maxw.gt.ncore) then
nblockcore=0
else
nblockcore=mod(ncore,maxw)
endif
endif
c nblock
if(ii*maxw.le.nbe+ncore) then
nblock=maxw
else
nblock=mod(max(nal+ncore,nbe+ncore),maxw)
endif
c ***read the integrals needed***
call dfillzero(vint2aa,maxw*(nbasis+ncore)**3)
C frozen core
rewind(scrfile1)
C RHF and ROHF cases
if(scftype.ne.'uhf ') then
1046 read(scrfile1,*) sum2,i,j,k,l
if(l.eq.0) goto 1046
if((i.gt.(ii-1)*maxw).and.(i.le.(ii-1)*maxw+nblockal)) then
vint2aa(l,j,k,i-(ii-1)*maxw)=sum2
vint2aa(k,j,l,i-(ii-1)*maxw)=sum2
endif
if((j.gt.(ii-1)*maxw).and.(j.le.(ii-1)*maxw+nblockal)) then
vint2aa(k,i,l,j-(ii-1)*maxw)=sum2
vint2aa(l,i,k,j-(ii-1)*maxw)=sum2
endif
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblockal)) then
vint2aa(j,l,i,k-(ii-1)*maxw)=sum2
vint2aa(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblockal)) then
vint2aa(i,k,j,l-(ii-1)*maxw)=sum2
vint2aa(j,k,i,l-(ii-1)*maxw)=sum2
endif
1047 read(scrfile1,*) sum2,i,j,k,l
if((i.gt.(ii-1)*maxw).and.(i.le.(ii-1)*maxw+nblockal)) then
vint2aa(l,j,k,i-(ii-1)*maxw)=sum2
vint2aa(k,j,l,i-(ii-1)*maxw)=sum2
endif
if((j.gt.(ii-1)*maxw).and.(j.le.(ii-1)*maxw+nblockal)) then
vint2aa(k,i,l,j-(ii-1)*maxw)=sum2
vint2aa(l,i,k,j-(ii-1)*maxw)=sum2
endif
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblockal)) then
vint2aa(j,l,i,k-(ii-1)*maxw)=sum2
vint2aa(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblockal)) then
vint2aa(i,k,j,l-(ii-1)*maxw)=sum2
vint2aa(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1047
else
C UHF case
C aaaa list
1048 read(scrfile1,*) sum2,i,j,k,l
if(i.ne.0) goto 1048
1049 read(scrfile1,*) sum2,i,j,k,l
if((i.gt.(ii-1)*maxw).and.(i.le.(ii-1)*maxw+nblockal)) then
vint2aa(l,j,k,i-(ii-1)*maxw)=sum2
vint2aa(k,j,l,i-(ii-1)*maxw)=sum2
endif
if((j.gt.(ii-1)*maxw).and.(j.le.(ii-1)*maxw+nblockal)) then
vint2aa(k,i,l,j-(ii-1)*maxw)=sum2
vint2aa(l,i,k,j-(ii-1)*maxw)=sum2
endif
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblockal)) then
vint2aa(j,l,i,k-(ii-1)*maxw)=sum2
vint2aa(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblockal)) then
vint2aa(i,k,j,l-(ii-1)*maxw)=sum2
vint2aa(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1049
call dfillzero(vint2bb,maxw*(nbasis+ncore)**3)
call dfillzero(vint2ab,maxw*(nbasis+ncore)**3)
C bbbb list
1050 read(scrfile1,*) sum2,i,j,k,l
if(i.ne.0) goto 1050
1051 read(scrfile1,*) sum2,i,j,k,l
if((i.gt.(ii-1)*maxw).and.(i.le.(ii-1)*maxw+nblockbe)) then
vint2bb(l,j,k,i-(ii-1)*maxw)=sum2
vint2bb(k,j,l,i-(ii-1)*maxw)=sum2
endif
if((j.gt.(ii-1)*maxw).and.(j.le.(ii-1)*maxw+nblockbe)) then
vint2bb(k,i,l,j-(ii-1)*maxw)=sum2
vint2bb(l,i,k,j-(ii-1)*maxw)=sum2
endif
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblockbe)) then
vint2bb(j,l,i,k-(ii-1)*maxw)=sum2
vint2bb(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblockbe)) then
vint2bb(i,k,j,l-(ii-1)*maxw)=sum2
vint2bb(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1051
c abab list
1052 read(scrfile1,*) sum2,k,l,i,j
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblockbe)) then
vint2ab(j,l,i,k-(ii-1)*maxw)=sum2
vint2ab(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblockbe)) then
vint2ab(i,k,j,l-(ii-1)*maxw)=sum2
vint2ab(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1052
endif
nuc=sum2
c BUILD THE MATRIX
if(scftype.ne.'rohf ') then
do a=1,nvirtal
do l=1,nal+ncore
sum=0.d0
do i=1,nal
do f=1,nblockcore
sum=sum+1.d0/(fockmoa(f+(ii-1)*maxw,f+(ii-1)*maxw)
$-fockmoa(i+ncore,i+ncore))*
$xintfa(i,f+(ii-1)*maxw)*(4.d0*fd*vint2aa(a+nal+ncore,i+ncore,l,f)-
$vint2aa(i+ncore,a+nal+ncore,l,f)-vint2aa(i+ncore,l,a+nal+ncore,f))
enddo
enddo
C UHF case
if(scftype.eq.'uhf ') then
do i=1,nbe
do f=1,nblockcore
sum=sum+2.d0/(fockmob(f+(ii-1)*maxw,f+(ii-1)*maxw)
$-fockmob(i+ncore,i+ncore))*
$xintfb(i,f+(ii-1)*maxw)*vint2ab(a+nal+ncore,i+ncore,l,f)
enddo
enddo
endif
xinta(a,l)=xinta(a,l)+sum
enddo
enddo
if(scftype.eq.'uhf ') then
do a=1,nvirtbe
do l=1,nblockbe
sum=0.d0
do i=1,nbe
do f=1,ncore
sum=sum+1.d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*
$xintfb(i,f)*(4.d0*fd*vint2bb(i+ncore,a+nbe+ncore,f,l)-
$vint2bb(a+nbe+ncore,i+ncore,f,l)-vint2bb(i+ncore,f,a+nbe+ncore,l))
enddo
enddo
do i=1,nal
do f=1,ncore
sum=sum+2.d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*
$xintfa(i,f)*vint2ab(i+ncore,a+nbe+ncore,f,l)
enddo
enddo
xintb(a,l+(ii-1)*maxw)=xintb(a,l+(ii-1)*maxw)+sum
enddo
enddo
endif
C ROHF case
else
c the last index is half occupied********
info=0
c the last index is either doubly or half occupied (or virtual)
if(((ii-1)*maxw.lt.nocc).and.((ii-1)*maxw+nblockal.gt.
$nocc)) then
nosblock=min(nos,(ii-1)*maxw+nblockal-nocc)
nosshift=nocc-(ii-1)*maxw
info=1
else
c the last index is half occupied
if(((ii-1)*maxw.ge.nocc).and.((ii-1)*maxw+nblockal.le.
$nocc+nos)) then
nosblock=nblockal
nosshift=0
info=1
else
c the last index is either half occupied or virtual
if(((ii-1)*maxw.lt.nocc+nos).and.((ii-1)*maxw+nblockal.gt.
$nocc+nos)) then
nosblock=nocc+nos-(ii-1)*maxw
nosshift=0
info=1
endif
endif
endif
if(nosblock.eq.0) info=0
c half occupied - doubly occupied
if(info.eq.1) then
do k=nosshift+1,nosshift+nosblock
do l=1,nocc
sum=0.d0
do i=1,nocc-ncore
do f=1,ncore
sum=sum+xintfrohf(i,f)/
$(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*
$(
$2.d0*vint2aa(i+ncore,l,f,k)-1.5d0*(vint2aa(f,i+ncore,l,k)+
$vint2aa(i+ncore,f,l,k))
$)
if(i+ncore.eq.l) then
do m=nocc+1,nocc+nos
sum=sum+xintfrohf(i,f)/
$(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*0.5d0*vint2aa(f,m,m,k)
enddo
endif
if(f.eq.l) then
do m=nocc+1,nocc+nos
sum=sum+xintfrohf(i,f)/
$(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*
$0.5d0*vint2aa(i+ncore,m,m,k)
enddo
endif
enddo
enddo
xintrohf(k-nocc+(ii-1)*maxw,l)
$=xintrohf(k-nocc+(ii-1)*maxw,l)+sum
enddo
enddo
endif
c virtual - doubly occupied
do k=nocc+nos+1,nocc+nos+nvirtal
do l=1,nocc
sum=0.d0
do i=1,nocc-ncore
do f=1,nblockcore
sum=sum+xintfrohf(i,f+(ii-1)*maxw)/
$(fockmoa(f+(ii-1)*maxw,f+(ii-1)*maxw)-fockmoa(i+ncore,i+ncore))*
$(
$4.d0*vint2aa(l,i+ncore,k,f)-(vint2aa(i+ncore,l,k,f)+
$vint2aa(i+ncore,k,l,f))
$)
if(i+ncore.eq.l) then
do m=nocc+1,nocc+nos
sum=sum+xintfrohf(i,f+(ii-1)*maxw)/
$(fockmoa(f+(ii-1)*maxw,f+(ii-1)*maxw)
$-fockmoa(i+ncore,i+ncore))*vint2aa(m,m,k,f)
enddo
endif
enddo
enddo
xintrohf(k-nocc,l)=xintrohf(k-nocc,l)+sum
sum=0.d0
do i=1,nocc-ncore
do f=1,ncore
if(f.eq.l.and.info.eq.1) then
do m=nosshift+1,nosshift+nosblock
sum=sum+xintfrohf(i,f)/
$(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))
$*vint2aa(i+ncore,k,m+(ii-1)*maxw,m)
enddo
endif
enddo
enddo
xintrohf(k-nocc,l)=xintrohf(k-nocc,l)+sum
enddo
enddo
c virtual - half occupied
do k=nocc+nos+1,nocc+nos+nvirtal
do l=nocc+1,nocc+nos
sum=0.d0
do i=1,nocc-ncore
do f=1,nblockcore
sum=sum+xintfrohf(i,f+(ii-1)*maxw)/
$(fockmoa(f+(ii-1)*maxw,f+(ii-1)*maxw)-fockmoa(i+ncore,i+ncore))*
$(
$2.d0*vint2aa(l,i+ncore,k,f)+0.5d0*(vint2aa(i+ncore,l,k,f)+
$vint2aa(i+ncore,k,l,f))
$)
if(i+ncore.eq.l) then
do m=nocc+1,nocc+nos
sum=sum+xintfrohf(i,f+(ii-1)*maxw)/
$(fockmoa(f+(ii-1)*maxw,f+(ii-1)*maxw)-fockmoa(i+ncore,i+ncore))
$*0.5d0*vint2aa(m,m,k,f)
enddo
endif
enddo
enddo
xintrohf(k-nocc,l)=xintrohf(k-nocc,l)+sum
sum=0.d0
do i=1,nocc-ncore
do f=1,ncore
if(f.eq.l.and.info.eq.1) then
do m=nosshift+1,nosshift+nosblock
sum=sum+xintfrohf(i,f)/
$(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*
$0.5d0*vint2aa(i+ncore,k,m+(ii-1)*maxw,m)
enddo
endif
enddo
enddo
xintrohf(k-nocc,l)=xintrohf(k-nocc,l)+sum
enddo
enddo
endif
c end of loop on the blocks
enddo
close(scrfile1)
endif
C Eliminate os-os block
if(scftype.eq.'rohf ') then
call dcopy((nvirt+nos)*(nocc+nos),xintrohf,1,xintrohf2,1)
r=0
do q=1,(nvirt+nos)*(nocc+nos)
r=r+1
e=mod(r-1,nvirt+nos)+1
m=(r-1)/(nos+nvirt)+1
if(m.gt.nocc.and.e.le.nos) then
do i=q,(nvirt+nos)*(nocc+nos)
xintrohf2_2(i)=xintrohf2_2(i+nos)
enddo
r=r+nos
endif
enddo
endif
C **********************************************************************
C Build the "G" hypermatrix of the Z-vector method
C **********************************************************************
C ROHF case
if(scftype.eq.'rohf ') then
C Eliminate os-os block
C Rows
r=0
do q=1,(nvirt+nos)*(nocc+nos)
r=r+1
e=mod(r-1,nvirt+nos)+1
m=(r-1)/(nos+nvirt)+1
if(m.gt.nocc.and.e.le.nos) then
do p=1,(nvirt+nos)*(nocc+nos)
do i=q,(nvirt+nos)*(nocc+nos)
gmtxrohf_2(i,p)=gmtxrohf_2(i+nos,p)
enddo
enddo
r=r+nos
endif
enddo
C Columns
r=0
do q=1,(nvirt+nos)*(nocc+nos)
r=r+1
e=mod(r-1,nvirt+nos)+1
m=(r-1)/(nos+nvirt)+1
if(m.gt.nocc.and.e.le.nos) then
do p=1,(nvirt+nos)*(nocc+nos)
do i=q,(nvirt+nos)*(nocc+nos)
gmtxrohf_2(p,i)=gmtxrohf_2(p,i+nos)
enddo
enddo
r=r+nos
endif
enddo
endif
C *********************************************************************
C Solve the Z-vector equation
C *********************************************************************
info=0
C RHF case
if(scftype.eq.'rhf ') call dgesv(nvirtal*(nal+ncore),1,gmtx,
$nvirtal*(nal+ncore),ipiv,xinta,nvirtal*(nal+ncore),info)
C ROHF case
if(scftype.eq.'rohf ') then
call dgesv(
$(nvirt+nos)*(nocc+nos)-nos*nos,
$1,gmtxrohf,(nvirt+nos)*(nocc+nos),ipivrohf2,xintrohf2,
$(nvirt+nos)*(nocc+nos),info)
C Fill the original X-intermediate with the values of the solution of
C the Z-vector eq.
call dfillzero(xintrohf,(nvirt+nos)*(nocc+nos))
r=0
do q=1,(nvirt+nos)*(nocc+nos)
e=mod(q-1,nvirt+nos)+1
m=(q-1)/(nos+nvirt)+1
if(m.gt.nocc.and.e.le.nos) then
r=r+1
xintrohf_2(q)=0.d0
else
xintrohf_2(q)=xintrohf2_2(q-r)
endif
enddo
endif
C UHF case
if(scftype.eq.'uhf ') call dgesv(nvirtal*(nal+ncore)+nvirtbe*
$(nbe+ncore),1,
$gmtxuhf,nvirtal*(nal+ncore)+nvirtbe*(nbe+ncore),ipivuhf,xinta,
$nvirtal*(nal+ncore)+nvirtbe*(nbe+ncore),info)
C *********************************************************************
C Build I intermediate
C *********************************************************************
c SCF grad issues...
call dfillzero(sumscf,(nbasis+ncore)**2)
c virtual-occupied flip
C RHF case
if(scftype.ne.'rohf ') then
do a=1,nvirtal
do i=1,nal+ncore
iprimea(a+nal+ncore,i)=iprimea(i,a+nal+ncore)-fockmoa(i,i)*
$xinta(a,i)
iprimea(i,a+nal+ncore)=iprimea(a+nal+ncore,i)
enddo
enddo
c UHF case
if(scftype.eq.'uhf ') then
do a=1,nvirtbe
do i=1,nbe+ncore
iprimeb(a+nbe+ncore,i)=iprimeb(i,a+nbe+ncore)-fockmob(i,i)
$*xintb(a,i)
iprimeb(i,a+nbe+ncore)=iprimeb(a+nbe+ncore,i)
enddo
enddo
endif
else
c ROHF case
nosshift=0
c active -- frozen flip
do i=ncore+1,nocc
do f=1,ncore
iprimea(i,f)=iprimea(f,i)
enddo
enddo
c half occ. -- doubly occ. flip
do h=nocc+1,nocc+nos
do l=1,nocc
iprimea(h,l)=iprimea(l,h)
enddo
enddo
c other flips
do p=1,nbasis+ncore
do q=1,nbasis+ncore
C virt-os block
if(p.gt.nocc+nos.and.q.gt.nocc.and.q.le.nocc+nos)
$iprimea(p,q)=iprimea(q,p)
c virt-occ block
if(p.gt.nocc+nos.and.q.le.nocc) iprimea(p,q)=iprimea(q,p)
enddo
enddo
endif
if(core.ne.'corr '.and.ncore.ne.0)
$ open(unit=scrfile1,file='PLUSCORE.55')
c loop on the blocks
do ii=1,max(nal+ncore,nbe+ncore)/maxw+1
c block sizes
c alpha
if(ii*maxw.le.nal+ncore) then
nblockal=maxw
else
if((ii-1)*maxw.gt.nal+ncore) then
nblockal=0
else
nblockal=mod(nal+ncore,maxw)
endif
endif
c beta
if(ii*maxw.le.nbe+ncore) then
nblockbe=maxw
else
if((ii-1)*maxw.gt.nbe+ncore) then
nblockbe=0
else
nblockbe=mod(nbe+ncore,maxw)
endif
endif
c core
if(ii*maxw.le.ncore) then
nblockcore=maxw
else
if((ii-1)*maxw.gt.ncore) then
nblockcore=0
else
nblockcore=mod(ncore,maxw)
endif
endif
c nblock
if(ii*maxw.le.nbe+ncore) then
nblock=maxw
else
nblock=mod(max(nal+ncore,nbe+ncore),maxw)
endif
c ***read the integrals needed***
call dfillzero(vint2aa,maxw*(nbasis+ncore)**3)
C frozen core
if(core.ne.'corr '.and.ncore.ne.0) then
rewind(scrfile1)
C RHF and ROHF cases
if(scftype.ne.'uhf ') then
1039 read(scrfile1,*) sum2,i,j,k,l
if(l.eq.0) goto 1039
if((i.gt.(ii-1)*maxw).and.(i.le.(ii-1)*maxw+nblockal)) then
vint2aa(l,j,k,i-(ii-1)*maxw)=sum2
vint2aa(k,j,l,i-(ii-1)*maxw)=sum2
endif
if((j.gt.(ii-1)*maxw).and.(j.le.(ii-1)*maxw+nblockal)) then
vint2aa(k,i,l,j-(ii-1)*maxw)=sum2
vint2aa(l,i,k,j-(ii-1)*maxw)=sum2
endif
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblockal)) then
vint2aa(j,l,i,k-(ii-1)*maxw)=sum2
vint2aa(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblockal)) then
vint2aa(i,k,j,l-(ii-1)*maxw)=sum2
vint2aa(j,k,i,l-(ii-1)*maxw)=sum2
endif
1040 read(scrfile1,*) sum2,i,j,k,l
if((i.gt.(ii-1)*maxw).and.(i.le.(ii-1)*maxw+nblockal)) then
vint2aa(l,j,k,i-(ii-1)*maxw)=sum2
vint2aa(k,j,l,i-(ii-1)*maxw)=sum2
endif
if((j.gt.(ii-1)*maxw).and.(j.le.(ii-1)*maxw+nblockal)) then
vint2aa(k,i,l,j-(ii-1)*maxw)=sum2
vint2aa(l,i,k,j-(ii-1)*maxw)=sum2
endif
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblockal)) then
vint2aa(j,l,i,k-(ii-1)*maxw)=sum2
vint2aa(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblockal)) then
vint2aa(i,k,j,l-(ii-1)*maxw)=sum2
vint2aa(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1040
else
C UHF case
C aaaa list
1041 read(scrfile1,*) sum2,i,j,k,l
if(i.ne.0) goto 1041
1042 read(scrfile1,*) sum2,i,j,k,l
if((i.gt.(ii-1)*maxw).and.(i.le.(ii-1)*maxw+nblockal)) then
vint2aa(l,j,k,i-(ii-1)*maxw)=sum2
vint2aa(k,j,l,i-(ii-1)*maxw)=sum2
endif
if((j.gt.(ii-1)*maxw).and.(j.le.(ii-1)*maxw+nblockal)) then
vint2aa(k,i,l,j-(ii-1)*maxw)=sum2
vint2aa(l,i,k,j-(ii-1)*maxw)=sum2
endif
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblockal)) then
vint2aa(j,l,i,k-(ii-1)*maxw)=sum2
vint2aa(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblockal)) then
vint2aa(i,k,j,l-(ii-1)*maxw)=sum2
vint2aa(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1042
call dfillzero(vint2bb,maxw*(nbasis+ncore)**3)
call dfillzero(vint2ab,maxw*(nbasis+ncore)**3)
C bbbb list
1043 read(scrfile1,*) sum2,i,j,k,l
if(i.ne.0) goto 1043
1044 read(scrfile1,*) sum2,i,j,k,l
if((i.gt.(ii-1)*maxw).and.(i.le.(ii-1)*maxw+nblockbe)) then
vint2bb(l,j,k,i-(ii-1)*maxw)=sum2
vint2bb(k,j,l,i-(ii-1)*maxw)=sum2
endif
if((j.gt.(ii-1)*maxw).and.(j.le.(ii-1)*maxw+nblockbe)) then
vint2bb(k,i,l,j-(ii-1)*maxw)=sum2
vint2bb(l,i,k,j-(ii-1)*maxw)=sum2
endif
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblockbe)) then
vint2bb(j,l,i,k-(ii-1)*maxw)=sum2
vint2bb(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblockbe)) then
vint2bb(i,k,j,l-(ii-1)*maxw)=sum2
vint2bb(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1044
c abab list
1045 read(scrfile1,*) sum2,k,l,i,j
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblockbe)) then
vint2ab(j,l,i,k-(ii-1)*maxw)=sum2
vint2ab(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblockbe)) then
vint2ab(i,k,j,l-(ii-1)*maxw)=sum2
vint2ab(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1045
endif
nuc=sum2
else
C correlated core
rewind(inp)
read(inp,*) i,j
read(inp,*) (i,j=1,nbasis+ncore)
read(inp,*) i
C aaaa list
1029 read(inp,*) sum2,i,j,k,l
if((i.gt.(ii-1)*maxw).and.(i.le.(ii-1)*maxw+nblockal)) then
vint2aa(l,j,k,i-(ii-1)*maxw)=sum2
vint2aa(k,j,l,i-(ii-1)*maxw)=sum2
endif
if((j.gt.(ii-1)*maxw).and.(j.le.(ii-1)*maxw+nblockal)) then
vint2aa(k,i,l,j-(ii-1)*maxw)=sum2
vint2aa(l,i,k,j-(ii-1)*maxw)=sum2
endif
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblockal)) then
vint2aa(j,l,i,k-(ii-1)*maxw)=sum2
vint2aa(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblockal)) then
vint2aa(i,k,j,l-(ii-1)*maxw)=sum2
vint2aa(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1029
c UHF case
if(scftype.eq.'uhf ') then
call dfillzero(vint2bb,maxw*(nbasis+ncore)**3)
call dfillzero(vint2ab,maxw*(nbasis+ncore)**3)
C bbbb list
1030 read(inp,*) sum2,i,j,k,l
if((i.gt.(ii-1)*maxw).and.(i.le.(ii-1)*maxw+nblockbe)) then
vint2bb(l,j,k,i-(ii-1)*maxw)=sum2
vint2bb(k,j,l,i-(ii-1)*maxw)=sum2
endif
if((j.gt.(ii-1)*maxw).and.(j.le.(ii-1)*maxw+nblockbe)) then
vint2bb(k,i,l,j-(ii-1)*maxw)=sum2
vint2bb(l,i,k,j-(ii-1)*maxw)=sum2
endif
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblockbe)) then
vint2bb(j,l,i,k-(ii-1)*maxw)=sum2
vint2bb(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblockbe)) then
vint2bb(i,k,j,l-(ii-1)*maxw)=sum2
vint2bb(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1030
c abab list
1031 read(inp,*) sum2,i,j,k,l
if((k.gt.(ii-1)*maxw).and.(k.le.(ii-1)*maxw+nblockbe)) then
vint2ab(j,l,i,k-(ii-1)*maxw)=sum2
vint2ab(i,l,j,k-(ii-1)*maxw)=sum2
endif
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblockbe)) then
vint2ab(i,k,j,l-(ii-1)*maxw)=sum2
vint2ab(j,k,i,l-(ii-1)*maxw)=sum2
endif
if(l.ne.0) goto 1031
endif
endif
c BUILD THE MATRIX
if(scftype.ne.'rohf ') then
do l=1,nal+ncore
do m=1,nal+ncore
sum=0.d0
do e=1,nvirtal
do n=1,nblockal
sum=sum+xinta(e,n+(ii-1)*maxw)*
$(fd*4.d0*vint2aa(l,e+nal+ncore,m,n)-
$vint2aa(e+nal+ncore,l,m,n)-vint2aa(l,m,e+nal+ncore,n))
enddo
enddo
C frozen core - J(m,l) term
if(core.ne.'corr '.and.ncore.ne.0) then
sum2=0.d0
do i=1,nal
do f=1,nblockcore
sum2=sum2+1.d0/(fockmoa(f+(ii-1)*maxw,f+(ii-1)*maxw)
$-fockmoa(i+ncore,i+ncore))*
$xintfa(i,f+(ii-1)*maxw)*(4.d0*fd*vint2aa(m,i+ncore,l,f)-
$vint2aa(i+ncore,m,l,f)-vint2aa(i+ncore,l,m,f))
enddo
enddo
sum=sum-sum2
endif
iprimea(l,m)=iprimea(l,m)-sum
enddo
enddo
C UHF case
if(scftype.eq.'uhf ') then
do l=1,nal+ncore
do m=1,nal+ncore
sum=0.d0
do e=1,nvirtbe
do n=1,nblockbe
sum=sum+xintb(e,n+(ii-1)*maxw)
$*2.d0*vint2ab(l,e+nbe+ncore,m,n)
enddo
enddo
C frozen core - J(m,l) term
if(core.ne.'corr '.and.ncore.ne.0) then
sum2=0.d0
do i=1,nbe
do f=1,nblockcore
sum2=sum2+1.d0/(fockmob(f+(ii-1)*maxw,f+(ii-1)*maxw)
$-fockmob(i+ncore,i+ncore))*
$xintfb(i,f+(ii-1)*maxw)*2.d0*vint2ab(l,i+ncore,m,f)
enddo
enddo
sum=sum-sum2
endif
iprimea(l,m)=iprimea(l,m)-sum
enddo
enddo
C Beta block
do l=1,nbe+ncore
do m=1,nblockbe
sum=0.d0
do e=1,nvirtbe
do n=1,nbe+ncore
sum=sum+xintb(e,n)
$*(2.d0*vint2bb(e+nbe+ncore,l,n,m)-
$vint2bb(l,e+nbe+ncore,n,m)-vint2bb(l,n,e+nbe+ncore,m))
enddo
enddo
do e=1,nvirtal
do n=1,nal+ncore
sum=sum+xinta(e,n)*2.d0*vint2ab(n,l,e+nal+ncore,m)
enddo
enddo
C frozen core - J(m,l) term
if(core.ne.'corr '.and.ncore.ne.0) then
sum2=0.d0
do i=1,nbe
do f=1,ncore
sum2=sum2+1.d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*
$xintfb(i,f)*(4.d0*fd*vint2bb(f,l,i+ncore,m)-
$vint2bb(l,f,i+ncore,m)-vint2bb(l,i+ncore,f,m))
enddo
enddo
do i=1,nal
do f=1,ncore
sum2=sum2+1.d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*
$xintfa(i,f)*2.d0*vint2ab(i+ncore,l,f,m)
enddo
enddo
sum=sum-sum2
endif
iprimeb(l,m+(ii-1)*maxw)=iprimeb(l,m+(ii-1)*maxw)-sum
enddo
enddo
endif
else
c ROHF case
c the last index is 'i'********
info2=0
c the last index is not virtual
if((ii-1)*maxw+nblock.le.nocc+nos) then
noccblock=maxw
info2=2
else
c the last index is either virtual or not
if(((ii-1)*maxw.lt.nocc+nos).and.((ii-1)*maxw+nblock.gt.
$nocc+nos)) then
noccblock=nocc+nos-(ii-1)*maxw
info2=2
endif
endif
if(noccblock.eq.0) info2=0
c the last index is half occupied********
info=0
c the last index is either doubly or half occupied (or virtual)
if(((ii-1)*maxw.lt.nocc).and.((ii-1)*maxw+nblockal.gt.
$nocc)) then
nosblock=min(nos,(ii-1)*maxw+nblockal-nocc)
nosshift=nocc-(ii-1)*maxw
info=1
else
c the last index is half occupied
if(((ii-1)*maxw.ge.nocc).and.((ii-1)*maxw+nblockal.le.
$nocc+nos)) then
nosblock=nblockal
nosshift=0
info=1
else
c the last index is either half occupied or virtual
if(((ii-1)*maxw.lt.nocc+nos).and.((ii-1)*maxw+nblockal.gt.
$nocc+nos)) then
nosblock=nocc+nos-(ii-1)*maxw
nosshift=0
info=1
endif
endif
endif
if(nosblock.eq.0) info=0
c the last index is doubly occupied********
info3=0
c the last index is doubly occuppied
if((ii-1)*maxw+nblockal.le.nocc) then
ndoublyoccblock=maxw
info3=3
else
c the last index is either doubly occuppied or not
if(((ii-1)*maxw.lt.nocc).and.((ii-1)*maxw+nblockal.gt.
$nocc)) then
ndoublyoccblock=nocc-(ii-1)*maxw
info3=3
endif
endif
do p=1,nbasis+ncore
do q=1,nbasis+ncore
C virt-os block
if(p.gt.nocc+nos.and.q.gt.nocc.and.q.le.nocc+nos
$.and.info.eq.1) then
sum=0.d0
do i=1,nocc+nos
do n=nosshift+1,nosshift+nosblock
sum=sum+xintrohf(p-nocc,i)
$*vint2aa(i,q,n+(ii-1)*maxw,n)
enddo
enddo
iprimea(p,q)=iprimea(p,q)+sum
endif
c virt-occ block
if(p.gt.nocc+nos.and.q.le.nocc.and.info.eq.1) then
sum=0.d0
do i=1,nocc+nos
do n=nosshift+1,nosshift+nosblock
sum=sum+xintrohf(p-nocc,i)*vint2aa(i,q,n+(ii-1)*maxw,n)
enddo
enddo
sum=sum*2.d0
iprimea(p,q)=iprimea(p,q)+sum
endif
c occ-occ block
c the last index is 'i'
if(p.le.nocc.and.q.le.nocc.and.info2.eq.2) then
sum=0.d0
do a=1,nos+nvirt
do i=1,noccblock
sum=sum-xintrohf(a,i+(ii-1)*maxw)
$*(2.d0*vint2aa(p,a+nocc,q,i)-vint2aa(p,q,a+nocc,i))
enddo
enddo
sum=sum*2.d0
iprimea(p,q)=iprimea(p,q)+sum
endif
c the last index is half occupied
if(p.le.nocc.and.q.le.nocc.and.info.eq.1) then
sum=0.d0
do a=1,nos+nvirt
do n=nosshift+1,nosshift+nosblock
sum=sum+xintrohf(a,p)*vint2aa(a+nocc,q,n+(ii-1)*maxw,n)
enddo
do i=1,nocc+nos
sum2=0.d0
do n=nosshift+1,nosshift+nosblock
if(i.eq.q) sum2=sum2+vint2aa(a+nocc,p,n+(ii-1)*maxw,n)
enddo
sum=sum-xintrohf(a,i)*sum2
enddo
enddo
sum=sum*2.d0
iprimea(p,q)=iprimea(p,q)+sum
endif
c os-os block
c the last index is half occupied
if(p.gt.nocc.and.p.le.nocc+nos.and.q.gt.nocc
$.and.q.le.nocc+nos.and.info.eq.1) then
sum=0.d0
do a=1,nos+nvirt
do i=1,nocc+nos
sum2=0.d0
do n=nosshift+1,nosshift+nosblock
if(a+nocc.eq.q) sum2=sum2+vint2aa(i,p,n+(ii-1)*maxw,n)
if(i.eq.q) sum2=sum2+vint2aa(a+nocc,p,n+(ii-1)*maxw,n)
enddo
sum2=sum2*0.5d0
sum=sum-xintrohf(a,i)*sum2
enddo
enddo
do i=1,nocc+nos
do n=nosshift+1,nosshift+nosblock
sum=sum+xintrohf(p-nocc,i)*vint2aa(i,q,n+(ii-1)*maxw,n)
$*0.5d0
enddo
enddo
do a=1,nos+nvirt
do n=nosshift+1,nosshift+nosblock
sum=sum+xintrohf(a,q)*vint2aa(a+nocc,p,n+(ii-1)*maxw,n)
$*0.5d0
enddo
enddo
sum=sum*2.d0
iprimea(p,q)=iprimea(p,q)+sum
endif
c the last index is 'i'
if(p.gt.nocc.and.p.le.nocc+nos.and.q.gt.nocc
$.and.q.le.nocc+nos.and.info2.eq.2) then
sum=0.d0
do a=1,nos+nvirt
do i=1,noccblock
c Occupation number
if(a.gt.nos) then
fa=0.d0
else
fa=1.d0
endif
if(i+(ii-1)*maxw.le.nocc) then
fi=2.d0
else
fi=1.d0
endif
c
sum=sum-xintrohf(a,i+(ii-1)*maxw)*(vint2aa(p,a+nocc,q,i)
$-0.5d0*(3.d0-fa-fi)*vint2aa(p,q,a+nocc,i))
enddo
enddo
sum=sum*2.d0
iprimea(p,q)=iprimea(p,q)+sum
endif
c os-occ block
c the last index is half occupied
if(q.le.nocc.and.p.gt.nocc.and.p.le.nocc+nos.and.info.eq.1)
$then
sum=0.d0
do a=1,nos+nvirt
do i=1,nocc+nos
sum2=0.d0
do n=nosshift+1,nosshift+nosblock
if(i.eq.q) sum2=sum2+vint2aa(a+nocc,p,n+(ii-1)*maxw,n)
if(a+nocc.eq.p) sum2=sum2+vint2aa(i,q,n+(ii-1)*maxw,n)
if(i.eq.p) sum2=sum2+vint2aa(a+nocc,q,n+(ii-1)*maxw,n)
enddo
sum2=sum2*0.5d0
sum=sum-xintrohf(a,i)*sum2
enddo
enddo
do i=1,nos+nocc
do n=nosshift+1,nosshift+nosblock
sum=sum+xintrohf(p-nocc,i)*vint2aa(i,q,n+(ii-1)*maxw,n)
enddo
enddo
do a=1,nos+nvirt
do n=nosshift+1,nosshift+nosblock
sum=sum+xintrohf(a,p)*vint2aa(a+nocc,q,n+(ii-1)*maxw,n)
$+xintrohf(a,q)*vint2aa(a+nocc,p,n+(ii-1)*maxw,n)*0.5d0
enddo
enddo
sum=sum*2.d0
iprimea(p,q)=iprimea(p,q)+sum
endif
c the last index is 'i'
if(q.le.nocc.and.p.gt.nocc.and.p.le.nocc+nos.and.info2.eq.2)
$then
sum=0.d0
do a=1,nos+nvirt
do i=1,noccblock
c Occupation number
if(a.gt.nos) then
fa=0.d0
else
fa=1.d0
endif
if(i+(ii-1)*maxw.le.nocc) then
fi=2.d0
else
fi=1.d0
endif
c
sum=sum-xintrohf(a,i+(ii-1)*maxw)
$*(2.d0*vint2aa(q,a+nocc,p,i)-0.5d0
$*(3.d0-fa-fi)*(vint2aa(q,p,a+nocc,i)+vint2aa(p,q,a+nocc,i)))
enddo
enddo
sum=sum*2.d0
iprimea(p,q)=iprimea(p,q)+sum
endif
enddo
enddo
c frozen core
if(core.ne.'corr '.and.ncore.ne.0) then
c doubly occ. -- doubly occ.
c the last index is doubly occ
if(info3.eq.3) then
do l=1,ndoublyoccblock
do m=1,nocc
sum=0.d0
do i=ncore+1,nocc
do f=1,ncore
sum=sum+xintfrohf(i-ncore,f)/(fockmoa(f,f)-fockmoa(i,i))*
$(2.d0*vint2aa(i,m,f,l)-vint2aa(f,i,m,l))
if(i.eq.m) then
do h=nocc+1,nocc+nos
sum=sum+xintfrohf(i-ncore,f)
$/(fockmoa(f,f)-fockmoa(i,i))*vint2aa(f,h,h,l)
enddo
endif
if(f.eq.m) then
do h=nocc+1,nocc+nos
sum=sum+xintfrohf(i-ncore,f)
$/(fockmoa(f,f)-fockmoa(i,i))*vint2aa(i,h,h,l)
enddo
endif
enddo
enddo
iprimea(l+(ii-1)*maxw,m)=iprimea(l+(ii-1)*maxw,m)+2.d0*sum
enddo
enddo
endif
c the last index is half occ
if(info.eq.1) then
c half occ. -- half occ.
do h=nocc+1,nocc+nos
do hh=nosshift+1,nosshift+nosblock
sum=0.d0
do i=ncore+1,nocc
do f=1,ncore
sum=sum+xintfrohf(i-ncore,f)/(fockmoa(f,f)-fockmoa(i,i))*
$(vint2aa(i,h,f,hh)+0.5d0*vint2aa(i,f,h,hh))
enddo
enddo
iprimea(hh+(ii-1)*maxw,h)=iprimea(hh+(ii-1)*maxw,h)+2.d0*sum
enddo
enddo
c doubly occ. -- half occ.
do l=1,nocc
do h=nosshift+1,nosshift+nosblock
sum=0.d0
do i=ncore+1,nocc
do f=1,ncore
sum=sum+xintfrohf(i-ncore,f)/(fockmoa(f,f)-fockmoa(i,i))*
$(2.d0*vint2aa(i,l,f,h)+0.5d0*(vint2aa(i,f,l,h)+vint2aa(f,i,l,h)))
if(i.eq.l) then
do hh=nocc+1,nocc+nos
sum=sum
$+xintfrohf(i-ncore,f)/(fockmoa(f,f)-fockmoa(i,i))*
$0.5d0*vint2aa(f,hh,hh,h)
enddo
endif
if(f.eq.l) then
do hh=nocc+1,nocc+nos
sum=sum
$+xintfrohf(i-ncore,f)/(fockmoa(f,f)-fockmoa(i,i))*
$0.5d0*vint2aa(i,hh,hh,h)
enddo
endif
enddo
enddo
iprimea(h+(ii-1)*maxw,l)=iprimea(h+(ii-1)*maxw,l)+2.d0*sum
enddo
enddo
c active -- half occ.
do i=ncore+1,nocc
do h=nosshift+1,nosshift+nosblock
sum=0.d0
do f=1,ncore
do hh=nocc+1,nocc+nos
sum=sum+xintfrohf(i-ncore,f)/(fockmoa(f,f)-fockmoa(i,i))*
$(-0.5d0*vint2aa(f,hh,hh,h))
enddo
enddo
iprimea(h+(ii-1)*maxw,i)=iprimea(h+(ii-1)*maxw,i)+2.d0*sum
enddo
enddo
c half occ. -- frozen
do h=nosshift+1,nosshift+nosblock
do f=1,ncore
sum=0.d0
do i=ncore+1,nocc
do hh=nocc+1,nocc+nos
sum=sum+xintfrohf(i-ncore,f)/(fockmoa(f,f)-fockmoa(i,i))*
$(-0.5d0*vint2aa(i,hh,hh,h))
enddo
enddo
iprimea(h+(ii-1)*maxw,f)=iprimea(h+(ii-1)*maxw,f)+2.d0*sum
enddo
enddo
c active -- doubly occ.
do i=ncore+1,nocc
do l=1,nocc
sum=0.d0
do f=1,ncore
do h=nosshift+1,nosshift+nosblock
sum=sum+xintfrohf(i-ncore,f)/(fockmoa(f,f)-fockmoa(i,i))*
$(-vint2aa(l,f,h+(ii-1)*maxw,h))
enddo
enddo
iprimea(i,l)=iprimea(i,l)+2.d0*sum
enddo
enddo
c doubly occ. -- frozen
do l=1,nocc
do f=1,ncore
sum=0.d0
do i=ncore+1,nocc
do h=nosshift+1,nosshift+nosblock
sum=sum+xintfrohf(i-ncore,f)/(fockmoa(f,f)-fockmoa(i,i))*
$(-vint2aa(i,l,h+(ii-1)*maxw,h))
enddo
enddo
iprimea(l,f)=iprimea(l,f)+2.d0*sum
enddo
enddo
endif
c end of frozen core
endif
c SCF gradient contrubution
c the last index is half occ
if(info.eq.1) then
do mu=1,nbasis+ncore
do nu=1,nbasis+ncore
sum=0.d0
do i=1,nocc+nos
do j=1,nocc+nos
C Occupation number
if(i.le.nocc) then
fi=2.d0
else
fi=1.d0
endif
if(j.le.nocc) then
fj=2.d0
else
fj=1.d0
endif
do m=nosshift+1,nosshift+nosblock
sum=sum-0.5d0*fi*fj*vint2aa(i,j,m+(ii-1)*maxw,m)
$*ca(mu,i)*ca(nu,j)
enddo
enddo
sumscf(mu,nu)=sumscf(mu,nu)+sum
enddo
enddo
enddo
endif
c end of ROHF
endif
c end of the loop on the blocks
enddo
if(core.ne.'corr '.and.ncore.ne.0) close(scrfile1)
c quantities with less than 4 indices
if(scftype.ne.'rohf ') then
C frozen core - K(i,f) term
if(core.ne.'corr '.and.ncore.ne.0) then
do i=1,nal
do f=1,ncore
iprimea(i+ncore,f)=iprimea(f,i+ncore)
iprimea(i+ncore,f)=iprimea(i+ncore,f)+1.d0*fockmoa(f,f)/
$(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)
iprimea(f,i+ncore)=iprimea(f,i+ncore)+1.d0*fockmoa(f,f)/
$(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)
enddo
enddo
endif
C UHF case
if(scftype.eq.'uhf ') then
C frozen core - K(i,f) term
if(core.ne.'corr '.and.ncore.ne.0) then
do i=1,nbe
do f=1,ncore
iprimeb(i+ncore,f)=iprimeb(f,i+ncore)
iprimeb(i+ncore,f)=iprimeb(i+ncore,f)+1.d0*fockmob(f,f)/
$(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)
iprimeb(f,i+ncore)=iprimeb(f,i+ncore)+1.d0*fockmob(f,f)/
$(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)
enddo
enddo
endif
endif
else
c ROHF case
do p=1,nbasis+ncore
do q=1,nbasis+ncore
C virt-os block
if(p.gt.nocc+nos.and.q.gt.nocc.and.q.le.nocc+nos)
$iprimea(p,q)=iprimea(p,q)-2.d0*fockmoa(q,q)*xintrohf(p-nocc,q)
c virt-occ block
if(p.gt.nocc+nos.and.q.le.nocc) iprimea(p,q)=iprimea(p,q)
$-2.d0*fockmoa(q,q)*xintrohf(p-nocc,q)
c os-os block
if(p.gt.nocc.and.p.le.nocc+nos.and.q.gt.nocc
$.and.q.le.nocc+nos) iprimea(p,q)=iprimea(p,q)
$-2.d0*fockmoa(q,q)*xintrohf(p-nocc,q)
c os-occ block
if(q.le.nocc.and.p.gt.nocc.and.p.le.nocc+nos) iprimea(p,q)
$=iprimea(p,q)-2.d0*fockmoa(q,q)*xintrohf(p-nocc,q)
enddo
enddo
c frozen core: active -- frozen
if(core.ne.'corr '.and.ncore.ne.0) then
do i=ncore+1,nocc
do f=1,ncore
iprimea(i,f)=iprimea(i,f)+2.d0*xintfrohf(i-ncore,f)*
$fockmoa(f,f)/(fockmoa(f,f)-fockmoa(i,i))
enddo
enddo
endif
endif
C End of correlated section
endif
if(trim(grdens).ne.'off'.and.dens.eq.1) goto 1234
C **********************************************************************
C AO transformations
C **********************************************************************
C Transform the two-electron densities into AO-basis
c reset the matrices
call dfillzero(gaoa,maxw*(nbasis+ncore)**3)
c call dfillzero(g2int,maxw*(nbasis+ncore)**3)
c call dfillzero(g2int2,maxw*(nbasis+ncore)**3)
if(scftype.eq.'uhf ') then
call dfillzero(gaob,maxw*(nbasis+ncore)**3)
call dfillzero(gaoab,maxw*(nbasis+ncore)**3)
endif
c open new AO density file
open(scrfile5 ,status='unknown',access='direct',
$recl=(nbasis+ncore)**2*ifltln,file='gaofile')
c open the MO density file
open(unit=densfile,file='CCDENSITIES')
c loop on the MO blocks
do ii=1,(nbasis+ncore)/maxw+1
c block size
if(ii*maxw.le.nbasis+ncore) then
nblock=maxw
else
nblock=mod(nbasis+ncore,maxw)
endif
if(nblock.gt.0) then
c Read MO densities
call dfillzero(g2aa,maxw*(nbasis+ncore)**3)
rewind(densfile)
c aaaa list
1053 read(densfile,*) sum2,i,k,j,l
if((i+ncore.gt.(ii-1)*maxw).and.
$(i+ncore.le.(ii-1)*maxw+nblock)) then
g2aa(l+ncore,j+ncore,k+ncore,i-(ii-1)*maxw+ncore)=sum2
g2aa(k+ncore,j+ncore,l+ncore,i-(ii-1)*maxw+ncore)=sum2
endif
if((j+ncore.gt.(ii-1)*maxw).and.
$(j+ncore.le.(ii-1)*maxw+nblock)) then
g2aa(k+ncore,i+ncore,l+ncore,j-(ii-1)*maxw+ncore)=sum2
g2aa(l+ncore,i+ncore,k+ncore,j-(ii-1)*maxw+ncore)=sum2
endif
if((k+ncore.gt.(ii-1)*maxw).and.
$(k+ncore.le.(ii-1)*maxw+nblock)) then
g2aa(j+ncore,l+ncore,i+ncore,k-(ii-1)*maxw+ncore)=sum2
g2aa(i+ncore,l+ncore,j+ncore,k-(ii-1)*maxw+ncore)=sum2
endif
if((l+ncore.gt.(ii-1)*maxw).and.
$(l+ncore.le.(ii-1)*maxw+nblock)) then
g2aa(i+ncore,k+ncore,j+ncore,l-(ii-1)*maxw+ncore)=sum2
g2aa(j+ncore,k+ncore,i+ncore,l-(ii-1)*maxw+ncore)=sum2
endif
if(l.ne.0) goto 1053
c build the frozen parts
if(core.ne.'corr '.and.ncore.ne.0) then
do f=1,ncore
if((scftype.ne.'uhf ').and.(f.gt.(ii-1)*maxw).and.
$(f.le.(ii-1)*maxw+nblock)) g2aa(f,f,f,f-(ii-1)*maxw)=2.d0
enddo
do f=1,ncore
do g=1,ncore
if(g.ne.f) then
if((f.gt.(ii-1)*maxw).and.(f.le.(ii-1)*maxw+nblock))
$g2aa(g,f,g,f-(ii-1)*maxw)=4.d0*fg
if((g.gt.(ii-1)*maxw).and.(g.le.(ii-1)*maxw+nblock)) then
g2aa(g,f,f,g-(ii-1)*maxw)=-1.d0*fd
g2aa(f,f,g,g-(ii-1)*maxw)=-1.d0*fd
endif
endif
enddo
enddo
do f=1,ncore
do p=ncore+1,ncore+nbasis
do q=ncore+1,ncore+nbasis
if((f.gt.(ii-1)*maxw).and.(f.le.(ii-1)*maxw+nblock)) then
g2aa(p,f,q,f-(ii-1)*maxw)=2.d0*da(p,q)*fd
g2aa(p,q,f,f-(ii-1)*maxw)=-0.5d0*da(p,q)
g2aa(f,q,p,f-(ii-1)*maxw)=-0.5d0*da(p,q)
endif
if((q.gt.(ii-1)*maxw).and.(q.le.(ii-1)*maxw+nblock)) then
g2aa(f,p,f,q-(ii-1)*maxw)=2.d0*da(p,q)*fd
g2aa(p,f,f,q-(ii-1)*maxw)=-0.5d0*da(p,q)
g2aa(f,f,p,q-(ii-1)*maxw)=-0.5d0*da(p,q)
endif
enddo
enddo
enddo
c Extra term in the gradient expression
if(scftype.ne.'rohf ') then
do l=1,nal+ncore
do i=1,nal
do f=1,ncore
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblock)) then
g2aa(i+ncore,l,f,l-(ii-1)*maxw)=
$g2aa(i+ncore,l,f,l-(ii-1)*maxw)-
$0.5d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)*4.d0*fd
g2aa(f,l,i+ncore,l-(ii-1)*maxw)=
$g2aa(f,l,i+ncore,l-(ii-1)*maxw)-
$0.5d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)*4.d0*fd
g2aa(i+ncore,f,l,l-(ii-1)*maxw)=
$g2aa(i+ncore,f,l,l-(ii-1)*maxw)+
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)*2.d0
g2aa(l,f,i+ncore,l-(ii-1)*maxw)=
$g2aa(l,f,i+ncore,l-(ii-1)*maxw)+
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)*2.d0
g2aa(l,i+ncore,f,l-(ii-1)*maxw)=
$g2aa(l,i+ncore,f,l-(ii-1)*maxw)+
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)*2.d0
g2aa(f,i+ncore,l,l-(ii-1)*maxw)=
$g2aa(f,i+ncore,l,l-(ii-1)*maxw)+
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)*2.d0
endif
if((f.gt.(ii-1)*maxw).and.(f.le.(ii-1)*maxw+nblock)) then
g2aa(l,i+ncore,l,f-(ii-1)*maxw)=
$g2aa(l,i+ncore,l,f-(ii-1)*maxw)-
$0.5d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)*4.d0*fd
g2aa(i+ncore,l,l,f-(ii-1)*maxw)=
$g2aa(i+ncore,l,l,f-(ii-1)*maxw)+
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)*2.d0
g2aa(l,l,i+ncore,f-(ii-1)*maxw)=
$g2aa(l,l,i+ncore,f-(ii-1)*maxw)+
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)*2.d0
endif
if((i+ncore.gt.(ii-1)*maxw).and.
$(i+ncore.le.(ii-1)*maxw+nblock)) then
g2aa(l,f,l,i+ncore-(ii-1)*maxw)=
$g2aa(l,f,l,i+ncore-(ii-1)*maxw)-
$0.5d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)*4.d0*fd
g2aa(l,l,f,i+ncore-(ii-1)*maxw)=
$g2aa(l,l,f,i+ncore-(ii-1)*maxw)+
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)*2.d0
g2aa(f,l,l,i+ncore-(ii-1)*maxw)=
$g2aa(f,l,l,i+ncore-(ii-1)*maxw)+
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)*2.d0
endif
enddo
enddo
enddo
C ROHF case
else
do l=1,nocc
do i=1,nocc-ncore
do f=1,ncore
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblock)) then
g2aa(i+ncore,l,f,l-(ii-1)*maxw)=
$g2aa(i+ncore,l,f,l-(ii-1)*maxw)-
$0.5d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)*4.d0
$*fd
g2aa(f,l,i+ncore,l-(ii-1)*maxw)=
$g2aa(f,l,i+ncore,l-(ii-1)*maxw)-
$0.5d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)*4.d0
$*fd
g2aa(i+ncore,f,l,l-(ii-1)*maxw)=
$g2aa(i+ncore,f,l,l-(ii-1)*maxw)+
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)*2.d0
g2aa(l,f,i+ncore,l-(ii-1)*maxw)=
$g2aa(l,f,i+ncore,l-(ii-1)*maxw)+
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)*2.d0
g2aa(l,i+ncore,f,l-(ii-1)*maxw)=
$g2aa(l,i+ncore,f,l-(ii-1)*maxw)+
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)*2.d0
g2aa(f,i+ncore,l,l-(ii-1)*maxw)=
$g2aa(f,i+ncore,l,l-(ii-1)*maxw)+
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)*2.d0
endif
if((f.gt.(ii-1)*maxw).and.(f.le.(ii-1)*maxw+nblock)) then
g2aa(l,i+ncore,l,f-(ii-1)*maxw)=
$g2aa(l,i+ncore,l,f-(ii-1)*maxw)-
$0.5d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)*4.d0
$*fd
g2aa(i+ncore,l,l,f-(ii-1)*maxw)=
$g2aa(i+ncore,l,l,f-(ii-1)*maxw)+
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)*2.d0
g2aa(l,l,i+ncore,f-(ii-1)*maxw)=
$g2aa(l,l,i+ncore,f-(ii-1)*maxw)+
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)*2.d0
endif
if((i+ncore.gt.(ii-1)*maxw).and.
$(i+ncore.le.(ii-1)*maxw+nblock)) then
g2aa(l,f,l,i+ncore-(ii-1)*maxw)=
$g2aa(l,f,l,i+ncore-(ii-1)*maxw)-
$0.5d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)*4.d0
$*fd
g2aa(l,l,f,i+ncore-(ii-1)*maxw)=
$g2aa(l,l,f,i+ncore-(ii-1)*maxw)+
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)*2.d0
g2aa(f,l,l,i+ncore-(ii-1)*maxw)=
$g2aa(f,l,l,i+ncore-(ii-1)*maxw)+
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)*2.d0
endif
enddo
enddo
enddo
do m=nocc+1,nocc+nos
do i=1,nocc-ncore
do f=1,ncore
if((m.gt.(ii-1)*maxw).and.(m.le.(ii-1)*maxw+nblock)) then
g2aa(i+ncore,m,f,m-(ii-1)*maxw)=
$g2aa(i+ncore,m,f,m-(ii-1)*maxw)-
$0.5d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)*2.d0
$*fd
g2aa(f,m,i+ncore,m-(ii-1)*maxw)=
$g2aa(f,m,i+ncore,m-(ii-1)*maxw)-
$0.5d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)*2.d0
$*fd
g2aa(i+ncore,f,m,m-(ii-1)*maxw)=
$g2aa(i+ncore,f,m,m-(ii-1)*maxw)-
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)
g2aa(m,f,i+ncore,m-(ii-1)*maxw)=
$g2aa(m,f,i+ncore,m-(ii-1)*maxw)-
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)
g2aa(m,i+ncore,f,m-(ii-1)*maxw)=
$g2aa(m,i+ncore,f,m-(ii-1)*maxw)-
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)
g2aa(f,i+ncore,m,m-(ii-1)*maxw)=
$g2aa(f,i+ncore,m,m-(ii-1)*maxw)-
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)
endif
if((f.gt.(ii-1)*maxw).and.(f.le.(ii-1)*maxw+nblock)) then
g2aa(m,i+ncore,m,f-(ii-1)*maxw)=
$g2aa(m,i+ncore,m,f-(ii-1)*maxw)-
$0.5d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)*2.d0
$*fd
g2aa(i+ncore,m,m,f-(ii-1)*maxw)=
$g2aa(i+ncore,m,m,f-(ii-1)*maxw)-
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)
g2aa(m,m,i+ncore,f-(ii-1)*maxw)=
$g2aa(m,m,i+ncore,f-(ii-1)*maxw)-
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)
endif
if((i+ncore.gt.(ii-1)*maxw).and.
$(i+ncore.le.(ii-1)*maxw+nblock)) then
g2aa(m,f,m,i+ncore-(ii-1)*maxw)=
$g2aa(m,f,m,i+ncore-(ii-1)*maxw)-
$0.5d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)*2.d0
$*fd
g2aa(m,m,f,i+ncore-(ii-1)*maxw)=
$g2aa(m,m,f,i+ncore-(ii-1)*maxw)-
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)
g2aa(f,m,m,i+ncore-(ii-1)*maxw)=
$g2aa(f,m,m,i+ncore-(ii-1)*maxw)-
$0.25d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfrohf(i,f)
endif
enddo
enddo
enddo
endif
endif
c UHF case
if(scftype.eq.'uhf ') then
call dfillzero(g2bb,maxw*(nbasis+ncore)**3)
call dfillzero(g2ab,maxw*(nbasis+ncore)**3)
c bbbb list
1054 read(densfile,*) sum2,i,k,j,l
if((i+ncore.gt.(ii-1)*maxw).and.
$(i+ncore.le.(ii-1)*maxw+nblock)) then
g2bb(l+ncore,j+ncore,k+ncore,i-(ii-1)*maxw+ncore)=sum2
g2bb(k+ncore,j+ncore,l+ncore,i-(ii-1)*maxw+ncore)=sum2
endif
if((j+ncore.gt.(ii-1)*maxw).and.
$(j+ncore.le.(ii-1)*maxw+nblock)) then
g2bb(k+ncore,i+ncore,l+ncore,j-(ii-1)*maxw+ncore)=sum2
g2bb(l+ncore,i+ncore,k+ncore,j-(ii-1)*maxw+ncore)=sum2
endif
if((k+ncore.gt.(ii-1)*maxw).and.
$(k+ncore.le.(ii-1)*maxw+nblock)) then
g2bb(j+ncore,l+ncore,i+ncore,k-(ii-1)*maxw+ncore)=sum2
g2bb(i+ncore,l+ncore,j+ncore,k-(ii-1)*maxw+ncore)=sum2
endif
if((l+ncore.gt.(ii-1)*maxw).and.
$(l+ncore.le.(ii-1)*maxw+nblock)) then
g2bb(i+ncore,k+ncore,j+ncore,l-(ii-1)*maxw+ncore)=sum2
g2bb(j+ncore,k+ncore,i+ncore,l-(ii-1)*maxw+ncore)=sum2
endif
if(l.ne.0) goto 1054
c build the frozen parts
if(core.ne.'corr '.and.ncore.ne.0) then
do f=1,ncore
do g=1,ncore
if(g.ne.f) then
if((f.gt.(ii-1)*maxw).and.(f.le.(ii-1)*maxw+nblock))
$g2bb(g,f,g,f-(ii-1)*maxw)=1.d0
if((g.gt.(ii-1)*maxw).and.(g.le.(ii-1)*maxw+nblock)) then
g2bb(g,f,f,g-(ii-1)*maxw)=-0.5d0
g2bb(f,f,g,g-(ii-1)*maxw)=-0.5d0
endif
endif
enddo
enddo
do f=1,ncore
do p=ncore+1,ncore+nbasis
do q=ncore+1,ncore+nbasis
if((f.gt.(ii-1)*maxw).and.(f.le.(ii-1)*maxw+nblock)) then
g2bb(p,f,q,f-(ii-1)*maxw)=db(p,q)
g2bb(p,q,f,f-(ii-1)*maxw)=-0.5d0*db(p,q)
g2bb(f,q,p,f-(ii-1)*maxw)=-0.5d0*db(p,q)
endif
if((q.gt.(ii-1)*maxw).and.(q.le.(ii-1)*maxw+nblock)) then
g2bb(f,p,f,q-(ii-1)*maxw)=db(p,q)
g2bb(p,f,f,q-(ii-1)*maxw)=-0.5d0*db(p,q)
g2bb(f,f,p,q-(ii-1)*maxw)=-0.5d0*db(p,q)
endif
enddo
enddo
enddo
c Extra term in the gradient expression
do l=1,nbe+ncore
do i=1,nbe
do f=1,ncore
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblock)) then
g2bb(i+ncore,l,f,l-(ii-1)*maxw)=
$g2bb(i+ncore,l,f,l-(ii-1)*maxw)-
$0.5d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)*4.d0*fd
g2bb(f,l,i+ncore,l-(ii-1)*maxw)=
$g2bb(f,l,i+ncore,l-(ii-1)*maxw)-
$0.5d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)*4.d0*fd
g2bb(i+ncore,f,l,l-(ii-1)*maxw)=
$g2bb(i+ncore,f,l,l-(ii-1)*maxw)+
$0.25d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)*2.d0
g2bb(l,f,i+ncore,l-(ii-1)*maxw)=
$g2bb(l,f,i+ncore,l-(ii-1)*maxw)+
$0.25d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)*2.d0
g2bb(l,i+ncore,f,l-(ii-1)*maxw)=
$g2bb(l,i+ncore,f,l-(ii-1)*maxw)+
$0.25d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)*2.d0
g2bb(f,i+ncore,l,l-(ii-1)*maxw)=
$g2bb(f,i+ncore,l,l-(ii-1)*maxw)+
$0.25d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)*2.d0
endif
if((f.gt.(ii-1)*maxw).and.(f.le.(ii-1)*maxw+nblock)) then
g2bb(l,i+ncore,l,f-(ii-1)*maxw)=
$g2bb(l,i+ncore,l,f-(ii-1)*maxw)-
$0.5d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)*4.d0*fd
g2bb(i+ncore,l,l,f-(ii-1)*maxw)=
$g2bb(i+ncore,l,l,f-(ii-1)*maxw)+
$0.25d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)*2.d0
g2bb(l,l,i+ncore,f-(ii-1)*maxw)=
$g2bb(l,l,i+ncore,f-(ii-1)*maxw)+
$0.25d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)*2.d0
endif
if((i+ncore.gt.(ii-1)*maxw).and.
$(i+ncore.le.(ii-1)*maxw+nblock)) then
g2bb(l,f,l,i+ncore-(ii-1)*maxw)=
$g2bb(l,f,l,i+ncore-(ii-1)*maxw)-
$0.5d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)*4.d0*fd
g2bb(l,l,f,i+ncore-(ii-1)*maxw)=
$g2bb(l,l,f,i+ncore-(ii-1)*maxw)+
$0.25d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)*2.d0
g2bb(f,l,l,i+ncore-(ii-1)*maxw)=
$g2bb(f,l,l,i+ncore-(ii-1)*maxw)+
$0.25d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)*2.d0
endif
enddo
enddo
enddo
endif
C abab list
1055 read(densfile,*) sum2,i,k,j,l
if((k+ncore.gt.(ii-1)*maxw).and.
$(k+ncore.le.(ii-1)*maxw+nblock)) then
g2ab(j+ncore,l+ncore,i+ncore,k-(ii-1)*maxw+ncore)=sum2
g2ab(i+ncore,l+ncore,j+ncore,k-(ii-1)*maxw+ncore)=sum2
endif
if((l+ncore.gt.(ii-1)*maxw).and.
$(l+ncore.le.(ii-1)*maxw+nblock)) then
g2ab(i+ncore,k+ncore,j+ncore,l-(ii-1)*maxw+ncore)=sum2
g2ab(j+ncore,k+ncore,i+ncore,l-(ii-1)*maxw+ncore)=sum2
endif
if(l.ne.0) goto 1055
call dscal(maxw*(nbasis+ncore)**3,fd,g2ab,1)
c build the frozen parts
if(core.ne.'corr '.and.ncore.ne.0) then
do f=1,ncore
if((f.gt.(ii-1)*maxw).and.(f.le.(ii-1)*maxw+nblock))
$g2ab(f,f,f,f-(ii-1)*maxw)=1.d0
enddo
do f=1,ncore
do g=1,ncore
if((g.ne.f).and.(f.gt.(ii-1)*maxw).and.
$(f.le.(ii-1)*maxw+nblock)) g2ab(g,f,g,f-(ii-1)*maxw)=1.d0
enddo
enddo
do f=1,ncore
do p=ncore+1,ncore+nbasis
do q=ncore+1,ncore+nbasis
if((f.gt.(ii-1)*maxw).and.(f.le.(ii-1)*maxw+nblock))
$g2ab(p,f,q,f-(ii-1)*maxw)=da(p,q)
if((q.gt.(ii-1)*maxw).and.(q.le.(ii-1)*maxw+nblock))
$g2ab(f,p,f,q-(ii-1)*maxw)=db(p,q)
enddo
enddo
enddo
c Extra term in the gradient expression
do l=1,nbe+ncore
do i=1,nal
do f=1,ncore
if((l.gt.(ii-1)*maxw).and.(l.le.(ii-1)*maxw+nblock)) then
g2ab(i+ncore,l,f,l-(ii-1)*maxw)=
$g2ab(i+ncore,l,f,l-(ii-1)*maxw)-
$0.5d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)*4.d0*fd
g2ab(f,l,i+ncore,l-(ii-1)*maxw)=
$g2ab(f,l,i+ncore,l-(ii-1)*maxw)-
$0.5d0/(fockmoa(f,f)-fockmoa(i+ncore,i+ncore))*xintfa(i,f)*4.d0*fd
endif
enddo
enddo
enddo
do l=1,nal+ncore
do i=1,nbe
do f=1,ncore
if((f.gt.(ii-1)*maxw).and.(f.le.(ii-1)*maxw+nblock))
$g2ab(l,i+ncore,l,f-(ii-1)*maxw)=g2ab(l,i+ncore,l,f-(ii-1)*maxw)-
$0.5d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)*4.d0*fd
if((i+ncore.gt.(ii-1)*maxw).and.
$(i+ncore.le.(ii-1)*maxw+nblock)) g2ab(l,f,l,i+ncore-(ii-1)*maxw)=
$g2ab(l,f,l,i+ncore-(ii-1)*maxw)-
$0.5d0/(fockmob(f,f)-fockmob(i+ncore,i+ncore))*xintfb(i,f)*4.d0*fd
enddo
enddo
enddo
endif
endif
c loop on the AO blocks
do jj=1,(nbasis+ncore)/maxw+1
c block size
if(jj*maxw.le.nbasis+ncore) then
nblockao=maxw
else
nblockao=mod(nbasis+ncore,maxw)
endif
if(nblockao.gt.0) then
c aaaa
if(ii.gt.1) call readgao(scrfile5,jj,gaoa,maxw,nbasis+ncore)
call transf2(g2aa,gaoa,g2int,g2int2,ca,ca,ii,jj,maxw,nblock,
$nblockao)
call writegao(scrfile5,jj,gaoa,maxw,nbasis+ncore)
C UHF case
if(scftype.eq.'uhf ') then
c bbbb
if(ii.gt.1) call readgao(scrfile5,jj+(nbasis+ncore)/maxw+1,
$gaob,maxw,nbasis+ncore)
call transf2(g2bb,gaob,g2int,g2int2,cb,cb,ii,jj,maxw,nblock,
$nblockao)
call writegao(scrfile5,jj+(nbasis+ncore)/maxw+1,gaob,maxw,
$nbasis+ncore)
c abab
if(ii.gt.1) call readgao(scrfile5,
$jj+2*((nbasis+ncore)/maxw+1),gaoab,maxw,nbasis+ncore)
call transf2(g2ab,gaoab,g2int,g2int2,ca,cb,ii,jj,maxw,nblock,
$nblockao)
call writegao(scrfile5,jj+2*((nbasis+ncore)/maxw+1),gaoab,
$maxw,nbasis+ncore)
endif
endif
c end of the loop on the AO blocks
enddo
endif
c end of the loop on the MO blocks
enddo
c close the MO density file
close(densfile)
C UHF case -- collect the densities in one matrix
if(scftype.eq.'uhf ') then
do jj=1,(nbasis+ncore)/maxw+1
c block size
if(jj*maxw.le.nbasis+ncore) then
nblock=maxw
else
nblock=mod(nbasis+ncore,maxw)
endif
if(nblock.gt.0) then
call readgao(scrfile5,jj,gaoa,maxw,nbasis+ncore)
call readgao(scrfile5,jj+(nbasis+ncore)/maxw+1,gaob,maxw,
$nbasis+ncore)
call readgao(scrfile5,jj+2*((nbasis+ncore)/maxw+1),gaoab,maxw,
$nbasis+ncore)
call daxpy(maxw*(nbasis+ncore)**3,1.d0,gaob,1,gaoa,1)
call daxpy(maxw*(nbasis+ncore)**3,1.d0,gaoab,1,gaoa,1)
call writegao(scrfile5,jj,gaoa,maxw,nbasis+ncore)
endif
enddo
c flipped abab blocks
do jj=1,(nbasis+ncore)/maxw+1
c block size
if(jj*maxw.le.nbasis+ncore) then
nblock=maxw
else
nblock=mod(nbasis+ncore,maxw)
endif
if(nblock.gt.0) then
call readgao(scrfile5,jj+2*((nbasis+ncore)/maxw+1),gaoab,maxw,
$nbasis+ncore)
do ii=1,(nbasis+ncore)/maxw+1
c block size
if(ii*maxw.le.nbasis+ncore) then
nblockao=maxw
else
nblockao=mod(nbasis+ncore,maxw)
endif
if(nblockao.gt.0) then
call readgao(scrfile5,ii,gaoa,maxw,nbasis+ncore)
do mu=1,nbasis+ncore
do nu=1,maxw
do lambda=1,nbasis+ncore
do sigma=1,maxw
gaoa(mu,nu+(jj-1)*maxw,lambda,sigma)=
$gaoa(mu,nu+(jj-1)*maxw,lambda,sigma)+
$gaoab(lambda,sigma+(ii-1)*maxw,mu,nu)
enddo
enddo
enddo
enddo
call writegao(scrfile5,ii,gaoa,maxw,nbasis+ncore)
endif
enddo
endif
enddo
endif
1234 continue
C Atomic decomposition of coupled-cluster energy
if(popul.eq.'deco ') then
call getvar('nuc ',nuc)
call dsymm('R','U',nbasis+ncore,nbasis+ncore,1.d0,da,
$nbasis+ncore,ca,nbasis+ncore,0.d0,da2,nbasis+ncore)
call dgemm('N','T',nbasis+ncore,nbasis+ncore,nbasis+ncore,1d0,
$da2,nbasis+ncore,ca,nbasis+ncore,0.d0,daoa,nbasis+ncore)
if(scftype.eq.'uhf ') then
call dsymm('R','U',nbasis+ncore,nbasis+ncore,1.d0,db,
$nbasis+ncore,cb,nbasis+ncore,0.d0,db2,nbasis+ncore)
call dgemm('N','T',nbasis+ncore,nbasis+ncore,nbasis+ncore,1d0,
$db2,nbasis+ncore,cb,nbasis+ncore,1.d0,daoa,nbasis+ncore)
endif
open(oeintfile,file='OEINT',form='UNFORMATTED')
read(oeintfile)
call roeint(r8heap,i4heap,ha,oeintfile,nbasis+ncore)
close(oeintfile)
call dcopy((nbasis+ncore)**2,ha,1,xdip,1)
if(scftype.eq.'rhf ')
$call daxpy((nbasis+ncore)**2,1.d0,focka,1,xdip,1)
open(tedatfile,file='TEDAT',form='unformatted')
read(tedatfile) npos,intpos
close(tedatfile)
write(iout,*)
write(iout,*)
$'Atomic decomposition of coupled-cluster energy [au]:'
write(iout,*)
$' Hartree-Fock Correlation Nuclear Total'
ii=0
do iatoms=1,natoms
do i=natrange(1,iatoms)+1,natrange(2,iatoms)
ii=ii+1
bfat(ii)=iatoms
enddo
enddo
deco=0.d0
do ii=1,(nbasis+ncore)/maxw+1
if(ii*maxw.le.nbasis+ncore) then
nblock=maxw
else
nblock=mod(nbasis+ncore,maxw)
endif
if(nblock.ne.0) then
call readgao(scrfile5,ii,gaoa,maxw,nbasis+ncore)
if(scftype.eq.'uhf ') then
call readgao(scrfile5,ii+(nbasis+ncore)/maxw+1,gaob,maxw,
$nbasis+ncore)
call readgao(scrfile5,ii+2*((nbasis+ncore)/maxw+1),gaoab,
$maxw,nbasis+ncore)
endif
call intopenrsq(teintf)
IPOS=0
MAXM=IFLTLN*MIN(MAXCOR-(IMEM-IMEM1),
$ 20*(nbasis+ncore)*((nbasis+ncore)+1)/2)
DO
NNINTEG=0
IPOSLO=IPOS+1
DO WHILE(IPOS.LT.NPOS.AND.12*NNINTEG.LT.MAXM)
IPOS=IPOS+1
NNINTEG=NNINTEG+INTPOS(3,IPOS)
ENDDO
IF(12*NNINTEG.GT.MAXM) THEN
NNINTEG=NNINTEG-INTPOS(3,IPOS)
IPOS=IPOS-1
ENDIF
CALL INTREADSQ(dcore(imem),NNINTEG,TEINTF)
NNINTEG=0
DO IP=IPOSLO,IPOS
ii2 =INTPOS(1,IP)
kk2 =INTPOS(2,IP)
NINTEG=INTPOS(3,IP)
call rhfdeco(nbasis+ncore,dcore(imem),ii2,kk2,ninteg,
$nninteg,deco,gaoa,ii,nblock,maxw,bfat)
NNINTEG=NNINTEG+6*NINTEG
ENDDO
IF(IPOS.EQ.NPOS) EXIT
ENDDO
call intclose(teintf)
endif
enddo
sum3=0.d0
sum4=0.d0
sum5=0.d0
do iatoms=1,natoms
sum2=0.5d0*
$ddot((nbasis+ncore)*(natrange(2,iatoms)-natrange(1,iatoms)),
$dscftao(1,natrange(1,iatoms)+1),1,xdip(1,natrange(1,iatoms)+1),1)
if(scftype.eq.'uhf ') then
sum2=sum2+0.5d0*(
$ ddot((nbasis+ncore)*(natrange(2,iatoms)-natrange(1,iatoms)),
$dscfaao(1,natrange(1,iatoms)+1),1,focka(1,natrange(1,iatoms)+1),1)
$+ddot((nbasis+ncore)*(natrange(2,iatoms)-natrange(1,iatoms)),
$dscfbao(1,natrange(1,iatoms)+1),1,fockb(1,natrange(1,iatoms)+1),1)
$)
endif
sum=nucrepat(iatoms,natoms,coord,atchg)
sum6=
$ddot((nbasis+ncore)*(natrange(2,iatoms)-natrange(1,iatoms)),
$daoa(1,natrange(1,iatoms)+1),1,ha(1,natrange(1,iatoms)+1),1)
sum6=sum6+deco(iatoms)
write(iout,'(i4,1x,a3,4f15.8)') iatoms,atsymbol(iatoms),
$sum2,sum6-sum2,sum,sum2+sum
sum3=sum3+sum2
sum4=sum4+sum
sum5=sum5+sum6-sum2
enddo
write(iout,"(2x,66('-'))")
write(iout,"(' Total ',4f15.8)") sum3,sum5,sum4,sum3+sum4+sum5
return
endif
C *********************************************************************
C Add the solution to the occupied-virtual block of the one-electron
C density matrix
c RHF case
if(scftype.ne.'rohf ') then
do a=1,nvirtal
do i=1,nal+ncore
da(a+nal+ncore,i)=da(a+nal+ncore,i)+xinta(a,i)
da(i,a+nal+ncore)=da(i,a+nal+ncore)+xinta(a,i)
enddo
enddo
C UHF case
if(scftype.eq.'uhf ') then
do a=1,nvirtbe
do i=1,nbe+ncore
db(a+nbe+ncore,i)=db(a+nbe+ncore,i)+xintb(a,i)
db(i,a+nbe+ncore)=db(i,a+nbe+ncore)+xintb(a,i)
enddo
enddo
endif
else
C ROHF case
C Add the solution to the OE-density...
do a=1,nvirt+nos
do i=1,nocc+nos
da(a+nocc,i)=da(a+nocc,i)+xintrohf(a,i)
da(i,a+nocc)=da(i,a+nocc)+xintrohf(a,i)
enddo
enddo
endif
c **********************************************************************
C Transform the correlated one-electron aa-density matrix into AO-basis
call dsymm('R','U',nbasis+ncore,nbasis+ncore,1.d0,da,nbasis+ncore,
$ca,nbasis+ncore,0.d0,da2,nbasis+ncore)
call dgemm('N','T',nbasis+ncore,nbasis+ncore,nbasis+ncore,1.d0,
$da2,nbasis+ncore,ca,nbasis+ncore,0.d0,daoa,nbasis+ncore)
c Transform the solution of the Z-vector equation (alpha) into AO-basis
C RHF and UHF cases
call dfillzero(xintaoa,(nbasis+ncore)**2)
if(scftype.ne.'rohf ') then
do i=1,nbasis+ncore
do j=1,nbasis+ncore
sum=0.d0
do k=1,nvirtal
do l=1,nal+ncore
sum=sum+ca(i,k+nal+ncore)*ca(j,l)*xinta(k,l)
enddo
enddo
xintaoa(i,j)=sum
enddo
enddo
else
C ROHF case
do i=1,nbasis
do j=1,nbasis
sum=0.d0
do k=1,nvirt+nos
do l=1,nocc+nos
sum=sum+ca(i,k+nocc)*ca(j,l)*xintrohf(k,l)
enddo
enddo
xintaoa(i,j)=sum
enddo
enddo
endif
C UHF case
C Transform the correlated one-electron bb-density matrix into AO-basis
if(scftype.eq.'uhf ') then
call dsymm('R','U',nbasis+ncore,nbasis+ncore,1.d0,db,
$nbasis+ncore,cb,nbasis+ncore,0.d0,db2,nbasis+ncore)
call dgemm('N','T',nbasis+ncore,nbasis+ncore,nbasis+ncore,1.d0,
$db2,nbasis+ncore,cb,nbasis+ncore,0.d0,daob,nbasis+ncore)
C Transform the solution of the Z-vector equation (beta) into AO-basis
call dfillzero(xintaob,(nbasis+ncore)**2)
do i=1,nbasis+ncore
do j=1,nbasis+ncore
sum=0.d0
do k=1,nvirtbe
do l=1,nbe+ncore
sum=sum+cb(i,k+nbe+ncore)*cb(j,l)*xintb(k,l)
enddo
enddo
xintaob(i,j)=sum
enddo
enddo
endif
c **********************************************************************
C Calculate dipole moments at SCF level
c **********************************************************************
write(iout,*)
write(iout,*) 'Calculation of SCF first-order properties...'
call calcmom(r8heap,i4heap,dscftao,xdip,nbasis+ncore,oeintfile,
$angtobohr,iout,1.d0,echesu*angtobohr,eref,dipx,dipy,dipz)
c **********************************************************************
C Calculate dipole moments at CC level
c **********************************************************************
if(trim(grdens).ne.'off'.and.scftype.eq.'uhf ') then
write(iout,*)
call symmat(daoa,nbasis+ncore)
call symmat(daob,nbasis+ncore)
c open(unit=scrfile1,file='SCFDENSITIES',form='unformatted')
c call woeint(r8heap,i4heap,daoa,scrfile1,0.d0,nbasis+ncore)
c call woeint(r8heap,i4heap,daob,scrfile1,0.d0,nbasis+ncore)
c close(scrfile1)
call griddens(nbasis+ncore,daoa,daob,scftype,dft)
endif
if(scftype.eq.'uhf ') call daxpy((nbasis+ncore)**2,1.d0,daob,1,
$daoa,1)
C Symmetrize correlated one-prticle density
call symmat(daoa,nbasis+ncore)
write(iout,*)
write(iout,*) 'Calculation of CC first-order properties...'
call calcmom(r8heap,i4heap,daoa,xdip,nbasis+ncore,oeintfile,
$angtobohr,iout,1.d0,echesu*angtobohr,ecc,dipx,dipy,dipz)
if(trim(grdens).ne.'off'.and.scftype.ne.'uhf ') then
write(iout,*)
c open(unit=scrfile1,file='SCFDENSITIES',form='unformatted')
c call woeint(r8heap,i4heap,daoa,scrfile1,0.d0,nbasis+ncore)
c close(scrfile1)
call griddens(nbasis+ncore,daoa,daob,scftype,dft)
endif
if(trim(grdens).ne.'off'.and.dens.eq.1) return
C SCF Cartesian gradient
C Redefine SCF one-electron density matrix in ROHF case
c if(scftype.eq.'rohf ') then
c call dfillzero(dscfa,nbasis**2)
c do i=1,nocc+nos
c dscfa(i,i)=1.d0
c enddo
C Transform the SCF one-electron density matrix into AO-basis
c call dfillzero(dscfaao,nbasis**2)
c do i=1,nbasis
c do j=1,nbasis
c sum=0.d0
c do k=1,nbasis
c do l=1,nbasis
c sum=sum+ca(i,k)*ca(j,l)*dscfa(k,l)
c enddo
c enddo
c dscfaao(i,j)=sum
c enddo
c enddo
c endif
c **********************************************************************
C Cartesian gradient
c **********************************************************************
c Reset the grad values
do xyz=1,3
do datoms=1,natoms
do i=1,5
grads(xyz,datoms,i)=0.d0
enddo
enddo
enddo
C **********************************************************************
C NEW VERSION
C **********************************************************************
write(iout,*)
write(iout,*)
lcatoms=1
do datoms=1,natoms
if(lcatoms(datoms).eq.1) then
lcatoms(datoms)=0
write(iout,"(' Calculation of Cartesian gradient for atom',
$i4)")datoms
c loop on the blocks
call dfillzero(tdera,3*(nbasis+ncore)**2)
call dfillzero(tdera2,3*(nbasis+ncore)**2)
if(scftype.ne.'rhf ') then
call dfillzero(tderb,3*(nbasis+ncore)**2)
call dfillzero(tderb2,3*(nbasis+ncore)**2)
endif
call dfillzero(tegrad,3)
call dfillzero(tegrad2,3)
do ii=1,(nbasis+ncore)/maxw+1
c block size
if(ii*maxw.le.nbasis+ncore) then
nblock=maxw
else
nblock=mod(nbasis+ncore,maxw)
endif
if(nblock.ne.0) then
c read AO-densities
call readgao(scrfile5,ii,gaoa,maxw,nbasis+ncore)
if(scftype.eq.'uhf ') then
call readgao(scrfile5,ii+(nbasis+ncore)/maxw+1,gaob,maxw,
$nbasis+ncore)
call readgao(scrfile5,ii+2*((nbasis+ncore)/maxw+1),gaoab,
$maxw,nbasis+ncore)
endif
call dfillzero(tdera,3*(nbasis+ncore)**2)
call dfillzero(tegrad,3)
if(scftype.ne.'rhf ') call dfillzero(tderb,3*(nbasis+
$ncore)**2)
call direct_fock_build(dscfaao,dscfbao,tdera,tderb,tderb,100,
$scftype,2,1,datoms,sder,hder,dcore,.false.,gaoa,gaob,gaoab,tegrad,
$i,i,i,i,i,1.d0,.false.,dcore,0,chfx,iout,varsfile,
$icore,dcore,nbset,oeintfile,nal+ncore,scrfile3,scrfile4,maxmem,
$imem,tedatfile,dfnbasis,nbasis+ncore,0,efield,minpfile,
$.false.,1,dcore,irecln,1,.false.,ii,nblock,maxw,1.d0,daoa,1,dcore,
$0,0,i,i,i,i,i,i,i,i,0.d0,0.d0,0.d0,.false.,grads,dcore,dcore,i,
$.false.,.false.,dcore,.false.,.false.,0,i,.false.,dcore,dcore,
$.false.,.false.,nal)
c collect the tdera, tderb and tegrad terms
call daxpy(3*(nbasis+ncore)**2,1.d0,tdera,1,tdera2,1)
if(scftype.ne.'rhf ')
$call daxpy(3*(nbasis+ncore)**2,1.d0,tderb,1,tderb2,1)
do xyz=1,3
tegrad2(xyz)=tegrad2(xyz)+tegrad(xyz)
enddo
endif
c end of loop on atoms
enddo
call dcopy(3*(nbasis+ncore)**2,tdera2,1,tdera,1)
if(scftype.ne.'rhf ')
$call dcopy(3*(nbasis+ncore)**2,tderb2,1,tderb,1)
do xyz=1,3
tegrad(xyz)=tegrad2(xyz)
enddo
do xyz=1,3
sum=0.d0
sum2=0.d0
sum3=0.d0
call mxsym(tdera(1,1,xyz),nbasis+ncore)
if(scftype.ne.'rhf ') call mxsym(tderb(1,1,xyz),
$nbasis+ncore)
C **********************************************************************
C Correlated gradient
C **********************************************************************
c Two-electron integral derivatives part 1
c if(ii.eq.1) sum2=4.d0*tegrad(xyz)
sum2=4.d0*tegrad(xyz)
do mu=1,nbasis+ncore
do nu=1,nbasis+ncore
c Transform I intermediate into AO-basis
sum1=0.d0
sum4=0.d0
c only once!
c if(ii.eq.1) then
do p=1,nbasis+ncore
do q=1,nbasis+ncore
sum1=sum1+iprimea(p,q)*ca(mu,p)*ca(nu,q)
sum4=sum4+ca(mu,p)*ca(nu,q)
if(scftype.eq.'uhf ')
C UHF case (iprime = iprimea + iprimeb)
$sum1=sum1+iprimeb(p,q)*cb(mu,p)*cb(nu,q)
enddo
enddo
c Core Hamiltonian derivative (UHF case: both daoa and daob are included in daoa at this point)
sum=sum+daoa(mu,nu)*hder(mu,nu,xyz)
c endif
c Two electron integral derivatives part 2 (in every loop)
if(scftype.ne.'rohf ') then
sum2=sum2+xintaoa(mu,nu)*tdera(mu,nu,xyz)*2.d0
if(scftype.eq.'uhf ')
$sum2=sum2+xintaob(mu,nu)*tderb(mu,nu,xyz)*2.d0
else
C ROHF case
do a=1,nvirt+nos
do i=1,nocc+nos
C Occupation number
if(a.gt.nos) then
fa=0.d0
else
fa=1.d0
endif
if(i.le.nocc) then
fi=2.d0
else
fi=1.d0
endif
sum2=sum2+2.d0*xintrohf(a,i)*ca(mu,a+nocc)*ca(nu,i)*(
$tdera(mu,nu,xyz)
$-(2.d0-fa-fi)*tderb(mu,nu,xyz)
$)
enddo
enddo
endif
c Reorthonormalization gradient -- only once!
c if(ii.eq.1) sum3=sum3+sum1*sder(mu,nu,xyz)
sum3=sum3+sum1*sder(mu,nu,xyz)
enddo
enddo
sum3=-sum3
C Vnn -- only once!
vnn=0.d0
c if(ii.eq.1) then
do iatoms=1,natoms
if (iatoms.ne.datoms)then
rab3=(dsqrt((coord(1,iatoms)-coord(1,datoms))**2+
$ (coord(2,iatoms)-coord(2,datoms))**2+
$ (coord(3,iatoms)-coord(3,datoms))**2))**3
C
vnn=vnn+dfloat(atnum(iatoms))*
$ (coord(xyz,iatoms)-coord(xyz,datoms))
$ /rab3
endif
enddo
if(qmmm.eq.'amber'.and.ncent.gt.natoms) then
do iatoms=natoms+1,ncent
if (iatoms.ne.datoms)then
rab3=(dsqrt((coord(1,iatoms)-coord(1,datoms))**2+
$ (coord(2,iatoms)-coord(2,datoms))**2+
$ (coord(3,iatoms)-coord(3,datoms))**2))**3
if(dabs(rab3).gt.1d-12) then
C
vnn=vnn+atchg(iatoms)*
$ (coord(xyz,iatoms)-coord(xyz,datoms))
$ /rab3
endif
endif
enddo
endif
vnn=vnn*atnum(datoms)
c endif
c grads(xyz,datoms,1)=grads(xyz,datoms,1)+sum
c grads(xyz,datoms,2)=grads(xyz,datoms,2)+sum2
c grads(xyz,datoms,3)=grads(xyz,datoms,3)+sum3
c grads(xyz,datoms,4)=grads(xyz,datoms,4)+vnn
c grads(xyz,datoms,5)=grads(xyz,datoms,5)+sum+sum2-sum3+vnn
grads(xyz,datoms,1)=sum
grads(xyz,datoms,2)=sum2
grads(xyz,datoms,3)=sum3
grads(xyz,datoms,4)=vnn
grads(xyz,datoms,5)=sum+sum2-sum3+vnn
enddo
C Symmetry equivalent atoms
do i=1,nir
iatoms=symat(datoms,i)
if(lcatoms(iatoms).eq.1) then
lcatoms(iatoms)=0
do j=1,5
call dgemv('n',3,3,1.d0,pggen(1,1,i),3,
$grads(1,datoms,j),1,0.d0,grads(1,iatoms,j),1)
enddo
call dgemv('n',3,3,1.d0,pggen(1,1,i),3,
$efield(1,datoms),1,0.d0,efield(1,iatoms),1)
endif
enddo
call timer
write(iout,*)
endif
enddo
c endif
c end of loop on atoms
c enddo
C **********************************************************************
C END OF NEW VERSION
C **********************************************************************
C QM/MM stuff (electric field calculation for external point charges)
if(natoms.ne.ncent) then
write(iout,'(a)')
$'Calculation of electric field at atomic centers...'
call direct_fock_build(dscfaao,dscfbao,tdera,tderb,tderb,100,
$scftype,2,1,natoms,sder,hder,dcore,.false.,gaoa,gaob,gaoab,tegrad,
$i,i,i,i,i,1.d0,.false.,dcore,0,chfx,iout,varsfile,
$icore,dcore,nbset,oeintfile,nal+ncore,scrfile3,scrfile4,maxmem,
$imem,tedatfile,dfnbasis,nbasis+ncore,1,efield,minpfile,
$.false.,1,dcore,irecln,1,.false.,ii,nblock,maxw,1.d0,daoa,1,dcore,
$0,0,i,i,i,i,i,i,i,i,0.d0,0.d0,0.d0,.false.,grads,dcore,dcore,i,
$.false.,.false.,dcore,.false.,.false.,0,i,.false.,dcore,dcore,
$.false.,.false.,nal)
endif
C Print electric field
write(iout,*)
write(iout,"(' Center Electric field [au]')")
write(iout,
$"(' x y z')")
do datoms=1,ncent
if(datoms.gt.natoms) atsymbol(datoms)=' '
write(iout,'(i5,1x,a3,3f17.10)') datoms,atsymbol(datoms),
$efield(1,datoms),efield(2,datoms),efield(3,datoms)
enddo
write(iout,*)
C
write(iout,*) 'Core Hamiltonian gradient [au]:'
do datoms=1,natoms
write(iout,'(i5,1x,a3,3f17.10)') datoms,atsymbol(datoms),
$ grads(1,datoms,1),grads(2,datoms,1),grads(3,datoms,1)
enddo
write(iout,*)
write(iout,*) 'Two-electron integral gradient [au]:'
do datoms=1,natoms
write(iout,'(i5,1x,a3,3f17.10)') datoms,atsymbol(datoms),
$ grads(1,datoms,2),grads(2,datoms,2),grads(3,datoms,2)
enddo
write(iout,*)
write(iout,*) 'Reorthonormalization gradient [au]:'
do datoms=1,natoms
write(iout,'(i5,1x,a3,3f17.10)') datoms,atsymbol(datoms),
$ grads(1,datoms,3),grads(2,datoms,3),grads(3,datoms,3)
enddo
write(iout,*)
write(iout,*) 'Nuclear repulsion gradient [au]:'
do datoms=1,natoms
write(iout,'(i5,1x,a3,3f17.10)') datoms,atsymbol(datoms),
$ grads(1,datoms,4),grads(2,datoms,4),grads(3,datoms,4)
enddo
write(iout,*)
write(iout,*) 'Molecular gradient [au]:'
open(scrfile6,file='GRAD',status='unknown',form='unformatted')
do datoms=1,natoms
write(iout,'(i5,1x,a3,3f17.10)') datoms,atsymbol(datoms),
$ grads(1,datoms,5),grads(2,datoms,5),grads(3,datoms,5)
write(scrfile6)
$ grads(1,datoms,5),grads(2,datoms,5),grads(3,datoms,5)
enddo
close(scrfile6)
if(qmmm.eq.'amber') then
call write_dat_file_for_amber(finalener,selfener,2,
$natoms,ncent,grads(1,1,5),efield,dipx,dipy,dipz)
endif
C
c do xyz=1,3
c sum=0.d0
c Hmunu + Qnumu
c do mu=1,nbasis
c do nu=1,nbasis
c sum1=0.d0
c do i=1,nocc
c sum1=sum1+fockmoa(i,i)*ca(i,mu)*ca(i,nu)
c enddo
c write(40,*)sder(mu,nu,xyz)
c sum=sum+dscfaao(mu,nu)*(hder(mu,nu,xyz)+tdera(mu,nu,xyz))-
c $sum1*sder(mu,nu,xyz)
c enddo
c enddo
c sum=2.d0*sum
c Vnn
c vnn=0.d0
c do iatoms=1,natoms
c if (iatoms.ne.datoms)then
c rab3=(dsqrt((coord(1,iatoms)-coord(1,datoms))**2+
c $ (coord(2,iatoms)-coord(2,datoms))**2+
c $ (coord(3,iatoms)-coord(3,datoms))**2))**3
c
c vnn=vnn+dfloat(atnum(iatoms))*
c $ (coord(xyz,iatoms)-coord(xyz,datoms))
c $ /rab3
c endif
c enddo
c vnn=vnn*atnum(datoms)
c sum=sum+vnn
c write(6,'(2f17.12)')sum,vnn
c enddo
c enddo
C
return
end
C
************************************************************************
subroutine mxsym(mx,lsize)
************************************************************************
* symmetrising a matrix
************************************************************************
implicit none
c arguments
integer lsize
real*8 mx(lsize,lsize)
c variables
integer i,j
real*8 temp
c
do i=1,lsize
do j=1,i-1
temp=(mx(i,j)+mx(j,i))/2.d0
mx(i,j)=temp
mx(j,i)=temp
enddo
enddo
c
return
end
c
************************************************************************
subroutine transf2(g,gao,gint,gint2,c1,c2,i,j,maxword,numblock,
$numblockao)
************************************************************************
* transforms two-electron densities into AO-basis
************************************************************************
#include "MRCCCOMMON"
c arguments
real*8 g(nbasis+ncore,nbasis+ncore,nbasis+ncore,maxword)
real*8 gao(nbasis+ncore,nbasis+ncore,nbasis+ncore,maxword)
real*8 gint(nbasis+ncore,nbasis+ncore,nbasis+ncore,maxword)
real*8 gint2(nbasis+ncore,nbasis+ncore,nbasis+ncore,maxword)
real*8 c1(nbasis+ncore,nbasis+ncore),c2(nbasis+ncore,nbasis+ncore)
c i: MO block index, j: AO block index
integer i,j,maxword,numblock,numblockao
c variables
integer p,q,r,z,mu,nu,lambda,sigma
real*8 sum
c
call dfillzero(gint,maxword*(nbasis+ncore)**3)
call dfillzero(gint2,maxword*(nbasis+ncore)**3)
if(i.eq.1) call dfillzero(gao,maxword*(nbasis+ncore)**3)
C
call dgemm('n','t',(nbasis+ncore)**3,numblockao,numblock,
$1.d0,g,(nbasis+ncore)**3,c2((j-1)*maxword+1,(i-1)*maxword+1),
$nbasis+ncore,0.d0,gint,(nbasis+ncore)**3)
do sigma=1,numblockao
call dgemm('n','t',(nbasis+ncore)**2,nbasis+ncore,nbasis+ncore,
$1.d0,gint(1,1,1,sigma),(nbasis+ncore)**2,c1,nbasis+ncore,0.d0,
$gint2(1,1,1,sigma),(nbasis+ncore)**2)
enddo
do sigma=1,numblockao
do lambda=1,nbasis+ncore
call dgemm('n','t',nbasis+ncore,nbasis+ncore,nbasis+ncore,
$1.d0,gint2(1,1,lambda,sigma),nbasis+ncore,c2,nbasis+ncore,0.d0,
$gint(1,1,lambda,sigma),nbasis+ncore)
enddo
enddo
do nu=1,nbasis+ncore
do lambda=1,nbasis+ncore
do sigma=1,numblockao
call dgemv('n',nbasis+ncore,nbasis+ncore,1.d0,c1,
$nbasis+ncore,gint(1,nu,lambda,sigma),1,1.d0,
$gao(1,lambda,nu,sigma),1) ! gao in (11|22) format !!!
enddo
enddo
enddo
C
return
end
C
************************************************************************
subroutine readgao(file,jj,gao,m,n)
************************************************************************
* Read from gaofile
************************************************************************
implicit none
integer file,jj,m,n,i,j
real*8 gao(n**2,n*m)
C
j=(jj-1)*m*n
do i=1,m*n
read(file,rec=j+i) gao(1:n**2,i)
enddo
C
return
end
C
************************************************************************
subroutine writegao(file,jj,gao,m,n)
************************************************************************
* Write to gaofile
************************************************************************
implicit none
integer file,jj,m,n,i,j
real*8 gao(n**2,n*m)
C
j=(jj-1)*m*n
do i=1,m*n
write(file,rec=j+i) gao(1:n**2,i)
enddo
C
return
end
C
************************************************************************
subroutine rhfdeco(nbasis,i2heap,i,k,ninteg,nninteg,deco,g,
$blockloop,blocksize,maxword,bfat)
************************************************************************
* RHF Fock-matrix construction
************************************************************************
implicit none
integer*2 i2heap(*),i,j,k,l,ssi(4)
integer nbasis,q,ninteg,nninteg,bfat(nbasis)
integer blockloop,blocksize,maxword,iatoms,jatoms,katoms,latoms
real*8 temp,g(nbasis,nbasis,nbasis,maxword),deco(*)
equivalence(temp,ssi)
C
DO Q=nninteg+1,nninteg+6*NINTEG,6
C GETTING ACTUAL INTEGRALS AND INDICES
ssi(1)=i2heap(q )
ssi(2)=i2heap(q+1)
ssi(3)=i2heap(q+2)
ssi(4)=i2heap(q+3)
j =i2heap(q+4)
l =i2heap(q+5)
c check if L is in the block -- if not: don't do anything
if((L.gt.(blockloop-1)*maxword).and.
$(L.le.(blockloop-1)*maxword+blocksize)) then
IF(I.EQ.J) TEMP=0.5D0*TEMP
IF(K.EQ.L) TEMP=0.5D0*TEMP
IF(I.EQ.K.AND.J.EQ.L) TEMP=0.5D0*TEMP
temp=temp*g(i,j,k,l-(blockloop-1)*maxword)
iatoms=bfat(i)
jatoms=bfat(j)
katoms=bfat(k)
latoms=bfat(l)
deco(iatoms)=deco(iatoms)+temp
deco(jatoms)=deco(jatoms)+temp
deco(katoms)=deco(katoms)+temp
deco(latoms)=deco(latoms)+temp
c write(6,"(4i3,3f14.10)")
c $I,K,J,L,TEMP,g(i,j,k,l-(blockloop-1)*maxword),deco
endif
ENDDO
C
return
end
C