Building and solving numerical problems
Systems are numerically solved by building and solving the appropriate problem type. Numerical solvers expect to receive functions taking a predefeined set of arguments and returning specific values. This format of argument and return value depends on the function and the problem. ModelingToolkit is capable of compiling and generating code for a variety of such numerical problems.
Dynamical systems
SciMLBase.ODEFunction
— TypeSciMLBase.ODEFunction(sys::System; kwargs...)
SciMLBase.ODEFunction{iip}(sys::System; kwargs...)
SciMLBase.ODEFunction{iip, specialize}(sys::System; kwargs...)
Create a SciMLBase.ODEFunction
from the given sys
. iip
is a boolean indicating whether the function should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.ODEFunction.
Keyword arguments
u0
: Theu0
vector for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_u0
.p
: The parameter object for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_p
.t
: The initial time for the corresponding problem, if available.eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.checkbounds
: Whether to enable bounds checking in the generated code.simplify
: Whether tosimplify
any symbolically computed jacobians/hessians/etc.cse
: Whether to enable Common Subexpression Elimination (CSE) on the generated code. This typically improves performance of the generated code but reduces readability.sparse
: Whether to generate jacobian/hessian/etc. functions that return/operate on sparse matrices. Also controls whether the mass matrix is sparse, wherever applicable.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.ODEFunction
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise. Constructing the expression does not support callbacksjac
: Whether to symbolically compute and generate code for the jacobian function.tgrad
: Whether to symbolically compute and generate code for thetgrad
function.sparsity
: Whether to provide symbolically compute and provide sparsity patterns for the jacobian/hessian/etc.
All other keyword arguments are forwarded to the SciMLBase.ODEFunction
struct constructor.
SciMLBase.ODEProblem
— TypeSciMLBase.SciMLBase.ODEProblem(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.SciMLBase.ODEProblem{iip}(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.SciMLBase.ODEProblem{iip, specialize}(sys::System, op, tspan::NTuple{2}; kwargs...)
Build a SciMLBase.ODEProblem
given a system sys
and operating point op
and timespan tspan
. iip
is a boolean indicating whether the problem should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.ODEFunction. The operating point should be an iterable collection of key-value pairs mapping variables/parameters in the system to the (initial) values they should take in SciMLBase.ODEProblem
. Any values not provided will fallback to the corresponding default (if present).
ModelingToolkit will build an initialization problem where all initial values for unknowns or observables of sys
(either explicitly provided or in defaults) will be constraints. To remove an initial condition in the defaults (without providing a replacement) give the corresponding variable a value of nothing
in the operating point. The initialization problem will also run parameter initialization. See the Initialization documentation for more information.
Keyword arguments
eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.check_initialization_units
: Enable or disable unit checks when constructing the initialization problem.tofloat
: Passed tovarmap_to_vars
when building the parameter vector of a non-split system.u0_eltype
: Theeltype
of theu0
vector. Ifnothing
, finds the promoted floating point type fromop
.u0_constructor
: A function to apply to theu0
value returned fromvarmap_to_vars
. to construct the finalu0
value.p_constructor
: A function to apply to each array buffer created when constructing the parameter object.warn_cyclic_dependency
: Whether to emit a warning listing out cycles in initial conditions provided for unknowns and parameters.circular_dependency_max_cycle_length
: Maximum length of cycle to check for. Only applicable ifwarn_cyclic_dependency == true
.circular_dependency_max_cycles
: Maximum number of cycles to check for. Only applicable ifwarn_cyclic_dependency == true
.substitution_limit
: The number times to substitute initial conditions into each other to attempt to arrive at a numeric value.callback
: An extra callback orCallbackSet
to add to the problem, in addition to the ones defined symbolically in the system.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.ODEProblem
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise. Constructing the expression does not support callbacks
All other keyword arguments are forwarded to the SciMLBase.ODEFunction constructor.
Extended docs
The following API is internal and may change or be removed without notice. Its usage is highly discouraged.
build_initializeprob
: Iffalse
, avoids building the initialization problem.check_length
: Whether to check the number of equations along with number of unknowns and length ofu0
vector for consistency. Iffalse
, do not check with equations. This is forwarded tocheck_eqs_u0
.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
SciMLBase.DAEFunction
— TypeSciMLBase.DAEFunction(sys::System; kwargs...)
SciMLBase.DAEFunction{iip}(sys::System; kwargs...)
SciMLBase.DAEFunction{iip, specialize}(sys::System; kwargs...)
Create a SciMLBase.DAEFunction
from the given sys
. iip
is a boolean indicating whether the function should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.DAEFunction.
Keyword arguments
u0
: Theu0
vector for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_u0
.p
: The parameter object for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_p
.t
: The initial time for the corresponding problem, if available.eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.checkbounds
: Whether to enable bounds checking in the generated code.simplify
: Whether tosimplify
any symbolically computed jacobians/hessians/etc.cse
: Whether to enable Common Subexpression Elimination (CSE) on the generated code. This typically improves performance of the generated code but reduces readability.sparse
: Whether to generate jacobian/hessian/etc. functions that return/operate on sparse matrices. Also controls whether the mass matrix is sparse, wherever applicable.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.DAEFunction
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise. Constructing the expression does not support callbacksjac
: Whether to symbolically compute and generate code for the jacobian function.tgrad
: Whether to symbolically compute and generate code for thetgrad
function.sparsity
: Whether to provide symbolically compute and provide sparsity patterns for the jacobian/hessian/etc.
All other keyword arguments are forwarded to the SciMLBase.DAEFunction
struct constructor.
SciMLBase.DAEProblem
— TypeSciMLBase.SciMLBase.DAEProblem(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.SciMLBase.DAEProblem{iip}(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.SciMLBase.DAEProblem{iip, specialize}(sys::System, op, tspan::NTuple{2}; kwargs...)
Build a SciMLBase.DAEProblem
given a system sys
and operating point op
and timespan tspan
. iip
is a boolean indicating whether the problem should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.DAEFunction. The operating point should be an iterable collection of key-value pairs mapping variables/parameters in the system to the (initial) values they should take in SciMLBase.DAEProblem
. Any values not provided will fallback to the corresponding default (if present).
ModelingToolkit will build an initialization problem where all initial values for unknowns or observables of sys
(either explicitly provided or in defaults) will be constraints. To remove an initial condition in the defaults (without providing a replacement) give the corresponding variable a value of nothing
in the operating point. The initialization problem will also run parameter initialization. See the Initialization documentation for more information.
Keyword arguments
eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.check_initialization_units
: Enable or disable unit checks when constructing the initialization problem.tofloat
: Passed tovarmap_to_vars
when building the parameter vector of a non-split system.u0_eltype
: Theeltype
of theu0
vector. Ifnothing
, finds the promoted floating point type fromop
.u0_constructor
: A function to apply to theu0
value returned fromvarmap_to_vars
. to construct the finalu0
value.p_constructor
: A function to apply to each array buffer created when constructing the parameter object.warn_cyclic_dependency
: Whether to emit a warning listing out cycles in initial conditions provided for unknowns and parameters.circular_dependency_max_cycle_length
: Maximum length of cycle to check for. Only applicable ifwarn_cyclic_dependency == true
.circular_dependency_max_cycles
: Maximum number of cycles to check for. Only applicable ifwarn_cyclic_dependency == true
.substitution_limit
: The number times to substitute initial conditions into each other to attempt to arrive at a numeric value.callback
: An extra callback orCallbackSet
to add to the problem, in addition to the ones defined symbolically in the system.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.DAEProblem
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise. Constructing the expression does not support callbacks
All other keyword arguments are forwarded to the SciMLBase.DAEFunction constructor.
Extended docs
The following API is internal and may change or be removed without notice. Its usage is highly discouraged.
build_initializeprob
: Iffalse
, avoids building the initialization problem.check_length
: Whether to check the number of equations along with number of unknowns and length ofu0
vector for consistency. Iffalse
, do not check with equations. This is forwarded tocheck_eqs_u0
.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
SciMLBase.SDEFunction
— TypeSciMLBase.SDEFunction(sys::System; kwargs...)
SciMLBase.SDEFunction{iip}(sys::System; kwargs...)
SciMLBase.SDEFunction{iip, specialize}(sys::System; kwargs...)
Create a SciMLBase.SDEFunction
from the given sys
. iip
is a boolean indicating whether the function should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.SDEFunction.
Keyword arguments
u0
: Theu0
vector for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_u0
.p
: The parameter object for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_p
.t
: The initial time for the corresponding problem, if available.eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.checkbounds
: Whether to enable bounds checking in the generated code.simplify
: Whether tosimplify
any symbolically computed jacobians/hessians/etc.cse
: Whether to enable Common Subexpression Elimination (CSE) on the generated code. This typically improves performance of the generated code but reduces readability.sparse
: Whether to generate jacobian/hessian/etc. functions that return/operate on sparse matrices. Also controls whether the mass matrix is sparse, wherever applicable.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.SDEFunction
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise. Constructing the expression does not support callbacksjac
: Whether to symbolically compute and generate code for the jacobian function.tgrad
: Whether to symbolically compute and generate code for thetgrad
function.sparsity
: Whether to provide symbolically compute and provide sparsity patterns for the jacobian/hessian/etc.
All other keyword arguments are forwarded to the SciMLBase.SDEFunction
struct constructor.
SciMLBase.SDEProblem
— TypeSciMLBase.SciMLBase.SDEProblem(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.SciMLBase.SDEProblem{iip}(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.SciMLBase.SDEProblem{iip, specialize}(sys::System, op, tspan::NTuple{2}; kwargs...)
Build a SciMLBase.SDEProblem
given a system sys
and operating point op
and timespan tspan
. iip
is a boolean indicating whether the problem should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.SDEFunction. The operating point should be an iterable collection of key-value pairs mapping variables/parameters in the system to the (initial) values they should take in SciMLBase.SDEProblem
. Any values not provided will fallback to the corresponding default (if present).
ModelingToolkit will build an initialization problem where all initial values for unknowns or observables of sys
(either explicitly provided or in defaults) will be constraints. To remove an initial condition in the defaults (without providing a replacement) give the corresponding variable a value of nothing
in the operating point. The initialization problem will also run parameter initialization. See the Initialization documentation for more information.
Keyword arguments
eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.check_initialization_units
: Enable or disable unit checks when constructing the initialization problem.tofloat
: Passed tovarmap_to_vars
when building the parameter vector of a non-split system.u0_eltype
: Theeltype
of theu0
vector. Ifnothing
, finds the promoted floating point type fromop
.u0_constructor
: A function to apply to theu0
value returned fromvarmap_to_vars
. to construct the finalu0
value.p_constructor
: A function to apply to each array buffer created when constructing the parameter object.warn_cyclic_dependency
: Whether to emit a warning listing out cycles in initial conditions provided for unknowns and parameters.circular_dependency_max_cycle_length
: Maximum length of cycle to check for. Only applicable ifwarn_cyclic_dependency == true
.circular_dependency_max_cycles
: Maximum number of cycles to check for. Only applicable ifwarn_cyclic_dependency == true
.substitution_limit
: The number times to substitute initial conditions into each other to attempt to arrive at a numeric value.callback
: An extra callback orCallbackSet
to add to the problem, in addition to the ones defined symbolically in the system.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.SDEProblem
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise. Constructing the expression does not support callbacks
All other keyword arguments are forwarded to the SciMLBase.SDEFunction constructor.
Extended docs
The following API is internal and may change or be removed without notice. Its usage is highly discouraged.
build_initializeprob
: Iffalse
, avoids building the initialization problem.check_length
: Whether to check the number of equations along with number of unknowns and length ofu0
vector for consistency. Iffalse
, do not check with equations. This is forwarded tocheck_eqs_u0
.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
SciMLBase.DDEFunction
— TypeSciMLBase.DDEFunction(sys::System; kwargs...)
SciMLBase.DDEFunction{iip}(sys::System; kwargs...)
SciMLBase.DDEFunction{iip, specialize}(sys::System; kwargs...)
Create a SciMLBase.DDEFunction
from the given sys
. iip
is a boolean indicating whether the function should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.DDEFunction.
Keyword arguments
u0
: Theu0
vector for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_u0
.p
: The parameter object for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_p
.t
: The initial time for the corresponding problem, if available.eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.checkbounds
: Whether to enable bounds checking in the generated code.simplify
: Whether tosimplify
any symbolically computed jacobians/hessians/etc.cse
: Whether to enable Common Subexpression Elimination (CSE) on the generated code. This typically improves performance of the generated code but reduces readability.sparse
: Whether to generate jacobian/hessian/etc. functions that return/operate on sparse matrices. Also controls whether the mass matrix is sparse, wherever applicable.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.DDEFunction
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise. Constructing the expression does not support callbacks
All other keyword arguments are forwarded to the SciMLBase.DDEFunction
struct constructor.
SciMLBase.DDEProblem
— TypeSciMLBase.SciMLBase.DDEProblem(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.SciMLBase.DDEProblem{iip}(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.SciMLBase.DDEProblem{iip, specialize}(sys::System, op, tspan::NTuple{2}; kwargs...)
Build a SciMLBase.DDEProblem
given a system sys
and operating point op
and timespan tspan
. iip
is a boolean indicating whether the problem should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.DDEFunction. The operating point should be an iterable collection of key-value pairs mapping variables/parameters in the system to the (initial) values they should take in SciMLBase.DDEProblem
. Any values not provided will fallback to the corresponding default (if present).
ModelingToolkit will build an initialization problem where all initial values for unknowns or observables of sys
(either explicitly provided or in defaults) will be constraints. To remove an initial condition in the defaults (without providing a replacement) give the corresponding variable a value of nothing
in the operating point. The initialization problem will also run parameter initialization. See the Initialization documentation for more information.
Keyword arguments
eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.check_initialization_units
: Enable or disable unit checks when constructing the initialization problem.tofloat
: Passed tovarmap_to_vars
when building the parameter vector of a non-split system.u0_eltype
: Theeltype
of theu0
vector. Ifnothing
, finds the promoted floating point type fromop
.u0_constructor
: A function to apply to theu0
value returned fromvarmap_to_vars
. to construct the finalu0
value.p_constructor
: A function to apply to each array buffer created when constructing the parameter object.warn_cyclic_dependency
: Whether to emit a warning listing out cycles in initial conditions provided for unknowns and parameters.circular_dependency_max_cycle_length
: Maximum length of cycle to check for. Only applicable ifwarn_cyclic_dependency == true
.circular_dependency_max_cycles
: Maximum number of cycles to check for. Only applicable ifwarn_cyclic_dependency == true
.substitution_limit
: The number times to substitute initial conditions into each other to attempt to arrive at a numeric value.callback
: An extra callback orCallbackSet
to add to the problem, in addition to the ones defined symbolically in the system.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.DDEProblem
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise. Constructing the expression does not support callbacks
All other keyword arguments are forwarded to the SciMLBase.DDEFunction constructor.
Extended docs
The following API is internal and may change or be removed without notice. Its usage is highly discouraged.
build_initializeprob
: Iffalse
, avoids building the initialization problem.check_length
: Whether to check the number of equations along with number of unknowns and length ofu0
vector for consistency. Iffalse
, do not check with equations. This is forwarded tocheck_eqs_u0
.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
SciMLBase.SDDEFunction
— TypeSciMLBase.SDDEFunction(sys::System; kwargs...)
SciMLBase.SDDEFunction{iip}(sys::System; kwargs...)
SciMLBase.SDDEFunction{iip, specialize}(sys::System; kwargs...)
Create a SciMLBase.SDDEFunction
from the given sys
. iip
is a boolean indicating whether the function should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.SDDEFunction.
Keyword arguments
u0
: Theu0
vector for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_u0
.p
: The parameter object for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_p
.t
: The initial time for the corresponding problem, if available.eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.checkbounds
: Whether to enable bounds checking in the generated code.simplify
: Whether tosimplify
any symbolically computed jacobians/hessians/etc.cse
: Whether to enable Common Subexpression Elimination (CSE) on the generated code. This typically improves performance of the generated code but reduces readability.sparse
: Whether to generate jacobian/hessian/etc. functions that return/operate on sparse matrices. Also controls whether the mass matrix is sparse, wherever applicable.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.SDDEFunction
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise. Constructing the expression does not support callbacks
All other keyword arguments are forwarded to the SciMLBase.SDDEFunction
struct constructor.
SciMLBase.SDDEProblem
— TypeSciMLBase.SciMLBase.SDDEProblem(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.SciMLBase.SDDEProblem{iip}(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.SciMLBase.SDDEProblem{iip, specialize}(sys::System, op, tspan::NTuple{2}; kwargs...)
Build a SciMLBase.SDDEProblem
given a system sys
and operating point op
and timespan tspan
. iip
is a boolean indicating whether the problem should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.SDDEFunction. The operating point should be an iterable collection of key-value pairs mapping variables/parameters in the system to the (initial) values they should take in SciMLBase.SDDEProblem
. Any values not provided will fallback to the corresponding default (if present).
ModelingToolkit will build an initialization problem where all initial values for unknowns or observables of sys
(either explicitly provided or in defaults) will be constraints. To remove an initial condition in the defaults (without providing a replacement) give the corresponding variable a value of nothing
in the operating point. The initialization problem will also run parameter initialization. See the Initialization documentation for more information.
Keyword arguments
eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.check_initialization_units
: Enable or disable unit checks when constructing the initialization problem.tofloat
: Passed tovarmap_to_vars
when building the parameter vector of a non-split system.u0_eltype
: Theeltype
of theu0
vector. Ifnothing
, finds the promoted floating point type fromop
.u0_constructor
: A function to apply to theu0
value returned fromvarmap_to_vars
. to construct the finalu0
value.p_constructor
: A function to apply to each array buffer created when constructing the parameter object.warn_cyclic_dependency
: Whether to emit a warning listing out cycles in initial conditions provided for unknowns and parameters.circular_dependency_max_cycle_length
: Maximum length of cycle to check for. Only applicable ifwarn_cyclic_dependency == true
.circular_dependency_max_cycles
: Maximum number of cycles to check for. Only applicable ifwarn_cyclic_dependency == true
.substitution_limit
: The number times to substitute initial conditions into each other to attempt to arrive at a numeric value.callback
: An extra callback orCallbackSet
to add to the problem, in addition to the ones defined symbolically in the system.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.SDDEProblem
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise. Constructing the expression does not support callbacks
All other keyword arguments are forwarded to the SciMLBase.SDDEFunction constructor.
Extended docs
The following API is internal and may change or be removed without notice. Its usage is highly discouraged.
build_initializeprob
: Iffalse
, avoids building the initialization problem.check_length
: Whether to check the number of equations along with number of unknowns and length ofu0
vector for consistency. Iffalse
, do not check with equations. This is forwarded tocheck_eqs_u0
.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
JumpProcesses.JumpProblem
— TypeSciMLBase.JumpProcesses.JumpProblem(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.JumpProcesses.JumpProblem{iip}(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.JumpProcesses.JumpProblem{iip, specialize}(sys::System, op, tspan::NTuple{2}; kwargs...)
Build a JumpProcesses.JumpProblem
given a system sys
and operating point op
and timespan tspan
. iip
is a boolean indicating whether the problem should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the inner SciMLFunction. The operating point should be an iterable collection of key-value pairs mapping variables/parameters in the system to the (initial) values they should take in JumpProcesses.JumpProblem
. Any values not provided will fallback to the corresponding default (if present).
ModelingToolkit will build an initialization problem where all initial values for unknowns or observables of sys
(either explicitly provided or in defaults) will be constraints. To remove an initial condition in the defaults (without providing a replacement) give the corresponding variable a value of nothing
in the operating point. The initialization problem will also run parameter initialization. See the Initialization documentation for more information.
Keyword arguments
eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.check_initialization_units
: Enable or disable unit checks when constructing the initialization problem.tofloat
: Passed tovarmap_to_vars
when building the parameter vector of a non-split system.u0_eltype
: Theeltype
of theu0
vector. Ifnothing
, finds the promoted floating point type fromop
.u0_constructor
: A function to apply to theu0
value returned fromvarmap_to_vars
. to construct the finalu0
value.p_constructor
: A function to apply to each array buffer created when constructing the parameter object.warn_cyclic_dependency
: Whether to emit a warning listing out cycles in initial conditions provided for unknowns and parameters.circular_dependency_max_cycle_length
: Maximum length of cycle to check for. Only applicable ifwarn_cyclic_dependency == true
.circular_dependency_max_cycles
: Maximum number of cycles to check for. Only applicable ifwarn_cyclic_dependency == true
.substitution_limit
: The number times to substitute initial conditions into each other to attempt to arrive at a numeric value.callback
: An extra callback orCallbackSet
to add to the problem, in addition to the ones defined symbolically in the system.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aJumpProcesses.JumpProblem
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise. Constructing the expression does not support callbacks
All other keyword arguments are forwarded to the inner SciMLFunction constructor.
Extended docs
The following API is internal and may change or be removed without notice. Its usage is highly discouraged.
build_initializeprob
: Iffalse
, avoids building the initialization problem.check_length
: Whether to check the number of equations along with number of unknowns and length ofu0
vector for consistency. Iffalse
, do not check with equations. This is forwarded tocheck_eqs_u0
.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
SciMLBase.BVProblem
— TypeSciMLBase.SciMLBase.BVProblem(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.SciMLBase.BVProblem{iip}(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.SciMLBase.BVProblem{iip, specialize}(sys::System, op, tspan::NTuple{2}; kwargs...)
Build a SciMLBase.BVProblem
given a system sys
and operating point op
and timespan tspan
. iip
is a boolean indicating whether the problem should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.ODEFunction. The operating point should be an iterable collection of key-value pairs mapping variables/parameters in the system to the (initial) values they should take in SciMLBase.BVProblem
. Any values not provided will fallback to the corresponding default (if present).
ModelingToolkit will build an initialization problem where all initial values for unknowns or observables of sys
(either explicitly provided or in defaults) will be constraints. To remove an initial condition in the defaults (without providing a replacement) give the corresponding variable a value of nothing
in the operating point. The initialization problem will also run parameter initialization. See the Initialization documentation for more information.
Keyword arguments
eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.check_initialization_units
: Enable or disable unit checks when constructing the initialization problem.tofloat
: Passed tovarmap_to_vars
when building the parameter vector of a non-split system.u0_eltype
: Theeltype
of theu0
vector. Ifnothing
, finds the promoted floating point type fromop
.u0_constructor
: A function to apply to theu0
value returned fromvarmap_to_vars
. to construct the finalu0
value.p_constructor
: A function to apply to each array buffer created when constructing the parameter object.warn_cyclic_dependency
: Whether to emit a warning listing out cycles in initial conditions provided for unknowns and parameters.circular_dependency_max_cycle_length
: Maximum length of cycle to check for. Only applicable ifwarn_cyclic_dependency == true
.circular_dependency_max_cycles
: Maximum number of cycles to check for. Only applicable ifwarn_cyclic_dependency == true
.substitution_limit
: The number times to substitute initial conditions into each other to attempt to arrive at a numeric value.callback
: An extra callback orCallbackSet
to add to the problem, in addition to the ones defined symbolically in the system.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.BVProblem
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise. Constructing the expression does not support callbacks
All other keyword arguments are forwarded to the SciMLBase.ODEFunction constructor.
Extended docs
The following API is internal and may change or be removed without notice. Its usage is highly discouraged.
build_initializeprob
: Iffalse
, avoids building the initialization problem.check_length
: Whether to check the number of equations along with number of unknowns and length ofu0
vector for consistency. Iffalse
, do not check with equations. This is forwarded tocheck_eqs_u0
.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
SciMLBase.DiscreteProblem
— TypeSciMLBase.SciMLBase.DiscreteProblem(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.SciMLBase.DiscreteProblem{iip}(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.SciMLBase.DiscreteProblem{iip, specialize}(sys::System, op, tspan::NTuple{2}; kwargs...)
Build a SciMLBase.DiscreteProblem
given a system sys
and operating point op
and timespan tspan
. iip
is a boolean indicating whether the problem should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.DiscreteFunction. The operating point should be an iterable collection of key-value pairs mapping variables/parameters in the system to the (initial) values they should take in SciMLBase.DiscreteProblem
. Any values not provided will fallback to the corresponding default (if present).
ModelingToolkit will build an initialization problem where all initial values for unknowns or observables of sys
(either explicitly provided or in defaults) will be constraints. To remove an initial condition in the defaults (without providing a replacement) give the corresponding variable a value of nothing
in the operating point. The initialization problem will also run parameter initialization. See the Initialization documentation for more information.
Keyword arguments
eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.check_initialization_units
: Enable or disable unit checks when constructing the initialization problem.tofloat
: Passed tovarmap_to_vars
when building the parameter vector of a non-split system.u0_eltype
: Theeltype
of theu0
vector. Ifnothing
, finds the promoted floating point type fromop
.u0_constructor
: A function to apply to theu0
value returned fromvarmap_to_vars
. to construct the finalu0
value.p_constructor
: A function to apply to each array buffer created when constructing the parameter object.warn_cyclic_dependency
: Whether to emit a warning listing out cycles in initial conditions provided for unknowns and parameters.circular_dependency_max_cycle_length
: Maximum length of cycle to check for. Only applicable ifwarn_cyclic_dependency == true
.circular_dependency_max_cycles
: Maximum number of cycles to check for. Only applicable ifwarn_cyclic_dependency == true
.substitution_limit
: The number times to substitute initial conditions into each other to attempt to arrive at a numeric value.callback
: An extra callback orCallbackSet
to add to the problem, in addition to the ones defined symbolically in the system.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.DiscreteProblem
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise. Constructing the expression does not support callbacks
All other keyword arguments are forwarded to the SciMLBase.DiscreteFunction constructor.
Extended docs
The following API is internal and may change or be removed without notice. Its usage is highly discouraged.
build_initializeprob
: Iffalse
, avoids building the initialization problem.check_length
: Whether to check the number of equations along with number of unknowns and length ofu0
vector for consistency. Iffalse
, do not check with equations. This is forwarded tocheck_eqs_u0
.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
SciMLBase.ImplicitDiscreteProblem
— TypeSciMLBase.SciMLBase.ImplicitDiscreteProblem(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.SciMLBase.ImplicitDiscreteProblem{iip}(sys::System, op, tspan::NTuple{2}; kwargs...)
SciMLBase.SciMLBase.ImplicitDiscreteProblem{iip, specialize}(sys::System, op, tspan::NTuple{2}; kwargs...)
Build a SciMLBase.ImplicitDiscreteProblem
given a system sys
and operating point op
and timespan tspan
. iip
is a boolean indicating whether the problem should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.ImplicitDiscreteFunction. The operating point should be an iterable collection of key-value pairs mapping variables/parameters in the system to the (initial) values they should take in SciMLBase.ImplicitDiscreteProblem
. Any values not provided will fallback to the corresponding default (if present).
ModelingToolkit will build an initialization problem where all initial values for unknowns or observables of sys
(either explicitly provided or in defaults) will be constraints. To remove an initial condition in the defaults (without providing a replacement) give the corresponding variable a value of nothing
in the operating point. The initialization problem will also run parameter initialization. See the Initialization documentation for more information.
Keyword arguments
eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.check_initialization_units
: Enable or disable unit checks when constructing the initialization problem.tofloat
: Passed tovarmap_to_vars
when building the parameter vector of a non-split system.u0_eltype
: Theeltype
of theu0
vector. Ifnothing
, finds the promoted floating point type fromop
.u0_constructor
: A function to apply to theu0
value returned fromvarmap_to_vars
. to construct the finalu0
value.p_constructor
: A function to apply to each array buffer created when constructing the parameter object.warn_cyclic_dependency
: Whether to emit a warning listing out cycles in initial conditions provided for unknowns and parameters.circular_dependency_max_cycle_length
: Maximum length of cycle to check for. Only applicable ifwarn_cyclic_dependency == true
.circular_dependency_max_cycles
: Maximum number of cycles to check for. Only applicable ifwarn_cyclic_dependency == true
.substitution_limit
: The number times to substitute initial conditions into each other to attempt to arrive at a numeric value.callback
: An extra callback orCallbackSet
to add to the problem, in addition to the ones defined symbolically in the system.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.ImplicitDiscreteProblem
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise. Constructing the expression does not support callbacks
All other keyword arguments are forwarded to the SciMLBase.ImplicitDiscreteFunction constructor.
Extended docs
The following API is internal and may change or be removed without notice. Its usage is highly discouraged.
build_initializeprob
: Iffalse
, avoids building the initialization problem.check_length
: Whether to check the number of equations along with number of unknowns and length ofu0
vector for consistency. Iffalse
, do not check with equations. This is forwarded tocheck_eqs_u0
.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
Linear and Nonlinear systems
SciMLBase.NonlinearFunction
— TypeSciMLBase.NonlinearFunction(sys::System; kwargs...)
SciMLBase.NonlinearFunction{iip}(sys::System; kwargs...)
SciMLBase.NonlinearFunction{iip, specialize}(sys::System; kwargs...)
Create a SciMLBase.NonlinearFunction
from the given sys
. iip
is a boolean indicating whether the function should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.NonlinearFunction.
Keyword arguments
u0
: Theu0
vector for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_u0
.p
: The parameter object for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_p
.eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.checkbounds
: Whether to enable bounds checking in the generated code.simplify
: Whether tosimplify
any symbolically computed jacobians/hessians/etc.cse
: Whether to enable Common Subexpression Elimination (CSE) on the generated code. This typically improves performance of the generated code but reduces readability.sparse
: Whether to generate jacobian/hessian/etc. functions that return/operate on sparse matrices. Also controls whether the mass matrix is sparse, wherever applicable.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.NonlinearFunction
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise.resid_prototype
: The prototype of the residual functionf
for a problem involving a nonlinear solve where the residual andu0
have different sizes.jac
: Whether to symbolically compute and generate code for the jacobian function.sparsity
: Whether to provide symbolically compute and provide sparsity patterns for the jacobian/hessian/etc.
All other keyword arguments are forwarded to the SciMLBase.NonlinearFunction
struct constructor.
SciMLBase.NonlinearProblem
— TypeSciMLBase.SciMLBase.NonlinearProblem(sys::System, op; kwargs...)
SciMLBase.SciMLBase.NonlinearProblem{iip}(sys::System, op; kwargs...)
SciMLBase.SciMLBase.NonlinearProblem{iip, specialize}(sys::System, op; kwargs...)
Build a SciMLBase.NonlinearProblem
given a system sys
and operating point op
. iip
is a boolean indicating whether the problem should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.NonlinearFunction. The operating point should be an iterable collection of key-value pairs mapping variables/parameters in the system to the (initial) values they should take in SciMLBase.NonlinearProblem
. Any values not provided will fallback to the corresponding default (if present).
ModelingToolkit will build an initialization problem that will run parameter initialization. Since it does not solve for initial values of unknowns, observed equations will not be initialization constraints. If an initialization equation of the system must involve the initial value of an unknown x
, it must be used as Initial(x)
in the equation. For example, an equation to be used to solve for parameter p
in terms of unknowns x
and y
must be provided as Initial(x) + Initial(y) ~ p
instead of x + y ~ p
. See the Initialization documentation for more information.
Keyword arguments
eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.check_initialization_units
: Enable or disable unit checks when constructing the initialization problem.tofloat
: Passed tovarmap_to_vars
when building the parameter vector of a non-split system.u0_eltype
: Theeltype
of theu0
vector. Ifnothing
, finds the promoted floating point type fromop
.u0_constructor
: A function to apply to theu0
value returned fromvarmap_to_vars
. to construct the finalu0
value.p_constructor
: A function to apply to each array buffer created when constructing the parameter object.warn_cyclic_dependency
: Whether to emit a warning listing out cycles in initial conditions provided for unknowns and parameters.circular_dependency_max_cycle_length
: Maximum length of cycle to check for. Only applicable ifwarn_cyclic_dependency == true
.circular_dependency_max_cycles
: Maximum number of cycles to check for. Only applicable ifwarn_cyclic_dependency == true
.substitution_limit
: The number times to substitute initial conditions into each other to attempt to arrive at a numeric value.
check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.NonlinearProblem
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise.
All other keyword arguments are forwarded to the SciMLBase.NonlinearFunction constructor.
Extended docs
The following API is internal and may change or be removed without notice. Its usage is highly discouraged.
build_initializeprob
: Iffalse
, avoids building the initialization problem.check_length
: Whether to check the number of equations along with number of unknowns and length ofu0
vector for consistency. Iffalse
, do not check with equations. This is forwarded tocheck_eqs_u0
.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
SciMLBase.SCCNonlinearProblem
— TypeSciMLBase.SciMLBase.SCCNonlinearProblem(sys::System, op; kwargs...)
SciMLBase.SciMLBase.SCCNonlinearProblem{iip}(sys::System, op; kwargs...)
SciMLBase.SciMLBase.SCCNonlinearProblem{iip, specialize}(sys::System, op; kwargs...)
Build a SciMLBase.SCCNonlinearProblem
given a system sys
and operating point op
. iip
is a boolean indicating whether the problem should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.NonlinearFunction. The operating point should be an iterable collection of key-value pairs mapping variables/parameters in the system to the (initial) values they should take in SciMLBase.SCCNonlinearProblem
. Any values not provided will fallback to the corresponding default (if present).
ModelingToolkit will build an initialization problem that will run parameter initialization. Since it does not solve for initial values of unknowns, observed equations will not be initialization constraints. If an initialization equation of the system must involve the initial value of an unknown x
, it must be used as Initial(x)
in the equation. For example, an equation to be used to solve for parameter p
in terms of unknowns x
and y
must be provided as Initial(x) + Initial(y) ~ p
instead of x + y ~ p
. See the Initialization documentation for more information.
Keyword arguments
eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.check_initialization_units
: Enable or disable unit checks when constructing the initialization problem.tofloat
: Passed tovarmap_to_vars
when building the parameter vector of a non-split system.u0_eltype
: Theeltype
of theu0
vector. Ifnothing
, finds the promoted floating point type fromop
.u0_constructor
: A function to apply to theu0
value returned fromvarmap_to_vars
. to construct the finalu0
value.p_constructor
: A function to apply to each array buffer created when constructing the parameter object.warn_cyclic_dependency
: Whether to emit a warning listing out cycles in initial conditions provided for unknowns and parameters.circular_dependency_max_cycle_length
: Maximum length of cycle to check for. Only applicable ifwarn_cyclic_dependency == true
.circular_dependency_max_cycles
: Maximum number of cycles to check for. Only applicable ifwarn_cyclic_dependency == true
.substitution_limit
: The number times to substitute initial conditions into each other to attempt to arrive at a numeric value.
check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.SCCNonlinearProblem
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise.
All other keyword arguments are forwarded to the SciMLBase.NonlinearFunction constructor.
Extended docs
The following API is internal and may change or be removed without notice. Its usage is highly discouraged.
build_initializeprob
: Iffalse
, avoids building the initialization problem.check_length
: Whether to check the number of equations along with number of unknowns and length ofu0
vector for consistency. Iffalse
, do not check with equations. This is forwarded tocheck_eqs_u0
.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
SciMLBase.NonlinearLeastSquaresProblem
— TypeSciMLBase.SciMLBase.NonlinearLeastSquaresProblem(sys::System, op; kwargs...)
SciMLBase.SciMLBase.NonlinearLeastSquaresProblem{iip}(sys::System, op; kwargs...)
SciMLBase.SciMLBase.NonlinearLeastSquaresProblem{iip, specialize}(sys::System, op; kwargs...)
Build a SciMLBase.NonlinearLeastSquaresProblem
given a system sys
and operating point op
. iip
is a boolean indicating whether the problem should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.NonlinearFunction. The operating point should be an iterable collection of key-value pairs mapping variables/parameters in the system to the (initial) values they should take in SciMLBase.NonlinearLeastSquaresProblem
. Any values not provided will fallback to the corresponding default (if present).
ModelingToolkit will build an initialization problem that will run parameter initialization. Since it does not solve for initial values of unknowns, observed equations will not be initialization constraints. If an initialization equation of the system must involve the initial value of an unknown x
, it must be used as Initial(x)
in the equation. For example, an equation to be used to solve for parameter p
in terms of unknowns x
and y
must be provided as Initial(x) + Initial(y) ~ p
instead of x + y ~ p
. See the Initialization documentation for more information.
Keyword arguments
eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.check_initialization_units
: Enable or disable unit checks when constructing the initialization problem.tofloat
: Passed tovarmap_to_vars
when building the parameter vector of a non-split system.u0_eltype
: Theeltype
of theu0
vector. Ifnothing
, finds the promoted floating point type fromop
.u0_constructor
: A function to apply to theu0
value returned fromvarmap_to_vars
. to construct the finalu0
value.p_constructor
: A function to apply to each array buffer created when constructing the parameter object.warn_cyclic_dependency
: Whether to emit a warning listing out cycles in initial conditions provided for unknowns and parameters.circular_dependency_max_cycle_length
: Maximum length of cycle to check for. Only applicable ifwarn_cyclic_dependency == true
.circular_dependency_max_cycles
: Maximum number of cycles to check for. Only applicable ifwarn_cyclic_dependency == true
.substitution_limit
: The number times to substitute initial conditions into each other to attempt to arrive at a numeric value.
check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.NonlinearLeastSquaresProblem
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise.
All other keyword arguments are forwarded to the SciMLBase.NonlinearFunction constructor.
Extended docs
The following API is internal and may change or be removed without notice. Its usage is highly discouraged.
build_initializeprob
: Iffalse
, avoids building the initialization problem.check_length
: Whether to check the number of equations along with number of unknowns and length ofu0
vector for consistency. Iffalse
, do not check with equations. This is forwarded tocheck_eqs_u0
.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
SciMLBase.SteadyStateProblem
— TypeSciMLBase.SciMLBase.SteadyStateProblem(sys::System, op; kwargs...)
SciMLBase.SciMLBase.SteadyStateProblem{iip}(sys::System, op; kwargs...)
SciMLBase.SciMLBase.SteadyStateProblem{iip, specialize}(sys::System, op; kwargs...)
Build a SciMLBase.SteadyStateProblem
given a system sys
and operating point op
. iip
is a boolean indicating whether the problem should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.ODEFunction. The operating point should be an iterable collection of key-value pairs mapping variables/parameters in the system to the (initial) values they should take in SciMLBase.SteadyStateProblem
. Any values not provided will fallback to the corresponding default (if present).
ModelingToolkit will build an initialization problem that will run parameter initialization. Since it does not solve for initial values of unknowns, observed equations will not be initialization constraints. If an initialization equation of the system must involve the initial value of an unknown x
, it must be used as Initial(x)
in the equation. For example, an equation to be used to solve for parameter p
in terms of unknowns x
and y
must be provided as Initial(x) + Initial(y) ~ p
instead of x + y ~ p
. See the Initialization documentation for more information.
Keyword arguments
eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.check_initialization_units
: Enable or disable unit checks when constructing the initialization problem.tofloat
: Passed tovarmap_to_vars
when building the parameter vector of a non-split system.u0_eltype
: Theeltype
of theu0
vector. Ifnothing
, finds the promoted floating point type fromop
.u0_constructor
: A function to apply to theu0
value returned fromvarmap_to_vars
. to construct the finalu0
value.p_constructor
: A function to apply to each array buffer created when constructing the parameter object.warn_cyclic_dependency
: Whether to emit a warning listing out cycles in initial conditions provided for unknowns and parameters.circular_dependency_max_cycle_length
: Maximum length of cycle to check for. Only applicable ifwarn_cyclic_dependency == true
.circular_dependency_max_cycles
: Maximum number of cycles to check for. Only applicable ifwarn_cyclic_dependency == true
.substitution_limit
: The number times to substitute initial conditions into each other to attempt to arrive at a numeric value.
check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.SteadyStateProblem
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise.
All other keyword arguments are forwarded to the SciMLBase.ODEFunction constructor.
Extended docs
The following API is internal and may change or be removed without notice. Its usage is highly discouraged.
build_initializeprob
: Iffalse
, avoids building the initialization problem.check_length
: Whether to check the number of equations along with number of unknowns and length ofu0
vector for consistency. Iffalse
, do not check with equations. This is forwarded tocheck_eqs_u0
.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
SciMLBase.IntervalNonlinearFunction
— TypeSciMLBase.IntervalNonlinearFunction(sys::System; kwargs...)
SciMLBase.IntervalNonlinearFunction{iip}(sys::System; kwargs...)
SciMLBase.IntervalNonlinearFunction{iip, specialize}(sys::System; kwargs...)
Create a SciMLBase.IntervalNonlinearFunction
from the given sys
. iip
is a boolean indicating whether the function should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.IntervalNonlinearFunction.
Keyword arguments
u0
: Theu0
vector for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_u0
.p
: The parameter object for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_p
.eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.checkbounds
: Whether to enable bounds checking in the generated code.simplify
: Whether tosimplify
any symbolically computed jacobians/hessians/etc.cse
: Whether to enable Common Subexpression Elimination (CSE) on the generated code. This typically improves performance of the generated code but reduces readability.sparse
: Whether to generate jacobian/hessian/etc. functions that return/operate on sparse matrices. Also controls whether the mass matrix is sparse, wherever applicable.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.IntervalNonlinearFunction
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise.
All other keyword arguments are forwarded to the SciMLBase.IntervalNonlinearFunction
struct constructor.
SciMLBase.IntervalNonlinearProblem
— TypeSciMLBase.IntervalNonlinearProblem(sys::System, uspan::NTuple{2}, parammap = SciMLBase.NullParameters(); kwargs...)
Create an IntervalNonlinearProblem
from the given sys
. This is only valid for a system of nonlinear equations with a single equation and unknown. uspan
is the interval in which the root is to be found, and parammap
is an iterable collection of key-value pairs providing values for the parameters in the system.
ModelingToolkit will build an initialization problem that will run parameter initialization. Since it does not solve for initial values of unknowns, observed equations will not be initialization constraints. If an initialization equation of the system must involve the initial value of an unknown x
, it must be used as Initial(x)
in the equation. For example, an equation to be used to solve for parameter p
in terms of unknowns x
and y
must be provided as Initial(x) + Initial(y) ~ p
instead of x + y ~ p
. See the Initialization documentation for more information.
Keyword arguments
eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.check_initialization_units
: Enable or disable unit checks when constructing the initialization problem.tofloat
: Passed tovarmap_to_vars
when building the parameter vector of a non-split system.u0_eltype
: Theeltype
of theu0
vector. Ifnothing
, finds the promoted floating point type fromop
.u0_constructor
: A function to apply to theu0
value returned fromvarmap_to_vars
. to construct the finalu0
value.p_constructor
: A function to apply to each array buffer created when constructing the parameter object.warn_cyclic_dependency
: Whether to emit a warning listing out cycles in initial conditions provided for unknowns and parameters.circular_dependency_max_cycle_length
: Maximum length of cycle to check for. Only applicable ifwarn_cyclic_dependency == true
.circular_dependency_max_cycles
: Maximum number of cycles to check for. Only applicable ifwarn_cyclic_dependency == true
.substitution_limit
: The number times to substitute initial conditions into each other to attempt to arrive at a numeric value.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.IntervalNonlinearProblem
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise.
All other keyword arguments are forwarded to the IntervalNonlinearFunction
constructor.
Extended docs
The following API is internal and may change or be removed without notice. Its usage is highly discouraged.
build_initializeprob
: Iffalse
, avoids building the initialization problem.check_length
: Whether to check the number of equations along with number of unknowns and length ofu0
vector for consistency. Iffalse
, do not check with equations. This is forwarded tocheck_eqs_u0
.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
ModelingToolkit.HomotopyContinuationProblem
— TypeSciMLBase.ModelingToolkit.HomotopyContinuationProblem(sys::System, op; kwargs...)
SciMLBase.ModelingToolkit.HomotopyContinuationProblem{iip}(sys::System, op; kwargs...)
SciMLBase.ModelingToolkit.HomotopyContinuationProblem{iip, specialize}(sys::System, op; kwargs...)
Build a ModelingToolkit.HomotopyContinuationProblem
given a system sys
and operating point op
. iip
is a boolean indicating whether the problem should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.HomotopyNonlinearFunction. The operating point should be an iterable collection of key-value pairs mapping variables/parameters in the system to the (initial) values they should take in ModelingToolkit.HomotopyContinuationProblem
. Any values not provided will fallback to the corresponding default (if present).
Keyword arguments
eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.check_initialization_units
: Enable or disable unit checks when constructing the initialization problem.tofloat
: Passed tovarmap_to_vars
when building the parameter vector of a non-split system.u0_eltype
: Theeltype
of theu0
vector. Ifnothing
, finds the promoted floating point type fromop
.u0_constructor
: A function to apply to theu0
value returned fromvarmap_to_vars
. to construct the finalu0
value.p_constructor
: A function to apply to each array buffer created when constructing the parameter object.warn_cyclic_dependency
: Whether to emit a warning listing out cycles in initial conditions provided for unknowns and parameters.circular_dependency_max_cycle_length
: Maximum length of cycle to check for. Only applicable ifwarn_cyclic_dependency == true
.circular_dependency_max_cycles
: Maximum number of cycles to check for. Only applicable ifwarn_cyclic_dependency == true
.substitution_limit
: The number times to substitute initial conditions into each other to attempt to arrive at a numeric value.
check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aModelingToolkit.HomotopyContinuationProblem
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise.
All other keyword arguments are forwarded to the SciMLBase.HomotopyNonlinearFunction constructor.
Extended docs
The following API is internal and may change or be removed without notice. Its usage is highly discouraged.
build_initializeprob
: Iffalse
, avoids building the initialization problem.check_length
: Whether to check the number of equations along with number of unknowns and length ofu0
vector for consistency. Iffalse
, do not check with equations. This is forwarded tocheck_eqs_u0
.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
SciMLBase.HomotopyNonlinearFunction
— TypeSciMLBase.HomotopyNonlinearFunction(sys::System; kwargs...)
SciMLBase.HomotopyNonlinearFunction{iip}(sys::System; kwargs...)
SciMLBase.HomotopyNonlinearFunction{iip, specialize}(sys::System; kwargs...)
Create a HomotopyNonlinearFunction
from the given sys
. iip
is a boolean indicating whether the function should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the func.
Keyword arguments
u0
: Theu0
vector for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_u0
.p
: The parameter object for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_p
.eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.checkbounds
: Whether to enable bounds checking in the generated code.simplify
: Whether tosimplify
any symbolically computed jacobians/hessians/etc.cse
: Whether to enable Common Subexpression Elimination (CSE) on the generated code. This typically improves performance of the generated code but reduces readability.fraction_cancel_fn
: The function to use to simplify fractions in the polynomial expression. A more powerful function can increase processing time but be able to eliminate more rational functions, thus improving solve time. Should be a function that takes a symbolic expression containing zero or more fraction expressions and returns the simplified expression. While this defaults toSymbolicUtils.simplify_fractions
, a viable alternative isSymbolicUtils.quick_cancel
All keyword arguments are forwarded to the wrapped NonlinearFunction
constructor.
SciMLBase.LinearProblem
— TypeSciMLBase.LinearProblem(sys::System, op; kwargs...)
SciMLBase.LinearProblem{iip}(sys::System, op; kwargs...)
Build a LinearProblem
given a system sys
and operating point op
. iip
is a boolean indicating whether the problem should be in-place. The operating point should be an iterable collection of key-value pairs mapping variables/parameters in the system to the (initial) values they should take in LinearProblem
. Any values not provided will fallback to the corresponding default (if present).
Note that since u0
is optional for LinearProblem
, values of unknowns do not need to be specified in op
to create a LinearProblem
. In such a case, prob.u0
will be nothing
and attempting to symbolically index the problem with an unknown, observable, or expression depending on unknowns/observables will error.
Updating the parameters automatically updates the A
and b
arrays.
Keyword arguments
eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.check_initialization_units
: Enable or disable unit checks when constructing the initialization problem.tofloat
: Passed tovarmap_to_vars
when building the parameter vector of a non-split system.u0_eltype
: Theeltype
of theu0
vector. Ifnothing
, finds the promoted floating point type fromop
.u0_constructor
: A function to apply to theu0
value returned fromvarmap_to_vars
. to construct the finalu0
value.p_constructor
: A function to apply to each array buffer created when constructing the parameter object.warn_cyclic_dependency
: Whether to emit a warning listing out cycles in initial conditions provided for unknowns and parameters.circular_dependency_max_cycle_length
: Maximum length of cycle to check for. Only applicable ifwarn_cyclic_dependency == true
.circular_dependency_max_cycles
: Maximum number of cycles to check for. Only applicable ifwarn_cyclic_dependency == true
.substitution_limit
: The number times to substitute initial conditions into each other to attempt to arrive at a numeric value.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.LinearProblem
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise.
All other keyword arguments are forwarded to the func constructor.
Extended docs
The following API is internal and may change or be removed without notice. Its usage is highly discouraged.
build_initializeprob
: Iffalse
, avoids building the initialization problem.check_length
: Whether to check the number of equations along with number of unknowns and length ofu0
vector for consistency. Iffalse
, do not check with equations. This is forwarded tocheck_eqs_u0
.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
Optimization and optimal control
SciMLBase.OptimizationFunction
— TypeSciMLBase.OptimizationFunction(sys::System; kwargs...)
SciMLBase.OptimizationFunction{iip}(sys::System; kwargs...)
SciMLBase.OptimizationFunction{iip, specialize}(sys::System; kwargs...)
Create a SciMLBase.OptimizationFunction
from the given sys
. iip
is a boolean indicating whether the function should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.OptimizationFunction.
Keyword arguments
u0
: Theu0
vector for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_u0
.p
: The parameter object for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_p
.eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.checkbounds
: Whether to enable bounds checking in the generated code.simplify
: Whether tosimplify
any symbolically computed jacobians/hessians/etc.cse
: Whether to enable Common Subexpression Elimination (CSE) on the generated code. This typically improves performance of the generated code but reduces readability.sparse
: Whether to generate jacobian/hessian/etc. functions that return/operate on sparse matrices. Also controls whether the mass matrix is sparse, wherever applicable.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.OptimizationFunction
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise.jac
: Whether to symbolically compute and generate code for the jacobian function.grad
: Whether the symbolically compute and generate code for the gradient of the cost function with respect to unknowns.hess
: Whether to symbolically compute and generate code for the hessian function.cons_h
: Whether to symbolically compute and generate code for the hessian function of constraints. Since the constraint function is vector-valued, the hessian is a vector of hessian matrices.cons_j
: Whether to symbolically compute and generate code for the jacobian function of constraints.sparsity
: Whether to provide symbolically compute and provide sparsity patterns for the jacobian/hessian/etc.cons_sparse
: Identical to thesparse
keyword, but specifically for jacobian/hessian functions of the constraints.
All other keyword arguments are forwarded to the SciMLBase.OptimizationFunction
struct constructor.
SciMLBase.OptimizationProblem
— TypeSciMLBase.SciMLBase.OptimizationProblem(sys::System, op; kwargs...)
SciMLBase.SciMLBase.OptimizationProblem{iip}(sys::System, op; kwargs...)
SciMLBase.SciMLBase.OptimizationProblem{iip, specialize}(sys::System, op; kwargs...)
Build a SciMLBase.OptimizationProblem
given a system sys
and operating point op
. iip
is a boolean indicating whether the problem should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.OptimizationFunction. The operating point should be an iterable collection of key-value pairs mapping variables/parameters in the system to the (initial) values they should take in SciMLBase.OptimizationProblem
. Any values not provided will fallback to the corresponding default (if present).
ModelingToolkit will build an initialization problem that will run parameter initialization. Since it does not solve for initial values of unknowns, observed equations will not be initialization constraints. If an initialization equation of the system must involve the initial value of an unknown x
, it must be used as Initial(x)
in the equation. For example, an equation to be used to solve for parameter p
in terms of unknowns x
and y
must be provided as Initial(x) + Initial(y) ~ p
instead of x + y ~ p
. See the Initialization documentation for more information.
Keyword arguments
eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.check_initialization_units
: Enable or disable unit checks when constructing the initialization problem.tofloat
: Passed tovarmap_to_vars
when building the parameter vector of a non-split system.u0_eltype
: Theeltype
of theu0
vector. Ifnothing
, finds the promoted floating point type fromop
.u0_constructor
: A function to apply to theu0
value returned fromvarmap_to_vars
. to construct the finalu0
value.p_constructor
: A function to apply to each array buffer created when constructing the parameter object.warn_cyclic_dependency
: Whether to emit a warning listing out cycles in initial conditions provided for unknowns and parameters.circular_dependency_max_cycle_length
: Maximum length of cycle to check for. Only applicable ifwarn_cyclic_dependency == true
.circular_dependency_max_cycles
: Maximum number of cycles to check for. Only applicable ifwarn_cyclic_dependency == true
.substitution_limit
: The number times to substitute initial conditions into each other to attempt to arrive at a numeric value.
check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.OptimizationProblem
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise.
All other keyword arguments are forwarded to the SciMLBase.OptimizationFunction constructor.
Extended docs
The following API is internal and may change or be removed without notice. Its usage is highly discouraged.
build_initializeprob
: Iffalse
, avoids building the initialization problem.check_length
: Whether to check the number of equations along with number of unknowns and length ofu0
vector for consistency. Iffalse
, do not check with equations. This is forwarded tocheck_eqs_u0
.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
SciMLBase.ODEInputFunction
— TypeSciMLBase.ODEInputFunction(sys::System; kwargs...)
SciMLBase.ODEInputFunction{iip}(sys::System; kwargs...)
SciMLBase.ODEInputFunction{iip, specialize}(sys::System; kwargs...)
Create a SciMLBase.ODEInputFunction
from the given sys
. iip
is a boolean indicating whether the function should be in-place. specialization
is a SciMLBase.AbstractSpecalize
subtype indicating the level of specialization of the SciMLBase.ODEInputFunction.
Keyword arguments
u0
: Theu0
vector for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_u0
.p
: The parameter object for the corresponding problem, if available. Can be obtained usingModelingToolkit.get_p
.t
: The initial time for the corresponding problem, if available.eval_expression
: Whether to compile any functions viaeval
orRuntimeGeneratedFunctions
.eval_module
: Ifeval_expression == true
, the module toeval
into. Otherwise, the module in which to generate theRuntimeGeneratedFunction
.checkbounds
: Whether to enable bounds checking in the generated code.simplify
: Whether tosimplify
any symbolically computed jacobians/hessians/etc.cse
: Whether to enable Common Subexpression Elimination (CSE) on the generated code. This typically improves performance of the generated code but reduces readability.sparse
: Whether to generate jacobian/hessian/etc. functions that return/operate on sparse matrices. Also controls whether the mass matrix is sparse, wherever applicable.check_compatibility
: Whether to check if the given systemsys
contains all the information necessary to create aSciMLBase.ODEInputFunction
and no more. If disabled, assumes thatsys
at least contains the necessary information.expression
:Val{true}
to return anExpr
that constructs the corresponding problem instead of the problem itself.Val{false}
otherwise. Constructing the expression does not support callbacksinputs
: The variables in the input vector. The system must have been simplified usingmtkcompile
with these variables passed asinputs
.disturbance_inputs
: The disturbance input variables. The system must have been simplified usingmtkcompile
with these variables passed asdisturbance_inputs
.jac
: Whether to symbolically compute and generate code for the jacobian function.tgrad
: Whether to symbolically compute and generate code for thetgrad
function.controljac
: Whether to symbolically compute and generate code for the jacobian of the ODE with respect to the inputs.sparsity
: Whether to provide symbolically compute and provide sparsity patterns for the jacobian/hessian/etc.
All other keyword arguments are forwarded to the SciMLBase.ODEInputFunction
struct constructor.
ModelingToolkit.JuMPDynamicOptProblem
— FunctionJuMPDynamicOptProblem(sys::System, op, tspan; dt, steps, guesses, kwargs...)
Convert an System representing an optimal control system into a JuMP model for solving using optimization. Must provide either dt
, the timestep between collocation points (which, along with the timespan, determines the number of points), or directly provide the number of points as steps
.
To construct the problem, please load InfiniteOpt along with ModelingToolkit.
ModelingToolkit.InfiniteOptDynamicOptProblem
— FunctionInfiniteOptDynamicOptProblem(sys::System, op, tspan; dt)
Convert an System representing an optimal control system into a InfiniteOpt model for solving using optimization. Must provide dt
for determining the length of the interpolation arrays.
Related to JuMPDynamicOptProblem
, but directly adds the differential equations of the system as derivative constraints, rather than using a solver tableau.
To construct the problem, please load InfiniteOpt along with ModelingToolkit.
ModelingToolkit.CasADiDynamicOptProblem
— FunctionCasADiDynamicOptProblem(sys::System, op, tspan; dt, steps, guesses, kwargs...)
Convert an System representing an optimal control system into a CasADi model for solving using optimization. Must provide either dt
, the timestep between collocation points (which, along with the timespan, determines the number of points), or directly provide the number of points as steps
.
To construct the problem, please load CasADi along with ModelingToolkit.
Missing docstring for ModelingToolkit.DynamicOptSolution
. Check Documenter's build log for details.
The state vector and parameter object
Typically the unknowns of the system are present as a Vector
of the appropriate length in the numerical problem. The state vector can also be constructed manually without building a problem.
ModelingToolkit.get_u0
— Functionget_u0(
sys::ModelingToolkit.AbstractSystem,
varmap;
kwargs...
) -> Any
Return the u0
vector for the given system sys
and variable-value mapping varmap
. All keyword arguments are forwarded to varmap_to_vars
.
ModelingToolkit.varmap_to_vars
— Functionvarmap_to_vars(
varmap::AbstractDict,
vars::Vector;
tofloat,
use_union,
container_type,
buffer_eltype,
toterm,
check,
allow_symbolic,
is_initializeprob,
substitution_limit
) -> Any
Return an array of values where the i
th element corresponds to the value of vars[i]
in varmap
. Will mutate varmap
by symbolically substituting it into itself.
Keyword arguments:
container_type
: The type of the returned container.allow_symbolic
: Whether the returned container of values can have symbolic expressions.buffer_eltype
: Theeltype
of the returned container if!allow_symbolic
. IfNothing
, automatically promotes the values in the container to a commoneltype
.tofloat
: Whether to promote values to floating point numbers ifbuffer_eltype == Nothing
.use_union
: Whether to allow using aUnion
as theeltype
ifbuffer_eltype == Nothing
.toterm
: Thetoterm
function for canonicalizing keys ofvarmap
. A value ofnothing
disables this process.check
: Whether to check if all ofvars
are keys ofvarmap
.is_initializeprob
: Whether an initialization problem is being constructed. Used for better error messages.substitution_limit
: The maximum number of times to recursively substitutevarmap
into itself to get a numeric value for each variable invars
.
By default, the parameters of the system are stored in a custom data structure called MTKParameters
. The internals of this data structure are undocumented, and it should only be interacted with through defined public API. SymbolicIndexingInterface.jl contains functionality useful for this purpose.
ModelingToolkit.MTKParameters
— Typefunction MTKParameters(sys::AbstractSystem, p, u0 = Dict(); t0 = nothing)
Create an MTKParameters
object for the system sys
. p
(u0
) are symbolic maps from parameters (unknowns) to their values. The values can also be symbolic expressions, which are evaluated given the values of other parameters/unknowns. u0
is only required if the values of parameters depend on the unknowns. t0
is the initial time, for time- dependent systems. It is only required if the symbolic expressions also use the independent variable of the system.
This requires that complete
has been called on the system (usually via mtkcompile
or @mtkcompile
) and the keyword split = true
was passed (which is the default behavior).
ModelingToolkit.get_p
— Functionget_p(
sys::ModelingToolkit.AbstractSystem,
varmap;
split,
kwargs...
) -> Any
Return the p
object for the given system sys
and variable-value mapping varmap
. All keyword arguments are forwarded to MTKParameters
for split systems and varmap_to_vars
for non-split systems.
The following functions are useful when working with MTKParameters
objects, and especially the Tunables
portion. For more information about the "portions" of MTKParameters
, refer to the SciMLStructures.jl
documentation.
ModelingToolkit.reorder_dimension_by_tunables!
— Functionreorder_dimension_by_tunables!(dest::AbstractArray, sys::AbstractSystem, arr::AbstractArray, syms; dim = 1)
Assuming the order of values in dimension dim
of arr
correspond to the order of tunable parameters in the system, reorder them according to the order described in syms
. syms
must be a permutation of tunable_parameters(sys)
. The result is written to dest
. The size
of dest
and arr
must be equal. Return dest
.
See also: MTKParameters
, tunable_parameters
, reorder_dimension_by_tunables
.
ModelingToolkit.reorder_dimension_by_tunables
— Functionreorder_dimension_by_tunables(sys::AbstractSystem, arr::AbstractArray, syms; dim = 1)
Out-of-place version of reorder_dimension_by_tunables!
.
Initialization
ModelingToolkit.generate_initializesystem
— Functiongenerate_initializesystem(
sys::ModelingToolkit.AbstractSystem;
time_dependent_init,
kwargs...
) -> System
Generate the initialization system for sys
. The initialization system is a system of nonlinear equations that solve for the full set of initial conditions of sys
given specified constraints.
The initialization system can be of two types: time-dependent and time-independent. Time-dependent initialization systems solve for the initial values of unknowns as well as the values of solvable parameters of the system. Time-independent initialization systems only solve for solvable parameters of the system.
Keyword arguments
time_dependent_init
: Whether to create an initialization system for a time-dependent system. A time-dependent initialization requires a time-dependentsys
, but a time- independent initialization can be created regardless.op
: The operating point of user-specified initial conditions of variables insys
.initialization_eqs
: Additional initialization equations to use apart from those ininitialization_equations(sys)
.guesses
: Additional guesses to use apart from those inguesses(sys)
.default_dd_guess
: Default guess for dummy derivative variables in time-dependent initialization.algebraic_only
: Iffalse
, does not use initialization equations (provided via the keyword or part of the system) to construct initialization.check_defguess
: Whether to error when a variable does not have a default or guess despite ModelingToolkit expecting it to.name
: The name of the initialization system.
All other keyword arguments are forwarded to the System
constructor.
ModelingToolkit.InitializationProblem
— TypeInitializationProblem(sys::AbstractSystem, t, op = Dict(); kwargs...)
InitializationProblem{iip}(sys::AbstractSystem, t, op = Dict(); kwargs...)
InitializationProblem{iip, specialize}(sys::AbstractSystem, t, op = Dict(); kwargs...)
Generate a NonlinearProblem
, SCCNonlinearProblem
or NonlinearLeastSquaresProblem
to represent a consistent initialization of sys
given the initial time t
and operating point op
. The initial time can be nothing
for time-independent systems.
Keyword arguments
guesses
: The guesses for variables in the system, used as initial values for the initialization problem.warn_initialize_determined
: Warn if the initialization system is under/over-determined.initialization_eqs
: Extra equations to use in the initialization problem.fully_determined
: Override whether the initialization system is fully determined.use_scc
: Whether to useSCCNonlinearProblem
for initialization if the system is fully determined.time_dependent_init
: Whether to build a time-dependent initialization for the problem. A time-dependent initialization solves for a consistentu0
, whereas a time-independent one only runs parameter initialization.algebraic_only
: Whether to build the initialization problem using only algebraic equations.allow_incomplete
: Whether to allow incomplete initialization problems.
All other keyword arguments are forwarded to the wrapped nonlinear problem constructor.
Linear analysis
ModelingToolkit.linearization_function
— Functionlin_fun, simplified_sys = linearization_function(sys::AbstractSystem, inputs, outputs; simplify = false, initialize = true, initialization_solver_alg = TrustRegion(), kwargs...)
Return a function that linearizes the system sys
. The function linearize
provides a higher-level and easier to use interface.
lin_fun
is a function (variables, p, t) -> (; f_x, f_z, g_x, g_z, f_u, g_u, h_x, h_z, h_u)
, i.e., it returns a NamedTuple with the Jacobians of f,g,h
for the nonlinear sys
(technically for simplified_sys
) on the form
\[\begin{aligned} ẋ &= f(x, z, u) \\ 0 &= g(x, z, u) \\ y &= h(x, z, u) \end{aligned}\]
where x
are differential unknown variables, z
algebraic variables, u
inputs and y
outputs. To obtain a linear statespace representation, see linearize
. The input argument variables
is a vector defining the operating point, corresponding to unknowns(simplified_sys)
and p
is a vector corresponding to the parameters of simplified_sys
. Note: all variables in inputs
have been converted to parameters in simplified_sys
.
The simplified_sys
has undergone mtkcompile
and had any occurring input or output variables replaced with the variables provided in arguments inputs
and outputs
. The unknowns of this system also indicate the order of the unknowns that holds for the linearized matrices.
Arguments:
sys
: ASystem
of ODEs. This function will automatically apply simplification passes onsys
and return the resultingsimplified_sys
.inputs
: A vector of variables that indicate the inputs of the linearized input-output model.outputs
: A vector of variables that indicate the outputs of the linearized input-output model.simplify
: Apply simplification in tearing.initialize
: If true, a check is performed to ensure that the operating point is consistent (satisfies algebraic equations). If the op is not consistent, initialization is performed.initialization_solver_alg
: A NonlinearSolve algorithm to use for solving for a feasible set of state and algebraic variables that satisfies the specified operating point.autodiff
: AnADType
supported by DifferentiationInterface.jl to use for calculating the necessary jacobians. Defaults to usingAutoForwardDiff()
kwargs
: Are passed on tofind_solvables!
See also linearize
which provides a higher-level interface.
ModelingToolkit.LinearizationProblem
— Typemutable struct LinearizationProblem{F<:ModelingToolkit.LinearizationFunction, T}
A struct representing a linearization operation to be performed. Can be symbolically indexed to efficiently update the operating point for multiple linearizations in a loop. The value of the independent variable can be set by mutating the .t
field of this struct.
ModelingToolkit.linearize
— Function(; A, B, C, D), simplified_sys, extras = linearize(sys, inputs, outputs; t=0.0, op = Dict(), allow_input_derivatives = false, zero_dummy_der=false, kwargs...)
(; A, B, C, D), extras = linearize(simplified_sys, lin_fun; t=0.0, op = Dict(), allow_input_derivatives = false, zero_dummy_der=false)
Linearize sys
between inputs
and outputs
, both vectors of variables. Return a NamedTuple with the matrices of a linear statespace representation on the form
\[\begin{aligned} ẋ &= Ax + Bu\\ y &= Cx + Du \end{aligned}\]
The first signature automatically calls linearization_function
internally, while the second signature expects the outputs of linearization_function
as input.
op
denotes the operating point around which to linearize. If none is provided, the default values of sys
are used.
If allow_input_derivatives = false
, an error will be thrown if input derivatives ($u̇$) appear as inputs in the linearized equations. If input derivatives are allowed, the returned B
matrix will be of double width, corresponding to the input [u; u̇]
.
zero_dummy_der
can be set to automatically set the operating point to zero for all dummy derivatives.
The return value extras
is a NamedTuple (; x, p, t)
containing the result of the initialization problem that was solved to determine the operating point.
See also linearization_function
which provides a lower-level interface, linearize_symbolic
and ModelingToolkit.reorder_unknowns
.
See extended help for an example.
The implementation and notation follows that of "Linear Analysis Approach for Modelica Models", Allain et al. 2009
Extended help
This example builds the following feedback interconnection and linearizes it from the input of F
to the output of P
.
r ┌─────┐ ┌─────┐ ┌─────┐
───►│ ├──────►│ │ u │ │
│ F │ │ C ├────►│ P │ y
└─────┘ ┌►│ │ │ ├─┬─►
│ └─────┘ └─────┘ │
│ │
└─────────────────────┘
using ModelingToolkit
using ModelingToolkit: t_nounits as t, D_nounits as D
function plant(; name)
@variables x(t) = 1
@variables u(t)=0 y(t)=0
eqs = [D(x) ~ -x + u
y ~ x]
System(eqs, t; name = name)
end
function ref_filt(; name)
@variables x(t)=0 y(t)=0
@variables u(t)=0 [input = true]
eqs = [D(x) ~ -2 * x + u
y ~ x]
System(eqs, t, name = name)
end
function controller(kp; name)
@variables y(t)=0 r(t)=0 u(t)=0
@parameters kp = kp
eqs = [
u ~ kp * (r - y),
]
System(eqs, t; name = name)
end
@named f = ref_filt()
@named c = controller(1)
@named p = plant()
connections = [f.y ~ c.r # filtered reference to controller reference
c.u ~ p.u # controller output to plant input
p.y ~ c.y]
@named cl = System(connections, t, systems = [f, c, p])
lsys0, ssys = linearize(cl, [f.u], [p.x])
desired_order = [f.x, p.x]
lsys = ModelingToolkit.reorder_unknowns(lsys0, unknowns(ssys), desired_order)
@assert lsys.A == [-2 0; 1 -2]
@assert lsys.B == [1; 0;;]
@assert lsys.C == [0 1]
@assert lsys.D[] == 0
## Symbolic linearization
lsys_sym, _ = ModelingToolkit.linearize_symbolic(cl, [f.u], [p.x])
@assert substitute(lsys_sym.A, ModelingToolkit.defaults(cl)) == lsys.A
Missing docstring for CommonSolve.solve(::LinearizationProblem)
. Check Documenter's build log for details.
Missing docstring for linearize_symbolic
. Check Documenter's build log for details.
There are also utilities for manipulating the results of these analyses in a symbolic context.
ModelingToolkit.similarity_transform
— Function(; Ã, B̃, C̃, D̃) = similarity_transform(sys, T; unitary=false)
Perform a similarity transform T : Tx̃ = x
on linear system represented by matrices in NamedTuple sys
such that
à = T⁻¹AT
B̃ = T⁻¹ B
C̃ = CT
D̃ = D
If unitary=true
, T
is assumed unitary and the matrix adjoint is used instead of the inverse.
ModelingToolkit.reorder_unknowns
— Functionreorder_unknowns(sys::NamedTuple, old, new)
Permute the state representation of sys
obtained from linearize
so that the state unknown is changed from old
to new
Example:
lsys, ssys = linearize(pid, [reference.u, measurement.u], [ctr_output.u])
desired_order = [int.x, der.x] # Unknowns that are present in unknowns(ssys)
lsys = ModelingToolkit.reorder_unknowns(lsys, unknowns(ssys), desired_order)
See also ModelingToolkit.similarity_transform
Analysis point transformations
Linear analysis can also be done using analysis points to perform several common workflows.
ModelingToolkit.get_sensitivity_function
— Functionget_sensitivity_function(
sys::ModelingToolkit.AbstractSystem,
aps;
kwargs...
) -> Tuple{ModelingToolkit.LinearizationFunction{DI, AI, _A, P, _B, _C, J1, J2, J3, J4, IA, @NamedTuple{abstol::Float64, reltol::Float64, nlsolve_alg::NonlinearSolveFirstOrder.GeneralizedFirstOrderAlgorithm{Missing, NonlinearSolveFirstOrder.GenericTrustRegionScheme{NonlinearSolveFirstOrder.RadiusUpdateSchemes.__Simple, Rational{Int64}, Rational{Int64}, Rational{Int64}, Rational{Int64}, Rational{Int64}, Rational{Int64}, Rational{Int64}}, NonlinearSolveBase.Dogleg{NonlinearSolveBase.NewtonDescent{Nothing}, NonlinearSolveBase.SteepestDescent{Nothing}}, Nothing, Nothing, Nothing, Val{false}}}} where {DI<:AbstractVector{Int64}, AI<:AbstractVector{Int64}, _A, P<:ODEProblem, _B, _C, J1<:Union{Nothing, ModelingToolkit.PreparedJacobian{true, DifferentiationInterfaceForwardDiffExt.ForwardDiffTwoArgJacobianPrep{Nothing, C, Tuple{Nothing, Nothing}}, ModelingToolkit.var"#uff#957"{var"#1479#fun"}, _A, ADTypes.AutoForwardDiff{nothing, Nothing}} where {C<:(ForwardDiff.JacobianConfig{T, _A, _B, <:Tuple{Any, Any}} where {T<:(ForwardDiff.Tag{F} where F<:ModelingToolkit.var"#uff#957"), _A, _B}), var"#1479#fun", _A}}, J2<:Union{Nothing, ModelingToolkit.PreparedJacobian{true, DifferentiationInterfaceForwardDiffExt.ForwardDiffTwoArgJacobianPrep{Nothing, C, Tuple{Nothing, Nothing}}, _A, _B, ADTypes.AutoForwardDiff{nothing, Nothing}} where {C<:ForwardDiff.JacobianConfig, _A, _B}}, J3<:Union{Nothing, ModelingToolkit.PreparedJacobian{true, DifferentiationInterfaceForwardDiffExt.ForwardDiffTwoArgJacobianPrep{Nothing, C, Tuple{Nothing, Nothing, Nothing}}, ModelingToolkit.var"#pff#958"{var"#1480#fun", var"#1481#setter"}, _A, ADTypes.AutoForwardDiff{nothing, Nothing}} where {C<:(ForwardDiff.JacobianConfig{T, _A, _B, <:Tuple{Any, Any}} where {T<:(ForwardDiff.Tag{F} where F<:ModelingToolkit.var"#pff#958"), _A, _B}), var"#1480#fun", var"#1481#setter", _A}}, J4<:(ModelingToolkit.PreparedJacobian{true, DifferentiationInterfaceForwardDiffExt.ForwardDiffTwoArgJacobianPrep{Nothing, C, Tuple{Nothing, Nothing, Nothing}}, ModelingToolkit.var"#hpf#956"{fun, setter}, _A, ADTypes.AutoForwardDiff{nothing, Nothing}} where {C<:(ForwardDiff.JacobianConfig{T, _A, _B, <:Tuple{Any, Any}} where {T<:(ForwardDiff.Tag{F} where F<:ModelingToolkit.var"#hpf#956"), _A, _B}), fun, setter, _A}), IA<:Union{SciMLBase.NoInit, SciMLBase.OverrideInit{Nothing, Nothing, Nothing}}}, Any}
Return the sensitivity function for the analysis point(s) aps
, and the modified system simplified with the appropriate inputs and outputs.
Keyword Arguments
- `loop_openings`: A list of analysis points whose connections should be removed and
the outputs set to the input as a part of the linear analysis.
- `system_modifier`: A function taking the transformed system and applying any
additional transformations, returning the modified system. The modified system
is passed to `linearization_function`.
All other keyword arguments are forwarded to linearization_function
.
ModelingToolkit.get_sensitivity
— Function get_sensitivity(sys, ap::AnalysisPoint; kwargs)
get_sensitivity(sys, ap_name::Symbol; kwargs)
Compute the sensitivity function in analysis point ap
. The sensitivity function is obtained by introducing an infinitesimal perturbation d
at the input of ap
, linearizing the system and computing the transfer function between d
and the output of ap
.
Arguments:
kwargs
: Are sent toModelingToolkit.linearize
See also get_comp_sensitivity
, get_looptransfer
.
ModelingToolkit.get_comp_sensitivity_function
— Functionget_comp_sensitivity_function(
sys::ModelingToolkit.AbstractSystem,
aps;
kwargs...
) -> Tuple{ModelingToolkit.LinearizationFunction{DI, AI, _A, P, _B, _C, J1, J2, J3, J4, IA, @NamedTuple{abstol::Float64, reltol::Float64, nlsolve_alg::NonlinearSolveFirstOrder.GeneralizedFirstOrderAlgorithm{Missing, NonlinearSolveFirstOrder.GenericTrustRegionScheme{NonlinearSolveFirstOrder.RadiusUpdateSchemes.__Simple, Rational{Int64}, Rational{Int64}, Rational{Int64}, Rational{Int64}, Rational{Int64}, Rational{Int64}, Rational{Int64}}, NonlinearSolveBase.Dogleg{NonlinearSolveBase.NewtonDescent{Nothing}, NonlinearSolveBase.SteepestDescent{Nothing}}, Nothing, Nothing, Nothing, Val{false}}}} where {DI<:AbstractVector{Int64}, AI<:AbstractVector{Int64}, _A, P<:ODEProblem, _B, _C, J1<:Union{Nothing, ModelingToolkit.PreparedJacobian{true, DifferentiationInterfaceForwardDiffExt.ForwardDiffTwoArgJacobianPrep{Nothing, C, Tuple{Nothing, Nothing}}, ModelingToolkit.var"#uff#957"{var"#1479#fun"}, _A, ADTypes.AutoForwardDiff{nothing, Nothing}} where {C<:(ForwardDiff.JacobianConfig{T, _A, _B, <:Tuple{Any, Any}} where {T<:(ForwardDiff.Tag{F} where F<:ModelingToolkit.var"#uff#957"), _A, _B}), var"#1479#fun", _A}}, J2<:Union{Nothing, ModelingToolkit.PreparedJacobian{true, DifferentiationInterfaceForwardDiffExt.ForwardDiffTwoArgJacobianPrep{Nothing, C, Tuple{Nothing, Nothing}}, _A, _B, ADTypes.AutoForwardDiff{nothing, Nothing}} where {C<:ForwardDiff.JacobianConfig, _A, _B}}, J3<:Union{Nothing, ModelingToolkit.PreparedJacobian{true, DifferentiationInterfaceForwardDiffExt.ForwardDiffTwoArgJacobianPrep{Nothing, C, Tuple{Nothing, Nothing, Nothing}}, ModelingToolkit.var"#pff#958"{var"#1480#fun", var"#1481#setter"}, _A, ADTypes.AutoForwardDiff{nothing, Nothing}} where {C<:(ForwardDiff.JacobianConfig{T, _A, _B, <:Tuple{Any, Any}} where {T<:(ForwardDiff.Tag{F} where F<:ModelingToolkit.var"#pff#958"), _A, _B}), var"#1480#fun", var"#1481#setter", _A}}, J4<:(ModelingToolkit.PreparedJacobian{true, DifferentiationInterfaceForwardDiffExt.ForwardDiffTwoArgJacobianPrep{Nothing, C, Tuple{Nothing, Nothing, Nothing}}, ModelingToolkit.var"#hpf#956"{fun, setter}, _A, ADTypes.AutoForwardDiff{nothing, Nothing}} where {C<:(ForwardDiff.JacobianConfig{T, _A, _B, <:Tuple{Any, Any}} where {T<:(ForwardDiff.Tag{F} where F<:ModelingToolkit.var"#hpf#956"), _A, _B}), fun, setter, _A}), IA<:Union{SciMLBase.NoInit, SciMLBase.OverrideInit{Nothing, Nothing, Nothing}}}, Any}
Return the complementary sensitivity function for the analysis point(s) aps
, and the modified system simplified with the appropriate inputs and outputs.
Keyword Arguments
- `loop_openings`: A list of analysis points whose connections should be removed and
the outputs set to the input as a part of the linear analysis.
- `system_modifier`: A function taking the transformed system and applying any
additional transformations, returning the modified system. The modified system
is passed to `linearization_function`.
All other keyword arguments are forwarded to linearization_function
.
ModelingToolkit.get_comp_sensitivity
— Functionget_comp_sensitivity(sys, ap::AnalysisPoint; kwargs)
get_comp_sensitivity(sys, ap_name::Symbol; kwargs)
Compute the complementary sensitivity function in analysis point ap
. The complementary sensitivity function is obtained by introducing an infinitesimal perturbation d
at the output of ap
, linearizing the system and computing the transfer function between d
and the input of ap
.
Arguments:
kwargs
: Are sent toModelingToolkit.linearize
See also get_sensitivity
, get_looptransfer
.
ModelingToolkit.get_looptransfer_function
— Functionget_looptransfer_function(
sys::ModelingToolkit.AbstractSystem,
aps;
kwargs...
) -> Tuple{ModelingToolkit.LinearizationFunction{DI, AI, _A, P, _B, _C, J1, J2, J3, J4, IA, @NamedTuple{abstol::Float64, reltol::Float64, nlsolve_alg::NonlinearSolveFirstOrder.GeneralizedFirstOrderAlgorithm{Missing, NonlinearSolveFirstOrder.GenericTrustRegionScheme{NonlinearSolveFirstOrder.RadiusUpdateSchemes.__Simple, Rational{Int64}, Rational{Int64}, Rational{Int64}, Rational{Int64}, Rational{Int64}, Rational{Int64}, Rational{Int64}}, NonlinearSolveBase.Dogleg{NonlinearSolveBase.NewtonDescent{Nothing}, NonlinearSolveBase.SteepestDescent{Nothing}}, Nothing, Nothing, Nothing, Val{false}}}} where {DI<:AbstractVector{Int64}, AI<:AbstractVector{Int64}, _A, P<:ODEProblem, _B, _C, J1<:Union{Nothing, ModelingToolkit.PreparedJacobian{true, DifferentiationInterfaceForwardDiffExt.ForwardDiffTwoArgJacobianPrep{Nothing, C, Tuple{Nothing, Nothing}}, ModelingToolkit.var"#uff#957"{var"#1479#fun"}, _A, ADTypes.AutoForwardDiff{nothing, Nothing}} where {C<:(ForwardDiff.JacobianConfig{T, _A, _B, <:Tuple{Any, Any}} where {T<:(ForwardDiff.Tag{F} where F<:ModelingToolkit.var"#uff#957"), _A, _B}), var"#1479#fun", _A}}, J2<:Union{Nothing, ModelingToolkit.PreparedJacobian{true, DifferentiationInterfaceForwardDiffExt.ForwardDiffTwoArgJacobianPrep{Nothing, C, Tuple{Nothing, Nothing}}, _A, _B, ADTypes.AutoForwardDiff{nothing, Nothing}} where {C<:ForwardDiff.JacobianConfig, _A, _B}}, J3<:Union{Nothing, ModelingToolkit.PreparedJacobian{true, DifferentiationInterfaceForwardDiffExt.ForwardDiffTwoArgJacobianPrep{Nothing, C, Tuple{Nothing, Nothing, Nothing}}, ModelingToolkit.var"#pff#958"{var"#1480#fun", var"#1481#setter"}, _A, ADTypes.AutoForwardDiff{nothing, Nothing}} where {C<:(ForwardDiff.JacobianConfig{T, _A, _B, <:Tuple{Any, Any}} where {T<:(ForwardDiff.Tag{F} where F<:ModelingToolkit.var"#pff#958"), _A, _B}), var"#1480#fun", var"#1481#setter", _A}}, J4<:(ModelingToolkit.PreparedJacobian{true, DifferentiationInterfaceForwardDiffExt.ForwardDiffTwoArgJacobianPrep{Nothing, C, Tuple{Nothing, Nothing, Nothing}}, ModelingToolkit.var"#hpf#956"{fun, setter}, _A, ADTypes.AutoForwardDiff{nothing, Nothing}} where {C<:(ForwardDiff.JacobianConfig{T, _A, _B, <:Tuple{Any, Any}} where {T<:(ForwardDiff.Tag{F} where F<:ModelingToolkit.var"#hpf#956"), _A, _B}), fun, setter, _A}), IA<:Union{SciMLBase.NoInit, SciMLBase.OverrideInit{Nothing, Nothing, Nothing}}}, Any}
Return the loop-transfer function for the analysis point(s) aps
, and the modified system simplified with the appropriate inputs and outputs.
Keyword Arguments
- `loop_openings`: A list of analysis points whose connections should be removed and
the outputs set to the input as a part of the linear analysis.
- `system_modifier`: A function taking the transformed system and applying any
additional transformations, returning the modified system. The modified system
is passed to `linearization_function`.
All other keyword arguments are forwarded to linearization_function
.
ModelingToolkit.get_looptransfer
— Functionget_looptransfer(sys, ap::AnalysisPoint; kwargs)
get_looptransfer(sys, ap_name::Symbol; kwargs)
Compute the (linearized) loop-transfer function in analysis point ap
, from ap.out
to ap.in
.
Feedback loops often use negative feedback, and the computed loop-transfer function will in this case have the negative feedback included. Standard analysis tools often assume a loop-transfer function without the negative gain built in, and the result of this function may thus need negation before use.
Arguments:
kwargs
: Are sent toModelingToolkit.linearize
See also get_sensitivity
, get_comp_sensitivity
, open_loop
.
ModelingToolkit.open_loop
— Functionopen_loop(
sys,
ap::Union{Symbol, AnalysisPoint};
system_modifier
) -> Tuple{Any, Tuple{Any, Any}}
Apply LoopTransferTransform
to the analysis point ap
and return the result of apply_transformation
.
Keyword Arguments
system_modifier
: a function which takes the modified system and returns a new system with any required further modifications performed.