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

575 lines
20 KiB
Fortran

program orbloc
C
#include "MRCCCOMMON"
C
real*8 m,mc,mo,mv
real*8 itol,tresh_macro,tresh_micro
character*4 citol,loctypec,loctypeo,loctypev,localcc,ccprog
character*5 boys,gboys(2),scftype
character*5 orbloco5,orblocv5,orblocc5
character*13 rohftype
character*16 orblocc,orbloco,orblocv,c16
character*8 embed,corembed,qmreg,cscr,orblocguess,boysalg
character*6 core
integer mmbond,nc
logical lembed
logical locvirt,lcvs
integer icvsorb,nocc_cvs,iscr,ncore2
integer i,inatrange,is,iv,intalloc,dblalloc,natoms,iind,icoord
integer idistc,idisti,idte,ifock,iatmo,j,iatsymbol,printfile,nos
integer inatmo,imoat,inmoat,iatdom,imodom,inatdom,inmodom,iatchg
integer mnbasis,iam,is1,is1is,is2is,ict,io,iot,ip12,id,imoa,imob
integer iuboys
parameter(printfile=scrfile6)
equivalence(orblocc5,orblocc)
equivalence(orbloco5,orbloco)
equivalence(orblocv5,orblocv)
C
call mrccini
C
call getkey('boysalg',7,boysalg,8)
call getkey('orblocguess',11,orblocguess,8)
C Open files
open(oeintfile,file='OEINT',form='UNFORMATTED')
open(printfile,file='PRINT',form='UNFORMATTED')
if (orblocguess.eq.'read '.or.orblocguess.eq.'restart ') then
call ishell('cp MOCOEF.LOC MOCOEF')
endif
C Read variables
call getkey('orblocc',7,orblocc,16)
call getkey('orbloco',7,orbloco,16)
call getkey('orblocv',7,orblocv,16)
call getkey('scftype',7,scftype,5)
call getkey('rohftype',8,rohftype,13)
call getkey('localcc',7,localcc,4)
call getkey('embed',5,embed,8)
call getkey('corembed',8,corembed,8)
call getvar('nbasis ',nbasis)
call getvar('nal ',nal)
call getvar('nbe ',nbe)
call getvar('ncore ',ncore)
call getvar('nbf ',nbf)
call getkey('core',4,core,6)
if(core.eq.'corr ') ncore=0
lembed = embed.ne.'off '
if(trim(orbloco).eq.'fragment'.or.
c $ trim(orbloco).eq.'no'.or.
$ trim(orbloco).eq.'spade') ncore=0
if(scftype.eq.'rohf '.and.rohftype.eq.'semicanonical'
$ .and.localcc.eq.'off ') scftype='uhf '
nelec=nal+nbe
nvirtal=nbasis-nal
nvirtbe=nbasis-nbe
nal=nal-ncore
nbe=nbe-ncore
nocc=min(nal,nbe)
nvirt=min(nvirtal,nvirtbe)
nos=max(nal,nbe)-min(nal,nbe)
C
CCC start HB
nfroz=0
! In case of Huzinaga QM/MM + embedding, the frozen
! boundary orbitals should be unchanged during the 2nd
! SCF run, hence localization of these localized orbitals
! only results in unnecessary mixing frozen MOs.
open(minpfile,file='MINP')
call getkeym('qmreg',5,qmreg,8)
close(minpfile)
call getvar('mmbond ',mmbond)
if(qmreg.eq.'0 '.and.(embed.ne.'off '.or.
$corembed.ne.'off ')) then
do i=1,mmbond
write(cscr,'(i8)') i
open(scrfile1,file='MOCOEF.' // trim(adjustl(cscr)),
$form='unformatted')
read(scrfile1) nc
close(scrfile1)
nfroz=nfroz+nc
enddo
endif
CCC end HB
!!!!!!!!!!!!!!! HB
C Allocate memory
call memalloc
C Read MO coefficients
imoa=dblalloc(nbasis**2)
open(mocoeffile,file='MOCOEF',form='unformatted')
call readmo(dcore(imem),dcore(imem),dcore(imoa),mocoeffile,
$nbasis,nbasis)
if(scftype.eq.'uhf ') then
imob=dblalloc(nbasis**2) ! It must be right after moa
call readmo(dcore(imem),dcore(imem),dcore(imob),mocoeffile,
$nbasis,nbasis)
else
imob=imoa
endif
C Set tolerance of zero elements of MO coefficients
call getkey('scftol',6,citol,4)
read(citol,*) i
itol=dsqrt(10.d0**(-i))
tresh_macro=itol
tresh_micro=itol
C
C Reading keys
call getvar('nbset ',nbset)
mc=1.d0
mo=1.d0
mv=1.d0
loctypec='boys'
loctypeo='boys'
loctypev='boys'
locvirt=trim(orblocv).ne.'off'.and.
$ trim(orblocv).ne.'pao '
iuboys=1
if(trim(orblocc).ne.'off') then
if(trim(orblocc).eq.'cholesky') then
loctypec='chol'
else if(trim(orblocc).eq.'pm') then
loctypec='pm '
else if(trim(orblocc).eq.'boys') then
loctypec='boys'
else if(trim(orblocc).eq.'ibo') then
loctypec='ibo '
else if(trim(orblocc).eq.'fragment') then
loctypec='frag'
else if(trim(orblocc).eq.'no') then
loctypec='no '
else if(trim(orblocc).eq.'spade') then
loctypec='spad'
else if(orblocc5.eq.'gboys') then
call getkey('orblocc',7,gboys,10)
read(gboys(2),*) mc
endif
endif
if(trim(orbloco).ne.'off') then !HB
if(trim(orbloco).eq.'cholesky') then
loctypeo='chol'
else if(trim(orbloco).eq.'pm') then
loctypeo='pm '
else if(trim(orbloco).eq.'boys') then
loctypeo='boys'
else if(trim(orbloco).eq.'ibo') then
loctypeo='ibo '
else if(trim(orbloco).eq.'fragment') then
loctypeo='frag'
else if(trim(orbloco).eq.'no') then
loctypeo='no '
else if(trim(orbloco).eq.'spade') then
loctypeo='spad'
else if(orbloco5.eq.'gboys') then
call getkey('orbloco',7,gboys,10)
read(gboys(2),*) mo
endif
endif !HB
if(trim(orblocv).ne.'off') then
if(trim(orblocv).eq.'cholesky') then
loctypev='chol'
else if(trim(orblocv).eq.'pm') then
loctypev='pm '
else if(trim(orblocv).eq.'boys') then
loctypev='boys'
else if(trim(orblocv).eq.'ibo') then
loctypev='ibo '
else if(trim(orblocv).eq.'fragment') then
loctypev='frag'
else if(trim(orblocv).eq.'no') then
loctypev='no '
else if(trim(orblocv).eq.'spade') then
loctypev='spad'
else if(trim(orblocv).eq.'pao-subsys ') then
loctypev='pao'
else if(orblocv5.eq.'gboys') then
call getkey('orblocv',7,gboys,10)
read(gboys(2),*) mv
endif
endif
if(loctypec.eq.'boys'.or.loctypeo.eq.'boys'.or.
& loctypev.eq.'boys') iuboys=dblalloc(nocc*nocc)
call getvar('natoms ',natoms)
C Calculation of IAOs
if(trim(orblocc).eq.'ibo'.or.trim(orbloco).eq.'ibo'.or.
$ trim(orblocv).eq.'ibo ') then
write(iout,*)
write(iout,*) 'Constructing IAOs...'
mnbasis=nbf(4)
iam=dblalloc(nbasis**2)
is1=dblalloc(nbasis**2)
is1is=dblalloc(nbasis**2)
is2is=dblalloc(nbasis**2)
ict=dblalloc(nbasis**2)
io=dblalloc(nbasis**2)
iot=dblalloc(nbasis**2)
ip12=dblalloc(nbasis**2)
id=dblalloc(nbasis**2)
open(scrfile1,file='SCFDENSITIES',form='UNFORMATTED')
call roeint(dcore(imem),dcore(imem),dcore(id),scrfile1,nbasis)
close(scrfile1)
call roeint(dcore(imem),dcore(imem),dcore(is1),oeintfile,nbasis)
open(scrfile1,file='SROOT',form='unformatted')
read(scrfile1)
call roeint(dcore(imem),dcore(imem),dcore(is1is),scrfile1,
$nbasis)
close(scrfile1)
open(scrfile1,file='S12MAT',form='UNFORMATTED')
call roeint(dcore(imem),dcore(imem),dcore(is2is),scrfile1,
$mnbasis)
call rtdmx (dcore(imem),dcore(imem),dcore(iam),scrfile1,
$nbasis,mnbasis)
close(scrfile1)
call iao(nbasis,mnbasis,nelec/2,dcore(is1),dcore(is1is),
$dcore(is2is),dcore(iam),dcore(imoa),dcore(ict),dcore(imem),
$dcore(io),dcore(iot),dcore(ip12),dcore(id))
inatrange=intalloc(2*natoms*nbset)
iatsymbol=intalloc(natoms)
iatchg=dblalloc(natoms)
call iaoprt(nbasis,natoms,icore(inatrange),dcore(iatchg),
$icore(iatsymbol),dcore(is2is),dcore(iam),iout)
call dbldealloc(is1)
open(scrfile1,file='IAO',form='UNFORMATTED')
call wtdmx(dcore(imem),dcore(imem),dcore(iam),scrfile1,itol,
$nbasis,mnbasis)
close(scrfile1)
else
iam=imem
endif
C Localization of core orbitals
C
call getkey('ccprog',6,ccprog,4)
call getkey('cvs',3,c16,16)
lcvs=trim(c16).ne.'off'.and.trim(ccprog).eq.'cis'
if(lcvs) then
icvsorb=intalloc(2*nocc)
open(minpfile,file='MINP')
call embedat(nocc,minpfile,icore(icvsorb),iout,'cvs ',3)
close(minpfile)
call dcopy(nocc,icore(icvsorb),1,icore(icvsorb+nocc),1)
nocc_cvs=0
do i=1,nocc
if(icore(icvsorb-1+i).eq.1) then
nocc_cvs=nocc_cvs+1
icore(icvsorb-1+nocc_cvs)=i
endif
enddo
icore(icvsorb+nocc_cvs:icvsorb+nocc-1)=0
C
iscr=dblalloc(nbasis*(nocc+ncore))
call cvs_sort_loc(nbasis,nocc,nocc_cvs,dcore(imoa),dcore(iscr),
$icore(icvsorb),1,icore(icvsorb+nocc),ncore)
C
ncore2=ncore
ncore=nocc_cvs
nocc=nocc-nocc_cvs
C
endif
C
if((trim(orblocc).ne.'off'.and.ncore.gt.0).or.lcvs) then
write(*,"(/,' Localization of core orbitals...')")
if (orblocguess.eq.'read '.or.orblocguess.eq.'restart ')
$ write(iout,*)'Guess orbitals are read from the MOCOEF.LOC file'
if(loctypec.eq.'chol') then
write(*,"(' Localization type: Cholesky decomposition')")
else if(loctypec.eq.'pm ') then
write(*,"(' Localization type: Pipek-Mezey')")
else if(loctypec.eq.'ibo ') then
write(*,
$"(' Localization type: intrinsic bond orbitals (IBOs)')")
else if(loctypec.eq.'frag') then
write(*,"(' Localization type: Fragment')")
else if(loctypec.eq.'no ') then
write(*,"(' Localization type: Natural orbitals')")
else if(loctypec.eq.'spad') then
write(*,"(' Localization type: SPADE')")
else if(mc.eq.1.0d0) then
write(*,"(' Localization type: Boys')")
if(boysalg.eq.'jacobi ') then
write(*,"(' Localization algorithm: Jacobi')")
elseif(boysalg.eq.'newton ') then
write(*,"(' Localization algorithm: Newton')")
endif
else
write(*,"(' Localization type: Generalized Boys ',
& '(m = ',f5.1,')')") mo
endif
write(*,*)
call flush(6)
CCC start HB
call localize('c',ncore-nfroz,loctypec,mc,tresh_macro,
$tresh_micro,itol,natoms,dcore(iam),mnbasis,printfile,orblocguess,
$nos,dcore(imoa),nocc,lembed,dcore(iuboys),.false.)
if(scftype.eq.'uhf ') then
call localize('c',ncore-nfroz,loctypec,mc,tresh_macro,
$tresh_micro,itol,natoms,dcore(iam),mnbasis,printfile,orblocguess,
$nos,dcore(imob),nocc,lembed,dcore(iuboys),.true.)
endif
CCC end HB
endif
C Localization of occupied orbitals
if(trim(orbloco).ne.'off') then
if (nocc.gt.0)
$ write(*,"(/,' Localization of occupied orbitals...')")
if (orblocguess.eq.'read '.or.orblocguess.eq.'restart ')
$ write(iout,*)'Guess orbitals are read from the MOCOEF.LOC file'
if (loctypeo.eq.'chol') then
write(*,"(' Localization type: Cholesky decomposition')")
elseif (loctypeo.eq.'pm ') then
write(*,"(' Localization type: Pipek-Mezey')")
else if(loctypeo.eq.'ibo ') then
write(*,
$"(' Localization type: intrinsic bond orbitals (IBOs)')")
else if(loctypeo.eq.'frag') then
write(*,"(' Localization type: Fragment')")
else if(loctypeo.eq.'no ') then
write(*,"(' Localization type: Natural orbitals')")
else if(loctypeo.eq.'spad') then
write(*,"(' Localization type: SPADE')")
elseif (mo.eq.1.d0) then
write(*,"(' Localization type: Boys')")
if(boysalg.eq.'jacobi ') then
write(*,"(' Localization algorithm: Jacobi')")
elseif(boysalg.eq.'newton ') then
write(*,"(' Localization algorithm: Newton')")
endif
else
write(*,"(' Localization type: Generalized Boys ',
& '(m = ',f5.1,')')") mo
endif
write(*,*)
call flush(6)
if (nal.gt.0) then
if(scftype.eq.'uhf ') then
call localize('o',nal,loctypeo,mo,tresh_macro,tresh_micro,itol,
$natoms,dcore(iam),mnbasis,printfile,orblocguess,nos,dcore(imoa),
$nocc,lembed,dcore(iuboys),.false.)
endif
endif
if (nocc.gt.0) then
call localize('o',nocc,loctypeo,mo,tresh_macro,tresh_micro,itol,
$natoms,dcore(iam),mnbasis,printfile,orblocguess,nos,dcore(imob),
$nocc,lembed,dcore(iuboys),.true. .and. scftype.eq.'uhf ')
endif
C
if(lcvs) then
call cvs_sort_loc(nbasis,nocc+nocc_cvs,nocc_cvs,dcore(imoa),
$dcore(iscr),icore(icvsorb),2,icore(icvsorb+nocc+nocc_cvs),ncore2)
call dbldealloc(iscr)
call intdealloc(icvsorb)
endif
C
if (ncore-nfroz.le.0.and. ! core e only is allowed, eg. Na+
$ ((nocc.le.0.and.nos.le.0).or. ! exit if 0 core + 0 valence (RHF,ROHF)
$ (scftype.eq.'uhf '.and.nal.le.0)))then! exit if 0 core + 0 valence (UHF)
write(*,"(' Error: zero occupied orbitals in localization')")
call mrccend(1)
endif
if(scftype.eq.'rohf '.and.nos.gt.0) then
write(*,"(/,' Localization of singly-occupied orbitals...')")
write(*,*)
call flush(6)
call localize('s',nos,loctypeo,mo,tresh_macro,tresh_micro,itol,
$natoms,dcore(iam),mnbasis,printfile,orblocguess,nos,dcore(imoa),
$nocc,lembed,dcore(iuboys),.false.)
endif
endif
C Localization of virtual orbitals
if (locvirt.eqv..true.) then
write(*,"(/,' Localization of virtual orbitals...')")
if (orblocguess.eq.'read '.or.orblocguess.eq.'restart ')
$ write(iout,*)'Guess orbitals are read from the MOCOEF.LOC file'
if (loctypev.eq.'chol') then
write(*,"(' Localization type: Cholesky decomposition')")
elseif (loctypev.eq.'pm ') then
write(*,"(' Localization type: Pipek-Mezey')")
else if(loctypev.eq.'ibo ') then
write(*,
$"(' Localization type: intrinsic bond orbitals (IBOs)')")
else if(loctypev.eq.'frag') then
write(*,"(' Localization type: Fragment')")
else if(loctypev.eq.'no ') then
write(*,"(' Localization type: Natural orbitals')")
else if(loctypev.eq.'spad') then
write(*,"(' Localization type: SPADE')")
else if(loctypev.eq.'pao ') then
write(*,"(' Localization type: Projected Atomic Orbitals')")
write(*,*)
write(*,"('
$(using the projector of the full occupied subspace) ')")
elseif (mv.eq.1.d0) then
write(*,"(' Localization type: Boys')")
if(boysalg.eq.'jacobi ') then
write(*,"(' Localization algorithm: Jacobi')")
elseif(boysalg.eq.'newton ') then
write(*,"(' Localization algorithm: Newton')")
endif
else
write(*,"(' Localization type: Generalized Boys ',
& '(m = ',f5.1,')')") mv
endif
write(*,*)
call flush(6)
call localize('v',nvirt,loctypev,mv,tresh_macro,tresh_micro,
&itol,natoms,dcore(iam),mnbasis,printfile,orblocguess,nos,
$dcore(imoa),nal,lembed,dcore(iuboys),.false.)
if(scftype.eq.'uhf ') then
call localize('v',nvirtbe,loctypev,mv,tresh_macro,tresh_micro,
&itol,natoms,dcore(iam),mnbasis,printfile,orblocguess,nos,
$dcore(imob),nbe,lembed,dcore(iuboys),.true.)
endif
endif
write(*,"(/,' Localization succeeded.',/)")
C Save localized MOs
rewind(mocoeffile)
write(iout,
$ '(" File MOCOEF is overwritten with localized orbitals.")')
write(iout,*)
call wrtmo(dcore(imem),dcore(imem),dcore(imoa),mocoeffile,0.d0,
$nbasis,nbasis)
if(scftype.eq.'uhf ')
$call wrtmo(dcore(imem),dcore(imem),dcore(imob),mocoeffile,0.d0,
$nbasis,nbasis)
C%%% Memoriafoglalas!
C Write to MOLDEN file
call molden(dcore(imoa),dcore(imem),dcore(imem+2*nbasis),
$icore(imem+4*nbasis),icore(imem+5*nbasis),nbasis,moldenfile,
$scftype)
if(loctypev.eq.'pao') then
call edit_molden(scrfile1,scrfile2,nbasis,ncore,nocc,nvfroz)
endif
C Close files
close(oeintfile)
close(printfile)
close(mocoeffile)
C Timing
call timer
C
call mrccend(0)
end
C
************************************************************************
subroutine molden(moa,eigenvalue,occupation,perm,spin,nbasis,
$moldenfile,scftype)
************************************************************************
* Write molden file
************************************************************************
implicit none
integer nbasis,i,j,k,ii,jj,nn,perm(nbasis),moldenfile
real*8 eigenvalue(nbasis,2),occupation(nbasis,2)
real*8 moa(nbasis,nbasis,2)
character trash(7),spin(12*nbasis,2)
character*4 mark
character*5 scftype
character*6 mold
character*15 cscr
C
call getkey('molden',6,mold,6)
if(mold.eq.'off ') return
nn=1
if(scftype.eq.'uhf ') nn=2
C Write MO coefficients to moldenfile
open(unit=moldenfile,file='MOLDEN.perm')
read(moldenfile,*) (perm(i),i=1,nbasis)
close(moldenfile)
open(moldenfile,file='MOLDEN')
mark=' '
do while (mark.ne.'[MO]')
read(moldenfile,"(4a)")mark
enddo
do k=1,nn
do i=1,nbasis
read(moldenfile,"(5a,f14.4)")
& (trash(j),j=1,5),eigenvalue(i,k)
read(moldenfile,"(12a)")
& (spin(j,k),j=(i-1)*12+1,i*12)
read(moldenfile,"(7a,f14.4)")
& (trash(j),j=1,7),occupation(i,k)
do j=1,nbasis
read(moldenfile,*)
enddo
enddo
enddo
rewind(moldenfile)
mark='aaaa'
do while (mark.ne.'[MO]')
read(moldenfile,"(4a)")mark
enddo
do k=1,nn
do i=1,nbasis
write(moldenfile,"(' Ene=',f14.4)") eigenvalue(i,k)
write(moldenfile,"(12a)")
& (spin(j,k),j=(i-1)*12+1,i*12)
write(moldenfile,"(' Occup=',f14.4)") occupation(i,k)
do j=1,nbasis
write(moldenfile,"(i5,f18.10)") j,moa(perm(j),i,k)
enddo
enddo
enddo
close(moldenfile)
C
end
C
************************************************************************
subroutine iaoprt(nbasis,natoms,natrange,atchg,atsymbol,m1,m2,
$iout)
************************************************************************
* Print IAO charges
************************************************************************
implicit none
integer nbasis,natoms,iatoms,natrange(2,natoms,*),iout
real*8 atchg(natoms),ddot,m1(nbasis,nbasis),m2(nbasis,nbasis)
character*2 atsymbol(natoms)
C
call getvar('atchg ',atchg)
call getvar_c('atsymbol ',atsymbol)
call getvar('natrange ',natrange)
write(iout,*)
write(iout,*) ' IAO atomic charges'
do iatoms=1,natoms
write(iout,'(i4,1x,a3,3f12.6)') iatoms,atsymbol(iatoms),
$atchg(iatoms)-ddot(nbasis*(natrange(2,iatoms,4)-
$natrange(1,iatoms,4)),m2(1,natrange(1,iatoms,4)+1),1,
$m1(1,natrange(1,iatoms,4)+1),1)
enddo
C
return
end
C
************************************************************************
subroutine cvs_sort_loc(nbasis,nocc,nocc_cvs,cmo,scr,cvsorb,flag,
$cvsorb2,ncore)
************************************************************************
implicit none
integer nbasis,nocc,nocc_cvs,cvsorb(nocc),i,flag,cvsorb2(nocc),ii
integer jj,ncore
real*8 cmo(nbasis,ncore+nocc),scr(nbasis,ncore+nocc)
C
if(flag.eq.1) then
call dcopy(nbasis*(ncore+nocc),cmo,1,scr,1)
ii=0
jj=0
do i=1,nocc
if(cvsorb2(i).eq.1) then
ii=ii+1
cmo(1:nbasis,ii)=scr(1:nbasis,ncore+i)
else
jj=jj+1
cmo(1:nbasis,nocc_cvs+jj)=scr(1:nbasis,ncore+i)
endif
enddo
elseif(flag.eq.2) then
do i=1,nocc_cvs
scr(1:nbasis,ncore+cvsorb(i))=cmo(1:nbasis,i)
enddo
C
ii=0
do i=1,nocc
if(cvsorb2(i).eq.0) then
ii=ii+1
scr(1:nbasis,ncore+i)=cmo(1:nbasis,nocc_cvs+ii)
endif
enddo
C
call dcopy(nbasis*(ncore+nocc),scr,1,cmo,1)
endif
C
return
end