Skip to content

Alignment API Usage

Import package

from pyBioTools import Alignment
from pyBioTools.common import jhelp, head

Reads_index

jhelp(Alignment.Reads_index)

Reads_index (input_fn, skip_unmapped, skip_secondary, skip_supplementary, verbose, quiet, progress, kwargs)

Index reads found in a coordinated sorted bam file by read_id. The created index file can be used to randon access the alignment file per read_id


  • input_fn (required) [str]

Path to the bam file to index

  • skip_unmapped (default: False) [bool]

Do not include unmapped reads in index

  • skip_secondary (default: False) [bool]

Do not include secondary alignment in index

  • skip_supplementary (default: False) [bool]

Do not include supplementary alignment in index

  • verbose (default: False) [bool]

  • quiet (default: False) [bool]

  • progress (default: False) [bool]

  • kwargs

Basic usage

Alignment.Reads_index("./data/sample_1.bam", verbose=True)
head("./data/sample_1.bam.idx.gz")
## Running Alignment Reads_index ##
    Checking Bam file
    Parsing reads
    Read counts summary
     Reads retained
      total: 13,684
      primary: 10,584
      secondary: 1,496
      unmapped: 1,416
      supplementary: 188

read_id                        ref_id pointer    index 
4f20c07e-183a-4483-9313-220390 chr-I  31653888   0     
f5fc15c4-0bd5-4e56-9509-625574 chr-I  31681836   1     
4c144dc6-705a-4d0c-951f-02fd7f chr-I  31685499   2     
6449b796-5339-4d61-8cf4-14658f chr-I  1338834944 3     
f0bea8a5-a7bc-455c-b0e8-abc6ff chr-I  1338869081 4     
f50af728-e016-4c74-b48f-fe7459 chr-I  1338876255 5     
3213cb66-ac07-468a-862e-4e5fd1 chr-I  1338878547 6     
370dbed1-1e84-4b40-87ff-e6468b chr-I  1338879643 7     
9e21a7d4-1b3e-4967-856e-ca4e99 chr-I  1338889457 8     


Excluding reads from index

Alignment.Reads_index("./data/sample_1.bam", verbose=True, skip_secondary=True, skip_supplementary=True, skip_unmapped=True)
head("./data/sample_1.bam.idx.gz")
## Running Alignment Reads_index ##
    Checking Bam file
    Parsing reads
    Read counts summary
     Reads retained
      primary: 10,584
      total: 10,584
     Reads discarded
      total: 3,100
      secondary: 1,496
      unmapped: 1,416
      supplementary: 188

read_id                        ref_id pointer    index 
4f20c07e-183a-4483-9313-220390 chr-I  31653888   0     
f5fc15c4-0bd5-4e56-9509-625574 chr-I  31681836   1     
4c144dc6-705a-4d0c-951f-02fd7f chr-I  31685499   2     
6449b796-5339-4d61-8cf4-14658f chr-I  1338834944 3     
f0bea8a5-a7bc-455c-b0e8-abc6ff chr-I  1338869081 4     
f50af728-e016-4c74-b48f-fe7459 chr-I  1338876255 5     
3213cb66-ac07-468a-862e-4e5fd1 chr-I  1338878547 6     
370dbed1-1e84-4b40-87ff-e6468b chr-I  1338879643 7     
9e21a7d4-1b3e-4967-856e-ca4e99 chr-I  1338889457 8     


Reads_sample

jhelp(Alignment.Reads_sample)

Reads_sample (input_fn, output_folder, output_prefix, n_reads, n_samples, rand_seed, verbose, quiet, progress, kwargs)

Randomly sample n_reads reads from a bam file and write downsampled files in n_samples bam files. If the input bam file is not indexed by read_id index_reads is automatically called.


  • input_fn (required) [str]

Path to the indexed bam file

  • output_folder (default: ./) [str]

Path to a folder where to write sample files

  • output_prefix (default: out) [str]

Path to a folder where to write sample files

  • n_reads (default: 1000) [int]

Number of randomly selected reads in each sample

  • n_samples (default: 1) [int]

Number of samples to generate files for

  • rand_seed (default: 42) [int]

Seed to use for the pseudo randon generator. For non deterministic behaviour set to 0

  • verbose (default: False) [bool]

  • quiet (default: False) [bool]

  • progress (default: False) [bool]

  • kwargs

