codegen API¶
-
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.
-
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 equationsThe 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_side : SymPy 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].
coordinates : sequence 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.
speeds : sequence 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.
constants : sequence of SymPy Symbols
All of the constants present in the equations of motion. The order does not matter.
mass_matrix : sympy.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_derivatives : sympy.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.
specifieds : sequence of SymPy Functions
The specified exogenous inputs to the system. These should be functions of time and the order does not matter.
linear_sys_solver : string 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_type : string
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
ordictionary
. Note thatarray
is faster thandictionary
.specifieds_arg_type : string
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
, ordictionary
. The speed of each, from fast to slow, arearray
,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: rhs : function
A function which evaluates the derivaties of the states. See the function’s docstring for more details after generation.
Notes
The generated function still supports the pre-0.3.0 extra argument style, i.e. args = {‘constants’: ..., ‘specified’: ...}, but only if
constants_arg_type
andspecifieds_arg_type
are both set to None. This functionality is deprecated and will be removed in 0.4.0, so it’s best to adjust your code to support the new argument types. See the docstring for the generated function for more info on the new style of arguments.
-
class
pydy.codegen.cython_code.
CythonMatrixGenerator
(arguments, matrices, prefix='pydy_codegen')[source]¶ This class generates the Cython code for evaluating a sequence of matrices. It can compile the code and return a Python function.
-
compile
(tmp_dir=None, verbose=False)[source]¶ Returns a function which evaluates the matrices.
Parameters: tmp_dir : string
The path to an existing or non-existing directory where all of the generated files will be stored.
verbose : boolean
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_py : string
The text of the setup.py file used to compile the Cython extension.
cython_source : string
The text of the Cython pyx file which includes the wrapper function
eval
.c_header : string
The text of the C header file that exposes the evaluate function.
c_source : string
The text of the C source file containing the function that evaluates the matrices.
-
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.