def from_data_mapper_and_regularization(
        cls,
        visibilities: vis.Visibilities,
        noise_map: vis.VisibilitiesNoiseMap,
        transformer: trans.TransformerNUFFT,
        mapper: typing.Union[mappers.MapperRectangular, mappers.MapperVoronoi],
        regularization: reg.Regularization,
        settings=SettingsInversion(),
    ):

        transformed_mapping_matrix = transformer.transformed_mapping_matrix_from_mapping_matrix(
            mapping_matrix=mapper.mapping_matrix
        )

        data_vector = inversion_util.data_vector_via_transformed_mapping_matrix_from(
            transformed_mapping_matrix=transformed_mapping_matrix,
            visibilities=visibilities,
            noise_map=noise_map,
        )

        real_curvature_matrix = inversion_util.curvature_matrix_via_mapping_matrix_from(
            mapping_matrix=transformed_mapping_matrix.real, noise_map=noise_map.real
        )

        imag_curvature_matrix = inversion_util.curvature_matrix_via_mapping_matrix_from(
            mapping_matrix=transformed_mapping_matrix.imag, noise_map=noise_map.imag
        )

        regularization_matrix = regularization.regularization_matrix_from_mapper(
            mapper=mapper
        )

        curvature_matrix = np.add(real_curvature_matrix, imag_curvature_matrix)
        curvature_reg_matrix = np.add(curvature_matrix, regularization_matrix)

        try:
            values = np.linalg.solve(curvature_reg_matrix, data_vector)
        except np.linalg.LinAlgError:
            raise exc.InversionException()

        if settings.check_solution:
            if np.isclose(a=values[0], b=values[1], atol=1e-4).all():
                if np.isclose(a=values[0], b=values, atol=1e-4).all():
                    raise exc.InversionException()

        return InversionInterferometerMatrix(
            visibilities=visibilities,
            noise_map=noise_map,
            transformer=transformer,
            mapper=mapper,
            curvature_matrix=curvature_matrix,
            regularization=regularization,
            transformed_mapping_matrix=transformed_mapping_matrix,
            regularization_matrix=regularization_matrix,
            curvature_reg_matrix=curvature_reg_matrix,
            reconstruction=values,
            settings=settings,
        )
Exemple #2
0
    def from_data_mapper_and_regularization(cls, image, noise_map, convolver,
                                            mapper, regularization):

        blurred_mapping_matrix = convolver.convolve_mapping_matrix(
            mapping_matrix=mapper.mapping_matrix)

        data_vector = inversion_util.data_vector_from_blurred_mapping_matrix_and_data(
            blurred_mapping_matrix=blurred_mapping_matrix,
            image=image,
            noise_map=noise_map,
        )

        curvature_matrix = inversion_util.curvature_matrix_from_blurred_mapping_matrix(
            blurred_mapping_matrix=blurred_mapping_matrix, noise_map=noise_map)

        regularization_matrix = regularization.regularization_matrix_from_mapper(
            mapper=mapper)

        curvature_reg_matrix = np.add(curvature_matrix, regularization_matrix)

        try:
            values = np.linalg.solve(curvature_reg_matrix, data_vector)
        except np.linalg.LinAlgError:
            raise exc.InversionException()

        return InversionImaging(
            image=image,
            noise_map=noise_map,
            mapper=mapper,
            regularization=regularization,
            blurred_mapping_matrix=blurred_mapping_matrix,
            regularization_matrix=regularization_matrix,
            curvature_reg_matrix=curvature_reg_matrix,
            reconstruction=values,
        )