Basic usage

Alignment.Reads_sample("./data/sample_1.bam", "./output/sample_reads/", output_prefix="1K", n_reads=1000, n_samples=3)
## Running Alignment Reads_sample ##
    Checking Bam and index file
    Load index
    Write sample reads
    Indexing output bam file
    Indexing output bam file
    Indexing output bam file

References_sample

jhelp(Alignment.References_sample)

References_sample (input_fn, output_fn, selected_reads_fn, frac_reads, min_reads_ref, rand_seed, sorting_threads, verbose, quiet, progress, kwargs)

Randomly sample reads per references according to a fraction od the reads mapped to this reference for a one or several files and write selected reads in a new bam file


  • input_fn (required) [list(str)]

Bam file path or directory containing bam files or list of files, or regex or list of regex. It is quite flexible. All files need to be sorted and aligned to the same reference file.

  • output_fn (default: out.bam) [str]

Path to the output bam file (sorted and indexed)

  • selected_reads_fn (default: select_ref.txt) [str]

Path to the output text file containing all the read id selected

  • frac_reads (default: 0.5) [int]

Fraction of reads mapped to sample for each reference

  • min_reads_ref (default: 30) [int]

Minimal read coverage per file and reference before sampling

  • rand_seed (default: 42) [int]

Seed to use for the pseudo randon generator. For non deterministic behaviour set to None

  • sorting_threads (default: 4) [int]

Number of threads to use for bam file sorting

  • verbose (default: False) [bool]

  • quiet (default: False) [bool]

  • progress (default: False) [bool]

  • kwargs

Basic usage

Alignment.References_sample (
    input_fn = "./data/sample_*.bam",
    output_fn = "./output/sample_References_sample.bam",
    selected_reads_fn = "./output/sample_References_sample_refid.txt",
    frac_reads = 0.25,
    min_reads_ref = 100,
    progress = True)
## Running Alignment Ref_sample ##
## Index files ##
    Indexing alignment file ./data/sample_2.bam
    Reading : 13678 Reads [00:00, 20143.58 Reads/s]
    Indexing alignment file ./data/sample_1.bam
    Reading : 13684 Reads [00:00, 23549.99 Reads/s]
    Raw read counts summary
     primary reads: 21,185
     secondary reads: 2,966
     unmapped reads: 2,815
     supplementary reads: 396
## Randomly pick reads per references ##
## Sample reads and write to output file ##
    Writing selected reads for bam file ./data/sample_2.bam
    Writing : 100%|██████████| 2656/2656 [00:02<00:00, 1164.71 Reads/s]
    Writing selected reads for bam file ./data/sample_1.bam
    Writing : 100%|██████████| 2653/2653 [00:02<00:00, 1193.08 Reads/s]
    Sort BAM File
    Index sorted BAM File
    Selected read counts summary
     valid reads: 21,185
     valid sampled reads: 5,309
     valid references: 17

Filter

jhelp(Alignment.Filter)

Filter (input_fn, output_fn, selected_reads_fn, skip_unmapped, skip_secondary, skip_supplementary, index_reads, orientation, min_read_len, min_align_len, min_mapq, min_freq_identity, select_ref_fn, exclude_ref_fn, verbose, quiet, progress, kwargs)


  • input_fn (required) [str]

Path to the bam file to filter

  • output_fn (required) [str]

Path to the write filtered bam file

  • selected_reads_fn (default: None) [str]

Optional file where to write ids of selected reads

  • skip_unmapped (default: False) [bool]

Filter out unmapped reads

  • skip_secondary (default: False) [bool]

Filter out secondary alignment

  • skip_supplementary (default: False) [bool]

Filter out supplementary alignment

  • index_reads (default: False) [bool]

Index bam file with both pysam and pybiotools reads_index

  • orientation (default: .) [str]

Orientation of alignment on reference genome {"+","-" ,"."}

  • min_read_len (default: 0) [int]

Minimal query read length (basecalled length)

  • min_align_len (default: 0) [int]

Minimal query alignment length on reference

  • min_mapq (default: 0) [int]

Minimal mapping quality score (mapq)

  • min_freq_identity (default: 0) [float]

Minimal frequency of alignment identity [0 to 1]

  • select_ref_fn (default: None) [str]

