File: libcalignedsegment.pxd

package info (click to toggle)
python-pysam 0.23.3%2Bds-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 18,464 kB
  • sloc: ansic: 158,947; python: 8,806; makefile: 264; sh: 79
file content (116 lines) | stat: -rw-r--r-- 3,246 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
# cython: language_level=3
from pysam.libchtslib cimport *

cdef extern from "htslib_util.h":

    # add *nbytes* into the variable length data of *src* at *pos*
    bam1_t * pysam_bam_update(bam1_t * b,
           	         size_t nbytes_old,
                         size_t nbytes_new,
                         uint8_t * pos)

    # now: static
    int aux_type2size(int)

    char * pysam_bam_get_qname(bam1_t * b)
    uint32_t * pysam_bam_get_cigar(bam1_t * b)
    uint8_t * pysam_bam_get_seq(bam1_t * b)
    uint8_t * pysam_bam_get_qual(bam1_t * b)
    uint8_t * pysam_bam_get_aux(bam1_t * b)
    int pysam_bam_get_l_aux(bam1_t * b)
    char pysam_bam_seqi(uint8_t * s, int i)

    uint8_t pysam_get_qual(bam1_t * b)
    uint32_t pysam_get_n_cigar(bam1_t * b)
    void pysam_set_qual(bam1_t * b, uint8_t v)
    void pysam_set_n_cigar(bam1_t * b, uint32_t v)
    void pysam_update_flag(bam1_t * b, uint16_t v, uint16_t flag)


from pysam.libcalignmentfile cimport AlignmentFile, AlignmentHeader
ctypedef AlignmentFile AlignmentFile_t


cdef class _AlignedSegment_Cache:  # For internal use only
    cdef clear_query_sequences(self)
    cdef clear_query_qualities(self)

    cdef object query_sequence
    cdef object query_alignment_sequence
    cdef object query_qualities
    cdef object query_qualities_str
    cdef object query_alignment_qualities
    cdef object query_alignment_qualities_str


# Note: need to declare all C fields and methods here
cdef class AlignedSegment:

    # object that this AlignedSegment represents
    cdef bam1_t * _delegate

    # the header that a read is associated with
    cdef readonly AlignmentHeader header

    # caching of array properties for quick access
    cdef _AlignedSegment_Cache cache

    cdef object unused1
    cdef object unused2
    cdef object unused3

    # add an alignment tag with value to the AlignedSegment
    # an existing tag of the same name will be replaced.
    cpdef set_tag(self, tag, value, value_type=?, replace=?)

    # get an alignment tag from the AlignedSegment
    cpdef get_tag(self, tag, with_value_type=?)

    # return true if tag exists
    cpdef has_tag(self, tag)

    # returns a valid sam alignment string
    cpdef to_string(self)

    # returns a valid sam alignment string (deprecated)
    cpdef tostring(self, htsfile=*)


cdef class PileupColumn:
    cdef const bam_pileup1_t ** plp
    cdef int tid
    cdef int pos
    cdef int n_pu
    cdef AlignmentHeader header
    cdef uint32_t min_base_quality
    cdef kstring_t buf
    cdef char * reference_sequence

cdef class PileupRead:
    cdef int32_t  _qpos
    cdef AlignedSegment _alignment
    cdef int _indel
    cdef int _level
    cdef uint32_t _is_del
    cdef uint32_t _is_head
    cdef uint32_t _is_tail
    cdef uint32_t _is_refskip

# factory methods
cdef AlignedSegment makeAlignedSegment(
    bam1_t * src,
    AlignmentHeader header)

cdef PileupColumn makePileupColumn(
    const bam_pileup1_t ** plp,
    int tid,
    int pos,
    int n_pu,
    uint32_t min_base_quality,
    char * reference_sequence,
    AlignmentHeader header)

cdef PileupRead makePileupRead(const bam_pileup1_t * src,
		               AlignmentHeader header)

cdef uint32_t get_alignment_length(bam1_t * src)