[−][src]Trait rust_htslib::bam::Read
A trait for a BAM reader with a read method.
Required methods
pub fn read(&mut self, record: &mut Record) -> Option<Result<()>>[src]
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.
pub fn records(&mut self) -> Records<'_, Self>ⓘ[src]
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.
pub fn rc_records(&mut self) -> RcRecords<'_, Self>ⓘ[src]
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();
pub fn pileup(&mut self) -> Pileups<'_, Self>ⓘ[src]
Iterator over pileups.
pub fn htsfile(&self) -> *mut htsFile[src]
Return the htsFile struct
pub fn header(&self) -> &HeaderView[src]
Return the header.
pub fn set_thread_pool(&mut self, tpool: &ThreadPool) -> Result<()>[src]
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
pub fn seek(&mut self, offset: i64) -> Result<()>[src]
Seek to the given virtual offset in the file
pub fn tell(&self) -> i64[src]
Report the current virtual offset
pub fn set_threads(&mut self, n_threads: usize) -> Result<()>[src]
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.
Implementors
impl Read for IndexedReader[src]
pub fn read(&mut self, record: &mut Record) -> Option<Result<()>>[src]
pub fn records(&mut self) -> Records<'_, Self>ⓘ[src]
Iterator over the records of the fetched 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.
pub fn rc_records(&mut self) -> RcRecords<'_, Self>ⓘ[src]
pub fn pileup(&mut self) -> Pileups<'_, Self>ⓘ[src]
pub fn htsfile(&self) -> *mut htsFile[src]
pub fn header(&self) -> &HeaderView[src]
pub fn set_thread_pool(&mut self, tpool: &ThreadPool) -> Result<()>[src]
impl Read for Reader[src]
pub fn read(&mut self, record: &mut Record) -> Option<Result<()>>[src]
Read the next BAM record into the given Record.
Returns None if there are no more records.
This method is useful if you want to read records as fast as possible as the
Record can be reused. A more ergonomic approach is to use the records
iterator.
Errors
If there are any issues with reading the next record an error will be returned.
Examples
use rust_htslib::errors::Error; use rust_htslib::bam::{Read, Reader, Record}; let mut bam = Reader::from_path(&"test/test.bam")?; let mut record = Record::new(); // Print the TID of each record while let Some(r) = bam.read(&mut record) { r.expect("Failed to parse record"); println!("TID: {}", record.tid()) }
pub fn records(&mut self) -> Records<'_, Self>ⓘ[src]
Iterator over the records of the fetched 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.