kane.py 36 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859
  1. from sympy import zeros, Matrix, diff, eye, linear_eq_to_matrix
  2. from sympy.core.sorting import default_sort_key
  3. from sympy.physics.vector import (ReferenceFrame, dynamicsymbols,
  4. partial_velocity)
  5. from sympy.physics.mechanics.method import _Methods
  6. from sympy.physics.mechanics.particle import Particle
  7. from sympy.physics.mechanics.rigidbody import RigidBody
  8. from sympy.physics.mechanics.functions import (msubs, find_dynamicsymbols,
  9. _f_list_parser,
  10. _validate_coordinates,
  11. _parse_linear_solver)
  12. from sympy.physics.mechanics.linearize import Linearizer
  13. from sympy.utilities.iterables import iterable
  14. __all__ = ['KanesMethod']
  15. class KanesMethod(_Methods):
  16. r"""Kane's method object.
  17. Explanation
  18. ===========
  19. This object is used to do the "book-keeping" as you go through and form
  20. equations of motion in the way Kane presents in:
  21. Kane, T., Levinson, D. Dynamics Theory and Applications. 1985 McGraw-Hill
  22. The attributes are for equations in the form [M] udot = forcing.
  23. Attributes
  24. ==========
  25. q, u : Matrix
  26. Matrices of the generalized coordinates and speeds
  27. bodies : iterable
  28. Iterable of Particle and RigidBody objects in the system.
  29. loads : iterable
  30. Iterable of (Point, vector) or (ReferenceFrame, vector) tuples
  31. describing the forces on the system.
  32. auxiliary_eqs : Matrix
  33. If applicable, the set of auxiliary Kane's
  34. equations used to solve for non-contributing
  35. forces.
  36. mass_matrix : Matrix
  37. The system's dynamics mass matrix: [k_d; k_dnh]
  38. forcing : Matrix
  39. The system's dynamics forcing vector: -[f_d; f_dnh]
  40. mass_matrix_kin : Matrix
  41. The "mass matrix" for kinematic differential equations: k_kqdot
  42. forcing_kin : Matrix
  43. The forcing vector for kinematic differential equations: -(k_ku*u + f_k)
  44. mass_matrix_full : Matrix
  45. The "mass matrix" for the u's and q's with dynamics and kinematics
  46. forcing_full : Matrix
  47. The "forcing vector" for the u's and q's with dynamics and kinematics
  48. Parameters
  49. ==========
  50. frame : ReferenceFrame
  51. The inertial reference frame for the system.
  52. q_ind : iterable of dynamicsymbols
  53. Independent generalized coordinates.
  54. u_ind : iterable of dynamicsymbols
  55. Independent generalized speeds.
  56. kd_eqs : iterable of Expr, optional
  57. Kinematic differential equations, which linearly relate the generalized
  58. speeds to the time-derivatives of the generalized coordinates.
  59. q_dependent : iterable of dynamicsymbols, optional
  60. Dependent generalized coordinates.
  61. configuration_constraints : iterable of Expr, optional
  62. Constraints on the system's configuration, i.e. holonomic constraints.
  63. u_dependent : iterable of dynamicsymbols, optional
  64. Dependent generalized speeds.
  65. velocity_constraints : iterable of Expr, optional
  66. Constraints on the system's velocity, i.e. the combination of the
  67. nonholonomic constraints and the time-derivative of the holonomic
  68. constraints.
  69. acceleration_constraints : iterable of Expr, optional
  70. Constraints on the system's acceleration, by default these are the
  71. time-derivative of the velocity constraints.
  72. u_auxiliary : iterable of dynamicsymbols, optional
  73. Auxiliary generalized speeds.
  74. bodies : iterable of Particle and/or RigidBody, optional
  75. The particles and rigid bodies in the system.
  76. forcelist : iterable of tuple[Point | ReferenceFrame, Vector], optional
  77. Forces and torques applied on the system.
  78. explicit_kinematics : bool
  79. Boolean whether the mass matrices and forcing vectors should use the
  80. explicit form (default) or implicit form for kinematics.
  81. See the notes for more details.
  82. kd_eqs_solver : str, callable
  83. Method used to solve the kinematic differential equations. If a string
  84. is supplied, it should be a valid method that can be used with the
  85. :meth:`sympy.matrices.matrixbase.MatrixBase.solve`. If a callable is
  86. supplied, it should have the format ``f(A, rhs)``, where it solves the
  87. equations and returns the solution. The default utilizes LU solve. See
  88. the notes for more information.
  89. constraint_solver : str, callable
  90. Method used to solve the velocity constraints. If a string is
  91. supplied, it should be a valid method that can be used with the
  92. :meth:`sympy.matrices.matrixbase.MatrixBase.solve`. If a callable is
  93. supplied, it should have the format ``f(A, rhs)``, where it solves the
  94. equations and returns the solution. The default utilizes LU solve. See
  95. the notes for more information.
  96. Notes
  97. =====
  98. The mass matrices and forcing vectors related to kinematic equations
  99. are given in the explicit form by default. In other words, the kinematic
  100. mass matrix is $\mathbf{k_{k\dot{q}}} = \mathbf{I}$.
  101. In order to get the implicit form of those matrices/vectors, you can set the
  102. ``explicit_kinematics`` attribute to ``False``. So $\mathbf{k_{k\dot{q}}}$
  103. is not necessarily an identity matrix. This can provide more compact
  104. equations for non-simple kinematics.
  105. Two linear solvers can be supplied to ``KanesMethod``: one for solving the
  106. kinematic differential equations and one to solve the velocity constraints.
  107. Both of these sets of equations can be expressed as a linear system ``Ax = rhs``,
  108. which have to be solved in order to obtain the equations of motion.
  109. The default solver ``'LU'``, which stands for LU solve, results relatively low
  110. number of operations. The weakness of this method is that it can result in zero
  111. division errors.
  112. If zero divisions are encountered, a possible solver which may solve the problem
  113. is ``"CRAMER"``. This method uses Cramer's rule to solve the system. This method
  114. is slower and results in more operations than the default solver. However it only
  115. uses a single division by default per entry of the solution.
  116. While a valid list of solvers can be found at
  117. :meth:`sympy.matrices.matrixbase.MatrixBase.solve`, it is also possible to supply a
  118. `callable`. This way it is possible to use a different solver routine. If the
  119. kinematic differential equations are not too complex it can be worth it to simplify
  120. the solution by using ``lambda A, b: simplify(Matrix.LUsolve(A, b))``. Another
  121. option solver one may use is :func:`sympy.solvers.solveset.linsolve`. This can be
  122. done using `lambda A, b: tuple(linsolve((A, b)))[0]`, where we select the first
  123. solution as our system should have only one unique solution.
  124. Examples
  125. ========
  126. This is a simple example for a one degree of freedom translational
  127. spring-mass-damper.
  128. In this example, we first need to do the kinematics.
  129. This involves creating generalized speeds and coordinates and their
  130. derivatives.
  131. Then we create a point and set its velocity in a frame.
  132. >>> from sympy import symbols
  133. >>> from sympy.physics.mechanics import dynamicsymbols, ReferenceFrame
  134. >>> from sympy.physics.mechanics import Point, Particle, KanesMethod
  135. >>> q, u = dynamicsymbols('q u')
  136. >>> qd, ud = dynamicsymbols('q u', 1)
  137. >>> m, c, k = symbols('m c k')
  138. >>> N = ReferenceFrame('N')
  139. >>> P = Point('P')
  140. >>> P.set_vel(N, u * N.x)
  141. Next we need to arrange/store information in the way that KanesMethod
  142. requires. The kinematic differential equations should be an iterable of
  143. expressions. A list of forces/torques must be constructed, where each entry
  144. in the list is a (Point, Vector) or (ReferenceFrame, Vector) tuple, where
  145. the Vectors represent the Force or Torque.
  146. Next a particle needs to be created, and it needs to have a point and mass
  147. assigned to it.
  148. Finally, a list of all bodies and particles needs to be created.
  149. >>> kd = [qd - u]
  150. >>> FL = [(P, (-k * q - c * u) * N.x)]
  151. >>> pa = Particle('pa', P, m)
  152. >>> BL = [pa]
  153. Finally we can generate the equations of motion.
  154. First we create the KanesMethod object and supply an inertial frame,
  155. coordinates, generalized speeds, and the kinematic differential equations.
  156. Additional quantities such as configuration and motion constraints,
  157. dependent coordinates and speeds, and auxiliary speeds are also supplied
  158. here (see the online documentation).
  159. Next we form FR* and FR to complete: Fr + Fr* = 0.
  160. We have the equations of motion at this point.
  161. It makes sense to rearrange them though, so we calculate the mass matrix and
  162. the forcing terms, for E.o.M. in the form: [MM] udot = forcing, where MM is
  163. the mass matrix, udot is a vector of the time derivatives of the
  164. generalized speeds, and forcing is a vector representing "forcing" terms.
  165. >>> KM = KanesMethod(N, q_ind=[q], u_ind=[u], kd_eqs=kd)
  166. >>> (fr, frstar) = KM.kanes_equations(BL, FL)
  167. >>> MM = KM.mass_matrix
  168. >>> forcing = KM.forcing
  169. >>> rhs = MM.inv() * forcing
  170. >>> rhs
  171. Matrix([[(-c*u(t) - k*q(t))/m]])
  172. >>> KM.linearize(A_and_B=True)[0]
  173. Matrix([
  174. [ 0, 1],
  175. [-k/m, -c/m]])
  176. Please look at the documentation pages for more information on how to
  177. perform linearization and how to deal with dependent coordinates & speeds,
  178. and how do deal with bringing non-contributing forces into evidence.
  179. """
  180. def __init__(self, frame, q_ind, u_ind, kd_eqs=None, q_dependent=None,
  181. configuration_constraints=None, u_dependent=None,
  182. velocity_constraints=None, acceleration_constraints=None,
  183. u_auxiliary=None, bodies=None, forcelist=None,
  184. explicit_kinematics=True, kd_eqs_solver='LU',
  185. constraint_solver='LU'):
  186. """Please read the online documentation. """
  187. if not q_ind:
  188. q_ind = [dynamicsymbols('dummy_q')]
  189. kd_eqs = [dynamicsymbols('dummy_kd')]
  190. if not isinstance(frame, ReferenceFrame):
  191. raise TypeError('An inertial ReferenceFrame must be supplied')
  192. self._inertial = frame
  193. self._fr = None
  194. self._frstar = None
  195. self._forcelist = forcelist
  196. self._bodylist = bodies
  197. self.explicit_kinematics = explicit_kinematics
  198. self._constraint_solver = constraint_solver
  199. self._initialize_vectors(q_ind, q_dependent, u_ind, u_dependent,
  200. u_auxiliary)
  201. _validate_coordinates(self.q, self.u)
  202. self._initialize_kindiffeq_matrices(kd_eqs, kd_eqs_solver)
  203. self._initialize_constraint_matrices(
  204. configuration_constraints, velocity_constraints,
  205. acceleration_constraints, constraint_solver)
  206. def _initialize_vectors(self, q_ind, q_dep, u_ind, u_dep, u_aux):
  207. """Initialize the coordinate and speed vectors."""
  208. none_handler = lambda x: Matrix(x) if x else Matrix()
  209. # Initialize generalized coordinates
  210. q_dep = none_handler(q_dep)
  211. if not iterable(q_ind):
  212. raise TypeError('Generalized coordinates must be an iterable.')
  213. if not iterable(q_dep):
  214. raise TypeError('Dependent coordinates must be an iterable.')
  215. q_ind = Matrix(q_ind)
  216. self._qdep = q_dep
  217. self._q = Matrix([q_ind, q_dep])
  218. self._qdot = self.q.diff(dynamicsymbols._t)
  219. # Initialize generalized speeds
  220. u_dep = none_handler(u_dep)
  221. if not iterable(u_ind):
  222. raise TypeError('Generalized speeds must be an iterable.')
  223. if not iterable(u_dep):
  224. raise TypeError('Dependent speeds must be an iterable.')
  225. u_ind = Matrix(u_ind)
  226. self._udep = u_dep
  227. self._u = Matrix([u_ind, u_dep])
  228. self._udot = self.u.diff(dynamicsymbols._t)
  229. self._uaux = none_handler(u_aux)
  230. def _initialize_constraint_matrices(self, config, vel, acc, linear_solver='LU'):
  231. """Initializes constraint matrices."""
  232. linear_solver = _parse_linear_solver(linear_solver)
  233. # Define vector dimensions
  234. o = len(self.u)
  235. m = len(self._udep)
  236. p = o - m
  237. none_handler = lambda x: Matrix(x) if x else Matrix()
  238. # Initialize configuration constraints
  239. config = none_handler(config)
  240. if len(self._qdep) != len(config):
  241. raise ValueError('There must be an equal number of dependent '
  242. 'coordinates and configuration constraints.')
  243. self._f_h = none_handler(config)
  244. # Initialize velocity and acceleration constraints
  245. vel = none_handler(vel)
  246. acc = none_handler(acc)
  247. if len(vel) != m:
  248. raise ValueError('There must be an equal number of dependent '
  249. 'speeds and velocity constraints.')
  250. if acc and (len(acc) != m):
  251. raise ValueError('There must be an equal number of dependent '
  252. 'speeds and acceleration constraints.')
  253. if vel:
  254. # When calling kanes_equations, another class instance will be
  255. # created if auxiliary u's are present. In this case, the
  256. # computation of kinetic differential equation matrices will be
  257. # skipped as this was computed during the original KanesMethod
  258. # object, and the qd_u_map will not be available.
  259. if self._qdot_u_map is not None:
  260. vel = msubs(vel, self._qdot_u_map)
  261. self._k_nh, f_nh_neg = linear_eq_to_matrix(vel, self.u[:])
  262. self._f_nh = -f_nh_neg
  263. # If no acceleration constraints given, calculate them.
  264. if not acc:
  265. _f_dnh = (self._k_nh.diff(dynamicsymbols._t) * self.u +
  266. self._f_nh.diff(dynamicsymbols._t))
  267. if self._qdot_u_map is not None:
  268. _f_dnh = msubs(_f_dnh, self._qdot_u_map)
  269. self._f_dnh = _f_dnh
  270. self._k_dnh = self._k_nh
  271. else:
  272. if self._qdot_u_map is not None:
  273. acc = msubs(acc, self._qdot_u_map)
  274. self._k_dnh, f_dnh_neg = linear_eq_to_matrix(acc, self._udot[:])
  275. self._f_dnh = -f_dnh_neg
  276. # Form of non-holonomic constraints is B*u + C = 0.
  277. # We partition B into independent and dependent columns:
  278. # Ars is then -B_dep.inv() * B_ind, and it relates dependent speeds
  279. # to independent speeds as: udep = Ars*uind, neglecting the C term.
  280. B_ind = self._k_nh[:, :p]
  281. B_dep = self._k_nh[:, p:o]
  282. self._Ars = -linear_solver(B_dep, B_ind)
  283. else:
  284. self._f_nh = Matrix()
  285. self._k_nh = Matrix()
  286. self._f_dnh = Matrix()
  287. self._k_dnh = Matrix()
  288. self._Ars = Matrix()
  289. def _initialize_kindiffeq_matrices(self, kdeqs, linear_solver='LU'):
  290. """Initialize the kinematic differential equation matrices.
  291. Parameters
  292. ==========
  293. kdeqs : sequence of sympy expressions
  294. Kinematic differential equations in the form of f(u,q',q,t) where
  295. f() = 0. The equations have to be linear in the time-derivatives of
  296. the generalized coordinates and in the generalized speeds.
  297. """
  298. linear_solver = _parse_linear_solver(linear_solver)
  299. if kdeqs:
  300. if len(self.q) != len(kdeqs):
  301. raise ValueError('There must be an equal number of kinematic '
  302. 'differential equations and coordinates.')
  303. u = self.u
  304. qdot = self._qdot
  305. kdeqs = Matrix(kdeqs)
  306. u_zero = dict.fromkeys(u, 0)
  307. uaux_zero = dict.fromkeys(self._uaux, 0)
  308. qdot_zero = dict.fromkeys(qdot, 0)
  309. # Extract the linear coefficient matrices as per the following
  310. # equation:
  311. #
  312. # k_ku(q,t)*u(t) + k_kqdot(q,t)*q'(t) + f_k(q,t) = 0
  313. #
  314. k_ku = kdeqs.jacobian(u)
  315. k_kqdot = kdeqs.jacobian(qdot)
  316. f_k = kdeqs.xreplace(u_zero).xreplace(qdot_zero)
  317. # The kinematic differential equations should be linear in both q'
  318. # and u so check for u and q' in the components.
  319. dy_syms = find_dynamicsymbols(k_ku.row_join(k_kqdot).row_join(f_k))
  320. nonlin_vars = [vari for vari in u[:] + qdot[:] if vari in dy_syms]
  321. if nonlin_vars:
  322. msg = ('The provided kinematic differential equations are '
  323. 'nonlinear in {}. They must be linear in the '
  324. 'generalized speeds and derivatives of the generalized '
  325. 'coordinates.')
  326. raise ValueError(msg.format(nonlin_vars))
  327. self._f_k_implicit = f_k.xreplace(uaux_zero)
  328. self._k_ku_implicit = k_ku.xreplace(uaux_zero)
  329. self._k_kqdot_implicit = k_kqdot
  330. # Solve for q'(t) such that the coefficient matrices are now in
  331. # this form:
  332. #
  333. # k_kqdot^-1*k_ku*u(t) + I*q'(t) + k_kqdot^-1*f_k = 0
  334. #
  335. # NOTE : Solving the kinematic differential equations here is not
  336. # necessary and prevents the equations from being provided in fully
  337. # implicit form.
  338. f_k_explicit = linear_solver(k_kqdot, f_k)
  339. k_ku_explicit = linear_solver(k_kqdot, k_ku)
  340. self._qdot_u_map = dict(zip(qdot, -(k_ku_explicit*u + f_k_explicit)))
  341. self._f_k = f_k_explicit.xreplace(uaux_zero)
  342. self._k_ku = k_ku_explicit.xreplace(uaux_zero)
  343. self._k_kqdot = eye(len(qdot))
  344. else:
  345. self._qdot_u_map = None
  346. self._f_k_implicit = self._f_k = Matrix()
  347. self._k_ku_implicit = self._k_ku = Matrix()
  348. self._k_kqdot_implicit = self._k_kqdot = Matrix()
  349. def _form_fr(self, fl):
  350. """Form the generalized active force."""
  351. if fl is not None and (len(fl) == 0 or not iterable(fl)):
  352. raise ValueError('Force pairs must be supplied in an '
  353. 'non-empty iterable or None.')
  354. N = self._inertial
  355. # pull out relevant velocities for constructing partial velocities
  356. vel_list, f_list = _f_list_parser(fl, N)
  357. vel_list = [msubs(i, self._qdot_u_map) for i in vel_list]
  358. f_list = [msubs(i, self._qdot_u_map) for i in f_list]
  359. # Fill Fr with dot product of partial velocities and forces
  360. o = len(self.u)
  361. b = len(f_list)
  362. FR = zeros(o, 1)
  363. partials = partial_velocity(vel_list, self.u, N)
  364. for i in range(o):
  365. FR[i] = sum(partials[j][i].dot(f_list[j]) for j in range(b))
  366. # In case there are dependent speeds
  367. if self._udep:
  368. p = o - len(self._udep)
  369. FRtilde = FR[:p, 0]
  370. FRold = FR[p:o, 0]
  371. FRtilde += self._Ars.T * FRold
  372. FR = FRtilde
  373. self._forcelist = fl
  374. self._fr = FR
  375. return FR
  376. def _form_frstar(self, bl):
  377. """Form the generalized inertia force."""
  378. if not iterable(bl):
  379. raise TypeError('Bodies must be supplied in an iterable.')
  380. t = dynamicsymbols._t
  381. N = self._inertial
  382. # Dicts setting things to zero
  383. udot_zero = dict.fromkeys(self._udot, 0)
  384. uaux_zero = dict.fromkeys(self._uaux, 0)
  385. uauxdot = [diff(i, t) for i in self._uaux]
  386. uauxdot_zero = dict.fromkeys(uauxdot, 0)
  387. # Dictionary of q' and q'' to u and u'
  388. q_ddot_u_map = {k.diff(t): v.diff(t).xreplace(
  389. self._qdot_u_map) for (k, v) in self._qdot_u_map.items()}
  390. q_ddot_u_map.update(self._qdot_u_map)
  391. # Fill up the list of partials: format is a list with num elements
  392. # equal to number of entries in body list. Each of these elements is a
  393. # list - either of length 1 for the translational components of
  394. # particles or of length 2 for the translational and rotational
  395. # components of rigid bodies. The inner most list is the list of
  396. # partial velocities.
  397. def get_partial_velocity(body):
  398. if isinstance(body, RigidBody):
  399. vlist = [body.masscenter.vel(N), body.frame.ang_vel_in(N)]
  400. elif isinstance(body, Particle):
  401. vlist = [body.point.vel(N),]
  402. else:
  403. raise TypeError('The body list may only contain either '
  404. 'RigidBody or Particle as list elements.')
  405. v = [msubs(vel, self._qdot_u_map) for vel in vlist]
  406. return partial_velocity(v, self.u, N)
  407. partials = [get_partial_velocity(body) for body in bl]
  408. # Compute fr_star in two components:
  409. # fr_star = -(MM*u' + nonMM)
  410. o = len(self.u)
  411. MM = zeros(o, o)
  412. nonMM = zeros(o, 1)
  413. zero_uaux = lambda expr: msubs(expr, uaux_zero)
  414. zero_udot_uaux = lambda expr: msubs(msubs(expr, udot_zero), uaux_zero)
  415. for i, body in enumerate(bl):
  416. if isinstance(body, RigidBody):
  417. M = zero_uaux(body.mass)
  418. I = zero_uaux(body.central_inertia)
  419. vel = zero_uaux(body.masscenter.vel(N))
  420. omega = zero_uaux(body.frame.ang_vel_in(N))
  421. acc = zero_udot_uaux(body.masscenter.acc(N))
  422. inertial_force = (M.diff(t) * vel + M * acc)
  423. inertial_torque = zero_uaux((I.dt(body.frame).dot(omega)) +
  424. msubs(I.dot(body.frame.ang_acc_in(N)), udot_zero) +
  425. (omega.cross(I.dot(omega))))
  426. for j in range(o):
  427. tmp_vel = zero_uaux(partials[i][0][j])
  428. tmp_ang = zero_uaux(I.dot(partials[i][1][j]))
  429. for k in range(o):
  430. # translational
  431. MM[j, k] += M*tmp_vel.dot(partials[i][0][k])
  432. # rotational
  433. MM[j, k] += tmp_ang.dot(partials[i][1][k])
  434. nonMM[j] += inertial_force.dot(partials[i][0][j])
  435. nonMM[j] += inertial_torque.dot(partials[i][1][j])
  436. else:
  437. M = zero_uaux(body.mass)
  438. vel = zero_uaux(body.point.vel(N))
  439. acc = zero_udot_uaux(body.point.acc(N))
  440. inertial_force = (M.diff(t) * vel + M * acc)
  441. for j in range(o):
  442. temp = zero_uaux(partials[i][0][j])
  443. for k in range(o):
  444. MM[j, k] += M*temp.dot(partials[i][0][k])
  445. nonMM[j] += inertial_force.dot(partials[i][0][j])
  446. # Compose fr_star out of MM and nonMM
  447. MM = zero_uaux(msubs(MM, q_ddot_u_map))
  448. nonMM = msubs(msubs(nonMM, q_ddot_u_map),
  449. udot_zero, uauxdot_zero, uaux_zero)
  450. fr_star = -(MM * msubs(Matrix(self._udot), uauxdot_zero) + nonMM)
  451. # If there are dependent speeds, we need to find fr_star_tilde
  452. if self._udep:
  453. p = o - len(self._udep)
  454. fr_star_ind = fr_star[:p, 0]
  455. fr_star_dep = fr_star[p:o, 0]
  456. fr_star = fr_star_ind + (self._Ars.T * fr_star_dep)
  457. # Apply the same to MM
  458. MMi = MM[:p, :]
  459. MMd = MM[p:o, :]
  460. MM = MMi + (self._Ars.T * MMd)
  461. # Apply the same to nonMM
  462. nonMM = nonMM[:p, :] + (self._Ars.T * nonMM[p:o, :])
  463. self._bodylist = bl
  464. self._frstar = fr_star
  465. self._k_d = MM
  466. self._f_d = -(self._fr - nonMM)
  467. return fr_star
  468. def to_linearizer(self, linear_solver='LU'):
  469. """Returns an instance of the Linearizer class, initiated from the
  470. data in the KanesMethod class. This may be more desirable than using
  471. the linearize class method, as the Linearizer object will allow more
  472. efficient recalculation (i.e. about varying operating points).
  473. Parameters
  474. ==========
  475. linear_solver : str, callable
  476. Method used to solve the several symbolic linear systems of the
  477. form ``A*x=b`` in the linearization process. If a string is
  478. supplied, it should be a valid method that can be used with the
  479. :meth:`sympy.matrices.matrixbase.MatrixBase.solve`. If a callable is
  480. supplied, it should have the format ``x = f(A, b)``, where it
  481. solves the equations and returns the solution. The default is
  482. ``'LU'`` which corresponds to SymPy's ``A.LUsolve(b)``.
  483. ``LUsolve()`` is fast to compute but will often result in
  484. divide-by-zero and thus ``nan`` results.
  485. Returns
  486. =======
  487. Linearizer
  488. An instantiated
  489. :class:`sympy.physics.mechanics.linearize.Linearizer`.
  490. """
  491. if (self._fr is None) or (self._frstar is None):
  492. raise ValueError('Need to compute Fr, Fr* first.')
  493. # Get required equation components. The Kane's method class breaks
  494. # these into pieces. Need to reassemble
  495. f_c = self._f_h
  496. if self._f_nh and self._k_nh:
  497. f_v = self._f_nh + self._k_nh*Matrix(self.u)
  498. else:
  499. f_v = Matrix()
  500. if self._f_dnh and self._k_dnh:
  501. f_a = self._f_dnh + self._k_dnh*Matrix(self._udot)
  502. else:
  503. f_a = Matrix()
  504. # Dicts to sub to zero, for splitting up expressions
  505. u_zero = dict.fromkeys(self.u, 0)
  506. ud_zero = dict.fromkeys(self._udot, 0)
  507. qd_zero = dict.fromkeys(self._qdot, 0)
  508. qd_u_zero = dict.fromkeys(Matrix([self._qdot, self.u]), 0)
  509. # Break the kinematic differential eqs apart into f_0 and f_1
  510. f_0 = msubs(self._f_k, u_zero) + self._k_kqdot*Matrix(self._qdot)
  511. f_1 = msubs(self._f_k, qd_zero) + self._k_ku*Matrix(self.u)
  512. # Break the dynamic differential eqs into f_2 and f_3
  513. f_2 = msubs(self._frstar, qd_u_zero)
  514. f_3 = msubs(self._frstar, ud_zero) + self._fr
  515. f_4 = zeros(len(f_2), 1)
  516. # Get the required vector components
  517. q = self.q
  518. u = self.u
  519. if self._qdep:
  520. q_i = q[:-len(self._qdep)]
  521. else:
  522. q_i = q
  523. q_d = self._qdep
  524. if self._udep:
  525. u_i = u[:-len(self._udep)]
  526. else:
  527. u_i = u
  528. u_d = self._udep
  529. # Form dictionary to set auxiliary speeds & their derivatives to 0.
  530. uaux = self._uaux
  531. uauxdot = uaux.diff(dynamicsymbols._t)
  532. uaux_zero = dict.fromkeys(Matrix([uaux, uauxdot]), 0)
  533. # Checking for dynamic symbols outside the dynamic differential
  534. # equations; throws error if there is.
  535. sym_list = set(Matrix([q, self._qdot, u, self._udot, uaux, uauxdot]))
  536. if any(find_dynamicsymbols(i, sym_list) for i in [self._k_kqdot,
  537. self._k_ku, self._f_k, self._k_dnh, self._f_dnh, self._k_d]):
  538. raise ValueError('Cannot have dynamicsymbols outside dynamic \
  539. forcing vector.')
  540. # Find all other dynamic symbols, forming the forcing vector r.
  541. # Sort r to make it canonical.
  542. r = list(find_dynamicsymbols(msubs(self._f_d, uaux_zero), sym_list))
  543. r.sort(key=default_sort_key)
  544. # Check for any derivatives of variables in r that are also found in r.
  545. for i in r:
  546. if diff(i, dynamicsymbols._t) in r:
  547. raise ValueError('Cannot have derivatives of specified \
  548. quantities when linearizing forcing terms.')
  549. return Linearizer(f_0, f_1, f_2, f_3, f_4, f_c, f_v, f_a, q, u, q_i,
  550. q_d, u_i, u_d, r, linear_solver=linear_solver)
  551. # TODO : Remove `new_method` after 1.1 has been released.
  552. def linearize(self, *, new_method=None, linear_solver='LU', **kwargs):
  553. """ Linearize the equations of motion about a symbolic operating point.
  554. Parameters
  555. ==========
  556. new_method
  557. Deprecated, does nothing and will be removed.
  558. linear_solver : str, callable
  559. Method used to solve the several symbolic linear systems of the
  560. form ``A*x=b`` in the linearization process. If a string is
  561. supplied, it should be a valid method that can be used with the
  562. :meth:`sympy.matrices.matrixbase.MatrixBase.solve`. If a callable is
  563. supplied, it should have the format ``x = f(A, b)``, where it
  564. solves the equations and returns the solution. The default is
  565. ``'LU'`` which corresponds to SymPy's ``A.LUsolve(b)``.
  566. ``LUsolve()`` is fast to compute but will often result in
  567. divide-by-zero and thus ``nan`` results.
  568. **kwargs
  569. Extra keyword arguments are passed to
  570. :meth:`sympy.physics.mechanics.linearize.Linearizer.linearize`.
  571. Explanation
  572. ===========
  573. If kwarg A_and_B is False (default), returns M, A, B, r for the
  574. linearized form, M*[q', u']^T = A*[q_ind, u_ind]^T + B*r.
  575. If kwarg A_and_B is True, returns A, B, r for the linearized form
  576. dx = A*x + B*r, where x = [q_ind, u_ind]^T. Note that this is
  577. computationally intensive if there are many symbolic parameters. For
  578. this reason, it may be more desirable to use the default A_and_B=False,
  579. returning M, A, and B. Values may then be substituted in to these
  580. matrices, and the state space form found as
  581. A = P.T*M.inv()*A, B = P.T*M.inv()*B, where P = Linearizer.perm_mat.
  582. In both cases, r is found as all dynamicsymbols in the equations of
  583. motion that are not part of q, u, q', or u'. They are sorted in
  584. canonical form.
  585. The operating points may be also entered using the ``op_point`` kwarg.
  586. This takes a dictionary of {symbol: value}, or a an iterable of such
  587. dictionaries. The values may be numeric or symbolic. The more values
  588. you can specify beforehand, the faster this computation will run.
  589. For more documentation, please see the ``Linearizer`` class.
  590. """
  591. linearizer = self.to_linearizer(linear_solver=linear_solver)
  592. result = linearizer.linearize(**kwargs)
  593. return result + (linearizer.r,)
  594. def kanes_equations(self, bodies=None, loads=None):
  595. """ Method to form Kane's equations, Fr + Fr* = 0.
  596. Explanation
  597. ===========
  598. Returns (Fr, Fr*). In the case where auxiliary generalized speeds are
  599. present (say, s auxiliary speeds, o generalized speeds, and m motion
  600. constraints) the length of the returned vectors will be o - m + s in
  601. length. The first o - m equations will be the constrained Kane's
  602. equations, then the s auxiliary Kane's equations. These auxiliary
  603. equations can be accessed with the auxiliary_eqs property.
  604. Parameters
  605. ==========
  606. bodies : iterable
  607. An iterable of all RigidBody's and Particle's in the system.
  608. A system must have at least one body.
  609. loads : iterable
  610. Takes in an iterable of (Particle, Vector) or (ReferenceFrame, Vector)
  611. tuples which represent the force at a point or torque on a frame.
  612. Must be either a non-empty iterable of tuples or None which corresponds
  613. to a system with no constraints.
  614. """
  615. if bodies is None:
  616. bodies = self.bodies
  617. if loads is None and self._forcelist is not None:
  618. loads = self._forcelist
  619. if loads == []:
  620. loads = None
  621. if not self._k_kqdot:
  622. raise AttributeError('Create an instance of KanesMethod with '
  623. 'kinematic differential equations to use this method.')
  624. fr = self._form_fr(loads)
  625. frstar = self._form_frstar(bodies)
  626. if self._uaux:
  627. if not self._udep:
  628. km = KanesMethod(self._inertial, self.q, self._uaux,
  629. u_auxiliary=self._uaux, constraint_solver=self._constraint_solver)
  630. else:
  631. km = KanesMethod(self._inertial, self.q, self._uaux,
  632. u_auxiliary=self._uaux, u_dependent=self._udep,
  633. velocity_constraints=(self._k_nh * self.u +
  634. self._f_nh),
  635. acceleration_constraints=(self._k_dnh * self._udot +
  636. self._f_dnh),
  637. constraint_solver=self._constraint_solver
  638. )
  639. km._qdot_u_map = self._qdot_u_map
  640. self._km = km
  641. fraux = km._form_fr(loads)
  642. frstaraux = km._form_frstar(bodies)
  643. self._aux_eq = fraux + frstaraux
  644. self._fr = fr.col_join(fraux)
  645. self._frstar = frstar.col_join(frstaraux)
  646. return (self._fr, self._frstar)
  647. def _form_eoms(self):
  648. fr, frstar = self.kanes_equations(self.bodylist, self.forcelist)
  649. return fr + frstar
  650. def rhs(self, inv_method=None):
  651. """Returns the system's equations of motion in first order form. The
  652. output is the right hand side of::
  653. x' = |q'| =: f(q, u, r, p, t)
  654. |u'|
  655. The right hand side is what is needed by most numerical ODE
  656. integrators.
  657. Parameters
  658. ==========
  659. inv_method : str
  660. The specific sympy inverse matrix calculation method to use. For a
  661. list of valid methods, see
  662. :meth:`~sympy.matrices.matrixbase.MatrixBase.inv`
  663. """
  664. rhs = zeros(len(self.q) + len(self.u), 1)
  665. kdes = self.kindiffdict()
  666. for i, q_i in enumerate(self.q):
  667. rhs[i] = kdes[q_i.diff()]
  668. if inv_method is None:
  669. rhs[len(self.q):, 0] = self.mass_matrix.LUsolve(self.forcing)
  670. else:
  671. rhs[len(self.q):, 0] = (self.mass_matrix.inv(inv_method,
  672. try_block_diag=True) *
  673. self.forcing)
  674. return rhs
  675. def kindiffdict(self):
  676. """Returns a dictionary mapping q' to u."""
  677. if not self._qdot_u_map:
  678. raise AttributeError('Create an instance of KanesMethod with '
  679. 'kinematic differential equations to use this method.')
  680. return self._qdot_u_map
  681. @property
  682. def auxiliary_eqs(self):
  683. """A matrix containing the auxiliary equations."""
  684. if not self._fr or not self._frstar:
  685. raise ValueError('Need to compute Fr, Fr* first.')
  686. if not self._uaux:
  687. raise ValueError('No auxiliary speeds have been declared.')
  688. return self._aux_eq
  689. @property
  690. def mass_matrix_kin(self):
  691. r"""The kinematic "mass matrix" $\mathbf{k_{k\dot{q}}}$ of the system."""
  692. return self._k_kqdot if self.explicit_kinematics else self._k_kqdot_implicit
  693. @property
  694. def forcing_kin(self):
  695. """The kinematic "forcing vector" of the system."""
  696. if self.explicit_kinematics:
  697. return -(self._k_ku * Matrix(self.u) + self._f_k)
  698. else:
  699. return -(self._k_ku_implicit * Matrix(self.u) + self._f_k_implicit)
  700. @property
  701. def mass_matrix(self):
  702. """The mass matrix of the system."""
  703. if not self._fr or not self._frstar:
  704. raise ValueError('Need to compute Fr, Fr* first.')
  705. return Matrix([self._k_d, self._k_dnh])
  706. @property
  707. def forcing(self):
  708. """The forcing vector of the system."""
  709. if not self._fr or not self._frstar:
  710. raise ValueError('Need to compute Fr, Fr* first.')
  711. return -Matrix([self._f_d, self._f_dnh])
  712. @property
  713. def mass_matrix_full(self):
  714. """The mass matrix of the system, augmented by the kinematic
  715. differential equations in explicit or implicit form."""
  716. if not self._fr or not self._frstar:
  717. raise ValueError('Need to compute Fr, Fr* first.')
  718. o, n = len(self.u), len(self.q)
  719. return (self.mass_matrix_kin.row_join(zeros(n, o))).col_join(
  720. zeros(o, n).row_join(self.mass_matrix))
  721. @property
  722. def forcing_full(self):
  723. """The forcing vector of the system, augmented by the kinematic
  724. differential equations in explicit or implicit form."""
  725. return Matrix([self.forcing_kin, self.forcing])
  726. @property
  727. def q(self):
  728. return self._q
  729. @property
  730. def u(self):
  731. return self._u
  732. @property
  733. def bodylist(self):
  734. return self._bodylist
  735. @property
  736. def forcelist(self):
  737. return self._forcelist
  738. @property
  739. def bodies(self):
  740. return self._bodylist
  741. @property
  742. def loads(self):
  743. return self._forcelist