Skip to content

Misc

Settings

Disabling auto-save

By default, if your SpatialData object is stored on-disk, it will also store the new elements on disk.

You can disable this behavior as follow:

sopa.settings.auto_save_on_disk = False

Parallelization backends

Some methods (for instance sopa.segmentation.cellpose) may need a parallelization backend to run fast enough.

You can set it as below:

sopa.settings.parallelization_backend = "dask" # using dask
sopa.settings.parallelization_backend = None # no backend (i.e., sequential)

Warning

The dask backend is still experimental. You can add a comment to this issue to help us improve it.

You can also pass some kwargs to the dask Client:

sopa.settings.dask_client_kwargs["n_workers"] = 4

Otherwise, if you don't use the API, you can also set the SOPA_PARALLELIZATION_BACKEND env variable, e.g.:

export SOPA_PARALLELIZATION_BACKEND=dask

Xenium Explorer

sopa.io.explorer.write(path, sdata, image_key=None, shapes_key=None, points_key=None, gene_column=None, pixel_size=0.2125, layer=None, polygon_max_vertices=13, lazy=True, ram_threshold_gb=4, mode=None, save_h5ad=False, run_name=None)

Transform a SpatialData object into inputs for the Xenium Explorer. After running this function, double-click on the experiment.xenium file to open it.

Software download

Make sure you have the latest version of the Xenium Explorer

Note

This function will create up to 7 files, depending on the SpatialData object and the arguments:

  • experiment.xenium contains some experiment metadata. Double-click on this file to open the Xenium Explorer. This file can also be created with write_metadata.

  • morphology.ome.tif is the primary image. This file can also be created with write_image. Add more images with align.

  • analysis.zarr.zip contains the cells categories (or clusters), i.e. adata.obs. This file can also be created with write_cell_categories.

  • cell_feature_matrix.zarr.zip contains the cell-by-gene counts. This file can also be created with write_gene_counts.

  • cells.zarr.zip contains the cells polygon boundaries. This file can also be created with write_polygons.

  • transcripts.zarr.zip contains transcripts locations. This file can also be created with write_transcripts.

  • adata.h5ad is the AnnData object from the SpatialData. This is not used by the Explorer, but only saved for convenience.

Parameters:

Name Type Description Default
path str

Path to the directory where files will be saved.

required
sdata SpatialData

SpatialData object.

required
image_key str | None

Name of the image of interest (key of sdata.images).

None
shapes_key str | None

Name of the cell shapes (key of sdata.shapes).

None
points_key str | None

Name of the transcripts (key of sdata.points).

None
gene_column str | None

Column name of the points dataframe containing the gene names.

None
pixel_size float

Number of microns in a pixel. Invalid value can lead to inconsistent scales in the Explorer.

0.2125
layer str | None

Layer of the AnnData table where the gene counts are saved. If None, uses table.X.

None
polygon_max_vertices int

Maximum number of vertices for the cell polygons.

13
lazy bool

If True, will not load the full images in memory (except if the image memory is below ram_threshold_gb).

True
ram_threshold_gb int | None

Threshold (in gygabytes) from which image can be loaded in memory. If None, the image is never loaded in memory.

4
mode str

string that indicated which files should be created. "-ib" means everything except images and boundaries, while "+tocm" means only transcripts/observations/counts/metadata (each letter corresponds to one explorer file). By default, keeps everything.

None
save_h5ad bool

Whether to save the adata as h5ad in the explorer directory (for convenience only, since h5ad is faster to open than the original .zarr table)

False
run_name str | None

Name of the run displayed in the Xenium Explorer. If None, uses the image_key.

