test_kane3.py 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315
  1. from sympy.core.numbers import pi
  2. from sympy.core.symbol import symbols
  3. from sympy.functions.elementary.miscellaneous import sqrt
  4. from sympy.functions.elementary.trigonometric import acos, sin, cos
  5. from sympy.matrices.dense import Matrix
  6. from sympy.physics.mechanics import (ReferenceFrame, dynamicsymbols,
  7. KanesMethod, inertia, Point, RigidBody,
  8. dot)
  9. from sympy.testing.pytest import slow
  10. @slow
  11. def test_bicycle():
  12. # Code to get equations of motion for a bicycle modeled as in:
  13. # J.P Meijaard, Jim M Papadopoulos, Andy Ruina and A.L Schwab. Linearized
  14. # dynamics equations for the balance and steer of a bicycle: a benchmark
  15. # and review. Proceedings of The Royal Society (2007) 463, 1955-1982
  16. # doi: 10.1098/rspa.2007.1857
  17. # Note that this code has been crudely ported from Autolev, which is the
  18. # reason for some of the unusual naming conventions. It was purposefully as
  19. # similar as possible in order to aide debugging.
  20. # Declare Coordinates & Speeds
  21. # Simple definitions for qdots - qd = u
  22. # Speeds are:
  23. # - u1: yaw frame ang. rate
  24. # - u2: roll frame ang. rate
  25. # - u3: rear wheel frame ang. rate (spinning motion)
  26. # - u4: frame ang. rate (pitching motion)
  27. # - u5: steering frame ang. rate
  28. # - u6: front wheel ang. rate (spinning motion)
  29. # Wheel positions are ignorable coordinates, so they are not introduced.
  30. q1, q2, q4, q5 = dynamicsymbols('q1 q2 q4 q5')
  31. q1d, q2d, q4d, q5d = dynamicsymbols('q1 q2 q4 q5', 1)
  32. u1, u2, u3, u4, u5, u6 = dynamicsymbols('u1 u2 u3 u4 u5 u6')
  33. u1d, u2d, u3d, u4d, u5d, u6d = dynamicsymbols('u1 u2 u3 u4 u5 u6', 1)
  34. # Declare System's Parameters
  35. WFrad, WRrad, htangle, forkoffset = symbols('WFrad WRrad htangle forkoffset')
  36. forklength, framelength, forkcg1 = symbols('forklength framelength forkcg1')
  37. forkcg3, framecg1, framecg3, Iwr11 = symbols('forkcg3 framecg1 framecg3 Iwr11')
  38. Iwr22, Iwf11, Iwf22, Iframe11 = symbols('Iwr22 Iwf11 Iwf22 Iframe11')
  39. Iframe22, Iframe33, Iframe31, Ifork11 = symbols('Iframe22 Iframe33 Iframe31 Ifork11')
  40. Ifork22, Ifork33, Ifork31, g = symbols('Ifork22 Ifork33 Ifork31 g')
  41. mframe, mfork, mwf, mwr = symbols('mframe mfork mwf mwr')
  42. # Set up reference frames for the system
  43. # N - inertial
  44. # Y - yaw
  45. # R - roll
  46. # WR - rear wheel, rotation angle is ignorable coordinate so not oriented
  47. # Frame - bicycle frame
  48. # TempFrame - statically rotated frame for easier reference inertia definition
  49. # Fork - bicycle fork
  50. # TempFork - statically rotated frame for easier reference inertia definition
  51. # WF - front wheel, again posses a ignorable coordinate
  52. N = ReferenceFrame('N')
  53. Y = N.orientnew('Y', 'Axis', [q1, N.z])
  54. R = Y.orientnew('R', 'Axis', [q2, Y.x])
  55. Frame = R.orientnew('Frame', 'Axis', [q4 + htangle, R.y])
  56. WR = ReferenceFrame('WR')
  57. TempFrame = Frame.orientnew('TempFrame', 'Axis', [-htangle, Frame.y])
  58. Fork = Frame.orientnew('Fork', 'Axis', [q5, Frame.x])
  59. TempFork = Fork.orientnew('TempFork', 'Axis', [-htangle, Fork.y])
  60. WF = ReferenceFrame('WF')
  61. # Kinematics of the Bicycle First block of code is forming the positions of
  62. # the relevant points
  63. # rear wheel contact -> rear wheel mass center -> frame mass center +
  64. # frame/fork connection -> fork mass center + front wheel mass center ->
  65. # front wheel contact point
  66. WR_cont = Point('WR_cont')
  67. WR_mc = WR_cont.locatenew('WR_mc', WRrad * R.z)
  68. Steer = WR_mc.locatenew('Steer', framelength * Frame.z)
  69. Frame_mc = WR_mc.locatenew('Frame_mc', - framecg1 * Frame.x
  70. + framecg3 * Frame.z)
  71. Fork_mc = Steer.locatenew('Fork_mc', - forkcg1 * Fork.x
  72. + forkcg3 * Fork.z)
  73. WF_mc = Steer.locatenew('WF_mc', forklength * Fork.x + forkoffset * Fork.z)
  74. WF_cont = WF_mc.locatenew('WF_cont', WFrad * (dot(Fork.y, Y.z) * Fork.y -
  75. Y.z).normalize())
  76. # Set the angular velocity of each frame.
  77. # Angular accelerations end up being calculated automatically by
  78. # differentiating the angular velocities when first needed.
  79. # u1 is yaw rate
  80. # u2 is roll rate
  81. # u3 is rear wheel rate
  82. # u4 is frame pitch rate
  83. # u5 is fork steer rate
  84. # u6 is front wheel rate
  85. Y.set_ang_vel(N, u1 * Y.z)
  86. R.set_ang_vel(Y, u2 * R.x)
  87. WR.set_ang_vel(Frame, u3 * Frame.y)
  88. Frame.set_ang_vel(R, u4 * Frame.y)
  89. Fork.set_ang_vel(Frame, u5 * Fork.x)
  90. WF.set_ang_vel(Fork, u6 * Fork.y)
  91. # Form the velocities of the previously defined points, using the 2 - point
  92. # theorem (written out by hand here). Accelerations again are calculated
  93. # automatically when first needed.
  94. WR_cont.set_vel(N, 0)
  95. WR_mc.v2pt_theory(WR_cont, N, WR)
  96. Steer.v2pt_theory(WR_mc, N, Frame)
  97. Frame_mc.v2pt_theory(WR_mc, N, Frame)
  98. Fork_mc.v2pt_theory(Steer, N, Fork)
  99. WF_mc.v2pt_theory(Steer, N, Fork)
  100. WF_cont.v2pt_theory(WF_mc, N, WF)
  101. # Sets the inertias of each body. Uses the inertia frame to construct the
  102. # inertia dyadics. Wheel inertias are only defined by principle moments of
  103. # inertia, and are in fact constant in the frame and fork reference frames;
  104. # it is for this reason that the orientations of the wheels does not need
  105. # to be defined. The frame and fork inertias are defined in the 'Temp'
  106. # frames which are fixed to the appropriate body frames; this is to allow
  107. # easier input of the reference values of the benchmark paper. Note that
  108. # due to slightly different orientations, the products of inertia need to
  109. # have their signs flipped; this is done later when entering the numerical
  110. # value.
  111. Frame_I = (inertia(TempFrame, Iframe11, Iframe22, Iframe33, 0, 0, Iframe31), Frame_mc)
  112. Fork_I = (inertia(TempFork, Ifork11, Ifork22, Ifork33, 0, 0, Ifork31), Fork_mc)
  113. WR_I = (inertia(Frame, Iwr11, Iwr22, Iwr11), WR_mc)
  114. WF_I = (inertia(Fork, Iwf11, Iwf22, Iwf11), WF_mc)
  115. # Declaration of the RigidBody containers. ::
  116. BodyFrame = RigidBody('BodyFrame', Frame_mc, Frame, mframe, Frame_I)
  117. BodyFork = RigidBody('BodyFork', Fork_mc, Fork, mfork, Fork_I)
  118. BodyWR = RigidBody('BodyWR', WR_mc, WR, mwr, WR_I)
  119. BodyWF = RigidBody('BodyWF', WF_mc, WF, mwf, WF_I)
  120. # The kinematic differential equations; they are defined quite simply. Each
  121. # entry in this list is equal to zero.
  122. kd = [q1d - u1, q2d - u2, q4d - u4, q5d - u5]
  123. # The nonholonomic constraints are the velocity of the front wheel contact
  124. # point dotted into the X, Y, and Z directions; the yaw frame is used as it
  125. # is "closer" to the front wheel (1 less DCM connecting them). These
  126. # constraints force the velocity of the front wheel contact point to be 0
  127. # in the inertial frame; the X and Y direction constraints enforce a
  128. # "no-slip" condition, and the Z direction constraint forces the front
  129. # wheel contact point to not move away from the ground frame, essentially
  130. # replicating the holonomic constraint which does not allow the frame pitch
  131. # to change in an invalid fashion.
  132. conlist_speed = [WF_cont.vel(N) & Y.x, WF_cont.vel(N) & Y.y, WF_cont.vel(N) & Y.z]
  133. # The holonomic constraint is that the position from the rear wheel contact
  134. # point to the front wheel contact point when dotted into the
  135. # normal-to-ground plane direction must be zero; effectively that the front
  136. # and rear wheel contact points are always touching the ground plane. This
  137. # is actually not part of the dynamic equations, but instead is necessary
  138. # for the lineraization process.
  139. conlist_coord = [WF_cont.pos_from(WR_cont) & Y.z]
  140. # The force list; each body has the appropriate gravitational force applied
  141. # at its mass center.
  142. FL = [(Frame_mc, -mframe * g * Y.z),
  143. (Fork_mc, -mfork * g * Y.z),
  144. (WF_mc, -mwf * g * Y.z),
  145. (WR_mc, -mwr * g * Y.z)]
  146. BL = [BodyFrame, BodyFork, BodyWR, BodyWF]
  147. # The N frame is the inertial frame, coordinates are supplied in the order
  148. # of independent, dependent coordinates, as are the speeds. The kinematic
  149. # differential equation are also entered here. Here the dependent speeds
  150. # are specified, in the same order they were provided in earlier, along
  151. # with the non-holonomic constraints. The dependent coordinate is also
  152. # provided, with the holonomic constraint. Again, this is only provided
  153. # for the linearization process.
  154. KM = KanesMethod(N, q_ind=[q1, q2, q5],
  155. q_dependent=[q4], configuration_constraints=conlist_coord,
  156. u_ind=[u2, u3, u5],
  157. u_dependent=[u1, u4, u6], velocity_constraints=conlist_speed,
  158. kd_eqs=kd,
  159. constraint_solver="CRAMER")
  160. (fr, frstar) = KM.kanes_equations(BL, FL)
  161. # This is the start of entering in the numerical values from the benchmark
  162. # paper to validate the eigen values of the linearized equations from this
  163. # model to the reference eigen values. Look at the aforementioned paper for
  164. # more information. Some of these are intermediate values, used to
  165. # transform values from the paper into the coordinate systems used in this
  166. # model.
  167. PaperRadRear = 0.3
  168. PaperRadFront = 0.35
  169. HTA = (pi / 2 - pi / 10).evalf()
  170. TrailPaper = 0.08
  171. rake = (-(TrailPaper*sin(HTA)-(PaperRadFront*cos(HTA)))).evalf()
  172. PaperWb = 1.02
  173. PaperFrameCgX = 0.3
  174. PaperFrameCgZ = 0.9
  175. PaperForkCgX = 0.9
  176. PaperForkCgZ = 0.7
  177. FrameLength = (PaperWb*sin(HTA)-(rake-(PaperRadFront-PaperRadRear)*cos(HTA))).evalf()
  178. FrameCGNorm = ((PaperFrameCgZ - PaperRadRear-(PaperFrameCgX/sin(HTA))*cos(HTA))*sin(HTA)).evalf()
  179. FrameCGPar = (PaperFrameCgX / sin(HTA) + (PaperFrameCgZ - PaperRadRear - PaperFrameCgX / sin(HTA) * cos(HTA)) * cos(HTA)).evalf()
  180. tempa = (PaperForkCgZ - PaperRadFront)
  181. tempb = (PaperWb-PaperForkCgX)
  182. tempc = (sqrt(tempa**2+tempb**2)).evalf()
  183. PaperForkL = (PaperWb*cos(HTA)-(PaperRadFront-PaperRadRear)*sin(HTA)).evalf()
  184. ForkCGNorm = (rake+(tempc * sin(pi/2-HTA-acos(tempa/tempc)))).evalf()
  185. ForkCGPar = (tempc * cos((pi/2-HTA)-acos(tempa/tempc))-PaperForkL).evalf()
  186. # Here is the final assembly of the numerical values. The symbol 'v' is the
  187. # forward speed of the bicycle (a concept which only makes sense in the
  188. # upright, static equilibrium case?). These are in a dictionary which will
  189. # later be substituted in. Again the sign on the *product* of inertia
  190. # values is flipped here, due to different orientations of coordinate
  191. # systems.
  192. v = symbols('v')
  193. val_dict = {WFrad: PaperRadFront,
  194. WRrad: PaperRadRear,
  195. htangle: HTA,
  196. forkoffset: rake,
  197. forklength: PaperForkL,
  198. framelength: FrameLength,
  199. forkcg1: ForkCGPar,
  200. forkcg3: ForkCGNorm,
  201. framecg1: FrameCGNorm,
  202. framecg3: FrameCGPar,
  203. Iwr11: 0.0603,
  204. Iwr22: 0.12,
  205. Iwf11: 0.1405,
  206. Iwf22: 0.28,
  207. Ifork11: 0.05892,
  208. Ifork22: 0.06,
  209. Ifork33: 0.00708,
  210. Ifork31: 0.00756,
  211. Iframe11: 9.2,
  212. Iframe22: 11,
  213. Iframe33: 2.8,
  214. Iframe31: -2.4,
  215. mfork: 4,
  216. mframe: 85,
  217. mwf: 3,
  218. mwr: 2,
  219. g: 9.81,
  220. q1: 0,
  221. q2: 0,
  222. q4: 0,
  223. q5: 0,
  224. u1: 0,
  225. u2: 0,
  226. u3: v / PaperRadRear,
  227. u4: 0,
  228. u5: 0,
  229. u6: v / PaperRadFront}
  230. # Linearizes the forcing vector; the equations are set up as MM udot =
  231. # forcing, where MM is the mass matrix, udot is the vector representing the
  232. # time derivatives of the generalized speeds, and forcing is a vector which
  233. # contains both external forcing terms and internal forcing terms, such as
  234. # centripital or coriolis forces. This actually returns a matrix with as
  235. # many rows as *total* coordinates and speeds, but only as many columns as
  236. # independent coordinates and speeds.
  237. A, B, _ = KM.linearize(
  238. A_and_B=True,
  239. op_point={
  240. # Operating points for the accelerations are required for the
  241. # linearizer to eliminate u' terms showing up in the coefficient
  242. # matrices.
  243. u1.diff(): 0,
  244. u2.diff(): 0,
  245. u3.diff(): 0,
  246. u4.diff(): 0,
  247. u5.diff(): 0,
  248. u6.diff(): 0,
  249. u1: 0,
  250. u2: 0,
  251. u3: v / PaperRadRear,
  252. u4: 0,
  253. u5: 0,
  254. u6: v / PaperRadFront,
  255. q1: 0,
  256. q2: 0,
  257. q4: 0,
  258. q5: 0,
  259. },
  260. linear_solver="CRAMER",
  261. )
  262. # As mentioned above, the size of the linearized forcing terms is expanded
  263. # to include both q's and u's, so the mass matrix must have this done as
  264. # well. This will likely be changed to be part of the linearized process,
  265. # for future reference.
  266. A_s = A.xreplace(val_dict)
  267. B_s = B.xreplace(val_dict)
  268. A_s = A_s.evalf()
  269. B_s = B_s.evalf()
  270. # Finally, we construct an "A" matrix for the form xdot = A x (x being the
  271. # state vector, although in this case, the sizes are a little off). The
  272. # following line extracts only the minimum entries required for eigenvalue
  273. # analysis, which correspond to rows and columns for lean, steer, lean
  274. # rate, and steer rate.
  275. A = A_s.extract([1, 2, 3, 5], [1, 2, 3, 5])
  276. # Precomputed for comparison
  277. Res = Matrix([[ 0, 0, 1.0, 0],
  278. [ 0, 0, 0, 1.0],
  279. [9.48977444677355, -0.891197738059089*v**2 - 0.571523173729245, -0.105522449805691*v, -0.330515398992311*v],
  280. [11.7194768719633, -1.97171508499972*v**2 + 30.9087533932407, 3.67680523332152*v, -3.08486552743311*v]])
  281. # Actual eigenvalue comparison
  282. eps = 1.e-12
  283. for i in range(6):
  284. error = Res.subs(v, i) - A.subs(v, i)
  285. assert all(abs(x) < eps for x in error)