mirror of
https://code.it4i.cz/sccs/easyconfigs-it4i.git
synced 2025-04-17 20:20:49 +01:00
743 lines
25 KiB
Fortran
743 lines
25 KiB
Fortran
************************************************************************
|
|
subroutine localize(job,nmatrix,loctype,m,tresh_macro,tresh_micro,
|
|
&itol,natoms,am,mnbasis,printfile,orblocguess,nos,moa,no,lembed,
|
|
&uboys,use_beta)
|
|
************************************************************************
|
|
C Subroutine for localization of occupied/virtual orbitals
|
|
************************************************************************
|
|
#include "MRCCCOMMON"
|
|
real*8 m,tresh_macro,tresh_micro,itol,am(nbasis,mnbasis)
|
|
real*8 moa(nbasis,mnbasis)
|
|
character job
|
|
character*4 loctype,c4
|
|
character*8 orblocguess,boysalg !NP
|
|
integer nmatrix,nn,natoms,inatrange,is,iv,mnbasis,printfile,ne,i,j
|
|
integer dblalloc,intalloc,isi,iss,idens,idensa,idensb,verbosity,no
|
|
integer ixdip,iydip,izdip,ixxsec,iyysec,izzsec,ism,ipm,isd,nos
|
|
integer idensity,iumatrix
|
|
integer imocoef_chol,imocoef_lmo,imocoef_cmo
|
|
integer iworkp,iwork,ijpvt,iwread,iwtransf
|
|
integer ixi1,ixi2,ixi2prec,ikappa,ikappaold,iq,iomega1,iomega2
|
|
integer ixi1m,ixi2m,ikappam,ikq,iqk,ikx,iky,ikz,ixk,iyk,izk
|
|
integer ibasis_vec,itr_basis_vec,ireduced_hessian,max_macro_it
|
|
integer max_micro_it,work_mem
|
|
integer npreopt,imem_old,imem1
|
|
integer nfroz_embed
|
|
double precision uboys(nmatrix,nmatrix)
|
|
real*8 pao_subsys_tol
|
|
character*10 c10
|
|
character*16 c16
|
|
logical lembed
|
|
logical use_beta
|
|
C
|
|
integer ixdip_mo
|
|
integer iydip_mo
|
|
integer izdip_mo
|
|
integer ixdip_mo2
|
|
integer iydip_mo2
|
|
integer izdip_mo2
|
|
common/memcom/ imem1
|
|
|
|
imem_old=imem
|
|
call getkey('verbosity',9,c4,4)
|
|
read(c4,'(i4)') verbosity
|
|
call getkey('boysalg',7,boysalg,8)
|
|
C
|
|
nn=(nmatrix-1)*nmatrix/2
|
|
C Define pointers
|
|
ixdip=dblalloc(nmatrix**2)
|
|
iydip=dblalloc(nmatrix**2)
|
|
izdip=dblalloc(nmatrix**2)
|
|
C
|
|
ixxsec=dblalloc(nmatrix**2)
|
|
iyysec=dblalloc(nmatrix**2)
|
|
izzsec=dblalloc(nmatrix**2)
|
|
C
|
|
imocoef_chol=dblalloc(nmatrix*nbasis)
|
|
C
|
|
imocoef_cmo=dblalloc(nmatrix*nbasis)
|
|
C%%% Read MO coefficients from file (canonical) -> kesobb egybol
|
|
C suruseg!
|
|
CCC start HB
|
|
call read_mocoef(job,nmatrix,nbasis,ncore,no,nos,
|
|
$moa,dcore(imocoef_cmo),dcore(imem),dcore(imem),nfroz)
|
|
CCC end HB
|
|
C
|
|
if (orblocguess.eq.'cholesky ') then ! NP
|
|
idensity=dblalloc(nbasis**2)
|
|
C%%% Calculate density matrix
|
|
call set_density(nmatrix,nbasis,dcore(idensity),
|
|
&dcore(imocoef_cmo))
|
|
C
|
|
iworkp=dblalloc(nmatrix*nbasis)
|
|
iwork=dblalloc(nbasis)
|
|
ijpvt=intalloc(nbasis)
|
|
C Cholesky decomposition of the density matrix
|
|
call init_cholesky(nmatrix,nbasis,dcore(imocoef_chol),
|
|
&dcore(idensity),dcore(iworkp),dcore(iwork),icore(ijpvt))
|
|
call dbldealloc(idensity) ! NP
|
|
elseif (orblocguess.eq.'restart '.or.orblocguess.eq.'read ') !NP
|
|
$then
|
|
call dcopy(nmatrix*nbasis,dcore(imocoef_cmo),1,
|
|
$ dcore(imocoef_chol),1)
|
|
if (orblocguess.eq.'read ') go to 111
|
|
endif
|
|
C**********************************************************
|
|
C Pipek-Mezey localization
|
|
if(loctype.eq."pm ") then
|
|
inatrange=intalloc(2*natoms*nbset)
|
|
is=dblalloc(nbasis**2)
|
|
iv=dblalloc(nbasis**2)
|
|
call motransp(nmatrix,nbasis,dcore(imocoef_chol),
|
|
$dcore(imocoef_cmo),.true.)
|
|
rewind(oeintfile)
|
|
call roeint(dcore(imem),dcore(imem),dcore(is),oeintfile,nbasis)
|
|
call pml(nmatrix,nbasis,natoms,icore(inatrange),dcore(is),
|
|
$dcore(imocoef_cmo),dcore(iv),dcore(imem),.true.)
|
|
c call pmlvec(nmatrix,nbasis,natoms,icore(inatrange),dcore(is),
|
|
c $dcore(imocoef_cmo),dcore(iv),dcore(imem),.true.)
|
|
call motransp(nmatrix,nbasis,dcore(imocoef_chol),
|
|
$dcore(imocoef_cmo),.false.)
|
|
call dbldealloc(inatrange)
|
|
endif
|
|
C**********************************************************
|
|
C Fragment localization
|
|
if(loctype.eq.'frag'.or.loctype.eq.'no '.or.
|
|
$ loctype.eq.'spad'.or.loctype.eq.'pao') then
|
|
ne = 0
|
|
mnbasis=nbf(4)
|
|
inatrange=intalloc(2*natoms*nbset)
|
|
call getvar('natrange ',icore(inatrange))
|
|
if( loctype.eq.'pao') then
|
|
call getkey('pao_subsys_tol',14,c16,16)
|
|
read(c16,*) pao_subsys_tol
|
|
nfroz_embed = 0
|
|
if(lembed) call getvar('nfroz ',nfroz_embed)
|
|
endif
|
|
isi=dblalloc(nbasis**2)
|
|
iss=dblalloc(nbasis**2)
|
|
ism=dblalloc(nbasis**2)
|
|
idens=dblalloc(nbasis**2)
|
|
idensa=dblalloc(nbasis**2)
|
|
idensb=dblalloc(nbasis**2)
|
|
ipm=dblalloc(nbasis*mnbasis)
|
|
isd=dblalloc(mnbasis**2)
|
|
c if(job.eq.'v') then
|
|
call dsyrk('u','t',nbasis,nmatrix,2.d0,dcore(imocoef_cmo),
|
|
$nmatrix,0.d0,dcore(idens),nbasis)
|
|
call filllo(dcore(idens),nbasis)
|
|
c else
|
|
c open(scrfile1,file='SCFDENSITIES',form='UNFORMATTED')
|
|
c call roeint(dcore(imem),dcore(imem),dcore(idens),scrfile1,
|
|
c $nbasis)
|
|
c close(scrfile1)
|
|
c endif
|
|
C
|
|
c rewind(oeintfile)
|
|
c call roeint(dcore(imem),dcore(imem),dcore(ism),oeintfile,nbasis)
|
|
open(scrfile1,file='SROOT',form='unformatted')
|
|
call roeint(dcore(imem),dcore(imem),dcore(iss),scrfile1,nbasis)
|
|
call roeint(dcore(imem),dcore(imem),dcore(isi),scrfile1,nbasis)
|
|
close(scrfile1)
|
|
if(trim(loctype).eq.'no') then
|
|
call noloc(nbasis,dcore(idens),dcore(imem+natoms**2+nbasis),
|
|
$dcore(imem+natoms**2),dcore(isi),dcore(iss),dcore(imocoef_cmo),
|
|
$dcore(idensa),natoms,icore(inatrange),nmatrix,2.d0,verbosity,
|
|
$dcore(idensb),dcore(ism),mnbasis,dcore(ipm),dcore(isd),scrfile1,
|
|
$dcore(imem))
|
|
else if(trim(loctype).eq.'pao') then
|
|
call pao_loc(iintln,scrfile1,oeintfile,minpfile,mocoeffile,
|
|
$iout,lembed,verbosity,icore(inatrange),nbasis,nmatrix,no,ncore,
|
|
$natoms,nfroz_embed,nvfroz,pao_subsys_tol,
|
|
$dcore(imocoef_cmo),dcore(isi),dcore(idens),dcore(idensa),
|
|
$dcore(idensb),dcore(imem),dcore(imem+nbasis+1))
|
|
call motransp(nmatrix,nbasis,dcore(imocoef_chol),
|
|
$dcore(imocoef_cmo),.false.)
|
|
else if(trim(loctype).eq.'spad') then
|
|
block
|
|
integer :: iheap1a_natoms, iheap1b_natoms
|
|
integer :: rheap1_nmatrix
|
|
integer :: ncorenew ! dummy
|
|
iheap1a_natoms = intalloc(natoms)
|
|
iheap1b_natoms = intalloc(natoms)
|
|
rheap1_nmatrix = dblalloc( nmatrix + 1 )
|
|
call dcopy
|
|
$(nbasis*nmatrix,dcore(imocoef_cmo),1,dcore(idens),1)
|
|
call motransp
|
|
$(nmatrix,nbasis,dcore(idens),dcore(imocoef_cmo),.true.)
|
|
call spade
|
|
$( nbasis , nmatrix , ne , dcore( imocoef_cmo ) , dcore( idens ) ,
|
|
$ dcore( imem ), dcore( rheap1_nmatrix ) , dcore( iss ) , natoms ,
|
|
$ minpfile , icore(inatrange) , verbosity, ncorenew ,
|
|
$ icore( iheap1a_natoms ) , icore( iheap1b_natoms ), lembed , job )
|
|
call dbldealloc( rheap1_nmatrix )
|
|
call intdealloc( iheap1b_natoms )
|
|
call intdealloc( iheap1a_natoms )
|
|
end block
|
|
else
|
|
call fragloc(nbasis,dcore(idens),dcore(imem+nbasis),
|
|
$dcore(imem),dcore(isi),dcore(iss),dcore(imocoef_cmo),
|
|
$dcore(idensa),natoms,minpfile,icore(inatrange),nmatrix,ne,2.d0,i,
|
|
$verbosity,dcore(idensb),.true.,j,dcore(ism),mnbasis,dcore(ipm),
|
|
$dcore(isd),scrfile1)
|
|
endif
|
|
call motransp(nmatrix,nbasis,dcore(imocoef_chol),
|
|
$dcore(imocoef_cmo),.false.)
|
|
call intdealloc(inatrange)
|
|
open(scrfile1,file='VARS',form='unformatted',position='append') !HB
|
|
if(job.eq.'o'.or.job.eq.'c')
|
|
$ c10 = merge( 'nembedb ' , 'nembeda ' , use_beta )
|
|
if(job.eq.'v')
|
|
$ c10 = merge( 'nvembedb ' , 'nvembeda ' , use_beta )
|
|
if(job.eq.'s')
|
|
$ c10 = merge( 'nsembedb ' , 'nsembeda ' , use_beta )
|
|
write(scrfile1) c10,iintln,ne
|
|
close(scrfile1) !HB
|
|
endif
|
|
C**********************************************************
|
|
C Intrinsic bond orbitals (IBOs)
|
|
if(loctype.eq."ibo ") then
|
|
inatrange=intalloc(2*natoms*nbset)
|
|
is=dblalloc(nbasis**2)
|
|
call motransp(nmatrix,nbasis,dcore(imocoef_chol),
|
|
$dcore(imocoef_cmo),.true.)
|
|
rewind(oeintfile)
|
|
call roeint(dcore(imem),dcore(imem),dcore(is),oeintfile,nbasis)
|
|
call ibo(nmatrix,nbasis,natoms,icore(inatrange),dcore(is),
|
|
$dcore(imocoef_cmo),dcore(imem),.true.,am,mnbasis)
|
|
call motransp(nmatrix,nbasis,dcore(imocoef_chol),
|
|
$dcore(imocoef_cmo),.false.)
|
|
call dbldealloc(inatrange)
|
|
endif
|
|
call dbldealloc(imocoef_cmo)
|
|
C
|
|
C**********************************************************
|
|
C Boys Localization using Jacobi algorithm
|
|
if(loctype.eq."boys".and.boysalg.eq.'jacobi ') then
|
|
iwread=dblalloc(nbasis**2)
|
|
iwtransf=dblalloc(nmatrix*nbasis)
|
|
C Transform AO integrals with the Cholesky MO coefficients
|
|
call set_moints(nmatrix,nbasis,printfile,dcore(imocoef_chol),
|
|
& dcore(ixdip),dcore(iydip),dcore(izdip),dcore(imem),
|
|
& dcore(iwread),dcore(iwtransf),dcore(imem),dcore(imem),.true.)
|
|
call dbldealloc(iwread)
|
|
C
|
|
iwork=dblalloc(2*nmatrix**2)
|
|
C Take a Boys step for initial guess
|
|
call init_boys(1000,nmatrix,dcore(ixdip),dcore(iydip),
|
|
$dcore(izdip),uboys,dcore(iwork))
|
|
call dbldealloc(iwork)
|
|
endif
|
|
C**********************************************************
|
|
C Boys Localization using Newton algorithm
|
|
if(loctype.eq."boys".and.boysalg.eq.'newton ') then
|
|
uboys=0.0d0
|
|
do i=1,nmatrix
|
|
uboys(i,i)=1.0d0
|
|
enddo
|
|
if(nn.lt.1) then
|
|
call dbldealloc(imem_old)
|
|
return
|
|
endif
|
|
C Define pointers
|
|
ixdip=dblalloc(nmatrix**2)
|
|
iydip=dblalloc(nmatrix**2)
|
|
izdip=dblalloc(nmatrix**2)
|
|
iq=dblalloc(nmatrix**2)
|
|
C
|
|
i=dblalloc(nbasis**2-nmatrix**2)
|
|
iwread=dblalloc(nbasis**2)
|
|
iwtransf=dblalloc(nmatrix*nbasis)
|
|
C
|
|
C Transform AO integrals with the Cholesky MO coefficients
|
|
call set_moints(nmatrix,nbasis,printfile,dcore(imocoef_chol),
|
|
& dcore(ixdip),dcore(iydip),dcore(izdip),dcore(iq),
|
|
& dcore(iwread),dcore(iwtransf),dcore(imem),dcore(imem),.false.)
|
|
call dbldealloc(i)
|
|
C
|
|
iwork=dblalloc(2*nmatrix**2)
|
|
C Take a Boys step for initial guess
|
|
if(job.eq.'x') then
|
|
npreopt=0
|
|
else
|
|
npreopt=0
|
|
endif
|
|
if(npreopt.ne.0)
|
|
& write(*,"(' Preoptimizing orbitals using Jacobi algorithm')")
|
|
call init_boys(npreopt,nmatrix,dcore(ixdip),dcore(iydip),
|
|
$dcore(izdip),uboys,dcore(iwork))
|
|
call dbldealloc(iwork)
|
|
C
|
|
max_macro_it=500
|
|
ixi1=dblalloc(nn)
|
|
ixi2=dblalloc(nn)
|
|
ixi2prec=dblalloc(nn)
|
|
ikappa=dblalloc(nn)
|
|
ikappaold=dblalloc(nn)
|
|
iomega1=dblalloc(nn)
|
|
ixi1m=dblalloc(nmatrix**2)
|
|
ikq=dblalloc(nmatrix**2)
|
|
ikx=dblalloc(nmatrix**2)
|
|
iky=dblalloc(nmatrix**2)
|
|
ikz=dblalloc(nmatrix**2)
|
|
|
|
work_mem=maxcor-(imem-imem1)
|
|
max_micro_it=min(100,work_mem/(2*nn)-2)
|
|
if(max_micro_it.lt.30) then
|
|
write(*,'(1X,A)')
|
|
& 'Not enough memory for trust region localization'
|
|
call mrccend(1)
|
|
endif
|
|
i=imem
|
|
ibasis_vec=dblalloc(nn*max_micro_it)
|
|
itr_basis_vec=dblalloc(nn*max_micro_it)
|
|
ireduced_hessian=dblalloc(max_micro_it**2)
|
|
ikappam=dblalloc(max_micro_it)
|
|
iwork=dblalloc(max(nn,10*max_micro_it**2))
|
|
C Transform integrals
|
|
if(m.ne.1.0d0) then
|
|
write(*,'(1X,A)') 'Preoptimizing orbitals with m=1.0'
|
|
call transform(1.0d0,nmatrix,tresh_macro,tresh_micro,itol,
|
|
& dcore(ixdip),dcore(iydip),dcore(izdip),uboys,
|
|
& dcore(ixi1),dcore(ixi2),dcore(ixi2prec),dcore(ikappa),
|
|
& dcore(ikappaold),dcore(iq),dcore(iomega1),
|
|
& dcore(ixi1m),dcore(ikappam),dcore(ikq),
|
|
& dcore(ikx),dcore(iky),dcore(ikz),
|
|
& dcore(ibasis_vec),dcore(itr_basis_vec),
|
|
& dcore(ireduced_hessian),max_macro_it,max_micro_it,dcore(iwork),
|
|
& dcore(i),work_mem)
|
|
call transform(m,nmatrix,tresh_macro,tresh_micro,itol,
|
|
& dcore(ixdip),dcore(iydip),dcore(izdip),
|
|
& uboys,dcore(ixi1),
|
|
& dcore(ixi2),dcore(ixi2prec),dcore(ikappa),dcore(ikappaold),
|
|
& dcore(iq),dcore(iomega1),
|
|
& dcore(ixi1m),dcore(ikappam),dcore(ikq),
|
|
& dcore(ikx),dcore(iky),dcore(ikz),
|
|
& dcore(ibasis_vec),dcore(itr_basis_vec),
|
|
& dcore(ireduced_hessian),max_macro_it,max_micro_it,dcore(iwork),
|
|
& dcore(i),work_mem)
|
|
else
|
|
call transform(m,nmatrix,tresh_macro,tresh_micro,itol,
|
|
& dcore(ixdip),dcore(iydip),dcore(izdip),
|
|
& uboys,dcore(ixi1),
|
|
& dcore(ixi2),dcore(ixi2prec),dcore(ikappa),dcore(ikappaold),
|
|
& dcore(iq),dcore(iomega1),
|
|
& dcore(ixi1m),dcore(ikappam),dcore(ikq),
|
|
& dcore(ikx),dcore(iky),dcore(ikz),
|
|
& dcore(ibasis_vec),dcore(itr_basis_vec),
|
|
& dcore(ireduced_hessian),max_macro_it,max_micro_it,dcore(iwork),
|
|
& dcore(i),work_mem)
|
|
endif
|
|
call dbldealloc(ixi1)
|
|
endif
|
|
C
|
|
111 continue
|
|
imocoef_lmo=dblalloc(nmatrix*nbasis) ! NP
|
|
iwork=dblalloc(nmatrix**2)
|
|
C Transform MO coefficients
|
|
call transform_mocoef(job,loctype,nmatrix,nbasis,ncore,no,
|
|
&nos,dcore(imocoef_chol),uboys,uboys,
|
|
&dcore(imocoef_lmo),moa,dcore(iwork),nfroz,orblocguess)
|
|
call dbldealloc(imem_old)
|
|
C
|
|
end
|
|
C
|
|
************************************************************************
|
|
************************************************************************
|
|
subroutine read_mocoef(
|
|
&job,
|
|
&nmatrix,
|
|
&nbasis,
|
|
&ncore,
|
|
&nocc,
|
|
&nos,
|
|
&mocoeffull,
|
|
&mocoef_cmo,
|
|
&r8heap,
|
|
&i4heap,
|
|
&nfroz) ! HB
|
|
************************************************************************
|
|
C Subroutine to read molecular orbital coefficients from file
|
|
************************************************************************
|
|
implicit none
|
|
C
|
|
C Job control, if (job.eq.'o') then we read the occupied MO coefficients
|
|
C and if (job.eq.'v') then we read the virtual ones
|
|
character job
|
|
C
|
|
C Size of marices, number of basis functions, core, occupied and virtual
|
|
C orbitals
|
|
integer nmatrix,nbasis,ncore,nocc,nos
|
|
C
|
|
C Matrix of the coefficients (full)
|
|
real*8 mocoeffull(nbasis,*)
|
|
C
|
|
C Arrays for reading coefficients
|
|
real*8 r8heap(*)
|
|
integer*4 i4heap(*)
|
|
C
|
|
C Matrix of the coefficients (canonical)
|
|
real*8 mocoef_cmo(nmatrix,*)
|
|
C
|
|
integer i,j,ii,jj,nn,nfroz ! HB
|
|
C
|
|
C
|
|
C%%% Veglegesben majd csak a nemnulla elemeket taroljuk!
|
|
call dfillzero(mocoef_cmo,nbasis*nmatrix)
|
|
C
|
|
if(job.eq.'c') then
|
|
do i=1,nmatrix
|
|
do j=1,nbasis
|
|
CCC start HB
|
|
mocoef_cmo(i,j)=mocoeffull(j,i+nfroz)
|
|
CCC end HB
|
|
enddo
|
|
enddo
|
|
else if(job.eq.'o') then
|
|
do i=1,nmatrix
|
|
do j=1,nbasis
|
|
mocoef_cmo(i,j)=mocoeffull(j,i+ncore)
|
|
enddo
|
|
enddo
|
|
else if(job.eq.'s') then
|
|
do i=1,nmatrix
|
|
do j=1,nbasis
|
|
mocoef_cmo(i,j)=mocoeffull(j,i+ncore+nocc)
|
|
enddo
|
|
enddo
|
|
else if(job.eq.'v') then
|
|
do i=1,nmatrix
|
|
do j=1,nbasis
|
|
mocoef_cmo(i,j)=mocoeffull(j,i+ncore+nocc+nos)
|
|
enddo
|
|
enddo
|
|
elseif(job.eq.'x') then
|
|
do i=1,nmatrix
|
|
do j=1,nbasis
|
|
mocoef_cmo(i,j)=mocoeffull(j,i)
|
|
enddo
|
|
enddo
|
|
endif
|
|
C
|
|
end
|
|
C
|
|
C***********************************************************************
|
|
C***********************************************************************
|
|
C Subroutine to transform MO coefficients
|
|
C***********************************************************************
|
|
C***********************************************************************
|
|
C
|
|
subroutine transform_mocoef(
|
|
&job,
|
|
&loctype,
|
|
&nmatrix,
|
|
&nbasis,
|
|
&ncore,
|
|
&nocc,
|
|
&nos,
|
|
&mocoef_chol,
|
|
&uboysinit,
|
|
&uboys,
|
|
&mocoef_lmo,
|
|
&mocoeffull,
|
|
&work,
|
|
&nfroz,orblocguess) ! HB
|
|
C
|
|
implicit none
|
|
C
|
|
C Job control, if (job.eq.'o') then we transform the occupied MO
|
|
C coefficients and if (job.eq.'v') then we read the virtual ones
|
|
character job
|
|
C
|
|
C Type of localization
|
|
character*4 loctype
|
|
character*8 orblocguess
|
|
C
|
|
C Size of marices, number of basis functions, core, occupied and virtual
|
|
C orbitals
|
|
integer nmatrix,nbasis,ncore,nocc,nos,nfroz ! HB
|
|
C
|
|
C MO coefficients (Cholesky)
|
|
real*8 mocoef_chol
|
|
C
|
|
C Matrix of traditional Boys transformation (initial guess)
|
|
real*8 uboysinit
|
|
C
|
|
C Matrix of own Boys transformation
|
|
real*8 uboys(nbasis,nbasis)
|
|
C
|
|
C MO coefficients (localized)
|
|
real*8 mocoef_lmo(nmatrix,*)
|
|
C
|
|
C Full matrix of the MO coefficients
|
|
real*8 mocoeffull(nbasis,*)
|
|
C
|
|
C Work array for dgemm subroutine
|
|
real*8 work
|
|
C
|
|
integer i,j
|
|
C
|
|
C
|
|
if ((loctype.eq."chol".or.loctype.eq."pm ".or.
|
|
$loctype.eq."ibo ".or.loctype.eq.'frag'.or.loctype.eq.'no '.or.
|
|
$loctype.eq.'spad').or.(orblocguess.eq.'read ') .or.
|
|
$loctype.eq.'pao') then
|
|
call dcopy(nbasis*nmatrix,mocoef_chol,1,mocoef_lmo,1)
|
|
else
|
|
c write(6,*) 'uboys x uboysinit x mocoef_chol',job
|
|
C These two lines are required if transform is commented out
|
|
call dgemm('n','n',nmatrix,nbasis,nmatrix,1.d0,uboysinit,
|
|
&nmatrix,mocoef_chol,nmatrix,0.d0,mocoef_lmo,nmatrix)
|
|
C These four lines are required if transform is active
|
|
c call dgemm('n','n',nmatrix,nmatrix,nmatrix,1.d0,uboys,nmatrix,
|
|
c & uboysinit,nmatrix,0.d0,work,nmatrix)
|
|
c call dgemm('n','n',nmatrix,nbasis,nmatrix,1.d0,work,nmatrix,
|
|
c & mocoef_chol,nmatrix,0.d0,mocoef_lmo,nmatrix)
|
|
endif
|
|
if(job.eq.'c') then
|
|
do i=1,nmatrix
|
|
do j=1,nbasis
|
|
CCC start HB
|
|
mocoeffull(j,i+nfroz)=mocoef_lmo(i,j)
|
|
CCC end HB
|
|
enddo
|
|
enddo
|
|
else if(job.eq.'o') then
|
|
do i=1,nmatrix
|
|
do j=1,nbasis
|
|
mocoeffull(j,i+ncore)=mocoef_lmo(i,j)
|
|
enddo
|
|
enddo
|
|
else if(job.eq.'s') then
|
|
do i=1,nmatrix
|
|
do j=1,nbasis
|
|
mocoeffull(j,i+ncore+nocc)=mocoef_lmo(i,j)
|
|
enddo
|
|
enddo
|
|
else if(job.eq.'v') then
|
|
do i=1,nmatrix
|
|
do j=1,nbasis
|
|
mocoeffull(j,i+ncore+nocc+nos)=mocoef_lmo(i,j)
|
|
enddo
|
|
enddo
|
|
else if(job.eq.'x') then
|
|
do i=1,nmatrix
|
|
do j=1,nbasis
|
|
mocoeffull(j,i)=mocoef_lmo(i,j)
|
|
enddo
|
|
enddo
|
|
endif
|
|
C
|
|
end
|
|
C
|
|
C***********************************************************************
|
|
C***********************************************************************
|
|
C Subroutine to calculate MO integrals (diploe and second moment) by
|
|
C transforming AO integrals with canonical MO coefficients
|
|
C***********************************************************************
|
|
C***********************************************************************
|
|
C
|
|
subroutine set_moints(
|
|
&nmatrix,
|
|
&nbasis,
|
|
&printfile,
|
|
&mocoef_chol,
|
|
&xdip,
|
|
&ydip,
|
|
&zdip,
|
|
&qmat,
|
|
&wread,
|
|
&wtransf,
|
|
&i4heap,
|
|
&r8heap,
|
|
&only_dip)
|
|
C
|
|
implicit none
|
|
C
|
|
C Size of matrices, number of basis functions
|
|
integer nmatrix,nbasis
|
|
C
|
|
C Variable of file handling
|
|
integer printfile
|
|
C
|
|
C MO coefficients from Cholesky decomposition
|
|
real*8 mocoef_chol
|
|
C
|
|
C Dipole and second moment matrices
|
|
real*8 xdip
|
|
real*8 ydip
|
|
real*8 zdip
|
|
real*8 qmat
|
|
C
|
|
C Arrays for reading coefficients
|
|
integer*4 i4heap(*)
|
|
real*8 r8heap(*)
|
|
C
|
|
C Work array for read_aoints subroutine
|
|
real*8 wread
|
|
C
|
|
C Work array for transf_integ_init subroutine
|
|
real*8 wtransf
|
|
C
|
|
logical only_dip
|
|
C
|
|
C Read first and second line of OEINT file
|
|
rewind(printfile)
|
|
read(printfile)
|
|
C
|
|
C Read full AO integral matrices and then transform them into MO
|
|
C integral matrices with the MO coefficients from the Cholesky
|
|
C decomposition
|
|
call read_aoints(nbasis,printfile,wread,r8heap,i4heap)
|
|
call transf_integ_init(nmatrix,nbasis,mocoef_chol,wread,xdip,
|
|
& wtransf)
|
|
call read_aoints(nbasis,printfile,wread,r8heap,i4heap)
|
|
call transf_integ_init(nmatrix,nbasis,mocoef_chol,wread,ydip,
|
|
& wtransf)
|
|
call read_aoints(nbasis,printfile,wread,r8heap,i4heap)
|
|
call transf_integ_init(nmatrix,nbasis,mocoef_chol,wread,zdip,
|
|
& wtransf)
|
|
c Ha nincs newton akkor a többi nem kell!!!
|
|
c return
|
|
C
|
|
if(only_dip) return
|
|
call read_aoints(nbasis,printfile,qmat,r8heap,i4heap)
|
|
call dcopy(nbasis**2,qmat,1,wread,1)
|
|
call read_aoints(nbasis,printfile,qmat,r8heap,i4heap)
|
|
call daxpy(nbasis**2,1.0d0,qmat,1,wread,1)
|
|
call read_aoints(nbasis,printfile,qmat,r8heap,i4heap)
|
|
call daxpy(nbasis**2,1.0d0,qmat,1,wread,1)
|
|
call transf_integ_init(nmatrix,nbasis,mocoef_chol,wread,qmat,
|
|
& wtransf)
|
|
C
|
|
end
|
|
C
|
|
C***********************************************************************
|
|
C***********************************************************************
|
|
C Subroutine to read full AO dipole moment and second moment
|
|
C integral matrices
|
|
C***********************************************************************
|
|
C***********************************************************************
|
|
C
|
|
subroutine read_aoints(
|
|
&nbasis,
|
|
&oeintfile,
|
|
&x,
|
|
&r8heap,
|
|
&i4heap)
|
|
|
|
implicit none
|
|
C
|
|
C Size of basis
|
|
integer nbasis
|
|
C
|
|
C Variable of file handling
|
|
integer oeintfile
|
|
C
|
|
C Full matrix to transform (dipole or second moment)
|
|
real*8 x(nbasis,*)
|
|
C
|
|
C Scratch arrays to write file
|
|
integer*4 i4heap(*)
|
|
real*8 r8heap(*)
|
|
C
|
|
integer i,j,nn,ii,jj
|
|
C
|
|
C
|
|
C%%% Empty matrix -> veglegesben sparse!
|
|
call dfillzero(x,nbasis**2)
|
|
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)
|
|
x(i,j)=r8heap(ii*2-1)
|
|
x(j,i)=r8heap(ii*2-1)
|
|
enddo
|
|
C
|
|
end
|
|
C
|
|
C***********************************************************************
|
|
C***********************************************************************
|
|
C Subroutine to transoform a matrix (a dipole moment or a second
|
|
C moment integral) for the initial guess
|
|
C***********************************************************************
|
|
C***********************************************************************
|
|
C
|
|
subroutine transf_integ_init(
|
|
&nmatrix,
|
|
&nbasis,
|
|
&x,
|
|
&inm,
|
|
&outm,
|
|
&work)
|
|
C
|
|
implicit none
|
|
C
|
|
C Size of matrices, number of basis functions
|
|
integer nbasis,nmatrix
|
|
C
|
|
C Matrix to transform a matrix with
|
|
real*8 x
|
|
C
|
|
C The full matrix to be transformated
|
|
real*8 inm
|
|
C
|
|
C The transformated matrix
|
|
real*8 outm
|
|
C
|
|
C Work array for dgemm subroutine
|
|
real*8 work
|
|
C
|
|
C
|
|
call dgemm('n','n',nmatrix,nbasis,nbasis,1.d0,x,nmatrix,inm,
|
|
& nbasis,0.d0,work,nmatrix)
|
|
call dgemm('n','t',nmatrix,nmatrix,nbasis,1.d0,work,nmatrix,x,
|
|
& nmatrix,0.d0,outm,nmatrix)
|
|
C
|
|
end
|
|
C
|
|
************************************************************************
|
|
subroutine loc_guess(imem,dcore,icmo,is,iover,iiold,ialpha,
|
|
&nbasis,nocc,scrfile3,mo,iout)
|
|
************************************************************************
|
|
* Initial guess for localization
|
|
************************************************************************
|
|
implicit none
|
|
integer imem,icmo,is,iover,iiold,ialpha,nbasis,nocc,scrfile3,iout
|
|
integer ierr
|
|
double precision dcore(*),mo(nbasis,nocc)
|
|
logical oldmo_exists
|
|
|
|
if(ialpha.le.1) then
|
|
inquire(file='OLDMO',exist=oldmo_exists)
|
|
if(oldmo_exists) open(scrfile3,file='OLDMO',form='unformatted')
|
|
else
|
|
inquire(file='OLDMOb',exist=oldmo_exists)
|
|
if(oldmo_exists) open(scrfile3,file='OLDMOb',form='unformatted')
|
|
endif
|
|
if(.not. oldmo_exists) return
|
|
call rtdmx(dcore(imem+nbasis*nocc),dcore(imem+nbasis*nocc),
|
|
$dcore(imem),scrfile3,nbasis,nocc)
|
|
close(scrfile3)
|
|
|
|
call dsymm('l','l',nbasis,nocc,1.d0,dcore(is),nbasis,
|
|
$dcore(imem),nbasis,0.d0,dcore(iover),nbasis)
|
|
call dgemm('t','n',nocc,nocc,nbasis,1.d0,mo,nbasis,
|
|
$dcore(iover),nbasis,0.d0,dcore(icmo),nocc)
|
|
call dsyrk('u','t',nocc,nocc,1.d0,dcore(icmo),nocc,0.d0,
|
|
$dcore(imem),nocc) !Mocskos Intel compiler
|
|
call filllo(dcore(imem),nocc)
|
|
call invsqrt(dcore(imem),nocc,dcore(imem+nbasis*nocc),iout,
|
|
$ierr,0.9d0,1)
|
|
if(ierr.eq.0) then
|
|
call dsymm('r','l',nocc,nocc,1.d0,dcore(imem),nocc,
|
|
$dcore(icmo),nocc,0.d0,dcore(imem+nbasis*nocc),nocc)
|
|
call dgemm('n','n',nbasis,nocc,nocc,1.d0,mo,nbasis,
|
|
$dcore(imem+nbasis*nocc),nocc,0.d0,dcore(imem),nbasis)
|
|
call dcopy(nbasis*nocc,dcore(imem),1,mo,1)
|
|
call dcopy(nocc*nocc,dcore(imem+nbasis*nocc),1,dcore(is),1)
|
|
endif
|
|
|
|
end subroutine
|
|
|