codegen

Introduction

The pydy.codegen package contains various tools to generate numerical code from symbolic descriptions of the equations of motion of systems. It allows you to generate code using a variety of backends depending on your needs. The generated code can also be auto-wrapped for immediate use in a Python session or script. Each component of the code generators and wrappers are accessible so that you can use just the raw code or the wrapper versions.

We currently support three backends:

lambdify

This generates NumPy-aware Python code which is defined in a Python lambda function, using the sympy.utilities.lambdify module and is the default generator.

Theano

This generates Theano trees that are compiled into low level code, using the sympy.printers.theano_code module.

Cython

This generates C code that can be called from Python, using SymPy’s C code printer utilities and Cython.

On Windows

For the Cython backend to work on Windows you must install a suitable compiler. See this Cython wiki page for instructions on getting a compiler installed. The easiest solution is to use the Microsoft Visual C++ Compiler for Python.

Example Use

The simplest entry point to the code generation tools is through the System class.

>>> from pydy.models import multi_mass_spring_damper
>>> sys = multi_mass_spring_damper()
>>> type(sys)
<class 'pydy.system.System'>
>>> rhs = sys.generate_ode_function()
>>> help(rhs) # rhs is a function:
Returns the derivatives of the states, i.e. numerically evaluates the right
hand side of the first order differential equation.

x' = f(x, t, p)

Parameters
==========
x : ndarray, shape(2,)
    The state vector is ordered as such:
        - x0(t)
        - v0(t)
t : float
    The current time.
p : dictionary len(3) or ndarray shape(3,)
    Either a dictionary that maps the constants symbols to their numerical
    values or an array with the constants in the following order:
        - m0
        - c0
        - k0

Returns
=======
dx : ndarray, shape(2,)
    The derivative of the state vector.


>>> import numpy as np
>>> rhs(np.array([1.0, 2.0]), 0.0, np.array([1.0, 2.0, 3.0]))
array([ 2., -7.])

You can also use the functional interface to the code generation/wrapper classes:

>>> from numpy import array
>>> from pydy.models import multi_mass_spring_damper
>>> from pydy.codegen.ode_function_generators import generate_ode_function
>>> sys = multi_mass_spring_damper()
>>> sym_rhs = sys.eom_method.rhs()
>>> q = sys.coordinates
>>> u = sys.speeds
>>> p = sys.constants_symbols
>>> rhs = generate_ode_function(sym_rhs, q, u, p)
>>> rhs(array([1.0, 2.0]), 0.0, array([1.0, 2.0, 3.0]))
array([ 2., -7.])

Other backends can be used by passing in the generator keyword argument, e.g.:

>>> rhs = generate_ode_function(sym_rhs, q, u, p, generator='cython')
>>> rhs(array([1.0, 2.0]), 0.0, array([1.0, 2.0, 3.0]))
array([ 2., -7.])

The backends are implemented as subclasses of ODEFunctionGenerator. You can make use of the ODEFunctionGenerator classes directly:

>>> from pydy.codegen.ode_function_generators import LambdifyODEFunctionGenerator
>>> g = LambdifyODEFunctionGenerator(sym_rhs, q, u, p)
>>> rhs = g.generate()
>>> rhs(array([1.0, 2.0]), 0.0, array([1.0, 2.0, 3.0]))
array([ 2., -7.])

Furthermore, for direct control over evaluating matrices you can use the lamdify and theano_functions in SymPy or utilize the CythonMatrixGenerator class in PyDy. For example, this shows you how to generate C and Cython code to evaluate matrices:

>>> from pydy.codegen.cython_code import CythonMatrixGenerator
>>> sys = multi_mass_spring_damper()
>>> q = sys.coordinates
>>> u = sys.speeds
>>> p = sys.constants_symbols
>>> sym_rhs = sys.eom_method.rhs()
>>> g = CythonMatrixGenerator([q, u, p], [sym_rhs])
>>> setup_py, cython_src, c_header, c_src = g.doprint()
>>> print(setup_py)
#!/usr/bin/env python

from distutils.core import setup
from distutils.extension import Extension

