Kane’s Method & Lagrange’s Method (Docstrings)¶
KaneMethod¶
-
class
sympy.physics.mechanics.kane.
KanesMethod
(frame, q_ind, u_ind, kd_eqs=None, q_dependent=None, configuration_constraints=None, u_dependent=None, velocity_constraints=None, acceleration_constraints=None, u_auxiliary=None)[source]¶ 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(BL, FL) >>> MM = KM.mass_matrix >>> forcing = KM.forcing >>> rhs = MM.inv() * forcing >>> rhs Matrix([[(-c*u(t) - k*q(t))/m]]) >>> KM.linearize(A_and_B=True, new_method=True)[0] Matrix([ [ 0, 1], [-k/m, -c/m]])
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
q, u (Matrix) Matrices of the generalized coordinates and speeds bodylist (iterable) Iterable of Point and RigidBody objects in the system. forcelist (iterable) Iterable of (Point, vector) or (ReferenceFrame, vector) tuples describing the forces on the system. 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 -
auxiliary_eqs
¶ A matrix containing the auxiliary equations.
-
forcing
¶ The forcing vector of the system.
-
forcing_full
¶ The forcing vector of the system, augmented by the kinematic differential equations.
-
kanes_equations
(bodies, loads=None)[source]¶ 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: bodies : iterable
An iterable of all RigidBody’s and Particle’s in the system. A system must have at least one body.
loads : iterable
Takes in an iterable of (Particle, Vector) or (ReferenceFrame, Vector) tuples which represent the force at a point or torque on a frame. Must be either a non-empty iterable of tuples or None which corresponds to a system with no constraints.
-
linearize
(**kwargs)[source]¶ Linearize the equations of motion about a symbolic operating point.
If kwarg A_and_B is False (default), returns M, A, B, r for the linearized form, M*[q’, u’]^T = A*[q_ind, u_ind]^T + B*r.
If kwarg A_and_B is True, returns A, B, r for the linearized form dx = A*x + B*r, where x = [q_ind, u_ind]^T. Note that this is computationally intensive if there are many symbolic parameters. For this reason, it may be more desirable to use the default A_and_B=False, returning M, A, and B. Values may then be substituted in to these matrices, and the state space form found as A = P.T*M.inv()*A, B = P.T*M.inv()*B, where P = Linearizer.perm_mat.
In both cases, r is found as all dynamicsymbols in the equations of motion that are not part of q, u, q’, or u’. They are sorted in canonical form.
The operating points may be also entered using the
op_point
kwarg. This takes a dictionary of {symbol: value}, or a an iterable of such dictionaries. The values may be numberic or symbolic. The more values you can specify beforehand, the faster this computation will run.As part of the deprecation cycle, the new method will not be used unless the kwarg
new_method
is set to True. If the kwarg is missing, or set to false, the old linearization method will be used. After next release the need for this kwarg will be removed.For more documentation, please see the
Linearizer
class.
-
mass_matrix
¶ The mass matrix of the system.
-
mass_matrix_full
¶ The mass matrix of the system, augmented by the kinematic differential equations.
-
rhs
(inv_method=None)[source]¶ 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. For a list of valid methods, see
inv()
-
LagrangesMethod¶
-
class
sympy.physics.mechanics.lagrange.
LagrangesMethod
(Lagrangian, qs, coneqs=None, forcelist=None, frame=None, hol_coneqs=None, nonhol_coneqs=None)[source]¶ 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 the generalized coordinates, at the bare minimum. If there are any constraint equations, they can be supplied as keyword arguments. The Lagrange multipliers are automatically generated and are equal in number to the constraint equations. Similarly any non-conservative forces can be supplied in an iterable (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 their derivatives. 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, an iterable of nonconservative forces/torques must be constructed, where each item is a (Point, Vector) or (ReferenceFrame, Vector) tuple, with the Vectors representing the nonconservative forces or torques.
>>> Pa = Particle('Pa', P, m) >>> Pa.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 the Lagrangian, and the generalized coordinates. The constraint equations, the forcelist, and the inertial frame may also be provided, 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
q, u (Matrix) Matrices of the generalized coordinates and speeds forcelist (iterable) Iterable of (Point, vector) or (ReferenceFrame, vector) tuples describing the forces on the system. 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
()[source]¶ Method to form Lagrange’s equations of motion.
Returns a vector of equations of motion using Lagrange’s equations of the second kind.
-
linearize
(q_ind=None, qd_ind=None, q_dep=None, qd_dep=None, **kwargs)[source]¶ Linearize the equations of motion about a symbolic operating point.
If kwarg A_and_B is False (default), returns M, A, B, r for the linearized form, M*[q’, u’]^T = A*[q_ind, u_ind]^T + B*r.
If kwarg A_and_B is True, returns A, B, r for the linearized form dx = A*x + B*r, where x = [q_ind, u_ind]^T. Note that this is computationally intensive if there are many symbolic parameters. For this reason, it may be more desirable to use the default A_and_B=False, returning M, A, and B. Values may then be substituted in to these matrices, and the state space form found as A = P.T*M.inv()*A, B = P.T*M.inv()*B, where P = Linearizer.perm_mat.
In both cases, r is found as all dynamicsymbols in the equations of motion that are not part of q, u, q’, or u’. They are sorted in canonical form.
The operating points may be also entered using the
op_point
kwarg. This takes a dictionary of {symbol: value}, or a an iterable of such dictionaries. The values may be numberic or symbolic. The more values you can specify beforehand, the faster this computation will run.For more documentation, please see the
Linearizer
class.
-
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
(inv_method=None, **kwargs)[source]¶ Returns equations that can be solved numerically
Parameters: inv_method : str
The specific sympy inverse matrix calculation method to use. For a list of valid methods, see
inv()
-
solve_multipliers
(op_point=None, sol_type='dict')[source]¶ Solves for the values of the lagrange multipliers symbolically at the specified operating point
Parameters: op_point : dict or iterable of dicts, optional
Point at which to solve at. The operating point is specified as a dictionary or iterable of dictionaries of {symbol: value}. The value may be numeric or symbolic itself.
sol_type : str, optional
Solution return type. Valid options are: - ‘dict’: A dict of {symbol : value} (default) - ‘Matrix’: An ordered column matrix of the solution
-
to_linearizer
(q_ind=None, qd_ind=None, q_dep=None, qd_dep=None)[source]¶ Returns an instance of the Linearizer class, initiated from the data in the LagrangesMethod class. This may be more desirable than using the linearize class method, as the Linearizer object will allow more efficient recalculation (i.e. about varying operating points).
Parameters: q_ind, qd_ind : array_like, optional
The independent generalized coordinates and speeds.
q_dep, qd_dep : array_like, optional
The dependent generalized coordinates and speeds.
-