Struct rust_htslib::bam::record::Record

source ·
pub struct Record {
    pub inner: bam1_t,
    /* private fields */
}
Expand description

A BAM record.

Fields§

§inner: bam1_t

Implementations§

source§

impl Record

source

pub fn new() -> Self

Create an empty BAM record.

source

pub fn from_inner(from: *mut bam1_t) -> Self

source

pub fn from_sam(header_view: &HeaderView, sam: &[u8]) -> Result<Record>

source

pub fn set_header(&mut self, header: Rc<HeaderView>)

source

pub fn inner_mut(&mut self) -> &mut bam1_t

source

pub fn inner(&self) -> &bam1_t

source

pub fn tid(&self) -> i32

Get target id.

source

pub fn set_tid(&mut self, tid: i32)

Set target id.

source

pub fn pos(&self) -> i64

Get position (0-based).

source

pub fn set_pos(&mut self, pos: i64)

Set position (0-based).

source

pub fn bin(&self) -> u16

source

pub fn set_bin(&mut self, bin: u16)

source

pub fn mapq(&self) -> u8

Get MAPQ.

source

pub fn set_mapq(&mut self, mapq: u8)

Set MAPQ.

source

pub fn strand(&self) -> ReqStrand

Get strand information from record flags.

source

pub fn flags(&self) -> u16

Get raw flags.

source

pub fn set_flags(&mut self, flags: u16)

Set raw flags.

source

pub fn unset_flags(&mut self)

Unset all flags.

source

pub fn mtid(&self) -> i32

Get target id of mate.

source

pub fn set_mtid(&mut self, mtid: i32)

Set target id of mate.

source

pub fn mpos(&self) -> i64

Get mate position.

source

pub fn set_mpos(&mut self, mpos: i64)

Set mate position.

source

pub fn insert_size(&self) -> i64

Get insert size.

source

pub fn set_insert_size(&mut self, insert_size: i64)

Set insert size.

source

pub fn qname(&self) -> &[u8]

Get qname (read name). Complexity: O(1).

source

pub fn set_data(&mut self, new_data: &[u8])

Set the variable length data buffer

source

pub fn set( &mut self, qname: &[u8], cigar: Option<&CigarString>, seq: &[u8], qual: &[u8] )

Set variable length data (qname, cigar, seq, qual). The aux data is left unchanged. qual is Phred-scaled quality values, without any offset. NOTE: seq.len() must equal qual.len() or this method will panic. If you don’t have quality values use let quals = vec![ 255 as u8; seq.len()]; as a placeholder that will be recognized as missing QVs by samtools.

source

pub fn set_qname(&mut self, new_qname: &[u8])

Replace current qname with a new one.

source

pub fn cigar_len(&self) -> usize

source

pub fn raw_cigar(&self) -> &[u32]

Get reference to raw cigar string representation (as stored in BAM file). Usually, the method Record::cigar should be used instead.

source

pub fn cigar(&self) -> CigarStringView

Return unpacked cigar string. This will create a fresh copy the Cigar data.

source

pub fn cigar_cached(&self) -> Option<&CigarStringView>

source

pub fn cache_cigar(&mut self)

Decode the cigar string and cache it inside the Record

source

pub fn seq_len(&self) -> usize

source

pub fn seq(&self) -> Seq<'_>

Get read sequence. Complexity: O(1).

source

pub fn qual(&self) -> &[u8]

Get base qualities (PHRED-scaled probability that base is wrong). This does not entail any offsets, hence the qualities can be used directly without e.g. subtracting 33. Complexity: O(1).

source

pub fn aux(&self, tag: &[u8]) -> Result<Aux<'_>>

Look up an auxiliary field by its tag.

Only the first two bytes of a given tag are used for the look-up of a field. See Aux for more details.

source

pub fn aux_iter(&self) -> AuxIter<'_>

Returns an iterator over the auxiliary fields of the record.

When an error occurs, the Err variant will be returned and the iterator will not be able to advance anymore.

source

