API

direct_collocation.py

class opty.direct_collocation.ConstraintCollocator(equations_of_motion, state_symbols, num_collocation_nodes, node_time_interval, known_parameter_map={}, known_trajectory_map={}, instance_constraints=None, time_symbol=None, tmp_dir=None, integration_method='backward euler', parallel=False)[source]

This class is responsible for generating the constraint function and the sparse Jacobian of the constraint function using direct collocation methods for a non-linear programming problem where the essential constraints are defined from the equations of motion of the system.

Notes

  • N : number of collocation nodes

  • n : number of states

  • m : number of input trajectories

  • p : number of parameters

  • q : number of unknown input trajectories

  • r : number of unknown parameters

  • o : number of instance constraints

  • nN + qN + r : number of free variables

  • n(N - 1) + o : number of constraints

Instantiates a ConstraintCollocator object.

Parameters:
  • equations_of_motion (sympy.Matrix, shape(n, 1)) – A column matrix of SymPy expressions defining the right hand side of the equations of motion when the left hand side is zero, e.g. 0 = x’(t) - f(x(t), u(t), p) or 0 = f(x’(t), x(t), u(t), p). These should be in first order form but not necessairly explicit.

  • state_symbols (iterable) – An iterable containing all of the SymPy functions of time which represent the states in the equations of motion.

  • num_collocation_nodes (integer) – The number of collocation nodes, N. All known trajectory arrays should be of this length.

  • node_time_interval (float) – The time interval between collocation nodes.

  • known_parameter_map (dictionary, optional) – A dictionary that maps the SymPy symbols representing the known constant parameters to floats. Any parameters in the equations of motion not provided in this dictionary will become free optimization variables.

  • known_trajectory_map (dictionary, optional) – A dictionary that maps the non-state SymPy functions of time to ndarrays of floats of shape(N,). Any time varying parameters in the equations of motion not provided in this dictionary will become free trajectories optimization variables.

  • instance_constraints (iterable of SymPy expressions, optional) – These expressions are for constraints on the states at specific time points. They can be expressions with any state instance and any of the known parameters found in the equations of motion. All states should be evaluated at a specific instant of time. For example, the constraint x(0) = 5.0 would be specified as x(0) - 5.0 and the constraint x(0) = x(5.0) would be specified as x(0) - x(5.0). Unknown parameters and time varying parameters other than the states are currently not supported.

  • time_symbol (SymPy Symbol, optional) – The symbol representating time in the equations of motion. If not given, it is assumed to be the default stored in dynamicsymbols._t.

  • tmp_dir (string, optional) – If you want to see the generated Cython and C code for the constraint and constraint Jacobian evaluations, pass in a path to a directory here.

  • integration_method (string, optional) – The integration method to use, either backward euler or midpoint.

  • parallel (boolean, optional) – If true and openmp is installed, constraints and the Jacobian of the constraints will be executed across multiple threads. This is only useful when the equations of motion have an extremely large number of operations.

generate_constraint_function()[source]

Returns a function which evaluates the constraints given the array of free optimization variables.

generate_jacobian_function()[source]

Returns a function which evaluates the Jacobian of the constraints given the array of free optimization variables.

jacobian_indices()[source]

Returns the row and column indices for the non-zero values in the constraint Jacobian.

Returns:

  • jac_row_idxs (ndarray, shape(2*n + q + r,)) – The row indices for the non-zero values in the Jacobian.

  • jac_col_idxs (ndarray, shape(n,)) – The column indices for the non-zero values in the Jacobian.

class opty.direct_collocation.Problem(obj, obj_grad, *args, **kwargs)[source]

This class allows the user to instantiate a problem object with the essential data required to solve a direct collocation optinal control or parameter identification problem.

