Trait rust_htslib::bam::Read
source · pub trait Read: Sized {
// Required methods
fn read(&mut self, record: &mut Record) -> Option<Result<()>>;
fn records(&mut self) -> Records<'_, Self> ⓘ;
fn rc_records(&mut self) -> RcRecords<'_, Self> ⓘ;
fn pileup(&mut self) -> Pileups<'_, Self> ⓘ;
fn htsfile(&self) -> *mut htsFile;
fn header(&self) -> &HeaderView;
fn set_thread_pool(&mut self, tpool: &ThreadPool) -> Result<()>;
// Provided methods
fn seek(&mut self, offset: i64) -> Result<()> { ... }
fn tell(&self) -> i64 { ... }
fn set_threads(&mut self, n_threads: usize) -> Result<()> { ... }
}
Expand description
A trait for a BAM reader with a read method.
Required Methods§
sourcefn read(&mut self, record: &mut Record) -> Option<Result<()>>
fn read(&mut self, record: &mut Record) -> Option<Result<()>>
Read next BAM record into given record. Use this method in combination with a single allocated record to avoid the reallocations occurring with the iterator.
Arguments
record
- the record to be filled
Returns
Some(Ok(())) if the record was read and None if no more records to read
Example:
use rust_htslib::errors::Error;
use rust_htslib::bam::{Read, IndexedReader, Record};
let mut bam = IndexedReader::from_path(&"test/test.bam").unwrap();
bam.fetch((0, 1000, 2000)); // reads on tid 0, from 1000bp to 2000bp
let mut record = Record::new();
while let Some(result) = bam.read(&mut record) {
match result {
Ok(_) => {
println!("Read sequence: {:?}", record.seq().as_bytes());
},
Err(_) => panic!("BAM parsing failed...")
}
}
Consider using rc_records
instead.
sourcefn records(&mut self) -> Records<'_, Self> ⓘ
fn records(&mut self) -> Records<'_, Self> ⓘ
Iterator over the records of the seeked region.
Note that, while being convenient, this is less efficient than pre-allocating a
Record
and reading into it with the read
method, since every iteration involves
the allocation of a new Record
.
This is about 10% slower than record in micro benchmarks.
Consider using rc_records
instead.
sourcefn rc_records(&mut self) -> RcRecords<'_, Self> ⓘ
fn rc_records(&mut self) -> RcRecords<'_, Self> ⓘ
Records iterator using an Rc to avoid allocating a Record each turn.
This is about 1% slower than the read
based API in micro benchmarks,
but has nicer ergonomics (and might not actually be slower in your applications).
Example:
use rust_htslib::errors::Error;
use rust_htslib::bam::{Read, Reader, Record};
use rust_htslib::htslib; // for BAM_F*
let mut bam = Reader::from_path(&"test/test.bam").unwrap();
for read in
bam.rc_records()
.map(|x| x.expect("Failure parsing Bam file"))
.filter(|read|
read.flags()
& (htslib::BAM_FUNMAP
| htslib::BAM_FSECONDARY
| htslib::BAM_FQCFAIL
| htslib::BAM_FDUP) as u16
== 0
)
.filter(|read| !read.is_reverse()) {
println!("Found a forward read: {:?}", read.qname());
}
//or to add the read qnames into a Vec
let collected: Vec<_> = bam.rc_records().map(|read| read.unwrap().qname().to_vec()).collect();
sourcefn header(&self) -> &HeaderView
fn header(&self) -> &HeaderView
Return the header.
sourcefn set_thread_pool(&mut self, tpool: &ThreadPool) -> Result<()>
fn set_thread_pool(&mut self, tpool: &ThreadPool) -> Result<()>
Use a shared thread-pool for writing. This permits controlling the total
thread count when multiple readers and writers are working simultaneously.
A thread pool can be created with crate::tpool::ThreadPool::new(n_threads)
Arguments
tpool
- thread pool to use for compression work.
Provided Methods§
sourcefn set_threads(&mut self, n_threads: usize) -> Result<()>
fn set_threads(&mut self, n_threads: usize) -> Result<()>
Activate multi-threaded BAM read support in htslib. This should permit faster reading of large BAM files.
Setting nthreads
to 0
does not change the current state. Note that it is not
possible to set the number of background threads below 1
once it has been set.
Arguments
n_threads
- number of extra background writer threads to use, must be> 0
.