make_reanno: Makes annotation from R imported reannotation mapping

View source: R/make_reanno.R

make_reannoR Documentation

Makes annotation from R imported reannotation mapping

Description

make_reanno makes a reannotation (reanno) object.

Usage

make_reanno(reanno_path, PAC, mis_fasta_check = FALSE, output = "S4")

Arguments

reanno_path

Path to a directory where reannotation .Rdata files can be found.

PAC

PAC-list object containing an Anno data.frame with sequences as row names.

mis_fasta_check

Logical TRUE/FALSE if checking against anno_misX.fa should be done. The anno_misX.fa is the fasta file generated by the map_reanno function after completing the last mismatch cycle. This file contains all PAC sequences that failed to receive an alignment in any of fasta references and after all mismatch cycles. This file should be availble in the same folder as the output files generated in the map_reanno function (=reanno_path).

output

Character indicating output format. If type="S4" (default), then the reanno object is returned as an S4 object. If type="list", the reanno object will be returned as a list of tables (tibbles).

Details

Given the path to the R reannotation files (reanno_mis0/1/2/3.Rdata) generated by map_reanno, this function will summarize the reannotation files into one output that matches the order of sequences in a PAC object.

Value

Contains two slots/objects: 1. Overview: A table with a summary of the mapping. 2. Full_anno: Lists of tables holding the full imports per mismatch cycle (mis0, mis1 etc) and reference (mi0-ref1, mis0-ref2, mis1-ref1, mis1-ref2 etc). All table rows are ordered according to the order of the PAC originally used for the mapping. If mis_fasta_check is specified, the function will look for a anno_misX.fa (X = the file with the highest number) previously generated by the reannotation workflow. This file is used to double check so that no sequences are missing before and after reannotation.

See Also

http://bowtie-bio.sourceforge.net/index.shtml for information about Bowtie https://github.com/Danis102 for updates on the current package.

Other PAC reannotation: add_reanno(), as.reanno(), import_reanno(), make_conv(), map_reanno(), simplify_reanno()

Examples


######################################################### 
##### Simple example for reference mapping 
##### Please see manual for ?simply_reanno an ?add_reanno for more advanced 
##### mapping
closeAllConnections()
######################################################### 
##### Create an reanno object

### First, if you haven't already generated Bowtie indexes for the included
# fasta references you need to do so. If you are having problem see the small
# RNA guide (vignette) for more info.
                                                             
 ## tRNA:
 trna_file <- system.file("extdata/trna", "tRNA.fa", 
                          package = "seqpac", mustWork = TRUE)
 trna_dir <- dirname(trna_file)
 
 if(!sum(stringr::str_count(list.files(trna_dir), ".ebwt")) ==6){     
     Rbowtie::bowtie_build(trna_file, 
                           outdir=trna_dir, 
                           prefix="tRNA", force=TRUE)
                           }
 ## rRNA:
 rrna_file <- system.file("extdata/rrna", "rRNA.fa", 
                          package = "seqpac", mustWork = TRUE)
 rrna_dir <- dirname(rrna_file)
 
 if(!sum(stringr::str_count(list.files(rrna_dir), ".ebwt")) ==6){
     Rbowtie::bowtie_build(rrna_file, 
                           outdir=rrna_dir, 
                           prefix="rRNA", force=TRUE)
                           }
                           
                           
##  Then load a PAC-object and remove previous mapping from anno:
 load(system.file("extdata", "drosophila_sRNA_pac_filt_anno.Rdata", 
                   package = "seqpac", mustWork = TRUE))
 anno(pac) <- anno(pac)[,1, drop = FALSE]
 
 ref_paths <- list(trna= trna_file, rrna= rrna_file)

                                    
##  You may add an output path of your choice, but here we use a temp folder:
 output <- paste0(tempdir(),"/seqpac/test")


##  Then map the PAC-object against the fasta references. 
# Warning: if you use your own data, you may want to use override=FALSE, to avoid
# deleting previous mapping by mistake.

map_reanno(pac, ref_paths=ref_paths, output_path=output,
               type="internal", mismatches=0,  import="biotype", 
               threads=2, keep_temp=FALSE, override=TRUE)
 
##  Then import and generate a reanno-object of the temporary bowtie-files
reanno_biotype <- make_reanno(output, PAC=pac, mis_fasta_check = TRUE)
 
                                    
## Now make some search terms against reference names to create shorter names
# Theses can be used to create factors in downstream analysis
# One search hit (regular expressions) gives one new short name 

bio_search <- list(
              rrna=c("5S", "5.8S", "12S", "16S", "18S", "28S", "pre_S45"),
              trna =c("_tRNA", "mt:tRNA"))
 
 
## You can merge directly with your PAC-object by adding your original 
# PAC-object, that you used with map_reanno, to merge_pac option.

pac <- add_reanno(reanno_biotype, bio_search=bio_search, 
                       type="biotype", bio_perfect=FALSE, 
                       mismatches = 0, merge_pac=pac)

table(anno(pac)$mis0_bio) 

############################################################################ 
## Similarly, you can use Bowtie indexed genome fasta references
## But don't forget to use import="genome" for coordinate import
#    
# ref_paths <- list(genome="<path_to_bowtie_indexed_fasta>")
# output_genome <- "<your_path_to_output_folder>"

## We can actually map against several genomes simultaneously:
# (to illustrate the principle we run the mycoplasma genome twice)

 mycoplasma_file <- system.file("extdata/mycoplasma_genome",
                                "mycoplasma.fa",
                                package = "seqpac", mustWork = TRUE)
 mycoplasma_dir<- gsub("mycoplasma.fa", "", mycoplasma_file)
 
## Don't forget to create bowtie indexes
 if(!sum(stringr::str_count(list.files(mycoplasma_dir), ".ebwt")) ==6){
     Rbowtie::bowtie_build(mycoplasma_file, 
                           outdir=mycoplasma_dir, 
                           prefix="mycoplasma", force=TRUE)
                           }
                           
 ref_paths <- list(genome1=mycoplasma_file, genome2=mycoplasma_file)
 map_reanno(PAC=pac, ref_paths=ref_paths, output_path=output, 
            type="internal", mismatches=0, import="genome", 
            threads=2, keep_temp=TRUE, override=TRUE)
 reanno_genome <- make_reanno(output, PAC=pac, mis_fasta_check = TRUE)
 
 table(overview(reanno_genome)$Any) # Low Mycoplasma contamination
 
############################################################################


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