mirror of
https://code.it4i.cz/sccs/easyconfigs-it4i.git
synced 2025-04-16 19:50:50 +01:00
1819 lines
67 KiB
Fortran
1819 lines
67 KiB
Fortran
************************************************************************
|
|
subroutine locdrv(intcpu,intwal,orbcpu,orbwal,ovirtcpud,ovirtwald,
|
|
$scfcpu,scfwal,localcc,calctype)
|
|
************************************************************************
|
|
* Driver for local CC calculations
|
|
************************************************************************
|
|
#include "MRCCCOMMON"
|
|
c taken from MRCCCOMMON: inp,minpfile,ifcfile,scrfile1,scrfile2,iout,eref,nvirt
|
|
c +vars in MRCCCOMMON has to be set/initialized before locdrv call: nocc,nbasis
|
|
integer nbasisold,noccold,iocc,i,idegen,nocorr,nbcorr,n1,n2,n3
|
|
integer nom,nbm,nocorrmin,nocorrmax,nbcorrmin,nbcorrmax,maxmet
|
|
integer imet,nmet,ndom,idom,ncmo,cmo,iccsd,idomst,nbasiss,noccs
|
|
integer idumy,ncoi,nact,nmetact,ccmem(8),ipart,imp2,nob,noccolda
|
|
integer noccoldb,ncorep,numcalcEDs,methact_inited
|
|
parameter(maxmet=30)
|
|
real*8 ecor,totcpu(maxmet),totwal(maxmet),ecl(maxmet),escf
|
|
real*8 mp2corr,ovirtcpu,ovirtwal,lmp2ene,scfcpu,scfwal
|
|
real*8 ovirtcpud,ovirtwald,lmp2ened,intcpu,intwal,intcpud,intwald
|
|
real*8 times(6),edenergy(9)
|
|
real*8 cccpu,ccwal,orbcpu,orbwal,mulcpu,mulwal,eenv,pairpllno
|
|
real*8 tcpu,twal,mp2cpu,mp2wal,mp2fact
|
|
character*3 lmp2dens,half,csapprox
|
|
character*5 scftype
|
|
character*4 ccprog,localcc,accprog,eccprog,setccprog,talg
|
|
character*6 core,cscr,ccorr,totcorr,dir
|
|
character*128 label,meth(maxmet),beg,end,ccmo,methact(maxmet)
|
|
logical log1,lseqc,mp2env,hfenv,lrc,ccenv,lanycc,localcc15p
|
|
logical laptoldefault,cycle_this_domain,domain_is_cs
|
|
logical localccXp,mp2fact_halved,domainlist,cc_restart
|
|
character*16 naf,envcalc,cscr16,calctype,uenvcalc,alaptol,elaptol
|
|
character*8 naftyp,lccrest,corembed,ccsdalg,lccoporder,scfalg
|
|
character*8 ednumread
|
|
integer,allocatable::calc_the_ed(:),dlist(:)
|
|
c previously in MRCCCOMMON
|
|
integer nalcp,nbecp
|
|
c mpi
|
|
integer mpi_rank,mpi_size
|
|
logical group_master_thread,inbcast_group
|
|
mpi_rank=0
|
|
mpi_size=1
|
|
dir=' .'
|
|
group_master_thread=.true. ! ??? todo
|
|
inbcast_group=.true. ! ??? todo
|
|
C Type of correction
|
|
ccorr='MP2 '
|
|
half=' '
|
|
mp2fact=1.0d0 ! MP2 correction scale for all but CCSD, for that mp2fact=0.5
|
|
mp2fact_halved=.false.
|
|
lmp2ened=0.d0 ! initialization in case of ovirt call
|
|
C Execute mulli
|
|
if(localcc.eq.'2013') then
|
|
call runit('mulli',mulcpu,mulwal)
|
|
call ishell('rm -f 55.*')
|
|
elseif (localcc15p(localcc)) then
|
|
mulcpu=0.d0
|
|
mulwal=0.d0
|
|
endif
|
|
C
|
|
call getkey('core',4,core,6)
|
|
call getkey('ccprog',6,ccprog,4)
|
|
call getkey('lmp2dens',8,lmp2dens,3)
|
|
call getkey('scfalg',6,scfalg,8)
|
|
call getvar('lseqc ',lseqc)
|
|
call getkey('lccrest',7,lccrest,8)
|
|
call getvar('eref ',escf)
|
|
call getkey('naf_cor',7,naf,16)
|
|
call getkey('naftyp',6,naftyp,8)
|
|
call getkey('talg',4,talg,4)
|
|
call getkey('corembed',8,corembed,8)
|
|
call getkey('ccsdalg',7,ccsdalg,8)
|
|
call getkey('lccoporder',10,lccoporder,8)
|
|
call getvar('nal ',nalcp) !for oslcc
|
|
call getvar('nbe ',nbecp) !for oslcc
|
|
call getkey('scftype',7,scftype,5) !for oslcc
|
|
call getkey('csapprox',8,csapprox,3)
|
|
call getvar('ncore ',ncore)
|
|
ncorep=ncore
|
|
if (trim(core).eq.'corr') ncorep=0
|
|
nmetact=1 ! number of invoked methods (for the active region if corembed)
|
|
if (.not.(corembed.ne.'off ')) envcalc='none '
|
|
c {{{ read/initialize localcc.restart
|
|
if (lccrest.eq.'on '.or.lccoporder.eq.'lccfirst') then
|
|
open(unit=scrfile2,file='localcc.restart',status='old')
|
|
call localccrestart('read',scrfile2,calctype,meth,nmet,
|
|
$ cc_restart,nbasis,nocc,ncore,idomst,nmet,
|
|
$ nocorr,nbcorr,nocorrmin,nocorrmax,nbcorrmin,nbcorrmax,nact,
|
|
$ ecl,totcpu,
|
|
$ totwal,lmp2ene,mp2corr,pairpllno,eenv,
|
|
$ edenergy(1),edenergy(2),edenergy(3),edenergy(4),edenergy(5),
|
|
$ edenergy(6),edenergy(7),edenergy(8),edenergy(9),ccmem
|
|
$ ,intcpu,intwal,scfcpu,scfwal,orbcpu,orbwal,ovirtcpu,ovirtwal)
|
|
if (lccoporder.eq.'lccfirst') then ! do not overwrite timings with uninitilaized values of localcc.restart
|
|
ovirtcpu=ovirtcpud ; ovirtwal=ovirtwald
|
|
endif
|
|
idomst=idomst+1
|
|
write(iout,*)
|
|
write(iout,122)
|
|
write(iout,122)
|
|
write(iout,*)
|
|
else
|
|
idomst=1
|
|
edenergy(1:9)=0.d0
|
|
c read corrections from the ED calculations
|
|
if(localcc15p(localcc)) then
|
|
open(unit=ifcfile,file='iface',status='old')
|
|
read(ifcfile,*) lmp2ene,mp2corr,eenv
|
|
close(ifcfile)
|
|
pairpllno=mp2corr
|
|
if (naftyp.eq.'jpqpilot'.or.naftyp.eq.'jpq ') then
|
|
mp2corr=lmp2ene-eenv ! second term of mp2corr is read from the result file of ccsd for all idom values
|
|
endif
|
|
else
|
|
lmp2ene=lmp2ened
|
|
mp2corr=0.d0
|
|
eenv=0.d0
|
|
pairpllno=mp2corr
|
|
endif
|
|
c inititalize localcc.restart for restarting with the first LIS
|
|
open(unit=scrfile2,file='localcc.restart')
|
|
rewind(scrfile2)
|
|
nmet=min(maxmet,30)
|
|
ecl(1:nmet)=0.d0
|
|
totcpu(1:nmet)=0.d0
|
|
totwal(1:nmet)=0.d0
|
|
ovirtcpu=ovirtcpud
|
|
ovirtwal=ovirtwald
|
|
ccmem(1:8)=0
|
|
call localccrestart('init',scrfile2,calctype,meth,nmet,
|
|
$ .false.,nbasis,nocc,ncore,0,nmet,
|
|
$ 0,0,nocc,0,nbasis,0,nocc,
|
|
$ ecl,totcpu,
|
|
$ totwal,lmp2ene,mp2corr,pairpllno,eenv,
|
|
$ 0.d0,0.d0,0.d0,0.d0,0.d0,0.d0,0.d0,0.d0,0.d0,ccmem,
|
|
$ intcpu,intwal,scfcpu,scfwal,orbcpu,orbwal,ovirtcpu,ovirtwal )
|
|
endif
|
|
c }}}
|
|
nbasiss=nbasis
|
|
noccs=nocc-ncorep
|
|
nbasisold=nbasis/2
|
|
noccolda=nalcp-ncorep
|
|
noccoldb=nbecp-ncorep
|
|
noccold=noccolda+noccoldb
|
|
nbasisold=nbasisold-ncorep
|
|
methact_inited=0
|
|
C Read the number of atomic domains
|
|
if(localcc15p(localcc)) then
|
|
ndom=nalcp-ncorep
|
|
else
|
|
open(unit=scrfile1,file='DOMAIN')
|
|
rewind(scrfile1)
|
|
read(scrfile1,*) ndom
|
|
close(scrfile1)
|
|
endif
|
|
if (idomst-1.eq.ndom) then
|
|
if (lccoporder.ne.'lccfirst')
|
|
$ write(iout,*) "Restarting local correlation computation",
|
|
$ " results for all domains are read from file: localcc.restart"
|
|
n1=1 ; n2=1 ; n3=1 ; label='ENERGY '
|
|
read(scrfile2,'(30a128)') (meth(imet),imet=1,nmet)
|
|
go to 101 ! all LIS calculations are completed, jump to final result print
|
|
endif
|
|
c
|
|
if (lccrest.ne.'on ') then
|
|
call lccorbnuminit(nocorr,nbcorr,nocorrmin,nocorrmax,nbcorrmin,
|
|
$ nbcorrmax,noccold,nbasisold)
|
|
elseif (lccrest.eq.'on ') then
|
|
write(iout,*) "Restarting local CC computation",
|
|
$ " with LIS:", idomst
|
|
write(iout,*) " results from previous domains are read from",
|
|
$ " file: localcc.restart"
|
|
endif
|
|
c
|
|
if (corembed.ne.'off ') then
|
|
call setenvcalc(envcalc,cscr16,minpfile,iout,lrc)
|
|
hfenv=trim(envcalc).eq.'scf'
|
|
mp2env=trim(envcalc).eq.'mp2'
|
|
ccenv=lanycc(envcalc)
|
|
if (ccenv) then
|
|
accprog=ccprog ! ccprog for the active orbitals
|
|
eccprog=setccprog(envcalc) ! ccprog for the environment orbitals
|
|
call getkey('laptol',6,alaptol,16) ! laptol for active orbitals
|
|
elaptol='1.d-2 ' ! default laptol for the environment orbitals warning make it consistend with the default in minp.f
|
|
laptoldefault=.false. ! is default set for laptol?
|
|
if (alaptol.eq.elaptol) laptoldefault=.true.
|
|
endif
|
|
nact=0
|
|
else
|
|
mp2env=.false.
|
|
hfenv=.false.
|
|
ccenv=.false.
|
|
lrc=.true.
|
|
nact=noccolda
|
|
laptoldefault=.true. ! not used
|
|
endif
|
|
c Domainlist keyword
|
|
domainlist=.false.
|
|
open(minpfile,file='MINP')
|
|
call getkeym('domainlist',10,ednumread,8)
|
|
if(ednumread.ne.' ') then
|
|
read(ednumread,"(i8)") numcalcEDs
|
|
if (numcalcEDs.le.ndom.and.numcalcEDs.ge.idomst) then
|
|
domainlist=.true.
|
|
allocate(calc_the_ed(ndom),dlist(numcalcEDs))
|
|
read(minpfile,*) (dlist(i),i=1,numcalcEDs)
|
|
calc_the_ed(1:ndom)=0
|
|
do i=1,numcalcEDs
|
|
calc_the_ed(dlist(i))=1
|
|
enddo
|
|
else
|
|
call unknown('domainlist',10)
|
|
endif ! numcalcEDs
|
|
endif !ednumread.ne.''
|
|
close(minpfile)
|
|
C Loop over domains
|
|
do idom=idomst,ndom
|
|
C Initialize files and variables
|
|
c {{{ localcc.eq.'2013'
|
|
if(localcc.eq.'2013') then
|
|
write(end,"(i9)") idom
|
|
end=adjustl(end)
|
|
beg='DOMAIN.' // end
|
|
beg='mv ' // trim(beg) // ' DOMAIN'
|
|
call ishell(beg)
|
|
beg='MOCOEF.' // end
|
|
beg='mv ' // trim(beg) // ' MOCOEF'
|
|
call ishell(beg)
|
|
beg='FOCK.' // end
|
|
beg='mv ' // trim(beg) // ' FOCK'
|
|
call ishell(beg)
|
|
open(unit=scrfile1,file='DOMAIN',form='unformatted')
|
|
rewind(scrfile1)
|
|
read(scrfile1) nbasis,nocc
|
|
read(scrfile1)
|
|
read(scrfile1) ncmo
|
|
c write(6,*) 'nbasis,nocc,ncmo',nbasis,nocc,ncmo
|
|
C Calculate integrals
|
|
if(trim(scfalg).ne.'disk'.or.ndom.gt.1.or..not.lseqc) then
|
|
call runit('integ',intcpud,intwald)
|
|
intcpu=intcpu+intcpud
|
|
intwal=intwal+intwald
|
|
endif
|
|
C Execute ovirt
|
|
call runit('ovirt',ovirtcpud,ovirtwald)
|
|
open(unit=ifcfile,file='iface',status='old')
|
|
read(ifcfile,*) lmp2ened
|
|
close(ifcfile)
|
|
lmp2ene=lmp2ene+lmp2ened
|
|
ovirtcpu=ovirtcpu+ovirtcpud
|
|
ovirtwal=ovirtwal+ovirtwald
|
|
else
|
|
ncmo=1
|
|
endif
|
|
c }}}
|
|
C Execute CC calculation for each occupied MO
|
|
do iocc=1,ncmo
|
|
ncoi=0 ! initilize domain as active
|
|
if(localcc.eq.'2013') then
|
|
read(scrfile1) (cmo,i=1,iocc)
|
|
backspace(scrfile1)
|
|
else
|
|
cmo=idom
|
|
endif
|
|
|
|
if (domainlist) then
|
|
if (calc_the_ed(idom).eq.0) then
|
|
write(iout,*)
|
|
write(iout,"(' CC is not calculated for occupied orbital ',
|
|
$ i7,' ...')") cmo
|
|
write(iout,*)'Due to domainlist keyword.'
|
|
cycle
|
|
endif
|
|
endif
|
|
|
|
domain_is_cs=.true.
|
|
write(ccmo,'(i6)') cmo
|
|
ccmo=adjustl(ccmo)
|
|
beg='55.' // trim(ccmo)
|
|
beg=adjustl(beg)
|
|
inquire(file=trim(beg),exist=log1)
|
|
if (log1) then
|
|
open(inp,status='old',file=trim(beg))
|
|
read(inp,*) nbasis,nocc,idumy,idegen,idumy,ncoi,nob
|
|
if (scftype.eq.'rohf ') then
|
|
domain_is_cs=.false.
|
|
!if (idegen.ne.0) then
|
|
!do i=1,6+nbasis**2
|
|
! read(inp,*)
|
|
!enddo
|
|
!read(inp,*) nob
|
|
if (2*nob.eq.nocc.and.trim(csapprox).eq.'on') then
|
|
!if (nob.eq.-1) then
|
|
! Closed-shell approximation in this domain
|
|
!nob=nocc/2
|
|
domain_is_cs=.true.
|
|
else if (trim(lccrest).eq.'on'.and.idom.eq.idomst) then
|
|
! Try to restart UCCSD calculation in domain
|
|
call changekey('ccsdrest',8,'ccsd',4)
|
|
endif
|
|
!endif !idegen.ne.0
|
|
else !scftype.ne.rohf
|
|
nob=nocc/2
|
|
endif !scftype
|
|
rewind(inp)
|
|
nvirt=2*nbasis-nocc
|
|
endif !log1
|
|
if (corembed.eq.'off ') ncoi=0 ! initilize domain as active
|
|
c
|
|
cycle_this_domain=(ncoi.eq.1.and.(mp2env.or.hfenv)) ! MP2/HF environment OR 0 virtual orbitals OR symmetry equivalent
|
|
$ .or.nvirt.eq.0.or.idegen.eq.0
|
|
c
|
|
c prepartion of CC calculation in the LIS
|
|
call preliscc(1,corembed,ncoi,inp,idegen,nact,ccenv,
|
|
$ccprog,accprog,eccprog,nbasis,calctype,envcalc,log1,
|
|
$cycle_this_domain,laptoldefault,elaptol,ifcfile,localcc,-1,
|
|
$domain_is_cs,nocc-nob,nob,iout)
|
|
c
|
|
if(log1.and.cycle_this_domain) then
|
|
write(iout,*)
|
|
write(iout,122)
|
|
write(iout,122)
|
|
write(iout,"(' CC is not calculated for occupied orbital ',
|
|
$ i7,' ...')") cmo
|
|
write(iout,*)
|
|
if (idegen.eq.0) then
|
|
write(iout,*)'Symmetry equivalent LIS is already completed'
|
|
else
|
|
write(iout,*) 'Inactive LIS: 0 correlation E contribution'
|
|
endif
|
|
write(iout,122)
|
|
cycle ! CC computation is not performed in the case of HF/MP2 environment model
|
|
endif
|
|
c
|
|
c perform CC calculation in the LIS
|
|
if(log1) then
|
|
if(ccprog.eq.'ccsd') then
|
|
call locfilecp(ccmo,'55 ',dir)
|
|
if(localcc.eq.'2013') then
|
|
call locfilemv(ccmo,'56 ')
|
|
call locfilemv(ccmo,'abcd.a ')
|
|
call locfilemv(ccmo,'abcd.s ')
|
|
call locfilemv(ccmo,'abci ')
|
|
call locfilemv(ccmo,'abij ')
|
|
call locfilemv(ccmo,'aibj ')
|
|
call locfilemv(ccmo,'aijk ')
|
|
call locfilemv(ccmo,'ijkl ')
|
|
else
|
|
call locfilecopy(0,cmo,naf,talg,nob,domain_is_cs,
|
|
$ mpi_rank,mpi_size,iout,inbcast_group)
|
|
endif
|
|
endif ! ccprog.eq.'ccsd'
|
|
end='cp ' // trim(beg) // ' fort.55'
|
|
call ishell(end)
|
|
open(inp,status='old',file='fort.55')
|
|
read(inp,*) nbasis,nocc,i,idegen
|
|
close(inp)
|
|
nbasis=2*nbasis
|
|
c nvirt=nbasis-nocc
|
|
C
|
|
write(iout,*)
|
|
write(iout,122)
|
|
write(iout,122)
|
|
write(iout,"(' Executing CC calculations for occupied orbital',
|
|
$i4,' ...')") cmo
|
|
write(iout,*)
|
|
write(iout,"(' Number of correlated electrons: ',i4,3x,
|
|
$'(',f5.1,'%)')") nocc,100.d0*dble(nocc)/dble(noccold)
|
|
write(iout,"(' Number of correlated basis functions:',i4,3x,
|
|
$'(',f5.1,'%)')") nbasis/2,100.d0*dble(nbasis/2)/dble(nbasisold)
|
|
write(iout,122)
|
|
c
|
|
c start the CC calculation in the LIS
|
|
call ddmrcc(ccprog,domain_is_cs,.true.)
|
|
c
|
|
c process the iface file, restore keywords, save to localcc.restart
|
|
c if(trim(scftype).ne.'rhf') nocc=nocc/2+nob
|
|
call postliscc(1,nocc,nbasis,meth,label,n1,n2,n3,eref,naftyp
|
|
$,idegen,methact_inited,corembed,ncoi,ifcfile,
|
|
$laptoldefault,alaptol,localcc,ccprog,maxmet, ! for localcc.restart
|
|
$scrfile2,calctype,methact,
|
|
$nmetact,nbasiss,noccs,ncore,idom,nmet,
|
|
$nocorr,nbcorr,nocorrmin,nocorrmax,nbcorrmin,nbcorrmax,nact,
|
|
$ecl,totcpu,totwal,lmp2ene,mp2corr,pairpllno,eenv,
|
|
$edenergy(1),edenergy(2),edenergy(3),edenergy(4),edenergy(5),
|
|
$edenergy(6),edenergy(7),edenergy(8),edenergy(9),ccmem,
|
|
$intcpu,intwal,scfcpu,scfwal,orbcpu,orbwal,ovirtcpu,ovirtwal,.true.
|
|
$,-1,iout,group_master_thread)
|
|
c
|
|
c cleanup files of LIS idom
|
|
if (localcc15p(localcc).and.ccprog.eq.'ccsd')
|
|
$ call locfiledel(cmo,naf,talg)
|
|
endif !log1
|
|
enddo ! idom
|
|
if(localcc.eq.'2013') close(scrfile1) ! close DOMAIN file
|
|
enddo !idom
|
|
c Deallocate variables used for domainlist keyword
|
|
if (domainlist) deallocate(calc_the_ed,dlist)
|
|
c
|
|
close(scrfile2) ! close localcc.restart file
|
|
c reset variables to their defaults in case of local embedding
|
|
if (corembed.ne.'off '.and.ccenv) then
|
|
ccprog=accprog
|
|
if (ccprog.eq.'mrcc')
|
|
$ call write56(nbasis,'none ',calctype,.true.,nocc/2,
|
|
$ nocc/2)
|
|
call changekey('calc',4,calctype,16)
|
|
nmet=nmetact
|
|
meth(1:nmet)=methact(1:nmet)
|
|
nact=noccolda
|
|
endif
|
|
c cleanup temporary files of localcc
|
|
end='rm -f DFINT_AI DFINT_AB DFINT_IJ 55'
|
|
if (naf.ne.'off ') end=trim(end) // ' ajb'
|
|
if (talg.eq.'lapl'.or.talg.eq.'topr'.or.talg.eq.'to '
|
|
$ .or.talg.eq.'lato') end=trim(end) // ' laplbas'
|
|
call ishell(end)
|
|
if (ccprog.eq.'ccsd'.and.ccsdalg.eq.'disk ') then
|
|
end='rm -f abcd.a abcd.s abci abij aibj ijkl aijk 56'
|
|
call ishell(end)
|
|
endif
|
|
101 continue
|
|
C Final results
|
|
write(iout,*)
|
|
write(iout,*)
|
|
write(iout,122)
|
|
write(iout,122)
|
|
write(iout,*)
|
|
write(iout,*) 'Final local CC results:'
|
|
write(iout,*)
|
|
|
|
if (nact.lt.1) nact=1
|
|
write(iout,"(' Minimum number of correlated electrons: ',
|
|
$i6,7x,'(',f5.1,'%)')") nocorrmin,
|
|
$100.d0*dble(nocorrmin)/dble(noccolda+noccoldb)
|
|
write(iout,"(' Maximum number of correlated electrons: ',
|
|
$i6,7x,'(',f5.1,'%)')") nocorrmax,
|
|
$100.d0*dble(nocorrmax)/dble(noccolda+noccoldb)
|
|
write(iout,"(' Average number of correlated electrons: ',
|
|
$f8.1,5x,'(',f5.1,'%)')") dble(nocorr)/dble(nact),
|
|
$100.d0*dble(nocorr)/dble(noccolda+noccoldb)/dble(nact)
|
|
write(iout,"(' Minimum number of correlated orbitals: ',
|
|
$i6,7x,'(',f5.1,'%)')") nbcorrmin,
|
|
$100.d0*dble(nbcorrmin)/dble(nbasisold)
|
|
write(iout,"(' Maximum number of correlated orbitals: ',
|
|
$i6,7x,'(',f5.1,'%)')") nbcorrmax,
|
|
$100.d0*dble(nbcorrmax)/dble(nbasisold)
|
|
write(iout,"(' Average number of correlated orbitals: ',
|
|
$f8.1,5x,'(',f5.1,'%)')") dble(nbcorr)/dble(nact),
|
|
$100.d0*dble(nbcorr)/dble(nbasisold)/dble(nact)
|
|
c
|
|
do imet=1,nmet
|
|
if(trim(meth(imet)).eq.'MP2') then
|
|
mp2cpu=totcpu(imet)
|
|
mp2wal=totwal(imet)
|
|
imp2=imet
|
|
elseif(trim(meth(imet)).eq.'CCSD') then
|
|
totcpu(imet)=totcpu(imet)+mp2cpu ! MP2 part of the LIS CCSD calculation is added to CCSD
|
|
totwal(imet)=totwal(imet)+mp2wal
|
|
cccpu=totcpu(imet)
|
|
ccwal=totwal(imet)
|
|
iccsd=imet
|
|
if (lccoporder.eq.'lccfirst') then
|
|
ovirtcpu=ovirtcpu-cccpu ! time of CCSD was part of time of drpa (ovirt...)
|
|
ovirtwal=ovirtwal-ccwal
|
|
endif
|
|
elseif(trim(meth(imet)).eq.'CCSD(T)') then
|
|
tcpu=totcpu(imet) ! time of (T) was part of time of drpa (ovirt...)
|
|
twal=totwal(imet)
|
|
ipart=imet
|
|
if (lccoporder.eq.'lccfirst') then
|
|
ovirtcpu=ovirtcpu-tcpu
|
|
ovirtwal=ovirtwal-twal
|
|
endif
|
|
endif
|
|
enddo
|
|
c
|
|
write(iout,*)
|
|
write(iout,"(' CPU time for integral calculations [min]: ',
|
|
$3x,f12.3)") intcpu/60.d0
|
|
write(iout,"(' Wall time for integral calculations [min]: ',
|
|
$3x,f12.3)") intwal/60.d0
|
|
write(iout,*)
|
|
write(iout,"(' CPU time for SCF calculation [min]: ',
|
|
$3x,f12.3)") scfcpu/60.d0
|
|
write(iout,"(' Wall time for SCF calculation [min]: ',
|
|
$3x,f12.3)") scfwal/60.d0
|
|
write(iout,*)
|
|
write(iout,"(' CPU time for orbital localization [min]: ',
|
|
$3x,f12.3)") orbcpu/60.d0
|
|
write(iout,"(' Wall time for orbital localization [min]: ',
|
|
$3x,f12.3)") orbwal/60.d0
|
|
if(localcc.eq.'2013') then
|
|
write(iout,*)
|
|
write(iout,"(' CPU time for domain construction [min]: ',
|
|
$3x,f12.3)") mulcpu/60.d0
|
|
write(iout,"(' Wall time for domain construction [min]: ',
|
|
$3x,f12.3)") mulwal/60.d0
|
|
endif
|
|
write(iout,*)
|
|
write(iout,"(' CPU time for integral transformations [min]: ',
|
|
$3x,f12.3)") ovirtcpu/60.d0
|
|
write(iout,"(' Wall time for integral transformations [min]: ',
|
|
$3x,f12.3)") ovirtwal/60.d0
|
|
write(iout,*)
|
|
beg=' Reference energy [au]:'
|
|
write(iout,7597) beg,escf
|
|
if(lmp2dens.eq.'on ') then
|
|
write(iout,*)
|
|
beg=' L' // trim(ccorr) // ' correlation energy [au]:'
|
|
write(iout,7597) beg,lmp2ene
|
|
beg=' Total L' // trim(ccorr) // ' energy [au]:'
|
|
write(iout,7597) beg,lmp2ene+escf
|
|
write(iout,*)
|
|
if(localcc.eq.'2013') then
|
|
do imet=1,nmet
|
|
if(trim(meth(imet)).eq.trim(ccorr))mp2corr=lmp2ene-ecl(imet)
|
|
enddo
|
|
endif
|
|
if (naftyp.eq.'jpqpilot'.or.naftyp.eq.'jpq ') then
|
|
beg=' Total MP2 correction for dropped NAFs [au]:'
|
|
write(iout,7597) beg,mp2corr-pairpllno
|
|
endif
|
|
beg=' Total ' // trim(ccorr) // ' correction [au]:'
|
|
write(iout,7597) beg,mp2corr
|
|
if (corembed.ne.'off ')
|
|
$ call setenvcalc(envcalc,cscr16,minpfile,iout,lrc)
|
|
if (trim(envcalc).eq.'mp2'.or.lanycc(envcalc)) then ! mp2env.or.ccenv
|
|
call uppercase(envcalc,uenvcalc,16)
|
|
beg=' Local ' // trim(uenvcalc) //
|
|
$ ' contribution of the environment [au]:'
|
|
write(iout,7598) beg,eenv
|
|
beg=' Total correction: ' // trim(ccorr) //
|
|
% ' + environment [au]:'
|
|
write(iout,7597) beg,eenv+mp2corr
|
|
totcorr='all '
|
|
else
|
|
totcorr=ccorr
|
|
endif
|
|
endif
|
|
write(iout,*)
|
|
C
|
|
open(unit=ifcfile,file='iface',status='replace')
|
|
write(ifcfile,"(68a)")
|
|
$'#property method sym st mul value' //
|
|
$' CPU(sec) Wall(sec)'
|
|
do imet=1,nmet
|
|
if(lmp2dens.eq.'off'.or.(trim(meth(imet)).ne.'MP2'.and.
|
|
$trim(meth(imet)).ne.trim(ccorr))) then
|
|
if(trim(meth(imet)).eq.'CCSD(T)') then
|
|
beg=' CPU time for (T) corrections [min]:'
|
|
else
|
|
beg=' CPU time for ' // meth(imet)
|
|
end=' calculations [min]:'
|
|
call trimloc(beg,end)
|
|
endif
|
|
write(iout,"(a54,f11.3)") beg,totcpu(imet)/60.d0
|
|
totcpu(imet)=totcpu(imet)+intcpu+scfcpu+orbcpu+mulcpu+ovirtcpu
|
|
if(trim(meth(imet)).eq.'CCSD(T)')
|
|
$totcpu(imet)=totcpu(imet)+cccpu
|
|
beg=' Total CPU time for ' // meth(imet)
|
|
end=' [min]:'
|
|
call trimloc(beg,end)
|
|
write(iout,"(a54,f11.3)") beg,totcpu(imet)/60.d0
|
|
if(trim(meth(imet)).eq.'CCSD(T)') then
|
|
beg=' Wall time for (T) corrections [min]:'
|
|
else
|
|
beg=' Wall time for ' // meth(imet)
|
|
end=' calculations [min]:'
|
|
call trimloc(beg,end)
|
|
endif
|
|
write(iout,"(a54,f11.3)") beg,totwal(imet)/60.d0
|
|
totwal(imet)=totwal(imet)+intwal+scfwal+orbwal+mulwal+ovirtwal
|
|
if(trim(meth(imet)).eq.'CCSD(T)')
|
|
$totwal(imet)=totwal(imet)+ccwal
|
|
beg=' Total wall time for ' // meth(imet)
|
|
end=' [min]:'
|
|
call trimloc(beg,end)
|
|
write(iout,"(a54,f11.3)") beg,totwal(imet)/60.d0
|
|
write(iout,*)
|
|
if (trim(meth(imet)).eq.'CCSD(T)') then
|
|
beg=' Total (T) correlation energy contribution [au]:'
|
|
write(iout,7597) beg,ecl(imet)-ecl(iccsd)
|
|
if (localccXp(localcc,2018).and.corembed.eq.'off ') then
|
|
beg=' Total (T) correlation energy+0.5 MP2 correction [au]:'
|
|
write(iout,7597) beg,ecl(imet)-ecl(iccsd)+mp2corr/2.d0
|
|
endif
|
|
elseif (trim(meth(imet)).eq.'CCSD'.and.localccXp(localcc,2018)
|
|
$ .and.trim(calctype).eq.'ccsd(t)') then
|
|
mp2fact=0.5d0
|
|
mp2fact_halved=.true.
|
|
if (totcorr.eq.ccorr) half='0.5'
|
|
endif
|
|
beg=' ' // meth(imet)
|
|
end=' correlation energy [au]:'
|
|
call trimloc(beg,end)
|
|
write(iout,7597) beg,ecl(imet)
|
|
beg=' Total ' // meth(imet)
|
|
end=' energy [au]:'
|
|
call trimloc(beg,end)
|
|
write(iout,7597) beg,ecl(imet)+escf
|
|
if(lmp2dens.eq.'on ') then
|
|
beg=' ' // meth(imet)
|
|
end=' correlation energy + '//
|
|
$trim(adjustl(half // ' ' // totcorr))// ' corrections [au]:'
|
|
call trimloc(beg,end)
|
|
write(iout,7597) beg,ecl(imet)+mp2fact*mp2corr+eenv
|
|
beg=' Total LNO-' // meth(imet)
|
|
end=' energy with '//trim(totcorr)//' corrections [au]:'
|
|
call trimloc(beg,end)
|
|
write(iout,7597) beg,ecl(imet)+escf+mp2fact*mp2corr+eenv
|
|
endif
|
|
write(iout,*)
|
|
flush(iout)
|
|
endif
|
|
C Write iface file
|
|
write(ifcfile,7596) label,meth(imet),n1,n2,n3,ecl(imet)+escf,
|
|
$totcpu(imet),totwal(imet)
|
|
if(lmp2dens.eq.'on ')
|
|
$write(ifcfile,7596) label,
|
|
$adjustl(trim(meth(imet)) // '+' // trim(ccorr)),n1,n2,n3,
|
|
$ecl(imet)+escf+mp2fact*mp2corr+eenv,totcpu(imet),totwal(imet)
|
|
if (mp2fact_halved) then
|
|
mp2fact=1.0d0
|
|
mp2fact_halved=.false.
|
|
if (totcorr.eq.ccorr) half=' '
|
|
endif
|
|
enddo
|
|
write(iout,122)
|
|
write(iout,122)
|
|
close(ifcfile)
|
|
c
|
|
122 format(1x,70('='))
|
|
7596 format(a8,1x,a15,1x,3(i2,1x),1pe23.15,2(1pe15.8))
|
|
7597 format(a54,f20.12,a21)
|
|
7598 format(a57,f17.12)
|
|
C
|
|
return
|
|
end
|
|
C
|
|
************************************************************************
|
|
subroutine runit(string,cput,walt)
|
|
************************************************************************
|
|
* Call external programs *
|
|
************************************************************************
|
|
#include "MRCCCOMMON"
|
|
c variables taken from mrcccommon: iout, outfilename, exitfile
|
|
c variables taken from mpif.h: MPI_COMM_WORLD
|
|
c implicit none
|
|
|
|
#if defined(MPI)
|
|
character(len = 4) :: localcc
|
|
character(len = 8) :: scfalg
|
|
integer :: mpi_rank, mpi_err
|
|
#endif
|
|
integer ierr
|
|
real*8 cput,walt
|
|
character*(*) string
|
|
character*16 devnul,calctype
|
|
logical lll,localccXp,lf12,dof12
|
|
common/output/ devnul
|
|
C
|
|
flush(iout)
|
|
|
|
lf12=.false. ! TODO resest mpi in scf in the case of F12
|
|
if(trim(string).ne.'minp') then
|
|
call getkey('calc',4,calctype,16)
|
|
lf12=dof12(calctype)
|
|
end if
|
|
|
|
#if defined(MPI)
|
|
call MPI_Comm_rank(MPI_COMM_WORLD, mpi_rank, mpi_err)
|
|
|
|
if(mpi_rank .ne. 0) close(iout)
|
|
#endif
|
|
|
|
call ishell('echo'//devnul)
|
|
call ishell(
|
|
$'echo " ************************ "`date +"%F %T"`"' //
|
|
$' *************************"'//trim(devnul))
|
|
call ishell('echo " Executing ' // string//'..."'//trim(devnul))
|
|
call ishell('echo'//devnul)
|
|
call flush(iout)
|
|
#if defined (G95)
|
|
if(string.eq.'goldstone') then
|
|
call ishell('export G95_IGNORE_ENDFILE=t; exec goldstone'
|
|
$// trim(devnul))
|
|
else
|
|
call ishell('export G95_IGNORE_ENDFILE=f; '// "exec "// string
|
|
$// trim(devnul))
|
|
endif
|
|
#else
|
|
|
|
#if defined(MPI)
|
|
c call MPI_Comm_size(MPI_COMM_WORLD, mpi_size, mpi_err)
|
|
call MPI_Barrier(MPI_COMM_WORLD, mpi_err)
|
|
|
|
localcc='off '
|
|
if (string.eq.'drpa') call getkey('localcc',7,localcc,4) ! cannot call getkey for string=minp
|
|
|
|
scfalg=' '
|
|
if(string.ne.'minp') then
|
|
inquire(file='KEYWD',exist=lll)
|
|
if(lll) call getkey('scfalg',6,scfalg,8)
|
|
endif
|
|
|
|
if(string .eq. 'ccsd' .or.
|
|
$ string .eq. 'mrcc' .or.
|
|
$ (string .eq. 'scf'.and.scfalg.ne.'disk '.and..not.lf12).or.
|
|
$ (string .eq. 'drpa' .and. localccxp(localcc,2022))) then
|
|
call run_mpi_program(string)
|
|
else if(string.eq.'xmrcc') then
|
|
call ishell( "exec xmrcc_mpi" // trim(devnul))
|
|
else
|
|
#endif
|
|
call ishell( "exec " // string // trim(devnul))
|
|
#if defined(MPI)
|
|
end if
|
|
|
|
call MPI_Barrier(MPI_COMM_WORLD, mpi_err)
|
|
|
|
if(mpi_rank .ne. 0)
|
|
$ open(unit = iout, file = outfilename, access = 'append')
|
|
#endif
|
|
|
|
#endif
|
|
open(unit=exitfile,file='EXIT',status='unknown')
|
|
read(exitfile,*) ierr,cput,walt
|
|
close(exitfile)
|
|
if(ierr.ne.0) then
|
|
write(iout,*)
|
|
write(iout,*) 'Fatal error in ',string,'.'
|
|
write(iout,*) 'Program will stop.'
|
|
call flush(iout)
|
|
call dmrccend(1)
|
|
endif
|
|
C
|
|
return
|
|
end
|
|
C
|
|
#if defined(MPI)
|
|
************************************************************************
|
|
subroutine run_mpi_program(string)
|
|
************************************************************************
|
|
#include "MRCCCOMMON"
|
|
c implicit none
|
|
c variables taken from mrcccommon:
|
|
c idfile,wdir,output_folder
|
|
c variables taken from mpif.h:
|
|
c MPI_COMM_WORLD,MPI_ARGV_NULL,MPI_INFO_NULL,MPI_ERRCODES_IGNORE,
|
|
c MPI_CHARACTER,MPI_ROOT,MPI_PROC_NULL,MPI_STATUS_IGNORE,MPI_MAX_PROCESSOR_NAME
|
|
character(len = *) :: string
|
|
|
|
logical :: flag
|
|
integer :: request, time
|
|
integer :: mpitasks, mpi_rank, mpi_err
|
|
c integer, allocatable :: error_codes(:)
|
|
c integer :: error_codes(1024)
|
|
character(len = 4) :: mpitasks_ch
|
|
|
|
character(len = 8) :: date
|
|
character(len = 10) :: times
|
|
character(len = 5) :: zone
|
|
integer :: dt(8)
|
|
character(len = MPI_MAX_PROCESSOR_NAME) :: pname
|
|
integer :: nlen
|
|
character(len = 64) :: id
|
|
integer :: children_comm
|
|
|
|
call MPI_Comm_rank(MPI_COMM_WORLD, mpi_rank, mpi_err)
|
|
|
|
call getkey('mpitasks', 8, mpitasks_ch, 4)
|
|
read(mpitasks_ch, '(i4)') mpitasks
|
|
|
|
call MPI_Get_processor_name(pname, nlen, mpi_err)
|
|
call date_and_time(date, times, zone, dt)
|
|
|
|
! create file for identification
|
|
open(idfile, file = "dmrcc.id")
|
|
id = date // times // zone // trim(pname)
|
|
write(idfile, *) trim (id); flush(idfile)
|
|
close(idfile)
|
|
|
|
! spawn mpi processes
|
|
c allocate(error_codes(mpi_size))
|
|
c allocate(error_codes(mpitasks))
|
|
call MPI_Comm_spawn(trim(string) // '_mpi', MPI_ARGV_NULL,
|
|
$ mpitasks,MPI_INFO_NULL, 0, MPI_COMM_WORLD, children_comm,
|
|
c $ error_codes, mpi_err)
|
|
$ MPI_ERRCODES_IGNORE, mpi_err)
|
|
c deallocate(error_codes)
|
|
|
|
! communicate working directory
|
|
if(mpi_rank .eq. 0) then
|
|
call MPI_Bcast(wdir, 512, MPI_CHARACTER, MPI_ROOT,
|
|
$ children_comm, mpi_err)
|
|
else
|
|
call MPI_Bcast(wdir, 512, MPI_CHARACTER, MPI_PROC_NULL,
|
|
$ children_comm, mpi_err)
|
|
end if
|
|
|
|
! communicate output directory
|
|
if(mpi_rank .eq. 0) then
|
|
call MPI_Bcast(output_folder, 32, MPI_CHARACTER, MPI_ROOT,
|
|
$ children_comm, mpi_err)
|
|
else
|
|
call MPI_Bcast(output_folder, 32, MPI_CHARACTER, MPI_PROC_NULL,
|
|
$ children_comm, mpi_err)
|
|
end if
|
|
|
|
! communicate id
|
|
if(mpi_rank .eq. 0) then
|
|
call MPI_Bcast(pname, MPI_MAX_PROCESSOR_NAME, MPI_CHARACTER,
|
|
$ MPI_ROOT, children_comm, mpi_err)
|
|
else
|
|
call MPI_Bcast(pname, MPI_MAX_PROCESSOR_NAME, MPI_CHARACTER,
|
|
$ MPI_PROC_NULL, children_comm, mpi_err)
|
|
end if
|
|
|
|
! wait for children to finish
|
|
call MPI_Ibarrier(children_comm, request, mpi_err)
|
|
flag = .false.
|
|
do while(.not. flag)
|
|
call sleep_f(1) ! required for OpenMPI to yield processor
|
|
call MPI_Test(request, flag, MPI_STATUS_IGNORE, mpi_err)
|
|
end do
|
|
|
|
call MPI_Comm_disconnect(children_comm, mpi_err)
|
|
c call sleep_f(2)
|
|
|
|
call delete_file('dmrcc.id', .true.)
|
|
|
|
end subroutine
|
|
#endif
|
|
C
|
|
************************************************************************
|
|
subroutine xmrccm(icyc,nactoo,nactvo,convero,almem,almemo,remem,
|
|
$lopt,llow,tomem,tomemo)
|
|
************************************************************************
|
|
* Driver for xmrcc *
|
|
************************************************************************
|
|
#include "MRCCCOMMON"
|
|
c variables taken from MRCCCOMMON:
|
|
c scrfile6, nacto,nactv,nocc,nvirt,conver,pert
|
|
integer icyc,nactoo,nactvo,convero
|
|
real*8 almem,almemo,remem,tomem,tomemo,cput,walt
|
|
logical lopt,llow
|
|
C
|
|
call runit('goldstone',cput,walt)
|
|
C
|
|
open(scrfile6,status='unknown',file='fort.24',form='unformatted')
|
|
rewind(scrfile6)
|
|
write(scrfile6) icyc,nactoo,nactvo,convero,almem,almemo,remem,
|
|
$lopt,llow,nacto,nactv,nocc,nvirt,conver,pert,tomem,tomemo
|
|
close(scrfile6)
|
|
C
|
|
call runit('xmrcc',cput,walt)
|
|
C
|
|
open(scrfile6,status='unknown',file='fort.24',form='unformatted')
|
|
rewind(scrfile6)
|
|
read(scrfile6) icyc,nactoo,nactvo,convero,almem,almemo,remem,
|
|
$lopt,llow,nacto,nactv,nocc,nvirt,conver,pert,tomem,tomemo
|
|
close(scrfile6)
|
|
C
|
|
return
|
|
end
|
|
C
|
|
************************************************************************
|
|
subroutine updatei(log,noccup,nbmax,nmax,lr)
|
|
************************************************************************
|
|
* Update input file *
|
|
************************************************************************
|
|
#include "MRCCCOMMON"
|
|
c implicit none
|
|
c variables taken from MRCCCOMMON:
|
|
c inpfile,ptfreq,op1,nsing,ntrip,rest,op,conver,iclsh,ihf,ndoub,
|
|
c nacto,nactv,ssmrcc,nbasis
|
|
|
|
integer i,j,jj,nline,iscr1(18),iscr2(2),noccup(*),nbmax,nmax
|
|
integer occv(nbmax),actv(nbmax),exv(nmax)
|
|
character line(1000)
|
|
real*8 eps
|
|
logical log,lr
|
|
C
|
|
open(inpfile,status='unknown',file='fort.56')
|
|
rewind(inpfile)
|
|
read(inpfile,*) iscr1,ptfreq,iscr2!,eps
|
|
nline=1
|
|
iscr1(1)=op1
|
|
iscr1(2)=nsing
|
|
iscr1(3)=ntrip
|
|
iscr1(4)=rest
|
|
if(lr) then
|
|
if(op1.ne.op) then
|
|
iscr1(1)=op1+1
|
|
iscr1(5)=3
|
|
else
|
|
iscr1(5)=1
|
|
endif
|
|
endif
|
|
iscr1(7)=conver
|
|
iscr1(10)=iclsh
|
|
iscr1(12)=ihf
|
|
iscr1(13)=ndoub
|
|
iscr1(14)=nacto
|
|
iscr1(15)=nactv
|
|
if(nactv*nacto.ne.0.and.ssmrcc.eq.0) iscr1(17)=iscr1(1)
|
|
read(inpfile,'(1000a1)',end=1001) line
|
|
nline=2
|
|
read(inpfile,*,end=1001) (occv(i),i=1,nbmax)
|
|
nline=3
|
|
read(inpfile,*,end=1001) (actv(i),i=1,nbmax)
|
|
nline=4
|
|
read(inpfile,*,end=1001) (exv(i),i=1,op)
|
|
nline=5
|
|
1001 continue
|
|
C
|
|
rewind(inpfile)
|
|
write(inpfile,"(18i8,f9.5,i7,i9,i7,1pe10.3)") iscr1,ptfreq,iscr2!,eps
|
|
jj=1
|
|
do j=1,1000
|
|
if(line(j).ne.' ') jj=j
|
|
enddo
|
|
write(inpfile,"(1000a1)") (line(j),j=1,jj)
|
|
if(log) then
|
|
write(inpfile,"(10000i3)") (noccup(j),j=1,nbmax)
|
|
else
|
|
if(nline.ge.3)write(inpfile,"(10000i3)") (occv(i),i=1,nbmax)
|
|
endif
|
|
if(nline.ge.4) write(inpfile,"(10000i3)") (actv(i),i=1,nbmax)
|
|
if(nline.ge.5) write(inpfile,"(10000i3)") (exv(i),i=1,op)
|
|
close(inpfile)
|
|
C
|
|
return
|
|
end
|
|
************************************************************************
|
|
subroutine memopt(icyc,nactoo,nactvo,convero,almem,almemo,remem,
|
|
$lopt,llow,llv,nact1,noc1,nacto1,tomem,tomemo,nbmax)
|
|
************************************************************************
|
|
* Optimize memory requirement *
|
|
************************************************************************
|
|
#include "MRCCCOMMON"
|
|
c variables taken from MRCCCOMMON: nbasis,op
|
|
integer i,icyc,nactoo,nactvo,convero,nn,nact1,noc1,nacto1,nbmax
|
|
real*8 almem,almemo,remem,tomem,tomemo
|
|
logical lopt,llow,log,llv
|
|
C
|
|
llv=.true.
|
|
if(almem.gt.remem.and.lopt.and.nact1.eq.0) then
|
|
nact1=nint(dble(noc1)/2.d0)
|
|
if(mod(nact1,2).eq.1) then
|
|
nact1=nact1+1
|
|
if(nact1.ge.noc1) nact1=max(0,nact1-2)
|
|
endif
|
|
if(nact1.gt.0) then
|
|
call updatei(.false.,i,nbmax,op,.false.)
|
|
call xmrccm(icyc,nactoo,nactvo,convero,almem,almemo,remem,
|
|
$lopt,llow,tomem,tomemo)
|
|
endif
|
|
endif
|
|
C
|
|
if(almem.gt.remem.and.lopt.and.nact1.ne.0) then
|
|
log=.true.
|
|
llow=.true.
|
|
nn=nact1
|
|
do while(llow.and.nact1.ge.4)
|
|
nact1=nact1-2
|
|
call updatei(.false.,i,nbmax,op,.false.)
|
|
call xmrccm(icyc,nactoo,nactvo,convero,almem,almemo,remem,
|
|
$lopt,llow,tomem,tomemo)
|
|
if(llow) log=.false.
|
|
if(almem.le.remem) then
|
|
llv=.false.
|
|
return
|
|
endif
|
|
enddo
|
|
if(log) then
|
|
llow=.true.
|
|
nact1=nn
|
|
do while(llow.and.nact1.le.noc1-4)
|
|
nact1=nact1+2
|
|
call updatei(.false.,i,nbmax,op,.false.)
|
|
call xmrccm(icyc,nactoo,nactvo,convero,almem,almemo,remem,
|
|
$lopt,llow,tomem,tomemo)
|
|
if(almem.le.remem) then
|
|
llv=.false.
|
|
return
|
|
endif
|
|
enddo
|
|
endif
|
|
nact1=nacto1
|
|
endif
|
|
C
|
|
return
|
|
end
|
|
C
|
|
************************************************************************
|
|
subroutine ddmrcc(ccprog,domain_is_cs,llg)
|
|
************************************************************************
|
|
* Driver for one mrcc calling
|
|
************************************************************************
|
|
#include "MRCCCOMMON"
|
|
c implicit none
|
|
c variables taken from MRCCCOMMON:
|
|
c ifcfile,inp,gfile,iout
|
|
c rest,ssmrcc,locno,maxmem,nvirt,nocc,pert,conver,nbasis,op,op1,nactv
|
|
c nacto,calc,nroot
|
|
integer icyc,nactoo,nactvo,convero,mm,nn,i,j,k,l,m,ilow,restold
|
|
integer nbmax
|
|
real*8 remem,almem,almemo,tomem,tomemo,cput,walt
|
|
logical lopt,llow,log,llv,log1,domain_is_cs,llg
|
|
character*4 ccprog,ovirt,localcc
|
|
character*5 scftype,dfintran
|
|
character*16 naf, eps
|
|
C
|
|
nbmax=nbasis/2
|
|
open(inp,status='unknown',file='fort.55')
|
|
rewind(inp)
|
|
read(inp,*,end=1006) nbmax
|
|
1006 close(inp)
|
|
C
|
|
if(llg) then
|
|
call getkey('naf_cor',7,naf,16)
|
|
call getkey('ovirt',5,ovirt,4)
|
|
call getkey('localcc',7,localcc,4)
|
|
call getkey('dfintran',8,dfintran,5)
|
|
call getkey('scftype',7,scftype,5)
|
|
else
|
|
naf=' '
|
|
ovirt=' '
|
|
localcc=' '
|
|
dfintran=' '
|
|
scftype=' '
|
|
endif
|
|
c call getkey('eps',3,eps,16)
|
|
c if((naf.eq.'off ' .and.
|
|
cc $ eps.eq.' ').or.ccprog.ne.'ccsd') then
|
|
c $ ovirt.ne.'mp2 ').or.ccprog.ne.'ccsd') then
|
|
if((.not.llg.or..not.(ccprog.eq.'ccsd' .and. localcc.eq.'off '
|
|
$.and.((trim(naf).ne.'off' .or. ovirt.eq.'mp2 ' .or.
|
|
$ ovirt.eq.'osv' .or. ovirt.eq.'ppl ') .and.
|
|
$dfintran.eq.'drpa '))).and.scftype.ne.'mcscf') then
|
|
C Initialize interface file
|
|
open(unit=ifcfile,file='iface',status='replace')
|
|
write(ifcfile,"(68a)")
|
|
$ '#property method sym st mul value' //
|
|
$ ' CPU(sec) Wall(sec)'
|
|
close(ifcfile)
|
|
end if
|
|
C Run Lori's ccsd
|
|
m=0
|
|
c ccprog=' '
|
|
if(llg) then
|
|
c call getkey('ccprog',6,ccprog,4)
|
|
if(ccprog.eq.'ccsd') then
|
|
if(scftype.eq.'uhf '.or.scftype.eq.'rohf ') then
|
|
if (domain_is_cs) then
|
|
call runit('ccsd',cput,walt)
|
|
else
|
|
call runit('uccsd',cput,walt)
|
|
endif
|
|
else
|
|
call runit('ccsd',cput,walt)
|
|
endif
|
|
return
|
|
else if(ccprog.eq.'cis ') then
|
|
call runit('cis',cput,walt)
|
|
return
|
|
endif
|
|
else
|
|
open(inp,status='unknown',file='fort.55')
|
|
rewind(inp)
|
|
read(inp,*,end=1005) m
|
|
1005 close(inp)
|
|
endif
|
|
C
|
|
lopt=.false.
|
|
icyc=0
|
|
almemo=2.d0**32.d0
|
|
tomemo=2.d0**32.d0
|
|
nactoo=0
|
|
nactvo=0
|
|
convero=0
|
|
almem=0.d0
|
|
remem=0.d0
|
|
tomem=0.d0
|
|
llow=.false.
|
|
restold=rest
|
|
C Check if input files have been changed
|
|
if(ssmrcc.gt.0.or.locno.ge.3) goto 1000
|
|
inquire(file='fort.14',exist=log1)
|
|
if(.not.log1) goto 1000
|
|
inquire(file='fort.10',exist=log1)
|
|
if(.not.log1) goto 1000
|
|
rewind(gfile)
|
|
read(gfile,*,end=1000) i,j,k,maxmem,l
|
|
close(gfile)
|
|
call ishell('cat fort.56 > fort.11')
|
|
call ishell("sed '/E/,$d' fort.55 >> fort.11")
|
|
call ishell(
|
|
$"sed '1,/inpfile/d' fort.14 | diff fort.11 - | wc -l > fort.99")
|
|
call ishell(
|
|
$"sed '1,/inpfile/d' fort.10 | diff fort.11 - | wc -l >> fort.99")
|
|
rewind(99)
|
|
read(99,*,end=1000) nn
|
|
read(99,*,end=1000) mm
|
|
close(99,status='delete')
|
|
if(m.gt.0.and.mm+nn.eq.0.and.j.eq.nvirt.and.k.eq.nocc) goto 1004
|
|
1000 continue
|
|
C Calculate memory requirement using the input data
|
|
call xmrccm(icyc,nactoo,nactvo,convero,almem,almemo,remem,lopt,
|
|
$llow,tomem,tomemo)
|
|
C Approximate methods
|
|
if(almem.gt.remem.and.pert.gt.0.and.conver.ne.2.and.locno.ne.3)
|
|
$then
|
|
conver=2
|
|
call updatei(.false.,i,nbmax,op,.false.)
|
|
call xmrccm(icyc,nactoo,nactvo,convero,almem,almemo,remem,lopt,
|
|
$llow,tomem,tomemo)
|
|
if(.not.llow) conver=0
|
|
endif
|
|
C Optimizing the number of active particles
|
|
call memopt(icyc,nactoo,nactvo,convero,almem,almemo,remem,lopt,
|
|
$llow,llv,nactv,nvirt,nactvo,tomem,tomemo,nbmax)
|
|
C Optimizing the number of active holes
|
|
if(llv.and.locno.ne.3) call memopt(icyc,nactoo,nactvo,convero,
|
|
$almem,almemo,remem,lopt,llow,llv,nacto,nocc,nactoo,tomem,tomemo,
|
|
$nbmax)
|
|
C Restore the inputs for the optimal case
|
|
if(lopt.and.(nactoo.ne.nacto.or.nactv.ne.nactvo.or.
|
|
$convero.ne.conver.or..not.llow)) then
|
|
nacto=nactoo
|
|
nactv=nactvo
|
|
conver=convero
|
|
call updatei(.false.,i,nbmax,op,.false.)
|
|
call xmrccm(icyc,nactoo,nactvo,convero,almem,almemo,remem,lopt,
|
|
$llow,tomem,tomemo)
|
|
endif
|
|
C
|
|
if(almem.gt.remem) then
|
|
write(iout,*)
|
|
$'Unable to optimize memory usage: Insufficient memory!'
|
|
call flush(iout)
|
|
call dmrccend(1)
|
|
else if(lopt.and.abs(rest).eq.1.and.nacto+nactv.gt.0) then
|
|
write(iout,*) 'The number of active orbitals has been changed!'
|
|
write(iout,*) 'It is dangerous to restart the program!'
|
|
call flush(iout)
|
|
endif
|
|
C Save input data
|
|
call ishell('echo "inpfile" >> fort.14')
|
|
call ishell('cat fort.56 >> fort.14')
|
|
call ishell("sed '/E/,$d' fort.55 >> fort.14")
|
|
call ishell('echo "inpfile" >> fort.10')
|
|
call ishell('cat fort.56 >> fort.10')
|
|
call ishell("sed '/E/,$d' fort.55 >> fort.10")
|
|
C Execute lower-level calculations if rest=2 or rest=4
|
|
if(abs(rest).eq.2.or.abs(rest).eq.4) then
|
|
ilow=2
|
|
if(calc.ge.2.and.calc.le.9) ilow=3
|
|
if(nroot.gt.1.and.abs(rest).ne.4) ilow=1
|
|
do op1=ilow,op
|
|
if(op1.eq.ilow) then
|
|
if(abs(rest).eq.4) then
|
|
rest=isign(3,rest)
|
|
else
|
|
rest=0
|
|
endif
|
|
else
|
|
rest=isign(2,rest)
|
|
endif
|
|
call updatei(.false.,i,nbmax,op,calc.eq.1.and.nroot.eq.1)
|
|
call runit('goldstone',cput,walt)
|
|
call runit('xmrcc',cput,walt)
|
|
if(op1.ne.op) call runit('mrcc',cput,walt)
|
|
enddo
|
|
endif
|
|
C Execute mrcc
|
|
1004 continue
|
|
C #if defined (MPI) || (defined (PGF) && defined (OMP))
|
|
#if defined (MPI)
|
|
if(ssmrcc.gt.0) then
|
|
call runit('mrcc',cput,walt) !!!szemet!!!
|
|
else
|
|
c write(iout,*)
|
|
c write(iout,*) 'Now launch mrcc in parallel mode!'
|
|
c write(iout,*)
|
|
call run_mpi_program('mrcc')
|
|
endif
|
|
#else
|
|
call runit('mrcc',cput,walt)
|
|
#endif
|
|
C Restore input file
|
|
if(rest.ne.restold) then
|
|
rest=restold
|
|
op1=op
|
|
call updatei(.false.,i,nbmax,op,.false.)
|
|
endif
|
|
C
|
|
return
|
|
end
|
|
C
|
|
c
|
|
c************************************************************************
|
|
c logical function lany_55i_exist(nocc)
|
|
c************************************************************************
|
|
c* Check if any of the 55.i file exists required for local CC runs
|
|
c************************************************************************
|
|
c implicit none
|
|
c integer nocc,i
|
|
c character(len=8) ilmoname
|
|
c character(len=11) fort55name
|
|
c logical fexist
|
|
c
|
|
c lany_55i_exist=.false.
|
|
c do i=1,nocc
|
|
c write(ilmoname,'(i6)') i
|
|
c fort55name=trim('55') // '.' // trim(ilmoname)
|
|
c inquire(file=trim(fort55name),exist=fexist)
|
|
c lany_55i_exist=lany_55i_exist.or.fexist
|
|
c enddo
|
|
cc
|
|
c end
|
|
cC
|
|
************************************************************************
|
|
logical function locdrvdone(mode)
|
|
************************************************************************
|
|
* communicate whether locdrv() call is completed in ldrpa
|
|
* mode=0 delete locdrvdone
|
|
* mode=1 set locdrvdone as completed
|
|
* mode=2 check whether locdrv() was called in ldrpa
|
|
************************************************************************
|
|
implicit none
|
|
integer mode,isitdone
|
|
logical fexist
|
|
|
|
locdrvdone=.false.
|
|
inquire(file='locdrvdone',exist=fexist)
|
|
if (mode.eq.0.and.fexist) call delete_file('locdrvdone',.true.)
|
|
if (mode.eq.1) then
|
|
open(unit=301,file='locdrvdone')
|
|
if (fexist) rewind(301)
|
|
write(301,*) 1
|
|
close(301)
|
|
elseif (mode.eq.2) then
|
|
isitdone=0
|
|
if (fexist) then
|
|
open(unit=301,file='locdrvdone')
|
|
rewind(301)
|
|
read(301,*,err=102) isitdone
|
|
close(301)
|
|
102 continue
|
|
else
|
|
isitdone=0
|
|
endif
|
|
if (isitdone.eq.1) locdrvdone=.true.
|
|
if (locdrvdone) call delete_file('locdrvdone',.true.)
|
|
endif
|
|
c
|
|
end function locdrvdone
|
|
C
|
|
************************************************************************
|
|
subroutine set_vars4locdrv(mode)
|
|
************************************************************************
|
|
* set/reset variables in MRCCCOMMON before/after locdrv() call
|
|
************************************************************************
|
|
#include "MRCCCOMMON"
|
|
integer mode
|
|
c
|
|
if (mode.eq.1) then ! set
|
|
nocc=2*nocc
|
|
nbasis=2*nbasis
|
|
elseif (mode.eq.2) then ! reset
|
|
nocc=nocc/2
|
|
nbasis=nbasis/2
|
|
endif
|
|
c
|
|
return
|
|
end
|
|
c
|
|
************************************************************************
|
|
subroutine localccrestart(mode,scrfile2,calctype,methact,nmetact,
|
|
$ cc_restart,nbasis,nocc,ncore,idom,nmet,
|
|
$ nocorr,nbcorr,nocorrmin,nocorrmax,nbcorrmin,nbcorrmax,nact,
|
|
$ ecl,totcpu,
|
|
$ totwal,lmp2ene,mp2corr,pairpllno,eenv,
|
|
$ emp2,epar,ecor,edeenv,eparenv,
|
|
$ esec_im1,epar_im1,ecor_im1,edeenv_im1,ccmem,
|
|
$ intcpu,intwal,scfcpu,scfwal,orbcpu,orbwal,ovirtcpu,ovirtwal)
|
|
************************************************************************
|
|
* read, initialize or write restart file localcc.restart for ground state LNO-CC
|
|
************************************************************************
|
|
implicit none
|
|
integer scrfile2,nbasis,nocc,ncore,idom,nmet,imet,nmetact,ccmem(8)
|
|
integer nocorr,nbcorr,nocorrmin,nocorrmax,nbcorrmin,nbcorrmax,nact
|
|
real*8 ecl(*),totcpu(*),totwal(*),epar_im1,ecor_im1,edeenv_im1
|
|
real*8 lmp2ene,mp2corr,pairpllno,eenv,esec_im1
|
|
real*8 intcpu,intwal,scfcpu,scfwal,orbcpu,orbwal,ovirtcpu,ovirtwal
|
|
real*8 emp2,epar,ecor,edeenv,eparenv
|
|
character*4 mode
|
|
character*16 calctype
|
|
character*128 methact(nmetact)
|
|
logical cc_restart
|
|
c
|
|
if (mode.eq.'save') then
|
|
c
|
|
write(scrfile2,'(a16,12i8,l3,4i16,1000e28.20)') calctype,
|
|
$ nbasis,nocc,ncore,idom,nmet, ! basiss in nbasis, noccs in nocc
|
|
$ nocorr,nbcorr,nocorrmin,nocorrmax,nbcorrmin,nbcorrmax,nact,
|
|
$ cc_restart,ccmem(1),ccmem(2),ccmem(3),ccmem(4),
|
|
$ (ecl(imet),imet=1,nmet),(totcpu(imet),imet=1,nmet),
|
|
$ (totwal(imet),imet=1,nmet),lmp2ene,mp2corr,pairpllno,eenv,
|
|
$ emp2,epar,ecor,edeenv,eparenv,
|
|
$ esec_im1,epar_im1,ecor_im1,edeenv_im1,
|
|
$ intcpu,intwal,scfcpu,scfwal,orbcpu,orbwal,ovirtcpu,ovirtwal
|
|
write(scrfile2,'(30a128)') (methact(imet),imet=1,nmetact)
|
|
c
|
|
elseif (mode.eq.'read') then
|
|
c
|
|
read(scrfile2,'(a16,12i8,l3,4i16,1000e28.20)') calctype,
|
|
$ nbasis,nocc,ncore,idom,nmet, ! idomst in idom
|
|
$ nocorr,nbcorr,nocorrmin,nocorrmax,nbcorrmin,nbcorrmax,nact,
|
|
$ cc_restart,ccmem(1),ccmem(2),ccmem(3),ccmem(4),
|
|
$ (ecl(imet),imet=1,nmet),(totcpu(imet),imet=1,nmet),
|
|
$ (totwal(imet),imet=1,nmet),lmp2ene,mp2corr,pairpllno,eenv,
|
|
$ emp2,epar,ecor,edeenv,eparenv,
|
|
$ esec_im1,epar_im1,ecor_im1,edeenv_im1,
|
|
$ intcpu,intwal,scfcpu,scfwal,orbcpu,orbwal,ovirtcpu,ovirtwal
|
|
c
|
|
elseif (mode.eq.'init') then
|
|
c
|
|
write(scrfile2,'(a16,12i8,l3,4i16,1000e28.20)') calctype,
|
|
$ nbasis,nocc,ncore,0,nmet,
|
|
$ 0,0,nocc,0,nbasis,0,nocc,cc_restart,ccmem(1),ccmem(2),ccmem(3),
|
|
$ ccmem(4),(ecl(imet),imet=1,nmet),(totcpu(imet),imet=1,nmet),
|
|
$ (totwal(imet),imet=1,nmet),lmp2ene,mp2corr,pairpllno,eenv,
|
|
$ emp2,epar,ecor,edeenv,eparenv,
|
|
$ esec_im1,epar_im1,ecor_im1,edeenv_im1,
|
|
$ intcpu,intwal,scfcpu,scfwal,orbcpu,orbwal,ovirtcpu,ovirtwal
|
|
c
|
|
else
|
|
write(*,*) 'unsupported mode in localccrestart'
|
|
call mrccend(1)
|
|
endif
|
|
return
|
|
end
|
|
************************************************************************
|
|
subroutine locfilemv(ccmo,fname)
|
|
************************************************************************
|
|
* Move files in the case of local CC calculations
|
|
************************************************************************
|
|
implicit none
|
|
character*9 fname
|
|
character*128 beg,ccmo
|
|
C
|
|
beg='mv ' // trim(fname) // '.' // trim(ccmo) // ' ' //trim(fname)
|
|
call ishell(beg)
|
|
C
|
|
return
|
|
end
|
|
C
|
|
************************************************************************
|
|
subroutine locfilecp(ccmo,fname,dir)
|
|
************************************************************************
|
|
* Copy files in the case of local CC calculations
|
|
************************************************************************
|
|
implicit none
|
|
character*9 fname
|
|
character*6 ccmo,dir
|
|
character*128 beg
|
|
C
|
|
beg='cp ' // trim(adjustl(dir)) // '/' // trim(fname) // '.' //
|
|
$ trim(ccmo) // ' ' //trim(fname)
|
|
call ishell(beg)
|
|
C
|
|
return
|
|
end
|
|
C
|
|
************************************************************************
|
|
subroutine locfilecopy(mode,cmo,naf,talg,nob,domain_is_cs,
|
|
$ mpi_rank,mpi_size,iout,inbcast_group)
|
|
************************************************************************
|
|
* copy file.i of LIS i to file readable by the CC programs
|
|
* mode = 1 : *55, DFINT*, (ajb,rLNO2cLNO,S_LNO)
|
|
* mode = 3 : cp files needed for domain i in MPI execution
|
|
* mode = 4 : cleanup in MPI execution
|
|
************************************************************************
|
|
implicit none
|
|
integer mode,cmo,nob,mpi_rank,system_error,system,iout,mpi_size
|
|
character*4 talg
|
|
character*5 scftype
|
|
character*6 dir
|
|
character*16 naf
|
|
character*128 ccmo,beg,end
|
|
logical domain_is_cs,inbcast_group
|
|
c mpi
|
|
character*512 wdir,tmp_folder
|
|
character*16 rank_ch
|
|
C
|
|
call getkey('scftype',7,scftype,5)
|
|
write(ccmo,'(i6)') cmo
|
|
ccmo=adjustl(ccmo)
|
|
dir=" ."
|
|
if (mpi_size.gt.1) dir=" .."
|
|
|
|
if (mode.lt.3) then
|
|
c
|
|
call locfilecp(ccmo,'DFINT_AI ',dir)
|
|
call locfilecp(ccmo,'DFINT_IJ ',dir)
|
|
call locfilecp(ccmo,'DFINT_AB ',dir)
|
|
if (naf.ne.'off ')
|
|
$ call locfilecp(ccmo,'ajb ',dir)
|
|
if (scftype.eq.'rohf '.and..not.domain_is_cs) then
|
|
if (nob.gt.0) then
|
|
call locfilecp(ccmo,'DFINT_AIb',dir)
|
|
call locfilecp(ccmo,'DFINT_IJb',dir)
|
|
call locfilecp(ccmo,'DFINT_ABb',dir)
|
|
endif !nob.gt.0
|
|
call locfilecp(ccmo,'rLNO2cLNO',dir)
|
|
call locfilecp(ccmo,'S_LNO ',dir)
|
|
endif !scftype.eq.rohf.and..not.domain_is_cs
|
|
if (talg.eq.'lapl'.or.talg.eq.'topr'.or.talg.eq.'to '
|
|
$ .or.talg.eq.'lato') call locfilecp(ccmo,'laplbas ',dir)
|
|
C
|
|
if (mode.eq.1) then
|
|
beg='55.' // trim(ccmo)
|
|
beg=adjustl(beg)
|
|
end='cp ' // trim(beg) // ' fort.55'
|
|
call ishell(end)
|
|
end='cp ' // trim(beg) // ' 55'
|
|
call ishell(end)
|
|
endif
|
|
#if defined(MPI)
|
|
elseif (mode.eq.3) then
|
|
call getcwd(wdir)
|
|
cdel write(iout,*) 'in localfilecp3',mpi_rank,trim(wdir) ; flush(6)
|
|
c write(rank_ch, '(i16)') mpi_rank
|
|
write(rank_ch,'(i16)') cmo
|
|
tmp_folder= trim(adjustl(wdir))// "/ccsdrun." //
|
|
$ trim(adjustl(rank_ch))
|
|
cdel write(iout,*) mpi_rank,trim(tmp_folder),mode
|
|
end='mkdir -p ' // trim(tmp_folder)
|
|
call ishell(end)
|
|
call chdir(trim(tmp_folder))
|
|
end='cp -u ../FOCK ../MOCOEF ../VARS ../KEYWD '//
|
|
$ '../MINP ../55.'//trim(ccmo) // " ."
|
|
call ishell(trim(end))
|
|
c end='rm ./mrcc.out.'//trim(adjustl(rank_ch))
|
|
c call ishell(trim(end))
|
|
c end='ln -s ../mrcc.out.'//trim(adjustl(rank_ch)) // " ."
|
|
c system_error=system('pwd')
|
|
c system_error=system('ls -lrt')
|
|
c write(iout,*) 'in localfilecp3',mpi_rank,trim(tmp_folder)
|
|
elseif (mode.eq.4) then
|
|
c system_error=system('pwd')
|
|
c system_error=system('ls -lrt')
|
|
c old end='rm 55 FOCK MOCOEF VARS KEYWD fort.55 MINP' ? why not fort.55
|
|
c end='rm 55 FOCK MOCOEF VARS KEYWD MINP'
|
|
c call ishell(trim(end))
|
|
call delete_file('55',inbcast_group)
|
|
call delete_file('FOCK',inbcast_group)
|
|
call delete_file('MOCOEF',inbcast_group)
|
|
call delete_file('VARS',inbcast_group)
|
|
call delete_file('KEYWD',inbcast_group)
|
|
call delete_file('MINP',inbcast_group)
|
|
c
|
|
c write(rank_ch, '(i16)') mpi_rank
|
|
write(rank_ch,'(i16)') cmo
|
|
c end='unlink ./mrcc.out.'//trim(adjustl(rank_ch))
|
|
c call ishell(trim(end))
|
|
c call getcwd(wdir)
|
|
cdel write(iout,*) 'in localfilecp4',mpi_rank,trim(wdir) ; flush(6)
|
|
call chdir('..')
|
|
call getcwd(wdir)
|
|
tmp_folder= trim(adjustl(wdir))// "/ccsdrun." //
|
|
$ trim(adjustl(rank_ch))
|
|
cdel write(iout,*) 'in localfilecp4',mpi_rank,trim(wdir) ; flush(6)
|
|
system_error = system('rmdir ' // trim(tmp_folder))
|
|
if(system_error .ne. 0) then
|
|
write(*,'(" Error: Couldn''t delete cc temp directory on ", a)')
|
|
$ trim(tmp_folder)
|
|
end if
|
|
#endif
|
|
endif
|
|
c
|
|
return
|
|
end
|
|
C
|
|
************************************************************************
|
|
subroutine locfiledel(cmo,naf,talg)
|
|
************************************************************************
|
|
* remove file.i of LIS i
|
|
************************************************************************
|
|
implicit none
|
|
integer cmo
|
|
character*4 talg
|
|
character*5 scftype
|
|
character*16 naf
|
|
character*128 ccmo,beg,end
|
|
logical log1
|
|
c
|
|
call getkey('scftype',7,scftype,5)
|
|
write(ccmo,'(i6)') cmo
|
|
ccmo=adjustl(ccmo)
|
|
c
|
|
beg='55.' // trim(ccmo)
|
|
end='rm -f ' // trim(beg)
|
|
inquire(file=trim(beg),exist=log1)
|
|
c write(*,*) '55.i delete', beg, log1
|
|
if (log1) call ishell(end)
|
|
beg='DFINT_AI.' // trim(ccmo)
|
|
end='rm -f ' // trim(beg)
|
|
inquire(file=trim(beg),exist=log1)
|
|
if (log1) call ishell(end)
|
|
beg='DFINT_IJ.' // trim(ccmo)
|
|
end='rm -f ' // trim(beg)
|
|
inquire(file=trim(beg),exist=log1)
|
|
if (log1) call ishell(end)
|
|
beg='DFINT_AB.' // trim(ccmo)
|
|
end='rm -f ' // trim(beg)
|
|
inquire(file=trim(beg),exist=log1)
|
|
if (log1) call ishell(end)
|
|
if (scftype.eq.'rohf ') then
|
|
beg='DFINT_AIb.' // trim(ccmo)
|
|
end='rm -f ' // trim(beg)
|
|
inquire(file=trim(beg),exist=log1)
|
|
if (log1) call ishell(end)
|
|
beg='DFINT_IJb.' // trim(ccmo)
|
|
end='rm -f ' // trim(beg)
|
|
inquire(file=trim(beg),exist=log1)
|
|
if (log1) call ishell(end)
|
|
beg='DFINT_ABb.' // trim(ccmo)
|
|
end='rm -f ' // trim(beg)
|
|
inquire(file=trim(beg),exist=log1)
|
|
if (log1) call ishell(end)
|
|
c
|
|
beg='rLNO2cLNO.' // trim(ccmo)
|
|
call delete_file(beg, .true.)
|
|
beg='S_LNO.' // trim(ccmo)
|
|
call delete_file(beg, .true.)
|
|
endif
|
|
if (naf.ne.'off ') then
|
|
beg='ajb.' // trim(ccmo)
|
|
end='rm -f ' // trim(beg)
|
|
inquire(file=trim(beg),exist=log1)
|
|
if (log1) call ishell(end)
|
|
endif
|
|
if (talg.eq.'lapl'.or.talg.eq.'topr'.or.talg.eq.'to '
|
|
$ .or.talg.eq.'lato') then
|
|
beg='laplbas.' // trim(ccmo)
|
|
end='rm -f ' // trim(beg)
|
|
inquire(file=trim(beg),exist=log1)
|
|
if (log1) call ishell(end)
|
|
endif
|
|
c
|
|
return
|
|
end
|
|
C
|
|
************************************************************************
|
|
subroutine lccorbnuminit(nocorr,nbcorr,nocorrmin,nocorrmax,
|
|
$nbcorrmin,nbcorrmax,nelectronold,nbasisold)
|
|
************************************************************************
|
|
* initialize orbital size counters for localcc
|
|
************************************************************************
|
|
implicit none
|
|
integer nocorr,nbcorr,nocorrmin,nocorrmax,nbcorrmin,nbcorrmax
|
|
integer nelectronold,nbasisold ! full molecule quantities
|
|
C
|
|
nocorr=0
|
|
nbcorr=0
|
|
nocorrmin=nelectronold
|
|
nocorrmax=0
|
|
nbcorrmin=nbasisold
|
|
nbcorrmax=0
|
|
C
|
|
return
|
|
end
|
|
C
|
|
************************************************************************
|
|
subroutine postliscc(mode,nocc,nbasis,meth,label,n1,n2,n3,eref,
|
|
$naftyp,idegen,methact_inited,corembed,ncoi,ifcfile,
|
|
$laptoldefault,alaptol,localcc,ccprog,maxmet, ! for localcc.restart
|
|
$scrfile2,calctype,methact,
|
|
$nmetact,nbasiss,noccs,ncore,idom,nmet,
|
|
$nocorr,nbcorr,nocorrmin,nocorrmax,nbcorrmin,nbcorrmax,nact,
|
|
$ecl,totcpu,totwal,lmp2ene,mp2corr,pairpllno,eenv,
|
|
$emp2,epar,edecor,edeenv,eparenv,
|
|
$esec_im1,epar_im1,ecor_im1,eenv_im1,ccmem,
|
|
$intcpu,intwal,scfcpu,scfwal,orbcpu,orbwal,ovirtcpu,ovirtwal,
|
|
$prtrestart,mpi_rank,iout,group_master_thread)
|
|
************************************************************************
|
|
* process iface file, restore keywords, save to localcc.restart
|
|
************************************************************************
|
|
implicit none
|
|
integer mode,idegen,nocorr,nbcorr,n1,n2,n3,iout
|
|
integer nocorrmin,nocorrmax,nbcorrmin,nbcorrmax
|
|
integer imet,nmet,idom,methact_inited,nbasiss,noccs,maxmet
|
|
integer ncoi,nact,nmetact,ccmem(8),j,mpi_rank
|
|
real*8 degen,ecor,totcpu(maxmet),totwal(maxmet),ecl(maxmet)
|
|
real*8 mp2corr,ovirtcpu,ovirtwal,cpu,wal,lmp2ene,scfcpu,scfwal
|
|
real*8 ovirtcpud,ovirtwald,intcpu,intwal,eenv_im1
|
|
real*8 orbcpu,orbwal,eenv,pairpllno,esec_im1,epar_im1,ecor_im1
|
|
real*8 emp2,epar,edecor,edeenv,eparenv
|
|
character*4 ccprog,localcc
|
|
character*8 naftyp,corembed,idom_str
|
|
character*16 calctype,alaptol,rank_ch
|
|
character*128 label,meth(maxmet),methact(maxmet),ifacename
|
|
logical laptoldefault,localcc15p,prtrestart,lll,uccsd_file
|
|
logical group_master_thread
|
|
c from mrcccommon
|
|
integer nocc,nbasis,ifcfile,scrfile2,ncore
|
|
real*8 eref
|
|
C
|
|
cmpi write(*,*)'in postlcc',mpi_rank,group_master_thread
|
|
cmpi write(*,'(a12)',advance = 'no')'ecl postlcc'
|
|
cmpi write(*,'(i4,10es14.6)') mpi_rank,ecl(1:3)
|
|
cmpi write(*,*)'mp2corr postlcc',mpi_rank,mp2corr
|
|
if (mode.eq.1.and.group_master_thread) then ! do postliscc and save
|
|
nocorr=nocorr+idegen*nocc
|
|
nbcorr=nbcorr+idegen*nbasis/2
|
|
if(nocc.lt.nocorrmin) nocorrmin=nocc
|
|
if(nocc.gt.nocorrmax) nocorrmax=nocc
|
|
if(nbasis/2.lt.nbcorrmin) nbcorrmin=nbasis/2
|
|
if(nbasis/2.gt.nbcorrmax) nbcorrmax=nbasis/2
|
|
c
|
|
c process iface file
|
|
degen=dble(idegen)
|
|
ifacename='iface'
|
|
c#if defined (MPI)
|
|
c if (localcc.ne.'off '.and.mpi_rank.ne.-1) then
|
|
c write(rank_ch, '(i8)') mpi_rank
|
|
c ifacename='iface.' // trim(adjustl(rank_ch))
|
|
c endif
|
|
c#endif
|
|
open(unit=ifcfile,file=ifacename,status='old')
|
|
read(ifcfile,*)
|
|
nmet=0
|
|
do imet=1,maxmet
|
|
read(ifcfile,7596,end=86) label,meth(imet),n1,n2,n3,eref,
|
|
$cpu,wal
|
|
read(ifcfile,7596,end=86) label,meth(imet),n1,n2,n3,ecor,
|
|
$cpu,wal
|
|
if ((naftyp.eq.'jpqpilot'.or.naftyp.eq.'jpq ')
|
|
$ .and.imet.eq.1) mp2corr=mp2corr-degen*ecor ! Local MP2 correction for Jpq NAF truncation
|
|
nmet=nmet+1
|
|
if(ccprog.ne.'ccsd') ecor=ecor-eref
|
|
if (ncoi.eq.0) ecl(imet)=ecl(imet)+degen*ecor ! if embedded or corembed=off
|
|
c write(6,*) imet,degen,meth(imet),' ZZ ',ecor,ecl(imet)
|
|
totcpu(imet)=totcpu(imet)+cpu
|
|
totwal(imet)=totwal(imet)+wal
|
|
enddo
|
|
86 continue
|
|
close(ifcfile,status='delete')
|
|
c inquire(file = ifcfile, exist = lll)
|
|
c
|
|
do j=1,4
|
|
ccmem(j)=max(ccmem(j),ccmem(j+4))
|
|
enddo
|
|
c
|
|
c restore keywords and preparation for localcc.restart
|
|
if (methact_inited.eq.0) then
|
|
nmetact=max(nmet,nmetact)
|
|
methact(1:nmet)=meth(1:nmet)
|
|
methact_inited=1
|
|
endif
|
|
if (corembed.ne.'off ') then
|
|
if (ncoi.eq.0) then
|
|
nmetact=max(nmet,nmetact)
|
|
methact(1:nmetact)=meth(1:nmetact)
|
|
endif
|
|
if (nmetact.ne.0) nmet=nmetact
|
|
if (ncoi.eq.1) then
|
|
eenv=eenv+degen*ecor ! if embedded or corembed=off
|
|
if (.not.laptoldefault) then ! assuming environment is always normal
|
|
call changekey('laptol',6,alaptol,16) ! reset laptol to the value used for the active orbs.
|
|
endif
|
|
endif
|
|
endif
|
|
write(idom_str,'(i8)') idom
|
|
inquire(file='UCCSDREST.' // adjustl(trim(idom_str)),
|
|
$ exist=uccsd_file)
|
|
if (uccsd_file) then
|
|
! UCCSD restart file exists for current domain, delete it
|
|
open(97,file='UCCSDREST.' // adjustl(trim(idom_str)))
|
|
close(97,status='delete')
|
|
endif !uccsd_file
|
|
call changekey('ccsdrest',8,'off',3) ! reset ccsdrest (it was set to 'ccsd' if domain CC was restarted)
|
|
elseif (mode.eq.1.and..not.group_master_thread) then ! do postliscc and save
|
|
if (corembed.ne.'off '.and.ncoi.eq.1.and.
|
|
$ .not.laptoldefault) then
|
|
call changekey('laptol',6,alaptol,16) ! reset laptol to the value used for the active orbs.
|
|
endif
|
|
endif ! mode.eq.1
|
|
cmpi write(*,'(a12)',advance = 'no')'ecl after p'
|
|
cmpi write(*,'(i4,10es14.6)') mpi_rank,ecl(1:3)
|
|
cmpi write(*,*)'mp2corr after',mpi_rank,mp2corr
|
|
c save data to restart file: localcc.restart
|
|
if (localcc15p(localcc).and.prtrestart) then
|
|
rewind(scrfile2)
|
|
call localccrestart('save',scrfile2,calctype,methact,
|
|
$ nmetact,.false.,nbasiss,noccs,ncore,idom,nmet,
|
|
$ nocorr,nbcorr,nocorrmin,nocorrmax,nbcorrmin,nbcorrmax,nact,
|
|
$ ecl,totcpu,
|
|
$ totwal,lmp2ene,mp2corr,pairpllno,eenv,
|
|
$ emp2,epar,edecor,edeenv,eparenv,
|
|
$ esec_im1,epar_im1,ecor_im1,eenv_im1,ccmem,
|
|
$ intcpu,intwal,scfcpu,scfwal,orbcpu,orbwal,ovirtcpu,ovirtwal)
|
|
endif
|
|
c
|
|
7596 format(a8,1x,a15,1x,3(i2,1x),1pe23.15,2(1pe15.8))
|
|
C
|
|
return
|
|
end
|
|
c
|
|
************************************************************************
|
|
subroutine preliscc(mode,corembed,ncoi,inp,idegen,nact,ccenv,
|
|
$ccprog,accprog,eccprog,nbasis,calctype,envcalc,log1,
|
|
$cycle_this_domain,laptoldefault,elaptol,ifcfile,localcc,mpi_rank,
|
|
$domain_is_cs,lis_nal,lis_nbe,iout)
|
|
************************************************************************
|
|
* prepartion of CC calculation in the LIS
|
|
************************************************************************
|
|
implicit none
|
|
integer mode,ncoi,inp,idegen,nact,nbasis,ifcfile,mpi_rank,lis_nal
|
|
integer lis_nbe,iout
|
|
character*4 ccprog,accprog,eccprog,localcc
|
|
character*8 corembed,scftype
|
|
character*16 calctype,envcalc,elaptol,rank_ch,liscalctype
|
|
character*128 ifacename
|
|
logical ccenv,log1,cycle_this_domain,laptoldefault,domain_is_cs
|
|
c
|
|
if (mode.eq.2) then
|
|
ifacename='iface'
|
|
c#if defined (MPI) ! not needed if all LISs have dedicated folders
|
|
c if (localcc.eq.'2022'.and.mpi_rank.ne.-1) then ! -1 is set in call from locdrv, trffirst
|
|
c write(rank_ch, '(i8)') mpi_rank
|
|
c ifacename='./iface.' // trim(adjustl(rank_ch))
|
|
c endif
|
|
c#endif
|
|
c write(iout,*) 'preliscc ifcfile',mpi_rank,ifacename
|
|
C Initialize interface file if called from edcomp of ldrpa
|
|
open(unit=ifcfile,file=ifacename,status='replace')
|
|
write(ifcfile,"(68a)")
|
|
$ '#property method sym st mul value' //
|
|
$ ' CPU(sec) Wall(sec)'
|
|
close(ifcfile)
|
|
endif
|
|
c
|
|
liscalctype = calctype
|
|
if (corembed.ne.'off ') then
|
|
if (mode.eq.1) then ! calling from locdrv of dmrcc
|
|
if (ncoi.eq.0) close(inp)
|
|
if (ncoi.eq.1.and..not.ccenv) close(inp,status='delete')
|
|
if (ncoi.eq.0.and.idegen.gt.0) nact=nact+idegen ! number of active (not environment) orbitals
|
|
endif
|
|
if (ccenv) then
|
|
if (ncoi.eq.0) then
|
|
ccprog=accprog
|
|
if (ccprog.eq.'mrcc') liscalctype = calctype
|
|
elseif (ncoi.eq.1) then
|
|
ccprog=eccprog
|
|
if (ccprog.eq.'mrcc') liscalctype = envcalc
|
|
endif ! ncoi.eq.0
|
|
endif ! ccenv
|
|
endif
|
|
if (ccprog.eq.'mrcc') then
|
|
call getkey('scftype',7,scftype,8)
|
|
if ((trim(corembed).ne.'off'.and.ccenv)
|
|
$ .or.trim(scftype).ne.'rhf') then
|
|
call write56(nbasis,'none ',liscalctype,domain_is_cs,
|
|
$ lis_nal,lis_nbe)
|
|
endif
|
|
endif !ccprog
|
|
c
|
|
if (log1.and..not.cycle_this_domain) then
|
|
if (ccenv) then
|
|
if (ncoi.eq.0) then
|
|
call changekey('calc',4,calctype,16)
|
|
elseif (ncoi.eq.1) then
|
|
call changekey('calc',4,envcalc,16)
|
|
if (.not.laptoldefault) then ! assuming environment is always normal
|
|
call changekey('laptol',6,elaptol,16) ! reset laptol to the value used for the environment orbs.
|
|
endif
|
|
endif
|
|
endif
|
|
endif
|
|
c
|
|
return
|
|
end
|
|
C
|
|
#ifdef MPI
|
|
subroutine split_communicator(communicator, n)
|
|
implicit none
|
|
include "mpif.h"
|
|
|
|
integer, intent(in) :: n
|
|
integer, intent(inout) :: communicator
|
|
|
|
integer :: rank, error, key, color
|
|
integer :: new_communicator
|
|
|
|
call MPI_Comm_rank(communicator, rank, error)
|
|
|
|
! split communicator after n ranks
|
|
if(rank < n) then
|
|
color = 0
|
|
else
|
|
color = 1
|
|
end if
|
|
key = rank
|
|
call MPI_Comm_split(communicator, color, key, new_communicator,
|
|
$ error)
|
|
|
|
! overwrite old communicator with the new one
|
|
!call MPI_Comm_free(communicator, error)
|
|
call MPI_Comm_dup(new_communicator, communicator, error)
|
|
!call MPI_Comm_free(new_communicator, error)
|
|
|
|
end subroutine
|
|
#endif
|