QuantizedSystemSolver.jl

This is an ODE solver that implements the quantized state system methods: An approach that builds the solution by updating the system variables independently.
Author mongibellili
Popularity
4 Stars
Updated Last
3 Months Ago
Started In
December 2023

QuantizedSystemSolver

Dev CI Coverage

The growing intricacy of contemporary engineering systems, typically reduced to differential equations with events, poses a difficulty in digitally simulating them using traditional numerical integration techniques. The Quantized State System (QSS) and the Linearly Implicit Quantized State System (LIQSS) are different methods for tackling such problems. The QuantizedSystemSolver aims to solve a set of Ordinary differential equations with a set of events. It implements the quantized state system methods: An approach that builds the solution by updating the system variables independently as opposed to classic integration methods that update all the system variables every step.

Installation

Run Julia, enter ] to bring up Julia's package manager, and add the QuantizedSystemSolver.jl package:

julia> ]
 add QuantizedSystemSolver

Example: Buck circuit

The Buck is a converter that decreases voltage and increases current with a greater power efficiency than linear regulators. After a mesh analysis we get the problem discribed below. buckcircuit

Problem

The NLodeProblem function takes the following user code:

    odeprob = NLodeProblem(quote
          name=(buck,)
          #parameters
          C = 1e-4; L = 1e-4; R = 10.0;U = 24.0; T = 1e-4; DC = 0.5; ROn = 1e-5;ROff = 1e5;
          #discrete and continous variables
          discrete = [1e5,1e-5,1e-4,0.0,0.0];u = [0.0,0.0]
          #rename for convenience
          rd=discrete[1];rs=discrete[2];nextT=discrete[3];lastT=discrete[4];diodeon=discrete[5]
          il=u[1] ;uc=u[2]
          #helper equations
          id=(il*rs-U)/(rd+rs) # diode's current
          #differential equations
          du[1] =(-id*rd-uc)/L
          du[2]=(il-uc/R)/C
          #events 
          if t-nextT>0.0 
            lastT=nextT;nextT=nextT+T;rs=ROn
          end
          if t-lastT-DC*T>0.0 
            rs=ROff
          end                          
          if diodeon*(id)+(1.0-diodeon)*(id*rd-0.6)>0
            rd=ROn;diodeon=1.0
          else
            rd=ROff;diodeon=0.0
          end     
    end)

The output is an object that subtypes the

abstract type NLODEProblem{PRTYPE,T,Z,Y,CS} end

Solve

The solve function takes the previous problem (NLODEProblem{PRTYPE,T,Z,Y,CS}) with a chosen algorithm (ALGORITHM{N,O}) and some simulation settings:

tspan = (0.0, 0.001)
sol= solve(odeprob,nmliqss2(),tspan,abstol=1e-4,reltol=1e-3)

It outputs a solution of type Sol{T,O}.

Query the solution

# returns a vector of the values of all variables  at time 0.0005 # feature not available in v1.0.1
sol(0.0005)


# The value of variable 2  at time 0.0005
sol(2,0.0005)
19.209921627620943
# Get the value of variable 2 at time 0.0005
value_at_time = sol(0.0005,idxs=2)# feature not available in v1.0.1
# The total number of steps to end the simulation
sol.totalSteps # feature available only in v1.0.1 ->sol.stats in later versions
498
# The number of simultaneous steps during the simulation
sol.simulStepCount # feature available only in v1.0.1 ->sol.stats in later versions
132
# The total number of events during the simulation
sol.evCount # feature available only in v1.0.1 ->sol.stats in later versions
53
# The actual data is stored in two vectors:
sol.savedTimes
sol.savedVars

Plot the solution

# If the user wants to perform other tasks with the plot:
plot_Sol(sol)
# If the user wants to save the plot to a file:
save_Sol(sol)

plot_buck_nmLiqss2_()0 0001 _ft_0 001_2024_7_2_16_43_54

Error Analysis

#Compute the error against an analytic solution
error=getError(sol::Sol{T,O},index::Int,f::Function)
#Compute the error against a reference solution
error=getAverageErrorByRefs(sol,solRef::Vector{Any})

Used By Packages

No packages found.