Solve

The solve function is the primary interface for solving non-linear ODE problems using various QSS (Quantized State Systems) algorithms. It does this by dispatching on the problem and algorithm types to select the right solver.

Based on the QSS Algorithm provided, the solve function either selects a basic QSS integration method or, in the case of LiQSS, constructs additional data structures needed for implicit integration. The method defaults to sparse handling as false (Val(false)), tolerances (abstol=1e-4 and reltol=1e-3), and a maximum number of iterations (maxiters=1e7). These parameters can be adjusted based on the problem's complexity and desired accuracy. When a user doesn't provide a solver explicitly, the solve function defaults to using the modified second-order implicit algorith (mLiQSS2).

Helper Functions:

createCommonData: Sets up the common data required by the QSS solver. It initializes all necessary vectors (x, q, tx, tq, nextStateTime, etc.) used in the integration process, and it pre-allocates a cache of Taylor series to avoid repeated memory allocation during integration (taylorOpsCache).

getClosure: Creates a closure for the Jacobian and the dependency matrices, allowing in a flexible in a way that enables easy extension. The closure is used to access the Jacobian matrix during the integration process.

createLiqssData: For LiQSS solvers, this function sets up additional data structures and auxiliary variables used for storing linear approximation coefficients.

Internals

QuantizedSystemSolver.custom_SolveMethod
custom_Solve(prob::NLODEProblem{F,PRTYPE,T,D,Z,CS},al::QSSAlgorithm{Solver, Order},::Val{M},finalTime::Float64,saveat::Float64,initialTime::Float64,abstol::Float64,reltol::Float64,maxErr::Float64,maxiters::Int) where{F,PRTYPE,T,D,Z,CS,Solver,Order,M}

calls the integrator to solve the nonlinear ODE problem.

Arguments

  • prob::NLODEProblem{F,PRTYPE,T,D,Z,CS}: The nonlinear ODE problem to solve.
  • al::QSSAlgorithm{Solver, Order}: The QSS algorithm to use for solving the problem.
  • ::Val{M}: A type parameter indicating which detection mechanism to use.
  • finalTime::Float64: The final time for the simulation.
  • saveat::Float64: The time interval at which to save the solution.
  • initialTime::Float64: The initial time for the simulation.
  • abstol::Float64: The absolute tolerance for the solver.
  • reltol::Float64: The relative tolerance for the solver.
  • maxErr::Float64: The maximum allowable error.
  • maxiters::Int: The maximum number of iterations.
source
QuantizedSystemSolver.createCommonDataMethod
createCommonData(prob::NLODEProblem{F,PRTYPE,T,D,Z,CS},::Val{Order},finalTime::Float64,saveat::Float64,initialTime::Float64,abstol::Float64,reltol::Float64,maxErr::Float64,maxiters::Int) where{F,PRTYPE,T,D,Z,CS,Order}

creates the necessary data for the simulation and stores it in a CommonQSS_Data struct.

Arguments

  • prob::NLODEProblem{F,PRTYPE,T,D,Z,CS}: The nonlinear ODE problem to solve.
  • ::Val{Order}: The order of the algorithm.
  • finalTime::Float64: The final time for the simulation.
  • saveat::Float64: The time interval at which to save the solution.
  • initialTime::Float64: The initial time for the simulation.
  • abstol::Float64: The absolute tolerance for the solver.
  • reltol::Float64: The relative tolerance for the solver.
  • maxErr::Float64: The maximum allowable error.
  • maxiters::Int: The maximum number of iterations.

Returns

  • A data structure containing common data required for the QSS algorithm.
source
QuantizedSystemSolver.createLiqssDataMethod
createLiqssData(::Val{M},::Val{T},::Val{Order}) where {T,Order}

Creates LIQSS-specific data required for solving an ODE problem.

Arguments

  • ::Val{M}: A type parameter indicating which detection mechanism to use.
  • ::Val{T}: The number of continuous variables.
  • ::Val{Order}: The order of the algorithm.

Returns

  • A data structure containing LIQSS-specific data required for liQSS algorithms.
source
QuantizedSystemSolver.LiQSS_DataType
LiQSS_Data{O,M}

helper datastructures needed only for implicit case The field variables are:

  • cd::Val{M} #
  • cacheA::MVector{1,Float64}
  • qaux::Vector{MVector{O,Float64}}
  • dxaux::Vector{MVector{O,Float64}}
source