kane.py 28 KB

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