Kane’s Method, Lagrange’s Method & Associated Functions (Docstrings)

Kane

class sympy.physics.mechanics.kane.KanesMethod(frame, q_ind, u_ind, kd_eqs=None, q_dependent=, []configuration_constraints=, []u_dependent=, []velocity_constraints=, []acceleration_constraints=None, u_auxiliary=[])

Kane’s method object.

This object is used to do the “book-keeping” as you go through and form equations of motion in the way Kane presents in: Kane, T., Levinson, D. Dynamics Theory and Applications. 1985 McGraw-Hill

The attributes are for equations in the form [M] udot = forcing.

Examples

This is a simple example for a one degree of freedom translational spring-mass-damper.

In this example, we first need to do the kinematics. This involves creating generalized speeds and coordinates and their derivatives. Then we create a point and set its velocity in a frame:

>>> from sympy import symbols
>>> from sympy.physics.mechanics import dynamicsymbols, ReferenceFrame
>>> from sympy.physics.mechanics import Point, Particle, KanesMethod
>>> q, u = dynamicsymbols('q u')
>>> qd, ud = dynamicsymbols('q u', 1)
>>> m, c, k = symbols('m c k')
>>> N = ReferenceFrame('N')
>>> P = Point('P')
>>> P.set_vel(N, u * N.x)

Next we need to arrange/store information in the way that KanesMethod requires. The kinematic differential equations need to be stored in a dict. A list of forces/torques must be constructed, where each entry in the list is a (Point, Vector) or (ReferenceFrame, Vector) tuple, where the Vectors represent the Force or Torque. Next a particle needs to be created, and it needs to have a point and mass assigned to it. Finally, a list of all bodies and particles needs to be created:

>>> kd = [qd - u]
>>> FL = [(P, (-k * q - c * u) * N.x)]
>>> pa = Particle('pa', P, m)
>>> BL = [pa]

Finally we can generate the equations of motion. First we create the KanesMethod object and supply an inertial frame, coordinates, generalized speeds, and the kinematic differential equations. Additional quantities such as configuration and motion constraints, dependent coordinates and speeds, and auxiliary speeds are also supplied here (see the online documentation). Next we form FR* and FR to complete: Fr + Fr* = 0. We have the equations of motion at this point. It makes sense to rearrnge them though, so we calculate the mass matrix and the forcing terms, for E.o.M. in the form: [MM] udot = forcing, where MM is the mass matrix, udot is a vector of the time derivatives of the generalized speeds, and forcing is a vector representing “forcing” terms:

>>> KM = KanesMethod(N, q_ind=[q], u_ind=[u], kd_eqs=kd)
>>> (fr, frstar) = KM.kanes_equations(FL, BL)
>>> MM = KM.mass_matrix
>>> forcing = KM.forcing
>>> rhs = MM.inv() * forcing
>>> rhs
Matrix([[(-c*u(t) - k*q(t))/m]])
>>> KM.linearize()[0]
Matrix([
[ 0,  1],
[-k, -c]])

Please look at the documentation pages for more information on how to perform linearization and how to deal with dependent coordinates & speeds, and how do deal with bringing non-contributing forces into evidence.

Attributes

auxiliary Matrix If applicable, the set of auxiliary Kane’s equations used to solve for non-contributing forces.
mass_matrix Matrix The system’s mass matrix
forcing Matrix The system’s forcing vector
mass_matrix_full Matrix The “mass matrix” for the u’s and q’s
forcing_full Matrix The “forcing vector” for the u’s and q’s
kanes_equations(FL, BL)

Method to form Kane’s equations, Fr + Fr* = 0.

Returns (Fr, Fr*). In the case where auxiliary generalized speeds are present (say, s auxiliary speeds, o generalized speeds, and m motion constraints) the length of the returned vectors will be o - m + s in length. The first o - m equations will be the constrained Kane’s equations, then the s auxiliary Kane’s equations. These auxiliary equations can be accessed with the auxiliary_eqs().

Parameters:

FL : list

Takes in a list of (Point, Vector) or (ReferenceFrame, Vector) tuples which represent the force at a point or torque on a frame.

BL : list

A list of all RigidBody’s and Particle’s in the system.

kindiffdict()

Returns the qdot’s in a dictionary.

linearize()

Method used to generate linearized equations.

Note that for linearization, it is assumed that time is not perturbed, but only coordinates and positions. The “forcing” vector’s jacobian is computed with respect to the state vector in the form [Qi, Qd, Ui, Ud]. This is the “f_lin_A” matrix.

It also finds any non-state dynamicsymbols and computes the jacobian of the “forcing” vector with respect to them. This is the “f_lin_B” matrix; if this is empty, an empty matrix is created.

Consider the following: If our equations are: [M]qudot = f, where [M] is the full mass matrix, qudot is a vector of the deriatives of the coordinates and speeds, and f in the full forcing vector, the linearization process is as follows: [M]qudot = [f_lin_A]qu + [f_lin_B]y, where qu is the state vector, f_lin_A is the jacobian of the full forcing vector with respect to the state vector, f_lin_B is the jacobian of the full forcing vector with respect to any non-speed/coordinate dynamicsymbols which show up in the full forcing vector, and y is a vector of those dynamic symbols (each column in f_lin_B corresponds to a row of the y vector, each of which is a non-speed/coordinate dynamicsymbol).

To get the traditional state-space A and B matrix, you need to multiply the f_lin_A and f_lin_B matrices by the inverse of the mass matrix. Caution needs to be taken when inverting large symbolic matrices; substituting in numerical values before inverting will work better.

A tuple of (f_lin_A, f_lin_B, other_dynamicsymbols) is returned.

rhs(inv_method=None)

