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

1061 lines
42 KiB
Fortran

C TODO: T3 memory, T3 closed shell, xmrcc
************************************************************************
subroutine mrcc_f12(nb,nal,nbe,scrfile1,ecc,nconf,trec,
$nstr,nmax,nir,isympv,isympo,nactm,strec,dcore,imem,iimem,tfile,
$scrfile7,job,op,mrop,nactv,nacto,nactva,nactvb,nactoa,nactob,
$nirmax,csympair,icore,isym)
************************************************************************
* Initialize the calculation of F12 contributions
************************************************************************
implicit none
integer nb,nal,nbe,nval,nvbe,scrfile1,nconf,trec,nstr,nmax,nir
integer isympv,isympo,nactm,strec,imem,ivv,icabij
integer tfile,scrfile7,icabijaa,icabijbb,icabijab,ivaa,ivbb,ivab
integer iuaia,iuaib,ivpia,ivpib,icaia,icaib,itaone,itbone,itatwo
integer itbtwo,itmtwo,inewta1mat,inewtb1mat,inewta2mat,inewtb2mat
integer inewtab2mat,intadd,iimem,itaaa,itbbb,itaab,itabb
integer op,mrop,nactv,nacto,nactva,nactvb,nactoa
integer nactob,nirmax,csympair,icore,isym
real*8 ecc,dcore(*)
character(len=4) job
character(len=5) scftype
logical lt3
C
c lt3=job.eq.'resi' !xmrcc: comment in lf12; mp2f12: t3f12
lt3=.false.
call getkey('scftype',7,scftype,5)
nval=nb-nal
nvbe=nb-nbe
C
ivpia=imem
imem=imem+(nal*nal)
iuaia=imem
imem=imem+(nval*nal)
icaia=imem
imem=imem+(nval*nal)
inewta1mat=imem
imem=imem+(nval*nal)
inewta2mat=imem
imem=imem+(nval*nval*nal*nal)
inewtab2mat=imem
imem=imem+(nval*nvbe*nal*nbe)
itaone=imem
imem=imem+(nval*nal)
itatwo=imem
imem=imem+(nval*nval*nal*nal)
itmtwo=imem
imem=imem+(nval*nvbe*nal*nbe)
if(lt3) then
itaaa=imem
imem=imem+(nval*nval*nval*nal*nal*nal)
itbbb=imem
imem=imem+(nvbe*nvbe*nvbe*nbe*nbe*nbe)
itaab=imem
imem=imem+(nval*nval*nvbe*nal*nal*nbe)
itabb=imem
imem=imem+(nval*nvbe*nvbe*nal*nbe*nbe)
else
itaaa=imem
itbbb=imem
itaab=imem
itabb=imem
endif
if(scftype.eq.'rhf ') then
ivv=imem
imem=imem+nb*nb*(nal+1)*nal/2
icabij=imem
imem=imem+nval*nval*(nal+1)*nal/2
ivaa=imem
ivab=imem
ivbb=imem
icabijaa=imem
icabijbb=imem
icabijab=imem
ivpib=ivpia
iuaib=iuaia
icaib=icaia
inewtb1mat=inewta1mat
inewtb2mat=inewta2mat
itbone=itaone
itbtwo=itatwo
else
ivv=imem
icabij=imem
ivpib=imem
imem=imem+(nbe*nbe)
iuaib=imem
imem=imem+(nvbe*nbe)
icaib=imem
imem=imem+(nvbe*nbe)
ivaa=imem
imem=imem+((nb-1)*nb/2*(nal-1)*nal/2)
ivab=imem
imem=imem+(nb*nb*nal*nbe)
ivbb=imem
imem=imem+((nb-1)*nb/2*(nbe-1)*nbe/2)
icabijaa=imem
imem=imem+((nval-1)*nval/2*(nal-1)*nal/2)
icabijbb=imem
imem=imem+((nvbe-1)*nvbe/2*(nbe-1)*nbe/2)
icabijab=imem
imem=imem+(nval*nvbe*nal*nbe)
inewtb1mat=imem
imem=imem+(nvbe*nbe)
inewtb2mat=imem
imem=imem+(nvbe*nvbe*nbe*nbe)
itbone=imem
imem=imem+(nvbe*nbe)
itbtwo=imem
imem=imem+(nvbe*nvbe*nbe*nbe)
endif
C
iimem=intadd(imem)
call f12cont(nb,nal,nbe,nval,nvbe,scrfile1,ecc,nconf,trec,nstr,
$nmax,dcore(imem),nir,isympv,isympo,nactm,strec,imem,tfile,
$scrfile7,dcore(icabijaa),dcore(icabijbb),dcore(icabijab),
$dcore(ivaa),dcore(ivbb),dcore(ivab),dcore(iuaia),dcore(iuaib),
$dcore(ivpia),dcore(ivpib),dcore(icaia),dcore(icaib),dcore(itaone),
$dcore(itbone),dcore(itatwo),dcore(itbtwo),dcore(itmtwo),
$dcore(inewta1mat),dcore(inewtb1mat),dcore(inewta2mat),
$dcore(inewtb2mat),dcore(inewtab2mat),job,scftype,dcore(ivv),
$dcore(icabij),dcore(itaaa),dcore(itbbb),dcore(itaab),dcore(itabb),
$lt3,op,mrop,nactv,nacto,nactva,nactvb,nactoa,nactob,nirmax,
$csympair,icore,isym)
imem=ivpia
iimem=intadd(imem)
C
return
end
C
************************************************************************
subroutine f12cont(nb,nal,nbe,nval,nvbe,scrfile1,ecc,nconf,trec,
$nstr,nmax,v,nir,isympv,isympo,nactm,strec,imem,tfile,scrfile7,
$cabijaa,cabijbb,cabijab,vaa,vbb,vab,uaia,uaib,vpia,vpib,caia,caib,
$taone,tbone,tatwo,tbtwo,tmtwo,newta1mat,newtb1mat,newta2mat,
$newtb2mat,newtab2mat,job,scftype,vv,cabij,taaa,tbbb,taab,tabb,lt3,
$op,mrop,nactv,nacto,nactva,nactvb,nactoa,nactob,nirmax,
$csympair,icore,isym)
************************************************************************
* Add F12 contribution to residual and calculate F12 contrib. to energy
************************************************************************
implicit none
integer nb,nal,nbe,nval,nvbe,scrfile1,i,j,k,l,ij,kl,a,b,isym
integer nconf,trec,nstr,nmax,nir,isympv,isympo,nactm,imem,icore
integer strec,tfile,scrfile7,op,mrop,nactv,nacto,nactva
integer nactvb,nactoa,nactob,nirmax,csympair
real*8 vaa((nb-1)*nb/2,(nal-1)*nal/2),vab(nb,nb,nal,nbe),v(*)
real*8 vbb((nb-1)*nb/2,(nbe-1)*nbe/2),ecabs,emp2f12,ecc
real*8 vpia(nal,nal),vpib(nbe,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),newtab2mat(nval,nvbe,nal,nbe)
real*8 newta1mat(nval,nal),newtb1mat(nvbe,nbe)
real*8 newta2mat(nval,nval,nal,nal),newtb2mat(nvbe,nvbe,nbe,nbe)
real*8 taone(nval,nal),tbone(nvbe,nbe),tmtwo(nval,nvbe,nal,nbe)
real*8 tatwo(nval,nval,nal,nal),tbtwo(nvbe,nvbe,nbe,nbe),tmp,ef12
real*8 vv(nb,nb,(nal+1)*nal/2),cabij(nval,nval,(nal+1)*nal/2)
real*8 taaa(nval,nval,nval,nal,nal,nal)
real*8 tbbb(nvbe,nvbe,nvbe,nbe,nbe,nbe)
real*8 taab(nval,nval,nvbe,nal,nal,nbe)
real*8 tabb(nval,nvbe,nvbe,nal,nbe,nbe)
character(len=4) job
character(len=5) scftype
logical lt3
C
open(scrfile7,file='F12INTE',form='unformatted')
read(scrfile7) ecabs,emp2f12!,ecoup
read(scrfile7) !tfact,tscalea,tscaleb
taone=0.d0
tatwo=0.d0
tmtwo=0.d0
if(scftype.eq.'rhf ') then
if(job.eq.'ener') then
read(scrfile7)
$((((cabij(a,b,j*(j-1)/2+i),a=1,nval),b=1,nval),i=1,j-1),
$((cabij(a,b,j*(j-1)/2+j),a=1,nval),b=1,nval),j=1,nal)
read(scrfile7)
$((((vv(nal+a,nal+b,j*(j-1)/2+i),a=1,nval),b=1,nval),i=1,j-1),
$((vv(nal+a,nal+b,j*(j-1)/2+j),a=1,nval),b=1,nval),j=1,nal)
read(scrfile7) uaia
endif
else
tbone=0.d0
tbtwo=0.d0
if(nal.gt.1.and.nval.gt.1) read(scrfile7) cabijaa
if(nbe.gt.1.and.nvbe.gt.1) read(scrfile7) cabijbb
if(nal.gt.0.and.nval.gt.0.and.nbe.gt.0.and.nvbe.gt.0)
$read(scrfile7) cabijab
if(nal.gt.1.and.nval.gt.1) read(scrfile7)
$(((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) read(scrfile7)
$(((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)
$read(scrfile7) (((vab(nal+1:nb,nbe+b,i,j),b=1,nvbe),i=1,nal),
$j=1,nbe)
if(nal.gt.0.and.nval.gt.0) read(scrfile7) uaia
if(nbe.gt.0.and.nvbe.gt.0) read(scrfile7) uaib
endif
close(scrfile7)
C Convert cluster amplitudes
call tconv(nconf,trec,nstr,nmax,v,nir,isympv,isympo,nactm,strec,
$taone,tbone,tatwo,tbtwo,tmtwo,taaa,tbbb,taab,tabb,.true.,tfile,
$.false.,nval,nvbe,nal,nbe,nir,op,mrop,nactv,nacto,nactva,nactvb,
$nactoa,nactob,nirmax,csympair,imem,icore,isym)
C F12 contribution to the residual
if(job.eq.'resi') then
newta1mat=0.d0
newta2mat=0.d0
newtab2mat=0.d0
open(scrfile7,file='F12INT1',form='unformatted')
if(scftype.eq.'rhf ') then
read(scrfile7) 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
read(scrfile7) vv(1:nal,1:nal,ij)
enddo
enddo
read(scrfile7) vpia
do j=1,nal
read(scrfile7) (vv(1:nal,nal+1:nb,j*(j-1)/2+i),i=1,j),
$(((vv(l,k,i*(i-1)/2+j),k=1,nal),l=nal+1,nb),i=j+1,nal)
vv(nal+1:nb,1:nal,j*(j-1)/2+j)=
$ transpose(vv(1:nal,nal+1:nb,j*(j-1)/2+j))
enddo
else
newtb1mat=0.d0
newtb2mat=0.d0
if(nal.gt.0.and.nval.gt.0) then
read(scrfile7) vpia
read(scrfile7) caia
endif
if(nbe.gt.0.and.nvbe.gt.0) then
read(scrfile7) vpib
read(scrfile7) caib
endif
close(scrfile7)
open(scrfile7,file='F12INT2',form='unformatted')
if(nal.gt.1.and.nval.gt.0) read(scrfile7)
$(((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) read(scrfile7)
$((vaa(kl,ij),kl=1,(nal-1)*nal/2),ij=1,(nal-1)*nal/2)
if(nbe.gt.1.and.nvbe.gt.0) read(scrfile7)
$(((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) read(scrfile7)
$((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) read(scrfile7)
$((((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) read(scrfile7)
$(((vab(1:nal,nbe+b,i,j),b=1,nvbe),i=1,nal),j=1,nbe)
if(nal.gt.0.and.nval.gt.0.and.nbe.gt.0.and.nvbe.gt.0)
$read(scrfile7) (((vab(1:nal,b,i,j),b=1,nbe),i=1,nal),j=1,nbe)
endif
close(scrfile7)
C Convert residual
call tconv(nconf,trec,nstr,nmax,v,nir,isympv,isympo,nactm,strec,
$newta1mat,newtb1mat,newta2mat,newtb2mat,newtab2mat,taaa,tbbb,taab,
$tabb,.true.,scrfile1,lt3,nval,nvbe,nal,nbe,nir,op,mrop,nactv,
$nacto,nactva,nactvb,nactoa,nactob,nirmax,csympair,imem,icore,
$isym)
C T1 equations
C aa
do i=1,nal
do a=1,nval
newta1mat(a,i)=newta1mat(a,i)+caia(a,i)
tmp=0.d0
do k=1,nal
tmp=tmp+vpia(k,i)*taone(a,k)
enddo
newta1mat(a,i)=newta1mat(a,i)-tmp
enddo
enddo
C bb
if(scftype.ne.'rhf ') then
do i=1,nbe
do a=1,nvbe
newtb1mat(a,i)=newtb1mat(a,i)+caib(a,i)
tmp=0.d0
do k=1,nbe
tmp=tmp+vpib(k,i)*tbone(a,k)
enddo
newtb1mat(a,i)=newtb1mat(a,i)-tmp
enddo
enddo
endif
C T2 equations
C abab
if(scftype.eq.'rhf ') then
do j=1,nbe
do i=1,j
do b=1,nvbe
do a=1,nval
tmp=cabij(a,b,j*(j-1)/2+i)
do k=1,nal
tmp=tmp-vv(k,nal+b,j*(j-1)/2+i)*taone(a,k)
$ -vv(nal+a,k,j*(j-1)/2+i)*taone(b,k)
$ -vpia(k,i)*tmtwo(b,a,j,k)
$ -vpia(k,j)*tmtwo(a,b,i,k)
do l=1,nal
tmp=tmp+vv(l,k,j*(j-1)/2+i)*
$(tmtwo(a,b,k,l)+taone(a,k)*taone(b,l))
enddo
enddo
newtab2mat(a,b,i,j)=newtab2mat(a,b,i,j)+tmp
enddo
enddo
enddo
enddo
do j=1,nbe
do i=j+1,nal
do b=1,nvbe
do a=1,nval
newtab2mat(a,b,i,j)=newtab2mat(b,a,j,i)
enddo
enddo
enddo
enddo
do j=1,nal
do i=1,j-1
do b=1,nval
do a=1,b-1
tmp=newtab2mat(a,b,i,j)-newtab2mat(b,a,i,j)
newta2mat(a,b,i,j)= tmp
newta2mat(b,a,i,j)=-tmp
newta2mat(a,b,j,i)=-tmp
newta2mat(b,a,j,i)= tmp
enddo
enddo
enddo
enddo
else
C abab
do j=1,nbe
do i=1,nal
do b=1,nvbe
do a=1,nval
tmp=cabijab(a,b,i,j)
do k=1,nal
tmp=tmp-vab(k,nbe+b,i,j)*taone(a,k)
$ -vpia(k,i)*tmtwo(a,b,k,j)
do l=1,nbe
tmp=tmp+vab(k,l,i,j)*
$(tmtwo(a,b,k,l)+taone(a,k)*tbone(b,l))
enddo
enddo
do k=1,nbe
tmp=tmp-vab(nal+a,k,i,j)*tbone(b,k)
$ -vpib(k,j)*tmtwo(a,b,i,k)
enddo
newtab2mat(a,b,i,j)=newtab2mat(a,b,i,j)+tmp
enddo
enddo
enddo
enddo
C aaaa
do j=1,nal
do i=1,j-1
do b=1,nval
do a=1,b-1
tmp=cabijaa((b-1)*(b-2)/2+a,(j-1)*(j-2)/2+i)
do k=1,nal
tmp=tmp+vaa((nal+a-1)*(nal+a-2)/2+k,(j-1)*(j-2)/2+i)*taone(b,k)
$ -vaa((nal+b-1)*(nal+b-2)/2+k,(j-1)*(j-2)/2+i)*taone(a,k)
$ -vpia(k,j)*tatwo(a,b,i,k)
$ +vpia(k,i)*tatwo(a,b,j,k)
do l=1,k-1
tmp=tmp-vaa((k-1)*(k-2)/2+l,(j-1)*(j-2)/2+i)*
$(tatwo(a,b,k,l)+taone(a,k)*taone(b,l)-taone(a,l)*taone(b,k))
enddo
enddo
newta2mat(a,b,i,j)=newta2mat(a,b,i,j)+tmp
newta2mat(b,a,i,j)=newta2mat(b,a,i,j)-tmp
newta2mat(a,b,j,i)=newta2mat(a,b,j,i)-tmp
newta2mat(b,a,j,i)=newta2mat(b,a,j,i)+tmp
enddo
enddo
enddo
enddo
C bbbb
do j=1,nbe
do i=1,j-1
do b=1,nvbe
do a=1,b-1
tmp=cabijbb((b-1)*(b-2)/2+a,(j-1)*(j-2)/2+i)
do k=1,nbe
tmp=tmp+vbb((nbe+a-1)*(nbe+a-2)/2+k,(j-1)*(j-2)/2+i)*tbone(b,k)
$ -vbb((nbe+b-1)*(nbe+b-2)/2+k,(j-1)*(j-2)/2+i)*tbone(a,k)
$ -vpib(k,j)*tbtwo(a,b,i,k)
$ +vpib(k,i)*tbtwo(a,b,j,k)
do l=1,k-1
tmp=tmp-vbb((k-1)*(k-2)/2+l,(j-1)*(j-2)/2+i)*
$(tbtwo(a,b,k,l)+tbone(a,k)*tbone(b,l)-tbone(a,l)*tbone(b,k))
enddo
enddo
newtb2mat(a,b,i,j)=newtb2mat(a,b,i,j)+tmp
newtb2mat(b,a,i,j)=newtb2mat(b,a,i,j)-tmp
newtb2mat(a,b,j,i)=newtb2mat(a,b,j,i)-tmp
newtb2mat(b,a,j,i)=newtb2mat(b,a,j,i)+tmp
enddo
enddo
enddo
enddo
endif
C T3 equations
if(lt3) then
call t3read(nal,nbe,nval,nvbe,taaa,tbbb,taab,tabb,scrfile7)
endif
C Convert residual back
call tconv(nconf,trec,nstr,nmax,v,nir,isympv,isympo,nactm,strec,
$newta1mat,newtb1mat,newta2mat,newtb2mat,newtab2mat,taaa,tbbb,taab,
$tabb,.false.,scrfile1,lt3,nval,nvbe,nal,nbe,nir,op,mrop,nactv,
$nacto,nactva,nactvb,nactoa,nactob,nirmax,csympair,imem,icore,
$isym)
else if(job.eq.'ener') then
C F12 contribution to the energy
ef12=0.d0
if(scftype.eq.'rhf ') then
do j=1,nal
do a=1,nval
ef12=ef12+uaia(a,j)*taone(a,j)
do b=1,nval
ef12=ef12+cabij(a,b,j*(j-1)/2+j)*tmtwo(a,b,j,j)
enddo
do i=1,j-1
do b=1,nval
ef12=ef12+cabij(a,b,j*(j-1)/2+i)*tmtwo(a,b,i,j)
enddo
enddo
enddo
enddo
do j=1,nal
do b=1,nval
do a=1,nval
ef12=ef12+vv(nal+a,nal+b,j*(j-1)/2+j)*taone(a,j)*taone(b,j)
enddo
enddo
do i=1,j-1
do b=1,nval
do a=1,nval
ef12=ef12+vv(nal+a,nal+b,j*(j-1)/2+i)*taone(a,i)*taone(b,j)
enddo
enddo
enddo
enddo
ef12=2.d0*ef12
else
do i=1,nal
do a=1,nval
ef12=ef12+uaia(a,i)*taone(a,i)
do j=1,i-1
do b=1,a-1
ef12=ef12+cabijaa((a-1)*(a-2)/2+b,(i-1)*(i-2)/2+j)*tatwo(a,b,i,j)
enddo
enddo
enddo
enddo
do i=1,nbe
do a=1,nvbe
ef12=ef12+uaib(a,i)*tbone(a,i)
do j=1,i-1
do b=1,a-1
ef12=ef12+cabijbb((a-1)*(a-2)/2+b,(i-1)*(i-2)/2+j)*tbtwo(a,b,i,j)
enddo
enddo
enddo
enddo
do i=1,nal
do a=1,nval
do j=1,nbe
do b=1,nvbe
ef12=ef12+cabijab(a,b,i,j)*tmtwo(a,b,i,j)
enddo
enddo
enddo
enddo
do j=1,nal
do i=1,j-1
do b=1,nval
do a=1,b-1
ef12=ef12+vaa((nal+b-1)*(nal+b-2)/2+nal+a,(j-1)*(j-2)/2+i)*
$(taone(a,i)*taone(b,j)-taone(b,i)*taone(a,j))
enddo
enddo
enddo
enddo
do j=1,nbe
do i=1,j-1
do b=1,nvbe
do a=1,b-1
ef12=ef12+vbb((nbe+b-1)*(nbe+b-2)/2+nbe+a,(j-1)*(j-2)/2+i)*
$(tbone(a,i)*tbone(b,j)-tbone(b,i)*tbone(a,j))
enddo
enddo
enddo
enddo
do i=1,nal
do a=1,nval
do j=1,nbe
do b=1,nvbe
ef12=ef12+vab(nal+a,nbe+b,i,j)*taone(a,i)*tbone(b,j)
enddo
enddo
enddo
enddo
endif
ecc=ecc+emp2f12+ef12
endif
C
return
end
C
************************************************************************
subroutine tconv(nconf,trec,nstr,nmax,v,nnir,isympv,isympo,nactm,
$irec,ta,tb,taa,tbb,tab,taaa,tbbb,taab,tabb,ltt,ifile,lt3,nvirtal,
$nvirtbe,nal,nbe,nir,op,mrop,nactv,nacto,nactva,nactvb,nactoa,
$nactob,nirmax,csympair,imem,icore,isym)
************************************************************************
* Converts string-based cluster amplitudes to many-index quantities
************************************************************************
implicit none
integer nnir,nmax,nactm,i,j,k,nex,ii,jj,nn,ifile,mex,ir,op,mrop,i1
integer nconf(0:nactm,0:nactm,0:nactm,0:nactm,0:nmax,0:nmax),imem
integer n1,n2,n3,n4,ire,ile,l,nvirtal,nvirtbe,nal,nbe,nir,icore(*)
integer trec(0:nactm,0:nactm,0:nactm,0:nactm,0:nmax,0:nmax),isym
integer iactv,iacto,iactva,iactoa,iactvb,iactob,nactv,nacto
integer nstr(nnir,0:nactm,0:nmax,4),nactva,nactvb,nactoa,nactob
integer isympv(0:nnir,nnir,0:nactm,0:nactm,0:nmax,0:nmax,2)
integer isympo(0:nnir,nnir,0:nactm,0:nactm,0:nmax,0:nmax,2)
integer isymv,isymo,isymva,isymvb,isymoa,isymob,irv,iro,dbladd
integer nn1,nn2,nn3,nn4,irec(nnir,0:nactm,0:nmax,4)
integer nampvirtal,nampvirtbe,nampoccal,nampoccbe,intadd
integer nirmax,csympair(nirmax,nirmax,2)
real*8 v(*),ta(nvirtal,nal),tb(nvirtbe,nbe)
real*8 taa(nvirtal,nvirtal,nal,nal),tbb(nvirtbe,nvirtbe,nbe,nbe)
real*8 tab(nvirtal,nvirtbe,nal,nbe)
real*8 taaa(nvirtal,nvirtal,nvirtal,nal,nal,nal)
real*8 tbbb(nvirtbe,nvirtbe,nvirtbe,nbe,nbe,nbe)
real*8 taab(nvirtal,nvirtal,nvirtbe,nal,nal,nbe)
real*8 tabb(nvirtal,nvirtbe,nvirtbe,nal,nbe,nbe)
logical ltt,lt3
C Loop over excitations
mex=2
if(lt3) mex=3
do nex=1,min(mex,op)
do iactv=max(0,nex-mrop),min(nactv,nex)
do iacto=max(0,nex-mrop),min(nacto,nex)
do i1=0,nex
do iactva=max(0,iactv-nactvb),min(nactva,i1,iactv)
iactvb=iactv-iactva
do iactoa=max(0,iacto-nactob),min(nactoa,i1,iacto)
iactob=iacto-iactoa
nn=nconf(iactva,iactvb,iactoa,iactob,i1,nex)
if(nn.gt.0) then
nampvirtal=i1
nampvirtbe=nex-nampvirtal
nampoccal=i1
nampoccbe=nex-nampoccal
if(ltt) call getlst(ifile,
$trec(iactva,iactvb,iactoa,iactob,i1,nex),v(1),nn)
ii=1
do ir=1,nir
isymv=csympair(isym,ir,1)
isymo=csympair(isym,ir,2)
do irv=1,isympv(0,isymv,iactva,iactvb,nampvirtal,
$nampvirtbe,1)
isymva=isympv(irv,isymv,iactva,iactvb,nampvirtal,nampvirtbe,1)
isymvb=isympv(irv,isymv,iactva,iactvb,nampvirtal,nampvirtbe,2)
n1=nstr(isymva,iactva,nampvirtal,1)
n2=nstr(isymvb,iactvb,nampvirtbe,2)
k=n1*n2
do iro=1,isympo(0,isymo,iactoa,iactob,nampoccal,
$nampoccbe,1)
isymoa=isympo(iro,isymo,iactoa,iactob,nampoccal,nampoccbe,1)
isymob=isympo(iro,isymo,iactoa,iactob,nampoccal,nampoccbe,2)
n3=nstr(isymoa,iactoa,nampoccal,3)
n4=nstr(isymob,iactob,nampoccbe,4)
call strread(nn1,nn2,nn3,nn4,n1,n2,n3,n4,nampvirtal,nampvirtbe,
$nampoccal,nampoccbe,isymva,isymvb,isymoa,isymob,iactva,iactvb,
$iactoa,iactob,nnir,nactm,nmax,irec,jj,imem+nn)
call tconv1(nampvirtal,nampvirtbe,nampoccal,nampoccbe,n1,n2,n3,n4,
$icore(nn1),icore(nn2),icore(nn3),icore(nn4),v(ii),ta,tb,taa,
$tbb,tab,taaa,tbbb,taab,tabb,nvirtal,nvirtbe,nal,nbe,nex,ltt)
ii=ii+n3*n4*k
enddo
enddo
enddo
if(.not.ltt) call putlst(ifile,
$trec(iactva,iactvb,iactoa,iactob,i1,nex),v(1),nn)
endif
enddo
enddo
enddo
enddo
enddo
enddo
C
return
end
C
************************************************************************
subroutine tconv1(nva,nvb,noa,nob,nvstral,nvstrbe,nostral,nostrbe,
$istrva,istrvb,istroa,istrob,tt,ta,tb,taa,tbb,tab,taaa,tbbb,
$taab,tabb,nvirtal,nvirtbe,nal,nbe,nex,ltt)
************************************************************************
* Converts string-based cluster amplitudes to many-index quantities
************************************************************************
implicit none
integer nva,nvb,noa,nob,nvstral,nvstrbe,nostral,nostrbe,ii
integer i,ivstral,ivstrbe,iostral,iostrbe
integer istrva(nva,nvstral),istrvb(nvb,nvstrbe)
integer istroa(noa,nostral),istrob(nob,nostrbe),p,q,ivb(nvb)
integer nvirtal,nvirtbe,nal,nbe,nex,ioa(noa),iob(nob),iva(nva)
real*8 tt(*),ta(nvirtal,nal),tb(nvirtbe,nbe),tmp
real*8 taa(nvirtal,nvirtal,nal,nal),tbb(nvirtbe,nvirtbe,nbe,nbe)
real*8 tab(nvirtal,nvirtbe,nal,nbe)
real*8 taaa(nvirtal,nvirtal,nvirtal,nal,nal,nal)
real*8 tbbb(nvirtbe,nvirtbe,nvirtbe,nbe,nbe,nbe)
real*8 taab(nvirtal,nvirtal,nvirtbe,nal,nal,nbe)
real*8 tabb(nvirtal,nvirtbe,nvirtbe,nal,nbe,nbe)
logical ltt
C
ii=0
do iostrbe=1,nostrbe
do iostral=1,nostral
do ivstrbe=1,nvstrbe
do ivstral=1,nvstral
ii=ii+1
do i=1,noa
ioa(i)=istroa(i,iostral)
enddo
do i=1,nob
iob(i)=istrob(i,iostrbe)
enddo
do i=1,nva
iva(i)=istrva(i,ivstral)
enddo
do i=1,nvb
ivb(i)=istrvb(i,ivstrbe)
enddo
if(ltt) then
tmp=tt(ii)
if(nex.eq.1) then
if(nva.eq.1) then
ta(iva(1),ioa(1))=tmp
else
tb(ivb(1),iob(1))=tmp
endif
else if(nex.eq.2) then
if(nva.eq.2) then
taa(iva(1),iva(2),ioa(1),ioa(2))= tmp
taa(iva(2),iva(1),ioa(1),ioa(2))=-tmp
taa(iva(1),iva(2),ioa(2),ioa(1))=-tmp
taa(iva(2),iva(1),ioa(2),ioa(1))= tmp
else if(nva.eq.0) then
tbb(ivb(1),ivb(2),iob(1),iob(2))= tmp
tbb(ivb(2),ivb(1),iob(1),iob(2))=-tmp
tbb(ivb(1),ivb(2),iob(2),iob(1))=-tmp
tbb(ivb(2),ivb(1),iob(2),iob(1))= tmp
else
tab(iva(1),ivb(1),ioa(1),iob(1))= tmp
endif
else
if(nva.eq.3) then !aaa
taaa(iva(1),iva(2),iva(3),ioa(1),ioa(2),ioa(3))= tmp
taaa(iva(1),iva(3),iva(2),ioa(1),ioa(2),ioa(3))=-tmp
taaa(iva(3),iva(1),iva(2),ioa(1),ioa(2),ioa(3))= tmp
taaa(iva(3),iva(2),iva(1),ioa(1),ioa(2),ioa(3))=-tmp
taaa(iva(2),iva(3),iva(1),ioa(1),ioa(2),ioa(3))= tmp
taaa(iva(2),iva(1),iva(3),ioa(1),ioa(2),ioa(3))=-tmp
taaa(iva(1),iva(2),iva(3),ioa(1),ioa(3),ioa(2))=-tmp
taaa(iva(1),iva(3),iva(2),ioa(1),ioa(3),ioa(2))= tmp
taaa(iva(3),iva(1),iva(2),ioa(1),ioa(3),ioa(2))=-tmp
taaa(iva(3),iva(2),iva(1),ioa(1),ioa(3),ioa(2))= tmp
taaa(iva(2),iva(3),iva(1),ioa(1),ioa(3),ioa(2))=-tmp
taaa(iva(2),iva(1),iva(3),ioa(1),ioa(3),ioa(2))= tmp
taaa(iva(1),iva(2),iva(3),ioa(3),ioa(1),ioa(2))= tmp
taaa(iva(1),iva(3),iva(2),ioa(3),ioa(1),ioa(2))=-tmp
taaa(iva(3),iva(1),iva(2),ioa(3),ioa(1),ioa(2))= tmp
taaa(iva(3),iva(2),iva(1),ioa(3),ioa(1),ioa(2))=-tmp
taaa(iva(2),iva(3),iva(1),ioa(3),ioa(1),ioa(2))= tmp
taaa(iva(2),iva(1),iva(3),ioa(3),ioa(1),ioa(2))=-tmp
taaa(iva(1),iva(2),iva(3),ioa(3),ioa(2),ioa(1))=-tmp
taaa(iva(1),iva(3),iva(2),ioa(3),ioa(2),ioa(1))= tmp
taaa(iva(3),iva(1),iva(2),ioa(3),ioa(2),ioa(1))=-tmp
taaa(iva(3),iva(2),iva(1),ioa(3),ioa(2),ioa(1))= tmp
taaa(iva(2),iva(3),iva(1),ioa(3),ioa(2),ioa(1))=-tmp
taaa(iva(2),iva(1),iva(3),ioa(3),ioa(2),ioa(1))= tmp
taaa(iva(1),iva(2),iva(3),ioa(2),ioa(3),ioa(1))= tmp
taaa(iva(1),iva(3),iva(2),ioa(2),ioa(3),ioa(1))=-tmp
taaa(iva(3),iva(1),iva(2),ioa(2),ioa(3),ioa(1))= tmp
taaa(iva(3),iva(2),iva(1),ioa(2),ioa(3),ioa(1))=-tmp
taaa(iva(2),iva(3),iva(1),ioa(2),ioa(3),ioa(1))= tmp
taaa(iva(2),iva(1),iva(3),ioa(2),ioa(3),ioa(1))=-tmp
taaa(iva(1),iva(2),iva(3),ioa(2),ioa(1),ioa(3))=-tmp
taaa(iva(1),iva(3),iva(2),ioa(2),ioa(1),ioa(3))= tmp
taaa(iva(3),iva(1),iva(2),ioa(2),ioa(1),ioa(3))=-tmp
taaa(iva(3),iva(2),iva(1),ioa(2),ioa(1),ioa(3))= tmp
taaa(iva(2),iva(3),iva(1),ioa(2),ioa(1),ioa(3))=-tmp
taaa(iva(2),iva(1),iva(3),ioa(2),ioa(1),ioa(3))= tmp
else if(nva.eq.0) then !bbb
tbbb(ivb(1),ivb(2),ivb(3),iob(1),iob(2),iob(3))= tmp
tbbb(ivb(1),ivb(3),ivb(2),iob(1),iob(2),iob(3))=-tmp
tbbb(ivb(3),ivb(1),ivb(2),iob(1),iob(2),iob(3))= tmp
tbbb(ivb(3),ivb(2),ivb(1),iob(1),iob(2),iob(3))=-tmp
tbbb(ivb(2),ivb(3),ivb(1),iob(1),iob(2),iob(3))= tmp
tbbb(ivb(2),ivb(1),ivb(3),iob(1),iob(2),iob(3))=-tmp
tbbb(ivb(1),ivb(2),ivb(3),iob(1),iob(3),iob(2))=-tmp
tbbb(ivb(1),ivb(3),ivb(2),iob(1),iob(3),iob(2))= tmp
tbbb(ivb(3),ivb(1),ivb(2),iob(1),iob(3),iob(2))=-tmp
tbbb(ivb(3),ivb(2),ivb(1),iob(1),iob(3),iob(2))= tmp
tbbb(ivb(2),ivb(3),ivb(1),iob(1),iob(3),iob(2))=-tmp
tbbb(ivb(2),ivb(1),ivb(3),iob(1),iob(3),iob(2))= tmp
tbbb(ivb(1),ivb(2),ivb(3),iob(3),iob(1),iob(2))= tmp
tbbb(ivb(1),ivb(3),ivb(2),iob(3),iob(1),iob(2))=-tmp
tbbb(ivb(3),ivb(1),ivb(2),iob(3),iob(1),iob(2))= tmp
tbbb(ivb(3),ivb(2),ivb(1),iob(3),iob(1),iob(2))=-tmp
tbbb(ivb(2),ivb(3),ivb(1),iob(3),iob(1),iob(2))= tmp
tbbb(ivb(2),ivb(1),ivb(3),iob(3),iob(1),iob(2))=-tmp
tbbb(ivb(1),ivb(2),ivb(3),iob(3),iob(2),iob(1))=-tmp
tbbb(ivb(1),ivb(3),ivb(2),iob(3),iob(2),iob(1))= tmp
tbbb(ivb(3),ivb(1),ivb(2),iob(3),iob(2),iob(1))=-tmp
tbbb(ivb(3),ivb(2),ivb(1),iob(3),iob(2),iob(1))= tmp
tbbb(ivb(2),ivb(3),ivb(1),iob(3),iob(2),iob(1))=-tmp
tbbb(ivb(2),ivb(1),ivb(3),iob(3),iob(2),iob(1))= tmp
tbbb(ivb(1),ivb(2),ivb(3),iob(2),iob(3),iob(1))= tmp
tbbb(ivb(1),ivb(3),ivb(2),iob(2),iob(3),iob(1))=-tmp
tbbb(ivb(3),ivb(1),ivb(2),iob(2),iob(3),iob(1))= tmp
tbbb(ivb(3),ivb(2),ivb(1),iob(2),iob(3),iob(1))=-tmp
tbbb(ivb(2),ivb(3),ivb(1),iob(2),iob(3),iob(1))= tmp
tbbb(ivb(2),ivb(1),ivb(3),iob(2),iob(3),iob(1))=-tmp
tbbb(ivb(1),ivb(2),ivb(3),iob(2),iob(1),iob(3))=-tmp
tbbb(ivb(1),ivb(3),ivb(2),iob(2),iob(1),iob(3))= tmp
tbbb(ivb(3),ivb(1),ivb(2),iob(2),iob(1),iob(3))=-tmp
tbbb(ivb(3),ivb(2),ivb(1),iob(2),iob(1),iob(3))= tmp
tbbb(ivb(2),ivb(3),ivb(1),iob(2),iob(1),iob(3))=-tmp
tbbb(ivb(2),ivb(1),ivb(3),iob(2),iob(1),iob(3))= tmp
else if(nva.eq.2) then !aab
taab(iva(1),iva(2),ivb(1),ioa(1),ioa(2),iob(1))= tmp
taab(iva(2),iva(1),ivb(1),ioa(1),ioa(2),iob(1))=-tmp
taab(iva(1),iva(2),ivb(1),ioa(2),ioa(1),iob(1))=-tmp
taab(iva(2),iva(1),ivb(1),ioa(2),ioa(1),iob(1))= tmp
else !abb
tabb(iva(1),ivb(1),ivb(2),ioa(1),iob(1),iob(2))= tmp
tabb(iva(1),ivb(2),ivb(1),ioa(1),iob(1),iob(2))=-tmp
tabb(iva(1),ivb(1),ivb(2),ioa(1),iob(2),iob(1))=-tmp
tabb(iva(1),ivb(2),ivb(1),ioa(1),iob(2),iob(1))= tmp
endif
endif
else
if(nex.eq.1) then
if(nva.eq.1) then
tt(ii)=ta(iva(1),ioa(1))
else
tt(ii)=tb(ivb(1),iob(1))
endif
else if(nex.eq.2) then
if(nva.eq.2) then
tt(ii)=taa(iva(1),iva(2),ioa(1),ioa(2))
else if(nva.eq.0) then
tt(ii)=tbb(ivb(1),ivb(2),iob(1),iob(2))
else
tt(ii)=tab(iva(1),ivb(1),ioa(1),iob(1))
endif
else
if(nva.eq.3) then !aaa
tt(ii)=taaa(iva(1),iva(2),iva(3),ioa(1),ioa(2),ioa(3))
else if(nva.eq.0) then !bbb
tt(ii)=tbbb(ivb(1),ivb(2),ivb(3),iob(1),iob(2),iob(3))
else if(nva.eq.2) then !aab
tt(ii)=taab(iva(1),iva(2),ivb(1),ioa(1),ioa(2),iob(1))
else !abb
tt(ii)=tabb(iva(1),ivb(1),ivb(2),ioa(1),iob(1),iob(2))
endif
endif
endif
enddo
enddo
enddo
enddo
C
return
end
C
************************************************************************
subroutine f12mname(method)
************************************************************************
* Generates the name of the explicitly correlated method
************************************************************************
implicit none
integer i,j
character(len=1) method(128)
C
i=1
do while ((ichar(method(i)).ge.65.and.ichar(method(i)).le.90)
$ .or.(ichar(method(i)).ge.48.and.ichar(method(i)).le.57))
i=i+1
enddo
if(method(i-1).eq.'C'.and.method(i).eq.'(') then
i=i+1
do while (method(i).ne.')')
i=i+1
enddo
endif
do j=128,i+6,-1
method(j)=method(j-6)
enddo
method(i )='('
method(i+1)='F'
method(i+2)='1'
method(i+3)='2'
method(i+4)='*'
method(i+5)=')'
i=i+6
if(method(i).eq.'('.or.method(i).eq.'[') then
i=i+1
do while (method(i).ne.')'.and.method(i).ne.']')
i=i+1
enddo
do j=128,i+1,-1
method(j)=method(j-1)
enddo
method(i+1)=method(i)
method(i)='+'
endif
C
return
end
C
************************************************************************
subroutine t3scale(iactva,iactvb,iactoa,iactob,nampvirtal,nex,tt,
$nmem,nnir,nmax,nactm,isympv,isympo,nstr,irec)
************************************************************************
* Scale amplitudes of triple and higher excitations
************************************************************************
#include "MRCCCOMMON"
integer nnir,nmax,nactm,i,j,k,nex,ii,jj,nn,ifile,nmem
integer n1,n2,n3,n4,ire,ile,l
integer iactv,iacto,iactva,iactoa,iactvb,iactob
integer nstr(nnir,0:nactm,0:nmax,4)
integer isympv(0:nnir,nnir,0:nactm,0:nactm,0:nmax,0:nmax,2)
integer isympo(0:nnir,nnir,0:nactm,0:nactm,0:nmax,0:nmax,2)
integer isymv,isymo,isymva,isymvb,isymoa,isymob,irv,iro,dbladd
integer nn1,nn2,nn3,nn4,irec(nnir,0:nactm,0:nmax,4)
integer nampvirtal,nampvirtbe,nampoccal,nampoccbe,intadd
real*8 tt(*),tfact,tscalea(nal),tscaleb(nbe)
C
open(scrfile7,file='F12INTE',form='unformatted')
read(scrfile7) !ecabs,emp2f12,ecoup
if(cs) then
read(scrfile7) tfact,tscalea
tscaleb=tscalea
else
read(scrfile7) tfact,tscalea,tscaleb
endif
close(scrfile7)
tscalea=tscalea*1.d0/dble(nex)
tscaleb=tscaleb*1.d0/dble(nex)
C
nampvirtbe=nex-nampvirtal
nampoccal=nampvirtal
nampoccbe=nex-nampoccal
ii=1
do ir=1,nir
isymv=csympair(isym,ir,1)
isymo=csympair(isym,ir,2)
do irv=1,isympv(0,isymv,iactva,iactvb,nampvirtal,
$nampvirtbe,1)
isymva=isympv(irv,isymv,iactva,iactvb,nampvirtal,nampvirtbe,1)
isymvb=isympv(irv,isymv,iactva,iactvb,nampvirtal,nampvirtbe,2)
n1=nstr(isymva,iactva,nampvirtal,1)
n2=nstr(isymvb,iactvb,nampvirtbe,2)
k=n1*n2
do iro=1,isympo(0,isymo,iactoa,iactob,nampoccal,
$nampoccbe,1)
isymoa=isympo(iro,isymo,iactoa,iactob,nampoccal,nampoccbe,1)
isymob=isympo(iro,isymo,iactoa,iactob,nampoccal,nampoccbe,2)
n3=nstr(isymoa,iactoa,nampoccal,3)
n4=nstr(isymob,iactob,nampoccbe,4)
call strread(nn1,nn2,nn3,nn4,n1,n2,n3,n4,nampvirtal,nampvirtbe,
$nampoccal,nampoccbe,isymva,isymvb,isymoa,isymob,iactva,iactvb,
$iactoa,iactob,nnir,nactm,nmax,irec,jj,imem+nmem)
call t3scale1(nampvirtal,nampvirtbe,nampoccal,nampoccbe,n1,n2,n3,
$n4,icore(nn3),icore(nn4),tt(ii),tscalea,tscaleb)
ii=ii+n3*n4*k
enddo
enddo
enddo
C
return
end
C
************************************************************************
subroutine t3scale1(nva,nvb,noa,nob,nvstral,nvstrbe,nostral,
$nostrbe,istroa,istrob,tt,tscalea,tscaleb)
************************************************************************
* Scale amplitudes of triple and higher excitations
************************************************************************
implicit none
integer nva,nvb,noa,nob,nvstral,nvstrbe,nostral,nostrbe,ii
integer i,ivstral,ivstrbe,iostral,iostrbe
integer istroa(noa,nostral),istrob(nob,nostrbe)
real*8 tt(*),tscoa,tscob,tscalea(*),tscaleb(*)
C
ii=0
do iostrbe=1,nostrbe
tscob=0.d0
do i=1,nob
tscob=tscob+tscaleb(istrob(i,iostrbe))
enddo
do iostral=1,nostral
tscoa=tscob
do i=1,noa
tscoa=tscoa+tscalea(istroa(i,iostral))
enddo
do ivstrbe=1,nvstrbe
do ivstral=1,nvstral
ii=ii+1
c write(6,"(2f9.5,100i2)") tt(ii),tscoa,
c $istroa(1:noa,iostral),istrob(1:nob,iostrbe)
c tt(ii)=tt(ii)*dsqrt(tscoa)
tt(ii)=tt(ii)*tscoa
enddo
enddo
enddo
enddo
C
return
end
C
************************************************************************
subroutine t3read(nal,nbe,nval,nvbe,taaa,tbbb,taab,tabb,scrfile7)
************************************************************************
* Read F12 T2 contribution to T3
************************************************************************
implicit none
integer nal,nbe,nval,nvbe,scrfile7
integer*2 a,b,c,i,j,k
real*8 taaa(nval,nval,nval,nal,nal,nal),tmp!,sum
real*8 tbbb(nvbe,nvbe,nvbe,nbe,nbe,nbe)
real*8 taab(nval,nval,nvbe,nal,nal,nbe)
real*8 tabb(nval,nvbe,nvbe,nal,nbe,nbe)
C
open(scrfile7,file='F12T3',form='unformatted')
c sum=0.d0
C aaa
do
read(scrfile7) tmp,a,b,c,i,j,k
if(a.eq.0) exit
taaa(a,b,c,i,j,k)=taaa(a,b,c,i,j,k)+tmp
taaa(a,c,b,i,j,k)=taaa(a,c,b,i,j,k)-tmp
taaa(c,a,b,i,j,k)=taaa(c,a,b,i,j,k)+tmp
taaa(c,b,a,i,j,k)=taaa(c,b,a,i,j,k)-tmp
taaa(b,c,a,i,j,k)=taaa(b,c,a,i,j,k)+tmp
taaa(b,a,c,i,j,k)=taaa(b,a,c,i,j,k)-tmp
taaa(a,b,c,i,k,j)=taaa(a,b,c,i,k,j)-tmp
taaa(a,c,b,i,k,j)=taaa(a,c,b,i,k,j)+tmp
taaa(c,a,b,i,k,j)=taaa(c,a,b,i,k,j)-tmp
taaa(c,b,a,i,k,j)=taaa(c,b,a,i,k,j)+tmp
taaa(b,c,a,i,k,j)=taaa(b,c,a,i,k,j)-tmp
taaa(b,a,c,i,k,j)=taaa(b,a,c,i,k,j)+tmp
taaa(a,b,c,k,i,j)=taaa(a,b,c,k,i,j)+tmp
taaa(a,c,b,k,i,j)=taaa(a,c,b,k,i,j)-tmp
taaa(c,a,b,k,i,j)=taaa(c,a,b,k,i,j)+tmp
taaa(c,b,a,k,i,j)=taaa(c,b,a,k,i,j)-tmp
taaa(b,c,a,k,i,j)=taaa(b,c,a,k,i,j)+tmp
taaa(b,a,c,k,i,j)=taaa(b,a,c,k,i,j)-tmp
taaa(a,b,c,k,j,i)=taaa(a,b,c,k,j,i)-tmp
taaa(a,c,b,k,j,i)=taaa(a,c,b,k,j,i)+tmp
taaa(c,a,b,k,j,i)=taaa(c,a,b,k,j,i)-tmp
taaa(c,b,a,k,j,i)=taaa(c,b,a,k,j,i)+tmp
taaa(b,c,a,k,j,i)=taaa(b,c,a,k,j,i)-tmp
taaa(b,a,c,k,j,i)=taaa(b,a,c,k,j,i)+tmp
taaa(a,b,c,j,k,i)=taaa(a,b,c,j,k,i)+tmp
taaa(a,c,b,j,k,i)=taaa(a,c,b,j,k,i)-tmp
taaa(c,a,b,j,k,i)=taaa(c,a,b,j,k,i)+tmp
taaa(c,b,a,j,k,i)=taaa(c,b,a,j,k,i)-tmp
taaa(b,c,a,j,k,i)=taaa(b,c,a,j,k,i)+tmp
taaa(b,a,c,j,k,i)=taaa(b,a,c,j,k,i)-tmp
taaa(a,b,c,j,i,k)=taaa(a,b,c,j,i,k)-tmp
taaa(a,c,b,j,i,k)=taaa(a,c,b,j,i,k)+tmp
taaa(c,a,b,j,i,k)=taaa(c,a,b,j,i,k)-tmp
taaa(c,b,a,j,i,k)=taaa(c,b,a,j,i,k)+tmp
taaa(b,c,a,j,i,k)=taaa(b,c,a,j,i,k)-tmp
taaa(b,a,c,j,i,k)=taaa(b,a,c,j,i,k)+tmp
c sum=sum+tmp
enddo
C bbb
do
read(scrfile7) tmp,a,b,c,i,j,k
if(a.eq.0) exit
tbbb(a,b,c,i,j,k)=tbbb(a,b,c,i,j,k)+tmp
tbbb(a,c,b,i,j,k)=tbbb(a,c,b,i,j,k)-tmp
tbbb(c,a,b,i,j,k)=tbbb(c,a,b,i,j,k)+tmp
tbbb(c,b,a,i,j,k)=tbbb(c,b,a,i,j,k)-tmp
tbbb(b,c,a,i,j,k)=tbbb(b,c,a,i,j,k)+tmp
tbbb(b,a,c,i,j,k)=tbbb(b,a,c,i,j,k)-tmp
tbbb(a,b,c,i,k,j)=tbbb(a,b,c,i,k,j)-tmp
tbbb(a,c,b,i,k,j)=tbbb(a,c,b,i,k,j)+tmp
tbbb(c,a,b,i,k,j)=tbbb(c,a,b,i,k,j)-tmp
tbbb(c,b,a,i,k,j)=tbbb(c,b,a,i,k,j)+tmp
tbbb(b,c,a,i,k,j)=tbbb(b,c,a,i,k,j)-tmp
tbbb(b,a,c,i,k,j)=tbbb(b,a,c,i,k,j)+tmp
tbbb(a,b,c,k,i,j)=tbbb(a,b,c,k,i,j)+tmp
tbbb(a,c,b,k,i,j)=tbbb(a,c,b,k,i,j)-tmp
tbbb(c,a,b,k,i,j)=tbbb(c,a,b,k,i,j)+tmp
tbbb(c,b,a,k,i,j)=tbbb(c,b,a,k,i,j)-tmp
tbbb(b,c,a,k,i,j)=tbbb(b,c,a,k,i,j)+tmp
tbbb(b,a,c,k,i,j)=tbbb(b,a,c,k,i,j)-tmp
tbbb(a,b,c,k,j,i)=tbbb(a,b,c,k,j,i)-tmp
tbbb(a,c,b,k,j,i)=tbbb(a,c,b,k,j,i)+tmp
tbbb(c,a,b,k,j,i)=tbbb(c,a,b,k,j,i)-tmp
tbbb(c,b,a,k,j,i)=tbbb(c,b,a,k,j,i)+tmp
tbbb(b,c,a,k,j,i)=tbbb(b,c,a,k,j,i)-tmp
tbbb(b,a,c,k,j,i)=tbbb(b,a,c,k,j,i)+tmp
tbbb(a,b,c,j,k,i)=tbbb(a,b,c,j,k,i)+tmp
tbbb(a,c,b,j,k,i)=tbbb(a,c,b,j,k,i)-tmp
tbbb(c,a,b,j,k,i)=tbbb(c,a,b,j,k,i)+tmp
tbbb(c,b,a,j,k,i)=tbbb(c,b,a,j,k,i)-tmp
tbbb(b,c,a,j,k,i)=tbbb(b,c,a,j,k,i)+tmp
tbbb(b,a,c,j,k,i)=tbbb(b,a,c,j,k,i)-tmp
tbbb(a,b,c,j,i,k)=tbbb(a,b,c,j,i,k)-tmp
tbbb(a,c,b,j,i,k)=tbbb(a,c,b,j,i,k)+tmp
tbbb(c,a,b,j,i,k)=tbbb(c,a,b,j,i,k)-tmp
tbbb(c,b,a,j,i,k)=tbbb(c,b,a,j,i,k)+tmp
tbbb(b,c,a,j,i,k)=tbbb(b,c,a,j,i,k)-tmp
tbbb(b,a,c,j,i,k)=tbbb(b,a,c,j,i,k)+tmp
c sum=sum+tmp
enddo
C aab
do
read(scrfile7) tmp,a,b,c,i,j,k
if(a.eq.0) exit
taab(a,b,c,i,j,k)=taab(a,b,c,i,j,k)+tmp
taab(b,a,c,i,j,k)=taab(b,a,c,i,j,k)-tmp
taab(a,b,c,j,i,k)=taab(a,b,c,j,i,k)-tmp
taab(b,a,c,j,i,k)=taab(b,a,c,j,i,k)+tmp
c sum=sum+tmp
enddo
C abb
do
read(scrfile7) tmp,a,b,c,i,j,k
if(a.eq.0) exit
tabb(a,b,c,i,j,k)=tabb(a,b,c,i,j,k)+tmp
tabb(a,c,b,i,j,k)=tabb(a,c,b,i,j,k)-tmp
tabb(a,b,c,i,k,j)=tabb(a,b,c,i,k,j)-tmp
tabb(a,c,b,i,k,j)=tabb(a,c,b,i,k,j)+tmp
c sum=sum+tmp
enddo
close(scrfile7)
c write(6,"(' Checksum ',f15.10)") sum
c write(6,*)
C
return
end
C