pub struct Reader<R> { /* private fields */ }
Expand description

A BAM reader.

A BAM file is an encoded and compressed version of a SAM file. While a SAM file has a header and a list of records, a BAM is comprised of three parts:

  1. a SAM header,
  2. a list of reference sequences, and
  3. a list of encoded SAM records.

The reader reads records sequentially but can use virtual positions to seek to offsets from the start of a seekable stream.

Examples

use noodles_bam as bam;

let mut reader = File::open("sample.bam").map(bam::Reader::new)?;
reader.read_header()?;
reader.read_reference_sequences()?;

for result in reader.records() {
    let record = result?;
    println!("{:?}", record);
}

Implementations

Returns a reference to the underlying reader.

Examples
use noodles_bam as bam;
let data = [];
let reader = bam::Reader::from(&data[..]);
assert!(reader.get_ref().is_empty());

Returns a mutable reference to the underlying reader.

Examples
use noodles_bam as bam;
let data = [];
let mut reader = bam::Reader::from(&data[..]);
assert!(reader.get_mut().is_empty());

Returns the underlying reader.

Examples
use noodles_bam as bam;
let data = [];
let reader = bam::Reader::from(&data[..]);
assert!(reader.into_inner().is_empty());

Reads the raw SAM header.

The BAM magic number is also checked.

The position of the stream is expected to be at the start.

This returns the raw SAM header as a String. It can subsequently be parsed as a noodles_sam::Header.

Examples
use noodles_bam as bam;
let mut reader = File::open("sample.bam").map(bam::Reader::new)?;
let header = reader.read_header()?;

Reads the binary reference sequences after the SAM header.

This is not the same as the @SQ records in the SAM header. A BAM has a list of reference sequences containing name and length tuples after the SAM header and before the list of records.

The position of the stream is expected to be directly after the header.

This returns a reference sequence dictionary (noodles_sam::header::ReferenceSequences), which can be used to build a minimal noodles_sam::Header if the SAM header is empty.

Examples
use noodles_bam as bam;
let mut reader = File::open("sample.bam").map(bam::Reader::new)?;
reader.read_header()?;
let reference_sequences = reader.read_reference_sequences()?;

Reads a single record.

The record block size (bs) is read from the underlying stream and bs bytes are read into an internal buffer. This buffer is used to populate the given record.

The stream is expected to be directly after the reference sequences or at the start of another record.

It is more ergonomic to read records using an iterator (see Self::records and Self::query), but using this method directly allows the reuse of a single Record buffer.

If successful, the record block size is returned. If a block size of 0 is returned, the stream reached EOF.

Examples
use noodles_bam as bam;
use noodles_sam::alignment::Record;

let mut reader = File::open("sample.bam").map(bam::Reader::new)?;
reader.read_header()?;
reader.read_reference_sequences()?;

let mut record = Record::default();
reader.read_record(&mut record)?;

Reads a single record without eagerly decoding its fields.

The record block size (bs) is read from the underlying stream and bs bytes are read into the lazy record’s buffer. No fields are decoded, meaning the record is not necessarily valid. However, the structure of the byte stream is guaranteed to be record-like.

The stream is expected to be directly after the reference sequences or at the start of another record.

If successful, the record block size is returned. If a block size of 0 is returned, the stream reached EOF.

Examples
use noodles_bam as bam;

let mut reader = File::open("sample.bam").map(bam::Reader::new)?;
reader.read_header()?;
reader.read_reference_sequences()?;

let mut record = bam::lazy::Record::default();
reader.read_lazy_record(&mut record)?;

Returns an iterator over records starting from the current stream position.

The stream is expected to be directly after the reference sequences or at the start of another record.

Examples
use noodles_bam as bam;

let mut reader = File::open("sample.bam").map(bam::Reader::new)?;
reader.read_header()?;
reader.read_reference_sequences()?;

for result in reader.records() {
    let record = result?;
    println!("{:?}", record);
}

Returns an iterator over lazy records.

The stream is expected to be directly after the reference sequences or at the start of another record.

Examples
use noodles_bam as bam;

let mut reader = File::open("sample.bam").map(bam::Reader::new)?;
reader.read_header()?;
reader.read_reference_sequences()?;

for result in reader.lazy_records() {
    let record = result?;
    // ...
}

Creates a BAM reader.

The given reader must be a raw BGZF stream, as the underlying reader wraps it in a decoder.

Examples
use noodles_bam as bam;
let data = [];
let reader = bam::Reader::new(&data[..]);

Returns the current virtual position of the underlying BGZF reader.

Examples
use noodles_bam as bam;

let data = Vec::new();
let reader = bam::Reader::new(&data[..]);
let virtual_position = reader.virtual_position();

assert_eq!(virtual_position.compressed(), 0);
assert_eq!(virtual_position.uncompressed(), 0);

Seeks the underlying BGZF reader to the given virtual position.

Virtual positions typically come from the associated BAM index file.

Examples
use noodles_bam as bam;
use noodles_bgzf as bgzf;
let mut reader = bam::Reader::new(Cursor::new(Vec::new()));
let virtual_position = bgzf::VirtualPosition::default();
reader.seek(virtual_position)?;

Returns an iterator over records that intersect the given region.

Examples
use noodles_bam::{self as bam, bai};
use noodles_core::Region;
use noodles_sam as sam;

let mut reader = File::open("sample.bam").map(bam::Reader::new)?;
let header: sam::Header = reader.read_header()?.parse()?;

let reference_sequences = header.reference_sequences();
let index = bai::read("sample.bam.bai")?;
let region = "sq0:8-13".parse()?;
let query = reader.query(reference_sequences, &index, &region)?;

for result in query {
    let record = result?;
    println!("{:?}", record);
}

Returns an iterator of unmapped records after querying for the unmapped region.

Examples
use noodles_bam::{self as bam, bai};

let mut reader = File::open("sample.bam").map(bam::Reader::new)?;
let index = bai::read("sample.bam.bai")?;
let query = reader.query_unmapped(&index)?;

for result in query {
    let record = result?;
    println!("{:?}", record);
}

Trait Implementations

Reads a SAM header.

Returns an iterator over records.

Converts to this type from the input type.

Auto Trait Implementations

Blanket Implementations

Gets the TypeId of self. Read more

Immutably borrows from an owned value. Read more

Mutably borrows from an owned value. Read more

Converts to this type from the input type.

Returns the argument unchanged.

Instruments this type with the provided Span, returning an Instrumented wrapper. Read more

Instruments this type with the current Span, returning an Instrumented wrapper. Read more

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

The type returned in the event of a conversion error.

Performs the conversion.

The type returned in the event of a conversion error.

Performs the conversion.

Attaches the provided Subscriber to this type, returning a WithDispatch wrapper. Read more

Attaches the current default Subscriber to this type, returning a WithDispatch wrapper. Read more