From 89d60ab74157969f2291f1eafda2b348f992d129 Mon Sep 17 00:00:00 2001 From: Philipp Pracht <42216327+pprcht@users.noreply.github.com> Date: Wed, 28 Feb 2024 13:54:43 +0000 Subject: [PATCH] Fix for making ALPB h2o compatible with tblite calculator (#271) * tblite apparently only knew "alpb water", but not "alpb h2o". fixed that --------- Signed-off-by: Philipp Pracht --- src/calculator/tblite_api.F90 | 17 ++++++++++++----- 1 file changed, 12 insertions(+), 5 deletions(-) diff --git a/src/calculator/tblite_api.F90 b/src/calculator/tblite_api.F90 index 09e4329b..e638d641 100644 --- a/src/calculator/tblite_api.F90 +++ b/src/calculator/tblite_api.F90 @@ -173,7 +173,7 @@ subroutine tblite_add_solv(mol,chrg,uhf,ctx,wfn,tbcalc,smodel,solvent) class(tblite_solvation_type),allocatable :: solv type(solvation_input),allocatable :: solv_inp type(solvent_data) :: solv_data - character(len=:),allocatable :: str + character(len=:),allocatable :: str,solvdum logical :: pr real(wp) :: etemp_au,energy @@ -182,6 +182,8 @@ subroutine tblite_add_solv(mol,chrg,uhf,ctx,wfn,tbcalc,smodel,solvent) return end if + !> special case: tblite dosen't know 'h2o', only 'water' + pr = (ctx%verbosity > 0) !>--- make an mctcmol object from mol @@ -189,7 +191,12 @@ subroutine tblite_add_solv(mol,chrg,uhf,ctx,wfn,tbcalc,smodel,solvent) if(pr) call ctx%message("tblite> setting up tblite implicit solvation") !>--- generat solvation parametrization - solv_data = get_solvent_data(solvent) + if(solvent == 'h2o')then !> special case: tblite doesn't know 'h2o', only 'water' ... + solvdum = 'water' + else + solvdum = solvent + endif + solv_data = get_solvent_data(solvdum) if (solv_data%eps <= 0.0_wp) then if(pr) call ctx%message("tblite> Unknown solvent!") return @@ -197,15 +204,15 @@ subroutine tblite_add_solv(mol,chrg,uhf,ctx,wfn,tbcalc,smodel,solvent) allocate (solv_inp) select case (trim(smodel)) case ('gbsa') - if(pr) call ctx%message("tblite> using GBSA/"//solvent) + if(pr) call ctx%message("tblite> using GBSA/"//solvdum) allocate (solv_inp%alpb) solv_inp%alpb = alpb_input(solv_data%eps,alpb=.false.) case ('cpcm') - if(pr) call ctx%message("tblite> using CPCM/"//solvent) + if(pr) call ctx%message("tblite> using CPCM/"//solvdum) allocate (solv_inp%cpcm) solv_inp%cpcm = cpcm_input(solv_data%eps) case ('alpb') - if(pr) call ctx%message("tblite> using ALPB/"//solvent) + if(pr) call ctx%message("tblite> using ALPB/"//solvdum) allocate (solv_inp%alpb) solv_inp%alpb = alpb_input(solv_data%eps,alpb=.true.) case default