Skip to content

Commit

Permalink
Shortcircuit Hessian evaluation for Gaussian output
Browse files Browse the repository at this point in the history
  • Loading branch information
Cyrille Lavigne committed Oct 4, 2020
1 parent 18f768c commit 4e73b63
Show file tree
Hide file tree
Showing 5 changed files with 29 additions and 15 deletions.
6 changes: 6 additions & 0 deletions src/hessian.F90
Original file line number Diff line number Diff line change
Expand Up @@ -504,6 +504,12 @@ subroutine numhess( &
write(env%unit, '(A)') "DFTB+ style hessian.out written"
end if

! If we are currently computing a Hessian for Gaussian to use, we return it
! right here and now and skip all the post-processing.
if (geometry_inputfile .eq. p_geo_gaussian) then
return
end if

! prepare all for diag
! copy
k=0
Expand Down
3 changes: 3 additions & 0 deletions src/io/reader/gaussian.f90
Original file line number Diff line number Diff line change
Expand Up @@ -29,6 +29,7 @@ module xtb_io_reader_gaussian

subroutine readMoleculeGaussianExternal(mol, unit, status, iomsg)
use xtb_mctc_accuracy, only : wp
use xtb_setmod, only : set_geopref
logical, parameter :: debug = .false.
class(TMolecule), intent(out) :: mol
integer, intent(in) :: unit
Expand All @@ -44,6 +45,8 @@ subroutine readMoleculeGaussianExternal(mol, unit, status, iomsg)

status = .false.

call set_geopref('gaussian')

read(unit, '(4i10)', iostat=err) n, mode, chrg, spin
if (err.ne.0) then
iomsg = "Could not read number of atoms, check format!"
Expand Down
23 changes: 12 additions & 11 deletions src/prog/main.F90
Original file line number Diff line number Diff line change
Expand Up @@ -265,7 +265,7 @@ subroutine xtbMain(env, argParser)
call close_file(ich)
end if
endif

!> efield read: gfnff only
call open_file(ich,'.EFIELD','r')
if (ich.ne.-1) then
Expand Down Expand Up @@ -321,7 +321,7 @@ subroutine xtbMain(env, argParser)
fname = 'caffeine'
call get_coffee(mol)
call generateFileMetaInfo(fname, directory, basename, extension)
else
else
call generateFileMetaInfo(fname, directory, basename, extension)
ftype = getFileType(basename, extension)
call open_file(ich, fname, 'r')
Expand Down Expand Up @@ -584,7 +584,7 @@ subroutine xtbMain(env, argParser)
end select
call env%checkpoint("Single point calculation terminated")

!> write 2d => 3d converted structure
!> write 2d => 3d converted structure
if (struc_conversion_done) then
call generateFileName(tmpname, 'gfnff_convert', extension, mol%ftype)
write(env%unit,'(10x,a,1x,a,/)') &
Expand All @@ -593,7 +593,7 @@ subroutine xtbMain(env, argParser)
call writeMolecule(mol, ich, energy=res%e_total, gnorm=res%gnorm)
call close_file(ich)
end if

! ========================================================================
!> determine kopt for bhess including final biased geometry optimization
if (runtyp.eq.p_run_bhess) then
Expand Down Expand Up @@ -857,7 +857,8 @@ subroutine xtbMain(env, argParser)
call main_geometry(iprop,mol)
endif

if ((runtyp.eq.p_run_hess).or.(runtyp.eq.p_run_ohess).or.(runtyp.eq.p_run_bhess)) then
if (((runtyp.eq.p_run_hess).or.(runtyp.eq.p_run_ohess).or.(runtyp.eq.p_run_bhess)) &
.and.(geometry_inputfile.ne.p_geo_gaussian)) then
call generic_header(iprop,'Frequency Printout',49,10)
call main_freq(iprop,mol,chk%wfn,fres)
endif
Expand Down Expand Up @@ -888,7 +889,7 @@ subroutine xtbMain(env, argParser)
type is(TGFFCalculator)
call write_energy_gff(env%unit,res,fres, &
& (runtyp.eq.p_run_hess).or.(runtyp.eq.p_run_ohess).or.(runtyp.eq.p_run_bhess))
end select
end select


! ------------------------------------------------------------------------
Expand Down Expand Up @@ -1278,10 +1279,10 @@ subroutine parseArguments(env, args, inputFile, paramFile, accuracy, lgrad, &
call set_exttyp('eht')
call env%warning("The use of '"//flag//"' is discouraged, " //&
& "please use '--gfn 0' next time", source)

case('--gfnff')
call set_exttyp('ff')

case('--gff')
call set_exttyp('ff')

Expand Down Expand Up @@ -1315,7 +1316,7 @@ subroutine parseArguments(env, args, inputFile, paramFile, accuracy, lgrad, &

case('--json')
call set_write(env,'json','true')

case('--ceasefiles')
restart = .false.
verbose=.false.
Expand All @@ -1326,7 +1327,7 @@ subroutine parseArguments(env, args, inputFile, paramFile, accuracy, lgrad, &
call set_opt(env, 'logfile', 'NUL')
#else
call set_opt(env, 'logfile', '/dev/null')
#endif
#endif

case('--orca')
call set_exttyp('orca')
Expand Down Expand Up @@ -1461,7 +1462,7 @@ subroutine parseArguments(env, args, inputFile, paramFile, accuracy, lgrad, &
if (allocated(sec)) then
call set_opt(env,'optlevel',sec)
endif

case('--bhess')
call set_runtyp('bhess')
call args%nextArg(sec)
Expand Down
3 changes: 3 additions & 0 deletions src/set_module.f90
Original file line number Diff line number Diff line change
Expand Up @@ -947,6 +947,9 @@ subroutine set_geopref(typ)
case('vasp','poscar')
geometry_inputfile = p_geo_poscar

case('EIn', 'gaussian')
geometry_inputfile = p_geo_gaussian

end select
set = .false.
end subroutine set_geopref
Expand Down
9 changes: 5 additions & 4 deletions src/setparam.f90
Original file line number Diff line number Diff line change
Expand Up @@ -59,10 +59,11 @@ module xtb_setparam
logical :: periodic = .false.

! Geometry input type
integer,parameter :: p_geo_coord = 1
integer,parameter :: p_geo_xmol = 2
integer,parameter :: p_geo_sdf = 3
integer,parameter :: p_geo_poscar = 4
integer,parameter :: p_geo_coord = 1
integer,parameter :: p_geo_xmol = 2
integer,parameter :: p_geo_sdf = 3
integer,parameter :: p_geo_poscar = 4
integer,parameter :: p_geo_gaussian = 5
integer :: geometry_inputfile = p_geo_coord

!! ------------------------------------------------------------------------
Expand Down

0 comments on commit 4e73b63

Please sign in to comment.