from Cython.Build import cythonize
import numpy

extension = Extension(name="pydy_codegen",
                      sources=["pydy_codegen.pyx",
                               "pydy_codegen_c.c"],
                      include_dirs=[numpy.get_include()])

setup(name="pydy_codegen",
      ext_modules=cythonize([extension]))

>>> print(cython_src)
import numpy as np
cimport numpy as np
cimport cython

cdef extern from "pydy_codegen_c.h":
    void evaluate(
                  double* input_0,
                  double* input_1,
                  double* input_2,
                  double* output_0
                 )

@cython.boundscheck(False)
@cython.wraparound(False)
def eval(
         np.ndarray[np.double_t, ndim=1, mode='c'] input_0,
         np.ndarray[np.double_t, ndim=1, mode='c'] input_1,
         np.ndarray[np.double_t, ndim=1, mode='c'] input_2,
         np.ndarray[np.double_t, ndim=1, mode='c'] output_0
        ):

    evaluate(
             <double*> input_0.data,
             <double*> input_1.data,
             <double*> input_2.data,
             <double*> output_0.data
            )

    return (
            output_0
           )

>>> print(c_src)
#include <math.h>
#include "pydy_codegen_c.h"

void evaluate(
              double input_0[1],
              double input_1[1],
              double input_2[3],
              double output_0[2]
             )
{

    double pydy_0 = input_1[0];

    output_0[0] = pydy_0;
    output_0[1] = (-input_2[1]*pydy_0 - input_2[2]*input_0[0])/input_2[0];

}

>>> print(c_header)
void evaluate(
              double input_0[1],
              double input_1[1],
              double input_2[3],
              double output_0[2]
             );
/*

input_0[1] : [x0(t)]
input_1[1] : [v0(t)]
input_2[3] : [m0, c0, k0]

*/

>>> rhs = g.compile()
>>> res = array([0.0, 0.0])
>>> rhs(array([1.0]), array([2.0]), array([1.0, 2.0, 3.0]), res)
array([ 2., -7.])

We also support generating Octave/Matlab code as shown below:

>>> from pydy.codegen.octave_code import OctaveMatrixGenerator
>>> sys = multi_mass_spring_damper()
>>> q = sys.coordinates
>>> u = sys.speeds
>>> p = sys.constants_symbols
>>> sym_rhs = sys.eom_method.rhs()
>>> g = OctaveMatrixGenerator([q + u, p], [sym_rhs])
>>> m_src = g.doprint()
>>> print(m_src)
function [output_1] = eval_mats(input_1, input_2)
% function [output_1] = eval_mats(input_1, input_2)
%
% input_1 : [x0(t), v0(t)]
% input_2 : [k0, m0, c0]

    pydy_0 = input_1(2);

    output_1 = [pydy_0; (-input_2(3).*pydy_0 - ...
    input_2(1).*input_1(1))./input_2(2)];

end

API

This module contains source code dedicated to generating C code from matrices generated from sympy.physics.mechanics.

class pydy.codegen.c_code.CMatrixGenerator(arguments, matrices, cse=True)[source]

This class generates C source files that simultaneously numerically evaluate any number of SymPy matrices.

doprint(prefix=None)[source]

Returns a string each for the header and the source files.

Parameters
prefixstring, optional

A prefix for the name of the header file. This will cause an include statement to to be added to the source.

write(prefix, path=None)[source]

Writes a header and source file to disk.

Parameters
prefixstring

Two files will be generated: <prefix>.c and <prefix>.h.

class pydy.codegen.cython_code.CythonMatrixGenerator(arguments, matrices, prefix='pydy_codegen', cse=True)[source]

This class generates the Cython code for evaluating a sequence of matrices. It can compile the code and return a Python function.

__init__(arguments, matrices, prefix='pydy_codegen', cse=True)[source]
Parameters
argumentssequences of sequences of SymPy Symbol or Function.

Each of the sequences will be converted to input arrays in the Cython function. All of the symbols/functions contained in matrices need to be in the sequences, but the sequences can also contain extra symbols/functions that are not contained in the matrices.

