Skip to content

Commit 57bb27b

Browse files
committed
Make -d, --max-depth 0 set the depth to unlimited
This is to behave the same way as `samtools mpileup` Resolves #2435
1 parent bb6d531 commit 57bb27b

File tree

3 files changed

+19
-6
lines changed

3 files changed

+19
-6
lines changed

NEWS

Lines changed: 4 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -26,6 +26,10 @@ Changes affecting specific commands:
2626
- Fix an error in parsing -i/-e command line options where the `qry:` and `gt:` prefix was
2727
not stripped (#2432)
2828

29+
* bcftools mpileup
30+
31+
- Make `-d, --max-depth 0` set the depth to unlimited (#2435)
32+
2933
* bcftools norm
3034

3135
- Make the -i/-e filtering option work for all options, such as line merging and

doc/bcftools.txt

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -2308,7 +2308,8 @@ those scenarios.
23082308
This behavior was problematic when working with a combination of
23092309
single- and multi-sample bams, therefore in *bcftools mpileup* the user
23102310
is given the full control (and responsibility), and an informative message
2311-
is printed instead [250]
2311+
is printed instead. Passing zero for this option sets it to the highest possible
2312+
value, effectively removing the depth limit [250]
23122313

23132314
*-E, --redo-BAQ*::
23142315
Recalculate BAQ on the fly, ignore existing BQ tags

mpileup.c

Lines changed: 13 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -930,11 +930,19 @@ static int mpileup(mplp_conf_t *conf)
930930
// init mpileup
931931
conf->iter = bam_mplp_init(conf->nfiles, mplp_func, (void**)conf->mplp_data);
932932
if ( conf->flag & MPLP_SMART_OVERLAPS ) bam_mplp_init_overlaps(conf->iter);
933-
fprintf(stderr, "[%s] maximum number of reads per input file set to -d %d\n", __func__, conf->max_depth);
934-
if ( (double)conf->max_depth * conf->nfiles > 1<<20)
935-
fprintf(stderr, "Warning: Potential memory hog, up to %.0fM reads in the pileup!\n", (double)conf->max_depth*conf->nfiles);
936-
if ( (double)conf->max_depth * conf->nfiles / nsmpl < 250 )
937-
fprintf(stderr, "Note: The maximum per-sample depth with -d %d is %.1fx\n", conf->max_depth,(double)conf->max_depth * conf->nfiles / nsmpl);
933+
if ( !conf->max_depth )
934+
{
935+
conf->max_depth = INT_MAX;
936+
fprintf(stderr, "[%s] Max depth set to maximum value (%d)\n", __func__, INT_MAX);
937+
}
938+
else
939+
{
940+
fprintf(stderr, "[%s] maximum number of reads per input file set to -d %d\n", __func__, conf->max_depth);
941+
if ( (double)conf->max_depth * conf->nfiles > 1<<20)
942+
fprintf(stderr, "Warning: Potential memory hog, up to %.0fM reads in the pileup!\n", (double)conf->max_depth*conf->nfiles);
943+
if ( (double)conf->max_depth * conf->nfiles / nsmpl < 250 )
944+
fprintf(stderr, "Note: The maximum per-sample depth with -d %d is %.1fx\n", conf->max_depth,(double)conf->max_depth * conf->nfiles / nsmpl);
945+
}
938946
bam_mplp_set_maxcnt(conf->iter, conf->max_depth);
939947
conf->max_indel_depth = conf->max_indel_depth * nsmpl;
940948
conf->bcf_rec = bcf_init1();

0 commit comments

Comments
 (0)