def calculate_ed_data(
        self,
        structure,
        reciprocal_radius,
        rotation=(0, 0, 0),
        with_direct_beam=True,
        max_excitation_error=1e-2,
        shape_factor_width=None,
        debye_waller_factors={},
    ):
        """Calculates the Electron Diffraction data for a structure.

        Parameters
        ----------
        structure : diffpy.structure.structure.Structure
            The structure for which to derive the diffraction pattern.
            Note that the structure must be rotated to the appropriate
            orientation and that testing is conducted on unit cells
            (rather than supercells).
        reciprocal_radius : float
            The maximum radius of the sphere of reciprocal space to
            sample, in reciprocal Angstroms.
        rotation : tuple
            Euler angles, in degrees, in the rzxz convention. Default is
            (0, 0, 0) which aligns 'z' with the electron beam.
        with_direct_beam : bool
            If True, the direct beam is included in the simulated
            diffraction pattern. If False, it is not.
        max_excitation_error : float
            The cut-off for geometric excitation error in the z-direction
            in units of reciprocal Angstroms. Spots with a larger distance
            from the Ewald sphere are removed from the pattern.
            Related to the extinction distance and roungly equal to 1/thickness.
        shape_factor_width : float
            Determines the width of the reciprocal rel-rod, for fine-grained
            control. If not set will be set equal to max_excitation_error.
        debye_waller_factors : dict of str:value pairs
            Maps element names to their temperature-dependent Debye-Waller factors.

        Returns
        -------
        diffsims.sims.diffraction_simulation.DiffractionSimulation
            The data associated with this structure and diffraction setup.
        """
        # Specify variables used in calculation
        wavelength = self.wavelength
        latt = structure.lattice

        # Obtain crystallographic reciprocal lattice points within `reciprocal_radius` and
        # g-vector magnitudes for intensity calculations.
        recip_latt = latt.reciprocal()
        g_indices, cartesian_coordinates, g_hkls = get_points_in_sphere(
            recip_latt, reciprocal_radius
        )

        ai, aj, ak = (
            np.deg2rad(rotation[0]),
            np.deg2rad(rotation[1]),
            np.deg2rad(rotation[2]),
        )
        R = euler2mat(ai, aj, ak, axes="rzxz")
        cartesian_coordinates = np.matmul(R, cartesian_coordinates.T).T

        # Identify the excitation errors of candidate points
        r_sphere = 1 / wavelength
        r_spot = np.sqrt(np.sum(np.square(cartesian_coordinates[:, :2]), axis=1))
        z_spot = cartesian_coordinates[:, 2]

        z_sphere = -np.sqrt(r_sphere ** 2 - r_spot ** 2) + r_sphere
        excitation_error = z_sphere - z_spot

        # determine the pre-selection reflections
        if self.precession_angle == 0:
            intersection = np.abs(excitation_error) < max_excitation_error
        else:
            # only consider points that intersect the ewald sphere at some point
            # the center point of the sphere
            P_z = r_sphere * np.cos(np.deg2rad(self.precession_angle))
            P_t = r_sphere * np.sin(np.deg2rad(self.precession_angle))
            # the extremes of the ewald sphere
            z_surf_up = P_z - np.sqrt(r_sphere ** 2 - (r_spot + P_t) ** 2)
            z_surf_do = P_z - np.sqrt(r_sphere ** 2 - (r_spot - P_t) ** 2)
            intersection = (z_spot - max_excitation_error <= z_surf_up) & (
                z_spot + max_excitation_error >= z_surf_do
            )

        # select these reflections
        intersection_coordinates = cartesian_coordinates[intersection]
        excitation_error = excitation_error[intersection]
        r_spot = r_spot[intersection]
        g_indices = g_indices[intersection]
        g_hkls = g_hkls[intersection]

        if shape_factor_width is None:
            shape_factor_width = max_excitation_error
        # select and evaluate shape factor model
        if self.precession_angle == 0:
            # calculate shape factor
            shape_factor = self.shape_factor_model(
                excitation_error, shape_factor_width, **self.shape_factor_kwargs
            )
        else:
            if self.approximate_precession:
                shape_factor = lorentzian_precession(
                    excitation_error,
                    shape_factor_width,
                    r_spot,
                    np.deg2rad(self.precession_angle),
                )
            else:
                shape_factor = _shape_factor_precession(
                    excitation_error,
                    r_spot,
                    np.deg2rad(self.precession_angle),
                    self.shape_factor_model,
                    shape_factor_width,
                    **self.shape_factor_kwargs,
                )
        # Calculate diffracted intensities based on a kinematical model.
        intensities = get_kinematical_intensities(
            structure,
            g_indices,
            g_hkls,
            prefactor=shape_factor,
            scattering_params=self.scattering_params,
            debye_waller_factors=debye_waller_factors,
        )

        # Threshold peaks included in simulation as factor of maximum intensity.
        peak_mask = intensities > np.max(intensities) * self.minimum_intensity
        intensities = intensities[peak_mask]
        intersection_coordinates = intersection_coordinates[peak_mask]
        g_indices = g_indices[peak_mask]

        return DiffractionSimulation(
            coordinates=intersection_coordinates,
            indices=g_indices,
            intensities=intensities,
            with_direct_beam=with_direct_beam,
        )
    def calculate_profile_data(
        self,
        structure,
        reciprocal_radius=1.0,
        minimum_intensity=1e-3,
        debye_waller_factors={},
    ):
        """Calculates a one dimensional diffraction profile for a
        structure.

        Parameters
        ----------
        structure : diffpy.structure.structure.Structure
            The structure for which to calculate the diffraction profile.
        reciprocal_radius : float
            The maximum radius of the sphere of reciprocal space to
            sample, in reciprocal angstroms.
        minimum_intensity : float
            The minimum intensity required for a diffraction peak to be
            considered real. Deals with numerical precision issues.
        debye_waller_factors : dict of str:value pairs
            Maps element names to their temperature-dependent Debye-Waller factors.

        Returns
        -------
        diffsims.sims.diffraction_simulation.ProfileSimulation
            The diffraction profile corresponding to this structure and
            experimental conditions.
        """
        wavelength = self.wavelength
        latt = structure.lattice

        # Obtain crystallographic reciprocal lattice points within range
        recip_latt = latt.reciprocal()
        spot_indices, _, spot_distances = get_points_in_sphere(
            recip_latt, reciprocal_radius
        )

        ##spot_indicies is a numpy.array of the hkls allowd in the recip radius
        g_indices, multiplicities, g_hkls = get_intensities_params(
            recip_latt, reciprocal_radius
        )

        i_hkl = get_kinematical_intensities(
            structure,
            g_indices,
            np.asarray(g_hkls),
            prefactor=multiplicities,
            scattering_params=self.scattering_params,
            debye_waller_factors=debye_waller_factors,
        )

        if is_lattice_hexagonal(latt):
            # Use Miller-Bravais indices for hexagonal lattices.
            g_indices = (
                g_indices[0],
                g_indices[1],
                -g_indices[0] - g_indices[1],
                g_indices[2],
            )

        hkls_labels = ["".join([str(int(x)) for x in xs]) for xs in g_indices]

        peaks = {}
        for l, i, g in zip(hkls_labels, i_hkl, g_hkls):
            peaks[l] = [i, g]

        # Scale intensities so that the max intensity is 100.

        max_intensity = max([v[0] for v in peaks.values()])
        x = []
        y = []
        hkls = []
        for k in peaks.keys():
            v = peaks[k]
            if v[0] / max_intensity * 100 > minimum_intensity and (k != "000"):
                x.append(v[1])
                y.append(v[0])
                hkls.append(k)

        y = np.asarray(y) / max(y) * 100

        return ProfileSimulation(x, y, hkls)
    def calculate_profile_data(self,
                               structure,
                               reciprocal_radius=1.0,
                               magnitude_tolerance=1e-5,
                               minimum_intensity=1e-3):
        """
        Calculates a one dimensional diffraction profile for a structure.

        Parameters
        ----------
        structure : Structure
            The structure for which to calculate the diffraction profile.
        reciprocal_radius : float
            The maximum radius of the sphere of reciprocal space to sample, in
            reciprocal angstroms.
        magnitude_tolerance : float
            The minimum difference between diffraction magnitudes in reciprocal
            angstroms for two peaks to be consdiered different.
        minimum_intensity : float
            The minimum intensity required for a diffraction peak to be
            considered real. Deals with numerical precision issues.

        Returns
        -------
        diffsims.ProfileSimulation
            The diffraction profile corresponding to this structure and
            experimental conditions.
        """
        max_r = reciprocal_radius
        wavelength = self.wavelength
        scattering_params = self.scattering_params

        latt = structure.lattice
        is_hex = is_lattice_hexagonal(latt)

        coeffs, fcoords, occus, dwfactors = get_vectorized_list_for_atomic_scattering_factors(
            structure, {}, scattering_params=scattering_params)

        # Obtain crystallographic reciprocal lattice points within range
        recip_latt = latt.reciprocal()
        spot_indicies, _, spot_distances = get_points_in_sphere(
            recip_latt, reciprocal_radius)

        peaks = {}
        mask = np.logical_not((np.any(spot_indicies, axis=1) == 0))

        for hkl, g_hkl in zip(spot_indicies[mask], spot_distances[mask]):
            # Force miller indices to be integers.
            hkl = [int(round(i)) for i in hkl]

            d_hkl = 1 / g_hkl

            # Bragg condition
            # theta = asin(wavelength * g_hkl / 2)

            # s = sin(theta) / wavelength = 1 / 2d = |ghkl| / 2 (d =
            # 1/|ghkl|)
            s = g_hkl / 2

            # Store s^2 since we are using it a few times.
            s2 = s**2

            # Vectorized computation of g.r for all fractional coords and
            # hkl.
            g_dot_r = np.dot(fcoords, np.transpose([hkl])).T[0]

            # Highly vectorized computation of atomic scattering factors.
            fs = np.sum(coeffs[:, :, 0] * np.exp(-coeffs[:, :, 1] * s2),
                        axis=1)

            dw_correction = np.exp(-dwfactors * s2)

            # Structure factor = sum of atomic scattering factors (with
            # position factor exp(2j * pi * g.r and occupancies).
            # Vectorized computation.
            f_hkl = np.sum(fs * occus * np.exp(2j * np.pi * g_dot_r) *
                           dw_correction)

            # Intensity for hkl is modulus square of structure factor.
            i_hkl = (f_hkl * f_hkl.conjugate()).real

            # two_theta = degrees(2 * theta)

            if is_hex:
                # Use Miller-Bravais indices for hexagonal lattices.
                hkl = (hkl[0], hkl[1], -hkl[0] - hkl[1], hkl[2])

            peaks[g_hkl] = [i_hkl, [tuple(hkl)], d_hkl]

        # Scale intensities so that the max intensity is 100.
        max_intensity = max([v[0] for v in peaks.values()])
        x = []
        y = []
        hkls = []
        d_hkls = []
        for k in sorted(peaks.keys()):
            v = peaks[k]
            fam = get_unique_families(v[1])
            if v[0] / max_intensity * 100 > minimum_intensity:
                x.append(k)
                y.append(v[0])
                hkls.append(fam)
                d_hkls.append(v[2])

        y = y / max(y) * 100

        return ProfileSimulation(x, y, hkls)
    def calculate_ed_data(self,
                          structure,
                          reciprocal_radius,
                          rotation=(0, 0, 0),
                          with_direct_beam=True):
        """Calculates the Electron Diffraction data for a structure.

        Parameters
        ----------
        structure : Structure
            The structure for which to derive the diffraction pattern. Note that
            the structure must be rotated to the appropriate orientation and
            that testing is conducted on unit cells (rather than supercells).
        reciprocal_radius : float
            The maximum radius of the sphere of reciprocal space to sample, in
            reciprocal angstroms.
        rotation : tuple
            Euler angles, in degrees, in the rzxz convention. Default is (0,0,0)
            which aligns 'z' with the electron beam
        with_direct_beam : bool
            If True, the direct beam is included in the simulated diffraction
            pattern. If False, it is not.

        Returns
        -------
        diffsims.DiffractionSimulation
            The data associated with this structure and diffraction setup.

        """
        # Specify variables used in calculation
        wavelength = self.wavelength
        max_excitation_error = self.max_excitation_error
        debye_waller_factors = self.debye_waller_factors
        latt = structure.lattice
        scattering_params = self.scattering_params

        # Obtain crystallographic reciprocal lattice points within `max_r` and
        # g-vector magnitudes for intensity calculations.
        recip_latt = latt.reciprocal()
        spot_indicies, cartesian_coordinates, spot_distances = get_points_in_sphere(
            recip_latt, reciprocal_radius)

        ai, aj, ak = np.deg2rad(rotation[0]), np.deg2rad(
            rotation[1]), np.deg2rad(rotation[2])
        R = euler2mat(ai, aj, ak, axes='rzxz')
        cartesian_coordinates = np.matmul(R, cartesian_coordinates.T).T

        # Identify points intersecting the Ewald sphere within maximum
        # excitation error and store the magnitude of their excitation error.
        r_sphere = 1 / wavelength
        r_spot = np.sqrt(
            np.sum(np.square(cartesian_coordinates[:, :2]), axis=1))
        z_sphere = -np.sqrt(r_sphere**2 - r_spot**2) + r_sphere
        proximity = np.absolute(z_sphere - cartesian_coordinates[:, 2])
        intersection = proximity < max_excitation_error
        # Mask parameters corresponding to excited reflections.
        intersection_coordinates = cartesian_coordinates[intersection]
        intersection_indices = spot_indicies[intersection]
        proximity = proximity[intersection]
        g_hkls = spot_distances[intersection]

        # Calculate diffracted intensities based on a kinematical model.
        intensities = get_kinematical_intensities(
            structure, intersection_indices, g_hkls, proximity,
            max_excitation_error, debye_waller_factors, scattering_params)

        # Threshold peaks included in simulation based on minimum intensity.
        peak_mask = intensities > 1e-20
        intensities = intensities[peak_mask]
        intersection_coordinates = intersection_coordinates[peak_mask]
        intersection_indices = intersection_indices[peak_mask]

        return DiffractionSimulation(coordinates=intersection_coordinates,
                                     indices=intersection_indices,
                                     intensities=intensities,
                                     with_direct_beam=with_direct_beam)
    def calculate_ed_data(self,
                          structure,
                          reciprocal_radius,
                          rotation=(0, 0, 0),
                          shape_factor_model="linear",
                          max_excitation_error=1e-2,
                          with_direct_beam=True,
                          **kwargs):
        """Calculates the Electron Diffraction data for a structure.

        Parameters
        ----------
        structure : Structure
            The structure for which to derive the diffraction pattern. Note that
            the structure must be rotated to the appropriate orientation and
            that testing is conducted on unit cells (rather than supercells).
        reciprocal_radius : float
            The maximum radius of the sphere of reciprocal space to sample, in
            reciprocal angstroms.
        rotation : tuple
            Euler angles, in degrees, in the rzxz convention. Default is (0,0,0)
            which aligns 'z' with the electron beam
        shape_factor_model : function or str
            a function that takes excitation_error and max_excitation_error (and potentially **kwargs) and returns an intensity
            scaling factor. The code provides "linear" and "binary" options accessed with by parsing the associated strings
        max_excitation_error : float
            the exctinction distance for reflections, in reciprocal Angstroms
        with_direct_beam : bool
            If True, the direct beam is included in the simulated diffraction
            pattern. If False, it is not.
        **kwargs :
            passed to shape_factor_model

        Returns
        -------
        diffsims.DiffractionSimulation
            The data associated with this structure and diffraction setup.

        """
        # Specify variables used in calculation
        wavelength = self.wavelength
        latt = structure.lattice

        # Obtain crystallographic reciprocal lattice points within `reciprocal_radius` and
        # g-vector magnitudes for intensity calculations.
        recip_latt = latt.reciprocal()
        spot_indices, cartesian_coordinates, spot_distances = get_points_in_sphere(
            recip_latt, reciprocal_radius)

        ai, aj, ak = (
            np.deg2rad(rotation[0]),
            np.deg2rad(rotation[1]),
            np.deg2rad(rotation[2]),
        )
        R = euler2mat(ai, aj, ak, axes="rzxz")
        cartesian_coordinates = np.matmul(R, cartesian_coordinates.T).T

        # Identify points intersecting the Ewald sphere within maximum
        # excitation error and store the magnitude of their excitation error.
        r_sphere = 1 / wavelength
        r_spot = np.sqrt(
            np.sum(np.square(cartesian_coordinates[:, :2]), axis=1))
        z_sphere = -np.sqrt(r_sphere**2 - r_spot**2) + r_sphere
        excitation_error = np.absolute(z_sphere - cartesian_coordinates[:, 2])
        intersection = excitation_error < max_excitation_error
        # Mask parameters corresponding to excited reflections.
        intersection_coordinates = cartesian_coordinates[intersection]
        g_indices = spot_indices[intersection]
        excitation_error = excitation_error[intersection]
        g_hkls = spot_distances[intersection]

        if shape_factor_model == "linear":
            shape_factor = 1 - (excitation_error / max_excitation_error)
        elif shape_factor_model == "binary":
            shape_factor = 1
        else:
            shape_factor = shape_factor_model(excitation_error,
                                              max_excitation_error, **kwargs)

        # Calculate diffracted intensities based on a kinematical model.
        intensities = get_kinematical_intensities(
            structure,
            g_indices,
            g_hkls,
            prefactor=shape_factor,
            scattering_params=self.scattering_params,
            debye_waller_factors=self.debye_waller_factors,
        )

        # Threshold peaks included in simulation based on minimum intensity.
        peak_mask = intensities > 1e-20
        intensities = intensities[peak_mask]
        intersection_coordinates = intersection_coordinates[peak_mask]
        g_indices = g_indices[peak_mask]

        return DiffractionSimulation(
            coordinates=intersection_coordinates,
            indices=g_indices,
            intensities=intensities,
            with_direct_beam=with_direct_beam,
        )