matricessequence of SymPy.Matrix

A sequence of the matrices that should be evaluated in the function. The expressions should contain only sympy.Symbol or sympy.Function that are functions of me.dynamicsymbols._t.

prefixstring, optional

The desired prefix for the generated files.

cseboolean

Find and replace common sub-expressions in matrices if True.

compile(tmp_dir=None, verbose=False)[source]

Returns a function which evaluates the matrices.

Parameters
tmp_dirstring

The path to an existing or non-existing directory where all of the generated files will be stored.

verboseboolean

If true the output of the completed compilation steps will be printed.

doprint()[source]

Returns the text of the four source files needed to compile Cython wrapper that evaluates the provided SymPy matrices.

Returns
setup_pystring

The text of the setup.py file used to compile the Cython extension.

cython_sourcestring

The text of the Cython pyx file which includes the wrapper function eval.

c_headerstring

The text of the C header file that exposes the evaluate function.

c_sourcestring

The text of the C source file containing the function that evaluates the matrices.

write(path=None)[source]

Writes the four source files needed to compile the Cython function to the current working directory.

Parameters
pathstring

The absolute or relative path to an existing directory to place the files instead of the cwd.

class pydy.codegen.matrix_generator.MatrixGenerator(arguments, matrices, cse=True)[source]

This abstract base class generates source files that simultaneously numerically evaluate any number of SymPy matrices.

__init__(arguments, matrices, cse=True)[source]
Parameters
argumentssequence of sequences of SymPy Symbol or Function

Each of the sequences will be converted to input arrays in the generated function. All of the symbols/functions contained in matrices need to be in the sequences, but the sequences can also contain extra symbols/functions that are not contained in the matrices.

matricessequence of SymPy.Matrix

A sequence of the matrices that should be evaluated in the function. The expressions should contain only sympy.Symbol or sympy.Function that are functions of me.dynamicsymbols._t.

cseboolean

Find and replace common sub-expressions in matrices if True.

comma_lists()[source]

Returns a string output for each of the sequences of SymPy arguments.

class pydy.codegen.octave_code.OctaveMatrixGenerator(arguments, matrices, cse=True)[source]

This class generates Octave/Matlab source files that simultaneously numerically evaluate any number of SymPy matrices.

doprint(prefix='eval_mats')[source]

Returns a string that implements the function.

Parameters
prefixstring, optional

The name of the Octave/Matlab function.

write(prefix='eval_mats', path=None)[source]

Writes the <prefix>.m file to disc at the give path location.

class pydy.codegen.ode_function_generators.CythonODEFunctionGenerator(*args, **kwargs)[source]
__init__(*args, **kwargs)[source]

Generates a numerical function which can evaluate the right hand side of the first order ordinary differential equations from a system described by one of the following three symbolic forms:

[1] x’ = F(x, t, r, p)

[2] M(x, p) x’ = F(x, t, r, p)

[3] M(q, p) u’ = F(q, u, t, r, p)

q’ = G(q, u, t, r, p)

where

x : states, i.e. [q, u] t : time r : specified (exogenous) inputs p : constants q : generalized coordinates u : generalized speeds M : mass matrix (full or minimum) F : right hand side (full or minimum) G : right hand side of the kinematical differential equations

The generated function is of the form F(x, t, p) or F(x, t, r, p) depending on whether the system has specified inputs or not.

Parameters
right_hand_sideSymPy Matrix, shape(n, 1)

A column vector containing the symbolic expressions for the right hand side of the ordinary differential equations. If the right hand side has been solved for symbolically then only F is required, see form [1]; if not then the mass matrix must also be supplied, see forms [2, 3].

coordinatessequence of SymPy Functions

The generalized coordinates. These must be ordered in the same order as the rows in M, F, and/or G and be functions of time.

speedssequence of SymPy Functions

The generalized speeds. These must be ordered in the same order as the rows in M, F, and/or G and be functions of time.

constantssequence of SymPy Symbols, optional

All of the constants present in the equations of motion. The order does not matter.

