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
107  TraceDetectionSinkInterface &sink) const
108 {
109 
110  // detect peak positions on close curve : a peak is an intensity value
111  // strictly greater than the two surrounding values. In case of
112  // equality (very rare, can happen with some old old spectrometers) we
113  // take the last equal point to be the peak
114  qDebug();
115  std::size_t size = xic.size();
116  if(size < 4)
117  {
118  qDebug() << QObject::tr(
119  "The original XIC is too small to detect peaks (%1)")
120  .arg(xic.size());
121  return;
122  }
123  if(size <= m_maxMin.getMaxMinHalfEdgeWindows())
124  return;
125  if(size <= m_minMax.getMinMaxHalfEdgeWindows())
126  return;
127 
128  Trace xic_minmax(xic); //"close" courbe du haut
130  {
131  qDebug() << "f";
132  m_smooth.filter(xic_minmax);
133  }
134 
135  qDebug();
136  Trace xic_maxmin(xic_minmax); //"open" courbe du bas
137 
138  qDebug();
139 
140  try
141  {
142  qDebug() << "f1";
143  m_minMax.filter(xic_minmax);
144  qDebug() << "f2";
145  m_maxMin.filter(xic_maxmin);
146  }
147  catch(const ExceptionOutOfRange &e)
148  {
149  qDebug() << QObject::tr(
150  "The original XIC is too small to detect peaks (%1) :\n%2")
151  .arg(xic.size())
152  .arg(e.qwhat());
153  return;
154  }
155  qDebug() << "a";
156 
157  std::vector<DataPoint>::const_iterator previous_rt = xic_minmax.begin();
158  std::vector<DataPoint>::const_iterator current_rt = (previous_rt + 1);
159  std::vector<DataPoint>::const_iterator last_rt = (xic_minmax.end() - 1);
160 
161  std::vector<DataPoint>::const_iterator current_rt_on_maxmin =
162  (xic_maxmin.begin() + 1);
163 
164  std::vector<DataPoint>::const_iterator xic_position = xic.begin();
165  qDebug() << "b";
166  while(current_rt != last_rt)
167  // for (unsigned int i = 1, count = 0; i < xic_minmax.size() - 1; )
168  {
169  // conditions to have a peak
170  if((previous_rt->y < current_rt->y) &&
171  (current_rt->y > m_detectionThresholdOnMinMax) &&
172  (current_rt_on_maxmin->y > m_detectionThresholdOnMaxMin))
173  {
174  // here we test the last condition to have a peak
175 
176  // no peak case
177  if(current_rt->y < (current_rt + 1)->y)
178  {
179  //++i;
180  previous_rt = current_rt;
181  current_rt++;
182  current_rt_on_maxmin++;
183  }
184  // there is a peak here ! case
185  else if(current_rt->y > (current_rt + 1)->y)
186  {
187  // peak.setMaxXicElement(*current_rt);
188 
189  // find left boundary
190  std::vector<DataPoint>::const_iterator it_left =
191  moveLowerYLeftDataPoint(xic_minmax, current_rt);
192  // walk back
193  it_left = moveLowerYRigthDataPoint(xic_minmax, it_left);
194  // peak.setLeftBoundary(*it_left);
195 
196  // find right boundary
197  std::vector<DataPoint>::const_iterator it_right =
198  moveLowerYRigthDataPoint(xic_minmax, current_rt);
199  // walk back
200  it_right = moveLowerYLeftDataPoint(xic_minmax, it_right);
201  // peak.setRightBoundary(*it_right);
202 
203  // integrate peak surface :
204  auto it =
205  findFirstEqualOrGreaterX(xic_position, xic.end(), it_left->x);
206  xic_position =
207  findFirstEqualOrGreaterX(it, xic.end(), it_right->x) + 1;
208  // peak.setArea(areaTrace(it, xic_position));
209 
210  // find the maximum :
211  // peak.setMaxXicElement(*maxYDataPoint(it, xic_position));
212 
213  // areaTrace()
214  TracePeak peak(it, xic_position);
215  sink.setTracePeak(peak);
216  // }
217  //++i;
218  previous_rt = current_rt;
219  current_rt++;
220  current_rt_on_maxmin++;
221  }
222  // equality case, skipping equal points
223  else
224  {
225  // while (v_minmax[i] == v_minmax[i + 1]) {
226  //++i;
227  current_rt++;
228  current_rt_on_maxmin++;
229 
230  //++count;
231  }
232  }
233  // no chance to have a peak at all, continue looping
234  else
235  {
236  //++i;
237  previous_rt = current_rt;
238  current_rt++;
239  current_rt_on_maxmin++;
240  }
241  } // end loop for peaks
242  qDebug();
243 }
244 } // 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:204
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
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 detect(const Trace &xic, TraceDetectionSinkInterface &sink) const override
void setDetectionThresholdOnMinmax(double detectionThresholdOnMinMax)
A simple container of DataPoint instances.
Definition: trace.h:132
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:182
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:164