Degradome-seq.snakefile 21 KB
Newer Older
1
2
3
"""
Snakefile to analyse Degradome-seq data.
"""
4
5
#TODO: compute TPM, and compare with RNA-seq to have "degradation efficiency" (like TE in Ribo-seq pipeline)
#TODO: bigwig files
6
7
8
9
10
11
12
13
import sys
major, minor = sys.version_info[:2]
if major < 3 or (major == 3 and minor < 6):
    sys.exit("Need at least python 3.6\n")

import os
OPJ = os.path.join
from glob import glob
14
15
from pickle import load
import numpy as np
16
17
import pandas as pd

18
19
20
21
from itertools import chain, repeat
from functools import reduce
from operator import or_ as union

22
from smincludes import rules as irules
Blaise Li's avatar
Blaise Li committed
23
from smwrappers import wrappers_dir
24
25
from libworkflows import (
    cleanup_and_backup,
26
    column_converter,
27
28
29
30
    ensure_relative,
    feature_orientation2stranded,
    get_chrom_sizes,
    read_feature_counts,
31
    sum_by_family,
32
33
    sum_feature_counts,
    wc_applied)
34
35
36
37
38
39
40
41
# http://sailfish.readthedocs.io/en/master/library_type.html
LIB_TYPE = config["lib_type"]
# key: library name
# value: dictionary
#   key: replicate number
#   value: path to raw data
lib2raw = config["lib2raw"]
LIBS = list(lib2raw.keys())
42
43
44
45
46
# Libraries for which we have matching RNA-seq data
# so that degradation efficiency can be computed
EFF_LIBS = list(config["transcriptome_TPM"].keys())
# What type of "efficency" are we computing? Degradation efficiency.
[this_TPM, ref_TPM, eff_name] = ["Degradome_TPM", "RNA_TPM", "DE"]
47
48
49
50
51
52
53
54
55
56
REPS = config["replicates"]
CONDITIONS = [{
    "lib" : lib,
    "rep" : rep} for rep in REPS for lib in LIBS]
# We use this for various things in order to have always the same library order:
COND_NAMES = ["_".join((
    cond["lib"],
    cond["rep"])) for cond in CONDITIONS]
# COND_COLUMNS = pd.DataFrame(CONDITIONS).assign(
#     cond_name=pd.Series(COND_NAMES).values).set_index("cond_name")
57
58
59
60
# Differential degradation
DD_COND_PAIRS = config.get("dd_cond_pairs", [])
DD_CONTRASTS = [f"{cond1}_vs_{cond2}" for [cond1, cond2] in DD_COND_PAIRS]
CONTRAST2PAIR = dict(zip(DD_CONTRASTS, DD_COND_PAIRS))
61

62
# TODO: have a distinct subset of biotypes for join_all_counts, that are guaranteed overlap-free.
63
BIOTYPES = config["biotypes"]
64
ORIENTATIONS = ["all", "fwd", "rev"]
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86

aligner = config["aligner"]
alignment_settings = {
    "bowtie2": "",
    "hisat2": "",
    "crac": "-k 20 --stranded --use-x-in-cigar"}

########################
# Genome configuration #
########################
genome_dict = config["genome_dict"]
genome = genome_dict["name"]
chrom_sizes = get_chrom_sizes(genome_dict["size"])
genomelen = sum(chrom_sizes.values())
genome_db = genome_dict["db"][aligner]
# bed file binning the genome in 10nt bins
genome_binned = genome_dict["binned"]
annot_dir = genome_dict["annot_dir"]
# TODO: figure out the difference between OPJ(convert_dir, "wormid2name.pickle") and genome_dict["converter"]
convert_dir = genome_dict["convert_dir"]
gene_lists_dir = genome_dict["gene_lists_dir"]
avail_id_lists = set(glob(OPJ(gene_lists_dir, "*_ids.txt")))
87
88
89
90
91
#output_dir = config["output_dir"]
#workdir: config["output_dir"]
output_dir = os.path.abspath(".")
log_dir = OPJ(f"logs_{genome}")
data_dir = OPJ("data")
92

93
counter = "feature_count"
94
counts_dir = OPJ(aligner, f"mapped_{genome}", counter)
95
96
# Limit risks of ambiguity by imposing replicates to be numbers
# and restricting possible forms of some other wildcards
97

