Package 'PICS'

Title: Probabilistic inference of ChIP-seq
Description: Probabilistic inference of ChIP-Seq using an empirical Bayes mixture model approach.
Authors: Xuekui Zhang <[email protected]>, Raphael Gottardo <[email protected]>
Maintainer: Renan Sauteraud <[email protected]>
License: Artistic-2.0
Version: 2.49.0
Built: 2024-07-03 06:03:56 UTC
Source: https://github.com/bioc/PICS

Help Index


Pre-process bam files

Description

Reads a bam file using Rsamtools and extract the reads for each chromosome.

Usage

bam2gr(bamFile, chr = NULL, PE = FALSE, verbose = FALSE)

Arguments

bamFile

A character. The name of the .bam file to read.

chr

A character. An optional character string. If specified, only the selected chromosome will be returned. Speed up the computation.

PE

A logical. Set to TRUE for paired-end sequencing data.

verbose

A logical. Print additional information about the data.

Value

A GRanges of all the reads for each chromosome.

Note

The user might encounter a memory allocation error when using bam files of bigger sizes. Splitting the file by chromosome before calling bam2gr will solve this issue.

For Paired-End data, non matched reads are discarded.

See Also

segmentPICS


Identify candidate regions

Description

Identify candidate regions from paired-end sequencing data.

Usage

candidate.region(PE.RD, islandDepth, min_cut, max_cut)

Arguments

PE.RD

A character. File including sequenced pairs.

islandDepth

A numeric. The minimum depth of candidate regions.

min_cut

A numeric. The minimum region length.

max_cut

A numeric. The maximum region length.


Export a PICS object to GRanges

Description

Export a PICS object to GRanges

Usage

makeGRangesOutput(
  obj,
  type = "fixed",
  length = 100,
  filter = list(delta = c(0, Inf), se = c(0, Inf), sigmaSqF = c(0, Inf), sigmaSqR = c(0,
    Inf), score = c(0, Inf))
)

Arguments

obj

A PICS object. The output of the PICS function.

type

A character. One of "fixed", "bed", "ci" or "wig".

length

A numeric. The length of the region around the center (mu). Only used when type = "fixed.

filter

A list. Additional filtering options.

Value

A GRanges object.

See Also

PICS


Estimation of binding site positions

Description

This object contains Estimation of binding site positions and has the following slots: segReadsList, dataType.

Usage

PICS(
  segReadsList,
  dataType = NULL,
  paraEM = NULL,
  paraPrior = NULL,
  nCores = 1
)

Arguments

segReadsList

This object contains segmentation of Genome

dataType

A character. The type of data you are processing: 'TF' for transcription factor.

paraEM

A list of parameters for the EM algorithm as returned by the setParaEM function. The default parameters should be good enough for most usages.

paraPrior

A list of parameters for the priors as returned by setParaPrior.

nCores

An integer. The number of cores that should be used in parallel by the function.

Value

An object of class picsList containing the estimated binding site positions.


PICS class

Description

This object is used to gather all parameters from fitting PICS to a single candidate region.

Usage

## S4 method for signature 'pics'
show(object)

## S4 method for signature 'pics'
minRange(x)

## S4 method for signature 'pics'
maxRange(x)

## S4 method for signature 'pics'
score(x)

## S4 method for signature 'pics'
scoreReverse(x)

## S4 method for signature 'pics'
scoreForward(x)

## S4 method for signature 'pics'
chromosome(x)

## S4 method for signature 'pics'
se(x)

## S4 method for signature 'pics'
seF(x)

## S4 method for signature 'pics'
seR(x)

## S4 method for signature 'pics'
sigmaSqF(x)

## S4 method for signature 'pics'
sigmaSqR(x)

## S4 method for signature 'pics'
delta(x)

## S4 method for signature 'pics'
mu(x)

## S4 method for signature 'pics'
w(x)

## S4 method for signature 'pics'
K(x)

## S4 method for signature 'pics'
code(x)

## S4 method for signature 'pics'
summary(object)

Arguments

object, x

A pics object.

