Propagating nonadiabatic systems forward in time remains challenging particularly at crossings
between the ground state () and the first excited state ().
For conical intersections, conventional electronic structure methods
(such as HF/CIS and DFT/TDDFT) recover the wrong dimensionality of the crossing seam.
The TDDFT-1D method address this limitation by combining DFT and configuration interaction
to achieve smooth crossings of and states
1293
J. Phys. Chem. Lett.
(2019),
10,
pp. 3426.
Link
while keeping the cost of the
calculation remains low (comaprable to a standard TDDFT calculation).
The idea behind the CIS-1D and the TDDFT-1D method is to diagonalize a configuration interaction Hamiltonian whose basis includes the reference state, all singly excited configurations, and one unique doubly excited configuration. This unique double (the 1D in CIS-1D/TDDFT-1D) in the configuration interaction Hamiltonian gives rise to the CIS-1D method if one starts from a Hartree-Fock reference, and the TDDFT-1D method if one starts from a Kohn-Sham reference.
The doubly excited configuration is formed by exciting a pair of electrons from the HOMO () to the LUMO (). To find the unique double, the cannonical molecular orbitals are rotated such that the energy of the doubly excited configuration () is minimized
55
J. Chem. Phys.
(2021),
155,
pp. 154105.
Link
.
From these optimized orbitals, the following configuration interaction Hamiltonian is constructed, as explained previosly, in the basis of the reference state, all singles , and one unique double :
(7.47) |
Diagonalizing the above Hamiltonian yields the CIS-1D or TDDFT-1D states. The excitation energies computed with TDDFT-1D are generally at the same level of accuracy as TDDFT with the Tamm-Dancoff approximation (TDA) at geometries far from crossing points. However, TDDFT-1D’s key advantage appears near ground and excited state crossing points, where it produces smooth potential energy surfaces.
To facilitate nonadiabatic dynamics simulations, analytical gradients and derivative couplings
54
J. Chem. Phys.
(2022),
157,
pp. 244110.
Link
are available for the TDDFT-1D method. Further, minimum energy crossing points (MECPs) can be located for TDDFT-1D states.
The rem variables for TDDFT-1D calculations are similar to those for CIS calculations. The following rem variables are used to control the TDDFT-1D calculation:
CIS1D_N_ROOTS
CIS1D_N_ROOTS
Sets the number of CIS-1D and TDDFT-1D states to calculate. The lowest CIS-1D/TDDFT-1D state is the ground state, which is followed by excited states.
TYPE:
INTEGER
DEFAULT:
0
Do not calculate any CIS-1D or TDDFT-1D states.
OPTIONS:
, Calculate the lowest CIS-1D or TDDFT-1D states.
RECOMMENDATION:
None
CIS1D_ED_CONVERGENCE
CIS1D_ED_CONVERGENCE
The first step in TDDFT-1D is to find the optimized orbitals for the double by minimizing the energy of the doubly excited configuration, . This variable controls the convergence criterion for the minimization of . The orbitals are optimized when the error is lower than .
TYPE:
INTEGER
DEFAULT:
7
OPTIONS:
Convergence achieved when the error is lower than .
RECOMMENDATION:
The convergence criterion for the roots of the CIS-1D and TDDFT-1D calculations is set by the CIS_CONVERGENCE rem variable, which is set to 9 by default for a CIS-1D/TDDFT-1D calculation. If CIS_CONVERGENCE is increased, then CIS1D_ED_CONVERGENCE should also be increased.
CIS1D_SCALE_GD
CIS1D_SCALE_GD
The coupling element between the ground and doubly excited configuration in the Hamiltonian in Eq. (7.47), () are scaled by this factor.
TYPE:
INTEGER
DEFAULT:
100
OPTIONS:
n
, The coupling element is scaled by .
RECOMMENDATION:
Since the KS reference is not a true wave function, there is no rigorous way to determine the coupling elements in the configuration interaction Hamiltonian. It sometimes becomes necessary to scale the coupling elements to obtain accurate potential energy surfaces, and the scaling factor has to be determined with some benchmarking
55
J. Chem. Phys.
(2021),
155,
pp. 154105.
Link
.
CIS1D_SCALE_SD
CIS1D_SCALE_SD
The coupling element between the singles and the double in the Hamiltonian in Eq. (7.47), () are scaled by this factor.
TYPE:
INTEGER
DEFAULT:
100
OPTIONS:
n
, The coupling element is scaled by .
RECOMMENDATION:
Same as CIS1D_SCALE_GD.
CIS1D_STATE_DERIV
CIS1D_STATE_DERIV
Selects the CIS-1D/TDDFT-1D state for which gradients are calculated. This is useful for jobs such as geometry optimizations.
TYPE:
INTEGER
DEFAULT:
-1
Does not select any state
OPTIONS:
, Selects the th CIS-1D/TDDFT-1D state.
RECOMMENDATION:
None
CIS1D_DER_NUMSTATE
CIS1D_DER_NUMSTATE
Determines the number of states for which derivative couplings are to be calculated. The states are specified in the $derivative_coupling section
TYPE:
INTEGER
DEFAULT:
0
OPTIONS:
0
Do not calculate derivative couplings.
Calculate pairs of derivative couplings.
RECOMMENDATION:
None
Example 7.23 TDDFT-1D calculation of water using the B3LYP functional for the lowest 5 states.
$molecule 0 1 H 0.96 0 0 O 0.000000000000000 0.000000000000000 0.000000000000000 H -0.240364803892264 0.929421734762983 0 $end $rem METHOD b3lyp BASIS cc-pvdz SYM_IGNORE true cis1d_n_roots 5 $end
Example 7.24 TDDFT-1D calculation of the gradient of the 1st excited state of water in the presence of external charges
$molecule 0 1 H 0.96 0 0 O 0.000000000000000 0.000000000000000 0.000000000000000 H -0.240364803892264 0.929421734762983 0 $end $rem METHOD b3lyp BASIS cc-pvdz SYM_IGNORE true cis1d_n_roots 4 cis1d_state_deriv 1 jobtype force $end $external_charges 0.0 0.0 -1.0 -1.0 0.0 0.0 1.2 1.0 0.0 0.0 -2.0 -1.0 0.0 0.0 2.2 1.0 $end
Example 7.25 Minimum energy crossing point (MECP) between the ground and the 1st excited state of ethylene using the B97X functional.
$rem geom_opt_driver optimize scf_convergence 10 thresh_diis_switch 9 cis_convergence 8 jobtype opt mecp_opt true mecp_methods branching_plane MECP_PROJ_HESS true mecp_state1 [0,0] mecp_state2 [0,1] unrestricted false calc_nac 1 scf_algorithm diis_gdm cis1d_n_roots 5 exchange wb97x basis 6-31G* symmetry_ignore true $end $molecule 0 1 C 0.0188526101 0.0181382820 -0.3533168752 C 1.3813178426 0.0038321438 -0.0594695222 H 2.0062266305 -0.8967843506 0.0891927721 H -0.3788613100 -0.7769905091 -1.0012547284 H -0.0243170648 -0.6123230262 0.6288689086 H 1.9213240937 0.9262667800 0.1943673040 $end