Exemple #3
0
    def from_data_mapper_and_regularization(cls, visibilities, noise_map,
                                            transformer, mapper,
                                            regularization):

        transformed_mapping_matrices = transformer.transformed_mapping_matrices_from_mapping_matrix(
            mapping_matrix=mapper.mapping_matrix)

        real_data_vector = inversion_util.data_vector_from_transformed_mapping_matrix_and_data(
            transformed_mapping_matrix=transformed_mapping_matrices[0],
            visibilities=visibilities[:, 0],
            noise_map=noise_map[:, 0],
        )

        imag_data_vector = inversion_util.data_vector_from_transformed_mapping_matrix_and_data(
            transformed_mapping_matrix=transformed_mapping_matrices[1],
            visibilities=visibilities[:, 1],
            noise_map=noise_map[:, 1],
        )

        real_curvature_matrix = inversion_util.curvature_matrix_from_transformed_mapping_matrix(
            transformed_mapping_matrix=transformed_mapping_matrices[0],
            noise_map=noise_map[:, 0],
        )

        imag_curvature_matrix = inversion_util.curvature_matrix_from_transformed_mapping_matrix(
            transformed_mapping_matrix=transformed_mapping_matrices[1],
            noise_map=noise_map[:, 1],
        )

        regularization_matrix = regularization.regularization_matrix_from_mapper(
            mapper=mapper)

        real_curvature_reg_matrix = np.add(real_curvature_matrix,
                                           regularization_matrix)
        imag_curvature_reg_matrix = np.add(imag_curvature_matrix,
                                           regularization_matrix)

        data_vector = np.add(real_data_vector, imag_data_vector)
        curvature_reg_matrix = np.add(real_curvature_reg_matrix,
                                      imag_curvature_reg_matrix)

        try:
            values = np.linalg.solve(curvature_reg_matrix, data_vector)
        except np.linalg.LinAlgError:
            raise exc.InversionException()

        return InversionInterferometer(
            visibilities=visibilities,
            noise_map=noise_map,
            mapper=mapper,
            regularization=regularization,
            transformed_mapping_matrices=transformed_mapping_matrices,
            regularization_matrix=regularization_matrix,
            curvature_reg_matrix=curvature_reg_matrix,
            reconstruction=values,
        )
    def from_total_pixels_grid_and_weight_map(
        cls,
        total_pixels,
        grid,
        weight_map,
        n_iter=1,
        max_iter=5,
        seed=None,
        stochastic=False,
    ):
        """Calculate a Grid2DSparse from a Grid2D and weight map.

        This is performed by running a KMeans clustering algorithm on the weight map, such that Grid2DSparse (y,x)
        coordinates cluster around the weight map values with higher values.

        This function is used in the `Inversion` package to set up the VoronoiMagnification Pixelization.

        Parameters
        -----------
        total_pixels : int
            The total number of pixels in the Grid2DSparse and input into the KMeans algortihm.
        grid : Grid2D
            The grid of (y,x) coordinates corresponding to the weight map.
        weight_map : np.ndarray
            The 2D array of weight values that the KMeans clustering algorithm adapts to to determine the Grid2DSparse.
        n_iter : int
            The number of times the KMeans algorithm is repeated.
        max_iter : int
            The maximum number of iterations in one run of the KMeans algorithm.
        seed : int or None
            The random number seed, which can be used to reproduce Grid2DSparse's for the same inputs.
        stochastic : bool
            If True, the random number seed is randommly chosen every time the function is called, ensuring every
            pixel-grid is randomly determined and thus stochastic.
        """

        if stochastic:
            seed = np.random.randint(low=1, high=2**31)

        if total_pixels > grid.shape[0]:
            raise exc.GridException

        kmeans = KMeans(n_clusters=total_pixels,
                        random_state=seed,
                        n_init=n_iter,
                        max_iter=max_iter)

        try:
            kmeans = kmeans.fit(X=grid.binned, sample_weight=weight_map)
        except ValueError or OverflowError:
            raise exc.InversionException()

        return Grid2DSparse(
            grid=kmeans.cluster_centers_,
            sparse_index_for_slim_index=kmeans.labels_.astype("int"),
        )
    def interpolated_values_from_shape_native(self, values, shape_native=None):

        if shape_native is not None:

            grid = grid_2d.Grid2D.bounding_box(
                bounding_box=self.mapper.source_pixelization_grid.extent,
                shape_native=shape_native,
                buffer_around_corners=False,
            )

        elif (
            conf.instance["general"]["inversion"]["interpolated_grid_shape"]
            in "image_grid"
        ):

            grid = self.mapper.source_grid_slim

        elif (
            conf.instance["general"]["inversion"]["interpolated_grid_shape"]
            in "source_grid"
        ):

            dimension = int(np.sqrt(self.mapper.pixels))
            shape_native = (dimension, dimension)

            grid = grid_2d.Grid2D.bounding_box(
                bounding_box=self.mapper.source_pixelization_grid.extent,
                shape_native=shape_native,
                buffer_around_corners=False,
            )

        else:

            raise exc.InversionException(
                "In the genenal.ini config file a valid option was not found for the"
                "interpolated_grid_shape. Must be {image_grid, source_grid}"
            )

        interpolated_reconstruction = griddata(
            points=self.mapper.source_pixelization_grid,
            values=values,
            xi=grid.binned.native,
            method="linear",
        )

        interpolated_reconstruction[np.isnan(interpolated_reconstruction)] = 0.0

        return array_2d.Array2D.manual(
            array=interpolated_reconstruction, pixel_scales=grid.pixel_scales
        )
