oracular (3) Bio::DB::HTS::Alignment.3pm.gz

Provided by: libbio-db-hts-perl_3.01-4build3_amd64 bug

LICENSE

       Copyright [2015-2018] EMBL-European Bioinformatics Institute

       Licensed under the Apache License, Version 2.0 (the "License"); you may not use this file except in
       compliance with the License.  You may obtain a copy of the License at

            http://www.apache.org/licenses/LICENSE-2.0

       Unless required by applicable law or agreed to in writing, software distributed under the License is
       distributed on an "AS IS" BASIS, WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
       See the License for the specific language governing permissions and limitations under the License.

NAME

       Bio::DB::HTS::Alignment -- The HTS alignment object

SYNOPSIS

        use Bio::DB::HTS;

        my $sam = Bio::DB::HTS->new(-fasta=>"data/ex1.fa",
                                    -bam  =>"data/ex1.bam");

        my @alignments = $sam->get_features_by_location(-seq_id => 'seq2',
                                                        -start  => 500,
                                                        -end    => 800);
        for my $a (@alignments) {
           my $seqid  = $a->seq_id;
           my $start  = $a->start;
           my $end    = $a->end;
           my $strand = $a->strand;
           my $ref_dna= $a->dna;

           my $query_start  = $a->query->start;
           my $query_end    = $a->query->end;
           my $query_strand = $a->query->strand;
           my $query_dna    = $a->query->dna;

           my $cigar     = $a->cigar_str;
           my @scores    = $a->qscore;     # per-base quality scores
           my $score     = $a->qstring;    # TAM-style quality string
           my $match_qual= $a->qual;       # quality of the match

           my $paired = $a->get_tag_values('PAIRED');
        }

DESCRIPTION

       The Bio::DB::HTS::Alignment and Bio::DB::HTS::AlignWrapper classes together represent an alignment
       between a sequence read (the "query") and a reference sequence (the "target"). Bio::DB::HTS::Alignment
       adheres strictly to the C-level BAM library's definition of a bam1_t* and is used in the Bio::DB::HTS
       low-level API The latter adds convenience methods that make it similar to a BioPerl Bio::SeqFeatureI
       object. This manual page describes both.

