@@ -364,6 +364,151 @@ readVisiumHD <- function(data_dir, bin_size = c(2L, 8L, 16L),
364364 sfes
365365}
366366
367+ # ' Flip the Y-axis of cell or nucleus segmentations to align with H&E image.
368+ # '
369+ # ' @param sf a sf object read from a `.geojson` file.
370+ # ' @param type "POINT" for cell centroid, or "POLYGON" for cell segmentation
371+ # ' mask. Default is "POINT".
372+ # ' @param img_height total length along the Y axis of the image. Obtained by
373+ # ' reading in `hires` or `lowre`s `.png` under `/spatial` folder with
374+ # ' `magick::image_read()`.
375+ # ' @param scalef scaling factor from `/spatial/scalefactors_json.json` file
376+ # '
377+ # ' @returns a sf object with Y-axis of the points or polygons flipped.
378+ # '
379+ # ' @examples
380+ # ' geo_data_flipped <- flip_sf_Y(sf = geo_data, type = "POLYGON",
381+ # ' img_height = 3886, scalef = 0.079)
382+ flip_sf_Y <- function (sf , type = " POINT" , img_height , scalef ){
383+ st_geometry(sf ) <- st_sfc(
384+ lapply(st_geometry(sf ), function (geom ) {
385+ coords <- st_coordinates(geom )
386+ coords [, 2 ] <- img_height / scalef - coords [, 2 ]
387+ if (type == " POINT" ){
388+ st_point(coords )
389+ }else { # type == "POLYGON"
390+ st_polygon(list (matrix (coords [, 1 : 2 ], ncol = 2 )))
391+ }
392+ }),
393+ crs = st_crs(sf )
394+ )
395+
396+ return (sf )
397+ }
398+
399+
400+ # ' Reader for Visium HD cell segmented output.
401+ # '
402+ # ' @param td path to unzipped Visium HD data download with `/segmented_outputs`
403+ # ' @param res "hires" or "lowres" for what .png to read in.
404+ # '
405+ # ' @returns a `SpatialFeatureExperiment` object.
406+ # '
407+ # ' @examples
408+ # ' vhdsfe <- readVisiumHDCellSeg(td = "~/Desktop", res = "hires")
409+ readVisiumHDCellSeg <- function (td , res ){
410+ if (res == " hires" ){
411+ png_name = " tissue_hires_image.png"
412+ scalef_name = " tissue_hires_scalef"
413+ }else { # res = "lowres"
414+ png_name = " tissue_lowres_image.png"
415+ scalef_name = " tissue_lowres_scalef"
416+ }
417+
418+ # # Coord ---------------------------------------------------------------
419+ # Read the GeoJSON file of cell segmentation
420+ geo_data <- st_read(file.path(td , " segmented_outputs" ,
421+ " cell_segmentations.geojson" ))
422+ scalef <- rjson :: fromJSON(file = file.path(td , " segmented_outputs" , " spatial" ,
423+ " scalefactors_json.json" ))
424+ image <- image_read(file.path(td , " segmented_outputs" ,
425+ " spatial" , png_name ))
426+ info <- image_info(image )
427+
428+ # Get centroid
429+ st_crs(geo_data ) <- NA
430+ rownames(geo_data ) <- geo_data $ cell_id
431+ centroids <- st_centroid(geo_data )
432+ centroids $ cell_id <- as.character(centroids $ cell_id )
433+
434+
435+ # # Countmat ------------------------------------------------------------
436+ countmat_file <- file.path(td , " segmented_outputs" ,
437+ " filtered_feature_cell_matrix.h5" )
438+ vhdcellsce <- DropletUtils :: read10xCounts(countmat_file , col.names = TRUE )
439+
440+ colnames(vhdcellsce ) <- sub(" cellid_0*([0-9]+)-.*" , " \\ 1" , vhdcellsce $ Barcode )
441+
442+ rownames(vhdcellsce ) <- rowData(vhdcellsce )$ Symbol
443+ rownames(vhdcellsce ) <- make.unique(rownames(vhdcellsce ))
444+
445+
446+ # # Matching coordinates with count matrix ------------------------------
447+ # Subset to common cells
448+ common_cells <- intersect(centroids $ cell_id , colnames(vhdcellsce ))
449+
450+ geo_data <- geo_data [geo_data $ cell_id %in% common_cells , ]
451+ centroids <- centroids [centroids $ cell_id %in% common_cells , ]
452+ vhdcellsce <- vhdcellsce [, colnames(vhdcellsce ) %in% common_cells ]
453+
454+ # # Sanity check on ordering of cell ids
455+ # all(as.character(centroids$cell_id) == colnames(vhdcellsce))
456+
457+ # Extract coordinates
458+ coords <- st_coordinates(centroids )
459+ colnames(coords ) <- c(" pxl_col_in_fullres" , " pxl_row_in_fullres" )
460+
461+
462+ # # Construct SPE -------------------------------------------------------
463+ vhd <- SpatialExperiment(
464+ assays = list (counts = as(counts(vhdcellsce ), " dgCMatrix" )),
465+ rowData = rowData(vhdcellsce ),
466+ colData = cbind(colData(vhdcellsce ), coords ),
467+ spatialCoordsNames = colnames(coords ),
468+ scaleFactors = scalef $ tissue_hires_scalef ,
469+ imageSources = file.path(td , " segmented_outputs" , " spatial" , png_name ),
470+ loadImage = TRUE
471+ )
472+ imgData(vhd )$ image_id <- res
473+
474+ # vhd
475+
476+ # # Coerce to SFE -------------------------------------------------------
477+ vhdsfe <- toSpatialFeatureExperiment(vhd )
478+ # vhdsfe
479+
480+ # # Add polygons to colGeometries
481+ colGeometries(vhdsfe )$ cellSeg <- geo_data
482+
483+
484+ # Flip Y in colGeometries -------------------------------------------------
485+ # Flip Centroid Y
486+ centroids <- colGeometries(vhdsfe )$ centroids
487+ colGeometries(vhdsfe )$ updatecentroids <-
488+ flip_sf_Y(centroids , type = " POINT" , img_height = info $ height ,
489+ scalef = imgData(vhdsfe )$ scaleFactor )
490+
491+ # Flip Cellseg Y
492+ cellSeg <- colGeometries(vhdsfe )$ cellSeg
493+ colGeometries(vhdsfe )$ updatecellSeg <-
494+ flip_sf_Y(cellSeg , type = " POLYGON" , img_height = info $ height ,
495+ scalef = imgData(vhdsfe )$ scaleFactor )
496+
497+ # Flip Centroid Y in `spatialCoords()`
498+ spatialCoords(vhdsfe )[, " pxl_row_in_fullres" ] <-
499+ info $ height / imgData(vhdsfe )$ scaleFactor -
500+ spatialCoords(vhdsfe )[, " pxl_row_in_fullres" ]
501+
502+ # Add high res coords to `colData()`
503+ vhdsfe [[paste0(" pxl_col_in_" , res )]] <-
504+ spatialCoords(vhdsfe )[, " pxl_col_in_fullres" ] * imgData(vhdsfe )$ scaleFactor
505+ vhdsfe [[paste0(" pxl_row_in_" , res )]] <-
506+ spatialCoords(vhdsfe )[, " pxl_row_in_fullres" ] * imgData(vhdsfe )$ scaleFactor
507+
508+ return (vhdsfe )
509+ }
510+
511+
367512# ' @importFrom sf st_nearest_feature st_distance
368513# ' @importFrom stats median
369514.pixel2micron <- function (df ) {
0 commit comments