File: libcalignedsegment.pxd

package info (click to toggle)
python-pysam 0.10.0%2Bds-2
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 14,196 kB
  • ctags: 10,087
  • sloc: ansic: 79,627; python: 8,569; sh: 282; makefile: 215; perl: 41
file content (91 lines) | stat: -rw-r--r-- 2,977 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
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)

    uint16_t pysam_get_bin(bam1_t * b)
    uint8_t pysam_get_qual(bam1_t * b)
    uint8_t pysam_get_l_qname(bam1_t * b)
    uint16_t pysam_get_flag(bam1_t * b)
    uint16_t pysam_get_n_cigar(bam1_t * b)
    void pysam_set_bin(bam1_t * b, uint16_t v)
    void pysam_set_qual(bam1_t * b, uint8_t v)
    void pysam_set_l_qname(bam1_t * b, uint8_t v)
    void pysam_set_flag(bam1_t * b, uint16_t v)
    void pysam_set_n_cigar(bam1_t * b, uint16_t v)
    void pysam_update_flag(bam1_t * b, uint16_t v, uint16_t flag)


from pysam.libcalignmentfile cimport AlignmentFile
ctypedef AlignmentFile AlignmentFile_t


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

    # object that this AlignedSegment represents
    cdef bam1_t * _delegate

    # the file from which this AlignedSegment originates (can be None)
    cdef AlignmentFile _alignment_file

    # caching of array properties for quick access
    cdef object cache_query_qualities
    cdef object cache_query_alignment_qualities
    cdef object cache_query_sequence
    cdef object cache_query_alignment_sequence

    # 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=?)

    # add an alignment tag with value to the AlignedSegment
    # an existing tag of the same name will be replaced.
    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 tostring(self, AlignmentFile_t handle)


cdef class PileupColumn:
    cdef bam_pileup1_t ** plp
    cdef int tid
    cdef int pos
    cdef int n_pu
    cdef AlignmentFile _alignment_file


cdef class PileupRead:
    cdef AlignedSegment _alignment
    cdef int32_t  _qpos
    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

# factor methods
cdef makeAlignedSegment(bam1_t * src, AlignmentFile alignment_file)
cdef makePileupColumn(bam_pileup1_t ** plp, int tid, int pos, int n_pu, AlignmentFile alignment_file)
cdef inline makePileupRead(bam_pileup1_t * src, AlignmentFile alignment_file)
cdef inline uint32_t get_alignment_length(bam1_t * src)