DynamicsMethods
NQCDynamics.DynamicsMethods — Module
This module contains functions and types necessary for performing nonadiabatic molecular dynamics.
Dynamics is performed using DifferentialEquations.jl. As such, this module is centered around the implementation of the functions necessary to integrate the dynamics.
For deterministic Hamiltonian methods, the central function is DynamicsMethods.motion!, which is the inplace form of the function to be integrated by DifferentialEquations.jl.
Further, methods that have discontinuities, such as surface hopping, use the callback interface provided by DifferentialEquations.jl.
NQCDynamics.DynamicsMethods.Method — Type
Each type of dynamics subtypes Method which is passed to the AbstractSimulation as a parameter to determine the type of dynamics desired.
NQCDynamics.DynamicsMethods.DynamicsVariables — Method
DynamicsVariables(::AbstractSimulation, args...)For each dynamics method this function is implemented to provide the variables for the dynamics in the appropriate format.
By default, DynamicsVariables is set up for the classical case and takes sim, v, r as arguments and returns a ComponentVector(v=v, r=r) which is used as a container for the velocities and positions during classical dynamics.
NQCDynamics.DynamicsMethods.create_problem — Method
Provides the DEProblem for each type of simulation.
NQCDynamics.DynamicsMethods.get_callbacks — Method
Select the default callbacks for this simulation type.
NQCDynamics.DynamicsMethods.motion! — Function
motion!(du, u, sim, t)As per DifferentialEquations.jl, this function is implemented for each method and defines the time derivatives of the DynamicalVariables.
We require that each implementation ensures du and u are subtypes of DynamicalVariables and sim subtypes AbstractSimulation.
NQCDynamics.DynamicsMethods.select_algorithm — Method
Choose a default algorithm for solving the differential equation.
ClassicalMethods
NQCDynamics.DynamicsMethods.ClassicalMethods.Classical — Type
Classical <: DynamicsMethods.MethodType for performing classical molecular dynamics.
sim = Simulation{Classical}(Atoms(:H), Harmonic())
# output
Simulation{Classical}:
Atoms{Float64}([:H], [1], [1837.4715941070515])
Harmonic{Float64, Float64, Float64}
m: Float64 1.0
ω: Float64 1.0
r₀: Float64 0.0
dofs: Int64 1NQCDynamics.DynamicsMethods.ClassicalMethods.Langevin — Type
Type for performing Langevin molecular dynamics.
using Unitful
sim = Simulation{Langevin}(Atoms(:H), Free(); γ=2.5, temperature=100u"K")
# output
Simulation{Langevin{Float64}}:
Atoms{Float64}([:H], [1], [1837.4715941070515])
Free(1)NQCDynamics.DynamicsMethods.ClassicalMethods.MDEF — Type
\[dr = v dt\\ dv = -\Delta U/M dt - \Gamma v dt + \sigma \sqrt{2\Gamma} dW\]
$\Gamma$ is the friction tensor with units of inverse time. For thermal dynamics we set $\sigma = \sqrt{kT / M}$, where $T$ is the electronic temperature.
This is integrated using the BAOAB algorithm where the friction "O" step is performed in the tensor's eigenbasis. See src/dynamics/mdef_baoab.jl for details.
NQCDynamics.DynamicsMethods.ClassicalMethods.ThermalLangevin — Type
Type for performing Langevin ring polymer molecular dynamics.
Currently there are separate types for classical and ring polymer versions of Langevin dynamics but they should be combined. The reason they are not at the moment is that they use different integration algorithms and require slightly different fields.
using Unitful
RingPolymerSimulation{ThermalLangevin}(Atoms(:H), Free(), 10; γ=0.1, temperature=25u"K")
# output
RingPolymerSimulation{ThermalLangevin{Float64}}:
Atoms{Float64}([:H], [1], [1837.4715941070515])
Free(1)
with 10 beads.NQCDynamics.DynamicsMethods.ClassicalMethods.acceleration! — Method
f1 in DifferentialEquations.jl docs.
NQCDynamics.DynamicsMethods.ClassicalMethods.friction! — Method
friction!(g, r, sim, t)Evaluates friction tensor
MappingVariableMethods
NQCDynamics.DynamicsMethods.MappingVariableMethods.NRPMD — Type
NRPMD{T} <: DynamicsMethods.MethodNonadiabatic ring polymer molecular dynamics Uses Meyer-Miller-Stock-Thoss mapping variables for electronic degrees of freedom and ring polymer formalism for nuclear degrees of freedom.
RingPolymerSimulation{NRPMD}(Atoms(:H), DoubleWell(), 10)
# output
RingPolymerSimulation{NRPMD{Float64}}:
Atoms{Float64}([:H], [1], [1837.4715941070515])
DoubleWell{Int64, Int64, Int64, Int64}
mass: Int64 1
ω: Int64 1
γ: Int64 1
Δ: Int64 1
with 10 beads.SurfaceHoppingMethods
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods — Module
SurfaceHoppingMethodsImplementation for surface hopping methods.
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.AdiabaticIESH — Type
IESH{T} <: SurfaceHoppingIndependent electron surface hopping.
References
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.BCME — Type
BCME{T} <: ClassicalMasterEquationExtension to CME that incorporates broadening in the potential energy surfaces.
Note that we do not rescale the velocity as this is not mentioned only in the 2016 paper, not any of the later ones, so I presume they later decided not to do it and instead keep it the same as the original CME.
- Dou, Subotnik, J. Chem. Phys. 144, 024116 (2016)
- Dou, Subotnik, J. Phys. Chem. A, 24, 757-771 (2020)
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.CME — Type
CME{T} <: ClassicalMasterEquationSimple surface hopping method for Newns-Anderson (Anderson-Holstein) models.
- Dou, Nitzan, Subotnik, J. Chem. Phys. 142, 084110 (2015)
- Dou, Subotnik, J. Phys. Chem. A, 24, 757-771 (2020)
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.DecoherenceCorrectionEDC — Type
DecoherenceCorrectionEDC{T}Energy decoherence correction of Granucci and Persico in J. Chem. Phys. 126, 134114 (2007).
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.FSSH — Type
FSSH{T} <: SurfaceHoppingType for fewest-switches surface hopping
Simulation{FSSH}(Atoms(:H), Free())
# output
Simulation{FSSH{Float64}}:
Atoms{Float64}([:H], [1], [1837.4715941070515])
Free(1)NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.SurfaceHopping — Type
Abstract type for all surface hopping methods.
Surface hopping methods follow the structure set out in this file. The nuclear and electronic variables are propagated by the motion! function. The surface hopping procedure is handled by the HoppingCallback which uses the functions check_hop! and execute_hop! as its condition and affect!.
To add a new surface hopping scheme, you must create a new struct and define methods for evaluate_hopping_probability!, select_new_state, and rescale_velocity!.
See fssh.jl for an example implementation.
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.SurfaceHoppingVariables — Type
Due to the need to multiple dispatch this type mirrors that of NamedArrayPartition from the great
pacakge RecursiveArrayTools but has been renamed to SurfaceHoppingVariables. This was taken from
version 3.33.0 so any future changes to NamedArrayPartition will not be mirrored here.NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.apply_decoherence_correction! — Method
Eq. 17 in J. Chem. Phys. 126, 134114 (2007)
Modify the wavefunction coefficients in ψ after a successful surface hop.
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.calculate_Akj — Method
Equation 17 in Shenvi, Roy, Tully 2009. Uses equations 19 and 20.
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.calculate_a — Method
calculate_a(sim::AbstractSimulation{<:SurfaceHopping}, coupling::AbstractMatrix)Equation 40 from [36].
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.calculate_b — Method
calculate_b(coupling::AbstractMatrix, velocity::AbstractMatrix)Equation 41 from [36].
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.compute_overlap! — Method
Equation 20 in Shenvi, Roy, Tully 2009.
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.evaluate_broadening — Method
evaluate_broadening(h, μ, β, Γ)Evaluate the convolution of the Fermi function with a Lorentzian.
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.evaluate_hopping_probability! — Method
evaluate_hopping_probability!(sim::Simulation{<:FSSH}, u, dt)Evaluates the probability of hopping from the current state to all other states
Implementation
σis Hermitan so the choiceσ[m,s]orσ[s,m]is irrelevant; we take the real part.- 'd' is skew-symmetric so here the indices are important.
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.evaluate_hopping_probability! — Method
Hopping probability according to equation 21 in Shenvi, Roy, Tully 2009.
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.evaluate_hopping_probability! — Method
This function should set the field sim.method.hopping_probability.
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.extract_nonadiabatic_coupling — Method
extract_nonadiabatic_coupling(coupling, new_state, old_state)Extract the nonadiabatic coupling vector between states new_state and old_state
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.frustrated_hop_invert_velocity! — Method
frustrated_hop_invert_velocity!(
sim::AbstractSimulation{<:SurfaceHopping}, velocity, d
)Measures the component of velocity along the nonadiabatic coupling vector and inverts that component.
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.iesh_diabatic_population — Method
Calculate the diabatic population in J. Chem. Theory Comput. 2022, 18, 4615−4626. Eqs. 12, 13 and 14 describe the steps to calculat it though the code is written to use matrix operations instead of summations to calculate the result.
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.perform_rescaling! — Method
perform_rescaling!(
sim::AbstractSimulation{<:SurfaceHopping}, velocity, velocity_rescale, d
)Equation 33 from [36].
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.rescale_velocity! — Method
rescale_velocity!(sim::AbstractSimulation{<:SurfaceHopping}, u)::BoolRescale the velocity in the direction of the nonadiabatic coupling.
References
[36]
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.select_new_state — Method
This function should return the desired state determined by the probability. Should return the original state if no hop is desired.
NQCDynamics.DynamicsMethods.SurfaceHoppingMethods.unpack_states — Method
unpack_states(sim)Get the two states that we are hopping between.
NQCDynamics.DynamicsUtils.acceleration! — Method
Set the acceleration due to the force from the currently occupied states. See Eq. 12 of Shenvi, Tully JCP 2009 paper.
NQCDynamics.DynamicsUtils.set_quantum_derivative! — Method
Propagation of electronic wave function happens according to Eq. (14) in the Shenvi, Tully paper (JCP 2009)
In IESH each electron is independent so we can loop through electrons and set the derivative one at a time, in the standard way for FSSH.
EhrenfestMethods
NQCDynamics.DynamicsMethods.EhrenfestMethods.AbstractEhrenfest — Type
Abstract type for Ehrenfest method.
NQCDynamics.DynamicsMethods.EhrenfestMethods.Ehrenfest — Type
Ehrenfest{T} <: AbstractEhrenfestEhrenfest molecular dynamics. Classical molecular dynamics where the force is derived by averaging contributions from multiple electronic states.
Simulation{Ehrenfest}(Atoms(:H), DoubleWell())
# output
Simulation{Ehrenfest{Float64}}:
Atoms{Float64}([:H], [1], [1837.4715941070515])
DoubleWell{Int64, Int64, Int64, Int64}
mass: Int64 1
ω: Int64 1
γ: Int64 1
Δ: Int64 1IntegrationAlgorithms
NQCDynamics.DynamicsMethods.IntegrationAlgorithms.MInt — Type
MInt <: OrdinaryDiffEqAlgorithmSecond order symplectic momentum integral algorithm.
Reference
NQCDynamics.DynamicsMethods.IntegrationAlgorithms.RingPolymerMInt — Type
RingPolymerMInt <: OrdinaryDiffEqAlgorithmSecond order symplectic momentum integral algorithm applied to NRPMD.
Reference
NQCDynamics.DynamicsMethods.IntegrationAlgorithms.get_C_propagator — Method
Get the C propagator for the mapping variables.
NQCDynamics.DynamicsMethods.IntegrationAlgorithms.get_D_propagator — Method
Get the D propagator for the mapping variables.
NQCDynamics.DynamicsMethods.IntegrationAlgorithms.get_gamma — Method
Get the Γ variable used to calculate the nuclear propagators.
NQCDynamics.DynamicsMethods.IntegrationAlgorithms.get_mapping_nuclear_force — Method
Get the force due to the mapping variables.
Equivalent to this but doesn't allocate: return 0.5 * (q'Eq + p'Ep) - q'Fp
NQCDynamics.DynamicsMethods.IntegrationAlgorithms.get_mapping_nuclear_force — Method
Get the force due to the mapping variables.
NQCDynamics.DynamicsMethods.IntegrationAlgorithms.get_xi — Method
Get the Ξ variable used to calculate the nuclear propagators.
NQCDynamics.DynamicsMethods.IntegrationAlgorithms.set_C_propagator! — Method
Get the C propagator for the mapping variables.
NQCDynamics.DynamicsMethods.IntegrationAlgorithms.set_D_propagator! — Method
Get the D propagator for the mapping variables.
NQCDynamics.DynamicsMethods.IntegrationAlgorithms.set_gamma! — Method
Get the Γ variable used to calculate the nuclear propagators.
NQCDynamics.DynamicsMethods.IntegrationAlgorithms.set_xi! — Method
Get the Ξ variable used to calculate the nuclear propagators.
StochasticDiffEq.alg_cache — Method
Insecting the inputs into a Cache structure.