Functions

  • show,pics-method: show method

  • minRange,pics-method: Get start of range

  • maxRange,pics-method: Get end of range

  • score,pics-method: Score accessor.

  • scoreReverse,pics-method: Reverse score accessor.

  • scoreForward,pics-method: Forward score accessor.

  • chromosome,pics-method: Chromosome accessor

  • se,pics-method: se accessor

  • seF,pics-method: Forward se accessor

  • seR,pics-method: Reverse se accessor

  • sigmaSqF,pics-method: sigmaSqF accessor

  • sigmaSqR,pics-method: sigmaSqR accessor

  • delta,pics-method: delta accessor

  • mu,pics-method: mu accessor

  • w,pics-method: w accessor

  • K,pics-method: K accessor

  • code,pics-method: Error code accessor

  • summary,pics-method: Summary of the object

Slots

estimates

A list containing all parameters estimates as well as standard errors.

Nmerged

The number of binding events that were merged; binding events that overlap are merged.

converge

A logical value indicating whether the EM algorithm has converged.

chr

The candidate region's chromosome.

score

Score of the binding event

scoreF

Forward score of the binding event.

scoreR

Reverse score of the binding event.

range

Genomic ranges.


Generics associated with pics classes

Description

These generics are used in methods of pics, picsError and picsList classes. See class help pages for method documentation.

Usage

minRange(x, ...)

maxRange(x, ...)

score(x, ...)

scoreReverse(x, ...)

scoreForward(x, ...)

chromosome(x, ...)

se(x, ...)

seF(x, ...)

seR(x, ...)

sigmaSqF(x, ...)

sigmaSqR(x, ...)

delta(x, ...)

mu(x, ...)

w(x, ...)

K(x, ...)

code(x, ...)

## S4 method for signature 'data.frame'
score(x)

## S4 method for signature 'data.frame'
scoreReverse(x)

## S4 method for signature 'data.frame'
scoreForward(x)

## S4 method for signature 'data.frame'
chromosome(x)

## S4 method for signature 'data.frame'
se(x)

## S4 method for signature 'data.frame'
seF(x)

## S4 method for signature 'data.frame'
seR(x)

## S4 method for signature 'data.frame'
sigmaSqF(x)

## S4 method for signature 'data.frame'
sigmaSqR(x)

## S4 method for signature 'data.frame'
delta(x)

## S4 method for signature 'data.frame'
mu(x)

Arguments

x, ...

Object and argumemts passed to the methods.

Functions

  • minRange: Start accessor

  • maxRange: End accessor

  • score: Score accessor

  • scoreReverse: Reverse score accessor

  • scoreForward: Forward score accesor

  • chromosome: Chromosome accessor

  • se: se accessor

  • seF: Forward se accessor

  • seR: Reverse se accessor

  • sigmaSqF: sigmaSqF accessor

  • sigmaSqR: sigmaSqR accessor

  • delta: delta accessor

  • mu: mu accessor

  • w: w accessor

  • K: K accessor

  • code: Return error codes

  • score,data.frame-method: Score accessor

  • scoreReverse,data.frame-method: Reverse score accessor

  • scoreForward,data.frame-method: Forward score accessorse

  • chromosome,data.frame-method: chromosome accessor

  • se,data.frame-method: se accessor

  • seF,data.frame-method: seF accessor

  • seR,data.frame-method: seR accessor

  • sigmaSqF,data.frame-method: sigmaSqF accessor

  • sigmaSqR,data.frame-method: sigmaSqR accessor

  • delta,data.frame-method: delta accessor

  • mu,data.frame-method: mu accessor


picsError class

Description

This class is used in cases when the algorithm does not converge.

Usage

## S4 method for signature 'picsError'
show(object)

## S4 method for signature 'picsError'
minRange(x)

## S4 method for signature 'picsError'
maxRange(x)

## S4 method for signature 'picsError'
score(x)

## S4 method for signature 'picsError'
scoreReverse(x)

## S4 method for signature 'picsError'
scoreForward(x)

## S4 method for signature 'picsError'
chromosome(x)

## S4 method for signature 'picsError'
se(x)

## S4 method for signature 'picsError'
seF(x)

## S4 method for signature 'picsError'
seR(x)

## S4 method for signature 'picsError'
sigmaSqF(x)

## S4 method for signature 'picsError'
sigmaSqR(x)