mass_matrixsympy.Matrix, shape(n, n), optional

This can be either the “full” mass matrix as in [2] or the “minimal” mass matrix as in [3]. The rows and columns must be ordered to match the order of the coordinates and speeds. In the case of the full mass matrix, the speeds should always be ordered before the speeds, i.e. x = [q, u].

coordinate_derivativessympy.Matrix, shape(m, 1), optional

If the “minimal” mass matrix, form [3], is supplied, then this column vector represents the right hand side of the kinematical differential equations.

specifiedssequence of SymPy Functions

The specified exogenous inputs to the system. These should be functions of time and the order does not matter.

linear_sys_solverstring or function

Specify either numpy or scipy to use the linear solvers provided in each package or supply a function that solves a linear system Ax=b with the call signature x = solve(A, b). For example, if you need to use custom kwargs for the SciPy solver, pass in a lambda function that wraps the solver and sets them.

constants_arg_typestring

The generated function accepts two different types of arguments for the numerical values of the constants: either a ndarray of the constants values in the correct order or a dictionary mapping the constants symbols to the numerical values. If None, this is determined inside of the generated function and can cause a significant slow down for performance critical code. If you know apriori what arg types you need to support choose either array or dictionary. Note that array is faster than dictionary.

specifieds_arg_typestring

The generated function accepts three different types of arguments for the numerical values of the specifieds: either a ndarray of the specifieds values in the correct order, a function that generates the correctly ordered ndarray, or a dictionary mapping the specifieds symbols or tuples of thereof to floats, ndarrays, or functions. If None, this is determined inside of the generated function and can cause a significant slow down for performance critical code. If you know apriori what arg types you want to support choose either array, function, or dictionary. The speed of each, from fast to slow, are array, function, dictionary, None.

class pydy.codegen.ode_function_generators.LambdifyODEFunctionGenerator(*args, **kwargs)[source]
__init__(*args, **kwargs)[source]

Generates a numerical function which can evaluate the right hand side of the first order ordinary differential equations from a system described by one of the following three symbolic forms:

[1] x’ = F(x, t, r, p)

[2] M(x, p) x’ = F(x, t, r, p)

[3] M(q, p) u’ = F(q, u, t, r, p)

q’ = G(q, u, t, r, p)

where

x : states, i.e. [q, u] t : time r : specified (exogenous) inputs p : constants q : generalized coordinates u : generalized speeds M : mass matrix (full or minimum) F : right hand side (full or minimum) G : right hand side of the kinematical differential equations

The generated function is of the form F(x, t, p) or F(x, t, r, p) depending on whether the system has specified inputs or not.

Parameters
right_hand_sideSymPy Matrix, shape(n, 1)

A column vector containing the symbolic expressions for the right hand side of the ordinary differential equations. If the right hand side has been solved for symbolically then only F is required, see form [1]; if not then the mass matrix must also be supplied, see forms [2, 3].

coordinatessequence of SymPy Functions

The generalized coordinates. These must be ordered in the same order as the rows in M, F, and/or G and be functions of time.

speedssequence of SymPy Functions

The generalized speeds. These must be ordered in the same order as the rows in M, F, and/or G and be functions of time.

constantssequence of SymPy Symbols, optional

All of the constants present in the equations of motion. The order does not matter.

mass_matrixsympy.Matrix, shape(n, n), optional

This can be either the “full” mass matrix as in [2] or the “minimal” mass matrix as in [3]. The rows and columns must be ordered to match the order of the coordinates and speeds. In the case of the full mass matrix, the speeds should always be ordered before the speeds, i.e. x = [q, u].

coordinate_derivativessympy.Matrix, shape(m, 1), optional

If the “minimal” mass matrix, form [3], is supplied, then this column vector represents the right hand side of the kinematical differential equations.

specifiedssequence of SymPy Functions

The specified exogenous inputs to the system. These should be functions of time and the order does not matter.

linear_sys_solverstring or function

Specify either numpy or scipy to use the linear solvers provided in each package or supply a function that solves a linear system Ax=b with the call signature x = solve(A, b). For example, if you need to use custom kwargs for the SciPy solver, pass in a lambda function that wraps the solver and sets them.

