Skip to content

Commit 5aa1015

Browse files
committed
WIP: Implement SHwIT method
1 parent 0750f7a commit 5aa1015

3 files changed

Lines changed: 47 additions & 26 deletions

File tree

src/sh_integ.F90

Lines changed: 15 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -10,6 +10,7 @@ module mod_sh_integ
1010
public :: sh_set_initialwf, sh_set_energy_shift
1111
public :: sh_select_integrator, sh_integrate_wf
1212
public :: sh_decoherence_correction, check_popsum, sh_TFS_transmat
13+
public :: shwit_switch
1314

1415
! Restart functionality
1516
public :: sh_write_phase_bin, sh_read_phase_bin
@@ -445,6 +446,20 @@ subroutine sh_set_initialwf(initial_state)
445446
cel_re(initial_state) = 1.0D0
446447
end subroutine sh_set_initialwf
447448

449+
! A horrible hack for SHwIT, we just switch
450+
! switch the coefficients between S1-S0
451+
! TODO: Do something less gross
452+
subroutine shwit_switch()
453+
real(DP) :: tmp
454+
tmp = cel_re(1)
455+
cel_re(1) = cel_re(2)
456+
cel_re(2) = tmp
457+
458+
tmp = cel_im(1)
459+
cel_im(1) = cel_im(2)
460+
cel_im(2) = tmp
461+
end subroutine shwit_switch
462+
448463
subroutine sh_set_energy_shift(potential_energy)
449464
real(DP), intent(in) :: potential_energy
450465

src/surfacehop.F90

Lines changed: 27 additions & 22 deletions
Original file line numberDiff line numberDiff line change
@@ -77,7 +77,12 @@ module mod_sh
7777
! Stop simulation if total energy drift since the beginning exceeds energydriftthr
7878
real(DP) :: energydriftthr = 1.0D0
7979

80-
! Special case for adiabatic TDDFT, terminate when close to S1-S0 crossing
80+
! Surface Hopping with Induced Transitions,
81+
! typically used with methods that can't describe S1S0 intersections,
82+
! such as TDDFT or ADC2.
83+
! In this case, we force a hop to S0 when the energy difference
84+
! between S1 and S0 drops below user defined threshold.
85+
logical :: SHwIT = .false.
8186
real(DP) :: dE_S0S1_thr = 0.0D0 !eV
8287
! NA Couplings
8388
real(DP), allocatable :: nacx(:, :, :), nacy(:, :, :), nacz(:, :, :)
@@ -100,7 +105,7 @@ module mod_sh
100105

101106
namelist /sh/ istate_init, nstate, substep, deltae, integ, couplings, nohop, phase, decoh_alpha, popthr, ignore_state, &
102107
nac_accu1, nac_accu2, popsumthr, energydifthr, energydriftthr, velocity_rescaling, revmom, &
103-
dE_S0S1_thr, correct_decoherence
108+
shwit, dE_S0S1_thr, correct_decoherence
104109
save
105110

106111
contains
@@ -252,7 +257,7 @@ subroutine check_sh_parameters()
252257
write (stdout, '(A)') 'Using approximate Baeck-An couplings.'
253258
case ('none')
254259
inac = 2
255-
write (stdout, '(A)') 'Ignoring nonadaiabatic couplings.'
260+
write (stdout, '(A)') 'Ignoring nonadiabatic couplings.'
256261
case default
257262
write (stderr, '(A)') 'Parameter "couplings" must be "analytic", "baeck-an" or "none".'
258263
error = .true.
@@ -279,7 +284,7 @@ subroutine check_sh_parameters()
279284
error = .true.
280285
end if
281286

282-
if (inac == 2 .and. nohop == 0) then
287+
if (inac == 2 .and. nohop == 0 .and. .not. shwit) then
283288
write (stdout, '(A)') 'WARNING: For simulations without couplings(="none") hopping probability cannot be determined.'
284289
nohop = 1
285290
end if
@@ -821,6 +826,12 @@ subroutine surfacehop(x, y, z, vx, vy, vz, vx_old, vy_old, vz_old, dt, eclas)
821826
!itp loop
822827
end do
823828