None
Source code in sopa/io/explorer/converter.py
def write(
    path: str,
    sdata: SpatialData,
    image_key: str | None = None,
    shapes_key: str | None = None,
    points_key: str | None = None,
    gene_column: str | None = None,
    pixel_size: float = 0.2125,
    layer: str | None = None,
    polygon_max_vertices: int = 13,
    lazy: bool = True,
    ram_threshold_gb: int | None = 4,
    mode: str = None,
    save_h5ad: bool = False,
    run_name: str | None = None,
) -> None:
    """
    Transform a SpatialData object into inputs for the Xenium Explorer.
    After running this function, double-click on the `experiment.xenium` file to open it.

    !!! note "Software download"
        Make sure you have the latest version of the [Xenium Explorer](https://www.10xgenomics.com/support/software/xenium-explorer)

    Note:
        This function will create up to 7 files, depending on the `SpatialData` object and the arguments:

        - `experiment.xenium` contains some experiment metadata. Double-click on this file to open the Xenium Explorer. This file can also be created with [`write_metadata`](./#sopa.io.explorer.write_metadata).

        - `morphology.ome.tif` is the primary image. This file can also be created with [`write_image`](./#sopa.io.explorer.write_image). Add more images with `align`.

        - `analysis.zarr.zip` contains the cells categories (or clusters), i.e. `adata.obs`. This file can also be created with [`write_cell_categories`](./#sopa.io.explorer.write_cell_categories).

        - `cell_feature_matrix.zarr.zip` contains the cell-by-gene counts. This file can also be created with [`write_gene_counts`](./#sopa.io.explorer.write_gene_counts).

        - `cells.zarr.zip` contains the cells polygon boundaries. This file can also be created with [`write_polygons`](./#sopa.io.explorer.write_polygons).

        - `transcripts.zarr.zip` contains transcripts locations. This file can also be created with [`write_transcripts`](./#sopa.io.explorer.write_transcripts).

        - `adata.h5ad` is the `AnnData` object from the `SpatialData`. This is **not** used by the Explorer, but only saved for convenience.

    Args:
        path: Path to the directory where files will be saved.
        sdata: SpatialData object.
        image_key: Name of the image of interest (key of `sdata.images`).
        shapes_key: Name of the cell shapes (key of `sdata.shapes`).
        points_key: Name of the transcripts (key of `sdata.points`).
        gene_column: Column name of the points dataframe containing the gene names.
        pixel_size: Number of microns in a pixel. Invalid value can lead to inconsistent scales in the Explorer.
        layer: Layer of the AnnData table where the gene counts are saved. If `None`, uses `table.X`.
        polygon_max_vertices: Maximum number of vertices for the cell polygons.
        lazy: If `True`, will not load the full images in memory (except if the image memory is below `ram_threshold_gb`).
        ram_threshold_gb: Threshold (in gygabytes) from which image can be loaded in memory. If `None`, the image is never loaded in memory.
        mode: string that indicated which files should be created. "-ib" means everything except images and boundaries, while "+tocm" means only transcripts/observations/counts/metadata (each letter corresponds to one explorer file). By default, keeps everything.
        save_h5ad: Whether to save the adata as h5ad in the explorer directory (for convenience only, since h5ad is faster to open than the original .zarr table)
        run_name: Name of the run displayed in the Xenium Explorer. If `None`, uses the `image_key`.
    """
    path: Path = Path(path)
    _check_explorer_directory(path)

    image_key, _ = get_spatial_image(sdata, key=image_key, return_key=True)

    ### Saving cell categories and gene counts
    if SopaKeys.TABLE in sdata.tables:
        adata = sdata.tables[SopaKeys.TABLE]

        shapes_key = adata.uns[ATTRS_KEY]["region"]
        shapes_key = shapes_key[0] if isinstance(shapes_key, list) else shapes_key

        geo_df = sdata[shapes_key]

        if _should_save(mode, "c"):
            write_gene_counts(path, adata, layer=layer)
        if _should_save(mode, "o"):
            write_cell_categories(path, adata)

    ### Saving cell boundaries
    if shapes_key is None:
        shapes_key, geo_df = get_boundaries(sdata, return_key=True, warn=True)
    else:
        geo_df = sdata[shapes_key]

    if _should_save(mode, "b") and geo_df is not None:
        geo_df = to_intrinsic(sdata, geo_df, image_key)

        if SopaKeys.TABLE in sdata.tables:
            geo_df = geo_df.loc[adata.obs[adata.uns[ATTRS_KEY]["instance_key"]]]

        assert all(isinstance(geom, Polygon) for geom in geo_df.geometry), "All geometries must be a `shapely.Polygon`"

        write_polygons(path, geo_df.geometry, polygon_max_vertices, pixel_size=pixel_size)

    ### Saving transcripts
    df = None
    if len(sdata.points):
        df = get_spatial_element(sdata.points, key=points_key or sdata.attrs.get(SopaAttrs.TRANSCRIPTS))

    if _should_save(mode, "t") and df is not None:
        gene_column = gene_column or get_feature_key(df)
        if gene_column is not None:
            df = to_intrinsic(sdata, df, image_key)
            write_transcripts(path, df, gene_column, pixel_size=pixel_size)
        else:
            log.warning("The argument 'gene_column' has to be provided to save the transcripts")

    ### Saving image
    if _should_save(mode, "i"):
        write_image(
            path,
            sdata[image_key],
            lazy=lazy,
            ram_threshold_gb=ram_threshold_gb,
            pixel_size=pixel_size,
        )

    ### Saving experiment.xenium file
    if _should_save(mode, "m"):
        write_metadata(path, run_name or image_key, shapes_key, _get_n_obs(sdata, geo_df), pixel_size)

    if save_h5ad and SopaKeys.TABLE in sdata.tables:
        sdata.tables[SopaKeys.TABLE].write_h5ad(path / FileNames.H5AD)

    log.info(f"Saved files in the following directory: {path}")
    log.info(f"You can open the experiment with 'open {path / FileNames.METADATA}'")

