Package 'spatialSimGP'

Title: Simulate Spatial Transcriptomics Data with the Mean-variance Relationship
Description: This packages simulates spatial transcriptomics data with the mean- variance relationship using a Gaussian Process model per gene.
Authors: Kinnary Shah [aut, cre] , Boyi Guo [aut] , Stephanie C. Hicks [aut]
Maintainer: Kinnary Shah <[email protected]>
License: MIT + file LICENSE
Version: 1.1.0
Built: 2024-10-31 05:28:04 UTC
Source: https://github.com/bioc/spatialSimGP

Help Index


Simulate SpatialExperiment object

Description

Simulate a SpatialExperiment object with spatially varying genes

Usage

spatial_simulate(
  n_genes,
  proportion,
  coords,
  range_sigma.sq,
  range_beta,
  length_scale,
  length_scale_option = "fixed"
)

Arguments

n_genes

an integer specifying the number of genes to simulate.

proportion

a numeric value specifying the proportion of genes that will have no spatially varying patterns.

coords

a matrix of coordinates.

range_sigma.sq

a numeric vector of length 2 specifying the range of the spatial variance parameter.

range_beta

a numeric vector of length 2 specifying the range of the mean expression value.

length_scale

if length_scale_option is "fixed", a numeric value specifying the length scale parameter. If length_scale_option is "unique", a numeric vector of length n_genes specifying the length scale parameter for each gene.

length_scale_option

a character string specifying the length scale option. Options are "fixed" for a single length scale for all genes or "unique" for a unique length scale for each gene.

Details

This function simulates a SpatialExperiment object with spatially varying genes. The function takes in the number of genes to simulate, the proportion of genes that will have no spatially varying patterns, a matrix of coordinates, the range of the spatial variance parameter, the range of the mean expression value, the length scale parameter, and the length scale option.

Value

A SpatialExperiment object with the simulated data.

Examples

library(STexampleData)

set.seed(1)
n_genes <- 1
proportion <- 0.5
range_sigma.sq <- c(0.2, 3)
range_beta <- c(0.5, 9)
length_scale <- 60

spe_demo <- Visium_mouseCoronal()
colData(spe_demo)$subset <- ifelse(
  colData(spe_demo)$array_row > 20 &
  colData(spe_demo)$array_row < 65 &
  colData(spe_demo)$array_col > 30 &
  colData(spe_demo)$array_col < 65,
  TRUE, FALSE
)
spe_demo <- spe_demo[, colData(spe_demo)$subset]
coords <- spatialCoords(spe_demo)

spe <- spatial_simulate(n_genes, proportion, coords, range_sigma.sq, range_beta,
   length_scale, length_scale_option = "fixed")