PyTrajectory Modules Reference

PyTrajectory is a Python library for the determination of the feed forward control to achieve a transition between desired states of a nonlinear control system.

system Module

class pytrajectory.system.ControlSystem(ff, a=0.0, b=1.0, xa=[], xb=[], ua=[], ub=[], constraints=None, **kwargs)[source]

Bases: object

Base class of the PyTrajectory project.

Parameters:
  • ff (callable) – Vector field (rhs) of the control system.
  • a (float) – Left border of the considered time interval.
  • b (float) – Right border of the considered time interval.
  • xa (list) – Boundary values at the left border.
  • xb (list) – Boundary values at the right border.
  • ua (list) – Boundary values of the input variables at left border.
  • ub (list) – Boundary values of the input variables at right border.
  • constraints (dict) – Box-constraints of the state variables.
  • kwargs
    key default value meaning
    sx 10 Initial number of spline parts for the system variables
    su 10 Initial number of spline parts for the input variables
    kx 2 Factor for raising the number of spline parts
    maxIt 10 Maximum number of iteration steps
    eps 1e-2 Tolerance for the solution of the initial value problem
    ierr 1e-1 Tolerance for the error on the whole interval
    tol 1e-5 Tolerance for the solver of the equation system
    dt_sim 1e-2 Sample time for integration (initial value problem)
    use_chains True Whether or not to use integrator chains
    sol_steps 100 Maximum number of iteration steps for the eqs solver
    first_guess None to initiate free parameters (might be useful: {‘seed’: value})
set_param(param='', value=None)[source]

Alters the value of the method parameters.

Parameters:
  • param (str) – The method parameter
  • value – The new value
unconstrain(constraints)[source]

This method is used to enable compliance with desired box constraints given by the user. It transforms the vectorfield by projecting the constrained state variables on new unconstrained ones.

Parameters:constraints (dict) – The box constraints for the state variables
constrain()[source]

This method is used to determine the solution of the original constrained state variables by creating a composition of the saturation functions and the calculated solution for the introduced unconstrained variables.

solve()[source]

This is the main loop.

While the desired accuracy has not been reached, the collocation system will be set up and solved with a iteratively raised number of spline parts.

Returns:
  • callable – Callable function for the system state.
  • callable – Callable function for the input variables.
simulate()[source]

This method is used to solve the resulting initial value problem after the computation of a solution for the input trajectories.

check_accuracy()[source]

Checks whether the desired accuracy for the boundary values was reached.

It calculates the difference between the solution of the simulation and the given boundary values at the right border and compares its maximum against the tolerance.

If set by the user it also calculates some kind of consistency error that shows how “well” the spline functions comply with the system dynamic given by the vector field.

plot()[source]

Plot the calculated trajectories and show interval error functions.

This method calculates the error functions and then calls the visualisation.plotsim function.

save(fname=None)[source]

Save data using the python module pickle.

a
b
class pytrajectory.system.DynamicalSystem(f_sym, a=0.0, b=1.0, xa=[], xb=[], ua=[], ub=[])[source]

Bases: object

Provides access to information about the dynamical system that is the object of the control process.

Parameters:
  • f_sym (callable) – The (symbolic) vector field of the dynamical system
  • b (a,) – The initial end final time of the control process
  • xb (xa,) – The initial and final conditions for the state variables
  • ub (ua,) – The initial and final conditions for the input variables

trajectories Module

class pytrajectory.trajectories.Trajectory(sys, **kwargs)[source]

Bases: object

This class handles the creation and managing of the spline functions that are intended to approximate the desired trajectories.

Parameters:sys (system.DynamicalSystem) – Instance of a dynamical system providing information like vector field function and boundary values
n_parts_x

Number of polynomial spline parts for system variables.

n_parts_u

Number of polynomial spline parts for input variables.

x(t)[source]

Returns the current system state.

Parameters:t (float) – The time point in (a,b) to evaluate the system at.
u(t)[source]

Returns the state of the input variables.

Parameters:t (float) – The time point in (a,b) to evaluate the input variables at.
dx(t)[source]

Returns the state of the 1st derivatives of the system variables.

Parameters:t (float) – The time point in (a,b) to evaluate the 1st derivatives at.
init_splines()[source]

This method is used to create the necessary spline function objects.

Parameters:boundary_values (dict) – Dictionary of boundary values for the state and input splines functions.
set_coeffs(sol)[source]

Set found numerical values for the independent parameters of each spline.