98
99
100
101
102
103
104
105
106
107
wildcard_constraints:
    lib="|".join(LIBS),
    rep="\d+",


rule all:
    input:
        #expand(
        #    OPJ(data_dir, "trimmed", "{lib}_{rep}.fastq.gz"),
        #    lib=LIBS, rep=REPS),
108
        # expand(
109
        #     OPJ(aligner, f"mapped_{genome}", "{lib}_{rep}_on_%s_sorted.bam" % genome),
110
111
112
113
114
        #     lib=LIBS, rep=REPS),
        # expand(
        #     OPJ(counts_dir,
        #         f"all_on_{genome}", "{biotype}_{orientation}_counts.txt"),
        #     biotype=BIOTYPES, orientation=ORIENTATIONS),
115
        expand(
116
            OPJ(counts_dir, "summaries",
117
118
                "all_on_%s_{orientation}_counts.txt" % genome),
            orientation=ORIENTATIONS),
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
        # expand(
        #     OPJ(counts_dir, f"all_on_{genome}",
        #         "{biotype}_{orientation}_TPM.txt"),
        #     biotype=["alltypes"], orientation=ORIENTATIONS),
        # expand(
        #     OPJ(counts_dir, "{lib}_mean_on_%s" % genome,
        #         "{lib}_mean_{biotype}_{orientation}_TPM.txt"),
        #     lib=LIBS, biotype=["alltypes"], orientation=ORIENTATIONS),
        # degradation_efficiency
        # expand(
        #     OPJ(counts_dir, f"all_on_{genome}", "{biotype}_{orientation}_%s.txt" % eff_name),
        #     lib=LIBS, biotype=["alltypes"], orientation=ORIENTATIONS),
        expand(
            OPJ(counts_dir, "diff_%s" % eff_name, "{contrast}", "{orientation}_{biotype}", "{contrast}_diff_%s.txt" % eff_name),
            contrast=DD_CONTRASTS, orientation=["fwd"], biotype=["alltypes"]),
134
135
136
137
138
139
140
141
142
143
144
145
146
147


include: ensure_relative(irules["link_raw_data"], workflow.basedir)


rule trim_reads:
    input:
        raw = rules.link_raw_data.output.link,
    output:
        trimmed = OPJ(data_dir, "trimmed", "{lib}_{rep}.fastq.gz"),
    params:
        trim3 = 4,
    shell:
        """
148
        seqtk trimfq -b {params.trim3} {input.raw} | gzip > {output.trimmed}
149
150
151
152
153
154
155
156
        """


rule map_on_genome:
    input:
        fastq = rules.trim_reads.output.trimmed,
    output:
        # temp because it uses a lot of space
157
158
        sam = temp(OPJ(aligner, f"mapped_{genome}", "{lib}_{rep}_on_%s.sam" % genome)),
        nomap_fastq = OPJ(aligner, f"not_mapped_{genome}", "{lib}_{rep}_unmapped_on_%s.fastq.gz" % genome),
159
160
161
162
163
164
165
166
167
168
169
170
171
172
    params:
        aligner = aligner,
        index = genome_db,
        settings = alignment_settings[aligner],
    message:
        "Mapping {wildcards.lib}_{wildcards.rep} on %s." % genome
    log:
        log = OPJ(log_dir, "map_on_genome_{lib}_{rep}.log"),
        err = OPJ(log_dir, "map_on_genome_{lib}_{rep}.err")
    resources:
        io=5
    threads:
        4
    wrapper:
Blaise Li's avatar
Blaise Li committed
173
        f"file://{wrappers_dir[0]}/map_on_genome"
174
175
176
177
178
179


rule sam2indexedbam:
    input:
        sam = rules.map_on_genome.output.sam,
    output:
180
181
        sorted_bam = protected(OPJ(aligner, f"mapped_{genome}", "{lib}_{rep}_on_%s_sorted.bam" % genome)),
        index = protected(OPJ(aligner, f"mapped_{genome}", "{lib}_{rep}_on_%s_sorted.bam.bai" % genome)),
182
183
184
185
186
187
188
    message:
        "Sorting and indexing sam file for {wildcards.lib}_{wildcards.rep}."
    log:
        log = OPJ(log_dir, "sam2indexedbam_{lib}_{rep}.log"),
        err = OPJ(log_dir, "sam2indexedbam_{lib}_{rep}.err"),
    threads:
        8
