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

1368 lines
51 KiB
Fortran
Executable File

************************************************************************
subroutine basoptdrv(gbasfile)
************************************************************************
* Driver for basis set optimization
************************************************************************
use optim
implicit none
integer npara,gbasfile,natoms,iatoms,natomsmax,nangmax,i,j,ii
integer nprimmax,ncontrmax,iang,icontr,iprim,istat
parameter(natomsmax=50,nangmax=20,nprimmax=50,ncontrmax=50)
integer nang(natomsmax),nprim(0:nangmax,natomsmax)
integer angt(0:nangmax,natomsmax)
integer ncontr(0:nangmax,natomsmax),ebf(0:nangmax,natomsmax)
integer ibopt(natomsmax) ! whether or not the given bond center should be optimized
integer iaopt(0:nangmax,natomsmax) ! whether or not the given angular momentum should be optimized
integer ieopt(nprimmax,0:nangmax,natomsmax) ! whether or not the given exponent should be optimized
integer igopt(nprimmax,ncontrmax,0:nangmax,natomsmax) ! whether or not the given contr. coeff. should be optimized
character(len=2) :: cbsm(natomsmax) ! whether or not the given bond center contains special marks(++,--)
character(len=2) :: casm(0:nangmax,natomsmax) ! whether or not the given angular momentum contains special marks(++,--)
character(len=2) :: cesm(nprimmax,0:nangmax,natomsmax) ! whether or not the given exponent contains special marks(++,--)
character(len=2) :: cgsm(nprimmax,ncontrmax,0:nangmax,natomsmax) ! whether or not the given contr. coeff. contains special marks(++,--)
real*8 gexp(nprimmax,0:nangmax,natomsmax),gbondc(natomsmax)
real*8 gcoef(nprimmax,ncontrmax,0:nangmax,natomsmax)
real*8 ell(0:nangmax,natomsmax),fgoal
logical lbond,lopt,lcv ! lcv is a logical variable for CV basis set optimizations
character*1 c1(8)
character*4 c4
character*8 c8,optalg
character*16 c16
character*64 line
equivalence(c1,c8)
character(len=20),dimension(natomsmax) :: strb ! temporary string for bond center position
character(len=20),dimension(nangmax) :: stra ! temporary string for angular momentum
character(len=20),dimension(nprimmax) :: strp ! temporary string for the primitive bf.'s exponents
character(len=20),dimension(ncontrmax) :: strc ! temporary string for the contraction coeff.
character(len=64),dimension(natomsmax) :: title ! to store the basis title
character(len=64),dimension(natomsmax) :: comment ! to store the comment line
character(len=2) :: mark = ' ' ! whether or not the user marked sg in the GENBAS file
character*24 basopt
integer :: n_ae,n_fc ! frozen all electron and frozen core orbitals
integer :: optmaxit,maxfun,plevel,lnum
real*8 :: opttol,steptol
logical :: bool ! whether or not the optimization succeeded
#if !defined (gfortran) && !defined (G95)
integer*4 system
external system
#endif
interface
subroutine read_GENBAS(gbasfile,lbond,title,comment,
$natomsmax,nangmax,nprimmax,ncontrmax,
$iatoms,npara,nang,nprim,angt,ncontr,
$ibopt,iaopt,ieopt,igopt,cbsm,casm,cesm,cgsm,gexp,gcoef,
$gbondc,ebf,ell)
implicit none
integer,intent(in) :: gbasfile
integer,intent(in) :: natomsmax,nangmax,nprimmax,ncontrmax!,nangmax
logical,intent(out) :: lbond
integer,intent(out) :: iatoms,npara
integer,intent(out) :: nang(natomsmax),nprim(0:nangmax,natomsmax)
integer,intent(out) :: angt(0:nangmax,natomsmax)
integer,intent(out) :: ebf(0:nangmax,natomsmax)
integer,intent(out) :: ncontr(0:nangmax,natomsmax)
integer,intent(out) :: ibopt(natomsmax) ! whether or not the given bond center should be optimized
integer,intent(out) :: iaopt(0:nangmax,natomsmax) ! whether or not the given angular momentum should be optimized
integer,intent(out) :: ieopt(nprimmax,0:nangmax,natomsmax) ! whether or not the given exponent should be optimized
integer,intent(out) :: igopt(nprimmax,ncontrmax,0:nangmax,
$ natomsmax) ! whether or not the given contr. coeff. should be optimized
character(len=2), intent(out) :: cbsm(natomsmax) ! whether or not the given bond center contains special marks(++,--)
character(len=2), intent(out) :: casm(0:nangmax,natomsmax) ! whether or not the given angular momentum contains special marks(++,--)
character(len=2), intent(out) :: cesm(nprimmax,0:nangmax,
$ natomsmax) ! whether or not the given exponent contains special marks(++,--)
character(len=2), intent(out) :: cgsm(nprimmax,ncontrmax,
$ 0:nangmax,natomsmax) ! whether or not the given contr. coeff. contains special marks(++,--)
double precision, intent(out) :: gexp(nprimmax,0:nangmax,
$ natomsmax)
double precision, intent(out) :: gcoef(nprimmax,ncontrmax,
$ 0:nangmax,natomsmax)
double precision, intent(out) :: gbondc(natomsmax)
double precision, intent(out) :: ell(0:nangmax,natomsmax)
character(len=64),intent(out) :: title(natomsmax) ! to store the basis title
character(len=64),intent(out) :: comment(natomsmax) ! to store the comment line
end subroutine read_GENBAS
end interface
! take the parameters from the MINP file
call getkey('basopt',6,basopt,24)
if(basopt.eq.'on ') then
fgoal = 0.d0
elseif(basopt.eq.'cv ') then
lcv = .True.
fgoal = 0.d0
call getcvkeys(n_ae,n_fc)
else
read(basopt,*) fgoal
endif
call getkey('verbosity',9,c4,4)
read(c4,'(i4)') plevel
call getkey('optalg',6,optalg,8)
call getkey('optmaxit',8,c8,8)
read(c8, '(i8)') optmaxit
maxfun = optmaxit*15 ! this is quite arbitrary but should be fine in most situations
call getkey('optetol',7,c16,16)
read(c16, *) opttol
call getkey('optstol',7,c16,16)
read(c16, *) steptol
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
istat=system('cp GENBAS GENBAS.init')
if(istat.ne.0) then
write(6,*) 'Error while copying GENBAS to GENBAS.init'
stop
endif
call read_GENBAS(gbasfile,lbond,title,comment,
$natomsmax,nangmax,nprimmax,ncontrmax,
$natoms,npara,nang,nprim,angt,ncontr,
$ibopt,iaopt,ieopt,igopt,cbsm,casm,cesm,cgsm,gexp,gcoef,
$gbondc,ebf,ell)
bool=.False.
if(optalg.eq.'simplex ') then
call basopt_smplx(npara,gbasfile,title,comment,
$ natoms,lcv,n_ae,n_fc,lbond,gbondc,nang,
$ natomsmax,nangmax,nprimmax,ncontrmax,angt,nprim,ncontr,
$ ibopt,ieopt,igopt,cbsm,casm,cesm,cgsm,gexp,gcoef,ebf,ell,
$ fgoal,opttol,steptol,optmaxit,maxfun,plevel,bool)
else
write(*,*) 'At the moment this should not happen
$- csonti from basopt.f'
!
! call basoptcore(npara,gbasfile,natoms,lbond,gbondc,nang,nangmax,
! $nprimmax,ncontrmax,angt,nprim,ncontr,igopt,gexp,gcoef,ebf,ell)
endif
if(bool) then
istat=system('mv GENBAS.tmp GENBAS.opt')
if(istat.ne.0) then
write(6,*) 'Error while renaming GENBAS.tmp to GENBAS.opt'
stop
endif
endif
return
end subroutine basoptdrv
************************************************************************
************************************************************************
subroutine getcvkeys(n_all,n_frozen)
************************************************************************
*
* n_frozen - number of frozen core orbitals (integer)
* n_all - number of frozen orbitals in "all electron" calculations
* for Cl n_all is equal to 1 (1s^2) and n_frozen is equal to 5 (1s^22s^2p^6)
* n_all < n_frozen
*
************************************************************************
use optim, only: str_cut,str_itest
implicit none
integer, intent(out):: n_all,n_frozen
character(len=24) :: basopt
integer :: check
integer, parameter :: minpfile = 68
character(len=512) :: line,stmp,s_ae,s_fc
character(len=2) :: sm
logical bool
open(minpfile, file='MINP',status='old')
call getkeym('basopt',6,basopt,24)
read(minpfile,200,iostat=check) line
if(check .eq. 0) then
stmp=adjustl(trim(line))
call str_cut(stmp,',',s_ae,s_fc)
call str_itest(s_ae,n_all,bool,sm)
call str_itest(s_fc,n_frozen,bool,sm)
else
write(*,*) 'Error: a line is expected right after basopt=cv!'
write(*,*) 'Example: 1,5'
write(*,*) 'Explanation:'
write(*,*) 'the 1s^2 electrons are frozen in the "all electron"'
write(*,*) 'calculation, and the 1s^2 2s^2 and 3p^6 electrons'
write(*,*) 'are frozen in the frozen core calculation.'
endif
close(minpfile)
200 format(a512)
end subroutine getcvkeys
************************************************************************
************************************************************************
subroutine read_GENBAS(gbasfile,lbond,title,comment,
$natomsmax,nangmax,nprimmax,ncontrmax,
$iatoms,npara,nang,nprim,angt,ncontr,
$ibopt,iaopt,ieopt,igopt,cbsm,casm,cesm,cgsm,gexp,gcoef,
$gbondc,ebf,ell)
use optim
implicit none
integer,intent(in) :: gbasfile
integer,intent(in) :: natomsmax,nangmax,nprimmax,ncontrmax!,nangmax
logical,intent(out) :: lbond
integer,intent(out) :: iatoms,npara
integer,intent(out) :: nang(natomsmax),nprim(0:nangmax,natomsmax)
integer,intent(out) :: angt(0:nangmax,natomsmax)
integer,intent(out) :: ebf(0:nangmax,natomsmax)
integer,intent(out) :: ncontr(0:nangmax,natomsmax)
integer,intent(out) :: ibopt(natomsmax) ! whether or not the given bond center should be optimized
integer,intent(out) :: iaopt(0:nangmax,natomsmax) ! whether or not the given angular momentum should be optimized
integer,intent(out) :: ieopt(nprimmax,0:nangmax,natomsmax) ! whether or not the given exponent should be optimized
integer,intent(out) :: igopt(nprimmax,ncontrmax,0:nangmax,
$ natomsmax) ! whether or not the given contr. coeff. should be optimized
character(len=2), intent(out) :: cbsm(natomsmax) ! whether or not the given bond center contains special marks(++,--)
character(len=2), intent(out) :: casm(0:nangmax,natomsmax) ! whether or not the given angular momentum contains special marks(++,--)
character(len=2), intent(out) :: cesm(nprimmax,0:nangmax,
$ natomsmax) ! whether or not the given exponent contains special marks(++,--)
character(len=2), intent(out) :: cgsm(nprimmax,ncontrmax,
$ 0:nangmax,natomsmax) ! whether or not the given contr. coeff. contains special marks(++,--)
double precision, intent(out) :: gexp(nprimmax,0:nangmax,
$ natomsmax)
double precision, intent(out) :: gcoef(nprimmax,ncontrmax,
$ 0:nangmax,natomsmax)
double precision, intent(out) :: gbondc(natomsmax)
double precision, intent(out) :: ell(0:nangmax,natomsmax)
character(len=64),intent(out) :: title(natomsmax) ! to store the basis title
character(len=64),intent(out) :: comment(natomsmax) ! to store the comment line
!locals
integer iang,icontr,iprim,check,natoms,i,j
character(len=20),dimension(natomsmax) :: strb ! temporary string for bond center position
character(len=20),dimension(0:nangmax) :: stra ! temporary string for angular momentum
character(len=20),dimension(nprimmax) :: strp ! temporary string for the primitive bf.'s exponents
character(len=20),dimension(ncontrmax) :: strc ! temporary string for the contraction coeff.
character(len=1),dimension(8) :: c1
character(len=2) :: mark = ' ' ! whether or not the user marked sg in the GENBAS file
character(len=8) :: c8
character(len=64) :: line
equivalence(c1,c8)
logical lopt
lbond = .False.
iatoms= 1
npara = 0
check = 0
do j=1,natomsmax
do i=0,nangmax
ebf(i,j)=0
enddo
enddo
open(gbasfile,file='GENBAS',status='old')
natoms= 1
read(gbasfile,"(a8)") c8
close(gbasfile)
call lowercase(c1,c1,8)
if(c1(1).eq.'b'.and.c1(3).ne.':') lbond = .True. ! check whether it is a standard GENBAS or bond function form
if(lbond) then
call bondfuncs()
else
call standard()
endif
************************************************************************
contains ! read_GENBAS
************************************************************************
************************************************************************
subroutine standard()
open(gbasfile,file='GENBAS',status='old')
do while(check.eq.0)
read(gbasfile,"(a64)",iostat=check) title(iatoms) ! title line
read(gbasfile,"(a64)") comment(iatoms) ! comment line
read(gbasfile,*) ! blank line
read(gbasfile,*) nang(iatoms) ! number of angular momentums
nang(iatoms)=nang(iatoms)-1 ! to coincide with the usual convention s->0, p->1, etc.
read(gbasfile,*) (stra(iang),iang=0,nang(iatoms)) ! angular momentum's type (s,p,d,etc)
do iang=0, nang(iatoms)
iaopt(iang,iatoms)=1
lopt = .True.
mark = ' '
call str_itest(stra(iang),angt(iang,iatoms),lopt,mark)
casm(iang,iatoms)=mark
if(.not.lopt) then
iaopt(iang,iatoms)=0
endif
enddo
read(gbasfile,*) (ncontr(iang,iatoms),iang=0,nang(iatoms)) ! num. of contracted basis funcs.
read(gbasfile,*) (nprim (iang,iatoms),iang=0,nang(iatoms)) ! num. of primitive basis funcs.
do iang=0,nang(iatoms)
read(gbasfile,*)
read(gbasfile,*)
$(strp(iprim),iprim=1,nprim(iang,iatoms))
do iprim=1,nprim(iang,iatoms)
if(iaopt(iang,iatoms).gt.0) then
ieopt(iprim,iang,iatoms)=1
lopt = .True.
else
ieopt(iprim,iang,iatoms)=0
lopt = .False.
endif
mark = ' '
call str_dtest(strp(iprim),gexp(iprim,iang,iatoms),lopt,
$ mark)
cesm(iprim,iang,iatoms) = mark
if(.not.lopt) then
ieopt(iprim,iang,iatoms)=0
else
ieopt(iprim,iang,iatoms)=1
npara=npara+1
endif
enddo
read(gbasfile,"(a64)") line ! blank line between exponents and contraction coeffs.
do iprim=1,nprim(iang,iatoms)
read(gbasfile,*)
$(strc(icontr),icontr=1,ncontr(iang,iatoms))
do icontr=1,ncontr(iang,iatoms)
if(iaopt(iang,iatoms).gt.0) then
lopt = .True.
else
lopt = .False.
endif
mark = ' '
call str_dtest(strc(icontr),gcoef(iprim,icontr,iang,iatoms),lopt,
$ mark)
cgsm(iprim,icontr,iang,iatoms) = mark
if(dabs(gcoef(iprim,icontr,iang,iatoms)).gt.1d-10.and.
$ dabs(1.d0-gcoef(iprim,icontr,iang,iatoms)).gt.1d-10.and.
$ lopt) then
npara=npara+1
igopt(iprim,icontr,iang,iatoms)=1
else
igopt(iprim,icontr,iang,iatoms)=0
endif
enddo
enddo
enddo ! iang
read(gbasfile,"(a64)",iostat=check) line ! blank line
do while (check.eq.0 .and. trim(line).eq.'')
read(gbasfile,"(a64)",iostat=check) line ! blank lines between basis sets
enddo
if(check.eq.0) then
iatoms = iatoms + 1
backspace(gbasfile)
endif
enddo ! iatoms
close(gbasfile)
end subroutine standard
************************************************************************
************************************************************************
subroutine bondfuncs()
open(gbasfile,file='GENBAS',status='old')
iatoms= 1
natoms= 1
read(gbasfile,"(a64)") title(iatoms) ! title line
read(gbasfile,"(a64)") comment(iatoms) ! comment line
read(gbasfile,*)
read(gbasfile,*) natoms
do while(iatoms .le. natoms)
read(gbasfile,*)
read(gbasfile,*) strb(iatoms)
npara=npara+1
ibopt(iatoms)=1
lopt = .True.
mark = ' '
call str_dtest(strb(iatoms),gbondc(iatoms),lopt,mark)
cbsm(iatoms) = mark
if(.not.lopt) then
npara=npara-1
ibopt(iatoms)=0
endif
read(gbasfile,*)
read(gbasfile,*) nang(iatoms)
nang(iatoms)=nang(iatoms)-1
read(gbasfile,*) (stra(iang),iang=0,nang(iatoms)) ! angular momentum's type (s,p,d,etc)
do iang=0, nang(iatoms)
iaopt(iang,iatoms)=1
lopt = .True.
mark = ' '
call str_itest(stra(iang),angt(iang,iatoms),lopt,mark)
casm(iang,iatoms)=mark
if(.not.lopt) then
iaopt(iang,iatoms)=0
endif
enddo
read(gbasfile,*) (ncontr(iang,iatoms),iang=0,nang(iatoms)) ! num. of contracted basis funcs.
read(gbasfile,*) (nprim (iang,iatoms),iang=0,nang(iatoms)) ! num. of primitive basis funcs.
do iang=0,nang(iatoms)
read(gbasfile,*)
read(gbasfile,*)
$(strp(iprim),iprim=1,nprim(iang,iatoms))
do iprim=1,nprim(iang,iatoms)
if(iaopt(iang,iatoms).gt.0) then
ieopt(iprim,iang,iatoms)=1
lopt = .True.
else
ieopt(iprim,iang,iatoms)=0
lopt = .False.
endif
mark = ' '
call str_dtest(strp(iprim),gexp(iprim,iang,iatoms),lopt,
$ mark)
cesm(iprim,iang,iatoms) = mark
if(.not.lopt) then
ieopt(iprim,iang,iatoms)=0
else
ieopt(iprim,iang,iatoms)=1
npara=npara+1
endif
enddo
read(gbasfile,"(a64)") line
ebf(iang,iatoms)=0
if(trim(line).ne.'') then
backspace(gbasfile)
read(gbasfile,*) ell(iang,iatoms)
ebf(iang,iatoms)=1
npara=npara+1
endif
do iprim=1,nprim(iang,iatoms)
read(gbasfile,*)
$(strc(icontr),icontr=1,ncontr(iang,iatoms))
do icontr=1,ncontr(iang,iatoms)
if(iaopt(iang,iatoms).gt.0) then
lopt = .True.
else
lopt = .False.
endif
mark = ' '
call str_dtest(strc(icontr),gcoef(iprim,icontr,iang,iatoms),lopt,
$ mark)
cgsm(iprim,icontr,iang,iatoms) = mark
if(dabs(gcoef(iprim,icontr,iang,iatoms)).gt.1d-10.and.
$ dabs(1.d0-gcoef(iprim,icontr,iang,iatoms)).gt.1d-10.and.
$ lopt) then
npara=npara+1
igopt(iprim,icontr,iang,iatoms)=1
else
igopt(iprim,icontr,iang,iatoms)=0
endif
enddo
enddo
enddo
iatoms = iatoms+1
enddo
iatoms = iatoms-1
close(gbasfile)
end subroutine bondfuncs
************************************************************************
end subroutine read_GENBAS
************************************************************************
************************************************************************
subroutine basopt_smplx(npara,gbasfile,title,comment,
$natoms,lcv,n_ae,n_fc,lbond,gbondc,nang,
$natomsmax,nangmax,nprimmax,ncontrmax,angt,nprim,ncontr,ibopt,
$ieopt,igopt,cbsm,casm,cesm,cgsm,gexp,gcoef,ebf,ell,
$fgoal,ftol,xtol,maxiter,maxfun,plevel,bool)
************************************************************************
use optim
implicit none
integer npara,ispec,ireac,gbasfile,i
integer ipara,jpara,iparat,nit,ipiv(npara),mnit,maxit,mmaxit
real*8 para(npara),fpm,grd,hss,moldener,oldener
real*8 paraold(npara),eps,grad(npara),htokj,conv,f00,fpp,fmm,fmp
real*8 fm(npara),fp(npara),hess(npara,npara),param(npara),maxstep
real*8 scr(3*npara**2),hess2(npara,npara),eig(npara)
real*8 paraold2(npara),mconv,gtol,mmaxstep,maxstepin,fgoal
integer*4 istat,isyev
!
integer natoms,iatoms,iang,iprim,icontr,nangmax,nprimmax,ncontrmax
integer natomsmax
integer nang(*),nprim(0:nangmax,*),ebf(0:nangmax,*)
integer ncontr(0:nangmax,*),angt(0:nangmax,*)
integer ibopt(natomsmax)
integer ieopt(nprimmax,0:nangmax,*)
integer igopt(nprimmax,ncontrmax,0:nangmax,*)
real*8 gexp(nprimmax,0:nangmax,*),gbondc(*),ell(0:nangmax,*)
real*8 gcoef(nprimmax,ncontrmax,0:nangmax,*)
logical lbond,lcv
character(len=2) :: cbsm(natomsmax) ! whether or not the given bond center contains special marks(++,--)
character(len=2) :: casm(0:nangmax,natomsmax) ! whether or not the given angular momentum contains special marks(++,--)
character(len=2) :: cesm(nprimmax,0:nangmax,natomsmax) ! whether or not the given exponent contains special marks(++,--)
character(len=2) :: cgsm(nprimmax,ncontrmax,0:nangmax,natomsmax) ! whether or not the given contr. coeff. contains special marks(++,--)
character(len=64) :: title(natomsmax) ! to store the basis title
character(len=64) :: comment(natomsmax) ! to store the comment line
! ***
! simplex parameters
! ftol -> convergence criterium for the function values
! xtol -> convergence criterium for the function arguments
! maxiter -> maximum number of iterations allowed
! maxfun -> maximum number of function evaluations allowed
! fopt -> minimal value of f2min
! niter -> number of iterations performed
! nfeval -> number of function evaluations performed
! plevel -> determine the amount of information to be printed (the larger the more)
! bool -> whether or not the optimization succeeded
! ***
integer :: maxiter
integer :: maxfun
integer :: plevel
double precision :: ftol
double precision :: xtol
logical :: bool ! whether or not the optimization succeeded
integer :: ndim,iter,nfeval,n_ae,n_fc
double precision :: fopt,finit!,goal= 0.0d0
logical, parameter :: logopt=.True. ! whether the parameter space's logarithm is taken or not
logical, dimension(npara) :: lns ! to take notes of the parameter's sign before log()
integer, parameter :: gbastmpf = 78
! end simplex parameters
#if !defined (gfortran) && !defined (G95)
integer*4 system
external system
#endif
! initial guess for the parameters to be optimized
ipara=0
do iatoms=1,natoms
if(lbond) then
if(ibopt(iatoms).gt.0) then
ipara=ipara+1
para(ipara)=gbondc(iatoms)
endif
endif
do iang=0,nang(iatoms)
do iprim=1,nprim(iang,iatoms)
if(ieopt(iprim,iang,iatoms).gt.0) then
ipara=ipara+1
para(ipara)=gexp(iprim,iang,iatoms)
endif
enddo
if(ebf(iang,iatoms).gt.0) then
ipara=ipara+1
para(ipara)=ell(iang,iatoms)
endif
do iprim=1,nprim(iang,iatoms)
do icontr=1,ncontr(iang,iatoms)
if(igopt(iprim,icontr,iang,iatoms).gt.0) then
ipara=ipara+1
para(ipara)=gcoef(iprim,icontr,iang,iatoms)
endif
enddo
enddo
enddo
enddo
! end initial guess for the parameters to be optimized
! taking the logarithm of the arguments if the optimization takes place
! in that space
if (logopt) call log_trans(para,lns)
!***********************************************************************
write(*,1003)
write(*,'(21x,a29)') 'Basis Set Optimization STARTS'
write(*,*) ''
finit = goalfunc(para)! initial function value
call simplex(goalfunc,para,ftol,xtol,maxiter,maxfun,
$ fopt,iter,nfeval,plevel,bool)
write(*,1003)
write(*,*) ''
write(*,'(a40,f17.8)') 'Initial function value:',finit
write(*,*) ''
write(*,'(a40,f17.8)') 'Optimal function value:',fopt
write(*,'(a42,d16.3)') 'Function tolerance:',ftol
write(*,'(a42,d16.3)') 'Argument tolerance:',xtol
write(*,'(a42,i8)') 'Number of iterations:',iter
write(*,'(a42,i8)') 'Number of function evaluations:',nfeval
write(*,*) ''
write(*,'(21x,a27)') 'Basis Set Optimization ENDS'
write(*,1003)
1003 format(1x,70('%'))
!***********************************************************************
! disk cleanup
!
c istat=system('rm -f COORD.xyz DFINV DFTGRID EXIT FOCK fort.*')
c istat=system('rm -f iface KEYWD MOCOEF MOLDEN OCCUP OEINT')
c istat=system('rm -f SCFDENSITIES SYMTRA TEDAT TEINT VARS')
c istat=system('rm -f PRINT')
c
return
************************************************************************
contains !basopt_smplx
************************************************************************
************************************************************************
double precision function goalfunc(para)
************************************************************************
implicit none
double precision, dimension(:), intent(inout) :: para
double precision :: res
res = mrccener(para)
if(dabs(fgoal).lt.1d-10) then
goalfunc = res
else
goalfunc = dabs(fgoal-res)
endif
end function goalfunc
************************************************************************
************************************************************************
double precision function mrccener(para)
************************************************************************
* it executes MRCC and returns an energy
************************************************************************
implicit none
double precision, dimension(:), intent(inout) :: para
double precision e_all,e_fc
character(len=6) stmp
************************************************************************
! taking the exponential of the arguments if the optimization takes place
! in the logarithm space of the arguments otherwise we can't get a meaningful GENBAS
if (logopt) call exp_trans(para,lns)
if(lbond) then ! Are there any bond functions?
call update_GBAS_bf(para)
else
call update_GBAS_std(para)
endif
if(lcv) then
write(stmp,'(I0)') n_ae
call changekey('core',4,stmp,6)
call spoint(.true.)
e_all = get_last_energy()
write(stmp,'(I0)') n_fc
call changekey('core',4,stmp,6)
call spoint(.true.)
e_fc = get_last_energy()
mrccener = e_all - e_fc
else
call spoint(.true.)
mrccener = get_last_energy()
endif
! transformation of the arguments back to the logarithm of the space
if (logopt) call log_trans(para,lns)
return
end function mrccener
************************************************************************
************************************************************************
subroutine log_trans(vec,lvec)
************************************************************************
* takes the logarithm of the vec's elements and initialize lvec's
* elements depending on their sign
************************************************************************
implicit none
double precision, dimension(:), intent(inout) :: vec
logical, dimension(:), intent(inout) :: lvec
! local
integer i,ndim
ndim = size(vec)
do i=1,ndim
if(vec(i).gt.0.d0) then
lvec(i)=.True.
vec(i)=log(vec(i))
else
lvec(i)=.False.
vec(i)=log(-vec(i))
endif
enddo
end subroutine log_trans
************************************************************************
************************************************************************
subroutine exp_trans(vec,lvec)
************************************************************************
* it takes the exponential of the vec's elements according to the
* states of lvec
************************************************************************
implicit none
double precision, dimension(:), intent(inout) :: vec
logical, dimension(:), intent(inout) :: lvec
! local
integer i,ndim
ndim = size(vec)
do i=1,ndim
if(lvec(i)) then
vec(i)=exp(vec(i))
else
vec(i)=-exp(vec(i))
endif
enddo
end subroutine exp_trans
************************************************************************
************************************************************************
subroutine update_GBAS_bf(para)
************************************************************************
* it updates the GENBAS file (bond function version)
*
* GENBAS.tmp needs special handling, nevertheless it can be implemented
* by a common routine (near future project)
************************************************************************
implicit none
double precision, dimension(:), intent(inout) :: para
! do i=1,npara ! test the parameters to be optimized
! write(*,*) i, para(i)
! enddo
ipara=0
open(gbasfile,file='GENBAS',status='old')
open(gbastmpf,file='GENBAS.tmp',status='unknown')
write(gbasfile,"(a64)") title(1) ! title line
write(gbasfile,"(a64)") comment(1) ! comment line
write(gbasfile,*) ! blank line
write(gbasfile,100) natoms ! number of bond centers/number of atoms
write(gbastmpf,"(a64)") title(1)
write(gbastmpf,"(a64)") comment(1)
write(gbastmpf,*)
write(gbastmpf,100) natoms
do iatoms=1,natoms
write(gbasfile,*)
write(gbastmpf,*)
if(ibopt(iatoms).gt.0) then
ipara=ipara+1
gbondc(iatoms)=para(ipara)
endif
***************|
* std. GENBAS V
write(gbasfile,101) gbondc(iatoms)
* GENBAS.tmp
if(cbsm(iatoms).ne.' ') then
write(gbastmpf,103,advance='no')
$ gbondc(iatoms),cbsm(iatoms)
else
write(gbastmpf,103,advance='no')
$ gbondc(iatoms)
endif
write(gbastmpf,*)
****************
write(gbasfile,*)
write(gbasfile,100) nang(iatoms)+1
write(gbastmpf,*)
write(gbastmpf,100) nang(iatoms)+1
***************|
* std. GENBAS V
write(gbasfile,100) (angt(iang,iatoms),iang=0,nang(iatoms))
* GENBAS.tmp
do iang=0,nang(iatoms)
if(casm(iang,iatoms).ne.' ') then
write(gbastmpf,102,advance='no')
$ angt(iang,iatoms),casm(iang,iatoms)
else
write(gbastmpf,100,advance='no')
$ angt(iang,iatoms)
endif
enddo
write(gbastmpf,*)
****************
write(gbasfile,100) (ncontr(iang,iatoms),iang=0,nang(iatoms))
write(gbasfile,100) (nprim (iang,iatoms),iang=0,nang(iatoms))
write(gbastmpf,100) (ncontr(iang,iatoms),iang=0,nang(iatoms))
write(gbastmpf,100) (nprim (iang,iatoms),iang=0,nang(iatoms))
do iang=0,nang(iatoms)
write(gbasfile,*)
write(gbastmpf,*)
do iprim=1,nprim(iang,iatoms)
if(ieopt(iprim,iang,iatoms).gt.0) then
ipara=ipara+1
gexp(iprim,iang,iatoms)=para(ipara)
endif
enddo
***************|
* std. GENBAS V
write(gbasfile,101)
$(gexp(iprim,iang,iatoms),iprim=1,nprim(iang,iatoms))
* GENBAS.tmp
do iprim=1,nprim(iang,iatoms)
if(cesm(iprim,iang,iatoms).ne.' ') then
write(gbastmpf,103,advance='no')
$ gexp(iprim,iang,iatoms),cesm(iprim,iang,iatoms)
else
write(gbastmpf,101,advance='no')
$ gexp(iprim,iang,iatoms)
endif
enddo
write(gbastmpf,*)
****************
! if(ebf(iang,iatoms).gt.0) then
! ipara=ipara+1
! write(gbasfile,101) para(ipara)
! write(gbastmpf,101) para(ipara)
! else
write(gbasfile,*)
write(gbastmpf,*)
! endif
do iprim=1,nprim(iang,iatoms)
do icontr=1,ncontr(iang,iatoms)
if(igopt(iprim,icontr,iang,iatoms).gt.0) then
ipara=ipara+1
gcoef(iprim,icontr,iang,iatoms)=para(ipara)
endif
enddo
***************|
* std. GENBAS V
write(gbasfile,101)
$(gcoef(iprim,icontr,iang,iatoms),icontr=1,ncontr(iang,iatoms))
* GENBAS.tmp
do icontr=1,ncontr(iang,iatoms)
if(cgsm(iprim,icontr,iang,iatoms).ne.' ') then
write(gbastmpf,103, advance='no')
$ gcoef(iprim,icontr,iang,iatoms),
$ cgsm(iprim,icontr,iang,iatoms)
else
write(gbastmpf,101, advance='no')
$ gcoef(iprim,icontr,iang,iatoms)
endif
enddo
write(gbastmpf,*)
****************
enddo
enddo
enddo
close(gbasfile)
close(gbastmpf)
100 format(1000i5)
101 format(1000f22.12)
102 format(i5,a2)
103 format(f22.12,a2)
! end update GENBAS
end subroutine update_GBAS_bf
************************************************************************
************************************************************************
subroutine update_GBAS_std(para)
************************************************************************
* it updates the GENBAS file (standard version)
*
* GENBAS.tmp needs special handling, nevertheless it can be implemented
* by a common routine (near future project)
************************************************************************
implicit none
double precision, dimension(:), intent(inout) :: para
! do i=1,npara ! test the parameters to be optimized
! write(*,*) i, para(i)
! enddo
open(gbasfile,file='GENBAS',status='old')
open(gbastmpf,file='GENBAS.tmp',status='unknown')
ipara=0
do iatoms=1,natoms
write(gbasfile,"(a64)") title(iatoms) ! title line
write(gbasfile,"(a64)") comment(iatoms) ! comment line
write(gbasfile,*) ! blank line
write(gbasfile,100) nang(iatoms)+1 ! number of angular momentums
write(gbastmpf,"(a64)") title(iatoms)
write(gbastmpf,"(a64)") comment(iatoms)
write(gbastmpf,*)
write(gbastmpf,100) nang(iatoms)+1
***************|
* std. GENBAS V
write(gbasfile,100) (angt(iang,iatoms),iang=0,nang(iatoms)) ! types of angular momentum, s->0, p->1, etc
* GENBAS.tmp
do iang=0,nang(iatoms)
if(casm(iang,iatoms).ne.' ') then
write(gbastmpf,102,advance='no')
$ angt(iang,iatoms),casm(iang,iatoms)
else
write(gbastmpf,100,advance='no')
$ angt(iang,iatoms)
endif
enddo
write(gbastmpf,*)
****************
write(gbasfile,100) (ncontr(iang,iatoms),iang=0,nang(iatoms)) ! number of contracted basis functions for the given ang. mom.
write(gbasfile,100) (nprim (iang,iatoms),iang=0,nang(iatoms)) ! number of primitive basis functions for the given ang. mom.
write(gbastmpf,100) (ncontr(iang,iatoms),iang=0,nang(iatoms))
write(gbastmpf,100) (nprim (iang,iatoms),iang=0,nang(iatoms))
do iang=0,nang(iatoms)
write(gbasfile,*)
write(gbastmpf,*)
do iprim=1,nprim(iang,iatoms)
if(ieopt(iprim,iang,iatoms).gt.0) then
ipara=ipara+1
gexp(iprim,iang,iatoms)=para(ipara)
endif
enddo
***************|
* std. GENBAS V
write(gbasfile,101)
$(gexp(iprim,iang,iatoms),iprim=1,nprim(iang,iatoms))
* GENBAS.tmp
do iprim=1,nprim(iang,iatoms)
if(cesm(iprim,iang,iatoms).ne.' ') then
write(gbastmpf,103,advance='no')
$ gexp(iprim,iang,iatoms),cesm(iprim,iang,iatoms)
else
write(gbastmpf,101,advance='no')
$ gexp(iprim,iang,iatoms)
endif
enddo
write(gbastmpf,*)
****************
write(gbasfile,*)
write(gbastmpf,*)
do iprim=1,nprim(iang,iatoms)
do icontr=1,ncontr(iang,iatoms)
if(igopt(iprim,icontr,iang,iatoms).gt.0) then
ipara=ipara+1
gcoef(iprim,icontr,iang,iatoms)=para(ipara)
endif
enddo
***************|
* std. GENBAS V
write(gbasfile,101)
$(gcoef(iprim,icontr,iang,iatoms),icontr=1,ncontr(iang,iatoms))
* GENBAS.tmp
do icontr=1,ncontr(iang,iatoms)
if(cgsm(iprim,icontr,iang,iatoms).ne.' ') then
write(gbastmpf,103, advance='no')
$ gcoef(iprim,icontr,iang,iatoms),
$ cgsm(iprim,icontr,iang,iatoms)
else
write(gbastmpf,101, advance='no')
$ gcoef(iprim,icontr,iang,iatoms)
endif
enddo
write(gbastmpf,*)
****************
enddo
enddo
write(gbasfile,*)
write(gbastmpf,*)
enddo
close(gbasfile)
close(gbastmpf)
100 format(1000i5)
101 format(1000f22.12)
102 format(i5,a2)
103 format(f22.12,a2)
! end update GENBAS
end subroutine update_GBAS_std
************************************************************************
************************************************************************
subroutine set_MINP(n)
************************************************************************
* it interchanges the appropriate MINP files for CV basis set optimization
*
* MINP.cv core electrons are correlated
* MINP.fc core electrons are frozen
*
************************************************************************
implicit none
integer, intent(in):: n
integer :: istat
if(n .eq. 0) then
istat=system('cp MINP.cv MINP')
elseif (n .eq. 1) then
istat=system('cp MINP.fc MINP')
else
write(6,*) 'Non valid integer in set_MINP (basopt.f)'
endif
if(istat .ne. 0) then
write(6,*) 'Error while copying MINP.cv/MINP.fc'
stop
endif
end subroutine set_MINP
************************************************************************
************************************************************************
C
end subroutine basopt_smplx
C
************************************************************************
************************************************************************
subroutine basoptcore(npara,gbasfile,natoms,lbond,gbondc,nang,
$nangmax,nprimmax,ncontrmax,angt,nprim,ncontr,igopt,gexp,gcoef,
$ebf,ell)
************************************************************************
implicit none
integer npara,ispec,ireac,gbasfile,i
integer ipara,jpara,iparat,nit,ipiv(npara),mnit,maxit,mmaxit
real*8 para(npara),temp,fpm,grd,hss,moldener,oldener
real*8 paraold(npara),eps,grad(npara),htokj,conv,f00,fpp,fmm,fmp
real*8 fm(npara),fp(npara),hess(npara,npara),param(npara),maxstep
real*8 scr(3*npara**2),hess2(npara,npara),eig(npara)
real*8 paraold2(npara),mconv,gtol,mmaxstep,maxstepin
parameter(eps=0.01d0,maxstepin=0.1d0,gtol=1e-4,maxit=100)
parameter(mmaxit=maxit)
integer*4 istat,isyev
C
integer natoms,iatoms,iang,iprim,icontr,nangmax,nprimmax,ncontrmax
integer nang(*),nprim(0:nangmax,*),ebf(0:nangmax,*)
integer ncontr(0:nangmax,*),angt(0:nangmax,*)
integer igopt(nprimmax,ncontrmax,0:nangmax,*)
real*8 gexp(nprimmax,0:nangmax,*),gbondc(*),ell(0:nangmax,*)
real*8 gcoef(nprimmax,ncontrmax,0:nangmax,*)
logical lbond
c#if !defined (gfortran) && !defined (G95)
c integer*4 system
c external system
c#endif
C Initialize parameters
ipara=0
do iatoms=1,natoms
if(lbond) then
ipara=ipara+1
para(ipara)=gbondc(iatoms)
endif
do iang=0,nang(iatoms)
do iprim=1,nprim(iang,iatoms)
ipara=ipara+1
para(ipara)=gexp(iprim,iang,iatoms)
enddo
if(ebf(iang,iatoms).gt.0) then
ipara=ipara+1
para(ipara)=ell(iang,iatoms)
endif
do iprim=1,nprim(iang,iatoms)
do icontr=1,ncontr(iang,iatoms)
if(igopt(iprim,icontr,iang,iatoms).gt.0) then
ipara=ipara+1
para(ipara)=gcoef(iprim,icontr,iang,iatoms)
endif
enddo
enddo
enddo
enddo
C Iteration
conv=1d30
nit=0
maxstep=maxstepin
do while((conv.gt.0.05d0.or.eig(1).le.0.d0).and.nit.lt.maxit)
nit=nit+1
paraold(1:npara)=para(1:npara)
C Calculate function
call calcener(npara,gbasfile,para,f00,natoms,lbond,gbondc,nang,
$nangmax,nprimmax,ncontrmax,angt,nprim,ncontr,igopt,gexp,gcoef,ebf)
if(nit.eq.1) oldener=f00
C Calculate gradient
do ipara=1,npara
param(1:npara)=para(1:npara)
param(ipara)=param(ipara)+eps
call calcener(npara,gbasfile,param,fp(ipara),natoms,lbond,
$gbondc,nang,nangmax,nprimmax,ncontrmax,angt,nprim,ncontr,igopt,
$gexp,gcoef,ebf)
param(1:npara)=para(1:npara)
param(ipara)=param(ipara)-eps
call calcener(npara,gbasfile,param,fm(ipara),natoms,lbond,
$gbondc,nang,nangmax,nprimmax,ncontrmax,angt,nprim,ncontr,igopt,
$gexp,gcoef,ebf)
grad(ipara)=(fp(ipara)-fm(ipara))/(2.d0*eps)
c write(6,*) 'grad',fp(ipara),fm(ipara),grad(ipara)
enddo
C Calculate Hessian
do ipara=1,npara
hess(ipara,ipara)=(fp(ipara)-2.d0*f00+fm(ipara))/eps**2
do jpara=1,ipara-1
param(1:npara)=para(1:npara)
param(ipara)=param(ipara)+eps
param(jpara)=param(jpara)+eps
call calcener(npara,gbasfile,param,fpp,natoms,lbond,gbondc,
$nang,nangmax,nprimmax,ncontrmax,angt,nprim,ncontr,igopt,gexp,
$gcoef,ebf)
param(1:npara)=para(1:npara)
param(ipara)=param(ipara)-eps
param(jpara)=param(jpara)-eps
call calcener(npara,gbasfile,param,fmm,natoms,lbond,gbondc,
$nang,nangmax,nprimmax,ncontrmax,angt,nprim,ncontr,igopt,gexp,
$gcoef,ebf)
param(1:npara)=para(1:npara)
param(ipara)=param(ipara)+eps
param(jpara)=param(jpara)-eps
call calcener(npara,gbasfile,param,fpm,natoms,lbond,gbondc,
$nang,nangmax,nprimmax,ncontrmax,angt,nprim,ncontr,igopt,gexp,
$gcoef,ebf)
param(1:npara)=para(1:npara)
param(ipara)=param(ipara)-eps
param(jpara)=param(jpara)+eps
call calcener(npara,gbasfile,param,fmp,natoms,lbond,gbondc,
$nang,nangmax,nprimmax,ncontrmax,angt,nprim,ncontr,igopt,gexp,
$gcoef,ebf)
hess(ipara,jpara)=(fpp+fmm-fmp-fpm)/(4.d0*eps**2)
hess(jpara,ipara)=hess(ipara,jpara)
enddo
enddo
C Diagonalize Hessian
hess2(1:npara,1:npara)=hess(1:npara,1:npara)
call dsyev('V','L',npara,hess2,npara,eig,scr,3*npara**2,isyev)
if(isyev.ne.0) then
write(6,*) 'Fatal error at diagonalization!'
stop
endif
write(6,"(' Hessian eigenvalues:')")
write(6,*) eig
C Line search
if(eig(1).le.0.d0) then
write(6,*) 'Negative eigenvalue of the Hessian:',eig(1)
mconv=1d30
mnit=0
mmaxstep=maxstep
moldener=1d30
do while(mconv.gt.0.05d0.and.mnit.lt.mmaxit)
mnit=mnit+1
paraold2(1:npara)=para(1:npara)
C Calculate function
call calcener(npara,gbasfile,para,f00,natoms,lbond,gbondc,
$nang,nangmax,nprimmax,ncontrmax,angt,nprim,ncontr,igopt,gexp,
$gcoef,ebf)
C Calculate gradient
param(1:npara)=para(1:npara)+eps*hess2(1:npara,1)
call calcener(npara,gbasfile,param,fpp,natoms,lbond,gbondc,
$nang,nangmax,nprimmax,ncontrmax,angt,nprim,ncontr,igopt,gexp,
$gcoef,ebf)
param(1:npara)=para(1:npara)-eps*hess2(1:npara,1)
call calcener(npara,gbasfile,param,fmm,natoms,lbond,gbondc,
$nang,nangmax,nprimmax,ncontrmax,angt,nprim,ncontr,igopt,gexp,
$gcoef,ebf)
grd=(fpp-fmm)/(2.d0*eps)
C Calculate Hessian
hss=(fpp-2.d0*f00+fmm)/eps**2
write(6,*) grd,hss
C Update parameters
if(hss.le.0.d0.or.dabs(grd/hss).gt.0.5d0) then
para(1:npara)=
$para(1:npara)-sign(mmaxstep,grd)*hess2(1:npara,1)
else
para(1:npara)=para(1:npara)-
$sign(min(mmaxstep,dabs(grd/hss)),grd/hss)*hess2(1:npara,1)
endif
C Calculate difference
mconv=0.d0
do ipara=1,npara
mconv=max(mconv,dabs(paraold2(ipara)-para(ipara)))
enddo
if(moldener.lt.f00) then
mmaxstep=mmaxstep/2.d0
para(1:npara)=paraold2(1:npara)
else
moldener=f00
endif
write(6,
$"(' Microiteration',i3,3x,'Convergence:',f10.6,3x,
$'Energy:',f12.6)") mnit,mconv,f00
enddo
else
C Newton step
hess2(1:npara,1:npara)=hess(1:npara,1:npara)
param(1:npara)=-grad(1:npara)
call dsysv('l',npara,1,hess2,npara,ipiv,param,npara,scr,
$3*npara**2,i)
if(i.ne.0) then
write(6,*) 'Fatal error at the inversion of Hessian!'
write(6,*) 'Error code:',i
write(6,*) 'Eigenvalues:',eig
stop
endif
C Update parameters
do ipara=1,npara
para(ipara)=para(ipara)+
$sign(1.d0,param(ipara))*min(maxstep,dabs(param(ipara)))
enddo
endif
C Calculate difference
conv=0.d0
do ipara=1,npara
conv=max(conv,dabs(paraold(ipara)-para(ipara)))
enddo
C Print info
call calcener(npara,gbasfile,para,f00,natoms,lbond,gbondc,nang,
$nangmax,nprimmax,ncontrmax,angt,nprim,ncontr,igopt,gexp,gcoef,ebf)
if(oldener.lt.f00) then
maxstep=maxstep/2.d0
para(1:npara)=paraold(1:npara)
f00=oldener
else
oldener=f00
endif
write(6,"(' Iteration',i3,3x,'Convergence:',f10.6,3x,'Energy:',
$f12.6)") nit,conv,f00
write(6,"(' Parameter Gradient')")
do ipara=1,npara
write(6,"(7x,f10.4,2x,f10.4)") para(ipara),grad(ipara)
enddo
enddo
c istat=system('rm -f COORD.xyz DFINV DFTGRID EXIT FOCK fort.*')
c istat=system('rm -f iface KEYWD MOCOEF MOLDEN OCCUP OEINT')
c istat=system('rm -f SCFDENSITIES SYMTRA TEDAT TEINT VARS')
c istat=system('rm -f PRINT')
C
return
end subroutine basoptcore
C
************************************************************************
subroutine calcener(npara,gbasfile,para,ener,natoms,lbond,gbondc,
$nang,nangmax,nprimmax,ncontrmax,angt,nprim,ncontr,igopt,gexp,
$gcoef,ebf)
************************************************************************
* Calculate energies
************************************************************************
use optim, only: get_last_energy
implicit none
integer npara,gbasfile,ipara,iparat,ispec,i
integer ireac
real*8 para(npara),ener
character*64 c64(2)
integer*4 istat
C
integer natoms,iatoms,iang,iprim,icontr,nangmax,nprimmax,ncontrmax
integer nang(*),nprim(0:nangmax,*),ebf(0:nangmax,*)
integer ncontr(0:nangmax,*),angt(0:nangmax,*)
integer igopt(nprimmax,ncontrmax,0:nangmax,*)
real*8 gexp(nprimmax,0:nangmax,*),gbondc(*)
real*8 gcoef(nprimmax,ncontrmax,0:nangmax,*)
logical lbond
c#if !defined (gfortran) && !defined (G95)
c integer*4 system
c external system
c#endif
C Update GENBAS file
open(gbasfile,file='GENBAS',status='old')
read(gbasfile,"(a64)") c64(1)
read(gbasfile,"(a64)") c64(2)
c64(1)=trim(adjustl(c64(1)))
c64(2)=trim(adjustl(c64(2)))
rewind(gbasfile)
write(gbasfile,"(a64)") c64(1)
write(gbasfile,"(a64)") c64(2)
if(lbond) then
write(gbasfile,*)
write(gbasfile,100) natoms
write(gbasfile,*)
endif
ipara=0
do iatoms=1,natoms
if(lbond) then
ipara=ipara+1
write(gbasfile,101) para(ipara)
endif
write(gbasfile,*)
write(gbasfile,100) nang(iatoms)+1
write(gbasfile,100) (angt(iang,iatoms),iang=0,nang(iatoms))
write(gbasfile,100) (ncontr(iang,iatoms),iang=0,nang(iatoms))
write(gbasfile,100) (nprim (iang,iatoms),iang=0,nang(iatoms))
do iang=0,nang(iatoms)
write(gbasfile,*)
write(gbasfile,101)
$(para(iprim),iprim=ipara+1,ipara+nprim(iang,iatoms))
ipara=ipara+nprim(iang,iatoms)
if(ebf(iang,iatoms).gt.0) then
ipara=ipara+1
write(gbasfile,101) para(ipara)
else
write(gbasfile,*)
endif
do iprim=1,nprim(iang,iatoms)
do icontr=1,ncontr(iang,iatoms)
if(igopt(iprim,icontr,iang,iatoms).gt.0) then
ipara=ipara+1
gcoef(iprim,icontr,iang,iatoms)=para(ipara)
endif
enddo
write(gbasfile,101)
$(gcoef(iprim,icontr,iang,iatoms),icontr=1,ncontr(iang,iatoms))
enddo
enddo
write(gbasfile,*)
enddo
close(gbasfile)
100 format(1000i5)
101 format(1000f22.12)
C Execute mrcc
c istat=system('dmrcc > /dev/null')
c if(istat.ne.0) then
c write(6,*) 'Fatal error in dmrcc!'
c stop
c endif
call spoint(.true.)
ener = get_last_energy()
write(6,*) ener
C
return
end subroutine calcener
C