Title: | Simulator for Spatially Resolved Transcriptomics |
---|---|
Description: | An independent, reproducible, and flexible Spatially Resolved Transcriptomics (SRT) simulation framework that can be used to facilitate the development of SRT analytical methods for a wide variety of SRT-specific analyses. It utilizes spatial localization information to simulate SRT expression count data in a reproducible and scalable fashion. Two major simulation schemes are implemented in 'SRTsim': reference-based and reference-free. |
Authors: | Jiaqiang Zhu [aut, ctb, cre] , Lulu Shang [aut] , Xiang Zhou [aut] |
Maintainer: | Jiaqiang Zhu <[email protected]> |
License: | GPL (>= 3) |
Version: | 0.99.7 |
Built: | 2024-11-20 06:38:50 UTC |
Source: | CRAN |
Summarize metrics for reference data and synthetic data
compareSRT(simsrt)
compareSRT(simsrt)
simsrt |
A SRTsim object |
Returns an object with summarized metrics
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) ## Compute metrics toySRT <- compareSRT(toySRT)
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) ## Compute metrics toySRT <- compareSRT(toySRT)
Convert continuous coordinate into integer, essential for BayesSpace to determine the neighborhood info
convert_grid(x)
convert_grid(x)
x |
A numeric vector of continuous coordinate |
Returns a numeric vector oof integer coordinate
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Create New Locations within Profile toySRT2 <- srtsim_newlocs(toySRT,new_loc_num=1000) ## Convert non-integer x-coordinates into an integer value newGrid_x <- convert_grid(simcolData(toySRT2)$x)
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Create New Locations within Profile toySRT2 <- srtsim_newlocs(toySRT,new_loc_num=1000) ## Convert non-integer x-coordinates into an integer value newGrid_x <- convert_grid(simcolData(toySRT2)$x)
Create simSRT object
createSRT(count_in, loc_in, refID = "ref1")
createSRT(count_in, loc_in, refID = "ref1")
count_in |
A gene expression count |
loc_in |
A location |
refID |
A |
Returns a spatialExperiment-based object
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) ## Explore the object toySRT
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) ## Explore the object toySRT
Access Model Fitting Parameters
EstParam(x)
EstParam(x)
x |
SRTsim object |
Returns a list of estimated parameters by fitting models
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) EstParam(toySRT)
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) EstParam(toySRT)
A data list containing the a gene expression matrix and a location dataframe
exampleLIBD
exampleLIBD
A data list
A sparse matrix with 80 rows and 3611 columns.
A data frame with 3611 rows and 6 columns.
created based on a SpatialLIBD data (SampleID: 151673) to serve as an example
data(exampleLIBD) #Lazy loading. Data becomes visible as soon as called
data(exampleLIBD) #Lazy loading. Data becomes visible as soon as called
fitting data with poisson through optim function
fit_pos_optim(x, maxiter = 500)
fit_pos_optim(x, maxiter = 500)
x |
A vector of count values to be fitted |
maxiter |
number of iteration |
Returns a vector with mean paramter lambda, loglikelihood value llk, convergence
Extracted summarized metrics for reference data and synthetic data
get_metrics_pd(simsrt, metric = "GeneMean")
get_metrics_pd(simsrt, metric = "GeneMean")
simsrt |
A SRTsim object |
metric |
Specification of metrics to be plotted. |
Returns a dataframe for ggplot
Summarize gene-wise summary metrics
get_stats_gene(mat, group, log_trans = TRUE)
get_stats_gene(mat, group, log_trans = TRUE)
mat |
A count matrix |
group |
A group label |
log_trans |
A logical constant indicating whether to log transform the gene mean and variance |
Returns a n by 5 dataframe with location metrics
Summarize location-wise summary metrics
get_stats_loc(mat, group, log_trans = TRUE)
get_stats_loc(mat, group, log_trans = TRUE)
mat |
A count matrix |
group |
A group label |
log_trans |
A logical constant indicating whether to log transform the libsize |
Returns a n by 3 dataframe with location metrics
Access User-Specified Parameters
metaParam(x)
metaParam(x)
x |
SRTsim object |
Returns a list of user-specified parameters
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) metaParam(toySRT)
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) metaParam(toySRT)
Access reference colData
refcolData(x)
refcolData(x)
x |
SRTsim object |
Returns the colData of reference data
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) refcolData(toySRT)
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) refcolData(toySRT)
Access reference count matrix
refCounts(x)
refCounts(x)
x |
SRTsim object |
Returns a reference count matrix
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) refCounts(toySRT)[1:3,1:3]
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) refCounts(toySRT)[1:3,1:3]
Access reference rowData
refrowData(x)
refrowData(x)
x |
SRTsim object |
Returns the rowData of reference data
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) refrowData(toySRT)
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) refrowData(toySRT)
ReSimulate Count Data with Parameters Specification from Shiny
reGenCountshiny(shinyOutput, NewSeed = NULL)
reGenCountshiny(shinyOutput, NewSeed = NULL)
shinyOutput |
A list of Shiny Output. Including a simCount, simInfo,simcountParam,simLocParam |
NewSeed |
A new seed for data generation. Useful when multiple replicates are needed. |
Returns a Count DataFrame
## Re-generate Count Data based on ShinyOutput Parameters, should be same as simCount in ShinyOutput cMat <- reGenCountshiny(toyShiny) ## Generate Count Data with A New Seed based on ShinyOutput Parameters cMat2 <- reGenCountshiny(toyShiny,NewSeed=2) ## Comparison across the output toyShiny$simCount[1:3,1:3] cMat[1:3,1:3] cMat2[1:3,1:3]
## Re-generate Count Data based on ShinyOutput Parameters, should be same as simCount in ShinyOutput cMat <- reGenCountshiny(toyShiny) ## Generate Count Data with A New Seed based on ShinyOutput Parameters cMat2 <- reGenCountshiny(toyShiny,NewSeed=2) ## Comparison across the output toyShiny$simCount[1:3,1:3] cMat[1:3,1:3] cMat2[1:3,1:3]
Create a SRTsim object from reference-free shinyoutput
Shiny2SRT(shinyOutput)
Shiny2SRT(shinyOutput)
shinyOutput |
A list of Shiny Output. Including a simCount, simInfo,simcountParam,simLocParam |
Returns a SRTsim object with user-specified parameters stored in metaParam slot.
shinySRT <- Shiny2SRT(toyShiny) ## Explore the new SRT object shinySRT@metaParam shinySRT@simCounts[1:3,1:3] shinySRT@simcolData
shinySRT <- Shiny2SRT(toyShiny) ## Explore the new SRT object shinySRT@metaParam shinySRT@simCounts[1:3,1:3] shinySRT@simcolData
Access synthetic colData
simcolData(x)
simcolData(x)
x |
SRTsim object |
Returns the colData of synthetic data
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) simcolData(toySRT)
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) simcolData(toySRT)
Access synthetic count matrix
simCounts(x)
simCounts(x)
x |
SRTsim object |
Returns a synthetic count matrix
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) simCounts(toySRT)[1:3,1:3]
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) simCounts(toySRT)[1:3,1:3]
Fit the marginal distributions for single gene
simNewLocs(newN, lay_out = c("grid", "random"), preLoc)
simNewLocs(newN, lay_out = c("grid", "random"), preLoc)
newN |
A integer specifying the number of spatial locations in the synthetic data |
lay_out |
A character string specifying arrangement of new generated spatial locations. Default is "grid" |
preLoc |
A data frame of shape n by 3 that x, y coodinates and domain label |
Returns a n by 2 dataframe with newly generated spatial locations
Access synthetic rowData
simrowData(x)
simrowData(x)
x |
SRTsim object |
Returns the rowData of synthetic data
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) simrowData(toySRT)
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) simrowData(toySRT)
Generate Data with Cell-Cell Interaction Under Reference-Free Mode
srtsim_cci_free( zero_prop_in = 0, disper_in = Inf, mu_in = 1, numGene = 1000, location_in, region_cell_map, fc = 3, LR_in, sim_seed = 1, numKNN = 4, numSingleCellType = 2000 )
srtsim_cci_free( zero_prop_in = 0, disper_in = Inf, mu_in = 1, numGene = 1000, location_in, region_cell_map, fc = 3, LR_in, sim_seed = 1, numKNN = 4, numSingleCellType = 2000 )
zero_prop_in |
A number specifying zero proportion for the count model, default is 0 |
disper_in |
A number specifying dispersion for the count model, default is Inf. Same as the size parameter in rnbinom. |
mu_in |
A number specifying mean for background count model, default is 1 |
numGene |
An integer specifying the number of genes in the synthetic data, default is 1000 |
location_in |
A dataframe with x, y, and region_label |
region_cell_map |
A dataframe specifying the cell type proportion in each region. Row: region,Column: cell type. |
fc |
A number specifying effect size for ligand-receptor pairs that mediate the cel-cell communication, default is 3 |
LR_in |
A dataframe specifying ligand and receptor pairs, containing four columns: protein_a, protein_b, celltypeA, and celltype B |
sim_seed |
A number for reproducible purpose |
numKNN |
A number specifying number of nearest neighbors with elevated gene expressin levels, default is 4 |
numSingleCellType |
A number specifying number of spots in the background pool. Gene expression count are then sampled from this background pool. |
Returns a SRTsim object with a newly generated count matrix and correspoding parameters
Generate Data with Cell-Cell Interaction Under Reference-Based Mode
srtsim_cci_ref( EstParam = NULL, numGene = 1000, location_in, region_cell_map, fc = 3, LR_in, sim_seed = 1, numKNN = 4, numSingleCellType = 2000 )
srtsim_cci_ref( EstParam = NULL, numGene = 1000, location_in, region_cell_map, fc = 3, LR_in, sim_seed = 1, numKNN = 4, numSingleCellType = 2000 )
EstParam |
A list of estimated parameters from srtsim_fit function, EstParam slot if the simSRT object. |
numGene |
An integer specifying the number of genes in the synthetic data, default is 1000 |
location_in |
A dataframe with x, y, and region_label |
region_cell_map |
A dataframe specifying the cell type proportion in each region. Row: region,Column: cell type. |
fc |
A number specifying effect size for ligand-receptor pairs that mediate the cel-cell communication, default is 3 |
LR_in |
A dataframe specifying ligand and receptor pairs, containing four columns: protein_a, protein_b, celltypeA, and celltype B |
sim_seed |
A number for reproducible purpose |
numKNN |
A number specifying number of nearest neighbors with elevated gene expressin levels, default is 4 |
numSingleCellType |
A number specifying number of spots in the background pool. Gene expression count are then sampled from this background pool. |
Returns a SRTsim object with a newly generated count matrix and correspoding parameters
Generate Data with Estimated Parameters
srtsim_count( simsrt, breaktie = "random", total_count_new = NULL, total_count_old = NULL, rrr = NULL, nn_num = 5, nn_func = c("mean", "median", "ransam"), numCores = 1, verbose = FALSE )
srtsim_count( simsrt, breaktie = "random", total_count_new = NULL, total_count_old = NULL, rrr = NULL, nn_num = 5, nn_func = c("mean", "median", "ransam"), numCores = 1, verbose = FALSE )
simsrt |
A object with estimated parameters from fitting step |
breaktie |
A character string specifying how ties are treated. Same as the "tie.method" in rank function |
total_count_new |
The (expected) total number of reads or UMIs in the simulated count matrix. |
total_count_old |
The total number of reads or UMIs in the original count matrix. |
rrr |
The ratio applies to the gene-specific mean estimate, used for the fixing average sequencing depth simulation. Default is null. Its specification will override the specification of total_count_new and total_count_old. |
nn_num |
A integer of nearest neighbors, default is 5. |
nn_func |
A character string specifying how the psedo-count to be generated. options include 'mean','median' and 'ransam'. |
numCores |
The number of cores to use |
verbose |
Whether to show running information for srtsim_count |
Returns a SRTsim object with a newly generated count matrix
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) ## Explore the synthetic count matrix simCounts(toySRT)[1:3,1:3]
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) ## Explore the synthetic count matrix simCounts(toySRT)[1:3,1:3]
Generate Data with Estimated Parameters For A New Designed Pattern
srtsim_count_affine( simsrt, reflabel, targetlabel, breaktie = "random", nn_func = c("mean", "median", "ransam"), nn_num = 5, local_sid = NULL, numCores = 1 )
srtsim_count_affine( simsrt, reflabel, targetlabel, breaktie = "random", nn_func = c("mean", "median", "ransam"), nn_num = 5, local_sid = NULL, numCores = 1 )
simsrt |
A SRTsim object with estimated parameters from fitting step |
reflabel |
A character vector specifying labels for reference regions |
targetlabel |
A character vector specifying labels for target regions |
breaktie |
A character string specifying how ties are treated. Same as the "tie.method" in rank function |
nn_func |
A character string specifying how the psedo-count to be generated. options include 'mean','median' and 'ransam'. |
nn_num |
A integer of nearest neighbors, default is 5. |
local_sid |
A numberic seed used locally for the affine transformation. Default is NULL. |
numCores |
A number of cores to use |
Returns a SRTsim object with a newly generated count matrix
## Prepare Data From LIBD Sample subinfo <- exampleLIBD$info[,c("imagecol","imagerow","layer")] colnames(subinfo) <- c("x","y","label") gns <- c("ENSG00000168314","ENSG00000183036", "ENSG00000132639" ) ## Create a simSRT Object with Three Genes For a Fast Example simSRT1 <- createSRT(count_in= exampleLIBD$count[gns,],loc_in =subinfo) ## Estimate model parameters for data generation: domain-specific simSRT1 <- srtsim_fit(simSRT1,sim_schem="domain") ## Define New Layer Structures simSRT1@refcolData$target_label <- "NL1" simSRT1@refcolData$target_label[simSRT1@refcolData$label %in% paste0("Layer",4:5)] <- "NL2" simSRT1@refcolData$target_label[simSRT1@refcolData$label %in% c("Layer6","WM")] <- "NL3" ## Perform Data Generation for New Defined Layer Structures ## Reference: WM --> NL3, Layer5--> NL2, Layer3 --> NL1 simSRT1 <- srtsim_count_affine(simSRT1, reflabel=c("Layer3","Layer5","WM"), targetlabel=c("NL1","NL2","NL3"), nn_func="ransam" ) ## Visualize the Expression Pattern for Gene of Interest visualize_gene(simsrt=simSRT1,plotgn = "ENSG00000168314",rev_y=TRUE,ptsizeCount=1)
## Prepare Data From LIBD Sample subinfo <- exampleLIBD$info[,c("imagecol","imagerow","layer")] colnames(subinfo) <- c("x","y","label") gns <- c("ENSG00000168314","ENSG00000183036", "ENSG00000132639" ) ## Create a simSRT Object with Three Genes For a Fast Example simSRT1 <- createSRT(count_in= exampleLIBD$count[gns,],loc_in =subinfo) ## Estimate model parameters for data generation: domain-specific simSRT1 <- srtsim_fit(simSRT1,sim_schem="domain") ## Define New Layer Structures simSRT1@refcolData$target_label <- "NL1" simSRT1@refcolData$target_label[simSRT1@refcolData$label %in% paste0("Layer",4:5)] <- "NL2" simSRT1@refcolData$target_label[simSRT1@refcolData$label %in% c("Layer6","WM")] <- "NL3" ## Perform Data Generation for New Defined Layer Structures ## Reference: WM --> NL3, Layer5--> NL2, Layer3 --> NL1 simSRT1 <- srtsim_count_affine(simSRT1, reflabel=c("Layer3","Layer5","WM"), targetlabel=c("NL1","NL2","NL3"), nn_func="ransam" ) ## Visualize the Expression Pattern for Gene of Interest visualize_gene(simsrt=simSRT1,plotgn = "ENSG00000168314",rev_y=TRUE,ptsizeCount=1)
Fit the marginal distributions for each row of a count matrix
srtsim_fit( simsrt, marginal = c("auto_choose", "zinb", "nb", "poisson", "zip"), sim_scheme = c("tissue", "domain"), min_nonzero_num = 2, maxiter = 500 )
srtsim_fit( simsrt, marginal = c("auto_choose", "zinb", "nb", "poisson", "zip"), sim_scheme = c("tissue", "domain"), min_nonzero_num = 2, maxiter = 500 )
simsrt |
A SRTsim object |
marginal |
Specification of the types of marginal distribution.Default value is 'auto_choose' which chooses between ZINB, NB, ZIP, and Poisson by a likelihood ratio test (lrt),AIC and whether there is underdispersion.'zinb' will fit the ZINB model. If there is underdispersion, it will fit the Poisson model. If there is no zero at all or an error occurs, it will fit an NB model instead.'nb' fits the NB model and chooses between NB and Poisson depending on whether there is underdispersion. 'poisson' simply fits the Poisson model.'zip' fits the ZIP model and chooses between ZIP and Poisson by a likelihood ratio test |
sim_scheme |
a character string specifying simulation scheme. "tissue" stands for tissue-based simulation; "domain" stands for domain-specific simulation. Default is "tissue". |
min_nonzero_num |
The minimum number of non-zero values required for a gene to be fitted. Default is 2. |
maxiter |
The number of iterations for the model-fitting. Default is 500. |
Returns an object with estimated parameters
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue")
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue")
Fit the marginal distributions for each row of a count matrix
srtsim_newlocs( simsrt, new_loc_num = NULL, loc_lay_out = c("grid", "random"), voting_nn = 3 )
srtsim_newlocs( simsrt, new_loc_num = NULL, loc_lay_out = c("grid", "random"), voting_nn = 3 )
simsrt |
A SRTsim object |
new_loc_num |
A integer specifying the number of spatial locations in the synthetic data |
loc_lay_out |
a character string specifying arrangement of new generated spatial locations. Default is "grid" |
voting_nn |
A integer of nearest neighbors used in label assignment for new generated locations. Default is 3. |
Returns a object with estimated parameters
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Create New Locations within Profile toySRT2 <- srtsim_newlocs(toySRT,new_loc_num=1000) ## Explore New Generated Locations simcolData(toySRT2)
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Create New Locations within Profile toySRT2 <- srtsim_newlocs(toySRT,new_loc_num=1000) ## Explore New Generated Locations simcolData(toySRT2)
Run the SRTsim Shiny Application
SRTsim_shiny()
SRTsim_shiny()
A list that contains a count matrix, a location dataframe, and all parameter specifications.
## Not run: ## Will Load an Interactive Session shinyOutput <- SRTsim_shiny() ## End(Not run)
## Not run: ## Will Load an Interactive Session shinyOutput <- SRTsim_shiny() ## End(Not run)
Subset SRT object based on domain labels of interest
subsetSRT(simsrt, sel_label = NULL)
subsetSRT(simsrt, sel_label = NULL)
simsrt |
A SRTsim object |
sel_label |
A vector of selected domain labels used for the data generation |
Returns a spatialExperiment-based object
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Only Keep the Spatial Locations labelled as "A" in the reference data subtoySRT <- subsetSRT(toySRT,sel_label = "A")
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Only Keep the Spatial Locations labelled as "A" in the reference data subtoySRT <- subsetSRT(toySRT,sel_label = "A")
A data list containing the a gene expression matrix and a location dataframe
toyData
toyData
A data list
A sparse matrix with 100 rows and 251 columns.
A data frame with 251 rows and 3 columns.
created based on a ST Human Breast Cancer data to serve as an example
data(toyData) #Lazy loading. Data becomes visible as soon as called
data(toyData) #Lazy loading. Data becomes visible as soon as called
A list of shiny output
toyShiny
toyShiny
A list of shiny output
A data frame with 150 rows and 980 columns.
A data frame with 980 rows and 4 columns: x, y, group, foldchange
A list of user-specified parameters for count generation
A list of user-specified parameters for pattern design
created based using the SRTsim_shiny()
data(toyShiny) #Lazy loading. Data becomes visible as soon as called
data(toyShiny) #Lazy loading. Data becomes visible as soon as called
Visualize expression pattern for the gene of interest in reference data and synthetic data
visualize_gene( simsrt, plotgn = NULL, ptsizeCount = 2, textsizeCount = 12, rev_y = FALSE, virOption = "D", virDirection = -1 )
visualize_gene( simsrt, plotgn = NULL, ptsizeCount = 2, textsizeCount = 12, rev_y = FALSE, virOption = "D", virDirection = -1 )
simsrt |
A SRTsim object |
plotgn |
A gene name selected for visualization |
ptsizeCount |
Specification of point size. Default is 2. |
textsizeCount |
Specification of axis font size. Default is 12. |
rev_y |
Logical indicating whether to reverse the y axis. Default is FALSE. Useful for Visualize the LIBD data. |
virOption |
Specification of |
virDirection |
Specification of |
Returns two expression plots for the gene of interest
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Create New Locations within Profile toySRT2 <- srtsim_newlocs(toySRT,new_loc_num=1000) ## Estimate model parameters for data generation toySRT2 <- srtsim_fit(toySRT2,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT2 <- srtsim_count(toySRT2,rrr=1) ## compare the expression pattern of HLA-B in synthetic data and reference data visualize_gene(simsrt=toySRT2,plotgn = "HLA-B")
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Create New Locations within Profile toySRT2 <- srtsim_newlocs(toySRT,new_loc_num=1000) ## Estimate model parameters for data generation toySRT2 <- srtsim_fit(toySRT2,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT2 <- srtsim_count(toySRT2,rrr=1) ## compare the expression pattern of HLA-B in synthetic data and reference data visualize_gene(simsrt=toySRT2,plotgn = "HLA-B")
Visualize summarized metrics for reference data and synthetic data
visualize_metrics( simsrt, metric_type = c("all", "genewise", "locwise", "GeneMean", "GeneVar", "GeneCV", "GeneZeroProp", "LocZeroProp", "LocLibSize"), colorpalette = "Set3", axistextsize = 12 )
visualize_metrics( simsrt, metric_type = c("all", "genewise", "locwise", "GeneMean", "GeneVar", "GeneCV", "GeneZeroProp", "LocZeroProp", "LocLibSize"), colorpalette = "Set3", axistextsize = 12 )
simsrt |
A SRTsim object |
metric_type |
Specification of metrics to be plotted. Default value is 'all', which will plot all six metrics: including four gene-wise metrics and two location-wise metrics. "genewise" will produce violin plots for all four gene-wise metrics; "locwise" will produce violin plots for all two location-wise metrics; "GeneMean", "GeneVar", "GeneCV", "GeneZeroProp", "LocZeroProp", and "LocLibSize" will produce single violin plot for the corresponding metric. |
colorpalette |
Specification of color palette to be passed to |
axistextsize |
Specification of axis font size. Default is 12. |
Returns a list of ggplots
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) ## Compute metrics toySRT <- compareSRT(toySRT) ## Visualize Metrics visualize_metrics(toySRT)
## Create a simSRT object toySRT <- createSRT(count_in=toyData$toyCount,loc_in = toyData$toyInfo) set.seed(1) ## Estimate model parameters for data generation toySRT <- srtsim_fit(toySRT,sim_schem="tissue") ## Generate synthetic data with estimated parameters toySRT <- srtsim_count(toySRT) ## Compute metrics toySRT <- compareSRT(toySRT) ## Visualize Metrics visualize_metrics(toySRT)