mirror of
https://code.it4i.cz/sccs/easyconfigs-it4i.git
synced 2025-04-16 19:50:50 +01:00
1445 lines
46 KiB
Fortran
Executable File
1445 lines
46 KiB
Fortran
Executable File
program mulli
|
|
C
|
|
implicit none
|
|
!#include "MRCCCOMMON"
|
|
C
|
|
real*8 m,mo,mv
|
|
real*8 itol,bpcompo,bpcompv
|
|
character*2 ch2
|
|
character*4 citol,loctypeo,loctypev,localcc
|
|
character*5 boys,gboys(2)
|
|
character*16 orblocv,c16
|
|
character*26 cscr,cscr2,cscr3
|
|
logical locvirt
|
|
integer i,ic,is,iv,intalloc,dblalloc,iind,icoord,nbset
|
|
integer idistc,idisti,ifock,iatmo,ivi,isa,iai,icmo,incmo,ibm
|
|
integer inatmo,imoat,inmoat,iatdom,imodom,inatdom,inmodom,isi,isb
|
|
integer ifa,isr,iatnum,incorb,icold,ismm,ibmm,isbb,ieval,isyevnek
|
|
integer rmemory,dfloatl,nbasis,moldenfile
|
|
integer*8 memrr
|
|
integer scrfile1,scrfile2,scrfile3,oeintfile,nvirt,mocoeffile
|
|
integer natoms,j,k,ii,jj,nn,mu,nu,nnatoms,iiatoms,jjatoms,ndom
|
|
integer iatoms,jatoms,katoms,idom,imo,nc,no,nv,kk,nelec,nnn,ncent
|
|
integer ncore,nocc,i1,i2,i3,iout,l,kkatoms,jdom,kdom,nlarge
|
|
integer nangmax,ncontrmax,nprimmax,nshmax
|
|
|
|
integer*4 err
|
|
real*8 sum,sum1,ddot,nuc,domrad,cl
|
|
logical log,llog
|
|
logical, allocatable :: logv(:)
|
|
character*2, allocatable :: atsymbol(:)
|
|
character*20, allocatable :: basname(:,:),ecpname(:)
|
|
parameter(dfloatl=8)
|
|
real*8,allocatable, target :: RLONG(:)
|
|
integer, pointer :: natrange(:,:),ind(:),atmo(:,:)
|
|
integer, pointer :: natmo(:),moat(:,:),ncorecp(:)
|
|
integer, pointer :: nmoat(:),atdom(:,:),modom(:,:),natdom(:)
|
|
integer, pointer :: nmodom(:),cmo(:,:),ncmo(:),atnum(:),ncorb(:)
|
|
integer, pointer :: perm2(:),perm(:),perm3(:),i8heap(:),atnumt(:)
|
|
integer*4, pointer :: i4heap(:)
|
|
real*8, pointer :: c(:,:),s(:,:),v(:,:),coord(:,:),distc(:,:)
|
|
real*8, pointer :: disti(:,:),atchg(:)
|
|
real*8, pointer :: fock(:,:),sa(:,:),vi(:),ai(:),bm(:),si(:,:)
|
|
real*8, pointer :: sb(:),r8heap(:)
|
|
real*8, pointer :: cold(:,:),smm(:,:),bmm(:,:),sbb(:,:),eval(:)
|
|
real*8, pointer :: syevnek(:)
|
|
real*8, pointer :: fa(:,:),bmmm(:,:),pao(:,:),pao_trunk(:,:)
|
|
real*8, pointer :: density(:),lbb(:,:),lmm(:,:),lmmm(:,:)
|
|
character*16 lline1
|
|
character*1 lline2(16)
|
|
equivalence(lline1,lline2)
|
|
interface
|
|
subroutine PROJECT(c,p,s,no,project_c,mo_ini,mo_end,scr,nb)
|
|
implicit none
|
|
integer i,j,k,no,nb,mo_ini,mo_end
|
|
real*8 c(nb,nb),s(nb,nb),p(nb,nb),scr(nb,nb),project_c(nb,nb)
|
|
end subroutine
|
|
subroutine DENSM(p,c,no,mo_ini,mo_end,nb,nlarge)
|
|
implicit none
|
|
integer i,j,k,nb,no,mo_ini,mo_end,nlarge
|
|
real*8 p(nb,nb),c(nb,nb)
|
|
end subroutine
|
|
subroutine MO_OVERL(c,s,ov,no,mo_ini,mo_end,scr,nb,nlarge)
|
|
implicit none
|
|
integer i,j,k,no,nb,mo_ini,mo_end,nlarge
|
|
real*8 c(nb,nb),s(nb,nb),ov(nb,nb),scr(nb,nb)
|
|
end subroutine
|
|
subroutine rpoint1d(egydim1,egydim2,dim1,memr)
|
|
implicit none
|
|
integer dim1
|
|
integer*8 memr
|
|
real*8,target :: egydim1(dim1)
|
|
real*8, pointer :: egydim2(:)
|
|
end subroutine
|
|
subroutine ipoint1d(egydim1,egydim2,dim1,memr)
|
|
implicit none
|
|
integer dim1
|
|
integer*8 memr
|
|
real*8,target :: egydim1(dim1)
|
|
integer, pointer :: egydim2(:)
|
|
end subroutine
|
|
subroutine ipoint1d_half(egydim1,egydim2,dim1,memr)
|
|
implicit none
|
|
integer dim1
|
|
integer*8 memr
|
|
real*8,target :: egydim1(dim1)
|
|
integer*4, pointer :: egydim2(:)
|
|
end subroutine
|
|
subroutine ipoint2d(egydim,ketdim,dim1,dim2,memr)
|
|
implicit none
|
|
integer dim1,dim2
|
|
integer*8 memr
|
|
real*8,target :: egydim(dim1,dim2)
|
|
integer, pointer :: ketdim(:,:)
|
|
end subroutine
|
|
subroutine rpoint2d(egydim,ketdim,dim1,dim2,memr)
|
|
implicit none
|
|
integer dim1,dim2
|
|
integer*8 memr
|
|
real*8,target :: egydim(dim1,dim2)
|
|
real*8, pointer :: ketdim(:,:)
|
|
end subroutine
|
|
end interface
|
|
C
|
|
call mrccini
|
|
C Open files
|
|
C
|
|
C Read fort.55
|
|
rewind(55)
|
|
read(55,*) nbasis,nelec,ncore
|
|
nocc=nelec/2-ncore
|
|
nvirt=nbasis-nocc-ncore
|
|
call getvar('natoms ',natoms)
|
|
call getvar('ncent ',ncent)
|
|
call getvar('nbset ',nbset)
|
|
|
|
c rewind(777)
|
|
c read(777,*)ncore,nocc,nbasis,natoms
|
|
c write(6,*) ncore,nocc,nbasis,natoms
|
|
c call flush(6)
|
|
nlarge=ncore+nocc+nbasis
|
|
oeintfile=91
|
|
mocoeffile=92
|
|
moldenfile=93
|
|
scrfile1=25
|
|
scrfile2=26
|
|
scrfile3=27
|
|
open(oeintfile,file='OEINT',form='UNFORMATTED')
|
|
open(mocoeffile,file='MOCOEF',form='UNFORMATTED')
|
|
open(moldenfile,file='MOLDEN')
|
|
C
|
|
C Read fort.55
|
|
c rewind(inp)
|
|
c read(inp,*) nbasis,nelec,ncore
|
|
c nocc=nelec/2-ncore
|
|
nvirt=nbasis-nocc-ncore
|
|
call getkey('mem',3,lline1,16)
|
|
if(lline1.ne.' ') then
|
|
i=1
|
|
do while(lline2(i).ne.' '.and.i.lt.16)
|
|
i=i+1
|
|
enddo
|
|
call lowercase(lline2(i-2),ch2,2)
|
|
if(ch2.ne.'mb'.and.ch2.ne.'gb') then
|
|
write(iout,*) 'Unknown memory unit!'
|
|
call mrccend(1)
|
|
endif
|
|
lline2(i-2)=' '
|
|
lline2(i-1)=' '
|
|
read(lline1,*) rmemory
|
|
else
|
|
ch2='mb'
|
|
rmemory=256
|
|
endif
|
|
if(ch2.eq.'gb') rmemory=1024*rmemory
|
|
|
|
rmemory=rmemory*1024*1024/dfloatl
|
|
C
|
|
c write(6,*) 'alloc elott ',rmemory,dfloatl
|
|
c call flush(6)
|
|
C Allocate memory
|
|
c call memalloc
|
|
allocate(RLONG(1:rmemory))
|
|
memrr=0
|
|
C
|
|
C Set tolerance of zero elements of MO coefficients
|
|
call getkey('scftol',6,citol,4)
|
|
read(citol,*) i
|
|
itol=dsqrt(10.d0**(-i))
|
|
C
|
|
call getkey('orblocv',7,orblocv,16)
|
|
C
|
|
C%%% Memoriafoglalas!
|
|
call getkey('localcc',7,localcc,4)
|
|
cDDDDDDDDDDDDDDDDDDDDDDDDDDDDDDDD
|
|
IF(localcc.ne.'off ')THEN
|
|
cDDDDDDDDDDDDDDDDDDDDDDDDDDDDDDDD
|
|
call getvar('natoms ',natoms)
|
|
allocate(logv(nbasis),atsymbol(natoms),basname(natoms,nbset))
|
|
allocate(ecpname(natoms),atchg(ncent),ncorecp(natoms))
|
|
call ipoint1d(RLONG(memrr+1),perm,nbasis,memrr)
|
|
call ipoint1d(RLONG(memrr+1),perm2,nbasis,memrr)
|
|
call ipoint1d(RLONG(memrr+1),perm3,nbasis,memrr)
|
|
call ipoint1d(RLONG(memrr+1),natmo,nbasis,memrr)
|
|
call ipoint1d(RLONG(memrr+1),nmoat,natoms,memrr)
|
|
call ipoint1d(RLONG(memrr+1),natdom,nocc,memrr)
|
|
call ipoint1d(RLONG(memrr+1),nmodom,nocc,memrr)
|
|
call ipoint1d(RLONG(memrr+1),atnum,natoms,memrr)
|
|
call ipoint1d(RLONG(memrr+1),atnumt,natoms,memrr)
|
|
call ipoint1d(RLONG(memrr+1),ncmo,nocc,memrr)
|
|
call ipoint1d(RLONG(memrr+1),ncorb,natoms,memrr)
|
|
call ipoint1d(RLONG(memrr+1),ind,nbasis,memrr)
|
|
|
|
call ipoint2d(RLONG(memrr+1),natrange,2,natoms*nbset,memrr)
|
|
call ipoint2d(RLONG(memrr+1),atmo,natoms,nbasis,memrr)
|
|
call ipoint2d(RLONG(memrr+1),moat,natoms,nbasis,memrr)
|
|
call ipoint2d(RLONG(memrr+1),atdom,natoms,nocc,memrr)
|
|
call ipoint2d(RLONG(memrr+1),modom,nbasis,nocc,memrr)
|
|
call ipoint2d(RLONG(memrr+1),cmo,nocc,nocc,memrr)
|
|
|
|
call ipoint1d(RLONG(memrr+1),i8heap,2*nbasis,memrr)
|
|
|
|
i=memrr
|
|
call rpoint1d(RLONG(memrr+1),r8heap,2*nbasis**2,memrr)
|
|
memrr=i
|
|
call ipoint1d_half(RLONG(memrr+1),I4HEAP,4*nbasis**2,memrr)
|
|
|
|
call rpoint1d(RLONG(memrr+1),syevnek,nlarge**2,memrr)
|
|
call rpoint1d(RLONG(memrr+1),bm,nbasis**2,memrr)
|
|
call rpoint1d(RLONG(memrr+1),vi,nbasis,memrr)
|
|
call rpoint1d(RLONG(memrr+1),ai,nlarge,memrr)
|
|
call rpoint1d(RLONG(memrr+1),sb,nbasis**2,memrr)
|
|
call rpoint1d(RLONG(memrr+1),eval,nbasis,memrr)
|
|
|
|
call rpoint2d(RLONG(memrr+1),c,nbasis,nbasis,memrr)
|
|
call rpoint2d(RLONG(memrr+1),s,nbasis,nbasis,memrr)
|
|
call rpoint2d(RLONG(memrr+1),v,nbasis,nbasis,memrr)
|
|
call rpoint2d(RLONG(memrr+1),coord,3,ncent,memrr)
|
|
call rpoint2d(RLONG(memrr+1),distc,natoms,natoms,memrr)
|
|
call rpoint2d(RLONG(memrr+1),disti,natoms,natoms,memrr)
|
|
call rpoint2d(RLONG(memrr+1),fock,nbasis,nbasis,memrr)
|
|
call rpoint2d(RLONG(memrr+1),sa,nbasis,nbasis,memrr)
|
|
call rpoint2d(RLONG(memrr+1),si,nbasis,nbasis,memrr)
|
|
call rpoint2d(RLONG(memrr+1),fa,nbasis,nbasis,memrr)
|
|
call rpoint2d(RLONG(memrr+1),cold,nbasis,2*nbasis,memrr)
|
|
call rpoint2d(RLONG(memrr+1),smm,nbasis,nbasis,memrr)
|
|
call rpoint2d(RLONG(memrr+1),bmm,nbasis,nbasis,memrr)
|
|
call rpoint2d(RLONG(memrr+1),lmm,nbasis,nbasis,memrr)
|
|
call rpoint2d(RLONG(memrr+1),sbb,nbasis,nbasis,memrr)
|
|
call rpoint2d(RLONG(memrr+1),lbb,nlarge,nlarge,memrr)
|
|
call rpoint2d(RLONG(memrr+1),bmmm,nbasis,nbasis,memrr)
|
|
call rpoint2d(RLONG(memrr+1),lmmm,nbasis,nbasis,memrr)
|
|
call rpoint2d(RLONG(memrr+1),pao,nbasis,nbasis,memrr)
|
|
call rpoint2d(RLONG(memrr+1),pao_trunk,nbasis,nbasis,memrr)
|
|
|
|
c write(6,*) 'rmemory,memrr',rmemory,memrr
|
|
C***********************************************************************
|
|
C Calculate Mulliken populations and define fragments
|
|
C***********************************************************************
|
|
|
|
iout=6
|
|
C Read MO coefficients from mocoeffile
|
|
call dfillzero(c,nbasis**2)
|
|
call dfillzero(v,nbasis**2)
|
|
rewind(mocoeffile)
|
|
read(mocoeffile) nn,(r8heap(jj),jj=1,2*nn)
|
|
do ii=1,nn
|
|
i=i4heap((ii-1)*4+3)
|
|
j=i4heap((ii-1)*4+4)
|
|
c(i,j)=r8heap(ii*2-1)
|
|
c write(6,*)i,j,c(i,j),' i,j,c(i,j) '
|
|
enddo
|
|
call dcopy(nbasis**2,c,1,cold,1)
|
|
C Read overlap integrals
|
|
call dfillzero(s,nbasis**2)
|
|
rewind(oeintfile)
|
|
read(oeintfile) nn,(r8heap(jj),jj=1,2*nn)
|
|
do ii=1,nn
|
|
i=i4heap((ii-1)*4+3)
|
|
j=i4heap((ii-1)*4+4)
|
|
s(i,j)=r8heap(ii*2-1)
|
|
s(j,i)=r8heap(ii*2-1)
|
|
enddo
|
|
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
|
|
do i=1,nbasis
|
|
perm2(i)=i
|
|
perm3(i)=i
|
|
enddo
|
|
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
|
|
|
|
C Read AO Fock-matrix
|
|
call dfillzero(fock,nbasis**2)
|
|
open(scrfile1,file='FOCK',form='UNFORMATTED')
|
|
rewind(scrfile1)
|
|
read(scrfile1) nn,(r8heap(jj),jj=1,2*nn)
|
|
do ii=1,nn
|
|
i=i4heap((ii-1)*4+3)
|
|
j=i4heap((ii-1)*4+4)
|
|
fock(i,j)=r8heap(ii*2-1)
|
|
fock(j,i)=r8heap(ii*2-1)
|
|
enddo
|
|
close(scrfile1)
|
|
call dgemm('t','n',nbasis,nbasis,nbasis,1.d0,c,nbasis,fock,nbasis,
|
|
$0.d0,r8heap(1),nbasis)
|
|
call dgemm('n','n',nbasis,nbasis,nbasis,1.d0,r8heap(1),
|
|
$nbasis,c,nbasis,0.d0,bmm,nbasis)
|
|
do i=1,nbasis
|
|
bm(i)=bmm(i,i)
|
|
c write(6,*) i,bmm(i,i),' fock ok'
|
|
perm(i)=i
|
|
enddo
|
|
do k=1,nocc+ncore-1
|
|
do i=1,nocc+ncore-k
|
|
if(bm(i+1).lt.bm(i))then
|
|
j=perm(i+1)
|
|
perm(i+1)=perm(i)
|
|
perm(i)=j
|
|
cl=bm(i+1)
|
|
bm(i+1)=bm(i)
|
|
bm(i)=cl
|
|
endif
|
|
enddo
|
|
enddo
|
|
do k=1,nbasis
|
|
c write(6,*) k,perm(k),' perm '
|
|
do i=1,nbasis
|
|
c cold(1:nbasis,k)=c(1:nbasis,perm(k))
|
|
cold(i,k)=c(i,perm(k))
|
|
enddo
|
|
enddo
|
|
do k=1,nbasis
|
|
do i=1,nbasis
|
|
c(i,k)=cold(i,k)
|
|
enddo
|
|
enddo
|
|
call dgemm('t','n',nbasis,nbasis,nbasis,1.d0,c,nbasis,fock,nbasis,
|
|
$0.d0,r8heap(1),nbasis)
|
|
call dgemm('n','n',nbasis,nbasis,nbasis,1.d0,r8heap(1),
|
|
$nbasis,c,nbasis,0.d0,bmm,nbasis)
|
|
do i=1,nbasis
|
|
bm(i)=bmm(i,i)
|
|
c write(6,*) i,bmm(i,i),' fock uj'
|
|
perm(i)=i
|
|
enddo
|
|
c call mrccend(1)
|
|
C Read variables from VARS file
|
|
call getvar('natrange ',natrange)
|
|
call getvar('atchg ',atchg)
|
|
call getvar('atnum ',atnum)
|
|
call getvar('atnumt ',atnumt)
|
|
call getvar('ncorecp ',ncorecp)
|
|
call getvar_c('basname ',basname)
|
|
call getvar_c('ecpname ',ecpname)
|
|
call getvar_c('atsymbol ',atsymbol)
|
|
call getvar('nuc ',nuc)
|
|
call getvar('ncorb ',ncorb)
|
|
call getvar('distc ',distc)
|
|
call getvar('coord ',coord)
|
|
call getvar('nangmax ',nangmax)
|
|
call getvar('ncontrmax ',ncontrmax)
|
|
call getvar('nprimmax ',nprimmax)
|
|
call getvar('nshmax ',nshmax)
|
|
c write(6,'("num of cores",70I3)') ncorb(1:natoms)
|
|
call getkey('domrad',6,c16,16)
|
|
read(c16,*) domrad
|
|
c write(6,*) ' domrad =',domrad
|
|
C Boughton-Pulay algorithm
|
|
C V = S * C
|
|
call dsymm('l','l',nbasis,nbasis,1.d0,s,nbasis,c,nbasis,
|
|
$0.d0,v,nbasis)
|
|
C Calculate Mulliken populations and assign MOs to atoms
|
|
do iatoms=1,natoms
|
|
c write(6,*)natrange(1,iatoms)+1,natrange(2,iatoms),'atran ',iatoms
|
|
do mu=natrange(1,iatoms)+1,natrange(2,iatoms)
|
|
eval(mu)=dfloat(iatoms)+dfloat(mu)/1000.d0
|
|
enddo
|
|
enddo
|
|
cscr='NON'
|
|
c call molden_printPAO(eval,pao,nbasis,nbasis,nbasis,nbasis,
|
|
c $ moldenfile,perm2,cscr)
|
|
cscr='AO'
|
|
smm=0.d0
|
|
do i=1,nbasis
|
|
smm(i,i)=1.d0
|
|
enddo
|
|
c call molden_printPAO(eval,smm,nbasis,nbasis,nbasis,nbasis,
|
|
c $ moldenfile,perm2,cscr)
|
|
i=ncore+nocc
|
|
c i=nbasis
|
|
if(orblocv.ne.'pao ') i=nbasis
|
|
write(iout,*) 'Executing the Boughton-Pulay algorithm...'
|
|
call getkey('bpcompo',7,c16,16)
|
|
read(c16,*) bpcompo
|
|
call getkey('bpcompv',7,c16,16)
|
|
read(c16,*) bpcompv
|
|
call bopu(nbasis,i,natoms,natrange,c,v,ind,r8heap,s,sa,vi,ai,
|
|
$atmo,nmoat,moat,natmo,bpcompo,bpcompv,ncore+nocc,natoms,'bopu',
|
|
$i8heap,natrange,nbasis,c,r8heap,0,0.d0,0.0d0,.false.,iout,0.05d0)
|
|
call timer
|
|
write(iout,*)
|
|
C Construct atomic domains for occupied LMOs
|
|
C Loop over occupied orbitals excluding core orbitals
|
|
write(iout,*) 'Constructing extended domains...'
|
|
do i=1,nocc
|
|
natdom(i)=natmo(ncore+i)
|
|
nmodom(i)=0
|
|
do iiatoms=1,natmo(ncore+i)
|
|
atdom(iiatoms,i)=atmo(iiatoms,ncore+i)
|
|
enddo
|
|
C Loop over atoms on which occupied orbital i is localzed
|
|
do iiatoms=1,natmo(ncore+i)
|
|
iatoms=atmo(iiatoms,ncore+i)
|
|
do jatoms=1,natoms
|
|
C Include all atoms in the domain whose distance from atom iiatoms is
|
|
C smaller than a threshold
|
|
c if(distc(iatoms,jatoms).lt.8.d0) then
|
|
if(distc(iatoms,jatoms).lt.domrad) then
|
|
log=.true.
|
|
do katoms=1,natdom(i)
|
|
if(jatoms.eq.atdom(katoms,i)) then
|
|
log=.false.
|
|
exit
|
|
endif
|
|
enddo
|
|
if(log) then
|
|
natdom(i)=natdom(i)+1
|
|
atdom(natdom(i),i)=jatoms
|
|
endif
|
|
endif
|
|
enddo
|
|
enddo
|
|
nn=natdom(i)
|
|
C Loop over atoms which are currently included in the domain
|
|
do iiatoms=1,nn
|
|
iatoms=atdom(iiatoms,i)
|
|
do jj=1,nmoat(iatoms)
|
|
j=moat(iatoms,jj)
|
|
C Extend the domain: for each orbital j localized on iiatoms the atoms
|
|
C on which orbital j is localized is also added to the domain
|
|
c if(j.ne.ncore+i) then !all orbitals
|
|
c if(j.ne.ncore+i.and.j.le.ncore+nocc) then !occupied only
|
|
c if(j.ne.ncore+i.and.j.gt.ncore.and.
|
|
if(j.ne.ncore+i.and.
|
|
$(j.le.ncore+nocc.or.orblocv.ne.'pao ')) then
|
|
if(j.gt.ncore)then
|
|
do jjatoms=1,natmo(j)
|
|
jatoms=atmo(jjatoms,j)
|
|
log=.true.
|
|
do katoms=1,natdom(i)
|
|
if(jatoms.eq.atdom(katoms,i)) then
|
|
log=.false.
|
|
exit
|
|
endif
|
|
enddo
|
|
if(log) then
|
|
natdom(i)=natdom(i)+1
|
|
atdom(natdom(i),i)=jatoms
|
|
endif
|
|
enddo
|
|
endif
|
|
log=.true.
|
|
do k=1,nmodom(i)
|
|
if(j.eq.modom(k,i)) then
|
|
log=.false.
|
|
exit
|
|
endif
|
|
enddo
|
|
if(log) then
|
|
nmodom(i)=nmodom(i)+1
|
|
modom(nmodom(i),i)=j
|
|
endif
|
|
endif
|
|
enddo
|
|
enddo
|
|
nmodom(i)=nmodom(i)+1
|
|
modom(nmodom(i),i)=ncore+i
|
|
C Add core orbitals localized on the atoms of the domain to the domain
|
|
do iiatoms=1,natdom(i)
|
|
iatoms=atdom(iiatoms,i)
|
|
do jj=1,nmoat(iatoms)
|
|
j=moat(iatoms,jj)
|
|
if(j.le.ncore)then
|
|
log=.true.
|
|
do k=1,nmodom(i)
|
|
if(j.eq.modom(k,i)) then
|
|
log=.false.
|
|
exit
|
|
endif
|
|
enddo
|
|
if(log) then
|
|
nmodom(i)=nmodom(i)+1
|
|
modom(nmodom(i),i)=j
|
|
endif
|
|
endif
|
|
enddo
|
|
enddo
|
|
C Now we check if particular domains are subsets of other domains,
|
|
C if so, then the smaller domains are dropped.
|
|
C First, we reorder atom and MO serial numbers into ascending order
|
|
call qsortint(i8heap,natdom(i),atdom(1,i))
|
|
do j=1,natdom(i)
|
|
i8heap(natdom(i)+j)=atdom(j,i)
|
|
enddo
|
|
do j=1,natdom(i)
|
|
atdom(j,i)=i8heap(natdom(i)+i8heap(j))
|
|
enddo
|
|
call qsortint(i8heap,nmodom(i),modom(1,i))
|
|
do j=1,nmodom(i)
|
|
i8heap(nmodom(i)+j)=modom(j,i)
|
|
enddo
|
|
do j=1,nmodom(i)
|
|
modom(j,i)=i8heap(nmodom(i)+i8heap(j))
|
|
enddo
|
|
ncmo(i)=1
|
|
cmo(i,1)=ncore+i
|
|
ind(i)=i
|
|
c write(6,"(a2,100i3)") 'MO',(cmo(i,j),j=1,ncmo(i))
|
|
c write(6,"('MOs in domain :',1000i4)")(modom(j,i),j=1,nmodom(i))
|
|
c write(6,"('Atoms in domain:',1000i4)")(atdom(j,i),j=1,natdom(i))
|
|
enddo
|
|
C Second, we order domains according to the number of atoms and MOs in
|
|
C the domain
|
|
ndom=nocc
|
|
log=.true.
|
|
do while(log)
|
|
log=.false.
|
|
do i=1,ndom-1
|
|
if(natdom(i).lt.natdom(i+1).or.(natdom(i).eq.natdom(i+1).and.
|
|
$nmodom(i).lt.nmodom(i+1))) then
|
|
log=.true.
|
|
do j=1,natdom(i+1)
|
|
k=atdom(j,i)
|
|
atdom(j,i)=atdom(j,i+1)
|
|
atdom(j,i+1)=k
|
|
enddo
|
|
do j=1,max(nmodom(i),nmodom(i+1))
|
|
k=modom(j,i)
|
|
modom(j,i)=modom(j,i+1)
|
|
modom(j,i+1)=k
|
|
enddo
|
|
k=natdom(i)
|
|
natdom(i)=natdom(i+1)
|
|
natdom(i+1)=k
|
|
k=nmodom(i)
|
|
nmodom(i)=nmodom(i+1)
|
|
nmodom(i+1)=k
|
|
k=cmo(i,1)
|
|
cmo(i,1)=cmo(i+1,1)
|
|
cmo(i+1,1)=k
|
|
endif
|
|
enddo
|
|
enddo
|
|
c write(6,*)
|
|
c write(6,*) 'Domains are ordered according to the num of atoms '
|
|
c do i=1,ndom
|
|
c write(6,"(a2,100i3)") 'MO',(cmo(i,j),j=1,ncmo(i))
|
|
c write(6,"('MOs in domain :',1000i4)")(modom(j,i),j=1,nmodom(i))
|
|
c write(6,"('Atoms in domain:',1000i4)")(atdom(j,i),j=1,natdom(i))
|
|
c enddo
|
|
C Merge domains: drop a domain if it is a subset of a bigger one.
|
|
do idom=2,ndom
|
|
do i=1,idom-1
|
|
if(ncmo(i).ne.0) then
|
|
log=.true.
|
|
do j=1,natdom(idom)
|
|
llog=.false.
|
|
do k=1,natdom(i)
|
|
if(atdom(j,idom).eq.atdom(k,i)) then
|
|
llog=.true.
|
|
exit
|
|
endif
|
|
enddo
|
|
log=log.and.llog
|
|
enddo
|
|
if(log) then
|
|
do j=1,nmodom(idom)
|
|
llog=.false.
|
|
do k=1,nmodom(i)
|
|
if(modom(j,idom).eq.modom(k,i)) then
|
|
llog=.true.
|
|
exit
|
|
endif
|
|
enddo
|
|
log=log.and.llog
|
|
enddo
|
|
if(log) then
|
|
do j=1,ncmo(idom)
|
|
cmo(i,ncmo(i)+j)=cmo(idom,j)
|
|
enddo
|
|
ncmo(i)=ncmo(i)+ncmo(idom)
|
|
ncmo(idom)=0
|
|
exit
|
|
endif
|
|
endif
|
|
endif
|
|
enddo
|
|
enddo
|
|
ii=0
|
|
do i=1,ndom
|
|
if(ncmo(i).gt.0) then
|
|
ncmo(i-ii)=ncmo(i)
|
|
natdom(i-ii)=natdom(i)
|
|
nmodom(i-ii)=nmodom(i)
|
|
do j=1,ncmo(i-ii)
|
|
cmo(i-ii,j)=cmo(i,j)
|
|
enddo
|
|
do j=1,natdom(i-ii)
|
|
atdom(j,i-ii)=atdom(j,i)
|
|
enddo
|
|
do j=1,nmodom(i-ii)
|
|
modom(j,i-ii)=modom(j,i)
|
|
enddo
|
|
else
|
|
ii=ii+1
|
|
endif
|
|
enddo
|
|
ndom=ndom-ii
|
|
do i=1,ndom
|
|
log=.true.
|
|
do while(log)
|
|
log=.false.
|
|
do j=1,ncmo(i)-1
|
|
if(cmo(i,j).gt.cmo(i,j+1)) then
|
|
log=.true.
|
|
k=cmo(i,j)
|
|
cmo(i,j)=cmo(i,j+1)
|
|
cmo(i,j+1)=k
|
|
endif
|
|
enddo
|
|
enddo
|
|
enddo
|
|
call timer
|
|
write(iout,*)
|
|
write(iout,"(' Number of extended domains:',i4)") ndom
|
|
do idom=1,ndom
|
|
write(iout,*)
|
|
write(iout,"(' Domain',i4)") idom
|
|
write(iout,"(' Number of central MOs: ',i4)") ncmo(idom)
|
|
write(iout,"(' Number of occupied MOs:',i4)") nmodom(idom)
|
|
write(iout,"(' Number of atoms: ',i4)") natdom(idom)
|
|
write(iout,"(' Central MOs:',1000000i5)")
|
|
$(cmo(idom,j),j=1,ncmo(idom))
|
|
write(iout,"(' MOs: ',1000000i5)")
|
|
$(modom(j,idom),j=1,nmodom(idom))
|
|
write(iout,"(' Atoms: ',10000000(i4,A2))")
|
|
$ (atdom(j,idom),atsymbol(atdom(j,idom)),j=1,natdom(idom))
|
|
|
|
c do j=1,nmodom(i)
|
|
c write(6,'("atoms of MO",I3," : ",200(I3,A2))')
|
|
c $ modom(j,i),
|
|
c $ (atmo(k,modom(j,i)),atsymbol(atmo(modom(j,i),k)),
|
|
c $ k=1,natmo(modom(j,i)))
|
|
c enddo
|
|
c enddo
|
|
c stop 'teszttt '
|
|
C Truncate orbitals in domains: drop the MO coeffs which belong to
|
|
C atoms outside the domain, and minimize the least-squares residuals
|
|
C of the original and the truncated MO a la Boughton-Pulay
|
|
C Loop over domains
|
|
write(iout,*) 'Constructing PAOs...'
|
|
c do idom=1,ndom
|
|
write(cscr3,"(i9)") idom
|
|
cscr3=adjustl(cscr3)
|
|
cscr=cscr3
|
|
cscr2='MOCOEF.' // cscr
|
|
open(unit=scrfile1,file=trim(cscr2),form='unformatted')
|
|
c write(6,*) 'scrfile1= ',trim(cscr2)
|
|
cscr2='DOMAIN.' // cscr
|
|
open(unit=scrfile2,file=trim(cscr2),form='unformatted')
|
|
cscr2='FOCK.' // cscr
|
|
open(unit=scrfile3,file=trim(cscr2),form='unformatted')
|
|
rewind(scrfile1)
|
|
rewind(scrfile2)
|
|
rewind(scrfile3)
|
|
C Construct overlap and Fock matrices: drop particular rows and columns
|
|
c write(6,*) 'Domain, atoms, # of orbitals'
|
|
ii=0
|
|
do iiatoms=1,natdom(idom)
|
|
iatoms=atdom(iiatoms,idom)
|
|
c write(6,*) idom,iatoms,natrange(2,iatoms)-natrange(1,iatoms)
|
|
do mu=natrange(1,iatoms)+1,natrange(2,iatoms)
|
|
ii=ii+1
|
|
perm(ii)=mu
|
|
jj=0
|
|
do jjatoms=1,natdom(idom)
|
|
jatoms=atdom(jjatoms,idom)
|
|
do nu=natrange(1,jatoms)+1,natrange(2,jatoms)
|
|
jj=jj+1
|
|
sa(ii,jj)=s(mu,nu)
|
|
fa(ii,jj)=fock(mu,nu)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
enddo
|
|
C Invert overlap matrix in domain
|
|
call syminv(sa,si,ii,nbasis,r8heap(1))
|
|
C Count the number of core and occupied orbitals
|
|
nc=0
|
|
no=0
|
|
do imo=1,nmodom(idom)
|
|
i=modom(imo,idom)
|
|
if(i.le.ncore) then
|
|
nc=nc+1
|
|
else if(i.le.ncore+nocc) then
|
|
no=no+1
|
|
endif
|
|
enddo
|
|
nv=nmodom(idom)-nc-no
|
|
i1=nc+no-ncmo(idom)
|
|
i2=0
|
|
i3=nc+no
|
|
C Loop over MOs in the domain
|
|
do imo=1,nmodom(idom)
|
|
i=modom(imo,idom)
|
|
C Calculate new MO coefficients
|
|
C First, truncate v to the domain
|
|
ii=0
|
|
do iiatoms=1,natdom(idom)
|
|
iatoms=atdom(iiatoms,idom)
|
|
do mu=natrange(1,iatoms)+1,natrange(2,iatoms)
|
|
ii=ii+1
|
|
vi(ii)=v(mu,i)
|
|
enddo
|
|
enddo
|
|
C Second, multiply truncated v by the inverse of S
|
|
call dsymv('l',ii,1.d0,si,nbasis,vi,1,0.d0,ai,1)
|
|
C Normalize new MO
|
|
call dscal(ii,1.d0/dsqrt(ddot(ii,ai,1,vi,1)),ai,1)
|
|
c write(6,"(10f10.6)") (ai(j),j=1,ii)
|
|
C Check if MO i is a central MO in domain idom
|
|
log=.false.
|
|
do j=1,ncmo(idom)
|
|
if(i.eq.cmo(idom,j)) then
|
|
log=.true.
|
|
exit
|
|
endif
|
|
enddo
|
|
C Pack MO to matrix c
|
|
C Order: core, non-central occupied, central occupied, virtual
|
|
if(log) then
|
|
i1=i1+1
|
|
i=i1
|
|
else if(i.le.ncore+nocc) then
|
|
i2=i2+1
|
|
i=i2
|
|
else
|
|
i3=i3+1
|
|
i=i3
|
|
endif
|
|
call dcopy(ii,ai,1,c(1,i),1)
|
|
enddo
|
|
C Orthogonalize occs for test purpose
|
|
! call MO_OVERL(c=c,s=sa,ov=bmm,no=ii,mo_ini=nc+1,mo_end=nc+no,
|
|
! $ scr=syevnek,nb=nbasis)
|
|
c si=0.d0
|
|
c do j=1,ii
|
|
c do k=1,ii
|
|
c do i=nc+1,nc+no
|
|
c si(j,i)=si(j,i)+sa(j,k)*c(k,i)
|
|
c enddo
|
|
c enddo
|
|
c enddo
|
|
c bmm=0.d0
|
|
! call
|
|
! $ dsyev('V','U',no,bmm(nc+1,nc+1),nbasis,ai,syevnek,
|
|
! $ 3*nbasis,i)
|
|
! if(i.ne.0) then
|
|
! write(iout,*)'Fatal error at the orthogonalization of PAOs!'
|
|
! call mrccend(1)
|
|
! endif
|
|
! do i=nc+1,nc+no
|
|
! do j=nc+1,nc+no
|
|
! bmm(j,i)=bmm(j,i)/dsqrt(ai(i-nc))
|
|
! write(6,*) 'ai small',ai(i-nc)
|
|
! enddo
|
|
! enddo
|
|
! do j=1,ii
|
|
! do i=nc+1,nc+no
|
|
! smm(j,i)=0.d0
|
|
! do l=nc+1,nc+no
|
|
! smm(j,i)=smm(j,i)+bmm(l,i)*c(j,l)
|
|
! enddo
|
|
! enddo
|
|
! enddo
|
|
! do j=1,ii
|
|
! do l=nc+1,nc+no
|
|
! c(j,l)=smm(j,l)
|
|
! enddo
|
|
! enddo
|
|
c
|
|
c write(6,*) 'Fock'
|
|
c do i=1,ii
|
|
c do j=1,nc+no
|
|
c smm(i,j)=0.d0
|
|
c do k=1,ii
|
|
c smm(i,j)=smm(i,j)+fa(i,k)*c(k,j)
|
|
c enddo
|
|
c enddo
|
|
c enddo
|
|
c do i=1,nc+no
|
|
c do j=1,nc+no
|
|
c bmm(i,j)=0.d0
|
|
c do k=1,ii
|
|
c bmm(i,j)=bmm(i,j)+smm(k,j)*c(k,i)
|
|
c enddo
|
|
c enddo
|
|
c enddo
|
|
c call dgemm('t','n',nc+no,ii,ii,1.d0,c,nbasis,fa,nbasis,
|
|
c $ 0.d0,r8heap(1),nc+no)
|
|
c call dgemm('n','n',nc+no,nc+no,ii,1.d0,r8heap(1),
|
|
c $ nc+no,c,nbasis,0.d0,bm,nc+no)
|
|
|
|
c do j=1,nc+no
|
|
c write(6,"(i3,f10.6,' fock kicsi diag')") j,bmm(j,j)
|
|
c enddo
|
|
|
|
c write(6,*) 'nc,no,ii= ',nc,no,ii
|
|
c Project out occupied MOs from the core orbitals
|
|
bmm=0.d0
|
|
call DENSM(p=bmm,c=c,no=ii,mo_ini=nc+1,mo_end=nc+no,nb=nbasis,
|
|
$ nlarge=nbasis)
|
|
call PROJECT(c=c,p=bmm,s=sa,no=ii,project_c=sbb,mo_ini=1,
|
|
$ mo_end=nc,scr=si,nb=nbasis)
|
|
c At this point the occs have been projected out from the cores.
|
|
do k=1,ii
|
|
do i=1,nc
|
|
c(k,i)=sbb(k,i)
|
|
enddo
|
|
enddo
|
|
c Symmetric orthogonlization of core orbitals, drop linearly dependent
|
|
call MO_OVERL(c=c,s=sa,ov=smm,no=ii,mo_ini=1,mo_end=nc,
|
|
$ scr=bmm,nb=nbasis,nlarge=nbasis)
|
|
c overlap matrix of projected cores: smm
|
|
|
|
c diagonalize the overlap matrix smm
|
|
c do i=1,nc
|
|
c do j=1,nc
|
|
c write(6,*) smm(i,j),i,j,' core ovl'
|
|
c enddo
|
|
c enddo
|
|
err=0
|
|
call dsyev('V','U',nc,smm,nbasis,ai,syevnek,nbasis**2,err)
|
|
if(err.ne.0) then
|
|
write(iout,*)'Fatal error at the orthogonalization of PAOs!'
|
|
call mrccend(1)
|
|
endif
|
|
c write(6,*)'eigenvalues of core overlap'
|
|
c do k=1,nc
|
|
c write(6,*) k,ai(k)
|
|
c enddo
|
|
|
|
k=0
|
|
do j=1,natdom(idom)
|
|
k=k+ncorb(atdom(j,idom))
|
|
enddo
|
|
k=nc-k
|
|
do i=k+1,nc
|
|
do j=1,nc
|
|
smm(j,i)=smm(j,i)/dsqrt(ai(i))
|
|
enddo
|
|
enddo
|
|
do j=1,ii
|
|
do i=k+1,nc
|
|
c(j,i-k)=0.d0
|
|
do l=1,nc
|
|
c(j,i-k)=c(j,i-k)+smm(l,i)*sbb(j,l)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
do i=nc+1,ii
|
|
do j=1,ii
|
|
c(j,i-k)=c(j,i)
|
|
enddo
|
|
enddo
|
|
nc=nc-k
|
|
c write(6,*) 'Num of linearly independent cores ',nc
|
|
c write(6,*) 'Num of occupied orbitals ',no
|
|
c write(6,*) 'Num of orbitals ',ii
|
|
! At this point the cores are roughly orthogonal to the truncated occupied block
|
|
! and form an orthonormal set
|
|
|
|
C Construct PAOs
|
|
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! pao block !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
|
|
if(orblocv.eq.'pao ') then
|
|
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! pao block !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
|
|
C Calculate density matrix
|
|
bmm=0.d0
|
|
call DENSM(p=bmm,c=cold,no=nbasis,
|
|
$ mo_ini=1,mo_end=nocc+ncore,nb=nbasis,nlarge=nbasis)
|
|
|
|
call MO_OVERL(c=s,s=bmm,ov=sbb,no=nbasis,mo_ini=1,
|
|
$ mo_end=nbasis,scr=bmmm,nb=nbasis,nlarge=nbasis)
|
|
! sbb=SPS
|
|
|
|
kk=0
|
|
bmm=0.d0
|
|
do iiatoms=1,natdom(idom)
|
|
iatoms=atdom(iiatoms,idom)
|
|
do mu=natrange(1,iatoms)+1,natrange(2,iatoms)
|
|
kk=kk+1
|
|
jj=0
|
|
do jjatoms=1,natdom(idom)
|
|
jatoms=atdom(jjatoms,idom)
|
|
do nu=natrange(1,jatoms)+1,natrange(2,jatoms)
|
|
jj=jj+1
|
|
bmm(kk,jj)=sbb(mu,nu)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
enddo
|
|
! bmm=SPS in Domain
|
|
|
|
do i=1,ii
|
|
do j=1,ii
|
|
lmmm(j,i)=sa(j,i)
|
|
enddo
|
|
enddo
|
|
|
|
call dsyev('V','U',ii,lmmm,nbasis,ai,syevnek,nbasis**2,err)
|
|
|
|
c lmm(1:ii,1:ii)=0.d0
|
|
do i=1,ii
|
|
do j=1,ii
|
|
lmm(j,i)=0.d0
|
|
enddo
|
|
enddo
|
|
|
|
do i=1,ii
|
|
cl=1.d0/dsqrt(ai(i))
|
|
do j=1,ii
|
|
do k=1,ii
|
|
bmmm(j,k)=lmmm(j,i)*lmmm(k,i)
|
|
enddo
|
|
enddo
|
|
do j=1,ii
|
|
do k=1,ii
|
|
c lmm(1:ii,1:ii)=lmm(1:ii,1:ii)+cl*bmmm(1:ii,1:ii)
|
|
lmm(k,j)=lmm(k,j)+cl*bmmm(k,j)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
! coeff matrix for lowdin basis: lmm
|
|
|
|
call MO_OVERL(c=lmm,s=bmm,ov=sbb,no=ii,mo_ini=1,
|
|
$ mo_end=ii,scr=bmmm,nb=nbasis,nlarge=nbasis)
|
|
! sbb=s^-0.5 SPS s^-0.5 = ortoganal P matrix on the lowdin basis
|
|
|
|
err=0
|
|
call dsyev('V','U',ii,sbb,nbasis,ai,syevnek, nbasis**2,err)
|
|
*
|
|
if(err.ne.0) then
|
|
write(iout,*)'Fatal error at the orth of PAOs!'
|
|
call mrccend(1)
|
|
endif
|
|
|
|
nv=0
|
|
do i=1,ii
|
|
if(dabs(ai(i)).lt.1.d-1)then
|
|
nv=nv+1
|
|
c write(6,*) ai(i),' kept se of pao S ',i
|
|
do k=1,ii
|
|
c lmmm(1:ii,nv)=sbb(1:ii,i)
|
|
lmmm(k,nv)=sbb(k,i)
|
|
enddo
|
|
else
|
|
c write(6,*) ai(i),' dropped se of pao S ',i
|
|
endif
|
|
enddo
|
|
!CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
|
|
c(1:ii,no+nc+1:no+nc+nv)=0.d0
|
|
do i=no+nc+1,no+nc+nv
|
|
do j=1,ii
|
|
do k=1,ii
|
|
c(j,i)=c(j,i)+lmmm(k,i-no-nc)*lmm(j,k)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
!XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX
|
|
c write(6,*) 'Num of linearly independent cores ',nc
|
|
c write(6,*) 'Num of occupied orbitals ',no
|
|
c write(6,*) 'Num of virtuals ',nc+no+nv
|
|
c write(6,*) 'Num of orbitals ',ii
|
|
|
|
call MO_OVERL(c=c,s=sa,ov=sbb,no=ii,mo_ini=1,mo_end=nc+no+nv,
|
|
$ scr=bmmm,nb=nbasis,nlarge=nbasis)
|
|
c write(6,*) 'Overlap ',nc+no+nv
|
|
c call flush(6)
|
|
c do i=1,nc+no+nv
|
|
c do j=1,i
|
|
c if(dabs(sbb(j,i)).gt.1.d-8)then
|
|
c write(6,"(2i6,4f12.8)") i,j,sbb(j,i)
|
|
c endif
|
|
c enddo
|
|
c enddo
|
|
!XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX
|
|
c cold(1:nbasis,nocc+ncore+1:nbasis)=0.d0
|
|
do j=nocc+ncore+1,nbasis
|
|
do i=1,nbasis
|
|
cold(i,j)=0.d0
|
|
enddo
|
|
enddo
|
|
|
|
|
|
do i=nc+no+1,nc+no+nv
|
|
jj=0
|
|
do jjatoms=1,natdom(idom)
|
|
jatoms=atdom(jjatoms,idom)
|
|
do nu=natrange(1,jatoms)+1,natrange(2,jatoms)
|
|
jj=jj+1
|
|
cold(nu,nocc+ncore+i)=c(jj,i)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
|
|
call MO_OVERL(c=cold,s=s,ov=sbb,no=nbasis,
|
|
$ mo_ini=1,mo_end=nocc+ncore+nv,scr=bmmm,nb=nbasis,nlarge=nbasis)
|
|
c write(6,*) 'Overlap nagy ',nocc+ncore+nv
|
|
c call flush(6)
|
|
c do i=1,nocc+ncore+nv
|
|
c do j=1,i
|
|
c if(dabs(sbb(j,i)).gt.1.d-8)then
|
|
c write(6,"(2i6,4f12.8)") i,j,sbb(j,i)
|
|
c endif
|
|
c enddo
|
|
c enddo
|
|
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! pao block !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
|
|
endif !if pao!!!!!
|
|
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! pao block !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
|
|
c write(6,*) nc,no,nv,'nc,no,nv'
|
|
write(iout,*) 'Calculating overlap and Fock-matrices...'
|
|
C Overlap in the new basis
|
|
if(orblocv.ne.'pao ')then
|
|
nc=ncore
|
|
no=nocc
|
|
nv=nbasis-ncore-nocc
|
|
ii=nc+no+nv
|
|
do j=1,nbasis
|
|
do i=1,nbasis
|
|
fa(i,j)=fock(i,j)
|
|
enddo
|
|
enddo
|
|
endif
|
|
C Fock-matrix in the new basis
|
|
call MO_OVERL(c=c,s=fa,ov=sbb,no=ii,mo_ini=1,mo_end=nc+no+nv,
|
|
$ scr=bmmm,nb=nbasis,nlarge=nbasis)
|
|
|
|
c do j=1,nc+no
|
|
c write(6,"(i3,f10.6,' fockdiag')") j,sbb(j,j)
|
|
c enddo
|
|
c write(6,*)'**********'
|
|
c do j=nc+no+1,nc+no+nv
|
|
c write(6,"(i3,f10.6,' fockdiag')") j,sbb(j,j)
|
|
c enddo
|
|
do j=nc+no+1,nc+no+nv
|
|
do k=nc+no+1,nc+no+nv
|
|
bmm(j-nc-no,k-nc-no)=sbb(j,k)
|
|
enddo
|
|
enddo
|
|
|
|
call dsyev('V','U',nv,bmm,nbasis,ai,syevnek,nbasis**2,i)
|
|
smm=0.d0
|
|
do i=1,ii
|
|
do j=nc+no+1,nc+no+nv
|
|
do k=nc+no+1,nc+no+nv
|
|
smm(i,j)=smm(i,j)+bmm(k-nc-no,j-nc-no)*c(i,k)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
jj=nc+no
|
|
do j=nc+no+1,nc+no+nv
|
|
do i=1,ii
|
|
c(i,j)=smm(i,j)
|
|
enddo
|
|
enddo
|
|
|
|
call MO_OVERL(c=c,s=fa,ov=sbb,no=ii,mo_ini=1,mo_end=nc+no+nv,
|
|
$ scr=bmmm,nb=nbasis,nlarge=nbasis)
|
|
|
|
do j=1,nc+no
|
|
c write(6,"(i3,f10.6,' eps')") j,sbb(j,j)
|
|
eval(j)=sbb(j,j)
|
|
enddo
|
|
c write(6,*)'**********'
|
|
do j=nc+no+1,nc+no+nv
|
|
c write(6,"(i3,f10.6,' eps')") j,sbb(j,j)
|
|
eval(j)=sbb(j,j)
|
|
enddo
|
|
c call molden_printPAO(eval,c,nc+no,nbasis,nc+no+nv,ii,moldenfile,
|
|
c $ perm,cscr3)
|
|
CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
|
|
c cold(1:nbasis,nocc+ncore+1:nbasis)=0.d0
|
|
do j=nocc+ncore+1,nbasis
|
|
do i=1,nbasis
|
|
cold(i,j)=0.d0
|
|
enddo
|
|
enddo
|
|
|
|
do i=nc+no+1,nc+no+nv
|
|
jj=0
|
|
do jjatoms=1,natdom(idom)
|
|
jatoms=atdom(jjatoms,idom)
|
|
do nu=natrange(1,jatoms)+1,natrange(2,jatoms)
|
|
jj=jj+1
|
|
cold(nu,nocc+ncore+i)=c(jj,i)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
|
|
call MO_OVERL(c=cold,s=s,ov=sbb,no=nbasis,
|
|
$ mo_ini=1,mo_end=nocc+ncore+nv,scr=bmmm,nb=nbasis,nlarge=nbasis)
|
|
c write(6,*) 'Overlap kanonikus ',nocc+ncore+nv
|
|
c call flush(6)
|
|
c do i=1,nocc+ncore+nv
|
|
c do j=1,i
|
|
c if(dabs(sbb(j,i)).gt.1.d-8)then
|
|
c write(6,"(2i6,4f12.8)") i,j,sbb(j,i)
|
|
c endif
|
|
c enddo
|
|
c enddo
|
|
|
|
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
|
|
C Save MO coefficients
|
|
itol=1.d-10
|
|
nnn=0
|
|
c write(6,*) 'nn0',nnn
|
|
do i=1,nc+no+nv
|
|
do j=1,ii
|
|
if(dabs(c(j,i)).gt.itol) then
|
|
c write(6,*) 'nnn',nnn
|
|
nnn=nnn+1
|
|
c write(6,*) 'nnn2',nnn
|
|
i4heap((nnn-1)*4+3)=j
|
|
i4heap((nnn-1)*4+4)=i
|
|
c if(nn.eq.2034.or.nn.eq.1)
|
|
c write(6,*)'r8heap',nnn,c(j,i),r8heap(nnn*2),j,i
|
|
r8heap(nnn*2-1)=c(j,i)
|
|
endif
|
|
enddo
|
|
enddo
|
|
c stop
|
|
c write(6,*) 'heap size ',2*nn
|
|
c write(6,*) 'write 880 !!! '
|
|
DO J=1,nnn
|
|
write(880,*) J,R8HEAP(2*j-1),R8HEAP(2*j),i4heap((j-1)*4+3),
|
|
$ i4heap((j-1)*4+4)
|
|
ENDDO
|
|
write(scrfile1) nnn,(r8heap(j),j=1,2*nnn)
|
|
C Save Fock-matrix
|
|
c write(6,*) 'F= '
|
|
nnn=0
|
|
do i=1,ii
|
|
do j=1,i
|
|
if(dabs(fa(i,j)).gt.itol) then
|
|
nnn=nnn+1
|
|
i4heap((nnn-1)*4+3)=i
|
|
i4heap((nnn-1)*4+4)=j
|
|
r8heap(nnn*2-1)=fa(i,j)
|
|
c write(6,*) fa(i,j),i,j
|
|
endif
|
|
enddo
|
|
enddo
|
|
write(scrfile3) nnn,(r8heap(j),j=1,2*nnn)
|
|
C # basis functions, # electrons, ncore, # MOs
|
|
write(scrfile2) ii,2*(nc+no),nc,nc+no+nv,nangmax,ncontrmax,
|
|
$nprimmax,nshmax
|
|
write(scrfile2) nuc
|
|
write(scrfile2) ncmo(idom)
|
|
write(scrfile2) (cmo(idom,j),j=1,ncmo(idom))
|
|
write(scrfile2) natdom(idom)
|
|
do iiatoms=1,natdom(idom)
|
|
iatoms=atdom(iiatoms,idom)
|
|
write(scrfile2) atnum(iatoms),atchg(iatoms),
|
|
$(coord(j,iatoms),j=1,3),(basname(iatoms,k),k=1,nbset),
|
|
$ecpname(iatoms),atsymbol(iatoms),atnumt(iatoms),ncorecp(iatoms)
|
|
enddo
|
|
C
|
|
close(scrfile1)
|
|
close(scrfile2)
|
|
close(scrfile3)
|
|
C Timings
|
|
call timer
|
|
enddo
|
|
C
|
|
open(unit=scrfile1,file='DOMAIN')
|
|
rewind(scrfile1)
|
|
write(scrfile1,*) ndom
|
|
close(scrfile1)
|
|
cDDDDDDDDDDDDDDDDDDDDDDDDDDDDDDDD
|
|
ENDIF !localcc.ne.'off '
|
|
cDDDDDDDDDDDDDDDDDDDDDDDDDDDDDDDD
|
|
|
|
close(oeintfile)
|
|
close(mocoeffile)
|
|
close(moldenfile)
|
|
|
|
call mrccend(0)
|
|
|
|
end
|
|
C
|
|
subroutine ipoint1d(egydim1,egydim2,dim1,memr)
|
|
implicit none
|
|
integer dim1
|
|
integer*8 memr
|
|
integer,target :: egydim1(dim1)
|
|
integer, pointer :: egydim2(:)
|
|
egydim2 => egydim1(1:dim1)
|
|
memr=memr+dim1
|
|
end
|
|
subroutine rpoint1d(egydim1,egydim2,dim1,memr)
|
|
implicit none
|
|
integer dim1
|
|
integer*8 memr
|
|
real*8,target :: egydim1(dim1)
|
|
real*8, pointer :: egydim2(:)
|
|
egydim2 => egydim1(1:dim1)
|
|
memr=memr+dim1
|
|
end
|
|
subroutine ipoint1d_half(egydim1,egydim2,dim1,memr)
|
|
implicit none
|
|
integer dim1
|
|
integer*8 memr
|
|
integer*4,target :: egydim1(dim1)
|
|
integer*4, pointer :: egydim2(:)
|
|
egydim2 => egydim1(1:dim1)
|
|
memr=memr+dim1/2+1
|
|
end
|
|
subroutine ipoint2d(egydim,ketdim,dim1,dim2,memr)
|
|
implicit none
|
|
integer dim1,dim2
|
|
integer*8 memr
|
|
integer,target :: egydim(dim1,dim2)
|
|
integer, pointer :: ketdim(:,:)
|
|
ketdim => egydim(1:dim1,1:dim2)
|
|
memr=memr+dim1*dim2
|
|
end
|
|
subroutine rpoint2d(egydim,ketdim,dim1,dim2,memr)
|
|
implicit none
|
|
integer dim1,dim2
|
|
integer*8 memr
|
|
real*8,target :: egydim(dim1,dim2)
|
|
real*8, pointer :: ketdim(:,:)
|
|
ketdim => egydim(1:dim1,1:dim2)
|
|
memr=memr+dim1*dim2
|
|
end
|
|
|
|
|
|
|
|
SUBROUTINE molden_printPAO(eigenvalue,c,ndocc,nbasis,nb,ii,
|
|
$ moldenfile,perm2,cscr)
|
|
implicit none
|
|
************************************************************************
|
|
* subroutine to print MO coefficients to molden format text file
|
|
************************************************************************
|
|
!#include "MRCCCOMMON"
|
|
#define _RHF_ 1
|
|
#define _UHF_ 2
|
|
#define _ROHF_ 3
|
|
c variables
|
|
c arguments
|
|
real*8 eigenvalue(nbasis) !the orbital energies
|
|
real*8 c(nbasis,nbasis) !AOMO coefficients matrix
|
|
integer ndocc,nbasis,nb,ii !number of doubly occupied spinorbitals
|
|
integer nsocc,moldenfile !number of singly occupied spinorbitals
|
|
integer noccsp(2) !number of occupied spinorbitals alpha,beta
|
|
integer calcway
|
|
c arrays
|
|
character trash(15)
|
|
character*3 moldenval
|
|
character*4 mark
|
|
character*26 cscr,cscr2
|
|
character*72 cscr3
|
|
integer perm(nbasis),perm2(nbasis),inv(nbasis)
|
|
c scalars
|
|
integer q,i,j
|
|
real*8 occupation
|
|
C program
|
|
c getting if molden input is desired
|
|
c calcway=1
|
|
c write(6,*) ' halott '
|
|
c call flush(6)
|
|
cscr2='MOLDEN.PAO.' // cscr
|
|
call getkey('molden',6,moldenval,3)
|
|
select case (moldenval)
|
|
case ('off')
|
|
return
|
|
case ('on ')
|
|
case (' ')
|
|
case default
|
|
call unknown('molden',3)
|
|
end select
|
|
C open file to get permutation vector
|
|
cscr3='cp MOLDEN MOLDEN.PAO.' // cscr
|
|
call system(trim(cscr3))
|
|
open(unit=moldenfile,file='MOLDEN.perm')
|
|
read(moldenfile,*) (perm(q),q=1,nbasis)
|
|
close(moldenfile)
|
|
do q=1,nbasis
|
|
inv(perm(q))=q
|
|
enddo
|
|
c open the file for appending MO coefficients
|
|
c write(6,*) ' barki '
|
|
c call flush(6)
|
|
open(unit=moldenfile,file=trim(cscr2),status='old')
|
|
rewind(moldenfile)
|
|
mark='aaaa'
|
|
do while (mark.ne.'[MO]')
|
|
read(moldenfile,"(4a)")mark
|
|
enddo
|
|
c select case(calcway)
|
|
c case(_RHF_)
|
|
c write(6,*) ' idodon '
|
|
c call flush(6)
|
|
do i=1,nb
|
|
if(i.le.ndocc) then
|
|
occupation=2.d0
|
|
else
|
|
occupation=0.d0
|
|
endif
|
|
write(moldenfile,"(' Ene=',f14.4)") min(eigenvalue(i),1.d-7)!NP
|
|
write(moldenfile,"(' Spin= Alpha')")
|
|
write(moldenfile,"(' Occup=',f14.4)") occupation
|
|
do j=1,ii
|
|
write(moldenfile,"(I4,f18.10)")
|
|
$ perm(perm2(j)),c(j,i)
|
|
enddo
|
|
enddo
|
|
c write(6,*) ' barki '
|
|
c call flush(6)
|
|
c case(_UHF_)
|
|
c do i=1,nbasis !alpha electrons
|
|
c if(i.le.noccsp(1)) then
|
|
c occupation=1.d0
|
|
c else
|
|
c occupation=0.d0
|
|
c endif
|
|
c write(moldenfile,"(' Ene=',f14.4)") min(eigenvalue(perm(i)),1.d-7)!NP
|
|
c write(moldenfile,"(' Spin= Alpha')")
|
|
c write(moldenfile,"(' Occup=',f14.4)") occupation
|
|
c do j=1,nbasis
|
|
c write(moldenfile,"(I4,f18.10)")
|
|
c & j,c((i-1)*nbasis+perm(j))
|
|
c enddo
|
|
c enddo
|
|
c do i=1,nbasis !beta electrons
|
|
c if(i.le.noccsp(2)) then
|
|
c occupation=1.d0
|
|
c else
|
|
c occupation=0.d0
|
|
c endif
|
|
c write(moldenfile,"(' Ene=',f14.4)") min(eigenvalue(perm(i)),1.d-7)!NP
|
|
c write(moldenfile,"(' Spin= Beta')")
|
|
c write(moldenfile,"(' Occup=',f14.4)") occupation
|
|
c do j=1,nbasis
|
|
c write(moldenfile,"(I4,f18.10)")
|
|
c & j,c(nbasis**2+(i-1)*nbasis+perm(j))
|
|
c enddo
|
|
c enddo
|
|
c case(_ROHF_)
|
|
c%%% itt az alpha es beta indexek majdnem teljesesn ertelmetlenek....????
|
|
c do i=1,nbasis !alpha electrons
|
|
c if(i.le.ndocc) then
|
|
c occupation=2.d0
|
|
c else if (i.gt.ndocc.and.i.le.ndocc+nsocc) then
|
|
c occupation=1.d0
|
|
c else
|
|
c occupation=0.d0
|
|
c endif
|
|
c write(moldenfile,"(' Ene=',f14.4)") min(eigenvalue(perm(i)),1.d-7)!NP
|
|
c write(moldenfile,"(' Spin= Alpha')")
|
|
c write(moldenfile,"(' Occup=',f14.4)") occupation
|
|
c do j=1,nbasis
|
|
c write(moldenfile,"(I4,f18.10)") j,c((i-1)*nbasis+perm(j))
|
|
c enddo
|
|
c enddo
|
|
c end select
|
|
c closing and exiting
|
|
close(moldenfile)
|
|
c
|
|
return
|
|
end
|
|
subroutine DENSM(bmm,c,ii,mo_ini,mo_end,nb,nlarge)
|
|
implicit none
|
|
integer i,j,k,nb,ii,mo_ini,mo_end,nlarge
|
|
real*8 bmm(nb,nb),c(nb,nlarge)
|
|
c bmm=0.d0
|
|
do j=1,ii
|
|
do k=1,ii
|
|
do i=mo_ini,mo_end
|
|
bmm(k,j)=bmm(k,j)+c(j,i)*c(k,i)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
end
|
|
subroutine PROJECT_B(c,p,s,ii,pc,mo_ini,mo_end,si,nb)
|
|
implicit none
|
|
integer i,j,k,ii,nb,mo_ini,mo_end
|
|
real*8 c(nb,nb),s(nb,nb),p(nb,nb),si(nb,nb),pc(nb,nb)
|
|
si=0.d0
|
|
do i=1,ii
|
|
do j=1,ii
|
|
do k=1,nb
|
|
si(j,i)=si(j,i)-s(j,k)*p(k,i)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
do i=1,ii
|
|
si(i,i)=si(i,i)+1.d0
|
|
enddo
|
|
pc=0.d0
|
|
do i=1,ii
|
|
do j=mo_ini,mo_end
|
|
do k=1,ii
|
|
pc(i,j)=pc(i,j)+si(k,i)*c(k,j)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
end subroutine
|
|
subroutine PROJECT(c,p,s,ii,pc,mo_ini,mo_end,si,nb)
|
|
implicit none
|
|
integer i,j,k,ii,nb,mo_ini,mo_end
|
|
real*8 c(nb,nb),s(nb,nb),p(nb,nb),si(nb,nb),pc(nb,nb)
|
|
si=0.d0
|
|
do i=1,ii
|
|
do j=1,ii
|
|
do k=1,ii
|
|
si(j,i)=si(j,i)-s(j,k)*p(k,i)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
do i=1,ii
|
|
si(i,i)=si(i,i)+1.d0
|
|
enddo
|
|
pc=0.d0
|
|
do i=1,ii
|
|
do j=mo_ini,mo_end
|
|
do k=1,ii
|
|
pc(i,j)=pc(i,j)+si(k,i)*c(k,j)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
end subroutine
|
|
subroutine MO_OVERL(c,s,ov,ii,mo_ini,mo_end,si,nb,nlarge)
|
|
implicit none
|
|
integer i,j,k,ii,nb,mo_ini,mo_end,nlarge
|
|
real*8 c(nb,nlarge),s(nb,nb)
|
|
real*8 ov(nlarge,nlarge),si(nb,nlarge)
|
|
si=0.d0
|
|
do j=1,ii
|
|
do k=1,ii
|
|
do i=mo_ini,mo_end
|
|
si(j,i)=si(j,i)+s(j,k)*c(k,i)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
ov=0.d0
|
|
do k=mo_ini,mo_end
|
|
do i=mo_ini,mo_end
|
|
do j=1,ii
|
|
ov(k,i)=ov(k,i)+si(j,i)*c(j,k)
|
|
enddo
|
|
enddo
|
|
enddo
|
|
end subroutine
|
|
C
|
|
|