DiabaticModels

NQCModels.DiabaticModels.BosonBathType
BosonBath(density::SpectralDensity, N::Integer)

Bosonic bath with given spectral density.

Useful for sampling the bath uncoupled from the spin for spin-boson dynamics.

source
NQCModels.DiabaticModels.DiabaticFrictionModelType
DiabaticFrictionModel <: 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.

source
NQCModels.DiabaticModels.DiabaticModelType
DiabaticModel <: Model

DiabaticModels 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

DiabaticModels 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
source
NQCModels.DiabaticModels.ErpenbeckThossType
struct 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)
source
NQCModels.DiabaticModels.FullGaussLegendreType
FullGaussLegendre{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.

source
NQCModels.DiabaticModels.LargeDiabaticModelType
LargeDiabaticModel <: 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.

source
NQCModels.DiabaticModels.MiaoSubotnikType
MiaoSubotnik{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)
source
NQCModels.DiabaticModels.ReferenceGaussLegendreType
ReferenceGaussLegendre{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.
source
NQCModels.DiabaticModels.ShenviGaussLegendreType
ShenviGaussLegendre{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.

source