Map (series) inversion¶
(by Dario Izzo)
In this notebook we show how to use the map inversion algorithm to locally invert systems of non linear equations.
Importing stuff¶
[1]:
from pyaudi import gdual_double as gdual
from pyaudi import sin, cos
from pyaudi import invert_map
The Double Pendulum end-point¶
Let us consider, in 2D, the end point position of a double pendulum.
x=l1sin(θ1)+l2sin(θ2)y=−l1cos(θ1)−l2cos(θ2)
We will develop a model for the inversion of this system of equations, so that θ1 and θ2 can be expressed as a function of x and y
[2]:
# We show an image for clarity
from IPython.display import Image
from IPython.core.display import HTML
Image(url= "http://www.physicsandbox.com/assets/images/pendulum.jpg")
[2]:

Define our variables¶
[3]:
# Lengths for the pendulum arms
L1 = 1.
L2 = 1.3
# Nominal values for thetas
th1n = 0.3
th2n = -0.2
# Gduals (high order)
th1 = gdual(th1n, "\\theta_1", 11)
th2 = gdual(th2n, "\\theta_2", 11)
# Equations
x = L1*sin(th1)+L2*sin(th2)
y = -L1*cos(th1)-L2*cos(th2)
[4]:
# Lets visualize the Taylor polinomial expansion for x
x
[4]:
1.27409dθ2+2.63265e−06dθ19−0.000410445dθ16+3.51104e−06dθ29−0.159223dθ13+0.0372501+0.000358709dθ26−0.000252795dθ27−0.212348dθ23−2.39332e−08dθ111+7.32937e−06dθ18+0.00796114dθ15−0.0107613dθ24−0.14776dθ12+0.0106174dθ25−8.14374e−08dθ110−6.40551e−06dθ28−0.000189551dθ17+0.0123133dθ14+0.955336dθ1+…+O(12)
[5]:
# Lets visualize the Taylor polinomial expansion for y
y
[5]:
−0.25827dθ2+8.14374e−07dθ19+0.00132686dθ16−7.11723e−07dθ29−0.0492534dθ13−2.22942+0.00176956dθ26+5.12441e−05dθ27+0.043045dθ23−7.4034e−09dθ111−2.36939e−05dθ18+0.00246267dθ15−0.0530869dθ24+0.477668dθ12−0.00215225dθ25+2.63265e−07dθ110−3.15994e−05dθ28−5.8635e−05dθ17−0.0398057dθ14+0.29552dθ1+…+O(12)
Compute the inverse map¶
[6]:
# And let us invert the relationship
res = invert_map([x,y])
th1_map = res[0]
th2_map = res[1]
[7]:
# Lets visualize the Taylor polinomial expansion for th1
th1_map
[7]:
−33.3889dp02dp12+2.97638dp02dp1−1.53985e+08dp02dp18+7.49742dp04dp1−1.18854dp08−450.043dp14−0.461981dp02−60579.9dp0dp16+13499.7dp07dp13−19.2519dp03dp12+4.30881e+06dp06dp15+1.03732dp03dp1−2.89428e+08dp03dp18+4.30674dp05dp1−0.913559dp09+7.34785e+07dp19−5141.9dp03dp14−0.0793933dp03−45079.9dp05dp14−47904.3dp16+…+O(12)
[8]:
# Lets visualize the Taylor polinomial expansion for th2
th2_map
[8]:
25.7695dp02dp12−2.34715dp02dp1+1.18575e+08dp02dp18−5.78329dp04dp1+0.917051dp08+347.27dp14+0.343531dp02+46667.1dp0dp16−10399dp07dp13+14.7909dp03dp12−3.3184e+06dp06dp15−0.872053dp03dp1+2.22845e+08dp03dp18−3.30875dp05dp1+0.704389dp09−5.65823e+07dp19+3962.12dp03dp14+0.00876767dp03+34726dp05dp14+36914.5dp16+…+O(12)
Compute θ from x¶
[9]:
# First we extract the x,y position around the nominal thetas
xn = x.constant_cf
yn = y.constant_cf
[10]:
print("x nominal is: ", xn)
print("y nominal is: ", yn)
x nominal is: 0.037250076627759976
y nominal is: -2.22942304031922
[11]:
# Lets assume some desired (close to nominal) values for the end point
xd = 0.04
yd = -2.21
# And compute the change with respect to the nominal position
dx = xd - xn
dy = yd - yn
# We now compute the thetas
th1d = th1n + th1_map.evaluate({"dp0": dx, "dp1": dy})
th2d = th2n + th2_map.evaluate({"dp0": dx, "dp1": dy})
# Let us check that indeed they are producing the desired end point position
xdi = L1*sin(th1d)+L2*sin(th2d)
ydi = -L1*cos(th1d)-L2*cos(th2d)
[12]:
print("Error in x: ", xdi-xd)
print("Error in y: ", ydi-yd)
Error in x: -9.838373171699999e-12
Error in y: 1.7152856912616699e-10
[13]:
ydi
[13]:
-2.2099999998284714
The Kepler’s Equation¶
Let’s consider the long standing problem of inverting Kepler’s equation:
M=E−esin(E)
and face it using the map inversion functionality of pyaudi
Define our variables¶
[14]:
# Nominal values (expansion points) for eccentric anomaly and eccentricity
En = 0
en = 0
# gduals
E = gdual(En, "E", 11)
e = gdual(en, "e", 11)
# The equation
M = E - e*sin(E)
[15]:
# Lets visualize the Taylor polinomial expansion for M (mean anomaly)
M
[15]:
dE−0.00833333dE5de+0.000198413dE7de−dEde−2.75573e−06dE9de+0.166667dE3de+O(12)
Invert the map¶
[16]:
# We need to have the same symbol set in the variables
e.extend_symbol_set(["de", "dE"])
# We may now create the inverse map
res = invert_map([M, e])
E_map = res[0]
e_map = res[1]
[17]:
# Lets visualize the Taylor polinomial expansion for E as a function of M (p0) and e (p1)
E_map
[17]:
−20dp03dp18+0.00833333dp05dp1−3.33333dp03dp14+2.8dp05dp14+0.133333dp05dp12+dp0dp110−1.66667dp03dp13+0.758333dp05dp13+0.000705467dp09dp12+dp0dp18−0.000198413dp07dp1+2.75573e−06dp09dp1−1.07937dp07dp14−0.0126984dp07dp12+dp0dp16+dp0dp12−14dp03dp17+8.05dp05dp15−0.162698dp07dp13−5.83333dp03dp15+…+O(12)
Compute E from M¶
[18]:
# We now investigate how accurate this map is around the nominal point
from mpl_toolkits.mplot3d import Axes3D
import matplotlib.pyplot as plt
from matplotlib import cm
from matplotlib.ticker import LinearLocator, FormatStrFormatter
import numpy as np
N = 50
M_grid = np.linspace(0,0.1,N)
e_grid = np.linspace(0,0.1,N)
# Make data.
X, Y = np.meshgrid(M_grid, e_grid)
Z = np.zeros([N,N])
for i in range(N):
for j in range(N):
res = E_map.evaluate({"dp0": X[i,j], "dp1": Y[i,j]}) + En
Z[i,j] = np.log10(np.abs(X[i,j] - res + Y[i,j]*sin(res)) + 1e-14)
# Prepare the plot
fig = plt.figure(figsize=(12, 9))
ax = fig.gca(projection='3d')
# Plot the surface.
surf = ax.plot_surface(X, Y, Z, cmap=cm.coolwarm,
linewidth=0, antialiased=False)
ax.set_xlabel('E')
ax.set_ylabel('e')
ax.set_zlabel('Error (log)')
plt.show()