High-level Bio::DB::HTS::Alignment methods

       These methods are provided by Bio::DB::HTS::Alignment, and are intended to be compatible with the
       Bio::SeqFeatureI interfaces. Note that these objects are not compatible with Bio::Align::AlignI, as the
       BAM API is fundamentally incompatible with the BioPerl API for alignments (the first deals with the
       alignment of a single read against the reference sequence, while the second deals with a multiple
       alignment).

       Note that the high-level API return Bio::DB::HTS::AlignWrapper objects except in the case of the callback
       to the fast_pileup() method. In this case only, the object returned by calling $pileup->b() is a
       Bio::DB::HTS::Alignment object for performance reasons.

       $seq_id = $align->seq_id
           Return the seq_id of the reference (target) sequence. This method is only available in the
           Bio::DB::HTS::AlignWrapper extension.

       $start = $align->start
           Return the start of the alignment in 1-based reference sequence coordinates.

       $end = $align->end
           Return the end of the alignment in 1-based reference sequence coordinates.

       $len = $align->length
           Return the length of the alignment on the reference sequence.

       $mseqid = $align->mate_seq_id
           Return the seq_id of the mate's reference (target) sequence. This method is only available in the
           Bio::DB::AlignWrapper extension.

       $mstart = $align->mate_start
           For paired reads, return the start of the mate's alignment in 1-based reference sequence coordinates.

       $mend = $align->mate_end
           For paired reads, return the end position of the mate's alignment in 1-based reference sequence
           coordinates.

       $mlen = $align->mate_len
           For mate-pairs, retrieve the length of the mate's alignment on the reference sequence.

       $strand = $align->strand
           Return the strand of the alignment as -1 for reversed, +1 for forward.

       $mstrand = $align->mstrand
           If the read has a mate pair, return the strand of the mate in the format -1 or +1.

       $ref_dna        = $align->dna
           Returns the reference sequence's DNA across the aligned region. If an MD tag is present in the
           alignment, it will be used preferentially to reconstruct the reference sequence. Otherwise the
           reference DNA access object passed to Bio::DB::HTS->new() will be used.

       $ref_dna        = $align->seq
           The reference sequence's DNA as a Bio::PrimarySeqI object (useful for passing to BioPerl functions
           and for calculating subsequences and reverse complements).

       $query = $align->query
           This method returns a Bio::DB::Alignment::Query object that can be used to retrieve information about
           the query sequence. The next few entries show how to use this object.

       $read_name = $align->query->name
           The name of the read.

       $q_start   = $align->query->start
           This returns the start position of the query (read) sequence in 1-based coordinates. It acts via a
           transient Bio::DB::HTS::Query object that is provided for Bio::Graphics compatibility (see
           Bio::Graphics).

       $q_end     = $align->query->end
           This returns the end position of the query sequence in 1-based coordinates.

       $q_len     = $align->query->length
           Return the length of the alignment on the read.

       $scores = $align->query->score
           Return an array reference containing the unpacked quality scores for each base of the query sequence.
           The length of this array reference will be equal to the length of the read.

       $read_dna = $align->query->dna
           The read's DNA string.

       $read_seq = $align->query->seq
           The read's DNA as a Bio::PrimarySeqI object.

       $target  = $align->target;
           The target() method is similar to query(), except that it follows Bio::AlignIO conventions for how to
           represent minus strand alignments. The object returned has start(), end(), qscore(), dna() and seq()
           methods, but for minus strand alignments the sequence will be represented as it appears on the
           reverse strand, rather than on the forward strand. This has the advantage of giving you the read as
           it came off the machine, before being reverse complemented for use in the SAM file.

       $query   = $align->hit
           The hit() method is identical to target() and returns information about the read. It is present for
           compatibility with some of the Bio::Graphics glyphs, which use hit() to represent the non-reference
           sequence in aligned sequences.

       $primary_id = $align->primary_id
           This method synthesizes a unique ID for the alignment which can be passed to
           $sam->get_feature_by_id() to retrieve the alignment at a later date.

       @tags = $align->get_all_tags
           Return all tag names known to this alignment. This includes SAM flags such as M_UNMAPPED, as well as
           auxiliary flags such as H0. The behavior of this method depends on the value of -expand_flags when
           the SAM object was created. If false (the default), then the standard SAM flags will be concatenated
           together into a single string and stored in a tag named 'FLAGS'. The format of this tag value is the
           list of one or more flag constants separated by the "|" character, as in:
           "PAIRED|MAP_PAIR|REVERSED|SECOND_MATE". If -expand_flags was true, then each flag becomes its own
           named tag, such as "MAP_PAIR".

       @values = $align->get_tag_values($tag)
           Given a tag name, such as 'PAIRED' or 'H0', return its value(s). -expand_flags must be true in order
           to use the standard SAM flag constants as tags. Otherwise, they can be fetched by asking for the
           "FLAGS" tag, or by using the low-level methods described below.

       $is_true = $align->has_tag($tag)
           Return true if the alignment has the indicated tag.

       $string = $align->cigar_str
           Return the CIGAR string for this alignment in conventional human readable format (e.g. "M34D1M1").

       $arrayref = $align->cigar_array
           Return a reference to an array representing the CIGAR string. This is an array of arrays, in which
           each subarray consists of a CIGAR operation and a count. Example:

            [ ['M',34], ['D',1], ['M',1] ]

       ($ref,$matches,$query) = $align->padded_alignment
           Return three strings that show the alignment between the reference sequence (the target) and the
           query. It will look like this:

            $ref     AGTGCCTTTGTTCA-----ACCCCCTTGCAACAACC
            $matches ||||||||||||||     |||||||||||||||||
            $query   AGTGCCTTTGTTCACATAGACCCCCTTGCAACAACC

       $str = $align->aux
           Returns the text version of the SAM tags, e.g.  "XT:A:M NM:i:2 SM:i:37 AM:i:37 XM:i:1 XO:i:1 XG:i:1
           MD:Z:6^C0A47"

       $str = $align->tam_line
           Returns the TAM (text) representation of the alignment (available in the high-level "AlignWrapper"
           interface only).

       $tag = $align->primary_tag
           This is provided for Bio::SeqFeatureI compatibility. Return the string "match".

       $tag = $align->source_tag
           This is provided for Bio::SeqFeatureI compatibility. Return the string "sam/bam".

       @parts = $align->get_SeqFeatures
           Return subfeatures of this alignment. If you have fetched a "read_pair" feature, this will be the two
           mate pair objects (both of type Bio::DB::HTS::AlignWrapper). If you have -split_splices set to true
           in the Bio::DB::HTS database, calling get_SeqFeatures() will return the components of split
           alignments. See "Bio::DB::HTS Constructor and basic accessors" in Bio::DB::HTS for an example of how
           to use this.