pub fn push_aux(&mut self, tag: &[u8], value: Aux<'_>) -> Result<()>

Add auxiliary data.

source

pub fn remove_aux(&mut self, tag: &[u8]) -> Result<()>

source

pub fn basemods_iter(&self) -> Result<BaseModificationsIter<'_>>

Access the base modifications associated with this Record through the MM tag. Example:

   use rust_htslib::bam::{Read, Reader, Record};
   let mut bam = Reader::from_path(&"test/base_mods/MM-orient.sam").unwrap();
   let mut mod_count = 0;
   for r in bam.records() {
       let record = r.unwrap();
       if let Ok(mods) = record.basemods_iter() {
           // print metadata for the modifications present in this record
           for mod_code in mods.recorded() {
               if let Ok(mod_metadata) = mods.query_type(*mod_code) {
                   println!("mod found with code {}/{} flags: [{} {} {}]",
                             mod_code, *mod_code as u8 as char,
                             mod_metadata.strand, mod_metadata.implicit, mod_metadata.canonical as u8 as char);
               }
           }

           // iterate over the modifications in this record
           // the modifications are returned as a tuple with the
           // position within SEQ and an hts_base_mod struct
           for res in mods {
               if let Ok( (position, m) ) = res {
                   println!("{} {},{}", position, m.modified_base as u8 as char, m.qual);
                   mod_count += 1;
               }
           }
       };
   }
   assert_eq!(mod_count, 14);
source

pub fn basemods_position_iter( &self ) -> Result<BaseModificationsPositionIter<'_>>

An iterator that returns all of the modifications for each position as a vector. This is useful for the case where multiple possible modifications can be annotated at a single position (for example a C could be 5-mC or 5-hmC)

source

pub fn read_pair_orientation(&self) -> SequenceReadPairOrientation

Infer read pair orientation from record. Returns SequenceReadPairOrientation::None if record is not paired, mates are not mapping to the same contig, or mates start at the same position.

source

pub fn is_paired(&self) -> bool

source

pub fn set_paired(&mut self)

source

pub fn unset_paired(&mut self)

source

pub fn is_proper_pair(&self) -> bool

source

pub fn set_proper_pair(&mut self)

source

pub fn unset_proper_pair(&mut self)

source

pub fn is_unmapped(&self) -> bool

source

pub fn set_unmapped(&mut self)

source

pub fn unset_unmapped(&mut self)

source

pub fn is_mate_unmapped(&self) -> bool

source

pub fn set_mate_unmapped(&mut self)

source

pub fn unset_mate_unmapped(&mut self)

source

pub fn is_reverse(&self) -> bool

source

pub fn set_reverse(&mut self)

source

pub fn unset_reverse(&mut self)

source

pub fn is_mate_reverse(&self) -> bool

source

pub fn set_mate_reverse(&mut self)

source

pub fn unset_mate_reverse(&mut self)

source

pub fn is_first_in_template(&self) -> bool

source

pub fn set_first_in_template(&mut self)

source

pub fn unset_first_in_template(&mut self)

source

pub fn is_last_in_template(&self) -> bool

source

pub fn set_last_in_template(&mut self)

source

pub fn unset_last_in_template(&mut self)

source

pub fn is_secondary(&self) -> bool

source

pub fn set_secondary(&mut self)

source

pub fn unset_secondary(&mut self)

source

pub fn is_quality_check_failed(&self) -> bool

source

pub fn set_quality_check_failed(&mut self)

source

pub fn unset_quality_check_failed(&mut self)

source

pub fn is_duplicate(&self) -> bool

source

pub fn set_duplicate(&mut self)

source

pub fn unset_duplicate(&mut self)

source

pub fn is_supplementary(&self) -> bool

source

pub fn set_supplementary(&mut self)

source

pub fn unset_supplementary(&mut self)

Trait Implementations§

source§

impl AbstractInterval for Record

source§

fn contig(&self) -> &str

Return contig name. Panics if record does not know its header (which happens if it has not been read from a file).

source§

fn range(&self) -> Range<Position>

Return genomic range covered by alignment. Panics if Record::cache_cigar() has not been called first or Record::pos() is less than zero.

source§

fn contains<L>(&self, locus: L) -> boolwhere L: AbstractLocus,

Return true if interval contains given locus.
source§

impl BamRecordExtensions for Record

source§

fn aligned_blocks(&self) -> IterAlignedBlocks

iterator over start and end positions of aligned gapless blocks Read more
source§

fn introns(&self) -> IterIntrons

This scans the CIGAR for reference skips and reports their positions. It does not inspect the reported regions for actual splice sites. pysam: get_introns
source§

fn aligned_block_pairs(&self) -> IterAlignedBlockPairs

Iter over <([read_start, read_stop], [genome_start, genome_stop]) blocks of continously aligned reads. Read more
source§

fn aligned_pairs(&self) -> IterAlignedPairs

iter aligned read and reference positions on a basepair level Read more
source§

fn aligned_pairs_full(&self) -> IterAlignedPairsFull

iter list of read and reference positions on a basepair level. Read more
source§

fn cigar_stats_nucleotides(&self) -> HashMap<Cigar, i32>

the number of nucleotides covered by each Cigar::* variant. Read more
source§

fn cigar_stats_blocks(&self) -> HashMap<Cigar, i32>

the number of occurrences of each each Cigar::* variant Read more
source§

fn reference_positions(&self) -> Box<dyn Iterator<Item = i64>>

iter over reference positions that this read aligns to Read more
source§

fn reference_positions_full(&self) -> Box<dyn Iterator<Item = Option<i64>>>

iter over reference positions that this read aligns to Read more
source§

fn reference_start(&self) -> i64

left most aligned reference position of the read on the reference genome.
source§

fn reference_end(&self) -> i64

right most aligned reference position of the read on the reference genome.
source§

fn seq_len_from_cigar(&self, include_hard_clip: bool) -> usize

infer the query length from the cigar string, optionally include hard clipped bases Read more
source§

impl Clone for Record

source§

fn clone(&self) -> Self

Returns a copy of the value. Read more
1.0.0 · source§

fn clone_from(&mut self, source: &Self)

Performs copy-assignment from source. Read more
source§

impl Debug for Record

source§

fn fmt(&self, fmt: &mut Formatter<'_>) -> Result<(), Error>

Formats the value using the given formatter. Read more
source§

impl Default for Record

source§

fn default() -> Self

Returns the “default value” for a type. Read more
source§

impl Drop for Record

source§

fn drop(&mut self)

Executes the destructor for this type. Read more
source§

impl PartialEq<Record> for Record

source§

fn eq(&self, other: &Record) -> bool

This method tests for self and other values to be equal, and is used by ==.
1.0.0 · source§

fn ne(&self, other: &Rhs) -> bool

This method tests for !=. The default implementation is almost always sufficient, and should not be overridden without very good reason.
source§

impl SequenceRead for Record

source§

fn name(&self) -> &[u8]

Read name.
source§

fn base(&self, i: usize) -> u8

Base at position i in the read.
source§

fn base_qual(&self, i: usize) -> u8

Base quality at position i in the read.
source§

fn len(&self) -> usize

Read length.
source§

fn is_empty(&self) -> bool

Return true if read is empty.
source§

impl Eq for Record

source§

impl Send for Record

source§

impl Sync for Record

Auto Trait Implementations§

Blanket Implementations§

source§

impl<T> Any for Twhere T: 'static + ?Sized,

source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
source§

impl<T> Borrow<T> for Twhere T: ?Sized,

source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
source§

impl<T> BorrowMut<T> for Twhere T: ?Sized,

source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
source§

impl<T> From<T> for T

source§

fn from(t: T) -> T

Returns the argument unchanged.

source§

impl<T, U> Into<U> for Twhere U: From<T>,

source§

fn into(self) -> U

Calls U::from(self).

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

source§

impl<T> ToOwned for Twhere T: Clone,

§

type Owned = T

The resulting type after obtaining ownership.
source§

fn to_owned(&self) -> T

Creates owned data from borrowed data, usually by cloning. Read more
source§

fn clone_into(&self, target: &mut T)

Uses borrowed data to replace owned data, usually by cloning. Read more
source§

impl<T, U> TryFrom<U> for Twhere U: Into<T>,

§

type Error = Infallible

The type returned in the event of a conversion error.
source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
source§

impl<T, U> TryInto<U> for Twhere U: TryFrom<T>,

§

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.