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

3041 lines
93 KiB
Fortran
Executable File

************************************************************************
subroutine jpq2pqrs(nocc,nvirt,dfnbasis,dcore,irecln,ibufln,ibuf,
$maxmem,imem,itol,ccprog)
************************************************************************
implicit none
integer ijai,dblalloc,iaibj,ijij,ijab,iabci,iabcd,msize
integer nocc,nvirt,dfnbasis,i,j,a,b,k,l,c,d,ii,itmp,itmp2
integer aa,recn,nbll,recnfrst,bmin,iname
integer irecln,ibufln,maxmem,imem
real*8 ibuf(ibufln),itol
integer ijdim,imin,ij,kl,kldim,jk,ab,nbuf,d1,abdim,d2,d4
integer nblo,nbllo,nblocko,odim,nb,nbl,ifrst,iblock
integer ad,cd,amin,cb,cdfrst,nblock,cddim,iisyev,abcim,abcdm
integer ijai2,ijab2,tot,d3,dblock,dlast,dfrst,abfrst
integer maxbatch,ai,nbllv2,nblv2,nblockv2,td,tb,tc,ta(nvirt)
integer p1,p2,nblv,nbllv,nblockv,ablock,cblock,cc,afrst,cfrst
parameter(msize=20) ! largest block of pqrs is 20 Mb
real*8 dcore(*),ibuf2(ibufln),temp,temp2,temp3
real*8,pointer :: jai(:,:,:),aibj(:,:,:,:),jij2(:,:,:),jab2(:,:)
real*8,pointer :: ijkl(:,:), jij(:,:),aijk(:,:,:),jab(:,:)
real*8,pointer :: abij(:,:),abci(:,:,:),abcd(:,:)
real*8,pointer :: jai2(:,:),abci2(:,:)
integer*4 isyev
equivalence(isyev,iisyev) !For Intel
logical jabrearr,tcl,tbl
character*4 ccprog
c db
integer pos,pos2,pos0
logical hrprt,dbprt
c db
c naftyp=jab
integer iw,ijii,dfnbasis_old,itr
real*8 naftol
real*8,pointer :: w(:,:),tr2naf(:,:)
character*16 cscr16,naf
character*8 naftyp
c
interface
c {{{ interfaces for pointers
subroutine rpoint2d(egydim,haromdim,dim1,dim2)
implicit none
integer dim1,dim2
real*8,target :: egydim(dim1,dim2)
real*8, pointer :: haromdim(:,:)
end subroutine
subroutine rpoint3d(egydim,haromdim,dim1,dim2,dim3)
implicit none
integer dim1,dim2,dim3
real*8,target :: egydim(dim1,dim2,dim3)
real*8, pointer :: haromdim(:,:,:)
end subroutine
subroutine rpoint4d(egydim,haromdim,dim1,dim2,dim3,dim4)
implicit none
integer dim1,dim2,dim3,dim4
real*8,target :: egydim(dim1,dim2,dim3,dim4)
real*8, pointer :: haromdim(:,:,:,:)
end subroutine
c}}}
end interface
c db
hrprt=.false.
dbprt=.false.
c db
jabrearr=.false.
c
c {{{ memory requirement
c memory requirment for ijkl, aibj, abij, aijk, Jai, Jij and Jab
itmp=5+max(nocc**2*nvirt**2,nocc**2*((nocc+1)/2)**2,
$ nocc**2*(nocc+1)/2*nvirt)+ ! full ijkl, abij, aijk aibj
$ 2*dfnbasis*max(nocc*nvirt,nocc*(nocc+1)/2)+ ! full Jai, Jab and Jij
$ dfnbasis*nvirt*(nvirt+1)/2
c memory requirment for the smallest block of abcd, abci
if (nvirt*(nvirt+1)/2*nvirt*nocc.lt.msize*1024**2) then
abcim=nvirt*(nvirt+1)/2*nvirt*nocc
else
abcim=msize*1024**2
endif
if (ccprog.eq.'ccsd') then
if (nvirt*(nvirt+1)/2*nvirt*(nvirt+1)/2.lt.msize*1024**2) then
abcdm=nvirt*(nvirt+1)/2*nvirt*(nvirt+1)/2
else
abcdm=msize*1024**2
endif
itmp2=dfnbasis*nvirt*(nvirt+1)/2+ max( ! Jab
$ dfnbasis*max(nvirt*(nvirt+1)/2,nocc*nvirt)+ ! Jai
c $ max(msize*1024**2,nvirt**2*(nvirt+1)/2)+dfnbasis*nvirt*nocc ! abci block and jai2
$ abcim+dfnbasis*nvirt*nocc ! abci block and jai2
c $ ,max( max(msize*1024**2,nvirt**2*(nvirt+1)/2), ! abcd block
$ ,max( abcdm, ! abcd block
$ dfnbasis*nvirt*(nvirt+1)/2)) ! jab2
elseif (ccprog.eq.'mrcc') then
itmp2=dfnbasis*nvirt*(nvirt+1)/2+ max( ! Jab
$ dfnbasis*max(nvirt*(nvirt+1)/2,nocc*nvirt)+ ! Jai
$ abcim+dfnbasis*nvirt*nocc ! abci block and jai2
$ ,nvirt**2) ! smallest block of (ab|cd)
endif
if (max(itmp,itmp2).ge.maxmem) then
write(0,*) 'Insufficient memory for integral assembly in ccsd.f'
write(0,*) 'Required:',8*int(max(itmp,itmp2)/1024**2), 'MB'
call mrccend(1)
endif
c }}}
c
c open files
if(ccprog.eq.'ccsd') then
open(556,file='56')
open(599,file='abcd.s',ACCESS='DIRECT',RECL=irecln)
open(600,file='abcd.a',ACCESS='DIRECT',RECL=irecln)
open(601,file='ijkl',ACCESS='DIRECT',RECL=irecln)
open(602,file='abci',ACCESS='DIRECT',RECL=irecln)
open(603,file='aijk',ACCESS='DIRECT',RECL=irecln)
open(604,file='abij',ACCESS='DIRECT',RECL=irecln)
open(605,file='aibj',ACCESS='DIRECT',RECL=irecln)
chr if (hrprt) then
chr open(502,file='abci.hr')
chr open(505,file='aibj.hr')
chr open(501,file='ijkl.hr')
chr open(503,file='aijk.hr')
chr open(504,file='abij.hr')
chr open(499,file='abcd.s.hr')
chr open(500,file='abcd.a.hr')
chr open(497,file='abcd.hr')
chr open(498,file='adcb.hr')
chr endif
endif
c allocate memory
ijab=dblalloc(dfnbasis*nvirt*(nvirt+1)/2)
ijai=dblalloc(dfnbasis*max(nvirt*nocc,nvirt*(nvirt+1)/2))
ijij=dblalloc(dfnbasis*max(nocc*nvirt,nocc*(nocc+1)/2)) ! for jij and jij2
c read J_ai from DFINT_AI in J(i,a,P) order
call rpoint3d(dcore(ijai),jai,nvirt,nocc,dfnbasis)
call rpoint3d(dcore(ijij),jij2,nocc,nvirt,dfnbasis)
open(111,file='DFINT_AI',form='unformatted')
read(111) jij2 ! read J_ia into jij first and transpose to jai
close(111)
c transponse: J_ia->J_ai
do i=1,nocc
do a=1,nvirt
call dcopy(dfnbasis,jij2(i,a,1:dfnbasis),1,jai(a,i,1:dfnbasis),1)
enddo
enddo
c read J_ij from DFINT_IJ in (P,ij) order, upper triangle format for ij indices
call rpoint2d(dcore(ijij),jij,dfnbasis,nocc*(nocc+1)/2)
open(111,file='DFINT_IJ',form='unformatted')
read(111) jij
close(111)
c read J_ab from DFINT_AB in (P,ab) order, upper triangle format for ab indices
call rpoint2d(dcore(ijab),jab,dfnbasis,nvirt*(nvirt+1)/2)
open(111,file='DFINT_AB',form='unformatted')
read(111) jab
close(111)
c db
c if (dbprt) then
c write(*,*) 'jai read in ccsd.f'
c write(*,*) jai
c write(*,*) 'jij read in ccsd.f'
c write(*,*) jij
c write(*,*) 'jab read in ccsd.f',dfnbasis,nvirt
c write(*,*) jab
c endif
c db
c {{{ pilot Jpq based NAF construction
call getkey('naftyp',6,naftyp,8)
if (naftyp.eq.'jpqpilot') then
write(*,*) 'SZEMET !!!!, all jpq are added to W'
c
call getkey('naf_cor',7,naf,16)
if(naf.ne.'off ') then
read(naf,*) naftol
else
write(6,*) 'Inconsistent naftyp and naf_cor!'
call mrccend(1)
endif
c
itr=dblalloc(dfnbasis**2)
iw=dblalloc(dfnbasis**2)
ijii=iw+dfnbasis**2 ! imem
call rpoint2d(dcore(itr),tr2naf,dfnbasis,dfnbasis)
call rpoint2d(dcore(iw),w,dfnbasis,dfnbasis)
c write(*,*) 'nvirt,dfnbasis,pp',nvirt,dfnbasis
c write(*,'(100000f14.8)') jab
C add Jab*Jab^T to W
call dsyrk('u','n',dfnbasis,nvirt*(nvirt+1)/2,1.d0,jab,dfnbasis,
$ 0.d0,w,dfnbasis)
c db
c if (dbprt) then
c write(*,*) 'Jab contribution to W',dfnbasis,nvirt
c write(*,'(100000f14.8)') w(1:dfnbasis,1:dfnbasis) ; flush(6)
c endif
c db
C add Jai^T*Jai to W
call dsyrk('u','t',dfnbasis,nvirt*nocc,1.d0,jai,nvirt*nocc,1.d0
$,w,dfnbasis)
C add Jij*Jij^T to W
call dsyrk('u','n',dfnbasis,nocc*(nocc+1)/2,1.d0,jij,dfnbasis,
$ 1.d0,w,dfnbasis)
c db
c if (dbprt) then
c write(*,*) 'W before diag'
c write(*,*) w(1:dfnbasis,1:dfnbasis) ; flush(6)
c endif
c db
c construct NAFs from W
call dsyev('V','U',dfnbasis,w,dfnbasis,dcore(ijii),
$dcore(ijii+dfnbasis),20*dfnbasis,isyev)
if(isyev.ne.0) then
write(6,*) 'Fatal error at the construction of NAFs jab!'
call mrccend(1)
endif
c db
c if (dbprt) then
c write(*,*) 'W eigenvals with jpq', naftol
c write(*,*) dcore(ijii:ijii+dfnbasis-1)
c write(*,*) 'NAF coeff matrix'
c write(*,*) w(1:dfnbasis,1:dfnbasis)
c endif
c db
if (dcore(ijii).lt.-max(0.001d0,naftol)) then
write(*,*) 'Warning: auxiliary basis might be linearly dependent'
write(*,*) 'All NAFs with imaginary singular values are dropped'
endif
i=0
do while((dabs(dcore(ijii+i)).lt.naftol.or.
$ dcore(ijii+i) .lt.-1.d-10).and.i.lt.dfnbasis) ! to exclude imag singular values
i=i+1
enddo
dfnbasis_old=dfnbasis
dfnbasis=dfnbasis-i
call dcopy(dfnbasis_old*dfnbasis,w(1,i+1),1,tr2naf,1)
write(cscr16,'(f16.2)') 100.d0*dfnbasis/dfnbasis_old
write(6,"(' Number of NAFs: ',i4,1x,a9)")
$dfnbasis,trim('(' // adjustl(cscr16)) // '%)'
call dbldealloc(iw)
write(*,*) 'SZEMET !!!!, rearrange memory if non pilot'
itmp=
$ dblalloc(dfnbasis*max(nvirt*(nvirt+1)/2,nocc*(nocc+1)/2,
$ nocc*nvirt))
c trf jai to the NAF basis
call dgemm('n','n',nvirt*nocc,dfnbasis,dfnbasis_old,1.d0,
$jai,nocc*nvirt,tr2naf,dfnbasis_old,0.d0,dcore(itmp),nvirt*nocc)
call dcopy(nocc*nvirt*dfnbasis,dcore(itmp),1,jai,1)
c db
c if (dbprt) then
c write(*,*) 'jai in the NAF basis in ccsd.f'
c write(*,*) jai
c write(*,*) dcore(itmp:itmp+nocc*nvirt*dfnbasis-1)
c endif
c db
c trf jij to the NAF basis
call dgemm('t','n',dfnbasis,nocc*(nocc+1)/2,dfnbasis_old,1.d0,
$tr2naf,dfnbasis_old,jij,dfnbasis_old,0.d0,dcore(itmp),dfnbasis)
call dcopy(nocc*(nocc+1)/2*dfnbasis,dcore(itmp),1,jij,1)
c db
c if (dbprt) then
c write(*,*) 'jij in the NAF basis in ccsd.f'
c write(*,*) jij
c write(*,*) dcore(itmp:itmp+nocc*(nocc+1)/2*dfnbasis-1)
c endif
c db
c trf jab to the NAF basis
call dgemm('t','n',dfnbasis,nvirt*(nvirt+1)/2,dfnbasis_old,1.d0,
$tr2naf,dfnbasis_old,jab,dfnbasis_old,0.d0,dcore(itmp),dfnbasis)
call dcopy(nvirt*(nvirt+1)/2*dfnbasis,dcore(itmp),1,jab,1)
c db
c if (dbprt) then
c write(*,*) 'jab in the NAF basis in ccsd.f'
c write(*,*) dcore(itmp:itmp+nvirt*(nvirt+1)/2*dfnbasis-1)
c write(*,*) jab
c endif
c db
call dbldealloc(itr)
endif ! if (naftyp.eq.'jpqpilot')
c }}}
c assemble (ai|bj) in one batch
d1=max(nocc**2*nvirt**2,nocc**2*nvirt*(nocc+1)/2,
$ nocc**2*((nocc+1)/2)**2)
iaibj=dblalloc(d1) ! for ijkl , aibj (> abij), aijk
call rpoint4d(dcore(iaibj),aibj,nvirt,nocc,nvirt,nocc)
call dsyrk('u','n',nvirt*nocc,dfnbasis,1.d0,jai,nvirt*nocc,0.d0,
$aibj,nvirt*nocc)
call filllo(aibj,nvirt*nocc)
c {{{ print <ab|ij> for ccsd.f
if(ccprog.eq.'ccsd') then
aa=0
recn=1
c write(*,*) 'abij'
do i=1,nocc
recnfrst=recn
do j=i,nocc
do a=1,nvirt
if (i.eq.j) then
bmin=a
else ! i>j
bmin=1
endif ! (i.eq.j)
do b=bmin,nvirt
chr if (hrprt) write(504,7000) aibj(b,j,a,i),
chr $ nocc+a,nocc+b,j,i
c write(*,7000) aibj(b,j,a,i), nocc+a,nocc+b,j,i
aa=aa+1
ibuf(aa)=aibj(b,j,a,i)
if (aa.eq.ibufln) then
write(604,rec=recn) ibuf(1:aa)
aa=0
recn=recn+1
endif ! aa.eq.ibufln
enddo ! a
enddo ! b
enddo ! j
if (aa.ne.0) then
write(604,rec=recn) ibuf(1:aa)
aa=0
recn=recn+1
endif
write(556,'(6I24)') iname('abij'),1,1,i,recnfrst,recn-1
enddo ! i
elseif(ccprog.eq.'mrcc') then
call aibj_to_fort55(nocc,nvirt,aibj,itol,155)
endif
c }}}
c assemble (ij|kl) in one batch
call rpoint2d(dcore(iaibj),ijkl,nocc*(nocc+1)/2,nocc*(nocc+1)/2)
ijdim=nocc*(nocc+1)/2
call dsyrk('u','t',ijdim,dfnbasis,1.d0,jij,dfnbasis,0.d0,
$ijkl,ijdim)
c {{{ print <ij|kl> for ccsd.f
if(ccprog.eq.'ccsd') then
c write(*,*) 'ijkl', nocc
aa=0
recn=1
do l=1,nocc
recnfrst=recn
do j=l,nocc
do k=l,nocc
kl=k*(k-1)/2+l
if (l.eq.j) then
imin=k
elseif (l.ne.j) then ! i>=j
imin=j
endif
do i=imin,nocc
if (i.ge.j) then
ij=i*(i-1)/2+j
elseif (i.lt.j) then
ij=j*(j-1)/2+i
endif
aa=aa+1
if (kl.le.ij) then
chr if (hrprt) write(501,7000) ijkl(kl,ij),i,k,j,l
c write(*,7000) ijkl(kl,ij),i,k,j,l
ibuf(aa)=ijkl(kl,ij)
elseif (kl.gt.ij) then
chr if (hrprt) write(501,7000) ijkl(ij,kl),i,k,j,l
c write(*,7000) ijkl(ij,kl),i,k,j,l
ibuf(aa)=ijkl(ij,kl)
endif
if (aa.eq.ibufln) then
write(601,rec=recn) ibuf(1:aa)
aa=0
recn=recn+1
endif ! aa.lt.ibufln
enddo ! i
enddo ! k
enddo ! j
c print last, unfull record and input file 56 for ccsd.f
if (aa.ne.0) then
write(601,rec=recn) ibuf(1:aa)
aa=0
recn=recn+1
endif
write(556,'(6I24)') iname('ijkl'),1,1,l,recnfrst,recn-1
enddo ! l
elseif(ccprog.eq.'mrcc') then
call ijkl_to_fort55(nocc,ijkl,itol,155)
endif
c }}}
c assemble (ai|jk) in one batch
call rpoint3d(dcore(iaibj),aijk,nvirt,nocc,nocc*(nocc+1)/2)
call dgemm('n','n',nvirt*nocc,ijdim,dfnbasis,1.d0,jai,
$nvirt*nocc,jij,dfnbasis,0.d0,aijk,nvirt*nocc)
c {{{ print <ai|jk> for ccsd.f
if(ccprog.eq.'ccsd') then
aa=0
recn=1
c write(*,*) 'aijk'
do k=1,nocc
recnfrst=recn
do i=1,nocc
do j=k,nocc
jk=j*(j-1)/2+k
do a=1,nvirt
chr if (hrprt) write(503,7000) aijk(a,i,jk),nocc+a,j,i,k
c write(*,7000) aijk(a,i,jk),nocc+a,j,i,k
aa=aa+1
ibuf(aa)=aijk(a,i,jk)
if (aa.eq.ibufln) then
write(603,rec=recn) ibuf(1:aa)
aa=0
recn=recn+1
endif ! aa.lt.ibufln
enddo ! a
enddo ! j
enddo ! i
if (aa.ne.0) then
write(603,rec=recn) ibuf(1:aa)
aa=0
recn=recn+1
endif
write(556,'(6I24)') iname('aijk'),1,1,k,recnfrst,recn-1
enddo ! k
elseif(ccprog.eq.'mrcc') then
call aijk_to_fort55(aijk,nocc,nvirt,itol,155)
endif
c }}}
c assemble (ab|ij) in one batch
call rpoint2d(dcore(iaibj),abij,nvirt*(nvirt+1)/2,nocc*(nocc+1)/2)
abdim=nvirt*(nvirt+1)/2
call dgemm('t','n',abdim,ijdim,dfnbasis,1.d0,jab,dfnbasis,
$jij,dfnbasis,0.d0,abij,abdim)
c {{{ print <ai|bj> for ccsd.f
if(ccprog.eq.'ccsd') then
nbuf=0
recn=1
do j=1,nocc
recnfrst=recn
do b=1,nvirt
do i=j,nocc
ij=i*(i-1)/2+j
do a=b,nvirt
ab=a*(a-1)/2+b
chr if (hrprt)write(505,7000) abij(ab,ij),nocc+a,i,nocc+b,j
nbuf=nbuf+1
ibuf(nbuf)=abij(ab,ij)
if (nbuf.eq.ibufln) then
write(605,rec=recn) ibuf(1:nbuf)
nbuf=0
recn=recn+1
endif ! nbuf.eq.ibufln
enddo ! a
enddo ! i
enddo ! b
if (nbuf.ne.0) then
write(605,rec=recn) ibuf(1:nbuf)
nbuf=0
recn=recn+1
endif
write(556,'(6I24)') iname('aibj'),1,1,j,recnfrst,recn-1
enddo ! j
elseif(ccprog.eq.'mrcc') then
call abij_to_fort55(nocc,nvirt,abij,itol,155)
endif
c }}}
c assemble (ab|ci)
call dbldealloc(ijij)
d2=nvirt*(nvirt+1)/2*nvirt*nocc ! full abci
d3=(nvirt*(nvirt+1)/2)**2 ! full abcd
maxbatch=max(msize*1024**2,nvirt**2*(nvirt+1)/2) ! 160 MB or one 1 column wide block of abcd or abci
c assemble (ab|ci)
if (d2.le.maxbatch) then ! if abci fits into 160 MB
c db write(*,*) 'print (ab|ci) from memory in one batch'
c assemble (ab|ci) in one batch
iabci=dblalloc(nvirt**2*nocc*(nvirt+1)/2) ! for abci
call rpoint3d(dcore(iabci),abci,nvirt*(nvirt+1)/2,nvirt,nocc)
call dgemm('t','t',abdim,nvirt*nocc,dfnbasis,1.d0,jab,dfnbasis,
$jai,nvirt*nocc,0.d0,abci,abdim)
c {{{ print <ab|ci> for ccsd.f
if(ccprog.eq.'ccsd') then
nbuf=0
recn=1
c write(*,*) 'abci'
do i=1,nocc
recnfrst=recn
do b=1,nvirt
do c=1,nvirt
do a=b,nvirt
ab=a*(a-1)/2+b
chr if (hrprt) write(502,7000) abci(ab,c,i),nocc+a,nocc+c,nocc+b,i
c write(*,7000) abci(ab,c,i),nocc+a,nocc+c,nocc+b,i
nbuf=nbuf+1
ibuf(nbuf)=abci(ab,c,i)
if (nbuf.eq.ibufln) then
write(602,rec=recn) ibuf(1:nbuf)
nbuf=0
recn=recn+1
endif ! nbuf.lt.ibufln
enddo ! a
enddo ! c
enddo ! b
if (nbuf.ne.0) then
write(602,rec=recn) ibuf(1:nbuf)
nbuf=0
recn=recn+1
endif
write(556,'(6I24)') iname('abci'),1,1,i,recnfrst,recn-1
enddo ! i
elseif(ccprog.eq.'mrcc') then
call abci_bld_to_fort55(nvirt,nocc,dfnbasis,abdim,jab,jai,
$ abci,itol,.false.,155)
endif
c }}}
else
c {{{ assemble (ab|ci) in blocks
c memory requirement, block size determination: abci blocks are smaller than maxbatch
nblo=0
itmp=dfnbasis*nvirt*(nocc+2*(nvirt+1)/2) ! fix: Jai2 and Jab + orig Jai
do while(nblo.lt.nocc)
nblo=nblo+1
if (itmp+nvirt**2*(nvirt+1)/2*nblo.gt.maxmem) then
nblo=nblo-1
exit
endif
if (nvirt**2*(nvirt+1)/2*nblo.gt.maxbatch) then ! abci block mem req exceeds maxbatch
nblo=nblo-1
exit
endif
enddo ! while
nbllo=mod(nocc,nblo)
nblocko=(nocc-nbllo)/nblo
if(nbllo.ne.0) nblocko=nblocko+1
if(nbllo.eq.0) nbllo=nblo
c db write(*,*) 'print (ab|ci) from memory in', nblocko ,'batch'
ijai2=dblalloc(dfnbasis*nvirt*nocc)
call rpoint2d(dcore(ijai2),jai2,dfnbasis,nocc*nvirt)
iabci=dblalloc(nvirt**2*nblo*(nvirt+1)/2) ! for the largest block of abci
call rpoint2d(dcore(iabci),abci2,nvirt*(nvirt+1)/2,nvirt*nblo)
c rearrange Jai
do i=1,nocc
do a=1,nvirt
ai=(i-1)*nvirt+a
call dcopy(dfnbasis,jai(a,i,1),nocc*nvirt,jai2(1,ai),1)
enddo
enddo
c repartition the memory space of jai -> jab2
call rpoint2d(dcore(ijai),jab2,dfnbasis,nvirt*(nvirt+1)/2)
c rearrange second dimension of jab from by row to by column -> jab2
d=1
tot=0
do b=1,nvirt
do a=b,nvirt
tot=tot+1
c=tot+d-1
cd=(2*nvirt-(d-2))*(d-1)/2+(c-(d-1))
call dcopy(dfnbasis,jab(1,a*(a-1)/2+b),1,
$ jab2(1,cd),1)
if (tot.eq.nvirt-d+1) then
d=d+1
tot=0
endif
enddo
enddo
if (d3.le.maxbatch.or.ccprog.eq.'mrcc') then ! if abcd fits into maxbatch
jabrearr=.false.
else
call dcopy(dfnbasis*nvirt*(nvirt+1)/2,jab2,1,jab,1)
jabrearr=.true.
endif
c loops for <ac|bi> print
nbuf=0
recn=1
do i=1,nocc
recnfrst=recn
if (mod(i,nblo).eq.0) then ! which block according to i
iblock=i/nblo
else
iblock=i/nblo+1
endif
c do assembly only in the beginning of iblock blocks of (ab|ci)
if ((mod(i,nblo).eq.1).or.(nblo.eq.1)) then
ifrst=(iblock-1)*nblo
cdfrst=nvirt*ifrst+1 ! cifrst
cddim=nvirt*nblo ! cidim
odim=nblo ! idim
if (iblock.eq.nblocko) then
cddim=nvirt*nbllo
odim=nbllo
endif
abdim=nvirt*(nvirt+1)/2
! intel 15.0.2 error
call dgemm('t','n',abdim,cddim,dfnbasis,1.d0,
$jab2,dfnbasis,jai2(1,cdfrst),dfnbasis,0.d0,abci2,abdim)
c call dgemm('t','n',abdim,cddim,dfnbasis,1.d0,
c $jab2(1:dfnbasis,1:abdim),dfnbasis,
c $jai2(1:dfnbasis,cdfrst:cdfrst+cddim-1),dfnbasis,0.d0,
c $abci2(1:abdim,1:cddim),abdim)
endif ! end of do assembly
if(ccprog.eq.'ccsd') then
do b=1,nvirt
do c=1,nvirt
cd=nvirt*(i-ifrst-1)+c !ci
do a=b,nvirt
ab=(2*nvirt-(b-2))*(b-1)/2+(a-(b-1))
chr if(hrprt) write(502,7000) abci2(ab,cd),nocc+a,nocc+c,nocc+b,i
nbuf=nbuf+1
ibuf(nbuf)=abci2(ab,cd)
if (nbuf.eq.ibufln) then
write(602,rec=recn) ibuf(1:nbuf)
nbuf=0
recn=recn+1
endif ! nbuf.lt.ibufln
enddo ! a
enddo ! c
enddo ! b
if (nbuf.ne.0) then
write(602,rec=recn) ibuf(1:nbuf)
nbuf=0
recn=recn+1
endif
write(556,'(6I24)') iname('abci'),1,1,i,recnfrst,recn-1
elseif(ccprog.eq.'mrcc') then
do b=1,nvirt
do c=1,nvirt
cd=nvirt*(i-ifrst-1)+c !ci
do a=b,nvirt
ab=(2*nvirt-(b-2))*(b-1)/2+(a-(b-1))
temp3=abci2(ab,cd)
if(dabs(temp3).gt.itol)
$write(155,7000) temp3,nocc+a,nocc+b,nocc+c,i
enddo ! a
enddo ! c
enddo ! b
endif !if(ccprog.eq.
enddo ! i
c }}}
endif ! assemble abci
call dbldealloc(ijai)
c assemble (ab|cd)
if (ccprog.eq.'ccsd') then
if (d3.le.maxbatch) then
c assemble (ab|cd) in one batch
c write(*,*) 'print (ab|cd) from memory in one batch' ,d3,maxbatch
iabcd=dblalloc((nvirt*(nvirt+1)/2)**2) ! for full abcd
call rpoint2d(dcore(iabcd),abcd,nvirt*(nvirt+1)/2,
$ nvirt*(nvirt+1)/2)
call dsyrk('u','t',abdim,dfnbasis,1.d0,jab,dfnbasis,0.d0,
$abcd,abdim)
c {{{ print <ab|cd> for ccsd.f
nbuf=0
recn=1
do d=1,nvirt
recnfrst=recn
do b=d,nvirt
do c=d,nvirt
cd=c*(c-1)/2+d
if (c.eq.d) then
amin=b
elseif (c.ne.d) then
amin=c
endif
if (c.ge.b) then
cb=c*(c-1)/2+b
else
cb=b*(b-1)/2+c
endif
do a=amin,nvirt
if (a.ge.b) then
ab=a*(a-1)/2+b
else
ab=b*(b-1)/2+a
endif
ad=a*(a-1)/2+d
if (cd.le.ab) then
temp=abcd(cd,ab) ! = (ab|cd) = <ac|bd>
elseif (cd.gt.ab) then
temp=abcd(ab,cd)
endif
if (ad.le.cb) then
temp2=abcd(ad,cb) ! = (ad|cb) = <ac|db>
elseif (ad.gt.cb) then
temp2=abcd(cb,ad)
endif
nbuf=nbuf+1
ibuf(nbuf)=temp+temp2
ibuf2(nbuf)=temp-temp2
chr if (hrprt) write(497,7000) temp,nocc+a,nocc+c,
chr $ nocc+b,nocc+d
chr if (hrprt) write(498,7000) temp2,nocc+a,nocc+c,
chr $ nocc+b,nocc+d
chr if (hrprt) write(499,7000) ibuf(nbuf),nocc+a,nocc+c,
chr $ nocc+b,nocc+d
chr if (hrprt) write(500,7000) ibuf2(nbuf),nocc+a,nocc+c,
chr $ nocc+b,nocc+d
if (nbuf.eq.ibufln) then
write(599,rec=recn) ibuf(1:nbuf)
write(600,rec=recn) ibuf2(1:nbuf)
nbuf=0
recn=recn+1
endif ! nbuf.lt.ibufln
enddo ! a
enddo ! c
enddo ! b
if (nbuf.ne.0) then
write(599,rec=recn) ibuf(1:nbuf)
write(600,rec=recn) ibuf2(1:nbuf)
nbuf=0
recn=recn+1
endif
write(556,'(6I24)') iname('abcd'),1,1,d+nocc,recnfrst,recn-1
enddo ! d
c }}}
else
c {{{ assemble (ab|cd) in blocks
c memory requirement, block size determination: abci blocks are smaller than maxbatch
nblv2=0
if (.not.jabrearr) then
itmp=nvirt*(nvirt+1)/2*(dfnbasis+max(dfnbasis,nvirt)) ! fix 1 Jab and (abcd or Jab2)
else
itmp=nvirt*(nvirt+1)/2*(dfnbasis+nvirt)
endif
do while(nblv2.lt.nvirt)
nblv2=nblv2+1
if (itmp+nvirt*(nvirt+1)/2*(2*nvirt-(nblv2-1))*nblv2/2
$ .gt.maxmem) then
nblv2=nblv2-1
exit
endif
if (nvirt*(nvirt+1)/2*(2*nvirt-(nblv2-1))*nblv2/2
$ .gt.maxbatch) then ! abcd block mem req exceeds maxbatch
nblv2=nblv2-1
exit
endif
enddo ! while
nbllv2=mod(nvirt,nblv2)
nblockv2=(nvirt-nbllv2)/nblv2
if(nbllv2.ne.0) then
nblockv2=nblockv2+1
else
nbllv2=nblv2
endif
c write(*,*) 'print (ab|cd) from memory in', nblockv2 ,'batch'
c rearrange second dimension of jab from by row to by column -> jab2
if (.not.jabrearr) then ! (if it is not done already for abci assembly)
ijab2=dblalloc(dfnbasis*nvirt*(nvirt+1)/2)
call rpoint2d(dcore(ijab2),jab2,dfnbasis,nvirt*(nvirt+1)/2)
d=1
tot=0
do b=1,nvirt
do a=b,nvirt
tot=tot+1
c=tot+d-1
cd=(2*nvirt-(d-2))*(d-1)/2+(c-(d-1))
call dcopy(dfnbasis,jab(1,a*(a-1)/2+b),1,jab2(1,cd),1)
if (tot.eq.nvirt-d+1) then
d=d+1
tot=0
endif
enddo
enddo
call dcopy(dfnbasis*nvirt*(nvirt+1)/2,jab2,1,jab,1)
jabrearr=.true.
call dbldealloc(ijab2)
endif ! jabrearr
iabcd=dblalloc((nvirt*(nvirt+1)/2)*(2*nvirt-(nblv2-1))*nblv2/2) ! largest block of abcd
call rpoint2d(dcore(iabcd),abcd,nvirt*(nvirt+1)/2,
$ (2*nvirt-(nblv2-1))*nblv2/2)
c {{{ print (ab|cd) for ccsd.f
c write(*,*) 'print (ab|cd) using blocks'
nbuf=0
recn=1
do d=1,nvirt
recnfrst=recn
if (mod(d,nblv2).eq.0) then ! which block according to d
dblock=d/nblv2
else
dblock=d/nblv2+1
endif
c do assembly only in the beginning of dblock blocks of (ab|cd)
if ((mod(d,nblv2).eq.1).or.(nblv2.eq.1)) then
dlast=dblock*nblv2
if (dblock.eq.nblockv2) dlast=nvirt
c write(cscr,'(i8)') dblock
c cscr=adjustl(cscr)
c write(iout,*) 'Assembly of <ab|cd> integrals of batch '
c $ // trim(cscr) // '...', 'orbitals',dlast
dfrst=(dblock-1)*nblv2
cdfrst=(2*nvirt-(dfrst-1))*dfrst/2+1
cddim=(2*nvirt-(2*dblock-1)*nblv2+1)*nblv2/2
if (dblock.eq.nblockv2) cddim=nbllv2*(nbllv2+1)/2
abdim=(nvirt-dblock*nblv2)*(nvirt-dblock*nblv2+1)/2
abfrst=cdfrst+cddim
call dsyrk('u','t',cddim,dfnbasis,1.d0,
$jab(1,cdfrst),dfnbasis,0.d0,abcd,nvirt*(nvirt+1)/2)
! intel error with ifort 15.0.2
! call dsyrk('u','t',cddim,dfnbasis,1.d0,
! $jab(1:dfnbasis,cdfrst:cdfrst+cddim-1),dfnbasis,0.d0,
! $abcd(1:cddim,1:cddim),cddim)
if (dblock.ne.nblockv2)
$ call dgemm('t','n',abdim,cddim,dfnbasis,1.d0,
$jab(1,abfrst),dfnbasis,jab(1,cdfrst),dfnbasis,0.d0,
$abcd(cddim+1,1),nvirt*(nvirt+1)/2)
c intel error with ifort 15.0.2
c call dgemm('t','n',abdim,cddim,dfnbasis,1.d0,
c $jab(1:dfnbasis,abfrst:abfrst+abdim-1),dfnbasis,
c $jab(1:dfnbasis,cdfrst:cdfrst+cddim-1),dfnbasis,0.d0,
c $abcd(cddim+1:cddim+abdim,1:cddim),abdim)
endif ! end of do assembly
td=(2*nvirt-dfrst-(d-2))*(d-dfrst-1)/2-(d-1)
do a=1,nvirt
ta(a)=(2*nvirt-dfrst-(a-2))*(a-dfrst-1)/2-(a-1) ! temp vector for combined index comp.
enddo
do b=d,nvirt
tb=(2*nvirt-dfrst-(b-2))*(b-dfrst-1)/2-(b-1)
tbl=b.le.dlast
do c=d,nvirt
tc=(2*nvirt-dfrst-(c-2))*(c-dfrst-1)/2-(c-1)
tcl=c.le.dlast
cd=td+c
if (c.ge.b) then
cb=tb+c
else
cb=tc+b
endif
if (c.eq.d) then
amin=b
else
amin=c
endif
do a=amin,nvirt
ad=td+a
if (a.ge.b) then
ab=tb+a
else
ab=ta(a)+b
endif
if (cd.le.ab) then
if (tbl.and.(b.le.a)) then
temp=abcd(cd,ab)
else
if ((a.le.dlast).and.(a.lt.b)) then
temp=abcd(cd,ab)
else
temp=abcd(ab,cd)
endif
endif
else
temp=abcd(ab,cd)
endif
c
if (ad.le.cb) then
if (tbl.and.(b.le.c)) then
temp2=abcd(ad,cb)
else
if (tcl.and.(c.lt.b)) then
temp2=abcd(ad,cb)
else
temp2=abcd(cb,ad)
endif
endif
else
temp2=abcd(cb,ad)
endif
nbuf=nbuf+1
ibuf(nbuf)=temp+temp2
ibuf2(nbuf)=temp-temp2
chr if (hrprt) then
chr write(497,7000) temp,nocc+a,nocc+c,
chr $ nocc+b,nocc+d
chr write(498,7000) temp2,nocc+a,nocc+c,
chr $ nocc+b,nocc+d
chr write(499,7000) ibuf(nbuf),nocc+a,nocc+c,
chr $ nocc+b,nocc+d
chr write(500,7000) ibuf2(nbuf),nocc+a,nocc+c,
chr $ nocc+b,nocc+d
chr endif
if (nbuf.eq.ibufln) then
write(599,rec=recn) ibuf(1:nbuf)
write(600,rec=recn) ibuf2(1:nbuf)
nbuf=0
recn=recn+1
endif ! nbuf.lt.ibufln
enddo ! a
enddo ! c
enddo ! b
if (nbuf.ne.0) then
write(599,rec=recn) ibuf(1:nbuf)
write(600,rec=recn) ibuf2(1:nbuf)
nbuf=0
recn=recn+1
endif
write(556,'(6I24)') iname('abcd'),1,1,d+nocc,recnfrst,recn-1
enddo ! d
c }}}
c }}}
endif
elseif (ccprog.eq.'mrcc') then
d4=dfnbasis*nvirt*(nvirt+1)/2+(nvirt*(nvirt+1)/2)**2 ! mem req for Jab + full (ab|cd)
if (d4.le.maxmem) then
c {{{ assemble (ab|cd) in one batch and save to fort.55
iabcd=dblalloc((nvirt*(nvirt+1)/2)**2) ! for full abcd
call rpoint2d(dcore(iabcd),abcd,nvirt*(nvirt+1)/2,
$ nvirt*(nvirt+1)/2)
call abcd_bld_to_fort55(nocc,nvirt,nvirt*(nvirt+1)/2,dfnbasis,
$ jab,abcd,itol,155)
c }}}
else
c {{{ assemble (ab|cd) in blocks
nblv=0
itmp=dfnbasis*nvirt*(nvirt+1)/2 ! full Jab
do while(nblv.lt.nvirt)
nblv=nblv+1
if (itmp+((2*nvirt-(nblv-1))*nblv/2)**2.gt.maxmem) then ! full Jab + largest block of (ab|cd)
c if (itmp+((2*nvirt-(nblv-1))*nblv/2)**2.gt.
c % itmp+((2*nvirt-(3-1))*3/2)**2) then ! db
nblv=nblv-1
exit
endif
enddo ! while
nbllv=mod(nvirt,nblv)
nblockv=(nvirt-nbllv)/nblv
if(nbllv.ne.0) then
nblockv=nblockv+1
else
nbllv=nblv
endif
c write(*,*) 'save (ab|cd) from memory in',nblockv,'batch for mrcc'
itmp=(2*nvirt-(nblv-1))*nblv/2 ! max block size for Jab
iabcd=dblalloc(itmp**2) ! all blocks of (ab|cd) fits into this array
call rpoint2d(dcore(iabcd),abcd,itmp,itmp)
c assemly symmetry inequivalent blocks of (ab|cd)
nbl=nblv
do ablock=1,nblockv
if (ablock.eq.nblockv) nbl=nbllv
afrst=(ablock-1)*nblv
p1=afrst*(afrst+1)/2+1
abdim=(2*afrst+nbl+1)*nbl/2
do cblock=1,ablock-1
cfrst=(cblock-1)*nblv
p2=cfrst*(cfrst+1)/2+1
cddim=(2*cfrst+nblv+1)*nblv/2
call dgemm('t','n',abdim,cddim,dfnbasis,1.d0,jab(1,p1),dfnbasis
$ ,jab(1,p2),dfnbasis,0.d0,abcd,itmp)
c save (ab|cd) to disk
do aa=1,nbl
a=afrst+aa
do b=1,a
ab=(afrst+a)*(aa-1)/2+b
do cc=1,nblv
c=cfrst+cc
do d=1,c
cd=(cfrst+c)*(cc-1)/2+d
temp3=abcd(ab,cd)
if(dabs(temp3).gt.itol)
$write(155,7000) temp3,nocc+a,nocc+b,nocc+c,nocc+d
enddo !d
enddo !cc
enddo !b
enddo !aa
enddo !cblock
c do last block, ablock=cblock
call dsyrk('u','t',abdim,dfnbasis,1.d0,jab(1,p1),dfnbasis,0.d0,
$abcd,itmp)
c save (ab|cd) to disk
do aa=1,nbl
a=afrst+aa
do b=1,a
ab=(afrst+a)*(aa-1)/2+b
do cc=1,nbl
c=afrst+cc
do d=1,c
cd=(afrst+c)*(cc-1)/2+d
if (cd.le.ab) then
temp3=abcd(cd,ab)
if(dabs(temp3).gt.itol)
$write(155,7000) temp3,nocc+a,nocc+b,nocc+c,nocc+d
endif
enddo !d
enddo !cc
enddo !b
enddo !aa
enddo !ablock
c }}}
endif !(d4.le.maxmem)
endif ! if (ccprog.eq.
write(556,'(6I24)') 0,0,0,0,0,0
c close files
if(ccprog.eq.'ccsd') then
close(556)
close(599)
close(600)
close(601)
close(602)
close(603)
close(604)
close(605)
chr if (hrprt) then
chr close(501)
chr close(502)
chr close(503)
chr close(504)
chr close(505)
chr close(497) ; close(498) ;close(499) ;close(500)
chr endif
endif
7000 format(e28.20,4i6)
c 7000 format(e14.6,4i6)
c stop 'SZEMET!!!!! itegrals are ready'
call dbldealloc(ijab)
c if (ccprog.eq.'ccsd') then
c write(6,*) 'Integral assembly finished'
c call timer ! if called makes an error in ccsd time measurement
c endif
return
end subroutine jpq2pqrs
subroutine ijkl_to_fort55(nocc,ijkl,itol,fort55)
implicit none
integer nocc,i,j,k,l,imin,ij,kl,fort55
real*8 ijkl(nocc*(nocc+1)/2,nocc*(nocc+1)/2),itol,temp3
do l=1,nocc
do j=l,nocc
do k=l,nocc
kl=k*(k-1)/2+l
if(l.eq.j) then
imin=k
elseif(l.ne.j) then ! i>=j
imin=j
endif
do i=imin,nocc
if(i.ge.j) then
ij=i*(i-1)/2+j
elseif(i.lt.j) then
ij=j*(j-1)/2+i
endif
if(kl.le.ij) then
temp3=ijkl(kl,ij)
if(dabs(temp3).gt.itol)
$ write(fort55,"(e28.20,4i6)") temp3,i,j,k,l
elseif(kl.gt.ij) then
temp3=ijkl(ij,kl)
if(dabs(temp3).gt.itol)
$ write(fort55,"(e28.20,4i6)") temp3,i,j,k,l
endif
enddo ! i
enddo ! k
enddo ! j
enddo ! l
return
end subroutine ijkl_to_fort55
subroutine aijk_to_fort55(aijk,nocc,nvirt,itol,fort55)
implicit none
integer nocc,nvirt,jk,j,k,i,a,fort55
real*8 temp3,aijk(nvirt,nocc,nocc*(nocc+1)/2),itol
do k=1,nocc
do i=1,nocc
do j=k,nocc
jk=j*(j-1)/2+k
do a=1,nvirt
temp3=aijk(a,i,jk)
if(dabs(temp3).gt.itol)
$ write(fort55,"(es28.20,4i6)") temp3,nocc+a,i,j,k
enddo ! a
enddo ! j
enddo ! i
enddo ! k
return
end subroutine aijk_to_fort55
subroutine aibj_to_fort55(nocc,nvirt,aibj,itol,fort55)
implicit none
integer nocc,nvirt,i,j,a,bmin,b,fort55
real*8 itol,aibj(nvirt,nocc,nvirt,nocc),temp3
do i=1,nocc
do j=i,nocc
do a=1,nvirt
if(i.eq.j) then
bmin=a
else ! i>j
bmin=1
endif ! (i.eq.j)
do b=bmin,nvirt
temp3=aibj(a,i,b,j)
if(dabs(temp3).gt.itol)
$ write(fort55,"(es28.20,4i6)") temp3,nocc+a,i,nocc+b,j
enddo ! a
enddo ! b
enddo ! j
enddo ! i
return
end subroutine aibj_to_fort55
subroutine abij_to_fort55(nocc,nvirt,abij,itol,fort55)
implicit none
integer nocc,nvirt,j,b,i,a,ij,ab,fort55
real*8 temp3,itol,abij(nvirt*(nvirt+1)/2,nvirt*(nvirt+1)/2)
do j=1,nocc
do b=1,nvirt
do i=j,nocc
ij=i*(i-1)/2+j
do a=b,nvirt
ab=a*(a-1)/2+b
temp3=abij(ab,ij)
if(dabs(temp3).gt.itol)
$ write(fort55,"(es28.20,4i6)") temp3,nocc+a,nocc+b,i,j
enddo ! a
enddo ! i
enddo ! b
enddo ! j
return
end subroutine abij_to_fort55
************************************************************************
subroutine jpq2pqrs_openshell(nal,nbe,nval,nvbe,dfnbasis,scr,
&iout,imem,imem1,maxcor,faadr,fbadr,fmadr,
&recaadr,recbadr,recmadr,recjab,irecln,nfa,nfb,
&nfm,localcc,ccsdalg,ccprog,itol,fort55)
************************************************************************
implicit none
integer nal,nbe,nval,nvbe,dfnbasis,iout,imem,imem1,maxcor,ijadim
integer ijklaa,i,j,k,l,a,ii,abadim,ncab,b,c,i0,ix,ijbdim,abbdim,ba
integer aiadim,aibdim,ijaab,ijbab,ijaai,ijbai,ijaij,ijbij,iscr,bi
integer jaab,jaij,d,kl,ij,bj,ai,ab,cd,ipuf,ik,jl,jk,il,aj,ak
integer ci,f0,nblc,fx,iabc,ac,iiscr,f,e,ae,be,irecln,recind
integer irec,nfa,nfb,nfm,start,ef,jj,efln,efmln,memneeded
integer recaadr(nfa*(nfa+1)/2),recbadr(nfb*(nfb+1)/2),fort55
integer recmadr(nfm*(nfm+1)/2),faadr(nfa+1),fbadr(nfb+1)
integer fmadr(nfm+1),ijbabblock,finblock,fln,recjab(nfm*(nfm+1)/2)
real*8 scr(*),itol
real*8,pointer::ijkl(:,:),aijk(:,:),abij(:,:),iabj(:,:),abci(:,:)
real*8,pointer::jaab2(:,:,:),jbai(:,:,:),iabc2(:,:,:,:)
real*8,pointer::jaai(:,:,:),jbab2(:,:,:),abef(:,:),jbab(:,:)
real*8,pointer::jbabblock(:,:,:),jbabtri(:,:)
logical error,iabc1blc
c local
character*4 localcc, ccprog
character*8 ccsdalg
c {{{ Interfaces for pointers
interface
subroutine rpoint2d(egydim,haromdim,dim1,dim2)
implicit none
integer dim1,dim2
real*8,target :: egydim(dim1,dim2)
real*8, pointer :: haromdim(:,:)
end subroutine
subroutine rpoint3d(egydim,haromdim,dim1,dim2,dim3)
implicit none
integer dim1,dim2,dim3
real*8,target :: egydim(dim1,dim2,dim3)
real*8, pointer :: haromdim(:,:,:)
end subroutine
subroutine rpoint4d(egydim,haromdim,dim1,dim2,dim3,dim4)
implicit none
integer dim1,dim2,dim3,dim4
real*8,target :: egydim(dim1,dim2,dim3,dim4)
real*8, pointer :: haromdim(:,:,:,:)
end subroutine
end interface
c }}}
error=.false.
c Partitioning memory
ijadim=nal*(nal+1)/2
ijbdim=nbe*(nbe+1)/2
abadim=nval*(nval+1)/2
abbdim=nvbe*(nvbe+1)/2
aiadim=nal*nval
aibdim=nbe*nvbe
ijaab=0
ijbab=ijaab+dfnbasis*abadim
ijaai=ijbab+dfnbasis*abbdim
ijbai=ijaai+dfnbasis*aiadim
ijaij=ijbai+dfnbasis*aibdim
ijbij=ijaij+dfnbasis*ijadim
iscr=ijbij+dfnbasis*ijbdim
c Memory check for ccprog.eq.'mrcc'
if (ccprog.eq.'mrcc') then
!Only need to check abci and abcd because everything else is the same as ccprog.ccsd
memneeded = iscr + max(
$max(abadim*aiadim,abbdim*aibdim,abadim*aibdim,abbdim*aiadim), !abci four spin cases
$max(abadim**2,abbdim**2,abadim*abbdim) !abcd three spin cases
$)
call memcheck(memneeded,maxcor-imem+imem1,iout,error)
if (error) then
write(iout,
$ "(' Insufficient memory for open-shell DF assembly.')")
call mrccend(1)
endif !error
endif !ccprog.eq.'mrcc'
c Reading J_ab
c alpha
open(16,file='DFINT_AB',form='unformatted')
read(16) scr(ijaab+1:ijbab)
close(16)
c beta
if (nbe.gt.0.and.nvbe.gt.0) then
open(16,file='DFINT_ABb',form='unformatted')
read(16) scr(ijbab+1:ijaai)
close(16)
endif
c Reading J_ai
c alpha
open(16,file='DFINT_AI',form='unformatted')
if (localcc.ne.'off ') then
read(16) scr(ijaai+1:ijbai)
call jaiinp(scr(ijaai+1),nal,nval,dfnbasis,scr(iscr+1))
else
read(16) scr(ijaai+1:ijbai)
endif
close(16)
c beta
if (nbe.gt.0.and.nvbe.gt.0) then
open(16,file='DFINT_AIb',form='unformatted')
if (localcc.ne.'off ') then
read(16) scr(ijbai+1:ijaij)
call jaiinp(scr(ijbai+1),nbe,nvbe,dfnbasis,scr(iscr+1))
else
read(16) scr(ijbai+1:ijaij)
endif
close(16)
endif
c Reading J_ij
c alpha
open(16,file='DFINT_IJ',form='unformatted')
read(16) scr(ijaij+1:ijbij)
close(16)
c beta
if (nbe.gt.0) then
open(16,file='DFINT_IJb',form='unformatted')
read(16) scr(ijbij+1:iscr)
close(16)
endif
if (ccsdalg.eq.'dfdirect'.and.ccprog.eq.'ccsd') then
c {{{ Create auxiliary file DFINT_ABbuccsd for ccsdalg=dfdirect
ijbabblock=iscr+1
call rpoint2d(scr(ijbab+1),jbab,dfnbasis,abbdim)
open(16,file='DFINT_ABbuccsd',access='direct',recl=irecln)
irec=1
recind=1
do i=1,nfm !ciklus az oszlopokra
fln=fmadr(i+1)-fmadr(i)
do j=1,i-1 !ciklus a negyzetes csempekre
call rpoint3d(scr(ijbabblock),jbabblock,dfnbasis,
& fmadr(j+1)-fmadr(j),fln)
do f=fmadr(i),fmadr(i+1)-1 !csempen beluli oszlopokra
finblock=f-fmadr(i)+1
do e=fmadr(j),fmadr(j+1)-1 !csempen beluli sorokra
jbabblock(1:dfnbasis,e-fmadr(j)+1,finblock)=
& jbab(1:dfnbasis,f*(f-1)/2+e)
enddo
enddo
call putlst(16,irec,jbabblock,
& (fmadr(j+1)-fmadr(j))*fln*dfnbasis) !size of the block (end-beg)
recjab(recind)=irec
recind=recind+1
if (mod((fmadr(j+1)-fmadr(j))*fln*dfnbasis,irecln/8).ne.0)
& then
irec=irec+(fmadr(j+1)-fmadr(j))*fln*dfnbasis/(irecln/8)+1
else
irec=irec+(fmadr(j+1)-fmadr(j))*fln*dfnbasis/(irecln/8)
endif
enddo !j=1,i-1
c Triangle tile
call rpoint2d(scr(ijbabblock),jbabtri,dfnbasis,fln*(fln+1)/2)
do f=fmadr(i),fmadr(i+1)-1
finblock=f-fmadr(i)+1
do e=fmadr(i),f
jbabtri(1:dfnbasis,finblock*(finblock-1)/2+e-fmadr(i)+1)=
& jbab(1:dfnbasis,f*(f-1)/2+e)
enddo
enddo
call putlst(16,irec,jbabtri,dfnbasis*fln*(fln+1)/2)
recjab(recind)=irec
recind=recind+1
if (mod(((fln+1)*fln/2)*dfnbasis,irecln/8).ne.0) then
irec=irec+((fln+1)*fln/2)*dfnbasis/(irecln/8)+1
else
irec=irec+((fln+1)*fln/2)*dfnbasis/(irecln/8)
endif
enddo !i=1,nfm
close(16)
c }}}
endif !ccsdalg.eq.dfdirect
c Assembling two-electron, four center integrals from DF three-center integrals
c Pure alpha integrals -------------------------------------------------
c (a<b|c<d)
if (nval.gt.1) then
if (ccprog.eq.'ccsd'.and.ccsdalg.ne.'dfdirect') then
call abcd_bld(nval,abadim,dfnbasis,ijaab,iscr,scr,irecln,
$ faadr,recaadr,nfa,'abcdt','abcda',error)
elseif (ccprog.eq.'mrcc') then
call abcd_bld_to_fort55(nal,nval,abadim,dfnbasis,scr(ijaab+1),
$ scr(iscr+1),itol,fort55)
endif !ccprog
endif !nval.gt.1.and.ccsdalg.ne.'dfdirect'
c (i<j|k<l)
if (nal.gt.1)
$ call ijkl_bld(ijadim,nal,dfnbasis,scr(ijaij+1),scr(iscr+1),
$ scr(iscr+ijadim**2+1),'ijklaa',ccprog,itol,fort55)
c (a,i|j<k)
if (nval.gt.0.and.nal.gt.1)
$ call aijk_bld(aiadim,ijadim,dfnbasis,nal,nval,scr(ijaai+1),
$ scr(ijaij+1),scr(iscr+1),scr(iscr+aiadim*ijadim+1),
$ 'aijkaa',ccprog,itol,fort55)
c (ai|bj) and (ab|ij)
if (nal.gt.0.and.nval.gt.0) then
ipuf=max(iscr+aiadim**2,iscr+abadim*ijadim)
call aibj_abij_bld(aiadim,abadim,ijadim,dfnbasis,scr(ijaai+1),
$ scr(ijaij+1),scr(ijaab+1),nal,nval,
$ scr(iscr+1),scr(iscr+1),scr(ipuf+1),ccprog,
$ itol,'abijaa','iabjaa',fort55)
endif !nal.gt.0.and.nval.gt.0
c (ab|ci)
if (ccprog.eq.'mrcc') then
if (nval.gt.1.and.nal.gt.0)
$ call abci_bld_to_fort55(nval,nal,dfnbasis,abadim,scr(ijaab+1),
$ scr(ijaai+1),scr(iscr+1),itol,.true.,
$ fort55)
write(fort55,"(e28.20,4i6)") 0.d0, 0, 0, 0, 0
endif !ccprog.eq.'mrcc'
c Pure beta integrals --------------------------------------------------
c <a<b||c<d>
if (nvbe.gt.1) then
if (ccprog.eq.'ccsd'.and.ccsdalg.ne.'dfdirect') then
call abcd_bld(nvbe,abbdim,dfnbasis,ijbab,iscr,scr,irecln,
$ fbadr,recbadr,nfb,'abcdt','abcdb',error)
elseif (ccprog.eq.'mrcc') then
call abcd_bld_to_fort55(nbe,nvbe,abbdim,dfnbasis,scr(ijbab+1),
$ scr(iscr+1),itol,fort55)
endif !ccprog
endif !nvbe.gt.1
c (i<j|k<l)
if (nbe.gt.1)
$ call ijkl_bld(ijbdim,nbe,dfnbasis,scr(ijbij+1),scr(iscr+1),
$ scr(iscr+ijbdim**2+1),'ijklbb',ccprog,itol,fort55)
c (a,i|j<k)
if (nvbe.gt.0.and.nbe.gt.1)
$ call aijk_bld(aibdim,ijbdim,dfnbasis,nbe,nvbe,scr(ijbai+1),
$ scr(ijbij+1),scr(iscr+1),scr(iscr+aibdim*ijbdim+1),
$ 'aijkbb',ccprog,itol,fort55)
c (ai|bj) and (ab|ij)
if (nbe.gt.0.and.nvbe.gt.0) then
ipuf=max(iscr+aibdim**2,iscr+abbdim*ijbdim)
call aibj_abij_bld(aibdim,abbdim,ijbdim,dfnbasis,scr(ijbai+1),
$ scr(ijbij+1),scr(ijbab+1),nbe,nvbe,
$ scr(iscr+1),scr(iscr+1),scr(ipuf+1),ccprog,
$ itol,'abijbb','iabjbb',fort55)
endif !nal.gt.0.and.nval.gt.0
c (ab|ci)
if (ccprog.eq.'mrcc') then
if (nvbe.gt.1.and.nbe.gt.0)
$ call abci_bld_to_fort55(nvbe,nbe,dfnbasis,abbdim,scr(ijbab+1),
$ scr(ijbai+1),scr(iscr+1),itol,.true.,
$ fort55)
write(fort55,"(e28.20,4i6)") 0.d0, 0, 0, 0, 0
endif !ccprog.eq.'mrcc'
c Mixed alpha-beta integrals -------------------------------------------
c <a<b||c<d>
if (nval.gt.0.and.nvbe.gt.0) then
if (ccprog.eq.'ccsd'.and.ccsdalg.ne.'dfdirect') then
c {{{ Build aBcD
recind=0
irec=1
do i=1,nfm !oszlopokra
open(16,file='abcdt',form='unformatted') !temp file (need less memory)
do f=fmadr(i),fmadr(i+1)-1 !egy oszlop beolvas
call dgemm('t','n',abadim,f,dfnbasis,1.d0,scr(ijaab+1),
& dfnbasis,scr(ijbab+(f*(f-1)/2)*dfnbasis+1),dfnbasis,0.d0,
& scr(iscr+1),abadim)
ii=0
call rpoint2d(scr(iscr+1),abef,abadim,f)
do e=1,nval
do b=1,f
do a=1,ncab(f,b,e,nval)
ii=ii+1
c calculating indeces
if (a.le.e) then
ae=e*(e-1)/2+a
else
ae=a*(a-1)/2+e
endif
c rearranging for printing
scr(iscr+abadim*f+ii)=abef(ae,b)
enddo !e
enddo !b
enddo !a
write(16) scr(iscr+abadim*f+1:iscr+abadim*f+ii)
enddo !f
open(17,file='abcdm',access='direct',recl=irecln)
rewind(16)
ii=0
do f=fmadr(i),fmadr(i+1)-1 !egy oszlop beolvas
read(16) scr(iscr+ii+1:iscr+ii+(nval**2*(2*f-1)+nval)/2)
ii=ii+(nval**2*(2*f-1)+nval)/2
enddo
c building the tiles
start=0
do j=1,i-1 !egesz, negyzet csempekre
jj=0
do ef=1,efmln(i,fmadr,nfm,nval)
do ab=1,efmln(j,fmadr,nfm,nval)
jj=jj+1
scr(iscr+ii+jj)=
& scr(iscr+start+(ef-1)*(fmadr(i)-1)*nval+(ef-1)*ef/2+ab)
enddo
enddo
start=start+efmln(j,fmadr,nfm,nval)
c writing tile
recind=recind+1
recmadr(recind)=irec
call putlst(17,irec,scr(iscr+ii+1),jj)
if (mod(jj,irecln/8).ne.0) then
irec=irec+jj/(irecln/8)+1
else
irec=irec+jj/(irecln/8)
endif
enddo
jj=0
do ef=1,efmln(i,fmadr,nfm,nval) !utolso, haromszog csempe
do ab=1,ef
jj=jj+1
scr(iscr+ii+jj)=
& scr(iscr+start+(ef-1)*(fmadr(i)-1)*nval+(ef-1)*ef/2+ab)
enddo
enddo
c writing tile
recind=recind+1
recmadr(recind)=irec
call putlst(17,irec,scr(iscr+ii+1),jj)
if (mod(jj,irecln/8).ne.0) then
irec=irec+jj/(irecln/8)+1
else
irec=irec+jj/(irecln/8)
endif
close(17)
close(16,status='delete')
enddo !nfm blokkokra
c }}}
elseif (ccprog.eq.'mrcc') then
call abcdm_bld_to_fort55(nal,nbe,nval,nvbe,abadim,abbdim,
$ dfnbasis,scr(ijaab+1),scr(ijbab+1),scr(iscr+1),itol,fort55)
endif !ccprog
endif !nal.gt.0.and.nvbe.gt.0.and.ccsdalg.ne.dfdirect
if (nal.gt.0.and.nbe.gt.0) then
c {{{ (i<j|k<l)
call dgemm('t','n',ijadim,ijbdim,dfnbasis,1.d0,scr(ijaij+1),
&dfnbasis,scr(ijbij+1),dfnbasis,0.d0,scr(iscr+1),ijadim)
call rpoint2d(scr(iscr+1),ijkl,ijadim,ijbdim)
ipuf=iscr+ijadim*ijbdim
c Placing (i,k|j,l) and (i,l|j,k) in <i,j||k,l>
ii=0
do l=1,nbe
do k=1,nal
do j=1,l
do i=1,ncab(l,j,k,nal)
ii=ii+1
c Coulomb indeces
if (i.le.k) then
ik=k*(k-1)/2+i
else
ik=i*(i-1)/2+k
endif
jl=l*(l-1)/2+j
if (ccprog.eq.'ccsd') then
scr(ipuf+ii)=ijkl(ik,jl)
elseif (ccprog.eq.'mrcc') then
write(fort55,"(es28.20,4i6)") ijkl(ik,jl),i,k,j,l
endif !ccprog
enddo
enddo
enddo
enddo
c Write <i,j||k,l>
if (ccprog.eq.'ccsd') then
open(16,file='ijklab',form='unformatted')
write(16) scr(ipuf+1:ipuf+ii)
close(16)
endif
c }}}
endif !nbe>0,nal>0
if (nval.gt.0.and.nal.gt.0.and.nbe.gt.0) then
c {{{ (aj|IK)
call dgemm('n','n',aiadim,ijbdim,dfnbasis,1.d0,scr(ijaai+1),
& aiadim,scr(ijbij+1),dfnbasis,0.d0,scr(iscr+1),aiadim)
call rpoint2d(scr(iscr+1),aijk,aiadim,ijbdim)
ipuf=iscr+aiadim*ijbdim
! Placing (a,j|i<k) in <a,i||j,k>
ii=0
do k=1,nbe
do j=1,nal
do i=1,nbe
do a=1,nval
ii=ii+1
! Coulomb indeces
aj=(j-1)*nval+a
if (i.le.k) then
ik=k*(k-1)/2+i
else
ik=i*(i-1)/2+k
endif
if (ccprog.eq.'ccsd') then
scr(ipuf+ii)=aijk(aj,ik)
elseif (ccprog.eq.'mrcc') then
write(fort55,"(es28.20,4i6)")aijk(aj,ik),nal+a,j,i,k
endif
enddo
enddo
enddo
enddo
! Writing aijk
if (ccprog.eq.'ccsd') then
open(16,file='aijkab',form='unformatted')
write(16) scr(ipuf+1:ipuf+ii)
close(16)
endif !ccprog ! }}}
endif !nval>0,nal>0,nbe>0
if (nvbe.gt.0.and.nal.gt.0.and.nbe.gt.0) then
c {{{ (ik|AJ)
call dgemm('t','t',ijadim,aibdim,dfnbasis,1.d0,scr(ijaij+1),
& dfnbasis,scr(ijbai+1),aibdim,0.d0,scr(iscr+1),ijadim)
call rpoint2d(scr(iscr+1),aijk,ijadim,aibdim)
ipuf=iscr+aibdim*ijadim
! Placing (a,j|i<k) in <a,i||j,k>
ii=0
do k=1,nal
do j=1,nbe
do i=1,nal
do a=1,nvbe
ii=ii+1
! Coulomb indeces
aj=(j-1)*nvbe+a
if (i.le.k) then
ik=k*(k-1)/2+i
else
ik=i*(i-1)/2+k
endif
if (ccprog.eq.'ccsd') then
scr(ipuf+ii)=aijk(ik,aj)
elseif (ccprog.eq.'mrcc') then
write(fort55,"(es28.20,4i6)") aijk(ik,aj),i,k,nbe+a,j
endif
enddo
enddo
enddo
enddo
! Writing aijk
if (ccprog.eq.'ccsd') then
open(16,file='aijkba',form='unformatted')
write(16) scr(ipuf+1:ipuf+ii)
close(16)
endif !ccprog !}}}
endif !nvbe>0,nal>0,nbe>0
if (nal.gt.0.and.nval.gt.0.and.nbe.gt.0.and.nvbe.gt.0) then
c {{{ (ai|BJ)
call dgemm('n','t',aiadim,aibdim,dfnbasis,1.d0,scr(ijaai+1),
& aiadim,scr(ijbai+1),aibdim,0.d0,scr(iscr+1),aiadim)
call rpoint2d(scr(iscr+1),abij,aiadim,aibdim)
ipuf=iscr+aiadim*aibdim
! Placing (a,i|b,j) in <a,b||i,j>
ii=0
do j=1,nbe
do i=1,nal
do b=1,nvbe
do a=1,nval
ii=ii+1
! Coulomb indeces
ai=(i-1)*nval+a
bj=(j-1)*nvbe+b
if (ccprog.eq.'ccsd') then
scr(ipuf+ii)=abij(ai,bj)
elseif (ccprog.eq.'mrcc') then
write(fort55,"(es28.20,4i6)") abij(ai,bj),nal+a,i,nbe+b,j
endif
enddo
enddo
enddo
enddo
! Writing abij
if (ccprog.eq.'ccsd') then
open(16,file='abijab',form='unformatted')
write(16) scr(ipuf+1:ipuf+ii)
close(16)
endif !ccprog !}}}
endif !nal,nval,nbe,nvbe>0
if (nval.gt.0.and.nbe.gt.0) then
c {{{ (ab|IJ)
call dgemm('t','n',abadim,ijbdim,dfnbasis,1.d0,scr(ijaab+1),
& dfnbasis,scr(ijbij+1),dfnbasis,0.d0,scr(iscr+1),abadim)
call rpoint2d(scr(iscr+1),iabj,abadim,ijbdim)
ipuf=iscr+abadim*ijbdim
! Placing (a,b|i,j) in <i,a||b,j>
ii=0
do j=1,nbe
do b=1,nval
do a=1,nval
if (b.le.a) then
ba=a*(a-1)/2+b
else
ba=b*(b-1)/2+a
endif
do i=1,nbe
ii=ii+1
! Exchange indeces
if (i.le.j) then
ij=j*(j-1)/2+i
else
ij=i*(i-1)/2+j
endif
if (ccprog.eq.'ccsd') then
scr(ipuf+ii)=-iabj(ba,ij)
elseif (ccprog.eq.'mrcc') then
write(fort55,"(es28.20,4i6)") iabj(ba,ij),nal+b,nal+a,i,j
endif
enddo
enddo
enddo
enddo
! Writing iabjbaab
if (ccprog.eq.'ccsd') then
open(16,file='iabjbaab',form='unformatted')
write(16) scr(ipuf+1:ipuf+ii)
close(16)
endif !ccprog !}}}
endif !nval,nbe>0
if (nal.gt.0.and.nvbe.gt.0) then
c {{{ (ij|AB)
call dgemm('t','n',ijadim,abbdim,dfnbasis,1.d0,scr(ijaij+1),
& dfnbasis,scr(ijbab+1),dfnbasis,0.d0,scr(iscr+1),ijadim)
call rpoint2d(scr(iscr+1),iabj,ijadim,abbdim)
ipuf=iscr+ijadim*abbdim
! Placing -(i,j|a,b) in <i,a||b,j>
ii=0
do j=1,nal
do b=1,nvbe
do a=1,nvbe
if (b.le.a) then
ba=a*(a-1)/2+b
else
ba=b*(b-1)/2+a
endif
do i=1,nal
ii=ii+1
! Exchange indeces
if (i.le.j) then
ij=j*(j-1)/2+i
else
ij=i*(i-1)/2+j
endif
if (ccprog.eq.'ccsd') then
scr(ipuf+ii)=-iabj(ij,ba)
elseif (ccprog.eq.'mrcc') then
write(fort55,"(es28.20,4i6)") iabj(ij,ba),i,j,nbe+b,nbe+a
endif
enddo
enddo
enddo
enddo
! Writing iabjabba
if (ccprog.eq.'ccsd') then
open(16,file='iabjabba',form='unformatted')
write(16) scr(ipuf+1:ipuf+ii)
close(16)
endif !ccprog !}}}
endif !nal,nvbe>0
if (ccprog.eq.'mrcc') then
if (nval.gt.1.and.nvbe.gt.0.and.nbe.gt.0)
$ call abcim_bld_to_fort55(nal,nbe,nval,nvbe,dfnbasis,abadim,
$ scr(ijaab+1),scr(ijbai+1),scr(iscr+1),itol,'ab',fort55)
if (nvbe.gt.1.and.nval.gt.0.and.nal.gt.0)
$ call abcim_bld_to_fort55(nbe,nal,nvbe,nval,dfnbasis,abbdim,
$ scr(ijbab+1),scr(ijaai+1),scr(iscr+1),itol,'ba',fort55)
write(fort55,"(e28.20,4i6)") 0.d0, 0, 0, 0, 0
endif !ccprog.eq.'mrcc'
if (ccprog.eq.'ccsd') then
c {{{ Assembling the (ab|ci) integrals for ccsd in a memory efficient way
iscr=ijaij !delete stuff starting from Jij (we don't need it for (ab|ci))
! <i,c||a<b>
if (max(abadim*aibdim+nval*nvbe*nval*nbe,
& abbdim*aiadim+nvbe*nvbe*nval*nal,
& abadim*aiadim+(nval*(nval-1)/2)*nval*nal,
& abbdim*aibdim+(nvbe*(nvbe-1)/2)*nvbe*nbe)
& +iscr.lt.maxcor-imem+imem1) then
iabc1blc=.true.
call rpoint3d(scr(ijaai+1),jaai,nval,nal,dfnbasis)
call rpoint3d(scr(ijbai+1),jbai,nvbe,nbe,dfnbasis)
else !1blc
iabc1blc=.false.
! Rearranging J_AIb in place
call jaiinp(scr(ijbai+1),nvbe,nbe,dfnbasis,scr(iscr+1))
call rpoint3d(scr(ijbai+1),jbai,nbe,nvbe,dfnbasis)
! Rearranging J_AIa in place
call jaiinp(scr(ijaai+1),nval,nal,dfnbasis,scr(iscr+1))
call rpoint3d(scr(ijaai+1),jaai,nal,nval,dfnbasis)
! Unboxing J_ABa
call jabunb(scr(ijaab+1),nval,dfnbasis,scr(iscr+1))
call rpoint3d(scr(iscr+1),jaab2,dfnbasis,nval,nval)
iiscr=iscr+dfnbasis*nval**2
endif !1blc
c {{{ b,a,a,b
! Building iabc
if (nbe.gt.0.and.nvbe.gt.0.and.nval.gt.0) then
if (iabc1blc) then
call dgemm('t','t',abadim,aibdim,dfnbasis,1.d0,scr(ijaab+1),
& dfnbasis,jbai,aibdim,0.d0,scr(iscr+1),abadim)
call rpoint2d(scr(iscr+1),abci,abadim,aibdim)
ipuf=iscr+aibdim*abadim
! Placing (a<c|b,i) and (a<b|c,i) in <i,a||b<c>
open(16,file='iabcba',form='unformatted')
do c=1,nvbe
do b=1,nval
ii=0
do a=1,nval
do i=1,nbe
ii=ii+1
! Exchange indeces
ci=(i-1)*nvbe+c
if (a.le.b) then
ab=b*(b-1)/2+a
else
ab=a*(a-1)/2+b
endif
scr(ipuf+ii)=abci(ab,ci)
enddo
enddo
! Writing <i,a||b<c>
write(16) scr(ipuf+1:ipuf+ii)
enddo
enddo
close(16)
else !1blc
open(16,file='iabcba',form='unformatted')
f0=0
do while (f0.lt.nvbe)
nblc=0
fx=0
do while (nblc*nval*nbe+fx*nbe*nval**2+iiscr.lt.
& maxcor-imem+imem1.and.fx+f0.le.nvbe)
fx=fx+1
nblc=nblc+nval
enddo
nblc=nblc-nval
fx=fx-1
call dgemm('n','n',nbe*fx,nval**2,dfnbasis,1.d0,
& jbai(1,f0+1,1),nbe*nvbe,jaab2(1,1,1),dfnbasis,0.d0,scr(iiscr+1),
& nbe*fx)
call rpoint4d(scr(iiscr+1),iabc2,nbe,fx,nval,nval)
iabc=iiscr+nbe*fx*nval**2
do c=1,fx
do b=1,nval
do a=1,nval
do i=1,nbe
ii=(c-1)*nval*nval*nbe+(b-1)*nval*nbe+(a-1)*nbe+i
scr(iabc+ii)=iabc2(i,c,a,b)
enddo
enddo
enddo
enddo
! Writing iabc block
do i=1,nblc
write(16) scr(iabc+(i-1)*nval*nbe+1:iabc+i*nval*nbe)
enddo
f0=f0+fx
enddo
close(16)
endif !1blc
endif !nbe>0,nvbe>0,nval>0
c }}}
! a,a,a,a
if (nal.gt.0.and.nval.gt.1) then
call iabc1(nal,nval,dfnbasis,abadim,aiadim,'iabcaa',jaai,
& ijaab,iscr,scr,maxcor,imem,imem1,iabc1blc,error)
endif !nal>0,nval>0
if (.not.iabc1blc) then
! Unboxing J_ABb
call jabunb(scr(ijbab+1),nvbe,dfnbasis,scr(iscr+1))
call rpoint3d(scr(iscr+1),jbab2,dfnbasis,nvbe,nvbe)
iiscr=iscr+dfnbasis*nvbe**2
endif
c {{{ a,b,a,b
! a,b,a,b
if (nal.gt.0.and.nvbe.gt.0.and.nval.gt.0) then
if (iabc1blc) then
call dgemm('t','t',abbdim,aiadim,dfnbasis,1.d0,scr(ijbab+1),
& dfnbasis,jaai,aiadim,0.d0,scr(iscr+1),abbdim)
call rpoint2d(scr(iscr+1),abci,abbdim,aiadim)
ipuf=iscr+aiadim*abbdim
! Placing (a<c|b,i) and (a<b|c,i) in <i,a||b<c>
open(16,file='iabcab',form='unformatted')
do c=1,nvbe
do b=1,nval
ii=0
do a=1,nvbe
do i=1,nal
ii=ii+1
! Coulomb indeces
bi=(i-1)*nval+b
if (a.le.c) then
ac=c*(c-1)/2+a
else
ac=a*(a-1)/2+c
endif
scr(ipuf+ii)=abci(ac,bi)
enddo
enddo
! Writing <i,a||b<c>
write(16) scr(ipuf+1:ipuf+ii)
enddo
enddo
close(16)
else !1blc
open(16,file='iabcab',form='unformatted')
f0=0
do while (f0.lt.nvbe)
nblc=0
fx=0
do while (nblc*nvbe*nal+fx*nvbe*nal*nval+iiscr.lt.
& maxcor-imem+imem1.and.fx+f0.le.nvbe)
fx=fx+1
nblc=nblc+nval
enddo
nblc=nblc-nval
fx=fx-1
call dgemm('n','n',nal*nval,fx*nvbe,dfnbasis,1.d0,jaai,
& nal*nval,jbab2(1,1,f0+1),dfnbasis,0.d0,scr(iiscr+1),nal*nval)
call rpoint4d(scr(iiscr+1),iabc2,nal,nval,nvbe,fx)
iabc=iiscr+nal*nval*nvbe*fx
do c=1,fx
do b=1,nval
do a=1,nvbe
do i=1,nal
ii=(c-1)*nval*nvbe*nal+(b-1)*nvbe*nal+(a-1)*nal+i
scr(iabc+ii)=iabc2(i,b,a,c)
enddo
enddo
enddo
enddo
! Writing iabc block
do i=1,nblc
write(16) scr(iabc+(i-1)*nvbe*nal+1:iabc+i*nvbe*nal)
enddo
f0=f0+fx
enddo
close(16)
endif !1blc
endif !nal>0,nval>0,nvbe>0
c }}}
! b,b,b,b
if (nbe.gt.0.and.nvbe.gt.1) then
call iabc1(nbe,nvbe,dfnbasis,abbdim,aibdim,'iabcbb',jbai,
& ijbab,iscr,scr,maxcor,imem,imem1,iabc1blc,error)
endif !nbe>0,nvbe>1
c}}}
endif !ccprog.eq.'ccsd'
if (error) then
write(iout,*)
write(iout,"(' Error during DF integral assembly!')")
call mrccend(1)
endif
return
end subroutine jpq2pqrs_openshell
************************************************************************
subroutine jpq2aibj_uccsd(nal,nbe,nval,nvbe,dfnbasis,scr,
&iout,imem,imem1,maxcor,faadr,fbadr,fmadr,
&recaadr,recbadr,recmadr,recjab,irecln,nfa,nfb,
&nfm,localcc,error,ccsdalg)
************************************************************************
implicit none
integer nal,nbe,nval,nvbe,dfnbasis,iout,imem,imem1,maxcor,ijadim
integer ijklaa,i,j,k,l,a,ii,abadim,ncab,b,c,i0,ix,ijbdim,abbdim,ba
integer aiadim,aibdim,ijaab,ijbab,ijaai,ijbai,ijaij,ijbij,iscr,bi
integer jaab,jaij,d,kl,ij,bj,ai,ab,cd,ipuf,ik,jl,jk,il,aj,ak
integer ci,f0,nblc,fx,iabc,ac,iiscr,f,e,ae,be,irecln,recind
integer irec,nfa,nfb,nfm,start,ef,jj,efln,efmln
integer recaadr(nfa*(nfa+1)/2),recbadr(nfb*(nfb+1)/2)
integer recmadr(nfm*(nfm+1)/2),faadr(nfa+1),fbadr(nfb+1)
integer fmadr(nfm+1),ijbabblock,finblock,fln,recjab(nfm*(nfm+1)/2)
real*8 scr(*)
real*8,pointer::ijkl(:,:),aijk(:,:),abij(:,:),iabj(:,:),abci(:,:)
real*8,pointer::jaab2(:,:,:),jbai(:,:,:),iabc2(:,:,:,:)
real*8,pointer::jaai(:,:,:),jbab2(:,:,:),abef(:,:),jbab(:,:)
real*8,pointer::jbabblock(:,:,:),jbabtri(:,:)
logical error,iabc1blc
c local
character*4 localcc
character*8 ccsdalg
c {{{ Interfaces for pointers
interface
subroutine rpoint2d(egydim,haromdim,dim1,dim2)
implicit none
integer dim1,dim2
real*8,target :: egydim(dim1,dim2)
real*8, pointer :: haromdim(:,:)
end subroutine
subroutine rpoint3d(egydim,haromdim,dim1,dim2,dim3)
implicit none
integer dim1,dim2,dim3
real*8,target :: egydim(dim1,dim2,dim3)
real*8, pointer :: haromdim(:,:,:)
end subroutine
subroutine rpoint4d(egydim,haromdim,dim1,dim2,dim3,dim4)
implicit none
integer dim1,dim2,dim3,dim4
real*8,target :: egydim(dim1,dim2,dim3,dim4)
real*8, pointer :: haromdim(:,:,:,:)
end subroutine
end interface
c }}}
c Partitioning memory
ijadim=nal*(nal+1)/2
ijbdim=nbe*(nbe+1)/2
abadim=nval*(nval+1)/2
abbdim=nvbe*(nvbe+1)/2
aiadim=nal*nval
aibdim=nbe*nvbe
ijaab=0
ijbab=ijaab+dfnbasis*abadim
iscr=ijbab+dfnbasis*abbdim
ijaai=ijbab+dfnbasis*abbdim
ijbai=ijaai+dfnbasis*aiadim
ijaij=ijbai+dfnbasis*aibdim
ijbij=ijaij+dfnbasis*ijadim
iscr=ijbij+dfnbasis*ijbdim
c Reading J_ai
c alpha
open(16,file='DFINT_AI_NONAF',form='unformatted')
if (localcc.ne.'off ') then
read(16) scr(ijaai+1:ijbai)
call jaiinp(scr(ijaai+1),nal,nval,dfnbasis,scr(iscr+1))
else
read(16) scr(ijaai+1:ijbai)
endif
close(16)
c beta
if (nbe.gt.0.and.nvbe.gt.0) then
open(16,file='DFINT_AIb_NONAF',form='unformatted')
if (localcc.ne.'off ') then
read(16) scr(ijbai+1:ijaij)
call jaiinp(scr(ijbai+1),nbe,nvbe,dfnbasis,scr(iscr+1))
else
read(16) scr(ijbai+1:ijaij)
endif
close(16)
endif
c Assembling the integral lists <1,2||1,2>
c Building (a,i|b,j)
c a,a,a,a
if (nal.gt.0.and.nval.gt.0) then
call dsyrk('u','n',aiadim,dfnbasis,1.d0,scr(ijaai+1),aiadim,0.d0,
&scr(iscr+1),aiadim)
call rpoint2d(scr(iscr+1),abij,aiadim,aiadim)
ipuf=max(iscr+aiadim**2,iscr+abadim*ijadim)
c Placing (a,i|b,j) and (a,j|b,i) in <a<b||i<j>
ii=0
do j=1,nal
do i=1,j-1
do b=1,nval
do a=1,b-1
ii=ii+1
c Coulomb indeces
ai=(i-1)*nval+a
bj=(j-1)*nval+b
c Exchange indeces
aj=(j-1)*nval+a
bi=(i-1)*nval+b
scr(ipuf+ii)=abij(ai,bj)-abij(bi,aj)
enddo
enddo
enddo
enddo
c Writing abij
open(16,file='abijaa_nonaf',form='unformatted')
write(16) scr(ipuf+1:ipuf+ii)
close(16)
c Placing (b,i|a,j) in <i,a||b,j>
ii=0
do j=1,nal
do b=1,nval
do a=1,nval
do i=1,nal
ii=ii+1
c Coulomb indeces
bj=(j-1)*nval+b
ai=(i-1)*nval+a
if (ai.le.bj) then
scr(ipuf+ii)=-abij(ai,bj)
else
scr(ipuf+ii)=-abij(bj,ai)
endif
enddo
enddo
enddo
enddo
endif !nval>0,nal>0
c b,b,b,b
if (nvbe.gt.0.and.nbe.gt.0) then
call dsyrk('u','n',aibdim,dfnbasis,1.d0,scr(ijbai+1),aibdim,0.d0,
&scr(iscr+1),aibdim)
call rpoint2d(scr(iscr+1),abij,aibdim,aibdim)
ipuf=max(iscr+aibdim**2,iscr+abbdim*ijbdim)
c Placing (a,i|b,j) and (a,j|b,i) in <a<b||i<j>
ii=0
do j=1,nbe
do i=1,j-1
do b=1,nvbe
do a=1,b-1
ii=ii+1
c Coulomb indeces
ai=(i-1)*nvbe+a
bj=(j-1)*nvbe+b
c Exchange indeces
aj=(j-1)*nvbe+a
bi=(i-1)*nvbe+b
scr(ipuf+ii)=abij(ai,bj)-abij(bi,aj)
enddo
enddo
enddo
enddo
c Writing aijk
open(16,file='abijbb_nonaf',form='unformatted')
write(16) scr(ipuf+1:ipuf+ii)
close(16)
endif !nvbe>0,nbe>0
c a,b,a,b
if (nal.gt.0.and.nval.gt.0.and.nbe.gt.0.and.nvbe.gt.0) then
call dgemm('n','t',aiadim,aibdim,dfnbasis,1.d0,scr(ijaai+1),
&aiadim,scr(ijbai+1),aibdim,0.d0,scr(iscr+1),aiadim)
call rpoint2d(scr(iscr+1),abij,aiadim,aibdim)
ipuf=iscr+aiadim*aibdim
c Placing (a,i|b,j) in <a,b||i,j>
ii=0
do j=1,nbe
do i=1,nal
do b=1,nvbe
do a=1,nval
ii=ii+1
c Coulomb indeces
ai=(i-1)*nval+a
bj=(j-1)*nvbe+b
scr(ipuf+ii)=abij(ai,bj)
enddo
enddo
enddo
enddo
c Writing abij
open(16,file='abijab_nonaf',form='unformatted')
write(16) scr(ipuf+1:ipuf+ii)
close(16)
endif !nal,nval,nbe,nvbe>0
return
end
************************************************************************
subroutine iabc1(nal,nval,dfnbasis,abadim,aiadim,iabcaa,jaai,
&ijaab,iscr,scr,maxcor,imem,imem1,iabc1blc,error)
************************************************************************
c UCCSD
c Writes iabcaa file for DF calculations.
implicit none
integer nal,nval,ijaab,iiscr,abadim,aiadim,dfnbasis,iscr
integer ipuf,iiscr2,i,c,b,a,ii,f0,fx,nblc,ac,bi,ci,ab,maxcor,imem
integer imem1,iabc
real*8 scr(*),jaai(nval*nal*dfnbasis)
real*8,pointer::abci(:,:),jab2(:,:,:),jai2(:,:,:)
logical error,iabc1blc
character*6 iabcaa
c {{{ Interfaces for pointers
interface
subroutine rpoint2d(egydim,haromdim,dim1,dim2)
implicit none
integer dim1,dim2
real*8,target :: egydim(dim1,dim2)
real*8, pointer :: haromdim(:,:)
end subroutine
subroutine rpoint3d(egydim,haromdim,dim1,dim2,dim3)
implicit none
integer dim1,dim2,dim3
real*8,target :: egydim(dim1,dim2,dim3)
real*8, pointer :: haromdim(:,:,:)
end subroutine
subroutine rpoint4d(egydim,haromdim,dim1,dim2,dim3,dim4)
implicit none
integer dim1,dim2,dim3,dim4
real*8,target :: egydim(dim1,dim2,dim3,dim4)
real*8, pointer :: haromdim(:,:,:,:)
end subroutine
end interface
c}}}
iiscr=iscr
c Building iabc
if (iabc1blc) then
call dgemm('t','t',abadim,aiadim,dfnbasis,1.d0,scr(ijaab+1),
&dfnbasis,jaai,aiadim,0.d0,scr(iiscr+1),abadim)
call rpoint2d(scr(iiscr+1),abci,abadim,aiadim)
ipuf=iiscr+aiadim*abadim
c Placing (a<c|b,i) and (a<b|c,i) in <i,a||b<c>
open(16,file=iabcaa,form='unformatted')
do c=1,nval
do b=1,c-1
ii=0
do a=1,nval
do i=1,nal
ii=ii+1
c Coulomb indeces
bi=(i-1)*nval+b
if (a.le.c) then
ac=c*(c-1)/2+a
else
ac=a*(a-1)/2+c
endif
c Exchange indeces
ci=(i-1)*nval+c
if (a.le.b) then
ab=b*(b-1)/2+a
else
ab=a*(a-1)/2+b
endif
scr(ipuf+ii)=abci(ab,ci)-abci(ac,bi)
enddo
enddo
c Writing <i,a||b<c>
write(16) scr(ipuf+1:ipuf+ii)
enddo
enddo
close(16)
else !1blc
c Unpacked Jab pointer
call rpoint3d(scr(iiscr+1),jab2,dfnbasis,nval,nval)
iiscr=iiscr+dfnbasis*nval**2
open(16,file=iabcaa,form='unformatted')
f0=0
do while (f0.lt.nval)
nblc=0
fx=0
do while (nblc*nval*nal+2*fx*nal*nval*(f0+fx-1)+iiscr.lt.
&maxcor-imem+imem1.and.fx+f0.le.nval)
fx=fx+1
nblc=nblc+f0+fx-1
enddo
nblc=nblc-f0-fx+1
fx=fx-1
iiscr2=iiscr+nal*(f0+fx-1)*fx*nval
iabc=iiscr2+nal*(f0+fx-1)*fx*nval
if (f0+fx-1.gt.0) then
call dgemm('n','n',nal*(f0+fx-1),fx*nval,dfnbasis,1.d0,jaai,
&nal*nval,jab2(1,1,f0+1),dfnbasis,0.d0,scr(iiscr+1),nal*(f0+fx-1))
endif
call dgemm('n','n',nal*fx,nval*(f0+fx-1),dfnbasis,1.d0,
&jaai(f0*nal+1),nal*nval,jab2(1,1,1),dfnbasis,0.d0,scr(iiscr2+1),
&nal*fx)
do c=1,fx
do b=1,f0+c-1
do a=1,nval
do i=1,nal
ii=(f0*(c-1)+(c-1)*(c-2)/2+b-1)*nval*nal+(a-1)*nal+i
scr(iabc+ii)=-scr(iiscr+((c-1)*nval+a-1)*(f0+fx-1)*nal+
&(b-1)*nal+i)+scr(iiscr2+((b-1)*nval+a-1)*fx*nal+(c-1)*nal+i)
enddo
enddo
enddo
enddo
c Writing iabc block
do i=1,nblc
write(16) scr(iabc+(i-1)*nval*nal+1:iabc+i*nval*nal)
enddo
f0=f0+fx
enddo
close(16)
endif !1blc
return
end subroutine iabc1
************************************************************************
************************************************************************
subroutine abcd_bld(nval,abadim,dfnbasis,ijaab,iscr,scr,irecln,
&faadr,recaadr,nfa,abcdtmp,abcda,error)
************************************************************************
c UCCSD
c Assembles and writes the abcd alfa and beta integs.
implicit none
integer nval,abadim,dfnbasis,ijaab,a,b,e,f,ii,jj,i,j,iscr,ncab,ae
integer be,irecln,recind,irec,nfa,faadr(nfa+1),start,ef,efln,ab
integer recaadr(nfa*(nfa+1)/2)
real*8 scr(*)
real*8,pointer::abef(:,:)
character*5 abcdtmp,abcda
logical error
c {{{ Interface for pointer
interface
subroutine rpoint2d(egydim,haromdim,dim1,dim2)
implicit none
integer dim1,dim2
real*8,target :: egydim(dim1,dim2)
real*8, pointer :: haromdim(:,:)
end subroutine
end interface
c}}}
recind=0
irec=1
do i=1,nfa !oszlopokra
open(16,file=abcdtmp,form='unformatted')
do f=faadr(i),faadr(i+1)-1 !egy oszlop beolvas
call dgemm('t','n',f*(f+1)/2,f,dfnbasis,1.d0,scr(ijaab+1),
& dfnbasis,scr(ijaab+(f*(f-1)/2)*dfnbasis+1),dfnbasis,0.d0,
& scr(iscr+1),f*(f+1)/2)
call rpoint2d(scr(iscr+1),abef,f*(f+1)/2,f)
c write(*,*) "newabef",abef
ii=0
do e=1,f-1
do b=1,f
do a=1,ncab(f,b,e,b-1)
ii=ii+1
c calculating indeces
if (b.eq.f) then ! I. eset (b=f,a<=e)
ae=e*(e-1)/2+a !coulomb
c exchange b=f tehat b>e
be=b*(b-1)/2+e
elseif (a.le.e) then ! III. eset (b<f,a<=e)
ae=e*(e-1)/2+a !coulomb
c exchange
if (b.le.e) then ! III.a eset (b<=e)
be=e*(e-1)/2+b
else ! III.b eset (b>e)
be=b*(b-1)/2+e
endif
else ! II. eset (b<f,a>e)
c coulomb
ae=a*(a-1)/2+e
c exchange b>a>e tehat b>e
be=b*(b-1)/2+e
endif
c rearranging for printing
scr(iscr+f*(f+1)*f/2+ii)=abef(ae,b)-abef(be,a)
c write(*,"('new',f14.10)")scr(iscr+f*(f+1)*f/2+ii)
enddo
enddo
enddo
c printing
write(16) scr(iscr+f*(f+1)*f/2+1:iscr+f*(f+1)*f/2+ii)
enddo !f
open(17,file=abcda,access='direct',recl=irecln)
rewind(16)
ii=0
do f=faadr(i),faadr(i+1)-1 !egy oszlop beolvas
read(16) scr(iscr+ii+1:iscr+ii+(f**3-3*f**2+4*f-2)/2)
ii=ii+(f**3-3*f**2+4*f-2)/2
enddo
c building the tiles
start=0
do j=1,i-1 !egesz, negyzet csempekre
jj=0
do ef=1,efln(i,faadr,nfa)
do ab=1,efln(j,faadr,nfa)
jj=jj+1
scr(iscr+ii+jj)=
&scr(iscr+start+(ef-1)*(faadr(i)-1)*(faadr(i)-2)/2+(ef-1)*ef/2+ab)
enddo
enddo
start=start+efln(j,faadr,nfa)
c writing tile
recind=recind+1
recaadr(recind)=irec
call putlst(17,irec,scr(iscr+ii+1),jj)
c write(*,"('new',f14.10)")scr(iscr+ii+1:iscr+ii+jj)
if (mod(jj,irecln/8).ne.0) then
irec=irec+jj/(irecln/8)+1
else
irec=irec+jj/(irecln/8)
endif
enddo
jj=0
do ef=1,efln(i,faadr,nfa) !utolso, haromszog csempe
do ab=1,ef
jj=jj+1
scr(iscr+ii+jj)=
&scr(iscr+start+(ef-1)*(faadr(i)-1)*(faadr(i)-2)/2+(ef-1)*ef/2+ab)
enddo
enddo
c writing tile
recind=recind+1
recaadr(recind)=irec
call putlst(17,irec,scr(iscr+ii+1),jj)
c write(*,"('new',f14.10)")scr(iscr+ii+1:iscr+ii+jj)
if (mod(jj,irecln/8).ne.0) then
irec=irec+jj/(irecln/8)+1
else
irec=irec+jj/(irecln/8)
endif
close(17)
close(16,status='delete')
enddo !nfa blokkokra
return
end subroutine abcd_bld
************************************************************************
************************************************************************
function efmln(i,fmadr,nf,nval)
************************************************************************
c UCCSD
c Returns the length of the ith ef block of abcd integral (ab).
implicit none
integer i,nf,fmadr(nf+1),efmln,nval
efmln=(fmadr(i+1)-fmadr(i))*nval
return
end function efmln
************************************************************************
************************************************************************
function ncab(outd,inb,outc,na)
************************************************************************
* outd - outer, slower index *
* inb - inner, slower index *
* outc - outer, faster index *
* na - max value of the fastest index *
* *
* Function for determining the range of the fastest index of (ab,cd) *
* type integrals. *
************************************************************************
implicit none
integer outd,inb,outc,na,ncab
if (inb.eq.outd) then
ncab=outc
else
ncab=na
endif
return
end function ncab
************************************************************************
************************************************************************
subroutine jaiinp(jai,na,ni,np,puf)
************************************************************************
c UCCSD
c Transposes Jai in place
c J(a,i,P)->J(i,a,P)
implicit none
integer na,ni,np,a,i,p
real*8 jai(na*ni,np),puf(ni*na)
do p=1,np
do i=1,ni
do a=1,na
puf((a-1)*ni+i)=jai((i-1)*na+a,p)
enddo
enddo
jai(1:ni*na,p)=puf(1:ni*na)
enddo
return
end subroutine jaiinp
************************************************************************
************************************************************************
subroutine jabunb(jab,na,np,jabu)
************************************************************************
c UCCSD
c Expands Jab index
c J(P,a<b)->J(P,a,b)
implicit none
integer na,b,a,np,i,p
real*8 jab(np,na*(na+1)/2),jabu(np,na*na)
do b=1,na
do a=1,b
c do p=1,np
c jabu(p,(b-1)*na+a)=jab(p,b*(b-1)/2+a)
jabu(1:np,(b-1)*na+a)=jab(1:np,b*(b-1)/2+a)
c if (a.ne.b) jabu(p,(a-1)*na+b)=jab(p,b*(b-1)/2+a)
if (a.ne.b) jabu(1:np,(a-1)*na+b)=jab(1:np,b*(b-1)/2+a)
c enddo
enddo
enddo
return
end subroutine jabunb
************************************************************************
************************************************************************
function efln(i,faadr,nf)
************************************************************************
c UCCSD
c Returns the length of the ith ef block. (aa,bb)
implicit none
integer i,nf,faadr(nf+1),efln
efln=(((faadr(i+1)-2)+(faadr(i)-1))*(faadr(i+1)-faadr(i)))/2
return
end function efln
************************************************************************
************************************************************************
subroutine ijkl_bld(ijadim,nal,dfnbasis,jaij,ijkl,work,intfile,
$ ccprog,itol,fort55)
c Builds and writes the (ij|kl) or <ik||jl> integrals
************************************************************************
implicit none
integer ijadim,dfnbasis,iscr,ipuf,ncab,i,ii,j,k,l,nal,jl,ik,jk,il
integer fort55
real*8 work(*),jaij(dfnbasis,ijadim),ijkl(ijadim,ijadim),itol
character*4 ccprog
character*6 intfile
call dsyrk('u','t',ijadim,dfnbasis,1.d0,jaij,dfnbasis,
&0.d0,ijkl,ijadim)
if (ccprog.eq.'ccsd') then
c Placing (i<k|j<l) and (i<l|j<k) in <i<j|k<l>
ii=0
do l=1,nal
do k=1,l-1
do j=1,l
do i=1,ncab(l,j,k,j-1)
ii=ii+1
c Coulomb indeces
jl=l*(l-1)/2+j
if (i.le.k) then
ik=k*(k-1)/2+i
else
ik=i*(i-1)/2+k
endif
c Exchange indeces
if (j.le.k) then
jk=k*(k-1)/2+j
else
jk=j*(j-1)/2+k
endif
il=l*(l-1)/2+i
if (ik.le.jl) then
if (il.le.jk) then
work(ii)=ijkl(ik,jl)-ijkl(il,jk)
else
work(ii)=ijkl(ik,jl)-ijkl(jk,il)
endif
else
if (il.le.jk) then
work(ii)=ijkl(jl,ik)-ijkl(il,jk)
else
work(ii)=ijkl(jl,ik)-ijkl(jk,il)
endif
endif
enddo
enddo
enddo
enddo
c Write <i<j||k<l>
open(16,file=intfile,form='unformatted')
write(16) work(1:ii)
close(16)
elseif (ccprog.eq.'mrcc') then
call ijkl_to_fort55(nal,ijkl,itol,fort55)
endif !ccprog
return
end subroutine ijkl_bld
************************************************************************
************************************************************************
subroutine aijk_bld(aiadim,ijadim,dfnbasis,nal,nval,jaai,jaij,
$ aijk,work,intfile,ccprog,itol,fort55)
c Builds and writes the (aj|ik) or <ai||jk> integrals
************************************************************************
implicit none
integer aiadim,ijadim,dfnbasis,nal,nval,i,j,k,a,ii,aj,ak,ik,ij
integer fort55
real*8 jaai(aiadim,dfnbasis),jaij(dfnbasis,aiadim)
real*8 aijk(aiadim,ijadim),work(*),itol
character*4 ccprog
character*6 intfile
call dgemm('n','n',aiadim,ijadim,dfnbasis,1.d0,jaai,
$aiadim,jaij,dfnbasis,0.d0,aijk,aiadim)
if (ccprog.eq.'ccsd') then
c Placing (a,j|i<k) and (a,k|i<j) in <a,i||j<k>
ii=0
do k=1,nal
do j=1,k-1
do i=1,nal
do a=1,nval
ii=ii+1
aj=(j-1)*nval+a
ak=(k-1)*nval+a
! Coulomb indeces
if (i.le.k) then
ik=k*(k-1)/2+i
else
ik=i*(i-1)/2+k
endif
! Exchange indeces
if (i.le.j) then
ij=j*(j-1)/2+i
else
ij=i*(i-1)/2+j
endif
work(ii)=aijk(aj,ik)-aijk(ak,ij)
enddo
enddo
enddo
enddo
! Writing aijk
open(16,file=intfile,form='unformatted')
write(16) work(1:ii)
close(16)
elseif (ccprog.eq.'mrcc') then
call aijk_to_fort55(aijk,nal,nval,itol,fort55)
endif !ccprog
return
end subroutine aijk_bld
************************************************************************
************************************************************************
subroutine aibj_abij_bld(aiadim,abadim,ijadim,dfnbasis,jaai,jaij,
$jaab,nal,nval,abij,iabj,work,ccprog,itol,abijfile,iabjfile,fort55)
c Builds and writes the (ai|bj) and (ab|ij) or <ab||ij> and <ia||bj>integrals
************************************************************************
implicit none
integer aiadim,abadim,ijadim,dfnbasis,nal,nval,i,j,a,b,ii
integer ba,ij,fort55
real*8 jaai(aiadim,dfnbasis),jaij(dfnbasis,ijadim)
real*8 iabj(abadim,ijadim),itol
real*8 jaab(dfnbasis,abadim),abij(nval,nal,nval,nal),work(*)
character*4 ccprog
character*6 abijfile,iabjfile
call dsyrk('u','n',aiadim,dfnbasis,1.d0,jaai,aiadim,0.d0,
&abij,aiadim)
if (ccprog.eq.'ccsd') then
! Placing (a,i|b,j) and (a,j|b,i) in <a<b||i<j>
ii=0
do j=1,nal
do i=1,j-1
do b=1,nval
do a=1,b-1
ii=ii+1
work(ii)=abij(a,i,b,j)-abij(b,i,a,j)
enddo
enddo
enddo
enddo
! Writing abij
open(16,file=abijfile,form='unformatted')
write(16) work(1:ii)
close(16)
elseif (ccprog.eq.'mrcc') then
call aibj_to_fort55(nal,nval,abij,itol,fort55)
endif !ccprog
if (ccprog.eq.'ccsd') then
! Placing (b,i|a,j) in <i,a||b,j>
ii=0
do j=1,nal
do b=1,nval
do a=1,nval
do i=1,nal
ii=ii+1
if ((i-1)*nval+a.le.(j-1)*nval+b) then
work(ii)=-abij(a,i,b,j)
else
work(ii)=-abij(b,j,a,i)
endif
enddo
enddo
enddo
enddo
endif !ccprog
! Building (a<b|i<j)
call dgemm('t','n',abadim,ijadim,dfnbasis,1.d0,jaab,
&dfnbasis,jaij,dfnbasis,0.d0,iabj,abadim)
if (ccprog.eq.'ccsd') then
! Placing (a<b|i<j) in <i,a||b,j>
ii=0
do j=1,nal
do b=1,nval
do a=1,nval
do i=1,nal
ii=ii+1
! Exchange indeces
if (i.le.j) then
ij=j*(j-1)/2+i
else
ij=i*(i-1)/2+j
endif
if (b.le.a) then
ba=a*(a-1)/2+b
else
ba=b*(b-1)/2+a
endif
work(ii)=work(ii)+iabj(ba,ij)
enddo
enddo
enddo
enddo
! Writing <i,a||b,j>
open(16,file=iabjfile,form='unformatted')
write(16) work(1:ii)
close(16)
elseif (ccprog.eq.'mrcc') then
call abij_to_fort55(nal,nval,iabj,itol,fort55)
endif !ccprog
return
end subroutine aibj_abij_bld
************************************************************************
************************************************************************
subroutine abci_bld_to_fort55(nvirt,nocc,dfnbasis,abdim,jab,jai,
$abci,itol,assemble,fort55)
c Builds and writes the (ab|ci) integrals to fort.55
************************************************************************
implicit none
integer abdim,nvirt,nocc,dfnbasis,i,b,c,a,ab,fort55
real*8 jab(dfnbasis,abdim),jai(nvirt,nocc,dfnbasis)
real*8 abci(abdim,nvirt,nocc),temp3,itol
logical assemble
if (assemble)
$ call dgemm('t','t',abdim,nvirt*nocc,dfnbasis,1.d0,jab,dfnbasis,
$ jai,nvirt*nocc,0.d0,abci,abdim)
do i=1,nocc
do b=1,nvirt
do c=1,nvirt
do a=b,nvirt
ab=a*(a-1)/2+b
temp3=abci(ab,c,i)
if(dabs(temp3).gt.itol)
$ write(fort55,"(es28.20,4i6)") temp3,nocc+a,nocc+b,nocc+c,i
enddo ! a
enddo ! c
enddo ! b
enddo ! i
return
end subroutine abci_bld_to_fort55
************************************************************************
************************************************************************
subroutine abcim_bld_to_fort55(nal,nbe,nval,nvbe,dfnbasis,abdim,
$ jaab,jbai,abci,itol,spincase,fort55)
c Builds and writes the (ab|CI) (or (AB|ci) integrals to fort.55
************************************************************************
implicit none
integer nbe,nval,nvbe,abdim,aidim,a,b,c,i,ab,nal,dfnbasis,fort55
real*8 jaab(dfnbasis,abdim),jbai(nvbe,nbe,dfnbasis)
real*8 abci(abdim,nvbe,nbe),temp3,itol
character*2 spincase
call dgemm('t','t',abdim,nvbe*nbe,dfnbasis,1.d0,jaab,dfnbasis,
$ jbai,nvbe*nbe,0.d0,abci,abdim)
do i=1,nbe
do b=1,nval
do c=1,nvbe
do a=b,nval
ab=a*(a-1)/2+b
temp3=abci(ab,c,i)
if(dabs(temp3).gt.itol) then
if (spincase.eq.'ab') then
write(fort55,"(es28.20,4i6)") temp3,nal+a,nal+b,nbe+c,i
elseif (spincase.eq.'ba') then
write(fort55,"(es28.20,4i6)") temp3,nbe+c,i,nal+a,nal+b
endif !spincase
endif !itol
enddo ! a
enddo ! c
enddo ! b
enddo ! i
return
end subroutine abcim_bld_to_fort55
************************************************************************
************************************************************************
subroutine abcd_bld_to_fort55(nocc,nvirt,abdim,dfnbasis,jab,abcd,
$ itol,fort55)
c Builds and writes the (ab|cd) integrals to fort.55
************************************************************************
implicit none
integer nvirt,abdim,dfnbasis,a,b,c,d,ab,cd,nocc,fort55
real*8 jab(abdim,dfnbasis),abcd(abdim,abdim),temp3,itol
call dsyrk('u','t',abdim,dfnbasis,1.d0,jab,dfnbasis,0.d0,
$ abcd,abdim)
do a=1,nvirt
do b=1,a
ab=a*(a-1)/2+b
do c=1,nvirt
do d=1,c
cd=c*(c-1)/2+d
if (cd.le.ab) then
temp3=abcd(cd,ab)
if(dabs(temp3).gt.itol)
$ write(fort55,"(es28.20,4i6)")temp3,nocc+a,nocc+b,nocc+c,nocc+d
endif
enddo
enddo
enddo
enddo
return
end subroutine abcd_bld_to_fort55
************************************************************************
************************************************************************
subroutine abcdm_bld_to_fort55(nal,nbe,nval,nvbe,abadim,abbdim,
$ dfnbasis,jaab,jbab,abcd,itol,fort55)
c Builds and writes the (ab|CD) integrals to fort.55
************************************************************************
implicit none
integer nal,nbe,nval,nvbe,abadim,abbdim,dfnbasis,a,b,c,d,ab,cd
integer fort55
real*8 jaab(dfnbasis,abadim),jbab(dfnbasis,abbdim)
real*8 abcd(abadim,abbdim),temp3,itol
call dgemm('t','n',abadim,abbdim,dfnbasis,1.d0,jaab,dfnbasis,jbab,
$ dfnbasis,0.d0,abcd,abadim)
cd=0
do d=1,nvbe
do c=1,d
cd=cd+1
ab=0
do b=1,nval
do a=1,b
ab=ab+1
temp3=abcd(ab,cd)
if(dabs(temp3).gt.itol)
$ write(fort55,"(es28.20,4i6)")temp3,nal+a,nal+b,nbe+c,nbe+d
enddo
enddo
enddo
enddo
return
end subroutine abcdm_bld_to_fort55
************************************************************************
c {{{ Routines for memory requirement management
subroutine memcheck(needed,available,iout,error)
c Checks if the available memory is sufficient, writes warning to iout
implicit none
integer needed,available,iout
real*8 adjav,adjreq
logical error
character*3 memreq,memav
if (needed.gt.available) then
c Converting memory into sensible units
call memconv(needed,adjreq,memreq)
call memconv(available,adjav,memav)
c Writing error message
write(iout,*)
write(iout,*) 'Insufficient memory!'
write(iout,"(' Requested: ',f6.2,3a)")
&adjreq,memreq
write(iout,"(' Available: ',f6.2,3a)")
&adjav,memav
error=.true.
endif
return
end
subroutine memconv(memword,memvalue,memunit)
c Converts memory value from words (8 bytes) to the appropriate unit (B,KB,MB,GB,TB)
implicit none
integer memword
real*8 memvalue
character*3 memunit
select case(memword*8)
case (1:512)
memvalue=dble(memword*8)
memunit=' B '
case (513:512000)
memvalue=dble(memword*8)/dble(1024)
memunit=' KB'
case(512001:524288000)
memvalue=dble(memword*8)/dble(1024**2)
memunit=' MB'
case(524288001:536870912000)
memvalue=dble(memword*8)/dble(1024**3)
memunit=' GB'
case(536870912001:)
memvalue=dble(memword*8)/dble(1024**4)
memunit=' TB'
end select
return
end subroutine
c}}}
c {{{ subroutine rpoint1d rpoint2d rpoint3d rpoint4d
subroutine rpoint1d(egydim1,egydim2,dim1)
implicit none
integer dim1
real*8,target :: egydim1(dim1)
real*8, pointer :: egydim2(:)
egydim2 => egydim1(1:dim1)
end subroutine
subroutine rpoint2d(egydim1,egydim2,dim1,dim2)
implicit none
integer dim1,dim2
real*8,target :: egydim1(dim1,dim2)
real*8, pointer :: egydim2(:,:)
egydim2 => egydim1(1:dim1,1:dim2)
end subroutine
subroutine rpoint3d(egydim1,egydim2,dim1,dim2,dim3)
implicit none
integer dim1,dim2,dim3
real*8,target :: egydim1(dim1,dim2,dim3)
real*8, pointer :: egydim2(:,:,:)
egydim2 => egydim1(1:dim1,1:dim2,1:dim3)
end subroutine
subroutine rpoint4d(egydim1,egydim2,dim1,dim2,dim3,dim4)
implicit none
integer dim1,dim2,dim3,dim4
real*8,target :: egydim1(dim1,dim2,dim3,dim4)
real*8, pointer :: egydim2(:,:,:,:)
egydim2 => egydim1(1:dim1,1:dim2,1:dim3,1:dim4)
end subroutine
subroutine ipoint1d(egydim1,egydim2,dim1)
implicit none
integer dim1
integer, target :: egydim1(dim1)
integer, pointer :: egydim2(:)
egydim2 => egydim1(1:dim1)
end subroutine
subroutine ipoint2d(egydim1,egydim2,dim1,dim2)
implicit none
integer dim1,dim2
integer, target :: egydim1(dim1,dim2)
integer, pointer :: egydim2(:,:)
egydim2 => egydim1(1:dim1,1:dim2)
end subroutine
c }}}