constants_arg_typestring

The generated function accepts two different types of arguments for the numerical values of the constants: either a ndarray of the constants values in the correct order or a dictionary mapping the constants symbols to the numerical values. If None, this is determined inside of the generated function and can cause a significant slow down for performance critical code. If you know apriori what arg types you need to support choose either array or dictionary. Note that array is faster than dictionary.

specifieds_arg_typestring

The generated function accepts three different types of arguments for the numerical values of the specifieds: either a ndarray of the specifieds values in the correct order, a function that generates the correctly ordered ndarray, or a dictionary mapping the specifieds symbols or tuples of thereof to floats, ndarrays, or functions. If None, this is determined inside of the generated function and can cause a significant slow down for performance critical code. If you know apriori what arg types you want to support choose either array, function, or dictionary. The speed of each, from fast to slow, are array, function, dictionary, None.

class pydy.codegen.ode_function_generators.ODEFunctionGenerator(right_hand_side, coordinates, speeds, constants=(), mass_matrix=None, coordinate_derivatives=None, specifieds=None, linear_sys_solver='numpy', constants_arg_type=None, specifieds_arg_type=None)[source]

This is an abstract base class for all of the generators. A subclass is expected to implement the methods necessary to evaluate the arrays needed to compute xdot for the three different system specification types.

__init__(right_hand_side, coordinates, speeds, constants=(), mass_matrix=None, coordinate_derivatives=None, specifieds=None, linear_sys_solver='numpy', constants_arg_type=None, specifieds_arg_type=None)[source]

Generates a numerical function which can evaluate the right hand side of the first order ordinary differential equations from a system described by one of the following three symbolic forms:

[1] x’ = F(x, t, r, p)

[2] M(x, p) x’ = F(x, t, r, p)

[3] M(q, p) u’ = F(q, u, t, r, p)

q’ = G(q, u, t, r, p)

where

x : states, i.e. [q, u] t : time r : specified (exogenous) inputs p : constants q : generalized coordinates u : generalized speeds M : mass matrix (full or minimum) F : right hand side (full or minimum) G : right hand side of the kinematical differential equations

The generated function is of the form F(x, t, p) or F(x, t, r, p) depending on whether the system has specified inputs or not.

Parameters
right_hand_sideSymPy Matrix, shape(n, 1)

A column vector containing the symbolic expressions for the right hand side of the ordinary differential equations. If the right hand side has been solved for symbolically then only F is required, see form [1]; if not then the mass matrix must also be supplied, see forms [2, 3].

coordinatessequence of SymPy Functions

The generalized coordinates. These must be ordered in the same order as the rows in M, F, and/or G and be functions of time.

speedssequence of SymPy Functions

The generalized speeds. These must be ordered in the same order as the rows in M, F, and/or G and be functions of time.

constantssequence of SymPy Symbols, optional

All of the constants present in the equations of motion. The order does not matter.

mass_matrixsympy.Matrix, shape(n, n), optional

This can be either the “full” mass matrix as in [2] or the “minimal” mass matrix as in [3]. The rows and columns must be ordered to match the order of the coordinates and speeds. In the case of the full mass matrix, the speeds should always be ordered before the speeds, i.e. x = [q, u].

coordinate_derivativessympy.Matrix, shape(m, 1), optional

If the “minimal” mass matrix, form [3], is supplied, then this column vector represents the right hand side of the kinematical differential equations.

specifiedssequence of SymPy Functions

The specified exogenous inputs to the system. These should be functions of time and the order does not matter.

linear_sys_solverstring or function

Specify either numpy or scipy to use the linear solvers provided in each package or supply a function that solves a linear system Ax=b with the call signature x = solve(A, b). For example, if you need to use custom kwargs for the SciPy solver, pass in a lambda function that wraps the solver and sets them.

constants_arg_typestring

