DEtime_infer: Inferring perturbation time point from biological time course...

View source: R/DEtime_infer.R

DEtime_inferR Documentation

Inferring perturbation time point from biological time course data

Description

This is the main function in DEtime Package, which applies a mixedGP kernel to time course data under control and perturbed conditions. It returns the posterior distribution of these predefined perturbation time candidates and relevant statistical estimations of the inferred perturbation time point.

Usage

DEtime_infer(ControlTimes, ControlData, PerturbedTimes, PerturbedData,
  TestTimes = NULL, gene_ID = NULL, bound.lengthscale = NULL)

Arguments

ControlTimes

Experimental time point at which time course biological data for the control case are measured, they have to be repeated if there are replicated measurements

ControlData

Time course data measured under control condition

PerturbedTimes

Experimental time point at which time course biological data for the perturbed case are measured, they have to be repeated if there are replicated measurements

PerturbedData

Time course data measured under perturbed condition

TestTimes

The predefined evenly spaced time points upon which perturbation will be evaluated. If undefined, TestTimes <- seq(min(c(ControlTimes, PerturbedTimes)), max(c(ControlTimes, PerturbedTimes)),length=50)

gene_ID

ID of these genes addressed in this study. If undefinied, numbers will be used instead

bound.lengthscale

bounds for the lengthscale used in the DEtime RBF kernel. When not provided,bound.lengthscale <- c(min(ControlTimes,PerturbedTimes), 4*max(c(ControlTimes,PerturbedTimes)))

Details

ControlTimes and PerturbedTimes can be ordered by either time series, for instance time1, time1, time2, time2, time3, time3 ... or replicate sequences, for instance: time1, time2, time3, time1, time2, time3. ControlData and PerturbedData are two matrices where each row represents the time course data for one particular gene under either control or perturbed condition. The orders of the ControlData and PeruturbedData have to match those of the ControlTimes and PerturbedTimes, respectively.

Value

The function will return a DEtimeOutput object which includes:

  1. result: the statistical estimation for the inferred perturbation time

  2. ____$MAP: maximum a posterior solution to the inferred perturbation time

  3. ____$mean: mean of the posterior distribution of the inferred perturbation time

  4. ____$median: median of the posterior distribution of the inferred perturbation time

  5. ____$ptl5: 5 percentile of the posterior distribution of the inferred perturbation time

  6. _____$ptl95: 95 percentile of the posterior distribution of the inferred perturbation time

  7. posterior: posterior distribution of the tested perturbation time points

  8. model: optimized GP model which will be used for later GP regression work

  9. best_param: optimized hyperparameter for the optimized GP model

  10. ControlTimes: original experimental time points for the control case which will be used for future print or plot functions

  11. ControlData: original measured time course data for the control case which will be used for future print or plot functions

  12. PerturbedTimes: original experimental time points for the perturbed case which will be used for future print or plot functions

  13. PerturbedData: original measured time course data for the control case which will be used for future print or plot functions

  14. TestTimes: tested perturbation time points

  15. gene_ID: the ID of genes for the data

Examples

### import simulated data
data(SimulatedData)
### start perturbation time inference
res <- DEtime_infer(ControlTimes = ControlTimes, ControlData = ControlData, 
PerturbedTimes = PerturbedTimes, PerturbedData = PerturbedData)

ManchesterBioinference/DEtime documentation built on Feb. 9, 2024, 12:10 p.m.