The API for building models with FunctionalModels. Includes basic types, models, and functions.
t
FunctionalModels.t
— ConstantIndependent variable
D
FunctionalModels.D
— FunctionDifferential(t)
der
FunctionalModels.der
— FunctionDifferential(t)
Unknown
FunctionalModels.Unknown
— FunctionUnknown(value = NaN; name = :u)
Unknown
is a helper to create a variable with a default value. The default value determines the type and shape of the result. It also adds metadata to variables so that variable names don't clash. The viewable variable name is based on a gensym
. name
is stored as metadata, and when equations are flattened with system
, variables are renamed to include subsystem names and variable base name.
For example, Unknown(name = :v)
may show as var"##v#1057"(t)
, but after flattening, it will show as something like ss₊c1₊v(t)
(ss
and c1
are subsystems).
Parameter
FunctionalModels.Parameter
— FunctionParameter(value = 0.0; name)
Parameter
is a helper to create a parameter with default values. The default value determines the type and shape of the result. It also adds metadata to variables so that variable names don't clash. The viewable variable name is based on a gensym
. name
is stored as metadata, and when equations are flattened with system
, variables are renamed to include subsystem names and the variable base name.
For example, Parameter(name = :R)
may show as var"##R#1057"
, but after flattening, it will show as something like ss₊r1₊R
(ss
and r1
are subsystems).
RefBranch
FunctionalModels.RefBranch
— TypeA special ModelType to specify branch flows into nodes. When the model is flattened, equations are created to zero out branch flows into nodes.
See also Branch.
RefBranch(n, i)
Arguments
n
: the reference node.i
: the flow variable that goes with this node.
References
This nodal description is based on work by David Broman. See the following:
- http://www.eecs.berkeley.edu/Pubs/TechRpts/2012/EECS-2012-173.pdf
- http://www.bromans.com/software/mkl/mkl-source-1.0.0.zip
- https://github.com/david-broman/modelyze
Modelyze has both RefBranch
and Branch
.
Examples
Here is an example of RefBranch used in the definition of a HeatCapacitor in the standard library. hp
is the reference node (a HeatPort aka Temperature), and Q_flow
is the flow variable.
function HeatCapacitor(hp::HeatPort, C::Signal)
Q_flow = HeatFlow(compatible_values(hp))
[
RefBranch(hp, Q_flow)
der(hp) ~ Q_flow ./ C
]
end
Here is the definition of SignalCurrent from the standard library a model that injects current (a flow variable) between two nodes:
function SignalCurrent(n1::ElectricalNode, n2::ElectricalNode, I::Signal)
[
RefBranch(n1, I)
RefBranch(n2, -I)
]
end
Branch
FunctionalModels.Branch
— FunctionA helper Model to connect a branch between two different nodes and specify potential between nodes and the flow between nodes.
See also RefBranch.
Branch(n1, n2, v, i)
Arguments
n1
: the positive reference node.n2
: the negative reference node.v
: the potential variable between nodes.i
: the flow variable between nodes.
Returns
- The model, consisting of a RefBranch entry for each node and an equation assigning
v
ton1 - n2
.
References
This nodal description is based on work by David Broman. See the following:
- http://www.eecs.berkeley.edu/Pubs/TechRpts/2012/EECS-2012-173.pdf
- http://www.bromans.com/software/mkl/mkl-source-1.0.0.zip
- https://github.com/david-broman/modelyze
Examples
Here is the definition of an electrical resistor in the standard library:
function Resistor(n1::ElectricalNode, n2::ElectricalNode, R::Signal)
i = Current(compatible_values(n1, n2))
v = Voltage(value(n1) - value(n2))
[
Branch(n1, n2, v, i)
v ~ R .* i
]
end
system
FunctionalModels.system
— Functionsystem
is the main elaboration/flattening function that returns an ODESystem
.
system(a)
Arguments
a
: the input model containing a nested vector of equations
Optional/Keyword Arguments
simplify = true
: whetherstructural_simplify
is used to simplify results
Returns
::ODESystem
: the flattened model
default_value
FunctionalModels.default_value
— FunctionThe default or starting value of a variable.
default_value(x)
Arguments
x
: the reference variable or numeric value.
compatible_values
FunctionalModels.compatible_values
— FunctionA helper functions to return the base value from a variable to use when creating other variables. It is especially useful for taking two model arguments and creating a new variable compatible with both arguments.
compatible_values(x,y)
compatible_values(x)
It's still somewhat broken but works for basic cases. No type promotion is currently done.
Arguments
x
,y
: objects or variables
Returns
The returned object has zeros of type and length common to both x
and y
.
Examples
a = Unknown(45.0 + 10im)
x = Unknown(compatible_values(a)) # Initialized to 0.0 + 0.0im.
a = Unknown()
b = Unknown([1., 0.])
y = Unknown(compatible_values(a,b)) # Initialized to [0.0, 0.0].
compatible_shape
FunctionalModels.compatible_shape
— FunctionA helper functions to return the base NaN value from a variable to use when creating other variables. It is especially useful for taking two model arguments and creating a new variable compatible with both arguments. This differs fron compatible_values
in that it returns values filled with NaNs to indicate a variable without a default value.
compatible_shape(x,y)
compatible_shape(x)
It's still somewhat broken but works for basic cases. No type promotion is currently done.
Arguments
x
,y
: objects or variables
Returns
The returned object has NaNs of type and length common to both x
and y
.
Examples
a = Unknown(45.0 + 10im)
x = Unknown(compatible_shape(a)) # Initialized to NaN + NaNim.
a = Unknown()
b = Unknown([1., 0.])
y = Unknown(compatible_shape(a,b)) # Initialized to [NaN, NaN].