## S4 method for signature 'picsError'
delta(x)

## S4 method for signature 'picsError'
mu(x)

## S4 method for signature 'picsError'
w(x)

## S4 method for signature 'picsError'
K(x)

## S4 method for signature 'picsError'
code(x)

Arguments

object, x

A picsError object.

Functions

  • show,picsError-method: show method

  • minRange,picsError-method: Get start of range

  • maxRange,picsError-method: Get end of range

  • score,picsError-method: Score accessor.

  • scoreReverse,picsError-method: Reverse score accessor.

  • scoreForward,picsError-method: Forward score accessor.

  • chromosome,picsError-method: Chromosome accessor

  • se,picsError-method: se accessor

  • seF,picsError-method: Forward se accessor

  • seR,picsError-method: Reverse se accessor

  • sigmaSqF,picsError-method: sigmaSqF accessor

  • sigmaSqR,picsError-method: sigmaSqR accessor

  • delta,picsError-method: delta accessor

  • mu,picsError-method: mu accessor

  • w,picsError-method: w accessor

  • K,picsError-method: K accessor

  • code,picsError-method: Error code accessor

Slots

errorCode

The error code for debugging.


Estimate the FDR

Description

Estimate the false detection rate for an object of class pics or picsList.

Usage

picsFDR(
  picsIP,
  picsCont,
  filter = list(delta = c(0, Inf), se = c(0, Inf), sigmaSqF = c(0, Inf), sigmaSqR = c(0,
    Inf))
)

Arguments

picsIP

An object of class pics or picsList containing the informations for the IP reads.

picsCont

An object of class pics or picsList containing the informations for the control reads

filter

A list of ranges for filtering regions based on PICS parameters. By default filter is set to NULL and all regions are used.

  • delta: Length of the binding sites

  • se: Standard error

  • sigmaSqF: Forward peak variance

  • sigmaSqR: Reverse peak variance

Value

A data.frame with the following columns: FDR, score, N

See Also

picsList pics


List of PICS objects

Description

List of PICS objects

Usage

## S4 method for signature 'picsList'
show(object)

## S4 method for signature 'picsList'
minRange(x)

## S4 method for signature 'picsList'
maxRange(x)

## S4 method for signature 'picsList'
score(x)

## S4 method for signature 'picsList'
scoreReverse(x)

## S4 method for signature 'picsList'
scoreForward(x)

## S4 method for signature 'picsList'
chromosome(x)

## S4 method for signature 'picsList'
se(x)

## S4 method for signature 'picsList'
seF(x)

## S4 method for signature 'picsList'
seR(x)

## S4 method for signature 'picsList'
sigmaSqF(x)

## S4 method for signature 'picsList'
sigmaSqR(x)

## S4 method for signature 'picsList'
delta(x)

## S4 method for signature 'picsList'
mu(x)

## S4 method for signature 'picsList'
w(x)

## S4 method for signature 'picsList'
K(x)

## S4 method for signature 'picsList'
code(x)

## S4 method for signature 'picsList'
length(x)

## S4 method for signature 'picsList'
summary(object)

## S4 method for signature 'picsList,ANY,ANY,ANY'
x[i, j, ..., drop = FALSE]

## S4 method for signature 'picsList,ANY,ANY'
x[[i, j, ..., exact = TRUE]]

Arguments

object, x

A pics object.

i, j, ..., drop, exact

Arguments passed to subset functions

