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

786 lines
28 KiB
Fortran

C***********************************************************************
subroutine bopu(nbasis,nmatrix,natoms,natrange,c,v,ind,r8heap,
$s,sa,vi,ai,atmo,nmoat,moat,natmo,bpcompo,bpcompv,nocc,natmaxo,
$meth,i8heap,mnatrange,mnbasis,cc,distc,ndao,ovltol,mocomp,pmul,
$iout,minlc)
C***********************************************************************
C Algorithms for assigning atoms to localized MOs
C***********************************************************************
implicit none
integer i,iatoms,jatoms,iiatoms,jjatoms,mu,nu,nnatoms,ii,jj
integer nmatrix,natoms,natrange(2,natoms),ind(natoms),nbasis
integer atmo(natoms,nmatrix),nmoat(natoms),moat(natoms,nmatrix)
integer natmo(nmatrix),natmaxo,natmax,i8heap(*),nocc
integer mnatrange(2,natoms),mnbasis,ndao,sadim
real*8 smm,c(nbasis,nmatrix),v(nbasis,nmatrix),r8heap(*),ddot,tol
real*8 s(nbasis,nbasis),sa(nbasis,nbasis),vi(nbasis),ai(nbasis)
real*8 sum1,bpcompo,bpcompv,mch,cc(mnbasis,nmatrix),norm !,smch
real*8 distc(natoms,natoms),r8dummy,ovltol
real*8,pointer :: psa(:,:)
character*4 meth,mode,citol
c internal variables !NP
logical lll,remove,distrearr
real*8 dmch(natoms+1),psmm,chgtol,mchtol,fact,smch,mchtol0,minlc
integer ind2(natoms+1),knatoms,ind3(natoms),iout
character*16 c16
real*8 itol !SB
real*8 mocomp(*) !HB
logical pmul !HB
C
interface
c {{{ interfaces for pointers
subroutine rpoint2d2(egydim,haromdim,dim1,dim2)
implicit none
integer dim1,dim2
real*8,target :: egydim(dim1,dim2)
real*8, pointer :: haromdim(:,:)
end subroutine
c}}}
end interface
c
r8dummy=0.d0
remove=.false.
distrearr=.false.
call getkey('itol',4,citol,4)
read(citol,*) i
itol=10.d0**(-i)
c if (meth.eq.'nbp2') then
c meth='nbp '
c mchtol0=0.5d0
c endif
if(meth.eq.'nbp '.or.meth.eq.'niao') then
remove=.true.
distrearr=.true.
call getvar('distc ',distc)
if (meth.eq.'nbp ') mode='bopu'
if (meth.eq.'niao') mode='iao '
fact=5.d0
else
mode=meth
endif
if(mode.ne.'lowd') then
write(iout,*)
write(iout,*) ' MO Natoms Completeness Atoms'
endif
call ifillzero(nmoat,natoms)
do i=1,nmatrix
if (mode.eq.'iao '.and.i.le.nocc) then
norm=1.d0 ! warning: need to compute actual norm if it is not 1.d0, e.g. bopu is called from other than ldrpa
else
norm=ddot(nbasis,c(1,i),1,v(1,i),1) ! norm of i'th MO, can be 2 if called from closed shell SCF
endif
tol=bpcompo
natmax=natmaxo
mchtol0=minlc*norm
mchtol=mchtol0
if(i.gt.nocc) then
tol=bpcompv
natmax=natmaxo !natmaxv
endif
if (remove) mchtol=min((norm-tol)*fact,mchtol0)
c if (dabs(tol-norm).lt.1.d-12) then
if (dabs(tol-norm).lt.max(1000.d0*itol,ovltol)) then
c add all atoms to the BP domain automatically
natmo(i)=natoms
do iatoms=1,natoms
atmo(iatoms,i)=iatoms
nmoat(iatoms)=nmoat(iatoms)+1
moat(iatoms,nmoat(iatoms))=i
enddo
if(i.eq.nocc+1) write(iout,*)
write(iout,"(2i5,4x,f9.6,4x,100i4)") i,natoms,tol,
$(atmo(iatoms,i),iatoms=1,natoms)
cycle
endif
call dfillzero(dmch,natoms)
c sum1=0.d0
if(mode.eq.'iao '.and.i.le.nocc) then
do iatoms=1,natoms
smm=0.d0
do mu=mnatrange(1,iatoms)+1,mnatrange(2,iatoms)
smm=smm+cc(mu,i)*cc(mu,i)
enddo
c smm=2.d0*smm !NP
r8heap(iatoms)=-dabs(smm) !minus for qsortd!
r8heap(natoms+iatoms)=smm
enddo
else
do iatoms=1,natoms
smm=0.d0
do mu=natrange(1,iatoms)+1,natrange(2,iatoms)
smm=smm+c(mu,i)*v(mu,i)
enddo
r8heap(iatoms)=-dabs(smm) !minus for qsortd!
r8heap(natoms+iatoms)=smm
enddo
endif
c write(iout,"(1x,a7,f10.6)") 'Charge:',smm
c write(iout,"(1000f10.6)") (r8heap(natoms+ii),ii=1,natoms)
call qsortd(ind,natoms,r8heap)
c write(iout,*)
c write(iout,"('srtd abs Mull chargeof',i6,22f10.6)") i,
c $ (r8heap(natoms+ind(ii)),ii=1,min(10,natoms)),
c $ (dble(ind(ii)),ii=1,min(10,natoms)),
c $-sum(r8heap(1:natoms)),sum(r8heap(natoms+1:2*natoms))
c
c {{{ mode.eq.'bopu'.or.mode.eq.'iao ').and.distrearr !NP
if((mode.eq.'bopu'.or.mode.eq.'iao ').and.distrearr) then
c rearrange atom order in ind in case of meaninglessly small charges based on distance from the most important atom, ind(1)
chgtol=min(mchtol,0.001d0*norm)
ii=ind(1)
sa(1:natoms,1)=r8heap(1:natoms)
sa(1:natoms,2)=r8heap(natoms+1:2*natoms)
i8heap(1:natoms)=ind(1:natoms)
knatoms=0
c collect distances for atoms with small charge
do jj=1,natoms
if (dabs(sa(ind(jj),1)).lt.chgtol) then
knatoms=knatoms+1
i8heap(natoms+knatoms)=ind(jj)
ai(knatoms)=distc(ii,ind(jj))
endif
enddo
call qsortd(i8heap(2*natoms+1),knatoms,ai)
nnatoms=natoms-knatoms ! num of atoms with high enough charge
do jj=1,knatoms
ind(nnatoms+jj)=i8heap(natoms+i8heap(2*natoms+jj))
enddo
c write(iout,"('uj ind',1000i5)") (ind(ii),ii=1,natoms)
r8heap(1:natoms)=sa(1:natoms,1)
r8heap(natoms+1:2*natoms)=sa(1:natoms,2)
endif
c }}}
smm=0.d0
if(mode.eq.'lowd'.or.mode.eq.'mull') then
C Domains are constructed on the basis of Lowdin/Mulliken atomic charges
nnatoms=0
do iatoms=1,min(natoms,natmax)
if(dabs(r8heap(ind(iatoms))).gt.tol) then
nnatoms=nnatoms+1
smm=smm+r8heap(natoms+ind(iatoms))
endif
enddo
if (nnatoms.eq.0) then
write(iout,*) 'warning: 0 atoms in LDF domain of LMO',i
do iatoms=1,min(natoms,natmax)
if(smm.lt.2.d0*tol) then ! initialize LDF to avoid 0 selected atoms
nnatoms=nnatoms+1
smm=smm+r8heap(natoms+ind(iatoms))
endif
enddo
endif
else
C Domains are constructed according to the Boughton-Pulay algorithm
nnatoms=0
mch=0.d0
c smch=0.d0
psmm=0.d0
do while((smm.lt.tol.or.mch.gt.mchtol).and.
$nnatoms.lt.min(natoms,natmax))
nnatoms=nnatoms+1
lll=.true.
sum1=r8heap(ind(nnatoms))
c smch=smch+r8heap(natoms+ind(nnatoms))
do while(nnatoms.lt.min(natoms,natmax).and.lll)
lll=.false.
if((dabs(sum1-r8heap(ind(nnatoms+1))).lt.1d-6).and.
$ .not.remove) then
lll=.true.
nnatoms=nnatoms+1
c smch=smch+r8heap(natoms+ind(nnatoms))
endif
enddo
c project mo onto the selected mo domain and obtain completeness measure
call bpcompleteness(nnatoms,ind,s,v,natrange,sa,vi,ai,smm,
$r8heap(2*natoms+1),ndao,ovltol,nbasis,i,.true.,r8dummy,sadim)
dmch(nnatoms)=-smm+psmm ! completeness increment due to nnatoms^th atom*(-1)
psmm=smm
cc write(iout,*) i,nnatoms,'add atom',ind(nnatoms)
cc write(iout,*) smm,dmch(nnatoms)
c write(iout,"(1000f10.6)") vi
c write(iout,"(1000f10.6)") ai
if(mode.ne.'iao '.and.i.le.nocc.and.
$nnatoms.lt.min(natoms,natmax)) mch=dabs(r8heap(ind(nnatoms+1)))
enddo
endif ! mode
c call qsortint(i8heap,nnatoms,ind)
if (remove) then
c {{{ remove atoms with unnecessary small contribution
call qsortd(ind2,nnatoms,dmch)
knatoms=0
smm=0.d0
psmm=0.d0
do while(smm.lt.tol.and.knatoms.le.nnatoms)
lll=.true.
do while((knatoms.lt.nnatoms).and.lll)
lll=.false.
if(psmm-dmch(ind2(knatoms+1)).lt.tol.or.knatoms.eq.0) then ! knatoms.eq.0 for the case when completeness > tol already with 1 atom
lll=.true.
psmm=psmm-dmch(ind2(knatoms+1))
knatoms=knatoms+1
cc write(iout,*) i,knatoms,'add in rmv',ind2(knatoms)
cc write(iout,*) psmm,dmch(ind2(knatoms))
endif
enddo !do while(lll
c project mo onto the selected mo domain and obtain completeness measure
do iiatoms=1,knatoms
ind3(iiatoms)=ind(ind2(iiatoms))
enddo
call bpcompleteness(knatoms,ind3,s,v,natrange,sa,vi,ai,smm,
$r8heap(2*natoms+1),ndao,ovltol,nbasis,i,.true.,r8dummy,sadim)
cc write(iout,*) i,knatoms,'bp in rmv',smm
if (knatoms.eq.nnatoms) exit
if (smm.lt.tol.and.knatoms.lt.nnatoms) then
psmm=psmm-dmch(ind2(knatoms+1))
knatoms=knatoms+1
c write(iout,'("rmv end",3i6,9f10.6)') i,knatoms,nnatoms,smm,psmm
endif
enddo !do while(smm.lt.tol
i8heap(1:natoms)=ind(1:natoms)
nnatoms=knatoms
do ii=1,knatoms
ind(ii)=i8heap(ind2(ii))
enddo
c write(iout,"(2i5,4x,f9.6,4x,1000i4)") i,nnatoms,smm,
c $(ind(iatoms),iatoms=1,nnatoms)
c add those atoms as well that have more than mchtol Mulliken charge
do ii=1,natoms
lll=.false.
mch=dabs(r8heap(natoms+ii))
c write(iout,"('ii, mch',i5,2f10.6)") ii,mch,mchtol
if (mch.ge.mchtol) then
lll=.true. ! atom ii can be added
do jj=1,nnatoms
if (ii.eq.ind(jj)) then
lll=.false. ! atom ii is aready added
c write(iout,"('atom is already added',2i5,f10.6)") ii,nnatoms,mch
exit
endif
enddo ! jj
if (lll) then
nnatoms=nnatoms+1
ind(nnatoms)=ii
cc write(iout,"('add atom ',2i5,f10.6)") ii,nnatoms,mch
endif ! lll
endif ! .ge.mchtol
enddo ! ii
c one final completeness computation
sadim=0
c number of AOs in the BP atomlist
do iiatoms=1,nnatoms
iatoms=ind(iiatoms)
sadim=sadim+natrange(2,iatoms)-natrange(1,iatoms)
enddo
call rpoint2d2(sa,psa,sadim,sadim)
ii=0
do iiatoms=1,nnatoms
iatoms=ind(iiatoms)
do mu=natrange(1,iatoms)+1,natrange(2,iatoms)
ii=ii+1
vi(ii)=v(mu,i)
jj=0
do jjatoms=1,nnatoms
jatoms=ind(jjatoms)
do nu=natrange(1,jatoms)+1,natrange(2,jatoms)
jj=jj+1
psa(ii,jj)=s(mu,nu)
enddo
enddo
enddo
enddo
if (ndao.eq.0) then
call syminvpd(psa,ii,ii) ! exact inverse of the AO overlap
elseif (ndao.gt.0) then ! approximate inverse, subspace of AO overlap eigenvectors with eigenvalues below ovltol are discarded
call genprojmx(r8dummy,r8dummy,ii,ii,r8dummy,ovltol,r8dummy
$ ,psa,.true.,r8heap(2*natoms+1),3)
endif
call dsymv('l',ii,1.d0,psa,ii,vi,1,0.d0,ai,1)
smm=ddot(ii,ai,1,vi,1)
c }}}
endif ! remove
call qsortint(i8heap,nnatoms,ind)
if(mode.ne.'lowd') then
c if(mode.eq.'mull') then
c mnorm=ddot(nbasis,c(1,i),1,v(1,i),1)
smm=(norm-dabs(norm-smm))/norm
c endif
if(i.eq.nocc+1) write(iout,*)
write(iout,"(2i5,4x,f9.6,4x,1000i5)") i,nnatoms,smm,
$(ind(i8heap(iatoms)),iatoms=1,nnatoms)
endif
natmo(i)=nnatoms
do iatoms=1,nnatoms
atmo(iatoms,i)=ind(i8heap(iatoms)) ! add ind(i8heap(iatoms)) to the atom list of MO i
nmoat(ind(i8heap(iatoms)))=nmoat(ind(i8heap(iatoms)))+1
moat(ind(i8heap(iatoms)),nmoat(ind(i8heap(iatoms))))=i ! add i to the MO list assigned to atom ind(i8heap(iatoms))
enddo
if((i.le.nbasis).and.(pmul.eqv..true.)) mocomp(i)=smm/norm !HB
enddo ! i=1,nmatrix
C
return
end
C
************************************************************************
subroutine emftdens(nbasis,natoms,minpfile,natrange,ind,iout,d,da,
$caa,nocc,mo,lmo,lab)
************************************************************************
* Construct density matrix for embedded subsystem
************************************************************************
implicit none
integer nbasis,natoms,minpfile,iout,i,iatoms,jatoms
integer natrange(2,natoms),ind(natoms),nocc
real*8 d(nbasis,nbasis),da(nbasis,nbasis),caa(nbasis,nocc)
real*8 mo(nocc,nbasis)
character*8 embed
logical lmo,lab
C Read atoms of embedded subsystem
open(minpfile,file='MINP')
call embedat(natoms,minpfile,ind,iout,'embed ',5)
close(minpfile)
C Construct densities
da=0.d0
caa=0.d0
if(lab) then
do iatoms=1,natoms
if(lmo.and.ind(iatoms).eq.1) then
do i=1,nocc
caa(natrange(1,iatoms)+1:natrange(2,iatoms),i)=
$ mo(i,natrange(1,iatoms)+1:natrange(2,iatoms))
enddo
endif
do jatoms=1,natoms
if(ind(iatoms).eq.1.or.ind(jatoms).eq.1) then
da(natrange(1,iatoms)+1:natrange(2,iatoms),
$ natrange(1,jatoms)+1:natrange(2,jatoms))=
$ d(natrange(1,iatoms)+1:natrange(2,iatoms),
$ natrange(1,jatoms)+1:natrange(2,jatoms))
endif
enddo
enddo
else
do iatoms=1,natoms
if(ind(iatoms).eq.1) then
if(lmo) then
do i=1,nocc
caa(natrange(1,iatoms)+1:natrange(2,iatoms),i)=
$ mo(i,natrange(1,iatoms)+1:natrange(2,iatoms))
enddo
endif
do jatoms=1,natoms
if(ind(jatoms).eq.1) then
da(natrange(1,iatoms)+1:natrange(2,iatoms),
$ natrange(1,jatoms)+1:natrange(2,jatoms))=
$ d(natrange(1,iatoms)+1:natrange(2,iatoms),
$ natrange(1,jatoms)+1:natrange(2,jatoms))
endif
enddo
endif
enddo
endif
c write(6,"(7f9.5)") da
C
return
end
C
************************************************************************
subroutine fdmdens(nbasis,natoms,natrange,ind,d,da)
************************************************************************
* Construct density matrix for embedded subsystem
************************************************************************
implicit none
integer nbasis,natoms,iatoms,jatoms
integer natrange(2,natoms),ind(natoms)
real*8 d(nbasis,nbasis),da(nbasis,nbasis)
C Construct densities
da=0.d0
c if(.true.) then
if(.false.) then
do iatoms=1,natoms
do jatoms=1,natoms
if(ind(iatoms).eq.1.or.ind(jatoms).eq.1) then
da(natrange(1,iatoms)+1:natrange(2,iatoms),
$ natrange(1,jatoms)+1:natrange(2,jatoms))=
$ d(natrange(1,iatoms)+1:natrange(2,iatoms),
$ natrange(1,jatoms)+1:natrange(2,jatoms))
endif
enddo
enddo
else
do iatoms=1,natoms
if(ind(iatoms).eq.1) then
do jatoms=1,natoms
if(ind(jatoms).eq.1) then
da(natrange(1,iatoms)+1:natrange(2,iatoms),
$ natrange(1,jatoms)+1:natrange(2,jatoms))=
$ d(natrange(1,iatoms)+1:natrange(2,iatoms),
$ natrange(1,jatoms)+1:natrange(2,jatoms))
endif
enddo
endif
enddo
endif
C
return
end
C
************************************************************************
subroutine envdens(nbasis,natoms,natrange,ind,d,da)
************************************************************************
* Construct density matrix for embedded subsystem
************************************************************************
implicit none
integer nbasis,natoms,iatoms,jatoms
integer natrange(2,natoms),ind(natoms)
real*8 d(nbasis,nbasis),da(nbasis,nbasis)
C Construct densities
do iatoms=1,natoms
if(ind(iatoms).eq.0) then
do jatoms=1,natoms
if(ind(jatoms).eq.0) then
da(natrange(1,iatoms)+1:natrange(2,iatoms),
$ natrange(1,jatoms)+1:natrange(2,jatoms))=
$ d(natrange(1,iatoms)+1:natrange(2,iatoms),
$ natrange(1,jatoms)+1:natrange(2,jatoms))
endif
enddo
endif
enddo
C
return
end
C
************************************************************************
subroutine rpoint2d2(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 locdom(nbasis,nocc,natoms,natrange,atmo,nmoat,moat,
$natmo,excrad,atd,ipra,dfnmobasis,moadd,hailen,dfnbasis,
$dfnmobasis_old,linv,ifile,step,ind,iout,scftype,lwdfn,mo,dfatdom,
$dfnatdom,dfatdom_old,dfatind,flag,ialpha,locfit,
$nmobasis,natdom,atdom,naoat,aoat,atind,scftol,lscrohfldf)
C***********************************************************************
C Construct a domain for each occupied LMO and truncate the LMO to it
C flag=2 SCF fitting, flag=3 corr. fitting !MD
C***********************************************************************
implicit none
integer i,j,iatoms,jatoms,iiatoms,jjatoms,nnatoms,jj,iout,n
integer nocc,natoms,natrange(2,natoms,*),nbasis,katoms,hailen
integer atmo(natoms,nocc),nmoat(natoms),moat(natoms,nocc),ialpha
integer natmo(nocc),atd(natoms),step,ind(natoms),locfit
integer moadd(nocc),dfnmobasis_old(nocc),natdommax,state_index
integer dfnmobasis(nocc),dfnbasis,ifile,natdomavg,nmobasis(nocc)
integer naoat(natoms,2),aoat(natoms,nocc),flag,natdom(nocc)
integer dfatdom(natoms,nocc),dfnatdom(nocc),nmodommax,nmodomavg
integer dfatind(nocc,natoms),dfatdom_old(natoms,nocc)
integer ndfatdommax,ndfatdomavg,atind(nocc,natoms)
integer atdom(natoms,nocc),nmodommin,natdommin
real*8 excrad,ipra(natoms,natoms),mo(nbasis,nocc),dfipra(natoms)
real*8 idf1,idf2,scftol
logical log,lll,linv,lwdfn,llo,lscrohfldf
character*5 scftype
C
write(iout,*)
write(iout,*) 'Constructing local fitting domains...'
call getvar('ipra ',ipra)
if(flag.eq.2) then
call getvar('dfipra ',dfipra)
else
call getvar('dfipra_cor',dfipra)
endif
C Construct local fitting domains
hailen=0
ndfatdommax=0
ndfatdomavg=0
do i=1,nocc
dfnatdom(i)=natmo(i)
dfatdom(1:dfnatdom(i),i)=atmo(1:dfnatdom(i),i)
do iiatoms=1,natmo(i)
iatoms=atmo(iiatoms,i)
call locdom1(i,iatoms,dfnatdom(i),dfatdom(1,i),dfipra,natoms,
$ipra,excrad,nmoat,moat)
enddo
call qsortint(ind,dfnatdom(i),dfatdom(1,i))
do iatoms=1,dfnatdom(i)
atd(iatoms)=dfatdom(ind(iatoms),i)
enddo
dfatdom(1:dfnatdom(i),i)=atd(1:dfnatdom(i))
dfatdom(dfnatdom(i)+1:natoms,i)=0
ndfatdommax=max(ndfatdommax,dfnatdom(i))
ndfatdomavg=ndfatdomavg+dfnatdom(i)
jj=1
do iiatoms=1,dfnatdom(i)
iatoms=dfatdom(iiatoms,i)
dfatind(i,iatoms)=jj
jj=jj+natrange(2,iatoms,flag)-natrange(1,iatoms,flag)
enddo
jj=jj-1
dfnmobasis(i)=jj
moadd(i)=hailen
hailen=hailen+jj*nbasis
enddo
c do i=1,nocc
c write(*,'(I3,A,500I4)') i,':',dfatdom(1:dfnatdom(i),i)
c enddo
C Construct local AO domains
if(locfit.eq.2) then
llo=.true.
hailen=0
c naoat=0
c natdom=0
natdommax=0
natdommin=natoms
natdomavg=0
nmodommax=0
nmodommin=natoms
nmodomavg=0
do i=1,nocc
c do iatoms=1,natoms
c if(0.5d0*sum(mo(natrange(1,iatoms,1)+1:
c $ natrange(2,iatoms,1),i)**2).gt.min(1.d-6,scftol)) then
c naoat(iatoms,1)=naoat(iatoms,1)+1
c aoat(iatoms,naoat(iatoms,1))=i
c natdom(i)=natdom(i)+1
c atdom(natdom(i),i)=iatoms
c endif
c enddo
nmodommax=max(nmodommax,natdom(i))
nmodommin=min(nmodommin,natdom(i))
nmodomavg=nmodomavg+natdom(i)
llo=llo.and.natdom(i).eq.natoms
enddo
naoat(1:natoms,2)=naoat(1:natoms,1)
do i=1,nocc
n=natdom(i)
do iiatoms=1,n
iatoms=atdom(iiatoms,i)
call locdom1(i,iatoms,natdom(i),atdom(1,i),dfipra,natoms,
$ipra,min(1.d-2,scftol*1d4),naoat,aoat)
enddo
call qsortint(ind,natdom(i),atdom(1,i))
do iatoms=1,natdom(i)
atd(iatoms)=atdom(ind(iatoms),i)
enddo
atdom(1:natdom(i),i)=atd(1:natdom(i))
natdommax=max(natdommax,natdom(i))
natdommin=min(natdommin,natdom(i))
natdomavg=natdomavg+natdom(i)
jj=1
do iiatoms=1,natdom(i)
iatoms=atdom(iiatoms,i)
atind(i,iatoms)=jj-2
jj=jj+natrange(2,iatoms,1)-natrange(1,iatoms,1)
enddo
jj=jj-1
nmobasis(i)=jj
llo=llo.and.nmobasis(i).eq.nbasis
moadd(i)=hailen
hailen=hailen+nmobasis(i)*dfnmobasis(i)
enddo
if(llo) locfit=0
endif
C Check if domains have been changed
if(ialpha.le.1.or.(lscrohfldf.and.ialpha.gt.1)) then
open(ifile,file='DFLTOC',form='unformatted')
else
open(ifile,file='DFLTOCb',form='unformatted')
endif
linv=.true.
lll=.false.
if(lwdfn.and.((step.gt.2.and.flag.eq.2).or.
$(step.ge.2.and.flag.eq.3))) then
rewind(ifile)
read(ifile) dfnmobasis_old
linv=.false.
do i=1,nocc
if(dfnmobasis(i).ne.dfnmobasis_old(i)) then
linv=.true.
exit
endif
enddo
c if(step.gt.1.and.linv) then
c do i=1,nocc
c if(dfnmobasis(i).lt.dfnmobasis_old(i)) then
c lll=.true.
c exit
c endif
c enddo
c endif
if(.not.linv) then
read(ifile) dfatdom_old
do i=1,nocc
do iatoms=1,dfnatdom(i)
if(dfatdom(iatoms,i).ne.dfatdom_old(iatoms,i)) then
linv=.true.
exit
endif
enddo
enddo
endif
endif
C Restore old domains
c if(lll.and.lwdfn) then
c rewind(ifile)
c read(ifile) dfnmobasis
c read(ifile) moadd
c linv=.false.
c endif
C Print results
if(linv) then
lwdfn=.true.
rewind(ifile)
write(ifile) dfnmobasis
write(ifile) dfatdom,dfnatdom
write(ifile) dfatind
write(iout,
$"(' Maximum number of atoms in fitting domains:',i5)") ndfatdommax
write(iout,
$"(' Average number of atoms in fitting domains:',f7.1)")
$dble(ndfatdomavg)/dble(nocc)
else
write(iout,*) 'Local fitting domains are unchanged.'
endif
if(locfit.eq.2) then
c write(iout,
c $"(' Mimimum number of atoms in MO domains: ',i5)") nmodommin
write(iout,
$"(' Maximum number of atoms in MO domains: ',i5)") nmodommax
write(iout,
$"(' Average number of atoms in MO domains: ',f7.1)")
$dble(nmodomavg)/dble(nocc)
c write(iout,
c $"(' Mimimum number of atoms in AO domains: ',i5)") natdommin
write(iout,
$"(' Maximum number of atoms in AO domains: ',i5)") natdommax
write(iout,
$"(' Average number of atoms in AO domains: ',f7.1)")
$dble(natdomavg)/dble(nocc)
endif
close(ifile)
c copy DFLTOC for direct_fock_build
if(lscrohfldf.and.ialpha.gt.1) then
open(ifile,file='DFLTOCb',form='unformatted')
rewind(ifile)
write(ifile) dfnmobasis
write(ifile) dfatdom,dfnatdom
write(ifile) dfatind
close(ifile)
endif
C
return
end
C
************************************************************************
subroutine locdom1(i,iatoms,natdom,atdom,dfipra,natoms,ipra,
$excrad,nmoat,moat)
************************************************************************
* Construct a domain for an occupied LMO
************************************************************************
implicit none
integer i,j,natoms,natdom,atdom(natoms)
integer iatoms,jatoms,katoms,nmoat(natoms),moat(natoms,*)
real*8 excrad,ipra(natoms,natoms),dfipra(natoms)
real*8 idf1,idf2
logical log,lll
C
idf1=dfipra(iatoms)
do jatoms=1,natoms
idf2=max(idf1,dfipra(jatoms))
if(ipra(jatoms,iatoms)*idf2.ge.excrad) then
log=.true.
do katoms=1,natdom
if(jatoms.eq.atdom(katoms)) then
log=.false.
exit
endif
enddo
if(log) then
natdom=natdom+1
atdom(natdom)=jatoms
lll=.true.
do j=1,nmoat(jatoms)
if(i.eq.moat(jatoms,j)) then
lll=.false.
exit
endif
enddo
if(lll) then
nmoat(jatoms)=nmoat(jatoms)+1
moat(jatoms,nmoat(jatoms))=i
endif
endif
endif
enddo
C
return
end
C
************************************************************************
subroutine bpcompleteness(natoms,atoms,s,v,natrange,sa,vi,ai,bpc,
$ scr,ndao,ovltol,nbasis,imo,project,cmo,sadim)
************************************************************************
* project mo onto the AO basis of the selected mo domain in list atoms and obtain completeness measure, bpc for MO i
* see J. Comp. Chem. 14, 736 (1993)
* if !project then calculate only the completeness
************************************************************************
implicit none
integer imo,natoms,atoms(natoms),mu,nbasis,nu,ndao
integer sadim,iiatoms,natrange(2,*),jjatoms,iatoms,jatoms,jj,ii
real*8 s(nbasis,nbasis),sa(nbasis,nbasis),vi(nbasis),ai(nbasis)
real*8 scr(*),r8dummy,ovltol,bpc,ddot,v(nbasis,*),cmo(nbasis,*)
real*8,pointer :: psa(:,:)
logical project
C
interface
c {{{ interfaces for pointers
subroutine rpoint2d2(egydim,haromdim,dim1,dim2)
implicit none
integer dim1,dim2
real*8,target :: egydim(dim1,dim2)
real*8, pointer :: haromdim(:,:)
end subroutine
c}}}
end interface
c
sadim=0
do iiatoms=1,natoms
iatoms=atoms(iiatoms)
sadim=sadim+natrange(2,iatoms)-natrange(1,iatoms) ! number of AOs in the BP atomlist
enddo
call rpoint2d2(sa,psa,sadim,sadim)
ii=0
do iiatoms=1,natoms
iatoms=atoms(iiatoms)
do mu=natrange(1,iatoms)+1,natrange(2,iatoms)
ii=ii+1
vi(ii)=v(mu,imo)
if (.not.project) ai(ii)=cmo(mu,imo)
jj=0
do jjatoms=1,natoms
jatoms=atoms(jjatoms)
do nu=natrange(1,jatoms)+1,natrange(2,jatoms)
jj=jj+1
psa(ii,jj)=s(mu,nu)
enddo
enddo
enddo
enddo
if (project) then
if (ndao.eq.0) then
call syminvpd(psa,ii,ii) ! exact inverse of the AO overlap
elseif (ndao.gt.0) then ! approximate inverse, subspace of AO overlap eigenvectors with eigenvalues below ovltol are discarded
call genprojmx(r8dummy,r8dummy,ii,ii,r8dummy,ovltol,
$ r8dummy,psa,.true.,scr,3)
endif
call dsymv('l',ii,1.d0,psa,ii,vi,1,0.d0,ai,1) ! projected mo coefficients are in ai
endif
bpc=ddot(ii,ai,1,vi,1)
C
return
end