test_differential_cellularity() takes as input A `tbl` (with at least three columns for sample, feature and transcript abundance) or `SummarizedExperiment` (more convenient if abstracted to tibble with library(tidySummarizedExperiment)) and returns a consistent object (to the input) with additional columns for the statistics from the hypothesis test.

test_differential_cellularity(
  .data,
  .formula,
  .sample = NULL,
  .transcript = NULL,
  .abundance = NULL,
  method = "cibersort",
  reference = X_cibersort,
  significance_threshold = 0.05,
  ...
)

# S4 method for class 'spec_tbl_df'
test_differential_cellularity(
  .data,
  .formula,
  .sample = NULL,
  .transcript = NULL,
  .abundance = NULL,
  method = "cibersort",
  reference = X_cibersort,
  significance_threshold = 0.05,
  ...
)

# S4 method for class 'tbl_df'
test_differential_cellularity(
  .data,
  .formula,
  .sample = NULL,
  .transcript = NULL,
  .abundance = NULL,
  method = "cibersort",
  reference = X_cibersort,
  significance_threshold = 0.05,
  ...
)

# S4 method for class 'tidybulk'
test_differential_cellularity(
  .data,
  .formula,
  .sample = NULL,
  .transcript = NULL,
  .abundance = NULL,
  method = "cibersort",
  reference = X_cibersort,
  significance_threshold = 0.05,
  ...
)

# S4 method for class 'SummarizedExperiment'
test_differential_cellularity(
  .data,
  .formula,
  .sample = NULL,
  .transcript = NULL,
  .abundance = NULL,
  method = "cibersort",
  reference = X_cibersort,
  significance_threshold = 0.05,
  ...
)

# S4 method for class 'RangedSummarizedExperiment'
test_differential_cellularity(
  .data,
  .formula,
  .sample = NULL,
  .transcript = NULL,
  .abundance = NULL,
  method = "cibersort",
  reference = X_cibersort,
  significance_threshold = 0.05,
  ...
)

Arguments

.data

A `tbl` (with at least three columns for sample, feature and transcript abundance) or `SummarizedExperiment` (more convenient if abstracted to tibble with library(tidySummarizedExperiment))

.formula

A formula representing the desired linear model. The formula can be of two forms: multivariable (recommended) or univariable Respectively: \"factor_of_interest ~ .\" or \". ~ factor_of_interest\". The dot represents cell-type proportions, and it is mandatory. If censored regression is desired (coxph) the formula should be of the form \"survival::Surv\(y, dead\) ~ .\"

.sample

The name of the sample column

.transcript

The name of the transcript/gene column

.abundance

The name of the transcript/gene abundance column

method

A string character. Either \"cibersort\", \"epic\" or \"llsr\". The regression method will be chosen based on being multivariable: lm or cox-regression (both on logit-transformed proportions); or univariable: beta or cox-regression (on logit-transformed proportions). See .formula for multi- or univariable choice.

reference

A data frame. The transcript/cell_type data frame of integer transcript abundance

significance_threshold

A real between 0 and 1 (usually 0.05).

...

Further parameters passed to the method deconvolve_cellularity

Value

A consistent object (to the input) with additional columns for the statistics from the hypothesis test (e.g., log fold change, p-value and false discovery rate).

A `SummarizedExperiment` object

A `SummarizedExperiment` object

Details

`r lifecycle::badge("maturing")`

This routine applies a deconvolution method (e.g., Cibersort; DOI: 10.1038/nmeth.3337) and passes the proportions inferred into a generalised linear model (DOI:dx.doi.org/10.1007/s11749-010-0189-z) or a cox regression model (ISBN: 978-1-4757-3294-8)

Underlying method for the generalised linear model: data |> deconvolve_cellularity( !!.sample, !!.transcript, !!.abundance, method=method, reference = reference, action="get", ... ) [..] betareg::betareg(.my_formula, .)

Underlying method for the cox regression: data |> deconvolve_cellularity( !!.sample, !!.transcript, !!.abundance, method=method, reference = reference, action="get", ... ) [..] mutate(.proportion_0_corrected = .proportion_0_corrected |> boot::logit()) survival::coxph(.my_formula, .)

Examples


 # Regular regression
  test_differential_cellularity(
   tidybulk::se_mini ,
      . ~ condition,
      cores = 1
  )
#> # A tibble: 22 × 7
#>    .cell_type                       cell_type_proportions `estimate_(Intercept)`
#>    <chr>                            <list>                                 <dbl>
#>  1 cibersort.B.cells.naive          <tibble [5 × 8]>                       -2.37
#>  2 cibersort.B.cells.memory         <tibble [5 × 8]>                       -2.75
#>  3 cibersort.Plasma.cells           <tibble [5 × 8]>                       NA   
#>  4 cibersort.T.cells.CD8            <tibble [5 × 8]>                       NA   
#>  5 cibersort.T.cells.CD4.naive      <tibble [5 × 8]>                       NA   
#>  6 cibersort.T.cells.CD4.memory.re… <tibble [5 × 8]>                       -3.81
#>  7 cibersort.T.cells.CD4.memory.ac… <tibble [5 × 8]>                       NA   
#>  8 cibersort.T.cells.follicular.he… <tibble [5 × 8]>                       -6.01
#>  9 cibersort.T.cells.regulatory..T… <tibble [5 × 8]>                       NA   
#> 10 cibersort.T.cells.gamma.delta    <tibble [5 × 8]>                       NA   
#> # ℹ 12 more rows
#> # ℹ 4 more variables: estimate_conditionTRUE <dbl>,
#> #   std.error_conditionTRUE <dbl>, statistic_conditionTRUE <dbl>,
#> #   p.value_conditionTRUE <dbl>

  # Cox regression - multiple

tidybulk::se_mini |>

  # Test
  test_differential_cellularity(
      survival::Surv(days, dead) ~ .,
      cores = 1
  )
#> Warning: Ran out of iterations and did not converge
#> # A tibble: 22 × 6
#>    .cell_type         cell_type_proportions estimate std.error statistic p.value
#>    <chr>              <list>                   <dbl>     <dbl>     <dbl>   <dbl>
#>  1 cibersort.B.cells… <tibble [5 × 7]>        -0.616     3148.  -1.96e-4   1.00 
#>  2 cibersort.B.cells… <tibble [5 × 7]>      -100.        2796.  -3.58e-2   0.971
#>  3 cibersort.Plasma.… <tibble [5 × 7]>         0            0  NaN       NaN    
#>  4 cibersort.T.cells… <tibble [5 × 7]>       118.        4318.   2.74e-2   0.978
#>  5 cibersort.T.cells… <tibble [5 × 7]>      -164.        3555.  -4.61e-2   0.963
#>  6 cibersort.T.cells… <tibble [5 × 7]>        -1.93      2430.  -7.96e-4   0.999
#>  7 cibersort.T.cells… <tibble [5 × 7]>         0.935    13450.   6.95e-5   1.00 
#>  8 cibersort.T.cells… <tibble [5 × 7]>        13.3       5004.   2.67e-3   0.998
#>  9 cibersort.T.cells… <tibble [5 × 7]>      -302.       49862.  -6.07e-3   0.995
#> 10 cibersort.T.cells… <tibble [5 × 7]>         0            0  NaN       NaN    
#> # ℹ 12 more rows