Belle II Software development
Splitter Class Reference

Class that allows to split runs into the intervals of intended properties given by the lossFunction. More...

#include <Splitter.h>

Public Member Functions

template<typename Evt >
std::vector< std::vector< std::map< ExpRun, std::pair< double, double > > > > getIntervals (const std::map< ExpRun, std::pair< double, double > > &runs, std::vector< Evt > evts, TString lossFunctionOuter, TString lossFunctionInner, double atomSize=3./60)
 Function to merge/divide runs into the calibration intervals of given characteristic length.
 

Static Public Member Functions

static std::pair< double, double > getStartEnd (std::vector< std::map< ExpRun, std::pair< double, double > > > res)
 Get the start/end time of the calibration interval (vector of the calib.
 
static std::vector< double > getBreaks (std::vector< std::map< ExpRun, std::pair< double, double > > > res)
 Get vector with breaks of the calib.
 
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.
 

Private Member Functions

std::vector< int > dynamicBreaks (const std::vector< Atom > &runs)
 Get optimal break points using algorithm based on dynamic programming.
 
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.
 
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 included in calibration interval have indices between s and e
 
template<typename Evt >
std::vector< AtomcreateAtoms (const std::vector< std::pair< double, double > > &atomsTimes, const std::vector< Evt > &evts)
 Get the vector with parameters of the calibration Atoms.
 
TF1 * toTF1 (TString LossString)
 Convert lossFunction from string to TF1.
 

Static Private Member Functions

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.
 

Private Attributes

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)
 

Detailed Description

Class that allows to split runs into the intervals of intended properties given by the lossFunction.

It typically avoids of having large time gaps within calibration intervals. The algorithm is based on minimizing the total loss function Notice, that the calibration interval is defined as a map, where for each run the included period is defined by the start and end time

Definition at line 108 of file Splitter.h.

Constructor & Destructor Documentation

◆ Splitter()

Splitter ( )
inline

Definition at line 111 of file Splitter.h.

111: lossFun(nullptr) {}
TF1 * lossFun
loss function used for clustering of Atoms
Definition: Splitter.h:297

Member Function Documentation

◆ createAtoms()

std::vector< Atom > createAtoms ( const std::vector< std::pair< double, double > > &  atomsTimes,
const std::vector< Evt > &  evts 
)
inlineprivate

Get the vector with parameters of the calibration Atoms.

Parameters
atomsTimesvector with start/end times of the atoms
evtsvector with events
Returns
A vector of atoms, each including start/end times and the number of events inside

Definition at line 260 of file Splitter.h.

261 {
262 std::vector<Atom> atoms(atomsTimes.size());
263
264 // Store start/end times to atoms
265 for (unsigned a = 0; a < atoms.size(); ++a) {
266 atoms[a].t1 = atomsTimes[a].first;
267 atoms[a].t2 = atomsTimes[a].second;
268 }
269
270
271 // count the number of events in each atom
272 for (unsigned i = 0; i < evts.size(); ++i) {
273 for (unsigned a = 0; a < atoms.size(); ++a)
274 if (atoms[a].t1 <= evts[i].t && evts[i].t < atoms[a].t2)
275 ++atoms[a].nEv;
276
277 }
278
279
280
281 return atoms;
282 }

◆ getBreaks()

static std::vector< double > getBreaks ( std::vector< std::map< ExpRun, std::pair< double, double > > >  res)
inlinestatic

Get vector with breaks of the calib.

interval

Parameters
resA single calibration interval
Returns
: Vector of times [hours] of the break points, i.e. the subintervals boundaries

Definition at line 131 of file Splitter.h.

132 {
133 std::vector<double> breaks;
134 for (int k = 0; k < int(res.size()) - 1; ++k) {
135 double e = res.at(k).rbegin()->second.second; //end of the previous interval
136 double s = res.at(k + 1).begin()->second.first; //start of the next interval
137 breaks.push_back((e + s) / 2.); //the break time is in between
138 }
139 return breaks;
140 }

◆ getIntervals()

std::vector< std::vector< std::map< ExpRun, std::pair< double, double > > > > getIntervals ( const std::map< ExpRun, std::pair< double, double > > &  runs,
std::vector< Evt >  evts,
TString  lossFunctionOuter,
TString  lossFunctionInner,
double  atomSize = 3. / 60 
)
inline

Function to merge/divide runs into the calibration intervals of given characteristic length.