This method is used to get the actual splines by using the numerical solutions to set up the coefficients of the polynomial spline parts of every created spline.

Parameters:sol (numpy.ndarray) – The solution vector for the free parameters, i.e. the independent coefficients.
save()[source]

collocation Module

class pytrajectory.collocation.CollocationSystem(sys, **kwargs)[source]

Bases: object

This class represents the collocation system that is used to determine a solution for the free parameters of the control system, i.e. the independent coefficients of the trajectory splines.

Parameters:sys (system.DynamicalSystem) – Instance of a dynamical system
build()[source]

This method is used to set up the equations for the collocation equation system and defines functions for the numerical evaluation of the system and its jacobian.

get_guess()[source]

This method is used to determine a starting value (guess) for the solver of the collocation equation system.

If it is the first iteration step, then a vector with the same length as the vector of the free parameters with arbitrary values is returned.

Else, for every variable a spline has been created for, the old spline of the iteration before and the new spline are evaluated at specific points and a equation system is solved which ensures that they are equal in these points.

The solution of this system is the new start value for the solver.

save()[source]
solve(G, DG, new_solver=True)[source]

This method is used to solve the collocation equation system.

Parameters:
  • G (callable) – Function that “evaluates” the equation system.
  • DG (callable) – Function for the jacobian.
  • new_solver (bool) – flag to determine whether a new solver instance should be initialized (default True)
class pytrajectory.collocation.Container(**kwargs)[source]

Bases: object

Simple data structure to store additional internal information for debugging and checking the algorithms. Some of the attributes might indeed be neccessary

pytrajectory.collocation.collocation_nodes(a, b, npts, coll_type)[source]

Create collocation points/nodes for the equation system.

Parameters:
  • a (float) – The left border of the considered interval.
  • b (float) – The right border of the considered interval.
  • npts (int) – The number of nodes.
  • coll_type (str) – Specifies how to generate the nodes.
Returns:

The collocation nodes.

Return type:

numpy.ndarray

splines Module

class pytrajectory.splines.Spline(a=0.0, b=1.0, n=5, bv={}, tag='', use_std_approach=False, **kwargs)[source]

Bases: object

This class provides a representation of a cubic spline function.

It simultaneously enables access to the spline function itself as well as to its derivatives up to the 3rd order. Furthermore it has its own method to ensure the steadiness and smoothness conditions of its polynomial parts in the joining points.

For more information see: Candidate Functions

Parameters:
  • a (float) – Left border of the spline interval.
  • b (float) – Right border of the spline interval.
  • n (int) – Number of polynomial parts the spline will be devided up into.
  • tag (str) – The ‘name’ of the spline object.
  • bv (dict) – Boundary values the spline function and/or its derivatives should satisfy.
  • use_std_approach (bool) – Whether to use the standard spline interpolation approach or the ones used in the project thesis
f(t)[source]

This is just a wrapper to evaluate the spline itself.

df(t)[source]

This is just a wrapper to evaluate the spline’s 1st derivative.

ddf(t)[source]

This is just a wrapper to evaluate the spline’s 2nd derivative.

dddf(t)[source]

This is just a wrapper to evaluate the spline’s 3rd derivative.

boundary_values
make_steady()[source]

Please see pytrajectory.splines.make_steady

differentiate(d=1, new_tag='')[source]

Returns the d-th derivative of this spline function object.

Parameters:d (int) – The derivation order.
get_dependence_vectors(points, d=0)[source]

This method yields a provisionally evaluation of the spline while there are no numerical values for its free parameters.

It returns a two vectors which reflect the dependence of the spline’s or its d-th derivative’s coefficients on its free parameters (independent coefficients).

Parameters:
  • points (float) – The points to evaluate the provisionally spline at.
  • d (int) – The derivation order.
set_coefficients(free_coeffs=None, coeffs=None)[source]

This function is used to set up numerical values either for all the spline’s coefficients or its independent ones.

Parameters:
  • free_coeffs (numpy.ndarray) – Array with numerical values for the free coefficients of the spline.
  • coeffs (numpy.ndarray) – Array with coefficients of the polynomial spline parts.
interpolate(fnc=None, m0=None, mn=None)[source]

Determines the spline’s coefficients such that it interpolates a given function.

save()[source]
plot(show=True, ret_array=False)[source]

Plots the spline function or returns an array with its values at some points of the spline interval.

Parameters:
  • show (bool) – Whethter to plot the spline’s curve or not.
  • ret_array (bool) – Wheter to return an array with values of the spline at points of the interval.