Beispiel #6
0
def test_get_points_in_sphere():
    latt = diffpy.structure.lattice.Lattice(0.5, 0.5, 0.5, 90, 90, 90)
    ind, cord, dist = get_points_in_sphere(latt, 0.6)
    assert len(ind) == len(cord)
    assert len(ind) == len(dist)
    assert len(dist) == 1 + 6
def _generate_lookup_table(recip_latt,
                           reciprocal_radius: float,
                           unique: bool = True):
    """Generate a look-up table with all combinations of indices,
    including their reciprocal distances and the angle between
    them.

    Parameters
    ----------
    recip_latt : :class:`diffpy.structure.lattice.Lattice`
        Reciprocal lattice
    reciprocal_radius : float
        The maximum g-vector magnitude to be included in the library.
    unique : bool
        Return a unique list of phase measurements

    Returns
    -------
    indices : np.array
        Nx2x3 numpy array containing the miller indices for
        reflection1, reflection2
    measurements : np.array
        Nx3 numpy array containing len1, len2, angle

    """
    miller_indices, coordinates, distances = get_points_in_sphere(
        recip_latt, reciprocal_radius)

    # Create pair_indices for selecting all point pair combinations
    num_indices = len(miller_indices)
    pair_a_indices, pair_b_indices = np.mgrid[:num_indices, :num_indices]

    # Only select one of the permutations and don't pair an index with
    # itself (select above diagonal)
    upper_indices = np.triu_indices(num_indices, 1)
    pair_a_indices = pair_a_indices[upper_indices].ravel()
    pair_b_indices = pair_b_indices[upper_indices].ravel()

    # Mask off origin (0, 0, 0)
    origin_index = num_indices // 2
    pair_a_indices = pair_a_indices[pair_a_indices != origin_index]
    pair_b_indices = pair_b_indices[pair_b_indices != origin_index]

    pair_indices = np.vstack([pair_a_indices, pair_b_indices])

    # Create library entries
    angles = get_angle_cartesian_vec(coordinates[pair_a_indices],
                                     coordinates[pair_b_indices])
    pair_distances = distances[pair_indices.T]
    # Ensure longest vector is first
    len_sort = np.fliplr(pair_distances.argsort(axis=1))
    # phase_index_pairs is a list of [hkl1, hkl2]
    phase_index_pairs = np.take_along_axis(miller_indices[pair_indices.T],
                                           len_sort[:, :, np.newaxis],
                                           axis=1)
    # phase_measurements is a list of [len1, len2, angle]
    phase_measurements = np.column_stack((np.take_along_axis(pair_distances,
                                                             len_sort,
                                                             axis=1), angles))

    if unique:
        # Only keep unique triplets
        measurements, measurement_indices = np.unique(phase_measurements,
                                                      axis=0,
                                                      return_index=True)
        indices = phase_index_pairs[measurement_indices]
    else:
        measurements = phase_measurements
        indices = phase_index_pairs

    return measurements, indices
