Kinetica.jl API

Conditions

ConditionSet

Kinetica.ConditionSetType
conditions = ConditionSet(Dict(
    :C1 => ConditionType1(...),
    :C2 => ConditionType2(...))
    [, ts_update=nothing]
))

Container for all conditions in a kinetic simulation.

Conditions can be static or variable, and variable conditions can be gradient-based or directly usable.

Contains fields for:

  • Symbolic representation of conditions (symbols)
  • Condition profile for each symbol (profiles)
  • Whether discrete rate constant updates are enabled for the conditions in this condition set (discrete_updates)
  • Discrete rate constant update timestep, is nothing if discrete_updates = false (ts_update)

Constructor separates condition profiles from their symbols and parses numeric profiles into StaticConditionProfiles. Registers AbstractVariableProfiles with Symbolics to allow for proper computation down the chain.

If ts_update is provided as a keyword argument, creates tstops arrays within each variable profile for use in discrete rate update simulations.

source
Kinetica.get_initial_conditionsMethod
get_initial_conditions(conditions::ConditionSet)

Extract initial values of conditions from conditions.

Returns an array of Pairs linking Symbols to initial values. For AbstractStaticProfiles, initial values are their static values. For AbstractVariableProfiles, initial values are their X_start values.

source
Kinetica.isstaticFunction
isstatic(cs::ConditionSet[, sym::Symbol])

Determines if condition profiles in a ConditionSet are static.

When a Symbol sym is provided, only checks if the profile linked to this Symbol is static. If no Symbol is provided, checks is all profiles are static.

isstatic(profile<:AbstractConditionProfile)

Determines if a given condition profile is static.

source
Kinetica.isvariableFunction
isvariable(cs::ConditionSet[, sym::Symbol])

Determines if condition profiles in a ConditionSet are variable.

When a Symbol sym is provided, only checks if the profile linked to this Symbol is variable. If no Symbol is provided, checks is all profiles are variable.

isvariable(profile<:AbstractConditionProfile)

Determines if a given condition profile is variable.

source
Kinetica.get_tstopsMethod
get_tstops(cs::ConditionSet)

Retrieves a sorted array of unique time stops from all condition profiles in cs.

Should be used for passing a unified set of time stops to a discrete rate constant update solver. Will throw an error if all condition profiles are static, as they have no tstops.

source
Kinetica.get_t_finalMethod
get_t_final(cs::ConditionSet)

Retrieves the last necessary time point needed to encompass all variable condition profiles in cs.

Will throw an error if all condition profiles are static, they have no set endpoints.

source
Kinetica.solve_variable_conditions!Function
solve_variable_conditions!(cs::ConditionSet, pars::ODESimulationParams[, reset=false, solver=OwrenZen5(), solve_kwargs])

Solves all variable condition profiles over the timespan in pars.tspan.

Places all condition profile solutions in their sol field. In the case of AbstractDirectProfiles, this creates a DiffEqArray to mimic the regular DiffEq solver interface.

If condition profiles already exist, they will not be overwritten unless reset=true.

The solver and solve_kwargs arguments are used when solving gradient profiles. The OwrenZen5 solver has shown to be a stable, accurate solver capable of handling sudden gradient changes, and is a sensible default. solve_kwargs is a Dict of keyword arguments that get passed to the solve call, with default state:

solve_kwargs=Dict{Symbol, Any}(
    :abstol => 1e-6,
    :reltol => 1e-4
)

These have been shown to be sensible defaults for most gradient profiles.

source

Static Condition Profiles

Variable Condition Profiles

Kinetica.create_discrete_tstops!Function
create_discrete_tstops!(profile<:AbstractVariableProfile, ts_update<:AbstractFloat)

Creates a custom array of time stops within profile.tstops.

This array contains a time stop every ts_update, but attempts to intellegently avoid unnecessary time stops in areas where they are not needed, i.e. when the given profile is stationary.

source

Directly Variable Condition Profiles

Kinetica.solve_variable_condition!Method
solve_variable_condition!(profile<:AbstractDirectProfile, pars::ODESimulationParams[, reset=false])

Generates a solution for the specified directly-variable condition profile.

For profiles with direct functions, this requires calculating values for the specified pars.tspan and wrapping them within a DiffEqArray for compatibility with other interfaces (plotting, interpolation, etc.).

Arguments sym, solver and solve_kwargs are provided for compatibility with unified callers, do nothing and should be ignored.

source
Kinetica.NullDirectProfileType
NullDirectProfile(; X_start, t_end)

Container for null direct profile data and condition function.

This condition profile should only be used for debugging, as it has a condition function which always returns the initial condition. If only this constant condition is required, StaticODESolve should always be used with a StaticConditionProfile instead of a VariableODESolve with this condition profile.

Contains fields for:

  • Condition function (f)
  • Initial value of condition (X_start)
  • Time to stop calculation (t_end)
  • Times for the ODE solver to ensure calculation at (tstops)
  • Profile solution, constructed by call to solve_variable_condition! (sol)
source
Kinetica.NullDirectProfileMethod
NullDirectProfile(; X_start, t_end)

Container for null direct profile data and condition function.

This condition profile should only be used for debugging, as it has a condition function which always returns the initial condition. If only this constant condition is required, StaticODESolve should always be used with a StaticConditionProfile instead of a VariableODESolve with this condition profile.

