1
0
Fork 0
mirror of https://github.com/MillironX/XAM.jl.git synced 2024-11-15 14:53:11 +00:00
XAM.jl/dev/man/hts-files/index.html
2020-08-02 04:19:00 +00:00

107 lines
42 KiB
HTML

<!DOCTYPE html>
<html lang="en"><head><meta charset="UTF-8"/><meta name="viewport" content="width=device-width, initial-scale=1.0"/><title>SAM and BAM · XAM.jl</title><link href="https://fonts.googleapis.com/css?family=Lato|Roboto+Mono" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/font-awesome/5.11.2/css/fontawesome.min.css" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/font-awesome/5.11.2/css/solid.min.css" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/font-awesome/5.11.2/css/brands.min.css" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/KaTeX/0.11.1/katex.min.css" rel="stylesheet" type="text/css"/><script>documenterBaseURL="../.."</script><script src="https://cdnjs.cloudflare.com/ajax/libs/require.js/2.3.6/require.min.js" data-main="../../assets/documenter.js"></script><script src="../../siteinfo.js"></script><script src="../../../versions.js"></script><link class="docs-theme-link" rel="stylesheet" type="text/css" href="../../assets/themes/documenter-dark.css" data-theme-name="documenter-dark"/><link class="docs-theme-link" rel="stylesheet" type="text/css" href="../../assets/themes/documenter-light.css" data-theme-name="documenter-light" data-theme-primary/><script src="../../assets/themeswap.js"></script></head><body><div id="documenter"><nav class="docs-sidebar"><a class="docs-logo" href="../../"><img src="../../assets/logo.svg" alt="XAM.jl logo"/></a><div class="docs-package-name"><span class="docs-autofit">XAM.jl</span></div><form class="docs-search" action="../../search/"><input class="docs-search-query" id="documenter-search-query" name="q" type="text" placeholder="Search docs"/></form><ul class="docs-menu"><li><a class="tocitem" href="../../">Home</a></li><li class="is-active"><a class="tocitem" href>SAM and BAM</a><ul class="internal"><li><a class="tocitem" href="#Introduction-1"><span>Introduction</span></a></li><li><a class="tocitem" href="#Reading-SAM-and-BAM-files-1"><span>Reading SAM and BAM files</span></a></li><li><a class="tocitem" href="#SAM-and-BAM-Headers-1"><span>SAM and BAM Headers</span></a></li><li><a class="tocitem" href="#SAM-and-BAM-Records-1"><span>SAM and BAM Records</span></a></li><li><a class="tocitem" href="#Accessing-auxiliary-data-1"><span>Accessing auxiliary data</span></a></li><li><a class="tocitem" href="#Getting-records-in-a-range-1"><span>Getting records in a range</span></a></li><li><a class="tocitem" href="#Getting-records-overlapping-genomic-features-1"><span>Getting records overlapping genomic features</span></a></li><li><a class="tocitem" href="#Writing-files-1"><span>Writing files</span></a></li></ul></li><li><a class="tocitem" href="../api/">API Reference</a></li></ul><div class="docs-version-selector field has-addons"><div class="control"><span class="docs-label button is-static is-size-7">Version</span></div><div class="docs-selector control is-expanded"><div class="select is-fullwidth is-size-7"><select id="documenter-version-selector"></select></div></div></div></nav><div class="docs-main"><header class="docs-navbar"><nav class="breadcrumb"><ul class="is-hidden-mobile"><li class="is-active"><a href>SAM and BAM</a></li></ul><ul class="is-hidden-tablet"><li class="is-active"><a href>SAM and BAM</a></li></ul></nav><div class="docs-right"><a class="docs-edit-link" href="https://github.com/BioJulia/XAM.jl/blob/develop/docs/src/man/hts-files.md" title="Edit on GitHub"><span class="docs-icon fab"></span><span class="docs-label is-hidden-touch">Edit on GitHub</span></a><a class="docs-settings-button fas fa-cog" id="documenter-settings-button" href="#" title="Settings"></a><a class="docs-sidebar-button fa fa-bars is-hidden-desktop" id="documenter-sidebar-button" href="#"></a></div></header><article class="content" id="documenter-page"><h1 id="SAM-and-BAM-1"><a class="docs-heading-anchor" href="#SAM-and-BAM-1">SAM and BAM</a><a class="docs-heading-anchor-permalink" href="#SAM-and-BAM-1" title="Permalink"></a></h1><h2 id="Introduction-1"><a class="docs-heading-anchor" href="#Introduction-1">Introduction</a><a class="docs-heading-anchor-permalink" href="#Introduction-1" title="Permalink"></a></h2><p>The <code>XAM</code> package offers high-performance tools for SAM and BAM file formats, which are the most popular file formats.</p><p>If you have questions about the SAM and BAM formats or any of the terminology used when discussing these formats, see the published <a href="https://samtools.github.io/hts-specs/SAMv1.pdf">specification</a>, which is maintained by the <a href="https://samtools.github.io/">samtools group</a>.</p><p>A very very simple SAM file looks like the following:</p><pre><code class="language-none">@HD VN:1.6 SO:coordinate
@SQ SN:ref LN:45
r001 99 ref 7 30 8M2I4M1D3M = 37 39 TTAGATAAAGGATACTG *
r002 0 ref 9 30 3S6M1P1I4M * 0 0 AAAAGATAAGGATA *
r003 0 ref 9 30 5S6M * 0 0 GCCTAAGCTAA * SA:Z:ref,29,-,6H5M,17,0;
r004 0 ref 16 30 6M14N5M * 0 0 ATAGCTTCAGC *
r003 2064 ref 29 17 6H5M * 0 0 TAGGC * SA:Z:ref,9,+,5S6M,30,1;
r001 147 ref 37 30 9M = 7 -39 CAGCGGCAT * NM:i:1</code></pre><p>Where the first two lines are part of the &quot;header&quot;, and the following lines are &quot;records&quot;. Each record describes how a read aligns to some reference sequence. Sometimes one record describes one read, but there are other cases like chimeric reads and split alignments, where multiple records apply to one read. In the example above, <code>r003</code> is a chimeric read, and <code>r004</code> is a split alignment, and <code>r001</code> are mate pair reads. Again, we refer you to the official <a href="https://samtools.github.io/hts-specs/SAMv1.pdf">specification</a> for more details.</p><p>A BAM file stores this same information but in a binary and compressible format that does not make for pretty printing here!</p><h2 id="Reading-SAM-and-BAM-files-1"><a class="docs-heading-anchor" href="#Reading-SAM-and-BAM-files-1">Reading SAM and BAM files</a><a class="docs-heading-anchor-permalink" href="#Reading-SAM-and-BAM-files-1" title="Permalink"></a></h2><p>A typical script iterating over all records in a file looks like below:</p><pre><code class="language-julia">using XAM
# Open a BAM file.
reader = open(BAM.Reader, &quot;data.bam&quot;)
# Iterate over BAM records.
for record in reader
# `record` is a BAM.Record object.
if BAM.ismapped(record)
# Print the mapped position.
println(BAM.refname(record), &#39;:&#39;, BAM.position(record))
end
end
# Close the BAM file.
close(reader)</code></pre><p>The size of a BAM file is often extremely large. The iterator interface demonstrated above allocates an object for each record and that may be a bottleneck of reading data from a BAM file. In-place reading reuses a pre-allocated object for every record and less memory allocation happens in reading:</p><pre><code class="language-julia">reader = open(BAM.Reader, &quot;data.bam&quot;)
record = BAM.Record()
while !eof(reader)
empty!(record)
read!(reader, record)
# do something
end</code></pre><h2 id="SAM-and-BAM-Headers-1"><a class="docs-heading-anchor" href="#SAM-and-BAM-Headers-1">SAM and BAM Headers</a><a class="docs-heading-anchor-permalink" href="#SAM-and-BAM-Headers-1" title="Permalink"></a></h2><p>Both <code>SAM.Reader</code> and <code>BAM.Reader</code> implement the <code>header</code> function, which returns a <code>SAM.Header</code> object. To extract certain information out of the headers, you can use the <code>find</code> method on the header to extract information according to SAM/BAM tag. Again we refer you to the <a href="https://samtools.github.io/hts-specs/SAMv1.pdf">specification</a> for full details of all the different tags that can occur in headers, and what they mean.</p><p>Below is an example of extracting all the info about the reference sequences from the BAM header. In SAM/BAM, any description of a reference sequence is stored in the header, under a tag denoted <code>SQ</code> (think <code>reference SeQuence</code>!).</p><pre><code class="language-jlcon">julia&gt; reader = open(SAM.Reader, &quot;data.sam&quot;);
julia&gt; find(header(reader), &quot;SQ&quot;)
7-element Array{Bio.Align.SAM.MetaInfo,1}:
Bio.Align.SAM.MetaInfo:
tag: SQ
value: SN=Chr1 LN=30427671
Bio.Align.SAM.MetaInfo:
tag: SQ
value: SN=Chr2 LN=19698289
Bio.Align.SAM.MetaInfo:
tag: SQ
value: SN=Chr3 LN=23459830
Bio.Align.SAM.MetaInfo:
tag: SQ
value: SN=Chr4 LN=18585056
Bio.Align.SAM.MetaInfo:
tag: SQ
value: SN=Chr5 LN=26975502
Bio.Align.SAM.MetaInfo:
tag: SQ
value: SN=chloroplast LN=154478
Bio.Align.SAM.MetaInfo:
tag: SQ
value: SN=mitochondria LN=366924
</code></pre><p>In the above we can see there were 7 sequences in the reference: 5 chromosomes, one chloroplast sequence, and one mitochondrial sequence.</p><h2 id="SAM-and-BAM-Records-1"><a class="docs-heading-anchor" href="#SAM-and-BAM-Records-1">SAM and BAM Records</a><a class="docs-heading-anchor-permalink" href="#SAM-and-BAM-Records-1" title="Permalink"></a></h2><h3 id="SAM.Record-1"><a class="docs-heading-anchor" href="#SAM.Record-1">SAM.Record</a><a class="docs-heading-anchor-permalink" href="#SAM.Record-1" title="Permalink"></a></h3><p>The <code>XAM</code> package supports the following accessors for <code>SAM.Record</code> types.</p><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.flag" href="#XAM.SAM.flag"><code>XAM.SAM.flag</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">flag(record::Record)::UInt16</code></pre><p>Get the bitwise flag of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L150-L154">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.ismapped" href="#XAM.SAM.ismapped"><code>XAM.SAM.ismapped</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">ismapped(record::Record)::Bool</code></pre><p>Test if <code>record</code> is mapped.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L164-L168">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.isprimary" href="#XAM.SAM.isprimary"><code>XAM.SAM.isprimary</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">isprimary(record::Record)::Bool</code></pre><p>Test if <code>record</code> is a primary line of the read.</p><p>This is equivalent to <code>flag(record) &amp; 0x900 == 0</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L173-L179">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.refname" href="#XAM.SAM.refname"><code>XAM.SAM.refname</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">refname(record::Record)::String</code></pre><p>Get the reference sequence name of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L184-L188">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.position" href="#XAM.SAM.position"><code>XAM.SAM.position</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">position(record::Record)::Int</code></pre><p>Get the 1-based leftmost mapping position of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L201-L205">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.rightposition" href="#XAM.SAM.rightposition"><code>XAM.SAM.rightposition</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">rightposition(record::Record)::Int</code></pre><p>Get the 1-based rightmost mapping position of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L219-L223">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.isnextmapped" href="#XAM.SAM.isnextmapped"><code>XAM.SAM.isnextmapped</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">isnextmapped(record::Record)::Bool</code></pre><p>Test if the mate/next read of <code>record</code> is mapped.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L232-L236">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.nextrefname" href="#XAM.SAM.nextrefname"><code>XAM.SAM.nextrefname</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">nextrefname(record::Record)::String</code></pre><p>Get the reference name of the mate/next read of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L241-L245">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.nextposition" href="#XAM.SAM.nextposition"><code>XAM.SAM.nextposition</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">nextposition(record::Record)::Int</code></pre><p>Get the position of the mate/next read of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L258-L262">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.mappingquality" href="#XAM.SAM.mappingquality"><code>XAM.SAM.mappingquality</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">mappingquality(record::Record)::UInt8</code></pre><p>Get the mapping quality of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L276-L280">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.cigar" href="#XAM.SAM.cigar"><code>XAM.SAM.cigar</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">cigar(record::Record)::String</code></pre><p>Get the CIGAR string of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L294-L298">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.alignment" href="#XAM.SAM.alignment"><code>XAM.SAM.alignment</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">alignment(record::Record)::BioAlignments.Alignment</code></pre><p>Get the alignment of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L311-L315">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.alignlength" href="#XAM.SAM.alignlength"><code>XAM.SAM.alignlength</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">alignlength(record::Record)::Int</code></pre><p>Get the alignment length of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L329-L333">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.tempname" href="#XAM.SAM.tempname"><code>XAM.SAM.tempname</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">tempname(record::Record)::String</code></pre><p>Get the query template name of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L355-L359">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.templength" href="#XAM.SAM.templength"><code>XAM.SAM.templength</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">templength(record::Record)::Int</code></pre><p>Get the template length of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L372-L376">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.sequence" href="#XAM.SAM.sequence"><code>XAM.SAM.sequence</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">sequence(record::Record)::BioSequences.LongDNASeq</code></pre><p>Get the segment sequence of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L390-L394">source</a></section><section><div><pre><code class="language-julia">sequence(::Type{String}, record::Record)::String</code></pre><p>Get the segment sequence of <code>record</code> as <code>String</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L409-L413">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.seqlength" href="#XAM.SAM.seqlength"><code>XAM.SAM.seqlength</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">seqlength(record::Record)::Int</code></pre><p>Get the sequence length of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L419-L423">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.quality" href="#XAM.SAM.quality"><code>XAM.SAM.quality</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">quality(record::Record)::Vector{UInt8}</code></pre><p>Get the Phred-scaled base quality of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L436-L440">source</a></section><section><div><pre><code class="language-julia">quality(::Type{String}, record::Record)::String</code></pre><p>Get the ASCII-encoded base quality of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L457-L461">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.SAM.auxdata" href="#XAM.SAM.auxdata"><code>XAM.SAM.auxdata</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">auxdata(record::Record)::Dict{String,Any}</code></pre><p>Get the auxiliary data (optional fields) of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/sam/record.jl#L467-L471">source</a></section></article><h3 id="BAM.Record-1"><a class="docs-heading-anchor" href="#BAM.Record-1">BAM.Record</a><a class="docs-heading-anchor-permalink" href="#BAM.Record-1" title="Permalink"></a></h3><p>The <code>XAM</code> package supports the following accessors for <code>BAM.Record</code> types.</p><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.flag" href="#XAM.BAM.flag"><code>XAM.BAM.flag</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">flag(record::Record)::UInt16</code></pre><p>Get the bitwise flag of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L131-L135">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.ismapped" href="#XAM.BAM.ismapped"><code>XAM.BAM.ismapped</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">ismapped(record::Record)::Bool</code></pre><p>Test if <code>record</code> is mapped.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L145-L149">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.isprimary" href="#XAM.BAM.isprimary"><code>XAM.BAM.isprimary</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">isprimary(record::Record)::Bool</code></pre><p>Test if <code>record</code> is a primary line of the read.</p><p>This is equivalent to <code>flag(record) &amp; 0x900 == 0</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L154-L160">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.refid" href="#XAM.BAM.refid"><code>XAM.BAM.refid</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">refid(record::Record)::Int</code></pre><p>Get the reference sequence ID of <code>record</code>.</p><p>The ID is 1-based (i.e. the first sequence is 1) and is 0 for a record without a mapping position.</p><p>See also: <code>BAM.rname</code></p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L176-L184">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.refname" href="#XAM.BAM.refname"><code>XAM.BAM.refname</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">refname(record::Record)::String</code></pre><p>Get the reference sequence name of <code>record</code>.</p><p>See also: <code>BAM.refid</code></p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L205-L211">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.reflen" href="#XAM.BAM.reflen"><code>XAM.BAM.reflen</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">reflen(record::Record)::Int</code></pre><p>Get the length of the reference sequence this record applies to.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L218-L222">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.position" href="#XAM.BAM.position"><code>XAM.BAM.position</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">position(record::Record)::Int</code></pre><p>Get the 1-based leftmost mapping position of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L233-L237">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.rightposition" href="#XAM.BAM.rightposition"><code>XAM.BAM.rightposition</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">rightposition(record::Record)::Int</code></pre><p>Get the 1-based rightmost mapping position of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L247-L251">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.isnextmapped" href="#XAM.BAM.isnextmapped"><code>XAM.BAM.isnextmapped</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">isnextmapped(record::Record)::Bool</code></pre><p>Test if the mate/next read of <code>record</code> is mapped.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L261-L265">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.nextrefid" href="#XAM.BAM.nextrefid"><code>XAM.BAM.nextrefid</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">nextrefid(record::Record)::Int</code></pre><p>Get the next/mate reference sequence ID of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L270-L274">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.nextrefname" href="#XAM.BAM.nextrefname"><code>XAM.BAM.nextrefname</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">nextrefname(record::Record)::String</code></pre><p>Get the reference name of the mate/next read of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L284-L288">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.nextposition" href="#XAM.BAM.nextposition"><code>XAM.BAM.nextposition</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">nextposition(record::Record)::Int</code></pre><p>Get the 1-based leftmost mapping position of the next/mate read of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L304-L308">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.mappingquality" href="#XAM.BAM.mappingquality"><code>XAM.BAM.mappingquality</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">mappingquality(record::Record)::UInt8</code></pre><p>Get the mapping quality of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L318-L322">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.cigar" href="#XAM.BAM.cigar"><code>XAM.BAM.cigar</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">cigar(record::Record)::String</code></pre><p>Get the CIGAR string of <code>record</code>.</p><p>Note that in the BAM specification, the field called <code>cigar</code> typically stores the cigar string of the record. However, this is not always true, sometimes the true cigar is very long, and due to some constraints of the BAM format, the actual cigar string is stored in an extra tag: <code>CG:B,I</code>, and the <code>cigar</code> field stores a pseudo-cigar string.</p><p>Calling this method with <code>checkCG</code> set to <code>true</code> (default) this method will always yield the true cigar string, because this is probably what you want the vast majority of the time.</p><p>If you have a record that stores the true cigar in a <code>CG:B,I</code> tag, but you still want to access the pseudo-cigar that is stored in the <code>cigar</code> field of the BAM record, then you can set checkCG to <code>false</code>.</p><p>See also <code>BAM.cigar_rle</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L347-L360">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.alignment" href="#XAM.BAM.alignment"><code>XAM.BAM.alignment</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">alignment(record::Record)::BioAlignments.Alignment</code></pre><p>Get the alignment of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L432-L436">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.alignlength" href="#XAM.BAM.alignlength"><code>XAM.BAM.alignlength</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">alignlength(record::Record)::Int</code></pre><p>Get the alignment length of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L465-L469">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.tempname" href="#XAM.BAM.tempname"><code>XAM.BAM.tempname</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">tempname(record::Record)::String</code></pre><p>Get the query template name of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L483-L487">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.templength" href="#XAM.BAM.templength"><code>XAM.BAM.templength</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">templength(record::Record)::Int</code></pre><p>Get the template length of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L498-L502">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.sequence" href="#XAM.BAM.sequence"><code>XAM.BAM.sequence</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">sequence(record::Record)::BioSequences.LongDNASeq</code></pre><p>Get the segment sequence of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L512-L516">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.seqlength" href="#XAM.BAM.seqlength"><code>XAM.BAM.seqlength</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">seqlength(record::Record)::Int</code></pre><p>Get the sequence length of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L534-L538">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.quality" href="#XAM.BAM.quality"><code>XAM.BAM.quality</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">quality(record::Record)</code></pre><p>Get the base quality of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L548-L552">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="XAM.BAM.auxdata" href="#XAM.BAM.auxdata"><code>XAM.BAM.auxdata</code></a><span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">auxdata(record::Record)::BAM.AuxData</code></pre><p>Get the auxiliary data of <code>record</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/BioJulia/XAM.jl/blob/20dbf800cc2f326d68c719597d38714e782debbf/src/bam/record.jl#L564-L568">source</a></section></article><h2 id="Accessing-auxiliary-data-1"><a class="docs-heading-anchor" href="#Accessing-auxiliary-data-1">Accessing auxiliary data</a><a class="docs-heading-anchor-permalink" href="#Accessing-auxiliary-data-1" title="Permalink"></a></h2><p>SAM and BAM records support the storing of optional data fields associated with tags.</p><p>Tagged auxiliary data follows a format of <code>TAG:TYPE:VALUE</code>. <code>TAG</code> is a two-letter string, and each tag can only appear once per record. <code>TYPE</code> is a single case-sensetive letter which defined the format of <code>VALUE</code>.</p><table><tr><th style="text-align: right">Type</th><th style="text-align: right">Description</th></tr><tr><td style="text-align: right">&#39;A&#39;</td><td style="text-align: right">Printable character</td></tr><tr><td style="text-align: right">&#39;i&#39;</td><td style="text-align: right">Signed integer</td></tr><tr><td style="text-align: right">&#39;f&#39;</td><td style="text-align: right">Single-precision floating number</td></tr><tr><td style="text-align: right">&#39;Z&#39;</td><td style="text-align: right">Printable string, including space</td></tr><tr><td style="text-align: right">&#39;H&#39;</td><td style="text-align: right">Byte array in Hex format</td></tr><tr><td style="text-align: right">&#39;B&#39;</td><td style="text-align: right">Integer of numeric array</td></tr></table><p>For more information about these tags and their types we refer you to the <a href="https://samtools.github.io/hts-specs/SAMv1.pdf">SAM/BAM specification</a> and the additional <a href="https://samtools.github.io/hts-specs/SAMtags.pdf">optional fields specification</a> document.</p><p>There are some tags that are reserved, predefined standard tags, for specific uses.</p><p>To access optional fields stored in tags, you use <code>getindex</code> indexing syntax on the record object. Note that accessing optional tag fields will result in type instability in Julia. This is because the type of the optional data is not known until run-time, as the tag is being read. This can have a significant impact on performance. To limit this, if the user knows the type of a value in advance, specifying it as a type annotation will alleviate the problem:</p><p>Below is an example of looping over records in a bam file and using indexing syntax to get the data stored in the &quot;NM&quot; tag. Note the <code>UInt8</code> type assertion to alleviate type instability.</p><pre><code class="language-julia">for record in open(BAM.Reader, &quot;data.bam&quot;)
nm = record[&quot;NM&quot;]::UInt8
# do something
end</code></pre><h2 id="Getting-records-in-a-range-1"><a class="docs-heading-anchor" href="#Getting-records-in-a-range-1">Getting records in a range</a><a class="docs-heading-anchor-permalink" href="#Getting-records-in-a-range-1" title="Permalink"></a></h2><p>The <code>XAM</code> package supports the BAI index to fetch records in a specific range from a BAM file. <a href="https://samtools.github.io/">Samtools</a> provides <code>index</code> subcommand to create an index file (.bai) from a sorted BAM file.</p><pre><code class="language-console">$ samtools index -b SRR1238088.sort.bam
$ ls SRR1238088.sort.bam*
SRR1238088.sort.bam SRR1238088.sort.bam.bai</code></pre><p>The method <code>eachoverlap(reader, chrom, range)</code> returns an iterator of BAM records overlapping the query interval:</p><pre><code class="language-julia">reader = open(BAM.Reader, &quot;SRR1238088.sort.bam&quot;, index=&quot;SRR1238088.sort.bam.bai&quot;)
for record in eachoverlap(reader, &quot;Chr2&quot;, 10000:11000)
# `record` is a BAM.Record object
# ...
end
close(reader)</code></pre><h2 id="Getting-records-overlapping-genomic-features-1"><a class="docs-heading-anchor" href="#Getting-records-overlapping-genomic-features-1">Getting records overlapping genomic features</a><a class="docs-heading-anchor-permalink" href="#Getting-records-overlapping-genomic-features-1" title="Permalink"></a></h2><p>The <code>eachoverlap</code> method also accepts the <code>Interval</code> type defined in <a href="https://github.com/BioJulia/GenomicFeatures.jl">GenomicFeatures.jl</a>.</p><p>This allows you to do things like first read in the genomic features from a GFF3 file, and then for each feature, iterate over all the BAM records that overlap with that feature.</p><pre><code class="language-julia">using GenomicFeatures
using GFF3
using XAM
# Load genomic features from a GFF3 file.
features = open(collect, GFF3.Reader, &quot;TAIR10_GFF3_genes.gff&quot;)
# Keep mRNA features.
filter!(x -&gt; GFF3.featuretype(x) == &quot;mRNA&quot;, features)
# Open a BAM file and iterate over records overlapping mRNA transcripts.
reader = open(BAM.Reader, &quot;SRR1238088.sort.bam&quot;, index = &quot;SRR1238088.sort.bam.bai&quot;)
for feature in features
for record in eachoverlap(reader, feature)
# `record` overlaps `feature`.
# ...
end
end
close(reader)</code></pre><h2 id="Writing-files-1"><a class="docs-heading-anchor" href="#Writing-files-1">Writing files</a><a class="docs-heading-anchor-permalink" href="#Writing-files-1" title="Permalink"></a></h2><p>In order to write a BAM or SAM file, you must first create a <code>SAM.Header</code>.</p><p>A <code>SAM.Header</code> is constructed from a vector of <code>SAM.MetaInfo</code> objects.</p><p>For example, to create the following simple header:</p><pre><code class="language-none">@HD VN:1.6 SO:coordinate
@SQ SN:ref LN:45</code></pre><pre><code class="language-julia">julia&gt; a = SAM.MetaInfo(&quot;HD&quot;, [&quot;VN&quot; =&gt; 1.6, &quot;SO&quot; =&gt; &quot;coordinate&quot;])
SAM.MetaInfo:
tag: HD
value: VN=1.6 SO=coordinate
julia&gt; b = SAM.MetaInfo(&quot;SQ&quot;, [&quot;SN&quot; =&gt; &quot;ref&quot;, &quot;LN&quot; =&gt; 45])
SAM.MetaInfo:
tag: SQ
value: SN=ref LN=45
julia&gt; h = SAM.Header([a, b])
SAM.Header(SAM.MetaInfo[SAM.MetaInfo:
tag: HD
value: VN=1.6 SO=coordinate, SAM.MetaInfo:
tag: SQ
value: SN=ref LN=45])
</code></pre><p>Then to create the writer for a SAM file, construct a <code>SAM.Writer</code> using the header and an <code>IO</code> type:</p><pre><code class="language-julia">julia&gt; samw = SAM.Writer(open(&quot;my-data.sam&quot;, &quot;w&quot;), h)
SAM.Writer(IOStream(&lt;file my-data.sam&gt;))
</code></pre><p>To make a BAM Writer is slightly different, as you need to use a specific stream type from the <a href="https://github.com/BioJulia/BGZFStreams.jl">https://github.com/BioJulia/BGZFStreams.jl</a> package:</p><pre><code class="language-julia">julia&gt; using BGZFStreams
julia&gt; bamw = BAM.Writer(BGZFStream(open(&quot;my-data.bam&quot;, &quot;w&quot;), &quot;w&quot;))
BAM.Writer(BGZFStreams.BGZFStream{IOStream}(&lt;mode=write&gt;))
</code></pre><p>Once you have a BAM or SAM writer, you can use the <code>write</code> method to write <code>BAM.Record</code>s or <code>SAM.Record</code>s to file:</p><pre><code class="language-julia">julia&gt; write(bamw, rec) # Here rec is a `BAM.Record`
330780</code></pre></article><nav class="docs-footer"><a class="docs-footer-prevpage" href="../../">« Home</a><a class="docs-footer-nextpage" href="../api/">API Reference »</a></nav></div><div class="modal" id="documenter-settings"><div class="modal-background"></div><div class="modal-card"><header class="modal-card-head"><p class="modal-card-title">Settings</p><button class="delete"></button></header><section class="modal-card-body"><p><label class="label">Theme</label><div class="select"><select id="documenter-themepicker"><option value="documenter-light">documenter-light</option><option value="documenter-dark">documenter-dark</option></select></div></p><hr/><p>This document was generated with <a href="https://github.com/JuliaDocs/Documenter.jl">Documenter.jl</a> on <span class="colophon-date" title="Sunday 2 August 2020 04:18">Sunday 2 August 2020</span>. Using Julia version 1.4.2.</p></section><footer class="modal-card-foot"></footer></div></div></div></body></html>