Exemple #6
0
    def interpolated_values_from_shape_2d(self, values, shape_2d=None):

        if shape_2d is not None:

            grid = grids.Grid.bounding_box(
                bounding_box=self.mapper.pixelization_grid.extent,
                shape_2d=shape_2d,
                buffer_around_corners=False,
            )

        elif (conf.instance.general.get("inversion", "interpolated_grid_shape",
                                        str) in "image_grid"):

            grid = self.mapper.grid

        elif (conf.instance.general.get("inversion", "interpolated_grid_shape",
                                        str) in "source_grid"):

            dimension = int(np.sqrt(self.mapper.pixels))
            shape_2d = (dimension, dimension)

            grid = grids.Grid.bounding_box(
                bounding_box=self.mapper.pixelization_grid.extent,
                shape_2d=shape_2d,
                buffer_around_corners=False,
            )

        else:

            raise exc.InversionException(
                "In the genenal.ini config file a valid option was not found for the"
                "interpolated_grid_shape. Must be {image_grid, source_grid}")

        interpolated_reconstruction = griddata(
            points=self.mapper.pixelization_grid,
            values=values,
            xi=grid.in_2d_binned,
            method="linear",
        )

        interpolated_reconstruction[np.isnan(
            interpolated_reconstruction)] = 0.0

        return arrays.Array.manual_2d(array=interpolated_reconstruction,
                                      pixel_scales=grid.pixel_scales)
Exemple #7
0
    def log_determinant_of_matrix_cholesky(matrix):
        """There are two terms in the inversion's Bayesian likelihood function which require the log determinant of \
        a matrix. These are (Nightingale & Dye 2015, Nightingale, Dye and Massey 2018):

        ln[det(F + H)] = ln[det(curvature_reg_matrix)]
        ln[det(H)]     = ln[det(regularization_matrix)]

        The curvature_reg_matrix is positive-definite, which means the above log determinants can be computed \
        efficiently (compared to using np.det) by using a Cholesky decomposition first and summing the log of each \
        diagonal term.

        Parameters
        -----------
        matrix : ndarray
            The positive-definite matrix the log determinant is computed for.
        """
        try:
            return 2.0 * np.sum(np.log(np.diag(np.linalg.cholesky(matrix))))
        except np.linalg.LinAlgError:
            raise exc.InversionException()
    def from_data_mapper_and_regularization(
        cls,
        image: array_2d.Array2D,
        noise_map: array_2d.Array2D,
        convolver: conv.Convolver,
        mapper: typing.Union[mappers.MapperRectangular, mappers.MapperVoronoi],
        regularization: reg.Regularization,
        settings=SettingsInversion(),
        preloads=pload.Preloads(),
    ):

        if preloads.blurred_mapping_matrix is None:

            blurred_mapping_matrix = convolver.convolve_mapping_matrix(
                mapping_matrix=mapper.mapping_matrix
            )

        else:

            blurred_mapping_matrix = preloads.blurred_mapping_matrix

        data_vector = inversion_util.data_vector_via_blurred_mapping_matrix_from(
            blurred_mapping_matrix=blurred_mapping_matrix,
            image=image,
            noise_map=noise_map,
        )

        if preloads.curvature_matrix_sparse_preload is None:

            curvature_matrix = inversion_util.curvature_matrix_via_mapping_matrix_from(
                mapping_matrix=blurred_mapping_matrix, noise_map=noise_map
            )

        else:

            curvature_matrix = inversion_util.curvature_matrix_via_sparse_preload_from(
                mapping_matrix=blurred_mapping_matrix,
                noise_map=noise_map,
                curvature_matrix_sparse_preload=preloads.curvature_matrix_sparse_preload.astype(
                    "int"
                ),
                curvature_matrix_preload_counts=preloads.curvature_matrix_preload_counts.astype(
                    "int"
                ),
            )

        regularization_matrix = regularization.regularization_matrix_from_mapper(
            mapper=mapper
        )

        curvature_reg_matrix = np.add(curvature_matrix, regularization_matrix)

        try:
            values = np.linalg.solve(curvature_reg_matrix, data_vector)
        except np.linalg.LinAlgError:
            raise exc.InversionException()

        if settings.check_solution:
            if np.isclose(a=values[0], b=values[1], atol=1e-4).all():
                if np.isclose(a=values[0], b=values, atol=1e-4).all():
                    raise exc.InversionException()

        return InversionImagingMatrix(
            image=image,
            noise_map=noise_map,
            convolver=convolver,
            mapper=mapper,
            regularization=regularization,
            blurred_mapping_matrix=blurred_mapping_matrix,
            curvature_matrix=curvature_matrix,
            regularization_matrix=regularization_matrix,
            curvature_reg_matrix=curvature_reg_matrix,
            reconstruction=values,
            settings=settings,
        )