Parameters:
  • obj (function) – Returns the value of the objective function given the free vector.

  • obj_grad (function) – Returns the gradient of the objective function given the free vector.

  • equations_of_motion (sympy.Matrix, shape(n, 1)) – A column matrix of SymPy expressions defining the right hand side of the equations of motion when the left hand side is zero, e.g. 0 = x’(t) - f(x(t), u(t), p) or 0 = f(x’(t), x(t), u(t), p). These should be in first order form but not necessairly explicit.

  • state_symbols (iterable) – An iterable containing all of the SymPy functions of time which represent the states in the equations of motion.

  • num_collocation_nodes (integer) – The number of collocation nodes, N. All known trajectory arrays should be of this length.

  • node_time_interval (float) – The time interval between collocation nodes.

  • known_parameter_map (dictionary, optional) – A dictionary that maps the SymPy symbols representing the known constant parameters to floats. Any parameters in the equations of motion not provided in this dictionary will become free optimization variables.

  • known_trajectory_map (dictionary, optional) – A dictionary that maps the non-state SymPy functions of time to ndarrays of floats of shape(N,). Any time varying parameters in the equations of motion not provided in this dictionary will become free trajectories optimization variables.

  • instance_constraints (iterable of SymPy expressions, optional) – These expressions are for constraints on the states at specific time points. They can be expressions with any state instance and any of the known parameters found in the equations of motion. All states should be evaluated at a specific instant of time. For example, the constraint x(0) = 5.0 would be specified as x(0) - 5.0 and the constraint x(0) = x(5.0) would be specified as x(0) - x(5.0). Unknown parameters and time varying parameters other than the states are currently not supported.

  • time_symbol (SymPy Symbol, optional) – The symbol representating time in the equations of motion. If not given, it is assumed to be the default stored in dynamicsymbols._t.

  • tmp_dir (string, optional) – If you want to see the generated Cython and C code for the constraint and constraint Jacobian evaluations, pass in a path to a directory here.

  • integration_method (string, optional) – The integration method to use, either backward euler or midpoint.

  • parallel (boolean, optional) – If true and openmp is installed, constraints and the Jacobian of the constraints will be executed across multiple threads. This is only useful when the equations of motion have an extremely large number of operations.

  • bounds (dictionary, optional) – This dictionary should contain a mapping from any of the symbolic states, unknown trajectories, or unknown parameters to a 2-tuple of floats, the first being the lower bound and the second the upper bound for that free variable, e.g. {x(t): (-1.0, 5.0)}.

addOption()

Add a keyword/value option pair to the problem.

Deprecated since version 1.0.0: addOption() will be removed in CyIpopt 1.1.0, it is replaced by add_option() because the latter complies with PEP8.

add_option()

Add a keyword/value option pair to the problem.

See the Ipopt documentaion for details on available options.

Parameters:
  • keyword (str) – Option name.

  • val (str, int or float) – Value of the option. The type of val should match the option definition as described in the Ipopt documentation.

close()

Deallocate memory resources used by the Ipopt package.

Called implicitly by the Problem class destructor.

constraints(free)[source]

Returns the value of the constraint functions given a solution to the problem.

Parameters:

free (ndarray, (n*N + q*N + r, )) – A solution to the optimization problem in the canonical form.

Returns:

constraints_val – The value of the constraint function.

Return type:

ndarray, shape(n*(N - 1) + o)

Notes

  • N : number of collocation nodes

  • n : number of unknown state trajectories

  • q : number of unknown input trajectories

  • r : number of unknown parameters

  • o : number of instance constraints

get_current_iterate()

Return the current iterate vectors during an Ipopt solve

The iterate contains vectors for primal variables, bound multipliers, constraint function values, and constraint multipliers. Here, the constraints are treated as a single function rather than separating equality and inequality constraints. This method can only be called during an intermediate callback.

Only supports Ipopt >=3.14.0

Parameters:

scaled (Bool) – Whether the scaled iterate vectors should be returned

Returns:

A dict containing the iterate vector with keys "x", "mult_x_L", "mult_x_U", "g", and "mult_g". If iterate vectors cannot be obtained, None is returned.

Return type:

dict or None

get_current_violations()

Return the current violation vectors during an Ipopt solve

Violations returned are primal variable bound violations, bound complementarities, the gradient of the Lagrangian, constraint violation, and constraint complementarity. Here, the constraints are treated as a single function rather than separating equality and inequality constraints. This method can only be called during an intermediate callback.

Only supports Ipopt >=3.14.0

Parameters:

scaled (Bool) – Whether to scale the returned violations

Returns:

A dict containing the violation vector with keys "x_L_violation", "x_U_violation", "compl_x_L", "compl_x_U", "grad_lag_x", "g_violation", and "compl_g". If violation vectors cannot be obtained, None is returned.

Return type:

dict or None

gradient(free)[source]

Returns the value of the gradient of the objective function given a solution to the problem.

Parameters:

free (ndarray, (n*N + q*N + r, )) – A solution to the optimization problem in the canonical form.

Returns:

gradient_val – The value of the gradient of the objective function.

Return type:

ndarray, shape(n*N + q*N + r, 1)

Notes

  • N : number of collocation nodes

  • n : number of unknown state trajectories

  • q : number of unknown input trajectories

  • r : number of unknown parameters

intermediate(*args)[source]

This method is called at every optimization iteration. Not for pubic use.

jacobian(free)[source]

Returns the non-zero values of the Jacobian of the constraint function.

Return type:

jac_vals = ndarray, shape()

jacobianstructure()[source]

Returns the sparsity structur of the Jacobian of the constraint function.

Returns:

  • jac_row_idxs (ndarray, shape(2*n + q + r,)) – The row indices for the non-zero values in the Jacobian.

  • jac_col_idxs (ndarray, shape(n,)) – The column indices for the non-zero values in the Jacobian.

