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

1774 lines
59 KiB
Fortran

module embed_grad
implicit none
integer nocc_a,nocc_b,nocc,nvirt ! number of occupied and virtual orbitals
integer ncore,nval ! number of core and valence orbitals
integer npos ! integral file stuff
integer imodip
integer iperm2ab,iperm2cv ! permutations to sort MOs to A/B and core-valence ordering
integer imo ! low level MOs (localized)
integer ifock2_ab ! low level Fock matrix
integer is ! AO overlap matrix
integer idiag ! diagonal for preconditioning
integer iumat,ieigval,ipiv ! other variables for preconditioning
double precision omega
double precision ll_chfx,ll_clrhfx,ll_csrhfx ! chfx for the low level calc
character*32 pcm,ll_dft ! pcm and low level dft
end module
!********************************************************************************
subroutine embed_grad_driver(scfdamp,itol,scftol)
!********************************************************************************
! Driver routine for solving the Z-vetor equations for embedded gradient
!********************************************************************************
use embed_grad, only: nocc_a,nocc_b,nocc,nvirt,ncore,nval,iperm2ab,iperm2cv, &
ll_dft,npos,ll_chfx,ll_clrhfx,ll_csrhfx,omega,pcm,ll_dft
use common_mod, only: dcore,icore,imem,nirmax,indexirrep_ptr,iintpos,nbasis, &
dfnbasis,iout,scrfile1
implicit none
integer i,imem_old,ialoc,itemp,ix,imo,ifock2_ab
integer ifock_embed,imo_embed,iy,ix_tilde,izmat,izloc
integer iamat,iatilde
integer offset(nirmax),sqroffset(nirmax),nc,ncorenew
data offset /1,1,1,1,1,1,1,1/
data sqroffset /1,1,1,1,1,1,1,1/
double precision exc
double precision devparr(2),itol,scftol
character*16 scfdamp
character*32 dft
integer dblalloc
double precision dnrm2
logical ll
imem_old=imem
devparr=1.0d0
! write(*,'(A,99I4)') 'perm ab: ',dcore(iperm2ab:iperm2ab+nocc-1)
! write(*,'(A,99I4)') 'perm cv: ',dcore(iperm2cv:iperm2cv+nocc-1)
ix_tilde=dblalloc(nocc_a*nocc_a) ! xtilde multiplier
iy=dblalloc(nocc_a*nocc_b) ! y multiplier
iamat=dblalloc(nbasis*nbasis) ! A matrix
imo_embed=dblalloc(nbasis*(nocc_a+nvirt))
imo=dblalloc(nbasis*nbasis)
ifock_embed=dblalloc(nbasis*nbasis)
! building xtilde and y multipliers
call read_from_disk('MOCOEF_AB ',dcore(imo),nbasis,.false.)
call read_from_disk('FOCK_EMBED ',dcore(ifock_embed),nbasis,.true.)
call read_subsys_mo(nbasis,nocc_a,nocc_b,nvirt,dcore(imo_embed))
call build_xtilde_multiplier(nbasis,nocc_a,dcore(imo_embed),&
dcore(ifock_embed),dcore(ix_tilde))
call build_y_multiplier(nbasis,nocc_a,nocc_b,dcore(ifock_embed),&
dcore(imo),dcore(imo_embed),dcore(iy))
call dbldealloc(ifock_embed)
! building the A matrix (constant terms in the Z-vector equations)
call dfillzero(dcore(iamat),nbasis*nbasis)
call build_amat_embed2(nbasis,nocc_a,nocc_b,nvirt,dcore(iy),dcore(imo),&
dcore(imo_embed),scfdamp,dcore(iamat))
! Solving Z-vector equations
write(iout,'(A)') ' Solving Z-vector equations'
! Solving the localization part
write(iout,'(A)') ' Solving localization part'
izloc=dblalloc(nocc*(nocc-1)/2)
ialoc=dblalloc(nbasis*nocc)
call dfillzero(dcore(izloc),nocc*(nocc-1)/2)
call dfillzero(dcore(ialoc),nocc*nbasis)
call solve_loc(dcore(iamat),dcore(izloc),dcore(ialoc))
call timer
write(iout,*)
! Solving Brillouin part
write(iout,'(A)') ' Solving Brillouin part'
izmat=dblalloc(nvirt*nocc)
call solve_bri(nbasis,dcore(izmat),dcore(iamat),dcore(ialoc),&
'FOCK2_AB ','MOCOEF_AB ',.true.)
call timer
write(iout,*)
itemp=dblalloc(nbasis*nbasis)
ifock2_ab=dblalloc(nbasis*nbasis) ! low level Fock with D^{AB}
call read_from_disk('FOCK2_AB ',dcore(ifock2_ab),nbasis,.true.)
call build_atilde_embed(nbasis,dfnbasis,nocc,nvirt,dcore(izmat), &
dcore(imo),dcore(ifock2_ab),npos,itemp,ll_chfx,ll_clrhfx,ll_csrhfx,&
omega,pcm,ll_dft)
iatilde=ifock2_ab
call dgemm('t','n',nbasis,nbasis,nbasis,1.0d0,dcore(imo),nbasis,&
dcore(itemp),nbasis,0.0d0,dcore(iatilde),nbasis)
ix=dblalloc(nbasis*nbasis)
call build_x_embed(nbasis,nocc,dcore(iamat),dcore(iatilde),dcore(ialoc),dcore(ix))
open(unit=scrfile1,file='LAGRANGE',form='unformatted')
if(ncore*(ncore-1)/2>0) write(scrfile1) dcore(izloc:izloc+ncore*(ncore-1)/2-1)
i=izloc+ncore*(ncore-1)/2
if(nval*(nval-1)/2>0) write(scrfile1) dcore(i:i+nval*(nval-1)-1)
write(scrfile1) dcore(izmat:izmat+nocc*nvirt-1)
i=izloc+ncore*(ncore-1)/2+nval*(nval-1)/2
if(ncore*nval>0) write(scrfile1) dcore(i:i+ncore*nval-1)
write(scrfile1) dcore(ix:ix+nbasis**2-1)
write(scrfile1) dcore(ix_tilde:ix_tilde+nocc_a**2-1)
write(scrfile1) dcore(iy:iy+nocc_a*nocc_b-1)
close(scrfile1)
write(iout,*) "Multipliers are written to the file 'LAGRANGE'"
call dbldealloc(imem_old)
end subroutine
!********************************************************************************
subroutine setup_mp2_vars(scftol,itol,scfdamp)
!********************************************************************************
! Initializing variables
!********************************************************************************
use embed_grad
use common_mod, only: dcore,icore,scfalg,indexirrep_ptr,iintpos,tedatfile, &
scrfile1,lsa,scftype,nbset,nbf,rs_ptr,imem,pepsilon,nir,inatrange,natoms, &
sqrsize,nfunc,nal,nbe,orbperir,fock_ptr,scfmaxit,sqrsize2,rmat_ptr,qmat_ptr,&
symtra_ptr,minpfile,nbasis
use mcscf, only: qscf_init
implicit none
double precision scftol,itol
character*16 scfdamp
integer n_intpos,i
integer dblalloc,intalloc
character*4 cscr4
character*8 qscf
symtra_ptr=imem
rs_ptr=imem
fock_ptr=imem
rmat_ptr=imem
qmat_ptr=imem
call getvar('nbasis ',nbasis)
call getvar('nbf ',nbf)
call getvar('nal ',nocc_a)
nir=1
lsa=.false.
scftype=1
sqrsize=nbasis*nbasis
sqrsize2=sqrsize
n_intpos=3*(nbasis+1)*nbasis/2
nfunc=0
nfunc(1)=nbasis
nocc_b=0
nal=nocc_a
nbe=nal
nocc=nocc_a+nocc_b
nvirt=nbasis-(nocc_a+nocc_b)
orbperir=0.0d0
orbperir(1)=nal
orbperir(1+nir)=nbe
indexirrep_ptr=intalloc(nbasis)
icore(indexirrep_ptr:indexirrep_ptr+nbasis-1)=1
iintpos=intalloc(0)
inatrange=intalloc(2*natoms*nbset)
call getvar('natrange ',icore(inatrange))
call getkey('scftol',6,cscr4,4)
read (cscr4, '(i4)') i
scftol=10.d0**dble(-i)
pepsilon=scftol
call getkey('itol',4,cscr4,4)
read(cscr4,*) i
itol=10.d0**(-i)
call getkey('scfmaxit',8,cscr4,4)
read (cscr4, '(i4)') scfmaxit
call getkey('pcm',3,pcm,32)
call getvar('omega ',omega)
call getkey('scfdamp',7,scfdamp,16)
! call getkey('scfalg',6,scfalg,8)
scfalg='direct '
open(tedatfile,file='TEDAT',form='UNFORMATTED')
! reading integral batch sizes from tedat file
if(trim(scfalg).eq.'disk') then
iintpos=intalloc(n_intpos)
rewind(tedatfile)
read(tedatfile) npos,icore(iintpos:iintpos+n_intpos-1)
endif
close(tedatfile)
call getkey('qscf',4,qscf,8)
call qscf_init(qscf)
nval=nocc_a+nocc_b-ncore
call getkey('dft',3,ll_dft,32)
call getvar('chfx ',ll_chfx)
call getvar('clrhfx ',ll_clrhfx)
call getvar('csrhfx ',ll_csrhfx)
end subroutine
!********************************************************************************
subroutine setup_embed_grad(scfdamp,itol,scftol)
!********************************************************************************
! Initializing variables
!********************************************************************************
use embed_grad
use common_mod, only: dcore,icore,scfalg,indexirrep_ptr,iintpos,tedatfile, &
scrfile1,lsa,scftype,nbset,nbf,rs_ptr,imem,pepsilon,nir,inatrange,natoms, &
sqrsize,nfunc,nal,nbe,orbperir,fock_ptr,scfmaxit,sqrsize2,rmat_ptr,qmat_ptr,&
symtra_ptr,minpfile,nbasis
use mcscf, only: qscf_init
implicit none
double precision scftol,itol
character*16 scfdamp
integer n_intpos,i,imo_loc,itemp,imo_ab
integer dblalloc,intalloc
character*4 cscr4
character*8 embed,qscf
symtra_ptr=imem
rs_ptr=imem
fock_ptr=imem
rmat_ptr=imem
qmat_ptr=imem
nir=1
lsa=.false.
scftype=1
sqrsize=nbasis*nbasis
sqrsize2=sqrsize
nfunc=0
nfunc(1)=nbasis
n_intpos=3*(nbasis+1)*nbasis/2
call getvar('natoms ',natoms)
call getvar('nbset ',nbset)
call getvar('nbf ',nbf)
call getvar('nal ',nocc_a)
call getvar('nfroz ',nocc_b)
! call getvar('ncore ',ncore)
open(scrfile1,file='NCORE')
read(scrfile1,'(I10)') ncore
close(scrfile1)
nal=nocc_a
nbe=nal
nocc_a=nocc_a-nocc_b
nocc=nocc_a+nocc_b
nvirt=nbasis-(nocc_a+nocc_b)
orbperir=0.0d0
orbperir(1)=nal
orbperir(1+nir)=nbe
indexirrep_ptr=intalloc(nbasis)
icore(indexirrep_ptr:indexirrep_ptr+nbasis-1)=1
iintpos=intalloc(0)
inatrange=intalloc(2*natoms*nbset)
call getvar('natrange ',icore(inatrange))
call getkey('scftol',6,cscr4,4)
read (cscr4, '(i4)') i
scftol=10.d0**dble(-i)
pepsilon=scftol
call getkey('itol',4,cscr4,4)
read(cscr4,*) i
itol=10.d0**(-i)
call getkey('scfmaxit',8,cscr4,4)
read (cscr4, '(i4)') scfmaxit
call getkey('pcm',3,pcm,32)
call getvar('omega ',omega)
call getkey('scfdamp',7,scfdamp,16)
call getkey('scfalg',6,scfalg,8)
open(tedatfile,file='TEDAT',form='UNFORMATTED')
! reading integral batch sizes from tedat file
if(trim(scfalg).eq.'disk') then
iintpos=intalloc(n_intpos)
rewind(tedatfile)
read(tedatfile) npos,icore(iintpos:iintpos+n_intpos-1)
endif
close(tedatfile)
open(scrfile1,file='CHFX')
read(scrfile1,'(ES30.16)') ll_chfx
read(scrfile1,'(2ES30.16)') ll_clrhfx,ll_csrhfx
close(scrfile1)
open(minpfile,file='MINP',status='old')
call getkeym('embed',5,embed,8)
read(minpfile,*)
read(minpfile,*) ll_dft
ll_dft=adjustl(ll_dft)
call lowercase(ll_dft,ll_dft,32)
call removed3(ll_dft)
if(trim(ll_dft).eq.'hf') ll_dft='off '
close(minpfile)
call getkey('qscf',4,qscf,8)
call qscf_init(qscf)
nval=nocc_a+nocc_b-ncore
write(*,*) 'nval: ',nval
iperm2ab=intalloc(nal)
iperm2cv=intalloc(nal)
imo_loc=dblalloc(nbasis*nbasis)
imo_ab=dblalloc(nbasis**2)
itemp=dblalloc(nbasis*nbasis)
call read_from_disk('MOCOEF.LOC ',dcore(imo_loc),nbasis,.false.)
call read_from_disk('OEINT ',dcore(itemp),nbasis,.true.)
call read_from_disk('MOCOEF_AB ',dcore(imo_ab),nbasis,.false.)
call get_perm(nbasis,nocc_a+nocc_b,dcore(imo_loc),dcore(imo_ab),&
dcore(itemp),icore(iperm2ab),icore(iperm2cv))
call dbldealloc(imo_loc)
end subroutine
!********************************************************************************
subroutine get_perm(nbasis,nocc,mo_cv,mo_ab,s,perm2ab,perm2cv)
!********************************************************************************
!********************************************************************************
use common_mod, only: dcore,imem
implicit none
integer nbasis,nocc
integer perm2ab(nocc),perm2cv(nocc)
double precision mo_cv(nbasis,nocc),mo_ab(nbasis,nocc)
double precision s(nbasis,nbasis)
integer i,j,ioverlap,itemp
integer dblalloc
ioverlap=dblalloc(nocc*nocc)
call get_mo_overlap(nbasis,nocc,mo_cv,mo_ab,s,dcore(ioverlap))
do i=1,nocc
j=1
do while(dabs(dcore(ioverlap+j-1+(i-1)*nocc)) < 0.1d0)
j=j+1
enddo
perm2ab(j)=i
perm2cv(i)=j
enddo
call dbldealloc(ioverlap)
end subroutine
!********************************************************************************
subroutine calc_diff_dens(nbasis,nocc_a,mo_a,mo_embed_a,diff_dens)
!********************************************************************************
! Calcualtes the difference of the high level and the low level SCF densities
!********************************************************************************
implicit none
integer nocc_a,nbasis
double precision mo_a(nbasis,nocc_a),mo_embed_a(nbasis,nocc_a)
double precision diff_dens(nbasis,nbasis)
call dsyrk('u','n',nbasis,nocc_a,0.5d0,mo_embed_a,nbasis,0.0d0,&
diff_dens,nbasis)
call dsyrk('u','n',nbasis,nocc_a,-0.5d0,mo_a,nbasis,1.0d0, &
diff_dens,nbasis)
call filllo(diff_dens,nbasis)
end subroutine
!********************************************************************************
subroutine build_y_multiplier(nbasis,nocc_a,nocc_b,embed_fock,mo_b,embed_mo,y)
!********************************************************************************
! Calculates the Lagrange multiplier of the orthogonality condition of the
! two subsystems
!********************************************************************************
use common_mod, only: dcore
implicit none
integer nbasis,nocc_a,nocc_b
double precision embed_fock(nbasis,nbasis)
double precision embed_mo(nbasis,nocc_a),mo_b(nbasis,nocc_b)
double precision y(nocc_a,nocc_b)
integer iw
integer dblalloc
iw=dblalloc(nbasis*nocc_b)
call dsymm('l','u',nbasis,nocc_b,1.0d0,embed_fock,nbasis,mo_b,nbasis,&
0.0d0,dcore(iw),nbasis)
call dgemm('t','n',nocc_a,nocc_b,nbasis,-4.0d0,embed_mo,nbasis, &
dcore(iw),nbasis,0.0d0,y,nocc_a)
call dbldealloc(iw)
end subroutine
!********************************************************************************
subroutine build_xtilde_multiplier(nbasis,nocc_a,mo_embed,fock_embed,x)
!********************************************************************************
! Calculates the Lagrange multiplier of the orthogonality condition of the
! embedded subsystems
!********************************************************************************
use common_mod, only: dcore
implicit none
integer nbasis,nocc_a
double precision mo_embed(nbasis,nocc_a),fock_embed(nbasis,nbasis)
double precision x(nocc_a,nocc_a)
integer iw
integer dblalloc
iw=dblalloc(nbasis*nocc_a)
call dsymm('l','u',nbasis,nocc_a,1.0d0,fock_embed,nbasis,mo_embed,nbasis,&
0.0d0,dcore(iw),nbasis)
call dgemm('t','n',nocc_a,nocc_a,nbasis,-2.0d0,mo_embed,nbasis,dcore(iw),&
nbasis,0.0d0,x,nocc_a)
call dbldealloc(iw)
end subroutine
!********************************************************************************
subroutine get_dipole_cv(nbasis,ncore,nval,dip,modip,mo_core,mo_val,w)
!********************************************************************************
! Get the AO and MO dipole moment integrals in case of core-valence separation
!********************************************************************************
implicit none
integer nbasis,ncore,nval
double precision mo_core(nbasis,ncore),mo_val(nbasis,nval)
double precision dip(nbasis,nbasis,3),modip(ncore**2*3+nval**2*3)
double precision w(*)
integer xyz
integer, parameter :: printfile=24
open(printfile,file='PRINT',form='UNFORMATTED')
read(printfile) !kinetic energy integrals
do xyz=1,3
call read_aoints(nbasis,printfile,dip(1,1,xyz),w,w)
if(ncore > 0) then
call dsymm('l','u',nbasis,ncore,1.0d0,dip(1,1,xyz),nbasis,mo_core,nbasis, &
0.0d0,w,nbasis)
call dgemm('t','n',ncore,ncore,nbasis,1.0d0,mo_core,nbasis,w,nbasis,0.0d0,&
modip(1+(xyz-1)*ncore**2),ncore)
endif
if(nval > 0) then
call dsymm('l','u',nbasis,nval,1.0d0,dip(1,1,xyz),nbasis,mo_val,nbasis, &
0.0d0,w,nbasis)
call dgemm('t','n',nval,nval,nbasis,1.0d0,mo_val,nbasis,w,nbasis,0.0d0, &
modip(1+3*ncore**2+(xyz-1)*nval**2),nval)
endif
enddo
close(printfile)
end subroutine
!********************************************************************************
subroutine aloc_loc(ncore,nval,nocc,nbasis,zloc,aloc)
!********************************************************************************
! Calcualting a(zloc) matrix for core and val orbitals (first index in AO basis)
! The occupied indices are separated to core-val orbitals.
! TODO delete!!
!********************************************************************************
use common_mod, only: dcore,imem,scfalg,dfnbasis
implicit none
integer ncore,nval,nocc,nbasis
double precision zloc(nocc*(nocc-1)/2),aloc(nbasis,nocc)
integer idip,imem_old,i,izfull,itemp,n,imo,imodip
integer dblalloc
imem_old=imem
n=max(ncore,nval)
imo=dblalloc(nbasis*nbasis)
call read_from_disk('MOCOEF.LOC ',dcore(imo),nbasis,.false.)
! Getting dipole integrals
idip=dblalloc(3*nbasis**2)
imodip=dblalloc(3*ncore**2+3*nval**2)
call get_dipole_cv(nbasis,ncore,nval,dcore(idip),dcore(imodip),&
dcore(imo),dcore(imo+nbasis*ncore),dcore(imem))
! Build a(zloc)
izfull=dblalloc(n*n)
itemp=dblalloc(max(nbasis**2,2*nbasis*n))
if(ncore*(ncore-1)/2 > 0) then
call build_aloc(nbasis,ncore,ncore,dcore(idip),dcore(imo),zloc,aloc,&
dcore(izfull),dcore(itemp),dcore(itemp+nbasis*n),dcore(imodip), &
'rhf ')
endif
if(nval*(nval-1)/2 > 0) then
i=1+ncore*(ncore-1)/2
call build_aloc(nbasis,nval,nval,dcore(idip),dcore(imo+nbasis*ncore),&
zloc(i),aloc(1,ncore+1),dcore(izfull),dcore(itemp),&
dcore(itemp+nbasis*n),dcore(imodip+3*ncore**2),'rhf ')
endif
call dbldealloc(imem_old)
end subroutine
!********************************************************************************
subroutine build_aloc_embed(ncore,nval,nocc,nbasis,zloc,aloc,perm2ab)
!********************************************************************************
! Calcualting a(zloc) matrix (first index in MO basis)
! The occupied indices are separated to A/B subsystems.
!********************************************************************************
use common_mod, only: dcore,imem,scfalg,dfnbasis
use embed_grad, only: ll_chfx,ll_csrhfx,ll_clrhfx,omega,pcm,ll_dft
implicit none
integer ncore,nval,nocc,nbasis,perm2ab(nocc)
double precision zloc(nocc*(nocc-1)/2),aloc(nbasis,nocc)
integer idip,imem_old,i,izfull,itemp
integer n,imo,imodip,is
character*16 fock_file
character*32 dft
integer dblalloc
imem_old=imem
n=max(ncore,nval)
fock_file='FOCK2_AB '
imo=dblalloc(nbasis*nbasis)
call read_from_disk('MOCOEF.LOC ',dcore(imo),nbasis,.false.)
! Getting dipole integrals
idip=dblalloc(3*nbasis**2)
imodip=dblalloc(3*ncore**2+3*nval**2)
call get_dipole_cv(nbasis,ncore,nval,dcore(idip),dcore(imodip),&
dcore(imo),dcore(imo+nbasis*ncore),dcore(imem))
! Build a(zloc)
izfull=dblalloc(n*n)
itemp=dblalloc(max(nbasis**2,2*nbasis*n))
if(ncore*(ncore-1)/2 > 0) then
call build_aloc(nbasis,ncore,ncore,dcore(idip),dcore(imo),zloc,aloc,&
dcore(izfull),dcore(itemp),dcore(itemp+nbasis*n),dcore(imodip), &
'rhf ')
endif
if(nval*(nval-1)/2 > 0) then
i=1+ncore*(ncore-1)/2
call build_aloc(nbasis,nval,nval,dcore(idip),dcore(imo+nbasis*ncore),&
zloc(i),aloc(1,ncore+1),dcore(izfull),dcore(itemp),&
dcore(itemp+nbasis*n),dcore(imodip+3*ncore**2),'rhf ')
endif
call dbldealloc(idip)
if(ncore*nval > 0) then
i=1+ncore*(ncore-1)/2+nval*(nval-1)/2
dft=ll_dft
if(trim(dft)=='user') dft='ll_user '
call aloc_core_val(nbasis,dfnbasis,ncore,nval,zloc(i),dcore(imo),&
scfalg,ll_chfx,ll_csrhfx,ll_clrhfx,omega,pcm,dft,aloc,fock_file)
endif
! transforming to MO basis
is=dblalloc(nbasis*nocc)
call read_from_disk('MOCOEF_AB ',dcore(imo),nbasis,.false.)
call dgemm('t','n',nbasis,nocc,nbasis,1.0d0,dcore(imo),nbasis,aloc,&
nbasis,0.0d0,dcore(is),nbasis)
! reordering MOs from core-valence to A-B subsystem
do i=1,nocc
aloc(:,perm2ab(i))=dcore(is+(i-1)*nbasis:is+i*nbasis-1)
enddo
call dbldealloc(imem_old)
end subroutine
!********************************************************************************
subroutine aloc_core_val(nbasis,dfnbasis,ncore,nval,zloc,mo,scfalg,chfx,csrhfx, &
clrhfx,omega,pcm,dft,aloc,fock_file)
!********************************************************************************
! Building the core-val separation contribution to the a(zloc) matrix
! First index is in AO basis
!********************************************************************************
use common_mod, only: dcore,imem,nirmax,orbperir
use mcscf, only: ldf,symm_dens,build_mmat_
use embed_grad, only: npos
implicit none
integer nbasis,ncore,nval,dfnbasis
double precision zloc(ncore,nval),mo(nbasis,ncore+nval)
double precision aloc(nbasis,ncore+nval)
double precision chfx,csrhfx,clrhfx,omega
character*8 scfalg
character*16 fock_file
integer nocca,noccb,ic,ic2(2),idens(2),im(2),itmp1,itmp2,imem_old,nn
integer offset(nirmax),sqroffset(nirmax),nocc,nc,ncorenew,i,j,lfin
double precision alpha,beta,delta,exc,devparr(2),x
character*16 scfdamp
character*32 pcm,dft,dvxc_dft
integer dblalloc
equivalence(ic2(1),ic)
data offset /1,1,1,1,1,1,1,1/
data sqroffset /1,1,1,1,1,1,1,1/
interface
subroutine setup_mmat_vars(nbasis,n1,n2,dens,c,z,mo1,mo2,scfalg)
implicit none
integer n1,n2,nbasis
double precision dens(nbasis,nbasis),c(2*nbasis*min(n1,n2))
double precision, target :: mo1(nbasis,n1),mo2(nbasis,n2)
double precision z(n1,n2)
character*8 scfalg
end subroutine
end interface
imem_old=imem
nocca=min(nval,ncore)
noccb=nocca
nocc=ncore+nval
ncorenew=0
scfdamp='off '
devparr=1.0d0
lfin=0
! Build modified MOs for build_mmat
idens(1)=dblalloc(nbasis**2)
idens(2)=idens(1)
ic=dblalloc(2*ncore*nbasis)
ic2(2)=ic
call setup_mmat_vars(nbasis,ncore,nval,dcore(idens(1)),dcore(ic),&
zloc,mo,mo(1,1+ncore),scfalg)
! Building M matrix
im(1)=dblalloc(nbasis**2)
im(2)=im(1)
! fock_build calculates the number of occupied orbitals from orbperir :(
nn=orbperir(1)
orbperir(1)=nocca
call build_mmat_(offset,sqroffset,idens,im,npos,zloc,mo,chfx,100,dfnbasis, &
ncore,ncore,lfin,exc,'off ',.false.,scfdamp,' ',0,nc,ncorenew,&
devparr,clrhfx,csrhfx,omega,.false.,pcm,nocca,noccb,ic2)
orbperir(1)=nn
! Transforming to MO basis
call dsymm('l','u',nbasis,nocc,-1.0d0,dcore(im(1)),nbasis,mo,nbasis,1.0d0,&
aloc,nbasis)
! Contribution from the Fock matrix
itmp1=dblalloc(nbasis*nbasis)
call read_from_disk(fock_file,dcore(im(1)),nbasis,.true.)
call dgemm('n','n',nbasis,nocc,nbasis,1.0d0,dcore(im(1)),nbasis,mo,nbasis, &
0.0d0,dcore(itmp1),nbasis)
call dgemm('n','t',nbasis,ncore,nval,1.0d0,dcore(itmp1+ncore*nbasis),nbasis,&
zloc,ncore,1.0d0,aloc,nbasis)
call dgemm('n','n',nbasis,nval,ncore,1.0d0,dcore(itmp1),nbasis,zloc,ncore, &
1.0d0,aloc(1,1+ncore),nbasis)
! DFT contribution
if(trim(dft)/='off') then
call dgemm('n','t',ncore,nbasis,nval,1.0d0,zloc,ncore,mo(1,1+ncore),&
nbasis,0.0d0,dcore(itmp1),ncore)
call dgemm('n','n',nbasis,nbasis,ncore,1.0d0,mo,nbasis,dcore(itmp1),&
ncore,0.0d0,dcore(im(1)),nbasis)
call symmat(dcore(im(1)),nbasis)
dvxc_dft=dft
! if(trim(dft)=='user') dvxc_dft='ll_user '
call amat_embed_dft2(nbasis,nocc,dcore(im(1)),mo,dcore(itmp1),dvxc_dft)
call dsymm('l','u',nbasis,nocc,1.0d0,dcore(itmp1),nbasis,mo,nbasis,&
1.0d0,aloc,nbasis)
endif
call dbldealloc(imem_old)
end subroutine
!********************************************************************************
subroutine build_amat_embed2(nbasis,nocc_a,nocc_b,nvirt,y,mo,mo_embed,scfdamp,amat)
!********************************************************************************
! Builds the A matrix for embedding calculaions (only SCF part)
!********************************************************************************
use common_mod, only: dcore,icore,imem,indexirrep_ptr,nirmax,iintpos
use embed_grad, only: npos,ll_chfx,ll_clrhfx,ll_csrhfx,omega,pcm,ll_dft
implicit none
integer nbasis,nocc_a,nocc_b,nvirt
double precision mo(nbasis,nbasis),mo_embed(nbasis,nocc_a+nvirt)
double precision amat(nbasis,nocc_a+nocc_b),y(nocc_a*nocc_b)
character*16 scfdamp
integer ig2_embed,idens_a,ifock2_ab,ifock2_a,if2_xc_a,is
integer idiff_dens,ia_ao
integer lfin,offset(nirmax),sqroffset(nirmax),nc,ncorenew
integer dblalloc
double precision exc,devparr(2)
character*32 dft
logical ll
data offset /1,1,1,1,1,1,1,1/
data sqroffset /1,1,1,1,1,1,1,1/
! 1st step: calculation of g2[D^Atilde]
dft='off '
lfin=0
devparr=1.0d0
ig2_embed=dblalloc(nbasis*nbasis)
idens_a=dblalloc(nbasis*nbasis)
call dsyrk('u','n',nbasis,nocc_a,2.0d0,mo_embed,nbasis,0.0d0,&
dcore(idens_a),nbasis)
call filllo(dcore(idens_a),nbasis)
call fock_build(dcore(idens_a),dcore(idens_a),dcore(ig2_embed), &
dcore(ig2_embed),dcore,icore(indexirrep_ptr),offset,sqroffset,npos,&
icore(iintpos),2,mo_embed,lfin,exc,dft,ll_chfx,.false.,'off ', &
scfdamp,' ',0,nc,ncorenew,devparr,dcore,0,ll_clrhfx,ll_csrhfx, &
omega,dcore,.false.,1,1,pcm,dcore,dcore,.false.)
call dbldealloc(idens_a)
! 2nd step: reading variables from disk
ifock2_ab=dblalloc(nbasis*nbasis) ! low level Fock with D^{AB}
ifock2_a=dblalloc(nbasis*nbasis) ! low level Fock with D^A
if2_xc_a=dblalloc(nbasis*nbasis) ! low level XC matrix with D^A
idiff_dens=dblalloc(nbasis*nbasis) ! D^Atilde - D^A
is=dblalloc(nbasis*nbasis)
call read_from_disk('FOCK2_AB ',dcore(ifock2_ab),nbasis,.true.)
call read_from_disk('FOCK2_A ',dcore(ifock2_a ),nbasis,.true.)
! inquire(file='F2_XC_A',exist=ll)
! if(ll) then
write(*,*) 'll_dft: '//ll_dft
if(trim(ll_dft).ne.'off') then
call read_from_disk('F2_XC_A ',dcore(if2_xc_a),nbasis,.true.)
else
call dfillzero(dcore(if2_xc_a),nbasis**2)
endif
call read_from_disk('OEINT ',dcore(is),nbasis,.true.)
call calc_diff_dens(nbasis,nocc_a,mo(1,1+nocc_b),mo_embed,dcore(idiff_dens))
! 3rd step: calculating A
ia_ao=dblalloc(nbasis*nbasis)
call build_amat_embed(nbasis,nocc_a,nocc_b,mo,dcore(ifock2_ab), &
dcore(ifock2_a),dcore(if2_xc_a),dcore(ig2_embed),dcore(idiff_dens),&
dcore(ia_ao),ll_dft)
! Contribution from the orthogonality of the two subsystems
call dsymm('l','u',nbasis,nocc_a,0.25d0,dcore(is),nbasis,mo_embed,nbasis,&
0.0d0,dcore(idens_a),nbasis)
call dgemm('n','n',nbasis,nocc_b,nocc_a,1.0d0,dcore(idens_a),nbasis,y, &
nocc_a,1.0d0,dcore(ia_ao),nbasis)
! transforming to MO basis
call dgemm('t','n',nbasis,nbasis,nbasis,4.0d0,mo,nbasis,dcore(ia_ao),&
nbasis,0.0d0,amat,nbasis)
call dbldealloc(ig2_embed)
end subroutine
!********************************************************************************
subroutine build_amat_embed(nbasis,nocc_a,nocc_b,mo_full,fock2_ab,fock2_a, &
f2_xc_a,g2_embed,diff_dens,amat,dft)
!********************************************************************************
! Calculates the constant terms in the derivative of the Langrangian
!********************************************************************************
use common_mod, only: dcore,imem,scrfile1
implicit none
integer nbasis,nocc_a,nocc_b
double precision mo_full(nbasis,nbasis) ! [mo_b,mo_a,mo_virt]
double precision fock2_ab(nbasis,nbasis),fock2_a(nbasis,nbasis)
double precision f2_xc_a(nbasis,nbasis),g2_embed(nbasis,nbasis)
double precision diff_dens(nbasis,nbasis)
double precision amat(nbasis,nbasis)
character*32 dft
integer imem_old,ifock1_a,ixc_ab,ixc_a,iha
integer nbasis2
character*32 dvxc_dft
logical ldft
integer dblalloc
imem_old=imem
nbasis2=nbasis*nbasis
ldft=trim(dft)/='off'
dvxc_dft=dft
if(trim(dft)=='user') dvxc_dft='ll_user '
call dfillzero(amat,nbasis2)
ixc_ab=dblalloc(nbasis2)
call dfillzero(dcore(ixc_ab),nbasis2)
if(ldft) then
ixc_a =dblalloc(nbasis2)
call dfillzero(dcore(ixc_a),nbasis2)
! d^2E^XC[D^AB]/d^2D^AB
call amat_embed_dft2(nbasis,nocc_a+nocc_b,diff_dens,mo_full, &
dcore(ixc_ab),dvxc_dft)
! d^2E^XC[D^A]/d^2D^A
call amat_embed_dft2(nbasis,nocc_a,diff_dens,mo_full(1,1+nocc_b),&
dcore(ixc_a),dvxc_dft)
call dscal(nbasis2,-1.0d0,dcore(ixc_a),1)
call daxpy(nbasis2,1.0d0,dcore(ixc_ab),1,dcore(ixc_a),1)
call dsymm('l','u',nbasis,nocc_a,1.0d0,dcore(ixc_a),nbasis,mo_full(1,1+nocc_b),&
nbasis,1.0d0,amat(1,1+nocc_b),nbasis)
call dbldealloc(ixc_a)
endif
! F_2[D^AB]+d^2E^XC[D^AB]/d^2D^AB*(D^Atilde-D^A)+G_2[D^Atilde]-G_2[D^A]
call daxpy(nbasis2,1.0d0,fock2_ab,1,dcore(ixc_ab),1)
call daxpy(nbasis2,1.0d0,g2_embed,1,dcore(ixc_ab),1)
! -G_2[D^A]=-(F_2[D^A]-h^A-F_2^XC[D^A])=-F_2[D^A]+h^A+F_2^XC[D^A]
call daxpy(nbasis2,-1.0d0,fock2_a,1,dcore(ixc_ab),1)
if(ldft) call daxpy(nbasis2, 1.0d0,f2_xc_a,1,dcore(ixc_ab),1)
iha=dblalloc(nbasis2)
open(unit=scrfile1,file='FOCK2_A',form='unformatted')
read(scrfile1)
call roeint(dcore(imem),dcore(imem),dcore(iha),scrfile1,nbasis)
close(scrfile1)
call daxpy(nbasis2,1.0d0,dcore(iha),1,dcore(ixc_ab),1)
call dsymm('l','u',nbasis,nocc_b,1.0d0,dcore(ixc_ab),nbasis,mo_full,&
nbasis,1.0d0,amat,nbasis)
call dbldealloc(imem_old)
end subroutine
!********************************************************************************
subroutine amat_embed_dft2(nbasis,nocc,diff_dens,mo,res,dft)
!********************************************************************************
! Contribution to the A matrix from the second derivative of the XC energy
!********************************************************************************
use common_mod, only: dcore,imem,iout,ifltln,minpfile,maxcor
implicit none
integer nbasis,nocc
double precision diff_dens(nbasis,nbasis),mo(nbasis,nocc)
double precision res(nbasis,nbasis)
character*32 dft
integer imo,itemp,imem_old
integer grfile,iscftype,multiplicity
parameter(grfile = 15)
double precision exc
character*4 mult
integer dblalloc
imem_old=imem
imo=dblalloc(nocc*nbasis)
iscftype=1
call getkey('mult', 4, mult, 4)
read(mult,*) multiplicity
call dvxc(nbasis,nocc,nocc,res,res,dcore(imo),dcore(imo),grfile,dcore,&
iout,exc,dft,minpfile,iscftype,ifltln,maxcor,imem,dcore,0,2,dcore, &
'vxcd',diff_dens,diff_dens,mo,mo,dcore,dcore,dcore,dcore,200,0,1,1, &
multiplicity)
call dscal(nbasis*nbasis,0.5d0,diff_dens,1) ! dvxc multiplies the density with 2
call dbldealloc(imem_old)
end subroutine
!********************************************************************************
subroutine amat_e12(nbasis,nocc_a,nocc_b,nvirt,fock2_ab,fock2_a,fock1_a,mo_a, &
mo_b,chat,amat)
!********************************************************************************
! E12 contribution to A
!********************************************************************************
use common_mod, only: dcore
implicit none
integer nbasis,nocc_a,nocc_b,nvirt
double precision mo_a(nbasis,nocc_a),mo_b(nbasis,nocc_b)
double precision fock2_ab(nbasis,nbasis),fock2_a(nbasis,nbasis)
double precision fock1_a(nbasis,nbasis),chat(nbasis,nocc_a+nvirt)
double precision amat(nbasis,nbasis)
integer ifock2_b
integer dblalloc
ifock2_b=dblalloc(nbasis*nbasis)
! F_2[D^{AB}]_{mu,p}
call dsymm('l','u',nbasis,nocc_b,1.0d0,fock2_ab,nbasis,mo_b,nbasis,&
1.0d0,amat,nbasis)
! ! F_1[D^Atilde]_{mu p}, p: subsystem MO or virtual
! call dsymm('l','u',nbasis,nocc_a+nvirt,1.0d0,fock1_a,nbasis,chat,nbasis,&
! 1.0d0,amat(1,1+nocc_b),nbasis)
! F_2[D^B]=F_2[D^{AB}]-F_2[D^A]
call dcopy(nbasis*nbasis,fock2_ab,1,dcore(ifock2_b),1)
call daxpy(nbasis*nbasis,-1.0d0,fock2_a,1,dcore(ifock2_b),1)
call dsymm('l','u',nbasis,nocc_a,1.0d0,dcore(ifock2_b),nbasis,mo_a,&
nbasis,1.0d0,amat(1,1+nocc_b),nbasis)
call dbldealloc(ifock2_b)
end subroutine
!********************************************************************************
subroutine get_mo_overlap(nbasis,nocc,mo_low,mo_high,s,umat)
!********************************************************************************
! Calculates MO overlap matrix
!********************************************************************************
use common_mod, only: dcore
implicit none
integer nbasis,nocc
double precision mo_low(nbasis,nocc),mo_high(nbasis,nocc)
double precision s(nbasis,nbasis),umat(nbasis,nbasis)
integer itemp
integer dblalloc
itemp=dblalloc(nocc*nbasis)
call dsymm('l','u',nbasis,nocc,1.0d0,s,nbasis,mo_high,nbasis,0.0d0, &
dcore(itemp),nbasis)
call dgemm('t','n',nocc,nocc,nbasis,1.0d0,mo_low,nbasis,dcore(itemp),&
nbasis,0.0d0,umat,nocc)
call dbldealloc(itemp)
end subroutine
!********************************************************************************
subroutine solve_loc(amat,zloc,aloc)
!********************************************************************************
! Solve Z-vector equations for zloc
!********************************************************************************
use common_mod, only: dcore,icore,imem,nbasis,scrfile1,iout,dfnbasis,scfalg
use embed_grad, only: nocc,ncore,nval,imodip,imo,iperm2ab,iperm2cv,idiag, &
omega,ll_chfx,ll_clrhfx,ll_csrhfx,pcm,ll_dft
implicit none
double precision amat(nbasis,nocc),zloc(ncore*(ncore-1)/2+nval*(nval-1)/2)
double precision aloc(nbasis,nocc)
integer imem_old,irhs,itemp,idipole,imo_can,ialoc_ao,ia_occ
integer ndim,maxit,i,j,n,info
double precision tol
character*4 c4
character*8 c8
character*16 vecfile,errfile,fock_file
character*32 dft
integer dblalloc,idamax
interface
subroutine multiply_loc(ndim,zloc,res)
integer ndim
double precision zloc(ndim),res(ndim)
end subroutine
subroutine precond_with_diag(ndim,x,y)
integer ndim
double precision x(ndim),y(ndim)
end subroutine
double precision function norm(ndim,x)
integer ndim
double precision x(ndim)
end function
end interface
imem_old=imem
vecfile='fort.18 '
errfile='fort.19 '
open(scrfile1,file=trim(vecfile),form='unformatted')
close(scrfile1,status='delete')
open(scrfile1,file=trim(errfile),form='unformatted')
close(scrfile1,status='delete')
call getkey('scfmaxit',8,c4,4)
read(c4,*) maxit
call getkey('ldfgrad_tol',11,c8,8)
read(c8,*) i
tol=10.d0**dble(-i)
fock_file='FOCK2_AB '
! sorting occ-occ block of the A matrix to core-val order
ia_occ=dblalloc(nocc*nocc)
itemp=dblalloc(nocc*nocc)
call dlacpy('f',nocc,nocc,amat,nbasis,dcore(ia_occ),nocc)
call sort_index(nocc,nocc,dcore(ia_occ),dcore(itemp),icore(iperm2cv))
call dbldealloc(itemp)
! Reading localized MOs and transform dipole moment integrals to MO basis
imo=dblalloc(nbasis*nbasis)
imodip=dblalloc(ncore*ncore*3+nval*nval*3)
idipole=dblalloc(3*nbasis**2)
call read_from_disk('MOCOEF.LOC ',dcore(imo),nbasis,.false.)
call get_dipole_cv(nbasis,ncore,nval,dcore(idipole),dcore(imodip),&
dcore(imo),dcore(imo+nbasis*ncore),dcore(imem))
call dbldealloc(idipole)
! Calculating r.h.s. of the equation
ndim=ncore*(ncore-1)/2+nval*(nval-1)/2
irhs=dblalloc(ndim)
call rhs_loc(nbasis,ncore,nval,nocc,ndim,icore(ia_occ),dcore(irhs))
! Prepare diagonal for preconditioner
n=ndim
idiag=dblalloc(n)
call build_diag(ncore,dcore(imodip),dcore(idiag))
call build_diag(nval,dcore(imodip+3*ncore**2),dcore(idiag+ncore*(ncore-1)/2))
! Solving the equation
call general_subspace_solver(ndim,zloc,dcore(irhs),maxit,tol,.false.,&
vecfile,errfile,norm,multiply_loc,precond_with_diag,.true.)
! Building the a(zloc) matrix for core and val part
ialoc_ao=dblalloc(nbasis*nocc)
call dfillzero(dcore(ialoc_ao),nbasis*nocc)
call aloc_loc(ncore,nval,nocc,nbasis,zloc,dcore(ialoc_ao))
call dgemm('t','n',nbasis,nocc,nbasis,1.0d0,dcore(imo),nbasis,dcore(ialoc_ao),&
nbasis,0.0d0,aloc,nbasis)
! call dbldealloc(ialoc_ao)
! Solving the core-val part
if(ncore*nval > 0) then
call dbldealloc(imodip)
imo_can=dblalloc(nbasis**2)
call read_from_disk('MOCOEF_CAN ',dcore(imo_can),nbasis,.false.)
call solve_core_val(nbasis,nocc,ncore,nval,zloc(1+ndim),dcore(ia_occ),&
dcore(imo),dcore(imo_can),aloc)
!updating the a(zloc) matrix
ialoc_ao=imo_can
call dfillzero(dcore(ialoc_ao),nbasis*nbasis)
dft=ll_dft
if(trim(dft)=='user') dft='ll_user '
call aloc_core_val(nbasis,dfnbasis,ncore,nval,zloc(1+ndim),dcore(imo),&
scfalg,ll_chfx,ll_csrhfx,ll_clrhfx,omega,pcm,dft,dcore(ialoc_ao), &
fock_file)
call dgemm('t','n',nbasis,nocc,nbasis,1.0d0,dcore(imo),nbasis,&
dcore(ialoc_ao),nbasis,1.0d0,aloc,nbasis)
endif
call sort_index(nbasis,nocc,aloc,dcore(imo),dcore(iperm2ab))
call dbldealloc(imem_old)
end subroutine
!********************************************************************************
subroutine solve_core_val(nbasis,nocc,ncore,nval,zloc,amat,mo_loc,mo_can,aloc)
!********************************************************************************
! Determining the Lagrange multiplier of the core-val separation condition
! amat stores the occ-occ part of the A matrix stored in core-val order
!********************************************************************************
use common_mod, only: dcore,imem,scrfile1
implicit none
integer nbasis,ncore,nval,nocc
double precision zloc(ncore,nval),amat(nocc,nocc)
double precision mo_loc(nbasis,nocc),mo_can(nbasis,nocc),aloc(nbasis,nocc)
integer irhs,iu,is,ieig,ifock,itemp,imem_old,i,j
integer dblalloc
imem_old=imem
! Read eigenvalues of the Fock matrix
ieig=dblalloc(nocc)
ifock=dblalloc(nbasis*nbasis)
call read_from_disk('FOCK2_AB ',dcore(ifock),nbasis,.true.)
call dgemm('n','n',nbasis,nocc,nbasis,1.0d0,dcore(ifock),nbasis,mo_can,&
nbasis,0.0d0,dcore(imem+nocc**2),nbasis)
call dgemm('t','n',nocc,nocc,nbasis,1.0d0,mo_can,nbasis,dcore(imem+nocc**2),&
nbasis,0.0d0,dcore(imem),nocc)
do i=1,nocc
dcore(ieig+i-1)=dcore(imem+i-1+(i-1)*nocc)
enddo
call dbldealloc(ifock)
! Building r.h.s. of the equations
irhs=dblalloc(ncore*nval)
do j=1,nval
do i=1,ncore
dcore(irhs+i-1+(j-1)*ncore)=-(&
amat(i,j+ncore)-amat(j+ncore,i)+aloc(i,j+ncore)-aloc(j+ncore,i))
enddo
enddo
! calculate U matrix which transforms local MOs to canonical MOs
iu=dblalloc(nocc*nocc)
is=dblalloc(nbasis*nbasis)
call read_from_disk('OEINT ',dcore(is),nbasis,.true.)
call get_mo_overlap(nbasis,nocc,mo_loc,mo_can,dcore(is),dcore(iu))
! Transform r.h.s. to canocial basis
call dgemm('n','n',ncore,nval,nval,1.0d0,dcore(irhs),ncore,&
dcore(iu+ncore+ncore*nocc),nocc,0.0d0,dcore(is),ncore)
call dgemm('t','n',ncore,nval,ncore,1.0d0,dcore(iu),nocc, &
dcore(is),ncore,0.0d0,dcore(irhs),ncore)
! Solving the equations in canonical basis
do j=1,nval
do i=1,ncore
dcore(irhs+i-1+(j-1)*ncore)=dcore(irhs+i-1+(j-1)*ncore)/&
(dcore(ieig+i-1)-dcore(ieig+j+ncore-1))
enddo
enddo
! transforming solution to local basis
call dgemm('n','t',ncore,nval,nval,1.0d0,dcore(irhs),ncore,&
dcore(iu+ncore+ncore*nocc),nocc,0.0d0,dcore(is),ncore)
call dgemm('n','n',ncore,nval,ncore,1.0d0,dcore(iu),nocc, &
dcore(is),ncore,0.0d0,zloc,ncore)
call dbldealloc(imem_old)
end subroutine
!********************************************************************************
subroutine multiply_loc(ndim,zloc,res)
!********************************************************************************
! Calculates the a_ij(zloc)-a_ji(zloc) vector (i<j)
!********************************************************************************
use common_mod, only: dcore,icore,imem,nbasis
use embed_grad, only: nocc,ncore,imodip
implicit none
integer ndim
double precision zloc(ndim),res(ndim)
integer n,nval
integer ialoc,izfull,itemp,imem_old
integer dblalloc
imem_old=imem
nval=nocc-ncore
n=max(ncore,nval)
ialoc=dblalloc(n*n)
izfull=dblalloc(n*n)
itemp=dblalloc(nocc*nocc)
! Build the a(zloc) matrix
if(ncore*(ncore-1)/2 > 0) then
call build_aloc_occ(ncore,dcore(imodip),zloc,dcore(ialoc),&
dcore(izfull),dcore(itemp))
call antisymm(ncore,dcore(ialoc),res)
endif
if(nval*(nval-1)/2 > 0) then
call build_aloc_occ(nval,dcore(imodip+3*ncore**2), &
zloc(1+ncore*(ncore-1)/2),dcore(ialoc),dcore(izfull),&
dcore(itemp))
call antisymm(nval,dcore(ialoc),res(1+ncore*(ncore-1)/2))
endif
call dbldealloc(imem_old)
end subroutine
!********************************************************************************
subroutine sort_index(nbasis,nocc,aloc1,aloc2,perm)
!********************************************************************************
! Sorts the indices of the aloc1 matrix according to perm
! The ith row/column becomes the perm(i)th row/column
!********************************************************************************
implicit none
integer nbasis,nocc,perm(nocc)
double precision aloc1(nbasis,nocc),aloc2(nbasis,nocc)
integer i
do i=1,nocc
aloc2(:,perm(i))=aloc1(:,i)
enddo
do i=1,nocc
aloc1(perm(i),:)=aloc2(i,:)
aloc1(nocc+1:nbasis,i)=aloc2(nocc+1:nbasis,i)
enddo
end subroutine
!********************************************************************************
subroutine rhs_loc(nbasis,ncore,nval,nocc,ndim,amat,rhs)
!********************************************************************************
! Builds the r.h.s. of the localization part of the Z-vector eqautions
! amat stores the occ-occ part of the A matrix stored in core-val order
!********************************************************************************
use common_mod, only: dcore,icore,imem
implicit none
integer nbasis,ncore,nval,ndim
double precision amat(nocc,nocc)
double precision rhs(ndim)
integer itemp,nocc,i,j,n
integer dblalloc
itemp=dblalloc(nocc*nocc)
if(ncore > 0) then
call dlacpy('f',ncore,ncore,amat,nocc,dcore(itemp),ncore)
call antisymm(ncore,dcore(itemp),rhs)
endif
if(nval > 0) then
call dlacpy('f',nval,nval,amat(1+ncore,1+ncore),nocc,&
dcore(itemp),nval)
call antisymm(nval,dcore(itemp),rhs(1+ncore*(ncore-1)/2))
endif
call dscal(ndim,-1.0d0,rhs,1)
call dbldealloc(itemp)
end subroutine
!********************************************************************************
subroutine precond_with_diag(ndim,x,y)
!********************************************************************************
! Preconditioning with the diagonal
!********************************************************************************
use embed_grad, only: idiag
use common_mod, only: dcore
implicit none
integer ndim,i
double precision x(ndim),y(ndim)
do i=1,ndim
y(i)=x(i)/dcore(idiag+i-1)
enddo
end subroutine
!********************************************************************************
subroutine solve_bri(nbasis,zmat,amat,aloc,fock_file,mocoef_file,lloc)
!********************************************************************************
! Solves the Brillouin part of the Z-vector equations
!********************************************************************************
use embed_grad, only: imo,ifock2_ab,nvirt,nocc,ieigval,iumat,nocc_b
use common_mod, only: dcore,imem,scrfile1,iout
implicit none
integer nbasis
double precision amat(nocc+nvirt,nocc+nvirt),aloc(nocc+nvirt,nocc)
double precision zmat(nvirt*nocc)
character*16 fock_file,mocoef_file
logical lloc ! true => aloc is nonzero
integer nf,i,a,maxit,irhs,imem_old,ndim,imofock,lwork,iwork
double precision tol,mem
character*4 c4
character*8 c8
character*16 vecfile,errfile
integer dblalloc,idamax
interface
subroutine multiply_bri(ndim,zloc,res)
integer ndim
double precision zloc(ndim),res(ndim)
end subroutine
subroutine precond_bri(ndim,x,y)
integer ndim
double precision x(ndim),y(ndim)
end subroutine
double precision function norm(ndim,x)
integer ndim
double precision x(ndim)
end function
end interface
imem_old=imem
vecfile='fort.18 '
errfile='fort.19 '
open(scrfile1,file=trim(vecfile),form='unformatted')
close(scrfile1,status='delete')
open(scrfile1,file=trim(errfile),form='unformatted')
close(scrfile1,status='delete')
call getkey('scfmaxit',8,c4,4)
read(c4,*) maxit
call getkey('ldfgrad_tol',11,c8,8)
read(c8,*) i
tol=10.d0**dble(-i)
nf=nocc+nvirt
ndim=nocc*nvirt
! Reading MOs and Fock matrix from disk
imo=dblalloc(nf*nbasis)
ifock2_ab=dblalloc(nbasis*nbasis)
if(nf==nbasis) then
call read_from_disk(mocoef_file,dcore(imo),nbasis,.false.)
call read_from_disk(fock_file,dcore(ifock2_ab),nbasis,.true.)
else
call read_from_disk(mocoef_file,dcore(ifock2_ab),nbasis,.false.)
call dcopy(nbasis*nf,dcore(ifock2_ab+nocc_b*nbasis),1,dcore(imo),1)
call read_from_disk(fock_file,dcore(ifock2_ab),nbasis,.true.)
endif
! Building r.h.s of the equation
irhs=dblalloc(ndim)
if(lloc) then
do i=1,nocc
do a=1,nvirt
dcore(irhs+a-1+(i-1)*nvirt)=&
-(amat(nocc+a,i)-amat(i,nocc+a)+aloc(nocc+a,i))
enddo
enddo
else
do i=1,nocc
do a=1,nvirt
dcore(irhs+a-1+(i-1)*nvirt)=amat(i,nocc+a)-amat(nocc+a,i)
enddo
enddo
endif
! Diagonalizing MO Fock matrix for the preconditioner
iumat=dblalloc(nocc**2)
ieigval=dblalloc(nbasis)
imofock=dblalloc(nf**2)
call dgemm('n','n',nbasis,nf,nbasis,1.0d0,dcore(ifock2_ab),nbasis, &
dcore(imo),nbasis,0.0d0,dcore(imem),nbasis)
call dgemm('t','n',nf,nf,nbasis,1.0d0,dcore(imo),nbasis,dcore(imem),&
nbasis,0.0d0,dcore(imofock),nbasis)
call dlacpy('f',nocc,nocc,dcore(imofock),nf,dcore(iumat),nocc)
call dsyev('v','u',nocc,dcore(iumat),nocc,dcore(ieigval),mem,-1,i)
lwork=int(mem)
iwork=dblalloc(lwork)
call dsyev('v','u',nocc,dcore(iumat),nocc,dcore(ieigval),dcore(iwork),lwork,i)
if(i/=0) then
write(iout,'(A)') 'Error at diagonalization!'
write(iout,'(A,I8)') 'Error code: ',i
call mrccend(1)
endif
! Virtual orbital energies
do i=1,nvirt
dcore(ieigval+nocc+i-1)=dcore(imofock+nocc+i-1+(nocc+i-1)*nbasis)
enddo
call dbldealloc(imofock)
! call dscal(nbasis*nocc,dsqrt(2.0d0),dcore(imo),1)
call general_subspace_solver(ndim,zmat,dcore(irhs),maxit,tol,.false.,vecfile,&
errfile,norm,multiply_bri,precond_bri,.true.)
call dbldealloc(imem_old)
end subroutine
!********************************************************************************
subroutine precond_bri(ndim,zmat,res)
!********************************************************************************
! Preconditioner for the Brillouin part of the Z-vector equations
!********************************************************************************
use embed_grad, only: nocc,nvirt,iumat,ieigval
use common_mod, only: dcore,imem
implicit none
integer ndim
double precision zmat(ndim),res(ndim)
integer iz_can,i,a
integer dblalloc
iz_can=dblalloc(nvirt*nocc)
call dgemm('n','n',nvirt,nocc,nocc,1.0d0,zmat,nvirt,dcore(iumat),nocc,0.0d0,&
dcore(iz_can),nvirt)
do i=1,nocc
do a=1,nvirt
dcore(iz_can+a-1+(i-1)*nvirt)=dcore(iz_can+a-1+(i-1)*nvirt)/&
(dcore(ieigval+nocc+a-1)-dcore(ieigval+i-1))
enddo
enddo
call dgemm('n','t',nvirt,nocc,nocc,1.0d0,dcore(iz_can),nvirt,dcore(iumat),&
nocc,0.0d0,res,nvirt)
end subroutine
!********************************************************************************
subroutine multiply_bri(ndim,zmat,res)
!********************************************************************************
! Calculates the Atilde(z) matrix.
!********************************************************************************
use embed_grad, only: imo,ifock2_ab,ll_chfx,ll_csrhfx,ll_clrhfx,omega,pcm,&
ll_dft,nocc,nvirt,npos
use common_mod, only: dcore,nbasis,dfnbasis
implicit none
integer ndim
double precision zmat(ndim),res(ndim)
integer iatilde,itemp,a,i,nf
integer dblalloc
nf=nocc+nvirt
itemp=dblalloc(nbasis*nbasis)
call build_atilde_embed(nbasis,dfnbasis,nocc,nvirt,zmat,dcore(imo),&
dcore(ifock2_ab),npos,itemp,ll_chfx,ll_csrhfx,ll_clrhfx,omega,pcm,ll_dft)
iatilde=dblalloc(nbasis*nbasis)
call dgemm('t','n',nf,nf,nbasis,1.0d0,dcore(imo),nbasis,&
dcore(itemp),nbasis,0.0d0,dcore(iatilde),nbasis)
do i=1,nocc
do a=1,nvirt
res(a+(i-1)*nvirt)=dcore(iatilde+nocc+a-1+(i-1)*nf)-&
dcore(iatilde+i-1+(nocc+a-1)*nf)
enddo
enddo
call dbldealloc(itemp)
end subroutine
!********************************************************************************
subroutine build_atilde_embed(nbasis,dfnbasis,nocc,nvirt,zmat,mo,fock2_ab,npos, &
iatilde,chfx,clrhfx,csrhfx,omega,pcm,dft)
!********************************************************************************
! Builds Atilde(z) matrix for embedding gradient
!********************************************************************************
use mcscf, only: build_mmat_
use common_mod, only: dcore,imem,nirmax,scfalg,orbperir
implicit none
integer nbasis,dfnbasis,nocc,nvirt,npos
integer iatilde
double precision zmat(nvirt,nocc),mo(nbasis,nbasis),fock2_ab(nbasis,nbasis)
double precision chfx,csrhfx,clrhfx,omega
character*32 pcm,dft
integer nbasis2,nc,ncorenew,lfin,nocca,noccb,nf,nn
integer offset(nirmax),sqroffset(nirmax)
integer idens(2),im(2),ic,ic2(2),iw,imem_old
integer dblalloc
double precision exc,devparr(2)
character*8 embed
character*16 scfdamp,calctype
character*32 dvxc_dft
data offset /1,1,1,1,1,1,1,1/
data sqroffset /1,1,1,1,1,1,1,1/
equivalence(ic2(1),ic)
interface
subroutine setup_mmat_vars(nbasis,n1,n2,dens,c,z,mo1,mo2,scfalg)
implicit none
integer n1,n2,nbasis
double precision dens(nbasis,nbasis),c(2*nbasis*min(n1,n2))
double precision, target :: mo1(nbasis,n1),mo2(nbasis,n2)
double precision z(n1,n2)
character*8 scfalg
end subroutine
end interface
imem_old=imem
call getkey('calc',4,calctype,16)
call getkey('embed',5,embed,8)
nf=nocc+nvirt
nbasis2=nbasis**2
scfdamp='off '
ncorenew=0
devparr=1.0d0
lfin=0
nocca=min(nocc,nvirt)
noccb=nocca
idens(1)=dblalloc(nbasis2)
idens(2)=idens(1)
ic=dblalloc(2*nocc*nbasis)
ic2(2)=ic
im=iatilde
call dfillzero(dcore(iatilde),nbasis**2)
! Contribution of the 2-electron integrals
! It is the same as the QSCF RHF Hessian M matrix (see Molecular electronic structure theory 10.8.64)
call setup_mmat_vars(nbasis,nvirt,nocc,dcore(idens(1)),dcore(ic),zmat,&
mo(1,nocc+1),mo,scfalg)
nn=orbperir(1)
orbperir(1)=nocca
call build_mmat_(offset,sqroffset,idens,im,npos,zmat,mo,chfx,100,dfnbasis,&
nocc,nocc,lfin,exc,'off ',.false.,scfdamp,' ',0,nc,ncorenew, &
devparr,clrhfx,csrhfx,omega,.false.,pcm,nocca,noccb,ic2)
orbperir(1)=nn
! call build_mmat(offset,sqroffset,idens,im,npos,zmat,mo,chfx,2,&
! dfnbasis,nocc,nocc,lfin,exc,'off ',.false.,scfdamp,' ',0, &
! nc,ncorenew,devparr,clrhfx,csrhfx,omega,.false.,pcm)
! -1 factor: build_mmat builds the -1*M matrix
call dsymm('l','u',nbasis,nocc,-1.0d0,dcore(iatilde),nbasis,mo,nbasis,0.0d0, &
dcore(idens(1)),nbasis)
call dfillzero(dcore(iatilde),nbasis2)
call dcopy(nbasis*nocc,dcore(idens(1)),1,dcore(iatilde),1)
! contribution from the Fock matrix
iw=dblalloc(nbasis**2)
call build_atilde_fock(nbasis,nf,nocc,nvirt,zmat,fock2_ab,mo,dcore(iatilde), &
dcore(iw))
! Contribution from the second derivative of the XC functional
if(trim(dft)/='off') then
dvxc_dft=dft
if(trim(dft)=='user'.and.embed.ne.'off ') &
dvxc_dft='ll_user '
! call dgemm('n','n',nbasis,nocc,nvirt,1.0d0,zmat,nvirt,mo(1,1+nocc),&
! nbasis,0.0d0,dcore(iw),nbasis)
! call dgemm('n','t',nbasis,nbasis,nocc,1.0d0,dcore(iw),nbasis,mo, &
! nbasis,0.0d0,dcore(idens(1)),nbasis)
! call symmat(dcore(idens(1)),nbasis)
! call amat_embed_dft2(nbasis,nocc,dcore(idens(1)),mo,dcore(iw),dvxc_dft)
! call dsymm('l','u',nbasis,nocc,1.0d0,dcore(iw),nbasis,mo,nbasis,&
! 1.0d0,dcore(iatilde),nbasis)
call dfillzero(dcore(iw),nbasis**2)
call build_atilde_dft(nbasis,nocc,nocc,zmat,zmat,mo,mo,'rhf ',&
dvxc_dft,dcore(iatilde))
endif
call dbldealloc(imem_old)
end subroutine
!********************************************************************************
subroutine build_x_embed(nbasis,nocc,amat,atilde,aloc,x)
!********************************************************************************
! Builds the x matrix (the Lagrange multipliers of the othogonality condition of
! the low level molecular orbitals)
!********************************************************************************
implicit none
integer nbasis,nocc
double precision amat(nbasis,nbasis),atilde(nbasis,nbasis),aloc(nbasis,nocc)
double precision x(nbasis,nbasis)
call dcopy(nbasis**2,atilde,1,x,1)
call daxpy(nbasis*nbasis,1.0d0,amat,1,x,1)
call daxpy(nbasis*nocc,1.0d0,aloc,1,x,1)
call dscal(nbasis**2,-0.5d0,x,1)
call symmat(x,nbasis)
end subroutine
!********************************************************************************
subroutine loc_dens_embed(nbasis,nocc,ncore,nval,zloc,dens)
!********************************************************************************
! Builds the density that is contracted with the localization constrain
!********************************************************************************
use common_mod, only: dcore
implicit none
integer nbasis,nocc,ncore,nval
double precision dens(nbasis,nbasis)
double precision zloc(ncore*(ncore-1)/2+nval*(nval-1)/2)
integer izfull,imo,i,itemp
integer dblalloc
izfull=dblalloc(nocc*nocc)
itemp=dblalloc(nbasis*nocc)
imo=dblalloc(nbasis**2)
call read_from_disk('MOCOEF.LOC ',dcore(imo),nbasis,.false.)
dens=0.0d0
if(ncore>0) then
call create_full(ncore,zloc,dcore(izfull))
call dgemm('n','n',nbasis,ncore,ncore,1.0d0,dcore(imo),nbasis,dcore(izfull),&
ncore,0.0d0,dcore(itemp),nbasis)
call dgemm('n','t',nbasis,nbasis,ncore,1.0d0,dcore(itemp),nbasis,dcore(imo),&
nbasis,1.0d0,dens,nbasis)
endif
if(nval>0) then
i=ncore*(ncore-1)/2
call create_full(ncore,zloc(i),dcore(izfull))
call dgemm('n','n',nbasis,nval,nval,1.0d0,dcore(imo+nbasis*ncore),nbasis, &
dcore(izfull),ncore,0.0d0,dcore(itemp),nbasis)
call dgemm('n','t',nbasis,nbasis,nval,1.0d0,dcore(itemp),nbasis,&
dcore(imo+nbasis*ncore),nbasis,1.0d0,dens,nbasis)
endif
call dbldealloc(izfull)
end subroutine
!********************************************************************************
subroutine read_subsys_mo(nbasis,nocc_a,nocc_b,nvirt,mo)
!********************************************************************************
! Reads subsystem MOs (occ. and virt.) calculated with the high level method
! The MOs of the environment is not included
!********************************************************************************
use common_mod, only: dcore,imem
implicit none
integer nbasis,nocc_a,nocc_b,nvirt
double precision mo(nbasis,nocc_a+nvirt)
integer imo_full
integer dblalloc
imo_full=dblalloc(nbasis*nbasis)
call read_from_disk('MOCOEF ',dcore(imo_full),nbasis,.false.)
call dcopy(nbasis*(nocc_a+nvirt),dcore(imo_full+nbasis*nocc_b),1,mo,1)
call dbldealloc(imo_full)
end subroutine
!********************************************************************************
subroutine setup_mmat_vars(nbasis,n1,n2,dens,c,z,mo1,mo2,scfalg)
!********************************************************************************
! Calculates some variables for the M matrix construction
! (n1,n2): the dimansions of the z vector
! dens,c : transformed density and mo coefficients
! mo1,mo2: the mo coefficients for the first and second index of the z vector
!********************************************************************************
use common_mod, only: imem,dcore
use mcscf, only: ldf,symm_dens
implicit none
integer n1,n2,nbasis
double precision dens(nbasis,nbasis),c(2*nbasis*min(n1,n2))
double precision, target :: mo1(nbasis,n1),mo2(nbasis,n2)
double precision z(n1,n2)
character*8 scfalg
integer imem_old,itmp1,itmp2,itmp3
integer ns,nl
integer dblalloc
double precision, pointer, CONTIGUOUS :: mos(:,:),mol(:,:)
double precision alpha,beta,delta
character*1 transp
imem_old=imem
! Selecting the smaller and the larger dimensions of z
if(n1<=n2) then
ns=n1
nl=n2
mos=>mo1
mol=>mo2
transp='t'
else
ns=n2
nl=n1
mos=>mo2
mol=>mo1
transp='n'
endif
! Build modified MOs for build_mmat
itmp1=dblalloc(nbasis*ns)
itmp2=dblalloc(nbasis*ns)
call dgemm('n',transp,nbasis,ns,nl,1.0d0,mol,nbasis,z,n1,&
0.0d0,dcore(itmp1),nbasis)
call dfillzero(dcore(itmp2),nbasis*ns)
call daxpy(nbasis*ns,-1.0d0,dcore(itmp1),1,dcore(itmp2),1)
call daxpy(nbasis*ns,1.0d0,mos,1,dcore(itmp2),1)
call daxpy(nbasis*ns,1.0d0,mos,1,dcore(itmp1),1)
if(trim(scfalg).ne.'disk') then
call dcopy(nbasis*ns,dcore(itmp1),1,c,1)
call dcopy(nbasis*ns,dcore(itmp2),1,c(1+nbasis*ns),1)
else
itmp3=dblalloc(2*nbasis*ns)
call dcopy(nbasis*ns,dcore(itmp1),1,dcore(itmp3),1)
call dcopy(nbasis*ns,dcore(itmp2),1,dcore(itmp3+nbasis*ns),1)
call motransp(2*ns,nbasis,c,dcore(itmp3),.false.)
call dbldealloc(itmp3)
endif
call dscal(2*ns*nbasis,1.0d0/dsqrt(2.0d0),c,1)
! Building modified density
if(ldf) then
alpha = 1.0d0
beta = -1.0d0
delta = -0.5d0
else
alpha = -4.0d0 ! multiplies the whole modified density
beta = 1.0d0 ! multiplies only the diagonal part
delta = 1.0d0 ! multiplies only the offdiagonal part
endif
call dgemm('n','t',n1,nbasis,n2,alpha,z,n1,mo2,nbasis,0.0d0,&
dcore(itmp1),n1)
call dgemm('n','n',nbasis,nbasis,n1,1.0d0,mo1,nbasis,dcore(itmp1),&
n1,0.0d0,dens,nbasis)
call symm_dens(dens,nbasis,beta,delta)
call dbldealloc(imem_old)
end subroutine
double precision function norm(ndim,x)
implicit none
integer ndim
double precision x(ndim),dnrm2
norm=dnrm2(ndim,x,1)/dsqrt(dble(ndim))
return
end function
!!********************************************************************************
!subroutine kaczmarz(nbasis,nocc,nvirt,mofock,x,rhs,row,row_norms,tol,maxit)
!!********************************************************************************
!!********************************************************************************
! implicit none
! integer nbasis,nocc,nvirt,maxit
! double precision mofock(nbasis,nbasis),x(nvirt,nocc),rhs(nvirt,nocc)
! double precision row_norms(nvirt,nocc),row(nvirt,nocc)
! double precision tol
!
! integer a,b,i,j,ndim,step
! double precision temp,max_coef
! logical conv
!
! double precision ddot
!
! ndim=nocc*nvirt
! row_norms=0.0d0
! x=0.0d0
!
! do i=1,nocc
! do a=1,nvirt
! row=0.0d0
! do j=1,nocc
! row(a,j)=row(a,j)+mofock(j,i)
! enddo
! do b=1,nvirt
! row(b,i)=row(b,i)-mofock(nocc+b,nocc+a)
! enddo
! row_norms(a,i)=ddot(ndim,row,1,row,1)
! enddo
! enddo
!
! step=0
! conv=.false.
! write(*,'(A)') 'Step Convergence'
! do while(.not. conv .and. step<maxit)
! conv=.true.
! step=step+1
! max_coef=0.0d0
! do i=1,nocc
! do a=1,nvirt
! row=0.0d0
! do j=1,nocc
! row(a,j)=row(a,j)+mofock(j,i)
! enddo
! do b=1,nvirt
! row(b,i)=row(b,i)-mofock(nocc+b,nocc+a)
! enddo
! temp=(rhs(a,i)-ddot(ndim,row,1,x,1))
! max_coef=max(temp,max_coef)
! x=x+temp/row_norms(a,i)*row
! enddo
! enddo
! conv=max_coef<tol
! write(*,'(I4,5X,ES10.3)') step,max_coef
! enddo
!
!end subroutine
!subroutine test_loc(nocc,modip,mat)
! implicit none
! integer nocc
! double precision mat(nocc,nocc,nocc,nocc),modip(nocc,nocc,3)
! integer i,j,k,l,xyz
!
! mat=0.0d0
! do j=1,nocc
! do i=1,nocc
! do l=1,nocc
! do k=1,nocc
! do xyz=1,3
! if(l==i) then
! mat(k,l,i,j)=mat(k,l,i,j)+2.0d0*modip(k,l,xyz)*modip(i,j,xyz)+&
! (modip(i,i,xyz)-modip(j,j,xyz))*modip(k,j,xyz)
! endif
! if(l==j) then
! mat(k,l,i,j)=mat(k,l,i,j)-2.0d0*modip(k,l,xyz)*modip(i,j,xyz)+&
! (modip(i,i,xyz)-modip(j,j,xyz))*modip(k,i,xyz)
! endif
! enddo
! enddo
! enddo
! enddo
! enddo
! write(*,*) 'mat'
! call prmx_ldf(mat,nocc**2,nocc**2,nocc**2)
!end subroutine