1
from sympy import evalf, symbols, zeros, pi, sin, cos, sqrt, acos, Matrix
2
from sympy.physics.mechanics import (ReferenceFrame, dynamicsymbols, inertia,
3
KanesMethod, RigidBody, Point, dot)
4
from sympy.utilities.pytest import slow, ON_TRAVIS, skip
9
skip("Too slow for travis.")
10
# Code to get equations of motion for a bicycle modeled as in:
11
# J.P Meijaard, Jim M Papadopoulos, Andy Ruina and A.L Schwab. Linearized
12
# dynamics equations for the balance and steer of a bicycle: a benchmark
13
# and review. Proceedings of The Royal Society (2007) 463, 1955-1982
14
# doi: 10.1098/rspa.2007.1857
16
# Note that this code has been crudely ported from Autolev, which is the
17
# reason for some of the unusual naming conventions. It was purposefully as
18
# similar as possible in order to aide debugging.
20
# Declare Coordinates & Speeds
21
# Simple definitions for qdots - qd = u
22
# Speeds are: yaw frame ang. rate, roll frame ang. rate, rear wheel frame
23
# ang. rate (spinning motion), frame ang. rate (pitching motion), steering
24
# frame ang. rate, and front wheel ang. rate (spinning motion).
25
# Wheel positions are ignorable coordinates, so they are not introduced.
26
q1, q2, q4, q5 = dynamicsymbols('q1 q2 q4 q5')
27
q1d, q2d, q4d, q5d = dynamicsymbols('q1 q2 q4 q5', 1)
28
u1, u2, u3, u4, u5, u6 = dynamicsymbols('u1 u2 u3 u4 u5 u6')
29
u1d, u2d, u3d, u4d, u5d, u6d = dynamicsymbols('u1 u2 u3 u4 u5 u6', 1)
31
# Declare System's Parameters
32
WFrad, WRrad, htangle, forkoffset = symbols('WFrad WRrad htangle forkoffset')
33
forklength, framelength, forkcg1 = symbols('forklength framelength forkcg1')
34
forkcg3, framecg1, framecg3, Iwr11 = symbols('forkcg3 framecg1 framecg3 Iwr11')
35
Iwr22, Iwf11, Iwf22, Iframe11 = symbols('Iwr22 Iwf11 Iwf22 Iframe11')
36
Iframe22, Iframe33, Iframe31, Ifork11 = symbols('Iframe22 Iframe33 Iframe31 Ifork11')
37
Ifork22, Ifork33, Ifork31, g = symbols('Ifork22 Ifork33 Ifork31 g')
38
mframe, mfork, mwf, mwr = symbols('mframe mfork mwf mwr')
40
# Set up reference frames for the system
44
# WR - rear wheel, rotation angle is ignorable coordinate so not oriented
45
# Frame - bicycle frame
46
# TempFrame - statically rotated frame for easier reference inertia definition
48
# TempFork - statically rotated frame for easier reference inertia definition
49
# WF - front wheel, again posses a ignorable coordinate
50
N = ReferenceFrame('N')
51
Y = N.orientnew('Y', 'Axis', [q1, N.z])
52
R = Y.orientnew('R', 'Axis', [q2, Y.x])
53
Frame = R.orientnew('Frame', 'Axis', [q4 + htangle, R.y])
54
WR = ReferenceFrame('WR')
55
TempFrame = Frame.orientnew('TempFrame', 'Axis', [-htangle, Frame.y])
56
Fork = Frame.orientnew('Fork', 'Axis', [q5, Frame.x])
57
TempFork = Fork.orientnew('TempFork', 'Axis', [-htangle, Fork.y])
58
WF = ReferenceFrame('WF')
60
# Kinematics of the Bicycle First block of code is forming the positions of
62
# rear wheel contact -> rear wheel mass center -> frame mass center +
63
# frame/fork connection -> fork mass center + front wheel mass center ->
64
# front wheel contact point
65
WR_cont = Point('WR_cont')
66
WR_mc = WR_cont.locatenew('WR_mc', WRrad * R.z)
67
Steer = WR_mc.locatenew('Steer', framelength * Frame.z)
68
Frame_mc = WR_mc.locatenew('Frame_mc', - framecg1 * Frame.x
70
Fork_mc = Steer.locatenew('Fork_mc', - forkcg1 * Fork.x
72
WF_mc = Steer.locatenew('WF_mc', forklength * Fork.x + forkoffset * Fork.z)
73
WF_cont = WF_mc.locatenew('WF_cont', WFrad * (dot(Fork.y, Y.z) * Fork.y -
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.
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)
92
# Form the velocities of the previously defined points, using the 2 - point
93
# theorem (written out by hand here). Accelerations again are calculated
94
# automatically when first needed.
96
WR_mc.v2pt_theory(WR_cont, N, WR)
97
Steer.v2pt_theory(WR_mc, N, Frame)
98
Frame_mc.v2pt_theory(WR_mc, N, Frame)
99
Fork_mc.v2pt_theory(Steer, N, Fork)
100
WF_mc.v2pt_theory(Steer, N, Fork)
101
WF_cont.v2pt_theory(WF_mc, N, WF)
103
# Sets the inertias of each body. Uses the inertia frame to construct the
104
# inertia dyadics. Wheel inertias are only defined by principle moments of
105
# inertia, and are in fact constant in the frame and fork reference frames;
106
# it is for this reason that the orientations of the wheels does not need
107
# to be defined. The frame and fork inertias are defined in the 'Temp'
108
# frames which are fixed to the appropriate body frames; this is to allow
109
# easier input of the reference values of the benchmark paper. Note that
110
# due to slightly different orientations, the products of inertia need to
111
# have their signs flipped; this is done later when entering the numerical
114
Frame_I = (inertia(TempFrame, Iframe11, Iframe22, Iframe33, 0, 0, Iframe31), Frame_mc)
115
Fork_I = (inertia(TempFork, Ifork11, Ifork22, Ifork33, 0, 0, Ifork31), Fork_mc)
116
WR_I = (inertia(Frame, Iwr11, Iwr22, Iwr11), WR_mc)
117
WF_I = (inertia(Fork, Iwf11, Iwf22, Iwf11), WF_mc)
119
# Declaration of the RigidBody containers. ::
121
BodyFrame = RigidBody('BodyFrame', Frame_mc, Frame, mframe, Frame_I)
122
BodyFork = RigidBody('BodyFork', Fork_mc, Fork, mfork, Fork_I)
123
BodyWR = RigidBody('BodyWR', WR_mc, WR, mwr, WR_I)
124
BodyWF = RigidBody('BodyWF', WF_mc, WF, mwf, WF_I)
127
# The kinematic differential equations; they are defined quite simply. Each
128
# entry in this list is equal to zero.
129
kd = [q1d - u1, q2d - u2, q4d - u4, q5d - u5]
131
# The nonholonomic constraints are the velocity of the front wheel contact
132
# point dotted into the X, Y, and Z directions; the yaw frame is used as it
133
# is "closer" to the front wheel (1 less DCM connecting them). These
134
# constraints force the velocity of the front wheel contact point to be 0
135
# in the inertial frame; the X and Y direction constraints enforce a
136
# "no-slip" condition, and the Z direction constraint forces the front
137
# wheel contact point to not move away from the ground frame, essentially
138
# replicating the holonomic constraint which does not allow the frame pitch
139
# to change in an invalid fashion.
141
conlist_speed = [WF_cont.vel(N) & Y.x, WF_cont.vel(N) & Y.y, WF_cont.vel(N) & Y.z]
143
# The holonomic constraint is that the position from the rear wheel contact
144
# point to the front wheel contact point when dotted into the
145
# normal-to-ground plane direction must be zero; effectively that the front
146
# and rear wheel contact points are always touching the ground plane. This
147
# is actually not part of the dynamic equations, but instead is necessary
148
# for the lineraization process.
150
conlist_coord = [WF_cont.pos_from(WR_cont) & Y.z]
152
# The force list; each body has the appropriate gravitational force applied
153
# at its mass center.
154
FL = [(Frame_mc, -mframe * g * Y.z),
155
(Fork_mc, -mfork * g * Y.z),
156
(WF_mc, -mwf * g * Y.z),
157
(WR_mc, -mwr * g * Y.z)]
158
BL = [BodyFrame, BodyFork, BodyWR, BodyWF]
161
# The N frame is the inertial frame, coordinates are supplied in the order
162
# of independent, dependent coordinates, as are the speeds. The kinematic
163
# differential equation are also entered here. Here the dependent speeds
164
# are specified, in the same order they were provided in earlier, along
165
# with the non-holonomic constraints. The dependent coordinate is also
166
# provided, with the holonomic constraint. Again, this is only provided
167
# for the linearization process.
169
KM = KanesMethod(N, q_ind=[q1, q2, q5],
170
q_dependent=[q4], configuration_constraints=conlist_coord,
172
u_dependent=[u1, u4, u6], velocity_constraints=conlist_speed,
174
(fr, frstar) = KM.kanes_equations(FL, BL)
176
# This is the start of entering in the numerical values from the benchmark
177
# paper to validate the eigen values of the linearized equations from this
178
# model to the reference eigen values. Look at the aforementioned paper for
179
# more information. Some of these are intermediate values, used to
180
# transform values from the paper into the coordinate systems used in this
184
HTA = evalf.N(pi / 2 - pi / 10)
186
rake = evalf.N(-(TrailPaper*sin(HTA)-(PaperRadFront*cos(HTA))))
192
FrameLength = evalf.N(PaperWb*sin(HTA)-(rake-(PaperRadFront-PaperRadRear)*cos(HTA)))
193
FrameCGNorm = evalf.N((PaperFrameCgZ - PaperRadRear-(PaperFrameCgX/sin(HTA))*cos(HTA))*sin(HTA))
194
FrameCGPar = evalf.N((PaperFrameCgX / sin(HTA) + (PaperFrameCgZ - PaperRadRear - PaperFrameCgX / sin(HTA) * cos(HTA)) * cos(HTA)))
195
tempa = evalf.N((PaperForkCgZ - PaperRadFront))
196
tempb = evalf.N((PaperWb-PaperForkCgX))
197
tempc = evalf.N(sqrt(tempa**2+tempb**2))
198
PaperForkL = evalf.N((PaperWb*cos(HTA)-(PaperRadFront-PaperRadRear)*sin(HTA)))
199
ForkCGNorm = evalf.N(rake+(tempc * sin(pi/2-HTA-acos(tempa/tempc))))
200
ForkCGPar = evalf.N(tempc * cos((pi/2-HTA)-acos(tempa/tempc))-PaperForkL)
202
# Here is the final assembly of the numerical values. The symbol 'v' is the
203
# forward speed of the bicycle (a concept which only makes sense in the
204
# upright, static equilibrium case?). These are in a dictionary which will
205
# later be substituted in. Again the sign on the *product* of inertia
206
# values is flipped here, due to different orientations of coordinate
209
val_dict = {WFrad: PaperRadFront,
213
forklength: PaperForkL,
214
framelength: FrameLength,
217
framecg1: FrameCGNorm,
218
framecg3: FrameCGPar,
242
u3: v / PaperRadRear,
245
u6: v / PaperRadFront}
247
# Linearizes the forcing vector; the equations are set up as MM udot =
248
# forcing, where MM is the mass matrix, udot is the vector representing the
249
# time derivatives of the generalized speeds, and forcing is a vector which
250
# contains both external forcing terms and internal forcing terms, such as
251
# centripital or coriolis forces. This actually returns a matrix with as
252
# many rows as *total* coordinates and speeds, but only as many columns as
253
# independent coordinates and speeds.
255
forcing_lin = KM.linearize()[0]
257
# As mentioned above, the size of the linearized forcing terms is expanded
258
# to include both q's and u's, so the mass matrix must have this done as
259
# well. This will likely be changed to be part of the linearized process,
260
# for future reference.
261
MM_full = KM.mass_matrix_full
263
MM_full_s = MM_full.subs(val_dict)
264
forcing_lin_s = forcing_lin.subs(KM.kindiffdict()).subs(val_dict)
267
MM_full_s = MM_full_s.evalf()
268
forcing_lin_s = forcing_lin_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
Amat = MM_full_s.inv() * forcing_lin_s
276
A = Amat.extract([1, 2, 4, 6], [1, 2, 3, 5])
278
# Precomputed for comparison
279
Res = Matrix([[ 0, 0, 1.0, 0],
281
[9.48977444677355, -0.891197738059089*v**2 - 0.571523173729245, -0.105522449805691*v, -0.330515398992311*v],
282
[11.7194768719633, -1.97171508499972*v**2 + 30.9087533932407, 3.67680523332152*v, -3.08486552743311*v]])
285
# Actual eigenvalue comparison
288
error = Res.subs(v, i) - A.subs(v, i)
289
assert all(abs(x) < eps for x in error)