189
190
    resources:
        mem_mb=4100
191
    wrapper:
Blaise Li's avatar
Blaise Li committed
192
        f"file://{wrappers_dir[0]}/sam2indexedbam"
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207


def biotype2annot(wildcards):
    if wildcards.biotype.endswith("_rmsk_families"):
        biotype = wildcards.biotype[:-9]
    else:
        biotype = wildcards.biotype
    return OPJ(annot_dir, f"{biotype}.gtf")


rule feature_count_reads:
    input:
        sorted_bam = rules.sam2indexedbam.output.sorted_bam,
        index = rules.sam2indexedbam.output.index,
    output:
208
209
        counts = OPJ(counts_dir, "{lib}_{rep}_on_%s" % genome, "{lib}_{rep}_{biotype}_{orientation}_counts.txt"),
        counts_converted = OPJ(counts_dir, "{lib}_{rep}_on_%s" % genome, "{lib}_{rep}_{biotype}_{orientation}_counts_gene_names.txt"),
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
    params:
        stranded = feature_orientation2stranded(LIB_TYPE),
        annot = biotype2annot,
        # pickled dictionary that associates gene ids to gene names
        converter = genome_dict["converter"]
    message:
        "Counting {wildcards.orientation} {wildcards.biotype} reads for {wildcards.lib}_{wildcards.rep} with featureCounts."
    benchmark:
        OPJ(log_dir, "feature_count_reads", "{lib}_{rep}_{biotype}_{orientation}_benchmark.txt")
    log:
        log = OPJ(log_dir, "feature_count_reads", "{lib}_{rep}_{biotype}_{orientation}.log"),
        err = OPJ(log_dir, "feature_count_reads", "{lib}_{rep}_{biotype}_{orientation}.err"),
    shell:
        """
        tmpdir=$(mktemp -dt "feature_{wildcards.lib}_{wildcards.rep}_{wildcards.biotype}_{wildcards.orientation}.XXXXXXXXXX")
        cmd="featureCounts -a {params.annot} -o {output.counts} -t transcript -g "gene_id" -O -M --primary -s {params.stranded} --fracOverlap 0 --tmpDir ${{tmpdir}} {input.sorted_bam}"
        featureCounts -v 2> {log.log}
        echo ${{cmd}} 1>> {log.log}
        eval ${{cmd}} 1>> {log.log} 2> {log.err} || error_exit "featureCounts failed"
        rm -rf ${{tmpdir}}
        cat {output.counts} | id2name.py {params.converter} > {output.counts_converted}
        """


rule summarize_counts:
    """For a given library, write a summary of the read counts for the various biotypes."""
    input:
237
238
239
        biotype_counts_files = expand(
            OPJ(counts_dir, "{{lib}}_{{rep}}_on_%s" % genome, "{{lib}}_{{rep}}_{biotype}_{{orientation}}_counts.txt"),
            biotype=BIOTYPES),
240
    output:
241
        summary = OPJ(counts_dir, "summaries", "{lib}_{rep}_on_%s_{orientation}_counts.txt" % genome)
242
243
244
245
246
247
248
249
250
251
252
253
254
    run:
        sum_counter = sum_feature_counts
        with open(output.summary, "w") as summary_file:
            header = "\t".join(BIOTYPES)
            summary_file.write("%s\n" % header)
            sums = "\t".join((str(sum_counter(counts_file)) for counts_file in input.biotype_counts_files))
            summary_file.write("%s\n" % sums)


rule gather_read_counts_summaries:
    """Gather read count summaries across libraries: lib_rep -> all."""
    input:
        summary_tables = expand(OPJ(
255
            counts_dir, "summaries",
256
257
258
            "{name}_on_%s_{{orientation}}_counts.txt" % genome), name=COND_NAMES),
    output:
        summary_table = OPJ(
259
            counts_dir, "summaries",
260
261
262
263
            "all_on_%s_{orientation}_counts.txt" % genome),
    run:
        #summary_files = (OPJ(
        summary_files = [OPJ(
264
            counts_dir, "summaries",
265
266
267
268
269
270
271
272
273
274
275
276
277
            f"{cond_name}_on_{genome}_{wildcards.orientation}_counts.txt") for cond_name in COND_NAMES]
        try:
            summaries = pd.concat((pd.read_table(summary_file).T.astype(int) for summary_file in summary_files), axis=1)
        except ValueError:
            warnings.warn("Some data may be missing. Using floats instead of ints.\n")
            summaries = pd.concat((pd.read_table(summary_file).T.astype(float) for summary_file in summary_files), axis=1)
        summaries.columns = COND_NAMES
        summaries.to_csv(output.summary_table, na_rep="NA", sep="\t")