The generated function accepts two different types of arguments for the numerical values of the constants: either a ndarray of the constants values in the correct order or a dictionary mapping the constants symbols to the numerical values. If None, this is determined inside of the generated function and can cause a significant slow down for performance critical code. If you know apriori what arg types you need to support choose either array or dictionary. Note that array is faster than dictionary.

specifieds_arg_typestring

The generated function accepts three different types of arguments for the numerical values of the specifieds: either a ndarray of the specifieds values in the correct order, a function that generates the correctly ordered ndarray, or a dictionary mapping the specifieds symbols or tuples of thereof to floats, ndarrays, or functions. If None, this is determined inside of the generated function and can cause a significant slow down for performance critical code. If you know apriori what arg types you want to support choose either array, function, or dictionary. The speed of each, from fast to slow, are array, function, dictionary, None.

define_inputs()[source]

Sets self.inputs to the list of sequences [q, u, p] or [q, u, r, p].

generate()[source]

Returns a function that evaluates the right hand side of the first order ordinary differential equations in one of two forms:

x’ = f(x, t, p)

or

x’ = f(x, t, r, p)

See the docstring of the generated function for more details.

static list_syms(indent, syms)[source]

Returns a string representation of a valid rst list of the symbols in the sequence syms and indents the list given the integer number of indentations.

class pydy.codegen.ode_function_generators.TheanoODEFunctionGenerator(*args, **kwargs)[source]
__init__(*args, **kwargs)[source]

Generates a numerical function which can evaluate the right hand side of the first order ordinary differential equations from a system described by one of the following three symbolic forms:

[1] x’ = F(x, t, r, p)

[2] M(x, p) x’ = F(x, t, r, p)

[3] M(q, p) u’ = F(q, u, t, r, p)

q’ = G(q, u, t, r, p)

where

x : states, i.e. [q, u] t : time r : specified (exogenous) inputs p : constants q : generalized coordinates u : generalized speeds M : mass matrix (full or minimum) F : right hand side (full or minimum) G : right hand side of the kinematical differential equations

The generated function is of the form F(x, t, p) or F(x, t, r, p) depending on whether the system has specified inputs or not.

Parameters
right_hand_sideSymPy Matrix, shape(n, 1)

A column vector containing the symbolic expressions for the right hand side of the ordinary differential equations. If the right hand side has been solved for symbolically then only F is required, see form [1]; if not then the mass matrix must also be supplied, see forms [2, 3].

coordinatessequence of SymPy Functions

The generalized coordinates. These must be ordered in the same order as the rows in M, F, and/or G and be functions of time.

speedssequence of SymPy Functions

The generalized speeds. These must be ordered in the same order as the rows in M, F, and/or G and be functions of time.

constantssequence of SymPy Symbols, optional

All of the constants present in the equations of motion. The order does not matter.

mass_matrixsympy.Matrix, shape(n, n), optional

This can be either the “full” mass matrix as in [2] or the “minimal” mass matrix as in [3]. The rows and columns must be ordered to match the order of the coordinates and speeds. In the case of the full mass matrix, the speeds should always be ordered before the speeds, i.e. x = [q, u].

coordinate_derivativessympy.Matrix, shape(m, 1), optional

If the “minimal” mass matrix, form [3], is supplied, then this column vector represents the right hand side of the kinematical differential equations.

specifiedssequence of SymPy Functions

The specified exogenous inputs to the system. These should be functions of time and the order does not matter.

linear_sys_solverstring or function

Specify either numpy or scipy to use the linear solvers provided in each package or supply a function that solves a linear system Ax=b with the call signature x = solve(A, b). For example, if you need to use custom kwargs for the SciPy solver, pass in a lambda function that wraps the solver and sets them.

constants_arg_typestring

The generated function accepts two different types of arguments for the numerical values of the constants: either a ndarray of the constants values in the correct order or a dictionary mapping the constants symbols to the numerical values. If None, this is determined inside of the generated function and can cause a significant slow down for performance critical code. If you know apriori what arg types you need to support choose either array or dictionary. Note that array is faster than dictionary.

specifieds_arg_typestring

