rust_htslib/
lib.rs

1// Copyright 2014-2021 Johannes Köster.
2// Licensed under the MIT license (http://opensource.org/licenses/MIT)
3// This file may not be copied, modified, or distributed
4// except according to those terms.
5
6//! Rust-Htslib provides a high level API to working with the common HTS file formats.
7//!
8//! Htslib itself is the *de facto* standard implementation for reading and writing files for
9//! HTS alignments (SAM and BAM) as well as variant calls in VCF and BCF format.
10//!
11//! For example, let's say that we use samtools to view the header of a test file:
12//!
13//! ```bash
14//! samtools view -H test/test.bam
15//! @SQ    SN:CHROMOSOME_I    LN:15072423
16//! @SQ    SN:CHROMOSOME_II    LN:15279345
17//! @SQ    SN:CHROMOSOME_III    LN:13783700
18//! @SQ    SN:CHROMOSOME_IV    LN:17493793
19//! @SQ    SN:CHROMOSOME_V    LN:20924149
20//! ```
21//!
22//! We can reproduce that with Rust-Htslib. Reading BAM files and printing the header
23//! to the the screen is as easy as
24//!
25//! ```
26//! use rust_htslib::{bam, bam::Read};
27//!
28//!
29//! let bam = bam::Reader::from_path(&"test/test.bam").unwrap();
30//! let header = bam::Header::from_template(bam.header());
31//!
32//! // print header records to the terminal, akin to samtool
33//! for (key, records) in header.to_hashmap() {
34//!     for record in records {
35//!          println!("@{}\tSN:{}\tLN:{}", key, record["SN"], record["LN"]);
36//!     }
37//! }
38//! ```
39//!
40//! which results in the following output, equivalent to samtools.
41//!
42//! ```bash
43//! @SQ    SN:CHROMOSOME_I    LN:15072423
44//! @SQ    SN:CHROMOSOME_II    LN:15279345
45//! @SQ    SN:CHROMOSOME_III    LN:13783700
46//! @SQ    SN:CHROMOSOME_IV    LN:17493793
47//! @SQ    SN:CHROMOSOME_V    LN:20924149
48//! ```
49//!
50//! We can also read directly from the BAM file and write to an output file
51//!
52//! ```
53//! use rust_htslib::{bam, bam::Read};
54//!
55//! let mut bam = bam::Reader::from_path(&"test/test.bam").unwrap();
56//! let header = bam::Header::from_template(bam.header());
57//! let mut out = bam::Writer::from_path(&"test/out.bam", &header, bam::Format::Bam).unwrap();
58//!
59//! // copy reverse reads to new BAM file
60//! for r in bam.records() {
61//!     let record = r.unwrap();
62//!     if record.is_reverse() {
63//!         out.write(&record).unwrap();
64//!     }
65//! }
66//! ```
67//!
68//! Pileups can be performed with
69//!
70//! ```
71//! use rust_htslib::{bam, bam::Read};
72//!
73//! let mut bam = bam::Reader::from_path(&"test/test.bam").unwrap();
74//!
75//! // pileup over all covered sites
76//! for p in bam.pileup() {
77//!     let pileup = p.unwrap();
78//!     println!("{}:{} depth {}", pileup.tid(), pileup.pos(), pileup.depth());
79//!
80//!     for alignment in pileup.alignments() {
81//!         if !alignment.is_del() && !alignment.is_refskip() {
82//!             println!("Base {}", alignment.record().seq()[alignment.qpos().unwrap()]);
83//!         }
84//!         // mark indel start
85//!         match alignment.indel() {
86//!             bam::pileup::Indel::Ins(len) => println!("Insertion of length {} between this and next position.", len),
87//!             bam::pileup::Indel::Del(len) => println!("Deletion of length {} between this and next position.", len),
88//!             bam::pileup::Indel::None => ()
89//!         }
90//!     }
91//! }
92//! ```
93//!
94//! In both cases, indexed BAM files can be seeked for specific regions using [`fetch`](bam/struct.IndexedReader.html#method.fetch), constraining either the record iterator or the pileups:
95//!
96//! ```
97//! use rust_htslib::{bam, bam::Read};
98//!
99//! let mut bam = bam::IndexedReader::from_path(&"test/test.bam").unwrap();
100//!
101//! bam.fetch(("CHROMOSOME_I", 0, 20)).unwrap();
102//! // afterwards, read or pileup in this region
103//! ```
104//!
105//! See
106//! * [`fetch`](bam/struct.IndexedReader.html#method.fetch)
107//! * [`records`](bam/struct.IndexedReader.html#method.records)
108//! * [`read`](bam/struct.IndexedReader.html#method.read)
109//! * [`pileup`](bam/struct.IndexedReader.html#method.pileup)
110
111#[macro_use]
112extern crate custom_derive;
113
114#[macro_use]
115extern crate newtype_derive;
116
117#[cfg(feature = "serde_feature")]
118extern crate serde;
119
120#[cfg(test)] // <-- not needed in examples + integration tests
121#[macro_use]
122extern crate pretty_assertions;
123#[cfg(all(test, feature = "serde_feature"))]
124extern crate serde_json;
125
126pub mod bam;
127pub mod bcf;
128pub mod bgzf;
129pub mod errors;
130pub mod faidx;
131pub mod htslib;
132pub mod tbx;
133pub mod tpool;
134pub mod utils;