def _is_conditionally_positive_definite(kernel, m): # Tests whether the kernel is conditionally positive definite of order m. # See chapter 7 of Fasshauer's "Meshfree Approximation Methods with # MATLAB". nx = 10 ntests = 100 for ndim in [1, 2, 3, 4, 5]: # Generate sample points with a Halton sequence to avoid samples that # are too close to eachother, which can make the matrix singular. seq = Halton(ndim, scramble=False, seed=np.random.RandomState()) for _ in range(ntests): x = 2 * seq.random(nx) - 1 A = _rbfinterp_pythran._kernel_matrix(x, kernel) P = _vandermonde(x, m - 1) Q, R = np.linalg.qr(P, mode='complete') # Q2 forms a basis spanning the space where P.T.dot(x) = 0. Project # A onto this space, and then see if it is positive definite using # the Cholesky decomposition. If not, then the kernel is not c.p.d. # of order m. Q2 = Q[:, P.shape[1]:] B = Q2.T.dot(A).dot(Q2) try: np.linalg.cholesky(B) except np.linalg.LinAlgError: return False return True
def test_chunking(self, monkeypatch): # If the observed data comes from a polynomial, then the interpolant # should be able to reproduce the polynomial exactly, provided that # `degree` is sufficiently high. rng = np.random.RandomState(0) seq = Halton(2, scramble=False, seed=rng) degree = 3 largeN = 1000 + 33 # this is large to check that chunking of the RBFInterpolator is tested x = seq.random(50) xitp = seq.random(largeN) P = _vandermonde(x, degree) Pitp = _vandermonde(xitp, degree) poly_coeffs = rng.normal(0.0, 1.0, P.shape[1]) y = P.dot(poly_coeffs) yitp1 = Pitp.dot(poly_coeffs) interp = self.build(x, y, degree=degree) ce_real = interp._chunk_evaluator def _chunk_evaluator(*args, **kwargs): kwargs.update(memory_budget=100) return ce_real(*args, **kwargs) monkeypatch.setattr(interp, '_chunk_evaluator', _chunk_evaluator) yitp2 = interp(xitp) assert_allclose(yitp1, yitp2, atol=1e-8)
def test_smoothing_misfit(self, kernel): # Make sure we can find a smoothing parameter for each kernel that # removes a sufficient amount of noise. rng = np.random.RandomState(0) seq = Halton(1, scramble=False, seed=rng) noise = 0.2 rmse_tol = 0.1 smoothing_range = 10**np.linspace(-4, 1, 20) x = 3 * seq.random(100) y = _1d_test_function(x) + rng.normal(0.0, noise, (100, )) ytrue = _1d_test_function(x) rmse_within_tol = False for smoothing in smoothing_range: ysmooth = self.build(x, y, epsilon=1.0, smoothing=smoothing, kernel=kernel)(x) rmse = np.sqrt(np.mean((ysmooth - ytrue)**2)) if rmse < rmse_tol: rmse_within_tol = True break assert rmse_within_tol
def test_extreme_domains(self, kernel): # Make sure the interpolant remains numerically stable for very # large/small domains. seq = Halton(2, scramble=False, seed=np.random.RandomState()) scale = 1e50 shift = 1e55 x = seq.random(100) y = _2d_test_function(x) xitp = seq.random(100) if kernel in _SCALE_INVARIANT: yitp1 = self.build(x, y, kernel=kernel)(xitp) yitp2 = self.build( x*scale + shift, y, kernel=kernel )(xitp*scale + shift) else: yitp1 = self.build(x, y, epsilon=5.0, kernel=kernel)(xitp) yitp2 = self.build( x*scale + shift, y, epsilon=5.0/scale, kernel=kernel )(xitp*scale + shift) assert_allclose(yitp1, yitp2, atol=1e-8)
def test_inconsistent_x_dimensions_error(self): # ValueError should be raised if the observation points and evaluation # points have a different number of dimensions. y = Halton(2, scramble=False, seed=np.random.RandomState()).random(10) d = _2d_test_function(y) x = Halton(1, scramble=False, seed=np.random.RandomState()).random(10) match = 'Expected the second axis of `x`' with pytest.raises(ValueError, match=match): self.build(y, d)(x)
def test_scale_invariance_2d(self, kernel): # Verify that the functions in _SCALE_INVARIANT are insensitive to the # shape parameter (when smoothing == 0) in 2d. seq = Halton(2, scramble=False, seed=np.random.RandomState()) x = seq.random(100) y = _2d_test_function(x) xitp = seq.random(100) yitp1 = self.build(x, y, epsilon=1.0, kernel=kernel)(xitp) yitp2 = self.build(x, y, epsilon=2.0, kernel=kernel)(xitp) assert_allclose(yitp1, yitp2, atol=1e-8)
def test_equivalent_to_rbf_interpolator(self): seq = Halton(1, scramble=False, seed=np.random.RandomState()) x = 3 * seq.random(50) xitp = 3 * seq.random(50) y = _1d_test_function(x) yitp1 = self.build(x, y)(xitp) yitp2 = RBFInterpolator(x, y)(xitp) assert_allclose(yitp1, yitp2, atol=1e-8)
def test_interpolation_misfit_2d(self, kernel): # Make sure that each kernel, with its default `degree` and an # appropriate `epsilon`, does a good job at interpolation in 2d. seq = Halton(2, scramble=False, seed=np.random.RandomState()) x = seq.random(100) xitp = seq.random(100) y = _2d_test_function(x) ytrue = _2d_test_function(xitp) yitp = self.build(x, y, epsilon=5.0, kernel=kernel)(xitp) mse = np.mean((yitp - ytrue)**2) assert mse < 2.0e-4
def test_pickleable(self): # Make sure we can pickle and unpickle the interpolant without any # changes in the behavior. seq = Halton(1, scramble=False, seed=np.random.RandomState()) x = 3 * seq.random(50) xitp = 3 * seq.random(50) y = _1d_test_function(x) interp = self.build(x, y) yitp1 = interp(xitp) yitp2 = pickle.loads(pickle.dumps(interp))(xitp) assert_array_equal(yitp1, yitp2)
def test_complex_data(self): # Interpolating complex input should be the same as interpolating the # real and complex components. seq = Halton(2, scramble=False, seed=np.random.RandomState()) x = seq.random(100) xitp = seq.random(100) y = _2d_test_function(x) + 1j * _2d_test_function(x[:, ::-1]) yitp1 = self.build(x, y)(xitp) yitp2 = self.build(x, y.real)(xitp) yitp3 = self.build(x, y.imag)(xitp) assert_allclose(yitp1.real, yitp2) assert_allclose(yitp1.imag, yitp3)
def test_vector_data(self): # Make sure interpolating a vector field is the same as interpolating # each component separately. seq = Halton(2, scramble=False, seed=np.random.RandomState()) x = seq.random(100) xitp = seq.random(100) y = np.array([_2d_test_function(x), _2d_test_function(x[:, ::-1])]).T yitp1 = self.build(x, y)(xitp) yitp2 = self.build(x, y[:, 0])(xitp) yitp3 = self.build(x, y[:, 1])(xitp) assert_allclose(yitp1[:, 0], yitp2) assert_allclose(yitp1[:, 1], yitp3)
def test_equivalent_to_rbf_interpolator(self): seq = Halton(2, scramble=False, seed=np.random.RandomState()) x = seq.random(100) xitp = seq.random(100) y = _2d_test_function(x) yitp1 = self.build(x, y)(xitp) yitp2 = [] tree = cKDTree(x) for xi in xitp: _, nbr = tree.query(xi, 20) yitp2.append(RBFInterpolator(x[nbr], y[nbr])(xi[None])[0]) assert_allclose(yitp1, yitp2, atol=1e-8)
def test_array_smoothing(self): # Test using an array for `smoothing` to give less weight to a known # outlier. rng = np.random.RandomState(0) seq = Halton(1, scramble=False, seed=rng) degree = 2 x = seq.random(50) P = _vandermonde(x, degree) poly_coeffs = rng.normal(0.0, 1.0, P.shape[1]) y = P.dot(poly_coeffs) y_with_outlier = np.copy(y) y_with_outlier[10] += 1.0 smoothing = np.zeros((50, )) smoothing[10] = 1000.0 yitp = self.build(x, y_with_outlier, smoothing=smoothing)(x) # Should be able to reproduce the uncorrupted data almost exactly. assert_allclose(yitp, y, atol=1e-4)
def test_smoothing_limit_2d(self): # For large smoothing parameters, the interpolant should approach a # least squares fit of a polynomial with the specified degree. seq = Halton(2, scramble=False, seed=np.random.RandomState()) degree = 3 smoothing = 1e8 x = seq.random(100) xitp = seq.random(100) y = _2d_test_function(x) yitp1 = self.build(x, y, degree=degree, smoothing=smoothing)(xitp) P = _vandermonde(x, degree) Pitp = _vandermonde(xitp, degree) yitp2 = Pitp.dot(np.linalg.lstsq(P, y, rcond=None)[0]) assert_allclose(yitp1, yitp2, atol=1e-8)
def test_polynomial_reproduction(self): # If the observed data comes from a polynomial, then the interpolant # should be able to reproduce the polynomial exactly, provided that # `degree` is sufficiently high. rng = np.random.RandomState(0) seq = Halton(2, scramble=False, seed=rng) degree = 3 x = seq.random(50) xitp = seq.random(50) P = _vandermonde(x, degree) Pitp = _vandermonde(xitp, degree) poly_coeffs = rng.normal(0.0, 1.0, P.shape[1]) y = P.dot(poly_coeffs) yitp1 = Pitp.dot(poly_coeffs) yitp2 = self.build(x, y, degree=degree)(xitp) assert_allclose(yitp1, yitp2, atol=1e-8)