29 #include <boost/math/special_functions/bessel.hpp>
33 #pragma clang diagnostic push
34 #pragma clang diagnostic ignored "-Wfloat-equal"
35 #pragma clang diagnostic ignored "-Wconversion"
36 #pragma clang diagnostic ignored "-Wshorten-64-to-32"
44 template <
typename PeakType>
64 typedef std::multimap<UInt, BoxElement>
Box;
130 (*trans_intens_)[i] = intens;
153 inline typename MSSpectrum::const_iterator
MZBegin(
const double mz)
const
160 inline typename MSSpectrum::const_iterator
MZEnd(
const double mz)
const
167 inline typename MSSpectrum::const_iterator
end()
const
174 inline typename MSSpectrum::const_iterator
begin()
const
231 const double ampl_cutoff,
const bool check_PPMs);
242 const UInt RT_votes_cutoff,
const Int front_bound = -1,
const Int end_bound = -1);
260 inline double getLinearInterpolation(
const typename MSSpectrum::const_iterator& left_iter,
const double mz_pos,
const typename MSSpectrum::const_iterator& right_iter)
262 return left_iter->getIntensity() + (right_iter->getIntensity() - left_iter->getIntensity()) / (right_iter->getMZ() - left_iter->getMZ()) * (mz_pos - left_iter->getMZ());
271 inline double getLinearInterpolation(
const double mz_a,
const double intens_a,
const double mz_pos,
const double mz_b,
const double intens_b)
273 return intens_a + (intens_b - intens_a) / (mz_b - mz_a) * (mz_pos - mz_a);
316 const double seed_mz,
const UInt c,
const double ampl_cutoff);
325 const double seed_mz,
const UInt c,
const double ampl_cutoff);
336 const UInt c,
const UInt scan_index,
const bool check_PPMs,
const double transintens,
const double prev_score);
346 const UInt c,
const UInt scan_index,
const bool check_PPMs,
const double transintens,
const double prev_score);
372 const double intens,
const double rt,
const UInt MZ_begin,
const UInt MZ_end,
const double ref_intens);
390 const double intens,
const double rt,
const UInt MZ_begin,
const UInt MZ_end);
405 const UInt scan_index,
const UInt c,
const bool check_PPMs);
412 const UInt scan_index,
const UInt c,
const bool check_PPMs);
426 double old_frac_mass = c_mass - (
Int)(c_mass);
428 double new_frac_mass = new_mass - (
Int)(new_mass);
430 if (new_frac_mass - old_frac_mass > 0.5)
435 if (new_frac_mass - old_frac_mass < -0.5)
446 inline double getPPMs_(
const double mass_a,
const double mass_b)
const
448 return fabs(mass_a - mass_b) / (0.5 * (mass_a + mass_b)) * 1e6;
470 template <
typename PeakType>
476 template <
typename PeakType>
482 template <
typename PeakType>
488 template <
typename PeakType>
494 template <
typename PeakType>
497 tmp_boxes_ =
new std::vector<std::multimap<double, Box> >(1);
501 max_num_peaks_per_pattern_ = 3;
506 template <
typename PeakType>
509 max_charge_ = max_charge;
510 max_scan_size_ = max_scan_size;
512 intenstype_ = intenstype;
513 tmp_boxes_ =
new std::vector<std::multimap<double, Box> >(max_charge);
514 if (max_scan_size <= 0)
523 Int size_estimate((
Int)ceil(max_scan_size_ / (max_mz - min_mz)));
525 psi_.reserve(to_reserve);
526 prod_.reserve(to_reserve);
527 xs_.reserve(to_reserve);
532 template <
typename PeakType>
538 template <
typename PeakType>
541 Int spec_size((
Int)c_ref.size());
545 double value, T_boundary_left, T_boundary_right, old, c_diff, current, old_pos, my_local_MZ, my_local_lambda, origin, c_mz;
547 for (
Int my_local_pos = 0; my_local_pos < spec_size; ++my_local_pos)
550 old = 0; old_pos = (my_local_pos - from_max_to_left_ - 1 >= 0) ? c_ref[my_local_pos - from_max_to_left_ - 1].getMZ() : c_ref[0].getMZ() - min_spacing_;
555 for (
Int current_conv_pos = std::max(0, my_local_pos - from_max_to_left_); c_diff < T_boundary_right; ++current_conv_pos)
557 if (current_conv_pos >= spec_size)
559 value += 0.5 * old * min_spacing_;
563 c_mz = c_ref[current_conv_pos].getMZ();
564 c_diff = c_mz + origin;
567 current = c_diff > T_boundary_left && c_diff <= T_boundary_right ?
IsotopeWavelet::getValueByLambda(my_local_lambda, c_diff * charge + 1.) * c_ref[current_conv_pos].getIntensity() : 0;
569 value += 0.5 * (current + old) * (c_mz - old_pos);
577 c_trans[my_local_pos].setIntensity(value);
581 template <
typename PeakType>
584 Int spec_size((
Int)c_ref.size());
588 double value, T_boundary_left, T_boundary_right, c_diff, current, my_local_MZ, my_local_lambda, origin, c_mz;
590 for (
Int my_local_pos = 0; my_local_pos < spec_size; ++my_local_pos)
599 for (
Int current_conv_pos = std::max(0, my_local_pos - from_max_to_left_); c_diff < T_boundary_right; ++current_conv_pos)
601 if (current_conv_pos >= spec_size)
606 c_mz = c_ref[current_conv_pos].getMZ();
607 c_diff = c_mz + origin;
610 current = c_diff > T_boundary_left && c_diff <= T_boundary_right ?
IsotopeWavelet::getValueByLambda(my_local_lambda, c_diff * charge + 1.) * c_ref[current_conv_pos].getIntensity() : 0;
615 c_trans[my_local_pos].setIntensity(value);
619 template <
typename PeakType>
622 data_length_ = (
UInt) c_ref.size();
623 computeMinSpacing(c_ref);
624 Int wavelet_length = 0, quarter_length = 0;
629 typename MSSpectrum::const_iterator start_iter, end_iter;
630 for (
UInt i = 0; i < data_length_; ++i)
633 start_iter = c_ref.
MZEnd(c_ref[i].getMZ());
634 end_iter = c_ref.
MZBegin(c_ref[i].getMZ() + c_mz_cutoff);
635 wavelet_length = std::max((
SignedSize) wavelet_length, distance(start_iter, end_iter) + 1);
637 quarter_length = std::max((
SignedSize) quarter_length, distance(end_iter, start_iter) + 1);
644 wavelet_length = (
UInt) ceil(max_mz_cutoff_ / min_spacing_);
649 if (wavelet_length > (
Int) c_ref.size())
651 std::cout <<
"Warning: the extremal length of the wavelet is larger (" << wavelet_length <<
") than the number of data points (" << c_ref.size() <<
"). This might (!) severely affect the transform." << std::endl;
652 std::cout <<
"Minimal spacing: " << min_spacing_ << std::endl;
653 std::cout <<
"Warning/Error generated at scan with RT " << c_ref.
getRT() <<
"." << std::endl;
657 from_max_to_left_ = max_index;
658 from_max_to_right_ = wavelet_length - 1 - from_max_to_left_;
661 template <
typename PeakType>
664 min_spacing_ = INT_MAX;
665 for (
UInt c_conv_pos = 1; c_conv_pos < c_ref.size(); ++c_conv_pos)
667 min_spacing_ = std::min(min_spacing_, c_ref[c_conv_pos].getMZ() - c_ref[c_conv_pos - 1].getMZ());
671 template <
typename PeakType>
673 const MSSpectrum& ref,
const UInt scan_index,
const UInt c,
const double ampl_cutoff,
const bool check_PPMs)
675 Size scan_size(candidates.size());
677 typename MSSpectrum::const_iterator iter_start, iter_end, seed_iter, iter2;
678 double mz_cutoff, seed_mz, c_av_intens = 0, c_score = 0, c_sd_intens = 0, threshold = 0, help_mz, share, share_pos, bwd, fwd;
679 UInt MZ_start, MZ_end;
682 diffed[0].setIntensity(0); diffed[scan_size - 1].setIntensity(0);
684 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
685 std::stringstream stream;
686 stream <<
"diffed_" << ref.
getRT() <<
"_" <<
c + 1 <<
".trans\0";
687 std::ofstream ofile(stream.str().c_str());
692 for (
UInt i = 0; i < scan_size - 2; ++i)
694 share = candidates[i + 1].getIntensity(); share_pos = candidates[i + 1].getMZ();
695 bwd = (share - candidates[i].getIntensity()) / (share_pos - candidates[i].getMZ());
696 fwd = (candidates[i + 2].getIntensity() - share) / (candidates[i + 2].getMZ() - share_pos);
698 if (!(bwd >= 0 && fwd <= 0) || share > ref[i + 1].getIntensity())
700 diffed[i + 1].setIntensity(0);
703 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
704 ofile << diffed[i + 1].getMZ() <<
"\t" << diffed[i + 1].getIntensity() << std::endl;
710 for (
UInt i = 0; i < scan_size - 2; ++i)
712 share = candidates[i + 1].getIntensity(); share_pos = candidates[i + 1].getMZ();
713 bwd = (share - candidates[i].getIntensity()) / (share_pos - candidates[i].getMZ());
714 fwd = (candidates[i + 2].getIntensity() - share) / (candidates[i + 2].getMZ() - share_pos);
716 if (!(bwd >= 0 && fwd <= 0))
718 diffed[i + 1].setIntensity(0);
721 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
722 ofile << diffed[i + 1].getMZ() <<
"\t" << diffed[i + 1].getIntensity() << std::endl;
726 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
735 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
736 std::stringstream stream2;
737 stream2 <<
"sorted_cpu_" << candidates.
getRT() <<
"_" <<
c + 1 <<
".trans\0";
738 std::ofstream ofile2(stream2.str().c_str());
739 for (iter = c_sorted_candidate.
end() - 1; iter != c_sorted_candidate.
begin(); --iter)
741 ofile2 << iter->getMZ() <<
"\t" << iter->getIntensity() << std::endl;
746 std::vector<UInt> processed(scan_size, 0);
754 c_av_intens = getAvIntens_(candidates);
755 c_sd_intens = getSdIntens_(candidates, c_av_intens);
756 threshold = ampl_cutoff * c_sd_intens + c_av_intens;
759 for (iter = c_sorted_candidate.
end() - 1; iter != c_sorted_candidate.
begin(); --iter)
761 if (iter->getIntensity() <= 0)
766 seed_mz = iter->getMZ();
767 seed_iter = ref.
MZBegin(seed_mz);
769 if (seed_iter == ref.end() || processed[distance(ref.begin(), seed_iter)])
780 iter_end = ref.
MZEnd(seed_iter, seed_mz + mz_cutoff / (
c + 1.), ref.end());
781 if (iter_end == ref.end())
786 MZ_start = distance(ref.begin(), iter_start);
787 MZ_end = distance(ref.begin(), iter_end);
789 memset(&(processed[MZ_start]), 1,
sizeof(
UInt) * (MZ_end - MZ_start + 1));
793 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
794 if (trunc(seed_mz) == 874)
795 std::cout << seed_mz <<
"\t" << c_score << std::endl;
798 if (c_score <= 0 && c_score != -1000)
805 push2TmpBox_(seed_mz, scan_index,
c, c_score, iter->getIntensity(), ref.
getRT(), MZ_start, MZ_end);
808 iter2 = candidates.
MZBegin(help_mz);
809 if (iter2 == candidates.end() || iter2 == candidates.begin())
817 if (iter2 != candidates.end())
819 push2TmpBox_(iter2->getMZ(), scan_index,
c, 0, getLinearInterpolation(iter2 - 1, help_mz, iter2), ref.
getRT(), MZ_start, MZ_end);
824 iter2 = candidates.
MZBegin(help_mz);
825 if (iter2 == candidates.end() || iter2 == candidates.begin())
833 if (iter2 != candidates.end())
835 push2TmpBox_(iter2->getMZ(), scan_index,
c, 0, getLinearInterpolation(iter2 - 1, help_mz, iter2), ref.
getRT(), MZ_start, MZ_end);
840 clusterSeeds_(candidates, ref, scan_index,
c, check_PPMs);
843 template <
typename PeakType>
845 const UInt peak_cutoff,
const double seed_mz,
const UInt c,
const double ampl_cutoff)
847 double c_score = 0, c_val;
848 typename MSSpectrum::const_iterator c_left_iter2, c_right_iter2;
849 Int signal_size((
Int)candidate.size());
854 Int p_h_ind = 1, end = 4 * (peak_cutoff - 1) - 1;
856 std::vector<double> positions(end);
857 for (
Int i = 0; i < end; ++i)
862 double l_score = 0, mid_val = 0;
863 Int start_index = distance(candidate.begin(), candidate.
MZBegin(positions[0])) - 1;
864 for (
Int v = 1; v <= end; ++v, ++p_h_ind)
868 if (start_index < signal_size - 1)
873 while (candidate[start_index].getMZ() < positions[v - 1]);
875 if (start_index <= 0 || start_index >= signal_size - 1)
880 c_left_iter2 = candidate.begin() + start_index - 1;
881 c_right_iter2 = c_left_iter2 + 1;
883 c_val = c_left_iter2->getIntensity() + (c_right_iter2->getIntensity() - c_left_iter2->getIntensity()) / (c_right_iter2->getMZ() - c_left_iter2->getMZ()) * (positions[v - 1] - c_left_iter2->getMZ());
885 if (v == (
int)(ceil(end / 2.)))
891 if (p_h_ind % 2 == 1)
894 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
895 if (trunc(seed_mz) == 874)
896 std::cout << -c_val << std::endl;
902 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
903 if (trunc(seed_mz) == 874)
904 std::cout << c_val << std::endl;
909 start_index = distance(candidate.begin(), c_left_iter2);
912 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
913 std::ofstream ofile_score(
"scores.dat", ios::app);
914 std::ofstream ofile_check_score(
"check_scores.dat", ios::app);
916 ofile_check_score.close();
919 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
920 if (trunc(seed_mz) == 874)
921 std::cout <<
"final_score: " << seed_mz <<
"\t" << c_score <<
"\t l_score: " << l_score <<
"\t" << c_score - l_score - mid_val <<
"\t" << c_score - mid_val <<
"\t" << ampl_cutoff << std::endl;
924 if (c_score - mid_val <= 0)
929 if (c_score - mid_val <= ampl_cutoff)
934 if (l_score <= 0 || c_score - l_score - mid_val <= 0)
942 template <
typename PeakType>
944 const UInt peak_cutoff,
const double seed_mz,
const UInt c,
const double ampl_cutoff)
946 double c_score = 0, c_val;
947 typename MSSpectrum::const_iterator c_left_iter2, c_right_iter2;
951 Int p_h_ind = 1, end = 4 * (peak_cutoff - 1) - 1;
953 std::vector<double> positions(end);
954 for (
Int i = 0; i < end; ++i)
959 double l_score = 0, mid_val = 0;
960 Int start_index = distance(candidate.
begin(), candidate.
MZBegin(positions[0])) - 1;
961 for (
Int v = 1; v <= end; ++v, ++p_h_ind)
965 if (start_index < signal_size - 1)
970 while (candidate.
getMZ(start_index) < positions[v - 1]);
972 if (start_index <= 0 || start_index >= signal_size - 1)
977 c_left_iter2 = candidate.
begin() + start_index - 1;
978 c_right_iter2 = c_left_iter2 + 1;
981 if (v == (
int)(ceil(end / 2.)))
987 if (p_h_ind % 2 == 1)
996 start_index = distance(candidate.
begin(), c_left_iter2);
999 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
1000 std::ofstream ofile_score(
"scores.dat", ios::app);
1001 std::ofstream ofile_check_score(
"check_scores.dat", ios::app);
1002 ofile_score << c_check_point <<
"\t" << c_score << std::endl;
1003 ofile_score.close();
1004 ofile_check_score.close();
1007 if (l_score <= 0 || c_score - l_score - mid_val <= 0 || c_score - mid_val <= ampl_cutoff)
1015 template <
typename PeakType>
1019 typename std::multimap<double, Box>::iterator iter;
1020 typename Box::iterator box_iter;
1021 std::vector<BoxElement> final_box;
1022 double c_mz, av_score = 0, av_mz = 0, av_intens = 0, av_abs_intens = 0, count = 0;
1023 double virtual_av_mz = 0, virtual_av_intens = 0, virtual_av_abs_intens = 0, virtual_count = 0;
1025 for (iter = tmp_boxes_->at(
c).begin(); iter != tmp_boxes_->at(
c).end(); ++iter)
1028 Box& c_box = iter->second;
1029 av_score = 0; av_mz = 0; av_intens = 0; av_abs_intens = 0; count = 0;
1030 virtual_av_mz = 0; virtual_av_intens = 0; virtual_av_abs_intens = 0; virtual_count = 0;
1033 for (box_iter = c_box.begin(); box_iter != c_box.end(); ++box_iter)
1035 if (box_iter->second.score == 0)
1040 c_mz = box_iter->second.mz;
1041 virtual_av_intens += box_iter->second.intens;
1042 virtual_av_abs_intens += fabs(box_iter->second.intens);
1043 virtual_av_mz += c_mz * fabs(box_iter->second.intens);
1048 c_mz = box_iter->second.mz;
1049 av_score += box_iter->second.score;
1050 av_intens += box_iter->second.intens;
1051 av_abs_intens += fabs(box_iter->second.intens);
1052 av_mz += c_mz * fabs(box_iter->second.intens);
1059 av_intens = virtual_av_intens / virtual_count;
1061 av_mz = virtual_av_mz / virtual_av_abs_intens;
1067 av_mz /= av_abs_intens;
1071 c_box_element.
mz = av_mz;
1072 c_box_element.
c =
c;
1073 c_box_element.
score = av_score;
1074 c_box_element.
intens = av_intens;
1076 c_box_element.
RT = c_box.begin()->second.RT;
1077 final_box.push_back(c_box_element);
1080 Size num_o_feature = final_box.size();
1081 if (num_o_feature == 0)
1083 tmp_boxes_->at(
c).clear();
1088 std::vector<double> bwd_diffs(num_o_feature, 0);
1091 for (
Size i = 1; i < num_o_feature; ++i)
1093 bwd_diffs[i] = (final_box[i].intens - final_box[i - 1].intens) / (final_box[i].mz - final_box[i - 1].mz);
1096 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
1097 std::ofstream ofile_bwd(
"bwd_cpu.dat");
1098 for (
Size i = 0; i < num_o_feature; ++i)
1100 ofile_bwd << final_box[i].mz <<
"\t" << bwd_diffs[i] << std::endl;
1106 for (
Size i = 0; i < num_o_feature - 1; ++i)
1108 while (i < num_o_feature - 2)
1110 if (final_box[i].score > 0 || final_box[i].score == -1000)
1115 if (bwd_diffs[i] > 0 && bwd_diffs[i + 1] < 0)
1117 checkPositionForPlausibility_(candidate, ref, final_box[i].mz, final_box[i].
c, scan_index, check_PPMs, final_box[i].intens, final_box[i].score);
1122 tmp_boxes_->at(
c).clear();
1125 template <
typename PeakType>
1128 double av_intens = 0;
1129 for (
UInt i = 0; i < scan.size(); ++i)
1131 if (scan[i].getIntensity() >= 0)
1133 av_intens += scan[i].getIntensity();
1136 return av_intens / (double)scan.size();
1139 template <
typename PeakType>
1142 double res = 0, intens;
1143 for (
UInt i = 0; i < scan.size(); ++i)
1145 if (scan[i].getIntensity() >= 0)
1147 intens = scan[i].getIntensity();
1148 res += (intens -
mean) * (intens -
mean);
1151 return sqrt(res / (
double)(scan.size() - 1));
1154 template <
typename PeakType>
1157 std::vector<double> diffs(scan.size() - 1, 0);
1158 for (
UInt i = 0; i < scan.size() - 1; ++i)
1160 diffs[i] = scan[i + 1].getMZ() - scan[i].getMZ();
1163 sort(diffs.begin(), diffs.end());
1164 double av_MZ_spacing = 0;
1165 for (
UInt i = 0; i < diffs.size() / 2; ++i)
1167 av_MZ_spacing += diffs[i];
1170 return av_MZ_spacing / (diffs.size() / 2);
1173 template <
typename PeakType>
1176 double av_intens = 0;
1177 for (
UInt i = 0; i < scan.
size(); ++i)
1184 return av_intens / (double)scan.
size();
1187 template <
typename PeakType>
1190 double res = 0, intens;
1191 for (
UInt i = 0; i < scan.
size(); ++i)
1196 res += (intens -
mean) * (intens -
mean);
1199 return sqrt(res / (
double)(scan.
size() - 1));
1202 template <
typename PeakType>
1204 const double score,
const double intens,
const double rt,
const UInt MZ_begin,
const UInt MZ_end,
double ref_intens)
1208 typename std::multimap<double, Box>::iterator upper_iter(open_boxes_.upper_bound(mz));
1209 typename std::multimap<double, Box>::iterator lower_iter(open_boxes_.lower_bound(mz));
1211 if (lower_iter != open_boxes_.end())
1214 if (mz != lower_iter->first && lower_iter != open_boxes_.begin())
1220 typename std::multimap<double, Box>::iterator insert_iter;
1221 bool create_new_box =
true;
1222 if (lower_iter == open_boxes_.end())
1227 if (!open_boxes_.empty())
1229 if (fabs((--lower_iter)->first - mz) < dist_constraint)
1231 create_new_box =
false;
1232 insert_iter = lower_iter;
1237 create_new_box =
true;
1242 if (upper_iter == open_boxes_.end() && fabs(lower_iter->first - mz) < dist_constraint)
1244 insert_iter = lower_iter;
1245 create_new_box =
false;
1249 create_new_box =
true;
1254 if (upper_iter != open_boxes_.end() && lower_iter != open_boxes_.end())
1259 double dist_lower = fabs(lower_iter->first - mz);
1260 double dist_upper = fabs(upper_iter->first - mz);
1261 dist_lower = (dist_lower < dist_constraint) ? dist_lower : INT_MAX;
1262 dist_upper = (dist_upper < dist_constraint) ? dist_upper : INT_MAX;
1264 if (dist_lower >= dist_constraint && dist_upper >= dist_constraint)
1266 create_new_box =
true;
1270 insert_iter = (dist_lower < dist_upper) ? lower_iter : upper_iter;
1271 create_new_box =
false;
1280 if (create_new_box ==
false)
1282 std::pair<UInt, BoxElement> help2(scan, element);
1283 insert_iter->second.insert(help2);
1286 Box replacement(insert_iter->second);
1290 double c_mz = insert_iter->first * (insert_iter->second.size() - 1) + mz;
1291 c_mz /= ((double) insert_iter->second.size());
1294 open_boxes_.erase(insert_iter);
1295 std::pair<double, std::multimap<UInt, BoxElement> > help3(c_mz, replacement);
1296 open_boxes_.insert(help3);
1300 std::pair<UInt, BoxElement> help2(scan, element);
1301 std::multimap<UInt, BoxElement> help3;
1302 help3.insert(help2);
1303 std::pair<double, std::multimap<UInt, BoxElement> > help4(mz, help3);
1304 open_boxes_.insert(help4);
1308 template <
typename PeakType>
1310 const double score,
const double intens,
const double rt,
const UInt MZ_begin,
const UInt MZ_end)
1314 std::multimap<double, Box>& tmp_box(tmp_boxes_->at(
c));
1315 typename std::multimap<double, Box>::iterator upper_iter(tmp_box.upper_bound(mz));
1316 typename std::multimap<double, Box>::iterator lower_iter(tmp_box.lower_bound(mz));
1318 if (lower_iter != tmp_box.end())
1321 if (mz != lower_iter->first && lower_iter != tmp_box.begin())
1327 typename std::multimap<double, Box>::iterator insert_iter;
1328 bool create_new_box =
true;
1329 if (lower_iter == tmp_box.end())
1334 if (!tmp_box.empty())
1336 if (fabs((--lower_iter)->first - mz) < dist_constraint)
1338 create_new_box =
false;
1339 insert_iter = lower_iter;
1344 create_new_box =
true;
1349 if (upper_iter == tmp_box.end() && fabs(lower_iter->first - mz) < dist_constraint)
1351 insert_iter = lower_iter;
1352 create_new_box =
false;
1356 create_new_box =
true;
1361 if (upper_iter != tmp_box.end() && lower_iter != tmp_box.end())
1364 double dist_lower = fabs(lower_iter->first - mz);
1365 double dist_upper = fabs(upper_iter->first - mz);
1366 dist_lower = (dist_lower < dist_constraint) ? dist_lower : INT_MAX;
1367 dist_upper = (dist_upper < dist_constraint) ? dist_upper : INT_MAX;
1369 if (dist_lower >= dist_constraint && dist_upper >= dist_constraint)
1371 create_new_box =
true;
1375 insert_iter = (dist_lower < dist_upper) ? lower_iter : upper_iter;
1376 create_new_box =
false;
1384 if (create_new_box ==
false)
1386 std::pair<UInt, BoxElement> help2(scan, element);
1387 insert_iter->second.insert(help2);
1390 Box replacement(insert_iter->second);
1394 double c_mz = insert_iter->first * (insert_iter->second.size() - 1) + mz;
1395 c_mz /= ((double) insert_iter->second.size());
1398 tmp_box.erase(insert_iter);
1399 std::pair<double, std::multimap<UInt, BoxElement> > help3(c_mz, replacement);
1400 tmp_box.insert(help3);
1404 std::pair<UInt, BoxElement> help2(scan, element);
1405 std::multimap<UInt, BoxElement> help3;
1406 help3.insert(help2);
1408 std::pair<double, std::multimap<UInt, BoxElement> > help4(mz, help3);
1409 tmp_box.insert(help4);
1413 template <
typename PeakType>
1415 const UInt RT_votes_cutoff,
const Int front_bound,
const Int end_bound)
1417 typename std::multimap<double, Box>::iterator iter, iter2;
1419 if ((
Int)scan_index == end_bound && end_bound != (
Int)map.
size() - 1)
1421 for (iter = open_boxes_.begin(); iter != open_boxes_.end(); ++iter)
1423 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
1424 std::cout <<
"LOW THREAD insert in end_box " << iter->first << std::endl;
1425 typename Box::iterator dings;
1426 for (dings = iter->second.begin(); dings != iter->second.end(); ++dings)
1427 std::cout << map[dings->first].getRT() <<
"\t" << dings->second.c + 1 << std::endl;
1429 end_boxes_.insert(*iter);
1431 open_boxes_.
clear();
1435 for (iter = open_boxes_.begin(); iter != open_boxes_.end(); )
1437 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
1438 if (front_bound > 0)
1440 std::cout <<
"HIGH THREAD open box. " << iter->first <<
"\t current scan index " << scan_index <<
"\t" << ((iter->second.begin()))->first <<
"\t of last scan " << map.
size() - 1 <<
"\t" << front_bound << std::endl;
1445 UInt lastScan = (--(iter->second.end()))->first;
1446 if (scan_index - lastScan > RT_interleave + 1 || scan_index == map.
size() - 1)
1448 if (iter->second.begin()->first - front_bound <= RT_interleave + 1 && front_bound > 0)
1452 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
1453 std::cout <<
"HIGH THREAD insert in front_box " << iter->first << std::endl;
1455 front_boxes_.insert(*iter);
1456 open_boxes_.erase(iter);
1466 if (iter->second.size() >= RT_votes_cutoff)
1470 closed_boxes_.insert(*(--iter));
1472 open_boxes_.erase(iter);
1482 template <
typename PeakType>
1485 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
1486 std::cout <<
"**** CHECKING FOR BOX EXTENSIONS ****" << std::endl;
1490 typename Box::const_iterator iter;
1491 std::vector<double> elution_profile(box.size());
1493 for (iter = box.begin(); iter != box.end(); ++iter, ++index)
1495 for (
Size i = iter->second.MZ_begin; i != iter->second.MZ_end; ++i)
1497 elution_profile[index] += map[iter->second.RT_index][i].getIntensity();
1499 elution_profile[index] /= iter->second.MZ_end - iter->second.MZ_begin + 1.;
1503 Int max_index = INT_MIN;
1504 for (
Size i = 0; i < elution_profile.size(); ++i)
1506 if (elution_profile[i] > max)
1509 max = elution_profile[i];
1513 Int max_extension = (
Int)(elution_profile.size()) - 2 * max_index;
1515 double av_elution = 0;
1516 for (
Size i = 0; i < elution_profile.size(); ++i)
1518 av_elution += elution_profile[i];
1520 av_elution /= (double)elution_profile.
size();
1522 double sd_elution = 0;
1523 for (
Size i = 0; i < elution_profile.size(); ++i)
1525 sd_elution += (av_elution - elution_profile[i]) * (av_elution - elution_profile[i]);
1527 sd_elution /= (double)(elution_profile.size() - 1);
1528 sd_elution = sqrt(sd_elution);
1532 for (iter = box.begin(); iter != box.end(); ++iter, ++index)
1534 av_mz += iter->second.mz;
1535 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
1536 std::cout << iter->second.RT <<
"\t" << iter->second.mz <<
"\t" << iter->second.c + 1 << std::endl;
1539 av_mz /= (double)box.size();
1543 if ((
Int)(box.begin()->second.RT_index) - 1 < 0)
1548 UInt pre_index = box.begin()->second.RT_index - 1;
1549 typename MSSpectrum::const_iterator c_iter = map[pre_index].MZBegin(av_mz);
1550 double pre_elution = 0;
1552 double mz_start = map[pre_index + 1][box.begin()->second.MZ_begin].getMZ();
1553 double mz_end = map[pre_index + 1][box.begin()->second.MZ_end].getMZ();
1555 typename MSSpectrum::const_iterator mz_start_iter = map[pre_index].MZBegin(mz_start), mz_end_iter = map[pre_index].MZBegin(mz_end);
1556 for (
typename MSSpectrum::const_iterator mz_iter = mz_start_iter; mz_iter != mz_end_iter; ++mz_iter)
1558 pre_elution += mz_iter->getIntensity();
1563 if (pre_elution <= av_elution - 2 * sd_elution)
1568 Int first_index = box.
begin()->second.RT_index;
1569 for (
Int i = 1; i < max_extension; ++i)
1571 Int c_index = first_index - i;
1578 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
1579 std::cout << box.begin()->second.RT <<
"\t" << av_mz <<
"\t" << box.begin()->second.c + 1 <<
"\t" <<
" extending the box " << std::endl;
1582 push2Box_(av_mz, c_index, box.begin()->second.c, box.begin()->second.score, c_iter->getIntensity(),
1583 map[c_index].getRT(), box.
begin()->second.MZ_begin, box.begin()->second.MZ_end);
1587 template <
typename PeakType>
1591 typename std::multimap<double, Box>::iterator iter;
1592 typename Box::iterator box_iter;
1593 std::vector<BoxElement> final_box;
1594 double c_mz, av_score = 0, av_mz = 0, av_intens = 0, av_abs_intens = 0, count = 0;
1595 double virtual_av_mz = 0, virtual_av_intens = 0, virtual_av_abs_intens = 0, virtual_count = 0;
1597 for (iter = tmp_boxes_->at(
c).begin(); iter != tmp_boxes_->at(
c).end(); ++iter)
1599 Box& c_box = iter->second;
1600 av_score = 0; av_mz = 0; av_intens = 0; av_abs_intens = 0; count = 0;
1601 virtual_av_mz = 0; virtual_av_intens = 0; virtual_av_abs_intens = 0; virtual_count = 0;
1603 for (box_iter = c_box.begin(); box_iter != c_box.end(); ++box_iter)
1605 if (box_iter->second.score == 0)
1610 c_mz = box_iter->second.mz;
1611 virtual_av_intens += box_iter->second.intens;
1612 virtual_av_abs_intens += fabs(box_iter->second.intens);
1613 virtual_av_mz += c_mz * fabs(box_iter->second.intens);
1618 c_mz = box_iter->second.mz;
1619 av_score += box_iter->second.score;
1620 av_intens += box_iter->second.intens;
1621 av_abs_intens += fabs(box_iter->second.intens);
1622 av_mz += c_mz * fabs(box_iter->second.intens);
1629 av_intens = virtual_av_intens / virtual_count;
1631 av_mz = virtual_av_mz / virtual_av_abs_intens;
1637 av_mz /= av_abs_intens;
1641 c_box_element.
mz = av_mz;
1642 c_box_element.
c =
c;
1643 c_box_element.
score = av_score;
1644 c_box_element.
intens = av_intens;
1646 c_box_element.
RT = c_box.begin()->second.RT;
1648 final_box.push_back(c_box_element);
1651 UInt num_o_feature = final_box.size();
1652 if (num_o_feature == 0)
1654 tmp_boxes_->at(
c).clear();
1659 std::vector<double> bwd_diffs(num_o_feature, 0);
1662 for (
UInt i = 1; i < num_o_feature; ++i)
1664 bwd_diffs[i] = (final_box[i].intens - final_box[i - 1].intens) / (final_box[i].mz - final_box[i - 1].mz);
1667 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
1668 std::ofstream ofile_bwd(
"bwd_gpu.dat");
1669 for (
UInt i = 0; i < num_o_feature; ++i)
1671 ofile_bwd << final_box[i].mz <<
"\t" << bwd_diffs[i] << std::endl;
1676 for (
UInt i = 0; i < num_o_feature - 1; ++i)
1678 while (i < num_o_feature - 2)
1680 if (final_box[i].score > 0 || final_box[i].score == -1000)
1685 if (bwd_diffs[i] > 0 && bwd_diffs[i + 1] < 0)
1687 checkPositionForPlausibility_(candidates, ref, final_box[i].mz, final_box[i].
c, scan_index, check_PPMs, final_box[i].intens, final_box[i].score);
1692 tmp_boxes_->at(
c).clear();
1695 template <
typename PeakType>
1699 typename std::multimap<double, Box>::iterator iter;
1700 typename Box::iterator box_iter;
1701 UInt best_charge_index;
double best_charge_score, c_mz, c_RT;
UInt c_charge;
1702 double av_intens = 0, av_ref_intens = 0, av_score = 0, av_mz = 0, av_RT = 0, mz_cutoff, sum_of_ref_intenses_g;
1704 for (iter = closed_boxes_.begin(); iter != closed_boxes_.end(); ++iter)
1706 sum_of_ref_intenses_g = 0;
1707 Box& c_box = iter->second;
1708 std::vector<double> charge_votes(max_charge_, 0), charge_binary_votes(max_charge_, 0);
1709 bool restart =
false;
1713 for (box_iter = c_box.begin(); box_iter != c_box.end(); ++box_iter)
1715 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
1716 if (trunc(box_iter->second.mz) == 874)
1717 std::cout << box_iter->second.c <<
"\t" << box_iter->second.intens <<
"\t" << box_iter->second.score << std::endl;
1720 if (box_iter->second.score == -1000)
1726 charge_votes[box_iter->second.c] += box_iter->second.intens;
1727 ++charge_binary_votes[box_iter->second.c];
1736 best_charge_index = 0; best_charge_score = 0;
1737 for (
UInt i = 0; i < max_charge_; ++i)
1739 if (charge_votes[i] > best_charge_score)
1741 best_charge_index = i;
1742 best_charge_score = charge_votes[i];
1747 if (charge_binary_votes[best_charge_index] < RT_votes_cutoff && RT_votes_cutoff <= map.
size())
1752 c_charge = best_charge_index + 1;
1754 av_intens = 0; av_ref_intens = 0; av_score = 0; av_mz = 0; av_RT = 0;
1756 std::vector<DPosition<2> > point_set;
1757 double sum_of_ref_intenses_l;
1758 for (box_iter = c_box.begin(); box_iter != c_box.end(); ++box_iter)
1760 sum_of_ref_intenses_l = 0;
1761 c_mz = box_iter->second.mz;
1762 c_RT = box_iter->second.RT;
1768 point_set.push_back(
DPosition<2>(c_RT, c_mz + mz_cutoff / (
double)c_charge));
1770 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
1771 std::cout <<
"Intenstype: " << intenstype_ << std::endl;
1773 if (intenstype_ ==
"ref")
1776 const MSSpectrum& c_spec(map[box_iter->second.RT_index]);
1778 for (
unsigned int i = 0; i < mz_cutoff; ++i)
1784 while (h_iter != c_spec.begin())
1787 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
1788 if (trunc(c_mz) == 874)
1790 std::cout <<
"cmz: " << c_mz + i *
Constants::IW_NEUTRON_MASS / c_charge <<
"\t" << hc_iter->getMZ() <<
"\t" << hc_iter->getIntensity() <<
"\t" << h_iter->getMZ() <<
"\t" << h_iter->getIntensity() << std::endl;
1795 if (h_iter->getIntensity() > hc_iter->getIntensity() || (h_iter->getIntensity() == hc_iter->getIntensity() && hc_iter->getIntensity() == 0))
1805 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
1806 if (trunc(c_mz) == 874)
1811 sum_of_ref_intenses_l += hc_iter->getIntensity();
1812 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
1813 if (trunc(c_mz) == 874)
1815 std::cout << sum_of_ref_intenses_l <<
"********" << std::endl;
1821 if (best_charge_index == box_iter->second.c)
1823 av_score += box_iter->second.score;
1824 av_intens += box_iter->second.intens;
1825 av_ref_intens += box_iter->second.ref_intens;
1826 sum_of_ref_intenses_g += sum_of_ref_intenses_l;
1827 av_mz += c_mz * box_iter->second.intens;
1833 av_ref_intens /= (double)charge_binary_votes[best_charge_index];
1834 av_score /= (double)charge_binary_votes[best_charge_index];
1835 av_RT /= (double)c_box.size();
1837 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
1838 if (trunc(av_mz) == 874)
1839 std::cout << av_mz <<
"\t" << best_charge_index <<
"\t" << best_charge_score << std::endl;
1846 c_feature.
setConvexHulls(std::vector<ConvexHull2D>(1, c_conv_hull));
1849 if (intenstype_ ==
"corrected")
1852 av_intens /= exp(-2 * lambda) * boost::math::cyl_bessel_i(0.0, 2.0 * lambda);
1854 if (intenstype_ ==
"ref")
1856 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
1857 if (trunc(c_mz) == 874)
1859 std::cout << sum_of_ref_intenses_g <<
"####" << std::endl;
1863 av_intens = sum_of_ref_intenses_g;
1866 c_feature.
setMZ(av_mz);
1868 c_feature.
setRT(av_RT);
1876 template <
typename PeakType>
1878 const MSSpectrum& ref,
const double seed_mz,
const UInt c,
const UInt scan_index,
const bool check_PPMs,
const double transintens,
const double prev_score)
1880 typename MSSpectrum::const_iterator iter, ref_iter;
1884 iter = candidate.
MZBegin(seed_mz);
1886 if (iter == candidate.begin() || iter == candidate.end())
1891 std::pair<double, double> reals;
1892 ref_iter = ref.
MZBegin(seed_mz);
1894 double real_mz, real_intens;
1901 typename MSSpectrum::const_iterator h_iter = ref_iter, hc_iter = ref_iter;
1902 while (h_iter != ref.begin())
1905 if (h_iter->getIntensity() > hc_iter->getIntensity() || (h_iter->getIntensity() == hc_iter->getIntensity() && hc_iter->getIntensity() == 0))
1919 reals = checkPPMTheoModel_(ref, h_iter->getMZ(),
c);
1920 real_mz = reals.first; real_intens = reals.second;
1921 if (real_mz <= 0 || real_intens <= 0)
1925 real_mz = h_iter->getMZ();
1926 real_intens = h_iter->getIntensity();
1931 reals = std::pair<double, double>(seed_mz, ref_iter->getIntensity());
1932 real_mz = reals.first; real_intens = reals.second;
1934 if (real_mz <= 0 || real_intens <= 0)
1936 typename MSSpectrum::const_iterator h_iter = ref_iter, hc_iter = ref_iter;
1937 while (h_iter != ref.begin())
1940 if (h_iter->getIntensity() > hc_iter->getIntensity() || (h_iter->getIntensity() == hc_iter->getIntensity() && hc_iter->getIntensity() == 0))
1954 real_mz = h_iter->getMZ(); real_intens = h_iter->getIntensity();
1955 if (real_mz <= 0 || real_intens <= 0)
1959 real_mz = h_iter->getMZ();
1960 real_intens = h_iter->getIntensity();
1964 double c_score = scoreThis_(candidate, peak_cutoff, real_mz,
c, 0);
1973 typename MSSpectrum::const_iterator real_r_MZ_iter = ref.
MZBegin(real_l_MZ_iter, real_mz + mz_cutoff / (
c + 1.), ref.end());
1974 if (real_r_MZ_iter == ref.end())
1980 UInt real_mz_begin = distance(ref.begin(), real_l_MZ_iter);
1981 UInt real_mz_end = distance(ref.begin(), real_r_MZ_iter);
1983 if (prev_score == -1000)
1985 push2Box_(real_mz, scan_index,
c, prev_score, transintens, ref.
getRT(), real_mz_begin, real_mz_end, real_intens);
1989 push2Box_(real_mz, scan_index,
c, c_score, transintens, ref.
getRT(), real_mz_begin, real_mz_end, real_intens);
1994 template <
typename PeakType>
1996 const MSSpectrum& ref,
const double seed_mz,
const UInt c,
const UInt scan_index,
const bool check_PPMs,
const double transintens,
const double prev_score)
1998 typename MSSpectrum::const_iterator iter, ref_iter;
2002 iter = candidate.
MZBegin(seed_mz);
2004 if (iter == candidate.
begin() || iter == candidate.
end())
2009 std::pair<double, double> reals;
2010 ref_iter = ref.
MZBegin(seed_mz);
2012 double real_mz, real_intens;
2015 reals = checkPPMTheoModel_(ref, iter->getMZ(),
c);
2016 real_mz = reals.first; real_intens = reals.second;
2019 auto h_iter = ref_iter, hc_iter = ref_iter;
2020 while (h_iter != ref.begin())
2023 if (h_iter->getIntensity() > hc_iter->getIntensity() || (h_iter->getIntensity() == hc_iter->getIntensity() && hc_iter->getIntensity() == 0))
2038 reals = checkPPMTheoModel_(ref, h_iter->getMZ(),
c);
2039 real_mz = reals.first; real_intens = reals.second;
2041 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
2042 std::cout <<
"Plausibility check old_mz: " << iter->getMZ() <<
"\t" << real_mz << std::endl;
2045 if (real_mz <= 0 || real_intens <= 0)
2049 real_mz = h_iter->getMZ();
2050 real_intens = h_iter->getIntensity();
2055 reals = std::pair<double, double>(seed_mz, ref_iter->getIntensity());
2056 real_mz = reals.first; real_intens = reals.second;
2058 if (real_mz <= 0 || real_intens <= 0)
2060 auto h_iter = ref_iter, hc_iter = ref_iter;
2061 while (h_iter != ref.begin())
2064 if (h_iter->getIntensity() > hc_iter->getIntensity() || (h_iter->getIntensity() == hc_iter->getIntensity() && hc_iter->getIntensity() == 0))
2078 real_mz = h_iter->getMZ(); real_intens = h_iter->getIntensity();
2079 if (real_mz <= 0 || real_intens <= 0)
2083 real_mz = h_iter->getMZ();
2084 real_intens = h_iter->getIntensity();
2088 double c_score = scoreThis_(candidate, peak_cutoff, real_mz,
c, 0);
2097 typename MSSpectrum::const_iterator real_r_MZ_iter = ref.
MZBegin(real_l_MZ_iter, real_mz + mz_cutoff / (
c + 1.), ref.end());
2098 if (real_r_MZ_iter == ref.end())
2104 UInt real_mz_begin = distance(ref.begin(), real_l_MZ_iter);
2105 UInt real_mz_end = distance(ref.begin(), real_r_MZ_iter);
2107 if (prev_score == -1000)
2109 push2Box_(real_mz, scan_index,
c, prev_score, transintens, ref.
getRT(), real_mz_begin, real_mz_end, real_intens);
2113 push2Box_(real_mz, scan_index,
c, c_score, transintens, ref.
getRT(), real_mz_begin, real_mz_end, real_intens);
2119 template <
typename PeakType>
2123 double ppms = getPPMs_(peptideMassRule_(mass), mass);
2126 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
2127 std::cout << ::std::setprecision(8) << std::fixed << c_mz <<
"\t =(" <<
"ISO_WAVE" <<
")> " <<
"REJECT \t" << ppms <<
" (rule: " << peptideMassRule_(mass) <<
" got: " << mass <<
")" << std::endl;
2129 return std::pair<double, double>(-1, -1);
2132 #ifdef OPENMS_DEBUG_ISOTOPE_WAVELET
2133 std::cout << ::std::setprecision(8) << std::fixed << c_mz <<
"\t =(" <<
"ISO_WAVE" <<
")> " <<
"ACCEPT \t" << ppms <<
" (rule: " << peptideMassRule_(mass) <<
" got: " << mass <<
")" << std::endl;
2135 return std::pair<double, double>(c_mz, ref.
MZBegin(c_mz)->getIntensity());
2140 #pragma clang diagnostic pop
void setCharge(const ChargeType &ch)
Set charge state.
Mutable iterator for the ConstRefVector.
Definition: ConstRefVector.h:196
This vector holds pointer to the elements of another container.
Definition: ConstRefVector.h:44
Iterator begin()
See std::vector documentation.
Definition: ConstRefVector.h:371
void sortByIntensity(bool reverse=false)
Sorting.
Definition: ConstRefVector.h:693
Iterator end()
See std::vector documentation.
Definition: ConstRefVector.h:377
A 2-dimensional hull representation in [counter]clockwise direction - depending on axis labelling.
Definition: ConvexHull2D.h:47
void addPoints(const PointArrayType &points)
void push_back(const VectorElement &f)
Definition: ExposedVector.h:161
A container for features.
Definition: FeatureMap.h:80
An LC-MS feature.
Definition: Feature.h:46
void setConvexHulls(const std::vector< ConvexHull2D > &hulls)
Set the convex hulls of single mass traces.
void setOverallQuality(QualityType q)
Set the overall quality.
static UInt getNumPeakCutOff(const double mass, const UInt z)
static IsotopeWavelet * init(const double max_m, const UInt max_charge)
static double getLambdaL(const double m)
Returns the mass-parameter lambda (linear fit).
static double getValueByLambda(const double lambda, const double tz1)
Returns the value of the isotope wavelet at position t via a fast table lookup. Usually,...
static UInt getMzPeakCutOffAtMonoPos(const double mass, const UInt z)
In-Memory representation of a mass spectrometry run.
Definition: MSExperiment.h:46
Iterator begin()
Definition: MSExperiment.h:156
Size size() const
The number of spectra.
Definition: MSExperiment.h:121
void clear(bool clear_meta_data)
Clears all data and meta data.
The representation of a 1D spectrum.
Definition: MSSpectrum.h:44
Iterator MZBegin(CoordinateType mz)
Binary search for peak range begin.
Iterator MZEnd(CoordinateType mz)
Binary search for peak range end (returns the past-the-end iterator)
A 2-dimensional raw data point or peak.
Definition: Peak2D.h:29
CoordinateType getMZ() const
Returns the m/z coordinate (index 1)
Definition: Peak2D.h:172
void setMZ(CoordinateType coordinate)
Mutable access to the m/z coordinate (index 1)
Definition: Peak2D.h:178
void setRT(CoordinateType coordinate)
Mutable access to the RT coordinate (index 0)
Definition: Peak2D.h:190
IntensityType getIntensity() const
Definition: Peak2D.h:142
void setIntensity(IntensityType intensity)
Sets data point intensity (height)
Definition: Peak2D.h:148
A more convenient string class.
Definition: String.h:34
int Int
Signed integer type.
Definition: Types.h:76
unsigned int UInt
Unsigned integer type.
Definition: Types.h:68
ptrdiff_t SignedSize
Signed Size type e.g. used as pointer difference.
Definition: Types.h:108
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:101
static double mean(IteratorType begin, IteratorType end)
Calculates the mean of a range of values.
Definition: StatisticFunctions.h:94
const double IW_PROTON_MASS
Definition: IsotopeWaveletConstants.h:42
const double PEPTIDE_MASS_RULE_FACTOR
Definition: IsotopeWaveletConstants.h:23
const double IW_NEUTRON_MASS
Definition: IsotopeWaveletConstants.h:28
const unsigned int DEFAULT_NUM_OF_INTERPOLATION_POINTS
Definition: IsotopeWaveletConstants.h:17
const double c
Definition: Constants.h:188
const double PEPTIDE_MASS_RULE_THEO_PPM_BOUND
Definition: IsotopeWaveletConstants.h:25
const double PEPTIDE_MASS_RULE_BOUND
Definition: IsotopeWaveletConstants.h:24
const double IW_HALF_NEUTRON_MASS
Definition: IsotopeWaveletConstants.h:29
const double IW_QUARTER_NEUTRON_MASS
Definition: IsotopeWaveletConstants.h:30
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:22
bool intensityComparator(const PeakType &a, const PeakType &b)
Definition: IsotopeWaveletTransform.h:471
bool intensityPointerComparator(PeakType *a, PeakType *b)
Definition: IsotopeWaveletTransform.h:483
bool positionComparator(const PeakType &a, const PeakType &b)
Definition: IsotopeWaveletTransform.h:489
bool intensityAscendingComparator(const PeakType &a, const PeakType &b)
Definition: IsotopeWaveletTransform.h:477