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

5175 lines
210 KiB
Fortran
Executable File

program mp2f12
#include "MRCCCOMMON"
integer nbfc,ifocka,ifockb,iexcha,iexchb,icmoa,icmob,ieva,ievb
integer dblalloc,verbosity,ihpja,ihpjb
character(len=4) c4
character(len=5) scftype
character(len=6) core
C
call mrccini
write(iout,*) 'MP2-F12 calculation'
write(iout,*)
call memalloc
C
call getvar('nbf ',nbf)
call getvar('nal ',nal)
call getvar('nbe ',nbe)
call getvar('ncore ',ncore)
call getkey('core',4,core,6)
call getkey('scftype',7,scftype,5)
call getkey('verbosity',9,c4,4)
read(c4,"(i4)") verbosity
if(core.eq.'corr ') ncore=0
nocc=nal+nbe
nvirtal=nbf(5)-nal
nvirtbe=nbf(5)-nbe
nvirt=nvirtal+nvirtbe
nbasis=nbf(5)
nbfc=nbf(7)
dfnbasis=nbf(3)
C
ifocka=dblalloc(nbfc*nbfc)
iexcha=dblalloc(nbfc*nbfc)
ihpja =dblalloc(nbfc*nbfc)
icmoa =dblalloc(nbfc*nbfc)
ieva =dblalloc(nbfc)
if(scftype.ne.'rhf ') then
ifockb=dblalloc(nbfc*nbfc)
iexchb=dblalloc(nbfc*nbfc)
ihpjb =dblalloc(nbfc*nbfc)
icmob =dblalloc(nbfc*nbfc)
ievb =dblalloc(nbfc)
else
ifockb=ifocka
iexchb=iexcha
ihpjb =ihpja
icmob =icmoa
ievb =ieva
endif
C
call mp2f12core(nbasis,nbf(6),nbfc,scrfile1,nocc,nvirt,nal,nbe,
$nvirtal,nvirtbe,ncore,scftype,iout,dcore(ifocka),dcore(ifockb),
$dcore(iexcha),dcore(iexchb),dcore(icmoa),dcore(icmob),
$dcore(ieva),dcore(ievb),dcore,imem,verbosity,diisfile,errfile,
$ifltln,dcore(ihpja),dcore(ihpjb),ifcfile,varsfile,icore,dfnbasis,
$nbf,maxcor,scrfile2,inp,iimem,nirmax)
C
call mrccend(0)
end
C
************************************************************************
subroutine mp2f12core(nbasis,nbfa,nbfc,scrfile1,nocc,nvirt,nal,
$nbe,nvirtal,nvirtbe,ncore,scftype,iout,focka,fockb,excha,exchb,
$cmoa,cmob,eva,evb,dcore,imem,verbosity,diisfile,errfile,ifltln,
$hpja,hpjb,ifcfile,varsfile,icore,dfnbasis,nbf,maxcor,scrfile2,inp,
$iimem,nirmax)
************************************************************************
* Driver for MP2-F12 calculations
************************************************************************
implicit none
integer nbfc,scrfile1,nbasis,nbfa,i,j,nbfan,nbfcn,iout,nalc,nbec
integer nocc,nvirt,nal,nbe,nvirtal,nvirtbe,ncore,imem,dblalloc,nn
integer nbf(*),iepaa,iepbb,iepab,ifcfile,if12ij,idaa,idbb,inp,nir
integer verbosity,diisfile,errfile,ifltln,varsfile,imem1,iimem
integer icf12aa,icf12bb,icf12ab,icf12ba,nb,itscalea,itscaleb
integer icore(*),isinv,is,is12,iconv,intalloc,dfnbasis,maxcor
integer icf122aa,icf122bb,icf122ab,icf122ba,scrfile2,nirmax
integer ivaa,ivbb,ivab,nacc,nbcc
integer ivpia,ivpib,icaia,icaib,iuaia,iuaib
integer icf12r12aa,icf12r12bb,icf12r12ab,icf12r12ba
real*8 focka(nbfc*nbfc),fockb(nbfc*nbfc),eva(nbfc),evb(nbfc),eref
real*8 excha(nbfc*nbfc),exchb(nbfc*nbfc),dcore(*),emp2f12,ecabs
real*8 hpja(nbfc*nbfc),hpjb(nbfc*nbfc),cctol,tfact,ecoup,nabtol
real*8 cmoa(nbfc,nbfc),cmob(nbfc,nbfc),naftol,eps,ef12old,emp2
real*8 ovltol
character(len=4) c4,ovirt,ccprog
character(len=5) scftype
character(len=16) c1,nab,naf_f12
character(len=32) calc
C Variables for integral calculation
integer natoms,nangmax,ncontrmax,nprimmax,ncartmax,nsphermax
integer nmboys,nbfshmax,inang,incontr,inprim,igexp,igcoef,icoord
integer ictostr,icf,iboysval,iindarr,igcn,ipre,inshrange,inangmin
integer iatnum,inzipr,inzjpr,inzkpr,inzlpr,inzint,ithad,ithcf2
integer iscoord,irqqij,irqqkl,ihrec,nbset,inbf,ispctostr,iatchg
integer ncent,iebf,necpatoms,nbfatmax,ispre,idfipre,itcdfint
integer idfrqq,inatrange,idfipra,ilogkc,igck,ikp,icprea
integer ivifaa,ivifbb,ivifab,ivifba,itcf12int,icpreb,iint1
integer ivintaa,ivintbb,ivintab,ivintba,iwaa,iwbb,iwab,iwba
real*8 itol
logical cartg,lcc,lfno,lhocc
common/memcom/ imem1
C
call getkey('ccprog',6,ccprog,4)
call getvar('eref ',eref)
call getkey('calc',4,calc,32)
call getkey('ovltol',6,c1,16)
read(c1,*) ovltol
call getkey('ovirt',5,ovirt,4)
nalc=nal-ncore
nbec=nbe-ncore
lcc=trim(calc).ne.'mp2-f12'
lfno=trim(ovirt).eq.'mp2'.and.lcc.and.nalc+nbec.gt.1
c lhocc=lcc.and.ccprog.eq.'mrcc' !mrccszemet
c $.and.trim(calc).ne.'ccsd-f12'.and.trim(calc).ne.'ccsd(t)-f12'
lhocc=.false.
call getkey('cctol',5,c4,4)
read(c4,*) i
cctol=10.d0**(-i)
C nbasis -size of AO basis
C nbfa - size of CABS basis
C nbfc - size of merged AO + CABS basis
C nbfan - size of non-redundant CABS basis
C nbfcn - size of non-redundant AO + CABS basis
C Construct CABS MOs
isinv=dblalloc(nbasis*nbasis)
is=dblalloc(nbfc*nbfc)
is12=dblalloc(nbasis*nbfa)
iconv=intalloc(nbfc)
call getvar('conv ',icore(iconv))
open(scrfile1,file='S12MATCABS',form='unformatted')
call rtdmx(dcore(imem),dcore(imem),dcore(is12),scrfile1,nbasis,
$nbfa)
call roeint(dcore(imem),dcore(imem),dcore(is),scrfile1,nbfc)
close(scrfile1)
open(scrfile1,file='SROOT_AO',form='unformatted')
read(scrfile1)
call roeint(dcore(imem),dcore(imem),dcore(isinv),scrfile1,nbasis)
close(scrfile1)
call dsymm('l','l',nbasis,nbasis,1.d0,dcore(isinv),nbasis,
$dcore(isinv),nbasis,0.d0,dcore(imem),nbasis)
call dgemm('n','n',nbasis,nbfa,nbasis,1.d0,dcore(imem),nbasis,
$dcore(is12),nbasis,0.d0,dcore(imem+nbasis**2),nbasis)
call dcopy(nbasis*nbfa,dcore(imem+nbasis**2),1,dcore(is12),1)
call aoproj(nbfc,nbfa,nbasis,icore(iconv),cmoa,dcore(is12))
call dsymm('l','u',nbfc,nbfa,1.d0,dcore(is),nbfc,cmoa(1,nbasis+1),
$nbfc,0.d0,dcore(imem+nbfc*nbfc),nbfc)
call dgemm('t','n',nbfa,nbfa,nbfc,1.d0,cmoa(1,nbasis+1),nbfc,
$dcore(imem+nbfc*nbfc),nbfc,0.d0,dcore(imem),nbfa)
call dsyev('V','U',nbfa,dcore(imem),nbfa,eva,
$dcore(imem+nbfa*nbfa),3*nbfa**2,i)
nbfan=0
do while(eva(nbfan+1).lt.ovltol.and.(nbfan+1).le.nbfa)
nbfan=nbfan+1
enddo
do j=nbfan+1,nbfa
call dscal(nbfa,1.d0/dsqrt(eva(j)),dcore(imem+(j-1)*nbfa),1)
enddo
call dgemm('n','n',nbfc,nbfa-nbfan,nbfa,1.d0,cmoa(1,nbasis+1),
$nbfc,dcore(imem+nbfan*nbfa),nbfa,0.d0,dcore(imem+nbfa*nbfa),nbfc)
nbfan=nbfa-nbfan
nbfcn=nbfc-(nbfa-nbfan)
call dcopy(nbfc*nbfan,dcore(imem+nbfa*nbfa),1,cmoa(1,nbasis+1),1)
if(scftype.ne.'rhf ') cmob=cmoa
C Read HF MO coefficients
open(scrfile1,file='MOCOEF',form='unformatted')
call moread(nbfc,nbasis,scrfile1,icore(iconv),cmoa,dcore(imem),
$dcore(imem+nbasis**2))
if(scftype.ne.'rhf ') call moread(nbfc,nbasis,scrfile1,
$icore(iconv),cmob,dcore(imem),dcore(imem+nbasis**2))
close(scrfile1)
call dbldealloc(isinv)
C Read Fock-matrix
open(scrfile1,file='FOCK',form='unformatted')
call roeint(dcore(imem),dcore(imem),focka,scrfile1,nbfc)
if(scftype.ne.'rhf ') then
call roeint(dcore(imem),dcore(imem),fockb,scrfile1,nbfc)
read(scrfile1)
endif
read(scrfile1)
call roeint(dcore(imem),dcore(imem),excha,scrfile1,nbfc)
if(scftype.ne.'rhf ')
$call roeint(dcore(imem),dcore(imem),exchb,scrfile1,nbfc)
close(scrfile1)
C Canonicalize CABS virtuals
call cancabs(nbfc,nbfan,nbasis,focka,cmoa,dcore(imem+nbfan*nbfan),
$eva,dcore(imem),iout)
if(scftype.ne.'rhf ')
$call cancabs(nbfc,nbfan,nbasis,fockb,cmob,dcore(imem+nbfan*nbfan),
$evb,dcore(imem),iout)
C Determine the symmetry of the orbitals
if(lhocc) then
call getvar('nir ',nir)
ikp=intalloc(nbfc*nir)
iint1=dblalloc(nbfcn*nbfc*nir)
igck=dblalloc(nbfc*nbfc)
icprea=dblalloc(nbfc*nbfc)
ilogkc=dblalloc(nir*nbfcn)
call symorb(nbfc,nbfcn,cmoa,cmob,icore(ikp),scrfile1,scftype,
$dcore(iint1),dcore(imem),nir,inp,nal,nbe,nirmax,dcore(igck),
$dcore(ilogkc),cctol,dcore(icprea))
call dbldealloc(ikp)
endif
C Transform Fock- and exchange matrix to CABS MO basis
call transf(nbfc,nbfcn,focka,cmoa,dcore(imem),eva,.true.)
call transf(nbfc,nbfcn,excha,cmoa,dcore(imem),eva,.false.)
call dcopy(nbfcn*nbfcn,focka,1,hpja,1)
call daxpy(nbfcn*nbfcn,-1.d0,excha,1,hpja,1)
if(scftype.ne.'rhf ') then
call transf(nbfc,nbfcn,fockb,cmob,dcore(imem),evb,.true.)
call transf(nbfc,nbfcn,exchb,cmob,dcore(imem),evb,.false.)
call dcopy(nbfcn*nbfcn,fockb,1,hpjb,1)
call daxpy(nbfcn*nbfcn,-1.d0,exchb,1,hpjb,1)
endif
if(verbosity.ge.3) then
write(iout,*)
write(iout,"(' Number of alpha electrons: ',14x,i6)") nal
write(iout,"(' Number of beta electrons: ',14x,i6)") nbe
write(iout,"(' Number of core orbitals: ',14x,i6)") ncore
write(iout,"(' Number of AOs: ',14x,i6)") nbasis
write(iout,"(' Number of auxiliary functions:',14x,i6)") nbfa
write(iout,
$"(' Number of non-redundant auxiliary functions:',i6)") nbfan
write(iout,
$"(' Number of non-redundant CABS functions: ',i6)") nbfcn
endif
C Calculate CABS singles correction
write(iout,*)
write(iout,*) 'Calculating CABS singles correction...'
ecabs=0d0
call cabssing(nbfcn,nbfan,nal,nvirtal,ecabs,iout,focka,eva,
$eva(nal+1),diisfile,errfile,ifltln,cctol,dcore(imem),
$dcore(imem+nbfan*nal),dcore(imem+2*nbfan*nal),verbosity)
if(scftype.eq.'rhf ') then
ecabs=2.d0*ecabs
else
call cabssing(nbfcn,nbfan,nbe,nvirtbe,ecabs,iout,fockb,evb,
$evb(nbe+1),diisfile,errfile,ifltln,cctol,dcore(imem),
$dcore(imem+nbfan*nbe),dcore(imem+2*nbfan*nbe),verbosity)
endif
call timer
c write(iout,*)
c write(iout,"(' Reference energy [au]: 'f24.12)") eref
c write(iout,"(' CABS singles correction [au]: 'f24.12)") ecabs
c write(iout,"(' Corrected reference energy [au]: 'f24.12)")
c $eref+ecabs
c call mrccend(0)
C Allocate memory for integrals, fitting coefficients, pair energies
nb=nbasis-ncore
nn=nalc
nacc=nalc
nbcc=nbec
if(lcc) nn=nb
if(lhocc) then
nn=nbfcn-ncore
nacc=nb
nbcc=nb
endif
itscalea=dblalloc(nalc*2)
ivifaa =dblalloc(nn *nbfcn*dfnbasis)
icf12aa =dblalloc(nacc*nbfcn*dfnbasis)
iepaa=dblalloc(nalc*(nalc+1)/2)
if(scftype.ne.'rhf ') then
itscaleb=dblalloc(nbec*2)
ivifbb =dblalloc(nn *nbfcn*dfnbasis)
icf12bb =dblalloc(nacc*nbfcn*dfnbasis)
ivifab =dblalloc(nn *nbfcn*dfnbasis)
icf12ab =dblalloc(nacc*nbfcn*dfnbasis)
ivifba =dblalloc(nn *nbfcn*dfnbasis)
icf12ba =dblalloc(nacc*nbfcn*dfnbasis)
iepbb=dblalloc(nbec*(nbec+1)/2)
iepab=dblalloc(nalc*nbec)
else
itscaleb=itscalea
ivifbb= ivifaa
icf12bb= icf12aa
ivifab= ivifaa
icf12ab= icf12aa
ivifba= ivifaa
icf12ba= icf12aa
iepbb=iepaa
iepab=iepaa
endif
if(lcc) then
icf12r12aa=dblalloc(max(nb*nalc*dfnbasis,nbfan*nvirtbe*nalc))
if(scftype.ne.'rhf ') then
icf12r12bb=dblalloc(max(nb*nbec*dfnbasis,nbfan*nvirtbe*nalc))
icf12r12ab=dblalloc(max(nb*nalc*dfnbasis,nbfan*nvirtbe*nalc))
icf12r12ba=dblalloc(max(nb*nbec*dfnbasis,nbfan*nvirtbe*nalc))
else
icf12r12bb=icf12r12aa
icf12r12ab=icf12r12aa
icf12r12ba=icf12r12aa
endif
else
icf12r12aa=imem
icf12r12bb=imem
icf12r12ab=imem
icf12r12ba=imem
endif
if(lfno) then
idaa=dblalloc(nvirtal*nvirtal)
call dfillzero(dcore(idaa),nvirtal*nvirtal)
if(scftype.ne.'rhf ') then
idbb=dblalloc(nvirtbe*nvirtbe)
call dfillzero(dcore(idbb),nvirtbe*nvirtbe)
else
idbb=idaa
endif
else
idaa=imem
idbb=imem
endif ! Caution: do not allocate memory after this point!
C Read data for integral calculation
open(varsfile,file='VARS',form='UNFORMATTED')
call tedatr(varsfile,natoms,nangmax,ncontrmax,nprimmax,
$ncartmax,cartg,nsphermax,nmboys,i,itol,nbfshmax,
$inang,incontr,inprim,igexp,igcoef,icoord,ictostr,icf,iboysval,
$iindarr,igcn,ipre,inshrange,inangmin,iatnum,icore,dcore,inzipr,
$inzjpr,inzkpr,inzlpr,inzint,ithad,ithcf2,iscoord,irqqij,irqqkl,
$ihrec,iout,0,nbset,inbf,ispctostr,iatchg,ncent,iebf,necpatoms,
$i,i,i,i,i,i,.false.,nbfatmax,ispre)
close(varsfile)
C Allocate memory for integral calculation
idfipre=dblalloc((nangmax+1)*natoms)
itcdfint=dblalloc((dfnbasis+1)*dfnbasis/2)
itcf12int=dblalloc(dfnbasis*dfnbasis)
idfrqq=dblalloc((4*nangmax+1)*nprimmax*(nangmax+1)*natoms)
inatrange=intalloc(2*natoms*nbset)
call getvar('natrange ',icore(inatrange))
idfipra=dblalloc(natoms)
ilogkc=intalloc((nangmax+1)*natoms)
igck=dblalloc(nprimmax*ncontrmax*(nangmax+1)*natoms)
ikp=intalloc(nprimmax*(nangmax+1)*natoms)
icprea=dblalloc(natoms*(nangmax+1))
if(scftype.ne.'rhf ') then
icpreb=dblalloc(natoms*(nangmax+1))
else
icpreb=icprea
endif
C Calculate and transform three-center integrals for MP2-F12
icf122aa=dblalloc(nalc*nbfcn*dfnbasis)
if(scftype.ne.'rhf ') then
icf122bb=dblalloc(nbec*nbfcn*dfnbasis)
icf122ab=dblalloc(nalc*nbfcn*dfnbasis)
icf122ba=dblalloc(nbec*nbfcn*dfnbasis)
else
icf122bb=icf122aa
icf122ab=icf122aa
icf122ba=icf122aa
endif
call mp2int(nbf(1),nalc,nbec,ncore,nb,dfnbasis,nbfc,nbfcn,cmoa,
$cmob,dcore(iepaa),dcore(iepbb),dcore(iepab),focka,fockb,hpja,hpjb,
$dcore,imem,natoms,iout,imem1,maxcor,nangmax,ncontrmax,nprimmax,
$ncartmax,cartg,nsphermax,nmboys,itol,nbfshmax,icore(inang),
$icore(incontr),icore(inprim),dcore(igexp),dcore(igcoef),
$dcore(icoord),dcore(ictostr),dcore(icf),dcore(iboysval),
$icore(iindarr),icore(igcn),dcore(ipre),icore(inshrange),
$icore(ithad),dcore(ithcf2),dcore(iscoord),dcore(irqqij),
$dcore(irqqkl),dcore(ihrec),nbset,dcore(ispctostr),dcore(idfipre),
$icore,dcore(itcdfint),dcore(idfrqq),dcore(ispre),icore(inatrange),
$dcore(idfipra),icore(ilogkc),icore(ikp),dcore(igck),dcore(icprea),
$dcore(icpreb),dcore(itcf12int),dcore(ivifaa),dcore(ivifbb),
$dcore(ivifab),dcore(icf12aa),dcore(icf12bb),dcore(icf12ab),
$dcore(icf122aa),dcore(icf122bb),dcore(icf122ab),dcore(ivifba),
$dcore(icf12ba),dcore(icf122ba),scrfile2,scftype,lcc,
$dcore(icf12r12aa),dcore(icf12r12bb),dcore(icf12r12ab),
$dcore(icf12r12ba),verbosity,'mp2 ')
if(lhocc) then
call dbldealloc(icf122aa)
else
call dbldealloc(inang)
endif
C Calculate MP2-F12 pair energies
call mp2pair(nbfcn,ncore,nal,nbe,nalc,nbec,nbasis,dfnbasis,nbfan,
$iout,scftype,verbosity,eref,emp2f12,ecabs,ifcfile,dcore,imem,
$focka,fockb,eva,evb,eva(nal+1),evb(nbe+1),nvirtal,nvirtbe,excha,
$exchb,dcore(ivifaa),dcore(ivifbb),dcore(ivifab),dcore(ivifba),
$dcore(icf12aa),dcore(icf12bb),dcore(icf12ab),dcore(icf12ba),
$dcore(iepaa),dcore(iepbb),dcore(iepab),dcore(itscalea),
$dcore(itscaleb),tfact,cctol,ecoup,dcore(idaa),dcore(idbb),lfno,
$ef12old,emp2)
if(trim(calc).eq.'mp2-f12') return
C Construct frozen natural orbitals
if(lfno) then
call getkey('lnoepsv',7,c1,16)
read(c1,*) eps
eps=0.5d0*eps ! ???????????
ivaa=dblalloc(nbfcn*nbfcn)
ivpia=dblalloc((nvirtal+nbfan)**2)
icaia=dblalloc(nbfcn*nbfcn)
if(scftype.eq.'rhf ') then
ivbb=ivaa
ivpib=ivpia
else
ivbb=dblalloc(nbfcn*nbfcn)
ivpib=dblalloc((nvirtbe+nbfan)**2)
endif
ivab=dblalloc(nbfcn)
if12ij=dblalloc(max(nb*nbfcn*dfnbasis,3*nvirtbe*nvirtbe,
$nbfcn*nbfcn,3*nbfcn))
call cfno(nb,dfnbasis,nalc,nbec,ncore,focka,fockb,nbfcn,nbfan,
$scftype,dcore(ivifaa),dcore(ivifbb),dcore(ivifab),dcore(ivifba),
$dcore(icf12aa),dcore(icf12bb),dcore(icf12ab),dcore(icf12ba),
$dcore(icf12r12aa),dcore(icf12r12bb),dcore(icf12r12ab),
$dcore(icf12r12ba),iout,verbosity,dcore(idaa),dcore(idbb),
$dcore(ivab),dcore(icaia),dcore(ivaa),dcore(ivbb),dcore(if12ij),
$eps,nvirtal,nvirtbe,dcore(ivpia),dcore(ivpib),scrfile1,nbasis,
$eva,evb)
call cecoup(nbfcn,ncore,nal,nbe,nalc,nbec,nbasis,dfnbasis,nbfan,
$iout,scftype,verbosity,eref,emp2f12,dcore,imem,focka,fockb,eva,
$evb,eva(nal+1),evb(nbe+1),nvirtal,nvirtbe,dcore(ivifaa),
$dcore(ivifbb),dcore(ivifab),dcore(ivifba),dcore(icf12aa),
$dcore(icf12bb),dcore(icf12ab),dcore(icf12ba),dcore(itscalea),
$dcore(itscaleb),tfact,cctol,ecoup,ef12old,emp2)
if(.not.lhocc) call dbldealloc(idaa)
endif
C Construct natural auxiliary basis
call getkey('nab',3,nab,16)
if(trim(nab).ne.'off') then
read(nab,*) nabtol
ivaa=dblalloc(nbfan*nbfan)
ivpia=dblalloc(nbfan*nbfan)
icaia=dblalloc(nbfcn*nbfcn)
if(scftype.eq.'rhf ') then
ivbb=ivaa
icaib=icaia
else
ivbb=dblalloc(nbfan*nbfan)
icaib=dblalloc(nbfcn*nbfcn)
endif
ivab=dblalloc(nbfan)
if12ij=dblalloc(max(nb*nbfcn*dfnbasis,nbfan*nbfan,3*nbfan))
call cnab(nb,dfnbasis,nalc,nbec,ncore,focka,fockb,nbfcn,nbfan,
$scftype,dcore(ivifaa),dcore(ivifbb),dcore(ivifab),dcore(ivifba),
$dcore(icf12aa),dcore(icf12bb),dcore(icf12ab),dcore(icf12ba),iout,
$verbosity,dcore(ivaa),dcore(ivbb),dcore(ivab),dcore(ivpia),
$dcore(icaia),dcore(icaib),dcore(if12ij),nabtol)
call dbldealloc(ivaa)
endif
C Construct natural auxiliary functions
call getkey('naf_f12',7,naf_f12,16)
if(trim(naf_f12).ne.'off') then
read(naf_f12,*) naftol
ivaa=dblalloc(dfnbasis*dfnbasis)
ivab=dblalloc(dfnbasis)
if12ij=dblalloc(max(nbfcn*dfnbasis,3*dfnbasis))
call cnaf(nb,dfnbasis,nalc,nbec,nbfcn,scftype,dcore(ivifaa),
$dcore(ivifbb),dcore(ivifab),dcore(ivifba),dcore(icf12aa),
$dcore(icf12bb),dcore(icf12ab),dcore(icf12ba),dcore(icf12r12aa),
$dcore(icf12r12bb),dcore(icf12r12ab),dcore(icf12r12ba),iout,
$verbosity,dcore(ivaa),dcore(ivab),dcore(if12ij),naftol)
call dbldealloc(ivaa)
endif
C Calculate CC intermediates
if12ij=dblalloc(max(nbfcn*nbfcn,2*nbfan*nvirtbe))
ivpia=dblalloc(nb*nalc)
icaia=dblalloc(nvirtal*nalc)
iuaia=dblalloc(nvirtal*nalc)
if(scftype.eq.'rhf ') then
ivaa=dblalloc(nb*nb*(nalc+1)*nalc/2)
ivbb=ivaa
ivab=ivaa
ivpib=ivpia
icaib=icaia
iuaib=iuaia
else
ivaa=dblalloc((nb-1)*nb/2*(nalc-1)*nalc/2)
ivbb=dblalloc((nb-1)*nb/2*(nbec-1)*nbec/2)
ivab=dblalloc(nb*nb*nalc*nbec)
ivpib=dblalloc(nb*nbec)
icaib=dblalloc(nvirtbe*nbec)
iuaib=dblalloc(nvirtbe*nbec)
endif
call ccimed(nb,dfnbasis,nalc,nbec,nvirtal,nvirtbe,ncore,focka,
$fockb,nbfcn,nbfan,emp2f12,ecabs,scrfile1,scftype,dcore(ivaa),
$dcore(ivaa),dcore(ivbb),dcore(ivab),dcore(ivpia),dcore(ivpib),
$dcore(icaia),dcore(icaib),dcore(iuaia),dcore(iuaib),dcore(ivifaa),
$dcore(ivifbb),dcore(ivifab),dcore(ivifba),dcore(icf12aa),
$dcore(icf12bb),dcore(icf12ab),dcore(icf12ba),iout,dcore,imem,
$imem1,maxcor,dcore(if12ij),dcore(if12ij),dcore(if12ij),
$dcore(if12ij+nbfan*nvirtbe),verbosity,dcore(icf12r12aa),
$dcore(icf12r12bb),dcore(icf12r12ab),dcore(icf12r12ba),
$dcore(itscalea),dcore(itscaleb),tfact,dcore(iepaa),dcore(iepbb),
$dcore(iepab),ecoup)
call dbldealloc(if12ij)
C Write integral list for mrcc
if(ccprog.eq.'mrcc') then
if(lhocc)
$call mp2int(nbf(1),nacc,nbcc,ncore,nn,dfnbasis,nbfc,nbfcn,cmoa,
$cmob,dcore(iepaa),dcore(iepbb),dcore(iepab),focka,fockb,hpja,hpjb,
$dcore,imem,natoms,iout,imem1,maxcor,nangmax,ncontrmax,nprimmax,
$ncartmax,cartg,nsphermax,nmboys,itol,nbfshmax,icore(inang),
$icore(incontr),icore(inprim),dcore(igexp),dcore(igcoef),
$dcore(icoord),dcore(ictostr),dcore(icf),dcore(iboysval),
$icore(iindarr),icore(igcn),dcore(ipre),icore(inshrange),
$icore(ithad),dcore(ithcf2),dcore(iscoord),dcore(irqqij),
$dcore(irqqkl),dcore(ihrec),nbset,dcore(ispctostr),dcore(idfipre),
$icore,dcore(itcdfint),dcore(idfrqq),dcore(ispre),icore(inatrange),
$dcore(idfipra),icore(ilogkc),icore(ikp),dcore(igck),dcore(icprea),
$dcore(icpreb),dcore(itcf12int),dcore(ivifaa),dcore(ivifbb),
$dcore(ivifab),dcore(icf12aa),dcore(icf12bb),dcore(icf12ab),
$dcore(icf122aa),dcore(icf122bb),dcore(icf122ab),dcore(ivifba),
$dcore(icf12ba),dcore(icf122ba),scrfile2,scftype,lcc,
$dcore(icf12r12aa),dcore(icf12r12bb),dcore(icf12r12ab),
$dcore(icf12r12ba),verbosity,'mrcc')
call write55(inp,nn,dfnbasis,nalc,nbec,ncore,nbfcn,itol,
$dcore(ivifaa),dcore(ivifbb),dcore(ivifab),dcore(ivifba),
$dcore(icf12aa),dcore(icf12bb),dcore(icf12ab),dcore(icf12ba),focka,
$fockb,scftype,icore(iimem),dcore(imem+nn),dcore(imem+nn),lhocc,nb)
endif
C Calculate F12 T2 contribution to T3
c if(trim(calc).ne.'ccsd-f12'.and.trim(calc).ne.'ccsd(t)-f12') then
c if(lhocc) then
c ivintaa=dblalloc(nbfan*nalc*nvirtal*nvirtal)
c ivintbb=dblalloc(nbfan*nbec*nvirtbe*nvirtbe)
c ivintab=dblalloc(nbfan*nbec*nvirtal*nvirtbe)
c ivintba=dblalloc(nalc*nbfan*nvirtal*nvirtbe)
c iwaa=dblalloc(nvirtal*nbfan*nalc*nalc)
c iwbb=dblalloc(nvirtbe*nbfan*nbec*nbec)
c iwab=dblalloc(nbfan*nvirtbe*nalc*nbec)
c iwba=dblalloc(nvirtal*nbfan*nalc*nbec)
c call t3f12(nalc,nbec,nvirtal,nvirtbe,ncore,nb,nbfcn,nbfan,
c $dfnbasis,dcore(ivifaa),dcore(ivifbb),dcore(ivifab),dcore(ivifba),
c $dcore(icf12aa),dcore(icf12bb),dcore(icf12ab),dcore(icf12ba),
c $dcore(ivintaa),dcore(ivintbb),dcore(ivintab),dcore(ivintba),
c $dcore(iwaa),dcore(iwbb),dcore(iwab),dcore(iwba),itol,scrfile1)
c call dbldealloc(ivintaa)
c endif
return
end
C
************************************************************************
subroutine mp2pair(nbfcn,ncore,nal,nbe,nalc,nbec,nbasis,dfnbasis,
$nbfan,iout,scftype,verbosity,eref,ef12,ecabs,ifcfile,dcore,imem,
$focka,fockb,eoa,eob,eva,evb,nvirtal,nvirtbe,excha,exchb,vifaa,
$vifbb,vifab,vifba,cf12aa,cf12bb,cf12ab,cf12ba,epaa,epbb,epab,
$tscalea,tscaleb,tfact,cctol,ecoup,daa,dbb,lfno,ef12old,emp2)
************************************************************************
* Calculate MP2-F12 energy
************************************************************************
implicit none
integer nbfcn,nocc,nvirt,i,ncore,nal,nbe,nbasis,dfnbasis,nbfan
integer verbosity,nvirtal,nvirtbe,ifcfile,iout,nalc,nbec,ivifaaf
integer iscr,irsa,if12ij,if12sy,ivintpq,ivintao,ivintoa,icf12aaf
integer dblalloc,imem,irsb,icf12bbf,icf12abf,icf12baf,ivifbbf
integer ivifabf,ivifbaf
real*8 focka(nbfcn,nbfcn),fockb(nbfcn,nbfcn),eoa(*),eob(*),eva(*)
real*8 excha(nbfcn,nbfcn),exchb(nbfcn,nbfcn),dcore(*),evb(*),tfact
real*8 vifaa,vifbb,vifab,cf12aa,cf12bb,cf12ab,epab(*),ecoup,cctol
real*8 emp2,ecabs,ef12,et1,eref,vifba,cf12ba,epaa(*),epbb(*)
real*8 tscalea(nalc,2),tscaleb(nbec,2),ef12old
real*8 daa(nvirtal,nvirtal),dbb(nvirtbe,nvirtbe)
character(len=5) scftype
logical lfno
C Calculate MP2 singles contribution
et1=0.d0
tscalea=0.d0
tscaleb=0.d0
call mp2sing(nbfcn,ncore,nal,nvirtal,focka,eoa,eva,et1,
$tscalea(1,2))
if(scftype.eq.'rhf ') then
et1=2.d0*et1
call dscal(nalc,2.d0,tscalea(1,2),1)
else
call mp2sing(nbfcn,ncore,nbe,nvirtbe,fockb,eob,evb,et1,
$tscaleb(1,2))
endif
emp2=et1
C Calculate F12 contribution
write(iout,*)
write(iout,*) 'Calculating pair energies...'
iscr =dblalloc(nbfcn*nbfcn)
irsa =dblalloc(nbfcn*nbfcn)
if12ij =dblalloc(nbfcn*nbfcn)
if12sy =dblalloc(nbfcn*nbfcn)
ivintpq =dblalloc(nbasis*nbasis)
ivintao =dblalloc(nbfan*nal)
ivintoa =dblalloc(nal*nbfan)
icf12aaf=dblalloc(nbfcn*dfnbasis*nalc)
ivifaaf =dblalloc(nbfcn*dfnbasis*nalc)
ef12=0.d0
ecoup=0.d0
if(verbosity.gt.2) write(iout,*)
$' MOs MP2 pair energy [au] F12 pair energy [au]'
if(scftype.eq.'rhf ') then
C Closed-shell
call epair_st(nal,nalc,nbfcn,nbfan,nvirtal,ncore,nbasis,
$dfnbasis,focka,excha,eoa,eva,iout,verbosity,ef12,emp2,ecoup,
$vifaa,cf12aa,epaa,dcore(irsa),dcore(if12ij),dcore(if12sy),
$dcore(iscr),dcore(ivintpq),dcore(ivintao),dcore(ivintoa),
$dcore(icf12aaf),dcore(ivifaaf),tscalea,cctol,daa,lfno)
else
C alpha-alpha
call epair_ss(nal,nalc,nbfcn,nbfan,nvirtal,ncore,nbasis,
$dfnbasis,focka,excha,eoa,eva,iout,verbosity,ef12,emp2,ecoup,
$vifaa,cf12aa,epaa,dcore(irsa),dcore(if12ij),dcore(if12sy),
$dcore(iscr),dcore(ivintpq),dcore(ivintao),dcore(ivintoa),
$dcore(icf12aaf),dcore(ivifaaf),'a',tscalea,cctol,daa,lfno)
C beta-beta
call epair_ss(nbe,nbec,nbfcn,nbfan,nvirtbe,ncore,nbasis,
$dfnbasis,fockb,exchb,eob,evb,iout,verbosity,ef12,emp2,ecoup,
$vifbb,cf12bb,epbb,dcore(irsa),dcore(if12ij),dcore(if12sy),
$dcore(iscr),dcore(ivintpq),dcore(ivintao),dcore(ivintoa),
$dcore(icf12aaf),dcore(ivifaaf),'b',tscaleb,cctol,dbb,lfno)
C alpha-beta
irsb =dblalloc(nbfcn*nbfcn)
icf12bbf=dblalloc(nbfcn*dfnbasis*nbec)
icf12abf=dblalloc(nbfcn*dfnbasis*nalc)
icf12baf=dblalloc(nbfcn*dfnbasis*nbec)
ivifbbf =dblalloc(nbfcn*dfnbasis*nbec)
ivifabf =dblalloc(nbfcn*dfnbasis*nalc)
ivifbaf =dblalloc(nbfcn*dfnbasis*nbec)
call epair_os(nal,nbe,nalc,nbec,nbfcn,nbfan,nvirtal,nvirtbe,
$ncore,nbasis,dfnbasis,focka,fockb,excha,exchb,eoa,eva,eob,evb,
$iout,verbosity,ef12,emp2,ecoup,vifaa,vifbb,vifab,vifba,cf12aa,
$cf12bb,cf12ab,cf12ba,epab,dcore(irsa),dcore(irsb),dcore(if12ij),
$dcore(if12sy),dcore(iscr),dcore(ivintpq),dcore(ivintao),
$dcore(ivintoa),dcore(icf12aaf),dcore(icf12bbf),dcore(icf12abf),
$dcore(icf12baf),dcore(ivifaaf),dcore(ivifbbf),dcore(ivifabf),
$dcore(ivifbaf),tscalea,tscaleb,cctol,daa,dbb,lfno)
endif
call timer
C Calculating scaling parameters for (T+)
call mtscale(nalc,nbec,ncore,iout,verbosity,tscalea,tscaleb,lfno,
$scftype,cctol)
C Print results
tfact=(emp2+ef12)/emp2
write(iout,*)
write(iout,"(' Reference energy [au]: 'f24.12)") eref
write(iout,"(' CABS singles correction [au]: 'f24.12)") ecabs
write(iout,"(' Corrected reference energy [au]: 'f24.12)")
$eref+ecabs
if(scftype.eq.'rohf ')
$write(iout,"(' MP2 singles contribution [au]: 'f24.12)") et1
write(iout,"(' MP2 correlation energy [au]: 'f24.12)") emp2
write(iout,"(' MP2 energy [au]: 'f24.12)")
$eref+emp2
write(iout,"(' F12 contribution [au]: 'f24.12)") ef12
write(iout,"(' Coupling term [au]: 'f24.12)") ecoup
write(iout,"(' MP2-F12 correlation energy [au]: 'f24.12)")
$emp2+ef12
write(iout,"(' MP2-F12 energy [au]: 'f24.12)")
$eref+ecabs+emp2+ef12
C
write(iout,*)
open(ifcfile,status='unknown',file='iface',position='append')
call prtenergc('MP2-F12',eref+ecabs+emp2+ef12,eref+ecabs,1)
close(ifcfile)
ef12old=emp2+ef12
ef12=ef12-ecoup
call dbldealloc(iscr)
C
return
end
C
************************************************************************
subroutine epair_st(nal,nalc,nbfcn,nbfan,nvirtal,ncore,nbasis,
$dfnbasis,fa,exa,eoa,eva,iout,verbosity,ef12,emp2,ecoup,vifaa,
$cf12aa,epaa,rsij,f12ij,f12sy,scr,vintpq,vintao,vintoa,cf12aaf,
$vifaaf,tscalea,cctol,daa,lfno)
************************************************************************
* Calculate MP2-F12 pair energies for RHF
************************************************************************
implicit none
integer nbfcn,nbfan,nal,nalc,nvirtal,i,a,j,b,ncore,ij
integer o,p,q,nbasis,dfnbasis,iout,verbosity
real*8 eoa(*),eva(*),fa(nbfcn,nbfcn),ddot,eij,eijb,emp2,tmp,ef12
real*8 exa(nbfcn,nbfcn),e12,ecoup,eco,epaa(*),tscalea(nalc,2)
real*8 vifaa(nbfcn,dfnbasis,nalc),cf12aa(nbfcn,dfnbasis,nalc)
real*8 rsij(nbfcn,nbfcn),f12ij(nbfcn,nbfcn),f12sy(nbfcn,nbfcn)
real*8 scr(nbfcn,nbfcn),vintpq(nbasis,nbasis),vintao(nbfan,nal)
real*8 vintoa(nal,nbfan),cf12aaf(nbfcn,dfnbasis,nalc),cctol
real*8 vifaaf(nbfcn,dfnbasis,nalc),daa(nvirtal,nvirtal)
logical lfno
C Transform fitting coefficients with Fock
call dsymm('r','l',nbfcn*dfnbasis,nalc,1d0,fa(ncore+1,ncore+1),
$nbfcn,cf12aa,nbfcn*dfnbasis,0d0,cf12aaf,nbfcn*dfnbasis)
call dsymm('r','l',nbfcn*dfnbasis,nalc,1d0,fa(ncore+1,ncore+1),
$nbfcn,vifaa,nbfcn*dfnbasis,0d0,vifaaf,nbfcn*dfnbasis)
C Loop over pairs
ij=0
do i=1,nalc
do j=1,i
C Intermediate B1, B2; V, F12/r12 contribution; X, F12^2 contribution
ij=ij+1
e12=epaa(ij)
C Intermediate B, term 3
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,1d0, vifaa(1,1,i),
$nbfcn,cf12aa(1,1,j),nbfcn,0d0,f12ij,nbfcn)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,1d0,cf12aa(1,1,i),
$nbfcn, vifaa(1,1,j),nbfcn,1d0,f12ij,nbfcn)
f12sy=transpose(f12ij)
call daxpy(nbfcn*nbfcn,7d0,f12ij,1,f12sy,1)
call dgemm('n','n',nbfcn,nbfcn,nbfcn,1d0,exa,nbfcn,f12sy,
$nbfcn,0d0,scr,nbfcn)
call dgemm('n','n',nbfcn,nbfcn,nbfcn,1d0,f12sy,nbfcn,exa,
$nbfcn,1d0,scr,nbfcn)
C Intermediate B, term 4
call dgemm('n','n',nbfcn,nal,nbfcn,-1d0,fa,nbfcn,f12sy,
$nbfcn,1d0,scr,nbfcn)
call dgemm('n','n',nal,nbfcn,nbfcn,-1d0,f12sy,nbfcn,fa,
$nbfcn,1d0,scr,nbfcn)
C Intermediate B, term 5
call dgemm('n','n',nal,nbfan,nal,1d0,fa,nbfcn,
$f12sy(1,nbasis+1),nbfcn,1d0,scr(1,nbasis+1),nbfcn)
call dgemm('n','n',nbfan,nal,nal,1d0,f12sy(nbasis+1,1),
$nbfcn,fa,nbfcn,1d0,scr(nbasis+1,1),nbfcn)
C Intermediate B, term 6
call dgemm('n','n',nbasis,nvirtal,nbasis,-1d0,fa,nbfcn,
$f12sy(1,nal+1),nbfcn,1d0,scr(1,nal+1),nbfcn)
call dgemm('n','n',nvirtal,nbasis,nbasis,-1d0,
$f12sy(nal+1,1),nbfcn,fa,nbfcn,1d0,scr(nal+1,1),nbfcn)
C Intermediate B, term 7
call dgemm('n','n',nal,nbfan,nbfcn,-2d0,fa,nbfcn,
$f12sy(1,nbasis+1),nbfcn,1d0,scr(1,nbasis+1),nbfcn)
call dgemm('n','n',nbfan,nal,nbfcn,-2d0,f12sy(nbasis+1,1),
$nbfcn,fa,nbfcn,1d0,scr(nbasis+1,1),nbfcn)
C Intermediate B, term 8
rsij=0.d0
call dgemm('n','n',nvirtal,nbasis,nbfan,1d0,f12ij(nal+1,
$nbasis+1),nbfcn,fa(nbasis+1,1),nbfcn,0d0,rsij(nal+1,1),nbfcn)
call dgemm('n','n',nbasis,nvirtal,nbfan,1d0,fa(1,nbasis+1),
$nbfcn,f12ij(nbasis+1,nal+1),nbfcn,1d0,rsij(1,nal+1),nbfcn)
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a)
do a=nal+1,nbasis
scr(nal+1:nbasis,a)=scr(nal+1:nbasis,a)-
$14d0*rsij(nal+1:nbasis,a)-2d0*rsij(a,nal+1:nbasis)
enddo
C$OMP END PARALLEL DO
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(o)
do o=1,nal
scr(nal+1:nbasis,o)=scr(nal+1:nbasis,o)-
$14d0*rsij(nal+1:nbasis,o)-2d0*rsij(o,nal+1:nbasis)
enddo
C$OMP END PARALLEL DO
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,o)
do a=nal+1,nbasis
scr(1:nal,a)=scr(1:nal,a)-14d0*rsij(1:nal,a)-2d0*rsij(a,1:nal)
enddo
C$OMP END PARALLEL DO
e12=e12+ddot(nbfcn*nbfcn,f12ij,1,scr,1)
C Intermediate X, F12 contribution
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0, vifaaf(1,1,i),
$nbfcn,cf12aa (1,1,j),nbfcn,0d0,vintpq,nbasis)
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0,cf12aaf(1,1,i),
$nbfcn, vifaa (1,1,j),nbfcn,1d0,vintpq,nbasis)
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0, vifaa (1,1,i),
$nbfcn,cf12aaf(1,1,j),nbfcn,1d0,vintpq,nbasis)
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0,cf12aa (1,1,i),
$nbfcn, vifaaf(1,1,j),nbfcn,1d0,vintpq,nbasis)
call dgemm('n','t',nbfan,nal,dfnbasis,1d0,
$ vifaaf(nbasis+1,1,i),nbfcn,cf12aa (1,1,j),nbfcn,0d0,vintao,nbfan)
call dgemm('n','t',nbfan,nal,dfnbasis,1d0,
$cf12aaf(nbasis+1,1,i),nbfcn, vifaa (1,1,j),nbfcn,1d0,vintao,nbfan)
call dgemm('n','t',nbfan,nal,dfnbasis,1d0,
$ vifaa (nbasis+1,1,i),nbfcn,cf12aaf(1,1,j),nbfcn,1d0,vintao,nbfan)
call dgemm('n','t',nbfan,nal,dfnbasis,1d0,
$cf12aa (nbasis+1,1,i),nbfcn, vifaaf(1,1,j),nbfcn,1d0,vintao,nbfan)
call dgemm('n','t',nal,nbfan,dfnbasis,1d0, vifaaf(1,1,i),
$nbfcn,cf12aa (nbasis+1,1,j),nbfcn,0d0,vintoa,nal)
call dgemm('n','t',nal,nbfan,dfnbasis,1d0,cf12aaf(1,1,i),
$nbfcn, vifaa (nbasis+1,1,j),nbfcn,1d0,vintoa,nal)
call dgemm('n','t',nal,nbfan,dfnbasis,1d0, vifaa (1,1,i),
$nbfcn,cf12aaf(nbasis+1,1,j),nbfcn,1d0,vintoa,nal)
call dgemm('n','t',nal,nbfan,dfnbasis,1d0,cf12aa (1,1,i),
$nbfcn, vifaaf(nbasis+1,1,j),nbfcn,1d0,vintoa,nal)
tmp=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(p,q) REDUCTION(+:tmp)
do p=1,nbasis
do q=1,nbasis
tmp=tmp+f12sy(p,q)*vintpq(p,q)
enddo
enddo
C$OMP END PARALLEL DO
e12=e12+tmp
tmp=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,o) REDUCTION(+:tmp)
do o=1,nal
do a=nbasis+1,nbfcn
tmp=tmp+f12sy(o,a)*vintoa(o,a-nbasis)
$ +f12sy(a,o)*vintao(a-nbasis,o)
enddo
enddo
C$OMP END PARALLEL DO
e12=e12+tmp
C Intermediate V, F12 contribution
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0,vifaa(1,1,i),
$nbfcn,vifaa(1,1,j),nbfcn,0d0,vintpq,nbasis)
call dgemm('n','t',nbfan,nal,dfnbasis,1d0,vifaa(nbasis+1,1,i),
$nbfcn,vifaa(1,1,j),nbfcn,0d0,vintao,nbfan)
call dgemm('n','t',nal,nbfan,dfnbasis,1d0,vifaa(1,1,i),
$nbfcn,vifaa(nbasis+1,1,j),nbfcn,0d0,vintoa,nal)
tmp=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(p,q) REDUCTION(+:tmp)
do p=1,nbasis
do q=1,nbasis
tmp=tmp-(5d0*f12ij(p,q)-f12ij(q,p))*vintpq(p,q)
enddo
enddo
C$OMP END PARALLEL DO
e12=e12+8d0*tmp
tmp=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,o) REDUCTION(+:tmp)
do o=1,nal
do a=nbasis+1,nbfcn
tmp=tmp-(5d0*f12ij(a,o)-f12ij(o,a))*vintao(a-nbasis,o)
$ -(5d0*f12ij(o,a)-f12ij(a,o))*vintoa(o,a-nbasis)
enddo
enddo
C$OMP END PARALLEL DO
e12=e12+8d0*tmp
C Intermediate C
eij=eoa(ncore+i)+eoa(ncore+j)
eco=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,b,eijb) REDUCTION(+:eco)
do b=nal+1,nbasis
eijb=eij-eva(b-nal)
do a=nal+1,nbasis
eco=eco+(40d0*vintpq(a,b)-8d0*vintpq(b,a)
$+7d0*rsij(a,b)+rsij(b,a))*rsij(a,b)/(eijb-eva(a-nal))
enddo
enddo
C$OMP END PARALLEL DO
e12=e12+eco
eco=eco/16d0
e12=e12/16d0
C
if(i.eq.j) then
eco=0.5d0*eco
e12=0.5d0*e12
endif
ecoup=ecoup+eco
ef12=ef12+e12
C Calculate MP2 energy
tmp=0.d0
C$OMP PARALLEL DO SCHEDULE(DYNAMIC) DEFAULT(SHARED) PRIVATE(a,b,eijb)
C$OMP& REDUCTION(+:tmp)
do b=nal+1,nal+nvirtal
eijb=eij-eva(b-nal)
tmp=tmp+0.5d0*vintpq(b,b)**2/(eijb-eva(b-nal))
do a=nal+1,b-1
tmp=tmp+(
$(vintpq(a,b)-vintpq(b,a))**2+vintpq(a,b)*vintpq(b,a))/
$(eijb-eva(a-nal))
enddo
enddo
C$OMP END PARALLEL DO
tmp=4d0*tmp
if(i.eq.j) tmp=0.5d0*tmp
emp2=emp2+tmp
if(verbosity.gt.2) write(iout,"(2(i4,a1),f20.12,2x,f20.12)")
$ncore+i,' ',ncore+j,' ',tmp,e12
if(e12.ge.cctol) write(iout,*)
$'Warning! Non-negative pair energy!'
tscalea(i,1)=tscalea(i,1)+0.5d0*e12
tscalea(i,2)=tscalea(i,2)+0.5d0*tmp
tscalea(j,1)=tscalea(j,1)+0.5d0*e12
tscalea(j,2)=tscalea(j,2)+0.5d0*tmp
epaa(ij)=0.d0
if(dabs(tmp).gt.cctol) epaa(ij)=(tmp+e12)/tmp
C Calculate MP2 density matrix
if(lfno) then
C$OMP PARALLEL DO SCHEDULE(DYNAMIC) DEFAULT(SHARED) PRIVATE(a,b,eijb)
do b=nal+1,nal+nvirtal
eijb=eij-eva(b-nal)
do a=nal+1,nal+nvirtal
vintpq(a,b)=vintpq(a,b)/(eijb-eva(a-nal))
enddo
enddo
C$OMP END PARALLEL DO
if(i.eq.j) then
call dsyrk('U','N',nvirtal,nvirtal,1d0,
$vintpq(nal+1,nal+1),nbasis,1d0,daa,nvirtal)
else
call dsyrk('U','N',nvirtal,nvirtal,2d0,
$vintpq(nal+1,nal+1),nbasis,1d0,daa,nvirtal)
call dsyrk('U','T',nvirtal,nvirtal,2d0,
$vintpq(nal+1,nal+1),nbasis,1d0,daa,nvirtal)
call dgemm('n','n',nvirtal,nvirtal,nvirtal,-1d0,
$vintpq(nal+1,nal+1),nbasis,vintpq(nal+1,nal+1),nbasis,1d0,daa,
$nvirtal)
call dgemm('t','t',nvirtal,nvirtal,nvirtal,-1d0,
$vintpq(nal+1,nal+1),nbasis,vintpq(nal+1,nal+1),nbasis,1d0,daa,
$nvirtal)
endif
endif
enddo
enddo
C
return
end
C
************************************************************************
subroutine epair_os(nal,nbe,nalc,nbec,nbfcn,nbfan,nvirtal,nvirtbe,
$ncore,nbasis,dfnbasis,fa,fb,exa,exb,eoa,eva,eob,evb,iout,
$verbosity,ef12,emp2,ecoup,vifaa,vifbb,vifab,vifba,cf12aa,cf12bb,
$cf12ab,cf12ba,epab,rsa,rsb,f12ij,f12sy,scr,vintpq,vintao,vintoa,
$cf12aaf,cf12bbf,cf12abf,cf12baf,vifaaf,vifbbf,vifabf,vifbaf,
$tscalea,tscaleb,cctol,daa,dbb,lfno)
************************************************************************
* Calculate MP2-F12 opposite-spin pair energies
************************************************************************
implicit none
integer nbfcn,nbfan,nal,nbe,nalc,nbec,nvirtal,nvirtbe,ncore,k,l
integer o,s,p,q,u,w,nbasis,c,m,n,iout,verbosity,i,a,j,b,dfnbasis
real*8 eoa(*),eva(*),eob(*),evb(*),fa(nbfcn,nbfcn),fb(nbfcn,nbfcn)
real*8 exa(nbfcn,nbfcn),exb(nbfcn,nbfcn),emp2,tmp,ef12,eij,eijb
real*8 rsa(nbfcn,nbfcn),rsb(nbfcn,nbfcn),epab(nalc,nbec),cctol
real*8 vifaa(nbfcn,dfnbasis,nalc),cf12aa(nbfcn,dfnbasis,nalc)
real*8 vifbb(nbfcn,dfnbasis,nbec),cf12bb(nbfcn,dfnbasis,nbec)
real*8 vifab(nbfcn,dfnbasis,nalc),cf12ab(nbfcn,dfnbasis,nalc)
real*8 vifba(nbfcn,dfnbasis,nbec),cf12ba(nbfcn,dfnbasis,nbec)
real*8 f12ij(nbfcn,nbfcn),f12sy(nbfcn,nbfcn),e12,ecoup,eco,ddot
real*8 scr(nbfcn,nbfcn),vintpq(nbasis,nbasis),vintao(nbfan,nbe)
real*8 vintoa(nal,nbfan),tscalea(nalc,2),tscaleb(nbec,2)
real*8 vifaaf(nbfcn,dfnbasis,nalc),cf12aaf(nbfcn,dfnbasis,nalc)
real*8 vifbbf(nbfcn,dfnbasis,nbec),cf12bbf(nbfcn,dfnbasis,nbec)
real*8 vifabf(nbfcn,dfnbasis,nalc),cf12abf(nbfcn,dfnbasis,nalc)
real*8 vifbaf(nbfcn,dfnbasis,nbec),cf12baf(nbfcn,dfnbasis,nbec)
real*8 daa(nvirtal,nvirtal),dbb(nvirtbe,nvirtbe)
logical lfno
C Transform fitting coefficients with Fock
call dsymm('r','l',nbfcn*dfnbasis,nalc,1d0,fa(ncore+1,ncore+1),
$nbfcn,cf12aa,nbfcn*dfnbasis,0d0,cf12aaf,nbfcn*dfnbasis)
call dsymm('r','l',nbfcn*dfnbasis,nbec,1d0,fb(ncore+1,ncore+1),
$nbfcn,cf12bb,nbfcn*dfnbasis,0d0,cf12bbf,nbfcn*dfnbasis)
call dsymm('r','l',nbfcn*dfnbasis,nalc,1d0,fa(ncore+1,ncore+1),
$nbfcn,cf12ab,nbfcn*dfnbasis,0d0,cf12abf,nbfcn*dfnbasis)
call dsymm('r','l',nbfcn*dfnbasis,nbec,1d0,fb(ncore+1,ncore+1),
$nbfcn,cf12ba,nbfcn*dfnbasis,0d0,cf12baf,nbfcn*dfnbasis)
call dsymm('r','l',nbfcn*dfnbasis,nalc,1d0,fa(ncore+1,ncore+1),
$nbfcn,vifaa,nbfcn*dfnbasis,0d0,vifaaf,nbfcn*dfnbasis)
call dsymm('r','l',nbfcn*dfnbasis,nbec,1d0,fb(ncore+1,ncore+1),
$nbfcn,vifbb,nbfcn*dfnbasis,0d0,vifbbf,nbfcn*dfnbasis)
call dsymm('r','l',nbfcn*dfnbasis,nalc,1d0,fa(ncore+1,ncore+1),
$nbfcn,vifab,nbfcn*dfnbasis,0d0,vifabf,nbfcn*dfnbasis)
call dsymm('r','l',nbfcn*dfnbasis,nbec,1d0,fb(ncore+1,ncore+1),
$nbfcn,vifba,nbfcn*dfnbasis,0d0,vifbaf,nbfcn*dfnbasis)
C
do i=1,nalc
do j=1,nbec
C Intermediate B1, B2; V, F12/r12 contribution; X, F12^2 contribution
e12=epab(i,j)
C Intermediate B, term 3
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,3d0/8d0, vifaa(1,1,i),
$nbfcn,cf12bb(1,1,j),nbfcn,0d0,f12sy,nbfcn)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,3d0/8d0,cf12aa(1,1,i),
$nbfcn, vifbb(1,1,j),nbfcn,1d0,f12sy,nbfcn)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,1d0/8d0, vifba(1,1,j),
$nbfcn,cf12ab(1,1,i),nbfcn,1d0,f12sy,nbfcn)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,1d0/8d0,cf12ba(1,1,j),
$nbfcn, vifab(1,1,i),nbfcn,1d0,f12sy,nbfcn)
call dgemm('n','n',nbfcn,nbfcn,nbfcn,1d0,exa,nbfcn,f12sy,
$nbfcn,0d0,scr,nbfcn)
call dgemm('n','n',nbfcn,nbfcn,nbfcn,1d0,f12sy,nbfcn,exb,
$nbfcn,1d0,scr,nbfcn)
C Intermediate B, term 4
call dgemm('n','n',nbfcn,nbe,nbfcn,-1d0,fa,nbfcn,f12sy,
$nbfcn,1d0,scr,nbfcn)
call dgemm('n','n',nal,nbfcn,nbfcn,-1d0,f12sy,nbfcn,fb,
$nbfcn,1d0,scr,nbfcn)
C Intermediate B, term 5
call dgemm('n','n',nal,nbfan,nal,1d0,fa,nbfcn,
$f12sy(1,nbasis+1),nbfcn,1d0,scr(1,nbasis+1),nbfcn)
call dgemm('n','n',nbfan,nbe,nbe,1d0,f12sy(nbasis+1,1),
$nbfcn,fb,nbfcn,1d0,scr(nbasis+1,1),nbfcn)
C Intermediate B, term 6
call dgemm('n','n',nbasis,nvirtbe,nbasis,-1d0,fa,nbfcn,
$f12sy(1,nbe+1),nbfcn,1d0,scr(1,nbe+1),nbfcn)
call dgemm('n','n',nvirtal,nbasis,nbasis,-1d0,
$f12sy(nal+1,1),nbfcn,fb,nbfcn,1d0,scr(nal+1,1),nbfcn)
C Intermediate B, term 7
call dgemm('n','n',nal,nbfan,nbfcn,-2d0,fa,nbfcn,
$f12sy(1,nbasis+1),nbfcn,1d0,scr(1,nbasis+1),nbfcn)
call dgemm('n','n',nbfan,nbe,nbfcn,-2d0,f12sy(nbasis+1,1),
$nbfcn,fb,nbfcn,1d0,scr(nbasis+1,1),nbfcn)
C Intermediate B, term 8
call dgemm('n','n',nvirtal,nbasis,nbfan,1d0,f12sy(nal+1,
$nbasis+1),nbfcn,fb(nbasis+1,1),nbfcn,0d0,rsa(nal+1,1),nbfcn)
call dgemm('n','n',nbasis,nvirtbe,nbfan,1d0,fa(1,nbasis+1),
$nbfcn,f12sy(nbasis+1,nbe+1),nbfcn,0d0,rsb(1,nbe+1),nbfcn)
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(q)
do q=1,nbasis
scr(nal+1:nal+nbasis,q)=
$ scr(nal+1:nal+nbasis,q)-2d0*rsa(nal+1:nal+nbasis,q)
enddo
C$OMP END PARALLEL DO
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(b)
do b=nbe+1,nbasis
scr(1:nbasis,b)=scr(1:nbasis,b)-2d0*rsb(1:nbasis,b)
enddo
C$OMP END PARALLEL DO
e12=e12+ddot(nbfcn*nbfcn,f12sy,1,scr,1)
C Intermediate X, F12 contribution
call dgemm('n','t',nbasis,nbasis,dfnbasis,3d0/8d0,
$ vifaaf(1,1,i),nbfcn,cf12bb (1,1,j),nbfcn,0d0,vintpq,nbasis)
call dgemm('n','t',nbasis,nbasis,dfnbasis,3d0/8d0,
$cf12aaf(1,1,i),nbfcn, vifbb (1,1,j),nbfcn,1d0,vintpq,nbasis)
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0/8d0,
$ vifba (1,1,j),nbfcn,cf12abf(1,1,i),nbfcn,1d0,vintpq,nbasis)
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0/8d0,
$cf12ba (1,1,j),nbfcn, vifabf(1,1,i),nbfcn,1d0,vintpq,nbasis)
call dgemm('n','t',nbasis,nbasis,dfnbasis,3d0/8d0,
$ vifaa (1,1,i),nbfcn,cf12bbf(1,1,j),nbfcn,1d0,vintpq,nbasis)
call dgemm('n','t',nbasis,nbasis,dfnbasis,3d0/8d0,
$cf12aa (1,1,i),nbfcn, vifbbf(1,1,j),nbfcn,1d0,vintpq,nbasis)
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0/8d0,
$ vifbaf(1,1,j),nbfcn,cf12ab (1,1,i),nbfcn,1d0,vintpq,nbasis)
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0/8d0,
$cf12baf(1,1,j),nbfcn, vifab (1,1,i),nbfcn,1d0,vintpq,nbasis)
call dgemm('n','t',nbfan,nbe,dfnbasis,3d0/8d0,
$ vifaaf(nbasis+1,1,i),nbfcn,cf12bb (1,1,j),nbfcn,0d0,vintao,nbfan)
call dgemm('n','t',nbfan,nbe,dfnbasis,3d0/8d0,
$cf12aaf(nbasis+1,1,i),nbfcn, vifbb (1,1,j),nbfcn,1d0,vintao,nbfan)
call dgemm('n','t',nbfan,nbe,dfnbasis,1d0/8d0,
$ vifba (nbasis+1,1,j),nbfcn,cf12abf(1,1,i),nbfcn,1d0,vintao,nbfan)
call dgemm('n','t',nbfan,nbe,dfnbasis,1d0/8d0,
$cf12ba (nbasis+1,1,j),nbfcn, vifabf(1,1,i),nbfcn,1d0,vintao,nbfan)
call dgemm('n','t',nbfan,nbe,dfnbasis,3d0/8d0,
$ vifaa (nbasis+1,1,i),nbfcn,cf12bbf(1,1,j),nbfcn,1d0,vintao,nbfan)
call dgemm('n','t',nbfan,nbe,dfnbasis,3d0/8d0,
$cf12aa (nbasis+1,1,i),nbfcn, vifbbf(1,1,j),nbfcn,1d0,vintao,nbfan)
call dgemm('n','t',nbfan,nbe,dfnbasis,1d0/8d0,
$ vifbaf(nbasis+1,1,j),nbfcn,cf12ab (1,1,i),nbfcn,1d0,vintao,nbfan)
call dgemm('n','t',nbfan,nbe,dfnbasis,1d0/8d0,
$cf12baf(nbasis+1,1,j),nbfcn, vifab (1,1,i),nbfcn,1d0,vintao,nbfan)
call dgemm('n','t',nal,nbfan,dfnbasis,3d0/8d0,
$ vifaaf(1,1,i),nbfcn,cf12bb (nbasis+1,1,j),nbfcn,0d0,vintoa,nal)
call dgemm('n','t',nal,nbfan,dfnbasis,3d0/8d0,
$cf12aaf(1,1,i),nbfcn, vifbb (nbasis+1,1,j),nbfcn,1d0,vintoa,nal)
call dgemm('n','t',nal,nbfan,dfnbasis,1d0/8d0,
$ vifba (1,1,j),nbfcn,cf12abf(nbasis+1,1,i),nbfcn,1d0,vintoa,nal)
call dgemm('n','t',nal,nbfan,dfnbasis,1d0/8d0,
$cf12ba (1,1,j),nbfcn, vifabf(nbasis+1,1,i),nbfcn,1d0,vintoa,nal)
call dgemm('n','t',nal,nbfan,dfnbasis,3d0/8d0,
$ vifaa (1,1,i),nbfcn,cf12bbf(nbasis+1,1,j),nbfcn,1d0,vintoa,nal)
call dgemm('n','t',nal,nbfan,dfnbasis,3d0/8d0,
$cf12aa (1,1,i),nbfcn, vifbbf(nbasis+1,1,j),nbfcn,1d0,vintoa,nal)
call dgemm('n','t',nal,nbfan,dfnbasis,1d0/8d0,
$ vifbaf(1,1,j),nbfcn,cf12ab (nbasis+1,1,i),nbfcn,1d0,vintoa,nal)
call dgemm('n','t',nal,nbfan,dfnbasis,1d0/8d0,
$cf12baf(1,1,j),nbfcn, vifab (nbasis+1,1,i),nbfcn,1d0,vintoa,nal)
tmp=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(p,q) REDUCTION(+:tmp)
do p=1,nbasis
do q=1,nbasis
tmp=tmp+f12sy(p,q)*vintpq(p,q)
enddo
enddo
C$OMP END PARALLEL DO
e12=e12+tmp
tmp=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,o) REDUCTION(+:tmp)
do o=1,nbe
do a=nbasis+1,nbfcn
tmp=tmp+f12sy(a,o)*vintao(a-nbasis,o)
enddo
enddo
C$OMP END PARALLEL DO
e12=e12+tmp
tmp=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,o) REDUCTION(+:tmp)
do a=nbasis+1,nbfcn
do o=1,nal
tmp=tmp+f12sy(o,a)*vintoa(o,a-nbasis)
enddo
enddo
C$OMP END PARALLEL DO
e12=e12+tmp
C Intermediate V, F12 contribution
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0,vifaa(1,1,i),
$nbfcn,vifbb(1,1,j),nbfcn,0d0,vintpq,nbasis)
call dgemm('n','t',nbfan,nbe,dfnbasis,1d0,vifaa(nbasis+1,1,i),
$nbfcn,vifbb(1,1,j),nbfcn,0d0,vintao,nbfan)
call dgemm('n','t',nal,nbfan,dfnbasis,1d0,vifaa(1,1,i),
$nbfcn,vifbb(nbasis+1,1,j),nbfcn,0d0,vintoa,nal)
tmp=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(p,q) REDUCTION(+:tmp)
do p=1,nbasis
do q=1,nbasis
tmp=tmp-f12sy(p,q)*vintpq(p,q)
enddo
enddo
C$OMP END PARALLEL DO
e12=e12+2.d0*tmp
tmp=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,o) REDUCTION(+:tmp)
do a=nbasis+1,nbfcn
do o=1,nal
tmp=tmp-f12sy(o,a)*vintoa(o,a-nbasis)
enddo
enddo
C$OMP END PARALLEL DO
e12=e12+2.d0*tmp
tmp=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,o) REDUCTION(+:tmp)
do a=nbasis+1,nbfcn
do o=1,nbe
tmp=tmp-f12sy(a,o)*vintao(a-nbasis,o)
enddo
enddo
C$OMP END PARALLEL DO
e12=e12+2.d0*tmp
C Intermediate C
eij=eoa(ncore+i)+eob(ncore+j)
eco=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,b,eijb,tmp) REDUCTION(+:eco)
do b=nbe+1,nbasis
eijb=eij-evb(b-nbe)
do a=nal+1,nbasis
tmp=rsb(a,b)+rsa(a,b)
eco=eco+(vintpq(a,b)+0.5d0*tmp)*tmp/(eijb-eva(a-nal))
enddo
enddo
C$OMP END PARALLEL DO
eco=2d0*eco
ecoup=ecoup+eco
e12=e12+eco
ef12=ef12+e12
C Calculate MP2 energy
tmp=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,b,eijb) REDUCTION(+:tmp)
do b=nbe+1,nbasis
eijb=eij-evb(b-nbe)
do a=nal+1,nbasis
tmp=tmp+vintpq(a,b)**2/(eijb-eva(a-nal))
enddo
enddo
C$OMP END PARALLEL DO
emp2=emp2+tmp
if(verbosity.gt.2) write(iout,"(2(i4,a1),f20.12,2x,f20.12)")
$ncore+i,'a',ncore+j,'b',tmp,e12
if(e12.ge.cctol) write(iout,*)
$'Warning! Non-negative pair energy!'
tscalea(i,1)=tscalea(i,1)+0.5d0*e12
tscalea(i,2)=tscalea(i,2)+0.5d0*tmp
tscaleb(j,1)=tscaleb(j,1)+0.5d0*e12
tscaleb(j,2)=tscaleb(j,2)+0.5d0*tmp
epab(i,j)=0.d0
if(dabs(tmp).gt.cctol) epab(i,j)=(tmp+e12)/tmp
C Calculate MP2 density matrix
if(lfno) then
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,b,eijb)
do b=nbe+1,nbasis
eijb=eij-evb(b-nbe)
do a=nal+1,nbasis
vintpq(a,b)=vintpq(a,b)/(eijb-eva(a-nal))
enddo
enddo
C$OMP END PARALLEL DO
call dsyrk('U','N',nvirtal,nvirtbe,1d0,vintpq(nal+1,nbe+1),
$nbasis,1d0,daa,nvirtal)
call dsyrk('U','T',nvirtbe,nvirtal,1d0,vintpq(nal+1,nbe+1),
$nbasis,1d0,dbb,nvirtbe)
endif
enddo
enddo
C
return
end
C
************************************************************************
subroutine epair_ss(nal,nalc,nbfcn,nbfan,nvirtal,ncore,nbasis,
$dfnbasis,fa,exa,eoa,eva,iout,verbosity,ef12,emp2,ecoup,vifaa,
$cf12aa,epaa,rsij,f12ij,f12sy,scr,vintpq,vintao,vintoa,cf12aaf,
$vifaaf,spl,tscalea,cctol,daa,lfno)
************************************************************************
* Calculate MP2-F12 same-spin pair energies
************************************************************************
implicit none
integer nbfcn,nbfan,nal,nalc,nvirtal,i,a,j,b,ncore,dfnbasis,ij
integer o,p,q,nbasis,iout,verbosity
real*8 emp2,tmp,ef12,ddot,temp,eijb,ecoup,tscalea(nalc,2),cctol
real*8 eoa(*),eva(*),fa(nbfcn,nbfcn),exa(nbfcn,nbfcn),eij,e12,eco
real*8 vifaa(nbfcn,dfnbasis,nalc),cf12aa(nbfcn,dfnbasis,nalc)
real*8 rsij(nbfcn,nbfcn),f12ij(nbfcn,nbfcn),f12sy(nbfcn,nbfcn)
real*8 scr(nbfcn,nbfcn),vintpq(nbasis,nbasis),vintao(nbfan,nal)
real*8 vintoa(nal,nbfan),cf12aaf(nbfcn,dfnbasis,nalc),epaa(*)
real*8 vifaaf(nbfcn,dfnbasis,nalc),daa(nvirtal,nvirtal)
character(len=1) spl
logical lfno
C Transform fitting coefficients with Fock
call dsymm('r','l',nbfcn*dfnbasis,nalc,1d0,fa(ncore+1,ncore+1),
$nbfcn,cf12aa,nbfcn*dfnbasis,0d0,cf12aaf,nbfcn*dfnbasis)
call dsymm('r','l',nbfcn*dfnbasis,nalc,1d0,fa(ncore+1,ncore+1),
$nbfcn,vifaa,nbfcn*dfnbasis,0d0,vifaaf,nbfcn*dfnbasis)
C Loop over pairs
ij=0
do i=1,nalc
do j=1,i-1
C Intermediate B1, B2; V, F12/r12 contribution; X, F12^2 contribution
ij=ij+1
e12=epaa(ij)
C Intermediate B, term 3
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,1d0, vifaa(1,1,i),
$nbfcn,cf12aa(1,1,j),nbfcn,0d0,f12ij,nbfcn)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,1d0,cf12aa(1,1,i),
$nbfcn, vifaa(1,1,j),nbfcn,1d0,f12ij,nbfcn)
f12sy=-transpose(f12ij)
call daxpy(nbfcn*nbfcn,1d0,f12ij,1,f12sy,1)
call dgemm('n','n',nbfcn,nbfcn,nbfcn,1d0,exa,nbfcn,f12sy,
$nbfcn,0d0,scr,nbfcn)
C Intermediate B, term 4
call dgemm('n','n',nbfcn,nal,nbfcn,-1d0,fa,nbfcn,f12sy,
$nbfcn,1d0,scr,nbfcn)
C Intermediate B, term 5
call dgemm('n','n',nal,nbfan,nal,1d0,fa,nbfcn,
$f12sy(1,nbasis+1),nbfcn,1d0,scr(1,nbasis+1),nbfcn)
C Intermediate B, term 6
call dgemm('n','n',nbasis,nvirtal,nbasis,-1d0,fa,nbfcn,
$f12sy(1,nal+1),nbfcn,1d0,scr(1,nal+1),nbfcn)
C Intermediate B, term 7
call dgemm('n','n',nal,nbfan,nbfcn,-2d0,fa,nbfcn,
$f12sy(1,nbasis+1),nbfcn,1d0,scr(1,nbasis+1),nbfcn)
C Intermediate B, term 8
rsij=0.d0
call dgemm('n','n',nbasis,nvirtal,nbfan,1d0,fa(1,nbasis+1),
$nbfcn,f12sy(nbasis+1,nal+1),nbfcn,0d0,rsij(1,nal+1),nbfcn)
call daxpy(nbfcn*nvirtal,-2d0,rsij(1,nal+1),1,scr(1,nal+1),1)
e12=e12+ddot(nbfcn*nbfcn,f12sy,1,scr,1)
C Intermediate X, F12 contribution
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0, vifaaf(1,1,i),
$nbfcn,cf12aa (1,1,j),nbfcn,0d0,vintpq,nbasis)
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0,cf12aaf(1,1,i),
$nbfcn, vifaa (1,1,j),nbfcn,1d0,vintpq,nbasis)
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0, vifaa (1,1,i),
$nbfcn,cf12aaf(1,1,j),nbfcn,1d0,vintpq,nbasis)
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0,cf12aa (1,1,i),
$nbfcn, vifaaf(1,1,j),nbfcn,1d0,vintpq,nbasis)
tmp=0.d0
C$OMP PARALLEL DO SCHEDULE(DYNAMIC) DEFAULT(SHARED) PRIVATE(p,q)
C$OMP& REDUCTION(+:tmp)
do p=1,nbasis
do q=1,p-1
tmp=tmp+f12sy(p,q)*(vintpq(p,q)-vintpq(q,p))
enddo
enddo
C$OMP END PARALLEL DO
e12=e12+tmp
call dgemm('n','t',nbfan,nal,dfnbasis, 1d0,
$ vifaa (nbasis+1,1,i),nbfcn,cf12aaf(1,1,j),nbfcn,0d0,vintao,nbfan)
call dgemm('n','t',nbfan,nal,dfnbasis, 1d0,
$cf12aa (nbasis+1,1,i),nbfcn, vifaaf(1,1,j),nbfcn,1d0,vintao,nbfan)
call dgemm('n','t',nbfan,nal,dfnbasis,-1d0,
$cf12aaf(nbasis+1,1,j),nbfcn, vifaa (1,1,i),nbfcn,1d0,vintao,nbfan)
call dgemm('n','t',nbfan,nal,dfnbasis,-1d0,
$ vifaaf(nbasis+1,1,j),nbfcn,cf12aa (1,1,i),nbfcn,1d0,vintao,nbfan)
tmp=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,o) REDUCTION(+:tmp)
do o=1,nal
do a=nbasis+1,nbfcn
tmp=tmp+f12sy(a,o)*vintao(a-nbasis,o)
enddo
enddo
C$OMP END PARALLEL DO
e12=e12+tmp
call dgemm('n','t',nal,nbfan,dfnbasis, 1d0, vifaaf(1,1,i),
$nbfcn,cf12aa (nbasis+1,1,j),nbfcn,0d0,vintoa,nal)
call dgemm('n','t',nal,nbfan,dfnbasis, 1d0,cf12aaf(1,1,i),
$nbfcn, vifaa (nbasis+1,1,j),nbfcn,1d0,vintoa,nal)
call dgemm('n','t',nal,nbfan,dfnbasis,-1d0,
$cf12aa (1,1,j),nbfcn, vifaaf(nbasis+1,1,i),nbfcn,1d0,vintoa,nal)
call dgemm('n','t',nal,nbfan,dfnbasis,-1d0,
$ vifaa (1,1,j),nbfcn,cf12aaf(nbasis+1,1,i),nbfcn,1d0,vintoa,nal)
tmp=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,o) REDUCTION(+:tmp)
do o=1,nal
do a=nbasis+1,nbfcn
tmp=tmp+f12sy(o,a)*vintoa(o,a-nbasis)
enddo
enddo
C$OMP END PARALLEL DO
e12=e12+tmp
C Intermediate V, F12 contribution
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0,vifaa(1,1,i),
$nbfcn,vifaa(1,1,j),nbfcn,0d0,vintpq,nbasis)
tmp=0.d0
C$OMP PARALLEL DO SCHEDULE(DYNAMIC) DEFAULT(SHARED) PRIVATE(p,q,temp)
C$OMP& REDUCTION(+:tmp)
do q=1,nbasis
do p=1,q-1
temp=vintpq(p,q)-vintpq(q,p)
vintpq(p,q)=temp
tmp=tmp-f12sy(p,q)*temp
enddo
enddo
C$OMP END PARALLEL DO
e12=e12+8d0*tmp
call dgemm('n','t',nbfan,nal,dfnbasis, 1d0,
$vifaa(nbasis+1,1,i),nbfcn,vifaa(1,1,j),nbfcn,0d0,vintao,nbfan)
call dgemm('n','t',nbfan,nal,dfnbasis,-1d0,
$vifaa(nbasis+1,1,j),nbfcn,vifaa(1,1,i),nbfcn,1d0,vintao,nbfan)
tmp=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,o) REDUCTION(+:tmp)
do o=1,nal
do a=nbasis+1,nbfcn
tmp=tmp-f12sy(a,o)*vintao(a-nbasis,o)
enddo
enddo
C$OMP END PARALLEL DO
e12=e12+8d0*tmp
C Intermediate C
eij=eoa(ncore+i)+eoa(ncore+j)
eco=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,b,eijb,tmp) REDUCTION(+:eco)
do b=nal+1,nal+nvirtal
eijb=eij-eva(b-nal)
do a=nal+1,b-1
tmp=0.25d0*(rsij(a,b)-rsij(b,a))
eco=eco+(2d0*vintpq(a,b)+tmp)*tmp/(eijb-eva(a-nal))
enddo
enddo
C$OMP END PARALLEL DO
ecoup=ecoup+eco
e12=e12/16d0+eco
ef12=ef12+e12
C Calculate MP2 energy
tmp=0.d0
C$OMP PARALLEL DO SCHEDULE(DYNAMIC) DEFAULT(SHARED) PRIVATE(a,b,eijb)
C$OMP& REDUCTION(+:tmp)
do b=nal+1,nal+nvirtal
eijb=eij-eva(b-nal)
do a=nal+1,b-1
tmp=tmp+vintpq(a,b)**2/(eijb-eva(a-nal))
enddo
enddo
C$OMP END PARALLEL DO
emp2=emp2+tmp
if(verbosity.gt.2) write(iout,"(2(i4,a1),f20.12,2x,f20.12)")
$ncore+j,spl,ncore+i,spl,tmp,e12
if(e12.ge.cctol) write(iout,*)
$'Warning! Non-negative pair energy!'
tscalea(i,1)=tscalea(i,1)+0.5d0*e12
tscalea(i,2)=tscalea(i,2)+0.5d0*tmp
tscalea(j,1)=tscalea(j,1)+0.5d0*e12
tscalea(j,2)=tscalea(j,2)+0.5d0*tmp
epaa(ij)=0.d0
if(dabs(tmp).gt.cctol) epaa(ij)=(tmp+e12)/tmp
C Calculate MP2 density matrix
if(lfno) then
C$OMP PARALLEL DO SCHEDULE(DYNAMIC) DEFAULT(SHARED) PRIVATE(a,b,eijb)
do b=nal+1,nal+nvirtal
vintpq(b,b)=0.d0
eijb=eij-eva(b-nal)
do a=nal+1,b-1
vintpq(a,b)=vintpq(a,b)/(eijb-eva(a-nal))
vintpq(b,a)=-vintpq(a,b)
enddo
enddo
C$OMP END PARALLEL DO
call dsyrk('U','N',nvirtal,nvirtal,1d0,vintpq(nal+1,nal+1),
$nbasis,1d0,daa,nvirtal)
endif
enddo
enddo
C
return
end
C
************************************************************************
subroutine mp2sing(nbfcn,ncore,nocc,nvirt,f,epso,epsv,et1,tscale)
************************************************************************
* Calculate MP2 singles contribution
************************************************************************
implicit none
integer nbfcn,ncore,nocc,nvirt,a,i
real*8 f(nbfcn,nbfcn),epso(*),epsv(*),et1,et2,tmp
real*8 tscale(nocc-ncore)
C
et2=0.d0
C$OMP PARALLEL DO
C$OMP& DEFAULT(PRIVATE) SHARED(nocc,nvirt,ncore,f,epsv,epso,tscale)
C$OMP& REDUCTION(+:et2)
do i=ncore+1,nocc
tmp=0.d0
do a=nocc+1,nocc+nvirt
tmp=tmp+f(a,i)**2/(epso(i)-epsv(a-nocc))
enddo
et2=et2+tmp
tscale(i-ncore)=tscale(i-ncore)+tmp
enddo
C$OMP END PARALLEL DO
et1=et1+et2
C
return
end
C
************************************************************************
subroutine cabssing(nbfcn,nbfan,nocc,nvirt,ecabs,iout,f,epso,epsv,
$diisfile,errfile,ifltln,tol,ft,t0,t1,verbosity)
************************************************************************
* Calculate CABS singles correction
************************************************************************
implicit none
integer nbfcn,nbfan,nocc,nvirt,nit,a,b,c,i,iout,maxit,diisfile
integer errfile,ifltln,verbosity
parameter(maxit=30)
real*8 epso(*),epsv(*),f(nbfcn,nbfcn),dnrm2,tmp,ecabs,norm,tol
real*8 bmat(maxit,maxit),invbmat(maxit+1,maxit+1),ddot,eold,enew
real*8 ft(nbfan,nocc),t1(nbfan,nocc),t0(nbfan,nocc)
C Construct f tilde
C$OMP PARALLEL DO
C$OMP& DEFAULT(PRIVATE)
C$OMP& SHARED(nocc,nvirt,nbfcn,f,epsv,epso,ft,t0)
do b=nocc+nvirt+1,nbfcn
do i=1,nocc
tmp=f(b,i)
do a=nocc+1,nocc+nvirt
tmp=tmp-f(b,a)*f(a,i)/(epsv(a-nocc)-epso(i))
enddo
ft(b-nocc-nvirt,i)=tmp
t0(b-nocc-nvirt,i)=tmp/(epso(i)-epsv(b-nocc))
enddo
enddo
C$OMP END PARALLEL DO
C Calculate singles amplitudes
nit=0
norm=1d0
enew=ddot(nbfan*nocc,ft,1,t0,1)
eold=enew+1.d0
if(verbosity.ge.3) write(iout,*) 'Iter. Convergence'
do while(nit.lt.maxit.and.
$(norm.gt.tol.or.dabs(enew-eold).gt.tol))
nit=nit+1
C$OMP PARALLEL DO
C$OMP& DEFAULT(PRIVATE)
C$OMP& SHARED(nocc,nvirt,nbfan,f,epsv,epso,ft,t0,t1)
do b=1,nbfan
do i=1,nocc
tmp=-ft(b,i)
do c=1,nbfan
do a=nocc+1,nocc+nvirt
tmp=tmp+f(b+nocc+nvirt,a)*f(c+nocc+nvirt,a)*t0(c,i)
$/(epsv(a-nocc)-epso(i))
enddo
enddo
t1(b,i)=tmp/(epsv(b+nvirt)-epso(i))
enddo
enddo
C$OMP END PARALLEL DO
t0=t0-t1
norm=dnrm2(nbfan*nocc,t0,1)
call diis(nit,nbfan*nocc,t1,t0,maxit,diisfile,errfile,ifltln,
$bmat,invbmat)
t0=t1
eold=enew
enew=ddot(nbfan*nocc,ft,1,t0,1)
if(verbosity.ge.3) write(iout,"(i4,f15.8)") nit,norm
enddo
if(norm.gt.tol) then
write(iout,*) 'Warning! Convergence has not been achieved!'
write(iout,"(' Residual: ',f15.8)") norm
write(iout,"(' Threshold: ',f15.8)") tol
endif
ecabs=ecabs+enew
C
return
end
C
************************************************************************
subroutine transf(nbfc,nbfcn,fock,cmo,scr,ev,lll)
************************************************************************
* Transform one-electron quantity to CABS MO basis
************************************************************************
implicit none
integer nbfc,nbfcn,i
real*8 fock(*),cmo(*),scr(*),ev(*)
logical lll
C
call dsymm('l','l',nbfc,nbfcn,1d0,fock,nbfc,cmo,nbfc,0d0,scr,nbfc)
call dgemm('t','n',nbfcn,nbfcn,nbfc,1.d0,cmo,nbfc,scr,nbfc,
$0.d0,fock,nbfcn)
if(lll) then
do i=1,nbfcn
ev(i)=fock((i-1)*nbfcn+i)
enddo
endif
C
return
end
C
************************************************************************
subroutine cancabs(nbfc,nbfan,nbasis,fock,cmo,scr,ev,um,iout)
************************************************************************
* Canonicalize CABS virtuals
************************************************************************
implicit none
integer nbfc,nbfan,nbasis,iisyev,iout
real*8 fock(*),cmo(nbfc,nbfc),scr(*),ev(nbfan),um(nbfan,nbfan)
integer*4 isyev
equivalence(isyev,iisyev) !For Intel
C
call dsymm('l','l',nbfc,nbfan,1d0,fock,nbfc,cmo(1,nbasis+1),nbfc,
$0d0,scr,nbfc)
call dgemm('t','n',nbfan,nbfan,nbfc,1.d0,cmo(1,nbasis+1),nbfc,scr,
$nbfc,0.d0,um,nbfan)
call dsyev('V','U',nbfan,um,nbfan,ev,scr,3*nbfan**2,isyev)
c write(6,"(7f9.5)") ev
if(isyev.ne.0) then
write(iout,*)
$'Fatal error at the canonicalization of CABS virtuals!'
call mrccend(1)
endif
call dgemm('n','n',nbfc,nbfan,nbfan,1.d0,cmo(1,nbasis+1),nbfc,um,
$nbfan,0.d0,scr,nbfc)
call dcopy(nbfc*nbfan,scr,1,cmo(1,nbasis+1),1)
C
return
end
C
************************************************************************
subroutine moread(nbfc,nbasis,scrfile1,conv,cmo,mo,scr)
************************************************************************
* Read MOs and convert them to CABS basis MOs
************************************************************************
implicit none
integer nbfc,nbasis,conv(nbfc),scrfile1,i
real*8 cmo(nbfc,nbfc),mo(nbasis,nbasis),scr(*)
call readmo(scr,scr,mo,scrfile1,nbasis,nbasis)
do i=1,nbfc
if(conv(i).gt.0) cmo(i,1:nbasis)=mo(conv(i),1:nbasis)
enddo
C
return
end
C
************************************************************************
subroutine aoproj(nbfc,nbfa,nbasis,conv,cmoa,s12)
************************************************************************
* Project out AO basis from CABS
************************************************************************
implicit none
integer nbfc,nbfa,nbasis,i,j,conv(nbfc)
real*8 cmoa(nbfc,nbfc),s12(nbasis,nbfa)
C
cmoa=0.d0
do i=1,nbfc
if(conv(i).lt.0) then
cmoa(i,nbasis+iabs(conv(i)))=1.d0
do j=1,nbfc
if(conv(j).gt.0) cmoa(j,nbasis+iabs(conv(i)))=
$cmoa(j,nbasis+iabs(conv(i)))-s12(conv(j),iabs(conv(i)))
enddo
endif
enddo
C
return
end
C
************************************************************************
subroutine ccimed(nb,dfnbasis,nal,nbe,nval,nvbe,ncore,focka,fockb,
$nbfcn,nbfan,emp2f12,ecabs,scrfile1,scftype,vv,vaa,vbb,vab,vpia,
$vpib,caia,caib,uaia,uaib,vifaa,vifbb,vifab,vifba,cf12aa,cf12bb,
$cf12ab,cf12ba,iout,dcore,imem,imem1,maxcor,f12ij,f12pq,f12ao,
$f12ai,verbosity,cf12r12aa,cf12r12bb,cf12r12ab,cf12r12ba,tscalea,
$tscaleb,tfact,epaa,epbb,epab,ecoup)
************************************************************************
* Calculate intermediates for CCSD calculation
************************************************************************
implicit none
integer nb,i,j,k,l,p,q,r,s,a,b,nal,nbe,ncore,o,rs,nbfcn,iout,nnb
integer nval,nvbe,c,kl,ij,jk,ab,dfnbasis,scrfile1,pq,imem,dblalloc
integer icp,icm,ivs,iva,imn,rsdims,slo,sup,imem1,maxcor,maxmem
integer verbosity,nbfan,rsdima,aodim,rsdim,ijs,ija
integer icabij,icabijaa,icabijbb,icabijab
real*8 focka(nbfcn,nbfcn),fockb(nbfcn,nbfcn),emp2f12,ecabs,tmp
real*8 vifaa(nbfcn,dfnbasis,nb),cf12aa(nbfcn,dfnbasis,nal),tfact
real*8 vifbb(nbfcn,dfnbasis,nb),cf12bb(nbfcn,dfnbasis,nbe),ecoup
real*8 vifab(nbfcn,dfnbasis,nb),cf12ab(nbfcn,dfnbasis,nal)
real*8 vifba(nbfcn,dfnbasis,nb),cf12ba(nbfcn,dfnbasis,nbe)
real*8 cf12r12aa(nb,dfnbasis,nal),cf12r12ab(nb,dfnbasis,nal)
real*8 cf12r12bb(nb,dfnbasis,nbe),cf12r12ba(nb,dfnbasis,nbe)
real*8 vaa((nb-1)*nb/2,(nal-1)*nal/2),vab(nb,nb,nal,nbe),dcore(*)
real*8 vbb((nb-1)*nb/2,(nbe-1)*nbe/2),vv(nb,nb,(nal+1)*nal/2)
real*8 vpia(nb,nal),vpib(nb,nbe),uaia(nval,nal),uaib(nvbe,nbe)
real*8 caia(nval,nal),caib(nvbe,nbe),f12pq(ncore+nb,ncore+nb)
real*8 f12ij(nbfcn,nbfcn),f12ao(nbfan,nval),f12ai(nbfan,nval)
real*8 tscalea(nal),tscaleb(nbe),epaa(nal*(nal+1)/2)
real*8 epbb(nbe*(nbe+1)/2),epab(nal,nbe)
character(len=5) scftype
character(len=8) c8
C Calculate F12 intermediates
write(iout,*)
write(iout,*) 'Calculating CCSD F12 intermediates...'
open(scrfile1,file='F12INTE',form='unformatted')
write(scrfile1) ecabs,emp2f12,ecoup
nnb=ncore+nb
if(scftype.eq.'rhf ') then
C Restricted orbitals
C Intermediate V
C F12/r12 contribution
if(verbosity.ge.3) write(iout,*)
$'Intermediate V, F12/r12 contribution...'
call f12r12cs(nbfcn,nb,dfnbasis,nal,ncore,vifaa,cf12r12aa,vv,
$f12ij)
if(verbosity.ge.3) call timer
if(verbosity.ge.3) write(iout,*)
$'Intermediate V, F12 contribution...'
C Calculate PPL-like contribution to intermediate V, term 1
C Build (anti)symmetrized F12 integrals
icp=dblalloc((nnb+1)*nnb/2*(nal+1)*nal/2)
icm=dblalloc((nnb-1)*nnb/2*(nal-1)*nal/2)
ij=icm-1
pq=icp-1
do j=1,nal
do i=1,j-1
call dgemm('n','t',nnb,nnb,dfnbasis,1d0, vifaa(1,1,i),
$nbfcn,cf12aa(1,1,j),nbfcn,0d0,f12pq,nnb)
call dgemm('n','t',nnb,nnb,dfnbasis,1d0,cf12aa(1,1,i),
$nbfcn, vifaa(1,1,j),nbfcn,1d0,f12pq,nnb)
call dscal(nnb*nnb,1d0/4d0,f12pq,1)
do q=1,nnb
do p=1,q
pq=pq+1
dcore(pq)= f12pq(p,q)+f12pq(q,p)
enddo
dcore(pq)=dcore(pq)/2d0
enddo
do q=1,nnb
do p=1,q-1
ij=ij+1
dcore(ij)=(f12pq(p,q)-f12pq(q,p))/2d0
enddo
enddo
enddo
call dsyr2k('U','N',nnb,dfnbasis,0.5d0,vifaa(1,1,j),nbfcn,
$cf12aa(1,1,j),nbfcn,0.d0,f12pq,nnb)
do q=1,nnb
do p=1,q
pq=pq+1
dcore(pq)=f12pq(p,q)
enddo
dcore(pq)=dcore(pq)/2d0
enddo
enddo
C Loop over blocks
maxmem=maxcor-(imem-imem1)
slo=1
pq=0
do sup=1,nb
rsdims=(sup-slo+1)*(slo+sup)/2
rsdima=(sup-slo+1)*(slo+sup-2)/2
if(sup.eq.nb.or.(rsdims+sup+1)*((nnb+1)*nnb/2+(nal+1)*nal/2)+
$(rsdima+sup)*(nnb-1)*nnb/2.gt.maxmem) then
if(verbosity.ge.3) then
pq=pq+1
write(c8,'(i8)') pq
write(iout,*) 'Term 1, block '//trim(adjustl(c8))//'...'
endif
ivs=dblalloc((nnb+1)*nnb/2*rsdims)
iva=dblalloc((nnb-1)*nnb/2*rsdima)
imn=dblalloc((nal+1)*nal/2*rsdims)
call vterm1(nbfcn,nb,nnb,dfnbasis,nal,vifaa,dcore(icp),
$dcore(icm),dcore(ivs),dcore(iva),vv,f12pq,dcore(imn),dcore(imn),
$rsdims,rsdima,slo,sup)
call dbldealloc(ivs)
slo=sup+1
if(verbosity.ge.3) call timer
endif
enddo
call dbldealloc(icp)
C Calculate PPL-like contribution to intermediate V, term 2
C Build (anti)symmetrized F12 integrals
icp=dblalloc(nbfan*(ncore+nal)*(nal+1)*nal/2)
icm=dblalloc(nbfan*(ncore+nal)*(nal-1)*nal/2)
aodim=nbfan*(ncore+nal)
ijs=icp
ija=icm
do j=1,nal
do i=1,j-1
call dgemm('n','t',nbfan,ncore+nal,dfnbasis,0.25d0,
$ vifaa(nnb+1,1,i),nbfcn,cf12aa(1,1,j),nbfcn,0d0,dcore(ijs),nbfan)
call dgemm('n','t',nbfan,ncore+nal,dfnbasis,0.25d0,
$cf12aa(nnb+1,1,i),nbfcn, vifaa(1,1,j),nbfcn,1d0,dcore(ijs),nbfan)
call dgemm('n','t',nbfan,ncore+nal,dfnbasis,1d0,
$ vifaa(nnb+1,1,j),nbfcn,cf12aa(1,1,i),nbfcn,0d0,f12ao,nbfan)
call dgemm('n','t',nbfan,ncore+nal,dfnbasis,1d0,
$cf12aa(nnb+1,1,j),nbfcn, vifaa(1,1,i),nbfcn,1d0,f12ao,nbfan)
call dcopy(aodim,dcore(ijs),1,dcore(ija),1)
call dscal(aodim,0.5d0,dcore(ija),1)
call daxpy(aodim,-0.125d0,f12ao,1,dcore(ija),1)
call daxpy(aodim,0.25d0,f12ao,1,dcore(ijs),1)
ijs=ijs+aodim
ija=ija+aodim
enddo
call dgemm('n','t',nbfan,ncore+nal,dfnbasis,0.5d0,
$ vifaa(nnb+1,1,j),nbfcn,cf12aa(1,1,j),nbfcn,0d0,dcore(ijs),nbfan)
call dgemm('n','t',nbfan,ncore+nal,dfnbasis,0.5d0,
$cf12aa(nnb+1,1,j),nbfcn, vifaa(1,1,j),nbfcn,1d0,dcore(ijs),nbfan)
ijs=ijs+aodim
enddo
C Loop over blocks
maxmem=maxcor-(imem-imem1)
slo=1
pq=0
do sup=1,nb
rsdims=(sup-slo+1)*(slo+sup)/2
rsdima=(sup-slo+1)*(slo+sup-2)/2
if(sup.eq.nb.or.(rsdims+sup+1)*(aodim+(nal+1)*nal/2)+
$(rsdima+sup)*aodim.gt.maxmem) then
if(verbosity.ge.3) then
pq=pq+1
write(c8,'(i8)') pq
write(iout,*) 'Term 2, block '//trim(adjustl(c8))//'...'
endif
ivs=dblalloc(aodim*rsdims)
iva=dblalloc(aodim*rsdima)
imn=dblalloc((nal+1)*nal/2*rsdims)
call vterm2(nbfcn,nb,nnb,dfnbasis,nal,vifaa,dcore(icp),
$dcore(icm),dcore(ivs),dcore(iva),vv,f12ao,dcore(imn),dcore(imn),
$rsdims,rsdima,slo,sup,ncore,nbfan,aodim)
call dbldealloc(ivs)
slo=sup+1
if(verbosity.ge.3) call timer
endif
enddo
call dbldealloc(icp)
C Intermediate Vpi
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(p,i,k,tmp)
do p=1,nb
do i=1,nal
tmp=vv(p,i,i*(i-1)/2+i)
do k=1,i-1
tmp=tmp+2d0*vv(k,p,i*(i-1)/2+k)-vv(p,k,i*(i-1)/2+k)
enddo
do k=i+1,nal
tmp=tmp+2d0*vv(p,k,k*(k-1)/2+i)-vv(k,p,k*(k-1)/2+i)
enddo
vpia(p,i)=tmp
enddo
enddo
C$OMP END PARALLEL DO
C Intermediates Uai and Cai
if(verbosity.ge.3) write(iout,*) 'Intermediate U...'
iva=dblalloc(nbfan*nal)
uaia(:,:)=vpia(nal+1:nal+nval,:)
caia=0d0
do j=1,nal
do i=1,j
call dgemm('n','t',nbfan,nval,dfnbasis,1d0,
$ vifaa(ncore+nb+1,1,i),nbfcn,cf12aa(ncore+nal+1,1,j),nbfcn,0d0,
$f12ai,nbfan)
call dgemm('n','t',nbfan,nval,dfnbasis,1d0,
$cf12aa(ncore+nb+1,1,i),nbfcn, vifaa(ncore+nal+1,1,j),nbfcn,1d0,
$f12ai,nbfan)
call dgemm('n','t',nbfan,nval,dfnbasis,-1d0,
$ vifaa(ncore+nb+1,1,j),nbfcn,cf12aa(ncore+nal+1,1,i),nbfcn,0d0,
$f12ao,nbfan)
call dgemm('n','t',nbfan,nval,dfnbasis,-1d0,
$cf12aa(ncore+nb+1,1,j),nbfcn, vifaa(ncore+nal+1,1,i),nbfcn,1d0,
$f12ao,nbfan)
call daxpy(nbfan*nval,5d0,f12ai,1,f12ao,1)
call dgemv('t',nbfan,nval,1d0/8d0,f12ao,nbfan,
$focka(ncore+nb+1,ncore+i),1,1d0,caia(1,j),1)
call dgemm('n','t',nbfan,nal,dfnbasis,1d0,
$vifaa(ncore+nb+1,1,i),nbfcn,vifaa(ncore+1,1,j),nbfcn,0d0,
$dcore(iva),nbfan)
call dgemm('t','n',nval,nal,nbfan,-1d0/4d0,f12ao,nbfan,
$dcore(iva),nbfan,1d0,uaia,nval)
if(i.ne.j) then
call daxpy(nbfan*nval,-5d0,f12ai,1,f12ao,1)
call dscal(nbfan*nval,-5d0,f12ao,1)
call daxpy(nbfan*nval,-1d0,f12ai,1,f12ao,1)
call dgemv('t',nbfan,nval,1d0/8d0,f12ao,nbfan,
$focka(ncore+nb+1,ncore+j),1,1d0,caia(1,i),1)
call dgemm('n','t',nbfan,nal,dfnbasis,1d0,
$vifaa(ncore+nb+1,1,j),nbfcn,vifaa(ncore+1,1,i),nbfcn,0d0,
$dcore(iva),nbfan)
call dgemm('t','n',nval,nal,nbfan,-1d0/4d0,f12ao,nbfan,
$dcore(iva),nbfan,1d0,uaia,nval)
endif
enddo
enddo
call daxpy(nval*nal,1d0,uaia,1,caia,1)
call dbldealloc(iva)
if(verbosity.ge.3) call timer
C Intermediate Cabij
if(verbosity.ge.3) write(iout,*) 'Intermediate C...'
icabij=dblalloc(nval*nval*(nal+1)*nal/2)
ivs=dblalloc(nbfan*nval*nal)
iva=dblalloc(nbfan*nval*nal)
imn=dblalloc(nbfan*nval*nal)
call cabijcalc(nb,ncore,nal,nval,nbfcn,nbfan,dfnbasis,vifaa,
$cf12aa,cf12r12aa,dcore(ivs),dcore(iva),dcore(imn),f12ao,focka,
$dcore(icabij),vv)
call dbldealloc(ivs)
if(verbosity.ge.3) call timer
C Save intermediates
if(verbosity.ge.3) write(iout,*) 'Saving intermediates...'
call saveimcs(nb,nal,nval,scrfile1,dcore(icabij),vv,uaia,vpia,
$caia,tscalea,tfact,epaa)
call dbldealloc(icabij)
else
C Unrestricted orbitals
C Intermediate V
if(verbosity.ge.3) write(iout,*)
$'Intermediate V, F12/r12 contribution...'
call f12r12os(nbfcn,nb,dfnbasis,nal,nbe,ncore,vifaa,vifbb,vifab,
$vifba,cf12r12aa,cf12r12bb,cf12r12ab,cf12r12ba,vaa,vbb,vab,f12ij)
if(verbosity.ge.3) call timer
if(verbosity.ge.3) write(iout,*)
$'Intermediate V, F12 contribution, alpha-alpha block...'
call vaacalc(nb,nnb,ncore,nal,nbfcn,nbfan,dfnbasis,dcore,
$imem,imem1,maxcor,vifaa,cf12aa,f12pq,vaa)
if(verbosity.ge.3) call timer
if(nbe.gt.0) then
if(verbosity.ge.3) write(iout,*)
$'Intermediate V, F12 contribution, beta-beta block...'
call vaacalc(nb,nnb,ncore,nbe,nbfcn,nbfan,dfnbasis,dcore,
$imem,imem1,maxcor,vifbb,cf12bb,f12pq,vbb)
if(verbosity.ge.3) call timer
if(verbosity.ge.3) write(iout,*)
$'Intermediate V, F12 contribution, alpha-beta block...'
call vabcalc(nb,nnb,ncore,nal,nbe,nbfcn,nbfan,dfnbasis,
$dcore,imem,imem1,maxcor,vifaa,vifbb,vifab,vifba,cf12aa,cf12bb,
$cf12ab,cf12ba,vab)
endif
if(verbosity.ge.3) call timer
C Intermediate Vpi
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(p,i,k,tmp)
do p=1,nb
do i=1,nal
tmp=0.d0
do k=1,min(i-1,p-1)
tmp=tmp+vaa((p-1)*(p-2)/2+k,(i-1)*(i-2)/2+k)
enddo
do k=i+1,min(nal,p-1)
tmp=tmp-vaa((p-1)*(p-2)/2+k,(k-1)*(k-2)/2+i)
enddo
do k=max(1,p+1),i-1
tmp=tmp-vaa((k-1)*(k-2)/2+p,(i-1)*(i-2)/2+k)
enddo
do k=max(i+1,p+1),nal
tmp=tmp+vaa((k-1)*(k-2)/2+p,(k-1)*(k-2)/2+i)
enddo
do k=1,nbe
tmp=tmp+vab(p,k,i,k)
enddo
vpia(p,i)=tmp
enddo
enddo
C$OMP END PARALLEL DO
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(p,i,k,tmp)
do p=1,nb
do i=1,nbe
tmp=0.d0
do k=1,min(i-1,p-1)
tmp=tmp+vbb((p-1)*(p-2)/2+k,(i-1)*(i-2)/2+k)
enddo
do k=i+1,min(nbe,p-1)
tmp=tmp-vbb((p-1)*(p-2)/2+k,(k-1)*(k-2)/2+i)
enddo
do k=max(1,p+1),i-1
tmp=tmp-vbb((k-1)*(k-2)/2+p,(i-1)*(i-2)/2+k)
enddo
do k=max(i+1,p+1),nbe
tmp=tmp+vbb((k-1)*(k-2)/2+p,(k-1)*(k-2)/2+i)
enddo
do k=1,nal
tmp=tmp+vab(k,p,k,i)
enddo
vpib(p,i)=tmp
enddo
enddo
C$OMP END PARALLEL DO
C Intermediates Uai and Cai
if(verbosity.ge.3) write(iout,*) 'Intermediate U...'
call uaiaacalc(nb,ncore,nal,nval,nbfcn,nbfan,dfnbasis,vifaa,
$cf12aa,focka,f12ao,f12ai,vpia,uaia,caia)
call uaiaacalc(nb,ncore,nbe,nvbe,nbfcn,nbfan,dfnbasis,vifbb,
$cf12bb,fockb,f12ao,f12ai,vpib,uaib,caib)
C Alpha-beta
do j=1,nbe
do i=1,nal
call dgemm('n','t',nval,nbfan,dfnbasis,3d0/8d0,
$ vifaa(ncore+nal+1,1,i),nbfcn,cf12bb(ncore+nb+1,1,j),nbfcn,0d0,
$f12ao,nval)
call dgemm('n','t',nval,nbfan,dfnbasis,3d0/8d0,
$cf12aa(ncore+nal+1,1,i),nbfcn, vifbb(ncore+nb+1,1,j),nbfcn,1d0,
$f12ao,nval)
call dgemm('n','t',nval,nbfan,dfnbasis,1d0/8d0,
$ vifba(ncore+nal+1,1,j),nbfcn,cf12ab(ncore+nb+1,1,i),nbfcn,1d0,
$f12ao,nval)
call dgemm('n','t',nval,nbfan,dfnbasis,1d0/8d0,
$cf12ba(ncore+nal+1,1,j),nbfcn, vifab(ncore+nb+1,1,i),nbfcn,1d0,
$f12ao,nval)
call dgemv('n',nval,nbfan,1d0,f12ao,nval,
$fockb(ncore+nb+1,ncore+j),1,1d0,caia(1,i),1)
call dgemm('n','t',nal,nbfan,dfnbasis,1d0,
$vifaa(ncore+1,1,i),nbfcn,vifbb(ncore+nb+1,1,j),nbfcn,0d0,f12ai,
$nal)
call dgemm('n','t',nval,nal,nbfan,-2d0,f12ao,nval,f12ai,
$nal,1d0,uaia,nval)
enddo
enddo
C
do j=1,nbe
do i=1,nal
call dgemm('n','t',nbfan,nvbe,dfnbasis,3d0/8d0,
$ vifaa(ncore+nb+1,1,i),nbfcn,cf12bb(ncore+nbe+1,1,j),nbfcn,0d0,
$f12ao,nbfan)
call dgemm('n','t',nbfan,nvbe,dfnbasis,3d0/8d0,
$cf12aa(ncore+nb+1,1,i),nbfcn, vifbb(ncore+nbe+1,1,j),nbfcn,1d0,
$f12ao,nbfan)
call dgemm('n','t',nbfan,nvbe,dfnbasis,1d0/8d0,
$ vifba(ncore+nb+1,1,j),nbfcn,cf12ab(ncore+nbe+1,1,i),nbfcn,1d0,
$f12ao,nbfan)
call dgemm('n','t',nbfan,nvbe,dfnbasis,1d0/8d0,
$cf12ba(ncore+nb+1,1,j),nbfcn, vifab(ncore+nbe+1,1,i),nbfcn,1d0,
$f12ao,nbfan)
call dgemv('t',nbfan,nvbe,1d0,f12ao,nbfan,
$focka(ncore+nb+1,ncore+i),1,1d0,caib(1,j),1)
call dgemm('n','t',nbfan,nbe,dfnbasis,1d0,
$vifaa(ncore+nb+1,1,i),nbfcn,vifbb(ncore+1,1,j),nbfcn,0d0,f12ai,
$nbfan)
call dgemm('t','n',nvbe,nbe,nbfan,-2d0,f12ao,nbfan,f12ai,
$nbfan,1d0,uaib,nvbe)
enddo
enddo
call daxpy(nval*nal,1d0,uaia,1,caia,1)
call daxpy(nvbe*nbe,1d0,uaib,1,caib,1)
if(verbosity.ge.3) call timer
C Intermediate Cabij
if(verbosity.ge.3) write(iout,*) 'Intermediate C, alpha...'
icabijaa=dblalloc((nval-1)*nval/2*(nal-1)*nal/2)
icabijbb=dblalloc((nvbe-1)*nvbe/2*(nbe-1)*nbe/2)
icabijab=dblalloc(nval*nvbe*nal*nbe)
call cabijal(nb,ncore,nal,nbe,nval,nvbe,nbfcn,nbfan,dfnbasis,
$vifaa,vifbb,vifab,vifba,cf12aa,cf12bb,cf12ab,cf12ba,cf12r12aa,
$cf12r12aa,cf12r12bb,cf12r12ab,cf12r12ba,cf12r12ba,cf12r12ba,f12ao,
$focka,fockb,dcore(icabijaa),dcore(icabijbb),dcore(icabijab),f12ij,
$f12ij,vaa,vbb,vab)
if(verbosity.ge.3) call timer
if(verbosity.ge.3) write(iout,*) 'Intermediate C, beta...'
call cabijbe(nb,ncore,nal,nbe,nval,nvbe,nbfcn,nbfan,dfnbasis,
$vifaa,vifbb,vifab,vifba,cf12aa,cf12bb,cf12ab,cf12ba,cf12r12aa,
$cf12r12aa,cf12r12bb,cf12r12ab,cf12r12ba,cf12r12ba,cf12r12ba,f12ao,
$focka,fockb,dcore(icabijaa),dcore(icabijbb),dcore(icabijab),f12ij,
$f12ij)
if(verbosity.ge.3) call timer
C Save intermediates
if(verbosity.ge.3) write(iout,*) 'Saving intermediates...'
call saveimos(nb,nal,nbe,nval,nvbe,scrfile1,dcore(icabijaa),
$dcore(icabijbb),dcore(icabijab),vaa,vbb,vab,uaia,uaib,vpia,vpib,
$caia,caib,tscalea,tscaleb,tfact,epaa,epbb,epab)
call dbldealloc(icabijaa)
endif
close(scrfile1)
call timer
C
return
end
C
************************************************************************
subroutine mp2int(nbf,nalc,nbec,ncore,nb,dfnbasis,nbfc,nbfcn,cmoa,
$cmob,epaa,epbb,epab,focka,fockb,hpja,hpjb,dcore,imem,natoms,iout,
$imem1,maxcor,nangmax,ncontrmax,nprimmax,ncartmax,cartg,nsphermax,
$nmboys,itol,nbfshmax,nang,ncontr,nprim,gexp,gcoef,coord,ctostr,cf,
$boysval,indarr,gcn,pre,nshrange,thad,thcf2,scoord,rqqij,rqqkl,
$hrec,nbset,spctostr,dfipre,icore,tcdfint,dfrqq,spre,natrange,
$dfipra,logkc,kp,gck,cprea,cpreb,tcf12int,vifaa,vifbb,vifab,cf12aa,
$cf12bb,cf12ab,cf122aa,cf122bb,cf122ab,vifba,cf12ba,cf122ba,
$scrfile2,scftype,lcc,cf12r12aa,cf12r12bb,cf12r12ab,cf12r12ba,
$verbosity,route)
************************************************************************
c Calculate and transform two- and three-center integrals
************************************************************************
implicit none
integer nbf,nalc,nbec,ncore,dfnbasis,imem,natoms,iout,imem1,i,nb
integer dblalloc,intalloc,icore(*),nbfc,nbfcn,scrfile2,p,q,r,s
integer nangmax,ncontrmax,nprimmax,ncartmax,nsphermax,nmboys,nbset
integer nbfshmax,nang(natoms,nbset),ncontr(0:nangmax,natoms,nbset)
integer nprim(0:nangmax,natoms,nbset),nangmin(natoms,nbset),iisyev
integer indarr(natoms,0:nangmax,ncontrmax,nsphermax,nbset),logkc
integer gcn(2,ncontrmax,0:nangmax,natoms,nbset),thad,kp,maxcor
integer nshrange(2,0:nangmax,natoms,nbset),ivhaia,ivhaib,ihai
integer natrange(2,natoms,nbset),if12ij,iscr,nalx,nbex,verbosity
real*8 dcore(*),tcdfint((dfnbasis+1)*dfnbasis/2),dfipra(natoms)
real*8 cmoa(nbfc,nbfcn),cmob(nbfc,nbfcn),cprea(natoms,0:nangmax)
real*8 itol,gexp(nprimmax,0:nangmax,natoms,nbset),ctostr,cf,gck
real*8 gcoef(nprimmax,ncontrmax,0:nangmax,natoms,nbset)
real*8 rqqij,rqqkl,hrec,spctostr,boysval,pre,thcf2,scoord,tcmax
real*8 coord(3,natoms),dfrqq,spre,dfipre(natoms,0:nangmax)
real*8 vifaa(dfnbasis,nbfcn,nbfcn),vifbb(dfnbasis,nbfcn,nbfcn)
real*8 vifab(dfnbasis,nbfcn,nbfcn),vifba(dfnbasis,nbfcn,nbfcn)
real*8 cpreb(natoms,0:nangmax),tcf12int(dfnbasis,dfnbasis)
real*8 cf12aa(nbfcn,dfnbasis,*),cf12bb(nbfcn,dfnbasis,*)
real*8 cf12ab(nbfcn,dfnbasis,*),cf12ba(nbfcn,dfnbasis,*)
real*8 cf122aa(nbfcn,dfnbasis,*),cf122bb(nbfcn,dfnbasis,*)
real*8 cf122ab(nbfcn,dfnbasis,*),cf122ba(nbfcn,dfnbasis,*)
real*8 epaa(*),epbb(*),epab(*)
real*8 hpja(nbfcn,nbfcn),hpjb(nbfcn,nbfcn)
real*8 focka(nbfcn,nbfcn),fockb(nbfcn,nbfcn)
real*8 cf12r12aa,cf12r12bb,cf12r12ab,cf12r12ba
character(len=4) route
character(len=5) scftype
logical cartg,lcc
integer*4 isyev
equivalence(isyev,iisyev) !For Intel
interface
subroutine boys
end
subroutine boysf12
end
subroutine boysf122
end
subroutine boysdelf122
end
subroutine boysf12r12
end
end interface
C MO coefficient prescreening
if(lcc) then
nalx=nb
nbex=nb
else
nalx=nalc
nbex=nbec
endif
cprea=0.d0
call dfcpre(cmoa(1,ncore+1),cprea,natoms,nangmax,ncontrmax,
$nsphermax,cartg,indarr,nang,ncontr,nbfc,nalx,.true.)
if(scftype.ne.'rhf ') then
cpreb=0.d0
call dfcpre(cmob(1,ncore+1),cpreb,natoms,nangmax,ncontrmax,
$nsphermax,cartg,indarr,nang,ncontr,nbfc,nbex,.true.)
endif
C Two-center Coulomb integrals
call df2int(natoms,nangmax,ncontrmax,nprimmax,nang,ncontr,nprim,
$gexp,coord,ncartmax,ctostr,cartg,nsphermax,cf,boysval,nmboys,
$dcore,imem,0,0.00001d0*itol,nshrange(1,0,1,3),iout,imem1,maxcor,
$thad,thcf2,scoord,rqqij,rqqkl,0,1,nang(1,3),ncontr(0,1,3),
$nprim(0,1,3),gexp(1,0,1,3),gcn(1,1,0,1,3),dfrqq,dfnbasis,
$gcoef(1,1,0,1,3),indarr(1,0,1,1,3),tcdfint,dfipre,i,i,.true.,
$scrfile2,tcdfint,spctostr,1,dfipra,tcmax,.false.,tcdfint,
$0.d0,boys)
iisyev=0
call dpptrf('L',dfnbasis,tcdfint,isyev)
if(isyev.ne.0) then
write(iout,*) 'Fatal error at the Cholesky decomposition!'
call mrccend(1)
endif
C Three-center Coulomb integrals + first half transformation
ihai=dblalloc(dfnbasis*max(4,nalc,nbec)*nbfc)
ivhaia=dblalloc(nbfc*dfnbasis*nalx)
if(scftype.ne.'rhf ') then
ivhaib=dblalloc(nbfc*dfnbasis*nbex)
else
ivhaib=ivhaia
endif
if(lcc) then
i=dblalloc(dfnbasis*nb*nbfc) !haixx
call dbldealloc(i)
endif
write(iout,*)
write(iout,*) 'Calculating Coulomb integrals...'
call gintf(natoms,nangmax,ncontrmax,nprimmax,nang,ncontr,nprim,
$gexp,gcoef,coord,ncartmax,ctostr,cartg,nsphermax,cf,boysval,
$nmboys,dcore,imem,indarr,nbfc,gcn,pre,itol,nbfshmax,nshrange,iout,
$imem1,maxcor,thad,thcf2,scoord,rqqij,rqqkl,hrec,tcdfint,dfipre,
$vifaa,vifab,cmoa,cmob,nbfcn,nalc,nalx,ncore,dfnbasis,spctostr,
$cprea,natrange,spre,nbset,dfipra,tcmax,logkc,kp,gck,nangmin,icore,
$scftype,dcore(ivhaia),dcore(ivhaia),dcore(imem),lcc,verbosity)
if(scftype.ne.'rhf ') then
C Construct beta list
call gintf(natoms,nangmax,ncontrmax,nprimmax,nang,ncontr,nprim,
$gexp,gcoef,coord,ncartmax,ctostr,cartg,nsphermax,cf,boysval,
$nmboys,dcore,imem,indarr,nbfc,gcn,pre,itol,nbfshmax,nshrange,iout,
$imem1,maxcor,thad,thcf2,scoord,rqqij,rqqkl,hrec,tcdfint,dfipre,
$vifbb,vifba,cmob,cmoa,nbfcn,nbec,nbex,ncore,dfnbasis,spctostr,
$cpreb,natrange,spre,nbset,dfipra,tcmax,logkc,kp,gck,nangmin,icore,
$scftype,dcore(ivhaib),dcore(ivhaib),dcore(imem),lcc,verbosity)
endif
call timer
C (DelF12)^2 integrals
call f12init(nmboys)
if(route.eq.'mp2 ') then
write(iout,*)
write(iout,*) 'Calculating (DelF12)^2 integrals...'
call fintf(natoms,nangmax,ncontrmax,nprimmax,nang,ncontr,nprim,
$gexp,gcoef,coord,ncartmax,ctostr,cartg,nsphermax,cf,boysval,
$nmboys,dcore,imem,indarr,nbfc,gcn,pre,itol,nbfshmax,nshrange,iout,
$imem1,maxcor,thad,thcf2,scoord,rqqij,rqqkl,hrec,tcdfint,dfipre,
$dcore(ivhaia),dcore(ivhaib),cmoa,cmob,nbfcn,nalc,nbec,ncore,
$nb,dfnbasis,spctostr,cprea,cpreb,natrange,spre,nbset,dfipra,tcmax,
$logkc,kp,gck,nangmin,icore,scftype,dcore(ihai),tcf12int,dfrqq,
$cf122aa,cf122bb,cf122ab,cf122ba,boysdelf122,'delf122 ',lcc,
$cf12r12aa,cf12r12bb,cf12r12ab,cf12r12ba,verbosity,route)
call timer
endif
C F12/r12 integrals
if(route.eq.'mp2 '.or.route.eq.'xsp ') then
write(iout,*)
write(iout,*) 'Calculating F12/r12 integrals...'
call fintf(natoms,nangmax,ncontrmax,nprimmax,nang,ncontr,nprim,
$gexp,gcoef,coord,ncartmax,ctostr,cartg,nsphermax,cf,boysval,
$nmboys,dcore,imem,indarr,nbfc,gcn,pre,itol,nbfshmax,nshrange,iout,
$imem1,maxcor,thad,thcf2,scoord,rqqij,rqqkl,hrec,tcdfint,dfipre,
$dcore(ivhaia),dcore(ivhaib),cmoa,cmob,nbfcn,nalc,nbec,ncore,
$nb,dfnbasis,spctostr,cprea,cpreb,natrange,spre,nbset,dfipra,tcmax,
$logkc,kp,gck,nangmin,icore,scftype,dcore(ihai),tcf12int,dfrqq,
$cf12aa,cf12bb,cf12ab,cf12ba,boysf12r12,'f12r12 ',lcc,
$cf12r12aa,cf12r12bb,cf12r12ab,cf12r12ba,verbosity,route)
C Calculate (DelF12)^2 and F12/r12 contributions to pair enegies
if(route.eq.'mp2 ')
$call vb1(cf122aa,cf122bb,cf122ab,cf122ba,cf12aa,cf12bb,cf12ab,
$cf12ba,nalc,nbec,dfnbasis,epaa,epbb,epab,scftype,vifaa,vifbb,
$vifab,vifba,ncore,nbfcn,dcore(ihai),dcore(ihai+dfnbasis),
$dcore(ihai+2*dfnbasis),dcore(ihai+3*dfnbasis))
call timer
C F12^2 integrals
write(iout,*)
write(iout,*) 'Calculating F12^2 integrals...'
call fintf(natoms,nangmax,ncontrmax,nprimmax,nang,ncontr,nprim,
$gexp,gcoef,coord,ncartmax,ctostr,cartg,nsphermax,cf,boysval,
$nmboys,dcore,imem,indarr,nbfc,gcn,pre,itol,nbfshmax,nshrange,iout,
$imem1,maxcor,thad,thcf2,scoord,rqqij,rqqkl,hrec,tcdfint,dfipre,
$dcore(ivhaia),dcore(ivhaib),cmoa,cmob,nbfcn,nalc,nbec,ncore,
$nb,dfnbasis,spctostr,cprea,cpreb,natrange,spre,nbset,dfipra,tcmax,
$logkc,kp,gck,nangmin,icore,scftype,dcore(ihai),tcf12int,dfrqq,
$cf122aa,cf122bb,cf122ab,cf122ba,boysf122,'f122 ',lcc,
$cf12r12aa,cf12r12bb,cf12r12ab,cf12r12ba,verbosity,route)
if(route.eq.'mp2 ') then
C Calculate F12^2 contribution to pair enegies
iscr =dblalloc(nbfcn*nbfcn)
if12ij=dblalloc(nbfcn*nbfcn)
if(scftype.eq.'rhf ') then
C Closed-shell
call epair_stf122(nalc,nbfcn,ncore,dfnbasis,focka,hpja,
$vifaa,cf122aa,epaa,dcore(if12ij),dcore(iscr))
else
C alpha-alpha
call epair_ssf122(nalc,nbfcn,ncore,dfnbasis,focka,hpja,
$vifaa,cf122aa,epaa,dcore(if12ij),dcore(iscr))
C beta-beta
call epair_ssf122(nbec,nbfcn,ncore,dfnbasis,fockb,hpjb,
$vifbb,cf122bb,epbb,dcore(if12ij),dcore(iscr))
C alpha-beta
call epair_osf122(nalc,nbec,nbfcn,ncore,dfnbasis,focka,
$fockb,hpja,hpjb,vifaa,vifbb,vifab,vifba,cf122aa,cf122bb,cf122ab,
$cf122ba,epab,dcore(if12ij))
endif
call dbldealloc(iscr)
endif
call timer
endif
C F12 integrals
write(iout,*)
write(iout,*) 'Calculating F12 integrals...'
call fintf(natoms,nangmax,ncontrmax,nprimmax,nang,ncontr,nprim,
$gexp,gcoef,coord,ncartmax,ctostr,cartg,nsphermax,cf,boysval,
$nmboys,dcore,imem,indarr,nbfc,gcn,pre,itol,nbfshmax,nshrange,iout,
$imem1,maxcor,thad,thcf2,scoord,rqqij,rqqkl,hrec,tcdfint,dfipre,
$dcore(ivhaia),dcore(ivhaib),cmoa,cmob,nbfcn,nalc,nbec,ncore,
$nb,dfnbasis,spctostr,cprea,cpreb,natrange,spre,nbset,dfipra,tcmax,
$logkc,kp,gck,nangmin,icore,scftype,dcore(ihai),tcf12int,dfrqq,
$cf12aa,cf12bb,cf12ab,cf12ba,boysf12,'f12 ',lcc,
$cf12r12aa,cf12r12bb,cf12r12ab,cf12r12ba,verbosity,route)
call timer
C
call dbldealloc(ihai)
C
return
end
C
************************************************************************
subroutine gintf(natoms,nangmax,ncontrmax,nprimmax,nang,ncontr,
$nprim,gexp,gcoef,coord,ncartmax,ctostr,cartg,nsphermax,cf,boysval,
$nmboys,dcore,imem,indarr,nbfc,gcn,pre,itol,nbfshmax,nshrange,iout,
$imem1,maxcor,thad,thcf2,scoord,rqqij,rqqkl,hrec,tcdfint,dfipre,
$vifaa,vifab,cmoa,cmob,nbfcn,nalc,nalx,ncore,dfnbasis,spctostr,
$cpre,natrange,spre,nbset,dfipra,tcmax,logkc,kp,gck,nangmin,icore,
$scftype,hai,haix,haixx,lcc,verbosity)
************************************************************************
c Calculate and transform three-center Coulomb integrals
************************************************************************
implicit none
integer dfnbasis,imem,natoms,iout,imem1,maxcor,i,nalc,ncore,nalx
integer icore(*),nbfc,nbfcn,nn,p,q,r,s,natrange(2,natoms,nbset)
integer nangmax,ncontrmax,nprimmax,ncartmax,nsphermax,nmboys,nbset
integer nbfshmax,nang(natoms,nbset),ncontr(0:nangmax,natoms,nbset)
integer nprim(0:nangmax,natoms,nbset),nangmin(natoms,nbset),iisyev
integer indarr(natoms,0:nangmax,ncontrmax,nsphermax,nbset),logkc
integer gcn(2,ncontrmax,0:nangmax,natoms,nbset),thad,kp,verbosity
integer nshrange(2,0:nangmax,natoms,nbset)
real*8 dcore(*),tcdfint((dfnbasis+1)*dfnbasis/2)
real*8 cmoa(nbfc,nbfcn),cmob(nbfc,nbfcn),cpre(natoms,0:nangmax)
real*8 itol,gexp(nprimmax,0:nangmax,natoms,nbset),ctostr,cf,gck
real*8 gcoef(nprimmax,ncontrmax,0:nangmax,natoms,nbset)
real*8 rqqij,rqqkl,hrec,spctostr,boysval,pre,thcf2,scoord,tcmax
real*8 coord(3,natoms),spre,dfipre(natoms,0:nangmax)
real*8 dfipra(natoms),hai(dfnbasis,nalx,nbfc)
real*8 vifaa(nbfcn,dfnbasis,nalx),vifab(nbfcn,dfnbasis,nalx)
real*8 haix(dfnbasis,nalc,nbfc),haixx(dfnbasis,nalx,nbfc)
character(len=5) scftype
logical cartg,lcc
integer*4 isyev
equivalence(isyev,iisyev) !For Intel
interface
subroutine boys
end
end interface
C Three-center Coulomb integrals + first half transformation
call df3int(natoms,nangmax,ncontrmax,nprimmax,nang,ncontr,nprim,
$gexp,gcoef,coord,ncartmax,ctostr,cartg,nsphermax,cf,boysval,
$nmboys,dcore,imem,indarr,i,nbfc,gcn,dcore(imem),pre,5,itol,
$nbfshmax,nshrange,i,iout,dcore,dcore,dcore,dcore,' ',imem1,
$maxcor,thad,thcf2,scoord,rqqij,rqqkl,hrec,0,1,nang(1,3),
$ncontr(0,1,3),nprim(0,1,3),gexp(1,0,1,3),gcn(1,1,0,1,3),dfnbasis,
$gcoef(1,1,0,1,3),indarr(1,0,1,1,3),dcore(imem),tcdfint,dfipre,
$dcore,hai,cmoa(1,ncore+1),nalx,dcore,dfnbasis*nalx*nbfc,i,i,i,i,i,
$i,i,itol,0,spctostr,1,nalx,cpre,.true.,.true.,0,1.d0,.false.,0,0,
$i,.false.,0.d0,dcore,0,nalx,i,i,i,i,i,i,dcore,i,i,i,i,7,dcore,1,
$.false.,dcore,.true.,nalx,natrange(1,1,3),i,i,i,i,i,i,i,i,i,
$dcore,i,dcore,dcore,dcore,dcore,dcore,spre,dfipra,tcmax,logkc,kp,
$gck,dcore,0,0,i,i,i,i,i,i,natrange,i,i,i,0,icore,0.d0,boys,1,
$dcore,i,i,i,.false.,dcore,i,dcore,i,dcore,dcore,i,nangmin,
$verbosity.ge.3,i,i,dcore,dcore,dcore,dcore,dcore,dcore,dcore,
$dcore,dcore,dcore,dcore,'off ',i,.false.,dcore,.false.,
$.false.,dcore,dcore,.false.,0,.false.,dcore(imem),i,nalx)
C Fitting
nn=nalx*nbfc
call dtptrs('L','N','N',dfnbasis,nn,tcdfint,hai,dfnbasis,isyev)
if(isyev.ne.0) then
write(iout,*)
$'Fatal error at the fitting of the Coulomb integrals!'
call mrccend(1)
endif
C Second half transformation
call dgemm('t','t',nbfcn,dfnbasis*nalx,nbfc,1.d0,cmoa,nbfc,hai,
$dfnbasis*nalx,0.d0,vifaa,nbfcn)
if(scftype.ne.'rhf ')
$call dgemm('t','t',nbfcn,dfnbasis*nalx,nbfc,1.d0,cmob,nbfc,hai,
$dfnbasis*nalx,0.d0,vifab,nbfcn)
C Reorder hai in the case of CC
if(lcc) then
call dcopy(dfnbasis*nalx*nbfc,hai,1,haixx,1)
do i=1,nalc
haix(:,i,:)=haixx(:,i,:)
enddo
endif
C
return
end
C
************************************************************************
subroutine fintf(natoms,nangmax,ncontrmax,nprimmax,nang,ncontr,
$nprim,gexp,gcoef,coord,ncartmax,ctostr,cartg,nsphermax,cf,boysval,
$nmboys,dcore,imem,indarr,nbfc,gcn,pre,itol,nbfshmax,nshrange,iout,
$imem1,maxcor,thad,thcf2,scoord,rqqij,rqqkl,hrec,tcdfint,dfipre,
$vhaia,vhaib,cmoa,cmob,nbfcn,nalc,nbec,ncore,nb,dfnbasis,spctostr,
$cprea,cpreb,natrange,spre,nbset,dfipra,tcmax,logkc,kp,gck,nangmin,
$icore,scftype,hai,tcf12int,dfrqq,cf12aa,cf12bb,cf12ab,cf12ba,boys,
$intype,lcc,cf12r12aa,cf12r12bb,cf12r12ab,cf12r12ba,verbosity,
$route)
************************************************************************
c Calculate and transform three-center Coulomb integrals
************************************************************************
implicit none
integer dfnbasis,imem,natoms,iout,imem1,maxcor,i,nb,verbosity
integer icore(*),nbfc,nbfcn,nalc,nbec,ncore
integer nangmax,ncontrmax,nprimmax,ncartmax,nsphermax,nmboys,nbset
integer nbfshmax,nang(natoms,nbset),ncontr(0:nangmax,natoms,nbset)
integer nprim(0:nangmax,natoms,nbset),nangmin(natoms,nbset),iisyev
integer indarr(natoms,0:nangmax,ncontrmax,nsphermax,nbset),logkc
integer gcn(2,ncontrmax,0:nangmax,natoms,nbset),thad,kp
integer nshrange(2,0:nangmax,natoms,nbset)
integer natrange(2,natoms,nbset)
real*8 dcore(*),tcdfint((dfnbasis+1)*dfnbasis/2)
real*8 cmoa(nbfc,nbfcn),cmob(nbfc,nbfcn),cprea(natoms,0:nangmax)
real*8 itol,gexp(nprimmax,0:nangmax,natoms,nbset),ctostr,cf,gck
real*8 gcoef(nprimmax,ncontrmax,0:nangmax,natoms,nbset)
real*8 rqqij,rqqkl,hrec,spctostr,boysval,pre,thcf2,scoord,tcmax
real*8 coord(3,natoms),dfrqq,spre,dfipre(natoms,0:nangmax)
real*8 dfipra(natoms),hai(dfnbasis,nalc,nbfc)
real*8 cf12aa(nbfcn,dfnbasis,*),vhaia(dfnbasis,nalc,nbfc)
real*8 cf12bb(nbfcn,dfnbasis,*),vhaib(dfnbasis,nalc,nbfc)
real*8 cf12ab(nbfcn,dfnbasis,*),cf12ba(nbfcn,dfnbasis,*)
real*8 tcf12int(dfnbasis,dfnbasis),cpreb(natoms,0:nangmax)
real*8 cf12r12aa,cf12r12bb,cf12r12ab,cf12r12ba
character(len=4) route
character(len=5) scftype
character(len=8) intype
logical cartg,lcc
integer*4 isyev
equivalence(isyev,iisyev) !For Intel
interface
subroutine boys
end
end interface
C Two-center integrals
call df2int(natoms,nangmax,ncontrmax,nprimmax,nang,ncontr,nprim,
$gexp,coord,ncartmax,ctostr,cartg,nsphermax,cf,boysval,nmboys,
$dcore,imem,0,0.00001d0*itol,nshrange(1,0,1,3),iout,imem1,maxcor,
$thad,thcf2,scoord,rqqij,rqqkl,0,1,nang(1,3),ncontr(0,1,3),
$nprim(0,1,3),gexp(1,0,1,3),gcn(1,1,0,1,3),dfrqq,dfnbasis,
$gcoef(1,1,0,1,3),indarr(1,0,1,1,3),tcf12int,dfipre,i,i,.true.,
$i,tcf12int,spctostr,1,dfipra,tcmax,.true.,tcf12int,1.d0,boys)
call dtptrs('L','N','N',dfnbasis,dfnbasis,tcdfint,tcf12int,
$dfnbasis,isyev)
if(isyev.ne.0) then
write(iout,*) 'Fatal error at the fitting of the F12 integrals!'
call mrccend(1)
endif
C Construct aa fitting coefficients
call fintfa(natoms,nangmax,ncontrmax,nprimmax,nang,ncontr,
$nprim,gexp,gcoef,coord,ncartmax,ctostr,cartg,nsphermax,cf,boysval,
$nmboys,dcore,imem,indarr,nbfc,gcn,pre,itol,nbfshmax,nshrange,iout,
$imem1,maxcor,thad,thcf2,scoord,rqqij,rqqkl,hrec,tcdfint,dfipre,
$nbfcn,nalc,nbec,ncore,nb,dfnbasis,spctostr,natrange,spre,nbset,
$dfipra,tcmax,logkc,kp,gck,nangmin,icore,hai,boys,tcf12int,vhaia,
$cf12aa,cf12ab,cmoa,cmob,cprea,intype,scftype,cf12r12aa,cf12r12ab,
$lcc,verbosity,route)
if(scftype.eq.'rhf '.or.nbec.eq.0) return
C Construct bb fitting coefficients
call fintfa(natoms,nangmax,ncontrmax,nprimmax,nang,ncontr,
$nprim,gexp,gcoef,coord,ncartmax,ctostr,cartg,nsphermax,cf,boysval,
$nmboys,dcore,imem,indarr,nbfc,gcn,pre,itol,nbfshmax,nshrange,iout,
$imem1,maxcor,thad,thcf2,scoord,rqqij,rqqkl,hrec,tcdfint,dfipre,
$nbfcn,nbec,nalc,ncore,nb,dfnbasis,spctostr,natrange,spre,nbset,
$dfipra,tcmax,logkc,kp,gck,nangmin,icore,hai,boys,tcf12int,vhaib,
$cf12bb,cf12ba,cmob,cmoa,cpreb,intype,scftype,cf12r12bb,cf12r12ba,
$lcc,verbosity,route)
C
return
end
C
************************************************************************
subroutine fintfa(natoms,nangmax,ncontrmax,nprimmax,nang,ncontr,
$nprim,gexp,gcoef,coord,ncartmax,ctostr,cartg,nsphermax,cf,boysval,
$nmboys,dcore,imem,indarr,nbfc,gcn,pre,itol,nbfshmax,nshrange,iout,
$imem1,maxcor,thad,thcf2,scoord,rqqij,rqqkl,hrec,tcdfint,dfipre,
$nbfcn,nalc,nbec,ncore,nb,dfnbasis,spctostr,natrange,spre,nbset,
$dfipra,tcmax,logkc,kp,gck,nangmin,icore,hai,boys,tcf12int,vhaia,
$cf12aa,cf12ab,cmoa,cmob,cpre,intype,scftype,cf12r12aa,cf12r12ab,
$lcc,verbosity,route)
************************************************************************
c Calculate and transform three-center Coulomb integrals
************************************************************************
implicit none
integer dfnbasis,imem,natoms,iout,imem1,maxcor,i,nb,verbosity
integer icore(*),nbfc,nbfcn,nn,nalc,nbec,ncore
integer nangmax,ncontrmax,nprimmax,ncartmax,nsphermax,nmboys,nbset
integer nbfshmax,nang(natoms,nbset),ncontr(0:nangmax,natoms,nbset)
integer nprim(0:nangmax,natoms,nbset),nangmin(natoms,nbset),iisyev
integer indarr(natoms,0:nangmax,ncontrmax,nsphermax,nbset),logkc
integer gcn(2,ncontrmax,0:nangmax,natoms,nbset),thad,kp
integer nshrange(2,0:nangmax,natoms,nbset)
integer natrange(2,natoms,nbset)
real*8 dcore(*),tcdfint((dfnbasis+1)*dfnbasis/2)
real*8 cmoa(nbfc,nbfcn),cmob(nbfc,nbfcn),cpre(natoms,0:nangmax)
real*8 itol,gexp(nprimmax,0:nangmax,natoms,nbset),ctostr,cf,gck
real*8 gcoef(nprimmax,ncontrmax,0:nangmax,natoms,nbset)
real*8 rqqij,rqqkl,hrec,spctostr,boysval,pre,thcf2,scoord,tcmax
real*8 coord(3,natoms),spre,dfipre(natoms,0:nangmax)
real*8 dfipra(natoms),hai(dfnbasis,nalc,nbfc)
real*8 vhaia(dfnbasis,nalc,nbfc),cf12aa(nbfcn,dfnbasis,nalc)
real*8 tcf12int(dfnbasis,dfnbasis),cf12ab(nbfcn,dfnbasis,nalc)
real*8 cf12r12aa,cf12r12ab
logical cartg,lcc
character(len=4) route
character(len=5) scftype
character(len=8) intype
integer*4 isyev
equivalence(isyev,iisyev) !For Intel
interface
subroutine boys
end
end interface
C Three-center integrals + first half transformation
call df3int(natoms,nangmax,ncontrmax,nprimmax,nang,ncontr,nprim,
$gexp,gcoef,coord,ncartmax,ctostr,cartg,nsphermax,cf,boysval,
$nmboys,dcore,imem,indarr,i,nbfc,gcn,dcore(imem),pre,5,itol,
$nbfshmax,nshrange,i,iout,dcore,dcore,dcore,dcore,' ',imem1,
$maxcor,thad,thcf2,scoord,rqqij,rqqkl,hrec,0,1,nang(1,3),
$ncontr(0,1,3),nprim(0,1,3),gexp(1,0,1,3),gcn(1,1,0,1,3),dfnbasis,
$gcoef(1,1,0,1,3),indarr(1,0,1,1,3),dcore(imem),tcdfint,dfipre,
$dcore,hai,cmoa(1,ncore+1),nalc,dcore,dfnbasis*nalc*nbfc,i,i,i,i,i,
$i,i,itol,0,spctostr,1,nalc,cpre,.true.,.true.,0,1.d0,.false.,0,0,
$i,.false.,0.d0,dcore,0,nalc,i,i,i,i,i,i,dcore,i,i,i,i,7,dcore,1,
$.false.,dcore,.true.,nalc,natrange(1,1,3),i,i,i,i,i,i,i,i,i,
$dcore,i,dcore,dcore,dcore,dcore,dcore,spre,dfipra,tcmax,logkc,kp,
$gck,dcore,0,0,i,i,i,i,i,i,natrange,i,i,i,0,icore,1.d0,boys,1,
$dcore,i,i,i,.false.,dcore,i,dcore,i,dcore,dcore,i,nangmin,
$verbosity.ge.3,i,i,dcore,dcore,dcore,dcore,dcore,dcore,dcore,
$dcore,dcore,dcore,dcore,'off ',i,.false.,dcore,.false.,
$.false.,dcore,dcore,.false.,0,.false.,dcore(imem),i,nalc)
C Fitting
nn=nalc*nbfc
call dgemm('t','n',dfnbasis,nn,dfnbasis,-0.5d0,tcf12int,dfnbasis,
$vhaia,dfnbasis,1.d0,hai,dfnbasis)
call dtptrs('L','N','N',dfnbasis,nn,tcdfint,hai,dfnbasis,isyev)
if(isyev.ne.0) then
write(iout,*) 'Fatal error at the fitting of the F12 integrals!'
call mrccend(1)
endif
C Second half transformation
if(intype.eq.'f12 '.or.intype.eq.'f122 '.or.
$route.eq.'xsp ') then
if(intype.eq.'f12r12 ') then
call dgemm('t','t',nbfcn,dfnbasis*nalc,nbfc,1.d0,cmoa,nbfc,
$hai,dfnbasis*nalc,0.d0,cf12r12aa,nbfcn)
else
call dgemm('t','t',nbfcn,dfnbasis*nalc,nbfc,1.d0,cmoa,nbfc,
$hai,dfnbasis*nalc,0.d0,cf12aa,nbfcn)
endif
else
call dgemm('n','n',dfnbasis*nalc,nalc,nbfc,1.d0,hai,
$dfnbasis*nalc,cmoa(1,ncore+1),nbfc,0.d0,cf12aa,dfnbasis*nalc)
if(intype.eq.'f12r12 '.and.lcc)
$ call dgemm('t','t',nb,dfnbasis*nalc,nbfc,1.d0,cmoa(1,ncore+1),
$nbfc,hai,dfnbasis*nalc,0.d0,cf12r12aa,nb)
endif
if(scftype.eq.'rhf ') return
if(intype.eq.'f12 '.or.intype.eq.'f122 '.or.
$route.eq.'xsp ') then
if(intype.eq.'f12r12 ') then
call dgemm('t','t',nbfcn,dfnbasis*nalc,nbfc,1.d0,cmob,nbfc,
$hai,dfnbasis*nalc,0.d0,cf12r12ab,nbfcn)
else
call dgemm('t','t',nbfcn,dfnbasis*nalc,nbfc,1.d0,cmob,nbfc,
$hai,dfnbasis*nalc,0.d0,cf12ab,nbfcn)
endif
else
call dgemm('n','n',dfnbasis*nalc,nbec,nbfc,1.d0,hai,
$dfnbasis*nalc,cmob(1,ncore+1),nbfc,0.d0,cf12ab,dfnbasis*nalc)
if(intype.eq.'f12r12 '.and.lcc)
$ call dgemm('t','t',nb,dfnbasis*nalc,nbfc,1.d0,cmob(1,ncore+1),
$nbfc,hai,dfnbasis*nalc,0.d0,cf12r12ab,nb)
endif
C
return
end
C
************************************************************************
subroutine vb1(cdelf122aa,cdelf122bb,cdelf122ab,cdelf122ba,
$cf12r12aa,cf12r12bb,cf12r12ab,cf12r12ba,nalc,nbec,dfnbasis,epaa,
$epbb,epab,scftype,vifaa,vifbb,vifab,vifba,ncore,nbfcn,vifii,vifjj,
$vifij,vifji)
************************************************************************
c Calculate intermediate B1 + intermediate V, F12/r12 contribution
************************************************************************
implicit none
integer nalc,nbec,dfnbasis,i,j,ij,nbfc,nbfcn,ncore
real*8 cdelf122aa(dfnbasis,nalc,nalc),epaa(*),epbb(*),ddot
real*8 cdelf122bb(dfnbasis,nbec,nbec),vifaa(nbfcn,dfnbasis,nalc)
real*8 cdelf122ab(dfnbasis,nalc,nbec),vifbb(nbfcn,dfnbasis,nbec)
real*8 cf12r12aa(dfnbasis,nalc,nalc),cf12r12bb(dfnbasis,nbec,nbec)
real*8 cf12r12ab(dfnbasis,nalc,nbec),vifii(dfnbasis)
real*8 vifab(nbfcn,dfnbasis,nalc),vifjj(dfnbasis),vifij(dfnbasis)
real*8 vifba(nbfcn,dfnbasis,nbec),cdelf122ba(dfnbasis,nbec,nalc)
real*8 cf12r12ba(dfnbasis,nbec,nalc),epab(nalc,nbec)
real*8 vifji(dfnbasis)
character(len=5) scftype
C
if(scftype.eq.'rhf ') then
ij=0
do i=1,nalc
vifii=vifaa(ncore+i,1:dfnbasis,i)
do j=1,i
ij=ij+1
vifjj=vifaa(ncore+j,1:dfnbasis,j)
vifij=vifaa(ncore+i,1:dfnbasis,j)
epaa(ij)=
C Intermediate B, term 1 (= tau)
$ 7d0*(ddot(dfnbasis,vifii ,1,cdelf122aa(1,j,j),1)+
$ ddot(dfnbasis,cdelf122aa(1,i,i),1,vifjj ,1))
$ +ddot(dfnbasis,vifij ,1,cdelf122aa(1,i,j),1)
$ +ddot(dfnbasis,cdelf122aa(1,i,j),1,vifij ,1)
C Intermediate V, F12/r12 contribution
$+8d0*(5d0*(ddot(dfnbasis,vifii ,1,cf12r12aa(1,j,j),1)+
$ ddot(dfnbasis,cf12r12aa(1,i,i),1,vifjj ,1))
$ -ddot(dfnbasis,vifij ,1,cf12r12aa(1,i,j),1)
$ -ddot(dfnbasis,cf12r12aa(1,i,j),1,vifij ,1))
enddo
enddo
else
C alpha-alpha
ij=0
do i=1,nalc
vifii=vifaa(ncore+i,1:dfnbasis,i)
do j=1,i-1
vifjj=vifaa(ncore+j,1:dfnbasis,j)
vifij=vifaa(ncore+i,1:dfnbasis,j)
ij=ij+1
epaa(ij)=
C Intermediate B, term 1 (= tau)
$ ddot(dfnbasis,vifii ,1,cdelf122aa(1,j,j),1)
$ +ddot(dfnbasis,cdelf122aa(1,i,i),1,vifjj ,1)
$ -ddot(dfnbasis,vifij ,1,cdelf122aa(1,i,j),1)
$ -ddot(dfnbasis,cdelf122aa(1,i,j),1,vifij ,1)
C Intermediate V, F12/r12 contribution
$ +8d0*(ddot(dfnbasis,vifii ,1,cf12r12aa(1,j,j),1)
$ +ddot(dfnbasis,cf12r12aa(1,i,i),1,vifjj ,1)
$ -ddot(dfnbasis,vifij ,1,cf12r12aa(1,i,j),1)
$ -ddot(dfnbasis,cf12r12aa(1,i,j),1,vifij ,1))
enddo
enddo
C beta-beta
ij=0
do i=1,nbec
vifii=vifbb(ncore+i,1:dfnbasis,i)
do j=1,i-1
vifjj=vifbb(ncore+j,1:dfnbasis,j)
vifij=vifbb(ncore+i,1:dfnbasis,j)
ij=ij+1
epbb(ij)=
C Intermediate B, term 1 (= tau)
$ ddot(dfnbasis,vifii ,1,cdelf122bb(1,j,j),1)
$ +ddot(dfnbasis,cdelf122bb(1,i,i),1,vifjj ,1)
$ -ddot(dfnbasis,vifij ,1,cdelf122bb(1,i,j),1)
$ -ddot(dfnbasis,cdelf122bb(1,i,j),1,vifij ,1)
C Intermediate V, F12/r12 contribution
$ +8d0*(ddot(dfnbasis,vifii ,1,cf12r12bb(1,j,j),1)
$ +ddot(dfnbasis,cf12r12bb(1,i,i),1,vifjj ,1)
$ -ddot(dfnbasis,vifij ,1,cf12r12bb(1,i,j),1)
$ -ddot(dfnbasis,cf12r12bb(1,i,j),1,vifij ,1))
enddo
enddo
C alpha-beta
do i=1,nalc
vifii=vifaa(ncore+i,1:dfnbasis,i)
do j=1,nbec
vifjj=vifbb(ncore+j,1:dfnbasis,j)
vifij=vifba(ncore+i,1:dfnbasis,j)
vifji=vifab(ncore+j,1:dfnbasis,i)
epab(i,j)=
C Intermediate B, term 1 (= tau)
$ 5d0/32d0*(ddot(dfnbasis,vifii ,1,cdelf122bb(1,j,j),1)
$ +ddot(dfnbasis,cdelf122aa(1,i,i),1,vifjj ,1))
$+3d0/32d0*(ddot(dfnbasis,vifij ,1,cdelf122ab(1,i,j),1)
$ +ddot(dfnbasis,cdelf122ba(1,j,i),1,vifji ,1))
C Intermediate V, F12/r12 contribution
$ +3d0/4d0*(ddot(dfnbasis,vifii ,1,cf12r12bb(1,j,j),1)
$ +ddot(dfnbasis,cf12r12aa(1,i,i),1,vifjj ,1))
$ +1d0/4d0*(ddot(dfnbasis,vifij ,1,cf12r12ab(1,i,j),1)
$ +ddot(dfnbasis,cf12r12ba(1,j,i),1,vifji ,1))
enddo
enddo
endif
C
return
end
C
************************************************************************
subroutine epair_stf122(nalc,nbfcn,ncore,dfnbasis,fa,hpja,vifaa,
$cf122aa,epaa,f12ij,scr)
************************************************************************
* Calculate F12^2 contribution to MP2-F12 pair energies for RHF
************************************************************************
implicit none
integer nbfcn,nalc,i,j,ncore,ij,dfnbasis
real*8 fa(nbfcn,nbfcn),hpja(nbfcn,nbfcn),e12,epaa(*),ddot
real*8 vifaa(nbfcn,dfnbasis,nalc),cf122aa(nbfcn,dfnbasis,nalc)
real*8 f12ij(nbfcn,nbfcn),scr(nbfcn,nbfcn)
C Loop over pairs
ij=0
do i=1,nalc
do j=1,i
ij=ij+1
C Intermediate B, term 2 (= pi)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,1d0, vifaa(1,1,i),
$nbfcn,cf122aa(1,1,j),nbfcn,0d0,scr,nbfcn)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,1d0,cf122aa(1,1,i),
$nbfcn, vifaa(1,1,j),nbfcn,1d0,scr,nbfcn)
f12ij=transpose(scr)
call daxpy(nbfcn*nbfcn,7d0,scr,1,f12ij,1)
e12= ddot(nbfcn,f12ij(1,ncore+j),1 ,hpja(1,ncore+i),1)
e12=e12+ddot(nbfcn,f12ij(ncore+i,1),nbfcn,hpja(1,ncore+j),1)
C Intermediate X, F12^2 contribution
e12=e12-ddot(nalc,f12ij(ncore+1,ncore+j),1 ,
$fa(ncore+1,ncore+i),1)
e12=e12-ddot(nalc,f12ij(ncore+i,ncore+1),nbfcn,
$fa(ncore+1,ncore+j),1)
epaa(ij)=epaa(ij)+e12
enddo
enddo
C
return
end
C
************************************************************************
subroutine epair_ssf122(nalc,nbfcn,ncore,dfnbasis,fa,hpja,vifaa,
$cf122aa,epaa,f12ij,scr)
************************************************************************
* Calculate same-spin F12^2 contribution to MP2-F12 pair energies
************************************************************************
implicit none
integer nbfcn,nalc,i,j,ncore,ij,dfnbasis
real*8 fa(nbfcn,nbfcn),hpja(nbfcn,nbfcn),e12,epaa(*),ddot
real*8 vifaa(nbfcn,dfnbasis,nalc),cf122aa(nbfcn,dfnbasis,nalc)
real*8 f12ij(nbfcn,nbfcn),scr(nbfcn,nbfcn)
C Loop over pairs
ij=0
do i=1,nalc
do j=1,i-1
C Intermediate B, term 2 (= pi)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,1d0, vifaa(1,1,i),
$nbfcn,cf122aa(1,1,j),nbfcn,0d0,scr,nbfcn)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,1d0,cf122aa(1,1,i),
$nbfcn, vifaa(1,1,j),nbfcn,1d0,scr,nbfcn)
f12ij=-transpose(scr)
call daxpy(nbfcn*nbfcn,1d0,scr,1,f12ij,1)
e12= ddot(nbfcn,f12ij(1,ncore+j),1 ,hpja(1,ncore+i),1)
e12=e12+ddot(nbfcn,f12ij(ncore+i,1),nbfcn,hpja(1,ncore+j),1)
C Intermediate X, F12^2 contribution
e12=e12-ddot(nalc,f12ij(ncore+1,ncore+j),1 ,
$fa(ncore+1,ncore+i),1)
e12=e12-ddot(nalc,f12ij(ncore+i,ncore+1),nbfcn,
$fa(ncore+1,ncore+j),1)
ij=ij+1
epaa(ij)=epaa(ij)+e12
enddo
enddo
C
return
end
C
************************************************************************
subroutine epair_osf122(nalc,nbec,nbfcn,ncore,dfnbasis,fa,fb,hpja,
$hpjb,vifaa,vifbb,vifab,vifba,cf122aa,cf122bb,cf122ab,cf122ba,epab,
$f12ij)
************************************************************************
* Calculate opposite-spin F12^2 contribution to MP2-F12 pair energies
************************************************************************
implicit none
integer nbfcn,nalc,nbec,ncore,i,j,dfnbasis
real*8 fa(nbfcn,nbfcn),fb(nbfcn,nbfcn),epab(nalc,nbec)
real*8 hpja(nbfcn,nbfcn),hpjb(nbfcn,nbfcn),e12,ddot
real*8 vifaa(nbfcn,dfnbasis,nalc),cf122aa(nbfcn,dfnbasis,nalc)
real*8 vifbb(nbfcn,dfnbasis,nbec),cf122bb(nbfcn,dfnbasis,nbec)
real*8 vifab(nbfcn,dfnbasis,nalc),cf122ab(nbfcn,dfnbasis,nalc)
real*8 vifba(nbfcn,dfnbasis,nbec),cf122ba(nbfcn,dfnbasis,nbec)
real*8 f12ij(nbfcn,nbfcn)
C
do i=1,nalc
do j=1,nbec
C Intermediate B, term 2 (= pi)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,5d0/32d0,
$ vifaa(1,1,i),nbfcn,cf122bb(1,1,j),nbfcn,0d0,f12ij,nbfcn)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,5d0/32d0,
$cf122aa(1,1,i),nbfcn, vifbb(1,1,j),nbfcn,1d0,f12ij,nbfcn)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,3d0/32d0,
$ vifba(1,1,j),nbfcn,cf122ab(1,1,i),nbfcn,1d0,f12ij,nbfcn)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,3d0/32d0,
$cf122ba(1,1,j),nbfcn, vifab(1,1,i),nbfcn,1d0,f12ij,nbfcn)
e12= ddot(nbfcn,f12ij(1,ncore+j),1 ,hpja(1,ncore+i),1)
e12=e12+ddot(nbfcn,f12ij(ncore+i,1),nbfcn,hpjb(1,ncore+j),1)
C Intermediate X, F12^2 contribution
e12=e12-ddot(nalc,f12ij(ncore+1,ncore+j),1 ,
$fa(ncore+1,ncore+i),1)
e12=e12-ddot(nbec,f12ij(ncore+i,ncore+1),nbfcn,
$fb(ncore+1,ncore+j),1)
epab(i,j)=epab(i,j)+e12
enddo
enddo
C
return
end
C
************************************************************************
subroutine f12r12cs(nbfcn,nb,dfnbasis,nal,ncore,vifaa,cf12r12aa,
$vv,f12ij)
************************************************************************
* Assemble F12/r12 integrals for intermediate V, closed-shell
************************************************************************
implicit none
integer nbfcn,nb,dfnbasis,nal,ncore,i,j,ij
real*8 vifaa(nbfcn,dfnbasis,nal),cf12r12aa(nb,dfnbasis,nal)
real*8 vv(nb,nb,(nal+1)*nal/2),f12ij(nb,nb)
C
ij=0
do j=1,nal
do i=1,j
ij=ij+1
call dgemm('n','t',nb,nb,dfnbasis,1d0,vifaa(ncore+1,1,i),
$nbfcn,cf12r12aa(1,1,j) ,nb ,0d0,f12ij,nb)
call dgemm('n','t',nb,nb,dfnbasis,1d0, cf12r12aa(1,1,i),
$nb ,vifaa(ncore+1,1,j),nbfcn,1d0,f12ij,nb)
call dscal(nb*nb,1d0/8d0,f12ij,1)
vv(:,:,ij)=transpose(f12ij)
call daxpy(nb*nb,3d0,f12ij,1,vv(1,1,ij),1)
enddo
enddo
C
return
end
C
************************************************************************
subroutine f12r12os(nbfcn,nb,dfnbasis,nal,nbe,ncore,vifaa,vifbb,
$vifab,vifba,cf12r12aa,cf12r12bb,cf12r12ab,cf12r12ba,vaa,vbb,vab,
$f12ij)
************************************************************************
* Assemble F12/r12 integrals for intermediate V, open-shell
************************************************************************
implicit none
integer nbfcn,nb,dfnbasis,nal,nbe,ncore,i,j,ij,p,q,pq
real*8 vifaa(nbfcn,dfnbasis,nal),cf12r12aa(nb,dfnbasis,nal)
real*8 vifbb(nbfcn,dfnbasis,nbe),cf12r12bb(nb,dfnbasis,nbe)
real*8 vifab(nbfcn,dfnbasis,nal),cf12r12ab(nb,dfnbasis,nal)
real*8 vifba(nbfcn,dfnbasis,nbe),cf12r12ba(nb,dfnbasis,nbe)
real*8 vaa((nb-1)*nb/2,(nal-1)*nal/2),vab(nb,nb,nal,nbe)
real*8 vbb((nb-1)*nb/2,(nbe-1)*nbe/2),f12ij(nb,nb)
C
ij=0
do j=1,nal
do i=1,j-1
ij=ij+1
call dgemm('n','t',nb,nb,dfnbasis,1d0,vifaa(ncore+1,1,i),
$nbfcn,cf12r12aa(1,1,j) ,nb ,0d0,f12ij,nb)
call dgemm('n','t',nb,nb,dfnbasis,1d0, cf12r12aa(1,1,i),
$nb ,vifaa(ncore+1,1,j),nbfcn,1d0,f12ij,nb)
pq=0
do q=1,nb
do p=1,q-1
pq=pq+1
vaa(pq,ij)=(f12ij(p,q)-f12ij(q,p))/4d0
enddo
enddo
enddo
enddo
ij=0
do j=1,nbe
do i=1,j-1
ij=ij+1
call dgemm('n','t',nb,nb,dfnbasis,1d0,vifbb(ncore+1,1,i),
$nbfcn,cf12r12bb(1,1,j) ,nb ,0d0,f12ij,nb)
call dgemm('n','t',nb,nb,dfnbasis,1d0, cf12r12bb(1,1,i),
$nb ,vifbb(ncore+1,1,j),nbfcn,1d0,f12ij,nb)
pq=0
do q=1,nb
do p=1,q-1
pq=pq+1
vbb(pq,ij)=(f12ij(p,q)-f12ij(q,p))/4d0
enddo
enddo
enddo
enddo
do i=1,nal
do j=1,nbe
call dgemm('n','t',nb,nb,dfnbasis,3d0/8d0,vifaa(ncore+1,1,i),
$nbfcn,cf12r12bb(1,1,j) ,nb ,0d0,vab(1,1,i,j),nb)
call dgemm('n','t',nb,nb,dfnbasis,3d0/8d0, cf12r12aa(1,1,i),
$nb ,vifbb(ncore+1,1,j),nbfcn,1d0,vab(1,1,i,j),nb)
call dgemm('n','t',nb,nb,dfnbasis,1d0/8d0,vifba(ncore+1,1,j),
$nbfcn,cf12r12ab(1,1,i), nb ,1d0,vab(1,1,i,j),nb)
call dgemm('n','t',nb,nb,dfnbasis,1d0/8d0, cf12r12ba(1,1,j),
$nb ,vifab(ncore+1,1,i),nbfcn,1d0,vab(1,1,i,j),nb)
enddo
enddo
C
return
end
C
************************************************************************
subroutine vterm1(nbfcn,nb,nnb,dfnbasis,nal,vifaa,cp,cm,vs,
$va,vv,f12pq,mns,mna,rsdims,rsdima,slo,sup)
************************************************************************
* Calculate PPL-like contribution to intermediate V, term 1
************************************************************************
implicit none
integer nbfcn,nb,nnb,dfnbasis,nal,i,j,ij,p,q,pq,r,s,rs,rsdims,slo
integer sup,rsdima,rss,rsa,ijs,ija
real*8 vifaa(nbfcn,dfnbasis,nb),tmp
real*8 cp((nnb+1)*nnb/2*(nal+1)*nal/2),vv(nb,nb,(nal+1)*nal/2)
real*8 cm((nnb-1)*nnb/2*(nal-1)*nal/2),f12pq(nnb,nnb)
real*8 vs((nnb+1)*nnb/2,rsdims)
real*8 va((nnb-1)*nnb/2,rsdima)
real*8 mns(rsdims,(nal+1)*nal/2),mna(rsdima,(nal+1)*nal/2)
C Build (anti)symmetrized Coulomb integrals
rss=0
rsa=0
do s=slo,sup
do r=1,s-1
rss=rss+1
rsa=rsa+1
call dgemm('n','t',nnb,nnb,dfnbasis,1d0,vifaa(1,1,r),
$nbfcn,vifaa(1,1,s),nbfcn,0d0,f12pq,nnb)
pq=0
do q=1,nnb
do p=1,q
pq=pq+1
vs(pq,rss)=f12pq(p,q)+f12pq(q,p)
enddo
enddo
pq=0
do q=1,nnb
do p=1,q-1
pq=pq+1
va(pq,rsa)=f12pq(p,q)-f12pq(q,p)
enddo
enddo
enddo
rss=rss+1
call dsyrk('U','N',nnb,dfnbasis,2d0,vifaa(1,1,s),nbfcn,0d0,
$f12pq,nnb)
pq=0
do q=1,nnb
do p=1,q
pq=pq+1
vs(pq,rss)=f12pq(p,q)
enddo
enddo
enddo
C Term 2
call vppl(nb,nal,cp,cm,vs,va,vv,mns,mna,rsdims,rsdima,slo,
$sup,(nnb+1)*nnb/2,(nnb-1)*nnb/2)
C
return
end
C
************************************************************************
subroutine vterm2(nbfcn,nb,nnb,dfnbasis,nal,vifaa,cp,cm,vs,
$va,vv,f12ao,mns,mna,rsdims,rsdima,slo,sup,ncore,nbfan,aodim)
************************************************************************
* Calculate PPL-like contribution to intermediate V, term 2
************************************************************************
implicit none
integer nbfcn,nb,nnb,dfnbasis,nal,i,j,ij,p,q,pq,r,s,rs,rsdims,slo
integer sup,rsdima,rss,rsa,ijs,ija,a,o,ncore,nbfan,aodim
real*8 vifaa(nbfcn,dfnbasis,nb),tmp
real*8 cp(nbfan,ncore+nal,(nal+1)*nal/2),vv(nb,nb,(nal+1)*nal/2)
real*8 cm(nbfan,ncore+nal,(nal-1)*nal/2),f12ao(nbfan,ncore+nal)
real*8 vs(nbfan,ncore+nal,rsdims)
real*8 va(nbfan,ncore+nal,rsdima)
real*8 mns(rsdims,(nal+1)*nal/2),mna(rsdima,(nal+1)*nal/2)
C Build (anti)symmetrized Coulomb integrals
rss=0
rsa=0
do s=slo,sup
do r=1,s-1
rss=rss+1
rsa=rsa+1
call dgemm('n','t',nbfan,ncore+nal,dfnbasis,1d0,
$vifaa(nnb+1,1,r),nbfcn,vifaa(1,1,s),nbfcn,0d0,vs(1,1,rss),nbfan)
call dgemm('n','t',nbfan,ncore+nal,dfnbasis,1d0,
$vifaa(nnb+1,1,s),nbfcn,vifaa(1,1,r),nbfcn,0d0,f12ao ,nbfan)
call dcopy(aodim, vs(1,1,rss),1,va(1,1,rsa),1)
call daxpy(aodim, 1d0,f12ao,1,vs(1,1,rss),1)
call daxpy(aodim, -1d0,f12ao,1,va(1,1,rsa),1)
enddo
rss=rss+1
call dgemm('n','t',nbfan,ncore+nal,dfnbasis,2d0,
$vifaa(nnb+1,1,s),nbfcn,vifaa(1,1,s),nbfcn,0d0,vs(1,1,rss),nbfan)
enddo
C Term 2
call vppl(nb,nal,cp,cm,vs,va,vv,mns,mna,rsdims,rsdima,slo,
$sup,aodim,aodim)
C
return
end
C
************************************************************************
subroutine vppl(nb,nal,cp,cm,vs,va,vv,mns,mna,rsdims,rsdima,slo,
$sup,aodims,aodima)
************************************************************************
* Calculate PPL-like contribution to intermediate V
************************************************************************
implicit none
integer nb,nal,i,j,ij,r,s,rs,rsdims,slo,sup,rsdima,ijs,ija,aodims
integer aodima
real*8 tmp,vv(nb,nb,(nal+1)*nal/2)
real*8 cp(aodims,(nal+1)*nal/2),vs(aodims,rsdims)
real*8 cm(aodima,(nal-1)*nal/2),va(aodima,rsdima)
real*8 mns(rsdims,(nal+1)*nal/2),mna(rsdima,(nal+1)*nal/2)
C
call dgemm('t','n',rsdims,(nal+1)*nal/2,aodims,1d0,vs,aodims,cp,
$aodims,0d0,mns,rsdims)
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(ij,rs,r,s,tmp)
do ij=1,(nal+1)*nal/2
rs=0
do s=slo,sup
do r=1,s-1
rs=rs+1
tmp=mns(rs,ij)
vv(r,s,ij)=vv(r,s,ij)-tmp
vv(s,r,ij)=vv(s,r,ij)-tmp
enddo
rs=rs+1
tmp=mns(rs,ij)
vv(s,s,ij)=vv(s,s,ij)-tmp
enddo
enddo
C$OMP END PARALLEL DO
call dgemm('t','n',rsdima,(nal-1)*nal/2,aodima,1d0,va,aodima,cm,
$aodima,0d0,mna,rsdima)
C$OMP PARALLEL DO SCHEDULE(DYNAMIC) DEFAULT(SHARED)
C$OMP& PRIVATE(i,j,ijs,ija,rs,r,s,tmp)
do j=1,nal
ijs=(j-1)*j/2
ija=(j-1)*(j-2)/2
do i=1,j-1
ijs=ijs+1
ija=ija+1
rs=0
do s=slo,sup
do r=1,s-1
rs=rs+1
tmp=mna(rs,ija)
vv(r,s,ijs)=vv(r,s,ijs)-tmp
vv(s,r,ijs)=vv(s,r,ijs)+tmp
enddo
enddo
enddo
enddo
C$OMP END PARALLEL DO
C
return
end
C
************************************************************************
subroutine cabijcalc(nb,ncore,nal,nval,nbfcn,nbfan,dfnbasis,vifaa,
$cf12aa,f12caj,f12acj,f12abk,vinabk,f12ao,focka,cabij,vv)
************************************************************************
* Calculate Cabij intermediate, closed-shell
************************************************************************
implicit none
integer nb,ncore,nal,nval,nbfcn,nbfan,dfnbasis,i,j,k,ij,jk,kl,a,b
real*8 focka(nbfcn,nbfcn),f12caj(nbfan,nval,nal),f12ao(nbfan,nval)
real*8 vifaa(nbfcn,dfnbasis,nb),cf12aa(nbfcn,dfnbasis,nal)
real*8 f12acj(nbfan,nval,nal),cabij(nval,nval,(nal+1)*nal/2)
real*8 f12abk(nbfan,nval,nal),vinabk(nval,nbfan,nal)
real*8 vv(nb,nb,(nal+1)*nal/2)
C
cabij(:,:,:)=vv(nal+1:nal+nval,nal+1:nal+nval,:)
jk=0
do k=1,nal
do j=1,nal
call dgemm('n','t',nbfan,nval,dfnbasis,1d0,
$ vifaa(ncore+nb+1,1,j),nbfcn,cf12aa(ncore+nal+1,1,k),nbfcn,0d0,
$f12caj(1,1,j),nbfan)
call dgemm('n','t',nbfan,nval,dfnbasis,1d0,
$cf12aa(ncore+nb+1,1,j),nbfcn, vifaa(ncore+nal+1,1,k),nbfcn,1d0,
$f12caj(1,1,j),nbfan)
enddo
call dscal(nbfan*nval*nal,1d0/8d0,f12caj,1)
do j=1,nal
call dgemm('n','t',nbfan,nval,dfnbasis,1d0,
$ vifaa(ncore+nb+1,1,k),nbfcn,cf12aa(ncore+nal+1,1,j),nbfcn,0d0,
$f12acj(1,1,j),nbfan)
call dgemm('n','t',nbfan,nval,dfnbasis,1d0,
$cf12aa(ncore+nb+1,1,k),nbfcn, vifaa(ncore+nal+1,1,j),nbfcn,1d0,
$f12acj(1,1,j),nbfan)
enddo
call dscal(nbfan*nval*nal,1d0/8d0,f12acj,1)
call dcopy(nbfan*nval*nal,f12acj,1,f12abk,1)
call daxpy(nbfan*nval*nal,3d0,f12caj,1,f12abk,1)
do a=1,nval
call dgemm('n','t',nbfan,nal,dfnbasis,1d0,
$vifaa(ncore+nb+1,1,nal+a),nbfcn,vifaa(ncore+1,1,k),nbfcn,0d0,
$f12ao,nbfan)
call dcopy(nbfan*nal,f12ao,1,vinabk(a,1,1),nval)
enddo
ij=0
kl=jk
do j=1,nal
if(j.le.k) then
jk=jk+1
call dgemm('n','n',nval,nval,nbfan,1d0,
$focka(ncore+nal+1,ncore+nb+1),nbfcn,f12abk(1,1,j),nbfan,1d0,
$cabij(1,1,jk),nval)
endif
do i=1,j
ij=ij+1
call dgemm('n','n',nval,nval,nbfan,-1d0,
$vinabk(1,1,j),nval,f12abk(1,1,i),nbfan,1d0,cabij(1,1,ij),nval)
call dgemm('t','t',nval,nval,nbfan,-1d0,
$f12abk(1,1,j),nbfan,vinabk(1,1,i),nval,1d0,cabij(1,1,ij),nval)
enddo
enddo
call dcopy(nbfan*nval*nal,f12caj,1,f12abk,1)
call daxpy(nbfan*nval*nal,3d0,f12acj,1,f12abk,1)
ij=0
jk=kl
do j=1,nal
if(j.le.k) then
jk=jk+1
call dgemm('t','n',nval,nval,nbfan,1d0,
$f12abk(1,1,j),nbfan,focka(ncore+nb+1,ncore+nal+1),nbfcn,1d0,
$cabij(1,1,jk),nval)
endif
do i=1,j
ij=ij+1
call dgemm('n','n',nval,nval,nbfan,-1d0,
$vinabk(1,1,i),nval,f12abk(1,1,j),nbfan,1d0,cabij(1,1,ij),nval)
call dgemm('t','t',nval,nval,nbfan,-1d0,
$f12abk(1,1,i),nbfan,vinabk(1,1,j),nval,1d0,cabij(1,1,ij),nval)
enddo
enddo
call dcopy(nbfan*nval*nal,f12caj,1,f12abk,1)
f12abk=-f12abk
call daxpy(nbfan*nval*nal,5d0,f12acj,1,f12abk,1)
ij=0
do j=1,nal
call dgemm('n','t',nval,nbfan,dfnbasis,1d0,
$vifaa(ncore+nal+1,1,j),nbfcn,vifaa(ncore+nb+1,1,k),nbfcn,0d0,
$vinabk(1,1,j),nval)
do i=1,j
ij=ij+1
call dgemm('n','n',nval,nval,nbfan,1d0,
$vinabk(1,1,i),nval,f12abk(1,1,j),nbfan,1d0,cabij(1,1,ij),nval)
call dgemm('t','t',nval,nval,nbfan,1d0,
$f12abk(1,1,i),nbfan,vinabk(1,1,j),nval,1d0,cabij(1,1,ij),nval)
enddo
enddo
enddo
C
return
end
C
************************************************************************
subroutine vaacalc(nb,nnb,ncore,nal,nbfcn,nbfan,dfnbasis,dcore,
$imem,imem1,maxcor,vifaa,cf12aa,f12pq,vaa)
************************************************************************
* Calculate V intermediate, same-spin block
************************************************************************
implicit none
integer nb,ncore,nal,nbfcn,nbfan,dfnbasis,imem,imem1,maxcor,s,r,rs
integer i,j,ij,p,q,pq,aodim,ijdim,rsdim,icm,iva,dblalloc,slo,sup
integer maxmem,nnb
real*8 vifaa(nbfcn,dfnbasis,nb),cf12aa(nbfcn,dfnbasis,nal)
real*8 dcore(*),f12pq(nnb,nnb),vaa((nb-1)*nb/2,(nal-1)*nal/2)
C Calculate PPL-like contribution to intermediate V, term 1
C Build antisymmetrized F12 integrals
aodim=(nnb-1)*nnb/2
ijdim=(nal-1)*nal/2
icm=dblalloc(aodim*ijdim)
pq=icm-1
do j=1,nal
do i=1,j-1
call dgemm('n','t',nnb,nnb,dfnbasis,1d0, vifaa(1,1,i),
$nbfcn,cf12aa(1,1,j),nbfcn,0d0,f12pq,nnb)
call dgemm('n','t',nnb,nnb,dfnbasis,1d0,cf12aa(1,1,i),
$nbfcn, vifaa(1,1,j),nbfcn,1d0,f12pq,nnb)
do q=1,nnb
do p=1,q-1
pq=pq+1
dcore(pq)=(f12pq(p,q)-f12pq(q,p))/4d0
enddo
enddo
enddo
enddo
C Loop over blocks
maxmem=maxcor-(imem-imem1)
slo=1
rs=1
do sup=2,nb
rsdim=(sup-slo+1)*(slo+sup-2)/2
if(sup.eq.nb.or.(rsdim+sup)*aodim.gt.maxmem) then
iva=dblalloc(aodim*rsdim)
pq=iva-1
do s=slo,sup
do r=1,s-1
call dgemm('n','t',nnb,nnb,dfnbasis,1d0,vifaa(1,1,r),
$nbfcn,vifaa(1,1,s),nbfcn,0d0,f12pq,nnb)
do q=1,nnb
do p=1,q-1
pq=pq+1
dcore(pq)=f12pq(p,q)-f12pq(q,p)
enddo
enddo
enddo
enddo
call dgemm('t','n',rsdim,ijdim,aodim,-1d0,dcore(iva),
$aodim,dcore(icm),aodim,1d0,vaa(rs,1),(nb-1)*nb/2)
call dbldealloc(iva)
slo=sup+1
rs=rs+rsdim
endif
enddo
call dbldealloc(icm)
C Calculate PPL-like contribution to intermediate V, term 2
C Build antisymmetrized F12 integrals
aodim=nbfan*(ncore+nal)
icm=dblalloc(aodim*ijdim)
ij=icm
do j=1,nal
do i=1,j-1
call dgemm('n','t',nbfan,ncore+nal,dfnbasis, 1d0,
$ vifaa(nnb+1,1,i),nbfcn,cf12aa(1,1,j),nbfcn,0d0,dcore(ij),nbfan)
call dgemm('n','t',nbfan,ncore+nal,dfnbasis, 1d0,
$cf12aa(nnb+1,1,i),nbfcn, vifaa(1,1,j),nbfcn,1d0,dcore(ij),nbfan)
call dgemm('n','t',nbfan,ncore+nal,dfnbasis,-1d0,
$ vifaa(nnb+1,1,j),nbfcn,cf12aa(1,1,i),nbfcn,1d0,dcore(ij),nbfan)
call dgemm('n','t',nbfan,ncore+nal,dfnbasis,-1d0,
$cf12aa(nnb+1,1,j),nbfcn, vifaa(1,1,i),nbfcn,1d0,dcore(ij),nbfan)
call dscal(aodim,0.25d0,dcore(ij),1)
ij=ij+aodim
enddo
enddo
C Loop over blocks
maxmem=maxcor-(imem-imem1)
slo=1
rs=1
do sup=2,nb
rsdim=(sup-slo+1)*(slo+sup-2)/2
if(sup.eq.nb.or.(rsdim+sup)*aodim.gt.maxmem) then
iva=dblalloc(aodim*rsdim)
pq=iva
do s=slo,sup
do r=1,s-1
call dgemm('n','t',nbfan,ncore+nal,dfnbasis, 1d0,
$vifaa(nnb+1,1,r),nbfcn,vifaa(1,1,s),nbfcn,0d0,dcore(pq),nbfan)
call dgemm('n','t',nbfan,ncore+nal,dfnbasis,-1d0,
$vifaa(nnb+1,1,s),nbfcn,vifaa(1,1,r),nbfcn,1d0,dcore(pq),nbfan)
pq=pq+aodim
enddo
enddo
call dgemm('t','n',rsdim,ijdim,aodim,-1d0,dcore(iva),
$aodim,dcore(icm),aodim,1d0,vaa(rs,1),(nb-1)*nb/2)
call dbldealloc(iva)
slo=sup+1
rs=rs+rsdim
endif
enddo
call dbldealloc(icm)
C
return
end
C
************************************************************************
subroutine vabcalc(nb,nnb,ncore,nal,nbe,nbfcn,nbfan,dfnbasis,
$dcore,imem,imem1,maxcor,vifaa,vifbb,vifab,vifba,cf12aa,cf12bb,
$cf12ab,cf12ba,vab)
************************************************************************
* Calculate V intermediate, opposite-spin block
************************************************************************
implicit none
integer nb,ncore,nal,nbfcn,nbfan,dfnbasis,imem,imem1,maxcor,s,r,rs
integer i,j,ij,p,q,pq,aodim,ijdim,rsdim,icm,iva,dblalloc,slo,sup
integer maxmem,nnb,nbe
real*8 vifaa(nbfcn,dfnbasis,nb),cf12aa(nbfcn,dfnbasis,nal)
real*8 vifbb(nbfcn,dfnbasis,nb),cf12bb(nbfcn,dfnbasis,nbe)
real*8 vifab(nbfcn,dfnbasis,nb),cf12ab(nbfcn,dfnbasis,nal)
real*8 vifba(nbfcn,dfnbasis,nb),cf12ba(nbfcn,dfnbasis,nbe)
real*8 dcore(*),vab(nb,nb,nal,nbe)
C Calculate PPL-like contribution to intermediate V, term 1
C Build F12 integrals
aodim=nnb*nnb
ijdim=nal*nbe
icm=dblalloc(aodim*ijdim)
pq=icm
do j=1,nbe
do i=1,nal
call dgemm('n','t',nnb,nnb,dfnbasis,3d0/8d0, vifaa(1,1,i),
$nbfcn,cf12bb(1,1,j),nbfcn,0d0,dcore(pq),nnb)
call dgemm('n','t',nnb,nnb,dfnbasis,3d0/8d0,cf12aa(1,1,i),
$nbfcn, vifbb(1,1,j),nbfcn,1d0,dcore(pq),nnb)
call dgemm('n','t',nnb,nnb,dfnbasis,1d0/8d0, vifba(1,1,j),
$nbfcn,cf12ab(1,1,i),nbfcn,1d0,dcore(pq),nnb)
call dgemm('n','t',nnb,nnb,dfnbasis,1d0/8d0,cf12ba(1,1,j),
$nbfcn, vifab(1,1,i),nbfcn,1d0,dcore(pq),nnb)
pq=pq+aodim
enddo
enddo
C Loop over blocks
maxmem=maxcor-(imem-imem1)
slo=1
rs=1
do sup=1,nb
rsdim=(sup-slo+1)*nb
if(sup.eq.nb.or.(rsdim+nb)*aodim.gt.maxmem) then
iva=dblalloc(aodim*rsdim)
pq=iva
do s=slo,sup
do r=1,nb
call dgemm('n','t',nnb,nnb,dfnbasis,1d0,vifaa(1,1,r),
$nbfcn,vifbb(1,1,s),nbfcn,0d0,dcore(pq),nnb)
pq=pq+aodim
enddo
enddo
call dgemm('t','n',rsdim,ijdim,aodim,-1d0,dcore(iva),
$aodim,dcore(icm),aodim,1d0,vab(1,slo,1,1),nb*nb)
call dbldealloc(iva)
slo=sup+1
rs=rs+rsdim
endif
enddo
call dbldealloc(icm)
C Calculate PPL-like contribution to intermediate V, term 2
C Build F12 integrals
aodim=nbfan*(ncore+nbe)
icm=dblalloc(aodim*ijdim)
pq=icm
do j=1,nbe
do i=1,nal
call dgemm('n','t',nbfan,ncore+nbe,dfnbasis,3d0/8d0,
$ vifaa(nnb+1,1,i),nbfcn,cf12bb(1,1,j),nbfcn,0d0,dcore(pq),nbfan)
call dgemm('n','t',nbfan,ncore+nbe,dfnbasis,3d0/8d0,
$cf12aa(nnb+1,1,i),nbfcn, vifbb(1,1,j),nbfcn,1d0,dcore(pq),nbfan)
call dgemm('n','t',nbfan,ncore+nbe,dfnbasis,1d0/8d0,
$ vifba(nnb+1,1,j),nbfcn,cf12ab(1,1,i),nbfcn,1d0,dcore(pq),nbfan)
call dgemm('n','t',nbfan,ncore+nbe,dfnbasis,1d0/8d0,
$cf12ba(nnb+1,1,j),nbfcn, vifab(1,1,i),nbfcn,1d0,dcore(pq),nbfan)
pq=pq+aodim
enddo
enddo
C Loop over blocks
maxmem=maxcor-(imem-imem1)
slo=1
rs=1
do sup=1,nb
rsdim=(sup-slo+1)*nb
if(sup.eq.nb.or.(rsdim+nb)*aodim.gt.maxmem) then
iva=dblalloc(aodim*rsdim)
pq=iva
do s=slo,sup
do r=1,nb
call dgemm('n','t',nbfan,ncore+nbe,dfnbasis,1d0,
$vifaa(nnb+1,1,r),nbfcn,vifbb(1,1,s),nbfcn,0d0,dcore(pq),nbfan)
pq=pq+aodim
enddo
enddo
call dgemm('t','n',rsdim,ijdim,aodim,-1d0,dcore(iva),
$aodim,dcore(icm),aodim,1d0,vab(1,slo,1,1),nb*nb)
call dbldealloc(iva)
slo=sup+1
rs=rs+rsdim
endif
enddo
call dbldealloc(icm)
C Calculate PPL-like contribution to intermediate V, term 3
C Build F12 integrals
aodim=nbfan*(ncore+nal)
icm=dblalloc(aodim*ijdim)
pq=icm
do j=1,nbe
do i=1,nal
call dgemm('n','t',ncore+nal,nbfan,dfnbasis,3d0/8d0,
$ vifaa(1,1,i),nbfcn,cf12bb(nnb+1,1,j),nbfcn,0d0,dcore(pq),
$ncore+nal)
call dgemm('n','t',ncore+nal,nbfan,dfnbasis,3d0/8d0,
$cf12aa(1,1,i),nbfcn, vifbb(nnb+1,1,j),nbfcn,1d0,dcore(pq),
$ncore+nal)
call dgemm('n','t',ncore+nal,nbfan,dfnbasis,1d0/8d0,
$ vifba(1,1,j),nbfcn,cf12ab(nnb+1,1,i),nbfcn,1d0,dcore(pq),
$ncore+nal)
call dgemm('n','t',ncore+nal,nbfan,dfnbasis,1d0/8d0,
$cf12ba(1,1,j),nbfcn, vifab(nnb+1,1,i),nbfcn,1d0,dcore(pq),
$ncore+nal)
pq=pq+aodim
enddo
enddo
C Loop over blocks
maxmem=maxcor-(imem-imem1)
slo=1
rs=1
do sup=1,nb
rsdim=(sup-slo+1)*nb
if(sup.eq.nb.or.(rsdim+nb)*aodim.gt.maxmem) then
iva=dblalloc(aodim*rsdim)
pq=iva
do s=slo,sup
do r=1,nb
call dgemm('n','t',ncore+nal,nbfan,dfnbasis,1d0,
$vifaa(1,1,r),nbfcn,vifbb(nnb+1,1,s),nbfcn,0d0,dcore(pq),ncore+nal)
pq=pq+aodim
enddo
enddo
call dgemm('t','n',rsdim,ijdim,aodim,-1d0,dcore(iva),
$aodim,dcore(icm),aodim,1d0,vab(1,slo,1,1),nb*nb)
call dbldealloc(iva)
slo=sup+1
rs=rs+rsdim
endif
enddo
call dbldealloc(icm)
C
return
end
C
************************************************************************
subroutine uaiaacalc(nb,ncore,nal,nval,nbfcn,nbfan,dfnbasis,vifaa,
$cf12aa,focka,f12ao,f12ai,vpia,uaia,caia)
************************************************************************
* Calculate Uai and Cai intermediates, same-spin block
************************************************************************
implicit none
integer nb,ncore,nal,nval,nbfcn,nbfan,dfnbasis,i,j
real*8 vpia(nb,nal),uaia(nval,nal),caia(nval,nal)
real*8 f12ao(nbfan,nval),f12ai(nbfan,nval),focka(nbfcn,nbfcn)
real*8 vifaa(nbfcn,dfnbasis,nb),cf12aa(nbfcn,dfnbasis,nal)
C
uaia(:,:)=vpia(nal+1:nal+nval,:)
caia=0d0
do j=1,nal
do i=1,j-1
call dgemm('n','t',nbfan,nval,dfnbasis,1d0,
$ vifaa(ncore+nb+1,1,i),nbfcn,cf12aa(ncore+nal+1,1,j),nbfcn,0d0,
$f12ao,nbfan)
call dgemm('n','t',nbfan,nval,dfnbasis,1d0,
$cf12aa(ncore+nb+1,1,i),nbfcn, vifaa(ncore+nal+1,1,j),nbfcn,1d0,
$f12ao,nbfan)
call dgemm('n','t',nbfan,nval,dfnbasis,-1d0,
$ vifaa(ncore+nb+1,1,j),nbfcn,cf12aa(ncore+nal+1,1,i),nbfcn,1d0,
$f12ao,nbfan)
call dgemm('n','t',nbfan,nval,dfnbasis,-1d0,
$cf12aa(ncore+nb+1,1,j),nbfcn, vifaa(ncore+nal+1,1,i),nbfcn,1d0,
$f12ao,nbfan)
call dscal(nbfan*nval,0.25d0,f12ao,1)
call dgemv('t',nbfan,nval, 1d0,f12ao,nbfan,
$focka(ncore+nb+1,ncore+i),1,1d0,caia(1,j),1)
call dgemv('t',nbfan,nval,-1d0,f12ao,nbfan,
$focka(ncore+nb+1,ncore+j),1,1d0,caia(1,i),1)
call dgemm('n','t',nbfan,nal,dfnbasis, 1d0,
$vifaa(ncore+nb+1,1,i),nbfcn,vifaa(ncore+1,1,j),nbfcn,0d0,f12ai,
$nbfan)
call dgemm('n','t',nbfan,nal,dfnbasis,-1d0,
$vifaa(ncore+nb+1,1,j),nbfcn,vifaa(ncore+1,1,i),nbfcn,1d0,f12ai,
$nbfan)
call dgemm('t','n',nval,nal,nbfan,-2d0,f12ao,nbfan,f12ai,
$nbfan,1d0,uaia,nval)
enddo
enddo
C
return
end
C
************************************************************************
subroutine cabijal(nb,ncore,nal,nbe,nval,nvbe,nbfcn,nbfan,
$dfnbasis,vifaa,vifbb,vifab,vifba,cf12aa,cf12bb,cf12ab,cf12ba,
$f12caj,f12abb,f12acj,f12abk,vinabk,vinabb,vinaab,f12ao,focka,
$fockb,cabijaa,cabijbb,cabijab,f12pq,f12rs,vaa,vbb,vab)
************************************************************************
* Calculate Cabij intermediate, open-shell, alpha loop out
************************************************************************
implicit none
integer nb,ncore,nal,nbe,nval,nvbe,nbfcn,nbfan,dfnbasis
integer i,j,k,ij,jk,kl,a,b,c,ab
real*8 focka(nbfcn,nbfcn),fockb(nbfcn,nbfcn),f12ao(nbfan,nval),tmp
real*8 vifaa(nbfcn,dfnbasis,nb),cf12aa(nbfcn,dfnbasis,nal)
real*8 vifbb(nbfcn,dfnbasis,nb),cf12bb(nbfcn,dfnbasis,nbe)
real*8 vifab(nbfcn,dfnbasis,nb),cf12ab(nbfcn,dfnbasis,nal)
real*8 vifba(nbfcn,dfnbasis,nb),cf12ba(nbfcn,dfnbasis,nbe)
real*8 f12caj(nbfan,nval,nal),f12acj(nbfan,nval,nal)
real*8 f12abk(nbfan,nval,nal),vinabk(nval,nbfan,nal)
real*8 vinabb(nvbe,nbfan,nbe),vinaab(nvbe,nbfan,nal)
real*8 cabijaa((nval-1)*nval/2,(nal-1)*nal/2),f12pq(nval,nval)
real*8 cabijbb((nvbe-1)*nvbe/2,(nbe-1)*nbe/2),f12rs(nvbe,nvbe)
real*8 cabijab(nval,nvbe,nal,nbe),f12abb(nbfan,nvbe,nbe)
real*8 vaa((nb-1)*nb/2,(nal-1)*nal/2),vab(nb,nb,nal,nbe)
real*8 vbb((nb-1)*nb/2,(nbe-1)*nbe/2)
C Initialize C
ab=0
do b=1,nval
do a=1,b-1
ab=ab+1
cabijaa(ab,:)=vaa((nal+b-1)*(nal+b-2)/2+nal+a,:)
enddo
enddo
ab=0
do b=1,nvbe
do a=1,b-1
ab=ab+1
cabijbb(ab,:)=vbb((nbe+b-1)*(nbe+b-2)/2+nbe+a,:)
enddo
enddo
cabijab(:,:,:,:)=vab(nal+1:nal+nval,nbe+1:nbe+nvbe,:,:)
C Alpha loop
jk=0
do k=1,nal
call cabijss(k,jk,nb,ncore,nal,nval,nbfcn,nbfan,dfnbasis,vifaa,
$cf12aa,f12caj,f12acj,f12abk,vinabk,f12ao,focka,cabijaa,f12pq)
do j=1,nbe
call dgemm('n','t',nbfan,nvbe,dfnbasis,3d0/8d0,
$ vifaa(ncore+nb+1,1,k),nbfcn,cf12bb(ncore+nbe+1,1,j),nbfcn,0d0,
$f12abb(1,1,j),nbfan)
call dgemm('n','t',nbfan,nvbe,dfnbasis,3d0/8d0,
$cf12aa(ncore+nb+1,1,k),nbfcn, vifbb(ncore+nbe+1,1,j),nbfcn,1d0,
$f12abb(1,1,j),nbfan)
call dgemm('n','t',nbfan,nvbe,dfnbasis,1d0/8d0,
$ vifba(ncore+nb+1,1,j),nbfcn,cf12ab(ncore+nbe+1,1,k),nbfcn,1d0,
$f12abb(1,1,j),nbfan)
call dgemm('n','t',nbfan,nvbe,dfnbasis,1d0/8d0,
$cf12ba(ncore+nb+1,1,j),nbfcn, vifab(ncore+nbe+1,1,k),nbfcn,1d0,
$f12abb(1,1,j),nbfan)
do i=1,nal
call dgemm('n','n',nval,nvbe,nbfan,1d0,
$vinabk(1,1,i),nval,f12abb(1,1,j),nbfan,1d0,cabijab(1,1,i,j),nval)
enddo
enddo
do a=1,nvbe
call dgemm('n','t',nbfan,nbe,dfnbasis, 1d0,
$vifaa(ncore+nb+1,1,k),nbfcn,vifbb(ncore+1,1,nbe+a),nbfcn,0d0,
$f12ao,nbfan)
call dcopy(nbfan*nbe,f12ao,1,vinabb(a,1,1),nvbe)
enddo
ij=0
do j=1,nbe
do i=1,j-1
ij=ij+1
call dgemm('n','n',nvbe,nvbe,nbfan, 1d0,
$vinabb(1,1,i),nvbe,f12abb(1,1,j),nbfan,0d0,f12rs,nvbe)
call dgemm('n','n',nvbe,nvbe,nbfan,-1d0,
$vinabb(1,1,j),nvbe,f12abb(1,1,i),nbfan,1d0,f12rs,nvbe)
ab=0
do b=1,nvbe
do a=1,b-1
ab=ab+1
cabijbb(ab,ij)=cabijbb(ab,ij)+f12rs(a,b)-f12rs(b,a)
enddo
enddo
enddo
enddo
do j=1,nbe
do i=1,nal
call dgemm('t','t',nval,nvbe,nbfan,1d0,
$f12abk(1,1,i),nbfan,vinabb(1,1,j),nvbe,1d0,cabijab(1,1,i,j),nval)
enddo
enddo
do a=1,nvbe
call dgemm('n','t',nbfan,nal,dfnbasis, 1d0,
$vifbb(ncore+nb+1,1,nbe+a),nbfcn,vifaa(ncore+1,1,k),nbfcn,0d0,
$f12ao,nbfan)
call dcopy(nbfan*nal,f12ao,1,vinaab(a,1,1),nvbe)
enddo
do j=1,nbe
call dgemm('n','t',nval,nbfan,dfnbasis,3d0/8d0,
$ vifaa(ncore+nal+1,1,k),nbfcn,cf12bb(ncore+nb+1,1,j),nbfcn,0d0,
$f12ao,nval)
call dgemm('n','t',nval,nbfan,dfnbasis,3d0/8d0,
$cf12aa(ncore+nal+1,1,k),nbfcn, vifbb(ncore+nb+1,1,j),nbfcn,1d0,
$f12ao,nval)
call dgemm('n','t',nval,nbfan,dfnbasis,1d0/8d0,
$ vifba(ncore+nal+1,1,j),nbfcn,cf12ab(ncore+nb+1,1,k),nbfcn,1d0,
$f12ao,nval)
call dgemm('n','t',nval,nbfan,dfnbasis,1d0/8d0,
$cf12ba(ncore+nal+1,1,j),nbfcn, vifab(ncore+nb+1,1,k),nbfcn,1d0,
$f12ao,nval)
do i=1,nal
call dgemm('n','t',nval,nvbe,nbfan,-1d0,
$f12ao,nval,vinaab(1,1,i),nvbe,1d0,cabijab(1,1,i,j),nval)
enddo
enddo
enddo
C
return
end
C
************************************************************************
subroutine cabijss(k,jk,nb,ncore,nal,nval,nbfcn,nbfan,dfnbasis,
$vifaa,cf12aa,f12caj,f12acj,f12abk,vinabk,f12ao,focka,cabijaa,
$f12pq)
************************************************************************
* Calculate Cabij intermediate, same-spin contribution
************************************************************************
implicit none
integer nb,ncore,nal,nbe,nval,nvbe,nbfcn,nbfan,dfnbasis
integer i,j,k,ij,jk,kl,a,b,c,ab
real*8 focka(nbfcn,nbfcn),f12ao(nbfan,nval),f12pq(nval,nval),tmp
real*8 vifaa(nbfcn,dfnbasis,nb),cf12aa(nbfcn,dfnbasis,nal)
real*8 f12caj(nbfan,nval,nal),f12acj(nbfan,nval,nal)
real*8 f12abk(nbfan,nval,nal),vinabk(nval,nbfan,nal)
real*8 cabijaa((nval-1)*nval/2,(nal-1)*nal/2)
C
do j=1,nal
call dgemm('n','t',nbfan,nval,dfnbasis,1d0,
$ vifaa(ncore+nb+1,1,j),nbfcn,cf12aa(ncore+nal+1,1,k),nbfcn,0d0,
$f12caj(1,1,j),nbfan)
call dgemm('n','t',nbfan,nval,dfnbasis,1d0,
$cf12aa(ncore+nb+1,1,j),nbfcn, vifaa(ncore+nal+1,1,k),nbfcn,1d0,
$f12caj(1,1,j),nbfan)
enddo
call dscal(nbfan*nval*nal,1d0/4d0,f12caj,1)
do j=1,nal
call dgemm('n','t',nbfan,nval,dfnbasis,1d0,
$ vifaa(ncore+nb+1,1,k),nbfcn,cf12aa(ncore+nal+1,1,j),nbfcn,0d0,
$f12acj(1,1,j),nbfan)
call dgemm('n','t',nbfan,nval,dfnbasis,1d0,
$cf12aa(ncore+nb+1,1,k),nbfcn, vifaa(ncore+nal+1,1,j),nbfcn,1d0,
$f12acj(1,1,j),nbfan)
enddo
call dscal(nbfan*nval*nal,1d0/4d0,f12acj,1)
call dcopy(nbfan*nval*nal,f12acj,1,f12abk,1)
f12abk=-f12abk
call daxpy(nbfan*nval*nal,1d0,f12caj,1,f12abk,1)
do j=1,k-1
jk=jk+1
call dgemm('n','n',nval,nval,nbfan, 1d0,
$focka(ncore+nal+1,ncore+nb+1),nbfcn,f12abk(1,1,j),nbfan,0d0,
$f12pq,nval)
ab=0
do b=1,nval
do a=1,b-1
ab=ab+1
cabijaa(ab,jk)=cabijaa(ab,jk)+f12pq(a,b)-f12pq(b,a)
enddo
enddo
enddo
call dcopy(nbfan*nval*nal,f12caj,1,f12abk,1)
f12abk=-f12abk
call daxpy(nbfan*nval*nal,1d0,f12acj,1,f12abk,1)
do a=1,nval
call dgemm('n','t',nbfan,nal,dfnbasis, 1d0,
$vifaa(ncore+nb+1,1,k),nbfcn,vifaa(ncore+1,1,nal+a),nbfcn,0d0,
$f12ao,nbfan)
call dgemm('n','t',nbfan,nal,dfnbasis,-1d0,
$vifaa(ncore+nb+1,1,nal+a),nbfcn,vifaa(ncore+1,1,k),nbfcn,1d0,
$f12ao,nbfan)
call dcopy(nbfan*nal,f12ao,1,vinabk(a,1,1),nval)
enddo
ij=0
do j=1,nal
do i=1,j-1
ij=ij+1
call dgemm('n','n',nval,nval,nbfan, 1d0,
$vinabk(1,1,i),nval,f12abk(1,1,j),nbfan,0d0,f12pq,nval)
call dgemm('n','n',nval,nval,nbfan,-1d0,
$vinabk(1,1,j),nval,f12abk(1,1,i),nbfan,1d0,f12pq,nval)
ab=0
do b=1,nval
do a=1,b-1
ab=ab+1
cabijaa(ab,ij)=cabijaa(ab,ij)+f12pq(a,b)-f12pq(b,a)
enddo
enddo
enddo
enddo
C
return
end
C
************************************************************************
subroutine cabijbe(nb,ncore,nal,nbe,nval,nvbe,nbfcn,nbfan,
$dfnbasis,vifaa,vifbb,vifab,vifba,cf12aa,cf12bb,cf12ab,cf12ba,
$f12caj,f12abb,f12acj,f12abk,vinabk,vinabb,vinaab,f12ao,focka,
$fockb,cabijaa,cabijbb,cabijab,f12pq,f12rs)
************************************************************************
* Calculate Cabij intermediate, open-shell, beta loop out
************************************************************************
implicit none
integer nb,ncore,nal,nbe,nval,nvbe,nbfcn,nbfan,dfnbasis
integer i,j,k,ij,jk,kl,a,b,c,ab
real*8 focka(nbfcn,nbfcn),fockb(nbfcn,nbfcn),f12ao(nbfan,nval),tmp
real*8 vifaa(nbfcn,dfnbasis,nb),cf12aa(nbfcn,dfnbasis,nal)
real*8 vifbb(nbfcn,dfnbasis,nb),cf12bb(nbfcn,dfnbasis,nbe)
real*8 vifab(nbfcn,dfnbasis,nb),cf12ab(nbfcn,dfnbasis,nal)
real*8 vifba(nbfcn,dfnbasis,nb),cf12ba(nbfcn,dfnbasis,nbe)
real*8 f12caj(nbfan,nvbe,nbe),f12acj(nbfan,nvbe,nbe)
real*8 f12abk(nbfan,nvbe,nbe),vinabk(nvbe,nbfan,nbe)
real*8 vinabb(nval,nbfan,nal),vinaab(nval,nbfan,nal)
real*8 cabijaa((nval-1)*nval/2,(nal-1)*nal/2),f12pq(nvbe,nvbe)
real*8 cabijbb((nvbe-1)*nvbe/2,(nbe-1)*nbe/2),f12rs(nval,nval)
real*8 cabijab(nval,nvbe,nal,nbe),f12abb(nval,nbfan,nal)
C Beta loop
jk=0
do k=1,nbe
call cabijss(k,jk,nb,ncore,nbe,nvbe,nbfcn,nbfan,dfnbasis,vifbb,
$cf12bb,f12caj,f12acj,f12abk,vinabk,f12ao,fockb,cabijbb,f12pq)
do j=1,nal
call dgemm('n','t',nval,nbfan,dfnbasis,3d0/8d0,
$ vifaa(ncore+nal+1,1,j),nbfcn,cf12bb(ncore+nb+1,1,k),nbfcn,0d0,
$f12abb(1,1,j),nval)
call dgemm('n','t',nval,nbfan,dfnbasis,3d0/8d0,
$cf12aa(ncore+nal+1,1,j),nbfcn, vifbb(ncore+nb+1,1,k),nbfcn,1d0,
$f12abb(1,1,j),nval)
call dgemm('n','t',nval,nbfan,dfnbasis,1d0/8d0,
$ vifba(ncore+nal+1,1,k),nbfcn,cf12ab(ncore+nb+1,1,j),nbfcn,1d0,
$f12abb(1,1,j),nval)
call dgemm('n','t',nval,nbfan,dfnbasis,1d0/8d0,
$cf12ba(ncore+nal+1,1,k),nbfcn, vifab(ncore+nb+1,1,j),nbfcn,1d0,
$f12abb(1,1,j),nval)
enddo
do j=1,nal
call dgemm('n','n',nval,nvbe,nbfan,1d0,f12abb(1,1,j),nval,
$fockb(ncore+nb+1,ncore+nbe+1),nbfcn,1d0,cabijab(1,1,j,k),nval)
do i=1,nal
call dgemm('n','t',nval,nvbe,nbfan,1d0,
$f12abb(1,1,i),nval,vinabk(1,1,j),nvbe,1d0,cabijab(1,1,i,j),nval)
enddo
enddo
do a=1,nval
call dgemm('n','t',nbfan,nal,dfnbasis, 1d0,
$vifbb(ncore+nb+1,1,k),nbfcn,vifaa(ncore+1,1,nal+a),nbfcn,0d0,
$f12ao,nbfan)
call dcopy(nbfan*nal,f12ao,1,vinabb(a,1,1),nval)
enddo
ij=0
do j=1,nal
do i=1,j-1
ij=ij+1
call dgemm('n','t',nval,nval,nbfan, 1d0,
$vinabb(1,1,i),nval,f12abb(1,1,j),nval,0d0,f12rs,nval)
call dgemm('n','t',nval,nval,nbfan,-1d0,
$vinabb(1,1,j),nval,f12abb(1,1,i),nval,1d0,f12rs,nval)
ab=0
do b=1,nval
do a=1,b-1
ab=ab+1
cabijaa(ab,ij)=cabijaa(ab,ij)+f12rs(a,b)-f12rs(b,a)
enddo
enddo
enddo
enddo
do j=1,nbe
do i=1,nal
call dgemm('n','n',nval,nvbe,nbfan,1d0,
$vinabb(1,1,i),nval,f12abk(1,1,j),nbfan,1d0,cabijab(1,1,i,j),nval)
enddo
enddo
do a=1,nval
call dgemm('n','t',nbfan,nal,dfnbasis, 1d0,
$vifaa(ncore+nb+1,1,nal+a),nbfcn,vifbb(ncore+1,1,k),nbfcn,0d0,
$f12ao,nbfan)
call dcopy(nbfan*nal,f12ao,1,vinaab(a,1,1),nval)
enddo
do i=1,nal
call dgemm('n','t',nbfan,nvbe,dfnbasis,3d0/8d0,
$ vifaa(ncore+nb+1,1,i),nbfcn,cf12bb(ncore+nbe+1,1,k),nbfcn,0d0,
$f12caj(1,1,i),nbfan)
call dgemm('n','t',nbfan,nvbe,dfnbasis,3d0/8d0,
$cf12aa(ncore+nb+1,1,i),nbfcn, vifbb(ncore+nbe+1,1,k),nbfcn,1d0,
$f12caj(1,1,i),nbfan)
call dgemm('n','t',nbfan,nvbe,dfnbasis,1d0/8d0,
$ vifba(ncore+nb+1,1,k),nbfcn,cf12ab(ncore+nbe+1,1,i),nbfcn,1d0,
$f12caj(1,1,i),nbfan)
call dgemm('n','t',nbfan,nvbe,dfnbasis,1d0/8d0,
$cf12ba(ncore+nb+1,1,k),nbfcn, vifab(ncore+nbe+1,1,i),nbfcn,1d0,
$f12caj(1,1,i),nbfan)
call dgemm('n','n',nval,nvbe,nbfan,1d0,
$focka(ncore+nal+1,ncore+nb+1),nbfcn,f12caj(1,1,i),nbfan,1d0,
$cabijab(1,1,i,k),nval)
enddo
do j=1,nbe
do i=1,nal
call dgemm('n','n',nval,nvbe,nbfan,-1d0,
$vinaab(1,1,j),nval,f12caj(1,1,i),nbfan,1d0,cabijab(1,1,i,j),nval)
enddo
enddo
enddo
C
return
end
C
************************************************************************
subroutine saveimos(nb,nal,nbe,nval,nvbe,scrfile1,cabijaa,
$cabijbb,cabijab,vaa,vbb,vab,uaia,uaib,vpia,vpib,caia,caib,tscalea,
$tscaleb,tfact,epaa,epbb,epab)
************************************************************************
* Save intermediates, open-shell
************************************************************************
implicit none
integer nb,nal,nbe,nval,nvbe,scrfile1,i,j,k,ij,kl,a,b
real*8 vaa((nb-1)*nb/2,(nal-1)*nal/2),vab(nb,nb,nal,nbe)
real*8 vbb((nb-1)*nb/2,(nbe-1)*nbe/2)
real*8 vpia(nb,nal),vpib(nb,nbe),uaia(nval,nal),uaib(nvbe,nbe)
real*8 caia(nval,nal),cabijaa((nval-1)*nval/2,(nal-1)*nal/2)
real*8 caib(nvbe,nbe),cabijbb((nvbe-1)*nvbe/2,(nbe-1)*nbe/2)
real*8 cabijab(nval,nvbe,nal,nbe),tscalea(nal),tscaleb(nbe),tfact
real*8 epaa(nal*(nal-1)/2),epbb(nbe*(nbe-1)/2),epab(nal,nbe)
C
write(scrfile1) tfact,tscalea,tscaleb,epaa,epbb,epab
if(nal.gt.1.and.nval.gt.1) write(scrfile1) cabijaa
if(nbe.gt.1.and.nvbe.gt.1) write(scrfile1) cabijbb
if(nal.gt.0.and.nval.gt.0.and.nbe.gt.0.and.nvbe.gt.0)
$write(scrfile1) cabijab
if(nal.gt.1.and.nval.gt.1) write(scrfile1)
$(((vaa((b-1)*(b-2)/2+a,ij),a=nal+1,b-1),b=nal+1,nb),
$ij=1,(nal-1)*nal/2)
if(nbe.gt.1.and.nvbe.gt.1) write(scrfile1)
$(((vbb((b-1)*(b-2)/2+a,ij),a=nbe+1,b-1),b=nbe+1,nb),
$ij=1,(nbe-1)*nbe/2)
if(nal.gt.0.and.nval.gt.0.and.nbe.gt.0.and.nvbe.gt.0)
$write(scrfile1) (((vab(nal+1:nb,nbe+b,i,j),b=1,nvbe),i=1,nal),
$j=1,nbe) ! before: vab(nal+1:nb,nbe+1:nb,:,:)
if(nal.gt.0.and.nval.gt.0) write(scrfile1) uaia
if(nbe.gt.0.and.nvbe.gt.0) write(scrfile1) uaib
close(scrfile1)
open(scrfile1,file='F12INT1',form='unformatted')
if(nal.gt.0.and.nval.gt.0) then
write(scrfile1) vpia(1:nal,:)
write(scrfile1) caia
endif
if(nbe.gt.0.and.nvbe.gt.0) then
write(scrfile1) vpib(1:nbe,:)
write(scrfile1) caib
endif
close(scrfile1)
open(scrfile1,file='F12INT2',form='unformatted')
if(nal.gt.1.and.nval.gt.0) write(scrfile1)
$(((vaa((a-1)*(a-2)/2+k,ij),k=1,nal),a=nal+1,nb),
$ij=1,(nal-1)*nal/2)
if(nal.gt.1.and.nval.gt.1) write(scrfile1)
$((vaa(kl,ij),kl=1,(nal-1)*nal/2),ij=1,(nal-1)*nal/2)
if(nbe.gt.1.and.nvbe.gt.0) write(scrfile1)
$(((vbb((a-1)*(a-2)/2+k,ij),k=1,nbe),a=nbe+1,nb),
$ij=1,(nbe-1)*nbe/2)
if(nbe.gt.1.and.nvbe.gt.1) write(scrfile1)
$((vbb(kl,ij),kl=1,(nbe-1)*nbe/2),ij=1,(nbe-1)*nbe/2)
if(nal.gt.0.and.nval.gt.0.and.nbe.gt.0) write(scrfile1)
$((((vab(nal+a,b,i,j),b=1,nbe),a=1,nval),i=1,nal),j=1,nbe)
if(nal.gt.0.and.nvbe.gt.0.and.nbe.gt.0) write(scrfile1)
$(((vab(1:nal,nbe+b,i,j),b=1,nvbe),i=1,nal),j=1,nbe) ! before vab(1:nal,nbe+1:nb,:,:)
if(nal.gt.0.and.nval.gt.0.and.nbe.gt.0.and.nvbe.gt.0)
$write(scrfile1) (((vab(1:nal,b,i,j),b=1,nbe),i=1,nal),j=1,nbe) ! before vab(1:nal,1:nbe,:,:)
C
return
end
C
************************************************************************
subroutine saveimcs(nb,nal,nval,scrfile1,cabij,vv,uaia,vpia,caia,
$tscale,tfact,epaa)
************************************************************************
* Save intermediates, closed-shell
************************************************************************
implicit none
integer nb,nal,nval,scrfile1,i,j,k,ij,a,b
real*8 vv(nb,nb,(nal+1)*nal/2),cabij(nval,nval,(nal+1)*nal/2)
real*8 vpia(nb,nal),uaia(nval,nal),caia(nval,nal),tscale(nal)
real*8 tfact,epaa(nal*(nal+1)/2)
C
write(scrfile1) tfact,tscale,epaa
write(scrfile1)
$((((2.d0*cabij(a,b,j*(j-1)/2+i)-cabij(b,a,j*(j-1)/2+i),
$a=1,nval),b=1,nval),i=1,j-1),
$((0.5d0*cabij(a,b,j*(j-1)/2+j),a=1,nval),b=1,nval),j=1,nal)
write(scrfile1)
$((((2.d0*vv(nal+a,nal+b,j*(j-1)/2+i)-vv(nal+b,nal+a,j*(j-1)/2+i),
$a=1,nval),b=1,nval),i=1,j-1),
$((0.5d0*vv(nal+a,nal+b,j*(j-1)/2+j),a=1,nval),b=1,nval),j=1,nal)
write(scrfile1) uaia
close(scrfile1)
open(scrfile1,file='F12INT1',form='unformatted')
write(scrfile1) caia,
$((((cabij(a,b,j*(j-1)/2+i),a=1,nval),b=1,nval),i=1,j),j=1,nal)
ij=0
do j=1,nal
do i=1,j
ij=ij+1
write(scrfile1) transpose(vv(1:nal,1:nal,ij))
enddo
enddo
write(scrfile1) vpia(1:nal,:)
do j=1,nal
write(scrfile1) (vv(1:nal,nal+1:nb,j*(j-1)/2+i),i=1,j),
$ (transpose(vv(nal+1:nb,1:nal,i*(i-1)/2+j)),i=j+1,nal)
enddo
C
return
end
C
************************************************************************
subroutine cnab(nb,dfnbasis,nal,nbe,ncore,focka,fockb,nbfcn,nbfan,
$scftype,vifaa,vifbb,vifab,vifba,cf12aa,cf12bb,cf12ab,cf12ba,iout,
$verbosity,wnaba,wnabb,eig,um,fockan,fockbn,scr,nabtol)
************************************************************************
* Construct the natural auxiliary basis
************************************************************************
implicit none
integer nb,dfnbasis,nal,nbe,ncore,nbfcn,nbfan,iout,verbosity
integer iisyev,nnab,nbasis
real*8 focka(nbfcn,nbfcn),fockb(nbfcn,nbfcn),um(nbfan,nbfan)
real*8 vifaa(nbfcn,dfnbasis,nb),cf12aa(nbfcn,dfnbasis,nal),scr(*)
real*8 vifbb(nbfcn,dfnbasis,nb),cf12bb(nbfcn,dfnbasis,nbe),nabtol
real*8 vifab(nbfcn,dfnbasis,nb),cf12ab(nbfcn,dfnbasis,nal)
real*8 vifba(nbfcn,dfnbasis,nb),cf12ba(nbfcn,dfnbasis,nbe)
real*8 eig(nbfan),wnaba(nbfan,nbfan),wnabb(nbfan,nbfan)
real*8 fockan(nbfcn,nbfcn),fockbn(nbfcn,nbfcn)
character(len=5) scftype
character(len=8) c8
character(len=16) c16
integer*4 isyev
equivalence(isyev,iisyev)
C
write(iout,*)
write(iout,*) 'Constructing the natural auxiliary basis...'
nbasis=ncore+nb
nnab=0
C SVD of (p'q|1/r12|Q) integrals
call svdf12(nbasis,nbfcn,nbfan,nnab,vifaa,wnaba,eig,scr,iout,
$nabtol,dfnbasis*nb)
if(scftype.ne.'rhf ')
$call svdf12(nbasis,nbfcn,nbfan,nnab,vifbb,wnabb,eig,scr,iout,
$nabtol,dfnbasis*nb)
C SVD of (p'i|f12|Q) integrals
c call svdf12(nbasis,nbfcn,nbfan,nnab,cf12aa,wnaba,eig,scr,iout,
c $nabtol,dfnbasis*nal)
c if(scftype.ne.'rhf ')
c $call svdf12(nbasis,nbfcn,nbfan,nnab,cf12bb,wnabb,eig,scr,iout,
c $nabtol,dfnbasis*nbe)
c nnab=80 !szemet
write(c8,'(i8)') nnab
write(c16,'(f16.1)') 100.d0*dble(nnab)/dble(nbfan)
write(iout,*) 'Number of natural auxiliary basis functions: ' //
$trim(adjustl(c8)) // trim(' (' // adjustl(c16)) // '%)'
C Canonicalize NAB virtuals
call cannab(nbasis,nbfcn,nbfan,nnab,focka,wnaba,um,scr,eig,iout)
if(scftype.ne.'rhf ')
$call cannab(nbasis,nbfcn,nbfan,nnab,fockb,wnabb,um,scr,eig,iout)
C Transform Fock-matrix to NAB basis
call trfnab(nbasis,nbfcn,nbfan,nnab,focka,fockan,wnaba,scr)
call trinab(nbasis,nbfcn,nbfan,nnab,dfnbasis*nb,vifaa,wnaba,scr)
call trinab(nbasis,nbfcn,nbfan,nnab,dfnbasis*nal,cf12aa,wnaba,scr)
if(scftype.ne.'rhf ') then
call trfnab(nbasis,nbfcn,nbfan,nnab,fockb,fockbn,wnabb,scr)
call trinab(nbasis,nbfcn,nbfan,nnab,dfnbasis*nb,vifbb,wnabb,scr)
call trinab(nbasis,nbfcn,nbfan,nnab,dfnbasis*nb,vifab,wnabb,scr)
call trinab(nbasis,nbfcn,nbfan,nnab,dfnbasis*nb,vifba,wnaba,scr)
call trinab(nbasis,nbfcn,nbfan,nnab,dfnbasis*nbe,cf12bb,wnabb,scr)
call trinab(nbasis,nbfcn,nbfan,nnab,dfnbasis*nal,cf12ab,wnabb,scr)
call trinab(nbasis,nbfcn,nbfan,nnab,dfnbasis*nbe,cf12ba,wnaba,scr)
endif
nbfan=nnab
nbfcn=nbasis+nnab
call timer
C
return
end
C
************************************************************************
subroutine trfnab(nbasis,nbfcn,nbfan,nnab,fock,fockn,wnab,scr)
************************************************************************
* Transform Fock-matrix to NAB basis
************************************************************************
implicit none
integer nbasis,nbfcn,nbfan,nnab
real*8 fock(nbfcn,nbfcn),fockn(nbasis+nnab,nbasis+nnab),scr(*)
real*8 wnab(nbfan,nnab)
C
fockn(1:nbasis,1:nbasis)=fock(1:nbasis,1:nbasis)
call dsymm('l','l',nbfan,nnab,1d0,fock(nbasis+1,nbasis+1),nbfcn,
$wnab,nbfan,0d0,scr,nbfan)
call dgemm('t','n',nnab,nnab,nbfan,1d0,wnab,nbfan,scr,nbfan,0d0,
$fockn(nbasis+1,nbasis+1),nbasis+nnab)
call dgemm('n','n',nbasis,nnab,nbfan,1d0,fock(1,nbasis+1),nbfcn,
$wnab,nbfan,0d0,fockn(1,nbasis+1),nbasis+nnab)
call dgemm('t','n',nnab,nbasis,nbfan,1d0,wnab,nbfan,
$fock(nbasis+1,1),nbfcn,0d0,fockn(nbasis+1,1),nbasis+nnab)
call dcopy((nbasis+nnab)*(nbasis+nnab),fockn,1,fock,1)
C
return
end
C
************************************************************************
subroutine trinab(nbasis,nbfcn,nbfan,nnab,nr,vv,wnab,scr)
************************************************************************
* Transform intermediates to NAB basis
************************************************************************
implicit none
integer nbasis,nbfcn,nbfan,nnab,nr
real*8 vv(nbfcn,nr),wnab(nbfan,nnab),scr(nbasis+nnab,nr)
C
scr(1:nbasis,1:nr)=vv(1:nbasis,1:nr)
call dgemm('t','n',nnab,nr,nbfan,1d0,wnab,nbfan,
$vv(nbasis+1,1),nbfcn,0d0,scr(nbasis+1,1),nbasis+nnab)
call dcopy((nbasis+nnab)*nr,scr,1,vv,1)
C
return
end
C
************************************************************************
subroutine cannab(nbasis,nbfcn,nbfan,nnab,fock,wnab,um,scr,eig,
$iout)
************************************************************************
* Canonicalize NAB virtuals
************************************************************************
implicit none
integer nbasis,nbfcn,nbfan,nnab,iout,iisyev
real*8 fock(nbfcn,nbfcn),wnab,um,scr,eig
integer*4 isyev
equivalence(isyev,iisyev)
C
call dsymm('l','l',nbfan,nnab,1d0,fock(nbasis+1,nbasis+1),nbfcn,
$wnab,nbfan,0d0,scr,nbfan)
call dgemm('t','n',nnab,nnab,nbfan,1.d0,wnab,nbfan,scr,nbfan,0.d0,
$um,nnab)
call dsyev('V','U',nnab,um,nnab,eig,scr,3*nnab**2,isyev)
if(isyev.ne.0) then
write(iout,*)
$'Fatal error at the canonicalization of CABS virtuals!'
call mrccend(1)
endif
call dgemm('n','n',nbfan,nnab,nnab,1.d0,wnab,nbfan,um,nnab,0.d0,
$scr,nbfan)
call dcopy(nbfan*nnab,scr,1,wnab,1)
C
return
end
C
************************************************************************
subroutine svdf12(nbasis,nbfcn,nbfan,nnab,cf12aa,wnaba,eig,scr,
$iout,nabtol,nr)
************************************************************************
* SVD of (pi|f12|Q) integrals
************************************************************************
implicit none
integer nbasis,nbfcn,nbfan,nnab,iout,iisyev,nr
real*8 cf12aa(nbfcn,nr),eig(nbfan),wnaba(nbfan,nbfan)
real*8 scr(*),nabtol
integer*4 isyev
equivalence(isyev,iisyev)
C
call dsyrk('U','N',nbfan,nr,-1d0,cf12aa(nbasis+1,1),
$nbfcn,0.d0,wnaba,nbfan)
call dsyev('V','U',nbfan,wnaba,nbfan,eig,scr,
$max(nbfan*nbfan,3*nbfan),isyev)
if(isyev.ne.0) then
write(iout,*) 'Fatal error at the construction of NAFs!'
call mrccend(1)
endif
c write(6,"(7f9.5)") -eig
do while(nnab.lt.nbfan.and.dabs(eig(nnab+1)).ge.nabtol)
nnab=nnab+1
enddo
C
return
end
C
************************************************************************
subroutine cnaf(nb,dfnbasis,nal,nbe,nbfcn,scftype,vifaa,vifbb,
$vifab,vifba,cf12aa,cf12bb,cf12ab,cf12ba,cf12r12aa,cf12r12bb,
$cf12r12ab,cf12r12ba,iout,verbosity,wnaf,eig,scr,naftol)
************************************************************************
* Construct the natural auxiliary functions
************************************************************************
implicit none
integer nb,dfnbasis,nal,nbe,nbfcn,iout,verbosity,i,iisyev,dfnb
real*8 vifaa(nbfcn,dfnbasis,nb),cf12aa(nbfcn,dfnbasis,nal),scr(*)
real*8 vifbb(nbfcn,dfnbasis,nb),cf12bb(nbfcn,dfnbasis,nbe),naftol
real*8 vifab(nbfcn,dfnbasis,nb),cf12ab(nbfcn,dfnbasis,nal)
real*8 vifba(nbfcn,dfnbasis,nb),cf12ba(nbfcn,dfnbasis,nbe)
real*8 cf12r12aa(nb,dfnbasis,nal),cf12r12ab(nb,dfnbasis,nal)
real*8 cf12r12bb(nb,dfnbasis,nbe),cf12r12ba(nb,dfnbasis,nbe)
real*8 eig(dfnbasis),wnaf(dfnbasis,dfnbasis)
character(len=5) scftype
character(len=8) c8,nafalg
character(len=16) c16
integer*4 isyev
equivalence(isyev,iisyev)
C
write(iout,*)
write(iout,*) 'Constructing natural auxiliary functions...'
C Build matrix W
wnaf=0.d0
do i=1,nb
call dsyrk('U','T',dfnbasis,nbfcn,-1d0,vifaa(1,1,i),nbfcn,1d0,
$wnaf,dfnbasis)
enddo
call getkey('nafalg',6,nafalg,8)
if(scftype.ne.'rhf '.and.trim(nafalg).eq.'albe') then
do i=1,nb
call dsyrk('U','T',dfnbasis,nbfcn,-1d0,vifbb(1,1,i),nbfcn,1d0,
$wnaf,dfnbasis)
call dsyrk('U','T',dfnbasis,nbfcn,-1d0,vifab(1,1,i),nbfcn,1d0,
$wnaf,dfnbasis)
call dsyrk('U','T',dfnbasis,nbfcn,-1d0,vifba(1,1,i),nbfcn,1d0,
$wnaf,dfnbasis)
enddo
call dscal(dfnbasis*dfnbasis,0.25d0,wnaf,1)
endif
C Diagonalize W
call dsyev('V','U',dfnbasis,wnaf,dfnbasis,eig,scr,
$max(nbfcn*dfnbasis,3*dfnbasis),isyev)
if(isyev.ne.0) then
write(iout,*) 'Fatal error at the construction of NAFs!'
call mrccend(1)
endif
c write(6,"(7f9.5)") -eig
dfnb=0
do while(dfnb.lt.dfnbasis.and.dabs(eig(dfnb+1)).ge.naftol)
dfnb=dfnb+1
enddo
write(c8,'(i8)') dfnb
write(c16,'(f16.1)') 100.d0*dble(dfnb)/dble(dfnbasis)
write(iout,*) 'Number of natural auxiliary functions: ' //
$trim(adjustl(c8)) // trim(' (' // adjustl(c16)) // '%)'
C Transform intermediates to NAF basis
call trinaf(nbfcn,dfnbasis,nb ,dfnb,vifaa ,vifaa ,wnaf,scr)
call trinaf(nbfcn,dfnbasis,nal,dfnb,cf12aa ,cf12aa ,wnaf,scr)
call trinaf(nb ,dfnbasis,nal,dfnb,cf12r12aa,cf12r12aa,wnaf,scr)
if(scftype.ne.'rhf ') then
call trinaf(nbfcn,dfnbasis,nb ,dfnb,vifbb ,vifbb ,wnaf,scr)
call trinaf(nbfcn,dfnbasis,nb ,dfnb,vifab ,vifab ,wnaf,scr)
call trinaf(nbfcn,dfnbasis,nb ,dfnb,vifba ,vifba ,wnaf,scr)
call trinaf(nbfcn,dfnbasis,nbe,dfnb,cf12bb ,cf12bb ,wnaf,scr)
call trinaf(nbfcn,dfnbasis,nal,dfnb,cf12ab ,cf12ab ,wnaf,scr)
call trinaf(nbfcn,dfnbasis,nbe,dfnb,cf12ba ,cf12ba ,wnaf,scr)
call trinaf(nb ,dfnbasis,nbe,dfnb,cf12r12bb,cf12r12bb,wnaf,scr)
call trinaf(nb ,dfnbasis,nal,dfnb,cf12r12ab,cf12r12ab,wnaf,scr)
call trinaf(nb ,dfnbasis,nbe,dfnb,cf12r12ba,cf12r12ba,wnaf,scr)
endif
dfnbasis=dfnb
call timer
C
return
end
C
************************************************************************
subroutine trinaf(nbfcn,dfnbasis,nb,dfnb,vvo,vvn,wnaf,scr)
************************************************************************
* Transform intermediates to NAF basis
************************************************************************
implicit none
integer nbfcn,dfnbasis,nb,dfnb,i
real*8 vvo(nbfcn,dfnbasis,nb),vvn(nbfcn,dfnb,nb)
real*8 scr(nbfcn,dfnbasis),wnaf(dfnbasis,dfnb)
C
do i=1,nb
call dcopy(nbfcn*dfnbasis,vvo(1,1,i),1,scr,1)
call dgemm('n','n',nbfcn,dfnb,dfnbasis,1d0,scr,nbfcn,
$wnaf,dfnbasis,0d0,vvn(1,1,i),nbfcn)
enddo
C
return
end
C
************************************************************************
subroutine diagd(nvirtal,nva,iout,daa,eig,scr,eps)
************************************************************************
* Diagonalize MP2 density matrix
************************************************************************
implicit none
integer iisyev,nvirtal,nva,iout
real*8 eig(nvirtal),daa(nvirtal,nvirtal),scr(*),eps
integer*4 isyev
equivalence(isyev,iisyev)
C
daa=-daa
call dsyev('V','U',nvirtal,daa,nvirtal,eig,scr,
$max(nvirtal*nvirtal,3*nvirtal),isyev)
if(isyev.ne.0) then
write(iout,*) 'Fatal error at the construction of NOs!'
call mrccend(1)
endif
c write(6,"(7f9.5)") -eig
nva=0
do while(nva.lt.nvirtal.and.dabs(eig(nva+1)).ge.eps)
nva=nva+1
enddo
C
return
end
C
************************************************************************
subroutine cfno(nb,dfnbasis,nal,nbe,ncore,focka,fockb,nbfcn,nbfan,
$scftype,vifaa,vifbb,vifab,vifba,cf12aa,cf12bb,cf12ab,cf12ba,
$cf12r12aa,cf12r12bb,cf12r12ab,cf12r12ba,iout,verbosity,daa,dbb,
$eig,um,fockan,fockbn,scr,eps,nvirtal,nvirtbe,uma,umb,scrfile1,
$nbasis,eva,evb)
************************************************************************
* Construct the frozen natural orbitals
************************************************************************
implicit none
integer nb,dfnbasis,nal,nbe,ncore,nbfcn,nbfan,iout,verbosity,i,j
integer iisyev,nnb,nbasis,nvirtal,nvirtbe,nva,nvb,nnbfan,scrfile1
real*8 focka(nbfcn,nbfcn),fockb(nbfcn,nbfcn),um(nbfcn,nbfcn)
real*8 vifaa(nbfcn,dfnbasis,nb),cf12aa(nbfcn,dfnbasis,nal)
real*8 vifbb(nbfcn,dfnbasis,nb),cf12bb(nbfcn,dfnbasis,nbe),eps
real*8 vifab(nbfcn,dfnbasis,nb),cf12ab(nbfcn,dfnbasis,nal)
real*8 vifba(nbfcn,dfnbasis,nb),cf12ba(nbfcn,dfnbasis,nbe)
real*8 cf12r12aa(nb,dfnbasis,nal),cf12r12ab(nb,dfnbasis,nal)
real*8 cf12r12bb(nb,dfnbasis,nbe),cf12r12ba(nb,dfnbasis,nbe)
real*8 eig(nbfcn),daa(nvirtal,nvirtal),dbb(nvirtbe,nvirtbe)
real*8 fockan(nbfcn,nbfcn),fockbn(nbfcn,nbfcn),eva(*),evb(*)
real*8 uma(nvirtal+nbfan,nvirtal+nbfan),scr(nbfcn,dfnbasis,nal)
real*8 umb(nvirtbe+nbfan,nvirtbe+nbfan)
character(len=5) scftype
character(len=8) c8
character(len=16) c16
C
write(iout,*)
write(iout,*) 'Constructing natural orbitals...'
call diagd(nvirtal,nva,iout,daa,eig,scr,eps)
if(scftype.ne.'rhf ') then
call diagd(nvirtbe,nvb,iout,dbb,eig,scr,eps)
else
nvb=nva
endif
nnb=max(nal+nva,nbe+nvb)
nva=nnb-nal
nvb=nnb-nbe
nnbfan=nbfan+nb-nnb
write(c8,'(i8)') nva
write(c16,'(f16.2)') 100.d0*dble(nva)/dble(nvirtal)
write(iout,*) 'Number of virtual natural orbitals: ' //
$trim(adjustl(c8)) // trim(' (' // adjustl(c16)) // '%)'
C Canonicalize virtuals NOs
call cannab(ncore+nal,nbfcn,nvirtal,nva,focka,daa,um,scr,eig,iout)
if(scftype.ne.'rhf ')
$call cannab(ncore+nbe,nbfcn,nvirtbe,nvb,fockb,dbb,um,scr,eig,iout)
C Canonicalize new CABS virtuals -- not needed for CC!
call cannca(nbasis,nbfcn,nbfan,nnbfan,focka,daa,um,scr,eig,iout,
$nvirtal,ncore+nal,nva,uma)
if(scftype.ne.'rhf ')
$call cannca(nbasis,nbfcn,nbfan,nnbfan,fockb,dbb,um,scr,eig,iout,
$nvirtbe,ncore+nbe,nvb,umb)
C Transform Fock-matrix to the new virtual + CABS virtual basis
call trfnab(ncore+nal,nbfcn,nvirtal+nbfan,nvirtal+nbfan,focka,
$fockan,uma,scr)
call trinca(nbfcn,dfnbasis,nb,nvirtal,nvirtal,nbfan,ncore,nal,nal,
$nva,vifaa,uma,daa,scr)
call trinab(ncore+nal,nbfcn,nvirtal+nbfan,nvirtal+nbfan,
$dfnbasis*nal,cf12aa,uma,scr)
call trinab(nal,nb,nvirtal,nva,dfnbasis*nal,cf12r12aa,daa,scr)
do i=ncore+nal+1,nbfcn
eva(i)=focka(i,i)
enddo
if(scftype.ne.'rhf ') then
call trfnab(ncore+nbe,nbfcn,nvirtbe+nbfan,nvirtbe+nbfan,fockb,
$fockbn,umb,scr)
call trinca(nbfcn,dfnbasis,nb,nvirtbe,nvirtbe,nbfan,ncore,nbe,nbe,
$nvb,vifbb,umb,dbb,scr)
call trinca(nbfcn,dfnbasis,nb,nvirtal,nvirtbe,nbfan,ncore,nal,nbe,
$nva,vifab,umb,daa,scr)
call trinca(nbfcn,dfnbasis,nb,nvirtbe,nvirtal,nbfan,ncore,nbe,nal,
$nvb,vifba,uma,dbb,scr)
call trinab(ncore+nbe,nbfcn,nvirtbe+nbfan,nvirtbe+nbfan,
$dfnbasis*nbe,cf12bb,umb,scr)
call trinab(ncore+nbe,nbfcn,nvirtbe+nbfan,nvirtbe+nbfan,
$dfnbasis*nal,cf12ab,umb,scr)
call trinab(ncore+nal,nbfcn,nvirtal+nbfan,nvirtal+nbfan,
$dfnbasis*nbe,cf12ba,uma,scr)
call trinab(nbe,nb,nvirtbe,nvb,dfnbasis*nbe,cf12r12bb,dbb,scr)
call trinab(nbe,nb,nvirtbe,nvb,dfnbasis*nal,cf12r12ab,dbb,scr)
call trinab(nal,nb,nvirtal,nva,dfnbasis*nbe,cf12r12ba,daa,scr)
do i=ncore+nbe+1,nbfcn
evb(i)=fockb(i,i)
enddo
endif
C MO coefficients in the new basis
c call trfmo(nbasis,ncore,nal,nbe,nvirtal,nvirtbe,scrfile1,fockan,
c $fockbn,daa,dbb,scr,scftype)
C Save MO to NO transformation matrix - szemet
open(scrfile1,file='MONO',form='unformatted')
write(scrfile1) nva
write(scrfile1) ((daa(i,j),i=1,nvirtal),j=nva+1,nvirtal),
$ ((daa(i,j),i=1,nvirtal),j=1,nva)
write(scrfile1) ((focka(i,j),i=ncore+1,nbasis),j=ncore+1,nbasis)
if(scftype.ne.'rhf ') then
write(scrfile1) nvb
write(scrfile1) ((dbb(i,j),i=1,nvirtbe),j=nvb+1,nvirtbe),
$ ((dbb(i,j),i=1,nvirtbe),j=1,nvb)
write(scrfile1) ((fockb(i,j),i=ncore+1,nbasis),j=ncore+1,nbasis)
endif
close(scrfile1)
C
nb=nnb
nbasis=ncore+nb
nbfan=nnbfan
nvirtal=nva
nvirtbe=nvb
call timer
C
return
end
C
************************************************************************
subroutine cannca(nbasis,nbfcn,nbfan,nnbfan,fock,daa,um,scr,eig,
$iout,nvirtal,nal,nva,uma)
************************************************************************
* Canonicalize new CABS virtuals
************************************************************************
implicit none
integer nbasis,nbfcn,nbfan,nnbfan,iout,iisyev,nn,nvirtal,nva,nal,i
real*8 fock(nbfcn,nbfcn),daa(nvirtal,nvirtal),um(nnbfan,nnbfan)
real*8 scr,eig(nnbfan),uma(nvirtal+nbfan,nvirtal+nbfan)
integer*4 isyev
equivalence(isyev,iisyev)
C
nn=nnbfan-nbfan
um(nn+1:nnbfan,nn+1:nnbfan)=fock(nbasis+1:nbfcn,nbasis+1:nbfcn)
call dsymm('l','l',nvirtal,nn,1d0,fock(nal+1,nal+1),nbfcn,
$daa(1,nva+1),nvirtal,0d0,scr,nvirtal)
call dgemm('t','n',nn,nn,nvirtal,1.d0,daa(1,nva+1),nvirtal,scr,
$nvirtal,0.d0,um,nnbfan)
call dgemm('t','n',nn,nbfan,nvirtal,1.d0,daa(1,nva+1),nvirtal,
$fock(nal+1,nbasis+1),nbfcn,0.d0,um(1,nn+1),nnbfan)
call dsyev('V','U',nnbfan,um,nnbfan,eig,scr,
$max(nnbfan*nnbfan,3*nnbfan),isyev)
if(isyev.ne.0) then
write(iout,*)
$'Fatal error at the canonicalization of CABS virtuals!'
call mrccend(1)
endif
uma=0.d0
uma(1:nvirtal,1:nva)=daa(1:nvirtal,1:nva)
call dgemm('n','n',nvirtal,nnbfan,nn,1.d0,daa(1,nva+1),
$nvirtal,um,nnbfan,0.d0,uma(1,nva+1),nvirtal+nbfan)
uma(nvirtal+1:nvirtal+nbfan,nva+1:nvirtal+nbfan)=
$um(nn+1:nnbfan,1:nnbfan)
C
return
end
C
************************************************************************
subroutine trinca(nbfcn,dfnbasis,nb,nvirtal,nvirtbe,nbfan,ncore,
$nal,nbe,nva,vif,umb,daa,scr)
************************************************************************
* Transform intermediate vif to the new vrtual + virtual CABS basis
************************************************************************
implicit none
integer nbfcn,dfnbasis,nb,nvirtal,nvirtbe,nbfan,ncore,nal,nbe,nva
real*8 vif(nbfcn,dfnbasis,nb),scr(nbfcn,dfnbasis,nb)
real*8 umb(nvirtbe+nbfan,nvirtbe+nbfan),daa(nvirtal,nvirtal)
C
scr(1:ncore+nbe,:,:)=vif(1:ncore+nbe,:,:)
call dgemm('t','n',nvirtbe+nbfan,dfnbasis*nb,nvirtbe+nbfan,1d0,
$umb,nvirtbe+nbfan,vif(ncore+nbe+1,1,1),nbfcn,0d0,
$scr(ncore+nbe+1,1,1),nbfcn)
vif(:,:,1:nal)=scr(:,:,1:nal)
call dgemm('n','n',nbfcn*dfnbasis,nva,nvirtal,1d0,scr(1,1,nal+1),
$nbfcn*dfnbasis,daa,nvirtal,0d0,vif(1,1,nal+1),nbfcn*dfnbasis)
C
return
end
C
************************************************************************
subroutine trfmo(nbasis,ncore,nal,nbe,nvirtal,nvirtbe,scrfile1,
$cmo,mo,daa,dbb,scr,scftype)
************************************************************************
* Construct MO coefficients in NO basis
************************************************************************
implicit none
integer nbasis,ncore,nal,nbe,nvirtal,nvirtbe,scrfile1
real*8 mo(nbasis,nbasis),cmo(nbasis,nbasis),scr(*)
real*8 daa(nvirtal,nvirtal),dbb(nvirtbe,nvirtbe)
character(len=5) scftype
C
open(scrfile1,file='MOCOEF',form='unformatted')
call readmo(scr,scr,mo,scrfile1,nbasis,nbasis)
call dcopy(nbasis*(ncore+nal),mo,1,cmo,1)
call dgemm('n','n',nbasis,nvirtal,nvirtal,1d0,mo(1,ncore+nal+1),
$nbasis,daa,nvirtal,0d0,cmo(1,ncore+nal+1),nbasis)
if(scftype.ne.'rhf ') then
call readmo(scr,scr,mo,scrfile1,nbasis,nbasis)
rewind(scrfile1)
call wrtmo(scr,scr,cmo,scrfile1,0.d0,nbasis,nbasis)
call dcopy(nbasis*(ncore+nbe),mo,1,cmo,1)
call dgemm('n','n',nbasis,nvirtbe,nvirtbe,1d0,mo(1,ncore+nbe+1),
$nbasis,dbb,nvirtbe,0d0,cmo(1,ncore+nbe+1),nbasis)
else
rewind(scrfile1)
endif
call wrtmo(scr,scr,cmo,scrfile1,0.d0,nbasis,nbasis)
close(scrfile1)
C
return
end
C
************************************************************************
subroutine cecoup(nbfcn,ncore,nal,nbe,nalc,nbec,nbasis,dfnbasis,
$nbfan,iout,scftype,verbosity,eref,ef12,dcore,imem,focka,fockb,eoa,
$eob,eva,evb,nvirtal,nvirtbe,vifaa,vifbb,vifab,vifba,cf12aa,cf12bb,
$cf12ab,cf12ba,tscalea,tscaleb,tfact,cctol,ecoup,ef12old,emp2)
************************************************************************
* Calculate MP2-F12 energy
************************************************************************
implicit none
integer nbfcn,nocc,nvirt,i,ncore,nal,nbe,nbasis,dfnbasis,nbfan
integer verbosity,nvirtal,nvirtbe,iout,nalc,nbec
integer irsa,if12ij,if12sy,ivintpq,dblalloc,imem,irsb
real*8 focka(nbfcn,nbfcn),fockb(nbfcn,nbfcn),eoa(*),eob(*),eva(*)
real*8 dcore(*),evb(*),tfact,ecoup,cctol,emp2,ef12,et1,eref
real*8 vifaa,vifbb,vifab,vifba,cf12aa,cf12bb,cf12ab,cf12ba,emp2old
real*8 tscalea(nalc,2),tscaleb(nbec,2),ecoupold,ef12old
character(len=5) scftype
C Calculate MP2 singles contribution
write(iout,*)
write(iout,*)
$'Calculating MP2 energy and coupling term in NO basis...'
ecoupold=ecoup
et1=0.d0
tscalea(:,2)=0.d0
tscaleb(:,2)=0.d0
call mp2sing(nbfcn,ncore,nal,nvirtal,focka,eoa,eva,et1,
$tscalea(1,2))
if(scftype.eq.'rhf ') then
et1=2.d0*et1
call dscal(nalc,2.d0,tscalea(1,2),1)
else
call mp2sing(nbfcn,ncore,nbe,nvirtbe,fockb,eob,evb,et1,
$tscaleb(1,2))
endif
emp2old=emp2
emp2=et1
C Calculate F12 contribution
irsa =dblalloc(nbfcn*nbfcn)
if12ij =dblalloc(nbfcn*nbfcn)
if12sy =dblalloc(nbfcn*nbfcn)
ivintpq =dblalloc(nbasis*nbasis)
ecoup=0.d0
if(scftype.eq.'rhf ') then
C Closed-shell
call ecoup_st(nal,nalc,nbfcn,nbfan,nvirtal,ncore,nbasis,
$dfnbasis,focka,eoa,eva,iout,verbosity,emp2,ecoup,vifaa,cf12aa,
$dcore(irsa),dcore(if12ij),dcore(ivintpq),tscalea,cctol)
else
C alpha-alpha
call ecoup_ss(nal,nalc,nbfcn,nbfan,nvirtal,ncore,nbasis,
$dfnbasis,focka,eoa,eva,iout,verbosity,emp2,ecoup,vifaa,cf12aa,
$dcore(irsa),dcore(if12ij),dcore(if12sy),dcore(ivintpq),
$tscalea,cctol)
C beta-beta
call ecoup_ss(nbe,nbec,nbfcn,nbfan,nvirtbe,ncore,nbasis,
$dfnbasis,fockb,eob,evb,iout,verbosity,emp2,ecoup,vifbb,cf12bb,
$dcore(irsa),dcore(if12ij),dcore(if12sy),dcore(ivintpq),
$tscaleb,cctol)
C alpha-beta
irsb=dblalloc(nbfcn*nbfcn)
call ecoup_os(nal,nbe,nalc,nbec,nbfcn,nbfan,nvirtal,nvirtbe,
$ncore,nbasis,dfnbasis,focka,fockb,eoa,eva,eob,evb,iout,verbosity,
$emp2,ecoup,vifaa,vifbb,vifab,vifba,cf12aa,cf12bb,cf12ab,cf12ba,
$dcore(irsa),dcore(irsb),dcore(if12sy),dcore(ivintpq),tscalea,
$tscaleb,cctol)
endif
C
write(iout,"(' MP2 correlation energy [au]: 'f24.12)") emp2
write(iout,"(' MP2 energy [au]: 'f24.12)")
$eref+emp2
write(iout,"(' MP2 correction [au]: 'f24.12)")
$emp2old-emp2
write(iout,"(' Coupling term [au]: 'f24.12)") ecoup
write(iout,"(' Coupling term correction [au]: 'f24.12)")
$ecoupold-ecoup
tfact=ef12old/emp2
ef12=ef12+ecoupold-ecoup
call timer
C Calculating scaling parameters for (T+)
call mtscale(nalc,nbec,ncore,iout,verbosity,tscalea,tscaleb,
$.false.,scftype,cctol)
call dbldealloc(irsa)
C
return
end
C
************************************************************************
subroutine ecoup_st(nal,nalc,nbfcn,nbfan,nvirtal,ncore,nbasis,
$dfnbasis,fa,eoa,eva,iout,verbosity,emp2,ecoup,vifaa,cf12aa,
$rsij,f12ij,vintpq,tscalea,cctol)
************************************************************************
* Calculate coupling term for RHF
************************************************************************
implicit none
integer nbfcn,nbfan,nal,nalc,nvirtal,i,a,j,b,ncore,ij
integer o,p,q,nbasis,dfnbasis,iout,verbosity
real*8 eoa(*),eva(*),fa(nbfcn,nbfcn),eij,eijb,emp2,tmp,cctol
real*8 e12,ecoup,eco,tscalea(nalc,2)
real*8 vifaa(nbfcn,dfnbasis,nalc),cf12aa(nbfcn,dfnbasis,nalc)
real*8 rsij(nbfcn,nbfcn),f12ij(nbfcn,nbfcn),vintpq(nbasis,nbasis)
C Loop over pairs
ij=0
do i=1,nalc
do j=1,i
ij=ij+1
C Coupling term
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,1d0, vifaa(1,1,i),
$nbfcn,cf12aa(1,1,j),nbfcn,0d0,f12ij,nbfcn)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,1d0,cf12aa(1,1,i),
$nbfcn, vifaa(1,1,j),nbfcn,1d0,f12ij,nbfcn)
rsij=0.d0
call dgemm('n','n',nvirtal,nbasis,nbfan,1d0,f12ij(nal+1,
$nbasis+1),nbfcn,fa(nbasis+1,1),nbfcn,0d0,rsij(nal+1,1),nbfcn)
call dgemm('n','n',nbasis,nvirtal,nbfan,1d0,fa(1,nbasis+1),
$nbfcn,f12ij(nbasis+1,nal+1),nbfcn,1d0,rsij(1,nal+1),nbfcn)
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0,vifaa(1,1,i),
$nbfcn,vifaa(1,1,j),nbfcn,0d0,vintpq,nbasis)
C Intermediate C
eij=eoa(ncore+i)+eoa(ncore+j)
eco=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,b,eijb) REDUCTION(+:eco)
do b=nal+1,nbasis
eijb=eij-eva(b-nal)
do a=nal+1,nbasis
eco=eco+(40d0*vintpq(a,b)-8d0*vintpq(b,a)
$+7d0*rsij(a,b)+rsij(b,a))*rsij(a,b)/(eijb-eva(a-nal))
enddo
enddo
C$OMP END PARALLEL DO
eco=eco/16d0
if(i.eq.j) eco=0.5d0*eco
ecoup=ecoup+eco
C Calculate MP2 energy
tmp=0.d0
C$OMP PARALLEL DO SCHEDULE(DYNAMIC) DEFAULT(SHARED) PRIVATE(a,b,eijb)
C$OMP& REDUCTION(+:tmp)
do b=nal+1,nal+nvirtal
eijb=eij-eva(b-nal)
tmp=tmp+0.5d0*vintpq(b,b)**2/(eijb-eva(b-nal))
do a=nal+1,b-1
tmp=tmp+(
$(vintpq(a,b)-vintpq(b,a))**2+vintpq(a,b)*vintpq(b,a))/
$(eijb-eva(a-nal))
enddo
enddo
C$OMP END PARALLEL DO
tmp=4d0*tmp
if(i.eq.j) tmp=0.5d0*tmp
emp2=emp2+tmp
tscalea(i,2)=tscalea(i,2)+0.5d0*tmp
tscalea(j,2)=tscalea(j,2)+0.5d0*tmp
enddo
enddo
C
return
end
C
************************************************************************
subroutine ecoup_ss(nal,nalc,nbfcn,nbfan,nvirtal,ncore,nbasis,
$dfnbasis,fa,eoa,eva,iout,verbosity,emp2,ecoup,vifaa,cf12aa,
$rsij,f12ij,f12sy,vintpq,tscalea,cctol)
************************************************************************
* Calculate coupling term for same-spin pairs
************************************************************************
implicit none
integer nbfcn,nbfan,nal,nalc,nvirtal,i,a,j,b,ncore,dfnbasis,ij
integer o,p,q,nbasis,iout,verbosity
real*8 emp2,tmp,eijb,ecoup,tscalea(nalc,2),cctol
real*8 eoa(*),eva(*),fa(nbfcn,nbfcn),eij,e12,eco
real*8 vifaa(nbfcn,dfnbasis,nalc),cf12aa(nbfcn,dfnbasis,nalc)
real*8 rsij(nbfcn,nbfcn),f12ij(nbfcn,nbfcn),f12sy(nbfcn,nbfcn)
real*8 vintpq(nbasis,nbasis)
C Loop over pairs
ij=0
do i=1,nalc
do j=1,i-1
ij=ij+1
C Coupling term
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,1d0, vifaa(1,1,i),
$nbfcn,cf12aa(1,1,j),nbfcn,0d0,f12ij,nbfcn)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,1d0,cf12aa(1,1,i),
$nbfcn, vifaa(1,1,j),nbfcn,1d0,f12ij,nbfcn)
f12sy=-transpose(f12ij)
call daxpy(nbfcn*nbfcn,1d0,f12ij,1,f12sy,1)
call dgemm('n','n',nbasis,nvirtal,nbfan,1d0,fa(1,nbasis+1),
$nbfcn,f12sy(nbasis+1,nal+1),nbfcn,0d0,rsij(1,nal+1),nbfcn)
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0,vifaa(1,1,i),
$nbfcn,vifaa(1,1,j),nbfcn,0d0,vintpq,nbasis)
C$OMP PARALLEL DO SCHEDULE(DYNAMIC) DEFAULT(SHARED) PRIVATE(p,q,tmp)
do q=1,nbasis
do p=1,q-1
tmp=vintpq(p,q)-vintpq(q,p)
vintpq(p,q)=tmp
enddo
enddo
C$OMP END PARALLEL DO
eij=eoa(ncore+i)+eoa(ncore+j)
eco=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,b,eijb,tmp) REDUCTION(+:eco)
do b=nal+1,nal+nvirtal
eijb=eij-eva(b-nal)
do a=nal+1,b-1
tmp=0.25d0*(rsij(a,b)-rsij(b,a))
eco=eco+(2d0*vintpq(a,b)+tmp)*tmp/(eijb-eva(a-nal))
enddo
enddo
C$OMP END PARALLEL DO
ecoup=ecoup+eco
C Calculate MP2 energy
tmp=0.d0
C$OMP PARALLEL DO SCHEDULE(DYNAMIC) DEFAULT(SHARED) PRIVATE(a,b,eijb)
C$OMP& REDUCTION(+:tmp)
do b=nal+1,nal+nvirtal
eijb=eij-eva(b-nal)
do a=nal+1,b-1
tmp=tmp+vintpq(a,b)**2/(eijb-eva(a-nal))
enddo
enddo
C$OMP END PARALLEL DO
emp2=emp2+tmp
tscalea(i,2)=tscalea(i,2)+0.5d0*tmp
tscalea(j,2)=tscalea(j,2)+0.5d0*tmp
enddo
enddo
C
return
end
C
************************************************************************
subroutine ecoup_os(nal,nbe,nalc,nbec,nbfcn,nbfan,nvirtal,nvirtbe,
$ncore,nbasis,dfnbasis,fa,fb,eoa,eva,eob,evb,iout,verbosity,emp2,
$ecoup,vifaa,vifbb,vifab,vifba,cf12aa,cf12bb,cf12ab,cf12ba,
$rsa,rsb,f12sy,vintpq,tscalea,tscaleb,cctol)
************************************************************************
* Calculate coupling term for opposite-spin pairs
************************************************************************
implicit none
integer nbfcn,nbfan,nal,nbe,nalc,nbec,nvirtal,nvirtbe,ncore,k,l
integer o,s,p,q,u,w,nbasis,c,m,n,iout,verbosity,i,a,j,b,dfnbasis
real*8 eoa(*),eva(*),eob(*),evb(*),fa(nbfcn,nbfcn),fb(nbfcn,nbfcn)
real*8 emp2,tmp,eij,eijb
real*8 rsa(nbfcn,nbfcn),rsb(nbfcn,nbfcn),cctol
real*8 vifaa(nbfcn,dfnbasis,nalc),cf12aa(nbfcn,dfnbasis,nalc)
real*8 vifbb(nbfcn,dfnbasis,nbec),cf12bb(nbfcn,dfnbasis,nbec)
real*8 vifab(nbfcn,dfnbasis,nalc),cf12ab(nbfcn,dfnbasis,nalc)
real*8 vifba(nbfcn,dfnbasis,nbec),cf12ba(nbfcn,dfnbasis,nbec)
real*8 f12sy(nbfcn,nbfcn),ecoup,eco
real*8 vintpq(nbasis,nbasis),tscalea(nalc,2),tscaleb(nbec,2)
C
do i=1,nalc
do j=1,nbec
C Coupling term
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,3d0/8d0, vifaa(1,1,i),
$nbfcn,cf12bb(1,1,j),nbfcn,0d0,f12sy,nbfcn)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,3d0/8d0,cf12aa(1,1,i),
$nbfcn, vifbb(1,1,j),nbfcn,1d0,f12sy,nbfcn)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,1d0/8d0, vifba(1,1,j),
$nbfcn,cf12ab(1,1,i),nbfcn,1d0,f12sy,nbfcn)
call dgemm('n','t',nbfcn,nbfcn,dfnbasis,1d0/8d0,cf12ba(1,1,j),
$nbfcn, vifab(1,1,i),nbfcn,1d0,f12sy,nbfcn)
call dgemm('n','n',nvirtal,nbasis,nbfan,1d0,f12sy(nal+1,
$nbasis+1),nbfcn,fb(nbasis+1,1),nbfcn,0d0,rsa(nal+1,1),nbfcn)
call dgemm('n','n',nbasis,nvirtbe,nbfan,1d0,fa(1,nbasis+1),
$nbfcn,f12sy(nbasis+1,nbe+1),nbfcn,0d0,rsb(1,nbe+1),nbfcn)
call dgemm('n','t',nbasis,nbasis,dfnbasis,1d0,vifaa(1,1,i),
$nbfcn,vifbb(1,1,j),nbfcn,0d0,vintpq,nbasis)
eij=eoa(ncore+i)+eob(ncore+j)
eco=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,b,eijb,tmp) REDUCTION(+:eco)
do b=nbe+1,nbasis
eijb=eij-evb(b-nbe)
do a=nal+1,nbasis
tmp=rsb(a,b)+rsa(a,b)
eco=eco+(vintpq(a,b)+0.5d0*tmp)*tmp/(eijb-eva(a-nal))
enddo
enddo
C$OMP END PARALLEL DO
eco=2d0*eco
ecoup=ecoup+eco
C Calculate MP2 energy
tmp=0.d0
C$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(a,b,eijb) REDUCTION(+:tmp)
do b=nbe+1,nbasis
eijb=eij-evb(b-nbe)
do a=nal+1,nbasis
tmp=tmp+vintpq(a,b)**2/(eijb-eva(a-nal))
enddo
enddo
C$OMP END PARALLEL DO
emp2=emp2+tmp
tscalea(i,2)=tscalea(i,2)+0.5d0*tmp
tscaleb(j,2)=tscaleb(j,2)+0.5d0*tmp
enddo
enddo
C
return
end
C
************************************************************************
subroutine mtscale(nalc,nbec,ncore,iout,verbosity,tscalea,tscaleb,
$lfno,scftype,cctol)
************************************************************************
* Calculating scaling parameters for (T+)
************************************************************************
implicit none
integer nalc,nbec,ncore,iout,verbosity,i
real*8 tscalea(nalc,2),tscaleb(nbec,2),cctol
character(len=5) scftype
logical lfno
C Calculating scaling parameters for (T+)
if(verbosity.gt.2) then
write(iout,*)
write(iout,*) 'Contribution of orbitals to correlation energy:'
write(iout,*)
$' MO MP2 [au] F12 [au] ' //
$'(MP2+F12)/MP2'
if(scftype.eq.'rhf ') then
do i=1,nalc
if(dabs(tscalea(i,2)).gt.cctol)
$write(iout,"(i4,1x,f20.12,2x,f20.12,2x,f20.12)") ncore+i,
$tscalea(i,2),tscalea(i,1),(tscalea(i,1)+tscalea(i,2))/tscalea(i,2)
enddo
else
do i=1,nalc
if(dabs(tscalea(i,2)).gt.cctol)
$write(iout,"(i4,a1,f20.12,2x,f20.12,2x,f20.12)") ncore+i,'a',
$tscalea(i,2),tscalea(i,1),(tscalea(i,1)+tscalea(i,2))/tscalea(i,2)
enddo
do i=1,nbec
if(dabs(tscaleb(i,2)).gt.cctol)
$write(iout,"(i4,a1,f20.12,2x,f20.12,2x,f20.12)") ncore+i,'b',
$tscaleb(i,2),tscaleb(i,1),(tscaleb(i,1)+tscaleb(i,2))/tscaleb(i,2)
enddo
endif
endif
if(.not.lfno) then
do i=1,nalc
if(dabs(tscalea(i,2)).gt.cctol) then
tscalea(i,1)=(tscalea(i,1)+tscalea(i,2))/tscalea(i,2)
else
tscalea(i,1)=0.d0
endif
enddo
if(scftype.ne.'rhf ') then
do i=1,nbec
if(dabs(tscaleb(i,2)).gt.cctol) then
tscaleb(i,1)=(tscaleb(i,1)+tscaleb(i,2))/tscaleb(i,2)
else
tscaleb(i,1)=0.d0
endif
enddo
endif
endif
C
return
end
C
************************************************************************
subroutine t3f12(nal,nbe,nval,nvbe,ncore,nb,nbfcn,nbfan,dfnbasis,
$vifaa,vifbb,vifab,vifba,cf12aa,cf12bb,cf12ab,cf12ba,vintaa,vintbb,
$vintab,vintba,waa,wbb,wab,wba,itol,scrfile1)
************************************************************************
* Calculating F12 T2 contribution to T3
************************************************************************
implicit none
integer nal,nbe,nval,nvbe,ncore,nb,nbfcn,nbfan,dfnbasis,scrfile1
integer*2 a,b,c,d,i,j,k
real*8 vifaa(nbfcn,dfnbasis,nb),cf12aa(nbfcn,dfnbasis,nal),tmp!,sum
real*8 vifbb(nbfcn,dfnbasis,nb),cf12bb(nbfcn,dfnbasis,nbe),itol
real*8 vifab(nbfcn,dfnbasis,nb),cf12ab(nbfcn,dfnbasis,nal)
real*8 vifba(nbfcn,dfnbasis,nb),cf12ba(nbfcn,dfnbasis,nbe)
real*8 vintaa(nbfan,nal,nval,nval),vintbb(nbfan,nbe,nvbe,nvbe)
real*8 vintab(nbfan,nbe,nval,nvbe),vintba(nal,nbfan,nval,nvbe)
real*8 waa(nval,nbfan,nal,nal),wbb(nvbe,nbfan,nbe,nbe)
real*8 wab(nbfan,nvbe,nal,nbe),wba(nval,nbfan,nal,nbe)
C Assemble Coulomb integrals
do b=1,nval
do a=1,nval
call dgemm('n','t',nbfan,nal,dfnbasis, 1.d0,
$vifaa(ncore+nb+1,1,nal+a),nbfcn,vifaa(ncore+1,1,nal+b),nbfcn,0.d0,
$vintaa(1,1,a,b),nbfan)
call dgemm('n','t',nbfan,nal,dfnbasis,-1.d0,
$vifaa(ncore+nb+1,1,nal+b),nbfcn,vifaa(ncore+1,1,nal+a),nbfcn,1.d0,
$vintaa(1,1,a,b),nbfan)
enddo
enddo
do b=1,nvbe
do a=1,nvbe
call dgemm('n','t',nbfan,nbe,dfnbasis, 1.d0,
$vifbb(ncore+nb+1,1,nbe+a),nbfcn,vifbb(ncore+1,1,nbe+b),nbfcn,0.d0,
$vintbb(1,1,a,b),nbfan)
call dgemm('n','t',nbfan,nbe,dfnbasis,-1.d0,
$vifbb(ncore+nb+1,1,nbe+b),nbfcn,vifbb(ncore+1,1,nbe+a),nbfcn,1.d0,
$vintbb(1,1,a,b),nbfan)
enddo
enddo
do b=1,nvbe
do a=1,nval
call dgemm('n','t',nbfan,nbe,dfnbasis, 1.d0,
$vifaa(ncore+nb+1,1,nal+a),nbfcn,vifbb(ncore+1,1,nbe+b),nbfcn,0.d0,
$vintab(1,1,a,b),nbfan)
enddo
enddo
do b=1,nvbe
do a=1,nval
call dgemm('n','t',nal,nbfan,dfnbasis, 1.d0,
$vifaa(ncore+1,1,nal+a),nbfcn,vifbb(ncore+nb+1,1,nbe+b),nbfcn,0.d0,
$vintba(1,1,a,b),nal)
enddo
enddo
C Assemble F12 integrals
do j=1,nal
do i=1,nal
call dgemm('n','t',nval,nbfan,dfnbasis, 0.25d0,
$ vifaa(ncore+nal+1,1,i),nbfcn,cf12aa(ncore+nb+1,1,j),nbfcn,0.d0,
$waa(1,1,i,j),nval)
call dgemm('n','t',nval,nbfan,dfnbasis, 0.25d0,
$cf12aa(ncore+nal+1,1,i),nbfcn, vifaa(ncore+nb+1,1,j),nbfcn,1.d0,
$waa(1,1,i,j),nval)
call dgemm('n','t',nval,nbfan,dfnbasis,-0.25d0,
$ vifaa(ncore+nal+1,1,j),nbfcn,cf12aa(ncore+nb+1,1,i),nbfcn,1.d0,
$waa(1,1,i,j),nval)
call dgemm('n','t',nval,nbfan,dfnbasis,-0.25d0,
$cf12aa(ncore+nal+1,1,j),nbfcn, vifaa(ncore+nb+1,1,i),nbfcn,1.d0,
$waa(1,1,i,j),nval)
enddo
enddo
do j=1,nbe
do i=1,nbe
call dgemm('n','t',nvbe,nbfan,dfnbasis, 0.25d0,
$ vifbb(ncore+nbe+1,1,i),nbfcn,cf12bb(ncore+nb+1,1,j),nbfcn,0.d0,
$wbb(1,1,i,j),nvbe)
call dgemm('n','t',nvbe,nbfan,dfnbasis, 0.25d0,
$cf12bb(ncore+nbe+1,1,i),nbfcn, vifbb(ncore+nb+1,1,j),nbfcn,1.d0,
$wbb(1,1,i,j),nvbe)
call dgemm('n','t',nvbe,nbfan,dfnbasis,-0.25d0,
$ vifbb(ncore+nbe+1,1,j),nbfcn,cf12bb(ncore+nb+1,1,i),nbfcn,1.d0,
$wbb(1,1,i,j),nvbe)
call dgemm('n','t',nvbe,nbfan,dfnbasis,-0.25d0,
$cf12bb(ncore+nbe+1,1,j),nbfcn, vifbb(ncore+nb+1,1,i),nbfcn,1.d0,
$wbb(1,1,i,j),nvbe)
enddo
enddo
do j=1,nbe
do i=1,nal
call dgemm('n','t',nbfan,nvbe,dfnbasis,3d0/8d0,
$ vifaa(ncore+nb+1,1,i),nbfcn,cf12bb(ncore+nbe+1,1,j),nbfcn,0d0,
$wab(1,1,i,j),nbfan)
call dgemm('n','t',nbfan,nvbe,dfnbasis,3d0/8d0,
$cf12aa(ncore+nb+1,1,i),nbfcn, vifbb(ncore+nbe+1,1,j),nbfcn,1d0,
$wab(1,1,i,j),nbfan)
call dgemm('n','t',nbfan,nvbe,dfnbasis,1d0/8d0,
$ vifba(ncore+nb+1,1,j),nbfcn,cf12ab(ncore+nbe+1,1,i),nbfcn,1d0,
$wab(1,1,i,j),nbfan)
call dgemm('n','t',nbfan,nvbe,dfnbasis,1d0/8d0,
$cf12ba(ncore+nb+1,1,j),nbfcn, vifab(ncore+nbe+1,1,i),nbfcn,1d0,
$wab(1,1,i,j),nbfan)
enddo
enddo
do j=1,nbe
do i=1,nal
call dgemm('n','t',nval,nbfan,dfnbasis,3d0/8d0,
$ vifaa(ncore+nal+1,1,i),nbfcn,cf12bb(ncore+nb+1,1,j),nbfcn,0d0,
$wba(1,1,i,j),nval)
call dgemm('n','t',nval,nbfan,dfnbasis,3d0/8d0,
$cf12aa(ncore+nal+1,1,i),nbfcn, vifbb(ncore+nb+1,1,j),nbfcn,1d0,
$wba(1,1,i,j),nval)
call dgemm('n','t',nval,nbfan,dfnbasis,1d0/8d0,
$ vifba(ncore+nal+1,1,j),nbfcn,cf12ab(ncore+nb+1,1,i),nbfcn,1d0,
$wba(1,1,i,j),nval)
call dgemm('n','t',nval,nbfan,dfnbasis,1d0/8d0,
$cf12ba(ncore+nal+1,1,j),nbfcn, vifab(ncore+nb+1,1,i),nbfcn,1d0,
$wba(1,1,i,j),nval)
enddo
enddo
C
open(scrfile1,file='F12T3',form='unformatted')!,position='append')
c sum=0.d0
c a,a,a
do k=1,nal
do j=1,k-1
do i=1,j-1
do c=1,nval
do b=1,c-1
do a=1,b-1
tmp=0.d0
do d=1,nbfan
tmp=tmp+waa(a,d,j,k)*vintaa(d,i,b,c)
& -waa(b,d,j,k)*vintaa(d,i,a,c)
& -waa(c,d,j,k)*vintaa(d,i,b,a)
& -waa(a,d,i,k)*vintaa(d,j,b,c)
& +waa(b,d,i,k)*vintaa(d,j,a,c)
& +waa(c,d,i,k)*vintaa(d,j,b,a)
& -waa(a,d,j,i)*vintaa(d,k,b,c)
& +waa(b,d,j,i)*vintaa(d,k,a,c)
& +waa(c,d,j,i)*vintaa(d,k,b,a)
enddo
c sum=sum+tmp
c sum=sum+dabs(tmp)
if(dabs(tmp).gt.itol)write(scrfile1) tmp,a,b,c,i,j,k
enddo
enddo
enddo
enddo
enddo
enddo
write(scrfile1) 0.d0,0,0,0,0,0,0
c b,b,b
do k=1,nbe
do j=1,k-1
do i=1,j-1
do c=1,nvbe
do b=1,c-1
do a=1,b-1
tmp=0.d0
do d=1,nbfan
tmp=tmp+wbb(a,d,j,k)*vintbb(d,i,b,c)
& -wbb(b,d,j,k)*vintbb(d,i,a,c)
& -wbb(c,d,j,k)*vintbb(d,i,b,a)
& -wbb(a,d,i,k)*vintbb(d,j,b,c)
& +wbb(b,d,i,k)*vintbb(d,j,a,c)
& +wbb(c,d,i,k)*vintbb(d,j,b,a)
& -wbb(a,d,j,i)*vintbb(d,k,b,c)
& +wbb(b,d,j,i)*vintbb(d,k,a,c)
& +wbb(c,d,j,i)*vintbb(d,k,b,a)
enddo
c sum=sum+tmp
c sum=sum+dabs(tmp)
if(dabs(tmp).gt.itol)write(scrfile1) tmp,a,b,c,i,j,k
enddo
enddo
enddo
enddo
enddo
enddo
write(scrfile1) 0.d0,0,0,0,0,0,0
c a,a,b
do k=1,nbe
do j=1,nal
do i=1,j-1
do c=1,nvbe
do b=1,nval
do a=1,b-1
tmp=0.d0
do d=1,nbfan
tmp=tmp+wab(d,c,j,k)*vintaa(d,i,b,a)
& -wab(d,c,i,k)*vintaa(d,j,b,a)
& -waa(a,d,j,i)*vintab(d,k,b,c)
& +waa(b,d,j,i)*vintab(d,k,a,c)
& -wba(a,d,j,k)*vintba(i,d,b,c)
& +wba(b,d,j,k)*vintba(i,d,a,c)
& +wba(a,d,i,k)*vintba(j,d,b,c)
& -wba(b,d,i,k)*vintba(j,d,a,c)
enddo
c sum=sum+tmp
c sum=sum+dabs(tmp)
if(dabs(tmp).gt.itol)write(scrfile1) tmp,a,b,c,i,j,k
enddo
enddo
enddo
enddo
enddo
enddo
write(scrfile1) 0.d0,0,0,0,0,0,0
c a,b,b
do k=1,nbe
do j=1,k-1
do i=1,nal
do c=1,nvbe
do b=1,c-1
do a=1,nval
tmp=0.d0
do d=1,nbfan
tmp=tmp-wab(d,b,i,k)*vintab(d,j,a,c)
& +wab(d,c,i,k)*vintab(d,j,a,b)
& +wab(d,b,i,j)*vintab(d,k,a,c)
& -wab(d,c,i,j)*vintab(d,k,a,b)
& +wbb(b,d,j,k)*vintba(i,d,a,c)
& -wbb(c,d,j,k)*vintba(i,d,a,b)
& -wba(a,d,i,k)*vintbb(d,j,b,c)
& +wba(a,d,i,j)*vintbb(d,k,b,c)
enddo
c sum=sum+tmp
c sum=sum+dabs(tmp)
if(dabs(tmp).gt.itol)write(scrfile1) tmp,a,b,c,i,j,k
enddo
enddo
enddo
enddo
enddo
enddo
write(scrfile1) 0.d0,0,0,0,0,0,0
c write(6,"(' Checksum ',f15.10)") sum
c write(6,*)
close(scrfile1)
C
return
end
C
************************************************************************
subroutine write55(inp,nb,dfnbasis,nalc,nbec,ncore,nbfcn,itol,
$vifaa,vifbb,vifab,vifba,cf12aa,cf12bb,cf12ab,cf12ba,focka,fockb,
$scftype,orbsym,vint,fint,lhocc,nbasis)
************************************************************************
* Write fort.55 for mrcc
************************************************************************
implicit none
integer inp,nb,dfnbasis,nalc,nbec,ncore,nbfcn,orbsym(nb),p,q,r,s
integer i18(18),nbasis,i,nval,nvbe
real*8 vifaa(nbfcn,dfnbasis,nb),cf12aa(nbfcn,dfnbasis,nb),nuc,tmp
real*8 vifbb(nbfcn,dfnbasis,nb),cf12bb(nbfcn,dfnbasis,nb)
real*8 vifab(nbfcn,dfnbasis,nb),cf12ab(nbfcn,dfnbasis,nb)
real*8 vifba(nbfcn,dfnbasis,nb),cf12ba(nbfcn,dfnbasis,nb)
real*8 focka(nbfcn,nbfcn),fockb(nbfcn,nbfcn),itol,vint(nb,nb)
real*8 fint(nbfcn,nbfcn)
character(len=5) scftype
character(len=138) c138
logical lhocc
C Print Coulomb integrals (fort.55)
open(inp,status='unknown',file='fort.55')
rewind(inp)
read(inp,*) p,p,p,p
read(inp,*) (p,q=1,ncore),(orbsym(p),p=1,nb)
rewind(inp)
write(inp,*) nb,nalc+nbec,0
write(inp,"(10000000i2)") (orbsym(p),p=1,nb)
write(inp,"(i7)") 150000
C Alpha-alpha
do s=1,nb
do q=1,nb
call dgemm('n','t',nb,nb,dfnbasis,1.d0,vifaa(ncore+1,1,q),
$nbfcn,vifaa(ncore+1,1,s),nbfcn,0.d0,vint,nb)
do r=1,s
do p=1,min(q,(s-1)*s/2+r-(q-1)*q/2)
tmp=vint(p,r)
if(dabs(tmp).gt.itol) write(inp,"(e28.20,4i6)") tmp,p,q,r,s
enddo
enddo
enddo
enddo
if(scftype.ne.'rhf ') then
write(inp,"(e28.20,4i6)") 0.d0,0,0,0,0
C Beta-beta
do s=1,nb
do q=1,nb
call dgemm('n','t',nb,nb,dfnbasis,1.d0,vifbb(ncore+1,1,q),
$nbfcn,vifbb(ncore+1,1,s),nbfcn,0.d0,vint,nb)
do r=1,s
do p=1,min(q,(s-1)*s/2+r-(q-1)*q/2)
tmp=vint(p,r)
if(dabs(tmp).gt.itol) write(inp,"(e28.20,4i6)") tmp,p,q,r,s
enddo
enddo
enddo
enddo
write(inp,"(e28.20,4i6)") 0.d0,0,0,0,0
C Alpha-beta
do s=1,nb
do q=1,nb
call dgemm('n','t',nb,nb,dfnbasis,1.d0,vifaa(ncore+1,1,q),
$nbfcn,vifbb(ncore+1,1,s),nbfcn,0.d0,vint,nb)
do r=1,s
do p=1,q
tmp=vint(p,r)
if(dabs(tmp).gt.itol) write(inp,"(e28.20,4i6)") tmp,p,q,r,s
enddo
enddo
enddo
enddo
write(inp,"(e28.20,4i6)") 0.d0,0,0,0,0
endif
C Print Fock matrix
do p=1,nb
do q=1,p
tmp=focka(ncore+p,ncore+q)
if(dabs(tmp).gt.itol) write(inp,"(e28.20,4i6)") tmp,p,q,0,0
enddo
enddo
if(scftype.ne.'rhf ') then
write(inp,"(e28.20,4i6)") 0.d0,0,0,0,0
do p=1,nb
do q=1,p
tmp=fockb(ncore+p,ncore+q)
if(dabs(tmp).gt.itol) write(inp,"(e28.20,4i6)") tmp,p,q,0,0
enddo
enddo
write(inp,"(e28.20,4i6)") 0.d0,0,0,0,0
endif
C Print nuclear repulsion
call getvar('nuc ',nuc)
write(inp,"(e28.20,4i6)") nuc,0,0,0,0
close(inp)
C Write fort.56
if(.not.lhocc) return
open(inp,status='unknown',file='fort.56')
read(inp,*) i18,nuc,p,q,r,tmp
read(inp,"(a138)") c138
i18(17)=i18(1)
rewind(inp)
write(inp,"(i4,17i7,f9.5,i7,i9,i5,1pe10.3)") i18,nuc,p,q,r,tmp
write(inp,"(a138)") c138
write(inp,"(1000000i2)") (2,i=1,nbec),(1,i=1,nalc-nbec),
$(0,i=1,nb-nalc)
write(inp,"(1000000i2)") (1,i=1,nbasis),(0,i=nbasis+1,nb)
if(i18(5).eq.1) then
write(inp,"(1000000i2)") 0,2,(2,i=1,i18(1)-2)
else
write(inp,"(1000000i2)") 0,2,(2,i=1,i18(1)-3),0
endif
close(inp)
C Print F12 integrals (fort.54)
open(inp,status='unknown',file='fort.54')
rewind(inp)
write(inp,*) nb,nalc+nbec,0
write(inp,"(10000000i2)") (orbsym(p),p=1,nb)
write(inp,"(i7)") 150000
C Alpha-alpha
nval=nbfcn-ncore-nalc
do s=1,nbasis
do q=1,nbasis
call dgemm('n','t',nval,nval,dfnbasis,3d0/8d0,
$ vifaa(ncore+nalc+1,1,q),nbfcn,cf12aa(ncore+nalc+1,1,s),nbfcn,0d0,
$fint(ncore+nalc+1,ncore+nalc+1),nbfcn)
call dgemm('n','t',nval,nval,dfnbasis,3d0/8d0,
$cf12aa(ncore+nalc+1,1,q),nbfcn, vifaa(ncore+nalc+1,1,s),nbfcn,1d0,
$fint(ncore+nalc+1,ncore+nalc+1),nbfcn)
call dgemm('n','t',nval,nval,dfnbasis,1d0/8d0,
$ vifaa(ncore+nalc+1,1,s),nbfcn,cf12aa(ncore+nalc+1,1,q),nbfcn,1d0,
$fint(ncore+nalc+1,ncore+nalc+1),nbfcn)
call dgemm('n','t',nval,nval,dfnbasis,1d0/8d0,
$cf12aa(ncore+nalc+1,1,s),nbfcn, vifaa(ncore+nalc+1,1,q),nbfcn,1d0,
$fint(ncore+nalc+1,ncore+nalc+1),nbfcn)
do r=nalc+1,nb
do p=nalc+1,nb
tmp=fint(ncore+p,ncore+r)
if(dabs(tmp).gt.itol) write(inp,"(e28.20,4i6)") tmp,p,q,r,s
enddo
enddo
enddo
enddo
if(scftype.ne.'rhf ') then
write(inp,"(e28.20,4i6)") 0.d0,0,0,0,0
C Beta-beta
nvbe=nbfcn-ncore-nbec
do s=1,nbasis
do q=1,nbasis
call dgemm('n','t',nvbe,nvbe,dfnbasis,3d0/8d0,
$ vifbb(ncore+nbec+1,1,q),nbfcn,cf12bb(ncore+nbec+1,1,s),nbfcn,0d0,
$fint(ncore+nbec+1,ncore+nbec+1),nbfcn)
call dgemm('n','t',nvbe,nvbe,dfnbasis,3d0/8d0,
$cf12bb(ncore+nbec+1,1,q),nbfcn, vifbb(ncore+nbec+1,1,s),nbfcn,1d0,
$fint(ncore+nbec+1,ncore+nbec+1),nbfcn)
call dgemm('n','t',nvbe,nvbe,dfnbasis,1d0/8d0,
$ vifbb(ncore+nbec+1,1,s),nbfcn,cf12bb(ncore+nbec+1,1,q),nbfcn,1d0,
$fint(ncore+nbec+1,ncore+nbec+1),nbfcn)
call dgemm('n','t',nvbe,nvbe,dfnbasis,1d0/8d0,
$cf12bb(ncore+nbec+1,1,s),nbfcn, vifbb(ncore+nbec+1,1,q),nbfcn,1d0,
$fint(ncore+nbec+1,ncore+nbec+1),nbfcn)
do r=nbec+1,nb
do p=nbec+1,nb
tmp=fint(ncore+p,ncore+r)
if(dabs(tmp).gt.itol) write(inp,"(e28.20,4i6)") tmp,p,q,r,s
enddo
enddo
enddo
enddo
write(inp,"(e28.20,4i6)") 0.d0,0,0,0,0
C Alpha-beta
do s=1,nbasis
do q=1,nbasis
call dgemm('n','t',nval,nvbe,dfnbasis,3d0/8d0,
$ vifaa(ncore+nalc+1,1,q),nbfcn,cf12bb(ncore+nbec+1,1,s),nbfcn,0d0,
$fint(ncore+nalc+1,ncore+nbec+1),nbfcn)
call dgemm('n','t',nval,nvbe,dfnbasis,3d0/8d0,
$cf12aa(ncore+nalc+1,1,q),nbfcn, vifbb(ncore+nbec+1,1,s),nbfcn,1d0,
$fint(ncore+nalc+1,ncore+nbec+1),nbfcn)
call dgemm('n','t',nval,nvbe,dfnbasis,1d0/8d0,
$ vifba(ncore+nalc+1,1,s),nbfcn,cf12ab(ncore+nbec+1,1,q),nbfcn,1d0,
$fint(ncore+nalc+1,ncore+nbec+1),nbfcn)
call dgemm('n','t',nval,nvbe,dfnbasis,1d0/8d0,
$cf12ba(ncore+nalc+1,1,s),nbfcn, vifab(ncore+nbec+1,1,q),nbfcn,1d0,
$fint(ncore+nalc+1,ncore+nbec+1),nbfcn)
do r=nbec+1,nb
do p=nalc+1,nb
tmp=fint(ncore+p,ncore+r)
if(dabs(tmp).gt.itol) write(inp,"(e28.20,4i6)") tmp,p,q,r,s
enddo
enddo
enddo
enddo
write(inp,"(e28.20,4i6)") 0.d0,0,0,0,0
C
write(inp,"(e28.20,4i6)") 0.d0,0,0,0,0
endif
C
write(inp,"(e28.20,4i6)") 0.d0,0,0,0,0
write(inp,"(e28.20,4i6)") 0.d0,0,0,0,0
close(inp)
C
return
end
C
************************************************************************
subroutine symorb(nbfc,nbfcn,cmoa,cmob,bftrans,scrfile1,scftype,
$tmo,scr,nir,inp,nal,nbe,nirmax,smat,norm,cctol,proj)
************************************************************************
* Determine the symmetry of the orbitals
************************************************************************
implicit none
integer nbfc,nbfcn,scrfile1,nir,inp,i,j,nal,nbe,ir!,nfunc(8)
integer orbsyma(nbfcn),orbsymb(nbfcn),bftrans(nbfc,nir)
integer nirmax,chartab(nirmax,nirmax)
real*8 cmoa(nbfc,nbfcn),cmob(nbfc,nbfcn),smat(nbfc,nbfc),proj
real*8 tmo(nbfc,nbfcn),scr(nbfc,nbfcn),norm(nir,nbfcn),cctol
character(len=5) scftype
C Read overlap integrals
open(scrfile1,file='OEINT',form='unformatted')
rewind(scrfile1)
call roeint(scr,scr,smat,scrfile1,nbfc)
close(scrfile1)
C Read symmetry transformation matrix
open(unit=scrfile1,file='SYMTRA',status='old',form='unformatted')
read(scrfile1)
read(scrfile1) bftrans
close(scrfile1)
call getvar('chartab ',chartab)
c call getvar('nfunc ',nfunc)
c write(6,"(10i3)") nfunc(1:nir)
C Determine irrep
call getirrep(nbfc,nbfcn,bftrans,tmo,scr,nir,nirmax,chartab,
$smat,norm,cctol,proj,orbsyma,cmoa)
if(scftype.ne.'rhf ') then
call getirrep(nbfc,nbfcn,bftrans,tmo,scr,nir,nirmax,chartab,
$smat,norm,cctol,proj,orbsymb,cmob)
call dcopy(nbfc*nbfcn,cmob,1,tmo,1)
c write(6,"(10000000i2)") (orbsyma(i),i=1,nbfcn)
c write(6,"(10000000i2)") (orbsymb(i),i=1,nbfcn)
do i=1,nbfcn
j=1
do while(orbsymb(j).ne.orbsyma(i))
j=j+1
enddo
orbsymb(j)=0
call dcopy(nbfc,tmo(1,j),1,cmob(1,i),1)
enddo
endif
C Save irreps
open(inp,status='unknown',file='fort.55')
rewind(inp)
write(inp,*) nbfcn,nal+nbe,0,0
write(inp,"(10000000i2)") (orbsyma(i),i=1,nbfcn)
close(inp)
C
return
end
C
************************************************************************
subroutine getirrep(nbfc,nbfcn,bftrans,tmo,scr,nir,nirmax,chartab,
$smat,norm,cctol,proj,orbsyma,cmoa)
************************************************************************
* Determine the irrep of an oribtal
************************************************************************
implicit none
integer nbfc,nbfcn,nir,inp,i,j,n,ii,ir,bftrans(nbfc,nir),iop,mu
integer orbsyma(nbfcn)!,nfund(8)
integer nirmax,chartab(nirmax,nirmax)
real*8 cmoa(nbfc,nbfcn),cmob(nbfc,nbfcn),smat(nbfc,nbfc),ddot,tmp
real*8 tmo(nbfc,nbfcn,nir),scr(nbfc,nbfcn),norm(nir,nbfcn),cctol
real*8 proj(nbfc,nbfc)
C
c nfund=0
tmo=0.d0
do ir=1,nir
do i=1,nbfcn
do iop=1,nir
do mu=1,nbfc
tmo(iabs(bftrans(mu,iop)),i,ir)=tmo(iabs(bftrans(mu,iop)),i,ir)+
$dble(isign(1,bftrans(mu,iop))*chartab(iop,ir))*cmoa(mu,i)
enddo
enddo
enddo
call dscal(nbfc*nbfcn,1.d0/dble(nir),tmo(1,1,ir),1)
call dsymm('l','u',nbfc,nbfcn,1.d0,smat,nbfc,tmo(1,1,ir),nbfc,
$0.d0,scr,nbfc)
do i=1,nbfcn
norm(ir,i)=ddot(nbfc,tmo(1,i,ir),1,scr(1,i),1)
if(norm(ir,i).gt.cctol)
$ call dscal(nbfc,1.d0/dsqrt(norm(ir,i)),tmo(1,i,ir),1)
enddo
c write(6,*) 'ir=',ir
c write(6,"(7f9.5)") norm(ir,:)
enddo
C
proj=0.d0
do i=1,nbfc
proj(i,i)=1.d0
enddo
C
do i=1,nbfcn
4561 ir=1
do ii=2,nir
if(norm(ii,i).gt.norm(ir,i)) ir=ii
enddo
orbsyma(i)=ir
if(dabs(1.d0-norm(ir,i)).gt.cctol) then
c write(6,"(a4,i3,7f9.5)") 'norm',i,norm(ir,i)
norm(ir,i)=0.d0
call dgemv('n',nbfc,nbfc,1.d0,proj,nbfc,tmo(1,i,ir),1,0.d0,
$scr,1)
call dsymv('u',nbfc,1.d0,smat,nbfc,scr,1,0.d0,scr(1,2),1)
tmp=ddot(nbfc,scr,1,scr(1,2),1)
c write(6,"(a4,i3,7f14.5)") 'tmp ',i,tmp
c if(tmp.lt.cctol) write(6,*) tmp,ir
if(tmp.lt.cctol) goto 4561
call dscal(nbfc,1.d0/dsqrt(tmp),scr,1)
call dcopy(nbfc,scr,1,cmoa(1,i),1)
endif
call dsymv('u',nbfc,1.d0,smat,nbfc,cmoa(1,i),1,0.d0,scr,1)
call dger(nbfc,nbfc,-1.d0,cmoa(1,i),1,scr,1,proj,nbfc)
c nfund(ir)=nfund(ir)+1
enddo
c write(6,"(10i3)") nfund(1:nir)
c call dsymm('l','u',nbfc,nbfcn,1.d0,smat,nbfc,cmoa,nbfc,0.d0,scr,
c $nbfc)
c call dgemm('t','n',nbfcn,nbfcn,nbfc,1.d0,cmoa,nbfc,
c $scr,nbfc,0.d0,tmo,nbfc)
c do i=1,nbfcn
c do j=1,nbfcn
c if(dabs(tmo(i,j,1)).gt.cctol) write(6,*) i,j,tmo(i,j,1)
c enddo
c enddo
C
return
end
C