Contains fields for:

  • Condition function (f)
  • Initial value of condition (X_start)
  • Time to stop calculation (t_end)
  • Times for the ODE solver to ensure calculation at (tstops)
  • Profile solution, constructed by call to solve_variable_condition! (sol)
source
Kinetica.LinearDirectProfileType
LinearDirectProfile(; rate, X_start, X_end)

Container for linear condition ramp profile data and condition function.

This condition profile represents a linear condition increase/decrease from X_start to X_end. Determines the simulation end time from the provided conditions and rate, then constructs the condition function (which is a linear y = mx + c function).

Contains fields for:

  • Condition function (f)
  • Rate of change of condition (rate)
  • Initial value of condition (X_start)
  • Final value of condition (X_end)
  • Time to stop calculation (t_end)
  • Times for the ODE solver to ensure calculation at (tstops)
  • Profile solution, constructed by call to solve_variable_condition! (sol)
source
Kinetica.LinearDirectProfileMethod
LinearDirectProfile(; rate, X_start, X_end)

Container for linear condition ramp profile data and condition function.

This condition profile represents a linear condition increase/decrease from X_start to X_end. Determines the simulation end time from the provided conditions and rate, then constructs the condition function (which is a linear y = mx + c function).

Contains fields for:

  • Condition function (f)
  • Rate of change of condition (rate)
  • Initial value of condition (X_start)
  • Final value of condition (X_end)
  • Time to stop calculation (t_end)
  • Times for the ODE solver to ensure calculation at (tstops)
  • Profile solution, constructed by call to solve_variable_condition! (sol)
source

Gradient-Variable Condition Profiles

Kinetica.solve_variable_condition!Method
solve_variable_condition!(profile<:AbstractGradientProfile, pars::ODESimulationParams[, sym=nothing, reset=false, solver, solve_kwargs])

Generates a solution for the specified gradient-variable condition profile.

For gradient-based profiles, this requires constructing an ODEProblem around their MTK-derived symbolic gradient expressions and solving over the timespan in pars.

The ODE solver and the arguments passed to the solve() call can be controlled with the solver and solve_kwargs arguments respectively. If sym is passed a Symbol, this will bind the solution result to that symbol in the underlying ODESolution.

source
Kinetica.LinearGradientProfileType
LinearGradientProfile(; rate, X_start, X_end)

Container for linear condition ramp profile data and condition gradient function.

This condition profile represents a linear condition increase/decrease from X_start to X_end. Determines the simulation end time from the provided conditions and gradient, then constructs the condition gradient function (which returns rate for every timestep).

Contains fields for:

  • Condition gradient function (grad)
  • Rate of change of condition (rate)
  • Initial value of condition (X_start)
  • Final value of condition (X_end)
  • Time to stop calculation (t_end)
  • Times for the ODE solver to ensure calculation at (tstops)
  • Profile solution, constructed by call to solve_variable_condition! (sol)
source
Kinetica.LinearGradientProfileMethod
LinearGradientProfile(; rate, X_start, X_end)

Container for linear condition ramp profile data and condition gradient function.

This condition profile represents a linear condition increase/decrease from X_start to X_end. Determines the simulation end time from the provided conditions and gradient, then constructs the condition gradient function (which returns rate for every timestep).

Contains fields for:

  • Condition gradient function (grad)
  • Rate of change of condition (rate)
  • Initial value of condition (X_start)
  • Final value of condition (X_end)
  • Time to stop calculation (t_end)
  • Times for the ODE solver to ensure calculation at (tstops)
  • Profile solution, constructed by call to solve_variable_condition! (sol)
source
Kinetica.DoubleRampGradientProfileType
DoubleRampGradientProfile(; X_start, t_start_plateau, rate1, X_mid, t_mid_plateau, rate2, X_end, t_end_plateau[, t_blend])

Container for double condition ramp profile data and condition gradient function.

This condition profile represents two condition ramps with adjustable condition plateaus before, after and in between the ramps, i.e.

              ------   X_mid
      rate1  /      \
            /        \  rate2
X_start ----          \
                       ----- X_end

The profile starts at X_start and maintains that value for t_start_plateau. The condition then ramps with gradient rate1 to condition value X_mid. This value is maintained for t_mid_plateau. The condition then ramps with gradient rate2 to condition value X_end. This value is maintained for t_end_plateau until the calculated time t_end.

To smooth out gradient discontinuities, a blending time t_blend can be passed to linearly interploate between plateaus and ramps, forming a smooth function of time. Larger values of t_blend yield smoother functions but decrease accuracy of the ramps, so should be used carefully.

source
Kinetica.DoubleRampGradientProfileMethod
DoubleRampGradientProfile(; X_start, t_start_plateau, rate1, X_mid, t_mid_plateau, rate2, X_end, t_end_plateau[, t_blend])

Container for double condition ramp profile data and condition gradient function.

This condition profile represents two condition ramps with adjustable condition plateaus before, after and in between the ramps, i.e.

              ------   X_mid
      rate1  /      \
            /        \  rate2
X_start ----          \
                       ----- X_end

The profile starts at X_start and maintains that value for t_start_plateau. The condition then ramps with gradient rate1 to condition value X_mid. This value is maintained for t_mid_plateau. The condition then ramps with gradient rate2 to condition value X_end. This value is maintained for t_end_plateau until the calculated time t_end.

To smooth out gradient discontinuities, a blending time t_blend can be passed to linearly interploate between plateaus and ramps, forming a smooth function of time. Larger values of t_blend yield smoother functions but decrease accuracy of the ramps, so should be used carefully.

source