sopa.io.explorer.align(sdata, image, transformation_matrix_path, key_added=None, image_key=None, overwrite=False)

Add an image to the SpatialData object after alignment with the Xenium Explorer.

Parameters:

Name Type Description Default
sdata SpatialData

A SpatialData object

required
image DataArray

A DataArray object. Note that image.name is used as the key for the aligned image.

required
transformation_matrix_path str

Path to the .csv transformation matrix exported from the Xenium Explorer

required
key_added str | None

Optional name to add to the new image. If None, will use image.name.

None
image_key str

Optional name of the image on which it has been aligned. Required if multiple images in the SpatialData object.

None
overwrite bool

Whether to overwrite the image, if already existing.

False
Source code in sopa/io/explorer/images.py
def align(
    sdata: SpatialData,
    image: DataArray,
    transformation_matrix_path: str,
    key_added: str | None = None,
    image_key: str = None,
    overwrite: bool = False,
):
    """Add an image to the `SpatialData` object after alignment with the Xenium Explorer.

    Args:
        sdata: A `SpatialData` object
        image: A `DataArray` object. Note that `image.name` is used as the key for the aligned image.
        transformation_matrix_path: Path to the `.csv` transformation matrix exported from the Xenium Explorer
        key_added: Optional name to add to the new image. If `None`, will use `image.name`.
        image_key: Optional name of the image on which it has been aligned. Required if multiple images in the `SpatialData` object.
        overwrite: Whether to overwrite the image, if already existing.
    """
    key_added = key_added or image.name

    assert key_added is not None, "The image has no name, use the `key_added` argument to provide one"
    assert key_added not in sdata, f"Image '{key_added}' already exists in the `SpatialData` object"

    to_pixel = Sequence(
        [
            Affine(
                np.genfromtxt(transformation_matrix_path, delimiter=","),
                input_axes=("x", "y"),
                output_axes=("x", "y"),
            )
        ]
    )

    default_image = get_spatial_image(sdata, image_key)

    original_transformations = get_transformation(default_image, get_all=True)
    transformations = {cs: to_pixel.compose_with(t) for cs, t in original_transformations.items()}

    set_transformation(image, transformations, set_all=True)

    add_spatial_element(sdata, key_added, image, overwrite=overwrite)
    log.info(f"Added image '{key_added}' (aligned with the Xenium Explorer)")