Beispiel #8
0
    def calculate_ed_data(
        self,
        structure,
        reciprocal_radius,
        rotation=(0, 0, 0),
        with_direct_beam=True,
        max_excitation_error=1e-2,
        debye_waller_factors={},
    ):
        """Calculates the Electron Diffraction data for a structure.

        Parameters
        ----------
        structure : diffpy.structure.structure.Structure
            The structure for which to derive the diffraction pattern.
            Note that the structure must be rotated to the appropriate
            orientation and that testing is conducted on unit cells
            (rather than supercells).
        reciprocal_radius : float
            The maximum radius of the sphere of reciprocal space to
            sample, in reciprocal Angstroms.
        rotation : tuple
            Euler angles, in degrees, in the rzxz convention. Default is
            (0, 0, 0) which aligns 'z' with the electron beam.
        with_direct_beam : bool
            If True, the direct beam is included in the simulated
            diffraction pattern. If False, it is not.
        max_excitation_error : float
            The extinction distance for reflections, in reciprocal
            Angstroms. Roughly equal to 1/thickness.
        debye_waller_factors : dict of str:value pairs
            Maps element names to their temperature-dependent Debye-Waller factors.

        Returns
        -------
        diffsims.sims.diffraction_simulation.DiffractionSimulation
            The data associated with this structure and diffraction setup.
        """
        # Specify variables used in calculation
        wavelength = self.wavelength
        latt = structure.lattice

        # Obtain crystallographic reciprocal lattice points within `reciprocal_radius` and
        # g-vector magnitudes for intensity calculations.
        recip_latt = latt.reciprocal()
        g_indices, cartesian_coordinates, g_hkls = get_points_in_sphere(
            recip_latt, reciprocal_radius)

        ai, aj, ak = (
            np.deg2rad(rotation[0]),
            np.deg2rad(rotation[1]),
            np.deg2rad(rotation[2]),
        )
        R = euler2mat(ai, aj, ak, axes="rzxz")
        cartesian_coordinates = np.matmul(R, cartesian_coordinates.T).T

        # Identify the excitation errors of candidate points
        r_sphere = 1 / wavelength
        r_spot = np.sqrt(
            np.sum(np.square(cartesian_coordinates[:, :2]), axis=1))
        z_spot = cartesian_coordinates[:, 2]

        z_sphere = -np.sqrt(r_sphere**2 - r_spot**2) + r_sphere
        excitation_error = z_sphere - z_spot

        if self.precession_angle == 0:
            # Mask parameters corresponding to excited reflections.
            intersection = np.abs(excitation_error) < max_excitation_error
            intersection_coordinates = cartesian_coordinates[intersection]
            excitation_error = excitation_error[intersection]
            r_spot = r_spot[intersection]
            g_indices = g_indices[intersection]
            g_hkls = g_hkls[intersection]

            # calculate shape factor
            shape_factor = self.shape_factor_model(excitation_error,
                                                   max_excitation_error,
                                                   **self.shape_factor_kwargs)
        else:
            intersection_coordinates = cartesian_coordinates  #for naming simplicity

            if self.approximate_precession:
                shape_factor = lorentzian_precession(
                    excitation_error,
                    max_excitation_error,
                    r_spot,
                    np.deg2rad(self.precession_angle),
                )
            else:
                shape_factor = _shape_factor_precession(
                    excitation_error,
                    r_spot,
                    np.deg2rad(self.precession_angle),
                    self.shape_factor_model,
                    max_excitation_error,
                    **self.shape_factor_kwargs,
                )
        # Calculate diffracted intensities based on a kinematical model.
        intensities = get_kinematical_intensities(
            structure,
            g_indices,
            g_hkls,
            prefactor=shape_factor,
            scattering_params=self.scattering_params,
            debye_waller_factors=debye_waller_factors,
        )

        # Threshold peaks included in simulation based on minimum intensity.
        peak_mask = intensities > np.max(intensities) * self.minimum_intensity
        intensities = intensities[peak_mask]
        intersection_coordinates = intersection_coordinates[peak_mask]
        g_indices = g_indices[peak_mask]

        return DiffractionSimulation(
            coordinates=intersection_coordinates,
            indices=g_indices,
            intensities=intensities,
            with_direct_beam=with_direct_beam,
        )