The generated function accepts three different types of arguments for the numerical values of the specifieds: either a ndarray of the specifieds values in the correct order, a function that generates the correctly ordered ndarray, or a dictionary mapping the specifieds symbols or tuples of thereof to floats, ndarrays, or functions. If None, this is determined inside of the generated function and can cause a significant slow down for performance critical code. If you know apriori what arg types you want to support choose either array, function, or dictionary. The speed of each, from fast to slow, are array, function, dictionary, None.

define_inputs()[source]

Sets self.inputs to the list of sequences [q, u, p] or [q, u, r, p].

pydy.codegen.ode_function_generators.generate_ode_function(*args, **kwargs)[source]

Generates a numerical function which can evaluate the right hand side of the first order ordinary differential equations from a system described by one of the following three symbolic forms:

[1] x’ = F(x, t, r, p)

[2] M(x, p) x’ = F(x, t, r, p)

[3] M(q, p) u’ = F(q, u, t, r, p)

q’ = G(q, u, t, r, p)

where

x : states, i.e. [q, u] t : time r : specified (exogenous) inputs p : constants q : generalized coordinates u : generalized speeds M : mass matrix (full or minimum) F : right hand side (full or minimum) G : right hand side of the kinematical differential equations

The generated function is of the form F(x, t, p) or F(x, t, r, p) depending on whether the system has specified inputs or not.

Parameters
right_hand_sideSymPy Matrix, shape(n, 1)

A column vector containing the symbolic expressions for the right hand side of the ordinary differential equations. If the right hand side has been solved for symbolically then only F is required, see form [1]; if not then the mass matrix must also be supplied, see forms [2, 3].

coordinatessequence of SymPy Functions

The generalized coordinates. These must be ordered in the same order as the rows in M, F, and/or G and be functions of time.

speedssequence of SymPy Functions

The generalized speeds. These must be ordered in the same order as the rows in M, F, and/or G and be functions of time.

constantssequence of SymPy Symbols, optional

All of the constants present in the equations of motion. The order does not matter.

mass_matrixsympy.Matrix, shape(n, n), optional

This can be either the “full” mass matrix as in [2] or the “minimal” mass matrix as in [3]. The rows and columns must be ordered to match the order of the coordinates and speeds. In the case of the full mass matrix, the speeds should always be ordered before the speeds, i.e. x = [q, u].

coordinate_derivativessympy.Matrix, shape(m, 1), optional

If the “minimal” mass matrix, form [3], is supplied, then this column vector represents the right hand side of the kinematical differential equations.

specifiedssequence of SymPy Functions

The specified exogenous inputs to the system. These should be functions of time and the order does not matter.

linear_sys_solverstring or function

Specify either numpy or scipy to use the linear solvers provided in each package or supply a function that solves a linear system Ax=b with the call signature x = solve(A, b). For example, if you need to use custom kwargs for the SciPy solver, pass in a lambda function that wraps the solver and sets them.

constants_arg_typestring

The generated function accepts two different types of arguments for the numerical values of the constants: either a ndarray of the constants values in the correct order or a dictionary mapping the constants symbols to the numerical values. If None, this is determined inside of the generated function and can cause a significant slow down for performance critical code. If you know apriori what arg types you need to support choose either array or dictionary. Note that array is faster than dictionary.

specifieds_arg_typestring

The generated function accepts three different types of arguments for the numerical values of the specifieds: either a ndarray of the specifieds values in the correct order, a function that generates the correctly ordered ndarray, or a dictionary mapping the specifieds symbols or tuples of thereof to floats, ndarrays, or functions. If None, this is determined inside of the generated function and can cause a significant slow down for performance critical code. If you know apriori what arg types you want to support choose either array, function, or dictionary. The speed of each, from fast to slow, are array, function, dictionary, None.

generator : string or and ODEFunctionGenerator, optional

The method used for generating the numeric right hand side. The string options are {‘lambdify’|’theano’|’cython’} with ‘lambdify’ being the default. You can also pass in a custom subclass of ODEFunctionGenerator.

Returns
rhsfunction

A function which evaluates the derivaties of the states. See the function’s docstring for more details after generation.