View source: R/runCallsPipeline.R
generate_calls_workflow | R Documentation |
Main function running the workflow that generates present/absent calls from a file, a data.frame, or objects of the classe UserMetadata (please choose only 1 out of the 3). This workflow is highly tunable by editing default values of the slots of S4 objects. For more information on how to tune the workflow please have a look at the vignette and the documentation of the classes KallistoMetadata, AbundanceMetadata, UserMetadata and BgeeMetadata
generate_calls_workflow(
abundanceMetadata = new("KallistoMetadata"),
bgeeMetadata = new("BgeeMetadata"),
userMetadata = NULL,
userDataFrame = NULL,
userFile = NULL,
checkTxVersion = FALSE
)
abundanceMetadata |
A Class AbundanceMetadata object (optional) allowing to tune your gene quantification abundance analyze |
bgeeMetadata |
A Class BgeeMetadata object (optional) allowing to choose the version of reference intergenic sequences |
userMetadata |
A Class UserMetadata object (optional). generate present/absent calls using slots of the UserMetadata class. |
userDataFrame |
a data.frame comtaining all information to generate present/absent calls. Each line of this data.frame will generate calls for one RNA-Seq library. This data.frame must contains between 4 and 8 columns :
|
userFile |
path to a tsv file containing between 4 and 8 columns. these columns are the same than for userDataFrame (see above). a template of this file is available at the root of the package and accessible with the command system.file('userMetadataTemplate.tsv', package = 'BgeeCall') |
checkTxVersion |
boolean used to define if BgeeCall check rather transcript version should be removed. Default value is FALSE |
paths to the 5 results files (see vignette for more details)
Julien Wollbrett
AbundanceMetadata, KallistoMetadata, BgeeMetadata, UserMetadata
## Not run:
# import gene annotation and transcriptome from AnnotationHub
library(AnnotationHub)
ah <- AnnotationHub()
ah_resources <- query(ah, c('Ensembl', 'Caenorhabditis elegans', '84'))
annotation_object <- ah_resources[['AH50789']]
transcriptome_object <- rtracklayer::import.2bit(ah_resources[['AH50453']])
# instanciate BgeeCall object
# add annotation and transcriptome in the user_BgeeCall object
# it is possible to import them using an S4 object (GRanges, DNAStringSet)
# or a file (gtf, fasta) with methods setAnnotationFromFile() and
# setTranscriptomeFromFile()
user_BgeeCall <- setAnnotationFromObject(user_BgeeCall,
annotation_object,
'WBcel235_84')
user_BgeeCall <- setTranscriptomeFromObject(user_BgeeCall,
transcriptome_object,
'WBcel235')
# provide path to the directory of your RNA-Seq library
user_BgeeCall <- setRNASeqLibPath(user_BgeeCall,
system.file('extdata', 'SRX099901_subset',
package = 'BgeeCall'))
# run the full BgeeCall workflow
calls_output <- generate_calls_workflow(
userMetadata = user_BgeeCall)
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.