sopa.io.explorer.add_explorer_selection(sdata, path, key_added='explorer_selection', shapes_key=None, image_key=None, pixel_size=0.2125)

After saving a selection on the Xenium Explorer, it will add all polygons inside sdata.shapes[shapes_key]

Parameters:

Name Type Description Default
sdata SpatialData

A SpatialData object

required
path str

The path to the coordinates.csv selection file

required
key_added str

The name to provide to the selection as shapes

'explorer_selection'
shapes_key str | None

Deprecated. Use key_added instead.

None
image_key str | None

The original image name

None
pixel_size float

Number of microns in a pixel. It must be the same value as the one used in sopa.io.write

0.2125
Source code in sopa/io/explorer/utils.py
def add_explorer_selection(
    sdata: SpatialData,
    path: str,
    key_added: str = "explorer_selection",
    shapes_key: str | None = None,
    image_key: str | None = None,
    pixel_size: float = 0.2125,
):
    """After saving a selection on the Xenium Explorer, it will add all polygons inside `sdata.shapes[shapes_key]`

    Args:
        sdata: A `SpatialData` object
        path: The path to the `coordinates.csv` selection file
        key_added: The name to provide to the selection as shapes
        shapes_key: Deprecated. Use `key_added` instead.
        image_key: The original image name
        pixel_size: Number of microns in a pixel. It must be the same value as the one used in `sopa.io.write`
    """
    if shapes_key is not None:
        log.warning(
            "The `shapes_key` argument is deprecated and will be removed in sopa==2.1.0. Use `key_added` instead."
        )
        key_added = shapes_key

    polys = xenium_explorer_selection(path, pixel_size=pixel_size, return_list=True)
    image = get_spatial_element(sdata.images, key=image_key or sdata.attrs.get(SopaAttrs.CELL_SEGMENTATION))

    transformations = get_transformation(image, get_all=True).copy()

    geo_df = ShapesModel.parse(gpd.GeoDataFrame(geometry=polys), transformations=transformations)
    add_spatial_element(sdata, key_added, geo_df)

sopa.io.explorer.int_cell_id(explorer_cell_id)

Transforms an alphabetical cell id from the Xenium Explorer to an integer ID

E.g., int_cell_id('aaaachba-1') = 10000

Source code in sopa/io/explorer/utils.py
def int_cell_id(explorer_cell_id: str) -> int:
    """Transforms an alphabetical cell id from the Xenium Explorer to an integer ID

    E.g., int_cell_id('aaaachba-1') = 10000"""
    code = explorer_cell_id[:-2] if explorer_cell_id[-2] == "-" else explorer_cell_id
    coefs = [ord(c) - 97 for c in code][::-1]
    return sum(value * 16**i for i, value in enumerate(coefs))

sopa.io.explorer.str_cell_id(cell_id)

Transforms an integer cell ID into an Xenium Explorer alphabetical cell id

E.g., str_cell_id(10000) = 'aaaachba-1'

Source code in sopa/io/explorer/utils.py
def str_cell_id(cell_id: int) -> str:
    """Transforms an integer cell ID into an Xenium Explorer alphabetical cell id

    E.g., str_cell_id(10000) = 'aaaachba-1'"""
    coefs = []
    for _ in range(8):
        cell_id, coef = divmod(cell_id, 16)
        coefs.append(coef)
    return "".join([chr(97 + coef) for coef in coefs][::-1]) + "-1"

sopa.io.explorer.write_image(path, image, lazy=True, tile_width=TILE_SIZE, n_subscales=5, pixel_size=0.2125, ram_threshold_gb=4, is_dir=True)

Convert an image into a morphology.ome.tif file that can be read by the Xenium Explorer

Parameters:

Name Type Description Default
path str

Path to the Xenium Explorer directory where the image will be written

required
image DataTree | DataArray | ndarray

Image of shape (C, Y, X)

required
lazy bool

If False, the image will not be read in-memory (except if the image size is below ram_threshold_gb). If True, all the images levels are always loaded in-memory.