objective(free)[source]

Returns the value of the objective function given a solution to the problem.

Parameters:

free (ndarray, (n*N + q*N + r, )) – A solution to the optimization problem in the canonical form.

Returns:

obj_val – The value of the objective function.

Return type:

float

Notes

  • N : number of collocation nodes

  • n : number of unknown state trajectories

  • q : number of unknown input trajectories

  • r : number of unknown parameters

plot_constraint_violations(vector)[source]

Returns an axis with the state constraint violations plotted versus node number and the instance constraints as a bar graph.

Parameters:

vector (ndarray, (n*N + q*N + r, )) – The initial guess, solution, or any other vector that is in the canonical form.

Returns:

axes – A matplotlib axes with the constraint violations plotted.

Return type:

ndarray of AxesSubplot

Notes

  • N : number of collocation nodes

  • n : number of unknown state trajectories

  • q : number of unknown input trajectories

  • r : number of unknown parameters

plot_objective_value()[source]

Returns an axis with the objective value plotted versus the optimization iteration. solve() must be run first.

plot_trajectories(vector, axes=None)[source]

Returns the axes for two plots. The first plot displays the state trajectories versuse time and the second plot displays the input trjaectories versus time.

Parameters:
  • vector (ndarray, (n*N + q*N + r, )) – The initial guess, solution, or any other vector that is in the canonical form.

  • axes (ndarray of AxesSubplot, shape(n + m, )) – An array of matplotlib axes to plot to.

Returns:

axes – A matplotlib axes with the state and input trajectories plotted.

Return type:

ndarray of AxesSubplot

Notes

  • N : number of collocation nodes

  • n : number of unknown state trajectories

  • m : number of input trajectories

  • q : number of unknown input trajectories

  • r : number of unknown parameters

setProblemScaling()

Optional function for setting scaling parameters for the problem.

Deprecated since version 1.0.0: setProblemScaling() will be removed in CyIpopt 1.1.0, it is replaced by set_problem_scaling() because the latter complies with PEP8.

set_problem_scaling()

Optional function for setting scaling parameters for the problem.

To use the scaling parameters set the option nlp_scaling_method to user-scaling.

Parameters:
  • obj_scaling (float) – Determines, how Ipopt should internally scale the objective function. For example, if this number is chosen to be 10, then Ipopt solves internally an optimization problem that has 10 times the value of the original objective. In particular, if this value is negative, then Ipopt will maximize the objective function instead of minimizing it.

  • x_scaling (array-like, shape(n, )) – The scaling factors for the variables. If None, no scaling is done.

  • g_scaling (array-like, shape(m, )) – The scaling factors for the constrains. If None, no scaling is done.

solve()

Returns the optimal solution and an info dictionary.

Solves the posed optimization problem starting at point x.

Parameters:

x (array-like, shape(n, )) – Initial guess.

Returns:

  • x (array, shape(n, )) – Optimal solution.

  • info (dictionary) –

    x: ndarray, shape(n, )

    optimal solution

    g: ndarray, shape(m, )

    constraints at the optimal solution

    obj_val: float

    objective value at optimal solution

    mult_g: ndarray, shape(m, )

    final values of the constraint multipliers

    mult_x_L: ndarray, shape(n, )

    bound multipliers at the solution

    mult_x_U: ndarray, shape(n, )

    bound multipliers at the solution

    status: integer

    gives the status of the algorithm

    status_msg: string

    gives the status of the algorithm as a message

parameter_identification.py

opty.parameter_identification.objective_function(free, num_dis_points, num_states, dis_period, time_measured, y_measured)[source]

Returns the norm of the difference in the measured and simulated output.

Parameters:
  • free (ndarray, shape(n*N + q*N + r,)) – The flattened state array with n states at N time points, q input trajectories and N time points, and the r free model constants.

  • num_dis_points (integer) – The number of model discretization points.

  • num_states (integer) – The number of system states.

  • dis_period (float) – The discretization time interval.

  • y_measured (ndarray, shape(M, o)) – The measured trajectories of the o output variables at each sampled time instance.

Returns:

cost – The cost value.

Return type:

float

Notes

This assumes that the states are ordered:

[coord1, …, coordn, speed1, …, speedn]

y_measured is interpolated at the discretization time points and compared to the model output at the discretization time points.

opty.parameter_identification.objective_function_gradient(free, num_dis_points, num_states, dis_period, time_measured, y_measured)[source]

Returns the gradient of the objective function with respect to the free parameters.

Parameters:
  • free (ndarray, shape(n*N + q*N + r,)) – The flattened state array with n states at N time points, q input trajectories and N time points, and the r free model constants.

  • num_dis_points (integer) – The number of model discretization points.

  • num_states (integer) – The number of system states.

  • dis_period (float) – The discretization time interval.

  • y_measured (ndarray, shape(M, o)) – The measured trajectories of the o output variables at each sampled time instance.