File containing a list of references on which the reads have to be mapped.

  • exclude_ref_fn (default: None) [str]

File containing a list of references on which the reads should not be mapped.

  • verbose (default: False) [bool]

  • quiet (default: False) [bool]

  • progress (default: False) [bool]

  • kwargs

Basic usage

Filter all non primary reads

Alignment.Filter(
    "./data/sample_1.bam", 
    "./output/sample_1_filter.bam",
    skip_unmapped = True,
    skip_supplementary = True,
    skip_secondary = True,
    progress=True,
    verbose=True)

head("./output/sample_1_filter.bam")
## Running Alignment Filter ##
    Checking input bam file
    [DEBUG]: Make output directories
    Parsing reads
    Reading : 13684 Reads [00:05, 2633.81 Reads/s]
    Read counts summary
     Reads retained
      primary: 10,584
      total: 10,584
     Reads discarded
      total: 3,100
      secondary: 1,496
      unmapped: 1,416
      supplementary: 188

4f20c07e-183a-4483-9313-220390 16 chr-I 2    60 328S11M5I14M1D10M4I30M2D52M1D6 * 0 0 GGCACACCCACACCACACCCACACCCACAC ##$'&$(.-,.+./(+,/844-,-;7005, NM:i:1988 ms:i:18286 AS:i:18594 nn:i:0 tp:A:P cm:i:921  s1:i:6929 s2:i:609  de:f:0.085  SA:Z:chr-VIII,562456,+,14512S1 
f5fc15c4-0bd5-4e56-9509-625574 16 chr-I 16   60 263S37M1D18M1I9M2I19M1D10M3D4M * 0 0 CTGGGACACACCACACCCACACCACACCCA "$#"%+*-.1)''+)(*,$+-462/..;75 NM:i:183  ms:i:2086  AS:i:2086  nn:i:0 tp:A:P cm:i:103  s1:i:758  s2:i:357  de:f:0.0852 SA:Z:chr-VIII,562456,+,1742S12 
4c144dc6-705a-4d0c-951f-02fd7f 16 chr-I 542  60 8S25M1I49M1D49M1I14M1I12M1D6M1 * 0 0 TTTGAGGTACGGCACTTGCCTCAGCGGTCT "$$)*'+(.,+&%(&(10&("'*%$..,1- NM:i:221  ms:i:1888  AS:i:1888  nn:i:0 tp:A:P cm:i:80   s1:i:588  s2:i:370  de:f:0.1074 
6449b796-5339-4d61-8cf4-14658f 0  chr-I 1007 60 116S17M1I51M1D21M1I3M1I41M1D26 * 0 0 TCAATTGTACTCGTTAGTTCAAGATGGTGT ''%$+'(*&,+%&$*%).)#%'&'-(1(*, NM:i:2273 ms:i:24152 AS:i:24491 nn:i:0 tp:A:P cm:i:1195 s1:i:9277 s2:i:1734 de:f:0.0799 
f0bea8a5-a7bc-455c-b0e8-abc6ff 0  chr-I 1331 60 130S14M2D5M1D10M3D5M1I19M3D2M3 * 0 0 TCAGTGATAGCCGTTCAGTTGGTGTTAGCA 25+/-%$$$"$%/.7&%#&&(/-2-:(++, NM:i:459  ms:i:4616  AS:i:4616  nn:i:0 tp:A:P cm:i:213  s1:i:1646 s2:i:247  de:f:0.0948 
f50af728-e016-4c74-b48f-fe7459 0  chr-I 1708 60 57S16M2D10M2D35M2D12M1D2M1I34M * 0 0 ATATTGCTTCGTTCCGGTTATCTGGAAACG #$#(&',<3:79++($"$#"#&$$#()'(% NM:i:99   ms:i:1754  AS:i:1754  nn:i:0 tp:A:P cm:i:95   s1:i:698  s2:i:383  de:f:0.0686 
3213cb66-ac07-468a-862e-4e5fd1 0  chr-I 4220 60 103S12M3D6M2D56M1I33M2D19M1I35 * 0 0 ATCGGTGGCGCTTCGTTCAATGGGTGTTTA %*)&)''),#%#)%*,2)$##%&%%&),%" NM:i:30   ms:i:744   AS:i:744   nn:i:0 tp:A:P cm:i:40   s1:i:310  s2:i:0    de:f:0.0495 
370dbed1-1e84-4b40-87ff-e6468b 16 chr-I 4327 60 19S21M1D33M2D32M1I2M2D19M2D39M * 0 0 TTTCCAATAGTTTAAGGTTACTACATGACA '''((,&#(''&)*,))*%0%-*298522* NM:i:449  ms:i:8048  AS:i:8048  nn:i:0 tp:A:P cm:i:459  s1:i:3247 s2:i:0    de:f:0.0594 
9e21a7d4-1b3e-4967-856e-ca4e99 0  chr-I 4506 48 106S20M1I3M1D7M1D9M2D11M1D4M1I * 0 0 ACCGTGCTGGTTCGGATCAGCAAAGAAGTT "%$#%$%#'%)(+%'"#""#%%)'%(*),7 NM:i:256  ms:i:336   AS:i:336   nn:i:0 tp:A:P cm:i:10   s1:i:73   s2:i:0    de:f:0.1846 
b6a48773-8ef9-4491-b554-40d18f 16 chr-I 4949 60 242S10M1I30M1D24M1D18M2I22M2D5 * 0 0 TAAACTTGGGATGCAGAGCCTAAAAATGTT "),(##"#""""$$'('$#%$'+8.*'&)- NM:i:254  ms:i:524   AS:i:524   nn:i:0 tp:A:P cm:i:15   s1:i:119  s2:i:0    de:f:0.1917 


Multi criteria filtering

Remove unmapped, short reads and alignments, reads mapped on the minus strand, low mapq and low identity

Alignment.Filter(
    "./data/sample_1.bam", 
    "./output/sample_1_filter.bam",
    skip_unmapped = True,
    min_read_len=300,
    min_align_len=300,
    orientation = "+",
    min_mapq = 10,
    min_freq_identity=0.8,
    progress=True,
    verbose=True)

head("./output/sample_1_filter.bam") 
## Running Alignment Filter ##
    Checking input bam file
    [DEBUG]: Make output directories
    Parsing reads
    Reading : 13684 Reads [00:02, 5320.77 Reads/s]
    Read counts summary
     Reads discarded
      total: 9,211
      wrong_orientation: 6,168
      unmapped: 1,416
      low_mapping_quality: 973
      low_identity: 514
      short_alignment: 99
      short_read: 41
     Reads retained
      total: 4,473
      primary: 4,422
      supplementary: 51

6449b796-5339-4d61-8cf4-14658f 0 chr-I 1007  60 116S17M1I51M1D21M1I3M1I41M1D26 * 0 0 TCAATTGTACTCGTTAGTTCAAGATGGTGT ''%$+'(*&,+%&$*%).)#%'&'-(1(*, NM:i:2273 ms:i:24152 AS:i:24491 nn:i:0 tp:A:P cm:i:1195 s1:i:9277 s2:i:1734 de:f:0.0799 
f0bea8a5-a7bc-455c-b0e8-abc6ff 0 chr-I 1331  60 130S14M2D5M1D10M3D5M1I19M3D2M3 * 0 0 TCAGTGATAGCCGTTCAGTTGGTGTTAGCA 25+/-%$$$"$%/.7&%#&&(/-2-:(++, NM:i:459  ms:i:4616  AS:i:4616  nn:i:0 tp:A:P cm:i:213  s1:i:1646 s2:i:247  de:f:0.0948 
f50af728-e016-4c74-b48f-fe7459 0 chr-I 1708  60 57S16M2D10M2D35M2D12M1D2M1I34M * 0 0 ATATTGCTTCGTTCCGGTTATCTGGAAACG #$#(&',<3:79++($"$#"#&$$#()'(% NM:i:99   ms:i:1754  AS:i:1754  nn:i:0 tp:A:P cm:i:95   s1:i:698  s2:i:383  de:f:0.0686 
3213cb66-ac07-468a-862e-4e5fd1 0 chr-I 4220  60 103S12M3D6M2D56M1I33M2D19M1I35 * 0 0 ATCGGTGGCGCTTCGTTCAATGGGTGTTTA %*)&)''),#%#)%*,2)$##%&%%&),%" NM:i:30   ms:i:744   AS:i:744   nn:i:0 tp:A:P cm:i:40   s1:i:310  s2:i:0    de:f:0.0495 
a93f0bb3-2ac6-4139-9eea-618d3d 0 chr-I 5549  60 9S9M1I14M1I4M3D3M1I9M1D27M1D5M * 0 0 ACAATGGTTTTTTATAGAGACATATTTTTT ""'"""""%&%&&-0'+''-+)+)2-;84, NM:i:209  ms:i:2678  AS:i:2678  nn:i:0 tp:A:P cm:i:127  s1:i:1031 s2:i:0    de:f:0.078  
68ef24e6-a3c8-4c7e-a10e-b9db60 0 chr-I 6611  60 92S38M2D12M1D25M1D22M1D17M2D4M * 0 0 TCAGTACTTCGTTCAGTAACCAAGAAAGTT &"##&$*)*%(+22$###"$#*(%*.-1/< NM:i:64   ms:i:1366  AS:i:1366  nn:i:0 tp:A:P cm:i:73   s1:i:526  s2:i:0    de:f:0.058  
4eb0b340-2c98-4d05-921e-25a33f 0 chr-I 7967  60 112S12M1I80M1D7M1D9M3D15M1I25M * 0 0 CTCAGTATGCTTCGTTCAGTCAGACTGGGT #$$%#'%-$12<7666.'$)*$&%$&).43 NM:i:87   ms:i:1128  AS:i:1128  nn:i:0 tp:A:P cm:i:56   s1:i:429  s2:i:0    de:f:0.0752 
1c657032-c062-4728-ae39-f455bd 0 chr-I 12969 60 98S6M1I12M1I25M1D14M3D16M1I92M * 0 0 ATCGGTATGCTTCGTTCAGGTGTTTAACCA ""%'&)(.$33<42%####/')+>=48/3) NM:i:40   ms:i:940   AS:i:940   nn:i:0 tp:A:P cm:i:55   s1:i:364  s2:i:129  de:f:0.0527 
3b2856a0-7f61-4829-ae09-d3ed6b 0 chr-I 15150 60 101S12M1I12M1I6M1D2M1D13M2D6M2 * 0 0 CTCATTGCTTCGTTTCGGTTGCCTCGATTC %,(%&&"%+)0021+###$#"&(,.55:>8 NM:i:112  ms:i:806   AS:i:806   nn:i:0 tp:A:P cm:i:23   s1:i:229  s2:i:136  de:f:0.1199 
54df109b-0ddf-48db-b7c0-8eecfe 0 chr-I 16913 60 46S13M2D63M1D43M2D30M2D8M2D16M * 0 0 TCGGTTTTCGACATCGTGAAACGCTTTCGC ""#""**-'&&&+-,'$$,/(++)-/+##% NM:i:52   ms:i:898   AS:i:898   nn:i:0 tp:A:P cm:i:41   s1:i:319  s2:i:145  de:f:0.0558 


Select specific reference

with open ("data/select_ref.txt", "w") as fp:
    for ref in ['chr-I', 'chr-II', 'chr-III', 'chr-IV', 'chr-V', 'chr-VI']:
        fp.write(f"{ref}\n")

Alignment.Filter(
    input_fn="./data/sample_1.bam",
    output_fn="./output/sample_1_filter.bam",
    select_ref_fn="data/select_ref.txt",
    index=True,
    progress=True,
    verbose=True)
## Running Alignment Filter ##
    Checking input bam file
    Loading selected reference file ids
    [DEBUG]: Found 6 references to select
    [DEBUG]: Make output directories
    Parsing reads
    Reading : 13684 Reads [00:02, 6754.57 Reads/s]
    Read counts summary
     Reads retained
      total: 4,426
      primary: 2,650
      unmapped: 1,416
      secondary: 322
      supplementary: 38
     Reads discarded
      invalid_reference: 9,258
      total: 9,258

To_fastq

jhelp(Alignment.To_fastq)

To_fastq (input_fn, output_r1_fn, output_r2_fn, ignore_paired_end, verbose, quiet, progress, kwargs)

Dump reads from an alignment file or set of alignment file(s) to a fastq or pair of fastq file(s). Only the primary alignment are kept and paired_end reads are assumed to be interleaved. Compatible with unmapped or unaligned alignment files as well as files without header.


  • input_fn (required) [list(str)]

Path (or list of paths) to input BAM/CRAM/SAM file(s)

  • output_r1_fn (required) [str]

Path to an output fastq file (for Read1 in paired_end mode of output_r2_fn is provided). Automatically gzipped if the .gz extension is found

  • output_r2_fn (default: None) [str]

Optional Path to an output fastq file. Automatically gzipped if the .gz extension is found

  • ignore_paired_end (default: False) [bool]

Ignore paired_end information and output everything in a single file.

  • verbose (default: False) [bool]

  • quiet (default: False) [bool]

  • progress (default: False) [bool]

  • kwargs

Single end read usage from bam files

Alignment.To_fastq(
    input_fn=["./data/sample_1.bam", "./data/sample_2.bam"],
    output_r1_fn="./output/sample_1-2_SE_from_bam.fastq.gz",
    verbose=True,
    progress=True)
## Running Alignment To_fastq ##
    [DEBUG]: Opening file ./output/sample_1-2_SE_from_bam.fastq.gz in writing mode
    Parsing reads
    Reading input file ./data/sample_1.bam
    Reading: 12000 Reads [00:16, 740.81 Reads/s] 
    [DEBUG]: Reached end of input file ./data/sample_1.bam
    Reading input file ./data/sample_2.bam
    Reading: 12000 Reads [00:19, 630.06 Reads/s] 
    [DEBUG]: Reached end of input file ./data/sample_2.bam
    [DEBUG]: Closing file:./output/sample_1-2_SE_from_bam.fastq.gz
    [DEBUG]: Sequences writen: 24000

Paired-end reads usage from unaligned CRAM files

Alignment.To_fastq(
    input_fn=["./data/sample_1_20k.cram", "./data/sample_2_20k.cram"],
    output_r1_fn="./output/sample_1-2_PE_from_CRAM_1.fastq.gz",
    output_r2_fn="./output/sample_1-2_PE_from_CRAM_2.fastq.gz",
    verbose=True,
    progress=True)
## Running Alignment To_fastq ##
    [DEBUG]: Opening file ./output/sample_1-2_PE_from_CRAM_1.fastq.gz in writing mode
    [DEBUG]: Opening file ./output/sample_1-2_PE_from_CRAM_2.fastq.gz in writing mode
    Parsing reads
    Reading input file ./data/sample_1_20k.cram
    Reading: 12000 Reads [00:03, 3318.62 Reads/s]
    [DEBUG]: Reached end of input file ./data/sample_1_20k.cram
    Reading input file ./data/sample_2_20k.cram
    Reading: 12000 Reads [00:03, 3385.10 Reads/s]
    [DEBUG]: Reached end of input file ./data/sample_2_20k.cram
    [DEBUG]: Closing file:./output/sample_1-2_PE_from_CRAM_1.fastq.gz
    [DEBUG]: Sequences writen: 24000
    [DEBUG]: Closing file:./output/sample_1-2_PE_from_CRAM_2.fastq.gz
    [DEBUG]: Sequences writen: 24000

Split

jhelp(Alignment.Split)

Split (input_fn, output_dir, n_files, output_fn_list, index, verbose, quiet, progress, kwargs)

Split reads in a bam file in N files. The input bam file has to be sorted by coordinates and indexed. The last file can contain a few extra reads.


  • input_fn (required) [str]

Path to the bam file to filter

  • output_dir (default: "") [str]

Path to the directory where to write split bam files. Files generated have the same basename as the source file and are suffixed with numbers starting from 0

  • n_files (default: 10) [int]

Number of file to split the original file into

  • output_fn_list (default: []) [list(str)]

As an alternative to output_dir and n_files one can instead give a list of output files. Reads will be automatically split between the files in the same order as given

  • index (default: False) [bool]

Index output BAM files

  • verbose (default: False) [bool]

  • quiet (default: False) [bool]

  • progress (default: False) [bool]

  • kwargs

Usage with number of output files to generate

Alignment.Split(
    input_fn="./data/sample_1.bam",
    output_dir="./output/split_bam",
    n_files= 4,
    verbose=True)

!ls -lh "./output/split_bam"
## Running Alignment Split ##
    Checking input bam file
    [DEBUG]: List of output files to generate:
    [DEBUG]: * ./output/split_bam/sample_1_0.bam
    [DEBUG]: * ./output/split_bam/sample_1_1.bam
    [DEBUG]: * ./output/split_bam/sample_1_2.bam
    [DEBUG]: * ./output/split_bam/sample_1_3.bam
    Parsing reads
    [DEBUG]: Counting reads
    [DEBUG]: Open ouput file './output/split_bam/sample_1_0.bam'
    [DEBUG]: Close output file './output/split_bam/sample_1_0.bam'
    [DEBUG]: Reads written: 3,421
    [DEBUG]: Open ouput file './output/split_bam/sample_1_1.bam'
    [DEBUG]: Close output file './output/split_bam/sample_1_1.bam'
    [DEBUG]: Reads written: 3,421
    [DEBUG]: Open ouput file './output/split_bam/sample_1_2.bam'
    [DEBUG]: Close output file './output/split_bam/sample_1_2.bam'
    [DEBUG]: Reads written: 3,421
    [DEBUG]: Open ouput file './output/split_bam/sample_1_3.bam'
    [DEBUG]: Reached end of input file
    [DEBUG]: Close output file './output/split_bam/sample_1_3.bam'
    [DEBUG]: Reads written: 3,421
    Read counts summary
     Reads from index: 13,684
     Reads writen: 13,684
     Reads per file: 3,421

total 38M
-rw-rw-r-- 1 aleg aleg  11M Jan 19 14:53 sample_1_0.bam
-rw-rw-r-- 1 aleg aleg  12M Jan 19 14:53 sample_1_1.bam
-rw-rw-r-- 1 aleg aleg 9.4M Jan 19 14:53 sample_1_2.bam
-rw-rw-r-- 1 aleg aleg 5.6M Jan 19 14:53 sample_1_3.bam

Usage with a predefined list of output files

Alignment.Split(
    input_fn="./data/sample_2.bam",
    output_fn_list=["./output/split_bam_2/f1.bam", "./output/split_bam_2/f4.bam", "./output/split_bam_2/f3.bam"],
    index=True,
    verbose=True)

!ls -lh "./output/split_bam_2"
## Running Alignment Split ##
    Checking input bam file
    [DEBUG]: List of output files to generate:
    [DEBUG]: * ./output/split_bam_2/f1.bam
    [DEBUG]: * ./output/split_bam_2/f4.bam
    [DEBUG]: * ./output/split_bam_2/f3.bam
    Parsing reads
    [DEBUG]: Counting reads
    [DEBUG]: Open ouput file './output/split_bam_2/f1.bam'
    [DEBUG]: Close output file './output/split_bam_2/f1.bam'
    [DEBUG]: Reads written: 4,559
    [DEBUG]: index output file './output/split_bam_2/f1.bam'
    [DEBUG]: Open ouput file './output/split_bam_2/f4.bam'
    [DEBUG]: Close output file './output/split_bam_2/f4.bam'
    [DEBUG]: Reads written: 4,559
    [DEBUG]: index output file './output/split_bam_2/f4.bam'
    [DEBUG]: Open ouput file './output/split_bam_2/f3.bam'
    [DEBUG]: Reached end of input file
    [DEBUG]: Close output file './output/split_bam_2/f3.bam'
    [DEBUG]: Reads written: 4,560
    [DEBUG]: index output file './output/split_bam_2/f3.bam'
    Read counts summary
     Reads from index: 13,678
     Reads writen: 13,678
     Reads per file: 4,559

total 55M
-rw-rw-r-- 1 aleg aleg  17M Jan 19 14:53 f1.bam
-rw-rw-r-- 1 aleg aleg  10K Jan 19 14:53 f1.bam.bai
-rw-rw-r-- 1 aleg aleg  12M Nov 11 21:31 f2.bam
-rw-rw-r-- 1 aleg aleg 7.4K Jan 30  2020 f2.bam.bai
-rw-rw-r-- 1 aleg aleg  11M Jan 19 14:53 f3.bam
-rw-rw-r-- 1 aleg aleg 5.2K Jan 19 14:53 f3.bam.bai
-rw-rw-r-- 1 aleg aleg  16M Jan 19 14:53 f4.bam
-rw-rw-r-- 1 aleg aleg 8.6K Jan 19 14:53 f4.bam.bai