pytrajectory.splines.get_spline_nodes(a=0.0, b=1.0, n=10, nodes_type='equidistant')[source]

Generates n spline nodes in the interval [a,b] of given type.

Parameters:
  • a (float) – Lower border of the considered interval.
  • b (float) – Upper border of the considered interval.
  • n (int) – Number of nodes to generate.
  • nodes_type (str) – How to generate the nodes.
pytrajectory.splines.differentiate(spline_fnc)[source]

Returns the derivative of a callable spline function.

Parameters:spline_fnc (callable) – The spline function to derivate.
pytrajectory.splines.make_steady(S)[source]

This method sets up and solves equations that satisfy boundary conditions and ensure steadiness and smoothness conditions of the spline S in every joining point.

Please see the documentation for more details: Candidate Functions

Parameters:S (Spline) – The spline function object for which to solve smoothness and boundary conditions.
pytrajectory.splines.get_smoothness_matrix(S, N1, N2)[source]

Returns the coefficient matrix and right hand site for the equation system that ensures the spline’s smoothness in its joining points and its compliance with the boundary conditions.

Parameters:
  • S (Spline) – The spline function object to get the matrix for.
  • N1 (int) – First dimension of the matrix.
  • N2 (int) – Second dimension of the matrix.
Returns:

  • array_like – The coefficient matrix for the equation system.
  • array_like – The right hand site of the equation system.

solver Module

class pytrajectory.solver.Solver(F, DF, x0, tol=1e-05, reltol=2e-05, maxIt=50, method='leven', mu=0.0001)[source]

This class provides solver for the collocation equation system.

Parameters:
  • F (callable) – The callable function that represents the equation system
  • DF (callable) – The function for the jacobian matrix of the eqs
  • x0 (numpy.ndarray) – The start value for the sover
  • tol (float) – The (absolute) tolerance of the solver
  • maxIt (int) – The maximum number of iterations of the solver
  • method (str) – The solver to use
solve()[source]

This is just a wrapper to call the chosen algorithm for solving the collocation equation system.

leven()[source]

This method is an implementation of the Levenberg-Marquardt-Method to solve nonlinear least squares problems.

For more information see: Levenberg-Marquardt Method

simulation Module

class pytrajectory.simulation.Simulator(ff, T, start, u, dt=0.01)[source]

Bases: object

This class simulates the initial value problem that results from solving the boundary value problem of the control system.

Parameters:
  • ff (callable) – Vectorfield of the control system.
  • T (float) – Simulation time.
  • u (callable) – Function of the input variables.
  • dt (float) – Time step.
rhs(t, x)[source]

Retruns the right hand side (vector field) of the ode system.

calcStep()[source]

Calculates one step of the simulation.

simulate()[source]

Starts the simulation

Returns:
Return type:List of numpy arrays with time steps and simulation data of system and input variables.

auxiliary Module

exception pytrajectory.auxiliary.NanError[source]

Bases: exceptions.ValueError

class pytrajectory.auxiliary.IntegChain(lst)[source]

Bases: object

This class provides a representation of an integrator chain.

For the elements (x_i)_{i=1,...,n} of the chain the relation \dot{x}_i = x_{i+1} applies.

Parameters:lst (list) – Ordered list of the integrator chain’s elements.
elements

tuple – Ordered list of all elements that are part of the integrator chain

upper

str – Upper end of the integrator chain

lower

str – Lower end of the integrator chain

elements

Return an ordered list of the integrator chain’s elements.

upper

Returns the upper end of the integrator chain, i.e. the element of which all others are derivatives of.

lower

Returns the lower end of the integrator chain, i.e. the element which has no derivative in the integrator chain.

pytrajectory.auxiliary.find_integrator_chains(dyn_sys)[source]

Searches for integrator chains in given vector field matrix fi, i.e. equations of the form \dot{x}_i = x_j.

Parameters:dyn_sys (pytrajectory.system.DynamicalSystem) – Instance of a dynamical system
Returns:
  • list – Found integrator chains.
  • list – Indices of the equations that have to be solved using collocation.
pytrajectory.auxiliary.sym2num_vectorfield(f_sym, x_sym, u_sym, vectorized=False, cse=False)[source]

This function takes a callable vector field of a control system that is to be evaluated with symbols for the state and input variables and returns a corresponding function that can be evaluated with numeric values for these variables.

