Skip to main content

rust_htslib/bam/
mod.rs

1// Copyright 2014 Christopher Schröder, 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//! Module for working with SAM, BAM, and CRAM files.
7
8pub mod buffer;
9pub mod ext;
10pub mod header;
11pub mod index;
12pub mod pileup;
13pub mod record;
14
15#[cfg(feature = "serde_feature")]
16pub mod record_serde;
17
18use std::ffi;
19use std::os::raw::c_char;
20use std::path::Path;
21use std::rc::Rc;
22use std::slice;
23use std::str;
24use std::sync::Arc;
25
26use url::Url;
27
28use crate::errors::{Error, Result};
29use crate::htslib;
30use crate::tpool::ThreadPool;
31use crate::utils::path_as_bytes;
32
33pub use crate::bam::buffer::RecordBuffer;
34pub use crate::bam::header::Header;
35pub use crate::bam::record::Record;
36use hts_sys::{hts_fmt_option, sam_fields};
37use std::convert::{TryFrom, TryInto};
38use std::mem::MaybeUninit;
39
40/// # Safety
41///
42/// Implementation for `Read::set_threads` and `Writer::set_threads`.
43unsafe fn set_threads(htsfile: *mut htslib::htsFile, n_threads: usize) -> Result<()> {
44    assert!(n_threads != 0, "n_threads must be > 0");
45
46    if htslib::hts_set_threads(htsfile, n_threads as i32) != 0 {
47        Err(Error::SetThreads)
48    } else {
49        Ok(())
50    }
51}
52
53unsafe fn set_thread_pool(htsfile: *mut htslib::htsFile, tpool: &ThreadPool) -> Result<()> {
54    let mut b = tpool.handle.borrow_mut();
55
56    if htslib::hts_set_thread_pool(htsfile, &mut b.inner as *mut _) != 0 {
57        Err(Error::ThreadPool)
58    } else {
59        Ok(())
60    }
61}
62
63/// # Safety
64///
65/// Set the reference FAI index path in a `htslib::htsFile` struct for reading CRAM format.
66pub unsafe fn set_fai_filename<P: AsRef<Path>>(
67    htsfile: *mut htslib::htsFile,
68    fasta_path: P,
69) -> Result<()> {
70    let path = if let Some(ext) = fasta_path.as_ref().extension() {
71        fasta_path
72            .as_ref()
73            .with_extension(format!("{}.fai", ext.to_str().unwrap()))
74    } else {
75        fasta_path.as_ref().with_extension(".fai")
76    };
77    let p: &Path = path.as_ref();
78    let c_str = ffi::CString::new(p.to_str().unwrap().as_bytes()).unwrap();
79    if htslib::hts_set_fai_filename(htsfile, c_str.as_ptr()) == 0 {
80        Ok(())
81    } else {
82        Err(Error::BamInvalidReferencePath { path: p.to_owned() })
83    }
84}
85
86/// A trait for a BAM reader with a read method.
87pub trait Read: Sized {
88    /// Read next BAM record into given record.
89    /// Use this method in combination with a single allocated record to avoid the reallocations
90    /// occurring with the iterator.
91    ///
92    /// # Arguments
93    ///
94    /// * `record` - the record to be filled
95    ///
96    /// # Returns
97    ///
98    /// Some(Ok(())) if the record was read and None if no more records to read
99    ///
100    /// Example:
101    /// ```
102    /// use rust_htslib::errors::Error;
103    /// use rust_htslib::bam::{Read, IndexedReader, Record};
104    ///
105    /// let mut bam = IndexedReader::from_path("test/test.bam").unwrap();
106    /// bam.fetch((0, 1000, 2000)); // reads on tid 0, from 1000bp to 2000bp
107    /// let mut record = Record::new();
108    /// while let Some(result) = bam.read(&mut record) {
109    ///     match result {
110    ///         Ok(_) => {
111    ///             println!("Read sequence: {:?}", record.seq().as_bytes());
112    ///         }
113    ///         Err(_) => panic!("BAM parsing failed...")
114    ///     }
115    /// }
116    /// ```
117    ///
118    /// Consider using [`rc_records`](#tymethod.rc_records) instead.
119    fn read(&mut self, record: &mut record::Record) -> Option<Result<()>>;
120
121    /// Iterator over the records of the seeked region.
122    /// Note that, while being convenient, this is less efficient than pre-allocating a
123    /// `Record` and reading into it with the `read` method, since every iteration involves
124    /// the allocation of a new `Record`.
125    ///
126    /// This is about 10% slower than record in micro benchmarks.
127    ///
128    /// Consider using [`rc_records`](#tymethod.rc_records) instead.
129    fn records(&mut self) -> Records<'_, Self>;
130
131    /// Records iterator using an Rc to avoid allocating a Record each turn.
132    /// This is about 1% slower than the [`read`](#tymethod.read) based API in micro benchmarks,
133    /// but has nicer ergonomics (and might not actually be slower in your applications).
134    ///
135    /// Example:
136    /// ```
137    /// use rust_htslib::errors::Error;
138    /// use rust_htslib::bam::{Read, Reader, Record};
139    /// use rust_htslib::htslib; // for BAM_F*
140    /// let mut bam = Reader::from_path("test/test.bam").unwrap();
141    ///
142    /// for read in
143    ///     bam.rc_records()
144    ///     .map(|x| x.expect("Failure parsing Bam file"))
145    ///     .filter(|read|
146    ///         read.flags()
147    ///          & (htslib::BAM_FUNMAP
148    ///              | htslib::BAM_FSECONDARY
149    ///              | htslib::BAM_FQCFAIL
150    ///              | htslib::BAM_FDUP) as u16
151    ///          == 0
152    ///     )
153    ///     .filter(|read| !read.is_reverse()) {
154    ///     println!("Found a forward read: {:?}", read.qname());
155    /// }
156    ///
157    /// //or to add the read qnames into a Vec
158    /// let collected: Vec<_> = bam.rc_records().map(|read| read.unwrap().qname().to_vec()).collect();
159    ///
160    ///
161    /// ```
162    fn rc_records(&mut self) -> RcRecords<'_, Self>;
163
164    /// Iterator over pileups.
165    fn pileup(&mut self) -> pileup::Pileups<'_, Self>;
166
167    /// Return the htsFile struct
168    fn htsfile(&self) -> *mut htslib::htsFile;
169
170    /// Return the header.
171    fn header(&self) -> &HeaderView;
172
173    /// Seek to the given virtual offset in the file
174    fn seek(&mut self, offset: i64) -> Result<()> {
175        let htsfile = unsafe { self.htsfile().as_ref() }.expect("bug: null pointer to htsFile");
176        let ret = match htsfile.format.format {
177            htslib::htsExactFormat_cram => unsafe {
178                i64::from(htslib::cram_seek(
179                    htsfile.fp.cram,
180                    offset as libc::off_t,
181                    libc::SEEK_SET,
182                ))
183            },
184            _ => unsafe { htslib::bgzf_seek(htsfile.fp.bgzf, offset, libc::SEEK_SET) },
185        };
186
187        if ret == 0 {
188            Ok(())
189        } else {
190            Err(Error::FileSeek)
191        }
192    }
193
194    /// Report the current virtual offset
195    fn tell(&self) -> i64 {
196        // this reimplements the bgzf_tell macro
197        let htsfile = unsafe { self.htsfile().as_ref() }.expect("bug: null pointer to htsFile");
198        let bgzf = unsafe { *htsfile.fp.bgzf };
199        (bgzf.block_address << 16) | (i64::from(bgzf.block_offset) & 0xFFFF)
200    }
201
202    /// Activate multi-threaded BAM read support in htslib. This should permit faster
203    /// reading of large BAM files.
204    ///
205    /// Setting `nthreads` to `0` does not change the current state.  Note that it is not
206    /// possible to set the number of background threads below `1` once it has been set.
207    ///
208    /// # Arguments
209    ///
210    /// * `n_threads` - number of extra background writer threads to use, must be `> 0`.
211    fn set_threads(&mut self, n_threads: usize) -> Result<()> {
212        unsafe { set_threads(self.htsfile(), n_threads) }
213    }
214
215    /// Use a shared thread-pool for writing. This permits controlling the total
216    /// thread count when multiple readers and writers are working simultaneously.
217    /// A thread pool can be created with `crate::tpool::ThreadPool::new(n_threads)`
218    ///
219    /// # Arguments
220    ///
221    /// * `tpool` - thread pool to use for compression work.
222    fn set_thread_pool(&mut self, tpool: &ThreadPool) -> Result<()>;
223
224    /// If the underlying file is in CRAM format, allows modifying CRAM options.
225    /// Note that this method does *not* check that the underlying file actually is in CRAM format.
226    ///
227    /// # Examples
228    ///
229    /// Set the required fields to RNAME and FLAG,
230    /// potentially allowing htslib to skip over the rest,
231    /// resulting in faster iteration:
232    /// ```
233    /// use rust_htslib::bam::{Read, Reader};
234    /// use hts_sys;
235    /// let mut cram = Reader::from_path("test/test_cram.cram").unwrap();
236    /// cram.set_cram_options(hts_sys::hts_fmt_option_CRAM_OPT_REQUIRED_FIELDS,
237    ///             hts_sys::sam_fields_SAM_RNAME | hts_sys::sam_fields_SAM_FLAG).unwrap();
238    /// ```
239    fn set_cram_options(&mut self, fmt_opt: hts_fmt_option, fields: sam_fields) -> Result<()> {
240        unsafe {
241            if hts_sys::hts_set_opt(self.htsfile(), fmt_opt, fields) != 0 {
242                Err(Error::HtsSetOpt)
243            } else {
244                Ok(())
245            }
246        }
247    }
248}
249
250/// A BAM reader.
251#[derive(Debug)]
252pub struct Reader {
253    htsfile: *mut htslib::htsFile,
254    header: Arc<HeaderView>,
255    tpool: Option<ThreadPool>,
256}
257
258unsafe impl Send for Reader {}
259
260impl Reader {
261    /// Create a new Reader from path.
262    ///
263    /// # Arguments
264    ///
265    /// * `path` - the path to open.
266    pub fn from_path<P: AsRef<Path>>(path: P) -> Result<Self> {
267        Self::new(&path_as_bytes(path, true)?)
268    }
269
270    /// Create a new Reader from STDIN.
271    pub fn from_stdin() -> Result<Self> {
272        Self::new(b"-")
273    }
274
275    /// Create a new Reader from URL.
276    pub fn from_url(url: &Url) -> Result<Self> {
277        Self::new(url.as_str().as_bytes())
278    }
279
280    /// Create a new Reader.
281    ///
282    /// # Arguments
283    ///
284    /// * `path` - the path to open. Use "-" for stdin.
285    fn new(path: &[u8]) -> Result<Self> {
286        let htsfile = hts_open(path, b"r")?;
287
288        let header = unsafe { htslib::sam_hdr_read(htsfile) };
289        if header.is_null() {
290            return Err(Error::BamOpen {
291                target: String::from_utf8_lossy(path).to_string(),
292            });
293        }
294
295        // Invalidate the `text` representation of the header
296        unsafe {
297            let _ = htslib::sam_hdr_line_name(header, b"SQ".as_ptr().cast::<c_char>(), 0);
298        }
299
300        Ok(Reader {
301            htsfile,
302            header: Arc::new(HeaderView::new(header)),
303            tpool: None,
304        })
305    }
306
307    extern "C" fn pileup_read(
308        data: *mut ::std::os::raw::c_void,
309        record: *mut htslib::bam1_t,
310    ) -> i32 {
311        let mut _self = unsafe { (data as *mut Self).as_mut().unwrap() };
312        unsafe {
313            htslib::sam_read1(
314                _self.htsfile(),
315                _self.header().inner_ptr() as *mut hts_sys::sam_hdr_t,
316                record,
317            )
318        }
319    }
320
321    /// Iterator over the records between the (optional) virtual offsets `start` and `end`
322    ///
323    /// # Arguments
324    ///
325    /// * `start` - Optional starting virtual offset to seek to. Throws an error if it is not
326    ///   a valid virtual offset.
327    ///
328    /// * `end` - Read until the virtual offset is less than `end`
329    pub fn iter_chunk(&mut self, start: Option<i64>, end: Option<i64>) -> ChunkIterator<'_, Self> {
330        if let Some(pos) = start {
331            self.seek(pos)
332                .expect("Failed to seek to the starting position");
333        };
334
335        ChunkIterator { reader: self, end }
336    }
337
338    /// Set the reference path for reading CRAM files.
339    ///
340    /// # Arguments
341    ///
342    /// * `path` - path to the FASTA reference
343    pub fn set_reference<P: AsRef<Path>>(&mut self, path: P) -> Result<()> {
344        unsafe { set_fai_filename(self.htsfile, path) }
345    }
346}
347
348impl Read for Reader {
349    /// Read the next BAM record into the given `Record`.
350    /// Returns `None` if there are no more records.
351    ///
352    /// This method is useful if you want to read records as fast as possible as the
353    /// `Record` can be reused. A more ergonomic approach is to use the [records](Reader::records)
354    /// iterator.
355    ///
356    /// # Errors
357    /// If there are any issues with reading the next record an error will be returned.
358    ///
359    /// # Examples
360    ///
361    /// ```
362    /// use rust_htslib::errors::Error;
363    /// use rust_htslib::bam::{Read, Reader, Record};
364    ///
365    /// let mut bam = Reader::from_path("test/test.bam")?;
366    /// let mut record = Record::new();
367    ///
368    /// // Print the TID of each record
369    /// while let Some(r) = bam.read(&mut record) {
370    ///    r.expect("Failed to parse record");
371    ///    println!("TID: {}", record.tid())
372    /// }
373    /// # Ok::<(), Error>(())
374    /// ```
375    fn read(&mut self, record: &mut record::Record) -> Option<Result<()>> {
376        match unsafe {
377            htslib::sam_read1(
378                self.htsfile,
379                self.header().inner_ptr() as *mut hts_sys::sam_hdr_t,
380                record.inner_ptr_mut(),
381            )
382        } {
383            -1 => None,
384            -2 => Some(Err(Error::BamTruncatedRecord)),
385            -4 => Some(Err(Error::BamInvalidRecord)),
386            _ => {
387                record.set_header(Arc::clone(&self.header));
388
389                Some(Ok(()))
390            }
391        }
392    }
393
394    /// Iterator over the records of the fetched region.
395    /// Note that, while being convenient, this is less efficient than pre-allocating a
396    /// `Record` and reading into it with the `read` method, since every iteration involves
397    /// the allocation of a new `Record`.
398    fn records(&mut self) -> Records<'_, Self> {
399        Records { reader: self }
400    }
401
402    fn rc_records(&mut self) -> RcRecords<'_, Self> {
403        RcRecords {
404            reader: self,
405            record: Rc::new(record::Record::new()),
406        }
407    }
408
409    fn pileup(&mut self) -> pileup::Pileups<'_, Self> {
410        let _self = self as *const Self;
411        let itr = unsafe {
412            htslib::bam_plp_init(
413                Some(Reader::pileup_read),
414                _self as *mut ::std::os::raw::c_void,
415            )
416        };
417        pileup::Pileups::new(self, itr)
418    }
419
420    fn htsfile(&self) -> *mut htslib::htsFile {
421        self.htsfile
422    }
423
424    fn header(&self) -> &HeaderView {
425        &self.header
426    }
427
428    fn set_thread_pool(&mut self, tpool: &ThreadPool) -> Result<()> {
429        unsafe { set_thread_pool(self.htsfile(), tpool)? }
430        self.tpool = Some(tpool.clone());
431        Ok(())
432    }
433}
434
435impl Drop for Reader {
436    fn drop(&mut self) {
437        unsafe {
438            htslib::hts_close(self.htsfile);
439        }
440    }
441}
442
443/// Conversion type for start/stop coordinates
444/// only public because it's leaked by the conversions
445#[doc(hidden)]
446pub struct FetchCoordinate(i64);
447
448//the old sam spec
449impl From<i32> for FetchCoordinate {
450    fn from(coord: i32) -> FetchCoordinate {
451        FetchCoordinate(coord as i64)
452    }
453}
454
455// to support un-annotated literals (type interference fails on those)
456impl From<u32> for FetchCoordinate {
457    fn from(coord: u32) -> FetchCoordinate {
458        FetchCoordinate(coord as i64)
459    }
460}
461
462//the new sam spec
463impl From<i64> for FetchCoordinate {
464    fn from(coord: i64) -> FetchCoordinate {
465        FetchCoordinate(coord)
466    }
467}
468
469//what some of our header methods return
470impl From<u64> for FetchCoordinate {
471    fn from(coord: u64) -> FetchCoordinate {
472        FetchCoordinate(coord.try_into().expect("Coordinate exceeded 2^^63-1"))
473    }
474}
475
476/// Enum for [IndexdReader.fetch()](struct.IndexedReader.html#method.fetch) arguments.
477///
478/// tids may be converted From<>:
479/// * i32 (correct as per spec)
480/// * u32 (because of header.tid. Will panic if above 2^31-1).
481///
482///Coordinates may be (via FetchCoordinate)
483/// * i32 (as of the sam v1 spec)
484/// * i64 (as of the htslib 'large coordinate' extension (even though they are not supported in BAM)
485/// * u32 (because that's what rust literals will default to)
486/// * u64 (because of header.target_len(). Will panic if above 2^^63-1).
487#[derive(Debug)]
488pub enum FetchDefinition<'a> {
489    /// tid, start, stop,
490    Region(i32, i64, i64),
491    /// 'named-reference', start, stop tuple.
492    RegionString(&'a [u8], i64, i64),
493    ///complete reference. May be i32 or u32 (which panics if above 2^31-')
494    CompleteTid(i32),
495    ///complete reference by name (&[u8] or &str)
496    String(&'a [u8]),
497    /// Every read
498    All,
499    /// Only reads with the BAM flag BAM_FUNMAP (which might not be all reads with reference = -1)
500    Unmapped,
501}
502
503impl<'a, X: Into<FetchCoordinate>, Y: Into<FetchCoordinate>> From<(i32, X, Y)>
504    for FetchDefinition<'a>
505{
506    fn from(tup: (i32, X, Y)) -> FetchDefinition<'a> {
507        let start: FetchCoordinate = tup.1.into();
508        let stop: FetchCoordinate = tup.2.into();
509        FetchDefinition::Region(tup.0, start.0, stop.0)
510    }
511}
512
513impl<'a, X: Into<FetchCoordinate>, Y: Into<FetchCoordinate>> From<(u32, X, Y)>
514    for FetchDefinition<'a>
515{
516    fn from(tup: (u32, X, Y)) -> FetchDefinition<'a> {
517        let start: FetchCoordinate = tup.1.into();
518        let stop: FetchCoordinate = tup.2.into();
519        FetchDefinition::Region(
520            tup.0.try_into().expect("Tid exceeded 2^31-1"),
521            start.0,
522            stop.0,
523        )
524    }
525}
526
527//non tuple impls
528impl<'a> From<i32> for FetchDefinition<'a> {
529    fn from(tid: i32) -> FetchDefinition<'a> {
530        FetchDefinition::CompleteTid(tid)
531    }
532}
533
534impl<'a> From<u32> for FetchDefinition<'a> {
535    fn from(tid: u32) -> FetchDefinition<'a> {
536        let tid: i32 = tid.try_into().expect("tid exceeded 2^31-1");
537        FetchDefinition::CompleteTid(tid)
538    }
539}
540
541impl<'a> From<&'a str> for FetchDefinition<'a> {
542    fn from(s: &'a str) -> FetchDefinition<'a> {
543        FetchDefinition::String(s.as_bytes())
544    }
545}
546
547//also accept &[u8;n] literals
548impl<'a> From<&'a [u8]> for FetchDefinition<'a> {
549    fn from(s: &'a [u8]) -> FetchDefinition<'a> {
550        FetchDefinition::String(s)
551    }
552}
553
554//also accept &[u8;n] literals
555impl<'a, T: AsRef<[u8]>> From<&'a T> for FetchDefinition<'a> {
556    fn from(s: &'a T) -> FetchDefinition<'a> {
557        FetchDefinition::String(s.as_ref())
558    }
559}
560
561impl<'a, X: Into<FetchCoordinate>, Y: Into<FetchCoordinate>> From<(&'a str, X, Y)>
562    for FetchDefinition<'a>
563{
564    fn from(tup: (&'a str, X, Y)) -> FetchDefinition<'a> {
565        let start: FetchCoordinate = tup.1.into();
566        let stop: FetchCoordinate = tup.2.into();
567        FetchDefinition::RegionString(tup.0.as_bytes(), start.0, stop.0)
568    }
569}
570
571impl<'a, X: Into<FetchCoordinate>, Y: Into<FetchCoordinate>> From<(&'a [u8], X, Y)>
572    for FetchDefinition<'a>
573{
574    fn from(tup: (&'a [u8], X, Y)) -> FetchDefinition<'a> {
575        let start: FetchCoordinate = tup.1.into();
576        let stop: FetchCoordinate = tup.2.into();
577        FetchDefinition::RegionString(tup.0, start.0, stop.0)
578    }
579}
580
581//also accept &[u8;n] literals
582impl<'a, T: AsRef<[u8]>, X: Into<FetchCoordinate>, Y: Into<FetchCoordinate>> From<(&'a T, X, Y)>
583    for FetchDefinition<'a>
584{
585    fn from(tup: (&'a T, X, Y)) -> FetchDefinition<'a> {
586        let start: FetchCoordinate = tup.1.into();
587        let stop: FetchCoordinate = tup.2.into();
588        FetchDefinition::RegionString(tup.0.as_ref(), start.0, stop.0)
589    }
590}
591
592#[derive(Debug)]
593pub struct IndexedReader {
594    htsfile: *mut htslib::htsFile,
595    header: Arc<HeaderView>,
596    idx: Arc<IndexView>,
597    itr: Option<*mut htslib::hts_itr_t>,
598    tpool: Option<ThreadPool>,
599}
600
601unsafe impl Send for IndexedReader {}
602
603impl IndexedReader {
604    /// Create a new Reader from path.
605    ///
606    /// # Arguments
607    ///
608    /// * `path` - the path to open.
609    pub fn from_path<P: AsRef<Path>>(path: P) -> Result<Self> {
610        Self::new(&path_as_bytes(path, true)?)
611    }
612
613    pub fn from_path_and_index<P: AsRef<Path>>(path: P, index_path: P) -> Result<Self> {
614        Self::new_with_index_path(
615            &path_as_bytes(path, true)?,
616            &path_as_bytes(index_path, true)?,
617        )
618    }
619
620    pub fn from_url(url: &Url) -> Result<Self> {
621        Self::new(url.as_str().as_bytes())
622    }
623
624    /// Create a new Reader.
625    ///
626    /// # Arguments
627    ///
628    /// * `path` - the path. Use "-" for stdin.
629    fn new(path: &[u8]) -> Result<Self> {
630        let htsfile = hts_open(path, b"r")?;
631        let header = unsafe { htslib::sam_hdr_read(htsfile) };
632        let c_str = ffi::CString::new(path).unwrap();
633        let idx = unsafe { htslib::sam_index_load(htsfile, c_str.as_ptr()) };
634        if idx.is_null() {
635            Err(Error::BamInvalidIndex {
636                target: str::from_utf8(path).unwrap().to_owned(),
637            })
638        } else {
639            Ok(IndexedReader {
640                htsfile,
641                header: Arc::new(HeaderView::new(header)),
642                idx: Arc::new(IndexView::new(idx)),
643                itr: None,
644                tpool: None,
645            })
646        }
647    }
648    /// Create a new Reader.
649    ///
650    /// # Arguments
651    ///
652    /// * `path` - the path. Use "-" for stdin.
653    /// * `index_path` - the index path to use
654    fn new_with_index_path(path: &[u8], index_path: &[u8]) -> Result<Self> {
655        let htsfile = hts_open(path, b"r")?;
656        let header = unsafe { htslib::sam_hdr_read(htsfile) };
657        let c_str_path = ffi::CString::new(path).unwrap();
658        let c_str_index_path = ffi::CString::new(index_path).unwrap();
659        let idx = unsafe {
660            htslib::sam_index_load2(htsfile, c_str_path.as_ptr(), c_str_index_path.as_ptr())
661        };
662        if idx.is_null() {
663            Err(Error::BamInvalidIndex {
664                target: str::from_utf8(path).unwrap().to_owned(),
665            })
666        } else {
667            Ok(IndexedReader {
668                htsfile,
669                header: Arc::new(HeaderView::new(header)),
670                idx: Arc::new(IndexView::new(idx)),
671                itr: None,
672                tpool: None,
673            })
674        }
675    }
676
677    /// Define the region from which .read() or .records will retrieve reads.
678    ///
679    /// Both iterating (with [.records()](trait.Read.html#tymethod.records)) and looping without allocation (with [.read()](trait.Read.html#tymethod.read) are a two stage process:
680    /// 1. 'fetch' the region of interest
681    /// 2. iter/loop through the reads.
682    ///
683    /// Example:
684    /// ```
685    /// use rust_htslib::bam::{IndexedReader, Read};
686    /// let mut bam = IndexedReader::from_path("test/test.bam").unwrap();
687    /// bam.fetch(("chrX", 10000, 20000)); // coordinates 10000..20000 on reference named "chrX"
688    /// for read in bam.records() {
689    ///     println!("read name: {:?}", read.unwrap().qname());
690    /// }
691    /// ```
692    ///
693    /// The arguments may be anything that can be converted into a FetchDefinition
694    /// such as
695    ///
696    /// * fetch(tid: u32) -> fetch everything on this reference
697    /// * fetch(reference_name: &[u8] | &str) -> fetch everything on this reference
698    /// * fetch((tid: i32, start: i64, stop: i64)): -> fetch in this region on this tid
699    /// * fetch((reference_name: &[u8] | &str, start: i64, stop: i64) -> fetch in this region on this tid
700    /// * fetch(FetchDefinition::All) or fetch(".") -> Fetch overything
701    /// * fetch(FetchDefinition::Unmapped) or fetch("*") -> Fetch unmapped (as signified by the 'unmapped' flag in the BAM - might be unreliable with some aligners.
702    ///
703    /// The start / stop coordinates will take i64 (the correct type as of htslib's 'large
704    /// coordinates' expansion), i32, u32, and u64 (with a possible panic! if the coordinate
705    /// won't fit an i64).
706    ///
707    /// `start` and `stop` are zero-based. `start` is inclusive, `stop` is exclusive.
708    ///
709    /// This replaces the old fetch and fetch_str implementations.
710    pub fn fetch<'a, T: Into<FetchDefinition<'a>>>(&mut self, fetch_definition: T) -> Result<()> {
711        //this 'compile time redirect' safes us
712        //from monomorphing the 'meat' of the fetch function
713        self._inner_fetch(fetch_definition.into())
714    }
715
716    fn _inner_fetch(&mut self, fetch_definition: FetchDefinition) -> Result<()> {
717        match fetch_definition {
718            FetchDefinition::Region(tid, start, stop) => {
719                self._fetch_by_coord_tuple(tid, start, stop)
720            }
721            FetchDefinition::RegionString(s, start, stop) => {
722                let tid = self.header().tid(s);
723                match tid {
724                    Some(tid) => self._fetch_by_coord_tuple(tid as i32, start, stop),
725                    None => Err(Error::Fetch),
726                }
727            }
728            FetchDefinition::CompleteTid(tid) => {
729                let len = self.header().target_len(tid as u32);
730                match len {
731                    Some(len) => self._fetch_by_coord_tuple(tid, 0, len as i64),
732                    None => Err(Error::Fetch),
733                }
734            }
735            FetchDefinition::String(s) => {
736                // either a target-name or a samtools style definition
737                let tid = self.header().tid(s);
738                match tid {
739                    Some(tid) => {
740                        //'large position' spec says target len must will fit into an i64.
741                        let len: i64 = self.header.target_len(tid).unwrap().try_into().unwrap();
742                        self._fetch_by_coord_tuple(tid as i32, 0, len)
743                    }
744                    None => self._fetch_by_str(s),
745                }
746            }
747            FetchDefinition::All => self._fetch_by_str(b"."),
748            FetchDefinition::Unmapped => self._fetch_by_str(b"*"),
749        }
750    }
751
752    fn _fetch_by_coord_tuple(&mut self, tid: i32, beg: i64, end: i64) -> Result<()> {
753        if let Some(itr) = self.itr {
754            unsafe { htslib::hts_itr_destroy(itr) }
755        }
756        let itr = unsafe { htslib::sam_itr_queryi(self.index().inner_ptr(), tid, beg, end) };
757        if itr.is_null() {
758            self.itr = None;
759            Err(Error::Fetch)
760        } else {
761            self.itr = Some(itr);
762            Ok(())
763        }
764    }
765
766    fn _fetch_by_str(&mut self, region: &[u8]) -> Result<()> {
767        if let Some(itr) = self.itr {
768            unsafe { htslib::hts_itr_destroy(itr) }
769        }
770        let rstr = ffi::CString::new(region).unwrap();
771        let rptr = rstr.as_ptr();
772        let itr = unsafe {
773            htslib::sam_itr_querys(
774                self.index().inner_ptr(),
775                self.header().inner_ptr() as *mut hts_sys::sam_hdr_t,
776                rptr,
777            )
778        };
779        if itr.is_null() {
780            self.itr = None;
781            Err(Error::Fetch)
782        } else {
783            self.itr = Some(itr);
784            Ok(())
785        }
786    }
787
788    extern "C" fn pileup_read(
789        data: *mut ::std::os::raw::c_void,
790        record: *mut htslib::bam1_t,
791    ) -> i32 {
792        let _self = unsafe { (data as *mut Self).as_mut().unwrap() };
793        match _self.itr {
794            Some(itr) => itr_next(_self.htsfile, itr, record), // read fetched region
795            None => unsafe {
796                htslib::sam_read1(
797                    _self.htsfile,
798                    _self.header().inner_ptr() as *mut hts_sys::sam_hdr_t,
799                    record,
800                )
801            }, // ordinary reading
802        }
803    }
804
805    /// Set the reference path for reading CRAM files.
806    ///
807    /// # Arguments
808    ///
809    /// * `path` - path to the FASTA reference
810    pub fn set_reference<P: AsRef<Path>>(&mut self, path: P) -> Result<()> {
811        unsafe { set_fai_filename(self.htsfile, path) }
812    }
813
814    pub fn index(&self) -> &IndexView {
815        &self.idx
816    }
817
818    // Analogous to slow_idxstats in samtools, see
819    // https://github.com/samtools/samtools/blob/556c60fdff977c0e6cadc4c2581661f187098b4d/bam_index.c#L140-L199
820    unsafe fn slow_idxstats(&mut self) -> Result<Vec<(i64, u64, u64, u64)>> {
821        self.set_cram_options(
822            hts_sys::hts_fmt_option_CRAM_OPT_REQUIRED_FIELDS,
823            hts_sys::sam_fields_SAM_RNAME | hts_sys::sam_fields_SAM_FLAG,
824        )?;
825        let header = self.header();
826        let h = header.inner;
827        let mut ret;
828        let mut last_tid = -2;
829        let fp = self.htsfile();
830
831        let nref =
832            usize::try_from(hts_sys::sam_hdr_nref(h)).map_err(|_| Error::NoSequencesInReference)?;
833        if nref == 0 {
834            return Ok(vec![]);
835        }
836        let mut counts = vec![vec![0; 2]; nref + 1];
837        let mut bb: hts_sys::bam1_t = MaybeUninit::zeroed().assume_init();
838        let b = &mut bb as *mut hts_sys::bam1_t;
839        loop {
840            ret = hts_sys::sam_read1(fp, h, b);
841            if ret < 0 {
842                break;
843            }
844            let tid = (*b).core.tid;
845            if tid >= nref as i32 || tid < -1 {
846                return Err(Error::InvalidTid { tid });
847            }
848
849            if tid != last_tid {
850                if (last_tid >= -1) && (counts[tid as usize][0] + counts[tid as usize][1]) > 0 {
851                    return Err(Error::BamUnsorted);
852                }
853                last_tid = tid;
854            }
855
856            let idx = if ((*b).core.flag as u32 & hts_sys::BAM_FUNMAP) > 0 {
857                1
858            } else {
859                0
860            };
861            counts[(*b).core.tid as usize][idx] += 1;
862        }
863
864        if ret == -1 {
865            let res = (0..nref)
866                .map(|i| {
867                    (
868                        i as i64,
869                        header.target_len(i as u32).unwrap(),
870                        counts[i][0],
871                        counts[i][1],
872                    )
873                })
874                .chain([(-1, 0, counts[nref][0], counts[nref][1])])
875                .collect();
876            Ok(res)
877        } else {
878            Err(Error::SlowIdxStats)
879        }
880    }
881
882    /// Similar to samtools idxstats, this returns a vector of tuples
883    /// containing the target id, length, number of mapped reads, and number of unmapped reads.
884    /// The last entry in the vector corresponds to the unmapped reads for the entire file, with
885    /// the tid set to -1.
886    pub fn index_stats(&mut self) -> Result<Vec<(i64, u64, u64, u64)>> {
887        let header = self.header();
888        let index = self.index();
889        if index.inner_ptr().is_null() {
890            panic!("Index is null");
891        }
892        // the quick index stats method only works for BAM files, not SAM or CRAM
893        unsafe {
894            if (*self.htsfile()).format.format != htslib::htsExactFormat_bam {
895                return self.slow_idxstats();
896            }
897        }
898        Ok((0..header.target_count())
899            .map(|tid| {
900                let (mapped, unmapped) = index.number_mapped_unmapped(tid);
901                let tlen = header.target_len(tid).unwrap();
902                (tid as i64, tlen, mapped, unmapped)
903            })
904            .chain([(-1, 0, 0, index.number_unmapped())])
905            .collect::<_>())
906    }
907}
908
909#[derive(Debug)]
910pub struct IndexView {
911    inner: *mut hts_sys::hts_idx_t,
912}
913
914unsafe impl Send for IndexView {}
915unsafe impl Sync for IndexView {}
916
917impl IndexView {
918    fn new(hts_idx: *mut hts_sys::hts_idx_t) -> Self {
919        Self { inner: hts_idx }
920    }
921
922    #[inline]
923    pub fn inner(&self) -> &hts_sys::hts_idx_t {
924        unsafe { self.inner.as_ref().unwrap() }
925    }
926
927    #[inline]
928    // Pointer to inner hts_idx_t struct
929    pub fn inner_ptr(&self) -> *const hts_sys::hts_idx_t {
930        self.inner
931    }
932
933    #[inline]
934    pub fn inner_mut(&mut self) -> &mut hts_sys::hts_idx_t {
935        unsafe { self.inner.as_mut().unwrap() }
936    }
937
938    #[inline]
939    // Mutable pointer to hts_idx_t struct
940    pub fn inner_ptr_mut(&mut self) -> *mut hts_sys::hts_idx_t {
941        self.inner
942    }
943
944    /// Get the number of mapped and unmapped reads for a given target id
945    /// FIXME only valid for BAM, not SAM/CRAM
946    fn number_mapped_unmapped(&self, tid: u32) -> (u64, u64) {
947        let (mut mapped, mut unmapped) = (0, 0);
948        unsafe {
949            hts_sys::hts_idx_get_stat(self.inner, tid as i32, &mut mapped, &mut unmapped);
950        }
951        (mapped, unmapped)
952    }
953
954    /// Get the total number of unmapped reads in the file
955    /// FIXME only valid for BAM, not SAM/CRAM
956    fn number_unmapped(&self) -> u64 {
957        unsafe { hts_sys::hts_idx_get_n_no_coor(self.inner) }
958    }
959}
960
961impl Drop for IndexView {
962    fn drop(&mut self) {
963        unsafe {
964            htslib::hts_idx_destroy(self.inner);
965        }
966    }
967}
968
969impl Read for IndexedReader {
970    fn read(&mut self, record: &mut record::Record) -> Option<Result<()>> {
971        match self.itr {
972            Some(itr) => {
973                match itr_next(self.htsfile, itr, &mut record.inner as *mut htslib::bam1_t) {
974                    -1 => None,
975                    -2 => Some(Err(Error::BamTruncatedRecord)),
976                    -4 => Some(Err(Error::BamInvalidRecord)),
977                    _ => {
978                        record.set_header(Arc::clone(&self.header));
979
980                        Some(Ok(()))
981                    }
982                }
983            }
984            None => None,
985        }
986    }
987
988    /// Iterator over the records of the fetched region.
989    /// Note that, while being convenient, this is less efficient than pre-allocating a
990    /// `Record` and reading into it with the `read` method, since every iteration involves
991    /// the allocation of a new `Record`.
992    fn records(&mut self) -> Records<'_, Self> {
993        Records { reader: self }
994    }
995
996    fn rc_records(&mut self) -> RcRecords<'_, Self> {
997        RcRecords {
998            reader: self,
999            record: Rc::new(record::Record::new()),
1000        }
1001    }
1002
1003    fn pileup(&mut self) -> pileup::Pileups<'_, Self> {
1004        let _self = self as *const Self;
1005        let itr = unsafe {
1006            htslib::bam_plp_init(
1007                Some(IndexedReader::pileup_read),
1008                _self as *mut ::std::os::raw::c_void,
1009            )
1010        };
1011        pileup::Pileups::new(self, itr)
1012    }
1013
1014    fn htsfile(&self) -> *mut htslib::htsFile {
1015        self.htsfile
1016    }
1017
1018    fn header(&self) -> &HeaderView {
1019        &self.header
1020    }
1021
1022    fn set_thread_pool(&mut self, tpool: &ThreadPool) -> Result<()> {
1023        unsafe { set_thread_pool(self.htsfile(), tpool)? }
1024        self.tpool = Some(tpool.clone());
1025        Ok(())
1026    }
1027}
1028
1029impl Drop for IndexedReader {
1030    fn drop(&mut self) {
1031        unsafe {
1032            if self.itr.is_some() {
1033                htslib::hts_itr_destroy(self.itr.unwrap());
1034            }
1035            htslib::hts_close(self.htsfile);
1036        }
1037    }
1038}
1039
1040#[derive(Debug, Clone, Copy, PartialEq, Eq)]
1041pub enum Format {
1042    Sam,
1043    Bam,
1044    Cram,
1045}
1046
1047impl Format {
1048    fn write_mode(self) -> &'static [u8] {
1049        match self {
1050            Format::Sam => b"w",
1051            Format::Bam => b"wb",
1052            Format::Cram => b"wc",
1053        }
1054    }
1055}
1056
1057/// A BAM writer.
1058#[derive(Debug)]
1059pub struct Writer {
1060    f: *mut htslib::htsFile,
1061    header: Arc<HeaderView>,
1062    tpool: Option<ThreadPool>,
1063}
1064
1065unsafe impl Send for Writer {}
1066
1067impl Writer {
1068    /// Create a new SAM/BAM/CRAM file.
1069    ///
1070    /// # Arguments
1071    ///
1072    /// * `path` - the path.
1073    /// * `header` - header definition to use
1074    /// * `format` - the format to use (SAM/BAM/CRAM)
1075    pub fn from_path<P: AsRef<Path>>(
1076        path: P,
1077        header: &header::Header,
1078        format: Format,
1079    ) -> Result<Self> {
1080        Self::new(&path_as_bytes(path, false)?, format.write_mode(), header)
1081    }
1082
1083    /// Create a new SAM/BAM/CRAM file at STDOUT.
1084    ///
1085    /// # Arguments
1086    ///
1087    /// * `header` - header definition to use
1088    /// * `format` - the format to use (SAM/BAM/CRAM)
1089    pub fn from_stdout(header: &header::Header, format: Format) -> Result<Self> {
1090        Self::new(b"-", format.write_mode(), header)
1091    }
1092
1093    /// Create a new SAM/BAM/CRAM file.
1094    ///
1095    /// # Arguments
1096    ///
1097    /// * `path` - the path. Use "-" for stdout.
1098    /// * `mode` - write mode, refer to htslib::hts_open()
1099    /// * `header` - header definition to use
1100    fn new(path: &[u8], mode: &[u8], header: &header::Header) -> Result<Self> {
1101        let f = hts_open(path, mode)?;
1102
1103        // sam_hdr_parse does not populate the text and l_text fields of the header_record.
1104        // This causes non-SQ headers to be dropped in the output BAM file.
1105        // To avoid this, we copy the All header to a new C-string that is allocated with malloc,
1106        // and set this into header_record manually.
1107        let header_record = unsafe {
1108            let mut header_string = header.to_bytes();
1109            if !header_string.is_empty() && header_string[header_string.len() - 1] != b'\n' {
1110                header_string.push(b'\n');
1111            }
1112            let l_text = header_string.len();
1113            let text = ::libc::malloc(l_text + 1);
1114            libc::memset(text, 0, l_text + 1);
1115            libc::memcpy(
1116                text,
1117                header_string.as_ptr() as *const ::libc::c_void,
1118                header_string.len(),
1119            );
1120
1121            //println!("{}", str::from_utf8(&header_string).unwrap());
1122            let rec = htslib::sam_hdr_parse(l_text + 1, text as *const c_char);
1123
1124            (*rec).text = text as *mut c_char;
1125            (*rec).l_text = l_text;
1126            rec
1127        };
1128
1129        unsafe {
1130            htslib::sam_hdr_write(f, header_record);
1131        }
1132
1133        Ok(Writer {
1134            f,
1135            header: Arc::new(HeaderView::new(header_record)),
1136            tpool: None,
1137        })
1138    }
1139
1140    /// Activate multi-threaded BAM write support in htslib. This should permit faster
1141    /// writing of large BAM files.
1142    ///
1143    /// # Arguments
1144    ///
1145    /// * `n_threads` - number of extra background writer threads to use, must be `> 0`.
1146    pub fn set_threads(&mut self, n_threads: usize) -> Result<()> {
1147        unsafe { set_threads(self.f, n_threads) }
1148    }
1149
1150    /// Use a shared thread-pool for writing. This permits controlling the total
1151    /// thread count when multiple readers and writers are working simultaneously.
1152    /// A thread pool can be created with `crate::tpool::ThreadPool::new(n_threads)`
1153    ///
1154    /// # Arguments
1155    ///
1156    /// * `tpool` - thread pool to use for compression work.
1157    pub fn set_thread_pool(&mut self, tpool: &ThreadPool) -> Result<()> {
1158        unsafe { set_thread_pool(self.f, tpool)? }
1159        self.tpool = Some(tpool.clone());
1160        Ok(())
1161    }
1162
1163    /// Write record to BAM.
1164    ///
1165    /// # Arguments
1166    ///
1167    /// * `record` - the record to write
1168    pub fn write(&mut self, record: &record::Record) -> Result<()> {
1169        if unsafe { htslib::sam_write1(self.f, self.header.inner(), record.inner_ptr()) } == -1 {
1170            Err(Error::WriteRecord)
1171        } else {
1172            Ok(())
1173        }
1174    }
1175
1176    /// Return the header.
1177    pub fn header(&self) -> &HeaderView {
1178        &self.header
1179    }
1180
1181    /// Set the reference path for reading CRAM files.
1182    ///
1183    /// # Arguments
1184    ///
1185    /// * `path` - path to the FASTA reference
1186    pub fn set_reference<P: AsRef<Path>>(&mut self, path: P) -> Result<()> {
1187        unsafe { set_fai_filename(self.f, path) }
1188    }
1189
1190    /// Set the compression level for writing BAM/CRAM files.
1191    ///
1192    /// # Arguments
1193    ///
1194    /// * `compression_level` - `CompressionLevel` enum variant
1195    pub fn set_compression_level(&mut self, compression_level: CompressionLevel) -> Result<()> {
1196        let level = compression_level.convert()?;
1197        match unsafe {
1198            htslib::hts_set_opt(
1199                self.f,
1200                htslib::hts_fmt_option_HTS_OPT_COMPRESSION_LEVEL,
1201                level,
1202            )
1203        } {
1204            0 => Ok(()),
1205            _ => Err(Error::BamInvalidCompressionLevel { level }),
1206        }
1207    }
1208}
1209
1210/// Compression levels in BAM/CRAM files
1211///
1212/// * Uncompressed: No compression, zlib level 0
1213/// * Fastest: Lowest compression level, zlib level 1
1214/// * Maximum: Highest compression level, zlib level 9
1215/// * Level(i): Custom compression level in the range [0, 9]
1216#[derive(Debug, Clone, Copy)]
1217pub enum CompressionLevel {
1218    Uncompressed,
1219    Fastest,
1220    Maximum,
1221    Level(u32),
1222}
1223
1224impl CompressionLevel {
1225    // Convert and check the variants of the `CompressionLevel` enum to a numeric level
1226    fn convert(self) -> Result<u32> {
1227        match self {
1228            CompressionLevel::Uncompressed => Ok(0),
1229            CompressionLevel::Fastest => Ok(1),
1230            CompressionLevel::Maximum => Ok(9),
1231            CompressionLevel::Level(i @ 0..=9) => Ok(i),
1232            CompressionLevel::Level(i) => Err(Error::BamInvalidCompressionLevel { level: i }),
1233        }
1234    }
1235}
1236
1237impl Drop for Writer {
1238    fn drop(&mut self) {
1239        unsafe {
1240            htslib::hts_close(self.f);
1241        }
1242    }
1243}
1244
1245/// Iterator over the records of a BAM.
1246#[derive(Debug)]
1247pub struct Records<'a, R: Read> {
1248    reader: &'a mut R,
1249}
1250
1251impl<R: Read> Iterator for Records<'_, R> {
1252    type Item = Result<record::Record>;
1253
1254    fn next(&mut self) -> Option<Result<record::Record>> {
1255        let mut record = record::Record::new();
1256        match self.reader.read(&mut record) {
1257            None => None,
1258            Some(Ok(_)) => Some(Ok(record)),
1259            Some(Err(err)) => Some(Err(err)),
1260        }
1261    }
1262}
1263
1264/// Iterator over the records of a BAM, using an Rc.
1265///
1266/// See [rc_records](trait.Read.html#tymethod.rc_records).
1267#[derive(Debug)]
1268pub struct RcRecords<'a, R: Read> {
1269    reader: &'a mut R,
1270    record: Rc<record::Record>,
1271}
1272
1273impl<R: Read> Iterator for RcRecords<'_, R> {
1274    type Item = Result<Rc<record::Record>>;
1275
1276    fn next(&mut self) -> Option<Self::Item> {
1277        let record = match Rc::get_mut(&mut self.record) {
1278            //not make_mut, we don't need a clone
1279            Some(x) => x,
1280            None => {
1281                self.record = Rc::new(record::Record::new());
1282                Rc::get_mut(&mut self.record).unwrap()
1283            }
1284        };
1285
1286        match self.reader.read(record) {
1287            None => None,
1288            Some(Ok(_)) => Some(Ok(Rc::clone(&self.record))),
1289            Some(Err(err)) => Some(Err(err)),
1290        }
1291    }
1292}
1293
1294/// Iterator over the records of a BAM until the virtual offset is less than `end`
1295pub struct ChunkIterator<'a, R: Read> {
1296    reader: &'a mut R,
1297    end: Option<i64>,
1298}
1299
1300impl<R: Read> Iterator for ChunkIterator<'_, R> {
1301    type Item = Result<record::Record>;
1302    fn next(&mut self) -> Option<Result<record::Record>> {
1303        if let Some(pos) = self.end {
1304            if self.reader.tell() >= pos {
1305                return None;
1306            }
1307        }
1308        let mut record = record::Record::new();
1309        match self.reader.read(&mut record) {
1310            None => None,
1311            Some(Ok(_)) => Some(Ok(record)),
1312            Some(Err(err)) => Some(Err(err)),
1313        }
1314    }
1315}
1316
1317/// Wrapper for opening a BAM file.
1318fn hts_open(path: &[u8], mode: &[u8]) -> Result<*mut htslib::htsFile> {
1319    let cpath = ffi::CString::new(path).unwrap();
1320    let path = str::from_utf8(path).unwrap();
1321    let c_str = ffi::CString::new(mode).unwrap();
1322    let ret = unsafe { htslib::hts_open(cpath.as_ptr(), c_str.as_ptr()) };
1323    if ret.is_null() {
1324        Err(Error::BamOpen {
1325            target: path.to_owned(),
1326        })
1327    } else {
1328        if !mode.contains(&b'w') {
1329            unsafe {
1330                // Comparison against 'htsFormatCategory_sequence_data' doesn't handle text files correctly
1331                // hence the explicit checks against all supported exact formats
1332                if (*ret).format.format != htslib::htsExactFormat_sam
1333                    && (*ret).format.format != htslib::htsExactFormat_bam
1334                    && (*ret).format.format != htslib::htsExactFormat_cram
1335                {
1336                    return Err(Error::BamOpen {
1337                        target: path.to_owned(),
1338                    });
1339                }
1340            }
1341        }
1342        Ok(ret)
1343    }
1344}
1345
1346/// Wrapper for iterating an indexed BAM file.
1347fn itr_next(
1348    htsfile: *mut htslib::htsFile,
1349    itr: *mut htslib::hts_itr_t,
1350    record: *mut htslib::bam1_t,
1351) -> i32 {
1352    unsafe {
1353        htslib::hts_itr_next(
1354            (*htsfile).fp.bgzf,
1355            itr,
1356            record as *mut ::std::os::raw::c_void,
1357            htsfile as *mut ::std::os::raw::c_void,
1358        )
1359    }
1360}
1361
1362#[derive(Debug)]
1363pub struct HeaderView {
1364    inner: *mut htslib::bam_hdr_t,
1365}
1366
1367unsafe impl Send for HeaderView {}
1368unsafe impl Sync for HeaderView {}
1369
1370impl HeaderView {
1371    /// Create a new HeaderView from a pre-populated Header object
1372    pub fn from_header(header: &Header) -> Self {
1373        let mut header_string = header.to_bytes();
1374        if !header_string.is_empty() && header_string[header_string.len() - 1] != b'\n' {
1375            header_string.push(b'\n');
1376        }
1377        Self::from_bytes(&header_string)
1378    }
1379
1380    /// Create a new HeaderView from bytes
1381    pub fn from_bytes(header_string: &[u8]) -> Self {
1382        let header_record = unsafe {
1383            let l_text = header_string.len();
1384            let text = ::libc::malloc(l_text + 1);
1385            ::libc::memset(text, 0, l_text + 1);
1386            ::libc::memcpy(
1387                text,
1388                header_string.as_ptr() as *const ::libc::c_void,
1389                header_string.len(),
1390            );
1391
1392            let rec = htslib::sam_hdr_parse(l_text + 1, text as *const c_char);
1393            (*rec).text = text as *mut c_char;
1394            (*rec).l_text = l_text;
1395            rec
1396        };
1397
1398        HeaderView::new(header_record)
1399    }
1400
1401    /// Create a new HeaderView from the underlying Htslib type, and own it.
1402    fn new(inner: *mut htslib::bam_hdr_t) -> Self {
1403        HeaderView { inner }
1404    }
1405
1406    #[inline]
1407    pub fn inner(&self) -> &htslib::bam_hdr_t {
1408        unsafe { self.inner.as_ref().unwrap() }
1409    }
1410
1411    #[inline]
1412    // Pointer to inner bam_hdr_t struct
1413    pub fn inner_ptr(&self) -> *const htslib::bam_hdr_t {
1414        self.inner
1415    }
1416
1417    #[inline]
1418    pub fn inner_mut(&mut self) -> &mut htslib::bam_hdr_t {
1419        unsafe { self.inner.as_mut().unwrap() }
1420    }
1421
1422    #[inline]
1423    // Mutable pointer to bam_hdr_t struct
1424    pub fn inner_ptr_mut(&mut self) -> *mut htslib::bam_hdr_t {
1425        self.inner
1426    }
1427
1428    pub fn tid(&self, name: &[u8]) -> Option<u32> {
1429        let c_str = ffi::CString::new(name).expect("Expected valid name.");
1430        let tid = unsafe { htslib::sam_hdr_name2tid(self.inner, c_str.as_ptr()) };
1431        if tid < 0 {
1432            None
1433        } else {
1434            Some(tid as u32)
1435        }
1436    }
1437
1438    pub fn tid2name(&self, tid: u32) -> &[u8] {
1439        unsafe { ffi::CStr::from_ptr(htslib::sam_hdr_tid2name(self.inner, tid as i32)).to_bytes() }
1440    }
1441
1442    pub fn target_count(&self) -> u32 {
1443        self.inner().n_targets as u32
1444    }
1445
1446    pub fn target_names(&self) -> Vec<&[u8]> {
1447        let names = unsafe {
1448            slice::from_raw_parts(self.inner().target_name, self.target_count() as usize)
1449        };
1450        names
1451            .iter()
1452            .map(|name| unsafe { ffi::CStr::from_ptr(*name).to_bytes() })
1453            .collect()
1454    }
1455
1456    pub fn target_len(&self, tid: u32) -> Option<u64> {
1457        let inner = unsafe { *self.inner };
1458        if (tid as i32) < inner.n_targets {
1459            let l: &[u32] =
1460                unsafe { slice::from_raw_parts(inner.target_len, inner.n_targets as usize) };
1461            Some(l[tid as usize] as u64)
1462        } else {
1463            None
1464        }
1465    }
1466
1467    /// Retrieve the textual SAM header as bytes
1468    pub fn as_bytes(&self) -> &[u8] {
1469        unsafe {
1470            let rebuilt_hdr = htslib::sam_hdr_str(self.inner);
1471            if rebuilt_hdr.is_null() {
1472                return b"";
1473            }
1474            ffi::CStr::from_ptr(rebuilt_hdr).to_bytes()
1475        }
1476    }
1477}
1478
1479impl Clone for HeaderView {
1480    fn clone(&self) -> Self {
1481        HeaderView {
1482            inner: unsafe { htslib::sam_hdr_dup(self.inner) },
1483        }
1484    }
1485}
1486
1487impl Drop for HeaderView {
1488    fn drop(&mut self) {
1489        unsafe {
1490            htslib::sam_hdr_destroy(self.inner);
1491        }
1492    }
1493}
1494
1495#[cfg(test)]
1496mod tests {
1497    use super::header::HeaderRecord;
1498    use super::record::{Aux, Cigar, CigarString};
1499    use super::*;
1500    use std::collections::HashMap;
1501    use std::fs;
1502    use std::path::Path;
1503    use std::str;
1504
1505    type GoldType = (
1506        [&'static [u8]; 6],
1507        [u16; 6],
1508        [&'static [u8]; 6],
1509        [&'static [u8]; 6],
1510        [CigarString; 6],
1511    );
1512    fn gold() -> GoldType {
1513        let names = [
1514            &b"I"[..],
1515            &b"II.14978392"[..],
1516            &b"III"[..],
1517            &b"IV"[..],
1518            &b"V"[..],
1519            &b"VI"[..],
1520        ];
1521        let flags = [16u16, 16u16, 16u16, 16u16, 16u16, 2048u16];
1522        let seqs = [
1523            &b"CCTAGCCCTAACCCTAACCCTAACCCTAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCC\
1524TAAGCCTAAGCCTAAGCCTAA"[..],
1525            &b"CCTAGCCCTAACCCTAACCCTAACCCTAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCC\
1526TAAGCCTAAGCCTAAGCCTAA"[..],
1527            &b"CCTAGCCCTAACCCTAACCCTAACCCTAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCC\
1528TAAGCCTAAGCCTAAGCCTAA"[..],
1529            &b"CCTAGCCCTAACCCTAACCCTAACCCTAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCC\
1530TAAGCCTAAGCCTAAGCCTAA"[..],
1531            &b"CCTAGCCCTAACCCTAACCCTAACCCTAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCCTAAGCC\
1532TAAGCCTAAGCCTAAGCCTAA"[..],
1533            &b"ACTAAGCCTAAGCCTAAGCCTAAGCCAATTATCGATTTCTGAAAAAATTATCGAATTTTCTAGAAATTTTGCAAATTTT\
1534TTCATAAAATTATCGATTTTA"[..],
1535        ];
1536        let quals = [
1537            &b"#############################@B?8B?BA@@DDBCDDCBC@CDCDCCCCCCCCCCCCCCCCCCCCCCCCCCCC\
1538CCCCCCCCCCCCCCCCCCC"[..],
1539            &b"#############################@B?8B?BA@@DDBCDDCBC@CDCDCCCCCCCCCCCCCCCCCCCCCCCCCCCC\
1540CCCCCCCCCCCCCCCCCCC"[..],
1541            &b"#############################@B?8B?BA@@DDBCDDCBC@CDCDCCCCCCCCCCCCCCCCCCCCCCCCCCCC\
1542CCCCCCCCCCCCCCCCCCC"[..],
1543            &b"#############################@B?8B?BA@@DDBCDDCBC@CDCDCCCCCCCCCCCCCCCCCCCCCCCCCCCC\
1544CCCCCCCCCCCCCCCCCCC"[..],
1545            &b"#############################@B?8B?BA@@DDBCDDCBC@CDCDCCCCCCCCCCCCCCCCCCCCCCCCCCCC\
1546CCCCCCCCCCCCCCCCCCC"[..],
1547            &b"#############################@B?8B?BA@@DDBCDDCBC@CDCDCCCCCCCCCCCCCCCCCCCCCCCCCCCC\
1548CCCCCCCCCCCCCCCCCCC"[..],
1549        ];
1550        let cigars = [
1551            CigarString(vec![Cigar::Match(27), Cigar::Del(1), Cigar::Match(73)]),
1552            CigarString(vec![Cigar::Match(27), Cigar::Del(1), Cigar::Match(73)]),
1553            CigarString(vec![Cigar::Match(27), Cigar::Del(1), Cigar::Match(73)]),
1554            CigarString(vec![Cigar::Match(27), Cigar::Del(1), Cigar::Match(73)]),
1555            CigarString(vec![Cigar::Match(27), Cigar::Del(1), Cigar::Match(73)]),
1556            CigarString(vec![Cigar::Match(27), Cigar::Del(100000), Cigar::Match(73)]),
1557        ];
1558        (names, flags, seqs, quals, cigars)
1559    }
1560
1561    fn compare_inner_bam_cram_records(cram_records: &[Record], bam_records: &[Record]) {
1562        // Selectively compares bam1_t struct fields from BAM and CRAM
1563        for (c1, b1) in cram_records.iter().zip(bam_records.iter()) {
1564            // CRAM vs BAM l_data is off by 3, see: https://github.com/rust-bio/rust-htslib/pull/184#issuecomment-590133544
1565            // The rest of the fields should be identical:
1566            assert_eq!(c1.cigar(), b1.cigar());
1567            assert_eq!(c1.inner().core.pos, b1.inner().core.pos);
1568            assert_eq!(c1.inner().core.mpos, b1.inner().core.mpos);
1569            assert_eq!(c1.inner().core.mtid, b1.inner().core.mtid);
1570            assert_eq!(c1.inner().core.tid, b1.inner().core.tid);
1571            assert_eq!(c1.inner().core.bin, b1.inner().core.bin);
1572            assert_eq!(c1.inner().core.qual, b1.inner().core.qual);
1573            assert_eq!(c1.inner().core.l_extranul, b1.inner().core.l_extranul);
1574            assert_eq!(c1.inner().core.flag, b1.inner().core.flag);
1575            assert_eq!(c1.inner().core.l_qname, b1.inner().core.l_qname);
1576            assert_eq!(c1.inner().core.n_cigar, b1.inner().core.n_cigar);
1577            assert_eq!(c1.inner().core.l_qseq, b1.inner().core.l_qseq);
1578            assert_eq!(c1.inner().core.isize_, b1.inner().core.isize_);
1579            //... except m_data
1580        }
1581    }
1582
1583    #[test]
1584    fn test_read() {
1585        let (names, flags, seqs, quals, cigars) = gold();
1586        let mut bam = Reader::from_path(Path::new("test/test.bam")).expect("Error opening file.");
1587        let del_len = [1, 1, 1, 1, 1, 100000];
1588
1589        for (i, record) in bam.records().enumerate() {
1590            let rec = record.expect("Expected valid record");
1591            assert_eq!(rec.qname(), names[i]);
1592            assert_eq!(rec.flags(), flags[i]);
1593            assert_eq!(rec.seq().as_bytes(), seqs[i]);
1594
1595            let cigar = rec.cigar();
1596            assert_eq!(*cigar, cigars[i]);
1597
1598            let end_pos = cigar.end_pos();
1599            assert_eq!(end_pos, rec.pos() + 100 + del_len[i]);
1600            assert_eq!(
1601                cigar
1602                    .read_pos(end_pos as u32 - 10, false, false)
1603                    .unwrap()
1604                    .unwrap(),
1605                90
1606            );
1607            assert_eq!(
1608                cigar
1609                    .read_pos(rec.pos() as u32 + 20, false, false)
1610                    .unwrap()
1611                    .unwrap(),
1612                20
1613            );
1614            assert_eq!(cigar.read_pos(4000000, false, false).unwrap(), None);
1615            // fix qual offset
1616            let qual: Vec<u8> = quals[i].iter().map(|&q| q - 33).collect();
1617            assert_eq!(rec.qual(), &qual[..]);
1618        }
1619    }
1620
1621    #[test]
1622    fn test_seek() {
1623        let mut bam = Reader::from_path(Path::new("test/test.bam")).expect("Error opening file.");
1624
1625        let mut names_by_voffset = HashMap::new();
1626
1627        let mut offset = bam.tell();
1628        let mut rec = Record::new();
1629        while let Some(r) = bam.read(&mut rec) {
1630            r.expect("error reading bam");
1631            let qname = str::from_utf8(rec.qname()).unwrap().to_string();
1632            println!("{} {}", offset, qname);
1633            names_by_voffset.insert(offset, qname);
1634            offset = bam.tell();
1635        }
1636
1637        for (offset, qname) in names_by_voffset.iter() {
1638            println!("{} {}", offset, qname);
1639            bam.seek(*offset).unwrap();
1640            if let Some(r) = bam.read(&mut rec) {
1641                r.unwrap();
1642            };
1643            let rec_qname = str::from_utf8(rec.qname()).unwrap().to_string();
1644            assert_eq!(qname, &rec_qname);
1645        }
1646    }
1647
1648    #[test]
1649    fn test_read_sam_header() {
1650        let bam = Reader::from_path("test/test.bam").expect("Error opening file.");
1651
1652        let true_header = "@SQ\tSN:CHROMOSOME_I\tLN:15072423\n@SQ\tSN:CHROMOSOME_II\tLN:15279345\
1653             \n@SQ\tSN:CHROMOSOME_III\tLN:13783700\n@SQ\tSN:CHROMOSOME_IV\tLN:17493793\n@SQ\t\
1654             SN:CHROMOSOME_V\tLN:20924149\n"
1655            .to_string();
1656        let header_text = String::from_utf8(bam.header.as_bytes().to_owned()).unwrap();
1657        assert_eq!(header_text, true_header);
1658    }
1659
1660    #[test]
1661    fn test_read_against_sam() {
1662        let mut bam = Reader::from_path("./test/bam2sam_out.sam").unwrap();
1663        for read in bam.records() {
1664            let _read = read.unwrap();
1665        }
1666    }
1667
1668    fn _test_read_indexed_common(mut bam: IndexedReader) {
1669        let (names, flags, seqs, quals, cigars) = gold();
1670        let sq_1 = b"CHROMOSOME_I";
1671        let sq_2 = b"CHROMOSOME_II";
1672        let tid_1 = bam.header.tid(sq_1).expect("Expected tid.");
1673        let tid_2 = bam.header.tid(sq_2).expect("Expected tid.");
1674        assert!(bam.header.target_len(tid_1).expect("Expected target len.") == 15072423);
1675
1676        // fetch to position containing reads
1677        bam.fetch((tid_1, 0, 2))
1678            .expect("Expected successful fetch.");
1679        assert!(bam.records().count() == 6);
1680
1681        // compare reads
1682        bam.fetch((tid_1, 0, 2))
1683            .expect("Expected successful fetch.");
1684        for (i, record) in bam.records().enumerate() {
1685            let rec = record.expect("Expected valid record");
1686
1687            println!("{}", str::from_utf8(rec.qname()).unwrap());
1688            assert_eq!(rec.qname(), names[i]);
1689            assert_eq!(rec.flags(), flags[i]);
1690            assert_eq!(rec.seq().as_bytes(), seqs[i]);
1691            assert_eq!(*rec.cigar(), cigars[i]);
1692            // fix qual offset
1693            let qual: Vec<u8> = quals[i].iter().map(|&q| q - 33).collect();
1694            assert_eq!(rec.qual(), &qual[..]);
1695            assert_eq!(rec.aux(b"X"), Err(Error::BamAuxStringError));
1696            assert_eq!(rec.aux(b"NotAvailableAux"), Err(Error::BamAuxTagNotFound));
1697        }
1698
1699        // fetch to empty position
1700        bam.fetch((tid_2, 1, 1))
1701            .expect("Expected successful fetch.");
1702        assert!(bam.records().count() == 0);
1703
1704        // repeat with byte-string based fetch
1705
1706        // fetch to position containing reads
1707        // using coordinate-string chr:start-stop
1708        bam.fetch(format!("{}:{}-{}", str::from_utf8(sq_1).unwrap(), 0, 2).as_bytes())
1709            .expect("Expected successful fetch.");
1710        assert!(bam.records().count() == 6);
1711        // using &str and exercising some of the coordinate conversion funcs
1712        bam.fetch((str::from_utf8(sq_1).unwrap(), 0_u32, 2_u64))
1713            .expect("Expected successful fetch.");
1714        assert!(bam.records().count() == 6);
1715        // using a slice
1716        bam.fetch((&sq_1[..], 0, 2))
1717            .expect("Expected successful fetch.");
1718        assert!(bam.records().count() == 6);
1719        // using a literal
1720        bam.fetch((sq_1, 0, 2)).expect("Expected successful fetch.");
1721        assert!(bam.records().count() == 6);
1722
1723        // using a tid
1724        bam.fetch((0i32, 0u32, 2i64))
1725            .expect("Expected successful fetch.");
1726        assert!(bam.records().count() == 6);
1727        // using a tid:u32
1728        bam.fetch((0u32, 0u32, 2i64))
1729            .expect("Expected successful fetch.");
1730        assert!(bam.records().count() == 6);
1731
1732        // compare reads
1733        bam.fetch(format!("{}:{}-{}", str::from_utf8(sq_1).unwrap(), 0, 2).as_bytes())
1734            .expect("Expected successful fetch.");
1735        for (i, record) in bam.records().enumerate() {
1736            let rec = record.expect("Expected valid record");
1737
1738            println!("{}", str::from_utf8(rec.qname()).unwrap());
1739            assert_eq!(rec.qname(), names[i]);
1740            assert_eq!(rec.flags(), flags[i]);
1741            assert_eq!(rec.seq().as_bytes(), seqs[i]);
1742            assert_eq!(*rec.cigar(), cigars[i]);
1743            // fix qual offset
1744            let qual: Vec<u8> = quals[i].iter().map(|&q| q - 33).collect();
1745            assert_eq!(rec.qual(), &qual[..]);
1746            assert_eq!(rec.aux(b"NotAvailableAux"), Err(Error::BamAuxTagNotFound));
1747        }
1748
1749        // fetch to empty position
1750        bam.fetch(format!("{}:{}-{}", str::from_utf8(sq_2).unwrap(), 1, 1).as_bytes())
1751            .expect("Expected successful fetch.");
1752        assert!(bam.records().count() == 0);
1753
1754        //all on a tid
1755        bam.fetch(0).expect("Expected successful fetch.");
1756        assert!(bam.records().count() == 6);
1757        //all on a tid:u32
1758        bam.fetch(0u32).expect("Expected successful fetch.");
1759        assert!(bam.records().count() == 6);
1760
1761        //all on a tid - by &[u8]
1762        bam.fetch(sq_1).expect("Expected successful fetch.");
1763        assert!(bam.records().count() == 6);
1764        //all on a tid - by str
1765        bam.fetch(str::from_utf8(sq_1).unwrap())
1766            .expect("Expected successful fetch.");
1767        assert!(bam.records().count() == 6);
1768
1769        //all reads
1770        bam.fetch(FetchDefinition::All)
1771            .expect("Expected successful fetch.");
1772        assert!(bam.records().count() == 6);
1773
1774        //all reads
1775        bam.fetch(".").expect("Expected successful fetch.");
1776        assert!(bam.records().count() == 6);
1777
1778        //all unmapped
1779        bam.fetch(FetchDefinition::Unmapped)
1780            .expect("Expected successful fetch.");
1781        assert_eq!(bam.records().count(), 1); // expect one 'truncade record' Record.
1782
1783        bam.fetch("*").expect("Expected successful fetch.");
1784        assert_eq!(bam.records().count(), 1); // expect one 'truncade record' Record.
1785    }
1786
1787    #[test]
1788    fn test_read_indexed() {
1789        let bam = IndexedReader::from_path("test/test.bam").expect("Expected valid index.");
1790        _test_read_indexed_common(bam);
1791    }
1792
1793    #[test]
1794    fn test_read_indexed_different_index_name() {
1795        let bam = IndexedReader::from_path_and_index(
1796            &"test/test_different_index_name.bam",
1797            &"test/test.bam.bai",
1798        )
1799        .expect("Expected valid index.");
1800        _test_read_indexed_common(bam);
1801    }
1802
1803    #[test]
1804    fn test_set_record() {
1805        let (names, _, seqs, quals, cigars) = gold();
1806
1807        let mut rec = record::Record::new();
1808        rec.set_reverse();
1809        rec.set(names[0], Some(&cigars[0]), seqs[0], quals[0]);
1810        // note: this segfaults if you push_aux() before set()
1811        //       because set() obliterates aux
1812        rec.push_aux(b"NM", Aux::I32(15)).unwrap();
1813
1814        assert_eq!(rec.qname(), names[0]);
1815        assert_eq!(*rec.cigar(), cigars[0]);
1816        assert_eq!(rec.seq().as_bytes(), seqs[0]);
1817        assert_eq!(rec.qual(), quals[0]);
1818        assert!(rec.is_reverse());
1819        assert_eq!(rec.aux(b"NM").unwrap(), Aux::I32(15));
1820    }
1821
1822    #[test]
1823    fn test_set_repeated() {
1824        let mut rec = Record::new();
1825        rec.set(
1826            b"123",
1827            Some(&CigarString(vec![Cigar::Match(3)])),
1828            b"AAA",
1829            b"III",
1830        );
1831        rec.push_aux(b"AS", Aux::I32(12345)).unwrap();
1832        assert_eq!(rec.qname(), b"123");
1833        assert_eq!(rec.seq().as_bytes(), b"AAA");
1834        assert_eq!(rec.qual(), b"III");
1835        assert_eq!(rec.aux(b"AS").unwrap(), Aux::I32(12345));
1836
1837        rec.set(
1838            b"1234",
1839            Some(&CigarString(vec![Cigar::SoftClip(1), Cigar::Match(3)])),
1840            b"AAAA",
1841            b"IIII",
1842        );
1843        assert_eq!(rec.qname(), b"1234");
1844        assert_eq!(rec.seq().as_bytes(), b"AAAA");
1845        assert_eq!(rec.qual(), b"IIII");
1846        assert_eq!(rec.aux(b"AS").unwrap(), Aux::I32(12345));
1847
1848        rec.set(
1849            b"12",
1850            Some(&CigarString(vec![Cigar::Match(2)])),
1851            b"AA",
1852            b"II",
1853        );
1854        assert_eq!(rec.qname(), b"12");
1855        assert_eq!(rec.seq().as_bytes(), b"AA");
1856        assert_eq!(rec.qual(), b"II");
1857        assert_eq!(rec.aux(b"AS").unwrap(), Aux::I32(12345));
1858    }
1859
1860    #[test]
1861    fn test_set_qname() {
1862        let (names, _, seqs, quals, cigars) = gold();
1863
1864        assert!(names[0] != names[1]);
1865
1866        for i in 0..names.len() {
1867            let mut rec = record::Record::new();
1868            rec.set(names[i], Some(&cigars[i]), seqs[i], quals[i]);
1869            rec.push_aux(b"NM", Aux::I32(15)).unwrap();
1870
1871            assert_eq!(rec.qname(), names[i]);
1872            assert_eq!(*rec.cigar(), cigars[i]);
1873            assert_eq!(rec.seq().as_bytes(), seqs[i]);
1874            assert_eq!(rec.qual(), quals[i]);
1875            assert_eq!(rec.aux(b"NM").unwrap(), Aux::I32(15));
1876
1877            // Equal length qname
1878            assert!(rec.qname()[0] != b'X');
1879            rec.set_qname(b"X");
1880            assert_eq!(rec.qname(), b"X");
1881
1882            // Longer qname
1883            let mut longer_name = names[i].to_owned().clone();
1884            let extension = b"BuffaloBUffaloBUFFaloBUFFAloBUFFALoBUFFALO";
1885            longer_name.extend(extension.iter());
1886            rec.set_qname(&longer_name);
1887
1888            assert_eq!(rec.qname(), longer_name.as_slice());
1889            assert_eq!(*rec.cigar(), cigars[i]);
1890            assert_eq!(rec.seq().as_bytes(), seqs[i]);
1891            assert_eq!(rec.qual(), quals[i]);
1892            assert_eq!(rec.aux(b"NM").unwrap(), Aux::I32(15));
1893
1894            // Shorter qname
1895            let shorter_name = b"42";
1896            rec.set_qname(shorter_name);
1897
1898            assert_eq!(rec.qname(), shorter_name);
1899            assert_eq!(*rec.cigar(), cigars[i]);
1900            assert_eq!(rec.seq().as_bytes(), seqs[i]);
1901            assert_eq!(rec.qual(), quals[i]);
1902            assert_eq!(rec.aux(b"NM").unwrap(), Aux::I32(15));
1903
1904            // Zero-length qname
1905            rec.set_qname(b"");
1906
1907            assert_eq!(rec.qname(), b"");
1908            assert_eq!(*rec.cigar(), cigars[i]);
1909            assert_eq!(rec.seq().as_bytes(), seqs[i]);
1910            assert_eq!(rec.qual(), quals[i]);
1911            assert_eq!(rec.aux(b"NM").unwrap(), Aux::I32(15));
1912        }
1913    }
1914
1915    #[test]
1916    fn test_set_qname2() {
1917        let mut _header = Header::new();
1918        _header.push_record(
1919            HeaderRecord::new(b"SQ")
1920                .push_tag(b"SN", "1")
1921                .push_tag(b"LN", 10000000),
1922        );
1923        let header = HeaderView::from_header(&_header);
1924
1925        let line =
1926            b"blah1	0	1	1	255	1M	*	0	0	A	F	CB:Z:AAAA-1	UR:Z:AAAA	UB:Z:AAAA	GX:Z:G1	xf:i:1	fx:Z:G1\tli:i:0\ttf:Z:cC";
1927
1928        let mut rec = Record::from_sam(&header, line).unwrap();
1929        assert_eq!(rec.qname(), b"blah1");
1930        rec.set_qname(b"r0");
1931        assert_eq!(rec.qname(), b"r0");
1932    }
1933
1934    #[test]
1935    fn test_set_cigar() {
1936        let (names, _, seqs, quals, cigars) = gold();
1937
1938        assert!(names[0] != names[1]);
1939
1940        for i in 0..names.len() {
1941            let mut rec = record::Record::new();
1942            rec.set(names[i], Some(&cigars[i]), seqs[i], quals[i]);
1943            rec.push_aux(b"NM", Aux::I32(15)).unwrap();
1944
1945            assert_eq!(rec.qname(), names[i]);
1946            assert_eq!(*rec.cigar(), cigars[i]);
1947            assert_eq!(rec.seq().as_bytes(), seqs[i]);
1948            assert_eq!(rec.qual(), quals[i]);
1949            assert_eq!(rec.aux(b"NM").unwrap(), Aux::I32(15));
1950
1951            // boring cigar
1952            let new_cigar = CigarString(vec![Cigar::Match(rec.seq_len() as u32)]);
1953            assert_ne!(*rec.cigar(), new_cigar);
1954            rec.set_cigar(Some(&new_cigar));
1955            assert_eq!(*rec.cigar(), new_cigar);
1956
1957            assert_eq!(rec.qname(), names[i]);
1958            assert_eq!(rec.seq().as_bytes(), seqs[i]);
1959            assert_eq!(rec.qual(), quals[i]);
1960            assert_eq!(rec.aux(b"NM").unwrap(), Aux::I32(15));
1961
1962            // bizarre cigar
1963            let new_cigar = (0..rec.seq_len())
1964                .map(|i| {
1965                    if i % 2 == 0 {
1966                        Cigar::Match(1)
1967                    } else {
1968                        Cigar::Ins(1)
1969                    }
1970                })
1971                .collect::<Vec<_>>();
1972            let new_cigar = CigarString(new_cigar);
1973            assert_ne!(*rec.cigar(), new_cigar);
1974            rec.set_cigar(Some(&new_cigar));
1975            assert_eq!(*rec.cigar(), new_cigar);
1976
1977            assert_eq!(rec.qname(), names[i]);
1978            assert_eq!(rec.seq().as_bytes(), seqs[i]);
1979            assert_eq!(rec.qual(), quals[i]);
1980            assert_eq!(rec.aux(b"NM").unwrap(), Aux::I32(15));
1981
1982            // empty cigar
1983            let new_cigar = CigarString(Vec::new());
1984            assert_ne!(*rec.cigar(), new_cigar);
1985            rec.set_cigar(None);
1986            assert_eq!(*rec.cigar(), new_cigar);
1987
1988            assert_eq!(rec.qname(), names[i]);
1989            assert_eq!(rec.seq().as_bytes(), seqs[i]);
1990            assert_eq!(rec.qual(), quals[i]);
1991            assert_eq!(rec.aux(b"NM").unwrap(), Aux::I32(15));
1992        }
1993    }
1994
1995    #[test]
1996    fn test_remove_aux() {
1997        let mut bam = Reader::from_path(Path::new("test/test.bam")).expect("Error opening file.");
1998
1999        for record in bam.records() {
2000            let mut rec = record.expect("Expected valid record");
2001
2002            if rec.aux(b"XS").is_ok() {
2003                rec.remove_aux(b"XS").unwrap();
2004            }
2005
2006            if rec.aux(b"YT").is_ok() {
2007                rec.remove_aux(b"YT").unwrap();
2008            }
2009
2010            assert_eq!(rec.remove_aux(b"X"), Err(Error::BamAuxStringError));
2011            assert_eq!(rec.remove_aux(b"ab"), Err(Error::BamAuxTagNotFound));
2012
2013            assert_eq!(rec.aux(b"XS"), Err(Error::BamAuxTagNotFound));
2014            assert_eq!(rec.aux(b"YT"), Err(Error::BamAuxTagNotFound));
2015        }
2016    }
2017
2018    #[test]
2019    fn test_write() {
2020        let (names, _, seqs, quals, cigars) = gold();
2021
2022        let tmp = tempfile::Builder::new()
2023            .prefix("rust-htslib")
2024            .tempdir()
2025            .expect("Cannot create temp dir");
2026        let bampath = tmp.path().join("test.bam");
2027        println!("{:?}", bampath);
2028        {
2029            let mut bam = Writer::from_path(
2030                &bampath,
2031                Header::new().push_record(
2032                    HeaderRecord::new(b"SQ")
2033                        .push_tag(b"SN", "chr1")
2034                        .push_tag(b"LN", 15072423),
2035                ),
2036                Format::Bam,
2037            )
2038            .expect("Error opening file.");
2039
2040            for i in 0..names.len() {
2041                let mut rec = record::Record::new();
2042                rec.set(names[i], Some(&cigars[i]), seqs[i], quals[i]);
2043                rec.push_aux(b"NM", Aux::I32(15)).unwrap();
2044
2045                bam.write(&rec).expect("Failed to write record.");
2046            }
2047        }
2048
2049        {
2050            let mut bam = Reader::from_path(bampath).expect("Error opening file.");
2051
2052            for i in 0..names.len() {
2053                let mut rec = record::Record::new();
2054                if let Some(r) = bam.read(&mut rec) {
2055                    r.expect("Failed to read record.");
2056                };
2057
2058                assert_eq!(rec.qname(), names[i]);
2059                assert_eq!(*rec.cigar(), cigars[i]);
2060                assert_eq!(rec.seq().as_bytes(), seqs[i]);
2061                assert_eq!(rec.qual(), quals[i]);
2062                assert_eq!(rec.aux(b"NM").unwrap(), Aux::I32(15));
2063            }
2064        }
2065
2066        tmp.close().expect("Failed to delete temp dir");
2067    }
2068
2069    #[test]
2070    fn test_write_threaded() {
2071        let (names, _, seqs, quals, cigars) = gold();
2072
2073        let tmp = tempfile::Builder::new()
2074            .prefix("rust-htslib")
2075            .tempdir()
2076            .expect("Cannot create temp dir");
2077        let bampath = tmp.path().join("test.bam");
2078        println!("{:?}", bampath);
2079        {
2080            let mut bam = Writer::from_path(
2081                &bampath,
2082                Header::new().push_record(
2083                    HeaderRecord::new(b"SQ")
2084                        .push_tag(b"SN", "chr1")
2085                        .push_tag(b"LN", 15072423),
2086                ),
2087                Format::Bam,
2088            )
2089            .expect("Error opening file.");
2090            bam.set_threads(4).unwrap();
2091
2092            for i in 0..10000 {
2093                let mut rec = record::Record::new();
2094                let idx = i % names.len();
2095                rec.set(names[idx], Some(&cigars[idx]), seqs[idx], quals[idx]);
2096                rec.push_aux(b"NM", Aux::I32(15)).unwrap();
2097                rec.set_pos(i as i64);
2098
2099                bam.write(&rec).expect("Failed to write record.");
2100            }
2101        }
2102
2103        {
2104            let mut bam = Reader::from_path(bampath).expect("Error opening file.");
2105
2106            for (i, _rec) in bam.records().enumerate() {
2107                let idx = i % names.len();
2108
2109                let rec = _rec.expect("Failed to read record.");
2110
2111                assert_eq!(rec.pos(), i as i64);
2112                assert_eq!(rec.qname(), names[idx]);
2113                assert_eq!(*rec.cigar(), cigars[idx]);
2114                assert_eq!(rec.seq().as_bytes(), seqs[idx]);
2115                assert_eq!(rec.qual(), quals[idx]);
2116                assert_eq!(rec.aux(b"NM").unwrap(), Aux::I32(15));
2117            }
2118        }
2119
2120        tmp.close().expect("Failed to delete temp dir");
2121    }
2122
2123    #[test]
2124    fn test_write_shared_tpool() {
2125        let (names, _, seqs, quals, cigars) = gold();
2126
2127        let tmp = tempfile::Builder::new()
2128            .prefix("rust-htslib")
2129            .tempdir()
2130            .expect("Cannot create temp dir");
2131        let bampath1 = tmp.path().join("test1.bam");
2132        let bampath2 = tmp.path().join("test2.bam");
2133
2134        {
2135            let (mut bam1, mut bam2) = {
2136                let pool = crate::tpool::ThreadPool::new(4).unwrap();
2137
2138                let mut bam1 = Writer::from_path(
2139                    &bampath1,
2140                    Header::new().push_record(
2141                        HeaderRecord::new(b"SQ")
2142                            .push_tag(b"SN", "chr1")
2143                            .push_tag(b"LN", 15072423),
2144                    ),
2145                    Format::Bam,
2146                )
2147                .expect("Error opening file.");
2148
2149                let mut bam2 = Writer::from_path(
2150                    &bampath2,
2151                    Header::new().push_record(
2152                        HeaderRecord::new(b"SQ")
2153                            .push_tag(b"SN", "chr1")
2154                            .push_tag(b"LN", 15072423),
2155                    ),
2156                    Format::Bam,
2157                )
2158                .expect("Error opening file.");
2159
2160                bam1.set_thread_pool(&pool).unwrap();
2161                bam2.set_thread_pool(&pool).unwrap();
2162                (bam1, bam2)
2163            };
2164
2165            for i in 0..10000 {
2166                let mut rec = record::Record::new();
2167                let idx = i % names.len();
2168                rec.set(names[idx], Some(&cigars[idx]), seqs[idx], quals[idx]);
2169                rec.push_aux(b"NM", Aux::I32(15)).unwrap();
2170                rec.set_pos(i as i64);
2171
2172                bam1.write(&rec).expect("Failed to write record.");
2173                bam2.write(&rec).expect("Failed to write record.");
2174            }
2175        }
2176
2177        {
2178            let pool = crate::tpool::ThreadPool::new(2).unwrap();
2179
2180            for p in [bampath1, bampath2] {
2181                let mut bam = Reader::from_path(p).expect("Error opening file.");
2182                bam.set_thread_pool(&pool).unwrap();
2183
2184                for (i, _rec) in bam.iter_chunk(None, None).enumerate() {
2185                    let idx = i % names.len();
2186
2187                    let rec = _rec.expect("Failed to read record.");
2188
2189                    assert_eq!(rec.pos(), i as i64);
2190                    assert_eq!(rec.qname(), names[idx]);
2191                    assert_eq!(*rec.cigar(), cigars[idx]);
2192                    assert_eq!(rec.seq().as_bytes(), seqs[idx]);
2193                    assert_eq!(rec.qual(), quals[idx]);
2194                    assert_eq!(rec.aux(b"NM").unwrap(), Aux::I32(15));
2195                }
2196            }
2197        }
2198
2199        tmp.close().expect("Failed to delete temp dir");
2200    }
2201
2202    #[test]
2203    fn test_copy_template() {
2204        // Verify that BAM headers are transmitted correctly when using an existing BAM as a
2205        // template for headers.
2206
2207        let tmp = tempfile::Builder::new()
2208            .prefix("rust-htslib")
2209            .tempdir()
2210            .expect("Cannot create temp dir");
2211        let bampath = tmp.path().join("test.bam");
2212        println!("{:?}", bampath);
2213
2214        let mut input_bam = Reader::from_path("test/test.bam").expect("Error opening file.");
2215
2216        {
2217            let mut bam = Writer::from_path(
2218                &bampath,
2219                &Header::from_template(input_bam.header()),
2220                Format::Bam,
2221            )
2222            .expect("Error opening file.");
2223
2224            for rec in input_bam.records() {
2225                bam.write(&rec.unwrap()).expect("Failed to write record.");
2226            }
2227        }
2228
2229        {
2230            let copy_bam = Reader::from_path(bampath).expect("Error opening file.");
2231
2232            // Verify that the header came across correctly
2233            assert_eq!(input_bam.header().as_bytes(), copy_bam.header().as_bytes());
2234        }
2235
2236        tmp.close().expect("Failed to delete temp dir");
2237    }
2238
2239    #[test]
2240    fn test_pileup() {
2241        let (_, _, seqs, quals, _) = gold();
2242
2243        let mut bam = Reader::from_path("test/test.bam").expect("Error opening file.");
2244        let pileups = bam.pileup();
2245        for pileup in pileups.take(26) {
2246            let _pileup = pileup.expect("Expected successful pileup.");
2247            let pos = _pileup.pos() as usize;
2248            assert_eq!(_pileup.depth(), 6);
2249            assert!(_pileup.tid() == 0);
2250            for (i, a) in _pileup.alignments().enumerate() {
2251                assert_eq!(a.indel(), pileup::Indel::None);
2252                let qpos = a.qpos().unwrap();
2253                assert_eq!(qpos, pos - 1);
2254                assert_eq!(a.record().seq()[qpos], seqs[i][qpos]);
2255                assert_eq!(a.record().qual()[qpos], quals[i][qpos] - 33);
2256            }
2257        }
2258    }
2259
2260    #[test]
2261    fn test_idx_pileup() {
2262        let mut bam = IndexedReader::from_path("test/test.bam").expect("Error opening file.");
2263        // read without fetch
2264        for pileup in bam.pileup() {
2265            pileup.unwrap();
2266        }
2267        // go back again
2268        let tid = bam.header().tid(b"CHROMOSOME_I").unwrap();
2269        bam.fetch((tid, 0, 5)).unwrap();
2270        for p in bam.pileup() {
2271            println!("{}", p.unwrap().pos())
2272        }
2273    }
2274
2275    #[test]
2276    fn parse_from_sam() {
2277        use std::fs::File;
2278        use std::io::Read;
2279
2280        let bamfile = "./test/bam2sam_test.bam";
2281        let samfile = "./test/bam2sam_expected.sam";
2282
2283        // Load BAM file:
2284        let mut rdr = Reader::from_path(bamfile).unwrap();
2285        let bam_recs: Vec<Record> = rdr.records().map(|v| v.unwrap()).collect();
2286
2287        let mut sam = Vec::new();
2288        assert!(File::open(samfile).unwrap().read_to_end(&mut sam).is_ok());
2289
2290        let sam_recs: Vec<Record> = sam
2291            .split(|x| *x == b'\n')
2292            .filter(|x| !x.is_empty() && x[0] != b'@')
2293            .map(|line| Record::from_sam(rdr.header(), line).unwrap())
2294            .collect();
2295
2296        for (b1, s1) in bam_recs.iter().zip(sam_recs.iter()) {
2297            assert!(b1 == s1);
2298        }
2299    }
2300
2301    #[test]
2302    fn test_cigar_modes() {
2303        // test the cached and uncached ways of getting the cigar string.
2304
2305        let (_, _, _, _, cigars) = gold();
2306        let mut bam = Reader::from_path(Path::new("test/test.bam")).expect("Error opening file.");
2307
2308        for (i, record) in bam.records().enumerate() {
2309            let rec = record.expect("Expected valid record");
2310
2311            let cigar = rec.cigar();
2312            assert_eq!(*cigar, cigars[i]);
2313        }
2314
2315        for (i, record) in bam.records().enumerate() {
2316            let mut rec = record.expect("Expected valid record");
2317            rec.cache_cigar();
2318
2319            let cigar = rec.cigar_cached().unwrap();
2320            assert_eq!(**cigar, cigars[i]);
2321
2322            let cigar = rec.cigar();
2323            assert_eq!(*cigar, cigars[i]);
2324        }
2325    }
2326
2327    #[test]
2328    fn test_read_cram() {
2329        let cram_path = "./test/test_cram.cram";
2330        let bam_path = "./test/test_cram.bam";
2331        let ref_path = "./test/test_cram.fa";
2332
2333        // Load CRAM file, records
2334        let mut cram_reader = Reader::from_path(cram_path).unwrap();
2335        cram_reader.set_reference(ref_path).unwrap();
2336        let cram_records: Vec<Record> = cram_reader.records().map(|v| v.unwrap()).collect();
2337
2338        // Load BAM file, records
2339        let mut bam_reader = Reader::from_path(bam_path).unwrap();
2340        let bam_records: Vec<Record> = bam_reader.records().map(|v| v.unwrap()).collect();
2341
2342        compare_inner_bam_cram_records(&cram_records, &bam_records);
2343    }
2344
2345    #[test]
2346    fn test_write_cram() {
2347        // BAM file, records
2348        let bam_path = "./test/test_cram.bam";
2349        let ref_path = "./test/test_cram.fa";
2350        let mut bam_reader = Reader::from_path(bam_path).unwrap();
2351        let bam_records: Vec<Record> = bam_reader.records().map(|v| v.unwrap()).collect();
2352
2353        // New CRAM file
2354        let tmp = tempfile::Builder::new()
2355            .prefix("rust-htslib")
2356            .tempdir()
2357            .expect("Cannot create temp dir");
2358        let cram_path = tmp.path().join("test.cram");
2359
2360        // Write BAM records to new CRAM file
2361        {
2362            let mut header = Header::new();
2363            header.push_record(
2364                HeaderRecord::new(b"HD")
2365                    .push_tag(b"VN", "1.5")
2366                    .push_tag(b"SO", "coordinate"),
2367            );
2368            header.push_record(
2369                HeaderRecord::new(b"SQ")
2370                    .push_tag(b"SN", "chr1")
2371                    .push_tag(b"LN", 120)
2372                    .push_tag(b"M5", "20a9a0fb770814e6c5e49946750f9724")
2373                    .push_tag(b"UR", "test/test_cram.fa"),
2374            );
2375            header.push_record(
2376                HeaderRecord::new(b"SQ")
2377                    .push_tag(b"SN", "chr2")
2378                    .push_tag(b"LN", 120)
2379                    .push_tag(b"M5", "7a2006ccca94ea92b6dae5997e1b0d70")
2380                    .push_tag(b"UR", "test/test_cram.fa"),
2381            );
2382            header.push_record(
2383                HeaderRecord::new(b"SQ")
2384                    .push_tag(b"SN", "chr3")
2385                    .push_tag(b"LN", 120)
2386                    .push_tag(b"M5", "a66b336bfe3ee8801c744c9545c87e24")
2387                    .push_tag(b"UR", "test/test_cram.fa"),
2388            );
2389
2390            let mut cram_writer = Writer::from_path(&cram_path, &header, Format::Cram)
2391                .expect("Error opening CRAM file.");
2392            cram_writer.set_reference(ref_path).unwrap();
2393
2394            // Write BAM records to CRAM file
2395            for rec in bam_records.iter() {
2396                cram_writer
2397                    .write(rec)
2398                    .expect("Faied to write record to CRAM.");
2399            }
2400        }
2401
2402        // Compare written CRAM records with BAM records
2403        {
2404            // Load written CRAM file
2405            let mut cram_reader = Reader::from_path(cram_path).unwrap();
2406            cram_reader.set_reference(ref_path).unwrap();
2407            let cram_records: Vec<Record> = cram_reader.records().map(|v| v.unwrap()).collect();
2408
2409            // Compare CRAM records to BAM records
2410            compare_inner_bam_cram_records(&cram_records, &bam_records);
2411        }
2412
2413        tmp.close().expect("Failed to delete temp dir");
2414    }
2415
2416    #[test]
2417    fn test_compression_level_conversion() {
2418        // predefined compression levels
2419        assert_eq!(CompressionLevel::Uncompressed.convert().unwrap(), 0);
2420        assert_eq!(CompressionLevel::Fastest.convert().unwrap(), 1);
2421        assert_eq!(CompressionLevel::Maximum.convert().unwrap(), 9);
2422
2423        // numeric compression levels
2424        for level in 0..=9 {
2425            assert_eq!(CompressionLevel::Level(level).convert().unwrap(), level);
2426        }
2427        // invalid levels
2428        assert!(CompressionLevel::Level(10).convert().is_err());
2429    }
2430
2431    #[test]
2432    fn test_write_compression() {
2433        let tmp = tempfile::Builder::new()
2434            .prefix("rust-htslib")
2435            .tempdir()
2436            .expect("Cannot create temp dir");
2437        let input_bam_path = "test/test.bam";
2438
2439        // test levels with decreasing compression factor
2440        let levels_to_test = vec![
2441            CompressionLevel::Maximum,
2442            CompressionLevel::Level(6),
2443            CompressionLevel::Fastest,
2444            CompressionLevel::Uncompressed,
2445        ];
2446        let file_sizes: Vec<_> = levels_to_test
2447            .iter()
2448            .map(|level| {
2449                let output_bam_path = tmp.path().join("test.bam");
2450                {
2451                    let mut reader = Reader::from_path(input_bam_path).unwrap();
2452                    let header = Header::from_template(reader.header());
2453                    let mut writer =
2454                        Writer::from_path(&output_bam_path, &header, Format::Bam).unwrap();
2455                    writer.set_compression_level(*level).unwrap();
2456                    for record in reader.records() {
2457                        let r = record.unwrap();
2458                        writer.write(&r).unwrap();
2459                    }
2460                }
2461                fs::metadata(output_bam_path).unwrap().len()
2462            })
2463            .collect();
2464
2465        // check that out BAM file sizes are in decreasing order, in line with the expected compression factor
2466        println!("testing compression leves: {:?}", levels_to_test);
2467        println!("got compressed sizes: {:?}", file_sizes);
2468
2469        // libdeflate comes out with a slightly bigger file on Max compression
2470        // than on Level(6), so skip that check
2471        #[cfg(feature = "libdeflate")]
2472        assert!(file_sizes[1..].windows(2).all(|size| size[0] <= size[1]));
2473
2474        #[cfg(not(feature = "libdeflate"))]
2475        assert!(file_sizes.windows(2).all(|size| size[0] <= size[1]));
2476
2477        tmp.close().expect("Failed to delete temp dir");
2478    }
2479
2480    #[test]
2481    fn test_bam_fails_on_vcf() {
2482        let bam_path = "./test/test_left.vcf";
2483        let bam_reader = Reader::from_path(bam_path);
2484        assert!(bam_reader.is_err());
2485    }
2486
2487    #[test]
2488    fn test_indexde_bam_fails_on_vcf() {
2489        let bam_path = "./test/test_left.vcf";
2490        let bam_reader = IndexedReader::from_path(bam_path);
2491        assert!(bam_reader.is_err());
2492    }
2493
2494    #[test]
2495    fn test_bam_fails_on_toml() {
2496        let bam_path = "./Cargo.toml";
2497        let bam_reader = Reader::from_path(bam_path);
2498        assert!(bam_reader.is_err());
2499    }
2500
2501    #[test]
2502    fn test_sam_writer_example() {
2503        fn from_bam_with_filter<F>(bamfile: &str, samfile: &str, f: F) -> bool
2504        where
2505            F: Fn(&record::Record) -> Option<bool>,
2506        {
2507            let mut bam_reader = Reader::from_path(bamfile).unwrap(); // internal functions, just unwrap
2508            let header = header::Header::from_template(bam_reader.header());
2509            let mut sam_writer = Writer::from_path(samfile, &header, Format::Sam).unwrap();
2510            for record in bam_reader.records() {
2511                if record.is_err() {
2512                    return false;
2513                }
2514                let parsed = record.unwrap();
2515                match f(&parsed) {
2516                    None => return true,
2517                    Some(false) => {}
2518                    Some(true) => {
2519                        if sam_writer.write(&parsed).is_err() {
2520                            return false;
2521                        }
2522                    }
2523                }
2524            }
2525            true
2526        }
2527        use std::fs::File;
2528        use std::io::Read;
2529        let bamfile = "./test/bam2sam_test.bam";
2530        let samfile = "./test/bam2sam_out.sam";
2531        let expectedfile = "./test/bam2sam_expected.sam";
2532        let result = from_bam_with_filter(bamfile, samfile, |_| Some(true));
2533        assert!(result);
2534        let mut expected = Vec::new();
2535        let mut written = Vec::new();
2536        assert!(File::open(expectedfile)
2537            .unwrap()
2538            .read_to_end(&mut expected)
2539            .is_ok());
2540        assert!(File::open(samfile)
2541            .unwrap()
2542            .read_to_end(&mut written)
2543            .is_ok());
2544        assert_eq!(expected, written);
2545    }
2546
2547    // #[cfg(feature = "curl")]
2548    // #[test]
2549    // fn test_http_connect() {
2550    //     let url: Url = Url::parse(
2551    //         "https://raw.githubusercontent.com/brainstorm/tiny-test-data/master/wgs/mt.bam",
2552    //     )
2553    //     .unwrap();
2554    //     let r = Reader::from_url(&url);
2555    //     println!("{:#?}", r);
2556    //     let r = r.unwrap();
2557
2558    //     assert_eq!(r.header().target_names()[0], b"chr1");
2559    // }
2560
2561    #[test]
2562    fn test_rc_records() {
2563        let (names, flags, seqs, quals, cigars) = gold();
2564        let mut bam = Reader::from_path(Path::new("test/test.bam")).expect("Error opening file.");
2565        let del_len = [1, 1, 1, 1, 1, 100000];
2566
2567        for (i, record) in bam.rc_records().enumerate() {
2568            //let rec = record.expect("Expected valid record");
2569            let rec = record.unwrap();
2570            println!("{}", str::from_utf8(rec.qname()).ok().unwrap());
2571            assert_eq!(rec.qname(), names[i]);
2572            assert_eq!(rec.flags(), flags[i]);
2573            assert_eq!(rec.seq().as_bytes(), seqs[i]);
2574
2575            let cigar = rec.cigar();
2576            assert_eq!(*cigar, cigars[i]);
2577
2578            let end_pos = cigar.end_pos();
2579            assert_eq!(end_pos, rec.pos() + 100 + del_len[i]);
2580            assert_eq!(
2581                cigar
2582                    .read_pos(end_pos as u32 - 10, false, false)
2583                    .unwrap()
2584                    .unwrap(),
2585                90
2586            );
2587            assert_eq!(
2588                cigar
2589                    .read_pos(rec.pos() as u32 + 20, false, false)
2590                    .unwrap()
2591                    .unwrap(),
2592                20
2593            );
2594            assert_eq!(cigar.read_pos(4000000, false, false).unwrap(), None);
2595            // fix qual offset
2596            let qual: Vec<u8> = quals[i].iter().map(|&q| q - 33).collect();
2597            assert_eq!(rec.qual(), &qual[..]);
2598        }
2599    }
2600
2601    #[test]
2602    fn test_aux_arrays() {
2603        let bam_header = Header::new();
2604        let mut test_record = Record::from_sam(
2605            &HeaderView::from_header(&bam_header),
2606            "ali1\t4\t*\t0\t0\t*\t*\t0\t0\tACGT\tFFFF".as_bytes(),
2607        )
2608        .unwrap();
2609
2610        let array_i8: Vec<i8> = vec![i8::MIN, -1, 0, 1, i8::MAX];
2611        let array_u8: Vec<u8> = vec![u8::MIN, 0, 1, u8::MAX];
2612        let array_i16: Vec<i16> = vec![i16::MIN, -1, 0, 1, i16::MAX];
2613        let array_u16: Vec<u16> = vec![u16::MIN, 0, 1, u16::MAX];
2614        let array_i32: Vec<i32> = vec![i32::MIN, -1, 0, 1, i32::MAX];
2615        let array_u32: Vec<u32> = vec![u32::MIN, 0, 1, u32::MAX];
2616        let array_f32: Vec<f32> = vec![f32::MIN, 0.0, -0.0, 0.1, 0.99, f32::MAX];
2617
2618        test_record
2619            .push_aux(b"XA", Aux::ArrayI8((&array_i8).into()))
2620            .unwrap();
2621        test_record
2622            .push_aux(b"XB", Aux::ArrayU8((&array_u8).into()))
2623            .unwrap();
2624        test_record
2625            .push_aux(b"XC", Aux::ArrayI16((&array_i16).into()))
2626            .unwrap();
2627        test_record
2628            .push_aux(b"XD", Aux::ArrayU16((&array_u16).into()))
2629            .unwrap();
2630        test_record
2631            .push_aux(b"XE", Aux::ArrayI32((&array_i32).into()))
2632            .unwrap();
2633        test_record
2634            .push_aux(b"XF", Aux::ArrayU32((&array_u32).into()))
2635            .unwrap();
2636        test_record
2637            .push_aux(b"XG", Aux::ArrayFloat((&array_f32).into()))
2638            .unwrap();
2639
2640        {
2641            let tag = b"XA";
2642            if let Ok(Aux::ArrayI8(array)) = test_record.aux(tag) {
2643                // Retrieve aux array
2644                let aux_array_content = array.iter().collect::<Vec<_>>();
2645                assert_eq!(aux_array_content, array_i8);
2646
2647                // Copy the stored aux array to another record
2648                {
2649                    let mut copy_test_record = test_record.clone();
2650
2651                    // Pushing a field with an existing tag should fail
2652                    assert!(copy_test_record.push_aux(tag, Aux::I8(3)).is_err());
2653
2654                    // Remove aux array from target record
2655                    copy_test_record.remove_aux(tag).unwrap();
2656                    assert!(copy_test_record.aux(tag).is_err());
2657
2658                    // Copy array to target record
2659                    let src_aux = test_record.aux(tag).unwrap();
2660                    assert!(copy_test_record.push_aux(tag, src_aux).is_ok());
2661                    if let Ok(Aux::ArrayI8(array)) = copy_test_record.aux(tag) {
2662                        let aux_array_content_copied = array.iter().collect::<Vec<_>>();
2663                        assert_eq!(aux_array_content_copied, array_i8);
2664                    } else {
2665                        panic!("Aux tag not found");
2666                    }
2667                }
2668            } else {
2669                panic!("Aux tag not found");
2670            }
2671        }
2672
2673        {
2674            let tag = b"XB";
2675            if let Ok(Aux::ArrayU8(array)) = test_record.aux(tag) {
2676                // Retrieve aux array
2677                let aux_array_content = array.iter().collect::<Vec<_>>();
2678                assert_eq!(aux_array_content, array_u8);
2679
2680                // Copy the stored aux array to another record
2681                {
2682                    let mut copy_test_record = test_record.clone();
2683
2684                    // Pushing a field with an existing tag should fail
2685                    assert!(copy_test_record.push_aux(tag, Aux::U8(3)).is_err());
2686
2687                    // Remove aux array from target record
2688                    copy_test_record.remove_aux(tag).unwrap();
2689                    assert!(copy_test_record.aux(tag).is_err());
2690
2691                    // Copy array to target record
2692                    let src_aux = test_record.aux(tag).unwrap();
2693                    assert!(copy_test_record.push_aux(tag, src_aux).is_ok());
2694                    if let Ok(Aux::ArrayU8(array)) = copy_test_record.aux(tag) {
2695                        let aux_array_content_copied = array.iter().collect::<Vec<_>>();
2696                        assert_eq!(aux_array_content_copied, array_u8);
2697                    } else {
2698                        panic!("Aux tag not found");
2699                    }
2700                }
2701            } else {
2702                panic!("Aux tag not found");
2703            }
2704        }
2705
2706        {
2707            let tag = b"XC";
2708            if let Ok(Aux::ArrayI16(array)) = test_record.aux(tag) {
2709                // Retrieve aux array
2710                let aux_array_content = array.iter().collect::<Vec<_>>();
2711                assert_eq!(aux_array_content, array_i16);
2712
2713                // Copy the stored aux array to another record
2714                {
2715                    let mut copy_test_record = test_record.clone();
2716
2717                    // Pushing a field with an existing tag should fail
2718                    assert!(copy_test_record.push_aux(tag, Aux::I16(3)).is_err());
2719
2720                    // Remove aux array from target record
2721                    copy_test_record.remove_aux(tag).unwrap();
2722                    assert!(copy_test_record.aux(tag).is_err());
2723
2724                    // Copy array to target record
2725                    let src_aux = test_record.aux(tag).unwrap();
2726                    assert!(copy_test_record.push_aux(tag, src_aux).is_ok());
2727                    if let Ok(Aux::ArrayI16(array)) = copy_test_record.aux(tag) {
2728                        let aux_array_content_copied = array.iter().collect::<Vec<_>>();
2729                        assert_eq!(aux_array_content_copied, array_i16);
2730                    } else {
2731                        panic!("Aux tag not found");
2732                    }
2733                }
2734            } else {
2735                panic!("Aux tag not found");
2736            }
2737        }
2738
2739        {
2740            let tag = b"XD";
2741            if let Ok(Aux::ArrayU16(array)) = test_record.aux(tag) {
2742                // Retrieve aux array
2743                let aux_array_content = array.iter().collect::<Vec<_>>();
2744                assert_eq!(aux_array_content, array_u16);
2745
2746                // Copy the stored aux array to another record
2747                {
2748                    let mut copy_test_record = test_record.clone();
2749
2750                    // Pushing a field with an existing tag should fail
2751                    assert!(copy_test_record.push_aux(tag, Aux::U16(3)).is_err());
2752
2753                    // Remove aux array from target record
2754                    copy_test_record.remove_aux(tag).unwrap();
2755                    assert!(copy_test_record.aux(tag).is_err());
2756
2757                    // Copy array to target record
2758                    let src_aux = test_record.aux(tag).unwrap();
2759                    assert!(copy_test_record.push_aux(tag, src_aux).is_ok());
2760                    if let Ok(Aux::ArrayU16(array)) = copy_test_record.aux(tag) {
2761                        let aux_array_content_copied = array.iter().collect::<Vec<_>>();
2762                        assert_eq!(aux_array_content_copied, array_u16);
2763                    } else {
2764                        panic!("Aux tag not found");
2765                    }
2766                }
2767            } else {
2768                panic!("Aux tag not found");
2769            }
2770        }
2771
2772        {
2773            let tag = b"XE";
2774            if let Ok(Aux::ArrayI32(array)) = test_record.aux(tag) {
2775                // Retrieve aux array
2776                let aux_array_content = array.iter().collect::<Vec<_>>();
2777                assert_eq!(aux_array_content, array_i32);
2778
2779                // Copy the stored aux array to another record
2780                {
2781                    let mut copy_test_record = test_record.clone();
2782
2783                    // Pushing a field with an existing tag should fail
2784                    assert!(copy_test_record.push_aux(tag, Aux::I32(3)).is_err());
2785
2786                    // Remove aux array from target record
2787                    copy_test_record.remove_aux(tag).unwrap();
2788                    assert!(copy_test_record.aux(tag).is_err());
2789
2790                    // Copy array to target record
2791                    let src_aux = test_record.aux(tag).unwrap();
2792                    assert!(copy_test_record.push_aux(tag, src_aux).is_ok());
2793                    if let Ok(Aux::ArrayI32(array)) = copy_test_record.aux(tag) {
2794                        let aux_array_content_copied = array.iter().collect::<Vec<_>>();
2795                        assert_eq!(aux_array_content_copied, array_i32);
2796                    } else {
2797                        panic!("Aux tag not found");
2798                    }
2799                }
2800            } else {
2801                panic!("Aux tag not found");
2802            }
2803        }
2804
2805        {
2806            let tag = b"XF";
2807            if let Ok(Aux::ArrayU32(array)) = test_record.aux(tag) {
2808                // Retrieve aux array
2809                let aux_array_content = array.iter().collect::<Vec<_>>();
2810                assert_eq!(aux_array_content, array_u32);
2811
2812                // Copy the stored aux array to another record
2813                {
2814                    let mut copy_test_record = test_record.clone();
2815
2816                    // Pushing a field with an existing tag should fail
2817                    assert!(copy_test_record.push_aux(tag, Aux::U32(3)).is_err());
2818
2819                    // Remove aux array from target record
2820                    copy_test_record.remove_aux(tag).unwrap();
2821                    assert!(copy_test_record.aux(tag).is_err());
2822
2823                    // Copy array to target record
2824                    let src_aux = test_record.aux(tag).unwrap();
2825                    assert!(copy_test_record.push_aux(tag, src_aux).is_ok());
2826                    if let Ok(Aux::ArrayU32(array)) = copy_test_record.aux(tag) {
2827                        let aux_array_content_copied = array.iter().collect::<Vec<_>>();
2828                        assert_eq!(aux_array_content_copied, array_u32);
2829                    } else {
2830                        panic!("Aux tag not found");
2831                    }
2832                }
2833            } else {
2834                panic!("Aux tag not found");
2835            }
2836        }
2837
2838        {
2839            let tag = b"XG";
2840            if let Ok(Aux::ArrayFloat(array)) = test_record.aux(tag) {
2841                // Retrieve aux array
2842                let aux_array_content = array.iter().collect::<Vec<_>>();
2843                assert_eq!(aux_array_content, array_f32);
2844
2845                // Copy the stored aux array to another record
2846                {
2847                    let mut copy_test_record = test_record.clone();
2848
2849                    // Pushing a field with an existing tag should fail
2850                    assert!(copy_test_record.push_aux(tag, Aux::Float(3.0)).is_err());
2851
2852                    // Remove aux array from target record
2853                    copy_test_record.remove_aux(tag).unwrap();
2854                    assert!(copy_test_record.aux(tag).is_err());
2855
2856                    // Copy array to target record
2857                    let src_aux = test_record.aux(tag).unwrap();
2858                    assert!(copy_test_record.push_aux(tag, src_aux).is_ok());
2859                    if let Ok(Aux::ArrayFloat(array)) = copy_test_record.aux(tag) {
2860                        let aux_array_content_copied = array.iter().collect::<Vec<_>>();
2861                        assert_eq!(aux_array_content_copied, array_f32);
2862                    } else {
2863                        panic!("Aux tag not found");
2864                    }
2865                }
2866            } else {
2867                panic!("Aux tag not found");
2868            }
2869        }
2870
2871        // Test via `Iterator` impl
2872        for item in test_record.aux_iter() {
2873            match item.unwrap() {
2874                (b"XA", Aux::ArrayI8(array)) => {
2875                    assert_eq!(&array.iter().collect::<Vec<_>>(), &array_i8);
2876                }
2877                (b"XB", Aux::ArrayU8(array)) => {
2878                    assert_eq!(&array.iter().collect::<Vec<_>>(), &array_u8);
2879                }
2880                (b"XC", Aux::ArrayI16(array)) => {
2881                    assert_eq!(&array.iter().collect::<Vec<_>>(), &array_i16);
2882                }
2883                (b"XD", Aux::ArrayU16(array)) => {
2884                    assert_eq!(&array.iter().collect::<Vec<_>>(), &array_u16);
2885                }
2886                (b"XE", Aux::ArrayI32(array)) => {
2887                    assert_eq!(&array.iter().collect::<Vec<_>>(), &array_i32);
2888                }
2889                (b"XF", Aux::ArrayU32(array)) => {
2890                    assert_eq!(&array.iter().collect::<Vec<_>>(), &array_u32);
2891                }
2892                (b"XG", Aux::ArrayFloat(array)) => {
2893                    assert_eq!(&array.iter().collect::<Vec<_>>(), &array_f32);
2894                }
2895                _ => {
2896                    panic!();
2897                }
2898            }
2899        }
2900
2901        // Test via `PartialEq` impl
2902        assert_eq!(
2903            test_record.aux(b"XA").unwrap(),
2904            Aux::ArrayI8((&array_i8).into())
2905        );
2906        assert_eq!(
2907            test_record.aux(b"XB").unwrap(),
2908            Aux::ArrayU8((&array_u8).into())
2909        );
2910        assert_eq!(
2911            test_record.aux(b"XC").unwrap(),
2912            Aux::ArrayI16((&array_i16).into())
2913        );
2914        assert_eq!(
2915            test_record.aux(b"XD").unwrap(),
2916            Aux::ArrayU16((&array_u16).into())
2917        );
2918        assert_eq!(
2919            test_record.aux(b"XE").unwrap(),
2920            Aux::ArrayI32((&array_i32).into())
2921        );
2922        assert_eq!(
2923            test_record.aux(b"XF").unwrap(),
2924            Aux::ArrayU32((&array_u32).into())
2925        );
2926        assert_eq!(
2927            test_record.aux(b"XG").unwrap(),
2928            Aux::ArrayFloat((&array_f32).into())
2929        );
2930    }
2931
2932    #[test]
2933    fn test_aux_scalars() {
2934        let bam_header = Header::new();
2935        let mut test_record = Record::from_sam(
2936            &HeaderView::from_header(&bam_header),
2937            "ali1\t4\t*\t0\t0\t*\t*\t0\t0\tACGT\tFFFF".as_bytes(),
2938        )
2939        .unwrap();
2940
2941        test_record.push_aux(b"XA", Aux::I8(i8::MIN)).unwrap();
2942        test_record.push_aux(b"XB", Aux::I8(i8::MAX)).unwrap();
2943        test_record.push_aux(b"XC", Aux::U8(u8::MIN)).unwrap();
2944        test_record.push_aux(b"XD", Aux::U8(u8::MAX)).unwrap();
2945        test_record.push_aux(b"XE", Aux::I16(i16::MIN)).unwrap();
2946        test_record.push_aux(b"XF", Aux::I16(i16::MAX)).unwrap();
2947        test_record.push_aux(b"XG", Aux::U16(u16::MIN)).unwrap();
2948        test_record.push_aux(b"XH", Aux::U16(u16::MAX)).unwrap();
2949        test_record.push_aux(b"XI", Aux::I32(i32::MIN)).unwrap();
2950        test_record.push_aux(b"XJ", Aux::I32(i32::MAX)).unwrap();
2951        test_record.push_aux(b"XK", Aux::U32(u32::MIN)).unwrap();
2952        test_record.push_aux(b"XL", Aux::U32(u32::MAX)).unwrap();
2953        test_record
2954            .push_aux(b"XM", Aux::Float(std::f32::consts::PI))
2955            .unwrap();
2956        test_record
2957            .push_aux(b"XN", Aux::Double(std::f64::consts::PI))
2958            .unwrap();
2959        test_record
2960            .push_aux(b"XO", Aux::String("Test str"))
2961            .unwrap();
2962        test_record.push_aux(b"XP", Aux::I8(0)).unwrap();
2963
2964        let collected_aux_fields = test_record.aux_iter().collect::<Result<Vec<_>>>().unwrap();
2965        assert_eq!(
2966            collected_aux_fields,
2967            vec![
2968                (&b"XA"[..], Aux::I8(i8::MIN)),
2969                (&b"XB"[..], Aux::I8(i8::MAX)),
2970                (&b"XC"[..], Aux::U8(u8::MIN)),
2971                (&b"XD"[..], Aux::U8(u8::MAX)),
2972                (&b"XE"[..], Aux::I16(i16::MIN)),
2973                (&b"XF"[..], Aux::I16(i16::MAX)),
2974                (&b"XG"[..], Aux::U16(u16::MIN)),
2975                (&b"XH"[..], Aux::U16(u16::MAX)),
2976                (&b"XI"[..], Aux::I32(i32::MIN)),
2977                (&b"XJ"[..], Aux::I32(i32::MAX)),
2978                (&b"XK"[..], Aux::U32(u32::MIN)),
2979                (&b"XL"[..], Aux::U32(u32::MAX)),
2980                (&b"XM"[..], Aux::Float(std::f32::consts::PI)),
2981                (&b"XN"[..], Aux::Double(std::f64::consts::PI)),
2982                (&b"XO"[..], Aux::String("Test str")),
2983                (&b"XP"[..], Aux::I8(0)),
2984            ]
2985        );
2986    }
2987
2988    #[test]
2989    fn test_aux_array_partial_eq() {
2990        use record::AuxArray;
2991
2992        // Target types
2993        let one_data: Vec<i8> = vec![0, 1, 2, 3, 4, 5, 6];
2994        let one_aux_array = AuxArray::from(&one_data);
2995
2996        let two_data: Vec<i8> = vec![0, 1, 2, 3, 4, 5];
2997        let two_aux_array = AuxArray::from(&two_data);
2998
2999        assert_ne!(&one_data, &two_data);
3000        assert_ne!(&one_aux_array, &two_aux_array);
3001
3002        let one_aux = Aux::ArrayI8(one_aux_array);
3003        let two_aux = Aux::ArrayI8(two_aux_array);
3004        assert_ne!(&one_aux, &two_aux);
3005
3006        // Raw bytes
3007        let bam_header = Header::new();
3008        let mut test_record = Record::from_sam(
3009            &HeaderView::from_header(&bam_header),
3010            "ali1\t4\t*\t0\t0\t*\t*\t0\t0\tACGT\tFFFF".as_bytes(),
3011        )
3012        .unwrap();
3013
3014        test_record.push_aux(b"XA", one_aux).unwrap();
3015        test_record.push_aux(b"XB", two_aux).unwrap();
3016
3017        // RawLeBytes == RawLeBytes
3018        assert_eq!(
3019            test_record.aux(b"XA").unwrap(),
3020            test_record.aux(b"XA").unwrap()
3021        );
3022        // RawLeBytes != RawLeBytes
3023        assert_ne!(
3024            test_record.aux(b"XA").unwrap(),
3025            test_record.aux(b"XB").unwrap()
3026        );
3027
3028        // RawLeBytes == TargetType
3029        assert_eq!(
3030            test_record.aux(b"XA").unwrap(),
3031            Aux::ArrayI8((&one_data).into())
3032        );
3033        assert_eq!(
3034            test_record.aux(b"XB").unwrap(),
3035            Aux::ArrayI8((&two_data).into())
3036        );
3037        // RawLeBytes != TargetType
3038        assert_ne!(
3039            test_record.aux(b"XA").unwrap(),
3040            Aux::ArrayI8((&two_data).into())
3041        );
3042        assert_ne!(
3043            test_record.aux(b"XB").unwrap(),
3044            Aux::ArrayI8((&one_data).into())
3045        );
3046    }
3047
3048    /// Test if both text and binary representations of a BAM header are in sync (#156)
3049    #[test]
3050    fn test_bam_header_sync() {
3051        let reader = Reader::from_path("test/test_issue_156_no_text.bam").unwrap();
3052        let header_hashmap = Header::from_template(reader.header()).to_hashmap();
3053        let header_refseqs = header_hashmap.get("SQ").unwrap();
3054        assert_eq!(header_refseqs[0].get("SN").unwrap(), "ref_1",);
3055        assert_eq!(header_refseqs[0].get("LN").unwrap(), "10000000",);
3056    }
3057
3058    #[test]
3059    fn test_bam_new() {
3060        // Create the path to write the tmp test BAM
3061        let tmp = tempfile::Builder::new()
3062            .prefix("rust-htslib")
3063            .tempdir()
3064            .expect("Cannot create temp dir");
3065        let bampath = tmp.path().join("test.bam");
3066
3067        // write an unmapped BAM record (uBAM)
3068        {
3069            // Build the header
3070            let mut header = Header::new();
3071
3072            // Add the version
3073            header.push_record(
3074                HeaderRecord::new(b"HD")
3075                    .push_tag(b"VN", "1.6")
3076                    .push_tag(b"SO", "unsorted"),
3077            );
3078
3079            // Build the writer
3080            let mut writer = Writer::from_path(&bampath, &header, Format::Bam).unwrap();
3081
3082            // Build an empty record
3083            let record = Record::new();
3084
3085            // Write the record (this previously seg-faulted)
3086            assert!(writer.write(&record).is_ok());
3087        }
3088
3089        // Read the record
3090        {
3091            // Build th reader
3092            let mut reader = Reader::from_path(bampath).expect("Error opening file.");
3093
3094            // Read the record
3095            let mut rec = Record::new();
3096            match reader.read(&mut rec) {
3097                Some(r) => r.expect("Failed to read record."),
3098                None => panic!("No record read."),
3099            };
3100
3101            // Check a few things
3102            assert!(rec.is_unmapped());
3103            assert_eq!(rec.tid(), -1);
3104            assert_eq!(rec.pos(), -1);
3105            assert_eq!(rec.mtid(), -1);
3106            assert_eq!(rec.mpos(), -1);
3107        }
3108    }
3109
3110    #[test]
3111    fn test_idxstats_bam() {
3112        let mut reader = IndexedReader::from_path("test/test.bam").unwrap();
3113        let expected = vec![
3114            (0, 15072423, 6, 0),
3115            (1, 15279345, 0, 0),
3116            (2, 13783700, 0, 0),
3117            (3, 17493793, 0, 0),
3118            (4, 20924149, 0, 0),
3119            (-1, 0, 0, 0),
3120        ];
3121        let actual = reader.index_stats().unwrap();
3122        assert_eq!(expected, actual);
3123    }
3124
3125    #[test]
3126    fn test_number_mapped_and_unmapped_bam() {
3127        let reader = IndexedReader::from_path("test/test.bam").unwrap();
3128        let expected = (6, 0);
3129        let actual = reader.index().number_mapped_unmapped(0);
3130        assert_eq!(expected, actual);
3131    }
3132
3133    #[test]
3134    fn test_number_unmapped_global_bam() {
3135        let reader = IndexedReader::from_path("test/test_unmapped.bam").unwrap();
3136        let expected = 8;
3137        let actual = reader.index().number_unmapped();
3138        assert_eq!(expected, actual);
3139    }
3140
3141    #[test]
3142    fn test_idxstats_cram() {
3143        let mut reader = IndexedReader::from_path("test/test_cram.cram").unwrap();
3144        reader.set_reference("test/test_cram.fa").unwrap();
3145        let expected = vec![
3146            (0, 120, 2, 0),
3147            (1, 120, 2, 0),
3148            (2, 120, 2, 0),
3149            (-1, 0, 0, 0),
3150        ];
3151        let actual = reader.index_stats().unwrap();
3152        assert_eq!(expected, actual);
3153    }
3154
3155    #[test]
3156    fn test_slow_idxstats_cram() {
3157        let mut reader = IndexedReader::from_path("test/test_cram.cram").unwrap();
3158        reader.set_reference("test/test_cram.fa").unwrap();
3159        let expected = vec![
3160            (0, 120, 2, 0),
3161            (1, 120, 2, 0),
3162            (2, 120, 2, 0),
3163            (-1, 0, 0, 0),
3164        ];
3165        let actual = reader.index_stats().unwrap();
3166        assert_eq!(expected, actual);
3167    }
3168
3169    // #[test]
3170    // fn test_number_mapped_and_unmapped_cram() {
3171    //     let mut reader = IndexedReader::from_path("test/test_cram.cram").unwrap();
3172    //     reader.set_reference("test/test_cram.fa").unwrap();
3173    //     let expected = (2, 0);
3174    //     let actual = reader.index().number_mapped_unmapped(0);
3175    //     assert_eq!(expected, actual);
3176    // }
3177    //
3178    // #[test]
3179    // fn test_number_unmapped_global_cram() {
3180    //     let mut reader = IndexedReader::from_path("test/test_unmapped.cram").unwrap();
3181    //     let expected = 8;
3182    //     let actual = reader.index().number_unmapped();
3183    //     assert_eq!(expected, actual);
3184    // }
3185}