DiabaticModels
NQCModels.DiabaticModels
— ModuleDiabaticModels
Models defined within this module subtype the DiabaticModel
and provide potentials as Hermitian matrices and derivatives as arrays of Hermitian matrices.
NQCModels.DiabaticModels.AltDebyeSpectralDensity
— TypeAltDebyeSpectralDensity{T} <: SpectralDensity
Standard Debye spectral density but uses an alternative discretization scheme that requires a cutoff parameter ωᵐ
.
References
Najeh Rekik, Chang-Yu Hsieh, Holly Freedman, Gabriel Hanna, J. Chem. Phys. 138, 144106 (2013)
NQCModels.DiabaticModels.AnanthModelOne
— TypeAnanthModelOne(a=0.01, b=1.6, c=0.005, d=1.0)
Ananth's simple avoided crossing model (similar to Tully's first model) from J. Chem. Phys. 127, 084114 (2007).
NQCModels.DiabaticModels.AnanthModelTwo
— TypeAnanthModelTwo(a=0.04, b=0.01, c=0.005, d=1.0, e=0.7, f=1.6)
Ananth's asymmetric model from J. Chem. Phys. 127, 084114 (2007).
NQCModels.DiabaticModels.BosonBath
— TypeBosonBath(density::SpectralDensity, N::Integer)
Bosonic bath with given spectral density.
Useful for sampling the bath uncoupled from the spin for spin-boson dynamics.
NQCModels.DiabaticModels.DebyeSpectralDensity
— TypeDebyeSpectralDensity{T} <: SpectralDensity
Debye density as detailed in: Xin He, Jian Liu, J. Chem. Phys. 151, 024105 (2019)
NQCModels.DiabaticModels.DiabaticFrictionModel
— TypeDiabaticFrictionModel <: LargeDiabaticModel
These models are defined identically to the LargeDiabaticModel
but allocate extra temporary arrays when used with NQCDynamics.jl
.
This allows for the calculation of electronic friction internally from the diabatic potential after diagonalisation and calculation of nonadiabatic couplings.
NQCModels.DiabaticModels.DiabaticModel
— TypeDiabaticModel <: Model
DiabaticModel
s are used when a system has multiple electronic states that are presented in the diabatic representation. This is the case for the majority of model systems.
Implementation
DiabaticModel
s should implement:
potential(model, R)
derivative!(model, D, R)
nstates(model)
ndofs(model)
Example
In this example we create a simple 2 state, 1 dimensional diabatic model MyModel
. As noted above, we implement the 4 relevant functions then evaluate the potential. Since this is a 1D model the argument R
accepts a Real
value.
using StaticArrays: SMatrix
using LinearAlgebra: Hermitian
struct MyModel <: NQCModels.DiabaticModels.DiabaticModel end
NQCModels.nstates(::MyModel) = 2
NQCModels.ndofs(::MyModel) = 1
function NQCModels.potential(::MyModel, R::Real)
V11 = R
V22 = -R
V12 = 1
return Hermitian(SMatrix{2,2}(V11, V12, V12, V22))
end
function NQCModels.derivative!(::MyModel, D, R::Real)
return Hermitian(SMatrix{2,2}(1, 0, 0, 1))
end
model = MyModel()
NQCModels.potential(model, 10)
# output
2×2 Hermitian{Int64, SMatrix{2, 2, Int64, 4}}:
10 1
1 -10
NQCModels.DiabaticModels.DoubleWell
— TypeDoubleWell(mass=1, ω=1, γ=1, Δ=1)
Two state double well, also called the one-dimensional spin-boson model. See: J. Chem. Phys. 150, 244102 (2019)
NQCModels.DiabaticModels.ErpenbeckThoss
— Typestruct ErpenbeckThoss{T<:AbstractFloat} <: DiabaticModel
1D two-state diabatic system capable of modelling a molecule adsorbed on a metal surface or a single-molecule junction.
In the two references, all of the parameters are identical except for the particle mass m
and the vertical shift c
applied to the ϵ₀ state. Both references modify the shift to ensure the quantum ground-state has an energy of 0 eV. Note that the mass m
is specified in atomic mass units (amu) not atomic units. We calculate the offset automatically in the constructor from the Morse potential zero-point energy.
References
- PHYSICAL REVIEW B 97, 235452 (2018)
- J. Chem. Phys. 151, 191101 (2019)
NQCModels.DiabaticModels.FullGaussLegendre
— TypeFullGaussLegendre{T} <: WideBandBathDiscretisation
Use Gauss-Legendre quadrature to discretise the continuum across the entire band width. This is similar to the ShenviGaussLegendre except that splits the continuum at the Fermi level into two halves.
NQCModels.DiabaticModels.GatesHollowayElbow
— TypeGatesHollowayElbow()
Simple two state elbow potential from Gates and Holloway: Journal of Electron Spectroscopy and Related Phenomena, 64/65 (1993) 633-639
Has two diabatic states each comprised of the sum of a Morse and a repulsive potential. The coupling between them is an exponential function of z
(distance from the surface).
NQCModels.DiabaticModels.LargeDiabaticModel
— TypeLargeDiabaticModel <: DiabaticModel
Same as the DiabaticModels
but uses normal Julia arrays instead of StaticArrays and must implement the inplace potential!
rather than potential
. This is useful when nstates
is very large and StaticArrays are no longer efficient.
NQCModels.DiabaticModels.MiaoSubotnik
— TypeMiaoSubotnik{T<:AbstractFloat} <: DiabaticModel
Double well model with parameters matching those of Miao and Subotnik in the reference. This model should be paired with the AndersonHolstein
model to couple to the bath of metallic states.
References
- J. Chem. Phys. 150, 041711 (2019)
NQCModels.DiabaticModels.OhmicSpectralDensity
— TypeOhmicSpectralDensity{T} <: SpectralDensity
Ohmic density as detailed in: Xin He, Jian Liu, J. Chem. Phys. 151, 024105 (2019)
NQCModels.DiabaticModels.OuyangModelOne
— TypeOuyangModelOne(A=0.01, B=1.6, Γ=1e-4, N=10, ΔE=1.6e-2, D=1.0)
Model #1 from Ouyang and Subotnik. See also Ouyang's thesis.
NQCModels.DiabaticModels.ReferenceGaussLegendre
— TypeReferenceGaussLegendre{T}
Implementation translated from Fortran code used for simulations of Shenvi et al. in J. Chem. Phys. 130, 174107 (2009). Two differences from ShenviGaussLegendre:
- Position of minus sign in energy levels has been corrected.
- Division by sqrt(ΔE) in the coupling.
NQCModels.DiabaticModels.ShenviGaussLegendre
— TypeShenviGaussLegendre{T}
Defined as described by Shenvi et al. in J. Chem. Phys. 130, 174107 (2009). The position of the negative sign for the state energy level has been moved to ensure the states are sorted from lowest to highest.
NQCModels.DiabaticModels.SpinBoson
— TypeSpinBoson(density::SpectralDensity, N::Integer, ϵ, Δ)
Spin boson model with N
bosons with given spectral density.
References
NQCModels.DiabaticModels.ThreeStateMorse
— TypeThreeStateMorse()
Three state morse potential referred to as Model IA here: J. Chem. Phys. 150, 244102 (2019)
Models IB and IC retain the same functional form and need only a change of parameters.
NQCModels.DiabaticModels.TrapezoidalRule
— TypeTrapezoidalRule{B,T} <: WideBandBathDiscretisation
Discretise wide band continuum using trapezoidal rule. Leads to evenly spaced states and constant coupling.
NQCModels.DiabaticModels.TullyModelOne
— TypeTullyModelOne(a=0.01, b=1.6, c=0.005, d=1.0)
Tully's simple avoided crossing model from J. Chem. Phys. 93, 1061 (1990).
NQCModels.DiabaticModels.TullyModelThree
— TypeTullyModelThree(a=0.0006, b=0.1, c=0.9)
Tully's extended coupling with reflection model from J. Chem. Phys. 93, 1061 (1990).
NQCModels.DiabaticModels.TullyModelTwo
— TypeTullyModelTwo(a=0.1, b=0.28, c=0.015, d=0.06, e=0.05)
Tully's dual avoided crossing model from J. Chem. Phys. 93, 1061 (1990).
NQCModels.DiabaticModels.discretize
— MethodDiscretize a given spectral density for N oscillators. Returns frequencies and couplings.