rule gather_counts:
    """For a given biotype, gather counts from all libraries in one table."""
    input:
278
279
280
281
        counts_tables = expand(
            OPJ(counts_dir, "{lib}_{rep}_on_%s" % genome,
            "{lib}_{rep}_{{biotype}}_{{orientation}}_counts.txt"),
            lib=LIBS, rep=REPS),
282
    output:
283
284
        counts_table = OPJ(
            counts_dir, f"all_on_{genome}", "{biotype}_{orientation}_counts.txt"),
285
286
287
288
289
290
291
292
    wildcard_constraints:
        # Avoid ambiguity with join_all_counts
        biotype = "|".join(BIOTYPES)
    run:
        # Gathering the counts data
        ############################
        #counts_files = (OPJ(
        counts_files = [OPJ(
293
            counts_dir,
294
295
296
297
298
299
300
301
302
303
304
305
306
            f"{cond_name}_on_{genome}",
            f"{cond_name}_{wildcards.biotype}_{wildcards.orientation}_counts.txt") for cond_name in COND_NAMES]
            #"{biotype}_{orientation}_counts.txt".format(biotype=wildcards.biotype, orientation=wildcards.orientation)) for cond_name in COND_NAMES)
        counts_data = pd.concat(
            map(read_feature_counts, counts_files),
            axis=1).fillna(0).astype(int)
        counts_data.columns = COND_NAMES
        # Simple_repeat|Simple_repeat|(TTTTTTG)n:1
        # Simple_repeat|Simple_repeat|(TTTTTTG)n:2
        # Simple_repeat|Simple_repeat|(TTTTTTG)n:3
        # Simple_repeat|Simple_repeat|(TTTTTTG)n:4
        # -> Simple_repeat|Simple_repeat|(TTTTTTG)n
        if wildcards.biotype.endswith("_rmsk_families"):
307
            counts_data = sum_by_family(counts_data)
308
309
310
311
        counts_data.index.names = ["gene"]
        counts_data.to_csv(output.counts_table, na_rep="NA", sep="\t")


312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
def make_tag_association(dfs, tag):
    """Associates a tag "tag" to the union of the indices of dataframes *dfs*."""
    idx = reduce(union, (df.index for df in dfs))
    return pd.DataFrame(list(zip(idx, repeat(tag)))).set_index(0)


rule associate_biotype:
    """This rule uses the read count matrices to associate a biotype to each gene identifier."""
    input:
        counts_tables = expand(OPJ(
            counts_dir,
            "all_on_%s" % genome,
            "{biotype}_all_counts.txt"), biotype=[b_name for b_name in BIOTYPES if not b_name.startswith("protein_coding_")]),
    output:
        tags_table = OPJ(counts_dir, "all_on_%s" % genome, "id2tags.txt"),
    run:
        biotype2tags_tables = {}
        for biotype in BIOTYPES:
            if biotype.startswith("protein_coding_"):
                continue
            biotype2tags_tables[biotype] = make_tag_association(
                (pd.read_table(OPJ(
                    counts_dir,
                    f"all_on_{genome}",
                    f"{biotype}_all_counts.txt"), index_col="gene"),), biotype)
        tags_table = pd.concat(chain(biotype2tags_tables.values()))
        tags_table.index.names = ["gene"]
        tags_table.columns = ["biotype"]
        tags_table.to_csv(output.tags_table, sep="\t")


def add_tags_column(data, tags_table, tag_name):
    """Adds a column *tag_name* to *data* based on the DataFrame *tag_table*
    associating tags to row names."""
    # Why "inner" ?
    df = pd.concat((data, pd.read_table(tags_table, index_col=0)), join="inner", axis=1)
    df.columns = (*data.columns, tag_name)
    return df


