pub struct Record {
pub inner: bam1_t,
/* private fields */
}Expand description
A BAM record.
Fields
inner: bam1_tImplementations
sourceimpl Record
impl Record
pub fn from_inner(from: *mut bam1_t) -> Self
pub fn from_sam(header_view: &HeaderView, sam: &[u8]) -> Result<Record>
pub fn set_header(&mut self, header: Rc<HeaderView>)
pub fn inner_mut(&mut self) -> &mut bam1_t
pub fn inner(&self) -> &bam1_t
pub fn bin(&self) -> u16
pub fn set_bin(&mut self, bin: u16)
sourcepub fn unset_flags(&mut self)
pub fn unset_flags(&mut self)
Unset all flags.
sourcepub fn insert_size(&self) -> i64
pub fn insert_size(&self) -> i64
Get insert size.
sourcepub fn set_insert_size(&mut self, insert_size: i64)
pub fn set_insert_size(&mut self, insert_size: i64)
Set insert size.
sourcepub fn qname(&self) -> &[u8]ⓘNotable traits for &'_ mut [u8]impl<'_> Write for &'_ mut [u8]impl<'_> Read for &'_ [u8]
pub fn qname(&self) -> &[u8]ⓘNotable traits for &'_ mut [u8]impl<'_> Write for &'_ mut [u8]impl<'_> Read for &'_ [u8]
Get qname (read name). Complexity: O(1).
sourcepub fn set(
&mut self,
qname: &[u8],
cigar: Option<&CigarString>,
seq: &[u8],
qual: &[u8]
)
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.
pub fn cigar_len(&self) -> usize
sourcepub fn raw_cigar(&self) -> &[u32]
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.
sourcepub fn cigar(&self) -> CigarStringView
pub fn cigar(&self) -> CigarStringView
Return unpacked cigar string. This will create a fresh copy the Cigar data.
pub fn cigar_cached(&self) -> Option<&CigarStringView>
sourcepub fn cache_cigar(&mut self)
pub fn cache_cigar(&mut self)
Decode the cigar string and cache it inside the Record
pub fn seq_len(&self) -> usize
sourcepub fn qual(&self) -> &[u8]ⓘNotable traits for &'_ mut [u8]impl<'_> Write for &'_ mut [u8]impl<'_> Read for &'_ [u8]
pub fn qual(&self) -> &[u8]ⓘNotable traits for &'_ mut [u8]impl<'_> Write for &'_ mut [u8]impl<'_> Read for &'_ [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).
sourcepub fn aux(&self, tag: &[u8]) -> Result<Aux<'_>>
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.
sourcepub fn aux_iter(&self) -> AuxIter<'_>ⓘNotable traits for AuxIter<'a>impl<'a> Iterator for AuxIter<'a> type Item = Result<(&'a [u8], Aux<'a>)>;
pub fn aux_iter(&self) -> AuxIter<'_>ⓘNotable traits for AuxIter<'a>impl<'a> Iterator for AuxIter<'a> type Item = Result<(&'a [u8], Aux<'a>)>;
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.
pub fn remove_aux(&mut self, tag: &[u8]) -> Result<()>
sourcepub fn read_pair_orientation(&self) -> SequenceReadPairOrientation
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.
pub fn is_paired(&self) -> bool
pub fn set_paired(&mut self)
pub fn unset_paired(&mut self)
pub fn is_proper_pair(&self) -> bool
pub fn set_proper_pair(&mut self)
pub fn unset_proper_pair(&mut self)
pub fn is_unmapped(&self) -> bool
pub fn set_unmapped(&mut self)
pub fn unset_unmapped(&mut self)
pub fn is_mate_unmapped(&self) -> bool
pub fn set_mate_unmapped(&mut self)
pub fn unset_mate_unmapped(&mut self)
pub fn is_reverse(&self) -> bool
pub fn set_reverse(&mut self)
pub fn unset_reverse(&mut self)
pub fn is_mate_reverse(&self) -> bool
pub fn set_mate_reverse(&mut self)
pub fn unset_mate_reverse(&mut self)
pub fn is_first_in_template(&self) -> bool
pub fn set_first_in_template(&mut self)
pub fn unset_first_in_template(&mut self)
pub fn is_last_in_template(&self) -> bool
pub fn set_last_in_template(&mut self)
pub fn unset_last_in_template(&mut self)
pub fn is_secondary(&self) -> bool
pub fn set_secondary(&mut self)
pub fn unset_secondary(&mut self)
pub fn is_quality_check_failed(&self) -> bool
pub fn set_quality_check_failed(&mut self)
pub fn unset_quality_check_failed(&mut self)
pub fn is_duplicate(&self) -> bool
pub fn set_duplicate(&mut self)
pub fn unset_duplicate(&mut self)
pub fn is_supplementary(&self) -> bool
pub fn set_supplementary(&mut self)
pub fn unset_supplementary(&mut self)
Trait Implementations
sourceimpl AbstractInterval for Record
impl AbstractInterval for Record
sourcefn contig(&self) -> &str
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).
sourcefn range(&self) -> Range<Position>
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.
sourcefn contains<L>(&self, locus: L) -> bool where
L: AbstractLocus,
fn contains<L>(&self, locus: L) -> bool where
L: AbstractLocus,
Return true if interval contains given locus.
sourceimpl BamRecordExtensions for Record
impl BamRecordExtensions for Record
sourcefn aligned_blocks(&self) -> IterAlignedBlocksⓘNotable traits for IterAlignedBlocksimpl Iterator for IterAlignedBlocks type Item = [i64; 2];
fn aligned_blocks(&self) -> IterAlignedBlocksⓘNotable traits for IterAlignedBlocksimpl Iterator for IterAlignedBlocks type Item = [i64; 2];
iterator over start and end positions of aligned gapless blocks Read more
sourcefn introns(&self) -> IterIntronsⓘNotable traits for IterIntronsimpl Iterator for IterIntrons type Item = [i64; 2];
fn introns(&self) -> IterIntronsⓘNotable traits for IterIntronsimpl Iterator for IterIntrons type Item = [i64; 2];
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 Read more
sourcefn aligned_block_pairs(&self) -> IterAlignedBlockPairsⓘNotable traits for IterAlignedBlockPairsimpl Iterator for IterAlignedBlockPairs type Item = ([i64; 2], [i64; 2]);
fn aligned_block_pairs(&self) -> IterAlignedBlockPairsⓘNotable traits for IterAlignedBlockPairsimpl Iterator for IterAlignedBlockPairs type Item = ([i64; 2], [i64; 2]);
Iter over <([read_start, read_stop], [genome_start, genome_stop]) blocks of continously aligned reads. Read more
sourcefn aligned_pairs(&self) -> IterAlignedPairsⓘNotable traits for IterAlignedPairsimpl Iterator for IterAlignedPairs type Item = [i64; 2];
fn aligned_pairs(&self) -> IterAlignedPairsⓘNotable traits for IterAlignedPairsimpl Iterator for IterAlignedPairs type Item = [i64; 2];
iter aligned read and reference positions on a basepair level Read more
sourcefn aligned_pairs_full(&self) -> IterAlignedPairsFullⓘNotable traits for IterAlignedPairsFullimpl Iterator for IterAlignedPairsFull type Item = [Option<i64>; 2];
fn aligned_pairs_full(&self) -> IterAlignedPairsFullⓘNotable traits for IterAlignedPairsFullimpl Iterator for IterAlignedPairsFull type Item = [Option<i64>; 2];
iter list of read and reference positions on a basepair level. Read more
sourcefn cigar_stats_nucleotides(&self) -> HashMap<Cigar, i32>
fn cigar_stats_nucleotides(&self) -> HashMap<Cigar, i32>
the number of nucleotides covered by each Cigar::* variant. Read more
sourcefn cigar_stats_blocks(&self) -> HashMap<Cigar, i32>
fn cigar_stats_blocks(&self) -> HashMap<Cigar, i32>
the number of occurrences of each each Cigar::* variant Read more
sourcefn reference_positions(&self) -> Box<dyn Iterator<Item = i64>>
fn reference_positions(&self) -> Box<dyn Iterator<Item = i64>>
iter over reference positions that this read aligns to Read more
sourcefn reference_positions_full(&self) -> Box<dyn Iterator<Item = Option<i64>>>
fn reference_positions_full(&self) -> Box<dyn Iterator<Item = Option<i64>>>
iter over reference positions that this read aligns to Read more
sourcefn reference_start(&self) -> i64
fn reference_start(&self) -> i64
left most aligned reference position of the read on the reference genome.
sourcefn reference_end(&self) -> i64
fn reference_end(&self) -> i64
right most aligned reference position of the read on the reference genome.
sourcefn seq_len_from_cigar(&self, include_hard_clip: bool) -> usize
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
sourceimpl SequenceRead for Record
impl SequenceRead for Record
impl Eq for Record
impl Send for Record
impl Sync for Record
Auto Trait Implementations
Blanket Implementations
sourceimpl<T> BorrowMut<T> for T where
T: ?Sized,
impl<T> BorrowMut<T> for T where
T: ?Sized,
const: unstable · sourcefn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Mutably borrows from an owned value. Read more
sourceimpl<T> ToOwned for T where
T: Clone,
impl<T> ToOwned for T where
T: Clone,
type Owned = T
type Owned = T
The resulting type after obtaining ownership.
sourcefn clone_into(&self, target: &mut T)
fn clone_into(&self, target: &mut T)
toowned_clone_into)Uses borrowed data to replace owned data, usually by cloning. Read more