Tutorial for basic usageΒΆ

  • ../_images/sphx_glr_FOR_README_001.png
  • ../_images/sphx_glr_FOR_README_002.png
  • ../_images/sphx_glr_FOR_README_003.png
  • ../_images/sphx_glr_FOR_README_004.png

Out:

[[-3.14159265 -1.          2.        ]
 [-3.14159265 -0.5         2.        ]
 [-3.14159265  0.          2.        ]
 [-3.14159265  0.5         2.        ]
 [-3.14159265  1.          2.        ]
 [-2.74889357 -1.          1.27648914]
 [-2.74889357 -0.5         1.707008  ]
 [-2.74889357  0.          1.94618514]
 [-2.74889357  0.5         1.99402057]
 [-2.74889357  1.          1.85051428]
 [-2.35619449 -1.          0.66312624]
 [-2.35619449 -0.5         1.45862137]
 [-2.35619449  0.          1.90056311]
 [-2.35619449  0.5         1.98895146]
 [-2.35619449  1.          1.72378641]
 [-1.57079633 -1.          0.109375  ]
 [-1.57079633 -0.5         1.234375  ]
 [-1.57079633  0.          1.859375  ]
 [-1.57079633  0.5         1.984375  ]
 [-1.57079633  1.          1.609375  ]
 [-0.78539816 -1.          0.66312624]
 [-0.78539816 -0.5         1.45862137]
 [-0.78539816  0.          1.90056311]
 [-0.78539816  0.5         1.98895146]
 [-0.78539816  1.          1.72378641]
 [-0.39269908 -1.          1.27648914]
 [-0.39269908 -0.5         1.707008  ]
 [-0.39269908  0.          1.94618514]
 [-0.39269908  0.5         1.99402057]
 [-0.39269908  1.          1.85051428]
 [ 0.         -1.          2.        ]
 [ 0.         -0.5         2.        ]
 [ 0.          0.          2.        ]
 [ 0.          0.5         2.        ]
 [ 0.          1.          2.        ]
 [ 0.39269908 -1.          2.72351086]
 [ 0.39269908 -0.5         2.292992  ]
 [ 0.39269908  0.          2.05381486]
 [ 0.39269908  0.5         2.00597943]
 [ 0.39269908  1.          2.14948572]
 [ 0.78539816 -1.          3.33687376]
 [ 0.78539816 -0.5         2.54137863]
 [ 0.78539816  0.          2.09943689]
 [ 0.78539816  0.5         2.01104854]
 [ 0.78539816  1.          2.27621359]
 [ 1.57079633 -1.          3.890625  ]
 [ 1.57079633 -0.5         2.765625  ]
 [ 1.57079633  0.          2.140625  ]
 [ 1.57079633  0.5         2.015625  ]
 [ 1.57079633  1.          2.390625  ]
 [ 2.35619449 -1.          3.33687376]
 [ 2.35619449 -0.5         2.54137863]
 [ 2.35619449  0.          2.09943689]
 [ 2.35619449  0.5         2.01104854]
 [ 2.35619449  1.          2.27621359]
 [ 2.74889357 -1.          2.72351086]
 [ 2.74889357 -0.5         2.292992  ]
 [ 2.74889357  0.          2.05381486]
 [ 2.74889357  0.5         2.00597943]
 [ 2.74889357  1.          2.14948572]
 [ 3.14159265 -1.          2.        ]
 [ 3.14159265 -0.5         2.        ]
 [ 3.14159265  0.          2.        ]
 [ 3.14159265  0.5         2.        ]
 [ 3.14159265  1.          2.        ]]

Coefficients all same? True
Knots all same? True

Antiderivative of derivative:
 Coefficients differ by constant? True
Knots all same? True

Derivative of antiderivative:
 Coefficients the same? True
