hs-samtools-0.9.0.0: Read and write SAM, BAM, and CRAM files.
Copyright(c) Matthew Mosior 2023
LicenseBSD-style
Maintainermattm.github@gmail.com
Portabilityportable
Safe HaskellSafe-Inferred
LanguageHaskell2010

Data.SAM.Version1_6.Alignment.Base

Description

Description

This library enables the decoding/encoding of SAM, BAM and CRAM file formats.

Synopsis

SAM version 1.6 alignment mandatory and optional data types

data SAM_V1_6_Alignment Source #

Custom SAM (version 1.6) SAM_V1_6_Alignment data type.

See section 1.4 and 1.5 of the SAM v1.6 specification documentation.

Constructors

SAM_V1_6_Alignment 

Fields

  • sam_v1_6_alignment_qname :: ByteString

    Query template NAME. reads/segments having identical QNAME are regarded to come from the same template. A QNAME ‘*’ indicates the information is unavailable. In a SAM file, a read may occupy multiple alignment lines, when its alignment is chimeric or when multiple mappings are given.

  • sam_v1_6_alignment_flag :: Int

    Combination of bitwise FLAGs.

  • sam_v1_6_alignment_rname :: ByteString

    Reference sequence NAME of the alignment. If @SQ header lines are present, RNAME (if not ‘*’) must be present in one of the SQ-SN tag. An unmapped segment without coordinate has a ‘*’ at this field. However, an unmapped segment may also have an ordinary coordinate such that it can be placed at a desired position after sorting. If RNAME is ‘*’, no assumptions can be made about POS and CIGAR.

  • sam_v1_6_alignment_pos :: Integer

    1-based leftmost mapping POSition of the first CIGAR operation that “consumes” a reference base. The first base in a reference sequence has coordinate 1. POS is set as 0 for an unmapped read without coordinate. If POS is 0, no assumptions can be made about RNAME and CIGAR.

  • sam_v1_6_alignment_mapq :: Int

    MAPping Quality. It equals −10 log10 Pr{mapping position is wrong}, rounded to the nearest integer. A value 255 indicates that the mapping quality is not available.

  • sam_v1_6_alignment_cigar :: ByteString

    CIGAR string (set ‘*’ if unavailable).

  • sam_v1_6_alignment_rnext :: ByteString

    Reference sequence name of the primary alignment of the NEXT read in the template. For the last read, the next read is the first read in the template. If @SQ header lines are present, RNEXT (if not ‘*’ or ‘=’) must be present in one of the SQ-SN tag. This field is set as ‘*’ when the information is unavailable, and set as ‘=’ if RNEXT is identical RNAME. If not ‘=’ and the next read in the template has one primary mapping (see also bit 0x100 in FLAG), this field is identical to RNAME at the primary line of the next read. If RNEXT is ‘*’, no assumptions can be made on PNEXT and bit 0x20.

  • sam_v1_6_alignment_pnext :: Integer

    1-based Position of the primary alignment of the NEXT read in the template. Set as 0 when the information is unavailable. This field equals POS at the primary line of the next read. If PNEXT is 0, no assumptions can be made on RNEXT and bit 0x20.

  • sam_v1_6_alignment_tlen :: Integer

    signed observed Template LENgth. For primary reads where the primary alignments of all reads in the template are mapped to the same reference sequence, the absolute value of TLEN equals the distance between the mapped end of the template and the mapped start of the template, inclusively (i.e., end − start + 1). Note that mapped base is defined to be one that aligns to the reference as described by CIGAR, hence excludes soft-clipped bases. The TLEN field is positive for the leftmost segment of the template, negative for the rightmost, and the sign for any middle segment is undefined. If segments cover the same coordinates then the choice of which is leftmost and rightmost is arbitrary, but the two ends must still have differing signs. It is set as 0 for a single-segment template or when the information is unavailable (e.g., when the first or last segment of a multi-segment template is unmapped or when the two are mapped to different reference sequences).

  • sam_v1_6_alignment_seq :: ByteString

    segment SEQuence. This field can be a ‘*’ when the sequence is not stored. If not a ‘*’, the length of the sequence must equal the sum of lengths of MIS=X operations in CIGAR. An ‘=’ denotes the base is identical to the reference base. No assumptions can be made on the letter cases.

  • sam_v1_6_alignment_qual :: ByteString

    ASCII of base QUALity plus 33 (same as the quality string in the Sanger FASTQ format). A base quality is the phred-scaled base error probability which equals −10 log10 Pr{base is wrong}. This field can be a ‘*’ when quality is not stored. If not a ‘*’, SEQ must not be a ‘*’ and the length of the quality string ought to equal the length of SEQ.

  • sam_v1_6_alignment_aopt :: Maybe SAM_V1_6_Alignment_AOPT

    A - [!-~] - Printable characters.

  • sam_v1_6_alignment_iopt :: Maybe SAM_V1_6_Alignment_IOPT

    i - [-+]?[0-9]+ - Signed integer.

  • sam_v1_6_alignment_fopt :: Maybe SAM_V1_6_Alignment_FOPT

    f - [-+]?[0-9]*.?[0-9]+([eE][-+]?[0-9]+)? - Single-precision floating number.

  • sam_v1_6_alignment_zopt :: Maybe SAM_V1_6_Alignment_ZOPT

    Z - [ !-~]* - Printable string, including space.

  • sam_v1_6_alignment_hopt :: Maybe SAM_V1_6_Alignment_HOPT

    H - ([0-9A-F][0-9A-F])* - Byte array in the Hex format.

  • sam_v1_6_alignment_bopt :: Maybe SAM_V1_6_Alignment_BOPT

    B - [cCsSiIf]​(,[-+]?[0-9]*.?[0-9]+([eE][-+]?[0-9]+)?)* - Integer or numeric array.