Parameters
runsA map with key containing expNum+runNum and value start+end time in hours of the run
evtsA vector with all events
lossFunctionOuterA formula of the outer loss function (for calib. intervals)
lossFunctionInnerA formula of the inner loss function (for calib. subintervals)
atomSizelength of the small calibration interval - atom (in hours)
Returns
: Vector of the calib. intervals always containing vector of calib. subintervals. Each subinterval is defined as a map spanning in general over several runs

Definition at line 163 of file Splitter.h.

166 {
167 //sort events by time
168 std::sort(evts.begin(), evts.end(), [](const Evt & e1, const Evt & e2) { return (e1.t > e2.t); });
169
170 // Divide into small intervals
171 std::vector<std::pair<double, double>> smallRuns = splitToSmall(runs, atomSize);
172
173 std::vector<Atom> atoms = createAtoms(smallRuns, evts);
174
175
176 // Calculate breaks for the calib. intervals
177 lossFun = toTF1(lossFunctionOuter);
178 std::vector<int> breaksSize = dynamicBreaks(atoms);
179
180 std::vector<std::vector<std::map<ExpRun, std::pair<double, double>>>> splitsRun; //split intervals with runNumber info
181
182
183 // Calculate breaks for the calib. subintervals
184 lossFun = toTF1(lossFunctionInner);
185 std::vector<int> breaksVtx;
186 for (int i = 0; i < int(breaksSize.size()) + 1; ++i) { //loop over all calib. intervals
187 int s, e;
188 std::tie(s, e) = getStartEndIndexes(atoms.size(), breaksSize, i);
189
190 // Store only atoms lying in the current calib. interval
191 auto currVec = slice(atoms, s, e);
192
193 // Get optimal breaks for particular calib. interval
194 std::vector<int> breaksSmall = dynamicBreaks(currVec);
195
196 // Incorporate the runNumber information
197 auto splitSeg = breaks2intervalsSep(runs, currVec, breaksSmall);
198
199 // Put the vector of subintervals (representing division of the interval)
200 // into the vector
201 splitsRun.push_back(splitSeg);
202
203 if (s != 0) breaksVtx.push_back(s - 1);
204 for (auto b : breaksSmall)
205 breaksVtx.push_back(b + s);
206 if (e != int(breaksSize.size())) breaksVtx.push_back(e);
207 }
208
209 return splitsRun;
210 }
std::vector< Atom > createAtoms(const std::vector< std::pair< double, double > > &atomsTimes, const std::vector< Evt > &evts)
Get the vector with parameters of the calibration Atoms.
Definition: Splitter.h:260
TF1 * toTF1(TString LossString)
Convert lossFunction from string to TF1.
Definition: Splitter.h:286
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.
Definition: Splitter.cc:173
std::vector< Atom > slice(std::vector< Atom > vec, int s, int e)
Slice the vector to contain only elements with indexes s .. e (included)
Definition: Splitter.h:85
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.
Definition: Splitter.cc:280
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:83
std::vector< int > dynamicBreaks(const std::vector< Atom > &runs)
Get optimal break points using algorithm based on dynamic programming.
Definition: Splitter.cc:241

◆ getStartEnd()

static std::pair< double, double > getStartEnd ( std::vector< std::map< ExpRun, std::pair< double, double > > >  res)
inlinestatic

Get the start/end time of the calibration interval (vector of the calib.

subintervals).

Parameters
resA single calibration interval (for example from getIntervals function)
Returns
: A pair with start/end time of the whole calib. interval

Definition at line 120 of file Splitter.h.

121 {
122 return {res.front().begin()->second.first,
123 res.back().rbegin()->second.second};
124 }

◆ toTF1()

TF1 * toTF1 ( TString  LossString)
inlineprivate

Convert lossFunction from string to TF1.

Definition at line 286 of file Splitter.h.

287 {
288 LossString.ReplaceAll("rawTime", "[0]");
289 LossString.ReplaceAll("netTime", "[1]");
290 LossString.ReplaceAll("maxGap", "[2]");
291 LossString.ReplaceAll("nEv", "[3]");
292
293 return (new TF1("lossFun", LossString));
294 }

Member Data Documentation

◆ cache

std::vector<std::pair<double, std::vector<int> > > cache
private

cache used by the clustering algorithm (has to be reset every time)

Definition at line 301 of file Splitter.h.

◆ lossFun

TF1* lossFun
private

loss function used for clustering of Atoms

Definition at line 297 of file Splitter.h.


The documentation for this class was generated from the following files: