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 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754
|
<HTML>
<HEAD>
<TITLE>OfflinePrecursorIonSelection.h Source File</TITLE>
<LINK HREF="doxygen.css" REL="stylesheet" TYPE="text/css">
<LINK HREF="style_ini.css" REL="stylesheet" TYPE="text/css">
</HEAD>
<BODY BGCOLOR="#FFFFFF">
<A href="index.html">Home</A> ·
<A href="classes.html">Classes</A> ·
<A href="annotated.html">Annotated Classes</A> ·
<A href="modules.html">Modules</A> ·
<A href="functions_func.html">Members</A> ·
<A href="namespaces.html">Namespaces</A> ·
<A href="pages.html">Related Pages</A>
<HR style="height:1px; border:none; border-top:1px solid #c0c0c0;">
<!-- Generated by Doxygen 1.8.5 -->
<div id="nav-path" class="navpath">
<ul>
<li class="navelem"><a class="el" href="dir_e770f0cf77e550adde3f44739ef529fe.html">include</a></li><li class="navelem"><a class="el" href="dir_6a63c4937d4da007e55fff5dcf71e0f8.html">OpenMS</a></li><li class="navelem"><a class="el" href="dir_330fe73d44ce641f33dd6036f80161a5.html">ANALYSIS</a></li><li class="navelem"><a class="el" href="dir_432be27a4896f4299230a7289acf24ab.html">TARGETED</a></li> </ul>
</div>
</div><!-- top -->
<div class="header">
<div class="headertitle">
<div class="title">OfflinePrecursorIonSelection.h</div> </div>
</div><!--header-->
<div class="contents">
<a href="OfflinePrecursorIonSelection_8h.html">Go to the documentation of this file.</a><div class="fragment"><div class="line"><a name="l00001"></a><span class="lineno"> 1</span> <span class="comment">// --------------------------------------------------------------------------</span></div>
<div class="line"><a name="l00002"></a><span class="lineno"> 2</span> <span class="comment">// OpenMS -- Open-Source Mass Spectrometry</span></div>
<div class="line"><a name="l00003"></a><span class="lineno"> 3</span> <span class="comment">// --------------------------------------------------------------------------</span></div>
<div class="line"><a name="l00004"></a><span class="lineno"> 4</span> <span class="comment">// Copyright The OpenMS Team -- Eberhard Karls University Tuebingen,</span></div>
<div class="line"><a name="l00005"></a><span class="lineno"> 5</span> <span class="comment">// ETH Zurich, and Freie Universitaet Berlin 2002-2013.</span></div>
<div class="line"><a name="l00006"></a><span class="lineno"> 6</span> <span class="comment">//</span></div>
<div class="line"><a name="l00007"></a><span class="lineno"> 7</span> <span class="comment">// This software is released under a three-clause BSD license:</span></div>
<div class="line"><a name="l00008"></a><span class="lineno"> 8</span> <span class="comment">// * Redistributions of source code must retain the above copyright</span></div>
<div class="line"><a name="l00009"></a><span class="lineno"> 9</span> <span class="comment">// notice, this list of conditions and the following disclaimer.</span></div>
<div class="line"><a name="l00010"></a><span class="lineno"> 10</span> <span class="comment">// * Redistributions in binary form must reproduce the above copyright</span></div>
<div class="line"><a name="l00011"></a><span class="lineno"> 11</span> <span class="comment">// notice, this list of conditions and the following disclaimer in the</span></div>
<div class="line"><a name="l00012"></a><span class="lineno"> 12</span> <span class="comment">// documentation and/or other materials provided with the distribution.</span></div>
<div class="line"><a name="l00013"></a><span class="lineno"> 13</span> <span class="comment">// * Neither the name of any author or any participating institution</span></div>
<div class="line"><a name="l00014"></a><span class="lineno"> 14</span> <span class="comment">// may be used to endorse or promote products derived from this software</span></div>
<div class="line"><a name="l00015"></a><span class="lineno"> 15</span> <span class="comment">// without specific prior written permission.</span></div>
<div class="line"><a name="l00016"></a><span class="lineno"> 16</span> <span class="comment">// For a full list of authors, refer to the file AUTHORS.</span></div>
<div class="line"><a name="l00017"></a><span class="lineno"> 17</span> <span class="comment">// --------------------------------------------------------------------------</span></div>
<div class="line"><a name="l00018"></a><span class="lineno"> 18</span> <span class="comment">// THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"</span></div>
<div class="line"><a name="l00019"></a><span class="lineno"> 19</span> <span class="comment">// AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE</span></div>
<div class="line"><a name="l00020"></a><span class="lineno"> 20</span> <span class="comment">// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE</span></div>
<div class="line"><a name="l00021"></a><span class="lineno"> 21</span> <span class="comment">// ARE DISCLAIMED. IN NO EVENT SHALL ANY OF THE AUTHORS OR THE CONTRIBUTING</span></div>
<div class="line"><a name="l00022"></a><span class="lineno"> 22</span> <span class="comment">// INSTITUTIONS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,</span></div>
<div class="line"><a name="l00023"></a><span class="lineno"> 23</span> <span class="comment">// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,</span></div>
<div class="line"><a name="l00024"></a><span class="lineno"> 24</span> <span class="comment">// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS;</span></div>
<div class="line"><a name="l00025"></a><span class="lineno"> 25</span> <span class="comment">// OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY,</span></div>
<div class="line"><a name="l00026"></a><span class="lineno"> 26</span> <span class="comment">// WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR</span></div>
<div class="line"><a name="l00027"></a><span class="lineno"> 27</span> <span class="comment">// OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF</span></div>
<div class="line"><a name="l00028"></a><span class="lineno"> 28</span> <span class="comment">// ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.</span></div>
<div class="line"><a name="l00029"></a><span class="lineno"> 29</span> <span class="comment">//</span></div>
<div class="line"><a name="l00030"></a><span class="lineno"> 30</span> <span class="comment">// --------------------------------------------------------------------------</span></div>
<div class="line"><a name="l00031"></a><span class="lineno"> 31</span> <span class="comment">// $Maintainer: Alexandra Zerck $</span></div>
<div class="line"><a name="l00032"></a><span class="lineno"> 32</span> <span class="comment">// $Authors: Alexandra Zerck $</span></div>
<div class="line"><a name="l00033"></a><span class="lineno"> 33</span> <span class="comment">// --------------------------------------------------------------------------</span></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> <span class="preprocessor">#ifndef OPENMS_ANALYSIS_TARGETED_OFFLINEPRECURSORIONSELECTION_H</span></div>
<div class="line"><a name="l00036"></a><span class="lineno"> 36</span> <span class="preprocessor"></span><span class="preprocessor">#define OPENMS_ANALYSIS_TARGETED_OFFLINEPRECURSORIONSELECTION_H</span></div>
<div class="line"><a name="l00037"></a><span class="lineno"> 37</span> <span class="preprocessor"></span></div>
<div class="line"><a name="l00038"></a><span class="lineno"> 38</span> </div>
<div class="line"><a name="l00039"></a><span class="lineno"> 39</span> <span class="preprocessor">#include <<a class="code" href="FeatureMap_8h.html">OpenMS/KERNEL/FeatureMap.h</a>></span></div>
<div class="line"><a name="l00040"></a><span class="lineno"> 40</span> <span class="preprocessor">#include <<a class="code" href="MSExperiment_8h.html">OpenMS/KERNEL/MSExperiment.h</a>></span></div>
<div class="line"><a name="l00041"></a><span class="lineno"> 41</span> <span class="preprocessor">#include <<a class="code" href="IDMapper_8h.html">OpenMS/ANALYSIS/ID/IDMapper.h</a>></span></div>
<div class="line"><a name="l00042"></a><span class="lineno"> 42</span> <span class="preprocessor">#include <<a class="code" href="FeatureXMLFile_8h.html">OpenMS/FORMAT/FeatureXMLFile.h</a>></span></div>
<div class="line"><a name="l00043"></a><span class="lineno"> 43</span> <span class="preprocessor">#include <<a class="code" href="PSLPFormulation_8h.html">OpenMS/ANALYSIS/TARGETED/PSLPFormulation.h</a>></span></div>
<div class="line"><a name="l00044"></a><span class="lineno"> 44</span> </div>
<div class="line"><a name="l00045"></a><span class="lineno"> 45</span> <span class="keyword">namespace </span>OpenMS</div>
<div class="line"><a name="l00046"></a><span class="lineno"> 46</span> {</div>
<div class="line"><a name="l00047"></a><span class="lineno"> 47</span>  <span class="keyword">class </span>PeptideIdentification;</div>
<div class="line"><a name="l00048"></a><span class="lineno"> 48</span>  <span class="keyword">class </span>ProteinIdentification;</div>
<div class="line"><a name="l00049"></a><span class="lineno"> 49</span>  <span class="keyword">class </span>String;</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> </div>
<div class="line"><a name="l00061"></a><span class="lineno"><a class="line" href="classOpenMS_1_1OfflinePrecursorIonSelection.html"> 61</a></span>  <span class="keyword">class </span>OPENMS_DLLAPI <a class="code" href="classOpenMS_1_1OfflinePrecursorIonSelection.html">OfflinePrecursorIonSelection</a> :</div>
<div class="line"><a name="l00062"></a><span class="lineno"> 62</span>  <span class="keyword">public</span> <a class="code" href="classOpenMS_1_1DefaultParamHandler.html">DefaultParamHandler</a></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> <span class="keyword">public</span>:</div>
<div class="line"><a name="l00065"></a><span class="lineno"><a class="line" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#aff2cbb90bd3ebac2208d25616c0ecf7f"> 65</a></span>  <span class="keyword">typedef</span> <a class="code" href="structOpenMS_1_1PSLPFormulation_1_1IndexTriple.html">PSLPFormulation::IndexTriple</a> <a class="code" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#aff2cbb90bd3ebac2208d25616c0ecf7f">IndexTriple</a>;</div>
<div class="line"><a name="l00066"></a><span class="lineno"> 66</span> </div>
<div class="line"><a name="l00067"></a><span class="lineno"> 67</span>  <a class="code" href="classOpenMS_1_1OfflinePrecursorIonSelection.html">OfflinePrecursorIonSelection</a>();</div>
<div class="line"><a name="l00068"></a><span class="lineno"> 68</span>  <span class="keyword">virtual</span> ~<a class="code" href="classOpenMS_1_1OfflinePrecursorIonSelection.html">OfflinePrecursorIonSelection</a>();</div>
<div class="line"><a name="l00069"></a><span class="lineno"> 69</span> </div>
<div class="line"><a name="l00079"></a><span class="lineno"> 79</span>  <span class="keyword">template</span> <<span class="keyword">typename</span> InputPeakType></div>
<div class="line"><a name="l00080"></a><span class="lineno"> 80</span>  <span class="keywordtype">void</span> makePrecursorSelectionForKnownLCMSMap(<span class="keyword">const</span> <a class="code" href="classOpenMS_1_1FeatureMap.html">FeatureMap<></a>& features,</div>
<div class="line"><a name="l00081"></a><span class="lineno"> 81</span>  <span class="keyword">const</span> <a class="code" href="classOpenMS_1_1MSExperiment.html">MSExperiment<InputPeakType></a>& experiment,</div>
<div class="line"><a name="l00082"></a><span class="lineno"> 82</span>  <a class="code" href="classOpenMS_1_1MSExperiment.html">MSExperiment<InputPeakType></a>& ms2,</div>
<div class="line"><a name="l00083"></a><span class="lineno"> 83</span>  std::set<Int>& charges_set,</div>
<div class="line"><a name="l00084"></a><span class="lineno"> 84</span>  <span class="keywordtype">bool</span> feature_based);</div>
<div class="line"><a name="l00085"></a><span class="lineno"> 85</span> </div>
<div class="line"><a name="l00093"></a><span class="lineno"> 93</span>  <span class="keyword">template</span> <<span class="keyword">typename</span> InputPeakType></div>
<div class="line"><a name="l00094"></a><span class="lineno"> 94</span>  <span class="keywordtype">void</span> getMassRanges(<span class="keyword">const</span> <a class="code" href="classOpenMS_1_1FeatureMap.html">FeatureMap<></a>& features,</div>
<div class="line"><a name="l00095"></a><span class="lineno"> 95</span>  <span class="keyword">const</span> <a class="code" href="classOpenMS_1_1MSExperiment.html">MSExperiment<InputPeakType></a>& experiment,</div>
<div class="line"><a name="l00096"></a><span class="lineno"> 96</span>  std::vector<std::vector<std::pair<Size, Size> > >& indices);</div>
<div class="line"><a name="l00097"></a><span class="lineno"> 97</span> </div>
<div class="line"><a name="l00098"></a><span class="lineno"> 98</span>  <span class="keywordtype">void</span> createProteinSequenceBasedLPInclusionList(<a class="code" href="classOpenMS_1_1String.html">String</a> include, <a class="code" href="classOpenMS_1_1String.html">String</a> rt_model_file, <a class="code" href="classOpenMS_1_1String.html">String</a> pt_model_file, <a class="code" href="classOpenMS_1_1FeatureMap.html">FeatureMap<></a>& precursors);</div>
<div class="line"><a name="l00099"></a><span class="lineno"> 99</span> </div>
<div class="line"><a name="l00100"></a><span class="lineno"><a class="line" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#a2b206fdbe8649b90b854129f37fc594a"> 100</a></span>  <span class="keywordtype">void</span> <a class="code" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#a2b206fdbe8649b90b854129f37fc594a">setLPSolver</a>(<a class="code" href="classOpenMS_1_1LPWrapper.html#ab659f8690ca779b8c3d73fe524a39a00">LPWrapper::SOLVER</a> solver)</div>
<div class="line"><a name="l00101"></a><span class="lineno"> 101</span>  {</div>
<div class="line"><a name="l00102"></a><span class="lineno"> 102</span>  solver_ = solver;</div>
<div class="line"><a name="l00103"></a><span class="lineno"> 103</span>  std::cout << <span class="stringliteral">" LPSolver set to "</span> << solver_ << std::endl;</div>
<div class="line"><a name="l00104"></a><span class="lineno"> 104</span>  }</div>
<div class="line"><a name="l00105"></a><span class="lineno"> 105</span> </div>
<div class="line"><a name="l00106"></a><span class="lineno"><a class="line" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#ac945a0966726a3a644d194b81c28d960"> 106</a></span>  <a class="code" href="classOpenMS_1_1LPWrapper.html#ab659f8690ca779b8c3d73fe524a39a00">LPWrapper::SOLVER</a> <a class="code" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#ac945a0966726a3a644d194b81c28d960">getLPSolver</a>()</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>  <span class="keywordflow">return</span> solver_;</div>
<div class="line"><a name="l00109"></a><span class="lineno"> 109</span>  }</div>
<div class="line"><a name="l00110"></a><span class="lineno"> 110</span> </div>
<div class="line"><a name="l00111"></a><span class="lineno"> 111</span> <span class="keyword">private</span>:</div>
<div class="line"><a name="l00115"></a><span class="lineno"> 115</span>  <span class="keyword">template</span> <<span class="keyword">typename</span> InputPeakType></div>
<div class="line"><a name="l00116"></a><span class="lineno"> 116</span>  <span class="keywordtype">void</span> calculateXICs_(<span class="keyword">const</span> <a class="code" href="classOpenMS_1_1FeatureMap.html">FeatureMap<></a>& features,</div>
<div class="line"><a name="l00117"></a><span class="lineno"> 117</span>  <span class="keyword">const</span> std::vector<std::vector<std::pair<Size, Size> > >& mass_ranges,</div>
<div class="line"><a name="l00118"></a><span class="lineno"> 118</span>  <span class="keyword">const</span> <a class="code" href="classOpenMS_1_1MSExperiment.html">MSExperiment<InputPeakType></a>& experiment,</div>
<div class="line"><a name="l00119"></a><span class="lineno"> 119</span>  <span class="keyword">const</span> std::set<Int>& charges_set,</div>
<div class="line"><a name="l00120"></a><span class="lineno"> 120</span>  std::vector<std::vector<std::pair<Size, DoubleReal> > >& xics);</div>
<div class="line"><a name="l00121"></a><span class="lineno"> 121</span> </div>
<div class="line"><a name="l00125"></a><span class="lineno"> 125</span>  <span class="keyword">template</span> <<span class="keyword">typename</span> InputPeakType></div>
<div class="line"><a name="l00126"></a><span class="lineno"> 126</span>  <span class="keywordtype">void</span> checkMassRanges_(std::vector<std::vector<std::pair<Size, Size> > >& mass_ranges,</div>
<div class="line"><a name="l00127"></a><span class="lineno"> 127</span>  <span class="keyword">const</span> <a class="code" href="classOpenMS_1_1MSExperiment.html">MSExperiment<InputPeakType></a>& experiment);</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>  <span class="keyword">template</span> <<span class="keyword">typename</span> T></div>
<div class="line"><a name="l00130"></a><span class="lineno"> 130</span>  <span class="keywordtype">void</span> updateExclusionList_(std::vector<std::pair<T, Size> >& exclusion_list);</div>
<div class="line"><a name="l00131"></a><span class="lineno"> 131</span> </div>
<div class="line"><a name="l00132"></a><span class="lineno"> 132</span>  <span class="keywordtype">void</span> updateExclusionList_(std::map<std::pair<DoubleReal, DoubleReal>, <a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a>, <a class="code" href="structOpenMS_1_1PairComparatorSecondElement.html">PairComparatorSecondElement</a><std::pair<DoubleReal, DoubleReal> > >& exclusion_list);</div>
<div class="line"><a name="l00133"></a><span class="lineno"> 133</span> </div>
<div class="line"><a name="l00134"></a><span class="lineno"><a class="line" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#afb8f1202d548f01394d85fbeb51f9adf"> 134</a></span>  <a class="code" href="classOpenMS_1_1LPWrapper.html#ab659f8690ca779b8c3d73fe524a39a00">LPWrapper::SOLVER</a> <a class="code" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#afb8f1202d548f01394d85fbeb51f9adf">solver_</a>;</div>
<div class="line"><a name="l00135"></a><span class="lineno"> 135</span>  };</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>  <span class="keyword">template</span> <<span class="keyword">typename</span> InputPeakType></div>
<div class="line"><a name="l00138"></a><span class="lineno"><a class="line" href="namespaceOpenMS.html#aa36b2ed6c06f4252dec15cd4f6d72635"> 138</a></span>  <span class="keywordtype">bool</span> <a class="code" href="namespaceOpenMS.html#aa36b2ed6c06f4252dec15cd4f6d72635">enclosesBoundingBox</a>(<span class="keyword">const</span> <a class="code" href="classOpenMS_1_1Feature.html">Feature</a>& f, <span class="keyword">typename</span> <a class="code" href="classdouble.html">MSExperiment<InputPeakType>::CoordinateType</a> rt, <span class="keyword">typename</span> <a class="code" href="classdouble.html">MSExperiment<InputPeakType>::CoordinateType</a> mz)</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>  <span class="keywordtype">bool</span> enclose_hit = <span class="keyword">false</span>;</div>
<div class="line"><a name="l00141"></a><span class="lineno"> 141</span>  <span class="keyword">const</span> std::vector<ConvexHull2D>& hulls = f.<a class="code" href="classOpenMS_1_1Feature.html#ad99bb8ba2103d6c52ee10ae6369b46e7">getConvexHulls</a>();</div>
<div class="line"><a name="l00142"></a><span class="lineno"> 142</span>  <span class="keywordflow">for</span> (<a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> i = 0; i < hulls.size(); ++i)</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">if</span> (hulls[i].getBoundingBox().encloses(rt, mz))</div>
<div class="line"><a name="l00145"></a><span class="lineno"> 145</span>  {</div>
<div class="line"><a name="l00146"></a><span class="lineno"> 146</span>  enclose_hit = <span class="keyword">true</span>;</div>
<div class="line"><a name="l00147"></a><span class="lineno"> 147</span>  <span class="keywordflow">return</span> enclose_hit;</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>  <span class="keywordflow">return</span> enclose_hit;</div>
<div class="line"><a name="l00151"></a><span class="lineno"> 151</span>  }</div>
<div class="line"><a name="l00152"></a><span class="lineno"> 152</span> </div>
<div class="line"><a name="l00153"></a><span class="lineno"> 153</span>  <span class="keyword">template</span> <<span class="keyword">typename</span> InputPeakType></div>
<div class="line"><a name="l00154"></a><span class="lineno"><a class="line" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#ae895e8f290c7f5869249e59a7b53c417"> 154</a></span>  <span class="keywordtype">void</span> <a class="code" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#ae895e8f290c7f5869249e59a7b53c417">OfflinePrecursorIonSelection::getMassRanges</a>(<span class="keyword">const</span> <a class="code" href="classOpenMS_1_1FeatureMap.html">FeatureMap<></a>& features,</div>
<div class="line"><a name="l00155"></a><span class="lineno"> 155</span>  <span class="keyword">const</span> <a class="code" href="classOpenMS_1_1MSExperiment.html">MSExperiment<InputPeakType></a>& experiment,</div>
<div class="line"><a name="l00156"></a><span class="lineno"> 156</span>  std::vector<std::vector<std::pair<Size, Size> > >& indices)</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>  <span class="keywordflow">if</span> (experiment.<a class="code" href="classOpenMS_1_1MSExperiment.html#ac6e61de369e994009e36f344f99c15ad">empty</a>())</div>
<div class="line"><a name="l00159"></a><span class="lineno"> 159</span>  <span class="keywordflow">throw</span> <a class="code" href="classOpenMS_1_1Exception_1_1InvalidSize.html">Exception::InvalidSize</a>(__FILE__, __LINE__, __PRETTY_FUNCTION__, 0);</div>
<div class="line"><a name="l00160"></a><span class="lineno"> 160</span>  <span class="keywordflow">for</span> (<a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> f = 0; f < features.size(); ++f)</div>
<div class="line"><a name="l00161"></a><span class="lineno"> 161</span>  {</div>
<div class="line"><a name="l00162"></a><span class="lineno"> 162</span>  std::vector<std::pair<Size, Size> > vec;</div>
<div class="line"><a name="l00163"></a><span class="lineno"> 163</span> </div>
<div class="line"><a name="l00164"></a><span class="lineno"> 164</span>  <span class="keywordflow">for</span> (<a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> rt = 0; rt < experiment.<a class="code" href="classOpenMS_1_1MSExperiment.html#a9d78a687cf2a391198bb3cbc08bc06cb">size</a>(); ++rt)</div>
<div class="line"><a name="l00165"></a><span class="lineno"> 165</span>  {</div>
<div class="line"><a name="l00166"></a><span class="lineno"> 166</span>  <span class="comment">// is scan relevant?</span></div>
<div class="line"><a name="l00167"></a><span class="lineno"> 167</span>  <span class="keywordflow">if</span> (!enclosesBoundingBox<InputPeakType>(features[f], experiment[rt].getRT(), features[f].getMZ()))</div>
<div class="line"><a name="l00168"></a><span class="lineno"> 168</span>  <span class="keywordflow">continue</span>;</div>
<div class="line"><a name="l00169"></a><span class="lineno"> 169</span> </div>
<div class="line"><a name="l00170"></a><span class="lineno"> 170</span>  std::pair<Size, Size> start;</div>
<div class="line"><a name="l00171"></a><span class="lineno"> 171</span>  std::pair<Size, Size> end;</div>
<div class="line"><a name="l00172"></a><span class="lineno"> 172</span>  <span class="keywordtype">bool</span> start_found = <span class="keyword">false</span>;</div>
<div class="line"><a name="l00173"></a><span class="lineno"> 173</span>  <span class="keywordtype">bool</span> end_found = <span class="keyword">false</span>;</div>
<div class="line"><a name="l00174"></a><span class="lineno"> 174</span>  <span class="keyword">typename</span> <a class="code" href="classOpenMS_1_1MSSpectrum.html#a4895afe86c2aa0747ae0411c58c82114">MSSpectrum<InputPeakType>::ConstIterator</a> mz_iter = experiment[rt].MZBegin(features[f].getMZ());</div>
<div class="line"><a name="l00175"></a><span class="lineno"> 175</span>  <span class="keyword">typename</span> <a class="code" href="classOpenMS_1_1MSSpectrum.html#a4895afe86c2aa0747ae0411c58c82114">MSSpectrum<InputPeakType>::ConstIterator</a> mz_end = mz_iter;</div>
<div class="line"><a name="l00176"></a><span class="lineno"> 176</span>  <span class="keywordflow">if</span> (mz_iter == experiment[rt].end())</div>
<div class="line"><a name="l00177"></a><span class="lineno"> 177</span>  <span class="keywordflow">continue</span>;</div>
<div class="line"><a name="l00178"></a><span class="lineno"> 178</span>  <span class="comment">// check to the left</span></div>
<div class="line"><a name="l00179"></a><span class="lineno"> 179</span>  <span class="keywordflow">while</span> (enclosesBoundingBox<InputPeakType>(features[f], experiment[rt].getRT(), mz_iter->getMZ()))</div>
<div class="line"><a name="l00180"></a><span class="lineno"> 180</span>  {</div>
<div class="line"><a name="l00181"></a><span class="lineno"> 181</span>  start_found = <span class="keyword">true</span>;</div>
<div class="line"><a name="l00182"></a><span class="lineno"> 182</span>  start.first = rt;</div>
<div class="line"><a name="l00183"></a><span class="lineno"> 183</span>  start.second = distance(experiment[rt].begin(), mz_iter);</div>
<div class="line"><a name="l00184"></a><span class="lineno"> 184</span>  <span class="keywordflow">if</span> (mz_iter == experiment[rt].begin())</div>
<div class="line"><a name="l00185"></a><span class="lineno"> 185</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00186"></a><span class="lineno"> 186</span>  --mz_iter;</div>
<div class="line"><a name="l00187"></a><span class="lineno"> 187</span>  end_found = <span class="keyword">true</span>;</div>
<div class="line"><a name="l00188"></a><span class="lineno"> 188</span>  end.first = rt;</div>
<div class="line"><a name="l00189"></a><span class="lineno"> 189</span>  end.second = start.second;</div>
<div class="line"><a name="l00190"></a><span class="lineno"> 190</span>  }</div>
<div class="line"><a name="l00191"></a><span class="lineno"> 191</span>  <span class="comment">// and now to the right</span></div>
<div class="line"><a name="l00192"></a><span class="lineno"> 192</span>  <span class="keywordflow">while</span> (mz_end != experiment[rt].end() && enclosesBoundingBox<InputPeakType>(features[f], experiment[rt].getRT(), mz_end->getMZ()))</div>
<div class="line"><a name="l00193"></a><span class="lineno"> 193</span>  {</div>
<div class="line"><a name="l00194"></a><span class="lineno"> 194</span>  end_found = <span class="keyword">true</span>;</div>
<div class="line"><a name="l00195"></a><span class="lineno"> 195</span>  end.first = rt;</div>
<div class="line"><a name="l00196"></a><span class="lineno"> 196</span>  end.second = distance(experiment[rt].begin(), mz_end);</div>
<div class="line"><a name="l00197"></a><span class="lineno"> 197</span>  ++mz_end;</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> (start_found && end_found)</div>
<div class="line"><a name="l00200"></a><span class="lineno"> 200</span>  {</div>
<div class="line"><a name="l00201"></a><span class="lineno"> 201</span>  vec.push_back(start);</div>
<div class="line"><a name="l00202"></a><span class="lineno"> 202</span>  vec.push_back(end);</div>
<div class="line"><a name="l00203"></a><span class="lineno"> 203</span>  }</div>
<div class="line"><a name="l00204"></a><span class="lineno"> 204</span> <span class="preprocessor">#ifdef DEBUG_OPS</span></div>
<div class="line"><a name="l00205"></a><span class="lineno"> 205</span> <span class="preprocessor"></span> <span class="keywordflow">else</span> <span class="keywordflow">if</span> (start_found || end_found)</div>
<div class="line"><a name="l00206"></a><span class="lineno"> 206</span>  {</div>
<div class="line"><a name="l00207"></a><span class="lineno"> 207</span>  std::cout << <span class="stringliteral">"start "</span> << start_found << <span class="stringliteral">" end "</span> << end_found << std::endl;</div>
<div class="line"><a name="l00208"></a><span class="lineno"> 208</span>  std::cout << <span class="stringliteral">"feature: "</span> << f << <span class="stringliteral">" rt: "</span> << rt << std::endl;</div>
<div class="line"><a name="l00209"></a><span class="lineno"> 209</span>  }</div>
<div class="line"><a name="l00210"></a><span class="lineno"> 210</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00211"></a><span class="lineno"> 211</span> <span class="preprocessor"></span> }</div>
<div class="line"><a name="l00212"></a><span class="lineno"> 212</span> <span class="preprocessor">#ifdef DEBUG_OPS</span></div>
<div class="line"><a name="l00213"></a><span class="lineno"> 213</span> <span class="preprocessor"></span> <span class="keywordflow">if</span> (vec.size() > 0)</div>
<div class="line"><a name="l00214"></a><span class="lineno"> 214</span>  {</div>
<div class="line"><a name="l00215"></a><span class="lineno"> 215</span>  std::cout << vec.size() << <span class="stringliteral">" / 2 scans"</span> << std::endl;</div>
<div class="line"><a name="l00216"></a><span class="lineno"> 216</span>  <span class="keywordflow">for</span> (<a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> i = 0; i < vec.size(); i += 2)</div>
<div class="line"><a name="l00217"></a><span class="lineno"> 217</span>  {</div>
<div class="line"><a name="l00218"></a><span class="lineno"> 218</span>  std::cout << <span class="stringliteral">"Feature "</span> << f << <span class="stringliteral">" RT : "</span> << vec[i].first</div>
<div class="line"><a name="l00219"></a><span class="lineno"> 219</span>  << <span class="stringliteral">" MZ : "</span> << experiment[vec[i].first][vec[i].second].getMZ() << <span class="stringliteral">" "</span></div>
<div class="line"><a name="l00220"></a><span class="lineno"> 220</span>  << experiment[vec[i + 1].first][vec[i + 1].second].getMZ() << std::endl;</div>
<div class="line"><a name="l00221"></a><span class="lineno"> 221</span>  }</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> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00224"></a><span class="lineno"> 224</span> <span class="preprocessor"></span> <span class="keywordflow">if</span> (vec.empty())</div>
<div class="line"><a name="l00225"></a><span class="lineno"> 225</span>  {</div>
<div class="line"><a name="l00226"></a><span class="lineno"> 226</span> <span class="preprocessor">#ifdef DEBUG_OPS</span></div>
<div class="line"><a name="l00227"></a><span class="lineno"> 227</span> <span class="preprocessor"></span> std::cout << <span class="stringliteral">"According to the convex hulls no mass traces found for this feature->estimate!"</span></div>
<div class="line"><a name="l00228"></a><span class="lineno"> 228</span>  << features[f].getRT() << <span class="stringliteral">" "</span> << features[f].getMZ() << <span class="stringliteral">" "</span> << features[f].getCharge() << std::endl;</div>
<div class="line"><a name="l00229"></a><span class="lineno"> 229</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00230"></a><span class="lineno"> 230</span> <span class="preprocessor"></span> <span class="comment">// we estimate the convex hull</span></div>
<div class="line"><a name="l00231"></a><span class="lineno"> 231</span>  <span class="keyword">typename</span> <a class="code" href="classOpenMS_1_1MSExperiment.html#aee0a1a68d0fab63eefa6e14010b4a7f0">MSExperiment<InputPeakType>::ConstIterator</a> spec_iter = experiment.<a class="code" href="classOpenMS_1_1MSExperiment.html#adb2f0f3e10754ca4dc7be677f15689a3">RTBegin</a>(features[f].getRT());</div>
<div class="line"><a name="l00232"></a><span class="lineno"> 232</span>  <span class="keywordflow">if</span> (spec_iter == experiment.<a class="code" href="classOpenMS_1_1MSExperiment.html#ab45dae688fc5d8983727abffa4389003">end</a>())</div>
<div class="line"><a name="l00233"></a><span class="lineno"> 233</span>  --spec_iter;</div>
<div class="line"><a name="l00234"></a><span class="lineno"> 234</span> </div>
<div class="line"><a name="l00235"></a><span class="lineno"> 235</span>  <a class="code" href="classdouble.html">DoubleReal</a> dist1 = fabs(spec_iter->getRT() - features[f].getRT());</div>
<div class="line"><a name="l00236"></a><span class="lineno"> 236</span>  <a class="code" href="classdouble.html">DoubleReal</a> dist2 = std::numeric_limits<DoubleReal>::max();</div>
<div class="line"><a name="l00237"></a><span class="lineno"> 237</span>  <a class="code" href="classdouble.html">DoubleReal</a> dist3 = std::numeric_limits<DoubleReal>::max();</div>
<div class="line"><a name="l00238"></a><span class="lineno"> 238</span>  <span class="keywordflow">if</span> ((spec_iter + 1) != experiment.<a class="code" href="classOpenMS_1_1MSExperiment.html#ab45dae688fc5d8983727abffa4389003">end</a>())</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>  dist2 = fabs((spec_iter + 1)->getRT() - features[f].getRT());</div>
<div class="line"><a name="l00241"></a><span class="lineno"> 241</span>  }</div>
<div class="line"><a name="l00242"></a><span class="lineno"> 242</span>  <span class="keywordflow">if</span> (spec_iter != experiment.<a class="code" href="classOpenMS_1_1MSExperiment.html#a2387033802383edbdc95f9bbb12a707e">begin</a>())</div>
<div class="line"><a name="l00243"></a><span class="lineno"> 243</span>  {</div>
<div class="line"><a name="l00244"></a><span class="lineno"> 244</span>  dist3 = fabs((spec_iter - 1)->getRT() - features[f].getRT());</div>
<div class="line"><a name="l00245"></a><span class="lineno"> 245</span>  }</div>
<div class="line"><a name="l00246"></a><span class="lineno"> 246</span>  <span class="keywordflow">if</span> (dist3 <= dist1 && dist3 <= dist2)</div>
<div class="line"><a name="l00247"></a><span class="lineno"> 247</span>  {</div>
<div class="line"><a name="l00248"></a><span class="lineno"> 248</span>  --spec_iter;</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">else</span> <span class="keywordflow">if</span> (dist2 <= dist3 && dist2 <= dist1)</div>
<div class="line"><a name="l00251"></a><span class="lineno"> 251</span>  {</div>
<div class="line"><a name="l00252"></a><span class="lineno"> 252</span>  ++spec_iter;</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>  std::pair<Size, Size> start;</div>
<div class="line"><a name="l00255"></a><span class="lineno"> 255</span>  std::pair<Size, Size> end;</div>
<div class="line"><a name="l00256"></a><span class="lineno"> 256</span>  start.first = distance(experiment.<a class="code" href="classOpenMS_1_1MSExperiment.html#a2387033802383edbdc95f9bbb12a707e">begin</a>(), spec_iter);</div>
<div class="line"><a name="l00257"></a><span class="lineno"> 257</span>  end.first = start.first;</div>
<div class="line"><a name="l00258"></a><span class="lineno"> 258</span> </div>
<div class="line"><a name="l00259"></a><span class="lineno"> 259</span>  <span class="keyword">typename</span> <a class="code" href="classOpenMS_1_1MSSpectrum.html#a4895afe86c2aa0747ae0411c58c82114">MSSpectrum<InputPeakType>::ConstIterator</a> mz_iter = spec_iter->MZBegin(features[f].getMZ());</div>
<div class="line"><a name="l00260"></a><span class="lineno"> 260</span>  <span class="keywordflow">if</span> (spec_iter-><a class="code" href="classOpenMS_1_1MSExperiment.html#a2387033802383edbdc95f9bbb12a707e">begin</a>() == spec_iter-><a class="code" href="classOpenMS_1_1MSExperiment.html#ab45dae688fc5d8983727abffa4389003">end</a>())</div>
<div class="line"><a name="l00261"></a><span class="lineno"> 261</span>  {</div>
<div class="line"><a name="l00262"></a><span class="lineno"> 262</span>  indices.push_back(vec);</div>
<div class="line"><a name="l00263"></a><span class="lineno"> 263</span>  <span class="keywordflow">continue</span>;</div>
<div class="line"><a name="l00264"></a><span class="lineno"> 264</span>  }</div>
<div class="line"><a name="l00265"></a><span class="lineno"> 265</span>  <span class="keywordflow">if</span> (mz_iter == spec_iter-><a class="code" href="classOpenMS_1_1MSExperiment.html#ab45dae688fc5d8983727abffa4389003">end</a>() || (mz_iter->getMZ() > features[f].getMZ() && mz_iter != spec_iter-><a class="code" href="classOpenMS_1_1MSExperiment.html#a2387033802383edbdc95f9bbb12a707e">begin</a>()))</div>
<div class="line"><a name="l00266"></a><span class="lineno"> 266</span>  --mz_iter;</div>
<div class="line"><a name="l00267"></a><span class="lineno"> 267</span>  <span class="keywordflow">while</span> (mz_iter != spec_iter-><a class="code" href="classOpenMS_1_1MSExperiment.html#a2387033802383edbdc95f9bbb12a707e">begin</a>())</div>
<div class="line"><a name="l00268"></a><span class="lineno"> 268</span>  {</div>
<div class="line"><a name="l00269"></a><span class="lineno"> 269</span>  <span class="keywordflow">if</span> (fabs((mz_iter - 1)->getMZ() - features[f].getMZ()) < 0.5)</div>
<div class="line"><a name="l00270"></a><span class="lineno"> 270</span>  --mz_iter;</div>
<div class="line"><a name="l00271"></a><span class="lineno"> 271</span>  <span class="keywordflow">else</span></div>
<div class="line"><a name="l00272"></a><span class="lineno"> 272</span>  <span class="keywordflow">break</span>;</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>  start.second = distance(spec_iter-><a class="code" href="classOpenMS_1_1MSExperiment.html#a2387033802383edbdc95f9bbb12a707e">begin</a>(), mz_iter);</div>
<div class="line"><a name="l00275"></a><span class="lineno"> 275</span>  <span class="keyword">typename</span> <a class="code" href="classOpenMS_1_1MSSpectrum.html#a4895afe86c2aa0747ae0411c58c82114">MSSpectrum<InputPeakType>::ConstIterator</a> mz_end = mz_iter;</div>
<div class="line"><a name="l00276"></a><span class="lineno"> 276</span> <span class="preprocessor">#ifdef DEBUG_OPS</span></div>
<div class="line"><a name="l00277"></a><span class="lineno"> 277</span> <span class="preprocessor"></span> std::cout << features[f].getMZ() << <span class="stringliteral">" Start: "</span> << experiment[start.first].getRT() << <span class="stringliteral">" "</span> << experiment[start.first][start.second].getMZ();</div>
<div class="line"><a name="l00278"></a><span class="lineno"> 278</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00279"></a><span class="lineno"> 279</span> <span class="preprocessor"></span> <a class="code" href="group__Concept.html#ga7cc214a236ad3bb6ad435bdcf5262a3f">Int</a> charge = features[f].getCharge();</div>
<div class="line"><a name="l00280"></a><span class="lineno"> 280</span>  <span class="keywordflow">if</span> (charge == 0)</div>
<div class="line"><a name="l00281"></a><span class="lineno"> 281</span>  charge = 1;</div>
<div class="line"><a name="l00282"></a><span class="lineno"> 282</span>  <span class="keywordflow">while</span> (mz_end + 1 != spec_iter-><a class="code" href="classOpenMS_1_1MSExperiment.html#ab45dae688fc5d8983727abffa4389003">end</a>())</div>
<div class="line"><a name="l00283"></a><span class="lineno"> 283</span>  {</div>
<div class="line"><a name="l00284"></a><span class="lineno"> 284</span>  <span class="keywordflow">if</span> (fabs((mz_end + 1)->getMZ() - features[f].getMZ()) < 3.0 / (<a class="code" href="group__Concept.html#gace75bfb1aba684e874dffee13738bd15">DoubleReal</a>)charge)</div>
<div class="line"><a name="l00285"></a><span class="lineno"> 285</span>  ++mz_end;</div>
<div class="line"><a name="l00286"></a><span class="lineno"> 286</span>  <span class="keywordflow">else</span></div>
<div class="line"><a name="l00287"></a><span class="lineno"> 287</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00288"></a><span class="lineno"> 288</span>  }</div>
<div class="line"><a name="l00289"></a><span class="lineno"> 289</span>  end.second = distance(spec_iter-><a class="code" href="classOpenMS_1_1MSExperiment.html#a2387033802383edbdc95f9bbb12a707e">begin</a>(), mz_end);</div>
<div class="line"><a name="l00290"></a><span class="lineno"> 290</span> <span class="preprocessor">#ifdef DEBUG_OPS</span></div>
<div class="line"><a name="l00291"></a><span class="lineno"> 291</span> <span class="preprocessor"></span> std::cout << <span class="stringliteral">"\tEnd: "</span> << experiment[end.first].getRT() << <span class="stringliteral">" "</span> << experiment[end.first][end.second].getMZ() << std::endl;</div>
<div class="line"><a name="l00292"></a><span class="lineno"> 292</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00293"></a><span class="lineno"> 293</span> <span class="preprocessor"></span> vec.push_back(start);</div>
<div class="line"><a name="l00294"></a><span class="lineno"> 294</span>  vec.push_back(end);</div>
<div class="line"><a name="l00295"></a><span class="lineno"> 295</span>  }</div>
<div class="line"><a name="l00296"></a><span class="lineno"> 296</span> </div>
<div class="line"><a name="l00297"></a><span class="lineno"> 297</span>  indices.push_back(vec);</div>
<div class="line"><a name="l00298"></a><span class="lineno"> 298</span>  }</div>
<div class="line"><a name="l00299"></a><span class="lineno"> 299</span>  <span class="comment">// eliminate nearby peaks</span></div>
<div class="line"><a name="l00300"></a><span class="lineno"> 300</span>  <span class="keywordflow">if</span> (<a class="code" href="classOpenMS_1_1DefaultParamHandler.html#a28c73e623c63a4fe3bfceb1ae8274f39">param_</a>.<a class="code" href="classOpenMS_1_1Param.html#abd0f5571675ec99e1aa647b39e2a8711">getValue</a>(<span class="stringliteral">"exclude_overlapping_peaks"</span>) == <span class="stringliteral">"true"</span>)</div>
<div class="line"><a name="l00301"></a><span class="lineno"> 301</span>  <a class="code" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#acca2e18b840cc4b903f2958efca7aa1a">checkMassRanges_</a>(indices, experiment);</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> </div>
<div class="line"><a name="l00304"></a><span class="lineno"> 304</span>  <span class="keyword">template</span> <<span class="keyword">typename</span> InputPeakType></div>
<div class="line"><a name="l00305"></a><span class="lineno"><a class="line" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#a09cf655b910fa7166cd8b27d17e02453"> 305</a></span>  <span class="keywordtype">void</span> <a class="code" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#a09cf655b910fa7166cd8b27d17e02453">OfflinePrecursorIonSelection::calculateXICs_</a>(<span class="keyword">const</span> <a class="code" href="classOpenMS_1_1FeatureMap.html">FeatureMap<></a>& features,</div>
<div class="line"><a name="l00306"></a><span class="lineno"> 306</span>  <span class="keyword">const</span> std::vector<std::vector<std::pair<Size, Size> > >& mass_ranges,</div>
<div class="line"><a name="l00307"></a><span class="lineno"> 307</span>  <span class="keyword">const</span> <a class="code" href="classOpenMS_1_1MSExperiment.html">MSExperiment<InputPeakType></a>& experiment,</div>
<div class="line"><a name="l00308"></a><span class="lineno"> 308</span>  <span class="keyword">const</span> std::set<Int>& charges_set,</div>
<div class="line"><a name="l00309"></a><span class="lineno"> 309</span>  std::vector<std::vector<std::pair<Size, DoubleReal> > >& xics)</div>
<div class="line"><a name="l00310"></a><span class="lineno"> 310</span>  {</div>
<div class="line"><a name="l00311"></a><span class="lineno"> 311</span>  xics.clear();</div>
<div class="line"><a name="l00312"></a><span class="lineno"> 312</span>  xics.resize(experiment.<a class="code" href="classOpenMS_1_1MSExperiment.html#a9d78a687cf2a391198bb3cbc08bc06cb">size</a>());</div>
<div class="line"><a name="l00313"></a><span class="lineno"> 313</span>  <span class="comment">// for each feature</span></div>
<div class="line"><a name="l00314"></a><span class="lineno"> 314</span>  <span class="keywordflow">for</span> (<a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> f = 0; f < mass_ranges.size(); ++f)</div>
<div class="line"><a name="l00315"></a><span class="lineno"> 315</span>  {</div>
<div class="line"><a name="l00316"></a><span class="lineno"> 316</span>  <span class="comment">// is charge valid</span></div>
<div class="line"><a name="l00317"></a><span class="lineno"> 317</span>  <span class="keywordflow">if</span> (charges_set.count(features[f].getCharge()) < 1)</div>
<div class="line"><a name="l00318"></a><span class="lineno"> 318</span>  {</div>
<div class="line"><a name="l00319"></a><span class="lineno"> 319</span>  <span class="keywordflow">continue</span>;</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>  <span class="comment">// go through all scans where the feature occurs</span></div>
<div class="line"><a name="l00322"></a><span class="lineno"> 322</span>  <span class="keywordflow">for</span> (<a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> s = 0; s < mass_ranges[f].size(); s += 2)</div>
<div class="line"><a name="l00323"></a><span class="lineno"> 323</span>  {</div>
<div class="line"><a name="l00324"></a><span class="lineno"> 324</span>  <span class="comment">// sum intensity over all raw datapoints belonging to the feature in the current scan</span></div>
<div class="line"><a name="l00325"></a><span class="lineno"> 325</span>  <a class="code" href="classdouble.html">DoubleReal</a> weight = 0.;</div>
<div class="line"><a name="l00326"></a><span class="lineno"> 326</span>  <span class="keywordflow">for</span> (<a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> j = mass_ranges[f][s].second; j <= mass_ranges[f][s + 1].second; ++j)</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>  weight += experiment[mass_ranges[f][s].first][j].getIntensity();</div>
<div class="line"><a name="l00329"></a><span class="lineno"> 329</span>  }</div>
<div class="line"><a name="l00330"></a><span class="lineno"> 330</span>  <span class="comment">// enter xic in the vector for scan</span></div>
<div class="line"><a name="l00331"></a><span class="lineno"> 331</span>  xics[mass_ranges[f][s].first].push_back(std::make_pair(f, weight));</div>
<div class="line"><a name="l00332"></a><span class="lineno"> 332</span>  }</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> </div>
<div class="line"><a name="l00335"></a><span class="lineno"> 335</span>  <span class="keywordflow">for</span> (<a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> s = 0; s < xics.size(); ++s)</div>
<div class="line"><a name="l00336"></a><span class="lineno"> 336</span>  {</div>
<div class="line"><a name="l00337"></a><span class="lineno"> 337</span>  sort(xics[s].begin(), xics[s].end(), <a class="code" href="structOpenMS_1_1PairComparatorSecondElement.html">PairComparatorSecondElement</a><std::pair<Size, DoubleReal> >());</div>
<div class="line"><a name="l00338"></a><span class="lineno"> 338</span>  }</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> </div>
<div class="line"><a name="l00341"></a><span class="lineno"> 341</span>  <span class="keyword">template</span> <<span class="keyword">typename</span> InputPeakType></div>
<div class="line"><a name="l00342"></a><span class="lineno"><a class="line" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#afc3ec85793adf0cdaacc5803e59f6c17"> 342</a></span>  <span class="keywordtype">void</span> <a class="code" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#afc3ec85793adf0cdaacc5803e59f6c17">OfflinePrecursorIonSelection::makePrecursorSelectionForKnownLCMSMap</a>(<span class="keyword">const</span> <a class="code" href="classOpenMS_1_1FeatureMap.html">FeatureMap<></a>& features,</div>
<div class="line"><a name="l00343"></a><span class="lineno"> 343</span>  <span class="keyword">const</span> <a class="code" href="classOpenMS_1_1MSExperiment.html">MSExperiment<InputPeakType></a>& experiment,</div>
<div class="line"><a name="l00344"></a><span class="lineno"> 344</span>  <a class="code" href="classOpenMS_1_1MSExperiment.html">MSExperiment<InputPeakType></a>& ms2,</div>
<div class="line"><a name="l00345"></a><span class="lineno"> 345</span>  std::set<Int>& charges_set,</div>
<div class="line"><a name="l00346"></a><span class="lineno"> 346</span>  <span class="keywordtype">bool</span> feature_based)</div>
<div class="line"><a name="l00347"></a><span class="lineno"> 347</span>  {</div>
<div class="line"><a name="l00348"></a><span class="lineno"> 348</span> </div>
<div class="line"><a name="l00349"></a><span class="lineno"> 349</span>  <span class="keyword">const</span> <a class="code" href="classdouble.html">DoubleReal</a> window = <a class="code" href="classOpenMS_1_1DefaultParamHandler.html#a28c73e623c63a4fe3bfceb1ae8274f39">param_</a>.<a class="code" href="classOpenMS_1_1Param.html#abd0f5571675ec99e1aa647b39e2a8711">getValue</a>(<span class="stringliteral">"selection_window"</span>);</div>
<div class="line"><a name="l00350"></a><span class="lineno"> 350</span>  <span class="keyword">const</span> <a class="code" href="classdouble.html">DoubleReal</a> excl_window = <a class="code" href="classOpenMS_1_1DefaultParamHandler.html#a28c73e623c63a4fe3bfceb1ae8274f39">param_</a>.<a class="code" href="classOpenMS_1_1Param.html#abd0f5571675ec99e1aa647b39e2a8711">getValue</a>(<span class="stringliteral">"min_peak_distance"</span>);</div>
<div class="line"><a name="l00351"></a><span class="lineno"> 351</span> </div>
<div class="line"><a name="l00352"></a><span class="lineno"> 352</span>  <span class="comment">// get the mass ranges for each features for each scan it occurs in</span></div>
<div class="line"><a name="l00353"></a><span class="lineno"> 353</span>  std::vector<std::vector<std::pair<Size, Size> > > indices;</div>
<div class="line"><a name="l00354"></a><span class="lineno"> 354</span>  <a class="code" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#ae895e8f290c7f5869249e59a7b53c417">getMassRanges</a>(features, experiment, indices);</div>
<div class="line"><a name="l00355"></a><span class="lineno"> 355</span>  <a class="code" href="classdouble.html">DoubleReal</a> rt_dist = 0.;</div>
<div class="line"><a name="l00356"></a><span class="lineno"> 356</span>  <span class="keywordflow">if</span> (experiment.<a class="code" href="classOpenMS_1_1MSExperiment.html#a9d78a687cf2a391198bb3cbc08bc06cb">size</a>() > 1)</div>
<div class="line"><a name="l00357"></a><span class="lineno"> 357</span>  {</div>
<div class="line"><a name="l00358"></a><span class="lineno"> 358</span>  rt_dist = experiment[1].getRT() - experiment[0].getRT();</div>
<div class="line"><a name="l00359"></a><span class="lineno"> 359</span>  }</div>
<div class="line"><a name="l00360"></a><span class="lineno"> 360</span> </div>
<div class="line"><a name="l00361"></a><span class="lineno"> 361</span>  <span class="comment">// feature based selection (e.g. with LC-MALDI)</span></div>
<div class="line"><a name="l00362"></a><span class="lineno"> 362</span>  <span class="keywordflow">if</span> (feature_based)</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">// create ILP</span></div>
<div class="line"><a name="l00365"></a><span class="lineno"> 365</span>  <a class="code" href="classOpenMS_1_1PSLPFormulation.html">PSLPFormulation</a> ilp_wrapper;</div>
<div class="line"><a name="l00366"></a><span class="lineno"> 366</span> </div>
<div class="line"><a name="l00367"></a><span class="lineno"> 367</span>  std::vector<IndexTriple> variable_indices;</div>
<div class="line"><a name="l00368"></a><span class="lineno"> 368</span>  std::vector<int> solution_indices;</div>
<div class="line"><a name="l00369"></a><span class="lineno"> 369</span>  ilp_wrapper.<a class="code" href="classOpenMS_1_1PSLPFormulation.html#a2a93e2928cdbb7704c7d49a33f13b091">createAndSolveILPForKnownLCMSMapFeatureBased</a>(features, experiment, variable_indices,</div>
<div class="line"><a name="l00370"></a><span class="lineno"> 370</span>  indices, charges_set,</div>
<div class="line"><a name="l00371"></a><span class="lineno"> 371</span>  <a class="code" href="classOpenMS_1_1DefaultParamHandler.html#a28c73e623c63a4fe3bfceb1ae8274f39">param_</a>.<a class="code" href="classOpenMS_1_1Param.html#abd0f5571675ec99e1aa647b39e2a8711">getValue</a>(<span class="stringliteral">"ms2_spectra_per_rt_bin"</span>),</div>
<div class="line"><a name="l00372"></a><span class="lineno"> 372</span>  solution_indices);</div>
<div class="line"><a name="l00373"></a><span class="lineno"> 373</span> </div>
<div class="line"><a name="l00374"></a><span class="lineno"> 374</span>  sort(variable_indices.begin(), variable_indices.end(), <a class="code" href="structOpenMS_1_1PSLPFormulation_1_1IndexLess.html">PSLPFormulation::IndexLess</a>());</div>
<div class="line"><a name="l00375"></a><span class="lineno"> 375</span> <span class="preprocessor">#ifdef DEBUG_OPS</span></div>
<div class="line"><a name="l00376"></a><span class="lineno"> 376</span> <span class="preprocessor"></span> std::cout << <span class="stringliteral">"best_solution "</span> << std::endl;</div>
<div class="line"><a name="l00377"></a><span class="lineno"> 377</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00378"></a><span class="lineno"> 378</span> <span class="preprocessor"></span> <span class="comment">// print best solution</span></div>
<div class="line"><a name="l00379"></a><span class="lineno"> 379</span>  <span class="comment">// create inclusion list</span></div>
<div class="line"><a name="l00380"></a><span class="lineno"> 380</span>  <span class="keywordflow">for</span> (<a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> i = 0; i < solution_indices.size(); ++i)</div>
<div class="line"><a name="l00381"></a><span class="lineno"> 381</span>  {</div>
<div class="line"><a name="l00382"></a><span class="lineno"> 382</span>  <a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> feature_index = variable_indices[solution_indices[i]].feature;</div>
<div class="line"><a name="l00383"></a><span class="lineno"> 383</span>  <a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> feature_scan_idx = variable_indices[solution_indices[i]].scan;</div>
<div class="line"><a name="l00384"></a><span class="lineno"> 384</span>  <span class="keyword">typename</span> <a class="code" href="classOpenMS_1_1MSExperiment.html#aee0a1a68d0fab63eefa6e14010b4a7f0">MSExperiment<InputPeakType>::ConstIterator</a> scan = experiment.<a class="code" href="classOpenMS_1_1MSExperiment.html#a2387033802383edbdc95f9bbb12a707e">begin</a>() + feature_scan_idx;</div>
<div class="line"><a name="l00385"></a><span class="lineno"> 385</span>  <span class="keyword">typename</span> <a class="code" href="classOpenMS_1_1MSSpectrum.html">MSExperiment<InputPeakType>::SpectrumType</a> ms2_spec;</div>
<div class="line"><a name="l00386"></a><span class="lineno"> 386</span>  <a class="code" href="classOpenMS_1_1Precursor.html">Precursor</a> p;</div>
<div class="line"><a name="l00387"></a><span class="lineno"> 387</span>  std::vector<Precursor> pcs;</div>
<div class="line"><a name="l00388"></a><span class="lineno"> 388</span>  p.<a class="code" href="classOpenMS_1_1Peak1D.html#ab875c2479c6cd2a0189e4684bfa4f97d">setIntensity</a>(features[feature_index].getIntensity());</div>
<div class="line"><a name="l00389"></a><span class="lineno"> 389</span>  p.<a class="code" href="classOpenMS_1_1Peak1D.html#ab9c4ce5d172a75d6de2c6f5f1847c9ab">setMZ</a>(features[feature_index].getMZ());</div>
<div class="line"><a name="l00390"></a><span class="lineno"> 390</span>  p.<a class="code" href="classOpenMS_1_1Precursor.html#aca40a4630f5f8a0568899b74a62dcd4b">setCharge</a>(features[feature_index].getCharge());</div>
<div class="line"><a name="l00391"></a><span class="lineno"> 391</span>  pcs.push_back(p);</div>
<div class="line"><a name="l00392"></a><span class="lineno"> 392</span>  ms2_spec.<a class="code" href="classOpenMS_1_1SpectrumSettings.html#ad85e4e6fe1ec77097d9108f38a8e57ae">setPrecursors</a>(pcs);</div>
<div class="line"><a name="l00393"></a><span class="lineno"> 393</span>  ms2_spec.<a class="code" href="classOpenMS_1_1MSSpectrum.html#aac8c50f318154e9e5b7b01aee806e17f">setRT</a>(scan->getRT() + rt_dist / 2.0);</div>
<div class="line"><a name="l00394"></a><span class="lineno"> 394</span>  ms2_spec.<a class="code" href="classOpenMS_1_1MSSpectrum.html#a24f081fe74bbbf773533ad19e96d8358">setMSLevel</a>(2);</div>
<div class="line"><a name="l00395"></a><span class="lineno"> 395</span>  <span class="comment">// link ms2 spectrum with features overlapping its precursor</span></div>
<div class="line"><a name="l00396"></a><span class="lineno"> 396</span>  <span class="comment">// Warning: this depends on the current order of features in the map</span></div>
<div class="line"><a name="l00397"></a><span class="lineno"> 397</span>  <span class="comment">// Attention: make sure to name ALL features that overlap, not only one!</span></div>
<div class="line"><a name="l00398"></a><span class="lineno"> 398</span>  ms2_spec.<a class="code" href="classOpenMS_1_1MetaInfoInterface.html#a135a7e818125b31198a8a65d9ac7a3d1">setMetaValue</a>(<span class="stringliteral">"parent_feature_ids"</span>, <a class="code" href="classOpenMS_1_1IntList.html#ac244d9a257aab1a24c1a622ccacab745">IntList::create</a>(<a class="code" href="classOpenMS_1_1String.html">String</a>(feature_index)));</div>
<div class="line"><a name="l00399"></a><span class="lineno"> 399</span>  ms2.<a class="code" href="classOpenMS_1_1MSExperiment.html#ae9dafbb74a1a4c430ceef03e25ab8e10">addSpectrum</a>(ms2_spec);</div>
<div class="line"><a name="l00400"></a><span class="lineno"> 400</span>  std::cout << <span class="stringliteral">" MS2 spectra generated at: "</span> << scan->getRT() << <span class="stringliteral">" x "</span> << p.<a class="code" href="classOpenMS_1_1Peak1D.html#a6fefef82aae7a3208e71648e763d5ea4">getMZ</a>() << <span class="stringliteral">"\n"</span>;</div>
<div class="line"><a name="l00401"></a><span class="lineno"> 401</span> </div>
<div class="line"><a name="l00402"></a><span class="lineno"> 402</span>  }</div>
<div class="line"><a name="l00403"></a><span class="lineno"> 403</span> <span class="preprocessor">#ifdef DEBUG_OPS</span></div>
<div class="line"><a name="l00404"></a><span class="lineno"> 404</span> <span class="preprocessor"></span> std::cout << solution_indices.size() << <span class="stringliteral">" out of "</span> << features.size()</div>
<div class="line"><a name="l00405"></a><span class="lineno"> 405</span>  << <span class="stringliteral">" precursors are in best solution.\n"</span>;</div>
<div class="line"><a name="l00406"></a><span class="lineno"> 406</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00407"></a><span class="lineno"> 407</span> <span class="preprocessor"></span> }</div>
<div class="line"><a name="l00408"></a><span class="lineno"> 408</span>  <span class="keywordflow">else</span> <span class="comment">// scan based selection (take the x highest signals for each spectrum)</span></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> <span class="preprocessor">#ifdef DEBUG_OPS</span></div>
<div class="line"><a name="l00411"></a><span class="lineno"> 411</span> <span class="preprocessor"></span> std::cout << <span class="stringliteral">"scan based precursor selection"</span> << std::endl;</div>
<div class="line"><a name="l00412"></a><span class="lineno"> 412</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00413"></a><span class="lineno"> 413</span> <span class="preprocessor"></span> <span class="comment">// if the highest signals for each scan shall be selected we don't need an ILP formulation</span></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>  <span class="comment">//cache the values for each feature</span></div>
<div class="line"><a name="l00416"></a><span class="lineno"> 416</span>  std::vector<DoubleList> feature_elution_bounds;</div>
<div class="line"><a name="l00417"></a><span class="lineno"> 417</span>  std::vector<DoubleList> elution_profile_intensities;</div>
<div class="line"><a name="l00418"></a><span class="lineno"> 418</span>  std::vector<DoubleList> isotope_intensities;</div>
<div class="line"><a name="l00419"></a><span class="lineno"> 419</span> </div>
<div class="line"><a name="l00420"></a><span class="lineno"> 420</span>  <span class="keywordtype">bool</span> meta_values_present = <span class="keyword">false</span>;</div>
<div class="line"><a name="l00421"></a><span class="lineno"> 421</span> </div>
<div class="line"><a name="l00422"></a><span class="lineno"> 422</span>  <span class="keywordflow">if</span> (!features.empty() &&</div>
<div class="line"><a name="l00423"></a><span class="lineno"> 423</span>  features[0].metaValueExists(<span class="stringliteral">"elution_profile_bounds"</span>) &&</div>
<div class="line"><a name="l00424"></a><span class="lineno"> 424</span>  features[0].metaValueExists(<span class="stringliteral">"elution_profile_intensities"</span>) &&</div>
<div class="line"><a name="l00425"></a><span class="lineno"> 425</span>  features[0].metaValueExists(<span class="stringliteral">"isotope_intensities"</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="keywordflow">for</span> (<a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> feat = 0; feat < features.size(); ++feat)</div>
<div class="line"><a name="l00428"></a><span class="lineno"> 428</span>  {</div>
<div class="line"><a name="l00429"></a><span class="lineno"> 429</span>  feature_elution_bounds.push_back(features[feat].getMetaValue(<span class="stringliteral">"elution_profile_bounds"</span>));</div>
<div class="line"><a name="l00430"></a><span class="lineno"> 430</span>  elution_profile_intensities.push_back(features[feat].getMetaValue(<span class="stringliteral">"elution_profile_intensities"</span>));</div>
<div class="line"><a name="l00431"></a><span class="lineno"> 431</span>  isotope_intensities.push_back(features[feat].getMetaValue(<span class="stringliteral">"isotope_intensities"</span>));</div>
<div class="line"><a name="l00432"></a><span class="lineno"> 432</span>  }</div>
<div class="line"><a name="l00433"></a><span class="lineno"> 433</span>  meta_values_present = <span class="keyword">true</span>;</div>
<div class="line"><a name="l00434"></a><span class="lineno"> 434</span>  }</div>
<div class="line"><a name="l00435"></a><span class="lineno"> 435</span> </div>
<div class="line"><a name="l00436"></a><span class="lineno"> 436</span>  <span class="comment">//for each feature cache for which scans it has to be considered</span></div>
<div class="line"><a name="l00437"></a><span class="lineno"> 437</span>  std::vector<std::vector<Size> > scan_features(experiment.<a class="code" href="classOpenMS_1_1MSExperiment.html#a9d78a687cf2a391198bb3cbc08bc06cb">size</a>());</div>
<div class="line"><a name="l00438"></a><span class="lineno"> 438</span> </div>
<div class="line"><a name="l00439"></a><span class="lineno"> 439</span>  <span class="keywordflow">for</span> (<a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> feat = 0; feat < features.size(); ++feat)</div>
<div class="line"><a name="l00440"></a><span class="lineno"> 440</span>  {</div>
<div class="line"><a name="l00441"></a><span class="lineno"> 441</span>  <span class="keywordflow">if</span> (charges_set.count(features[feat].getCharge()))</div>
<div class="line"><a name="l00442"></a><span class="lineno"> 442</span>  {</div>
<div class="line"><a name="l00443"></a><span class="lineno"> 443</span>  <a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> lower_rt = features[feat].getConvexHull().getBoundingBox().minX();</div>
<div class="line"><a name="l00444"></a><span class="lineno"> 444</span>  <a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> upper_rt = features[feat].getConvexHull().getBoundingBox().maxX();</div>
<div class="line"><a name="l00445"></a><span class="lineno"> 445</span>  <span class="keyword">typename</span> <a class="code" href="classOpenMS_1_1MSExperiment.html#aee0a1a68d0fab63eefa6e14010b4a7f0">MSExperiment<InputPeakType>::ConstIterator</a> it;</div>
<div class="line"><a name="l00446"></a><span class="lineno"> 446</span>  <span class="keywordflow">for</span> (it = experiment.<a class="code" href="classOpenMS_1_1MSExperiment.html#adb2f0f3e10754ca4dc7be677f15689a3">RTBegin</a>(lower_rt); it != experiment.<a class="code" href="classOpenMS_1_1MSExperiment.html#aecc2b466c4b527df363f5b1287a5e646">RTEnd</a>(upper_rt); ++it)</div>
<div class="line"><a name="l00447"></a><span class="lineno"> 447</span>  {</div>
<div class="line"><a name="l00448"></a><span class="lineno"> 448</span>  scan_features[it - experiment.<a class="code" href="classOpenMS_1_1MSExperiment.html#a2387033802383edbdc95f9bbb12a707e">begin</a>()].push_back(feat);</div>
<div class="line"><a name="l00449"></a><span class="lineno"> 449</span>  }</div>
<div class="line"><a name="l00450"></a><span class="lineno"> 450</span>  }</div>
<div class="line"><a name="l00451"></a><span class="lineno"> 451</span>  }</div>
<div class="line"><a name="l00452"></a><span class="lineno"> 452</span> </div>
<div class="line"><a name="l00453"></a><span class="lineno"> 453</span>  <span class="keywordtype">bool</span> dynamic_exclusion = <a class="code" href="classOpenMS_1_1DefaultParamHandler.html#a28c73e623c63a4fe3bfceb1ae8274f39">param_</a>.<a class="code" href="classOpenMS_1_1Param.html#abd0f5571675ec99e1aa647b39e2a8711">getValue</a>(<span class="stringliteral">"Exclusion:use_dynamic_exclusion"</span>) == <span class="stringliteral">"true"</span> ? <span class="keyword">true</span> : <span class="keyword">false</span>;</div>
<div class="line"><a name="l00454"></a><span class="lineno"> 454</span>  <span class="keyword">typedef</span> std::map<std::pair<DoubleReal, DoubleReal>, <a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a>, <a class="code" href="structOpenMS_1_1PairComparatorSecondElement.html">PairComparatorSecondElement<std::pair<DoubleReal, DoubleReal></a> > > ExclusionListType;</div>
<div class="line"><a name="l00455"></a><span class="lineno"> 455</span>  ExclusionListType exclusion_list;</div>
<div class="line"><a name="l00456"></a><span class="lineno"> 456</span>  Size exclusion_specs = (<a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a>)(floor((<a class="code" href="classdouble.html">DoubleReal</a>)<a class="code" href="classOpenMS_1_1DefaultParamHandler.html#a28c73e623c63a4fe3bfceb1ae8274f39">param_</a>.<a class="code" href="classOpenMS_1_1Param.html#abd0f5571675ec99e1aa647b39e2a8711">getValue</a>(<span class="stringliteral">"Exclusion:exclusion_time"</span>) / (<a class="code" href="group__Concept.html#gace75bfb1aba684e874dffee13738bd15">DoubleReal</a>) rt_dist));</div>
<div class="line"><a name="l00457"></a><span class="lineno"> 457</span>  <span class="keywordflow">if</span> (!dynamic_exclusion)</div>
<div class="line"><a name="l00458"></a><span class="lineno"> 458</span>  {</div>
<div class="line"><a name="l00459"></a><span class="lineno"> 459</span>  <span class="comment">//if the dynamic exclusion if not active we use the eclusion list to guarantee no two peaks within min_peak_distance are selected for single scan</span></div>
<div class="line"><a name="l00460"></a><span class="lineno"> 460</span>  exclusion_specs = 0;</div>
<div class="line"><a name="l00461"></a><span class="lineno"> 461</span>  }</div>
<div class="line"><a name="l00462"></a><span class="lineno"> 462</span> </div>
<div class="line"><a name="l00463"></a><span class="lineno"> 463</span>  <span class="comment">//cache bounding boxes of features and mass traces (mass trace bb are also widened for effective discovery of enclosing peaks in intervalls)</span></div>
<div class="line"><a name="l00464"></a><span class="lineno"> 464</span>  std::map<Size, typename OpenMS::DBoundingBox<2> > bounding_boxes_f;</div>
<div class="line"><a name="l00465"></a><span class="lineno"> 465</span>  std::map<std::pair<Size, Size>, <span class="keyword">typename</span> <a class="code" href="classOpenMS_1_1DBoundingBox.html">OpenMS::DBoundingBox<2></a> > bounding_boxes;</div>
<div class="line"><a name="l00466"></a><span class="lineno"> 466</span>  <span class="keywordflow">for</span> (Size feature_num = 0; feature_num < features.size(); ++feature_num)</div>
<div class="line"><a name="l00467"></a><span class="lineno"> 467</span>  {</div>
<div class="line"><a name="l00468"></a><span class="lineno"> 468</span>  <span class="keywordflow">if</span> (charges_set.count(features[feature_num].getCharge()))</div>
<div class="line"><a name="l00469"></a><span class="lineno"> 469</span>  {</div>
<div class="line"><a name="l00470"></a><span class="lineno"> 470</span>  bounding_boxes_f.insert(std::make_pair(feature_num, features[feature_num].getConvexHull().getBoundingBox()));</div>
<div class="line"><a name="l00471"></a><span class="lineno"> 471</span>  <span class="keyword">const</span> std::vector<ConvexHull2D> mass_traces = features[feature_num].getConvexHulls();</div>
<div class="line"><a name="l00472"></a><span class="lineno"> 472</span>  <span class="keywordflow">for</span> (Size mass_trace_num = 0; mass_trace_num < mass_traces.size(); ++mass_trace_num)</div>
<div class="line"><a name="l00473"></a><span class="lineno"> 473</span>  {</div>
<div class="line"><a name="l00474"></a><span class="lineno"> 474</span>  <span class="keyword">typename</span> <a class="code" href="classOpenMS_1_1DBoundingBox.html">OpenMS::DBoundingBox<2></a> tmp_bbox = mass_traces[mass_trace_num].getBoundingBox();</div>
<div class="line"><a name="l00475"></a><span class="lineno"> 475</span>  tmp_bbox.<a class="code" href="classOpenMS_1_1Internal_1_1DIntervalBase.html#a3faa28974e310e9743446f1c184496a6">setMinY</a>(tmp_bbox.<a class="code" href="classOpenMS_1_1Internal_1_1DIntervalBase.html#ad7c0ea2d4aceb7723ebcf33ae386634f">minY</a>() - window);</div>
<div class="line"><a name="l00476"></a><span class="lineno"> 476</span>  tmp_bbox.<a class="code" href="classOpenMS_1_1Internal_1_1DIntervalBase.html#a71b4e57973d9ad60975b774f4eeaef06">setMaxY</a>(tmp_bbox.<a class="code" href="classOpenMS_1_1Internal_1_1DIntervalBase.html#a1e8f4d69eddc6bd53fb8b1101e24f765">maxY</a>() + window);</div>
<div class="line"><a name="l00477"></a><span class="lineno"> 477</span>  bounding_boxes.insert(std::make_pair(std::make_pair(feature_num, mass_trace_num), tmp_bbox));</div>
<div class="line"><a name="l00478"></a><span class="lineno"> 478</span>  }</div>
<div class="line"><a name="l00479"></a><span class="lineno"> 479</span>  }</div>
<div class="line"><a name="l00480"></a><span class="lineno"> 480</span>  }</div>
<div class="line"><a name="l00481"></a><span class="lineno"> 481</span> </div>
<div class="line"><a name="l00482"></a><span class="lineno"> 482</span>  Size max_spec = (<a class="code" href="group__Concept.html#ga7cc214a236ad3bb6ad435bdcf5262a3f">Int</a>)<a class="code" href="classOpenMS_1_1DefaultParamHandler.html#a28c73e623c63a4fe3bfceb1ae8274f39">param_</a>.<a class="code" href="classOpenMS_1_1Param.html#abd0f5571675ec99e1aa647b39e2a8711">getValue</a>(<span class="stringliteral">"ms2_spectra_per_rt_bin"</span>);</div>
<div class="line"><a name="l00483"></a><span class="lineno"> 483</span>  <span class="comment">// get best x signals for each scan</span></div>
<div class="line"><a name="l00484"></a><span class="lineno"> 484</span>  <span class="keywordflow">for</span> (Size i = 0; i < experiment.<a class="code" href="classOpenMS_1_1MSExperiment.html#a9d78a687cf2a391198bb3cbc08bc06cb">size</a>(); ++i)</div>
<div class="line"><a name="l00485"></a><span class="lineno"> 485</span>  {</div>
<div class="line"><a name="l00486"></a><span class="lineno"> 486</span> <span class="preprocessor">#ifdef DEBUG_OPS</span></div>
<div class="line"><a name="l00487"></a><span class="lineno"> 487</span> <span class="preprocessor"></span> std::cout << <span class="stringliteral">"scan "</span> << experiment[i].getRT() << <span class="stringliteral">":"</span>;</div>
<div class="line"><a name="l00488"></a><span class="lineno"> 488</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00489"></a><span class="lineno"> 489</span> <span class="preprocessor"></span></div>
<div class="line"><a name="l00490"></a><span class="lineno"> 490</span>  <a class="code" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#ac76fec6f118131b8d46e5e2885ad7d50">updateExclusionList_</a>(exclusion_list);</div>
<div class="line"><a name="l00491"></a><span class="lineno"> 491</span>  <a class="code" href="classOpenMS_1_1MSSpectrum.html">MSSpectrum<InputPeakType></a> scan = experiment[i];</div>
<div class="line"><a name="l00492"></a><span class="lineno"> 492</span>  scan.<a class="code" href="classOpenMS_1_1MSSpectrum.html#a93f1c9f620403aec160d98e83e125502">sortByIntensity</a>(<span class="keyword">true</span>);</div>
<div class="line"><a name="l00493"></a><span class="lineno"> 493</span>  Size selected_peaks = 0, j = 0;</div>
<div class="line"><a name="l00494"></a><span class="lineno"> 494</span> </div>
<div class="line"><a name="l00495"></a><span class="lineno"> 495</span>  <span class="keywordflow">while</span> (selected_peaks < max_spec && j < scan.size())</div>
<div class="line"><a name="l00496"></a><span class="lineno"> 496</span>  {</div>
<div class="line"><a name="l00497"></a><span class="lineno"> 497</span>  <a class="code" href="classdouble.html">DoubleReal</a> peak_mz = scan[j].getMZ();</div>
<div class="line"><a name="l00498"></a><span class="lineno"> 498</span>  <a class="code" href="classdouble.html">DoubleReal</a> peak_rt = scan.<a class="code" href="classOpenMS_1_1MSSpectrum.html#a3da529bd3240fa0d7148484bbef0b9d7">getRT</a>();</div>
<div class="line"><a name="l00499"></a><span class="lineno"> 499</span> </div>
<div class="line"><a name="l00500"></a><span class="lineno"> 500</span>  ExclusionListType::iterator it_low = exclusion_list.lower_bound(std::make_pair(peak_mz, peak_mz));</div>
<div class="line"><a name="l00501"></a><span class="lineno"> 501</span>  <span class="keywordflow">if</span> (it_low != exclusion_list.end() && it_low->first.first <= peak_mz)</div>
<div class="line"><a name="l00502"></a><span class="lineno"> 502</span>  {</div>
<div class="line"><a name="l00503"></a><span class="lineno"> 503</span>  ++j;</div>
<div class="line"><a name="l00504"></a><span class="lineno"> 504</span>  <span class="keywordflow">continue</span>;</div>
<div class="line"><a name="l00505"></a><span class="lineno"> 505</span>  }</div>
<div class="line"><a name="l00506"></a><span class="lineno"> 506</span>  ++selected_peaks;</div>
<div class="line"><a name="l00507"></a><span class="lineno"> 507</span> </div>
<div class="line"><a name="l00508"></a><span class="lineno"> 508</span>  <span class="comment">//find all features (mass traces that are in the window around peak_mz)</span></div>
<div class="line"><a name="l00509"></a><span class="lineno"> 509</span>  <span class="keyword">typename</span> <a class="code" href="classOpenMS_1_1MSSpectrum.html">MSExperiment<InputPeakType>::SpectrumType</a> ms2_spec;</div>
<div class="line"><a name="l00510"></a><span class="lineno"> 510</span>  std::vector<Precursor> pcs;</div>
<div class="line"><a name="l00511"></a><span class="lineno"> 511</span>  std::set<std::pair<Size, Size> > selected_mt;</div>
<div class="line"><a name="l00512"></a><span class="lineno"> 512</span>  <a class="code" href="classOpenMS_1_1IntList.html">IntList</a> parent_feature_ids;</div>
<div class="line"><a name="l00513"></a><span class="lineno"> 513</span> </div>
<div class="line"><a name="l00514"></a><span class="lineno"> 514</span>  <a class="code" href="classdouble.html">DoubleReal</a> local_mz = peak_mz;</div>
<div class="line"><a name="l00515"></a><span class="lineno"> 515</span>  <span class="comment">//std::cerr<<"MZ pos: "<<local_mz<<std::endl;</span></div>
<div class="line"><a name="l00516"></a><span class="lineno"> 516</span>  <span class="keywordflow">for</span> (Size scan_feat_id = 0; scan_feat_id < scan_features[i].size(); ++scan_feat_id)</div>
<div class="line"><a name="l00517"></a><span class="lineno"> 517</span>  {</div>
<div class="line"><a name="l00518"></a><span class="lineno"> 518</span>  Size feature_num = scan_features[i][scan_feat_id];</div>
<div class="line"><a name="l00519"></a><span class="lineno"> 519</span>  <span class="keywordflow">if</span> (bounding_boxes_f[feature_num].encloses(peak_rt, local_mz))</div>
<div class="line"><a name="l00520"></a><span class="lineno"> 520</span>  {</div>
<div class="line"><a name="l00521"></a><span class="lineno"> 521</span>  <span class="comment">//find a mass trace enclosing the point</span></div>
<div class="line"><a name="l00522"></a><span class="lineno"> 522</span>  <a class="code" href="classdouble.html">DoubleReal</a> feature_intensity = 0;</div>
<div class="line"><a name="l00523"></a><span class="lineno"> 523</span>  <span class="keywordflow">for</span> (Size mass_trace_num = 0; mass_trace_num < features[feature_num].getConvexHulls().size(); ++mass_trace_num)</div>
<div class="line"><a name="l00524"></a><span class="lineno"> 524</span>  {</div>
<div class="line"><a name="l00525"></a><span class="lineno"> 525</span>  <span class="keywordflow">if</span> (bounding_boxes[std::make_pair(feature_num, mass_trace_num)].encloses(<a class="code" href="classOpenMS_1_1DPosition.html">DPosition<2></a>(peak_rt, local_mz)))</div>
<div class="line"><a name="l00526"></a><span class="lineno"> 526</span>  {</div>
<div class="line"><a name="l00527"></a><span class="lineno"> 527</span>  <a class="code" href="classdouble.html">DoubleReal</a> elu_factor = 1.0, iso_factor = 1.0;</div>
<div class="line"><a name="l00528"></a><span class="lineno"> 528</span>  <span class="comment">//get the intensity factor for the position in the elution profile</span></div>
<div class="line"><a name="l00529"></a><span class="lineno"> 529</span>  <span class="keywordflow">if</span> (meta_values_present)</div>
<div class="line"><a name="l00530"></a><span class="lineno"> 530</span>  {</div>
<div class="line"><a name="l00531"></a><span class="lineno"> 531</span>  <a class="code" href="classOpenMS_1_1DoubleList.html">DoubleList</a> xxx = elution_profile_intensities[feature_num];</div>
<div class="line"><a name="l00532"></a><span class="lineno"> 532</span>  <a class="code" href="classOpenMS_1_1DoubleList.html">DoubleList</a> yyy = feature_elution_bounds[feature_num];</div>
<div class="line"><a name="l00533"></a><span class="lineno"> 533</span> <span class="comment">// std::cout << "PEAKRT: " << peak_rt << std::endl;</span></div>
<div class="line"><a name="l00534"></a><span class="lineno"> 534</span> <span class="comment">// std::cout << "Max: " << yyy[3] << " vs. " << bounding_boxes_f[feature_num].maxX() << std::endl;</span></div>
<div class="line"><a name="l00535"></a><span class="lineno"> 535</span> <span class="comment">// std::cout << "Min: " << yyy[1] << " vs. " << bounding_boxes_f[feature_num].minX() << std::endl;</span></div>
<div class="line"><a name="l00536"></a><span class="lineno"> 536</span>  <a class="code" href="group__Conditions.html#ga45a74ea28109e7e1ed992fbd2ee97778">OPENMS_PRECONDITION</a>(i - yyy[0] < xxx.size(), <span class="stringliteral">"Tried to access invalid index for elution factor"</span>);</div>
<div class="line"><a name="l00537"></a><span class="lineno"> 537</span>  elu_factor = xxx[i - yyy[0]]; <span class="comment">// segfault here: "i-yyy[0]" yields invalid index</span></div>
<div class="line"><a name="l00538"></a><span class="lineno"> 538</span>  iso_factor = isotope_intensities[feature_num][mass_trace_num];</div>
<div class="line"><a name="l00539"></a><span class="lineno"> 539</span>  }</div>
<div class="line"><a name="l00540"></a><span class="lineno"> 540</span>  feature_intensity += features[feature_num].getIntensity() * iso_factor * elu_factor;</div>
<div class="line"><a name="l00541"></a><span class="lineno"> 541</span>  }</div>
<div class="line"><a name="l00542"></a><span class="lineno"> 542</span>  }</div>
<div class="line"><a name="l00543"></a><span class="lineno"> 543</span>  <a class="code" href="classOpenMS_1_1Precursor.html">Precursor</a> p;</div>
<div class="line"><a name="l00544"></a><span class="lineno"> 544</span>  p.<a class="code" href="classOpenMS_1_1Peak1D.html#ab875c2479c6cd2a0189e4684bfa4f97d">setIntensity</a>(feature_intensity);</div>
<div class="line"><a name="l00545"></a><span class="lineno"> 545</span>  p.<a class="code" href="classOpenMS_1_1Peak1D.html#ab9c4ce5d172a75d6de2c6f5f1847c9ab">setMZ</a>(features[feature_num].getMZ());</div>
<div class="line"><a name="l00546"></a><span class="lineno"> 546</span>  p.<a class="code" href="classOpenMS_1_1Precursor.html#aca40a4630f5f8a0568899b74a62dcd4b">setCharge</a>(features[feature_num].getCharge());</div>
<div class="line"><a name="l00547"></a><span class="lineno"> 547</span>  pcs.push_back(p);</div>
<div class="line"><a name="l00548"></a><span class="lineno"> 548</span>  parent_feature_ids.push_back((<a class="code" href="group__Concept.html#ga7cc214a236ad3bb6ad435bdcf5262a3f">Int</a>)feature_num);</div>
<div class="line"><a name="l00549"></a><span class="lineno"> 549</span>  }</div>
<div class="line"><a name="l00550"></a><span class="lineno"> 550</span>  }</div>
<div class="line"><a name="l00551"></a><span class="lineno"> 551</span> </div>
<div class="line"><a name="l00552"></a><span class="lineno"> 552</span>  <span class="keywordflow">if</span> (!pcs.empty())</div>
<div class="line"><a name="l00553"></a><span class="lineno"> 553</span>  {</div>
<div class="line"><a name="l00554"></a><span class="lineno"> 554</span>  <span class="comment">//std::cerr<<"scan "<<i<<" added spectrum for features: "<<parent_feature_ids<<std::endl;</span></div>
<div class="line"><a name="l00555"></a><span class="lineno"> 555</span>  ms2_spec.<a class="code" href="classOpenMS_1_1SpectrumSettings.html#ad85e4e6fe1ec77097d9108f38a8e57ae">setPrecursors</a>(pcs);</div>
<div class="line"><a name="l00556"></a><span class="lineno"> 556</span>  ms2_spec.<a class="code" href="classOpenMS_1_1MSSpectrum.html#a24f081fe74bbbf773533ad19e96d8358">setMSLevel</a>(2);</div>
<div class="line"><a name="l00557"></a><span class="lineno"> 557</span>  ms2_spec.<a class="code" href="classOpenMS_1_1MSSpectrum.html#aac8c50f318154e9e5b7b01aee806e17f">setRT</a>(experiment[i].getRT() + rt_dist / 2.0); <span class="comment">//(selected_peaks+1)*rt_dist/(max_spec+1) );</span></div>
<div class="line"><a name="l00558"></a><span class="lineno"> 558</span>  ms2_spec.<a class="code" href="classOpenMS_1_1MetaInfoInterface.html#a135a7e818125b31198a8a65d9ac7a3d1">setMetaValue</a>(<span class="stringliteral">"parent_feature_ids"</span>, parent_feature_ids);</div>
<div class="line"><a name="l00559"></a><span class="lineno"> 559</span>  ms2.<a class="code" href="classOpenMS_1_1MSExperiment.html#ae9dafbb74a1a4c430ceef03e25ab8e10">addSpectrum</a>(ms2_spec);</div>
<div class="line"><a name="l00560"></a><span class="lineno"> 560</span>  }</div>
<div class="line"><a name="l00561"></a><span class="lineno"> 561</span> </div>
<div class="line"><a name="l00562"></a><span class="lineno"> 562</span>  <span class="comment">//add m/z window to exclusion list</span></div>
<div class="line"><a name="l00563"></a><span class="lineno"> 563</span>  exclusion_list.insert(std::make_pair(std::make_pair(peak_mz - excl_window, peak_mz + excl_window), exclusion_specs + 1));</div>
<div class="line"><a name="l00564"></a><span class="lineno"> 564</span> </div>
<div class="line"><a name="l00565"></a><span class="lineno"> 565</span>  ++j;</div>
<div class="line"><a name="l00566"></a><span class="lineno"> 566</span>  }</div>
<div class="line"><a name="l00567"></a><span class="lineno"> 567</span>  }</div>
<div class="line"><a name="l00568"></a><span class="lineno"> 568</span>  }</div>
<div class="line"><a name="l00569"></a><span class="lineno"> 569</span>  }</div>
<div class="line"><a name="l00570"></a><span class="lineno"> 570</span> </div>
<div class="line"><a name="l00571"></a><span class="lineno"> 571</span>  <span class="keyword">template</span> <<span class="keyword">typename</span> InputPeakType></div>
<div class="line"><a name="l00572"></a><span class="lineno"><a class="line" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#acca2e18b840cc4b903f2958efca7aa1a"> 572</a></span>  <span class="keywordtype">void</span> <a class="code" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#acca2e18b840cc4b903f2958efca7aa1a">OfflinePrecursorIonSelection::checkMassRanges_</a>(std::vector<std::vector<std::pair<Size, Size> > >& mass_ranges,</div>
<div class="line"><a name="l00573"></a><span class="lineno"> 573</span>  <span class="keyword">const</span> <a class="code" href="classOpenMS_1_1MSExperiment.html">MSExperiment<InputPeakType></a>& experiment)</div>
<div class="line"><a name="l00574"></a><span class="lineno"> 574</span>  {</div>
<div class="line"><a name="l00575"></a><span class="lineno"> 575</span>  std::vector<std::vector<std::pair<Size, Size> > > checked_mass_ranges;</div>
<div class="line"><a name="l00576"></a><span class="lineno"> 576</span>  <a class="code" href="classdouble.html">DoubleReal</a> min_peak_distance = <a class="code" href="classOpenMS_1_1DefaultParamHandler.html#a28c73e623c63a4fe3bfceb1ae8274f39">param_</a>.<a class="code" href="classOpenMS_1_1Param.html#abd0f5571675ec99e1aa647b39e2a8711">getValue</a>(<span class="stringliteral">"min_peak_distance"</span>);</div>
<div class="line"><a name="l00577"></a><span class="lineno"> 577</span>  checked_mass_ranges.reserve(mass_ranges.size());</div>
<div class="line"><a name="l00578"></a><span class="lineno"> 578</span>  <span class="keywordflow">for</span> (<a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> f = 0; f < mass_ranges.size(); ++f)</div>
<div class="line"><a name="l00579"></a><span class="lineno"> 579</span>  {</div>
<div class="line"><a name="l00580"></a><span class="lineno"> 580</span>  std::vector<std::pair<Size, Size> > checked_mass_ranges_f;</div>
<div class="line"><a name="l00581"></a><span class="lineno"> 581</span>  <span class="keywordflow">for</span> (<a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> s_idx = 0; s_idx < mass_ranges[f].size(); s_idx += 2)</div>
<div class="line"><a name="l00582"></a><span class="lineno"> 582</span>  {</div>
<div class="line"><a name="l00583"></a><span class="lineno"> 583</span>  <a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> s = mass_ranges[f][s_idx].first;</div>
<div class="line"><a name="l00584"></a><span class="lineno"> 584</span>  <span class="keywordtype">bool</span> overlapping_features = <span class="keyword">false</span>;</div>
<div class="line"><a name="l00586"></a><span class="lineno"> 586</span>  <span class="comment">// check if other features overlap with this feature in the current scan</span></div>
<div class="line"><a name="l00588"></a><span class="lineno"> 588</span> <span class="comment"></span> <span class="keyword">const</span> InputPeakType& peak_left_border = experiment[s][mass_ranges[f][s_idx].second];</div>
<div class="line"><a name="l00589"></a><span class="lineno"> 589</span>  <span class="keyword">const</span> InputPeakType& peak_right_border = experiment[s][mass_ranges[f][s_idx + 1].second];</div>
<div class="line"><a name="l00590"></a><span class="lineno"> 590</span>  <span class="keywordflow">for</span> (<a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> fmr = 0; fmr < mass_ranges.size(); ++fmr)</div>
<div class="line"><a name="l00591"></a><span class="lineno"> 591</span>  {</div>
<div class="line"><a name="l00592"></a><span class="lineno"> 592</span>  <span class="keywordflow">if</span> (fmr == f)</div>
<div class="line"><a name="l00593"></a><span class="lineno"> 593</span>  <span class="keywordflow">continue</span>;</div>
<div class="line"><a name="l00594"></a><span class="lineno"> 594</span>  <span class="keywordflow">for</span> (<a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> mr = 0; mr < mass_ranges[fmr].size(); mr += 2)</div>
<div class="line"><a name="l00595"></a><span class="lineno"> 595</span>  {</div>
<div class="line"><a name="l00596"></a><span class="lineno"> 596</span>  <span class="keywordflow">if</span> (mass_ranges[fmr][mr].first == s) <span class="comment">// same spectrum</span></div>
<div class="line"><a name="l00597"></a><span class="lineno"> 597</span>  {</div>
<div class="line"><a name="l00598"></a><span class="lineno"> 598</span>  <span class="keyword">const</span> InputPeakType& tmp_peak_left = experiment[s][mass_ranges[fmr][mr].second];</div>
<div class="line"><a name="l00599"></a><span class="lineno"> 599</span>  <span class="keyword">const</span> InputPeakType& tmp_peak_right = experiment[s][mass_ranges[fmr][mr + 1].second];</div>
<div class="line"><a name="l00600"></a><span class="lineno"> 600</span> <span class="preprocessor">#ifdef DEBUG_OPS</span></div>
<div class="line"><a name="l00601"></a><span class="lineno"> 601</span> <span class="preprocessor"></span> std::cout << tmp_peak_left.getMZ() << <span class="stringliteral">" < "</span></div>
<div class="line"><a name="l00602"></a><span class="lineno"> 602</span>  << peak_left_border.getMZ() - min_peak_distance << <span class="stringliteral">" && "</span></div>
<div class="line"><a name="l00603"></a><span class="lineno"> 603</span>  << tmp_peak_right.getMZ() << <span class="stringliteral">" < "</span></div>
<div class="line"><a name="l00604"></a><span class="lineno"> 604</span>  << peak_left_border.getMZ() - min_peak_distance << <span class="stringliteral">" ? "</span></div>
<div class="line"><a name="l00605"></a><span class="lineno"> 605</span>  << (tmp_peak_left.getMZ() < peak_left_border.getMZ() - min_peak_distance &&</div>
<div class="line"><a name="l00606"></a><span class="lineno"> 606</span>  tmp_peak_right.getMZ() < peak_left_border.getMZ() - min_peak_distance)</div>
<div class="line"><a name="l00607"></a><span class="lineno"> 607</span>  << <span class="stringliteral">" || "</span></div>
<div class="line"><a name="l00608"></a><span class="lineno"> 608</span>  << tmp_peak_left.getMZ() << <span class="stringliteral">" > "</span></div>
<div class="line"><a name="l00609"></a><span class="lineno"> 609</span>  << peak_right_border.getMZ() + min_peak_distance << <span class="stringliteral">" && "</span></div>
<div class="line"><a name="l00610"></a><span class="lineno"> 610</span>  << tmp_peak_right.getMZ() << <span class="stringliteral">" > "</span></div>
<div class="line"><a name="l00611"></a><span class="lineno"> 611</span>  << peak_right_border.getMZ() + min_peak_distance << <span class="stringliteral">" ? "</span></div>
<div class="line"><a name="l00612"></a><span class="lineno"> 612</span>  << (tmp_peak_left.getMZ() > peak_right_border.getMZ() + min_peak_distance &&</div>
<div class="line"><a name="l00613"></a><span class="lineno"> 613</span>  tmp_peak_right.getMZ() > peak_right_border.getMZ() + min_peak_distance)</div>
<div class="line"><a name="l00614"></a><span class="lineno"> 614</span>  << std::endl;</div>
<div class="line"><a name="l00615"></a><span class="lineno"> 615</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00616"></a><span class="lineno"> 616</span> <span class="preprocessor"></span> <span class="comment">// all other features have to be either completely left or</span></div>
<div class="line"><a name="l00617"></a><span class="lineno"> 617</span>  <span class="comment">// right of the current feature</span></div>
<div class="line"><a name="l00618"></a><span class="lineno"> 618</span>  <span class="keywordflow">if</span> (!((tmp_peak_left.getMZ() < peak_left_border.getMZ() - min_peak_distance &&</div>
<div class="line"><a name="l00619"></a><span class="lineno"> 619</span>  tmp_peak_right.getMZ() < peak_left_border.getMZ() - min_peak_distance) ||</div>
<div class="line"><a name="l00620"></a><span class="lineno"> 620</span>  (tmp_peak_left.getMZ() > peak_right_border.getMZ() + min_peak_distance &&</div>
<div class="line"><a name="l00621"></a><span class="lineno"> 621</span>  tmp_peak_right.getMZ() > peak_right_border.getMZ() + min_peak_distance)))</div>
<div class="line"><a name="l00622"></a><span class="lineno"> 622</span>  {</div>
<div class="line"><a name="l00623"></a><span class="lineno"> 623</span> <span class="preprocessor">#ifdef DEBUG_OPS</span></div>
<div class="line"><a name="l00624"></a><span class="lineno"> 624</span> <span class="preprocessor"></span> std::cout << <span class="stringliteral">"found overlapping peak"</span> << std::endl;</div>
<div class="line"><a name="l00625"></a><span class="lineno"> 625</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00626"></a><span class="lineno"> 626</span> <span class="preprocessor"></span> overlapping_features = <span class="keyword">true</span>;</div>
<div class="line"><a name="l00627"></a><span class="lineno"> 627</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00628"></a><span class="lineno"> 628</span>  }</div>
<div class="line"><a name="l00629"></a><span class="lineno"> 629</span>  }</div>
<div class="line"><a name="l00630"></a><span class="lineno"> 630</span>  }</div>
<div class="line"><a name="l00631"></a><span class="lineno"> 631</span>  }</div>
<div class="line"><a name="l00632"></a><span class="lineno"> 632</span>  <span class="keywordflow">if</span> (!overlapping_features)</div>
<div class="line"><a name="l00633"></a><span class="lineno"> 633</span>  {</div>
<div class="line"><a name="l00634"></a><span class="lineno"> 634</span> <span class="preprocessor">#ifdef DEBUG_OPS</span></div>
<div class="line"><a name="l00635"></a><span class="lineno"> 635</span> <span class="preprocessor"></span> std::cout << <span class="stringliteral">"feature in spec ok"</span> << mass_ranges[f][s_idx].second << <span class="stringliteral">" in spec "</span></div>
<div class="line"><a name="l00636"></a><span class="lineno"> 636</span>  << mass_ranges[f][s_idx].first << std::endl;</div>
<div class="line"><a name="l00637"></a><span class="lineno"> 637</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00638"></a><span class="lineno"> 638</span> <span class="preprocessor"></span> checked_mass_ranges_f.insert(checked_mass_ranges_f.end(),</div>
<div class="line"><a name="l00639"></a><span class="lineno"> 639</span>  mass_ranges[f].begin() + s_idx,</div>
<div class="line"><a name="l00640"></a><span class="lineno"> 640</span>  mass_ranges[f].begin() + s_idx + 2);</div>
<div class="line"><a name="l00641"></a><span class="lineno"> 641</span>  }</div>
<div class="line"><a name="l00642"></a><span class="lineno"> 642</span>  }</div>
<div class="line"><a name="l00643"></a><span class="lineno"> 643</span>  checked_mass_ranges.push_back(checked_mass_ranges_f);</div>
<div class="line"><a name="l00644"></a><span class="lineno"> 644</span>  }</div>
<div class="line"><a name="l00645"></a><span class="lineno"> 645</span>  mass_ranges.swap(checked_mass_ranges);</div>
<div class="line"><a name="l00646"></a><span class="lineno"> 646</span>  }</div>
<div class="line"><a name="l00647"></a><span class="lineno"> 647</span> </div>
<div class="line"><a name="l00648"></a><span class="lineno"> 648</span>  <span class="keyword">template</span> <<span class="keyword">typename</span> T></div>
<div class="line"><a name="l00649"></a><span class="lineno"><a class="line" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#ac76fec6f118131b8d46e5e2885ad7d50"> 649</a></span>  <span class="keywordtype">void</span> <a class="code" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#ac76fec6f118131b8d46e5e2885ad7d50">OfflinePrecursorIonSelection::updateExclusionList_</a>(std::vector<std::pair<T, Size> >& exclusion_list)</div>
<div class="line"><a name="l00650"></a><span class="lineno"> 650</span>  {</div>
<div class="line"><a name="l00651"></a><span class="lineno"> 651</span>  <span class="keywordflow">for</span> (<a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a> i = 0; i < exclusion_list.size(); ++i)</div>
<div class="line"><a name="l00652"></a><span class="lineno"> 652</span>  {</div>
<div class="line"><a name="l00653"></a><span class="lineno"> 653</span>  <span class="keywordflow">if</span> (exclusion_list[i].second > 0)</div>
<div class="line"><a name="l00654"></a><span class="lineno"> 654</span>  --exclusion_list[i].second;</div>
<div class="line"><a name="l00655"></a><span class="lineno"> 655</span>  }</div>
<div class="line"><a name="l00656"></a><span class="lineno"> 656</span>  sort(exclusion_list.begin(), exclusion_list.end(), <a class="code" href="structOpenMS_1_1PairComparatorSecondElementMore.html">PairComparatorSecondElementMore<std::pair<T, Size></a> >());</div>
<div class="line"><a name="l00657"></a><span class="lineno"> 657</span>  <span class="keyword">typename</span> std::vector<std::pair<T, Size> >::iterator iter = exclusion_list.begin();</div>
<div class="line"><a name="l00658"></a><span class="lineno"> 658</span>  <span class="keywordflow">while</span> (iter != exclusion_list.end() && iter->second != 0)</div>
<div class="line"><a name="l00659"></a><span class="lineno"> 659</span>  ++iter;</div>
<div class="line"><a name="l00660"></a><span class="lineno"> 660</span>  exclusion_list.erase(iter, exclusion_list.end());</div>
<div class="line"><a name="l00661"></a><span class="lineno"> 661</span>  }</div>
<div class="line"><a name="l00662"></a><span class="lineno"> 662</span> </div>
<div class="line"><a name="l00663"></a><span class="lineno"><a class="line" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#a0a78199f092a6aa84227a6decad602ca"> 663</a></span>  <span class="keyword">inline</span> <span class="keywordtype">void</span> <a class="code" href="classOpenMS_1_1OfflinePrecursorIonSelection.html#ac76fec6f118131b8d46e5e2885ad7d50">OfflinePrecursorIonSelection::updateExclusionList_</a>(std::map<std::pair<DoubleReal, DoubleReal>, <a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a>, <a class="code" href="structOpenMS_1_1PairComparatorSecondElement.html">PairComparatorSecondElement</a><std::pair<DoubleReal, DoubleReal> > >& exclusion_list)</div>
<div class="line"><a name="l00664"></a><span class="lineno"> 664</span>  {</div>
<div class="line"><a name="l00665"></a><span class="lineno"> 665</span>  std::map<std::pair<DoubleReal, DoubleReal>, <a class="code" href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">Size</a>, <a class="code" href="structOpenMS_1_1PairComparatorSecondElement.html">PairComparatorSecondElement<std::pair<DoubleReal, DoubleReal></a> > >::iterator it;</div>
<div class="line"><a name="l00666"></a><span class="lineno"> 666</span> </div>
<div class="line"><a name="l00667"></a><span class="lineno"> 667</span>  it = exclusion_list.begin();</div>
<div class="line"><a name="l00668"></a><span class="lineno"> 668</span> </div>
<div class="line"><a name="l00669"></a><span class="lineno"> 669</span>  <span class="keywordflow">while</span> (it != exclusion_list.end())</div>
<div class="line"><a name="l00670"></a><span class="lineno"> 670</span>  {</div>
<div class="line"><a name="l00671"></a><span class="lineno"> 671</span>  <span class="keywordflow">if</span> ((it->second--) == 1)</div>
<div class="line"><a name="l00672"></a><span class="lineno"> 672</span>  {</div>
<div class="line"><a name="l00673"></a><span class="lineno"> 673</span>  exclusion_list.erase(it++);</div>
<div class="line"><a name="l00674"></a><span class="lineno"> 674</span>  }</div>
<div class="line"><a name="l00675"></a><span class="lineno"> 675</span>  <span class="keywordflow">else</span></div>
<div class="line"><a name="l00676"></a><span class="lineno"> 676</span>  {</div>
<div class="line"><a name="l00677"></a><span class="lineno"> 677</span>  ++it;</div>
<div class="line"><a name="l00678"></a><span class="lineno"> 678</span>  }</div>
<div class="line"><a name="l00679"></a><span class="lineno"> 679</span>  }</div>
<div class="line"><a name="l00680"></a><span class="lineno"> 680</span>  }</div>
<div class="line"><a name="l00681"></a><span class="lineno"> 681</span> </div>
<div class="line"><a name="l00682"></a><span class="lineno"> 682</span> }</div>
<div class="line"><a name="l00683"></a><span class="lineno"> 683</span> </div>
<div class="line"><a name="l00684"></a><span class="lineno"> 684</span> <span class="preprocessor">#endif // OPENMS_ANALYSIS_ID_OFFLINEPRECURSORIONSELECTION_H</span></div>
<div class="ttc" id="classOpenMS_1_1Internal_1_1DIntervalBase_html_a1e8f4d69eddc6bd53fb8b1101e24f765"><div class="ttname"><a href="classOpenMS_1_1Internal_1_1DIntervalBase.html#a1e8f4d69eddc6bd53fb8b1101e24f765">OpenMS::Internal::DIntervalBase::maxY</a></div><div class="ttdeci">CoordinateType maxY() const </div><div class="ttdoc">Accessor for max_ coordinate maximum. </div><div class="ttdef"><b>Definition:</b> DIntervalBase.h:258</div></div>
<div class="ttc" id="classOpenMS_1_1MetaInfoInterface_html_a135a7e818125b31198a8a65d9ac7a3d1"><div class="ttname"><a href="classOpenMS_1_1MetaInfoInterface.html#a135a7e818125b31198a8a65d9ac7a3d1">OpenMS::MetaInfoInterface::setMetaValue</a></div><div class="ttdeci">void setMetaValue(const String &name, const DataValue &value)</div><div class="ttdoc">sets the DataValue corresponding to a name </div></div>
<div class="ttc" id="FeatureXMLFile_8h_html"><div class="ttname"><a href="FeatureXMLFile_8h.html">FeatureXMLFile.h</a></div></div>
<div class="ttc" id="classOpenMS_1_1String_html"><div class="ttname"><a href="classOpenMS_1_1String.html">OpenMS::String</a></div><div class="ttdoc">A more convenient string class. </div><div class="ttdef"><b>Definition:</b> String.h:56</div></div>
<div class="ttc" id="classOpenMS_1_1Precursor_html"><div class="ttname"><a href="classOpenMS_1_1Precursor.html">OpenMS::Precursor</a></div><div class="ttdoc">Precursor meta information. </div><div class="ttdef"><b>Definition:</b> Precursor.h:56</div></div>
<div class="ttc" id="classOpenMS_1_1Internal_1_1DIntervalBase_html_a3faa28974e310e9743446f1c184496a6"><div class="ttname"><a href="classOpenMS_1_1Internal_1_1DIntervalBase.html#a3faa28974e310e9743446f1c184496a6">OpenMS::Internal::DIntervalBase::setMinY</a></div><div class="ttdeci">void setMinY(CoordinateType const c)</div><div class="ttdoc">Mutator for max_ coordinate of the smaller point. </div><div class="ttdef"><b>Definition:</b> DIntervalBase.h:271</div></div>
<div class="ttc" id="classOpenMS_1_1MSExperiment_html_a9d78a687cf2a391198bb3cbc08bc06cb"><div class="ttname"><a href="classOpenMS_1_1MSExperiment.html#a9d78a687cf2a391198bb3cbc08bc06cb">OpenMS::MSExperiment::size</a></div><div class="ttdeci">Size size() const </div><div class="ttdef"><b>Definition:</b> MSExperiment.h:117</div></div>
<div class="ttc" id="classOpenMS_1_1DefaultParamHandler_html_a28c73e623c63a4fe3bfceb1ae8274f39"><div class="ttname"><a href="classOpenMS_1_1DefaultParamHandler.html#a28c73e623c63a4fe3bfceb1ae8274f39">OpenMS::DefaultParamHandler::param_</a></div><div class="ttdeci">Param param_</div><div class="ttdoc">Container for current parameters. </div><div class="ttdef"><b>Definition:</b> DefaultParamHandler.h:148</div></div>
<div class="ttc" id="classOpenMS_1_1LPWrapper_html_ab659f8690ca779b8c3d73fe524a39a00"><div class="ttname"><a href="classOpenMS_1_1LPWrapper.html#ab659f8690ca779b8c3d73fe524a39a00">OpenMS::LPWrapper::SOLVER</a></div><div class="ttdeci">SOLVER</div><div class="ttdef"><b>Definition:</b> LPWrapper.h:124</div></div>
<div class="ttc" id="structOpenMS_1_1PSLPFormulation_1_1IndexLess_html"><div class="ttname"><a href="structOpenMS_1_1PSLPFormulation_1_1IndexLess.html">OpenMS::PSLPFormulation::IndexLess</a></div><div class="ttdef"><b>Definition:</b> PSLPFormulation.h:138</div></div>
<div class="ttc" id="classOpenMS_1_1FeatureMap_html"><div class="ttname"><a href="classOpenMS_1_1FeatureMap.html">OpenMS::FeatureMap<></a></div></div>
<div class="ttc" id="group__Conditions_html_ga45a74ea28109e7e1ed992fbd2ee97778"><div class="ttname"><a href="group__Conditions.html#ga45a74ea28109e7e1ed992fbd2ee97778">OPENMS_PRECONDITION</a></div><div class="ttdeci">#define OPENMS_PRECONDITION(condition, message)</div><div class="ttdoc">Precondition macro. </div><div class="ttdef"><b>Definition:</b> Macros.h:107</div></div>
<div class="ttc" id="classOpenMS_1_1Internal_1_1DIntervalBase_html_a71b4e57973d9ad60975b774f4eeaef06"><div class="ttname"><a href="classOpenMS_1_1Internal_1_1DIntervalBase.html#a71b4e57973d9ad60975b774f4eeaef06">OpenMS::Internal::DIntervalBase::setMaxY</a></div><div class="ttdeci">void setMaxY(CoordinateType const c)</div><div class="ttdoc">Mutator for max_ coordinate of the larger point. </div><div class="ttdef"><b>Definition:</b> DIntervalBase.h:285</div></div>
<div class="ttc" id="classOpenMS_1_1OfflinePrecursorIonSelection_html_ac945a0966726a3a644d194b81c28d960"><div class="ttname"><a href="classOpenMS_1_1OfflinePrecursorIonSelection.html#ac945a0966726a3a644d194b81c28d960">OpenMS::OfflinePrecursorIonSelection::getLPSolver</a></div><div class="ttdeci">LPWrapper::SOLVER getLPSolver()</div><div class="ttdef"><b>Definition:</b> OfflinePrecursorIonSelection.h:106</div></div>
<div class="ttc" id="classOpenMS_1_1Peak1D_html_a6fefef82aae7a3208e71648e763d5ea4"><div class="ttname"><a href="classOpenMS_1_1Peak1D.html#a6fefef82aae7a3208e71648e763d5ea4">OpenMS::Peak1D::getMZ</a></div><div class="ttdeci">CoordinateType getMZ() const </div><div class="ttdoc">Non-mutable access to m/z. </div><div class="ttdef"><b>Definition:</b> Peak1D.h:108</div></div>
<div class="ttc" id="classOpenMS_1_1MSSpectrum_html_a4895afe86c2aa0747ae0411c58c82114"><div class="ttname"><a href="classOpenMS_1_1MSSpectrum.html#a4895afe86c2aa0747ae0411c58c82114">OpenMS::MSSpectrum::ConstIterator</a></div><div class="ttdeci">ContainerType::const_iterator ConstIterator</div><div class="ttdoc">Non-mutable iterator. </div><div class="ttdef"><b>Definition:</b> MSSpectrum.h:125</div></div>
<div class="ttc" id="classOpenMS_1_1MSExperiment_html_a2387033802383edbdc95f9bbb12a707e"><div class="ttname"><a href="classOpenMS_1_1MSExperiment.html#a2387033802383edbdc95f9bbb12a707e">OpenMS::MSExperiment::begin</a></div><div class="ttdeci">Iterator begin()</div><div class="ttdef"><b>Definition:</b> MSExperiment.h:147</div></div>
<div class="ttc" id="namespaceOpenMS_html_aa36b2ed6c06f4252dec15cd4f6d72635"><div class="ttname"><a href="namespaceOpenMS.html#aa36b2ed6c06f4252dec15cd4f6d72635">OpenMS::enclosesBoundingBox</a></div><div class="ttdeci">bool enclosesBoundingBox(const Feature &f, typename MSExperiment< InputPeakType >::CoordinateType rt, typename MSExperiment< InputPeakType >::CoordinateType mz)</div><div class="ttdef"><b>Definition:</b> OfflinePrecursorIonSelection.h:138</div></div>
<div class="ttc" id="structOpenMS_1_1PSLPFormulation_1_1IndexTriple_html"><div class="ttname"><a href="structOpenMS_1_1PSLPFormulation_1_1IndexTriple.html">OpenMS::PSLPFormulation::IndexTriple</a></div><div class="ttdoc">Struct that holds the indices of the precursors in the feature map and the ilp formulation. </div><div class="ttdef"><b>Definition:</b> PSLPFormulation.h:66</div></div>
<div class="ttc" id="classdouble_html"><div class="ttname"><a href="classdouble.html">double</a></div></div>
<div class="ttc" id="classOpenMS_1_1Peak1D_html_ab9c4ce5d172a75d6de2c6f5f1847c9ab"><div class="ttname"><a href="classOpenMS_1_1Peak1D.html#ab9c4ce5d172a75d6de2c6f5f1847c9ab">OpenMS::Peak1D::setMZ</a></div><div class="ttdeci">void setMZ(CoordinateType mz)</div><div class="ttdoc">Mutable access to m/z. </div><div class="ttdef"><b>Definition:</b> Peak1D.h:114</div></div>
<div class="ttc" id="classOpenMS_1_1Peak1D_html_ab875c2479c6cd2a0189e4684bfa4f97d"><div class="ttname"><a href="classOpenMS_1_1Peak1D.html#ab875c2479c6cd2a0189e4684bfa4f97d">OpenMS::Peak1D::setIntensity</a></div><div class="ttdeci">void setIntensity(IntensityType intensity)</div><div class="ttdoc">Mutable access to the data point intensity (height) </div><div class="ttdef"><b>Definition:</b> Peak1D.h:105</div></div>
<div class="ttc" id="classOpenMS_1_1MSExperiment_html_aecc2b466c4b527df363f5b1287a5e646"><div class="ttname"><a href="classOpenMS_1_1MSExperiment.html#aecc2b466c4b527df363f5b1287a5e646">OpenMS::MSExperiment::RTEnd</a></div><div class="ttdeci">ConstIterator RTEnd(CoordinateType rt) const </div><div class="ttdoc">Fast search for spectrum range end (returns the past-the-end iterator) </div><div class="ttdef"><b>Definition:</b> MSExperiment.h:363</div></div>
<div class="ttc" id="structOpenMS_1_1PairComparatorSecondElementMore_html"><div class="ttname"><a href="structOpenMS_1_1PairComparatorSecondElementMore.html">OpenMS::PairComparatorSecondElementMore</a></div><div class="ttdoc">Class for comparison of std::pair using second ONLY e.g. for use with std::sort. </div><div class="ttdef"><b>Definition:</b> ComparatorUtils.h:368</div></div>
<div class="ttc" id="classOpenMS_1_1MSExperiment_html_ab45dae688fc5d8983727abffa4389003"><div class="ttname"><a href="classOpenMS_1_1MSExperiment.html#ab45dae688fc5d8983727abffa4389003">OpenMS::MSExperiment::end</a></div><div class="ttdeci">Iterator end()</div><div class="ttdef"><b>Definition:</b> MSExperiment.h:157</div></div>
<div class="ttc" id="classOpenMS_1_1MSExperiment_html_adb2f0f3e10754ca4dc7be677f15689a3"><div class="ttname"><a href="classOpenMS_1_1MSExperiment.html#adb2f0f3e10754ca4dc7be677f15689a3">OpenMS::MSExperiment::RTBegin</a></div><div class="ttdeci">ConstIterator RTBegin(CoordinateType rt) const </div><div class="ttdoc">Fast search for spectrum range begin. </div><div class="ttdef"><b>Definition:</b> MSExperiment.h:349</div></div>
<div class="ttc" id="MSExperiment_8h_html"><div class="ttname"><a href="MSExperiment_8h.html">MSExperiment.h</a></div></div>
<div class="ttc" id="classOpenMS_1_1Internal_1_1DIntervalBase_html_ad7c0ea2d4aceb7723ebcf33ae386634f"><div class="ttname"><a href="classOpenMS_1_1Internal_1_1DIntervalBase.html#ad7c0ea2d4aceb7723ebcf33ae386634f">OpenMS::Internal::DIntervalBase::minY</a></div><div class="ttdeci">CoordinateType minY() const </div><div class="ttdoc">Accessor for max_ coordinate minimum. </div><div class="ttdef"><b>Definition:</b> DIntervalBase.h:246</div></div>
<div class="ttc" id="classOpenMS_1_1MSSpectrum_html"><div class="ttname"><a href="classOpenMS_1_1MSSpectrum.html">OpenMS::MSSpectrum< PeakType ></a></div></div>
<div class="ttc" id="classOpenMS_1_1Param_html_abd0f5571675ec99e1aa647b39e2a8711"><div class="ttname"><a href="classOpenMS_1_1Param.html#abd0f5571675ec99e1aa647b39e2a8711">OpenMS::Param::getValue</a></div><div class="ttdeci">const DataValue & getValue(const String &key) const </div><div class="ttdoc">Returns a value of a parameter. </div></div>
<div class="ttc" id="classOpenMS_1_1OfflinePrecursorIonSelection_html_ac76fec6f118131b8d46e5e2885ad7d50"><div class="ttname"><a href="classOpenMS_1_1OfflinePrecursorIonSelection.html#ac76fec6f118131b8d46e5e2885ad7d50">OpenMS::OfflinePrecursorIonSelection::updateExclusionList_</a></div><div class="ttdeci">void updateExclusionList_(std::vector< std::pair< T, Size > > &exclusion_list)</div><div class="ttdef"><b>Definition:</b> OfflinePrecursorIonSelection.h:649</div></div>
<div class="ttc" id="classOpenMS_1_1DPosition_html"><div class="ttname"><a href="classOpenMS_1_1DPosition.html">OpenMS::DPosition< 2 ></a></div></div>
<div class="ttc" id="classOpenMS_1_1Feature_html_ad99bb8ba2103d6c52ee10ae6369b46e7"><div class="ttname"><a href="classOpenMS_1_1Feature.html#ad99bb8ba2103d6c52ee10ae6369b46e7">OpenMS::Feature::getConvexHulls</a></div><div class="ttdeci">const std::vector< ConvexHull2D > & getConvexHulls() const </div><div class="ttdoc">Non-mutable access to the convex hulls. </div></div>
<div class="ttc" id="classOpenMS_1_1MSSpectrum_html_a93f1c9f620403aec160d98e83e125502"><div class="ttname"><a href="classOpenMS_1_1MSSpectrum.html#a93f1c9f620403aec160d98e83e125502">OpenMS::MSSpectrum::sortByIntensity</a></div><div class="ttdeci">void sortByIntensity(bool reverse=false)</div><div class="ttdoc">Lexicographically sorts the peaks by their intensity. </div><div class="ttdef"><b>Definition:</b> MSSpectrum.h:314</div></div>
<div class="ttc" id="classOpenMS_1_1OfflinePrecursorIonSelection_html_afc3ec85793adf0cdaacc5803e59f6c17"><div class="ttname"><a href="classOpenMS_1_1OfflinePrecursorIonSelection.html#afc3ec85793adf0cdaacc5803e59f6c17">OpenMS::OfflinePrecursorIonSelection::makePrecursorSelectionForKnownLCMSMap</a></div><div class="ttdeci">void makePrecursorSelectionForKnownLCMSMap(const FeatureMap<> &features, const MSExperiment< InputPeakType > &experiment, MSExperiment< InputPeakType > &ms2, std::set< Int > &charges_set, bool feature_based)</div><div class="ttdoc">Makes the precursor selection for a given feature map, either feature or scan based. </div><div class="ttdef"><b>Definition:</b> OfflinePrecursorIonSelection.h:342</div></div>
<div class="ttc" id="classOpenMS_1_1OfflinePrecursorIonSelection_html_afb8f1202d548f01394d85fbeb51f9adf"><div class="ttname"><a href="classOpenMS_1_1OfflinePrecursorIonSelection.html#afb8f1202d548f01394d85fbeb51f9adf">OpenMS::OfflinePrecursorIonSelection::solver_</a></div><div class="ttdeci">LPWrapper::SOLVER solver_</div><div class="ttdef"><b>Definition:</b> OfflinePrecursorIonSelection.h:134</div></div>
<div class="ttc" id="classOpenMS_1_1MSSpectrum_html_a24f081fe74bbbf773533ad19e96d8358"><div class="ttname"><a href="classOpenMS_1_1MSSpectrum.html#a24f081fe74bbbf773533ad19e96d8358">OpenMS::MSSpectrum::setMSLevel</a></div><div class="ttdeci">void setMSLevel(UInt ms_level)</div><div class="ttdoc">Sets the MS level. </div><div class="ttdef"><b>Definition:</b> MSSpectrum.h:237</div></div>
<div class="ttc" id="classOpenMS_1_1Feature_html"><div class="ttname"><a href="classOpenMS_1_1Feature.html">OpenMS::Feature</a></div><div class="ttdoc">An LC-MS feature. </div><div class="ttdef"><b>Definition:</b> Feature.h:66</div></div>
<div class="ttc" id="classOpenMS_1_1OfflinePrecursorIonSelection_html_a2b206fdbe8649b90b854129f37fc594a"><div class="ttname"><a href="classOpenMS_1_1OfflinePrecursorIonSelection.html#a2b206fdbe8649b90b854129f37fc594a">OpenMS::OfflinePrecursorIonSelection::setLPSolver</a></div><div class="ttdeci">void setLPSolver(LPWrapper::SOLVER solver)</div><div class="ttdef"><b>Definition:</b> OfflinePrecursorIonSelection.h:100</div></div>
<div class="ttc" id="IDMapper_8h_html"><div class="ttname"><a href="IDMapper_8h.html">IDMapper.h</a></div></div>
<div class="ttc" id="classOpenMS_1_1MSExperiment_html_ae9dafbb74a1a4c430ceef03e25ab8e10"><div class="ttname"><a href="classOpenMS_1_1MSExperiment.html#ae9dafbb74a1a4c430ceef03e25ab8e10">OpenMS::MSExperiment::addSpectrum</a></div><div class="ttdeci">void addSpectrum(const MSSpectrum< PeakT > &spectrum)</div><div class="ttdoc">adds a spectra to the list </div><div class="ttdef"><b>Definition:</b> MSExperiment.h:738</div></div>
<div class="ttc" id="classOpenMS_1_1Exception_1_1InvalidSize_html"><div class="ttname"><a href="classOpenMS_1_1Exception_1_1InvalidSize.html">OpenMS::Exception::InvalidSize</a></div><div class="ttdoc">Invalid UInt exception. </div><div class="ttdef"><b>Definition:</b> Exception.h:304</div></div>
<div class="ttc" id="classOpenMS_1_1OfflinePrecursorIonSelection_html_acca2e18b840cc4b903f2958efca7aa1a"><div class="ttname"><a href="classOpenMS_1_1OfflinePrecursorIonSelection.html#acca2e18b840cc4b903f2958efca7aa1a">OpenMS::OfflinePrecursorIonSelection::checkMassRanges_</a></div><div class="ttdeci">void checkMassRanges_(std::vector< std::vector< std::pair< Size, Size > > > &mass_ranges, const MSExperiment< InputPeakType > &experiment)</div><div class="ttdoc">Eliminates overlapping peaks. </div><div class="ttdef"><b>Definition:</b> OfflinePrecursorIonSelection.h:572</div></div>
<div class="ttc" id="classOpenMS_1_1MSExperiment_html"><div class="ttname"><a href="classOpenMS_1_1MSExperiment.html">OpenMS::MSExperiment</a></div><div class="ttdoc">Representation of a mass spectrometry experiment. </div><div class="ttdef"><b>Definition:</b> MSExperiment.h:68</div></div>
<div class="ttc" id="classOpenMS_1_1MSExperiment_html_aee0a1a68d0fab63eefa6e14010b4a7f0"><div class="ttname"><a href="classOpenMS_1_1MSExperiment.html#aee0a1a68d0fab63eefa6e14010b4a7f0">OpenMS::MSExperiment::ConstIterator</a></div><div class="ttdeci">std::vector< SpectrumType >::const_iterator ConstIterator</div><div class="ttdoc">Non-mutable iterator. </div><div class="ttdef"><b>Definition:</b> MSExperiment.h:103</div></div>
<div class="ttc" id="classOpenMS_1_1OfflinePrecursorIonSelection_html"><div class="ttname"><a href="classOpenMS_1_1OfflinePrecursorIonSelection.html">OpenMS::OfflinePrecursorIonSelection</a></div><div class="ttdoc">Implements different algorithms for precursor ion selection. </div><div class="ttdef"><b>Definition:</b> OfflinePrecursorIonSelection.h:61</div></div>
<div class="ttc" id="classOpenMS_1_1SpectrumSettings_html_ad85e4e6fe1ec77097d9108f38a8e57ae"><div class="ttname"><a href="classOpenMS_1_1SpectrumSettings.html#ad85e4e6fe1ec77097d9108f38a8e57ae">OpenMS::SpectrumSettings::setPrecursors</a></div><div class="ttdeci">void setPrecursors(const std::vector< Precursor > &precursors)</div><div class="ttdoc">sets the precursors </div></div>
<div class="ttc" id="classOpenMS_1_1PSLPFormulation_html_a2a93e2928cdbb7704c7d49a33f13b091"><div class="ttname"><a href="classOpenMS_1_1PSLPFormulation.html#a2a93e2928cdbb7704c7d49a33f13b091">OpenMS::PSLPFormulation::createAndSolveILPForKnownLCMSMapFeatureBased</a></div><div class="ttdeci">void createAndSolveILPForKnownLCMSMapFeatureBased(const FeatureMap<> &features, const MSExperiment< InputPeakType > &experiment, std::vector< IndexTriple > &variable_indices, std::vector< std::vector< std::pair< Size, Size > > > &mass_ranges, std::set< Int > &charges_set, UInt ms2_spectra_per_rt_bin, std::vector< int > &solution_indices)</div><div class="ttdoc">Encode ILP formulation for a given LC-MS map, but unknown protein sample. </div><div class="ttdef"><b>Definition:</b> PSLPFormulation.h:289</div></div>
<div class="ttc" id="classOpenMS_1_1OfflinePrecursorIonSelection_html_ae895e8f290c7f5869249e59a7b53c417"><div class="ttname"><a href="classOpenMS_1_1OfflinePrecursorIonSelection.html#ae895e8f290c7f5869249e59a7b53c417">OpenMS::OfflinePrecursorIonSelection::getMassRanges</a></div><div class="ttdeci">void getMassRanges(const FeatureMap<> &features, const MSExperiment< InputPeakType > &experiment, std::vector< std::vector< std::pair< Size, Size > > > &indices)</div><div class="ttdoc">Calculates the mass ranges for each feature and stores them as indices of the raw data...</div><div class="ttdef"><b>Definition:</b> OfflinePrecursorIonSelection.h:154</div></div>
<div class="ttc" id="classOpenMS_1_1MSSpectrum_html_aac8c50f318154e9e5b7b01aee806e17f"><div class="ttname"><a href="classOpenMS_1_1MSSpectrum.html#aac8c50f318154e9e5b7b01aee806e17f">OpenMS::MSSpectrum::setRT</a></div><div class="ttdeci">void setRT(DoubleReal rt)</div><div class="ttdoc">Sets the absolute retention time (is seconds) </div><div class="ttdef"><b>Definition:</b> MSSpectrum.h:221</div></div>
<div class="ttc" id="group__Concept_html_gaf9ecec2d692138fab9167164a457cbd4"><div class="ttname"><a href="group__Concept.html#gaf9ecec2d692138fab9167164a457cbd4">OpenMS::Size</a></div><div class="ttdeci">size_t Size</div><div class="ttdoc">Size type e.g. used as variable which can hold result of size() </div><div class="ttdef"><b>Definition:</b> Types.h:144</div></div>
<div class="ttc" id="structOpenMS_1_1PairComparatorSecondElement_html"><div class="ttname"><a href="structOpenMS_1_1PairComparatorSecondElement.html">OpenMS::PairComparatorSecondElement</a></div><div class="ttdoc">Class for comparison of std::pair using second ONLY e.g. for use with std::sort. </div><div class="ttdef"><b>Definition:</b> ComparatorUtils.h:340</div></div>
<div class="ttc" id="classOpenMS_1_1MSExperiment_html_ac6e61de369e994009e36f344f99c15ad"><div class="ttname"><a href="classOpenMS_1_1MSExperiment.html#ac6e61de369e994009e36f344f99c15ad">OpenMS::MSExperiment::empty</a></div><div class="ttdeci">bool empty() const </div><div class="ttdef"><b>Definition:</b> MSExperiment.h:127</div></div>
<div class="ttc" id="classOpenMS_1_1IntList_html_ac244d9a257aab1a24c1a622ccacab745"><div class="ttname"><a href="classOpenMS_1_1IntList.html#ac244d9a257aab1a24c1a622ccacab745">OpenMS::IntList::create</a></div><div class="ttdeci">static IntList create(const String &list)</div><div class="ttdoc">Returns a list that is created by splitting the given comma-separated string (String are not trimmed!...</div></div>
<div class="ttc" id="classOpenMS_1_1Precursor_html_aca40a4630f5f8a0568899b74a62dcd4b"><div class="ttname"><a href="classOpenMS_1_1Precursor.html#aca40a4630f5f8a0568899b74a62dcd4b">OpenMS::Precursor::setCharge</a></div><div class="ttdeci">void setCharge(Int charge)</div><div class="ttdoc">Mutable access to the charge. </div></div>
<div class="ttc" id="classOpenMS_1_1DefaultParamHandler_html"><div class="ttname"><a href="classOpenMS_1_1DefaultParamHandler.html">OpenMS::DefaultParamHandler</a></div><div class="ttdoc">A base class for all classes handling default parameters. </div><div class="ttdef"><b>Definition:</b> DefaultParamHandler.h:90</div></div>
<div class="ttc" id="classOpenMS_1_1DBoundingBox_html"><div class="ttname"><a href="classOpenMS_1_1DBoundingBox.html">OpenMS::DBoundingBox</a></div><div class="ttdoc">A D-dimensional bounding box. </div><div class="ttdef"><b>Definition:</b> DBoundingBox.h:51</div></div>
<div class="ttc" id="classOpenMS_1_1OfflinePrecursorIonSelection_html_aff2cbb90bd3ebac2208d25616c0ecf7f"><div class="ttname"><a href="classOpenMS_1_1OfflinePrecursorIonSelection.html#aff2cbb90bd3ebac2208d25616c0ecf7f">OpenMS::OfflinePrecursorIonSelection::IndexTriple</a></div><div class="ttdeci">PSLPFormulation::IndexTriple IndexTriple</div><div class="ttdef"><b>Definition:</b> OfflinePrecursorIonSelection.h:65</div></div>
<div class="ttc" id="classOpenMS_1_1MSSpectrum_html_a3da529bd3240fa0d7148484bbef0b9d7"><div class="ttname"><a href="classOpenMS_1_1MSSpectrum.html#a3da529bd3240fa0d7148484bbef0b9d7">OpenMS::MSSpectrum::getRT</a></div><div class="ttdeci">DoubleReal getRT() const </div><div class="ttdef"><b>Definition:</b> MSSpectrum.h:215</div></div>
<div class="ttc" id="group__Concept_html_ga7cc214a236ad3bb6ad435bdcf5262a3f"><div class="ttname"><a href="group__Concept.html#ga7cc214a236ad3bb6ad435bdcf5262a3f">OpenMS::Int</a></div><div class="ttdeci">int Int</div><div class="ttdoc">Signed integer type. </div><div class="ttdef"><b>Definition:</b> Types.h:100</div></div>
<div class="ttc" id="PSLPFormulation_8h_html"><div class="ttname"><a href="PSLPFormulation_8h.html">PSLPFormulation.h</a></div></div>
<div class="ttc" id="FeatureMap_8h_html"><div class="ttname"><a href="FeatureMap_8h.html">FeatureMap.h</a></div></div>
<div class="ttc" id="classOpenMS_1_1DoubleList_html"><div class="ttname"><a href="classOpenMS_1_1DoubleList.html">OpenMS::DoubleList</a></div><div class="ttdoc">DoubleReal list. </div><div class="ttdef"><b>Definition:</b> DoubleList.h:56</div></div>
<div class="ttc" id="classOpenMS_1_1OfflinePrecursorIonSelection_html_a09cf655b910fa7166cd8b27d17e02453"><div class="ttname"><a href="classOpenMS_1_1OfflinePrecursorIonSelection.html#a09cf655b910fa7166cd8b27d17e02453">OpenMS::OfflinePrecursorIonSelection::calculateXICs_</a></div><div class="ttdeci">void calculateXICs_(const FeatureMap<> &features, const std::vector< std::vector< std::pair< Size, Size > > > &mass_ranges, const MSExperiment< InputPeakType > &experiment, const std::set< Int > &charges_set, std::vector< std::vector< std::pair< Size, DoubleReal > > > &xics)</div><div class="ttdoc">Calculate the sum of intensities of relevant features for each scan separately. </div><div class="ttdef"><b>Definition:</b> OfflinePrecursorIonSelection.h:305</div></div>
<div class="ttc" id="group__Concept_html_gace75bfb1aba684e874dffee13738bd15"><div class="ttname"><a href="group__Concept.html#gace75bfb1aba684e874dffee13738bd15">OpenMS::DoubleReal</a></div><div class="ttdeci">double DoubleReal</div><div class="ttdoc">Double-precision real type. </div><div class="ttdef"><b>Definition:</b> Types.h:118</div></div>
<div class="ttc" id="classOpenMS_1_1IntList_html"><div class="ttname"><a href="classOpenMS_1_1IntList.html">OpenMS::IntList</a></div><div class="ttdoc">Int list. </div><div class="ttdef"><b>Definition:</b> IntList.h:56</div></div>
<div class="ttc" id="classOpenMS_1_1PSLPFormulation_html"><div class="ttname"><a href="classOpenMS_1_1PSLPFormulation.html">OpenMS::PSLPFormulation</a></div><div class="ttdoc">Implements ILP formulation of precursor selection problems. </div><div class="ttdef"><b>Definition:</b> PSLPFormulation.h:51</div></div>
</div><!-- fragment --></div><!-- contents -->
<HR style="height:1px; border:none; border-top:1px solid #c0c0c0;">
<TABLE width="100%" border="0">
<TR>
<TD><font color="#c0c0c0">OpenMS / TOPP release 1.11.1</font></TD>
<TD align="right"><font color="#c0c0c0">Documentation generated on Thu Nov 14 2013 11:19:18 using doxygen 1.8.5</font></TD>
</TR>
</TABLE>
</BODY>
</HTML>
|