deterministic
This module is defined such that operation on ode are all gathered in one place. Future extension of operations should be added here
-
class
pygom.model.deterministic.
DeterministicOde
(state=None, param=None, derived_param=None, transition=None, birth_death=None, ode=None) This contains the interface and operation built above the already defined set of ode
Parameters: - state: list
A list of states (string)
- param: list
A list of the parameters (string)
- derived_param: list
A list of the derived parameters (tuple of (string,string))
- transition: list
A list of transition (
Transition
)- birth_death: list
A list of birth or death process (
Transition
)- ode: list
A list of ode (
Transition
)
-
adjoint
(state, t, state_param, func=None) Compute the adjoint given the adjoint vector, time, state variable and the objective function. Note that this function is very restrictive in the sense that the (original) state variable changes through time but this assumes it is a constant, i.e. we assume that the original system is linear.
Parameters: - state: array like
The current value of lambda, where lambda’s are the Lagrangian multipliers of the differential equation.
- t: double
The current time.
- state_param: array like
The state vector that is (or maybe) required to evaluate the jacobian of the original system
- func: callable
This should take inputs similar to an ode, i.e. of the form func(y,t). If j(y,t) is the cost function, then func is a function that calculates \(\partial j \over \partial x\).
Returns: numpy.ndarray
output of the same length as the ode
Notes
The size of lambda should be the same as the state. The integral should be starting from T, the final time of the original system and is integrated backwards (for stability).
-
adjoint_T
(t, state, state_param, func=None) Same as
adjoint()
but with t as first parameter
-
adjoint_interpolate
(state, t, interpolant, func=None) Compute the adjoint given the adjoint vector, time, the functions which was used to interpolate the state variable
Parameters: - state: array like
The current value of lambda, where lambda’s are the Lagrangian multipliers of the differential equation.
- t: double
The current time.
- interpolant: list
list of interpolating functions of the state
- func: callable
This should take inputs similar to an ode, i.e. of the form func(y,t). If j(y,t) is the cost function, then func is a function that calculates \(\partial j \over \partial x\).
Returns: numpy.ndarray
output of the same length as the ode
-
adjoint_interpolate_T
(t, state, interpolant, objInput=None) Same as
adjoint_interpolate()
but with t as first parameter
-
adjoint_interpolate_jacobian
(state, t, interpolant, func=None) Compute the Jacobian of the adjoint given the adjoint vector, time, function of the interpolation on the state variables and the objective function. This is simply the same as the negative Jacobian of the ode transposed.
Parameters: - state: array like
The current value of lambda, where lambda’s are the Lagrangian multipliers of the differential equation.
- t: double
The current time.
- interpolant: list
list of interpolating functions of the state
- func: callable
This should take inputs similar to an ode, i.e. of the form func(y,t). If j(y,t) is the cost function, then func is a function that calculates \(\partial j \over \partial x\).
Returns: numpy.ndarray
output of is a two dimensional array of size [number of state x number of state]
See also
adjoint_jacobian()
Notes
Same as
adjoint_jacobian()
but takes a list of interpolating function instead of a single (vector) value
-
adjoint_interpolate_jacobian_T
(t, state, interpolant, func=None) Same as
adjoint_interpolate_jacobian()
but with t as first parameter
-
adjoint_jacobian
(state, t, state_param, func=None) Compute the jacobian of the adjoint given the adjoint vector, time, state variable and the objective function. This is simply the same as the negative jacobian of the ode transposed.
Parameters: - state: array like
The current value of lambda, where lambda’s are the Lagrangian multipliers of the differential equation.
- t: double
The current time.
- state_param: array like
The state vector that is (or maybe) required to evaluate the jacobian of the original system
- func: callable
This should take inputs similar to an ode, i.e. of the form func(y,t). If j(y,t) is the cost function, then func is a function that calculates \(\partial j \over \partial x\).
Returns: numpy.ndarray
output of is a two dimensional array of size [number of state x number of state]
See also
adjoint()
Notes
It takes the same number of argument as the adjoint for simplicity when integrating.
-
adjoint_jacobian_T
(t, state, state_param, func=None) Same as
adjoint_jacobian_T()
but with t being the first parameter
-
diff_jacobian
(state, t) Evaluate the differential of jacobian given state and time
Parameters: - state: array like
The current numerical value for the states which can be
numpy.ndarray
orlist
- t: double
The current time
Returns: numpy.ndarray
Matrix of dimension [number of state x number of state]
-
diff_jacobian_T
(t, state) Same as
diff_jacobian()
but with t as first parameter
-
eval_diff_jacobian
(parameters=None, time=None, state=None) Evaluate the differential of the jacobian given parameters, state and time. An extension of
diff_jacobian()
but now also include the parameters.Parameters: - parameters: list
see
parameters()
- time: double
The current time
- state: array list
The current numerical value for the states which can be
numpy.ndarray
orlist
Returns: numpy.matrix
ormpmath.matrix
Matrix of dimension [number of state x number of state]
See also
jacobian()
Notes
Name and order of state and time are also different.
-
eval_forwardforward
(FF, S, state, t) Evaluate a single f(x) of the forward-forward sensitivities
Parameters: - FF: array like
this is in fact a 3rd order Tensor, aka 3d array
- S: array like
sensitivities in matrix form
- state: array like
the current state
- t: numeric
time
Returns: numpy.ndarray
f(x) of size [number of state * (number of parameters * number of parameters)]
-
eval_grad
(parameters=None, time=None, state=None) Evaluate the gradient given parameters, state and time. An extension of
grad()
but now also include the parameters.Parameters: - parameters: list
see
parameters()
- time: double
The current time
- state: array list
The current numerical value for the states which can be
numpy.ndarray
orlist
Returns: numpy.matrix
ormpmath.matrix
Matrix of dimension [number of state x number of state]
See also
grad()
Notes
Name and order of state and time are also different.
-
eval_grad_jacobian
(parameters=None, time=None, state=None) Evaluate the jacobian of the gradient given parameters, state and time. An extension of
grad_jacobian()
but now also include the parameters.Parameters: - parameters: list
see
parameters()
- time: double
The current time
- state: array list
The current numerical value for the states which can be
numpy.ndarray
orlist
Returns: numpy.matrix
ormpmath.matrix
Matrix of dimension [number of state x number of state]
See also
grad_jacobian()
,get_grad_jacobian_eqn()
Notes
Name and order of state and time are also different.
-
eval_hessian
(parameters=None, time=None, state=None) Evaluate the hessian given parameters, state and time. An extension of
hessian()
but now also include the parameters.Parameters: - parameters: list
see
parameters()
- time: double
The current time
- state: array list
The current numerical value for the states which can be
numpy.ndarray
orlist
Returns: - list
list of dimension number of state, each with matrix [number of parameters x number of parameters] in
sympy.matricies.matricies
See also
grad()
,eval_grad()
-
eval_jacobian
(parameters=None, time=None, state=None) Evaluate the jacobian given parameters, state and time. An extension of
jacobian()
but now also include the parameters.Parameters: - parameters: list
see
parameters()
- time: double
The current time
- state: array list
The current numerical value for the states which can be
numpy.ndarray
orlist
Returns: numpy.matrix
ormpmath.matrix
Matrix of dimension [number of state x number of state]
See also
jacobian()
Notes
Name and order of state and time are also different.
-
eval_ode
(parameters=None, time=None, state=None) Evaluate the ode given time, state and parameters. An extension of
ode()
but now also include the parameters.Parameters: - parameters: list
see
parameters()
- time: numeric
The current time
- state: array like
The current numerical value for the states which can be
numpy.ndarray
orlist
Returns: numpy.matrix
ormpmath.matrix
output of the same length as the ode.
See also
ode()
Notes
There are differences between the output of this function and
ode()
. Name and order of state and time are also different.
-
eval_sens_jacobian_state
(time=None, state=None, sens=None) Evaluate the jacobian of the sensitivities w.r.t the states given parameters, state and time. An extension of
sens_jacobian_state()
but now also include the parameters.Parameters: - parameters: list
see
parameters()
- time: double
The current time
- state: array list
The current numerical value for the states which can be
numpy.ndarray
orlist
Returns: numpy.matrix
ormpmath.matrix
Matrix of dimension [number of state x number of state]
See also
sens_jacobian_state()
Notes
Name and order of state and time are also different.
-
eval_sensitivity
(S, t, state, by_state=False) Evaluate the sensitivity given state and time
Parameters: - S: array like
Which should be
numpy.ndarray
. The starting sensitivity of size [number of state x number of parameters]. Which are normally zero or one, depending on whether the initial conditions are also variables.- t: double
The current time
- state: array like
The current numerical value for the states which can be
numpy.ndarray
orlist
- by_state: bool
how we want the output to be arranged. Default is True so that we have a block diagonal structure
Returns: See also
sensitivity()
Notes
It is different to
eval_ode()
andeval_jacobian()
in that the extra input argument is not a parameter
-
eval_sensitivityIV
(S, IV, t, state) Evaluate the sensitivity with initial values given state and time
Parameters: - S: array like
Which should be
numpy.ndarray
. The starting sensitivity of size [number of state x number of parameters]. Which are normally zero or one, depending on whether the initial conditions are also variables.- IV: array like
sensitivities for the initial values
- t: double
The current time
- state: array like
The current numerical value for the states which can be
numpy.ndarray
orlist
Returns: numpy.ndarray
\(f(s(x,\theta))\) and \(f(s(x_{0}))\)
See also
sensitivityIV()
Notes
It is different to
eval_ode()
andeval_jacobian()
in that the extra input argument is not a parameter.
-
forwardforward
(ff, t, state, s) Evaluate a single \(f(x)\) of the forward-forward sensitivities
Parameters: - ff: array like
the forward-forward sensitivities in vector form
- t: numeric
time
- state: array like
the current state
- s: array like
forward sensitivities in vector form
Returns: numpy.ndarray
\(f(x)\) of size [number of state * (number of parameters * number of parameters)]
-
forwardforward_T
(t, ff, s, state) Same as
forwardforward()
but with t as the first parameter
-
get_diff_jacobian_eqn
() Returns the jacobian differentiate w.r.t. states in algebraic form
Returns: - list
list of size (num of state,) each with
sympy.matrices.matrices
of dimension [number of state x number of state]
-
get_grad_eqn
() Return the gradient of the ode in algebraic form
Returns: sympy.matrices.matrices
A matrix of dimension [number of state x number of parameters]
-
get_grad_jacobian_eqn
() Return the jacobian of the gradient in algebraic form
Returns: sympy.matrices.matrices
A matrix of dimension [number of state * number of parameters x number of state]
See also
get_grad_eqn()
-
get_hessian_eqn
() Return the Hessian of the ode in algebraic form
Returns: - list
list of dimension number of state, each with matrix [number of parameters x number of parameters] in
sympy.matricies.matricies
Notes
We deliberately return a list instead of a 3d array of a tensor to avoid confusion
-
get_jacobian_eqn
() Returns the jacobian in algebraic form
Returns: sympy.matrices.matrices
A matrix of dimension [number of state x number of state]
-
get_ode_eqn
(param_sub=False) Find the algebraic equations of the ode system.
Returns: sympy.matrices.matrices
ode in matrix form
-
get_transition_graph
(file_name=None, show=True) Returns the transition graph using graphviz
Parameters: - file_name: str, optional
name of the output file, defaults to None
- show: bool, optional
If the graph should be plotted, defaults to True
Returns: graphviz.Digraph
-
grad
(state, time) Evaluate the gradient given state and time
Parameters: - state: array like
The current numerical value for the states which can be
numpy.ndarray
orlist
- t: numeric
The current time
Returns: numpy.ndarray
Matrix of dimension [number of state x number of parameters]
-
grad_T
(t, state) Same as
grad_T()
but with t as first parameter
-
grad_jacobian
(state, time) Evaluate the Jacobian of the gradient given state and time
Parameters: - state: array like
The current numerical value for the states which can be
numpy.ndarray
orlist
- t: numeric
The current time
Returns: numpy.ndarray
Matrix of dimension [number of state x number of parameters]
See also
grad()
-
grad_jacobianT
(t, state) Same as
grad_jacobian()
but with t as first parameter
-
hessian
(state, time) Evaluate the hessian given state and time
Parameters: - state: array like
The current numerical value for the states which can be
numpy.ndarray
orlist
- t: double
The current time
Returns: - list
list of dimension number of state, each with matrix [number of parameters x number of parameters] in
sympy.matricies.matricies
-
initial_state
Return the initial state values
-
initial_time
Return the initial time
-
initial_values
Returns the initial values, both time and state as a tuple (x0, t0)
-
integrate
(t, full_output=False) Integrate over a range of t when t is an array and a output at time t
Parameters: - t: array like
the range of time points which we want to see the result of
- full_output: bool
if we want additional information
-
integrate2
(t, full_output=False, method=None) Integrate over a range of t when t is an array and a output at time t. Select a suitable method to integrate when method is None.
Parameters: - t: array like
the range of time points which we want to see the result of
- full_output: bool
if we want additional information
- method: str, optional
the integration method. All those available in
ode
are allowed with ‘vode’ and ‘ivode’ representing the non-stiff and stiff version respectively. Defaults to None, which tries to choose the integration method via eigenvalue analysis (only one) using the initial conditions
-
is_stiff
(state=None, t=None) Test on the eigenvalues of the jacobian. We classify the problem as stiff if any of the eigenvalues are positive
Parameters: - state: array like
The current numerical value for the states which can be
numpy.ndarray
orlist
- t: double
The current time
Returns: numpy.ndarray
eigenvalues of the system given input
-
jacobian
(state, t) Evaluate the jacobian given state and time
Parameters: - state: array like
The current numerical value for the states which can be
numpy.ndarray
orlist
- t: double
The current time
Returns: numpy.ndarray
Matrix of dimension [number of state x number of state]
-
jacobian_T
(t, state) Same as
jacobian()
but with t as first parameter
-
jacobian_eigenvalue
(state=None, t=None) Find out the eigenvalues of the jacobian given state and time. If None is given, the initial values are used.
Parameters: - state: array like
The current numerical value for the states which can be
numpy.ndarray
orlist
- t: double
The current time
Returns: - bool
True if any eigenvalue is positive
-
linear_ode
() To check whether the input ode is linear
Returns: - bool
True if it is linear, False otherwise
-
ode
(state, t) Evaluate the ode given state and time
Parameters: - state: array like
The current numerical value for the states which can be
numpy.ndarray
orlist
- t: double
The current time
Returns: numpy.ndarray
output of the same length as the ode
-
ode_T
(t, state) Same as
ode()
but with t as the first parameter
-
ode_and_forwardforward
(state_param, t) Evaluate a single f(x) of the ode and the forward-forward sensitivities
Parameters: - state_param: array like
state and forward-forward sensitivities in vector form
- t: numeric
time
Returns: numpy.ndarray
same size as the state_param input
-
ode_and_forwardforward_T
(t, state_param) Same as
odeAndForwardForward()
but with time as the first input
-
ode_and_forwardforward_jacobian
(state_param, t) Return the jacobian after evaluation given the input of the state and the forward forward sensitivities
Parameters: - state_param: array like
state and forward-forward sensitivities in vector form
- t: numeric
time
Returns: numpy.ndarray
size of (a,a) where a is the length of the state_param input
-
ode_and_forwardforward_jacobian_T
(t, state_param) Same as
ode_and_forwardforward_jacobian()
but with t being the first parameters
-
ode_and_sensitivity
(state_param, t, by_state=False) Evaluate the sensitivity given state and time
Parameters: - state_param: array like
The current numerical value for the states as well as the sensitivities values all in one. We assume that the state values comes first.
- t: double
The current time
- by_state: bool
Whether the output vector should be arranged by state or by parameters. If False, then it means that the vector of output is arranged according to looping i,j from Sensitivity_{i,j} with i being the state and j the param. This is the preferred way because it leds to a block diagonal Jacobian
Returns: list
concatenation of 2 element. First contains the ode, second the sensitivity. Both are of type
numpy.ndarray
See also
sensitivity()
,ode()
-
ode_and_sensitivityIV
(state_param, t) Evaluate the sensitivity given state and time
Parameters: - state_param: array like
The current numerical value for the states as well as the sensitivities values all in one. We assume that the state values comes first.
- t: double
The current time
Returns: list
concatenation of 3 element. First contains the ode, second the sensitivity, then the sensitivity of the initial value. All of them are of type
numpy.ndarray
See also
sensitivity()
,ode()
-
ode_and_sensitivityIV_T
(t, state_param) Same as
ode_and_sensitivityIV()
but with t as first parameter
-
ode_and_sensitivityIV_jacobian
(state_param, t) Evaluate the sensitivity given state and time. Output a block diagonal sparse matrix as default.
Parameters: - state_param: array like
The current numerical value for the states as well as the sensitivities values all in one. We assume that the state values comes first.
- t: double
The current time
- byState: bool
How the output is arranged, according to the vector of output. It can be in terms of state or parameters, where by state means that the jacobian is a block diagonal matrix.
Returns: numpy.ndarray
output of a square matrix of size: number of ode + 1 times number of parameters
See also
ode_and_sensitivity()
-
ode_and_sensitivityIV_jacobian_T
(t, state_param) Same as
ode_and_sensitivityIV_jacobian()
but with t as first parameter
-
ode_and_sensitivity_T
(t, state_param, by_state=False) Same as
ode_and_sensitivity()
but with t as first parameter
-
ode_and_sensitivity_jacobian
(state_param, t, by_state=False) Evaluate the sensitivity given state and time. Output a block diagonal sparse matrix as default.
Parameters: - state_param: array like
The current numerical value for the states as well as the sensitivities values all in one. We assume that the state values comes first.
- t: double
The current time
- by_state: bool
How the output is arranged, according to the vector of output. It can be in terms of state or parameters, where by state means that the jacobian is a block diagonal matrix.
Returns: numpy.ndarray
output of a square matrix of size: number of ode + 1 times number of parameters
See also
ode_and_sensitivity()
-
ode_and_sensitivity_jacobian_T
(t, state_param, by_state=False) Same as
ode_and_sensitivity_jacobian()
but with t as first parameter
-
plot
() Plot the results of the integration
Notes
If we have 3 states or more, it will always be arrange such that it has 3 columns. Uses the operation from
odeutils
-
print_ode
(latex_output=False) Prints the ode in symbolic form onto the screen/console in actual symbols rather than the word of the symbol.
Parameters: - latex_output: bool, optional
Defaults to false which prints the equation in terms of symbols, if set to yes then the formula in terms of latex equations will be printed onto the screen.
-
sens_jacobian_state
(state_param, t) Evaluate the jacobian of the sensitivity w.r.t. the state given state and time
Parameters: - state_param: array like
The current numerical value for the states as well as the sensitivities, which can be
numpy.ndarray
orlist
- t: double
The current time
Returns: numpy.ndarray
Matrix of dimension [number of state * number of parameters x number of state]
-
sens_jacobian_state_T
(t, state) Same as
sens_jacobian_state_T()
but with t as first parameter
-
sensitivity
(sens, t, state, by_state=False) Evaluate the sensitivity given state and time. The default is to output the values by parameters, i.e. \(s_{i},\ldots,s_{i+n}\) are partial derivatives w.r.t. the states for \(i \in {1,1+p,1+2p,1+3p, \ldots, 1+(n-1)p}\). This is to take advantage of the fact that we have a block diagonal jacobian that was already evaluated
Parameters: - sens: array like
The starting sensitivity of size [number of state x number of parameters]. Which are normally zero or one, depending on whether the initial conditions are also variables.
- t: double
The current time
- state: array like
The current numerical value for the states which can be
numpy.ndarray
orlist
- by_state: bool
how we want the output to be arranged. Default is True so that we have a block diagonal structure
Returns:
-
sensitivityIV
(sensIV, t, state) Evaluate the sensitivity which include the initial values as our parameters given state and time. The default is to output the values by parameters, i.e. \(s_{i},\ldots,s_{i+n}\) are partial derivatives w.r.t. the states for \(i \in {1,1+p,1+2p,1+3p, \ldots, 1+(n-1)p}\). This is to take advantage of the fact that we have a block diagonal Jacobian that was already evaluated.
Parameters: - sensIV: array like
The starting sensitivity of size [number of state x number of parameters] + [number of state x number of state] for the initial condition. The latter is an identity matrix at time zero.
- t: double
The current time
- state: array like
The current numerical value for the states which can be
numpy.ndarray
orlist
Returns: numpy.ndarray
output of the same length as the ode
-
sensitivityIV_T
(t, sensIV, state) Same as
sensitivityIV()
but with t as first parameter
-
sensitivity_T
(t, sens, state, by_state=False) Same as
sensitivity()
but with t as first parameter