SAMRecord

Encapsulates a SAM/BAM/CRAM record, using the bam1_t type for memory efficiency, and the htslib helper functions for speed.

Constructors

this
this()

Construct blank SAMRecord with empty backing bam1_t

this
this(bam1_t* b)

Construct SAMRecord from supplied bam1_t

Destructor

~this
~this()

Undocumented in source.

Members

Functions

cigar
Cigar cigar()

Create cigar from bam1_t record

cigar
void cigar(Cigar cigar)

Assign a cigar string

flag
ushort flag()
void flag(ushort fl)

bitwise flag

getAlignedCoordinates
auto getAlignedCoordinates()

get aligned coordinates per each cigar op

getAlignedCoordinates
auto getAlignedCoordinates(long start, long end)

get aligned coordinates per each cigar op within range range is 0-based half open using chromosomal coordinates

getAlignedPairs
auto getAlignedPairs(long start, long end)

get a range of aligned read and reference positions this is meant to recreate functionality from pysam: https://pysam.readthedocs.io/en/latest/api.html#pysam.AlignedSegment.get_aligned_pairs range is 0-based half open using chromosomal coordinates

getAlignedPairs
auto getAlignedPairs()

get a range of aligned read and reference positions this is meant to recreate functionality from pysam: https://pysam.readthedocs.io/en/latest/api.html#pysam.AlignedSegment.get_aligned_pairs

insertSize
long insertSize()
void insertSize(long isize)

Presumably Insert size, but is undocumented. Per samtools source, is measured 5' to 5' https://github.com/samtools/samtools/blob/bd1a409aa750d25d70a093405d174db8709a3f2c/bam_mate.c#L320

isMapped
bool isMapped()

is read mapped?

isMateMapped
bool isMateMapped()

is read mate mapped?

isPaired
bool isPaired()

is read paired?

isReversed
bool isReversed()

is read reversed? bool bam_is_rev(bam1_t *b) { return ( ((*b).core.flag & BAM_FREVERSE) != 0 ); }

isSecondary
bool isSecondary()

is read secondary?

isSupplementary
bool isSupplementary()

is read supplementary?

l_extranul
int l_extranul()

number of EXTRA trailing \0 in qname; 0 <= l_extranul <= 3, see l_qname. Never set directly.

l_qname
int l_qname()

length of query name. Terminated by 1-4 \0, see l_extranul. Never set directly.

length
int length()

query (and quality string) length

matePos
long matePos()
void matePos(long mpos)

0-based leftmost coordinate of next read in template

mateReversed
bool mateReversed()

is mate reversed? bool bam_is_mrev(bam1_t *b) { return( ((*b).core.flag & BAM_FMREVERSE) != 0); }

mateTID
int mateTID()
void mateTID(int mtid)

chromosome ID of next read in template, defined by bam_hdr_t

opIndex
TagValue opIndex(string val)

Get aux tag from bam1_t record and return a TagValue

opIndexAssign
void opIndexAssign(T value, string index)

Assign a tag key value pair

opIndexAssign
void opIndexAssign(T[] value, string index)

Assign a tag key array pair

pos
long pos()
void pos(long pos)

0-based leftmost coordinate

q_scores
void q_scores(const(char)[] seq)

Add qscore sequence given that it is the same length as the bam sequence

qscores
const(char)[] qscores()

Return char array of the qscores see samtools/sam_view.c: get_quality

qscores
void qscores(ubyte[] seq)

Add qscore sequence given that it is the same length as the bam sequence

qual
ubyte qual()
void qual(ubyte q)

mapping quality

queryName
char[] queryName()

auto bam_get_qname(bam1_t *b) { return (cast(char*)(*b).data); }

queryName
void queryName(string qname)

Set query name

sequence
const(char)[] sequence()

Return char array of the sequence see samtools/sam_view.c: get_read

sequence
void sequence(const(char)[] seq)

Assigns sequence and resets quality score

strand
char strand()

Return read strand + or - (as char)

tid
int tid()
void tid(int tid)

chromosome ID, defined by bam_hdr_t

Structs

AlignedPair
struct AlignedPair(bool refSeq)

Undocumented in source.

Variables

b
bam1_t* b;

Meta