Returns:

gradient – The gradient of the cost function with respect to the free parameters.

Return type:

ndarray, shape(n*N + q*N + r,)

Warning

This is currently only valid if the model outputs (and measurements) are simply the states. The chain rule will be needed if the function output_equations() is more than a simple selection.

opty.parameter_identification.output_equations(x)[source]

Returns the outputs of the system. For now just the an array of the generalized coordinates.

Parameters:

x (ndarray, shape(N, n)) – The trajectories of the system states.

Returns:

y – The trajectories of the generalized coordinates.

Return type:

ndarray, shape(N, o)

Notes

The order of the states is assumed to be:

[coord_1, …, coord_{n/2}, speed_1, …, speed_{n/2}]

[q_1, …, q_{n/2}, u_1, …., u_{n/2}]

As this is what generate_ode_function creates.

utils.py

opty.utils.controllable(a, b)[source]

Returns true if the system is controllable and false if not.

Parameters:
  • a (array_like, shape(n,n)) – The state matrix.

  • b (array_like, shape(n,r)) – The input matrix.

Returns:

controllable

Return type:

boolean

opty.utils.f_minus_ma(mass_matrix, forcing_vector, states)[source]

Returns Fr + Fr* from the mass_matrix and forcing vector.

opty.utils.openmp_installed()[source]

Returns true if openmp is installed, false if not.

Modified from: https://stackoverflow.com/questions/16549893/programatically-testing-for-openmp-support-from-a-python-setup-script

opty.utils.parse_free(free, n, q, N)[source]

Parses the free parameters vector and returns it’s components.

Parameters:
  • free (ndarray, shape(n*N + q*N + r)) – The free parameters of the system.

  • n (integer) – The number of states.

  • q (integer) – The number of free specified inputs.

  • N (integer) – The number of time steps.

Returns:

  • states (ndarray, shape(n, N)) – The array of n states through N time steps.

  • specified_values (ndarray, shape(r, N) or shape(N,), or None) – The array of r specified inputs through N time steps.

  • constant_values (ndarray, shape(q,)) – The array of q constants.

opty.utils.state_derivatives(states)[source]

Returns functions of time which represent the time derivatives of the states.

opty.utils.substitute_matrix(matrix, row_idxs, col_idxs, sub_matrix)[source]

Returns the matrix with the values given by the row and column indices with those in the sub-matrix.

Parameters:
  • matrix (ndarray, shape(n, m)) – A matrix (i.e. 2D array).

  • row_idxs (array_like, shape(p<=n,)) – The row indices which designate which entries should be replaced by the sub matrix entries.

  • col_idxs (array_like, shape(q<=m,)) – The column indices which designate which entries should be replaced by the sub matrix entries.

  • sub_matrix (ndarray, shape(p, q)) – A matrix of values to substitute into the specified rows and columns.

Notes

This makes a copy of the sub_matrix, so if it is large it may be slower than a more optimal implementation.

Examples

>>> a = np.zeros((3, 4))
>>> sub = np.arange(4).reshape((2, 2))
>>> substitute_matrix(a, [1, 2], [0, 2], sub)
array([[ 0.,  0.,  0.,  0.],
       [ 0.,  0.,  1.,  0.],
       [ 2.,  0.,  3.,  0.]])
opty.utils.sum_of_sines(sigma, frequencies, time)[source]

Returns a sum of sines centered at zero along with its first and second derivatives.

Parameters:
  • sigma (float) – The desired standard deviation of the series.

  • frequencies (iterable of floats) – The frequencies of the sin curves to be included in the sum.

  • time (array_like, shape(n,)) – The montonically increasing time vector.

Returns:

  • sines (ndarray, shape(n,)) – A sum of sines.

  • sines_prime (ndarray, shape(n,)) – The first derivative of the sum of sines.

  • sines_double_prime (ndarray, shape(n,)) – The second derivative of the sum of sines.

opty.utils.ufuncify_matrix(args, expr, const=None, tmp_dir=None, parallel=False)[source]

Returns a function that evaluates a matrix of expressions in a tight loop.

Parameters:
  • args (iterable of sympy.Symbol) – A list of all symbols in expr in the desired order for the output function.

  • expr (sympy.Matrix) – A matrix of expressions.

  • const (tuple, optional) – This should include any of the symbols in args that should be constant with respect to the loop.

  • tmp_dir (string, optional) – The path to a directory in which to store the generated files. If None then the files will be not be retained after the function is compiled.

  • parallel (boolean, optional) – If True and openmp is installed, the generated code will be parallelized across threads. This is only useful when expr are extremely large.