Low-level Bio::DB::HTS::Alignment methods

       These methods are available to objects of type Bio::DB::HTS::Alignment as well as
       Bio::DB::HTS::AlignWrapper and closely mirror the native C API.

       $align = Bio::DB::HTS::Alignment->new
           Create a new, empty alignment object. This is usually needed when iterating through a HTS file using
           Bio::DB::HTS->read1().

       $tid = $align->tid( [$new_tid] )
           Return the target ID of the alignment. Optionally you may change the tid by providing it as an
           argument (currently this is the only field that you can change; the functionality was implemented as
           a proof of principle).

       $read_name = $align->qname
           Returns the name of the read.

       $pos = $align->pos
           0-based leftmost coordinate of the aligned sequence on the reference sequence.

       $end = $align->calend
           The 0-based rightmost coordinate of the aligned sequence on the reference sequence after taking
           alignment gaps into account.

       $len = $align->cigar2qlen
           The length of the query sequence calculated from the CIGAR string.

       $quality = $align->qual
           The quality score for the alignment as a whole.

       $flag = $align->flag
           The bitwise flag field (see the SAM documentation).

       $mtid = $align->mtid
           For paired reads, the target ID of the mate's alignemnt.

       $mpos = $align->mpos
           For paired reads, the 0-based leftmost coordinate of the mate's alignment on the reference sequence.

       $n_cigar = $align->n_cigar
           Number of CIGAR operations in this alignment.

       $length = $align->l_qseq
           The length of the query sequence (the read).

       $dna = $align->qseq
           The actual DNA sequence of the query. As in the SAM file, reads that are aligned to the minus strand
           of the reference are returned in reverse complemented form.

       $score_str = $align->_qscore
           A packed binary string containing the quality scores for each base of the read. It will be the same
           length as the DNA. You may unpack it using unpack('C*',$score_str), or use the high-level qscore()
           method.

       $score_arry = $align->qscore
       @score_arry = $align->qscore
           In a scalar context return an array reference containing the unpacked quality scores for each base of
           the query sequence. In a list context return a list of the scores. This array is in the same
           orientation as the reference sequence.

       $score_str = $align->qstring
           Returns the quality string in the same format used in the SAM (TAM) file.

       $length = $align->isize
           The calculated insert size for mapped paired reads.

       $length = $align->l_aux
           The length of the align "auxiliary" data.

       $value = $align->aux_get("tag")
           Given an auxiliary tag, such as "H0", return its value.

       @keys  = $align->aux_keys
           Return the list of auxiliary tags known to this alignment.

       $data = $align->data
           Return a packed string containing the alignment data (sequence, quality scores and cigar string).

       $length = $align->data_len
           Return the current length of the alignment data.

       $length = $align->m_data
           Return the maximum length of the alignment data.

       $is_paired = $align->paired
           Return true if the aligned read is part of a mate/read pair (regardless of whether the mate mapped).

       $is_proper = $align->proper_pair
           Return true if the aligned read is part of a mate/read pair and both partners mapped to the reference
           sequence.

       $is_unmapped = $align->unmapped
           Return true if the read failed to align.

       $mate_is_unmapped = $align->munmapped
           Return true if the read's mate failed to align.

       $reversed = $align->reversed
           Return true if the aligned read was reverse complemented prior to aligning.

       $mate_reversed = $align->mreversed
           Return true if the aligned read's mate was reverse complemented prior to aligning.

       $isize = $align->isize
           For mate-pairs, return the computed insert size.

       $arrayref = $align->cigar
           This returns the CIGAR data in its native BAM format. You will receive an arrayref in which each
           operation and count are packed together into an 8-bit structure. To decode each element you must use
           the following operations:

            use Bio::DB::HTS::Constants;
            my $c   = $align->cigar;
            my $op  = $c->[0] & BAM_CIGAR_MASK;
            my $len = $c->[0] >> BAM_CIGAR_SHIFT;

AUTHOR

       Rishi Nag <rishi@ebi.ac.uk<gt>

SEE ALSO

       Bio::Perl, Bio::DB::HTS, Bio::DB::HTS::Constants