def test_spherical_quadrature():
    """
    Testing spherical quadrature rule versus numerical integration.
    """

    b = 8  # 10

    # Create grids on the sphere
    x_gl = S2.meshgrid(b=b, grid_type='Gauss-Legendre')
    x_cc = S2.meshgrid(b=b, grid_type='Clenshaw-Curtis')
    x_soft = S2.meshgrid(b=b, grid_type='SOFT')
    x_gl = np.c_[x_gl[0][..., None], x_gl[1][..., None]]
    x_cc = np.c_[x_cc[0][..., None], x_cc[1][..., None]]
    x_soft = np.c_[x_soft[0][..., None], x_soft[1][..., None]]

    # Compute quadrature weights
    w_gl = S2.quadrature_weights(b=b, grid_type='Gauss-Legendre')
    w_cc = S2.quadrature_weights(b=b, grid_type='Clenshaw-Curtis')
    w_soft = S2.quadrature_weights(b=b, grid_type='SOFT')

    # Define a polynomial function, to be evaluated at one point or at an array of points
    def f1a(xs):
        xc = S2.change_coordinates(coords=xs, p_from='S', p_to='C')
        return xc[..., 0]**2 * xc[..., 1] - 1.4 * xc[..., 2] * xc[
            ..., 1]**3 + xc[..., 1] - xc[..., 2]**2 + 2.

    def f1(theta, phi):
        xs = np.array([theta, phi])
        return f1a(xs)

    # Obtain the "true" value of the integral of the function over the sphere, using scipy's numerical integration
    # routines
    i1 = S2.integrate(f1, normalize=False)

    # Compute the integral using the quadrature formulae
    # i1_gl_w = (w_gl * f1a(x_gl)).sum()
    i1_gl_w = S2.integrate_quad(f1a(x_gl),
                                grid_type='Gauss-Legendre',
                                normalize=False,
                                w=w_gl)
    print(i1_gl_w, i1, 'diff:', np.abs(i1_gl_w - i1))
    assert np.isclose(np.abs(i1_gl_w - i1), 0.0)

    # i1_cc_w = (w_cc * f1a(x_cc)).sum()
    i1_cc_w = S2.integrate_quad(f1a(x_cc),
                                grid_type='Clenshaw-Curtis',
                                normalize=False,
                                w=w_cc)
    print(i1_cc_w, i1, 'diff:', np.abs(i1_cc_w - i1))
    assert np.isclose(np.abs(i1_cc_w - i1), 0.0)

    i1_soft_w = (w_soft * f1a(x_soft)).sum()
    print(i1_soft_w, i1, 'diff:', np.abs(i1_soft_w - i1))
    print(i1_soft_w)
    print(i1)
Esempio n. 2
0
def naive_S2_conv(f1, f2, alpha, beta, gamma, g_parameterization='EA323'):
    """
    Compute int_S^2 f1(x) f2(g^{-1} x)* dx,
    where x = (theta, phi) is a point on the sphere S^2,
    and g = (alpha, beta, gamma) is a point in SO(3) in Euler angle parameterization

    :param f1, f2: functions to be convolved
    :param alpha, beta, gamma: the rotation at which to evaluate the result of convolution
    :return:
    """

    # This fails
    def integrand(theta, phi):
        g_inv = SO3.invert((alpha, beta, gamma),
                           parameterization=g_parameterization)
        g_inv_theta, g_inv_phi, _ = SO3.transform_r3(
            g=g_inv,
            x=(theta, phi, 1.),
            g_parameterization=g_parameterization,
            x_parameterization='S')
        return f1(theta, phi) * f2(g_inv_theta, g_inv_phi).conj()

    return S2.integrate(f=integrand, normalize=True)
Esempio n. 3
0
def check_orthogonality(L_max=3,
                        grid_type='Gauss-Legendre',
                        field='real',
                        normalization='quantum',
                        condon_shortley=True):

    theta, phi = S2.meshgrid(b=L_max + 1, grid_type=grid_type)
    w = S2.quadrature_weights(b=L_max + 1, grid_type=grid_type)

    for l in range(L_max):
        for m in range(-l, l + 1):
            for l2 in range(L_max):
                for m2 in range(-l2, l2 + 1):
                    Ylm = sh(l, m, theta, phi, field, normalization,
                             condon_shortley)
                    Ylm2 = sh(l2, m2, theta, phi, field, normalization,
                              condon_shortley)

                    dot_numerical = S2.integrate_quad(Ylm * Ylm2.conj(),
                                                      grid_type=grid_type,
                                                      normalize=False,
                                                      w=w)

                    dot_numerical2 = S2.integrate(
                        lambda t, p: sh(l, m, t, p, field, normalization, condon_shortley) * \
                                     sh(l2, m2, t, p, field, normalization, condon_shortley).conj(), normalize=False)

                    sqnorm_analytical = sh_squared_norm(l,
                                                        normalization,
                                                        normalized_haar=False)
                    dot_analytical = sqnorm_analytical * (l == l2 and m == m2)

                    print(l, m, l2, m2, field, normalization, condon_shortley,
                          dot_analytical, dot_numerical, dot_numerical2)
                    assert np.isclose(dot_numerical, dot_analytical)
                    assert np.isclose(dot_numerical2, dot_analytical)