17#include <TGraphErrors.h>
24#include <reconstruction/calibration/BeamSpotBoostInvMass/Splitter.h>
38 std::map<ExpRun, std::pair<double, double>>
filter(
const std::map<
ExpRun, std::pair<double, double>>& runs,
double cut,
39 std::map<
ExpRun, std::pair<double, double>>& runsRemoved)
41 std::map<ExpRun, std::pair<double, double>> runsCopy;
45 double d = I.second - I.first;
47 runsCopy[r.first] = r.second;
49 runsRemoved[r.first] = r.second;
57 void plotRuns(std::vector<std::pair<double, double>> runs)
59 TGraphErrors* gr =
new TGraphErrors();
63 double m = (r.first + r.second) / 2;
64 double e = (r.second - r.first) / 2;
66 gr->SetPoint(gr->GetN(), m, 1);
67 gr->SetPointError(gr->GetN() - 1, e, 0);
70 gStyle->SetEndErrorSize(6);
73 gr->SetMarkerSize(40);
75 gr->GetYaxis()->SetRangeUser(-10, 10);
76 gr->GetXaxis()->SetRangeUser(0, 256);
78 gr->GetXaxis()->SetTitle(
"time [hours]");
85 B2ASSERT(
"There must be at least one interval", nIntervals >= 1);
86 B2ASSERT(
"Interval index must be positive", indx >= 0);
87 B2ASSERT(
"Interval index must be smaller than #breaks", indx <
int(breaks.size()) + 1);
88 int s = (indx == 0) ? 0 : breaks[indx - 1] + 1;
89 int e = (indx == int(breaks.size())) ? nIntervals - 1 : breaks[indx];
94 void plotSRuns(std::vector<std::pair<double, double>> runs, std::vector<int> breaks,
int offset = 2)
96 TGraphErrors* gr =
new TGraphErrors();
98 for (
int i = 0; i < int(breaks.size()) + 1; ++i) {
101 double a = runs[s].first;
102 double b = runs[e].second;
104 double m = (a + b) / 2;
105 double err = (b - a) / 2;
107 gr->SetPoint(gr->GetN(), m, offset);
108 gr->SetPointError(gr->GetN() - 1, err, 0);
112 gr->SetLineColor(kRed);
113 gr->SetMarkerColor(kRed);
121 std::vector<double> dist;
122 for (
auto r : runs) {
123 double d = r.second - r.first;
127 sort(dist.begin(), dist.end());
139 double rawTime = vec[e].t2 - vec[s].t1;
143 for (
int i = s; i <= e - 1; ++i) {
144 double d = vec[i + 1].t1 - vec[i].t2;
145 maxGap = std::max(maxGap, d);
150 for (
int i = s; i <= e; ++i) {
151 netTime += vec[i].t2 - vec[i].t1;
156 for (
int i = s; i <= e; ++i) {
159 if (nEv == 0) nEv = 0.1;
164 lossFun->SetParameters(rawTime, netTime, maxGap, nEv);
165 double lossNew =
lossFun->Eval(0);
176 std::vector<std::pair<double, double>> smallRuns;
178 for (
auto r : runs) {
181 smallRuns.push_back(I);
185 double runTime = I.second - I.first;
186 int nSplits = runTime / intSize;
187 nSplits = std::max(1, nSplits);
189 for (
int i = 0; i < nSplits; ++i) {
190 double L = I.first + i * (runTime / nSplits);
191 double H = I.first + (i + 1) * (runTime / nSplits);
192 smallRuns.push_back({L, H});
207 if (
cache[e].first >= 0) {
208 breaks =
cache[e].second;
209 return cache[e].first;
213 std::vector<int> breaksOpt;
214 double minVal = 1e30;
216 for (
int i = -1; i <= e - 1; ++i) {
217 auto breaksNow = breaks;
222 double tot = r1 + r2;
227 breaksOpt = breaksNow;
232 breaksOpt.push_back(iMin);
236 cache[e] = std::make_pair(minVal, breaks);
244 cache.resize(runs.size());
245 for (
auto& c :
cache)
246 c = std::make_pair(-1.0, std::vector<int>({}));
249 std::vector<int> breaks;
266 for (
auto r : runs) {
267 if (r.second.first <= t && t < r.second.second) {
273 B2ASSERT(
"Exactly one interval should be found", nFound == 1);
274 B2ASSERT(
"Assert that something was found", rFound !=
ExpRun(-1, -1));
282 const std::vector<Atom>& currVec,
const std::vector<int>& breaks)
284 std::vector<std::map<ExpRun, std::pair<double, double>>> splitsNow(breaks.size() + 1);
285 for (
int i = 0; i < int(breaks.size()) + 1; ++i) {
290 for (
int k = s; k <= e; ++k) {
291 ExpRun r =
getRun(runsMap, (currVec[k].t1 + currVec[k].t2) / 2.);
292 if (splitsNow[i].count(r)) {
293 splitsNow[i].at(r).first = std::min(splitsNow[i].at(r).first, currVec[k].t1);
294 splitsNow[i].at(r).second = std::max(splitsNow[i].at(r).second, currVec[k].t2);
296 splitsNow[i][r].first = currVec[k].t1;
297 splitsNow[i][r].second = currVec[k].t2;
307 std::map<
ExpRun, std::pair<double, double>> I2)
309 std::map<ExpRun, std::pair<double, double>> I = I1;
312 if (I.count(run) == 0)
315 I.at(run) = std::make_pair(std::min(I1.at(run).first, I2.at(run).first), std::max(I1.at(run).second, I2.at(run).second));
TF1 * lossFun
loss function used for clustering of Atoms
std::vector< std::pair< double, std::vector< int > > > cache
cache used by the clustering algorithm (has to be reset every time)
static std::map< ExpRun, std::pair< double, double > > mergeIntervals(std::map< ExpRun, std::pair< double, double > > I1, std::map< ExpRun, std::pair< double, double > > I2)
Merge two subintervals into one subinterval.
void plotRuns(std::vector< std::pair< double, double > > runs)
plot runs on time axis
static std::vector< std::pair< double, double > > splitToSmall(std::map< ExpRun, std::pair< double, double > > runs, double intSize=1./60)
Split the runs into small calibration intervals (atoms) of a specified size.
double getMinLoss(const std::vector< Atom > &vec, int e, std::vector< int > &breaks)
Recursive function to evaluate minimal sum of the lossFuctions for the optimal clustering.
static ExpRun getRun(std::map< ExpRun, std::pair< double, double > > runs, double t)
Get exp number + run number from time.
std::vector< std::map< ExpRun, std::pair< double, double > > > breaks2intervalsSep(const std::map< ExpRun, std::pair< double, double > > &runsMap, const std::vector< Atom > &currVec, const std::vector< int > &breaks)
Get calibration intervals according to the indexes of the breaks.
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
std::map< ExpRun, std::pair< double, double > > filter(const std::map< ExpRun, std::pair< double, double > > &runs, double cut, std::map< ExpRun, std::pair< double, double > > &runsRemoved)
filter events to remove runs shorter than cut, it stores removed runs in runsRemoved
std::vector< int > dynamicBreaks(const std::vector< Atom > &runs)
Get optimal break points using algorithm based on dynamic programming.
void plotSRuns(std::vector< std::pair< double, double > > runs, std::vector< int > breaks, int offset=2)
plot clusters or runs on time axis
double lossFunction(const std::vector< Atom > &vec, int s, int e) const
lossFunction of the calibration interval consisting of several "atoms" stored in vector vec The atoms...
void printBySize(std::vector< std::pair< double, double > > runs)
print sorted lengths of the runs
Abstract base class for different kinds of events.
Struct containing exp number and run number.