rule join_all_counts:
    """Concatenate counts for all biotypes into "alltypes"."""
    input:
        counts_tables = expand(
            OPJ(counts_dir,
                f"all_on_{genome}", "{biotype}_{{orientation}}_counts.txt"),
            biotype=BIOTYPES),
    output:
        counts_table = OPJ(
            counts_dir,
            f"all_on_{genome}", "alltypes_{orientation}_counts.txt"),
    run:
        counts_data = pd.concat((pd.read_table(table, index_col="gene") for table in input.counts_tables)).drop_duplicates()
        assert len(counts_data.index.unique()) == len(counts_data.index), "Some genes appear several times in the counts table."
        counts_data.index.names = ["gene"]
        counts_data.to_csv(output.counts_table, sep="\t")


@wc_applied
def source_counts(wildcards):
    """Determines from which rule the gathered counts should be sourced."""
    if wildcards.biotype == "alltypes":
        return rules.join_all_counts.output.counts_table
    else:
        # "Directly" from the counts gathered across libraries
        return rules.gather_counts.output.counts_table


rule compute_RPK:
    """For a given biotype, compute the corresponding RPK value (reads per kilobase)."""
    input:
        counts_data = source_counts,
        #counts_data = rules.gather_counts.output.counts_table,
385
        #counts_data = OPJ(aligner, f"mapped_{genome}", "feature_count",
386
387
388
389
390
391
392
393
394
395
396
397
398
        #    f"all_on_{genome}", "{biotype}_{orientation}_counts.txt"),
    output:
        rpk_file = OPJ(counts_dir,
            f"all_on_{genome}", "{biotype}_{orientation}_RPK.txt"),
    params:
        feature_lengths_file = OPJ(annot_dir, "union_exon_lengths.txt"),
    # run:
    #     counts_data = pd.read_table(input.counts_data, index_col="gene")
    #     feature_lengths = pd.read_table(params.feature_lengths_file, index_col="gene")
    #     common = counts_data.index.intersection(feature_lengths.index)
    #     rpk = 1000 * counts_data.loc[common].div(feature_lengths.loc[common]["union_exon_len"], axis="index")
    #     rpk.to_csv(output.rpk_file, sep="\t")
    wrapper:
Blaise Li's avatar
Blaise Li committed
399
        f"file://{wrappers_dir[0]}/compute_RPK"
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420


# Compute TPM using total number of mappers divided by genome length
# (not sum of RPK across biotypes: some mappers may not be counted)
# No, doesn't work: mappers / genome length not really comparable
# Needs to be done on all_types
rule compute_TPM:
    """For a given biotype, compute the corresponding TPM value (reads per kilobase per million mappers)."""
    input:
        rpk_file = rules.compute_RPK.output.rpk_file
    output:
        tpm_file = OPJ(counts_dir,
            f"all_on_{genome}", "{biotype}_{orientation}_TPM.txt"),
    # The sum must be done over all counted features
    wildcard_constraints:
        biotype = "alltypes"
    # run:
    #     rpk = pd.read_table(input.rpk_file, index_col="gene")
    #     tpm = 1000000 * rpk / rpk.sum()
    #     tpm.to_csv(output.tpm_file, sep="\t")
    wrapper:
Blaise Li's avatar
Blaise Li committed
421
        f"file://{wrappers_dir[0]}/compute_TPM"
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520


# TODO: Is it better to compute the mean and then the fold of the means?
# Here, the matching between Degradome-seq and RNA-seq replicates is arbitrary.
# We don't even always have the same number of replicates,
# so we work with means across replicates.
rule compute_mean_TPM:
    input:
        all_tmp_file = rules.compute_TPM.output.tpm_file
    output:
        tpm_file = OPJ(counts_dir,
            "{lib}_mean_on_%s" % genome, "{lib}_mean_{biotype}_{orientation}_TPM.txt"),
    run:
        tpm = pd.read_table(
            input.all_tmp_file, index_col="gene",
            usecols=["gene", *[f"{wildcards.lib}_{rep}" for rep in REPS]])
        tpm_mean = tpm.mean(axis=1)
        # This is a Series
        tpm_mean.name = wildcards.lib
        tpm_mean.to_csv(output.tpm_file, sep="\t", header=True)


