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
|
<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
<html xmlns="http://www.w3.org/1999/xhtml">
<head>
<meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/>
<meta http-equiv="X-UA-Compatible" content="IE=9"/>
<meta name="generator" content="Doxygen 1.8.14"/>
<meta name="viewport" content="width=device-width, initial-scale=1"/>
<title>iVar: src/trim_primer_quality.cpp Source File</title>
<link href="tabs.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="jquery.js"></script>
<script type="text/javascript" src="dynsections.js"></script>
<link href="search/search.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="search/searchdata.js"></script>
<script type="text/javascript" src="search/search.js"></script>
<link href="doxygen.css" rel="stylesheet" type="text/css" />
</head>
<body>
<div id="top"><!-- do not remove this div, it is closed by doxygen! -->
<div id="titlearea">
<table cellspacing="0" cellpadding="0">
<tbody>
<tr style="height: 56px;">
<td id="projectalign" style="padding-left: 0.5em;">
<div id="projectname">iVar
</div>
</td>
</tr>
</tbody>
</table>
</div>
<!-- end header part -->
<!-- Generated by Doxygen 1.8.14 -->
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:cf05388f2679ee054f2beb29a391d25f4e673ac3&dn=gpl-2.0.txt GPL-v2 */
var searchBox = new SearchBox("searchBox", "search",false,'Search');
/* @license-end */
</script>
<script type="text/javascript" src="menudata.js"></script>
<script type="text/javascript" src="menu.js"></script>
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:cf05388f2679ee054f2beb29a391d25f4e673ac3&dn=gpl-2.0.txt GPL-v2 */
$(function() {
initMenu('',true,false,'search.php','Search');
$(document).ready(function() { init_search(); });
});
/* @license-end */</script>
<div id="main-nav"></div>
<!-- window showing the filter options -->
<div id="MSearchSelectWindow"
onmouseover="return searchBox.OnSearchSelectShow()"
onmouseout="return searchBox.OnSearchSelectHide()"
onkeydown="return searchBox.OnSearchSelectKey(event)">
</div>
<!-- iframe showing the search results (closed by default) -->
<div id="MSearchResultsWindow">
<iframe src="javascript:void(0)" frameborder="0"
name="MSearchResults" id="MSearchResults">
</iframe>
</div>
<div id="nav-path" class="navpath">
<ul>
<li class="navelem"><a class="el" href="dir_68267d1309a1af8e8297ef4c3efbcdba.html">src</a></li> </ul>
</div>
</div><!-- top -->
<div class="header">
<div class="headertitle">
<div class="title">trim_primer_quality.cpp</div> </div>
</div><!--header-->
<div class="contents">
<div class="fragment"><div class="line"><a name="l00001"></a><span class="lineno"> 1</span> <span class="preprocessor">#include "htslib/sam.h"</span></div><div class="line"><a name="l00002"></a><span class="lineno"> 2</span> <span class="preprocessor">#include "htslib/bgzf.h"</span></div><div class="line"><a name="l00003"></a><span class="lineno"> 3</span> </div><div class="line"><a name="l00004"></a><span class="lineno"> 4</span> <span class="preprocessor">#include <iostream></span></div><div class="line"><a name="l00005"></a><span class="lineno"> 5</span> <span class="preprocessor">#include <stdexcept></span></div><div class="line"><a name="l00006"></a><span class="lineno"> 6</span> <span class="preprocessor">#include <string></span></div><div class="line"><a name="l00007"></a><span class="lineno"> 7</span> <span class="preprocessor">#include <vector></span></div><div class="line"><a name="l00008"></a><span class="lineno"> 8</span> <span class="preprocessor">#include <fstream></span></div><div class="line"><a name="l00009"></a><span class="lineno"> 9</span> <span class="preprocessor">#include <sstream></span></div><div class="line"><a name="l00010"></a><span class="lineno"> 10</span> <span class="preprocessor">#include <cstring></span></div><div class="line"><a name="l00011"></a><span class="lineno"> 11</span> </div><div class="line"><a name="l00012"></a><span class="lineno"> 12</span> <span class="preprocessor">#include "trim_primer_quality.h"</span></div><div class="line"><a name="l00013"></a><span class="lineno"> 13</span> <span class="preprocessor">#include "primer_bed.h"</span></div><div class="line"><a name="l00014"></a><span class="lineno"> 14</span> </div><div class="line"><a name="l00015"></a><span class="lineno"> 15</span> int32_t get_pos_on_query(uint32_t *cigar, uint32_t ncigar, int32_t pos, int32_t ref_start){</div><div class="line"><a name="l00016"></a><span class="lineno"> 16</span>  <span class="keywordtype">int</span> cig;</div><div class="line"><a name="l00017"></a><span class="lineno"> 17</span>  int32_t n;</div><div class="line"><a name="l00018"></a><span class="lineno"> 18</span>  int32_t ql = 0, rl = ref_start;</div><div class="line"><a name="l00019"></a><span class="lineno"> 19</span>  <span class="keywordflow">for</span> (uint32_t i = 0; i < ncigar; ++i){</div><div class="line"><a name="l00020"></a><span class="lineno"> 20</span>  cig = bam_cigar_op(cigar[i]);</div><div class="line"><a name="l00021"></a><span class="lineno"> 21</span>  n = bam_cigar_oplen(cigar[i]);</div><div class="line"><a name="l00022"></a><span class="lineno"> 22</span>  <span class="keywordflow">if</span> (bam_cigar_type(cig) & 2) { <span class="comment">// Only reference consuming</span></div><div class="line"><a name="l00023"></a><span class="lineno"> 23</span>  <span class="keywordflow">if</span> (pos <= rl + n) {</div><div class="line"><a name="l00024"></a><span class="lineno"> 24</span>  <span class="keywordflow">if</span> (bam_cigar_type(cig) & 1) <span class="comment">// Only query consuming</span></div><div class="line"><a name="l00025"></a><span class="lineno"> 25</span>  ql += (pos - rl); <span class="comment">// n consumed reference, check if it consumes query too.</span></div><div class="line"><a name="l00026"></a><span class="lineno"> 26</span>  <span class="keywordflow">return</span> ql;</div><div class="line"><a name="l00027"></a><span class="lineno"> 27</span>  }</div><div class="line"><a name="l00028"></a><span class="lineno"> 28</span>  rl += n;</div><div class="line"><a name="l00029"></a><span class="lineno"> 29</span>  }</div><div class="line"><a name="l00030"></a><span class="lineno"> 30</span>  <span class="keywordflow">if</span> (bam_cigar_type(cig) & 1) <span class="comment">// Only query consuming</span></div><div class="line"><a name="l00031"></a><span class="lineno"> 31</span>  ql += n;</div><div class="line"><a name="l00032"></a><span class="lineno"> 32</span>  }</div><div class="line"><a name="l00033"></a><span class="lineno"> 33</span>  <span class="keywordflow">return</span> ql;</div><div class="line"><a name="l00034"></a><span class="lineno"> 34</span> }</div><div class="line"><a name="l00035"></a><span class="lineno"> 35</span> </div><div class="line"><a name="l00036"></a><span class="lineno"> 36</span> int32_t get_pos_on_reference(uint32_t *cigar, uint32_t ncigar, uint32_t pos, uint32_t ref_start){</div><div class="line"><a name="l00037"></a><span class="lineno"> 37</span>  <span class="keywordtype">int</span> cig;</div><div class="line"><a name="l00038"></a><span class="lineno"> 38</span>  int32_t n;</div><div class="line"><a name="l00039"></a><span class="lineno"> 39</span>  uint32_t ql = 0, rl = ref_start;</div><div class="line"><a name="l00040"></a><span class="lineno"> 40</span>  <span class="keywordflow">for</span> (uint32_t i = 0; i < ncigar; ++i){</div><div class="line"><a name="l00041"></a><span class="lineno"> 41</span>  cig = bam_cigar_op(cigar[i]);</div><div class="line"><a name="l00042"></a><span class="lineno"> 42</span>  n = bam_cigar_oplen(cigar[i]);</div><div class="line"><a name="l00043"></a><span class="lineno"> 43</span>  <span class="keywordflow">if</span> (bam_cigar_type(cig) & 1) { <span class="comment">// Only query consuming</span></div><div class="line"><a name="l00044"></a><span class="lineno"> 44</span>  <span class="keywordflow">if</span> (pos <= ql + n) {</div><div class="line"><a name="l00045"></a><span class="lineno"> 45</span>  <span class="keywordflow">if</span> (bam_cigar_type(cig) & 2) <span class="comment">// Only reference consuming</span></div><div class="line"><a name="l00046"></a><span class="lineno"> 46</span>  rl += (pos - ql); <span class="comment">// n consumed reference, check if it consumes query too.</span></div><div class="line"><a name="l00047"></a><span class="lineno"> 47</span>  <span class="keywordflow">return</span> rl;</div><div class="line"><a name="l00048"></a><span class="lineno"> 48</span>  }</div><div class="line"><a name="l00049"></a><span class="lineno"> 49</span>  ql += n;</div><div class="line"><a name="l00050"></a><span class="lineno"> 50</span>  }</div><div class="line"><a name="l00051"></a><span class="lineno"> 51</span>  <span class="keywordflow">if</span> (bam_cigar_type(cig) & 2) <span class="comment">// Only reference consuming</span></div><div class="line"><a name="l00052"></a><span class="lineno"> 52</span>  rl += n;</div><div class="line"><a name="l00053"></a><span class="lineno"> 53</span>  }</div><div class="line"><a name="l00054"></a><span class="lineno"> 54</span>  <span class="keywordflow">return</span> rl;</div><div class="line"><a name="l00055"></a><span class="lineno"> 55</span> }</div><div class="line"><a name="l00056"></a><span class="lineno"> 56</span> </div><div class="line"><a name="l00057"></a><span class="lineno"> 57</span> <span class="keywordtype">void</span> reverse_qual(uint8_t *q, <span class="keywordtype">int</span> l){</div><div class="line"><a name="l00058"></a><span class="lineno"> 58</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i = 0; i < l/2; ++i){</div><div class="line"><a name="l00059"></a><span class="lineno"> 59</span>  q[i]^=q[l-i-1];</div><div class="line"><a name="l00060"></a><span class="lineno"> 60</span>  q[l-i-1]^=q[i];</div><div class="line"><a name="l00061"></a><span class="lineno"> 61</span>  q[i]^=q[l-i-1];</div><div class="line"><a name="l00062"></a><span class="lineno"> 62</span>  }</div><div class="line"><a name="l00063"></a><span class="lineno"> 63</span> }</div><div class="line"><a name="l00064"></a><span class="lineno"> 64</span> </div><div class="line"><a name="l00065"></a><span class="lineno"> 65</span> <span class="keywordtype">void</span> reverse_cigar(uint32_t *cigar, <span class="keywordtype">int</span> l){</div><div class="line"><a name="l00066"></a><span class="lineno"> 66</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i = 0; i < l/2; ++i){</div><div class="line"><a name="l00067"></a><span class="lineno"> 67</span>  cigar[i]^=cigar[l-i-1];</div><div class="line"><a name="l00068"></a><span class="lineno"> 68</span>  cigar[l-i-1]^=cigar[i];</div><div class="line"><a name="l00069"></a><span class="lineno"> 69</span>  cigar[i]^=cigar[l-i-1];</div><div class="line"><a name="l00070"></a><span class="lineno"> 70</span>  }</div><div class="line"><a name="l00071"></a><span class="lineno"> 71</span> }</div><div class="line"><a name="l00072"></a><span class="lineno"> 72</span> </div><div class="line"><a name="l00073"></a><span class="lineno"> 73</span> <span class="keywordtype">double</span> mean_quality(uint8_t *a, <span class="keywordtype">int</span> s, <span class="keywordtype">int</span> e){</div><div class="line"><a name="l00074"></a><span class="lineno"> 74</span>  <span class="keywordtype">double</span> m = 0;</div><div class="line"><a name="l00075"></a><span class="lineno"> 75</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i = s; i < e; ++i){</div><div class="line"><a name="l00076"></a><span class="lineno"> 76</span>  m += (double)a[i];</div><div class="line"><a name="l00077"></a><span class="lineno"> 77</span>  }</div><div class="line"><a name="l00078"></a><span class="lineno"> 78</span>  m = m/(e-s);</div><div class="line"><a name="l00079"></a><span class="lineno"> 79</span>  <span class="keywordflow">return</span> m;</div><div class="line"><a name="l00080"></a><span class="lineno"> 80</span> }</div><div class="line"><a name="l00081"></a><span class="lineno"> 81</span> </div><div class="line"><a name="l00082"></a><span class="lineno"> 82</span> <a class="code" href="structcigar__.html">cigar_</a> quality_trim(bam1_t* r, uint8_t qual_threshold, uint8_t sliding_window){</div><div class="line"><a name="l00083"></a><span class="lineno"> 83</span>  uint32_t *ncigar = (uint32_t*) malloc(<span class="keyword">sizeof</span>(uint32_t) * (r->core.n_cigar + 1)), <span class="comment">// Maximum edit is one more element with soft mask</span></div><div class="line"><a name="l00084"></a><span class="lineno"> 84</span>  *cigar = bam_get_cigar(r);</div><div class="line"><a name="l00085"></a><span class="lineno"> 85</span>  uint8_t *qual = bam_get_qual(r);</div><div class="line"><a name="l00086"></a><span class="lineno"> 86</span>  int32_t start_pos;</div><div class="line"><a name="l00087"></a><span class="lineno"> 87</span>  <span class="keywordflow">if</span>(bam_is_rev(r)){</div><div class="line"><a name="l00088"></a><span class="lineno"> 88</span>  reverse_qual(qual, r->core.l_qseq);</div><div class="line"><a name="l00089"></a><span class="lineno"> 89</span>  }</div><div class="line"><a name="l00090"></a><span class="lineno"> 90</span>  <span class="keywordtype">double</span> m = 60;</div><div class="line"><a name="l00091"></a><span class="lineno"> 91</span>  <span class="keywordtype">int</span> del_len, cig, temp;</div><div class="line"><a name="l00092"></a><span class="lineno"> 92</span>  uint32_t i = 0, j = 0;</div><div class="line"><a name="l00093"></a><span class="lineno"> 93</span>  <span class="keywordflow">while</span>(m >= qual_threshold && (i < r->core.l_qseq)){</div><div class="line"><a name="l00094"></a><span class="lineno"> 94</span>  m = mean_quality(qual, i, i+sliding_window);</div><div class="line"><a name="l00095"></a><span class="lineno"> 95</span>  i++;</div><div class="line"><a name="l00096"></a><span class="lineno"> 96</span>  <span class="keywordflow">if</span>(i > r->core.l_qseq - sliding_window)</div><div class="line"><a name="l00097"></a><span class="lineno"> 97</span>  sliding_window--;</div><div class="line"><a name="l00098"></a><span class="lineno"> 98</span>  }</div><div class="line"><a name="l00099"></a><span class="lineno"> 99</span>  del_len = r->core.l_qseq - i;</div><div class="line"><a name="l00100"></a><span class="lineno"> 100</span>  start_pos = get_pos_on_reference(cigar, r->core.n_cigar, del_len, r->core.pos); <span class="comment">// For reverse reads need to set core->pos</span></div><div class="line"><a name="l00101"></a><span class="lineno"> 101</span>  <span class="keywordflow">if</span>(bam_is_rev(r) && start_pos <= r->core.pos) {</div><div class="line"><a name="l00102"></a><span class="lineno"> 102</span>  <span class="keywordflow">return</span> {</div><div class="line"><a name="l00103"></a><span class="lineno"> 103</span>  cigar,</div><div class="line"><a name="l00104"></a><span class="lineno"> 104</span>  r->core.n_cigar,</div><div class="line"><a name="l00105"></a><span class="lineno"> 105</span>  r->core.pos</div><div class="line"><a name="l00106"></a><span class="lineno"> 106</span>  };</div><div class="line"><a name="l00107"></a><span class="lineno"> 107</span>  }</div><div class="line"><a name="l00108"></a><span class="lineno"> 108</span>  int32_t n;</div><div class="line"><a name="l00109"></a><span class="lineno"> 109</span>  i = 0;</div><div class="line"><a name="l00110"></a><span class="lineno"> 110</span>  <span class="keywordflow">if</span>(bam_is_rev(r)){</div><div class="line"><a name="l00111"></a><span class="lineno"> 111</span>  reverse_cigar(cigar, r->core.n_cigar);</div><div class="line"><a name="l00112"></a><span class="lineno"> 112</span>  }</div><div class="line"><a name="l00113"></a><span class="lineno"> 113</span>  reverse_cigar(cigar, r->core.n_cigar); <span class="comment">// Reverse cigar and trim the beginning of read.</span></div><div class="line"><a name="l00114"></a><span class="lineno"> 114</span>  <span class="keywordflow">while</span>(i < r->core.n_cigar){</div><div class="line"><a name="l00115"></a><span class="lineno"> 115</span>  <span class="keywordflow">if</span> (del_len == 0){</div><div class="line"><a name="l00116"></a><span class="lineno"> 116</span>  ncigar[j] = cigar[i];</div><div class="line"><a name="l00117"></a><span class="lineno"> 117</span>  i++;</div><div class="line"><a name="l00118"></a><span class="lineno"> 118</span>  j++;</div><div class="line"><a name="l00119"></a><span class="lineno"> 119</span>  <span class="keywordflow">continue</span>;</div><div class="line"><a name="l00120"></a><span class="lineno"> 120</span>  }</div><div class="line"><a name="l00121"></a><span class="lineno"> 121</span>  cig = bam_cigar_op(cigar[i]);</div><div class="line"><a name="l00122"></a><span class="lineno"> 122</span>  n = bam_cigar_oplen(cigar[i]);</div><div class="line"><a name="l00123"></a><span class="lineno"> 123</span>  <span class="keywordflow">if</span> ((bam_cigar_type(cig) & 1)){ <span class="comment">// Consumes Query</span></div><div class="line"><a name="l00124"></a><span class="lineno"> 124</span>  <span class="keywordflow">if</span>(del_len >= n ){</div><div class="line"><a name="l00125"></a><span class="lineno"> 125</span>  ncigar[j] = bam_cigar_gen(n, BAM_CSOFT_CLIP);</div><div class="line"><a name="l00126"></a><span class="lineno"> 126</span>  } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (del_len < n){</div><div class="line"><a name="l00127"></a><span class="lineno"> 127</span>  ncigar[j] = bam_cigar_gen(del_len, BAM_CSOFT_CLIP);</div><div class="line"><a name="l00128"></a><span class="lineno"> 128</span>  }</div><div class="line"><a name="l00129"></a><span class="lineno"> 129</span>  j++;</div><div class="line"><a name="l00130"></a><span class="lineno"> 130</span>  temp = n;</div><div class="line"><a name="l00131"></a><span class="lineno"> 131</span>  n = std::max(n - del_len, 0);</div><div class="line"><a name="l00132"></a><span class="lineno"> 132</span>  del_len = std::max(del_len - temp, 0);</div><div class="line"><a name="l00133"></a><span class="lineno"> 133</span>  <span class="keywordflow">if</span>(n > 0){</div><div class="line"><a name="l00134"></a><span class="lineno"> 134</span>  ncigar[j] = bam_cigar_gen(n, cig);</div><div class="line"><a name="l00135"></a><span class="lineno"> 135</span>  j++;</div><div class="line"><a name="l00136"></a><span class="lineno"> 136</span>  }</div><div class="line"><a name="l00137"></a><span class="lineno"> 137</span>  }</div><div class="line"><a name="l00138"></a><span class="lineno"> 138</span>  i++;</div><div class="line"><a name="l00139"></a><span class="lineno"> 139</span>  }</div><div class="line"><a name="l00140"></a><span class="lineno"> 140</span>  reverse_cigar(ncigar, j); <span class="comment">// Reverse Back</span></div><div class="line"><a name="l00141"></a><span class="lineno"> 141</span>  <span class="keywordflow">if</span>(bam_is_rev(r)){</div><div class="line"><a name="l00142"></a><span class="lineno"> 142</span>  reverse_cigar(ncigar, j);</div><div class="line"><a name="l00143"></a><span class="lineno"> 143</span>  }</div><div class="line"><a name="l00144"></a><span class="lineno"> 144</span>  <span class="keywordflow">return</span> {</div><div class="line"><a name="l00145"></a><span class="lineno"> 145</span>  ncigar,</div><div class="line"><a name="l00146"></a><span class="lineno"> 146</span>  j,</div><div class="line"><a name="l00147"></a><span class="lineno"> 147</span>  start_pos</div><div class="line"><a name="l00148"></a><span class="lineno"> 148</span>  };</div><div class="line"><a name="l00149"></a><span class="lineno"> 149</span> }</div><div class="line"><a name="l00150"></a><span class="lineno"> 150</span> </div><div class="line"><a name="l00151"></a><span class="lineno"> 151</span> <span class="keywordtype">void</span> print_cigar(uint32_t *cigar, <span class="keywordtype">int</span> nlength){</div><div class="line"><a name="l00152"></a><span class="lineno"> 152</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i = 0; i < nlength; ++i){</div><div class="line"><a name="l00153"></a><span class="lineno"> 153</span>  std::cout << ((cigar[i]) & BAM_CIGAR_MASK);</div><div class="line"><a name="l00154"></a><span class="lineno"> 154</span>  std::cout << <span class="stringliteral">"-"</span> << ((cigar[i]) >> BAM_CIGAR_SHIFT) << <span class="stringliteral">" "</span>;</div><div class="line"><a name="l00155"></a><span class="lineno"> 155</span>  }</div><div class="line"><a name="l00156"></a><span class="lineno"> 156</span> }</div><div class="line"><a name="l00157"></a><span class="lineno"> 157</span> </div><div class="line"><a name="l00158"></a><span class="lineno"> 158</span> <a class="code" href="structcigar__.html">cigar_</a> primer_trim(bam1_t *r, int32_t new_pos){</div><div class="line"><a name="l00159"></a><span class="lineno"> 159</span>  uint32_t *ncigar = (uint32_t*) malloc(<span class="keyword">sizeof</span>(uint32_t) * (r->core.n_cigar + 1)), <span class="comment">// Maximum edit is one more element with soft mask</span></div><div class="line"><a name="l00160"></a><span class="lineno"> 160</span>  *cigar = bam_get_cigar(r);</div><div class="line"><a name="l00161"></a><span class="lineno"> 161</span>  uint32_t i = 0, j = 0;</div><div class="line"><a name="l00162"></a><span class="lineno"> 162</span>  <span class="keywordtype">int</span> del_len, cig, temp;</div><div class="line"><a name="l00163"></a><span class="lineno"> 163</span>  <span class="keywordflow">if</span> (bam_is_rev(r)){</div><div class="line"><a name="l00164"></a><span class="lineno"> 164</span>  del_len = bam_cigar2qlen(r->core.n_cigar, bam_get_cigar(r)) - get_pos_on_query(cigar, r->core.n_cigar, new_pos, r->core.pos) - 1;</div><div class="line"><a name="l00165"></a><span class="lineno"> 165</span>  <span class="comment">// print_cigar(cigar, r->core.n_cigar);</span></div><div class="line"><a name="l00166"></a><span class="lineno"> 166</span>  reverse_cigar(cigar, r->core.n_cigar);</div><div class="line"><a name="l00167"></a><span class="lineno"> 167</span>  <span class="comment">// print_cigar(cigar, r->core.n_cigar);</span></div><div class="line"><a name="l00168"></a><span class="lineno"> 168</span>  <span class="comment">// std::cout << "Del Length " << del_len << std::endl;</span></div><div class="line"><a name="l00169"></a><span class="lineno"> 169</span>  } <span class="keywordflow">else</span> {</div><div class="line"><a name="l00170"></a><span class="lineno"> 170</span>  del_len = get_pos_on_query(cigar, r->core.n_cigar, new_pos, r->core.pos);</div><div class="line"><a name="l00171"></a><span class="lineno"> 171</span>  }</div><div class="line"><a name="l00172"></a><span class="lineno"> 172</span>  int32_t n;</div><div class="line"><a name="l00173"></a><span class="lineno"> 173</span>  <span class="keywordflow">while</span>(i < r->core.n_cigar){</div><div class="line"><a name="l00174"></a><span class="lineno"> 174</span>  <span class="keywordflow">if</span> (del_len == 0){</div><div class="line"><a name="l00175"></a><span class="lineno"> 175</span>  ncigar[j] = cigar[i];</div><div class="line"><a name="l00176"></a><span class="lineno"> 176</span>  i++;</div><div class="line"><a name="l00177"></a><span class="lineno"> 177</span>  j++;</div><div class="line"><a name="l00178"></a><span class="lineno"> 178</span>  <span class="keywordflow">continue</span>;</div><div class="line"><a name="l00179"></a><span class="lineno"> 179</span>  }</div><div class="line"><a name="l00180"></a><span class="lineno"> 180</span>  cig = bam_cigar_op(cigar[i]);</div><div class="line"><a name="l00181"></a><span class="lineno"> 181</span>  n = bam_cigar_oplen(cigar[i]);</div><div class="line"><a name="l00182"></a><span class="lineno"> 182</span>  <span class="keywordflow">if</span> ((bam_cigar_type(cig) & 1)){ <span class="comment">// Consumes Query</span></div><div class="line"><a name="l00183"></a><span class="lineno"> 183</span>  <span class="keywordflow">if</span>(del_len >= n ){</div><div class="line"><a name="l00184"></a><span class="lineno"> 184</span>  ncigar[j] = bam_cigar_gen(n, BAM_CSOFT_CLIP);</div><div class="line"><a name="l00185"></a><span class="lineno"> 185</span>  } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (del_len < n){</div><div class="line"><a name="l00186"></a><span class="lineno"> 186</span>  ncigar[j] = bam_cigar_gen(del_len, BAM_CSOFT_CLIP);</div><div class="line"><a name="l00187"></a><span class="lineno"> 187</span>  }</div><div class="line"><a name="l00188"></a><span class="lineno"> 188</span>  j++;</div><div class="line"><a name="l00189"></a><span class="lineno"> 189</span>  temp = n;</div><div class="line"><a name="l00190"></a><span class="lineno"> 190</span>  n = std::max(n - del_len, 0);</div><div class="line"><a name="l00191"></a><span class="lineno"> 191</span>  del_len = std::max(del_len - temp, 0);</div><div class="line"><a name="l00192"></a><span class="lineno"> 192</span>  <span class="keywordflow">if</span>(n > 0){</div><div class="line"><a name="l00193"></a><span class="lineno"> 193</span>  ncigar[j] = bam_cigar_gen(n, cig);</div><div class="line"><a name="l00194"></a><span class="lineno"> 194</span>  j++;</div><div class="line"><a name="l00195"></a><span class="lineno"> 195</span>  }</div><div class="line"><a name="l00196"></a><span class="lineno"> 196</span>  }</div><div class="line"><a name="l00197"></a><span class="lineno"> 197</span>  i++;</div><div class="line"><a name="l00198"></a><span class="lineno"> 198</span>  }</div><div class="line"><a name="l00199"></a><span class="lineno"> 199</span>  <span class="keywordflow">if</span>(bam_is_rev(r)){</div><div class="line"><a name="l00200"></a><span class="lineno"> 200</span>  reverse_cigar(ncigar, j);</div><div class="line"><a name="l00201"></a><span class="lineno"> 201</span>  }</div><div class="line"><a name="l00202"></a><span class="lineno"> 202</span>  <a class="code" href="structcigar__.html">cigar_</a> t = {</div><div class="line"><a name="l00203"></a><span class="lineno"> 203</span>  ncigar,</div><div class="line"><a name="l00204"></a><span class="lineno"> 204</span>  j</div><div class="line"><a name="l00205"></a><span class="lineno"> 205</span>  };</div><div class="line"><a name="l00206"></a><span class="lineno"> 206</span>  <span class="keywordflow">return</span> t;</div><div class="line"><a name="l00207"></a><span class="lineno"> 207</span> }</div><div class="line"><a name="l00208"></a><span class="lineno"> 208</span> </div><div class="line"><a name="l00209"></a><span class="lineno"> 209</span> <span class="keywordtype">void</span> replace_cigar(bam1_t *b, <span class="keywordtype">int</span> n, uint32_t *cigar){</div><div class="line"><a name="l00210"></a><span class="lineno"> 210</span>  <span class="keywordflow">if</span> (n != b->core.n_cigar) {</div><div class="line"><a name="l00211"></a><span class="lineno"> 211</span>  <span class="keywordtype">int</span> o = b->core.l_qname + b->core.n_cigar * 4;</div><div class="line"><a name="l00212"></a><span class="lineno"> 212</span>  <span class="keywordflow">if</span> (b->l_data + (n - b->core.n_cigar) * 4 > b->m_data) {</div><div class="line"><a name="l00213"></a><span class="lineno"> 213</span>  b->m_data = b->l_data + (n - b->core.n_cigar) * 4;</div><div class="line"><a name="l00214"></a><span class="lineno"> 214</span>  kroundup32(b->m_data);</div><div class="line"><a name="l00215"></a><span class="lineno"> 215</span>  b->data = (uint8_t*)realloc(b->data, b->m_data);</div><div class="line"><a name="l00216"></a><span class="lineno"> 216</span>  }</div><div class="line"><a name="l00217"></a><span class="lineno"> 217</span>  memmove(b->data + b->core.l_qname + n * 4, b->data + o, b->l_data - o);</div><div class="line"><a name="l00218"></a><span class="lineno"> 218</span>  memcpy(b->data + b->core.l_qname, cigar, n * 4);</div><div class="line"><a name="l00219"></a><span class="lineno"> 219</span>  b->l_data += (n - b->core.n_cigar) * 4;</div><div class="line"><a name="l00220"></a><span class="lineno"> 220</span>  b->core.n_cigar = n;</div><div class="line"><a name="l00221"></a><span class="lineno"> 221</span>  } <span class="keywordflow">else</span> memcpy(b->data + b->core.l_qname, cigar, n * 4);</div><div class="line"><a name="l00222"></a><span class="lineno"> 222</span> }</div><div class="line"><a name="l00223"></a><span class="lineno"> 223</span> </div><div class="line"><a name="l00224"></a><span class="lineno"> 224</span> int16_t get_overlapping_primer_indice(bam1_t* r, std::vector<primer> primers){</div><div class="line"><a name="l00225"></a><span class="lineno"> 225</span>  uint32_t query_pos, start_pos, *cigar = bam_get_cigar(r);</div><div class="line"><a name="l00226"></a><span class="lineno"> 226</span>  <span class="keywordflow">if</span>(bam_is_rev(r)){</div><div class="line"><a name="l00227"></a><span class="lineno"> 227</span>  start_pos = bam_endpos(r)-1;</div><div class="line"><a name="l00228"></a><span class="lineno"> 228</span>  query_pos = start_pos + (bam_cigar2qlen(r->core.n_cigar, cigar) - get_pos_on_query(cigar, r->core.n_cigar, start_pos, r->core.pos)) - 1;</div><div class="line"><a name="l00229"></a><span class="lineno"> 229</span>  } <span class="keywordflow">else</span> {</div><div class="line"><a name="l00230"></a><span class="lineno"> 230</span>  start_pos = r->core.pos;</div><div class="line"><a name="l00231"></a><span class="lineno"> 231</span>  query_pos = start_pos - get_pos_on_query(cigar, r->core.n_cigar, start_pos, r->core.pos);</div><div class="line"><a name="l00232"></a><span class="lineno"> 232</span>  }</div><div class="line"><a name="l00233"></a><span class="lineno"> 233</span>  <span class="keywordflow">for</span>(std::vector<int>::size_type i = 0; i!=primers.size();i++){</div><div class="line"><a name="l00234"></a><span class="lineno"> 234</span>  <span class="keywordflow">if</span>(query_pos >= primers[i].get_start() && query_pos <= primers[i].get_end() && start_pos >= primers[i].get_start() && start_pos <= primers[i].get_end()) <span class="comment">// Change int to int32_t in primer_bed.cpp</span></div><div class="line"><a name="l00235"></a><span class="lineno"> 235</span>  <span class="keywordflow">return</span> i;</div><div class="line"><a name="l00236"></a><span class="lineno"> 236</span>  }</div><div class="line"><a name="l00237"></a><span class="lineno"> 237</span>  <span class="keywordflow">return</span> -1;</div><div class="line"><a name="l00238"></a><span class="lineno"> 238</span> }</div><div class="line"><a name="l00239"></a><span class="lineno"> 239</span> </div><div class="line"><a name="l00240"></a><span class="lineno"> 240</span> <a class="code" href="structcigar__.html">cigar_</a> remove_trailing_query_ref_consumption(uint32_t* cigar, uint32_t n){</div><div class="line"><a name="l00241"></a><span class="lineno"> 241</span>  <span class="keywordtype">int</span> i = 0, len = 0, cig, start_pos = 0;</div><div class="line"><a name="l00242"></a><span class="lineno"> 242</span>  <span class="keywordflow">while</span>(i < n){</div><div class="line"><a name="l00243"></a><span class="lineno"> 243</span>  cig = bam_cigar_op(cigar[i]);</div><div class="line"><a name="l00244"></a><span class="lineno"> 244</span>  len = bam_cigar_oplen(cigar[i]);</div><div class="line"><a name="l00245"></a><span class="lineno"> 245</span>  <span class="keywordflow">if</span>((bam_cigar_type(cig) & 2) && (bam_cigar_type(cig) & 1))</div><div class="line"><a name="l00246"></a><span class="lineno"> 246</span>  <span class="keywordflow">break</span>;</div><div class="line"><a name="l00247"></a><span class="lineno"> 247</span>  <span class="keywordflow">if</span>(bam_cigar_type(cig) & 1){ <span class="comment">// Consumes only query - insertion</span></div><div class="line"><a name="l00248"></a><span class="lineno"> 248</span>  cigar[i] = bam_cigar_gen(len, BAM_CSOFT_CLIP);</div><div class="line"><a name="l00249"></a><span class="lineno"> 249</span>  }</div><div class="line"><a name="l00250"></a><span class="lineno"> 250</span>  <span class="keywordflow">if</span>(bam_cigar_type(cig) & 2){ <span class="comment">// Consumes only reference - deletion</span></div><div class="line"><a name="l00251"></a><span class="lineno"> 251</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> j = i; j < n-1; ++j){</div><div class="line"><a name="l00252"></a><span class="lineno"> 252</span>  cigar[j] = cigar[j+1];</div><div class="line"><a name="l00253"></a><span class="lineno"> 253</span>  }</div><div class="line"><a name="l00254"></a><span class="lineno"> 254</span>  n--;</div><div class="line"><a name="l00255"></a><span class="lineno"> 255</span>  i--;</div><div class="line"><a name="l00256"></a><span class="lineno"> 256</span>  start_pos += len;</div><div class="line"><a name="l00257"></a><span class="lineno"> 257</span>  }</div><div class="line"><a name="l00258"></a><span class="lineno"> 258</span>  i++;</div><div class="line"><a name="l00259"></a><span class="lineno"> 259</span>  }</div><div class="line"><a name="l00260"></a><span class="lineno"> 260</span>  reverse_cigar(cigar, n);</div><div class="line"><a name="l00261"></a><span class="lineno"> 261</span>  i = 0, len = 0;</div><div class="line"><a name="l00262"></a><span class="lineno"> 262</span>  <span class="keywordflow">while</span>(i < n){</div><div class="line"><a name="l00263"></a><span class="lineno"> 263</span>  cig = bam_cigar_op(cigar[i]);</div><div class="line"><a name="l00264"></a><span class="lineno"> 264</span>  len = bam_cigar_oplen(cigar[i]);</div><div class="line"><a name="l00265"></a><span class="lineno"> 265</span>  <span class="keywordflow">if</span>((bam_cigar_type(cig) & 2) && (bam_cigar_type(cig) & 1))</div><div class="line"><a name="l00266"></a><span class="lineno"> 266</span>  <span class="keywordflow">break</span>;</div><div class="line"><a name="l00267"></a><span class="lineno"> 267</span>  <span class="keywordflow">if</span>(bam_cigar_type(cig) & 1){ <span class="comment">// Consumes only query - insertion</span></div><div class="line"><a name="l00268"></a><span class="lineno"> 268</span>  cigar[i] = bam_cigar_gen(len, BAM_CSOFT_CLIP);</div><div class="line"><a name="l00269"></a><span class="lineno"> 269</span>  }</div><div class="line"><a name="l00270"></a><span class="lineno"> 270</span>  <span class="keywordflow">if</span>(bam_cigar_type(cig) & 2){ <span class="comment">// Consumes only reference - deletion</span></div><div class="line"><a name="l00271"></a><span class="lineno"> 271</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> j = i; j < n-1; ++j){</div><div class="line"><a name="l00272"></a><span class="lineno"> 272</span>  cigar[j] = cigar[j+1];</div><div class="line"><a name="l00273"></a><span class="lineno"> 273</span>  }</div><div class="line"><a name="l00274"></a><span class="lineno"> 274</span>  n--;</div><div class="line"><a name="l00275"></a><span class="lineno"> 275</span>  i--;</div><div class="line"><a name="l00276"></a><span class="lineno"> 276</span>  <span class="comment">// start_pos += len;</span></div><div class="line"><a name="l00277"></a><span class="lineno"> 277</span>  }</div><div class="line"><a name="l00278"></a><span class="lineno"> 278</span>  i++;</div><div class="line"><a name="l00279"></a><span class="lineno"> 279</span>  }</div><div class="line"><a name="l00280"></a><span class="lineno"> 280</span>  reverse_cigar(cigar, n);</div><div class="line"><a name="l00281"></a><span class="lineno"> 281</span>  <span class="keywordflow">return</span> {</div><div class="line"><a name="l00282"></a><span class="lineno"> 282</span>  cigar,</div><div class="line"><a name="l00283"></a><span class="lineno"> 283</span>  n,</div><div class="line"><a name="l00284"></a><span class="lineno"> 284</span>  start_pos</div><div class="line"><a name="l00285"></a><span class="lineno"> 285</span>  };</div><div class="line"><a name="l00286"></a><span class="lineno"> 286</span> }</div><div class="line"><a name="l00287"></a><span class="lineno"> 287</span> </div><div class="line"><a name="l00288"></a><span class="lineno"> 288</span> <a class="code" href="structcigar__.html">cigar_</a> condense_cigar(uint32_t* cigar, uint32_t n){</div><div class="line"><a name="l00289"></a><span class="lineno"> 289</span>  <span class="keywordtype">int</span> i = 0, len = 0, cig, next_cig, start_pos = 0;</div><div class="line"><a name="l00290"></a><span class="lineno"> 290</span>  <a class="code" href="structcigar__.html">cigar_</a> t = remove_trailing_query_ref_consumption(cigar, n);</div><div class="line"><a name="l00291"></a><span class="lineno"> 291</span>  cigar = t.cigar;</div><div class="line"><a name="l00292"></a><span class="lineno"> 292</span>  n = t.nlength;</div><div class="line"><a name="l00293"></a><span class="lineno"> 293</span>  start_pos = t.start_pos;</div><div class="line"><a name="l00294"></a><span class="lineno"> 294</span>  <span class="keywordflow">while</span>(i< n -1){</div><div class="line"><a name="l00295"></a><span class="lineno"> 295</span>  cig = bam_cigar_op(cigar[i]);</div><div class="line"><a name="l00296"></a><span class="lineno"> 296</span>  next_cig = bam_cigar_op(cigar[i+1]);</div><div class="line"><a name="l00297"></a><span class="lineno"> 297</span>  <span class="keywordflow">if</span>(cig == next_cig){</div><div class="line"><a name="l00298"></a><span class="lineno"> 298</span>  len = bam_cigar_oplen(cigar[i])+bam_cigar_oplen(cigar[i+1]);</div><div class="line"><a name="l00299"></a><span class="lineno"> 299</span>  cigar[i] = bam_cigar_gen(len, bam_cigar_op(cigar[i]));</div><div class="line"><a name="l00300"></a><span class="lineno"> 300</span>  <span class="keywordflow">for</span>(<span class="keywordtype">int</span> j = i+1; j < n - 1; j++){</div><div class="line"><a name="l00301"></a><span class="lineno"> 301</span>  cigar[j] = cigar[j+1];</div><div class="line"><a name="l00302"></a><span class="lineno"> 302</span>  }</div><div class="line"><a name="l00303"></a><span class="lineno"> 303</span>  n--;</div><div class="line"><a name="l00304"></a><span class="lineno"> 304</span>  } <span class="keywordflow">else</span> {</div><div class="line"><a name="l00305"></a><span class="lineno"> 305</span>  i++;</div><div class="line"><a name="l00306"></a><span class="lineno"> 306</span>  }</div><div class="line"><a name="l00307"></a><span class="lineno"> 307</span>  }</div><div class="line"><a name="l00308"></a><span class="lineno"> 308</span>  <span class="keywordflow">return</span> {</div><div class="line"><a name="l00309"></a><span class="lineno"> 309</span>  cigar,</div><div class="line"><a name="l00310"></a><span class="lineno"> 310</span>  n,</div><div class="line"><a name="l00311"></a><span class="lineno"> 311</span>  start_pos</div><div class="line"><a name="l00312"></a><span class="lineno"> 312</span>  };</div><div class="line"><a name="l00313"></a><span class="lineno"> 313</span> }</div><div class="line"><a name="l00314"></a><span class="lineno"> 314</span> </div><div class="line"><a name="l00315"></a><span class="lineno"> 315</span> <span class="keywordtype">int</span> trim_bam_qual_primer(std::string bam, std::string bed, std::string bam_out, std::string region_, uint8_t min_qual, uint8_t sliding_window, <span class="keywordtype">int</span> min_length = 30){</div><div class="line"><a name="l00316"></a><span class="lineno"> 316</span>  std::vector<primer> primers = populate_from_file(bed);</div><div class="line"><a name="l00317"></a><span class="lineno"> 317</span>  <span class="keywordflow">if</span>(bam.empty()){</div><div class="line"><a name="l00318"></a><span class="lineno"> 318</span>  std::cout << <span class="stringliteral">"Bam file in empty."</span> << std::endl;</div><div class="line"><a name="l00319"></a><span class="lineno"> 319</span>  <span class="keywordflow">return</span> -1;</div><div class="line"><a name="l00320"></a><span class="lineno"> 320</span>  }</div><div class="line"><a name="l00321"></a><span class="lineno"> 321</span>  bam_out += <span class="stringliteral">".bam"</span>;</div><div class="line"><a name="l00322"></a><span class="lineno"> 322</span>  samFile *in = hts_open(bam.c_str(), <span class="stringliteral">"r"</span>);</div><div class="line"><a name="l00323"></a><span class="lineno"> 323</span>  BGZF *out = bgzf_open(bam_out.c_str(), <span class="stringliteral">"w"</span>);</div><div class="line"><a name="l00324"></a><span class="lineno"> 324</span>  <span class="keywordflow">if</span>(in == NULL) {</div><div class="line"><a name="l00325"></a><span class="lineno"> 325</span>  std::cout << (<span class="stringliteral">"Unable to open BAM/SAM file."</span>) << std::endl;</div><div class="line"><a name="l00326"></a><span class="lineno"> 326</span>  <span class="keywordflow">return</span> -1;</div><div class="line"><a name="l00327"></a><span class="lineno"> 327</span>  }</div><div class="line"><a name="l00328"></a><span class="lineno"> 328</span>  <span class="comment">//Load the index</span></div><div class="line"><a name="l00329"></a><span class="lineno"> 329</span>  hts_idx_t *idx = sam_index_load(in, bam.c_str());</div><div class="line"><a name="l00330"></a><span class="lineno"> 330</span>  <span class="keywordflow">if</span>(idx == NULL) {</div><div class="line"><a name="l00331"></a><span class="lineno"> 331</span>  std::cout << (<span class="stringliteral">"Unable to open BAM/SAM index."</span>) << std::endl; <span class="comment">// TODO: Generate index</span></div><div class="line"><a name="l00332"></a><span class="lineno"> 332</span>  <span class="keywordflow">return</span> -1;</div><div class="line"><a name="l00333"></a><span class="lineno"> 333</span>  }</div><div class="line"><a name="l00334"></a><span class="lineno"> 334</span>  <span class="comment">//Get the header</span></div><div class="line"><a name="l00335"></a><span class="lineno"> 335</span>  bam_hdr_t *header = sam_hdr_read(in);</div><div class="line"><a name="l00336"></a><span class="lineno"> 336</span>  <span class="keywordflow">if</span>(header == NULL) {</div><div class="line"><a name="l00337"></a><span class="lineno"> 337</span>  sam_close(in);</div><div class="line"><a name="l00338"></a><span class="lineno"> 338</span>  std::cout << <span class="stringliteral">"Unable to open BAM/SAM header."</span> << std::endl;</div><div class="line"><a name="l00339"></a><span class="lineno"> 339</span>  }</div><div class="line"><a name="l00340"></a><span class="lineno"> 340</span>  <span class="keywordflow">if</span>(bam_hdr_write(out, header) < 0){</div><div class="line"><a name="l00341"></a><span class="lineno"> 341</span>  std::cout << <span class="stringliteral">"Unable to write BAM header to path."</span> << std::endl;</div><div class="line"><a name="l00342"></a><span class="lineno"> 342</span>  sam_close(in);</div><div class="line"><a name="l00343"></a><span class="lineno"> 343</span>  <span class="keywordflow">return</span> -1;</div><div class="line"><a name="l00344"></a><span class="lineno"> 344</span>  }</div><div class="line"><a name="l00345"></a><span class="lineno"> 345</span>  <span class="keywordflow">if</span> (region_.empty()){</div><div class="line"><a name="l00346"></a><span class="lineno"> 346</span>  std::cout << <span class="stringliteral">"Number of references: "</span> << header->n_targets << std::endl;</div><div class="line"><a name="l00347"></a><span class="lineno"> 347</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i = 0; i < header->n_targets; ++i){</div><div class="line"><a name="l00348"></a><span class="lineno"> 348</span>  std::cout << <span class="stringliteral">"Reference Name: "</span> << header->target_name[i] << std::endl;</div><div class="line"><a name="l00349"></a><span class="lineno"> 349</span>  std::cout << <span class="stringliteral">"Reference Length: "</span> << header->target_len[i] << std::endl;</div><div class="line"><a name="l00350"></a><span class="lineno"> 350</span>  <span class="keywordflow">if</span>(i==0){</div><div class="line"><a name="l00351"></a><span class="lineno"> 351</span>  region_.assign(header->target_name[i]);</div><div class="line"><a name="l00352"></a><span class="lineno"> 352</span>  }</div><div class="line"><a name="l00353"></a><span class="lineno"> 353</span>  }</div><div class="line"><a name="l00354"></a><span class="lineno"> 354</span>  std::cout << <span class="stringliteral">"Using Region: "</span> << region_ << std::endl;</div><div class="line"><a name="l00355"></a><span class="lineno"> 355</span>  }</div><div class="line"><a name="l00356"></a><span class="lineno"> 356</span>  std::string hdr_text(header->text);</div><div class="line"><a name="l00357"></a><span class="lineno"> 357</span>  <span class="keywordflow">if</span> (hdr_text.find(std::string(<span class="stringliteral">"SO:coordinate"</span>))) {</div><div class="line"><a name="l00358"></a><span class="lineno"> 358</span>  std::cout << <span class="stringliteral">"Sorted By Coordinate"</span> << std::endl; <span class="comment">// Sort by coordinate</span></div><div class="line"><a name="l00359"></a><span class="lineno"> 359</span>  } <span class="keywordflow">if</span>(hdr_text.find(std::string(<span class="stringliteral">"SO:queryname"</span>))) {</div><div class="line"><a name="l00360"></a><span class="lineno"> 360</span>  std::cout << <span class="stringliteral">"Sorted By Query Name"</span> << std::endl; <span class="comment">// Sort by name</span></div><div class="line"><a name="l00361"></a><span class="lineno"> 361</span>  } <span class="keywordflow">else</span> {</div><div class="line"><a name="l00362"></a><span class="lineno"> 362</span>  std::cout << <span class="stringliteral">"Not sorted"</span> << std::endl;</div><div class="line"><a name="l00363"></a><span class="lineno"> 363</span>  }</div><div class="line"><a name="l00364"></a><span class="lineno"> 364</span>  <span class="comment">//Initialize iterator</span></div><div class="line"><a name="l00365"></a><span class="lineno"> 365</span>  hts_itr_t *iter = NULL;</div><div class="line"><a name="l00366"></a><span class="lineno"> 366</span>  <span class="comment">//Move the iterator to the region we are interested in</span></div><div class="line"><a name="l00367"></a><span class="lineno"> 367</span>  iter = sam_itr_querys(idx, header, region_.c_str());</div><div class="line"><a name="l00368"></a><span class="lineno"> 368</span>  <span class="keywordflow">if</span>(header == NULL || iter == NULL) {</div><div class="line"><a name="l00369"></a><span class="lineno"> 369</span>  sam_close(in);</div><div class="line"><a name="l00370"></a><span class="lineno"> 370</span>  std::cout << <span class="stringliteral">"Unable to iterate to region within BAM/SAM."</span> << std::endl;</div><div class="line"><a name="l00371"></a><span class="lineno"> 371</span>  <span class="keywordflow">return</span> -1;</div><div class="line"><a name="l00372"></a><span class="lineno"> 372</span>  }</div><div class="line"><a name="l00373"></a><span class="lineno"> 373</span>  <span class="comment">//Initiate the alignment record</span></div><div class="line"><a name="l00374"></a><span class="lineno"> 374</span>  bam1_t *aln = bam_init1();</div><div class="line"><a name="l00375"></a><span class="lineno"> 375</span>  <span class="keywordtype">int</span> ctr = 0;</div><div class="line"><a name="l00376"></a><span class="lineno"> 376</span>  <a class="code" href="structcigar__.html">cigar_</a> t;</div><div class="line"><a name="l00377"></a><span class="lineno"> 377</span>  int16_t *p = (int16_t*)malloc(<span class="keyword">sizeof</span>(int16_t));</div><div class="line"><a name="l00378"></a><span class="lineno"> 378</span>  uint32_t primer_trim_count = 0;</div><div class="line"><a name="l00379"></a><span class="lineno"> 379</span>  <span class="keywordflow">while</span>(sam_itr_next(in, iter, aln) >= 0) {</div><div class="line"><a name="l00380"></a><span class="lineno"> 380</span>  *p = get_overlapping_primer_indice(aln, primers);</div><div class="line"><a name="l00381"></a><span class="lineno"> 381</span>  <span class="keywordflow">if</span>(*p != -1){</div><div class="line"><a name="l00382"></a><span class="lineno"> 382</span>  primer_trim_count++;</div><div class="line"><a name="l00383"></a><span class="lineno"> 383</span>  <span class="keywordflow">if</span>(bam_is_rev(aln)){</div><div class="line"><a name="l00384"></a><span class="lineno"> 384</span>  t = primer_trim(aln, primers[*p].get_start() - 1);</div><div class="line"><a name="l00385"></a><span class="lineno"> 385</span>  } <span class="keywordflow">else</span> {</div><div class="line"><a name="l00386"></a><span class="lineno"> 386</span>  t = primer_trim(aln, primers[*p].get_end() + 1);</div><div class="line"><a name="l00387"></a><span class="lineno"> 387</span>  aln->core.pos = primers[*p].get_end() + 1;</div><div class="line"><a name="l00388"></a><span class="lineno"> 388</span>  }</div><div class="line"><a name="l00389"></a><span class="lineno"> 389</span>  replace_cigar(aln, t.nlength, t.cigar);</div><div class="line"><a name="l00390"></a><span class="lineno"> 390</span>  }</div><div class="line"><a name="l00391"></a><span class="lineno"> 391</span>  t = quality_trim(aln, min_qual, sliding_window); <span class="comment">// Quality Trimming</span></div><div class="line"><a name="l00392"></a><span class="lineno"> 392</span>  <span class="keywordflow">if</span>(bam_is_rev(aln))</div><div class="line"><a name="l00393"></a><span class="lineno"> 393</span>  aln->core.pos = t.start_pos;</div><div class="line"><a name="l00394"></a><span class="lineno"> 394</span>  t = condense_cigar(t.cigar, t.nlength);</div><div class="line"><a name="l00395"></a><span class="lineno"> 395</span>  aln->core.pos += t.start_pos;</div><div class="line"><a name="l00396"></a><span class="lineno"> 396</span>  replace_cigar(aln, t.nlength, t.cigar);</div><div class="line"><a name="l00397"></a><span class="lineno"> 397</span>  <span class="keywordflow">if</span>(bam_cigar2rlen(aln->core.n_cigar, bam_get_cigar(aln)) >= min_length){</div><div class="line"><a name="l00398"></a><span class="lineno"> 398</span>  <span class="keywordflow">if</span>(*p != -1)</div><div class="line"><a name="l00399"></a><span class="lineno"> 399</span>  bam_aux_append(aln, <span class="stringliteral">"xa"</span>, <span class="charliteral">'i'</span>, 4, (uint8_t*) p);</div><div class="line"><a name="l00400"></a><span class="lineno"> 400</span>  <span class="keywordflow">if</span>(bam_write1(out, aln) < 0){</div><div class="line"><a name="l00401"></a><span class="lineno"> 401</span>  std::cout << <span class="stringliteral">"Not able to write to BAM"</span> << std::endl;</div><div class="line"><a name="l00402"></a><span class="lineno"> 402</span>  hts_itr_destroy(iter);</div><div class="line"><a name="l00403"></a><span class="lineno"> 403</span>  hts_idx_destroy(idx);</div><div class="line"><a name="l00404"></a><span class="lineno"> 404</span>  bam_destroy1(aln);</div><div class="line"><a name="l00405"></a><span class="lineno"> 405</span>  bam_hdr_destroy(header);</div><div class="line"><a name="l00406"></a><span class="lineno"> 406</span>  sam_close(in);</div><div class="line"><a name="l00407"></a><span class="lineno"> 407</span>  bgzf_close(out);</div><div class="line"><a name="l00408"></a><span class="lineno"> 408</span>  <span class="keywordflow">return</span> -1;</div><div class="line"><a name="l00409"></a><span class="lineno"> 409</span>  };</div><div class="line"><a name="l00410"></a><span class="lineno"> 410</span>  }</div><div class="line"><a name="l00411"></a><span class="lineno"> 411</span>  ctr++;</div><div class="line"><a name="l00412"></a><span class="lineno"> 412</span>  <span class="keywordflow">if</span>(ctr % 1000000 == 0){</div><div class="line"><a name="l00413"></a><span class="lineno"> 413</span>  std::cout << <span class="stringliteral">"Processed "</span> << ctr << <span class="stringliteral">"reads ... "</span> << std::endl;</div><div class="line"><a name="l00414"></a><span class="lineno"> 414</span>  }</div><div class="line"><a name="l00415"></a><span class="lineno"> 415</span>  }</div><div class="line"><a name="l00416"></a><span class="lineno"> 416</span>  std::cout << <span class="stringliteral">"Results: "</span> << std::endl;</div><div class="line"><a name="l00417"></a><span class="lineno"> 417</span>  std::cout << <span class="stringliteral">"Trimmed primers from "</span> << primer_trim_count << <span class="stringliteral">" reads."</span> << std::endl;</div><div class="line"><a name="l00418"></a><span class="lineno"> 418</span>  hts_itr_destroy(iter);</div><div class="line"><a name="l00419"></a><span class="lineno"> 419</span>  hts_idx_destroy(idx);</div><div class="line"><a name="l00420"></a><span class="lineno"> 420</span>  bam_destroy1(aln);</div><div class="line"><a name="l00421"></a><span class="lineno"> 421</span>  bam_hdr_destroy(header);</div><div class="line"><a name="l00422"></a><span class="lineno"> 422</span>  sam_close(in);</div><div class="line"><a name="l00423"></a><span class="lineno"> 423</span>  bgzf_close(out);</div><div class="line"><a name="l00424"></a><span class="lineno"> 424</span>  <span class="keywordflow">return</span> 0;</div><div class="line"><a name="l00425"></a><span class="lineno"> 425</span> }</div><div class="line"><a name="l00426"></a><span class="lineno"> 426</span> </div><div class="line"><a name="l00427"></a><span class="lineno"> 427</span> <span class="comment">// int main_old(int argc, char* argv[]) {</span></div><div class="line"><a name="l00428"></a><span class="lineno"> 428</span> <span class="comment">// std::cout << "Path " << argv[1] <<std::endl;</span></div><div class="line"><a name="l00429"></a><span class="lineno"> 429</span> <span class="comment">// std::string bam = std::string(argv[1]);</span></div><div class="line"><a name="l00430"></a><span class="lineno"> 430</span> <span class="comment">// std::string region_;</span></div><div class="line"><a name="l00431"></a><span class="lineno"> 431</span> <span class="comment">// std::string bed = std::string(argv[2]);</span></div><div class="line"><a name="l00432"></a><span class="lineno"> 432</span> <span class="comment">// std::string bam_out = std::string(argv[3]);</span></div><div class="line"><a name="l00433"></a><span class="lineno"> 433</span> <span class="comment">// std::vector<primer> primers = populate_from_file(bed);</span></div><div class="line"><a name="l00434"></a><span class="lineno"> 434</span> <span class="comment">// if(argc > 4) {</span></div><div class="line"><a name="l00435"></a><span class="lineno"> 435</span> <span class="comment">// region_ = std::string(argv[4]);</span></div><div class="line"><a name="l00436"></a><span class="lineno"> 436</span> <span class="comment">// }</span></div><div class="line"><a name="l00437"></a><span class="lineno"> 437</span> <span class="comment">// if(!bam.empty()) {</span></div><div class="line"><a name="l00438"></a><span class="lineno"> 438</span> <span class="comment">// //open BAM for reading</span></div><div class="line"><a name="l00439"></a><span class="lineno"> 439</span> <span class="comment">// samFile *in = hts_open(bam.c_str(), "r");</span></div><div class="line"><a name="l00440"></a><span class="lineno"> 440</span> <span class="comment">// BGZF *out = bgzf_open(bam_out.c_str(), "w");</span></div><div class="line"><a name="l00441"></a><span class="lineno"> 441</span> <span class="comment">// if(in == NULL) {</span></div><div class="line"><a name="l00442"></a><span class="lineno"> 442</span> <span class="comment">// throw std::runtime_error("Unable to open BAM/SAM file.");</span></div><div class="line"><a name="l00443"></a><span class="lineno"> 443</span> <span class="comment">// }</span></div><div class="line"><a name="l00444"></a><span class="lineno"> 444</span> <span class="comment">// //Load the index</span></div><div class="line"><a name="l00445"></a><span class="lineno"> 445</span> <span class="comment">// hts_idx_t *idx = sam_index_load(in, bam.c_str());</span></div><div class="line"><a name="l00446"></a><span class="lineno"> 446</span> <span class="comment">// if(idx == NULL) {</span></div><div class="line"><a name="l00447"></a><span class="lineno"> 447</span> <span class="comment">// throw std::runtime_error("Unable to open BAM/SAM index."); // Generate index</span></div><div class="line"><a name="l00448"></a><span class="lineno"> 448</span> <span class="comment">// }</span></div><div class="line"><a name="l00449"></a><span class="lineno"> 449</span> <span class="comment">// //Get the header</span></div><div class="line"><a name="l00450"></a><span class="lineno"> 450</span> <span class="comment">// bam_hdr_t *header = sam_hdr_read(in);</span></div><div class="line"><a name="l00451"></a><span class="lineno"> 451</span> <span class="comment">// bam_hdr_write(out, header);</span></div><div class="line"><a name="l00452"></a><span class="lineno"> 452</span> <span class="comment">// if(header == NULL) {</span></div><div class="line"><a name="l00453"></a><span class="lineno"> 453</span> <span class="comment">// sam_close(in);</span></div><div class="line"><a name="l00454"></a><span class="lineno"> 454</span> <span class="comment">// throw std::runtime_error("Unable to open BAM header.");</span></div><div class="line"><a name="l00455"></a><span class="lineno"> 455</span> <span class="comment">// }</span></div><div class="line"><a name="l00456"></a><span class="lineno"> 456</span> <span class="comment">// if (region_.empty()){</span></div><div class="line"><a name="l00457"></a><span class="lineno"> 457</span> <span class="comment">// std::cout << "Number of references: " << header->n_targets << std::endl;</span></div><div class="line"><a name="l00458"></a><span class="lineno"> 458</span> <span class="comment">// for (int i = 0; i < header->n_targets; ++i){</span></div><div class="line"><a name="l00459"></a><span class="lineno"> 459</span> <span class="comment">// std::cout << "Reference Name: " << header->target_name[i] << std::endl;</span></div><div class="line"><a name="l00460"></a><span class="lineno"> 460</span> <span class="comment">// std::cout << "Reference Length: " << header->target_len[i] << std::endl;</span></div><div class="line"><a name="l00461"></a><span class="lineno"> 461</span> <span class="comment">// if(i==0){</span></div><div class="line"><a name="l00462"></a><span class="lineno"> 462</span> <span class="comment">// region_.assign(header->target_name[i]);</span></div><div class="line"><a name="l00463"></a><span class="lineno"> 463</span> <span class="comment">// }</span></div><div class="line"><a name="l00464"></a><span class="lineno"> 464</span> <span class="comment">// }</span></div><div class="line"><a name="l00465"></a><span class="lineno"> 465</span> <span class="comment">// std::cout << "Using Region: " << region_ << std::endl;</span></div><div class="line"><a name="l00466"></a><span class="lineno"> 466</span> <span class="comment">// }</span></div><div class="line"><a name="l00467"></a><span class="lineno"> 467</span> <span class="comment">// std::string hdr_text(header->text);</span></div><div class="line"><a name="l00468"></a><span class="lineno"> 468</span> <span class="comment">// if (hdr_text.find(std::string("SO:coordinate"))) {</span></div><div class="line"><a name="l00469"></a><span class="lineno"> 469</span> <span class="comment">// std::cout << "Sorted By Coordinate" << std::endl; // Sort by coordinate</span></div><div class="line"><a name="l00470"></a><span class="lineno"> 470</span> <span class="comment">// } if(hdr_text.find(std::string("SO:queryname"))) {</span></div><div class="line"><a name="l00471"></a><span class="lineno"> 471</span> <span class="comment">// std::cout << "Sorted By Query Name" << std::endl; // Sort by name</span></div><div class="line"><a name="l00472"></a><span class="lineno"> 472</span> <span class="comment">// } else {</span></div><div class="line"><a name="l00473"></a><span class="lineno"> 473</span> <span class="comment">// std::cout << "Not sorted" << std::endl;</span></div><div class="line"><a name="l00474"></a><span class="lineno"> 474</span> <span class="comment">// }</span></div><div class="line"><a name="l00475"></a><span class="lineno"> 475</span> <span class="comment">// //Initialize iterator</span></div><div class="line"><a name="l00476"></a><span class="lineno"> 476</span> <span class="comment">// hts_itr_t *iter = NULL;</span></div><div class="line"><a name="l00477"></a><span class="lineno"> 477</span> <span class="comment">// //Move the iterator to the region we are interested in</span></div><div class="line"><a name="l00478"></a><span class="lineno"> 478</span> <span class="comment">// iter = sam_itr_querys(idx, header, region_.c_str());</span></div><div class="line"><a name="l00479"></a><span class="lineno"> 479</span> <span class="comment">// if(header == NULL || iter == NULL) {</span></div><div class="line"><a name="l00480"></a><span class="lineno"> 480</span> <span class="comment">// sam_close(in);</span></div><div class="line"><a name="l00481"></a><span class="lineno"> 481</span> <span class="comment">// throw std::runtime_error("Unable to iterate to region within BAM.");</span></div><div class="line"><a name="l00482"></a><span class="lineno"> 482</span> <span class="comment">// }</span></div><div class="line"><a name="l00483"></a><span class="lineno"> 483</span> <span class="comment">// //Initiate the alignment record</span></div><div class="line"><a name="l00484"></a><span class="lineno"> 484</span> <span class="comment">// bam1_t *aln = bam_init1();</span></div><div class="line"><a name="l00485"></a><span class="lineno"> 485</span> <span class="comment">// int ctr = 0;</span></div><div class="line"><a name="l00486"></a><span class="lineno"> 486</span> <span class="comment">// cigar_ t;</span></div><div class="line"><a name="l00487"></a><span class="lineno"> 487</span> <span class="comment">// int16_t *p = (int16_t*)malloc(sizeof(int16_t));</span></div><div class="line"><a name="l00488"></a><span class="lineno"> 488</span> <span class="comment">// while(sam_itr_next(in, iter, aln) >= 0) {</span></div><div class="line"><a name="l00489"></a><span class="lineno"> 489</span> <span class="comment">// // std::cout << "Query Length: " << bam_cigar2qlen(aln->core.n_cigar, cigar) << std::endl;</span></div><div class="line"><a name="l00490"></a><span class="lineno"> 490</span> <span class="comment">// // std::cout << "Read Length: " << bam_cigar2rlen(aln->core.n_cigar, cigar) << std::endl;</span></div><div class="line"><a name="l00491"></a><span class="lineno"> 491</span> <span class="comment">// // std::cout << "Query Start: " << get_query_start(aln) << std::endl;</span></div><div class="line"><a name="l00492"></a><span class="lineno"> 492</span> <span class="comment">// // std::cout << "Query End: " << get_query_end(aln) << std::endl;</span></div><div class="line"><a name="l00493"></a><span class="lineno"> 493</span> <span class="comment">// *p = get_overlapping_primer_indice(aln, primers);</span></div><div class="line"><a name="l00494"></a><span class="lineno"> 494</span> <span class="comment">// if(*p == -1)</span></div><div class="line"><a name="l00495"></a><span class="lineno"> 495</span> <span class="comment">// continue;</span></div><div class="line"><a name="l00496"></a><span class="lineno"> 496</span> <span class="comment">// if(bam_is_rev(aln)){</span></div><div class="line"><a name="l00497"></a><span class="lineno"> 497</span> <span class="comment">// t = primer_trim(aln, primers[*p].get_start() - 1);</span></div><div class="line"><a name="l00498"></a><span class="lineno"> 498</span> <span class="comment">// } else {</span></div><div class="line"><a name="l00499"></a><span class="lineno"> 499</span> <span class="comment">// t = primer_trim(aln, primers[*p].get_end() + 1);</span></div><div class="line"><a name="l00500"></a><span class="lineno"> 500</span> <span class="comment">// aln->core.pos = primers[*p].get_end() + 1;</span></div><div class="line"><a name="l00501"></a><span class="lineno"> 501</span> <span class="comment">// }</span></div><div class="line"><a name="l00502"></a><span class="lineno"> 502</span> <span class="comment">// // std::cout << *p << std::endl;</span></div><div class="line"><a name="l00503"></a><span class="lineno"> 503</span> <span class="comment">// replace_cigar(aln, t.nlength, t.cigar);</span></div><div class="line"><a name="l00504"></a><span class="lineno"> 504</span> <span class="comment">// // Quality Trimming</span></div><div class="line"><a name="l00505"></a><span class="lineno"> 505</span> <span class="comment">// t = quality_trim(aln);</span></div><div class="line"><a name="l00506"></a><span class="lineno"> 506</span> <span class="comment">// if(bam_is_rev(aln))</span></div><div class="line"><a name="l00507"></a><span class="lineno"> 507</span> <span class="comment">// aln->core.pos = t.start_pos;</span></div><div class="line"><a name="l00508"></a><span class="lineno"> 508</span> <span class="comment">// // std::cout << "Old: " << t.nlength << std::endl;</span></div><div class="line"><a name="l00509"></a><span class="lineno"> 509</span> <span class="comment">// // print_cigar(t.cigar, t.nlength);</span></div><div class="line"><a name="l00510"></a><span class="lineno"> 510</span> <span class="comment">// // std::cout << std::endl;</span></div><div class="line"><a name="l00511"></a><span class="lineno"> 511</span> <span class="comment">// t = condense_cigar(t.cigar, t.nlength);</span></div><div class="line"><a name="l00512"></a><span class="lineno"> 512</span> <span class="comment">// aln->core.pos += t.start_pos;</span></div><div class="line"><a name="l00513"></a><span class="lineno"> 513</span> <span class="comment">// // std::cout << "New: " << t.nlength << std::endl;</span></div><div class="line"><a name="l00514"></a><span class="lineno"> 514</span> <span class="comment">// // print_cigar(t.cigar, t.nlength);</span></div><div class="line"><a name="l00515"></a><span class="lineno"> 515</span> <span class="comment">// // std::cout << std::endl;</span></div><div class="line"><a name="l00516"></a><span class="lineno"> 516</span> <span class="comment">// replace_cigar(aln, t.nlength, t.cigar);</span></div><div class="line"><a name="l00517"></a><span class="lineno"> 517</span> <span class="comment">// // std::cout << "Name: " << name << std::endl;</span></div><div class="line"><a name="l00518"></a><span class="lineno"> 518</span> <span class="comment">// // std::cout << "Seq: " << seq << "\tQual: " << qual;</span></div><div class="line"><a name="l00519"></a><span class="lineno"> 519</span> <span class="comment">// // std::cout << std::endl;</span></div><div class="line"><a name="l00520"></a><span class="lineno"> 520</span> <span class="comment">// // if (strcmp(bam_get_qname(aln), "M01244:143:000000000-BHWC7:1:1104:13925:8758") == 0){</span></div><div class="line"><a name="l00521"></a><span class="lineno"> 521</span> <span class="comment">// // std::cout << "Start Pos: " << aln->core.pos << std::endl;</span></div><div class="line"><a name="l00522"></a><span class="lineno"> 522</span> <span class="comment">// // uint8_t *q = bam_get_qual(aln);</span></div><div class="line"><a name="l00523"></a><span class="lineno"> 523</span> <span class="comment">// // for (int i = 0; i < aln->core.l_qseq -4; ++i){</span></div><div class="line"><a name="l00524"></a><span class="lineno"> 524</span> <span class="comment">// // std::cout << mean_quality(q, i, i+4) << " ";</span></div><div class="line"><a name="l00525"></a><span class="lineno"> 525</span> <span class="comment">// // }</span></div><div class="line"><a name="l00526"></a><span class="lineno"> 526</span> <span class="comment">// // std::cout << std::endl;</span></div><div class="line"><a name="l00527"></a><span class="lineno"> 527</span> <span class="comment">// // }</span></div><div class="line"><a name="l00528"></a><span class="lineno"> 528</span> <span class="comment">// int min_length = 30;</span></div><div class="line"><a name="l00529"></a><span class="lineno"> 529</span> <span class="comment">// if(bam_cigar2rlen(aln->core.n_cigar, bam_get_cigar(aln)) >= min_length){</span></div><div class="line"><a name="l00530"></a><span class="lineno"> 530</span> <span class="comment">// // bam1_t *b, const char tag[2], char type, int len, const uint8_t *data</span></div><div class="line"><a name="l00531"></a><span class="lineno"> 531</span> <span class="comment">// bam_aux_append(aln, "xa", 'i', 4, (uint8_t*) p);</span></div><div class="line"><a name="l00532"></a><span class="lineno"> 532</span> <span class="comment">// bam_write1(out, aln);</span></div><div class="line"><a name="l00533"></a><span class="lineno"> 533</span> <span class="comment">// }</span></div><div class="line"><a name="l00534"></a><span class="lineno"> 534</span> <span class="comment">// ctr++;</span></div><div class="line"><a name="l00535"></a><span class="lineno"> 535</span> <span class="comment">// if(ctr % 100000 == 0){</span></div><div class="line"><a name="l00536"></a><span class="lineno"> 536</span> <span class="comment">// std::cout << ctr << std::endl;</span></div><div class="line"><a name="l00537"></a><span class="lineno"> 537</span> <span class="comment">// }</span></div><div class="line"><a name="l00538"></a><span class="lineno"> 538</span> <span class="comment">// }</span></div><div class="line"><a name="l00539"></a><span class="lineno"> 539</span> <span class="comment">// hts_itr_destroy(iter);</span></div><div class="line"><a name="l00540"></a><span class="lineno"> 540</span> <span class="comment">// hts_idx_destroy(idx);</span></div><div class="line"><a name="l00541"></a><span class="lineno"> 541</span> <span class="comment">// bam_destroy1(aln);</span></div><div class="line"><a name="l00542"></a><span class="lineno"> 542</span> <span class="comment">// bam_hdr_destroy(header);</span></div><div class="line"><a name="l00543"></a><span class="lineno"> 543</span> <span class="comment">// sam_close(in);</span></div><div class="line"><a name="l00544"></a><span class="lineno"> 544</span> <span class="comment">// bgzf_close(out);</span></div><div class="line"><a name="l00545"></a><span class="lineno"> 545</span> <span class="comment">// }</span></div><div class="line"><a name="l00546"></a><span class="lineno"> 546</span> <span class="comment">// return 0;</span></div><div class="line"><a name="l00547"></a><span class="lineno"> 547</span> <span class="comment">// }</span></div><div class="ttc" id="structcigar___html"><div class="ttname"><a href="structcigar__.html">cigar_</a></div><div class="ttdef"><b>Definition:</b> <a href="trim__primer__quality_8h_source.html#l00010">trim_primer_quality.h:10</a></div></div>
</div><!-- fragment --></div><!-- contents -->
<!-- start footer part -->
<hr class="footer"/><address class="footer"><small>
Generated by  <a href="http://www.doxygen.org/index.html">
<img class="footer" src="doxygen.png" alt="doxygen"/>
</a> 1.8.14
</small></address>
</body>
</html>
|