DADA2_MAKE_TABLE

https://img.shields.io/github/issues-pr/snakemake/snakemake-wrappers/bio/dada2/make-table?label=version%20update%20pull%20requests

DADA2 Build a sequence - sample table from denoised samples using dada2 makeSequenceTable function. Optional parameters are documented in the manual and the function is introduced in the dedicated tutorial section.

Example

This wrapper can be used in the following way:

rule dada2_make_table_se:
    input:
    # Inferred composition
        expand("denoised/{sample}.1.RDS", sample=['a','b'])
    output:
        "results/dada2/seqTab-se.RDS"
    # Even though this is an R wrapper, use named arguments in Python syntax
    # here, to specify extra parameters. Python booleans (`arg1=True`, `arg2=False`)
    # and lists (`list_arg=[]`) are automatically converted to R.
    # For a named list as an extra named argument, use a python dict
    #   (`named_list={name1=arg1}`).
    params:
        names=['a','b'] # Sample names instead of paths
    log:
        "logs/dada2/make-table/make-table-se.log"
    threads: 1 # set desired number of threads here
    wrapper:
        "v3.8.0/bio/dada2/make-table"

rule dada2_make_table_pe:
    input:
    # Merged composition
        expand("merged/{sample}.RDS", sample=['a','b'])
    output:
        "results/dada2/seqTab-pe.RDS"
    # Even though this is an R wrapper, use named arguments in Python syntax
    # here, to specify extra parameters. Python booleans (`arg1=True`, `arg2=False`)
    # and lists (`list_arg=[]`) are automatically converted to R.
    # For a named list as an extra named argument, use a python dict
    #   (`named_list={name1=arg1}`).
    params:
        names=['a','b'], # Sample names instead of paths
        orderBy="nsamples" # Change the ordering of samples
    log:
        "logs/dada2/make-table/make-table-pe.log"
    threads: 1 # set desired number of threads here
    wrapper:
        "v3.8.0/bio/dada2/make-table"

Note that input, output and log file paths can be chosen freely.

When running with

snakemake --use-conda

the software dependencies will be automatically deployed into an isolated environment before execution.

Software dependencies

  • bioconductor-dada2=1.30.0

Input/Output

Input:

  • A list of RDS files with denoised samples (se), or denoised and merged samples (pe)

Output:

  • RDS file with the table

Params

  • names: A list of sample names instead of paths

  • params: Any other optional arguments for makeSequenceTable(), please provide them as python key=value pairs

Authors

  • Charlie Pauvert

Code

# __author__ = "Charlie Pauvert"
# __copyright__ = "Copyright 2020, Charlie Pauvert"
# __email__ = "cpauvert@protonmail.com"
# __license__ = "MIT"

# Snakemake wrapper for building a sequence - sample table from denoised samples using dada2 makeSequenceTable function.

# Sink the stderr and stdout to the snakemake log file
# https://stackoverflow.com/a/48173272
log.file<-file(snakemake@log[[1]],open="wt")
sink(log.file)
sink(log.file,type="message")

library(dada2)

# If names are provided use them
nm<-if(is.null(snakemake@params[["names"]])) NULL else snakemake@params[["names"]]

# From a list of n lists to one named list of n elements
smps<-setNames(
               object=unlist(snakemake@input),
               nm=nm
               )
# Read the RDS into the list
smps<-lapply(smps, readRDS)

# Prepare arguments (no matter the order)
args<-list( samples = smps)
# Check if extra params are passed (apart from [["names"]])
if(length(snakemake@params) > 1 ){
       # Keeping only the named elements of the list for do.call() (apart from [["names"]])
       extra<-snakemake@params[ names(snakemake@params) != "" & names(snakemake@params) != "names" ]
       # Add them to the list of arguments
       args<-c(args, extra)
} else{
    message("No optional parameters. Using default parameters from dada2::makeSequenceTable()")
}

# Make table
seqTab<-do.call(makeSequenceTable, args)

# Store the table as a RDS file
saveRDS(seqTab, snakemake@output[[1]],compress = T)

# Proper syntax to close the connection for the log file
# but could be optional for Snakemake wrapper
sink(type="message")
sink()