Instances

Instances details
Generic SAM_V1_6_Alignment Source # 
Instance details

Defined in Data.SAM.Version1_6.Alignment.Base

Associated Types

type Rep SAM_V1_6_Alignment :: Type -> Type #

Show SAM_V1_6_Alignment Source # 
Instance details

Defined in Data.SAM.Version1_6.Alignment.Base

Eq SAM_V1_6_Alignment Source # 
Instance details

Defined in Data.SAM.Version1_6.Alignment.Base

type Rep SAM_V1_6_Alignment Source # 
Instance details

Defined in Data.SAM.Version1_6.Alignment.Base

type Rep SAM_V1_6_Alignment = D1 ('MetaData "SAM_V1_6_Alignment" "Data.SAM.Version1_6.Alignment.Base" "hs-samtools-0.9.0.0-inplace" 'False) (C1 ('MetaCons "SAM_V1_6_Alignment" 'PrefixI 'True) ((((S1 ('MetaSel ('Just "sam_v1_6_alignment_qname") 'NoSourceUnpackedness 'NoSourceStrictness 'DecidedLazy) (Rec0 ByteString) :*: S1 ('MetaSel ('Just "sam_v1_6_alignment_flag") 'NoSourceUnpackedness 'NoSourceStrictness 'DecidedLazy) (Rec0 Int)) :*: (S1 ('MetaSel ('Just "sam_v1_6_alignment_rname") 'NoSourceUnpackedness 'NoSourceStrictness 'DecidedLazy) (Rec0 ByteString) :*: S1 ('MetaSel ('Just "sam_v1_6_alignment_pos") 'NoSourceUnpackedness 'NoSourceStrictness 'DecidedLazy) (Rec0 Integer))) :*: ((S1 ('MetaSel ('Just "sam_v1_6_alignment_mapq") 'NoSourceUnpackedness 'NoSourceStrictness 'DecidedLazy) (Rec0 Int) :*: S1 ('MetaSel ('Just "sam_v1_6_alignment_cigar") 'NoSourceUnpackedness 'NoSourceStrictness 'DecidedLazy) (Rec0 ByteString)) :*: (S1 ('MetaSel ('Just "sam_v1_6_alignment_rnext") 'NoSourceUnpackedness 'NoSourceStrictness 'DecidedLazy) (Rec0 ByteString) :*: S1 ('MetaSel ('Just "sam_v1_6_alignment_pnext") 'NoSourceUnpackedness 'NoSourceStrictness 'DecidedLazy) (Rec0 Integer)))) :*: (((S1 ('MetaSel ('Just "sam_v1_6_alignment_tlen") 'NoSourceUnpackedness 'NoSourceStrictness 'DecidedLazy) (Rec0 Integer) :*: S1 ('MetaSel ('Just "sam_v1_6_alignment_seq") 'NoSourceUnpackedness 'NoSourceStrictness 'DecidedLazy) (Rec0 ByteString)) :*: (S1 ('MetaSel ('Just "sam_v1_6_alignment_qual") 'NoSourceUnpackedness 'NoSourceStrictness 'DecidedLazy) (Rec0 ByteString) :*: S1 ('MetaSel ('Just "sam_v1_6_alignment_aopt") 'NoSourceUnpackedness 'NoSourceStrictness 'DecidedLazy) (Rec0 (Maybe SAM_V1_6_Alignment_AOPT)))) :*: ((S1 ('MetaSel ('Just "sam_v1_6_alignment_iopt") 'NoSourceUnpackedness 'NoSourceStrictness 'DecidedLazy) (Rec0 (Maybe SAM_V1_6_Alignment_IOPT)) :*: S1 ('MetaSel ('Just "sam_v1_6_alignment_fopt") 'NoSourceUnpackedness 'NoSourceStrictness 'DecidedLazy) (Rec0 (Maybe SAM_V1_6_Alignment_FOPT))) :*: (S1 ('MetaSel ('Just "sam_v1_6_alignment_zopt") 'NoSourceUnpackedness 'NoSourceStrictness 'DecidedLazy) (Rec0 (Maybe SAM_V1_6_Alignment_ZOPT)) :*: (S1 ('MetaSel ('Just "sam_v1_6_alignment_hopt") 'NoSourceUnpackedness 'NoSourceStrictness 'DecidedLazy) (Rec0 (Maybe SAM_V1_6_Alignment_HOPT)) :*: S1 ('MetaSel ('Just "sam_v1_6_alignment_bopt") 'NoSourceUnpackedness 'NoSourceStrictness 'DecidedLazy) (Rec0 (Maybe SAM_V1_6_Alignment_BOPT))))))))