Get top genes with largest/smallest gradients per cell type
Source:R/interGradientsDL.R
topGradientsCellType.Rd
Retrieve feature names with the largest/smallest gradients per cell
type. These genes can be used to visualize their spatial expression
in the ST data (plotGeneSpatial
function) or to plot the calculated
gradients as a heatmap (plotGradHeatmap
function).
Arguments
- object
SpatialDDLS
object with aDeconvDLModel
object containing gradients in theinterpret.gradients
slot.- method
Method gradients were calculated by. It can be either
'class'
(gradients of predicted classes w.r.t. inputs) or'loss'
(gradients of loss w.r.t. input features).- top.n.genes
Top n genes (positive and negative) taken per cell type.
Examples
# \donttest{
set.seed(123)
sce <- SingleCellExperiment::SingleCellExperiment(
assays = list(
counts = matrix(
rpois(30, lambda = 5), nrow = 15, ncol = 10,
dimnames = list(paste0("Gene", seq(15)), paste0("RHC", seq(10)))
)
),
colData = data.frame(
Cell_ID = paste0("RHC", seq(10)),
Cell_Type = sample(x = paste0("CellType", seq(2)), size = 10,
replace = TRUE)
),
rowData = data.frame(
Gene_ID = paste0("Gene", seq(15))
)
)
SDDLS <- createSpatialDDLSobject(
sc.data = sce,
sc.cell.ID.column = "Cell_ID",
sc.gene.ID.column = "Gene_ID",
sc.filt.genes.cluster = FALSE
)
#> === Spatial transcriptomics data not provided
#> === Processing single-cell data
#> - Filtering features:
#> - Selected features: 15
#> - Discarded features: 0
#>
#> === No mitochondrial genes were found by using ^mt- as regrex
#>
#> === Final number of dimensions for further analyses: 15
SDDLS <- genMixedCellProp(
object = SDDLS,
cell.ID.column = "Cell_ID",
cell.type.column = "Cell_Type",
num.sim.spots = 50,
train.freq.cells = 2/3,
train.freq.spots = 2/3,
verbose = TRUE
)
#>
#> === The number of mixed profiles that will be generated is equal to 50
#>
#> === Training set cells by type:
#> - CellType1: 4
#> - CellType2: 3
#> === Test set cells by type:
#> - CellType1: 2
#> - CellType2: 1
#> === Probability matrix for training data:
#> - Mixed spots: 34
#> - Cell types: 2
#> === Probability matrix for test data:
#> - Mixed spots: 16
#> - Cell types: 2
#> DONE
SDDLS <- simMixedProfiles(SDDLS)
#> === Setting parallel environment to 1 thread(s)
#>
#> === Generating train mixed profiles:
#>
#> === Generating test mixed profiles:
#>
#> DONE
SDDLS <- trainDeconvModel(
object = SDDLS,
batch.size = 12,
num.epochs = 5
)
#> === Training and test from stored data
#> Using only simulated mixed samples
#> Using only simulated mixed samples
#> Model: "SpatialDDLS"
#> _____________________________________________________________________
#> Layer (type) Output Shape Param #
#> =====================================================================
#> Dense1 (Dense) (None, 200) 3200
#> _____________________________________________________________________
#> BatchNormalization1 (BatchNorm (None, 200) 800
#> _____________________________________________________________________
#> Activation1 (Activation) (None, 200) 0
#> _____________________________________________________________________
#> Dropout1 (Dropout) (None, 200) 0
#> _____________________________________________________________________
#> Dense2 (Dense) (None, 200) 40200
#> _____________________________________________________________________
#> BatchNormalization2 (BatchNorm (None, 200) 800
#> _____________________________________________________________________
#> Activation2 (Activation) (None, 200) 0
#> _____________________________________________________________________
#> Dropout2 (Dropout) (None, 200) 0
#> _____________________________________________________________________
#> Dense3 (Dense) (None, 2) 402
#> _____________________________________________________________________
#> BatchNormalization3 (BatchNorm (None, 2) 8
#> _____________________________________________________________________
#> ActivationSoftmax (Activation) (None, 2) 0
#> =====================================================================
#> Total params: 45,410
#> Trainable params: 44,606
#> Non-trainable params: 804
#> _____________________________________________________________________
#>
#> === Training DNN with 34 samples:
#>
#> === Evaluating DNN in test data (16 samples)
#> - loss: NaN
#> - accuracy: 0.5
#> - mean_absolute_error: NaN
#> - categorical_accuracy: 0.5
#>
#> === Generating prediction results using test data
#> DONE
## calculating gradients
SDDLS <- interGradientsDL(SDDLS)
listGradients <- topGradientsCellType(SDDLS)
lapply(listGradients, head, n = 5)
#> $CellType1
#> $CellType1$Absolute
#> [1] "Gene12" "Gene4" "Gene13" "Gene2" "Gene11" "Gene14" "Gene8" "Gene3"
#> [9] "Gene5" "Gene6" "Gene10" "Gene9" "Gene7" "Gene15" "Gene1"
#>
#> $CellType1$Positive
#> [1] "Gene12" "Gene4" "Gene13" "Gene2" "Gene11" "Gene14" "Gene8" "Gene3"
#> [9] "Gene5" "Gene6" "Gene10" "Gene9" "Gene7" "Gene15" "Gene1"
#>
#> $CellType1$Negative
#> [1] "Gene12" "Gene4" "Gene13" "Gene2" "Gene11" "Gene14" "Gene8" "Gene3"
#> [9] "Gene5" "Gene6" "Gene10" "Gene9" "Gene7" "Gene15" "Gene1"
#>
#>
#> $CellType2
#> $CellType2$Absolute
#> [1] "Gene15" "Gene10" "Gene9" "Gene7" "Gene5" "Gene3" "Gene8" "Gene1"
#> [9] "Gene6" "Gene13" "Gene4" "Gene14" "Gene2" "Gene12" "Gene11"
#>
#> $CellType2$Positive
#> [1] "Gene15" "Gene10" "Gene9" "Gene7" "Gene5" "Gene3" "Gene8" "Gene1"
#> [9] "Gene6" "Gene13" "Gene4" "Gene14" "Gene2" "Gene12" "Gene11"
#>
#> $CellType2$Negative
#> [1] "Gene15" "Gene10" "Gene9" "Gene7" "Gene5" "Gene3" "Gene8" "Gene1"
#> [9] "Gene6" "Gene13" "Gene4" "Gene14" "Gene2" "Gene12" "Gene11"
#>
#>
# }