Beispiel #9
0
    def get_vector_library(self, reciprocal_radius):
        """Calculates a library of diffraction vectors and pairwise inter-vector
        angles for a library of crystal structures.

        Parameters
        ----------
        reciprocal_radius : float
            The maximum g-vector magnitude to be included in the library.

        Returns
        -------
        vector_library : :class:`DiffractionVectorLibrary`
            Mapping of phase identifier to phase information in dictionary
            format.
        """
        # Define DiffractionVectorLibrary object to contain results
        vector_library = DiffractionVectorLibrary()
        # Get structures from structure library
        structure_library = self.structures.struct_lib
        # Iterate through phases in library.
        for phase_name in structure_library.keys():
            # Get diffpy.structure object associated with phase
            structure = structure_library[phase_name][0]
            # Get reciprocal lattice points within reciprocal_radius
            recip_latt = structure.lattice.reciprocal()
            miller_indices, coordinates, distances = get_points_in_sphere(
                recip_latt, reciprocal_radius)

            # Create pair_indices for selecting all point pair combinations
            num_indices = len(miller_indices)
            pair_a_indices, pair_b_indices = np.mgrid[:num_indices, :
                                                      num_indices]

            # Only select one of the permutations and don't pair an index with
            # itself (select above diagonal)
            upper_indices = np.triu_indices(num_indices, 1)
            pair_a_indices = pair_a_indices[upper_indices].ravel()
            pair_b_indices = pair_b_indices[upper_indices].ravel()

            # Mask off origin (0, 0, 0)
            origin_index = num_indices // 2
            pair_a_indices = pair_a_indices[pair_a_indices != origin_index]
            pair_b_indices = pair_b_indices[pair_b_indices != origin_index]

            pair_indices = np.vstack([pair_a_indices, pair_b_indices])

            # Create library entries
            angles = get_angle_cartesian_vec(coordinates[pair_a_indices],
                                             coordinates[pair_b_indices])
            pair_distances = distances[pair_indices.T]
            # Ensure longest vector is first
            len_sort = np.fliplr(pair_distances.argsort(axis=1))
            # phase_index_pairs is a list of [hkl1, hkl2]
            phase_index_pairs = np.take_along_axis(
                miller_indices[pair_indices.T],
                len_sort[:, :, np.newaxis],
                axis=1)
            # phase_measurements is a list of [len1, len2, angle]
            phase_measurements = np.column_stack(
                (np.take_along_axis(pair_distances, len_sort, axis=1), angles))

            # Only keep unique triplets
            unique_measurements, unique_measurement_indices = np.unique(
                phase_measurements, axis=0, return_index=True)
            vector_library[phase_name] = {
                'indices': phase_index_pairs[unique_measurement_indices],
                'measurements': unique_measurements
            }

        # Pass attributes to diffraction library from structure library.
        vector_library.identifiers = self.structures.identifiers
        vector_library.structures = self.structures.structures
        vector_library.reciprocal_radius = reciprocal_radius

        return vector_library