Skip to content

Commit

Permalink
fix some memory leaks
Browse files Browse the repository at this point in the history
  • Loading branch information
kriemo committed Sep 15, 2023
1 parent 8ab8443 commit f7c91a6
Show file tree
Hide file tree
Showing 5 changed files with 41 additions and 30 deletions.
2 changes: 1 addition & 1 deletion R/differential_editing.R
Original file line number Diff line number Diff line change
Expand Up @@ -658,7 +658,7 @@ find_scde_sites <- function(
assay.type = "nAlt",
BPPARAM = BPPARAM)

grps <- unique(sce[[group]])
grps <- unique(as.character(sce[[group]]))
grp_pairs <- as.data.frame(t(utils::combn(grps, 2)))
colnames(grp_pairs) <- c("first", "second")

Expand Down
2 changes: 1 addition & 1 deletion R/pileup.R
Original file line number Diff line number Diff line change
Expand Up @@ -278,7 +278,7 @@ setup_valid_regions <- function(bam, chroms, region = NULL, fasta = NULL) {
msg <- paste(missing_chroms, collapse = "\n")
cli::cli_warn(c(
"the following chromosomes are not present in the fasta file:",
"msg"
"{msg}"
))
chroms_to_process <- setdiff(chroms_to_process, missing_chroms)
}
Expand Down
44 changes: 27 additions & 17 deletions src/plp.c
Original file line number Diff line number Diff line change
Expand Up @@ -631,7 +631,7 @@ static int run_pileup(char** cbampaths, char** cindexes,

const bam_pileup1_t** plp;
mplp_ref_t mp_ref = MPLP_REF_INIT;
bam_mplp_t iter;
bam_mplp_t iter = NULL;
bam_hdr_t* h = NULL; /* header of first file in input list */
char* ref;

Expand Down Expand Up @@ -664,33 +664,43 @@ static int run_pileup(char** cbampaths, char** cindexes,
data[i] = calloc(1, sizeof(mplp_aux_t));
data[i]->fp = sam_open(cbampaths[i], "rb");
if (!data[i]->fp) {
Rf_error("[raer internal] failed to open %s: %s\n",
REprintf("[raer internal] failed to open %s: %s\n",
cbampaths[i], strerror(errno));
ret = -1;
goto fail;
}
if (conf->fai_fname) {
if (hts_set_fai_filename(data[i]->fp, conf->fai_fname) != 0) {
Rf_error("[raer internal] failed to process %s: %s\n",
REprintf("[raer internal] failed to process %s: %s\n",
conf->fai_fname, strerror(errno));
ret = -1;
goto fail;
}
}
data[i]->conf = conf;
data[i]->ref = &mp_ref;

h_tmp = sam_hdr_read(data[i]->fp);
if (!h_tmp) {
Rf_error("[raer internal] fail to read the header of %s\n", cbampaths[i]);
REprintf("[raer internal] fail to read the header of %s\n", cbampaths[i]);
ret = -1;
goto fail;
}

if (conf->reg) {
hts_idx_t* idx = NULL;
idx = sam_index_load2(data[i]->fp, cbampaths[i], cindexes[i]) ;

if (idx == NULL) {
Rf_error("[raer internal] fail to load bamfile index for %s\n", cbampaths[i]);
REprintf("[raer internal] fail to load bamfile index for %s\n", cbampaths[i]);
ret = -1;
goto fail;
}

if ((data[i]->iter=sam_itr_querys(idx, h_tmp, conf->reg)) == 0) {
Rf_error("[raer internal] fail to parse region '%s' with %s\n", conf->reg, cbampaths[i]);
REprintf("[raer internal] fail to parse region '%s' with %s\n", conf->reg, cbampaths[i]);
ret = -1;
goto fail;
}

if (i == 0) {
Expand Down Expand Up @@ -891,11 +901,11 @@ static int run_pileup(char** cbampaths, char** cindexes,
fail:
if (ret >= 0) finish_PLP_DATA(pd);

bam_mplp_destroy(iter);
bam_hdr_destroy(h);
if(iter) bam_mplp_destroy(iter);
if(h) bam_hdr_destroy(h);

for (i = 0; i < conf->nbam; ++i) {
sam_close(data[i]->fp);
if(data[i]->fp) sam_close(data[i]->fp);
if (data[i]->iter) hts_itr_destroy(data[i]->iter);
free(data[i]);
clear_pcounts(&plpc[i]);
Expand All @@ -908,9 +918,9 @@ static int run_pileup(char** cbampaths, char** cindexes,

}

free(data);
free(plp);
free(n_plp);
if(data) free(data);
if(plp) free(plp);
if(n_plp) free(n_plp);
free(mp_ref.ref[0]);
free(mp_ref.ref[1]);

Expand All @@ -922,11 +932,11 @@ static int run_pileup(char** cbampaths, char** cindexes,
R_Free(pall);
}

R_Free(plpc);
if (pd->fps) R_Free(pd->fps);
R_Free(pd->pdat);

if (site_stats) R_Free(site_stats);
if(plpc) R_Free(plpc);
if(pd->fps) R_Free(pd->fps);
if(pd->pdat) R_Free(pd->pdat);
if(pd->sdat) R_Free(pd->sdat);
if(site_stats) R_Free(site_stats);

return ret;
}
Expand Down
19 changes: 10 additions & 9 deletions src/sc-plp.c
Original file line number Diff line number Diff line change
Expand Up @@ -68,7 +68,8 @@ static void clear_umi(umimap_t uhash) {
}

/*! @function
@abstract free keys, values and clear cbumi_map_t hashmap
@abstract free UMI values and clear cbumi_map_t hashmap. Keys are retained to
be freed at end of pileup loop.
*/
static void clear_cb_umiset(cbumi_map_t cbhash) {
khint_t k;
Expand Down Expand Up @@ -100,6 +101,7 @@ static void free_hashmaps(cbumi_map_t cbhash, str2intmap_t cbidx) {
if (cbhash) {
for (k = kh_begin(cbhash); k < kh_end(cbhash); ++k) {
if (!kh_exist(cbhash, k)) continue;
free((char*) kh_key(cbhash, k));
cdat = kh_value(cbhash, k);
clear_umi(cdat->umi);
kh_destroy(umimap, cdat->umi);
Expand Down Expand Up @@ -248,14 +250,13 @@ static int count_record(const bam_pileup1_t* p, sc_mplp_conf_t* conf, payload_t*

if (pld->strand != strand) return (1);

cb_cpy = strdup(cb);

k = kh_get(str2intmap, conf->cbidx, cb_cpy);
k = kh_get(str2intmap, conf->cbidx, cb);
if (k == kh_end(conf->cbidx)) {
free(cb_cpy);
return (0);
}

cb_cpy = strdup(cb);

k = kh_put(cbumimap, conf->cbmap, cb_cpy, &cret);
if (cret < 0) {
free(cb_cpy);
Expand Down Expand Up @@ -729,17 +730,17 @@ static int run_scpileup(sc_mplp_conf_t* conf, char* bamfn, char* index, char* ba

fail:
if (iter) bam_mplp_destroy(iter);
bam_hdr_destroy(h);
if(h) bam_hdr_destroy(h);

for (i = 0; i < nbam; ++i) {
sam_close(data[0]->fp);
if (data[0]->iter) hts_itr_destroy(data[0]->iter);
if (data[0]->idx) hts_idx_destroy(data[0]->idx);
free(data[0]);
}
free(data);
free(plp);
free(n_plp);
if(data) free(data);
if(plp) free(plp);
if(n_plp) free(n_plp);

return ret;
}
Expand Down
4 changes: 2 additions & 2 deletions tests/testthat/test_pileup_sites.R
Original file line number Diff line number Diff line change
Expand Up @@ -104,15 +104,15 @@ test_that("pileup regional query works", {
expect_equal(end(res), c(203, 204, 205))

# chr1 does not exist
expect_error(suppressWarnings(pileup_sites(bamfn, fafn, region = "chr1")))
expect_error(pileup_sites(bamfn, fafn, region = "chr1"))

res <- pileup_sites(bamfn, fafn, chrom = "SSR3")
expect_equal(nrow(res), 529)
})

test_that("incorrect regional query is caught", {
# will produce warning that chrHello is not in bamfile and an error
expect_error(suppressWarnings(pileup_sites(bamfn, fafn, region = "chrHello")))
expect_error(pileup_sites(bamfn, fafn, region = "chrHello"))
})

test_that("missing files are caught", {
Expand Down

0 comments on commit f7c91a6

Please sign in to comment.