edit

SAM

Description

SAM is a text-based file format for representing sequence alignments.

  • Reader type: SAM.Reader
  • Writer type: SAM.Writer
  • Element type: SAM.Record

This module provides 16-bit flags defined in the SAM specs:

Flag Bit Description
SAM.FLAG_PAIRED 0x0001 template having multiple segments in sequencing
SAM.FLAG_PROPER_PAIR 0x0002 each segment properly aligned according to the aligner
SAM.FLAG_UNMAP 0x0004 segment unmapped
SAM.FLAG_MUNMAP 0x0008 next segment in the template unmapped
SAM.FLAG_REVERSE 0x0010 SEQ being reverse complemented
SAM.FLAG_MREVERSE 0x0020 SEQ of the next segment in the template being reverse complemented
SAM.FLAG_READ1 0x0040 the first segment in the template
SAM.FLAG_READ2 0x0080 the last segment in the template
SAM.FLAG_SECONDARY 0x0100 secondary alignment
SAM.FLAG_QCFAIL 0x0200 not passing filters, such as platform/vendor quality controls
SAM.FLAG_DUP 0x0400 PCR or optical duplicate
SAM.FLAG_SUPPLEMENTARY 0x0800 supplementary alignment

Examples

TODO

Accessors

# BioAlignments.SAM.ReaderType.

SAM.Reader(input::IO)

Create a data reader of the SAM file format.

Arguments

  • input: data source

source

# BioAlignments.SAM.headerFunction.

header(reader::Reader)::Header

Get the header of reader.

source

# BioAlignments.SAM.HeaderType.

SAM.Header()

Create an empty header.

source

# Base.findMethod.

find(header::Header, key::AbstractString)::Vector{MetaInfo}

Find metainfo objects satisfying SAM.tag(metainfo) == key.

source

# BioAlignments.SAM.WriterType.

Writer(output::IO, header::Header=Header())

Create a data writer of the SAM file format.

Arguments

  • output: data sink
  • header=Header(): SAM header object

source

# BioAlignments.SAM.MetaInfoType.

MetaInfo(str::AbstractString)

Create a SAM metainfo from str.

Examples

julia> SAM.MetaInfo("@CO    some comment")
BioAlignments.SAM.MetaInfo:
    tag: CO
  value: some comment

julia> SAM.MetaInfo("@SQ    SN:chr1 LN:12345")
BioAlignments.SAM.MetaInfo:
    tag: SQ
  value: SN=chr1 LN=12345

source

MetaInfo(tag::AbstractString, value)

Create a SAM metainfo with tag and value.

tag is a two-byte ASCII string. If tag is "CO", value must be a string; otherwise, value is an iterable object with key and value pairs.

Examples

julia> SAM.MetaInfo("CO", "some comment")
BioAlignments.SAM.MetaInfo:
    tag: CO
  value: some comment

julia> string(ans)
"@CO    some comment"

julia> SAM.MetaInfo("SQ", ["SN" => "chr1", "LN" => 12345])
BioAlignments.SAM.MetaInfo:
    tag: SQ
  value: SN=chr1 LN=12345

julia> string(ans)
"@SQ    SN:chr1 LN:12345"

source

# BioAlignments.SAM.iscommentFunction.

iscomment(metainfo::MetaInfo)::Bool

Test if metainfo is a comment (i.e. its tag is "CO").

source

# BioAlignments.SAM.tagFunction.

tag(metainfo::MetaInfo)::String

Get the tag of metainfo.

source

# BioAlignments.SAM.valueFunction.

value(metainfo::MetaInfo)::String

Get the value of metainfo as a string.

source

# BioAlignments.SAM.keyvaluesFunction.

keyvalues(metainfo::MetaInfo)::Vector{Pair{String,String}}

Get the values of metainfo as string pairs.

source

# BioAlignments.SAM.RecordType.

SAM.Record()

Create an unfilled SAM record.

source

SAM.Record(data::Vector{UInt8})

Create a SAM record from data. This function verifies the format and indexes fields for accessors. Note that the ownership of data is transferred to a new record object.

source

SAM.Record(str::AbstractString)

Create a SAM record from str. This function verifies the format and indexes fields for accessors.

source

# BioAlignments.SAM.flagFunction.

flag(record::Record)::UInt16

Get the bitwise flag of record.

source

# BioAlignments.SAM.ismappedFunction.

ismapped(record::Record)::Bool

Test if record is mapped.

source

# BioAlignments.SAM.isprimaryFunction.

isprimary(record::Record)::Bool

Test if record is a primary line of the read.

This is equivalent to flag(record) & 0x900 == 0.

source

# BioAlignments.SAM.refnameFunction.

refname(record::Record)::String

Get the reference sequence name of record.

source

# BioAlignments.SAM.positionFunction.

position(record::Record)::Int

Get the 1-based leftmost mapping position of record.

source

# BioAlignments.SAM.rightpositionFunction.

rightposition(record::Record)::Int

Get the 1-based rightmost mapping position of record.

source

# BioAlignments.SAM.isnextmappedFunction.

isnextmapped(record::Record)::Bool

Test if the mate/next read of record is mapped.

source

# BioAlignments.SAM.nextrefnameFunction.

nextrefname(record::Record)::String

Get the reference name of the mate/next read of record.

source

# BioAlignments.SAM.nextpositionFunction.

nextposition(record::Record)::Int

Get the position of the mate/next read of record.

source

# BioAlignments.SAM.mappingqualityFunction.

mappingquality(record::Record)::UInt8

Get the mapping quality of record.

source

# BioAlignments.SAM.cigarFunction.

cigar(record::Record)::String

Get the CIGAR string of record.

source

# BioAlignments.SAM.alignmentFunction.

alignment(record::Record)::BioAlignments.Alignment

Get the alignment of record.

source

# BioAlignments.SAM.alignlengthFunction.

alignlength(record::Record)::Int

Get the alignment length of record.

source

# BioAlignments.SAM.tempnameFunction.

tempname(record::Record)::String

Get the query template name of record.

source

# BioAlignments.SAM.templengthFunction.

templength(record::Record)::Int

Get the template length of record.

source

# BioAlignments.SAM.sequenceFunction.

sequence(record::Record)::BioSequences.DNASequence

Get the segment sequence of record.

source

sequence(::Type{String}, record::Record)::String

Get the segment sequence of record as String.

source

# BioAlignments.SAM.seqlengthFunction.

seqlength(record::Record)::Int

Get the sequence length of record.

source

# BioAlignments.SAM.qualityFunction.

quality(record::Record)::Vector{UInt8}

Get the Phred-scaled base quality of record.

source

quality(::Type{String}, record::Record)::String

Get the ASCII-encoded base quality of record.

source

# BioAlignments.SAM.auxdataFunction.

auxdata(record::Record)::Dict{String,Any}

Get the auxiliary data (optional fields) of record.

source