Belle II Software  release-06-01-15
Splitter.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 
10 #include <vector>
11 #include <map>
12 #include <algorithm>
13 #include <iostream>
14 #include <utility>
15 #include <cmath>
16 #include <numeric>
17 #include <TGraphErrors.h>
18 #include <TStyle.h>
19 #include <TAxis.h>
20 
21 
22 //If compiled within BASF2
23 #ifdef _PACKAGE_
24 #include <tracking/calibration/Splitter.h>
25 #else
26 #include <Splitter.h>
27 #endif
28 
29 
30 using namespace std;
31 
32 namespace Belle2 {
39  //return filtered runs intervals (remove very short runs)
40  map<ExpRun, pair<double, double>> filter(const map<ExpRun, pair<double, double>>& runs, double cut,
41  map<ExpRun, pair<double, double>>& runsRemoved)
42  {
43  map<ExpRun, pair<double, double>> runsCopy;
44 
45  for (auto r : runs) {
46  auto& I = r.second;
47  double d = I.second - I.first;
48  if (d > cut)
49  runsCopy[r.first] = r.second;
50  else
51  runsRemoved[r.first] = r.second;
52  }
53 
54  return runsCopy;
55  }
56 
57 
59  void plotRuns(vector<pair<double, double>> runs)
60  {
61  TGraphErrors* gr = new TGraphErrors();
62 
63 
64  for (auto r : runs) {
65  double m = (r.first + r.second) / 2;
66  double e = (r.second - r.first) / 2;
67 
68  gr->SetPoint(gr->GetN(), m, 1);
69  gr->SetPointError(gr->GetN() - 1, e, 0);
70  }
71 
72  gStyle->SetEndErrorSize(6);
73 
74  gr->SetLineWidth(1);
75  gr->SetMarkerSize(40);
76  gr->Draw("ape");
77  gr->GetYaxis()->SetRangeUser(-10, 10);
78  gr->GetXaxis()->SetRangeUser(0, 256);
79 
80  gr->GetXaxis()->SetTitle("time [hours]");
81 
82  }
83 
85  pair<int, int> getStartEndIndexes(int nIntervals, vector<int> breaks, int indx)
86  {
87  B2ASSERT("There must be at least one interval", nIntervals >= 1);
88  B2ASSERT("Interval index must be positive", indx >= 0);
89  B2ASSERT("Interval index must be smaller than #breaks", indx < int(breaks.size()) + 1); //There is one more interval than #breaks
90  int s = (indx == 0) ? 0 : breaks[indx - 1] + 1;
91  int e = (indx == int(breaks.size())) ? nIntervals - 1 : breaks[indx];
92  return {s, e};
93  }
94 
96  void plotSRuns(vector<pair<double, double>> runs, vector<int> breaks, int offset = 2)
97  {
98  TGraphErrors* gr = new TGraphErrors();
99 
100  for (int i = 0; i < int(breaks.size()) + 1; ++i) {
101  int s, e;
102  tie(s, e) = getStartEndIndexes(runs.size(), breaks, i);
103  double a = runs[s].first;
104  double b = runs[e].second;
105 
106  double m = (a + b) / 2;
107  double err = (b - a) / 2;
108 
109  gr->SetPoint(gr->GetN(), m, offset);
110  gr->SetPointError(gr->GetN() - 1, err, 0);
111 
112  }
113 
114  gr->SetLineColor(kRed);
115  gr->SetMarkerColor(kRed);
116  gr->Draw("pe same");
117  }
118 
119 
121  void printBySize(vector<pair<double, double>> runs)
122  {
123  vector<double> dist;
124  for (auto r : runs) {
125  double d = r.second - r.first;
126  dist.push_back(d);
127  }
128 
129  sort(dist.begin(), dist.end());
130 
131  for (auto d : dist)
132  B2INFO(d);
133 
134  }
135 
137  double Splitter::lossFunction(const vector<Atom>& vec, int s, int e) const
138  {
139 
140  //raw time
141  double rawTime = vec[e].t2 - vec[s].t1;
142 
143  //max gap
144  double maxGap = 0;
145  for (int i = s; i <= e - 1; ++i) {
146  double d = vec[i + 1].t1 - vec[i].t2;
147  maxGap = max(maxGap, d);
148  }
149 
150  //net time
151  double netTime = 0;
152  for (int i = s; i <= e; ++i) {
153  netTime += vec[i].t2 - vec[i].t1;
154  }
155 
156  // Number of events
157  double nEv = 0;
158  for (int i = s; i <= e; ++i) {
159  nEv += vec[i].nEv;
160  }
161  if (nEv == 0) nEv = 0.1;
162 
163  //double loss = pow(rawTime - tBest, 2) + gapPenalty * pow(maxGap, 2);
164  //double loss = 1./nEv + timePenalty * pow(rawTime, 2);
165 
166  lossFun->SetParameters(rawTime, netTime, maxGap, nEv);
167  double lossNew = lossFun->Eval(0);
168 
169  return lossNew;
170  }
171 
172 
173 
174  // split to many small intervals (atoms)
175  vector<pair<double, double>> Splitter::splitToSmall(map<ExpRun, pair<double, double>> runs, double intSize)
176  {
177  // split into small intervals
178  vector<pair<double, double>> smallRuns;
179 
180  for (auto r : runs) {
181  auto& I = r.second;
182  if (intSize < 0) {
183  smallRuns.push_back(I);
184  continue;
185  }
186 
187  double runTime = I.second - I.first;
188  int nSplits = runTime / intSize; //1-m intervals
189  nSplits = max(1, nSplits); //at least 1 interval
190 
191  for (int i = 0; i < nSplits; ++i) {
192  double L = I.first + i * (runTime / nSplits);
193  double H = I.first + (i + 1) * (runTime / nSplits);
194  smallRuns.push_back({L, H});
195  }
196  }
197  return smallRuns;
198  }
199 
200 
201 
202 
203 
204 
205  // Get the optimal clustering of the atoms with indeces 0 .. e (recursive function with cache)
206  double Splitter::getMinLoss(const vector<Atom>& vec, int e, vector<int>& breaks)
207  {
208  // If entry in cache (speed up)
209  if (cache[e].first >= 0) {
210  breaks = cache[e].second;
211  return cache[e].first;
212  }
213 
214 
215  vector<int> breaksOpt;
216  double minVal = 1e30;
217  int iMin = -10;
218  for (int i = -1; i <= e - 1; ++i) {
219  auto breaksNow = breaks;
220  double r1 = 0;
221  if (i != -1)
222  r1 = getMinLoss(vec, i, breaksNow);
223  double r2 = lossFunction(vec, i + 1, e);
224  double tot = r1 + r2;
225 
226  if (tot < minVal) { //store minimum
227  minVal = tot;
228  iMin = i;
229  breaksOpt = breaksNow;
230  }
231  }
232 
233  if (iMin != -1)
234  breaksOpt.push_back(iMin);
235 
236 
237  breaks = breaksOpt;
238  cache[e] = make_pair(minVal, breaks); //store solution to cache
239  return minVal;
240  }
241 
242  //Get the indexes of the break-points
243  vector<int> Splitter::dynamicBreaks(const vector<Atom>& runs)
244  {
245  //reset cache
246  cache.resize(runs.size());
247  for (auto& c : cache)
248  c = make_pair(-1.0, vector<int>({}));
249 
250 
251  vector<int> breaks;
252  getMinLoss(runs, runs.size() - 1, breaks); //the minLoss (output) currently not used, only breaks
253 
254  return breaks;
255  }
256 
257 
258 
264  static ExpRun getRun(map<ExpRun, pair<double, double>> runs, double t)
265  {
266  ExpRun rFound(-1, -1);
267  int nFound = 0;
268  for (auto r : runs) { //Linear search over runs
269  if (r.second.first <= t && t < r.second.second) {
270  rFound = r.first;
271  ++nFound;
272  }
273  }
274 
275  B2ASSERT("Exactly one interval should be found", nFound == 1);
276  B2ASSERT("Assert that something was found", rFound != ExpRun(-1, -1));
277  return rFound;
278  }
279 
280 
281 
282  vector<map<ExpRun, pair<double, double>>> breaks2intervalsSep(const map<ExpRun, pair<double, double>>& runsMap,
283  const vector<Atom>& currVec, const vector<int>& breaks)
284  {
285  vector<map<ExpRun, pair<double, double>>> splitsNow(breaks.size() + 1);
286  for (int i = 0; i < int(breaks.size()) + 1; ++i) {
287  int s, e;
288  tie(s, e) = getStartEndIndexes(currVec.size(), breaks, i);
289 
290  // loop over atoms in single calib interval
291  for (int k = s; k <= e; ++k) {
292  ExpRun r = getRun(runsMap, (currVec[k].t1 + currVec[k].t2) / 2.); //runexp of the atom
293  if (splitsNow[i].count(r)) { //if already there
294  splitsNow[i].at(r).first = min(splitsNow[i].at(r).first, currVec[k].t1);
295  splitsNow[i].at(r).second = max(splitsNow[i].at(r).second, currVec[k].t2);
296  } else { //if new
297  splitsNow[i][r].first = currVec[k].t1;
298  splitsNow[i][r].second = currVec[k].t2;
299  }
300  }
301  }
302 
303  return splitsNow;
304  }
305 
306  //Merge two intervals together
307  map<ExpRun, pair<double, double>> Splitter::mergeIntervals(map<ExpRun, pair<double, double>> I1,
308  map<ExpRun, pair<double, double>> I2)
309  {
310  map<ExpRun, pair<double, double>> I = I1;
311  for (auto r : I2) {
312  ExpRun run = r.first;
313  if (I.count(run) == 0)
314  I[run] = r.second;
315  else {
316  I.at(run) = make_pair(min(I1.at(run).first, I2.at(run).first), max(I1.at(run).second, I2.at(run).second));
317  }
318  }
319  return I;
320  }
321 
322 
324 }
std::pair< int, int > getStartEndIndexes(int nIntervals, std::vector< int > breaks, int indx)
get the range of interval with nIntervals and breaks stored in a vector
Definition: Splitter.cc:85
void plotRuns(vector< pair< double, double >> runs)
plot runs on time axis
Definition: Splitter.cc:59
static ExpRun getRun(map< ExpRun, pair< double, double >> runs, double t)
Get exp number + run number from time.
Definition: Splitter.cc:264
void printBySize(vector< pair< double, double >> runs)
print sorted lenghts of the runs
Definition: Splitter.cc:121
void plotSRuns(vector< pair< double, double >> runs, vector< int > breaks, int offset=2)
plot clusters or runs on time axis
Definition: Splitter.cc:96
vector< map< ExpRun, pair< double, double > > > breaks2intervalsSep(const map< ExpRun, pair< double, double >> &runsMap, const vector< Atom > &currVec, const vector< int > &breaks)
Get calibration intervals according to the indexes of the breaks.
Definition: Splitter.cc:282
map< ExpRun, pair< double, double > > filter(const map< ExpRun, pair< double, double >> &runs, double cut, map< ExpRun, pair< double, double >> &runsRemoved)
filter events to remove runs shorter than cut, it stores removed runs in runsRemoved
Definition: Splitter.cc:40
Abstract base class for different kinds of events.
Struct containing exp number and run number.
Definition: Splitter.h:51