829+
if (SHwIT) then
830+
call shwit_check(en_array(1), en_array(2), dE_S0S1_thr, &
831+
vx, vy, vz, eclas)
832+
end if
833+
834+
824835
! set tocalc array for the next step
825836
call set_tocalc()
826837

@@ -1211,15 +1222,6 @@ subroutine check_energy(vx_old, vy_old, vz_old, vx, vy, vz)
12111222
real(DP), intent(in) :: vx_old(:, :), vy_old(:, :), vz_old(:, :)
12121223
real(DP) :: ekin, ekin_old, entot, entot_old
12131224

1214-
! Special case for running AIMD with TDDFT:
1215-
! End the simulation when S1-S0 energy difference drops
1216-
! below a certain small threshold.
1217-
if (dE_S0S1_thr > 0.0D0 .and. nstate >= 2) then
1218-
call check_S0S1_gap(en_S0=en_array(1), &
1219-
en_S1=en_array(2), &
1220-
threshold_ev=dE_S0S1_thr)
1221-
end if
1222-
12231225
ekin = ekin_v(vx, vy, vz)
12241226
ekin_old = ekin_v(vx_old, vy_old, vz_old)
12251227

@@ -1242,26 +1244,29 @@ subroutine check_energy(vx_old, vy_old, vz_old, vx, vy, vz)
12421244

12431245
end subroutine check_energy
12441246

1245-
subroutine check_S0S1_gap(en_S0, en_S1, threshold_ev)
1246-
use mod_general, only: STOP_SIMULATION
1247+
subroutine shwit_check(en_S0, en_S1, threshold_ev, vx, vy, vz, eclas)
12471248
use mod_const, only: AUtoEV
12481249
real(DP), intent(in) :: en_S0, en_S1
12491250
real(DP), intent(in) :: threshold_ev
1251+
real(DP), intent(inout) :: vx(:, :), vy(:, :), vz(:, :), eclas
12501252
real(DP) :: dE_S0S1
12511253

1254+
! We only hop from S1 to S0
1255+
if (istate /= 2) return
1256+
12521257
dE_S0S1 = (en_S1 - en_S0) * AUtoEV
12531258

12541259
if (dE_S0S1 < threshold_ev) then
1255-
write (stdout, *) 'S1 - S0 gap dropped below threshold!'
1260+
write (stdout, *) 'SHwIT: S1 - S0 gap dropped below threshold!'
12561261
write (stdout, *) dE_S0S1, ' < ', threshold_ev
1257-
! TODO: This is an unfortunate hack,
1258-
! but we can't directly stop here since we need the last step
1259-
! to be written to the output files.
1260-
write (stdout, *) 'Simulation will stop at the end of this step.'
1262+
write (stdout, *) "Let's jump to S0 and continue!"
12611263
! This global flag is checked in the main loop in abin.F90
1262-
STOP_SIMULATION = .true.
1264+
! STOP_SIMULATION = .true.
1265+
call try_hop_simple_rescale(vx, vy, vz, 2, 1, eclas)
1266+
! Horrible hack to exchange c_el coefficients between S1 - S0
1267+
call shwit_switch()
12631268
end if
1264-
end subroutine check_S0S1_gap
1269+
end subroutine shwit_check
12651270

12661271
subroutine check_energydrift(vx, vy, vz)
12671272
use mod_const, only: AUtoEV

tests/SH_S0S1/input.in

Lines changed: 5 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -19,17 +19,18 @@ inose=0,
1919
/
2020

2121
&sh
22+
shwit=.true.
2223
dE_S0S1_thr=9.0
23-
istate_init=3,
24+
istate_init=2,
2425
nstate=3,
2526
substep=100,
2627
deltae=100.,
2728
integ='butcher',
28-
couplings='analytic', ! non-adiabatic coupling terms 'analytic', 'baeck-an', 'none'
29+
couplings='none', ! non-adiabatic coupling terms 'analytic', 'baeck-an', 'none'
2930
nohop=0,
3031
decoh_alpha=0.1
3132
popthr=0.01
32-
energydifthr=0.50
33-
energydriftthr=0.50
33+
energydifthr=1.00
34+
energydriftthr=1.00
3435
correct_decoherence=.false.
3536
/

0 commit comments

Comments
 (0)