BamWriter

Class for outputting BAM.
Compresses BGZF blocks in parallel. Tries to write reads so that they don't cross BGZF block borders.
Usage is very simple, see example below.

Constructors

this
this(contrib.undead.stream.Stream stream, int compression_level, std.parallelism.TaskPool task_pool, size_t buffer_size)
this(string filename, int compression_level, std.parallelism.TaskPool task_pool)

Creates new BAM writer outputting to file or stream. Automatically writes BAM magic number (4 bytes).

Members

Functions

disableAutoIndexCreation
void disableAutoIndexCreation()

By default, the writer attempts to automatically create index when writing coordinate-sorted files. If this behaviour is not desired, it can be switched off before writing SAM header.

finish
void finish()

Flushes buffer and closes output stream. Adds BAM EOF block automatically.

flush
void flush()

Flushes current BGZF block.

setFilename
void setFilename(string output_filename)

Can be called right after the stream constructor, only once

writeByteArray
void writeByteArray(const(ubyte[]) array)
Undocumented in source. Be warned that the author may not have intended to support it.
writeInteger
void writeInteger(T integer)
Undocumented in source. Be warned that the author may not have intended to support it.
writeRecord
void writeRecord(R read)

Writes BAM read. Throws exception if read reference ID is out of range.

writeReferenceSequenceInfo
void writeReferenceSequenceInfo(const(bio.std.hts.bam.referenceinfo.ReferenceSequenceInfo)[] reference_sequences)

Writes reference sequence information. Should be called after dumping SAM header. Writer will store this array to use later for resolving read reference IDs to names.

writeSamHeader
void writeSamHeader(bio.std.hts.sam.header.SamHeader header)
void writeSamHeader(string header_text)

Writes SAM header. Should be called after construction.

writeString
void writeString(string str)
Undocumented in source. Be warned that the author may not have intended to support it.

Examples

import bio.std.hts.bam.writer, bio.std.hts.bam.reader;
...
auto src = new BamReader("in.bam");
auto dst = new BamWriter("out.bam", 9); // maximal compression
scope (exit) dst.finish();              // close the stream at exit
dst.writeSamHeader(src.header);         // copy header and reference sequence info
dst.writeReferenceSequenceInfo(src.reference_sequences);
foreach (read; src.reads) {
    if (read.mapping_quality > 10)      // skip low-quality reads
        dst.writeRecord(read);
}

Meta