Stationary State

HierarchicalEOM.jl implements two different ways to calculate stationary states of all Auxiliary Density Operators (ADOs).

To solve the stationary state of the reduced state and also all the ADOs, you only need to call SteadyState. Different methods are implemented with different input parameters of the function which makes it easy to switch between different methods. The output of the function SteadyState for each methods will always be in the type of the auxiliary density operators ADOs.

Solve with LinearSolve.jl

The first method is implemented by solving the linear problem

\[0=\hat{\mathcal{M}}\rho^{(m,n,p)}_{\textbf{j} \vert \textbf{q}}(t)\]

HierarchicalEOM.jl wraps some of the functions in LinearSolve.jl, which is a very rich numerical library for solving the linear problems and provides many solvers. It offers quite a few options for the user to tailor the solver to their specific needs. The default solver (and its corresponding settings) are chosen to suit commonly encountered problems and should work fine for most of the cases. If you require more specialized methods, such as the choice of algorithm, please refer to benchmark for LinearSolve solvers and also the documentation of LinearSolve.jl.

See the docstring of this method:

HierarchicalEOM.HeomAPI.SteadyStateMethod
SteadyState(M; solver, verbose, SOLVEROptions...)

Solve the steady state of the auxiliary density operators based on LinearSolve.jl (i.e., solving $x$ where $A \times x = b$).

Parameters

  • M::AbstractHEOMLSMatrix : the matrix given from HEOM model, where the parity should be EVEN.
  • solver : solver in package LinearSolve.jl. Default to UMFPACKFactorization().
  • verbose::Bool : To display verbose output and progress bar during the process or not. Defaults to true.
  • SOLVEROptions : extra options for solver

For more details about solvers and extra options, please refer to LinearSolve.jl

Returns

  • ::ADOs : The steady state of auxiliary density operators.
source
# the HEOMLS matrix
M::AbstractHEOMLSMatrix  
ados_steady = SteadyState(M)
Unphysical solution

This method does not require an initial condition $\rho^{(m,n,p)}_{\textbf{j} \vert \textbf{q}}(0)$. Although this method works for most of the cases, it does not guarantee that one can obtain a physical (or unique) solution. If there is any problem within the solution, please try the second method which solves with an initial condition, as shown below.

Solve with DifferentialEquations.jl

The second method is implemented by solving the ordinary differential equation (ODE) method : SteadyStateDiffEq.jl

\[\partial_{t}\rho^{(m,n,p)}_{\textbf{j} \vert \textbf{q}}(t)=\hat{\mathcal{M}}\rho^{(m,n,p)}_{\textbf{j} \vert \textbf{q}}(t)\]

until finding a stationary solution.

HierarchicalEOM.jl wraps some of the functions in DifferentialEquations.jl, which is a very rich numerical library for solving the differential equations and provides many ODE solvers. It offers quite a few options for the user to tailor the solver to their specific needs. The default solver (and its corresponding settings) are chosen to suit commonly encountered problems and should work fine for most of the cases. If you require more specialized methods, such as the choice of algorithm, please refer to the documentation of DifferentialEquations.jl.

Given the initial state as Density Operator (AbstractMatrix type)

See the docstring of this method:

HierarchicalEOM.HeomAPI.SteadyStateFunction
SteadyState(M, ρ0, tspan; solver, verbose, SOLVEROptions...)

Solve the steady state of the auxiliary density operators based on time evolution (ordinary differential equations; SteadyStateDiffEq.jl) with initial state is given in the type of density-matrix (ρ0).

Parameters

  • M::AbstractHEOMLSMatrix : the matrix given from HEOM model, where the parity should be EVEN.
  • ρ0 : system initial state (density matrix)
  • tspan::Number : the time limit to find stationary state. Default to Inf
  • solver : The ODE solvers in package DifferentialEquations.jl. Default to DP5().
  • reltol::Real : Relative tolerance in adaptive timestepping. Default to 1.0e-8.
  • abstol::Real : Absolute tolerance in adaptive timestepping. Default to 1.0e-10.
  • save_everystep::Bool : Saves the result at every step. Defaults to false.
  • verbose::Bool : To display verbose output and progress bar during the process or not. Defaults to true.
  • SOLVEROptions : extra options for solver

For more details about solvers, and extra options, please refer to DifferentialEquations.jl

Returns

  • ::ADOs : The steady state of auxiliary density operators.
source
# the HEOMLS matrix
M::AbstractHEOMLSMatrix  

# the initial state of the system density operator
ρ0::AbstractMatrix

ados_steady = SteadyState(M, ρ0)

Given the initial state as Auxiliary Density Operators

See the docstring of this method:

HierarchicalEOM.HeomAPI.SteadyStateFunction
SteadyState(M, ados, tspan; solver, verbose, SOLVEROptions...)

Solve the steady state of the auxiliary density operators based on time evolution (ordinary differential equations; SteadyStateDiffEq.jl) with initial state is given in the type of ADOs.

Parameters

  • M::AbstractHEOMLSMatrix : the matrix given from HEOM model, where the parity should be EVEN.
  • ados::ADOs : initial auxiliary density operators
  • tspan::Number : the time limit to find stationary state. Default to Inf
  • solver : The ODE solvers in package DifferentialEquations.jl. Default to DP5().
  • reltol::Real : Relative tolerance in adaptive timestepping. Default to 1.0e-8.
  • abstol::Real : Absolute tolerance in adaptive timestepping. Default to 1.0e-10.
  • save_everystep::Bool : Saves the result at every step. Defaults to false.
  • verbose::Bool : To display verbose output and progress bar during the process or not. Defaults to true.
  • SOLVEROptions : extra options for solver

For more details about solvers, and extra options, please refer to DifferentialEquations.jl

Returns

  • ::ADOs : The steady state of auxiliary density operators.
source
# the HEOMLS matrix
M::AbstractHEOMLSMatrix  

# the initial state of the ADOs
ados::ADOs

ados_steady = SteadyState(M, ados)