17 #include <TGraphErrors.h>
24 #include <tracking/calibration/Splitter.h>
40 map<ExpRun, pair<double, double>>
filter(
const map<
ExpRun, pair<double, double>>& runs,
double cut,
41 map<
ExpRun, pair<double, double>>& runsRemoved)
43 map<ExpRun, pair<double, double>> runsCopy;
47 double d = I.second - I.first;
49 runsCopy[r.first] = r.second;
51 runsRemoved[r.first] = r.second;
59 void plotRuns(vector<pair<double, double>> runs)
61 TGraphErrors* gr =
new TGraphErrors();
65 double m = (r.first + r.second) / 2;
66 double e = (r.second - r.first) / 2;
68 gr->SetPoint(gr->GetN(), m, 1);
69 gr->SetPointError(gr->GetN() - 1, e, 0);
72 gStyle->SetEndErrorSize(6);
75 gr->SetMarkerSize(40);
77 gr->GetYaxis()->SetRangeUser(-10, 10);
78 gr->GetXaxis()->SetRangeUser(0, 256);
80 gr->GetXaxis()->SetTitle(
"time [hours]");
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);
90 int s = (indx == 0) ? 0 : breaks[indx - 1] + 1;
91 int e = (indx == int(breaks.size())) ? nIntervals - 1 : breaks[indx];
96 void plotSRuns(vector<pair<double, double>> runs, vector<int> breaks,
int offset = 2)
98 TGraphErrors* gr =
new TGraphErrors();
100 for (
int i = 0; i < int(breaks.size()) + 1; ++i) {
103 double a = runs[s].first;
104 double b = runs[e].second;
106 double m = (a + b) / 2;
107 double err = (b - a) / 2;
109 gr->SetPoint(gr->GetN(), m, offset);
110 gr->SetPointError(gr->GetN() - 1, err, 0);
114 gr->SetLineColor(kRed);
115 gr->SetMarkerColor(kRed);
124 for (
auto r : runs) {
125 double d = r.second - r.first;
129 sort(dist.begin(), dist.end());
137 double Splitter::lossFunction(
const vector<Atom>& vec,
int s,
int e)
const
141 double rawTime = vec[e].t2 - vec[s].t1;
145 for (
int i = s; i <= e - 1; ++i) {
146 double d = vec[i + 1].t1 - vec[i].t2;
147 maxGap = max(maxGap, d);
152 for (
int i = s; i <= e; ++i) {
153 netTime += vec[i].t2 - vec[i].t1;
158 for (
int i = s; i <= e; ++i) {
161 if (nEv == 0) nEv = 0.1;
166 lossFun->SetParameters(rawTime, netTime, maxGap, nEv);
167 double lossNew = lossFun->Eval(0);
175 vector<pair<double, double>> Splitter::splitToSmall(map<
ExpRun, pair<double, double>> runs,
double intSize)
178 vector<pair<double, double>> smallRuns;
180 for (
auto r : runs) {
183 smallRuns.push_back(I);
187 double runTime = I.second - I.first;
188 int nSplits = runTime / intSize;
189 nSplits = max(1, nSplits);
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});
206 double Splitter::getMinLoss(
const vector<Atom>& vec,
int e, vector<int>& breaks)
209 if (cache[e].first >= 0) {
210 breaks = cache[e].second;
211 return cache[e].first;
215 vector<int> breaksOpt;
216 double minVal = 1e30;
218 for (
int i = -1; i <= e - 1; ++i) {
219 auto breaksNow = breaks;
222 r1 = getMinLoss(vec, i, breaksNow);
223 double r2 = lossFunction(vec, i + 1, e);
224 double tot = r1 + r2;
229 breaksOpt = breaksNow;
234 breaksOpt.push_back(iMin);
238 cache[e] = make_pair(minVal, breaks);
243 vector<int> Splitter::dynamicBreaks(
const vector<Atom>& runs)
246 cache.resize(runs.size());
247 for (
auto& c : cache)
248 c = make_pair(-1.0, vector<int>({}));
252 getMinLoss(runs, runs.size() - 1, breaks);
268 for (
auto r : runs) {
269 if (r.second.first <= t && t < r.second.second) {
275 B2ASSERT(
"Exactly one interval should be found", nFound == 1);
276 B2ASSERT(
"Assert that something was found", rFound !=
ExpRun(-1, -1));
283 const vector<Atom>& currVec,
const vector<int>& breaks)
285 vector<map<ExpRun, pair<double, double>>> splitsNow(breaks.size() + 1);
286 for (
int i = 0; i < int(breaks.size()) + 1; ++i) {
291 for (
int k = s; k <= e; ++k) {
292 ExpRun r =
getRun(runsMap, (currVec[k].t1 + currVec[k].t2) / 2.);
293 if (splitsNow[i].count(r)) {
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);
297 splitsNow[i][r].first = currVec[k].t1;
298 splitsNow[i][r].second = currVec[k].t2;
307 map<ExpRun, pair<double, double>> Splitter::mergeIntervals(map<
ExpRun, pair<double, double>> I1,
308 map<
ExpRun, pair<double, double>> I2)
310 map<ExpRun, pair<double, double>> I = I1;
313 if (I.count(run) == 0)
316 I.at(run) = make_pair(min(I1.at(run).first, I2.at(run).first), max(I1.at(run).second, I2.at(run).second));
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
void plotRuns(vector< pair< double, double >> runs)
plot runs on time axis
static ExpRun getRun(map< ExpRun, pair< double, double >> runs, double t)
Get exp number + run number from time.
void printBySize(vector< pair< double, double >> runs)
print sorted lenghts of the runs
void plotSRuns(vector< pair< double, double >> runs, vector< int > breaks, int offset=2)
plot clusters or runs on time axis
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.
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
Abstract base class for different kinds of events.
Struct containing exp number and run number.