Functions

  • show,picsList-method: show method

  • minRange,picsList-method: Get start of range

  • maxRange,picsList-method: Get end of range

  • score,picsList-method: Score accessor.

  • scoreReverse,picsList-method: Reverse score accessor.

  • scoreForward,picsList-method: Forward score accessor.

  • chromosome,picsList-method: Chromosome accessor

  • se,picsList-method: se accessor

  • seF,picsList-method: Forward se accessor

  • seR,picsList-method: Reverse se accessor

  • sigmaSqF,picsList-method: sigmaSqF accessor

  • sigmaSqR,picsList-method: sigmaSqR accessor

  • delta,picsList-method: delta accessor

  • mu,picsList-method: mu accessor

  • w,picsList-method: w accessor

  • K,picsList-method: K accessor

  • code,picsList-method: Error code accessor

  • length,picsList-method: Return the length of the object

  • summary,picsList-method: Summary of the object

  • [,picsList,ANY,ANY,ANY-method: Subset list

  • [[,picsList,ANY,ANY-method: Subset element


Plot PICS FDR

Description

This method plots a curve showing the FDR as a function of the PICS scores.

Usage

## S4 method for signature 'picsList,picsList'
plot(x, y, filter = NULL, h = 0.1, ...)

Arguments

x

A picsList object as returned by the function PICS run on the treatment data.

y

A picsList object as returned by the function PICS run on the control data.

filter

A list of ranges for filtering regions based on PICS parameters. By default filter is set to 'NULL' and all regions are used.

h

A value between 0 and 1, representing the desired FDR. This simply draws a horizontal line at the given value.

...

Further graphical parameters passed to the generic function plot.

See Also

PICS


Plot methods for PICS objects

Description

Methods to plot pics and segReads objects and derived classes.

Usage

## S4 method for signature 'pics,segReads'
plot(
  x,
  y,
  addKernel = FALSE,
  addNucleosome = FALSE,
  addSe = TRUE,
  main = NULL,
  ...
)

## S4 method for signature 'picsError,segReads'
plot(x, y, addKernel = FALSE, main = NULL, ...)

## S4 method for signature 'picsList,segReadsList'
plot(
  x,
  y,
  regionIndex = NULL,
  addKernel = FALSE,
  addNucleosome = FALSE,
  addSe = TRUE,
  main = NULL,
  ...
)

Arguments

x, y

Objects

addKernel

Add kernel density estimate to the plot.

addNucleosome

Add a nucleosome track to the plot.

addSe

Add standard error to the plot.

main

Main title.

...

Arguments to be passed to the plot method.

regionIndex

Add region index to the plot.

Functions

  • plot,pics,segReads-method: Plot method for pics and segReads

  • plot,picsError,segReads-method: Plot method for picsError and segReads

  • plot,picsList,segReadsList-method: Plot method for picsList and segReadsList


Segmentation of paired-end sequencing data

Description

These two functions are part of the segmentation step for paired-end sequencing data and are exported to be used in PING package.


Segment the genome into candidate regions

Description

Pre-process bidirectional aligned reads data from a single ChIP-Seq experiment to detect candidate regions with a minimum number of forward and reverse reads. These candidate regions will then be processed by PICS.

Usage

segmentPICS(
  data,
  dataC = NULL,
  map = NULL,
  minReads = 2,
  minReadsInRegion = 3,
  jitter = FALSE,
  dataType = "TF",
  maxLregion = 0,
  minLregion = 100
)

Arguments

data

A GRanges object containing the IP reads. See details for more information on how to set up the data.

dataC

A GRanges object containing the control reads. Set to NULL by default, i.e. no control.

map

A GRanges object containing the mappability profiles. Set to NULL by default, i.e. no profiles.

minReads

A numeric. The minimum number of F/R reads to be present in the sliding window.

minReadsInRegion

A numeric. The minimum number of F/R reads to be present in the region.

jitter

A logical value stating whether some noise should be added to the read locations. This is recommended if the read positions have lots of duplicates.

dataType

A character. Type of experiment. "TF" or "H".

maxLregion

A numeric. The maximum length.

minLregion

A numeric. The minimum length.

Value

An object of class segReadsList containing the results for all pre-processed regions.

References

X. Zhang, G. Robertson, M. Krzywinski, K. Ning, A. Droit, S. Jones, and R. Gottardo, “PICS: Probabilistic Inference for ChIP-seq” arXiv, 0903.3206, 2009.

See Also

segReadsList

Examples

# Read data
path<-system.file("extdata",package="PICS")
## Note that the col name for the chromosome needs to be space and not chr
dataIP <- read.table(file.path(path, "Treatment_tags_chr21_sort.bed"), header=TRUE,
                     colClasses = c("factor","integer","integer","factor"))
dataIP <- as(dataIP, "GRanges")

dataCont <- read.table(file.path(path, "Input_tags_chr21_sort.bed"), header=TRUE,
                       colClasses = c("factor","integer","integer","factor"))
dataCont <- as(dataCont, "GRanges")

map <- read.table(file.path(path, "mapProfileShort"), header=TRUE,
                  colClasses = c("factor","integer","integer","NULL"))
map <- as(map, "GRanges")
seg <- segmentPICS(dataIP, dataC = dataCont, map = map, minReads = 1)

Classes and functions to segment the genome in candidate regions

Description

Pre-process bidirectional aligned reads data from a single ChIP-Seq experiment to detect candidate regions with a minimum number of forward and reverse reads. These candidate regions will then be processed by PICS.

Usage

segReads(yF, yR, cF, cR, map, chr)

segReadsList(List, paraSW, N, Nc)

## S4 method for signature 'segReads'
show(object)

## S4 method for signature 'segReadsList'
show(object)

map(x, ...)

## S4 method for signature 'segReads'
map(x)

## S4 method for signature 'segReadsList'
map(x)

## S4 method for signature 'segReadsList'
length(x)

## S4 method for signature 'segReadsList'
summary(object)

## S4 method for signature 'segReads'
summary(object)

## S4 method for signature 'segReadsList,ANY,ANY,ANY'
x[i, j, ..., drop = FALSE]

## S4 method for signature 'segReadsList,ANY,ANY'
x[[i, j, ..., exact = TRUE]]

Arguments

yF

A numeric vector. Forward reads.

yR

A numeric vector. Reverse reads.

cF

A numeric vector. Forward reads for the controls.

cR

A numeric vector. Reverse reads for the controls.

map

A matrix. The mappability profile.

chr

A character. Chromosome name.

List

A list of segReads objects.

paraSW

A list of parameters for the genomic regions.

N

A numeric. The number of reads in the data.

Nc

A numeric. The number of reads in the control.

object, x

A segReads object.

i, j, ..., exact, drop

Additional arguments passed to subset methods.

Functions

  • segReads: segReads Constructor

  • segReadsList: segReadsList Constructor

  • show,segReads-method: show method

  • show,segReadsList-method: show method

  • map: map generic

  • map,segReads-method: map method

  • map,segReadsList-method: map method

  • length,segReadsList-method: Return length of segReadsList

  • summary,segReadsList-method: Summary method

  • summary,segReads-method: Summary method

  • [,segReadsList,ANY,ANY,ANY-method: Subset methods

  • [[,segReadsList,ANY,ANY-method: Subset methods

Note

segReads and segReadsList objects are not meant to be built via the constructors. The constructors are used in segmentPICS.


Perform genome segmentation depending

Description

Perform genome segmentation depending

Usage

segReadsGeneric(
  data,
  dataC = NULL,
  map = NULL,
  minReads = 2,
  minReadsInRegion = 3,
  jitter = FALSE,
  maxLregion = 0,
  minLregion = 100,
  step = 20,
  width = 250,
  package = "PICS"
)

Arguments

data

A GRanges object containing the IP reads. See details for more information on how to set up the data.

dataC

A GRanges object containing the control reads. Set to NULL by default, i.e. no control.

map

A GRanges object containing the mappability profiles. Set to NULL by default, i.e. no profiles.

minReads

A numeric. The minimum number of F/R reads to be present in the sliding window.

minReadsInRegion

A numeric. The minimum number of F/R reads to be present in the region.

jitter

A logical value stating whether some noise should be added to the read locations. This is recommended if the read positions have lots of duplicates.

maxLregion

A numeric. The maximum length.

minLregion

A numeric. The minimum length.

step

A numeric. The increment of the sliding window.

width

A numeric. The width of the region.

package

A character. "PICS" or "PING"


Class and methods for list of candidate regions from paired-end data

Description

Class and methods for list of candidate regions from paired-end data

Usage

segReadsListPE(List, paraSW, N, NFm, NRm, Nc, NcFm, NcRm)

## S4 method for signature 'segReadsListPE,ANY,ANY,ANY'
x[i, j, ..., drop = FALSE]

## S4 method for signature 'segReadsListPE,ANY,ANY'
x[[i, j, ..., exact = TRUE]]

Arguments

List

A list of segReadsPE objects.

paraSW

A list of parameters for the genomic regions.

N, NFm, NRm

Read counts in the data.

Nc, NcFm, NcRm

Read counts in the control.

x, i, j, ..., drop, exact

Arguments passed to subset methods

Functions

  • segReadsListPE: segReadsListPE constructor.

  • [,segReadsListPE,ANY,ANY,ANY-method: subset method

  • [[,segReadsListPE,ANY,ANY-method: subset method


Classe and methods for candidate regions from paired-end data

Description

A segReadsPE object represents a single candoidate region, including all its informative reads and mappability profile.

Usage

segReadsPE(yF, yR, yFm, yRm, cF, cR, cFm, cRm, map, chr)

Arguments

yF

A numeric vector. Forward reads.

yR

A numeric vector. Reverse reads.

yFm

A numeric vector. Forward reads on paired end.

yRm

A numeric vector. Reverse reads on paired end.

cF

A numeric vector. Forward reads for the controls.

cR

A numeric vector. Reverse reads for the controls.

cFm

A numeric vector. Forward reads on paired end for the controls.

cRm

A numeric vector. Reverse reads on paired end for the controls.

map

A matrix. The mappability profile.

chr

A character. Chromosome name.

Functions

  • segReadsPE: segReadsPE constructor.

Note

segReadsPE objects are not meant to be built via the constructors. The constructors are used in segmentPICS.


List of parameters for the EM algorithm that can be used as an argument of PICS.

Description

This function takes from 0 to 7 EM algorithm parameters as argument, check if they are valid and returns a list to be used in a call to PICS.

Usage

setParaEM(
  minK = 1,
  maxK = 15,
  tol = 1e-04,
  B = 100,
  mSelect = "BIC",
  mergePeaks = TRUE,
  mapCorrect = TRUE,
  dataType = NULL
)

Arguments

minK

An integer. The minimum number of binding events per region. If the value is 0, the minimum number is automatically calculated.

maxK

An integer.The maximum number of binding events per region. If the value is 0, the maximum number is automatically calculated.

tol

A numeric. The tolerance for the EM algorithm.

B

An integer. The maximum number of iterations to be used.

mSelect

A character. pecifying the information criteria to be used when selecting the number of binding events.

mergePeaks

A logical stating whether overlapping binding events should be picked.

mapCorrect

A logical stating whether mappability profiles should be incorporated in the estimation, i.e: missing reads estimated.

dataType

A character. If a dataType is set, the algorithm will use the default parameters for this type of data (all the previous arguments will be ignored).

Value

A list of parameters to be used in PICS.

See Also

PICS setParaPrior


List of parameters that can be used as an argument of PICS.

Description

This function takes from 0 to 6 parameters as argument, check if they are valid and returns a list to be used in a call to PICS.

Usage

setParaPrior(
  xi = 200,
  rho = 1,
  alpha = 20,
  beta = 40000,
  lambda = 0,
  dMu = 0,
  dataType = NULL,
  PExi = 0
)

Arguments

xi

An integer. The average DNA fragment size.

rho

An integer. A variance parameter for the average DNA fragment size distribution.

alpha

An integer. First hyperparameter of the inverse Gamma distribution for $sigma^2$ in the PICS model.

beta

An integer. Second hyperparameter of the inverse Gamma distribution for $sigma^2$ in the PICS model.

lambda

An integer. The precision of the prior for mu used for histone data.

dMu

An integer. Our best guess for the distance between two neighboring nucleosomes.

dataType

A character string. If a valid dataType is specified, use our suggested parameters. "MNase" or "sonicated"

PExi

A numeric. With paired end data, ‘xi’ can be calculated directly from the reads. If PExi is set, it will overwrite the xi determined by the dataType.

Value

A list of 6 parameters to be used in PICS.

See Also

setParaEM PICS

Examples

# set prior for PICS data
paraPrior <- setParaPrior()
# set prior for sonicated data using our selected default parameters
paraPrior <- setParaPrior(dataType="sonicated")

Summarize segmentList objects

Description

Summarize segmentList objects into a data.frame

Usage

summarySeg(seg)

Arguments

seg

A segmentList object as returned by segmentPICS.

Value

A data.frame. With

  • chr: Chromosome id

  • NF: Number of forward reads

  • NR: Number of reverse reads

  • L: Length of segment

  • min: Start location of segments

  • max: End location of segments