libpappsomspp
Library for mass spectrometry
tracedetectionzivy.cpp
Go to the documentation of this file.
1 
2 /*******************************************************************************
3  * Copyright (c) 2015 Olivier Langella <Olivier.Langella@moulon.inra.fr>.
4  *
5  * This file is part of the PAPPSOms++ library.
6  *
7  * PAPPSOms++ is free software: you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation, either version 3 of the License, or
10  * (at your option) any later version.
11  *
12  * PAPPSOms++ is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with PAPPSOms++. If not, see <http://www.gnu.org/licenses/>.
19  *
20  * Contributors:
21  * Olivier Langella <Olivier.Langella@moulon.inra.fr> - initial API and
22  *implementation
23  ******************************************************************************/
24 
25 #include "tracedetectionzivy.h"
26 #include "tracepeak.h"
27 #include "../../exception/exceptionoutofrange.h"
28 #include "../../exception/exceptionnotpossible.h"
29 
30 namespace pappso
31 {
33  unsigned int smoothing_half_window_length,
34  unsigned int minmax_half_window_length,
35  unsigned int maxmin_half_window_length,
36  pappso_double detection_threshold_on_minmax,
37  pappso_double detection_threshold_on_maxmin)
38  : m_smooth(smoothing_half_window_length),
39  m_minMax(minmax_half_window_length),
40  m_maxMin(maxmin_half_window_length)
41 {
42  m_detectionThresholdOnMaxMin = detection_threshold_on_maxmin;
43  m_detectionThresholdOnMinMax = detection_threshold_on_minmax;
44 }
46 {
47 }
48 
49 void
51 {
52  m_smooth = smooth;
53 }
54 void
56 {
57  m_minMax = minMax;
58 }
59 void
61 {
62  m_maxMin = maxMin;
63 }
64 
65 void
67  double detectionThresholdOnMinMax)
68 {
69  m_detectionThresholdOnMinMax = detectionThresholdOnMinMax;
70 }
71 void
73  double detectionThresholdOnMaxMin)
74 {
75  m_detectionThresholdOnMaxMin = detectionThresholdOnMaxMin;
76 }
77 unsigned int
79 {
81 };
82 unsigned int
84 {
86 };
87 
88 unsigned int
90 {
92 };
95 {
97 };
100 {
102 };
103 
104 
105 void
108  bool remove_peak_base) const
109 {
110 
111  // detect peak positions on close curve : a peak is an intensity value
112  // strictly greater than the two surrounding values. In case of
113  // equality (very rare, can happen with some old old spectrometers) we
114  // take the last equal point to be the peak
115  qDebug();
116  std::size_t size = xic.size();
117  if(size < 4)
118  {
119  qDebug() << QObject::tr(
120  "The original XIC is too small to detect peaks (%1)")
121  .arg(xic.size());
122  return;
123  }
124  if(size <= m_maxMin.getMaxMinHalfEdgeWindows())
125  return;
126  if(size <= m_minMax.getMinMaxHalfEdgeWindows())
127  return;
128 
129  Trace xic_minmax(xic); //"close" courbe du haut
131  {
132  qDebug() << "f";
133  m_smooth.filter(xic_minmax);
134  }
135 
136  qDebug();
137  Trace xic_maxmin(xic_minmax); //"open" courbe du bas
138 
139  qDebug();
140 
141  try
142  {
143  qDebug() << "f1";
144  m_minMax.filter(xic_minmax);
145  qDebug() << "f2";
146  m_maxMin.filter(xic_maxmin);
147  }
148  catch(const ExceptionOutOfRange &e)
149  {
150  qDebug() << QObject::tr(
151  "The original XIC is too small to detect peaks (%1) :\n%2")
152  .arg(xic.size())
153  .arg(e.qwhat());
154  return;
155  }
156  qDebug() << "a";
157 
158  std::vector<DataPoint>::const_iterator previous_rt = xic_minmax.begin();
159  std::vector<DataPoint>::const_iterator current_rt = (previous_rt + 1);
160  std::vector<DataPoint>::const_iterator last_rt = (xic_minmax.end() - 1);
161 
162  std::vector<DataPoint>::const_iterator current_rt_on_maxmin =
163  (xic_maxmin.begin() + 1);
164 
165  std::vector<DataPoint>::const_iterator xic_position = xic.begin();
166  qDebug() << "b";
167  while(current_rt != last_rt)
168  // for (unsigned int i = 1, count = 0; i < xic_minmax.size() - 1; )
169  {
170  // conditions to have a peak
171  if((previous_rt->y < current_rt->y) &&
172  (current_rt->y > m_detectionThresholdOnMinMax) &&
173  (current_rt_on_maxmin->y > m_detectionThresholdOnMaxMin))
174  {
175  // here we test the last condition to have a peak
176 
177  // no peak case
178  if(current_rt->y < (current_rt + 1)->y)
179  {
180  //++i;
181  previous_rt = current_rt;
182  current_rt++;
183  current_rt_on_maxmin++;
184  }
185  // there is a peak here ! case
186  else if(current_rt->y > (current_rt + 1)->y)
187  {
188  // peak.setMaxXicElement(*current_rt);
189 
190  // find left boundary
191  std::vector<DataPoint>::const_iterator it_left =
192  moveLowerYLeftDataPoint(xic_minmax, current_rt);
193  // walk back
194  it_left = moveLowerYRigthDataPoint(xic_minmax, it_left);
195  // peak.setLeftBoundary(*it_left);
196 
197  // find right boundary
198  std::vector<DataPoint>::const_iterator it_right =
199  moveLowerYRigthDataPoint(xic_minmax, current_rt);
200  // walk back
201  it_right = moveLowerYLeftDataPoint(xic_minmax, it_right);
202  // peak.setRightBoundary(*it_right);
203 
204  // integrate peak surface :
205  auto it =
206  findFirstEqualOrGreaterX(xic_position, xic.end(), it_left->x);
207  xic_position =
208  findFirstEqualOrGreaterX(it, xic.end(), it_right->x) + 1;
209  // peak.setArea(areaTrace(it, xic_position));
210 
211  // find the maximum :
212  // peak.setMaxXicElement(*maxYDataPoint(it, xic_position));
213 
214  // areaTrace()
215  TracePeak peak(it, xic_position, remove_peak_base);
216  sink.setTracePeak(peak);
217  // }
218  //++i;
219  previous_rt = current_rt;
220  current_rt++;
221  current_rt_on_maxmin++;
222  }
223  // equality case, skipping equal points
224  else
225  {
226  // while (v_minmax[i] == v_minmax[i + 1]) {
227  //++i;
228  current_rt++;
229  current_rt_on_maxmin++;
230 
231  //++count;
232  }
233  }
234  // no chance to have a peak at all, continue looping
235  else
236  {
237  //++i;
238  previous_rt = current_rt;
239  current_rt++;
240  current_rt_on_maxmin++;
241  }
242  } // end loop for peaks
243  qDebug();
244 }
245 } // namespace pappso
transform the trace with the maximum of the minimum equivalent of the erode filter for pictures
Definition: filtermorpho.h:138
Trace & filter(Trace &data_points) const override
std::size_t getMaxMinHalfEdgeWindows() const
mean filter apply mean of y values inside the window : this results in a kind of smoothing
Definition: filtermorpho.h:210
std::size_t getMeanHalfEdgeWindows() const
transform the trace with the minimum of the maximum equivalent of the dilate filter for pictures
Definition: filtermorpho.h:118
std::size_t getMinMaxHalfEdgeWindows() const
Trace & filter(Trace &data_points) const override
virtual Trace & filter(Trace &data_points) const override
virtual const QString & qwhat() const
virtual void setTracePeak(TracePeak &xic_peak)=0
void setFilterMorphoMinMax(const FilterMorphoMinMax &m_minMax)
void setDetectionThresholdOnMaxmin(double detectionThresholdOnMaxMin)
unsigned int getMinMaxHalfEdgeWindows() const
void detect(const Trace &xic, TraceDetectionSinkInterface &sink, bool remove_peak_base) const override
detect peaks on a trace
unsigned int getSmoothingHalfEdgeWindows() const
pappso_double getDetectionThresholdOnMinmax() const
pappso_double getDetectionThresholdOnMaxmin() const
void setFilterMorphoMaxMin(const FilterMorphoMaxMin &maxMin)
TraceDetectionZivy(unsigned int smoothing_half_window_length, unsigned int minmax_half_window_length, unsigned int maxmin_half_window_length, pappso_double detection_threshold_on_minmax, pappso_double detection_threshold_on_maxmin)
pappso_double m_detectionThresholdOnMaxMin
void setFilterMorphoMean(const FilterMorphoMean &smooth)
pappso_double m_detectionThresholdOnMinMax
unsigned int getMaxMinHalfEdgeWindows() const
void setDetectionThresholdOnMinmax(double detectionThresholdOnMinMax)
A simple container of DataPoint instances.
Definition: trace.h:147
tries to keep as much as possible monoisotopes, removing any possible C13 peaks and changes multichar...
Definition: aa.cpp:39
std::vector< DataPoint >::iterator findFirstEqualOrGreaterX(std::vector< DataPoint >::iterator begin, std::vector< DataPoint >::iterator end, const double &value)
find the first element in which X is equal or greater than the value searched important : it implies ...
Definition: trace.cpp:32
std::vector< DataPoint >::const_iterator moveLowerYLeftDataPoint(const Trace &trace, std::vector< DataPoint >::const_iterator begin)
Move left to the lower value.
Definition: trace.cpp:184
double pappso_double
A type definition for doubles.
Definition: types.h:48
std::vector< DataPoint >::const_iterator moveLowerYRigthDataPoint(const Trace &trace, std::vector< DataPoint >::const_iterator begin)
Move right to the lower value.
Definition: trace.cpp:166