True
tile_width int

Xenium tile width (do not update).

TILE_SIZE
n_subscales int

Number of sub-scales in the pyramidal image.

5
pixel_size float

Xenium pixel size (do not update).

0.2125
ram_threshold_gb int | None

If an image (of any level of the pyramid) is below this threshold, it will be loaded in-memory.

4
is_dir bool

If False, then path is a path to a single file, not to the Xenium Explorer directory.

True
Source code in sopa/io/explorer/images.py
def write_image(
    path: str,
    image: DataTree | DataArray | np.ndarray,
    lazy: bool = True,
    tile_width: int = TILE_SIZE,
    n_subscales: int = 5,
    pixel_size: float = 0.2125,
    ram_threshold_gb: int | None = 4,
    is_dir: bool = True,
):
    """Convert an image into a `morphology.ome.tif` file that can be read by the Xenium Explorer

    Args:
        path: Path to the Xenium Explorer directory where the image will be written
        image: Image of shape `(C, Y, X)`
        lazy: If `False`, the image will not be read in-memory (except if the image size is below `ram_threshold_gb`). If `True`, all the images levels are always loaded in-memory.
        tile_width: Xenium tile width (do not update).
        n_subscales: Number of sub-scales in the pyramidal image.
        pixel_size: Xenium pixel size (do not update).
        ram_threshold_gb: If an image (of any level of the pyramid) is below this threshold, it will be loaded in-memory.
        is_dir: If `False`, then `path` is a path to a single file, not to the Xenium Explorer directory.
    """
    path = explorer_file_path(path, FileNames.IMAGE, is_dir)

    if isinstance(image, np.ndarray):
        assert len(image.shape) == 3, "Can only write channels with shape (C,Y,X)"
        log.info(f"Converting image of shape {image.shape} into a DataArray (with dims: C,Y,X)")
        image = DataArray(image, dims=["c", "y", "x"], name="image")

    image = _to_xenium_explorer_multiscale(image, n_subscales)

    image_writer = MultiscaleImageWriter(image, pixel_size=pixel_size, tile_width=tile_width)
    image_writer.write(path, lazy=lazy, ram_threshold_gb=ram_threshold_gb)

sopa.io.explorer.save_column_csv(path, adata, key)

Save one column of the AnnData object as a CSV that can be open interactively in the explorer, under the "cell" panel.

Parameters:

Name Type Description Default
path str

Path where to write the CSV that will be open in the Xenium Explorer

required
adata AnnData

An AnnData object

required
key str

Key of adata.obs containing the column to convert

required
Source code in sopa/io/explorer/table.py
def save_column_csv(path: str, adata: AnnData, key: str):
    """Save one column of the AnnData object as a CSV that can be open interactively in the explorer, under the "cell" panel.

    Args:
        path: Path where to write the CSV that will be open in the Xenium Explorer
        adata: An `AnnData` object
        key: Key of `adata.obs` containing the column to convert
    """
    df = pd.DataFrame({"cell_id": adata.obs_names, "group": adata.obs[key].values})
    df.to_csv(path, index=None)

Report

sopa.io.write_report(path, sdata)

Create a HTML report (or web report) after running Sopa.

Note

This report is automatically generated based on a custom python-to-html engine

Parameters:

Name Type Description Default
path str

Path to the .html report that has to be created

required
sdata SpatialData

A SpatialData object, after running Sopa

required
Source code in sopa/io/report/generate.py
def write_report(path: str, sdata: SpatialData):
    """Create a HTML report (or web report) after running Sopa.

    Note:
        This report is automatically generated based on a custom python-to-html engine

    Args:
        path: Path to the `.html` report that has to be created
        sdata: A `SpatialData` object, after running Sopa
    """
    with warnings.catch_warnings():
        warnings.simplefilter(action="ignore", category=FutureWarning)

        sections = SectionBuilder(sdata).compute_sections()

        log.info(f"Writing report to {path}")
        Root(sections).write(path)