Title: | Generate synthetic nucleosome maps |
---|---|
Description: | This package can generate a synthetic map with reads covering the nucleosome regions as well as a synthetic map with forward and reverse reads emulating next-generation sequencing. The synthetic hybridization data of “Tiling Arrays” can also be generated. The user has choice between three different distributions for the read positioning: Normal, Student and Uniform. In addition, a visualization tool is provided to explore the synthetic nucleosome maps. |
Authors: | Rawane Samb [aut], Astrid Deschênes [cre, aut] , Pascal Belleau [aut] , Arnaud Droit [aut] |
Maintainer: | Astrid Deschênes <[email protected]> |
License: | Artistic-2.0 |
Version: | 1.35.0 |
Built: | 2024-10-31 06:07:17 UTC |
Source: | https://github.com/bioc/nucleoSim |
This package can generate a synthetic map with reads covering the nucleosome regions as well as a synthetic map with forward and reverse reads emulating next-generation sequencing. The user has choice between three different distributions for the read positioning: Normal, Student and Uniform.
Rawane Samb, Astrid Deschenes, Pascal Belleau and Arnaud Droit
Maintainer: Astrid Deschenes <[email protected]>
syntheticNucMapFromDist
to generate a synthetic
nucleosome map
Generate a graph for
a list marked as an syntheticNucMap
class
## S3 method for class 'syntheticNucMap' plot(x, ...)
## S3 method for class 'syntheticNucMap' plot(x, ...)
x |
a list marked as an |
... |
|
a graph of a synthetic nucleosome map
Astrid Deschenes, Rawane Samb
## Generate a synthetic map with 20 well-positioned nucleosomes, 5 fuzzy ## nucleosomes and 10 deleted nucleosomes using a Student distribution ## with a variance of 10 for the well-positioned nucleosomes, ## a variance of 20 for the fuzzy nucleosomes and a seed of 15 syntheticMap <- syntheticNucMapFromDist(wp.num = 30, wp.del = 10, wp.var = 10, fuz.num = 5, fuz.var = 20, rnd.seed = 15, distr = "Student") ## Create graph using the synthetic map plot(syntheticMap, xlab="Position", ylab="Coverage")
## Generate a synthetic map with 20 well-positioned nucleosomes, 5 fuzzy ## nucleosomes and 10 deleted nucleosomes using a Student distribution ## with a variance of 10 for the well-positioned nucleosomes, ## a variance of 20 for the fuzzy nucleosomes and a seed of 15 syntheticMap <- syntheticNucMapFromDist(wp.num = 30, wp.del = 10, wp.var = 10, fuz.num = 5, fuz.var = 20, rnd.seed = 15, distr = "Student") ## Create graph using the synthetic map plot(syntheticMap, xlab="Position", ylab="Coverage")
Generate a graph for
a list marked as an syntheticNucReads
class
## S3 method for class 'syntheticNucReads' plot(x, ...)
## S3 method for class 'syntheticNucReads' plot(x, ...)
x |
a list marked as a |
... |
|
a graph of a synthetic nucleosome map containing forward and reverse reads
Astrid Deschenes, Rawane Samb
## Generate a synthetic map with 20 well-positioned nucleosomes, 5 fuzzy ## nucleosomes and 10 deleted nucleosomes using a Student distribution ## with a variance of 10 for the well-positioned nucleosomes, ## a variance of 20 for the fuzzy nucleosomes and a seed of 15 syntheticNucSample <- syntheticNucReadsFromDist(wp.num = 30, wp.del = 10, wp.var = 10, fuz.num = 5, fuz.var = 20, rnd.seed = 15, distr = "Student", offset = 1000) ## Create graph using the synthetic map plot(syntheticNucSample, xlab="Position", ylab="Coverage")
## Generate a synthetic map with 20 well-positioned nucleosomes, 5 fuzzy ## nucleosomes and 10 deleted nucleosomes using a Student distribution ## with a variance of 10 for the well-positioned nucleosomes, ## a variance of 20 for the fuzzy nucleosomes and a seed of 15 syntheticNucSample <- syntheticNucReadsFromDist(wp.num = 30, wp.del = 10, wp.var = 10, fuz.num = 5, fuz.var = 20, rnd.seed = 15, distr = "Student", offset = 1000) ## Create graph using the synthetic map plot(syntheticNucSample, xlab="Position", ylab="Coverage")
Generate a synthetic nucleosome map, a map with complete sequences covering the nucleosome regions, using the distribution selected by the user. The distribution is used to assign the start position to the sequences associated with the nucleosomes. The user has choice between three different distributions: Normal, Student and Uniform.
The synthetic nucleosome map creation is separated into 3 steps :
1. Adding well-positioned nucleosomes following specified parameters. The nucleosomes are all positioned at equidistance. Assigning sequences of variable length to each nucleosome using a normal distribution and specified variance.
2. Deleting some well-positioned nucleosomes following specified parameters. Each nucleosome has an equal probability to be selected.
3. Adding fuzzy nucleosomes following an uniform distribution ad specified parameters. Assigning sequences of variable length to each nucleosome using the specified distribution and parameters. The sequence length is always following a normal distribution.
This function is a modified version of the syntheticNucMap() function from Bioconductor nucleR package (Flores and Orozco, 2011).
syntheticNucMapFromDist( wp.num, wp.del, wp.var, fuz.num, fuz.var, max.cover = 100, nuc.len = 147, len.var = 10, lin.len = 20, rnd.seed = NULL, as.ratio = FALSE, distr = c("Uniform", "Normal", "Student") )
syntheticNucMapFromDist( wp.num, wp.del, wp.var, fuz.num, fuz.var, max.cover = 100, nuc.len = 147, len.var = 10, lin.len = 20, rnd.seed = NULL, as.ratio = FALSE, distr = c("Uniform", "Normal", "Student") )
wp.num |
a non-negative |
wp.del |
a non-negative |
wp.var |
a non-negative |
fuz.num |
a non-negative |
fuz.var |
a non-negative |
max.cover |
a positive |
nuc.len |
a non-negative |
len.var |
a non-negative |
lin.len |
a non-negative |
rnd.seed |
a single value, interpreted as an |
as.ratio |
a |
distr |
the name of the distribution used to generate the nucleosome
map. The choices are : |
an list
of class
"syntheticNucMap" containing the
following elements:
call
the matched call.
wp.starts
a vector
of integer
, the start
positions of all well-positioned nucleosome regions. The central
position of the nucleosome is calculated as wp.starts + round(nuc.len/2).
wp.nreads
a vector
of integer
, the number of
sequences associated to each well-positioned nucleosome.
wp.reads
a IRanges
containing the well-positioned
nucleosome sequences.
fuz.starts
a vector
of integer
, the
start position of all the fuzzy nucleosomes.
fuz.nreads
a vector
of integer
, the number
of sequences associated to each fuzzy nucleosome.
fuz.reads
a IRanges
containing the fuzzy nucleosome
sequences.
syn.reads
a IRanges
containing all the synthetic
nucleosome sequences (from both fuzzy and well-positioned nucleosomes).
nuc.len
a numeric
the nucleosome length.
The following elements will be only returned if as.ratio=TRUE
:
ctr.reads
a IRanges
containing the naked DNA
(control) sequences.
syn.ratio
a Rle
containing the calculated ratio
between the nucleosome coverage and the control coverage.
Rawane Samb, Astrid Deschenes
## Generate a synthetic map with 20 well-positioned nucleosomes and 10 fuzzy ## nucleosomes using a Normal distribution with a variance of 30 for the ## well-positioned nucleosomes, a variance of 40 for the fuzzy nucleosomes ## and a seed of 15. syntheticNucMapFromDist(wp.num = 20, wp.del = 0, wp.var = 30, fuz.num = 10, fuz.var = 40, rnd.seed = 15, distr = "Normal") ## Same output but with ratio syntheticNucMapFromDist(wp.num = 20, wp.del = 0, wp.var = 30, fuz.num = 10, fuz.var = 40, rnd.seed = 15, as.ratio = TRUE, distr = "Normal")
## Generate a synthetic map with 20 well-positioned nucleosomes and 10 fuzzy ## nucleosomes using a Normal distribution with a variance of 30 for the ## well-positioned nucleosomes, a variance of 40 for the fuzzy nucleosomes ## and a seed of 15. syntheticNucMapFromDist(wp.num = 20, wp.del = 0, wp.var = 30, fuz.num = 10, fuz.var = 40, rnd.seed = 15, distr = "Normal") ## Same output but with ratio syntheticNucMapFromDist(wp.num = 20, wp.del = 0, wp.var = 30, fuz.num = 10, fuz.var = 40, rnd.seed = 15, as.ratio = TRUE, distr = "Normal")
Generate a synthetic nucleosome map, a map with forward and reverses reads (paired-end reads) covering the nucleosome regions, using the distribution selected by the user. The distribution is used to assign the start position to the forward reads associated with the nucleosomes. The user has choice between three different distributions: Normal, Student and Uniform. The final map is composed of paired-end reads.
#' The synthetic nucleosome map creation is separated into 3 steps :
1. Adding well-positioned nucleosomes following specified parameters. The nucleosomes are all positioned at equidistance. Assigning the starting positions of forward reads using the specified distribution and parameters. The distance between starting positions of paired-end reads is assigned using a normal distribution and specified variance.
2. Deleting some well-positioned nucleosomes following specified parameters. Each nucleosome has an equal probability to be selected.
3. Adding fuzzy nucleosomes following an uniform distribution and specified parameters. Assigning the starting positions of forward reads using the specified distribution and parameters. The distance between starting positions of paired-end reads is assigned using a normal distribution and specified variance.
This function has been largely inspired by the Generating synthetic maps section of the nucleR package (Flores et Orozco, 2011).
syntheticNucReadsFromDist( wp.num, wp.del, wp.var, fuz.num, fuz.var, max.cover = 100, nuc.len = 147, len.var = 10, lin.len = 20, read.len = 40, rnd.seed = NULL, distr = c("Uniform", "Normal", "Student"), offset )
syntheticNucReadsFromDist( wp.num, wp.del, wp.var, fuz.num, fuz.var, max.cover = 100, nuc.len = 147, len.var = 10, lin.len = 20, read.len = 40, rnd.seed = NULL, distr = c("Uniform", "Normal", "Student"), offset )
wp.num |
a non-negative |
wp.del |
a non-negative |
wp.var |
a non-negative |
fuz.num |
a non-negative |
fuz.var |
a non-negative |
max.cover |
a positive |
nuc.len |
a positive |
len.var |
a positive |
lin.len |
a non-negative |
read.len |
a positive |
rnd.seed |
a single value, interpreted as an |
distr |
the name of the distribution used to generate the nucleosome
map. The choices are : |
offset |
a non-negative |
an list
of class
"syntheticNucReads" containing the
following elements:
call
the matched call.
dataIP
a data.frame
with the chromosome name, the
starting and ending positions and the direction of all forward
and reverse reads for all well-positioned and fuzzy nucleosomes.
Paired-end reads are identified with an unique id.
wp
a data.frame
with the positions of all the
well-positioned nucleosomes, as well as the number of paired-reads
associated to each one.
fuz
a data.frame
with the positions of all the fuzzy
nucleosomes, as well as the number of paired-reads associated to each one.
paired
a data.frame
with the starting and ending
positions of the reads used to generate the paired-end reads. Paired-end
reads are identified with an unique id.
Pascal Belleau, Rawane Samb, Astrid Deschenes
## Generate a synthetic map with 20 well-positioned + 10 fuzzy nucleosomes ## using a Normal distribution with a variance of 30 for the well-positioned ## nucleosomes, a variance of 40 for the fuzzy nucleosomes and a seed of 15. ## Because of the fixed seed, each time is going to be run, the results ## are going to be the seed. res <- syntheticNucReadsFromDist(wp.num = 20, wp.del = 0, wp.var = 30, fuz.num = 10, fuz.var = 40, rnd.seed = 15, distr = "Normal", offset = 1000)
## Generate a synthetic map with 20 well-positioned + 10 fuzzy nucleosomes ## using a Normal distribution with a variance of 30 for the well-positioned ## nucleosomes, a variance of 40 for the fuzzy nucleosomes and a seed of 15. ## Because of the fixed seed, each time is going to be run, the results ## are going to be the seed. res <- syntheticNucReadsFromDist(wp.num = 20, wp.del = 0, wp.var = 30, fuz.num = 10, fuz.var = 40, rnd.seed = 15, distr = "Normal", offset = 1000)
Generate a synthetic nucleosome map using a synthetic nucleosome map.
This function is using a modified version of the syntheticNucMap() function from Bioconductor nucleR package (Flores and Orozco, 2011).
syntheticNucReadsFromMap(syntheticNucMap, read.len = 40, offset)
syntheticNucReadsFromMap(syntheticNucMap, read.len = 40, offset)
syntheticNucMap |
a |
read.len |
a positive |
offset |
a non-negatvie |
a list
of class
"syntheticNucReads" containing the
following elements:
call
the matched call.
dataIP
a data.frame
with the chromosome name, the
starting and ending positions and the direction of all forward
and reverse reads for all well-positioned and fuzzy nucleosomes.
Paired-end reads are identified with an unique id.
wp
a data.frame
with the positions of all the
well-positioned nucleosomes, as well as the number of paired-reads
associated to each one.
fuz
a data.frame
with the positions of all the fuzzy
nucleosomes, as well as the number of paired-reads associated to each one.
paired
a data.frame
with the starting and ending
positions of the reads used to generate the paired-end reads. Paired-end
reads are identified with an unique id.
Pascal Belleau, Rawane Samb, Astrid Deschenes
## Generate a synthetic map with 20 well-positioned + 10 fuzzy nucleosomes ## using a Normal distribution with a variance of 30 for the well-positioned ## nucleosomes, a variance of 40 for the fuzzy nucleosomes and a seed of 15 ## Because of the fixed seed, each time is going to be run, the results ## are going to be the seed syntheticMap <- syntheticNucMapFromDist(wp.num = 20, wp.del = 0, wp.var = 30, fuz.num = 10, fuz.var = 40, rnd.seed = 335, as.ratio = FALSE, distr = "Uniform") res <- nucleoSim:::syntheticNucReadsFromMap(syntheticMap, read.len = 45, offset = 1000)
## Generate a synthetic map with 20 well-positioned + 10 fuzzy nucleosomes ## using a Normal distribution with a variance of 30 for the well-positioned ## nucleosomes, a variance of 40 for the fuzzy nucleosomes and a seed of 15 ## Because of the fixed seed, each time is going to be run, the results ## are going to be the seed syntheticMap <- syntheticNucMapFromDist(wp.num = 20, wp.del = 0, wp.var = 30, fuz.num = 10, fuz.var = 40, rnd.seed = 335, as.ratio = FALSE, distr = "Uniform") res <- nucleoSim:::syntheticNucReadsFromMap(syntheticMap, read.len = 45, offset = 1000)