# TODO
rule compute_efficiency:
    """This is actually a TPM fold between library types."""
    input:
        mean_tpm = rules.compute_mean_TPM.output.tpm_file,
        tags_table = rules.associate_biotype.output.tags_table,
    output:
        eff_file = OPJ(counts_dir,
            "{lib}_mean_on_%s" % genome, "{lib}_{biotype}_{orientation}_%s.txt" % eff_name),
    run:
        mean_tpm = pd.read_table(input.mean_tpm, index_col="gene")
        transcriptome_tpm = pd.read_table(config["transcriptome_TPM"][wildcards.lib], index_col="gene")
        transcriptome_tpm.columns = mean_tpm.columns
        lfc = np.log2((mean_tpm + 1) / (transcriptome_tpm + 1))
        # TE_mut = (ribo_tpm_mut / transcritome_tpm_mut)
        # TE_WT = (ribo_tpm_WT / transcritome_tpm_WT)
        # log2(TE_mut / TE_WT) = log2(TE_mut) - log2(TE_WT)
        tpm_and_eff = pd.concat([mean_tpm, transcriptome_tpm, lfc], axis=1)
        tpm_and_eff.columns = [this_TPM, ref_TPM, eff_name]
        tpm_and_eff.index.name = "gene"
        # Converting gene IDs
        ######################
        with open(OPJ(convert_dir, "wormid2cosmid.pickle"), "rb") as dict_file:
            tpm_and_eff = tpm_and_eff.assign(cosmid=tpm_and_eff.apply(
                column_converter(load(dict_file)), axis=1))
        with open(OPJ(convert_dir, "wormid2name.pickle"), "rb") as dict_file:
            tpm_and_eff = tpm_and_eff.assign(name=tpm_and_eff.apply(
                column_converter(load(dict_file)), axis=1))
        tpm_and_eff = add_tags_column(tpm_and_eff, input.tags_table, "biotype")
        tpm_and_eff.to_csv(output.eff_file, sep="\t", na_rep="NA")


rule gather_efficiency:
    input:
        eff_files = expand(OPJ(
            counts_dir,
            "{lib}_mean_on_%s" % genome,
            "{lib}_{{biotype}}_{{orientation}}_%s.txt" % eff_name), lib=EFF_LIBS),
        #tags_table = rules.associate_biotype.output.tags_table,
    output:
        eff_file = OPJ(
            counts_dir,
            "all_on_%s" % genome, "{biotype}_{orientation}_%s.txt" % eff_name),
    run:
        eff_data = pd.concat(
            [pd.read_table(eff_file, index_col="gene", usecols=["gene", eff_name]) for eff_file in input.eff_files],
            axis=1)
        eff_data.columns = [f"{lib}_{eff_name}" for lib in EFF_LIBS]
        eff_data.index.name = "gene"
        # Converting gene IDs
        ######################
        with open(OPJ(convert_dir, "wormid2cosmid.pickle"), "rb") as dict_file:
            eff_data = eff_data.assign(cosmid=eff_data.apply(
                column_converter(load(dict_file)), axis=1))
        with open(OPJ(convert_dir, "wormid2name.pickle"), "rb") as dict_file:
            eff_data = eff_data.assign(name=eff_data.apply(
                column_converter(load(dict_file)), axis=1))
        #eff_data = add_tags_column(eff_data, input.tags_table, "biotype")
        eff_data.to_csv(output.eff_file, sep="\t", na_rep="NA")


rule compute_efficiency_difference:
    input:
        eff_file = rules.gather_efficiency.output.eff_file,
    output:
        diff_eff_file = OPJ(counts_dir, "diff_%s" % eff_name,
            "{contrast}", "{orientation}_{biotype}", "{contrast}_diff_%s.txt" % eff_name),
    run:
        eff_data = pd.read_table(input.eff_file, index_col="gene")
        (cond, ref) = CONTRAST2PAIR[wildcards.contrast]
        # TE_mut = (ribo_tpm_mut / transcritome_tpm_mut)
        # TE_WT = (ribo_tpm_WT / transcritome_tpm_WT)
        # log2(TE_mut / TE_WT) = log2(TE_mut) - log2(TE_WT)
        diff_eff = eff_data[f"{cond}_{eff_name}"] - eff_data[f"{ref}_{eff_name}"]
        diff_eff.to_csv(output.diff_eff_file, sep="\t", na_rep="NA")


521
522
onsuccess:
    print("Degradome-seq analysis finished.")
523
    cleanup_and_backup(output_dir, config, delete=True)
524
525
526
527
528
529

onerror:
    shell(f"rm -rf {output_dir}_err")
    shell(f"cp -rp {output_dir} {output_dir}_err")
    cleanup_and_backup(output_dir + "_err", config)
    print("Degradome-seq analysis failed.")