예제 #1
0
파일: test_utils.py 프로젝트: theochem/grid
    def test_generate_real_spherical_is_orthonormal(self):
        """Test generated real spherical harmonics is an orthonormal set."""
        atgrid = AngularGrid(degree=7)
        pts = atgrid.points
        wts = atgrid.weights
        r = np.linalg.norm(pts, axis=1)
        # polar
        phi = np.arccos(pts[:, 2] / r)
        # azimuthal
        theta = np.arctan2(pts[:, 1], pts[:, 0])
        # generate spherical harmonics
        sph_h = generate_real_spherical_harmonics(3, theta, phi)  # l_max = 3
        assert sph_h.shape == (16, 26)
        for _ in range(100):
            n1, n2 = np.random.randint(0, 16, 2)
            re = sum(sph_h[n1] * sph_h[n2] * wts)
            if n1 != n2:
                print(n1, n2, re)
                assert_almost_equal(re, 0)
            else:
                print(n1, n2, re)
                assert_almost_equal(re, 1)

        for i in range(10):
            sph_h = generate_real_spherical_harmonics(i, theta, phi)
            assert sph_h.shape == ((i + 1)**2, 26)
예제 #2
0
    def test_fitting_product_of_spherical_harmonic(self):
        r"""Test fitting the radial components of r**2 times spherical harmonic."""
        max_degree = 7  # Maximum degree
        rad = OneDGrid(np.linspace(0.0, 1.0, num=10), np.ones(10), (0, np.inf))
        atom_grid = AtomGrid(rad, degrees=[max_degree])
        max_degree = atom_grid.l_max
        spherical = atom_grid.convert_cartesian_to_spherical()

        spherical_harmonic = generate_real_spherical_harmonics(
            4,
            spherical[:, 1],
            spherical[:, 2]  # theta, phi points
        )
        # Test on the function r^2 * Y^1_3
        func_vals = spherical[:,
                              0]**2.0 * spherical_harmonic[(3 + 1)**2 + 1, :]
        # Fit radial components
        fit = atom_grid.radial_component_splines(func_vals)
        radial_pts = np.arange(0.0, 1.0, 0.01)
        i = 0
        for l_value in range(0, max_degree // 2):
            for m in ([0] + [x for x in range(1, l_value + 1)] +
                      [x for x in range(-l_value, 0)]):
                if i != (3 + 1)**2 + 1:
                    assert_almost_equal(fit[i](radial_pts), 0.0, decimal=8)
                else:
                    # Test that on the right spherical harmonic the function r* Y^1_3 projects
                    # to \rho^{1,3}(r) = r
                    assert_almost_equal(fit[i](radial_pts), radial_pts**2.0)
                i += 1
예제 #3
0
    def test_integrating_angular_components(self):
        """Test radial points that contain zero."""
        odg = OneDGrid(np.array([0.0, 1e-16, 1e-8, 1e-4, 1e-2]), np.ones(5),
                       (0, np.inf))
        atom_grid = AtomGrid(odg, degrees=[3])
        spherical = atom_grid.convert_cartesian_to_spherical()
        # Evaluate all spherical harmonics on the atomic grid points (r_i, theta_j, phi_j).
        spherical_harmonics = generate_real_spherical_harmonics(
            3,
            spherical[:, 1],
            spherical[:, 2]  # theta, phi points
        )
        # Convert to three-dimensional array (Degrees, Order, Points)
        spherical_array = np.zeros((3, 2 * 3 + 1, len(atom_grid.points)))
        spherical_array[0, 0, :] = spherical_harmonics[0, :]  # (l,m) = (0,0)
        spherical_array[1, 0, :] = spherical_harmonics[1, :]  # = (1, 0)
        spherical_array[1, 1, :] = spherical_harmonics[2, :]  # = (1, 1)
        spherical_array[1, 2, :] = spherical_harmonics[3, :]  # = (1, -1)
        spherical_array[2, 0, :] = spherical_harmonics[4, :]  # = (2, 0)
        spherical_array[2, 1, :] = spherical_harmonics[5, :]  # = (2, 2)
        spherical_array[2, 2, :] = spherical_harmonics[6, :]  # = (2, 1)
        spherical_array[2, 3, :] = spherical_harmonics[7, :]  # = (2, -2)
        spherical_array[2, 4, :] = spherical_harmonics[8, :]  # = (2, -1)

        integral = atom_grid.integrate_angular_coordinates(spherical_array)
        assert integral.shape == (3, 2 * 3 + 1, 5)
        # Assert that all spherical harmonics except when l=0,m=0 are all zero.
        assert_almost_equal(integral[0, 0, :], np.sqrt(4.0 * np.pi))
        assert_almost_equal(integral[0, 1:, :], 0.0)
        assert_almost_equal(integral[1:, :, :], 0.0)
예제 #4
0
파일: atomgrid.py 프로젝트: theochem/grid
        def interpolate_laplacian(points):
            r_pts, theta, phi = self.convert_cartesian_to_spherical(points).T

            # Evaluate the radial components for function f
            r_values_f = np.array([spline(r_pts) for spline in radial_comps_f])
            # Evaluate the second derivatives of splines of radial component of function rf
            r_values_rf = np.array(
                [spline(r_pts, 2) for spline in radial_comps_rf])

            r_sph_harm = generate_real_spherical_harmonics(
                self.l_max // 2, theta, phi)

            # Divide \rho_{lm}^{rf} by r  and \rho_{lm}^{rf} by r^2
            # l means the angular (l,m) variables and n represents points.
            with np.errstate(divide="ignore", invalid="ignore"):
                r_values_rf /= r_pts
                r_values_rf[:, r_pts < 1e-10] = 0.0
                r_values_f /= r_pts**2.0
                r_values_f[:, r_pts**2.0 < 1e-10] = 0.0

            # Multiply \rho_{lm}^f by l(l+1), note 2l+1 orders for each degree l
            degrees = np.hstack([[x * (x + 1)] * (2 * x + 1)
                                 for x in np.arange(0, self.l_max // 2 + 1)])
            second_component = np.einsum("ln,l->ln", r_values_f, degrees)
            # Compute \rho_{lm}^{rf}/r - \rho_{lm}^f l(l + 1) / r^2
            component = r_values_rf - second_component
            # Multiply by spherical harmonics and take the sum
            return np.einsum("ln, ln -> n", component, r_sph_harm)
예제 #5
0
파일: test_utils.py 프로젝트: theochem/grid
 def test_generate_real_sph_harms_integrates_correctly(self):
     """Test generated real spherical harmonics integrates correctly."""
     angular = AngularGrid(degree=7)
     pts = angular.points
     wts = angular.weights
     r = np.linalg.norm(pts, axis=1)
     # polar
     phi = np.arccos(pts[:, 2] / r)
     # azimuthal
     theta = np.arctan2(pts[:, 1], pts[:, 0])
     # generate spherical harmonics
     lmax = 3
     sph_h = generate_real_spherical_harmonics(lmax, theta,
                                               phi)  # l_max = 3
     # Test the shape matches (order, degree, number of points)
     assert sph_h.shape == (1 + 3 + 5 + 7, 26)
     counter = 0
     for l_value in range(0, lmax + 1):
         for m in ([0] + [x for x in range(1, l_value + 1)] +
                   [-x for x in range(1, l_value + 1)]):
             # print(l_value, m)
             re = sum(sph_h[counter, :] * wts)
             if l_value == 0:
                 assert_almost_equal(re, np.sqrt(4.0 * np.pi))
             else:
                 assert_almost_equal(re, 0)
             # no nan in the final result
             assert np.sum(np.isnan(re)) == 0
             counter += 1
예제 #6
0
파일: atomgrid.py 프로젝트: theochem/grid
    def radial_component_splines(self, func_vals: np.ndarray):
        r"""
        Return spline to interpolate radial components wrt to expansion in real spherical harmonics.

        For fixed r, a function :math:`f(r, \theta, \phi)` is projected onto the spherical
        harmonic expansion

        .. math::
            f(r, \theta, \phi) = \sum_{l=0}^\infty \sum_{m=-l}^l \rho^{lm}(r) Y^m_l(\theta, \phi)

        where :math:`Y^m_l` is the real Spherical harmonic of order :math:`l` and degree :math:`m`.
        The radial components :math:`\rho^{lm}(r)` are interpolated using a cubic spline and
        are calculated via integration:

        .. math::
            \rho^{lm}(r) = \int \int f(r, \theta, \phi) Y^m_l(\theta, \phi) \sin(\theta)
             d\theta d\phi.

        Parameters
        ----------
        func_vals : ndarray(N,)
            The function values evaluated on all :math:`N` points on the atomic grid.

        Returns
        -------
        list[scipy.PPoly]
            A list of size :math:`(l_{max}//2 + 1)**2` of  CubicSpline object for interpolating the
            coefficients :math:`\rho^{lm}(r)`. The input of spline is array
            of :math:`N` points on :math:`[0, \infty)` and the output is :\{math:`\rho^{lm}(r_i)\}`.
            The list starts with degrees :math:`l=0,\cdots l_{max}`, and the for each degree,
            the zeroth order spline is first, followed by positive orders then negative.

        """
        if func_vals.size != self.size:
            raise ValueError(
                "The size of values does not match with the size of grid\n"
                f"The size of value array: {func_vals.size}\n"
                f"The size of grid: {self.size}")
        if self._basis is None:
            theta, phi = self.convert_cartesian_to_spherical().T[1:]
            # Going up to `self.l_max // 2` is explained below.
            self._basis = generate_real_spherical_harmonics(
                self.l_max // 2, theta, phi)
        # Multiply spherical harmonic basis with the function values to project.
        values = np.einsum("ln,n->ln", self._basis, func_vals)
        radial_components = self.integrate_angular_coordinates(values)
        # each shell can only integrate upto shell_degree // 2, so if shell_degree < l_max,
        # the f_{lm} should be set to zero for l > shell_degree // 2. Instead, one could set
        # truncate the basis of a given shell.
        for i in range(self.n_shells):
            if self.degrees[i] != self.l_max:
                num_nonzero_sph = (self.degrees[i] // 2 + 1)**2
                radial_components[num_nonzero_sph:, i] = 0.0

        # Return a spline for each spherical harmonic with maximum degree `self.l_max // 2`.
        return [
            CubicSpline(x=self.rgrid.points, y=sph_val)
            for sph_val in radial_components
        ]
예제 #7
0
    def test_integration_of_spherical_harmonic_not_accurate_beyond_degree(self):
        r"""Test integration of spherical harmonic of degree higher than grid is not accurate."""
        grid = AngularGrid(degree=3)
        r = np.linalg.norm(grid.points, axis=1)
        phi = np.arccos(grid.points[:, 2] / r)
        theta = np.arctan2(grid.points[:, 1], grid.points[:, 0])

        sph_harm = generate_real_spherical_harmonics(l_max=6, theta=theta, phi=phi)
        # Check that l=4,m=0 gives inaccurate results
        assert np.abs(grid.integrate(sph_harm[(4) ** 2, :])) > 1e-8
        # Check that l=6,m=0 gives inaccurate results
        assert np.abs(grid.integrate(sph_harm[(6) ** 2, :])) > 1e-8
예제 #8
0
파일: test_utils.py 프로젝트: theochem/grid
def test_derivative_of_spherical_harmonics_with_finite_difference(
        numb_pts, max_degree):
    """Test the derivative value from spherical harmonics with finite difference."""
    theta = np.array([1e-5])
    theta = np.hstack(
        (theta, np.random.uniform(0.0, 2.0 * np.pi, size=(numb_pts, ))))
    phi = np.array([1e-5])
    phi = np.hstack((phi, np.random.uniform(0.0, np.pi, size=(numb_pts, ))))
    eps = 1e-7
    l_max = max_degree
    value = generate_real_spherical_harmonics(l_max, theta, phi)
    value_at_eps_theta = generate_real_spherical_harmonics(
        l_max, theta + eps, phi)
    value_at_eps_phi = generate_real_spherical_harmonics(
        l_max, theta, phi + eps)
    actual_answer = generate_derivative_real_spherical_harmonics(
        l_max, theta, phi)
    deriv_theta = (value_at_eps_theta - value) / eps
    deriv_phi = (value_at_eps_phi - value) / eps

    assert_almost_equal(actual_answer[0, :], deriv_theta, decimal=3)
    assert_almost_equal(actual_answer[1, :], deriv_phi, decimal=3)
예제 #9
0
파일: test_utils.py 프로젝트: theochem/grid
 def test_generate_real_spherical_is_accurate(self):
     r"""Test generated real spherical harmonic up to degree 3 is accurate."""
     numb_pts = 100
     pts = np.random.uniform(-1.0, 1.0, size=(numb_pts, 3))
     sph_pts = convert_cart_to_sph(pts)
     r, theta, phi = sph_pts[:, 0], sph_pts[:, 1], sph_pts[:, 2]
     sph_h = generate_real_spherical_harmonics(3, theta, phi)  # l_max = 3
     # Test l=0, m=0
     assert_allclose(sph_h[0, :],
                     np.ones(len(theta)) / np.sqrt(4.0 * np.pi))
     # Test l=1, m=0, obtained from wikipedia
     assert_allclose(sph_h[1, :],
                     np.sqrt(3.0 / (4.0 * np.pi)) * pts[:, 2] / r)
     # Test l=1, m=1, m=0
     assert_allclose(sph_h[2, :],
                     np.sqrt(3.0 / (4.0 * np.pi)) * pts[:, 0] / r)
     assert_allclose(sph_h[3, :],
                     np.sqrt(3.0 / (4.0 * np.pi)) * pts[:, 1] / r)
예제 #10
0
 def test_integration_of_spherical_harmonic_up_to_degree_three(self):
     r"""Test integration of spherical harmonic up to degree three is accurate."""
     degree = 3
     grid = AngularGrid(degree=100)
     # Concert to spherical coordinates from Cartesian.
     r = np.linalg.norm(grid.points, axis=1)
     phi = np.arccos(grid.points[:, 2] / r)
     theta = np.arctan2(grid.points[:, 1], grid.points[:, 0])
     # Generate All Spherical Harmonics Up To Degree = 3
     #   Returns a three dimensional array where [order m, degree l, points]
     sph_harm = generate_real_spherical_harmonics(degree, theta, phi)
     for l_deg in range(0, 4):
         for m_ord in range(-l_deg, l_deg):
             sph_harm_one = sph_harm[l_deg**2 : (l_deg + 1) ** 2, :]
             if l_deg == 0 and m_ord == 0:
                 actual = np.sqrt(4.0 * np.pi)
                 assert_equal(actual, grid.integrate(sph_harm_one[m_ord, :]))
             else:
                 assert_almost_equal(0.0, grid.integrate(sph_harm_one[m_ord, :]))
예제 #11
0
    def test_fitting_spherical_harmonics(self):
        r"""Test fitting the radial components of spherical harmonics is just a one, rest zeros."""
        max_degree = 10  # Maximum degree
        rad = OneDGrid(np.linspace(0.0, 1.0, num=10), np.ones(10), (0, np.inf))
        atom_grid = AtomGrid(rad, degrees=[max_degree])
        max_degree = atom_grid.l_max
        spherical = atom_grid.convert_cartesian_to_spherical()
        # Evaluate all spherical harmonics on the atomic grid points (r_i, theta_j, phi_j).
        spherical_harmonics = generate_real_spherical_harmonics(
            max_degree,
            spherical[:, 1],
            spherical[:, 2]  # theta, phi points
        )

        i = 0
        # Go through each spherical harmonic up to max_degree // 2 and check if projection
        # for its radial component is one and the rest are all zeros.
        for l_value in range(0, max_degree // 2):
            for m in ([0] + [x for x in range(1, l_value + 1)] +
                      [x for x in range(-l_value, 0)]):
                spherical_harm = spherical_harmonics[i, :]
                radial_components = atom_grid.radial_component_splines(
                    spherical_harm)
                assert len(radial_components) == (atom_grid.l_max // 2 +
                                                  1)**2.0

                radial_pts = np.arange(0.0, 1.0, 0.01)
                # Go Through each (l, m)
                for j, radial_comp in enumerate(radial_components):
                    # If the current index is j, then projection of spherical harmonic
                    # onto itself should be all ones, else they are all zero.
                    if i == j:
                        assert_almost_equal(radial_comp(radial_pts), 1.0)
                    else:
                        assert_almost_equal(radial_comp(radial_pts),
                                            0.0,
                                            decimal=8)
                i += 1
예제 #12
0
    def test_value_fitting(self):
        """Test spline projection the same as spherical harmonics."""
        odg = OneDGrid(np.arange(10) + 1, np.ones(10), (0, np.inf))
        rad = IdentityRTransform().transform_1d_grid(odg)
        atgrid = AtomGrid.from_pruned(rad, 1, sectors_r=[], sectors_degree=[7])
        values = self.helper_func_power(atgrid.points)
        spls = atgrid.radial_component_splines(values)
        assert len(spls) == 16

        for shell in range(1, 11):
            sh_grid = atgrid.get_shell_grid(shell - 1, r_sq=False)
            spherical_pts = atgrid.convert_cartesian_to_spherical(
                sh_grid.points)
            _, theta, phi = spherical_pts.T
            l_max = atgrid.l_max // 2
            r_sph = generate_real_spherical_harmonics(l_max, theta, phi)
            r_sph_proj = np.sum(
                r_sph * self.helper_func_power(sh_grid.points) *
                sh_grid.weights,
                axis=-1,
            )
            assert_allclose(r_sph_proj, [spl(shell) for spl in spls],
                            atol=1e-10)
예제 #13
0
 def test_orthogonality_of_spherical_harmonic_up_to_degree_three(self):
     r"""Test orthogonality of spherical harmonic up to degree 3 is accurate."""
     degree = 3
     grid = AngularGrid(degree=10)
     # Concert to spherical coordinates from Cartesian.
     r = np.linalg.norm(grid.points, axis=1)
     phi = np.arccos(grid.points[:, 2] / r)
     theta = np.arctan2(grid.points[:, 1], grid.points[:, 0])
     # Generate All Spherical Harmonics Up To Degree = 3
     #   Returns a three dimensional array where [order m, degree l, points]
     sph_harm = generate_real_spherical_harmonics(degree, theta, phi)
     for l_deg in range(0, 4):
         for m_ord in range(-l_deg, l_deg + 1):
             for l2 in range(0, 4):
                 for m2 in range(-l2, l2 + 1):
                     sph_harm_one = sph_harm[l_deg**2 : (l_deg + 1) ** 2, :]
                     sph_harm_two = sph_harm[l2**2 : (l2 + 1) ** 2, :]
                     integral = grid.integrate(
                         sph_harm_one[m_ord, :] * sph_harm_two[m2, :]
                     )
                     if l2 != l_deg or m2 != m_ord:
                         assert np.abs(integral) < 1e-8
                     else:
                         assert np.abs(integral - 1.0) < 1e-8
예제 #14
0
파일: atomgrid.py 프로젝트: theochem/grid
        def interpolate_low(points,
                            deriv=0,
                            deriv_spherical=False,
                            only_radial_derivs=False):
            r"""Construct a spline like callable for intepolation.

            Parameters
            ----------
            points : ndarray(N, 3)
                Cartesian coordinates of :math:`N` points to evaluate the splines on.
            deriv : int, optional
                If deriv is zero, then only returns function values. If it is one, then returns
                the first derivative of the interpolated function with respect to either Cartesian
                or spherical coordinates. Only higher-order derivatives (`deriv` =2,3) are supported
                for the derivatives wrt to radial components. `deriv=3` only returns a constant.
            deriv_spherical : bool
                If True, then returns the derivatives with respect to spherical coordinates
                :math:`(r, \theta, \phi)`. Default False.
            only_radial_derivs : bool
                If true, then the derivative wrt to radius :math:`r` is returned.

            Returns
            -------
            ndarray(M,...) :
                The interpolated function values or its derivatives with respect to Cartesian
                :math:`(x,y,z)` or if `deriv_spherical` then :math:`(r, \theta, \phi)` or
                if `only_radial_derivs` then derivative wrt to :math:`r` is only returned.

            """
            if deriv_spherical and only_radial_derivs:
                warnings.warn(
                    "Since `only_radial_derivs` is true, then only the derivative wrt to"
                    "radius is returned and `deriv_spherical` value is ignored."
                )
            r_pts, theta, phi = self.convert_cartesian_to_spherical(points).T
            r_values = np.array([spline(r_pts, deriv) for spline in splines])
            r_sph_harm = generate_real_spherical_harmonics(
                self.l_max // 2, theta, phi)

            # If theta, phi derivaitves are wanted and the derivative is first-order.
            if not only_radial_derivs and deriv == 1:
                # Calculate derivative of f with respect to radial, theta, phi
                # Get derivative of spherical harmonics first.
                radial_components = np.array(
                    [spline(r_pts, 0) for spline in splines])
                deriv_sph_harm = generate_derivative_real_spherical_harmonics(
                    self.l_max // 2, theta, phi)
                deriv_r = np.einsum("ij, ij -> j", r_values, r_sph_harm)
                deriv_theta = np.einsum("ij,ij->j", radial_components,
                                        deriv_sph_harm[0, :, :])
                deriv_phi = np.einsum("ij,ij->j", radial_components,
                                      deriv_sph_harm[1, :, :])

                # If deriv spherical is wanted, then return that.
                if deriv_spherical:
                    return np.hstack((deriv_r, deriv_theta, deriv_phi))

                # Convert derivative from spherical to Cartesian:
                derivs = np.zeros((len(r_pts), 3))
                # TODO: this could be vectorized properly with memory management.
                for i_pt in range(0, len(r_pts)):
                    radial_i, theta_i, phi_i = r_pts[i_pt], theta[i_pt], phi[
                        i_pt]
                    derivs[
                        i_pt] = convert_derivative_from_spherical_to_cartesian(
                            deriv_r[i_pt],
                            deriv_theta[i_pt],
                            deriv_phi[i_pt],
                            radial_i,
                            theta_i,
                            phi_i,
                        )
                return derivs
            elif not only_radial_derivs and deriv != 0:
                raise ValueError(
                    f"Higher order derivatives are only supported for derivatives"
                    f"with respect to the radius. Deriv is {deriv}.")

            return np.einsum("ij, ij -> j", r_values, r_sph_harm)