ReUseData
provides functionalities to construct workflow-based data
recipes for fully tracked and reproducible data processing. Evaluation
of data recipes generates curated data resources in their generic
formats (e.g., VCF, bed), as well as a YAML manifest file recording
the recipe parameters, data annotations, and data file paths for
subsequent reuse. The datasets are locally cached using a database
infrastructure, where updating and searching of specific data is made easy.
The data reusability is assured through cloud hosting and enhanced interoperability with downstream software tools or analysis workflows. The workflow strategy enables cross platform reproducibility of curated data resources.
if (!requireNamespace("BiocManager", quietly = TRUE))
install.packages("BiocManager")
BiocManager::install("ReUseData")
Use the development version:
BiocManager::install("ReUseData", version = "devel")
suppressPackageStartupMessages(library(Rcwl))
library(ReUseData)
ReUseData
core functions for data managementHere we introduce the core functions of ReUseData
for data
management and reuse: getData
for reproducible data generation,
dataUpdate
for syncing and updating data cache, and dataSearch
for
multi-keywords searching of dataset of interest.
First, we can construct data recipes by transforming shell or other ad
hoc data preprocessing scripts into workflow-based data recipes. Some
prebuilt data recipes for public data resources (e.g., downloading,
unzipping and indexing) are available for direct use through
recipeSearch
and recipeLoad
functions. Then we will assign values
to the input parameters and evaluate the recipe to generate data of
interest.
## set cache in tempdir for test
Sys.setenv(cachePath = file.path(tempdir(), "cache"))
recipeUpdate()
#> Updating recipes...
#> STAR_index.R added
#> bowtie2_index.R added
#> echo_out.R added
#> ensembl_liftover.R added
#> gcp_broad_gatk_hg19.R added
#> gcp_broad_gatk_hg38.R added
#> gcp_gatk_mutect2_b37.R added
#> gcp_gatk_mutect2_hg38.R added
#> gencode_annotation.R added
#> gencode_genome_grch38.R added
#> gencode_transcripts.R added
#> hisat2_index.R added
#> reference_genome.R added
#> salmon_index.R added
#> ucsc_database.R added
#>
#> recipeHub with 15 records
#> cache path: /tmp/RtmptrHkdM/cache/ReUseDataRecipe
#> # recipeSearch() to query specific recipes using multipe keywords
#> # recipeUpdate() to update the local recipe cache
#>
#> name
#> BFC1 | STAR_index
#> BFC2 | bowtie2_index
#> BFC3 | echo_out
#> BFC4 | ensembl_liftover
#> BFC5 | gcp_broad_gatk_hg19
#> ... ...
#> BFC11 | gencode_transcripts
#> BFC12 | hisat2_index
#> BFC13 | reference_genome
#> BFC14 | salmon_index
#> BFC15 | ucsc_database
recipeSearch("echo")
#> recipeHub with 1 records
#> cache path: /tmp/RtmptrHkdM/cache/ReUseDataRecipe
#> # recipeSearch() to query specific recipes using multipe keywords
#> # recipeUpdate() to update the local recipe cache
#>
#> name
#> BFC3 | echo_out
echo_out <- recipeLoad("echo_out")
#> Note: you need to assign a name for the recipe: rcpName <- recipeLoad('xx')
#> Data recipe loaded!
#> Use inputs() to check required input parameters before evaluation.
#> Check here: https://rcwl.org/dataRecipes/echo_out.html
#> for user instructions (e.g., eligible input values, data source, etc.)
inputs(echo_out)
#> inputs:
#> input (input) (string):
#> outfile (outfile) (string):
Users can then assign values to the input parameters, and evaluate the
recipe (getData
) to generate data of interest. Users need to specify
an output directory for all files (desired data file, intermediate
files that are internally generated as workflow scripts or annotation
files). Detailed notes for the data is encouraged which will be used
for keywords matching for later data search.
We can install cwltool first to make sure a cwl-runner is available.
invisible(Rcwl::install_cwltool())
echo_out$input <- "Hello World!"
echo_out$outfile <- "outfile"
outdir <- file.path(tempdir(), "SharedData")
res <- getData(echo_out,
outdir = outdir,
notes = c("echo", "hello", "world", "txt"))
#> }[1;30mINFO[0m Final process status is success
The file path to newly generated dataset can be easily retrieved. It
can also be retrieved using dataSearch()
functions with multiple
keywords. Before that, dataUpdate()
needs to be done.
res$output
#> [1] "/tmp/RtmptrHkdM/SharedData/outfile.txt"
There are some automatically generated files to help track the data
recipe evaluation, including *.sh
to record the original shell
script, *.cwl
file as the official workflow script which was
internally submitted for data recipe evaluation, *.yml
file as part
of CWL workflow evaluation, which also record data annotations, and
*.md5
checksum file to check/verify the integrity of generated data
file.
list.files(outdir, pattern = "echo")
#> [1] "echo_out_Hello_World!_outfile.cwl" "echo_out_Hello_World!_outfile.md5"
#> [3] "echo_out_Hello_World!_outfile.sh" "echo_out_Hello_World!_outfile.yml"
The *.yml
file contains information about recipe input parameters,
the file path to output file, the notes for the dataset, and
auto-added date for data generation time. A later data search using
dataSearch()
will refer to this file for keywords match.
readLines(res$yml)
#> [1] "input: Hello World!"
#> [2] "outfile: outfile"
#> [3] "# output: /tmp/RtmptrHkdM/SharedData/outfile.txt"
#> [4] "# notes: echo hello world txt"
#> [5] "# date: 2024-10-29 19:07:12.162606"
dataUpdate()
creates (if first time use), syncs and updates the local
cache for curated datasets. It finds and reads all the .yml
files
recursively in the provided data folder, creates a cache record for
each dataset that is associated (including newly generated ones with
getData()
), and updates the local cache for later data searching and
reuse.
IMPORTANT: It is recommended that users create a specified folder for
data archival (e.g., file/path/to/SharedData
) that other group
members have access to, and use sub-folders for different kinds of
datasets (e.g., those generated from same recipe).
(dh <- dataUpdate(dir = outdir))
#>
#> Updating data record...
#> outfile.txt added
#> dataHub with 1 records
#> cache path: /tmp/RtmptrHkdM/cache/ReUseData
#> # dataUpdate() to update the local data cache
#> # dataSearch() to query a specific dataset
#> # Additional information can be retrieved using:
#> # dataNames(), dataParams(), dataNotes(), dataPaths(), dataTag() or mcols()
#>
#> name Path
#> BFC1 | outfile.txt /tmp/RtmptrHkdM/SharedData/outfile.txt
dataUpdate
and dataSearch
return a dataHub
object with a list of
all available or matching datasets.
One can subset the list with [
and use getter functions to retrieve
the annotation information about the data, e.g., data names,
parameters values to the recipe, notes, tags, and the corresponding
yaml file.
dh[1]
#> dataHub with 1 records
#> cache path: /tmp/RtmptrHkdM/cache/ReUseData
#> # dataUpdate() to update the local data cache
#> # dataSearch() to query a specific dataset
#> # Additional information can be retrieved using:
#> # dataNames(), dataParams(), dataNotes(), dataPaths(), dataTag() or mcols()
#>
#> name Path
#> BFC1 | outfile.txt /tmp/RtmptrHkdM/SharedData/outfile.txt
## dh["BFC1"]
dh[dataNames(dh) == "outfile.txt"]
#> dataHub with 1 records
#> cache path: /tmp/RtmptrHkdM/cache/ReUseData
#> # dataUpdate() to update the local data cache
#> # dataSearch() to query a specific dataset
#> # Additional information can be retrieved using:
#> # dataNames(), dataParams(), dataNotes(), dataPaths(), dataTag() or mcols()
#>
#> name Path
#> BFC1 | outfile.txt /tmp/RtmptrHkdM/SharedData/outfile.txt
dataNames(dh)
#> [1] "outfile.txt"
dataParams(dh)
#> [1] "input: Hello World!; outfile: outfile"
dataNotes(dh)
#> [1] "echo hello world txt"
dataTags(dh)
#> [1] ""
dataYml(dh)
#> [1] "/tmp/RtmptrHkdM/SharedData/echo_out_Hello_World!_outfile.yml"
ReUseData
, as the name suggests, commits to promoting the data reuse.
Data can be prepared in standard input formats (toList
), e.g.,
YAML and JSON, to be easily integrated in workflow methods that are
locally or cloud-hosted.
(dh1 <- dataSearch(c("echo", "hello", "world")))
#> dataHub with 1 records
#> cache path: /tmp/RtmptrHkdM/cache/ReUseData
#> # dataUpdate() to update the local data cache
#> # dataSearch() to query a specific dataset
#> # Additional information can be retrieved using:
#> # dataNames(), dataParams(), dataNotes(), dataPaths(), dataTag() or mcols()
#>
#> name Path
#> BFC1 | outfile.txt /tmp/RtmptrHkdM/SharedData/outfile.txt
toList(dh1, listNames = c("input_file"))
#> $input_file
#> [1] "/tmp/RtmptrHkdM/SharedData/outfile.txt"
toList(dh1, format = "yaml", listNames = c("input_file"))
#> [1] "input_file: /tmp/RtmptrHkdM/SharedData/outfile.txt"
toList(dh1, format = "json", file = file.path(tempdir(), "data.json"))
#> File is saved as: "/tmp/RtmptrHkdM/data.json"
#> {
#> "outfile.txt": "/tmp/RtmptrHkdM/SharedData/outfile.txt"
#> }
Data can also be aggregated from different resources by tagging with specific software tools.
dataSearch()
#> dataHub with 1 records
#> cache path: /tmp/RtmptrHkdM/cache/ReUseData
#> # dataUpdate() to update the local data cache
#> # dataSearch() to query a specific dataset
#> # Additional information can be retrieved using:
#> # dataNames(), dataParams(), dataNotes(), dataPaths(), dataTag() or mcols()
#>
#> name Path
#> BFC1 | outfile.txt /tmp/RtmptrHkdM/SharedData/outfile.txt
dataTags(dh[1]) <- "#gatk"
dataSearch("#gatk")
#> dataHub with 1 records
#> cache path: /tmp/RtmptrHkdM/cache/ReUseData
#> # dataUpdate() to update the local data cache
#> # dataSearch() to query a specific dataset
#> # Additional information can be retrieved using:
#> # dataNames(), dataParams(), dataNotes(), dataPaths(), dataTag() or mcols()
#>
#> name Path
#> BFC1 | outfile.txt /tmp/RtmptrHkdM/SharedData/outfile.txt
The package can also be used to add annotation and notes to existing data resources or experiment data for management. Here we add exisiting “exp_data” to local data repository.
exp_data <- file.path(tempdir(), "exp_data")
dir.create(exp_data)
We first add notes to the data, and then update data repository with information from the new dataset.
annData(exp_data, notes = c("experiment data"))
#> meta.yml added
#> [1] "/tmp/RtmptrHkdM/exp_data/meta.yml"
dataUpdate(exp_data)
#>
#> Updating data record...
#> exp_data added
#> dataHub with 2 records
#> cache path: /tmp/RtmptrHkdM/cache/ReUseData
#> # dataUpdate() to update the local data cache
#> # dataSearch() to query a specific dataset
#> # Additional information can be retrieved using:
#> # dataNames(), dataParams(), dataNotes(), dataPaths(), dataTag() or mcols()
#>
#> name Path
#> BFC1 | outfile.txt /tmp/RtmptrHkdM/SharedData/outfile.txt
#> BFC2 | exp_data /tmp/RtmptrHkdM/exp_data
Now our data hub cached meta information from two different directories, one from data recipe and one from exisiting data. Data can be retrieved by keywords.
dataSearch("experiment")
#> dataHub with 1 records
#> cache path: /tmp/RtmptrHkdM/cache/ReUseData
#> # dataUpdate() to update the local data cache
#> # dataSearch() to query a specific dataset
#> # Additional information can be retrieved using:
#> # dataNames(), dataParams(), dataNotes(), dataPaths(), dataTag() or mcols()
#>
#> name Path
#> BFC2 | exp_data /tmp/RtmptrHkdM/exp_data
NOTE: if the argument cloud=TRUE
is enabled, dataUpdate()
will
also cache the pregenerated data sets (from evaluation of public
ReUseData recipes) that are available on ReUseData google bucket and
return in the dataHub
object that are fully searchable. Please see
the following section for details.
With the prebuilt data recipes for curation (e.g., downloading, unzipping, indexing) of commonly used public data resources we have pregenerated some data sets and put them on the cloud space for direct use.
Before searching, one need to use dataUpdate(cloud=TRUE)
to sync the
existing data sets on cloud, then dataSearch()
can be used to search
any available data set either in local cache and on the cloud.
gcpdir <- file.path(tempdir(), "gcpData")
dataUpdate(gcpdir, cloud=TRUE)
#>
#> Updating data record...
#> 1814f367c2e284_GRCh38.primary_assembly.genome.fa.1.bt2 added
#> 1814f33b748b38_GRCh38.primary_assembly.genome.fa.2.bt2 added
#> 1814f335e14b54_GRCh38.primary_assembly.genome.fa.3.bt2 added
#> 1814f31f589f34_GRCh38.primary_assembly.genome.fa.4.bt2 added
#> 1814f37ade4ef0_GRCh38.primary_assembly.genome.fa.rev.1.bt2 added
#> 1814f37ae46e5e_GRCh38.primary_assembly.genome.fa.rev.2.bt2 added
#> 1814f3778a717a_outfile.txt added
#> 1814f3288ffb62_GRCh37_to_GRCh38.chain added
#> 1814f3453c45ca_GRCh37_to_NCBI34.chain added
#> 1814f36701f2c8_GRCh37_to_NCBI35.chain added
#> 1814f354c82e4_GRCh37_to_NCBI36.chain added
#> 1814f328c51c7a_GRCh38_to_GRCh37.chain added
#> 1814f33936d4da_GRCh38_to_NCBI34.chain added
#> 1814f31b300b_GRCh38_to_NCBI35.chain added
#> 1814f32308c5ea_GRCh38_to_NCBI36.chain added
#> 1814f314ee1c95_NCBI34_to_GRCh37.chain added
#> 1814f3638a6d68_NCBI34_to_GRCh38.chain added
#> 1814f39c30db8_NCBI35_to_GRCh37.chain added
#> 1814f3edfe85a_NCBI35_to_GRCh38.chain added
#> 1814f368940ef2_NCBI36_to_GRCh37.chain added
#> 1814f34d0fdf13_NCBI36_to_GRCh38.chain added
#> 1814f350acf595_GRCm38_to_NCBIM36.chain added
#> 1814f34c7295ed_GRCm38_to_NCBIM37.chain added
#> 1814f31097d15c_NCBIM36_to_GRCm38.chain added
#> 1814f354b9efbf_NCBIM37_to_GRCm38.chain added
#> 1814f31988429b_1000G_omni2.5.b37.vcf.gz added
#> 1814f31e87ae2e_1000G_omni2.5.b37.vcf.gz.tbi added
#> 1814f324bc4aaf_Mills_and_1000G_gold_standard.indels.b37.vcf.gz added
#> 1814f362721cfa_Mills_and_1000G_gold_standard.indels.b37.vcf.gz.tbi added
#> 1814f369fa985f_1000G_omni2.5.hg38.vcf.gz added
#> 1814f34e602aa1_1000G_omni2.5.hg38.vcf.gz.tbi added
#> 1814f34a34ff7f_Mills_and_1000G_gold_standard.indels.hg38.vcf.gz added
#> 1814f3256f2397_Mills_and_1000G_gold_standard.indels.hg38.vcf.gz.tbi added
#> 1814f344175f5_af-only-gnomad.raw.sites.vcf added
#> 1814f3698d9eb3_af-only-gnomad.raw.sites.vcf.idx added
#> 1814f3204d7287_Mutect2-exome-panel.vcf.idx added
#> 1814f37f25e453_Mutect2-WGS-panel-b37.vcf added
#> 1814f36118102d_Mutect2-WGS-panel-b37.vcf.idx added
#> 1814f348dd6dea_small_exac_common_3.vcf added
#> 1814f344622a1d_small_exac_common_3.vcf.idx added
#> 1814f3481a02f5_1000g_pon.hg38.vcf.gz added
#> 1814f34e29f0ce_1000g_pon.hg38.vcf.gz.tbi added
#> 1814f36d274698_af-only-gnomad.hg38.vcf.gz added
#> 1814f3150d7d0_af-only-gnomad.hg38.vcf.gz.tbi added
#> 1814f34e4520d9_small_exac_common_3.hg38.vcf.gz added
#> 1814f310300c82_small_exac_common_3.hg38.vcf.gz.tbi added
#> 1814f3163ef465_gencode.v41.annotation.gtf added
#> 1814f331cf8e41_gencode.v42.annotation.gtf added
#> 1814f319f31a3a_gencode.vM30.annotation.gtf added
#> 1814f3251edcbf_gencode.vM31.annotation.gtf added
#> 1814f31a639d33_gencode.v41.transcripts.fa added
#> 1814f36702f94d_gencode.v41.transcripts.fa.fai added
#> 1814f375cbd254_gencode.v42.transcripts.fa added
#> 1814f366d63321_gencode.v42.transcripts.fa.fai added
#> 1814f3779acaa9_gencode.vM30.pc_transcripts.fa added
#> 1814f34a85c214_gencode.vM30.pc_transcripts.fa.fai added
#> 1814f35e75bc_gencode.vM31.pc_transcripts.fa added
#> 1814f3162278d8_gencode.vM31.pc_transcripts.fa.fai added
#> 1814f36f420cc3_GRCh38.primary_assembly.genome.fa.1.ht2 added
#> 1814f362d092b7_GRCh38.primary_assembly.genome.fa.2.ht2 added
#> 1814f31d1137_GRCh38.primary_assembly.genome.fa.3.ht2 added
#> 1814f33da23764_GRCh38.primary_assembly.genome.fa.4.ht2 added
#> 1814f32d059236_GRCh38.primary_assembly.genome.fa.5.ht2 added
#> 1814f3258c34cf_GRCh38.primary_assembly.genome.fa.6.ht2 added
#> 1814f341e3ad5a_GRCh38.primary_assembly.genome.fa.7.ht2 added
#> 1814f3169330e9_GRCh38.primary_assembly.genome.fa.8.ht2 added
#> 1814f345d9a756_GRCh38_full_analysis_set_plus_decoy_hla.fa.fai added
#> 1814f3410991ad_GRCh38_full_analysis_set_plus_decoy_hla.fa.amb added
#> 1814f377ab4116_GRCh38_full_analysis_set_plus_decoy_hla.fa.ann added
#> 1814f3eb71540_GRCh38_full_analysis_set_plus_decoy_hla.fa.bwt added
#> 1814f356bbbcb_GRCh38_full_analysis_set_plus_decoy_hla.fa.pac added
#> 1814f33fc5440b_GRCh38_full_analysis_set_plus_decoy_hla.fa.sa added
#> 1814f35ce1060e_GRCh38_full_analysis_set_plus_decoy_hla.fa added
#> 1814f372930263_GRCh38.primary_assembly.genome.fa.fai added
#> 1814f341161bdb_GRCh38.primary_assembly.genome.fa.amb added
#> 1814f32b2626e7_GRCh38.primary_assembly.genome.fa.ann added
#> 1814f32c30ee5_GRCh38.primary_assembly.genome.fa.bwt added
#> 1814f357551040_GRCh38.primary_assembly.genome.fa.pac added
#> 1814f35cf5b528_GRCh38.primary_assembly.genome.fa.sa added
#> 1814f31cb62920_GRCh38.primary_assembly.genome.fa added
#> 1814f37c73ecff_hs37d5.fa.fai added
#> 1814f37759525c_hs37d5.fa.amb added
#> 1814f33b9226d_hs37d5.fa.ann added
#> 1814f3723fbf54_hs37d5.fa.bwt added
#> 1814f35e2f857d_hs37d5.fa.pac added
#> 1814f37b53ed17_hs37d5.fa.sa added
#> 1814f33cc58168_hs37d5.fa added
#> 1814f35e8dfb39_complete_ref_lens.bin added
#> 1814f3117665ef_ctable.bin added
#> 1814f32c078e2b_ctg_offsets.bin added
#> 1814f3415e8df0_duplicate_clusters.tsv added
#> 1814f311937726_info.json added
#> 1814f369a9c590_mphf.bin added
#> 1814f36e642026_pos.bin added
#> 1814f3371fabf5_pre_indexing.log added
#> 1814f32b8d72ea_rank.bin added
#> 1814f34f7510f_ref_indexing.log added
#> 1814f37cf9534c_refAccumLengths.bin added
#> 1814f36c970497_reflengths.bin added
#> 1814f37ca29225_refseq.bin added
#> 1814f3bb0688c_seq.bin added
#> 1814f37202c062_versionInfo.json added
#> 1814f33c67d631_salmon_index added
#> 1814f368916e9b_chrLength.txt added
#> 1814f36495c2c5_chrName.txt added
#> 1814f37d7df20c_chrNameLength.txt added
#> 1814f313b79582_chrStart.txt added
#> 1814f36758d1ab_exonGeTrInfo.tab added
#> 1814f354d3024d_exonInfo.tab added
#> 1814f370ad4aab_geneInfo.tab added
#> 1814f340efacb_Genome added
#> 1814f35146ef4c_genomeParameters.txt added
#> 1814f368069d07_Log.out added
#> 1814f37c81d38_SA added
#> 1814f34386aea0_SAindex added
#> 1814f346362284_sjdbInfo.txt added
#> 1814f331c0a4f_sjdbList.fromGTF.out.tab added
#> 1814f34c3008_sjdbList.out.tab added
#> 1814f324c41dbd_transcriptInfo.tab added
#> 1814f31492703e_GRCh38.GENCODE.v42_100 added
#> 1814f32c53be34_knownGene_hg38.sql added
#> 1814f36622abae_knownGene_hg38.txt added
#> 1814f32625e765_refGene_hg38.sql added
#> 1814f315fd83c4_refGene_hg38.txt added
#> 1814f35486cbd4_knownGene_mm39.sql added
#> 1814f35d45935a_knownGene_mm39.txt added
#> 1814f3418af6ae_refGene_mm39.sql added
#> 1814f3597e1ce4_refGene_mm39.txt added
#> dataHub with 130 records
#> cache path: /tmp/RtmptrHkdM/cache/ReUseData
#> # dataUpdate() to update the local data cache
#> # dataSearch() to query a specific dataset
#> # Additional information can be retrieved using:
#> # dataNames(), dataParams(), dataNotes(), dataPaths(), dataTag() or mcols()
#>
#> name
#> BFC1 | outfile.txt
#> BFC2 | exp_data
#> BFC3 | GRCh38.primary_assembly.genome.fa.1.bt2
#> BFC4 | GRCh38.primary_assembly.genome.fa.2.bt2
#> BFC5 | GRCh38.primary_assembly.genome.fa.3.bt2
#> ... ...
#> BFC126 | refGene_hg38.txt
#> BFC127 | knownGene_mm39.sql
#> BFC128 | knownGene_mm39.txt
#> BFC129 | refGene_mm39.sql
#> BFC130 | refGene_mm39.txt
#> Path
#> BFC1 /tmp/RtmptrHkdM/SharedData/outfile.txt
#> BFC2 /tmp/RtmptrHkdM/exp_data
#> BFC3 https://storage.googleapis.com/reusedata/bowtie2_index/GRCh38.pr...
#> BFC4 https://storage.googleapis.com/reusedata/bowtie2_index/GRCh38.pr...
#> BFC5 https://storage.googleapis.com/reusedata/bowtie2_index/GRCh38.pr...
#> ... ...
#> BFC126 https://storage.googleapis.com/reusedata/ucsc_database/refGene_h...
#> BFC127 https://storage.googleapis.com/reusedata/ucsc_database/knownGene...
#> BFC128 https://storage.googleapis.com/reusedata/ucsc_database/knownGene...
#> BFC129 https://storage.googleapis.com/reusedata/ucsc_database/refGene_m...
#> BFC130 https://storage.googleapis.com/reusedata/ucsc_database/refGene_m...
If the data of interest already exist on the cloud, then
getCloudData
will directly download the data to your computer. Add
it to the local caching system using dataUpdate()
for later use.
(dh <- dataSearch(c("ensembl", "GRCh38")))
#> dataHub with 8 records
#> cache path: /tmp/RtmptrHkdM/cache/ReUseData
#> # dataUpdate() to update the local data cache
#> # dataSearch() to query a specific dataset
#> # Additional information can be retrieved using:
#> # dataNames(), dataParams(), dataNotes(), dataPaths(), dataTag() or mcols()
#>
#> name
#> BFC10 | GRCh37_to_GRCh38.chain
#> BFC14 | GRCh38_to_GRCh37.chain
#> BFC15 | GRCh38_to_NCBI34.chain
#> BFC16 | GRCh38_to_NCBI35.chain
#> BFC17 | GRCh38_to_NCBI36.chain
#> BFC19 | NCBI34_to_GRCh38.chain
#> BFC21 | NCBI35_to_GRCh38.chain
#> BFC23 | NCBI36_to_GRCh38.chain
#> Path
#> BFC10 https://storage.googleapis.com/reusedata/ensembl_liftover/GRCh37...
#> BFC14 https://storage.googleapis.com/reusedata/ensembl_liftover/GRCh38...
#> BFC15 https://storage.googleapis.com/reusedata/ensembl_liftover/GRCh38...
#> BFC16 https://storage.googleapis.com/reusedata/ensembl_liftover/GRCh38...
#> BFC17 https://storage.googleapis.com/reusedata/ensembl_liftover/GRCh38...
#> BFC19 https://storage.googleapis.com/reusedata/ensembl_liftover/NCBI34...
#> BFC21 https://storage.googleapis.com/reusedata/ensembl_liftover/NCBI35...
#> BFC23 https://storage.googleapis.com/reusedata/ensembl_liftover/NCBI36...
getCloudData(dh[1], outdir = gcpdir)
#> Data is downloaded:
#> /tmp/RtmptrHkdM/gcpData/GRCh37_to_GRCh38.chain
Now we create the data cache with only local data files, and we can see that the downloaded data is available.
dataUpdate(gcpdir) ## Update local data cache (without cloud data)
#>
#> Updating data record...
#> GRCh37_to_GRCh38.chain added
#> dataHub with 131 records
#> cache path: /tmp/RtmptrHkdM/cache/ReUseData
#> # dataUpdate() to update the local data cache
#> # dataSearch() to query a specific dataset
#> # Additional information can be retrieved using:
#> # dataNames(), dataParams(), dataNotes(), dataPaths(), dataTag() or mcols()
#>
#> name
#> BFC1 | outfile.txt
#> BFC2 | exp_data
#> BFC3 | GRCh38.primary_assembly.genome.fa.1.bt2
#> BFC4 | GRCh38.primary_assembly.genome.fa.2.bt2
#> BFC5 | GRCh38.primary_assembly.genome.fa.3.bt2
#> ... ...
#> BFC127 | knownGene_mm39.sql
#> BFC128 | knownGene_mm39.txt
#> BFC129 | refGene_mm39.sql
#> BFC130 | refGene_mm39.txt
#> BFC131 | GRCh37_to_GRCh38.chain
#> Path
#> BFC1 /tmp/RtmptrHkdM/SharedData/outfile.txt
#> BFC2 /tmp/RtmptrHkdM/exp_data
#> BFC3 https://storage.googleapis.com/reusedata/bowtie2_index/GRCh38.pr...
#> BFC4 https://storage.googleapis.com/reusedata/bowtie2_index/GRCh38.pr...
#> BFC5 https://storage.googleapis.com/reusedata/bowtie2_index/GRCh38.pr...
#> ... ...
#> BFC127 https://storage.googleapis.com/reusedata/ucsc_database/knownGene...
#> BFC128 https://storage.googleapis.com/reusedata/ucsc_database/knownGene...
#> BFC129 https://storage.googleapis.com/reusedata/ucsc_database/refGene_m...
#> BFC130 https://storage.googleapis.com/reusedata/ucsc_database/refGene_m...
#> BFC131 /tmp/RtmptrHkdM/gcpData/GRCh37_to_GRCh38.chain
dataSearch() ## data is available locally!!!
#> dataHub with 131 records
#> cache path: /tmp/RtmptrHkdM/cache/ReUseData
#> # dataUpdate() to update the local data cache
#> # dataSearch() to query a specific dataset
#> # Additional information can be retrieved using:
#> # dataNames(), dataParams(), dataNotes(), dataPaths(), dataTag() or mcols()
#>
#> name
#> BFC1 | outfile.txt
#> BFC2 | exp_data
#> BFC3 | GRCh38.primary_assembly.genome.fa.1.bt2
#> BFC4 | GRCh38.primary_assembly.genome.fa.2.bt2
#> BFC5 | GRCh38.primary_assembly.genome.fa.3.bt2
#> ... ...
#> BFC127 | knownGene_mm39.sql
#> BFC128 | knownGene_mm39.txt
#> BFC129 | refGene_mm39.sql
#> BFC130 | refGene_mm39.txt
#> BFC131 | GRCh37_to_GRCh38.chain
#> Path
#> BFC1 /tmp/RtmptrHkdM/SharedData/outfile.txt
#> BFC2 /tmp/RtmptrHkdM/exp_data
#> BFC3 https://storage.googleapis.com/reusedata/bowtie2_index/GRCh38.pr...
#> BFC4 https://storage.googleapis.com/reusedata/bowtie2_index/GRCh38.pr...
#> BFC5 https://storage.googleapis.com/reusedata/bowtie2_index/GRCh38.pr...
#> ... ...
#> BFC127 https://storage.googleapis.com/reusedata/ucsc_database/knownGene...
#> BFC128 https://storage.googleapis.com/reusedata/ucsc_database/knownGene...
#> BFC129 https://storage.googleapis.com/reusedata/ucsc_database/refGene_m...
#> BFC130 https://storage.googleapis.com/reusedata/ucsc_database/refGene_m...
#> BFC131 /tmp/RtmptrHkdM/gcpData/GRCh37_to_GRCh38.chain
The data supports user-friendly discovery and access through the
ReUseData
portal, where detailed instructions are provided for
straight-forward incorporation into data analysis pipelines run on
local computing nodes, web resources, and cloud computing platforms
(e.g., Terra, CGC).
Here we provide a function meta_data()
to create a data frame that
contains all information about the data sets in the specified file
path (recursively), including the annotation file ($yml
column),
parameter values for the recipe ($params
column), data file path
($output
column), keywords for data file (notes
columns), date of
data generation (date
column), and any tag if available (tag
column).
Use cleanup = TRUE
to cleanup any invalid or expired/older
intermediate files.
mt <- meta_data(outdir)
head(mt)
#> yml
#> 1 /tmp/RtmptrHkdM/SharedData/echo_out_Hello_World!_outfile.yml
#> params output
#> 1 input: Hello World!; outfile: outfile /tmp/RtmptrHkdM/SharedData/outfile.txt
#> notes date
#> 1 echo hello world txt 2024-10-29 19:07:12.162606
sessionInfo()
#> R Under development (unstable) (2024-10-21 r87258)
#> Platform: x86_64-pc-linux-gnu
#> Running under: Ubuntu 24.04.1 LTS
#>
#> Matrix products: default
#> BLAS: /home/biocbuild/bbs-3.21-bioc/R/lib/libRblas.so
#> LAPACK: /usr/lib/x86_64-linux-gnu/lapack/liblapack.so.3.12.0
#>
#> locale:
#> [1] LC_CTYPE=en_US.UTF-8 LC_NUMERIC=C
#> [3] LC_TIME=en_GB LC_COLLATE=C
#> [5] LC_MONETARY=en_US.UTF-8 LC_MESSAGES=en_US.UTF-8
#> [7] LC_PAPER=en_US.UTF-8 LC_NAME=C
#> [9] LC_ADDRESS=C LC_TELEPHONE=C
#> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C
#>
#> time zone: America/New_York
#> tzcode source: system (glibc)
#>
#> attached base packages:
#> [1] stats4 stats graphics grDevices utils datasets methods
#> [8] base
#>
#> other attached packages:
#> [1] ReUseData_1.7.0 Rcwl_1.23.0 S4Vectors_0.45.0
#> [4] BiocGenerics_0.53.0 yaml_2.3.10 BiocStyle_2.35.0
#>
#> loaded via a namespace (and not attached):
#> [1] tidyselect_1.2.1 dplyr_1.1.4 blob_1.2.4
#> [4] filelock_1.0.3 R.utils_2.12.3 fastmap_1.2.0
#> [7] BiocFileCache_2.15.0 promises_1.3.0 digest_0.6.37
#> [10] base64url_1.4 mime_0.12 lifecycle_1.0.4
#> [13] RSQLite_2.3.7 magrittr_2.0.3 compiler_4.5.0
#> [16] rlang_1.1.4 sass_0.4.9 progress_1.2.3
#> [19] tools_4.5.0 utf8_1.2.4 data.table_1.16.2
#> [22] knitr_1.48 prettyunits_1.2.0 brew_1.0-10
#> [25] htmlwidgets_1.6.4 bit_4.5.0 curl_5.2.3
#> [28] reticulate_1.39.0 RColorBrewer_1.1-3 batchtools_0.9.17
#> [31] BiocParallel_1.41.0 purrr_1.0.2 withr_3.0.2
#> [34] R.oo_1.26.0 grid_4.5.0 fansi_1.0.6
#> [37] git2r_0.35.0 xtable_1.8-4 debugme_1.2.0
#> [40] cli_3.6.3 rmarkdown_2.28 DiagrammeR_1.0.11
#> [43] crayon_1.5.3 generics_0.1.3 httr_1.4.7
#> [46] visNetwork_2.1.2 DBI_1.2.3 cachem_1.1.0
#> [49] parallel_4.5.0 BiocManager_1.30.25 basilisk_1.19.0
#> [52] vctrs_0.6.5 Matrix_1.7-1 jsonlite_1.8.9
#> [55] dir.expiry_1.15.0 bookdown_0.41 hms_1.1.3
#> [58] bit64_4.5.2 jquerylib_0.1.4 RcwlPipelines_1.23.0
#> [61] glue_1.8.0 codetools_0.2-20 stringi_1.8.4
#> [64] later_1.3.2 tibble_3.2.1 pillar_1.9.0
#> [67] basilisk.utils_1.19.0 rappdirs_0.3.3 htmltools_0.5.8.1
#> [70] R6_2.5.1 dbplyr_2.5.0 evaluate_1.0.1
#> [73] shiny_1.9.1 lattice_0.22-6 R.methodsS3_1.8.2
#> [76] png_0.1-8 backports_1.5.0 memoise_2.0.1
#> [79] httpuv_1.6.15 bslib_0.8.0 Rcpp_1.0.13
#> [82] checkmate_2.3.2 xfun_0.48 pkgconfig_2.0.3