Returns the system’s equations of motion in first order form.

The output of this will be the right hand side of:

[qdot, udot].T = f(q, u, t)

Or, the equations of motion in first order form. The right hand side is what is needed by most numerical ODE integrators.

Parameters:

inv_method : str

The specific sympy inverse matrix calculation method to use.

partial_velocity

sympy.physics.mechanics.functions.partial_velocity(vel_list, u_list, frame)

Returns a list of partial velocities.

For a list of velocity or angular velocity vectors the partial derivatives with respect to the supplied generalized speeds are computed, in the specified ReferenceFrame.

The output is a list of lists. The outer list has a number of elements equal to the number of supplied velocity vectors. The inner lists are, for each velocity vector, the partial derivatives of that velocity vector with respect to the generalized speeds supplied.

Parameters:

vel_list : list

List of velocities of Point’s and angular velocities of ReferenceFrame’s

u_list : list

List of independent generalized speeds.

frame : ReferenceFrame

The ReferenceFrame the partial derivatives are going to be taken in.

Examples

>>> from sympy.physics.mechanics import Point, ReferenceFrame
>>> from sympy.physics.mechanics import dynamicsymbols
>>> from sympy.physics.mechanics import partial_velocity
>>> u = dynamicsymbols('u')
>>> N = ReferenceFrame('N')
>>> P = Point('P')
>>> P.set_vel(N, u * N.x)
>>> vel_list = [P.vel(N)]
>>> u_list = [u]
>>> partial_velocity(vel_list, u_list, N)
[[N.x]]

LagrangesMethod

class sympy.physics.mechanics.lagrange.LagrangesMethod(Lagrangian, q_list, coneqs=None, forcelist=None, frame=None)

Lagrange’s method object.

This object generates the equations of motion in a two step procedure. The first step involves the initialization of LagrangesMethod by supplying the Lagrangian and a list of the generalized coordinates, at the bare minimum. If there are any constraint equations, they can be supplied as keyword arguments. The Lagrangian multipliers are automatically generated and are equal in number to the constraint equations.Similarly any non-conservative forces can be supplied in a list (as described below and also shown in the example) along with a ReferenceFrame. This is also discussed further in the __init__ method.

Examples

This is a simple example for a one degree of freedom translational spring-mass-damper.

In this example, we first need to do the kinematics.$ This involves creating generalized coordinates and its derivative. Then we create a point and set its velocity in a frame:

>>> from sympy.physics.mechanics import LagrangesMethod, Lagrangian
>>> from sympy.physics.mechanics import ReferenceFrame, Particle, Point
>>> from sympy.physics.mechanics import dynamicsymbols, kinetic_energy
>>> from sympy import symbols
>>> q = dynamicsymbols('q')
>>> qd = dynamicsymbols('q', 1)
>>> m, k, b = symbols('m k b')
>>> N = ReferenceFrame('N')
>>> P = Point('P')
>>> P.set_vel(N, qd * N.x)

We need to then prepare the information as required by LagrangesMethod to generate equations of motion. First we create the Particle, which has a point attached to it. Following this the lagrangian is created from the kinetic and potential energies. Then, a list of nonconservative forces/torques must be constructed, where each entry in is a (Point, Vector) or (ReferenceFrame, Vector) tuple, where the Vectors represent the nonconservative force or torque.

>>> Pa = Particle('Pa', P, m)
>>> Pa.set_potential_energy(k * q**2 / 2.0)
>>> L = Lagrangian(N, Pa)
>>> fl = [(P, -b * qd * N.x)]

Finally we can generate the equations of motion. First we create the LagrangesMethod object.To do this one must supply an the Lagrangian, the list of generalized coordinates. Also supplied are the constraint equations, the forcelist and the inertial frame, if relevant. Next we generate Lagrange’s equations of motion, such that: Lagrange’s equations of motion = 0. We have the equations of motion at this point.

>>> l = LagrangesMethod(L, [q], forcelist = fl, frame = N)
>>> print(l.form_lagranges_equations())
Matrix([[b*Derivative(q(t), t) + 1.0*k*q(t) + m*Derivative(q(t), t, t)]])

We can also solve for the states using the ‘rhs’ method.

>>> print(l.rhs())
Matrix([[Derivative(q(t), t)], [(-b*Derivative(q(t), t) - 1.0*k*q(t))/m]])

Please refer to the docstrings on each method for more details.

Attributes

mass_matrix Matrix The system’s mass matrix
forcing Matrix The system’s forcing vector
mass_matrix_full Matrix The “mass matrix” for the qdot’s, qdoubledot’s, and the lagrange multipliers (lam)
forcing_full Matrix The forcing vector for the qdot’s, qdoubledot’s and lagrange multipliers (lam)
forcing

Returns the forcing vector from ‘lagranges_equations’ method.

forcing_full

Augments qdots to the forcing vector above.

form_lagranges_equations()

Method to form Lagrange’s equations of motion.

Returns a vector of equations of motion using Lagrange’s equations of the second kind.

mass_matrix

Returns the mass matrix, which is augmented by the Lagrange multipliers, if necessary.

If the system is described by ‘n’ generalized coordinates and there are no constraint equations then an n X n matrix is returned.

If there are ‘n’ generalized coordinates and ‘m’ constraint equations have been supplied during initialization then an n X (n+m) matrix is returned. The (n + m - 1)th and (n + m)th columns contain the coefficients of the Lagrange multipliers.

mass_matrix_full

Augments the coefficients of qdots to the mass_matrix.

rhs(method='GE')

Returns equations that can be solved numerically

Parameters:

method : string

The method by which matrix inversion of mass_matrix_full must be performed such as Gauss Elimination or LU decomposition.