Knots all same? True

  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
 import ndsplines
 import numpy as np
 import matplotlib.pyplot as plt

 # generate grid of independent variables
 x = np.array([-1, -7/8, -3/4, -1/2, -1/4, -1/8, 0, 1/8, 1/4, 1/2, 3/4, 7/8, 1])*np.pi
 y = np.array([-1, -1/2, 0, 1/2, 1])
 meshx, meshy = np.meshgrid(x, y, indexing='ij')
 gridxy = np.stack((meshx, meshy), axis=-1)


 # generate denser grid of independent variables to interpolate
 sparse_dense = 2**7
 xx = np.concatenate([np.linspace(x[i], x[i+1], sparse_dense) for i in range(x.size-1)]) # np.linspace(x[0], x[-1], x.size*sparse_dense)
 yy = np.concatenate([np.linspace(y[i], y[i+1], sparse_dense) for i in range(y.size-1)]) # np.linspace(y[0], y[-1], y.size*sparse_dense)
 gridxxyy = np.stack(np.meshgrid(xx, yy, indexing='ij'), axis=-1)

 def plots(sparse_data, dense_data, ylabel='f(x,y)'):
     fig, axes = plt.subplots(1, 2, constrained_layout=True)
     for yidx in range(sparse_data.shape[1]):
         axes[0].plot(x, sparse_data[:, yidx], 'o', color='C%d'%yidx, label='y=%.2f'%y[yidx])
         axes[0].plot(xx, dense_data[:, np.clip(yidx*sparse_dense, 0, yy.size-1)], color='C%d'%yidx)# label='y=%.1f'%y[yidx])

     axes[0].legend()
     axes[0].set_xlabel('x')
     axes[0].set_ylabel(ylabel)
     for xidx in range(sparse_data.shape[0]//2):
         axes[1].plot(yy, dense_data[(xidx+3)*sparse_dense, :], '--', color='C%d'%xidx,)# label='x=%.1f'%x[xidx+3])
         axes[1].plot(y, sparse_data[xidx+3, :], 'o', color='C%d'%xidx, label='x=%.1f'%x[xidx+3],)

     axes[1].legend()
     axes[1].set_xlabel('y')
     plt.show()

 # evaluate a function to interpolate over input grid
 meshf = np.sin(meshx) * (meshy-3/8)**2 + 2

 # create the interpolating splane
 interp = ndsplines.make_interp_spline(gridxy, meshf)

 # evaluate spline over denser grid
 meshff = interp(gridxxyy)


 plots(meshf, meshff)


 ##

 # as subplots
 fig, axes = plt.subplots(1,2, constrained_layout=True)

 gridxxy = np.stack(np.meshgrid(xx, y, indexing='ij'), axis=-1)
 meshff = interp(gridxxy)

 for yidx in range(meshf.shape[1]):
     axes[0].plot(x, meshf[:, yidx], 'o', color='C%d'%yidx, label='y=%.1f'%y[yidx])
     axes[0].plot(xx, meshff[:, yidx], color='C%d'%yidx)
 axes[0].legend()
 axes[0].set_xlabel('$x$')
 axes[0].set_ylabel('$f(x,y)$')

 # y-dir plot
 gridxyy = np.stack(np.meshgrid(x, yy, indexing='ij'), axis=-1)

 meshff = interp(gridxyy)
 for xidx in range(meshf.shape[0]//2):
     axes[1].plot(yy, meshff[xidx*1+3, :], '--', color='C%d'%xidx, label='x=%.1f'%x[xidx*1+3])
     axes[1].plot(y, meshf[xidx*1+3, :], 'o', color='C%d'%xidx)

 axes[1].legend()
 axes[1].set_xlabel('$y$')
 # plt.ylabel(r'$\frac{\partial f(x,y)}{\partial y}$')
 plt.show()

 ##

 # we could also use tidy data format to make the grid

 tidy_data = np.dstack((gridxy, meshf)).reshape((-1,3))
 print(tidy_data)

 tidy_interp = ndsplines.make_interp_spline_from_tidy(tidy_data, [0,1], [2])

 print("\nCoefficients all same?", np.all(tidy_interp.coefficients == interp.coefficients))
 print("Knots all same?", np.all([np.all(knot0 == knot1) for knot0, knot1 in zip(tidy_interp.knots, interp.knots)]))

 # send to example of least squares
 ##
 # two ways to evaluate derivative - y direction

 deriv_interp = interp.derivative(1)
 deriv1 = deriv_interp(gridxy)
 deriv2 = interp(gridxxyy, nus=np.array([0,1]))

 plots(deriv1, deriv2, r'$\frac{\partial f(x,y)}{\partial y}$')

 ##
 # two ways to evaluate derivatives x-direction: create a derivative spline or call with nus:
 deriv_interp = interp.derivative(0)
 deriv1 = deriv_interp(gridxy)
 deriv2 = interp(gridxxyy, nus=np.array([1,0]))

 plots(deriv1, deriv2, r'$\frac{\partial f(x,y)}{\partial x}$')
 ##

 # Calculus demonstration
 interp1 = deriv_interp.antiderivative(0)
 coeff_diff = interp1.coefficients - interp.coefficients
 print("\nAntiderivative of derivative:\n","Coefficients differ by constant?", np.allclose(interp1.coefficients+2.0, interp.coefficients))
 print("Knots all same?", np.all([np.all(knot0 == knot1) for knot0, knot1 in zip(interp1.knots, interp.knots)]))

 antideriv_interp = interp.antiderivative(0)

 interp2 = antideriv_interp.derivative(0)
 print("\nDerivative of antiderivative:\n","Coefficients the same?", np.allclose(interp2.coefficients, interp.coefficients))
 print("Knots all same?", np.all([np.all(knot0 == knot1) for knot0, knot1 in zip(interp2.knots, interp.knots)]))

Total running time of the script: ( 0 minutes 11.524 seconds)

Gallery generated by Sphinx-Gallery