Parameters:
  • f_sym (callable or array_like) – The callable (“symbolic”) vector field of the control system.
  • x_sym (iterable) – The symbols for the state variables of the control system.
  • u_sym (iterable) – The symbols for the input variables of the control system.
  • vectorized (bool) – Whether or not to return a vectorized function.
  • cse (bool) – Whether or not to make use of common subexpressions in vector field
Returns:

The callable (“numeric”) vector field of the control system.

Return type:

callable

pytrajectory.auxiliary.check_expression(expr)[source]

Checks whether a given expression is a sympy epression or a list of sympy expressions.

Throws an exception if not.

pytrajectory.auxiliary.make_cse_eval_function(input_args, replacement_pairs, ret_filter=None, namespace=None)[source]

Returns a function that evaluates the replacement pairs created by the sympy cse.

Parameters:
  • input_args (iterable) – List of additional symbols that are necessary to evaluate the replacement pairs
  • replacement_pairs (iterable) – List of (Symbol, expression) pairs created from sympy cse
  • ret_filter (iterable) – List of sympy symbols of those replacements that should be returned from the created function (if None, all are returned)
  • namespace (dict) – A namespace in which to define the function
pytrajectory.auxiliary.cse_lambdify(args, expr, **kwargs)[source]

Wrapper for sympy.lambdify which makes use of common subexpressions.

pytrajectory.auxiliary.saturation_functions(y_fnc, dy_fnc, y0, y1)[source]

Creates callable saturation function and its first derivative to project the solution found for an unconstrained state variable back on the original constrained one.

For more information, please have a look at Handling constraints.

Parameters:
  • y_fnc (callable) – The calculated solution function for an unconstrained variable.
  • dy_fnc (callable) – The first derivative of the unconstrained solution function.
  • y0 (float) – Lower saturation limit.
  • y1 (float) – Upper saturation limit.
Returns:

  • callable – A callable of a saturation function applied to a calculated solution for an unconstrained state variable.
  • callable – A callable for the first derivative of a saturation function applied to a calculated solution for an unconstrained state variable.

pytrajectory.auxiliary.consistency_error(I, x_fnc, u_fnc, dx_fnc, ff_fnc, npts=500, return_error_array=False)[source]

Calculates an error that shows how “well” the spline functions comply with the system dynamic given by the vector field.

Parameters:
  • I (tuple) – The considered time interval.
  • x_fnc (callable) – A function for the state variables.
  • u_fnc (callable) – A function for the input variables.
  • dx_fnc (callable) – A function for the first derivatives of the state variables.
  • ff_fnc (callable) – A function for the vectorfield of the control system.
  • npts (int) – Number of point to determine the error at.
  • return_error_array (bool) – Whether or not to return the calculated errors (mainly for plotting).
Returns:

  • float – The maximum error between the systems dynamic and its approximation.
  • numpy.ndarray – An array with all errors calculated on the interval.

visualisation Module

pytrajectory.visualisation.plot_simulation(sim_data, H=[], fname=None)[source]

This method provides graphics for each system variable, manipulated variable and error function and plots the solution of the simulation.

Parameters:
  • sim_data (tuple) – Contains collocation points, and simulation results of system and input variables.
  • H (dict) – Dictionary of the callable error functions
  • fname (str) – If not None, plot will be saved as <fname>.png
class pytrajectory.visualisation.Animation(drawfnc, simdata, plotsys=[], plotinputs=[], rcParams=None)[source]

Provides animation capabilities.

Given a callable function that draws an image of the system state and smiulation data this class provides a method to created an animated representation of the system.

Parameters:
  • drawfnc (callable) – Function that returns an image of the current system state according to simdata
  • simdata (numpy.ndarray) – Array that contains simulation data (time, system states, input states)
  • plotsys (list) – List of tuples with indices and labels of system variables that will be plotted along the picture
  • plotinputs (list) – List of tuples with indices and labels of input variables that will be plotted along the picture
class Image[source]

This is just a container for the drawn system.

reset()[source]
Animation.get_axes()[source]
Animation.set_limits(ax='ax_img', xlim=(0, 1), ylim=(0, 1))[source]
Animation.set_label(ax='ax_img', label='')[source]
Animation.show(t=0.0, xlim=None, ylim=None, axes_callback=None, save_fname=None, show=True)[source]

Plots one frame of the system animation.

Parameters:t (float) – The time for which to plot the system
Animation.animate()[source]

Starts the animation of the system.

Animation.save(fname, fps=None, dpi=200)[source]

Saves the animation as a video file or animated gif.