PAC_covplot: Plot sequence coverage over a reference

View source: R/PAC_covplot.R

PAC_covplotR Documentation

Plot sequence coverage over a reference

Description

PAC_covplot Plotting sequences in a PAC object using an PAC mapping object.

Usage

PAC_covplot(
  PAC,
  map,
  summary_target = NULL,
  map_target = NULL,
  style = "line",
  xseq = TRUE,
  colors = NULL,
  check_overide = FALSE
)

Arguments

PAC

PAC-list object.

map

PAC_map object generated by PAC_mapper.

summary_target

List with 2 character vectors: 1st object: Character indicating the name of the target PAC summary object which should be used as input data for the plots. If left empty the 1st summary object will be used. 2nd object: Character vector indicating the names of the columns in the summary object from which individual graphs should be plotted. If left empty all summary columns are plotted in the same graph. Summaries are generated by the PAC_summary function or custom generated by the user and stored in the summary 'folder' of a PAC object.

map_target

(optional) Character vector. Important: This is similar to an anno_target, but instead extracts target references in the PAC_mapping object. Should contain search terms that can find unique strings in the reference names. The search terms are parsed to grepl("<search terms>, names(<PAC_mapper object>)). (default=NULL)

style

Character string. If style="line" (default) then graphs are plotted as simple path/line plots (ggplot2::geom_path). If style="solid" then the area under the line is filled with same color as the lines but with added transparency (ggplot2::geom_line + geom_ribbon, alpha=0.5).

xseq

Logical indicating whether the nucleotides should be plotted on the x-axis. Plotting long references with xseq=FALSE will increase script performance. (default=TRUE).

colors

Character vector indicating the rgb colors to be parsed to ggplot2 for plotting the coverage lines (default = c("black", "red", "grey", "blue"))

check_overide

Logical, whether all sequences names in PAC must be found in the PAC_map object (TRUE) or not (FALSE). If the PAC_map object was generated from the provided PAC objects, then all sequences should match. If the PAC list object was subsetted or the PAC_map is used with an independent PAC list object in which not all sequences are available, then this function will create an error, unless check_overide=TRUE. In such case, all missing sequences in the PAC object will automatically be assigned a 0-value in all summary_target columns.

Details

Given a PAC object and a map object generated by PAC_mapper this function will attempt to plot the sequence coverage over long references.

Value

Plot list with plots generated by ggplot2

See Also

https://github.com/Danis102 for updates on the current package.

Other PAC analysis: PAC_deseq(), PAC_filter(), PAC_filtsep(), PAC_gtf(), PAC_jitter(), PAC_mapper(), PAC_nbias(), PAC_norm(), PAC_pca(), PAC_pie(), PAC_saturation(), PAC_sizedist(), PAC_stackbar(), PAC_summary(), PAC_trna(), as.PAC(), filtsep_bin(), map_rangetype(), tRNA_class()

Examples


###########################################################
### Simple example of how to use PAC_mapper and PAC_covplot
# Note: More details, see vignette and manuals.)
# Also see: ?map_rangetype, ?tRNA_class or ?PAC_trna for more examples
# on how to use PAC_mapper.

## Load PAC-object, make summaries and extract rRNA and tRNA
load(system.file("extdata", "drosophila_sRNA_pac_filt_anno.Rdata", 
                   package = "seqpac", mustWork = TRUE))

pac <- PAC_summary(pac, norm = "cpm", type = "means", 
                   pheno_target=list("stage", unique(pheno(pac)$stage)))

pac_tRNA <- PAC_filter(pac, anno_target = list("Biotypes_mis0", "tRNA"))


## Give paths to a fasta reference (with or without bowtie index)
#  (Here we use an rRNA/tRNA fasta included in seqpac) 

ref_tRNA <- system.file("extdata/trna", "tRNA.fa", 
                         package = "seqpac", mustWork = TRUE)


## Map using PAC-mapper
      
map_tRNA <- PAC_mapper(pac_tRNA, mismatches=0, 
                        threads=1, ref=ref_tRNA, override=TRUE)

## Plot tRNA using xseq=TRUE gives you reference sequence as X-axis:
# (OBS! Long reference will not print well with xseq=TRUE)
cov_tRNA <- PAC_covplot(pac_tRNA, map_tRNA, 
                        summary_target = list("cpmMeans_stage"), 
                        xseq=TRUE, style="line", 
                        color=c("red", "black", "blue"))

cov_tRNA[[1]]

## Check which tRNA reached a decent number of unique fragments 
# (OBS! This is a very down sampled dataset)
logi_hi <- unlist(lapply(map_tRNA, function(x){nrow(x$Alignments) > 10 }))

table(logi_hi)  
names(map_tRNA[logi_hi])

targets <- c("Ala-AGC-1-1", "Lys-CTT-1-13","Ser-AGA-2-2")

cov_tRNA_sub <- PAC_covplot(pac_tRNA, map_tRNA, 
                        summary_target = list("cpmMeans_stage"),
                        map_target = targets,
                        xseq=TRUE, style="line", 
                        color=c("red", "black", "blue"))

cowplot::plot_grid(plotlist= cov_tRNA_sub) 


Danis102/seqpac documentation built on Aug. 26, 2023, 10:15 a.m.