Belle II Software development
reconstruction

Modules

 reconstruction data objects
 
 reconstruction modules
 

Classes

class  BeamSpotAlgorithm
 Class implementing BeamSpot calibration algorithm. More...
 
class  BoostVectorAlgorithm
 Class implementing BoostVector calibration algorithm. More...
 
struct  CalibPars
 The parameters related to single calibration interval. More...
 
struct  CalibrationData
 Parameters and data relevant for single calibration interval. More...
 
class  ChebFitter
 Unbinned Maximum Likelihood fitter with a possibility to use Chebyshev interpolation. More...
 
class  InvariantMassAlgorithm
 Class implementing InvariantMass calibration algorithm. More...
 
struct  Atom
 Very small (few mins) calibration interval which cannot be further divided : Atom. More...
 
struct  ExpRun
 Struct containing exp number and run number. More...
 
struct  ExpRunEvt
 struct with expNum, runNum, evtNum More...
 
class  Splitter
 Class that allows to split runs into the intervals of intended properties given by the lossFunction. More...
 
struct  Spline
 Spline structure for zero-order & linear splines. More...
 
class  CDCDedx1DCellAlgorithm
 A calibration algorithm for CDC dE/dx electron: 1D enta cleanup correction. More...
 
class  CDCDedx2DCellAlgorithm
 A calibration algorithm for CDC dE/dx electron 2D enta vs doca correction. More...
 
class  CDCDedxBadWireAlgorithm
 A calibration algorithm for CDC dE/dx to find the bad wires. More...
 
class  CDCDedxCosEdgeAlgorithm
 A calibration algorithm for CDC dE/dx electron cos(theta) dependence. More...
 
class  CDCDedxCosineAlgorithm
 A calibration algorithm for CDC dE/dx electron cos(theta) dependence. More...
 
class  CDCDedxInjectTimeAlgorithm
 A calibration algorithm for CDC dE/dx injection time (HER/LER) More...
 
class  CDCDedxMomentumAlgorithm
 A calibration algorithm for CDC dE/dx electron cos(theta) dependence. More...
 
class  CDCDedxRunGainAlgorithm
 A calibration algorithm for CDC dE/dx run gains. More...
 
class  CDCDedxWireGainAlgorithm
 A calibration algorithm for CDC dE/dx wire gains. More...
 
class  CDCDedx1DCell
 dE/dx wire gain calibration constants More...
 
class  CDCDedx2DCell
 dE/dx wire gain calibration constants More...
 
class  CDCDedxADCNonLinearity
 dE/dx electronic ADC non-linearity correction for highly ionising particles (used in offline hadron saturation calibration) parameters are for X vs Y relation and sep for inner and outer layer vector array 0,1 for inner and 2,3 for outer layers More...
 
class  CDCDedxBadWires
 dE/dx wire gain calibration constants More...
 
class  CDCDedxCosineCor
 dE/dx wire gain calibration constants More...
 
class  CDCDedxCosineEdge
 dE/dx special large cosine calibration to fix bending shoulder at large costh More...
 
class  CDCDedxDatabaseImporter
 dE/dx database importer. More...
 
class  CDCDedxHadronCor
 dE/dx hadron saturation parameterization constants More...
 
class  CDCDedxInjectionTime
 dE/dx injection time calibration constants More...
 
class  CDCDedxMeanPars
 dE/dx mean (curve versus beta-gamma) parameterization constants More...
 
class  CDCDedxMomentumCor
 dE/dx wire gain calibration constants More...
 
class  CDCDedxRunGain
 dE/dx run gain calibration constants More...
 
class  CDCDedxScaleFactor
 dE/dx run gain calibration constants More...
 
class  CDCDedxSigmaPars
 dE/dx sigma (versus beta-gamma) parameterization constants More...
 
class  CDCDedxWireGain
 dE/dx wire gain calibration constants More...
 
class  DedxPDFs
 dE/dx wire gain calibration constants More...
 
class  EventsOfDoomParameters
 DBObject containing parameters used in EventsOfDoomBuster module. More...
 
class  CDCDedxMeanPred
 Class to hold the prediction of mean as a function of beta-gamma (bg) More...
 
class  CDCDedxSigmaPred
 Class to hold the prediction of resolution depending dE/dx, nhit, and cos(theta) More...
 

Typedefs

typedef std::map< std::string, double > Pars
 values of parameters in ML fit
 
typedef std::map< std::string, std::pair< double, double > > Limits
 limits of parameters in ML fit
 

Functions

TMatrixDSym toTMatrixDSym (Eigen::MatrixXd mIn)
 Function that converts Eigen symmetric matrix to ROOT matrix.
 
B2Vector3D toB2Vector3 (Eigen::VectorXd vIn)
 Function that converts Eigen vector to ROOT vector.
 
int getID (const std::vector< double > &breaks, double t)
 get id of the time point t
 
void extrapolateCalibration (std::vector< CalibrationData > &calVec)
 Extrapolate calibration to intervals where it failed.
 
void addShortRun (std::vector< CalibrationData > &calVec, std::pair< ExpRun, std::pair< double, double > > shortRun)
 Extrapolate calibration to the very short runs which were filtered before.
 
double encodeNumber (double val, unsigned num)
 Encode integer num into double val such that val is nearly not changed (maximally by a relative shift 1e-6).
 
unsigned decodeNumber (double val)
 Decode the integer number encoded in val.
 
template<typename Evt >
void storePayloads (const std::vector< Evt > &evts, const std::vector< CalibrationData > &calVecConst, std::string objName, std::function< TObject *(Eigen::VectorXd, Eigen::MatrixXd, Eigen::MatrixXd) > getCalibObj)
 Store payloads to files.
 
void storePayloadsNoIntraRun (const std::vector< CalibrationData > &calVecConst, std::string objName, std::function< TObject *(Eigen::VectorXd, Eigen::MatrixXd, Eigen::MatrixXd) > getCalibObj)
 Store payloads to files, where calib data have no intra-run dependence.
 
template<typename Evt , typename Fun >
CalibrationData runAlgorithm (const std::vector< Evt > &evts, std::vector< std::map< ExpRun, std::pair< double, double > > > range, Fun runCalibAnalysis)
 run calibration algorithm for single calibration interval
 
template<typename Fun1 , typename Fun2 >
CalibrationAlgorithm::EResult runCalibration (TTree *tracks, const std::string &calibName, Fun1 GetEvents, Fun2 calibAnalysis, std::function< TObject *(Eigen::VectorXd, Eigen::MatrixXd, Eigen::MatrixXd)> calibObjCreator, TString m_lossFunctionOuter, TString m_lossFunctionInner)
 Run the the calibration over the whole event sample.
 
std::pair< double, double > getMinima (std::vector< std::vector< double > > vals, int i0, int j0)
 Get minimum inside and outside of the smaller window defined by i0, j0.
 
std::vector< double > getMinimum (std::function< double(double, double)> fun, double xMin, double xMax, double yMin, double yMax)
 Get minimum of 2D function in the rectangular domain defined by xMin,xMax & yMin,yMax.
 
Eigen::VectorXd getWeights (int Size)
 Get the vector of weights to calculate the integral over the Chebyshev nodes The nodes are by definition between 0 and 1, there are Size nodes To get their positions, use getNodes.
 
Eigen::VectorXd getNodes (int Size)
 Get the vector of positions of the Chebyshev nodes The nodes are by definition between 0 and 1, there are Size nodes For the corresponding weights use getWeights.
 
Eigen::VectorXd getPols (int Size, double x)
 Evaluate Chebyshev polynomials up to Size at point x It returns a vector of the P_i(x) for i=0..Size-1 The polynomial is defined for x between 0 and 1.
 
Eigen::VectorXd getPolsSum (int Size, Eigen::VectorXd x)
 Calculate the Chebyshev polynomials of order i=0..Size-1 at points given in vector x_j and sum it over point index j It returns sum_j P_i(x_j) for i=0..Size-1 The Chebyshev polynomials are defined for x between 0 and 1.
 
Eigen::MatrixXd getCoefs (int Size, bool isInverse=false)
 Transformation matrix between Cheb nodes and coefficients of the Cheb polynomials Notice, that there are two alternative ways defining polynomial interpolation:
 
double evalPol (const Eigen::VectorXd &polCoef, double x)
 Evaluate Cheb.
 
Eigen::MatrixXd getCoefsCheb (int oldSize)
 Transformation matrix between Cheb nodes and Cheb coefficients with better normalization of the borders.
 
Eigen::VectorXd interpol (const Eigen::VectorXd &xi, double x)
 Get Interpolation vector k_i for point x from the function values at points xi (polynomial interpolation) In the second step, the function value at x can be evaluated as sum_i vals_i k_i.
 
double interpol (Eigen::VectorXd xi, Eigen::VectorXd vals, double x)
 Get interpolated function value at point x when function values vals at points xi are provided.
 
bool operator!= (ExpRun a, ExpRun b)
 Not equal for ExpRun.
 
bool operator< (ExpRun a, ExpRun b)
 less than for ExpRun
 
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::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::vector< Atomslice (std::vector< Atom > vec, int s, int e)
 Slice the vector to contain only elements with indexes s .. e (included)
 
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.
 
template<typename Evt >
std::map< ExpRun, std::pair< double, double > > getRunInfo (const std::vector< Evt > &evts)
 Get the map of runs, where each run contains pair with start/end time [hours].
 
template<typename Evt >
ExpRunEvt getPosition (const std::vector< Evt > &events, double tEdge)
 Get the exp-run-evt number from the event time [hours].
 
template<typename Evt >
std::vector< ExpRunEvtconvertSplitPoints (const std::vector< Evt > &events, std::vector< double > splitPoints)
 Convert splitPoints [hours] to breakPoints in ExpRunEvt.
 
TString rn ()
 Get random string.
 
std::vector< std::vector< double > > merge (std::vector< std::vector< std::vector< double > > > toMerge)
 merge { vector<double> a, vector<double> b} into {a, b}
 
Eigen::VectorXd vec2vec (std::vector< double > vec)
 std vector -> ROOT vector
 
std::vector< double > vec2vec (Eigen::VectorXd v)
 ROOT vector -> std vector.
 
Eigen::MatrixXd vecs2mat (std::vector< std::vector< double > > vecs)
 merge columns (from std::vectors) into ROOT matrix
 
std::vector< double > getRangeLin (int nVals, double xMin, double xMax)
 Equidistant range between xMin and xMax for spline of the first order.
 
std::vector< double > getRangeZero (int nVals, double xMin, double xMax)
 Equidistant range between xMin and xMax for spline of the zero order.
 
std::vector< double > slice (std::vector< double > v, unsigned ind, unsigned n)
 put slice of original vector v[ind:ind+n] into new one, n is number of elements
 
double eval (const std::vector< double > &spl, const std::vector< double > &vals, double x)
 Evaluate spline (zero order or first order) in point x.
 
VectorXd getPolsSum (int Size, VectorXd x)
 Calculate the Chebyshev polynomials of order i=0..Size-1 at points given in vector x_j and sum it over point index j It returns sum_j P_i(x_j) for i=0..Size-1 The Chebyshev polynomials are defined for x between 0 and 1.
 
double evalPol (const VectorXd &polCoef, double x)
 Evaluate Cheb.
 
VectorXd interpol (const VectorXd &xi, double x)
 Get Interpolation vector k_i for point x from the function values at points xi (polynomial interpolation) In the second step, the function value at x can be evaluated as sum_i vals_i k_i.
 
double interpol (VectorXd xi, VectorXd vals, double x)
 Get interpolated function value at point x when function values vals at points xi are provided.
 
void plotRuns (std::vector< std::pair< double, double > > runs)
 plot runs on time axis
 
void plotSRuns (std::vector< std::pair< double, double > > runs, std::vector< int > breaks, int offset=2)
 plot clusters or runs on time axis
 
void printBySize (std::vector< std::pair< double, double > > runs)
 print sorted lengths of the runs
 
static ExpRun getRun (std::map< ExpRun, std::pair< double, double > > runs, double t)
 Get exp number + run number from time.
 
CDCDedxTrack const * getDedxFromParticle (Particle const *particle)
 CDC dEdx value from particle.
 
VXDDedxTrack const * getSVDDedxFromParticle (Particle const *particle)
 SVD dEdx value from particle.
 
Eigen::VectorXd getLogFunction (Pars pars) const
 Get the -2*log(p(x)) on the Cheb nodes.
 
void init (int Size, double xMin, double xMax)
 Initialize the fitter (the Chebyshev coefficients)
 
double getLogLikelihoodSlow (const Pars &pars) const
 Calculate log likelihood using exact formula.
 
double getLogLikelihoodFast (const Pars &pars) const
 Calculate log likelihood using approximation based on Chebyshev polynomials (typically faster)
 
double operator() (const double *par) const
 Evaluate the log likelihood.
 
Eigen::VectorXd getDataGrid () const
 Calculate Chebyshev coefficients for the data set.
 
std::pair< Eigen::VectorXd, Eigen::MatrixXd > getDataGridWithCov () const
 Calculate Chebyshev coefficients with the covariance (useful for toy studies)
 
std::pair< Pars, Eigen::MatrixXd > fitData (Pars pars, Limits limits, bool UseCheb=true)
 Fit the data with specified initial values of parameters and limits on them.
 
double getFunctionFast (const Pars &pars, double x)
 Evaluate the fitted function approximated with the Chebyshev polynomial, typically runs faster.
 
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
 
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.
 
std::vector< int > dynamicBreaks (const std::vector< Atom > &runs)
 Get optimal break points using algorithm based on dynamic programming.
 
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.
 

Detailed Description

Typedef Documentation

◆ Limits

typedef std::map<std::string, std::pair<double, double> > Limits

limits of parameters in ML fit

Definition at line 28 of file ChebFitter.h.

◆ Pars

typedef std::map<std::string, double> Pars

values of parameters in ML fit

Definition at line 25 of file ChebFitter.h.

Function Documentation

◆ addShortRun()

void addShortRun ( std::vector< CalibrationData > &  calVec,
std::pair< ExpRun, std::pair< double, double > >  shortRun 
)
inline

Extrapolate calibration to the very short runs which were filtered before.

Definition at line 151 of file calibTools.h.

152 {
153 double shortStart = shortRun.second.first;
154 double shortEnd = shortRun.second.second;
155
156 double distMin = 1e20;
157 int iMin = -1, jMin = -1;
158
159 for (unsigned i = 0; i < calVec.size(); ++i) {
160 if (calVec[i].isCalibrated == false)
161 continue;
162 for (unsigned j = 0; j < calVec[i].subIntervals.size(); ++j) {
163 for (auto I : calVec[i].subIntervals[j]) {
164 double s = I.second.first;
165 double e = I.second.second;
166
167 double dist1 = (s - shortEnd >= 0) ? (s - shortEnd) : 1e20;
168 double dist2 = (shortStart - e >= 0) ? (shortStart - e) : 1e20;
169 double dist = std::min(dist1, dist2);
170
171 if (dist < distMin) {
172 distMin = dist;
173 iMin = i;
174 jMin = j;
175 }
176 }
177 }
178 }
179
180 B2ASSERT("Must be found", iMin != -1 && jMin != -1);
181 calVec[iMin].subIntervals[jMin].insert(shortRun);
182 }

◆ breaks2intervalsSep()

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.

Parameters
runsMapmap defining the time ranges of the runs
currVecvector with time intervals of the atoms (small non-divisible time intervals)
breaksvector with integer indexes of the breaks
Returns
: a vector of the calib. intervals where each interval is a map with exp-run as a key and start- end-time as a value

Definition at line 280 of file Splitter.cc.

283 {
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) {
286 int s, e;
287 std::tie(s, e) = getStartEndIndexes(currVec.size(), breaks, i);
288
289 // loop over atoms in single calib interval
290 for (int k = s; k <= e; ++k) {
291 ExpRun r = getRun(runsMap, (currVec[k].t1 + currVec[k].t2) / 2.); //runexp of the atom
292 if (splitsNow[i].count(r)) { //if already there
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);
295 } else { //if new
296 splitsNow[i][r].first = currVec[k].t1;
297 splitsNow[i][r].second = currVec[k].t2;
298 }
299 }
300 }
301
302 return splitsNow;
303 }

◆ convertSplitPoints()

std::vector< ExpRunEvt > convertSplitPoints ( const std::vector< Evt > &  events,
std::vector< double >  splitPoints 
)

Convert splitPoints [hours] to breakPoints in ExpRunEvt.

Parameters
eventsvector of events
splitPointsthe vector containing times of the edges of the calibration intervals [hours]
Returns
a vector with calibration break-points in the exp-run-evt format

Definition at line 363 of file Splitter.h.

364 {
365
366 std::vector<ExpRunEvt> breakPos;
367 for (auto p : splitPoints) {
368 auto pos = getPosition(events, p);
369 breakPos.push_back(pos);
370 }
371 return breakPos;
372 }

◆ decodeNumber()

unsigned decodeNumber ( double  val)
inline

Decode the integer number encoded in val.

Definition at line 208 of file calibTools.h.

209 {
210 double factor = pow(FLT_RADIX, DBL_MANT_DIG);
211 static const long long fEnc = pow(2, 32); //32 binary digits for encoded number
212
213 int e;
214 double mantisa = std::frexp(val, &e);
215 long long mantisaI = mantisa * factor;
216
217 return (mantisaI % fEnc);
218 }

◆ dynamicBreaks()

std::vector< int > dynamicBreaks ( const std::vector< Atom > &  runs)
private

Get optimal break points using algorithm based on dynamic programming.

Parameters
runsVector of atoms, where each atom is an intervals in time
Returns
: Optimal indexes of the break points

Definition at line 241 of file Splitter.cc.

242 {
243 //reset cache
244 cache.resize(runs.size());
245 for (auto& c : cache)
246 c = std::make_pair(-1.0, std::vector<int>({}));
247
248
249 std::vector<int> breaks;
250 getMinLoss(runs, runs.size() - 1, breaks); //the minLoss (output) currently not used, only breaks
251
252 return breaks;
253 }
std::vector< std::pair< double, std::vector< int > > > cache
cache used by the clustering algorithm (has to be reset every time)
Definition: Splitter.h:301
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.
Definition: Splitter.cc:204

◆ encodeNumber()

double encodeNumber ( double  val,
unsigned  num 
)
inline

Encode integer num into double val such that val is nearly not changed (maximally by a relative shift 1e-6).

It is use to store time information to the payloads

Definition at line 186 of file calibTools.h.

187 {
188 double factor = pow(FLT_RADIX, DBL_MANT_DIG);
189 static const long long fEnc = pow(2, 32); //32 binary digits for encoded number
190
191 int e; //exponent of the number
192 double mantisa = std::frexp(val, &e);
193 long long mantisaI = mantisa * factor; //mantissa as integer
194
195 if (val != 0)
196 mantisaI = (mantisaI / fEnc) * fEnc + num; //adding encoded number to last digits of mantissa
197 else {
198 mantisaI = factor / 2 + num;
199 e = -100; //if the val is zero, ensure very small number by the exponent
200 }
201
202 double newVal = ldexp(mantisaI / factor, e);
203
204 return newVal;
205 }

◆ eval()

double eval ( const std::vector< double > &  spl,
const std::vector< double > &  vals,
double  x 
)
inline

Evaluate spline (zero order or first order) in point x.

Definition at line 115 of file tools.h.

116 {
117 int order = -1;
118 if (spl.size() == 0)
119 order = 0;
120 else if (spl.size() == vals.size() - 1)
121 order = 0;
122 else if (spl.size() == vals.size())
123 order = 1;
124 else {
125 B2FATAL("Unknown order of spline");
126 }
127 B2ASSERT("Spline order should be zero or one", order == 0 || order == 1);
128
129 if (order == 1) {
130 B2ASSERT("Linear spline only meaningful for two or more nodes", spl.size() >= 2);
131 B2ASSERT("As nodes as values in lin. spline", spl.size() == vals.size());
132
133 if (x <= spl[0]) return vals[0];
134 if (x >= spl.back()) return vals.back();
135
136 // binary search for position
137 int i1 = lower_bound(spl.begin(), spl.end(), x) - spl.begin() - 1;
138
139 if (!(spl[i1] <= x && x <= spl[i1 + 1])) {
140 B2FATAL("Wrong place founded : " << spl[i1] << " " << x << " " << spl[i1 + 1]);
141 }
142
143 // Linear interpolation between neighbouring nodes
144 double v = ((spl[i1 + 1] - x) * vals[i1] + (x - spl[i1]) * vals[i1 + 1]) / (spl[i1 + 1] - spl[i1]);
145 return v;
146 } else if (order == 0) { //zero order polynomial
147 B2ASSERT("#values vs #nodes in zero-order spline", spl.size() + 1 == vals.size());
148 if (vals.size() == 1) {
149 return vals[0];
150 } else {
151 double res = vals[0]; //by default value from lowest node
152 for (unsigned i = 0; i < spl.size(); ++i) {
153 if (spl[i] <= x) res = vals[i + 1];
154 else break;
155 }
156 return res;
157 }
158 }
159 return -99;
160 }

◆ evalPol() [1/2]

double evalPol ( const Eigen::VectorXd &  polCoef,
double  x 
)

Evaluate Cheb.

pol at point x when the coefficients of the expansion are provided

◆ evalPol() [2/2]

double evalPol ( const VectorXd &  polCoef,
double  x 
)

Evaluate Cheb.

pol at point x when the coefficients of the expansion are provided

Definition at line 176 of file nodes.cc.

177 {
178 VectorXd pols = getPols(polCoef.size(), x);
179
180 double s = pols.dot(polCoef);
181
182 return s;
183 }

◆ extrapolateCalibration()

void extrapolateCalibration ( std::vector< CalibrationData > &  calVec)
inline

Extrapolate calibration to intervals where it failed.

Definition at line 102 of file calibTools.h.

103 {
104 //put closest neighbor, where the statistic was low or algo failed
105 for (unsigned i = 0; i < calVec.size(); ++i) {
106 if (calVec[i].pars.cnt.size() != 0) continue;
107 const auto& r = calVec[i].subIntervals;
108 double Start, End;
109 std::tie(Start, End) = Splitter::getStartEnd(r);
110
111 Eigen::Vector3d ipNow;
112 Eigen::MatrixXd ipeNow;
113 Eigen::MatrixXd sizeMatNow;
114
115 double distMin = 1e20;
116 //Find the closest calibrated interval
117 for (unsigned j = 0; j < calVec.size(); ++j) {
118 if (calVec[j].isCalibrated == false) continue; //skip not-calibrated intervals
119 const auto& rJ = calVec[j].subIntervals;
120 for (unsigned jj = 0; jj < rJ.size(); ++jj) { //loop over subintervals
121 const auto& rNow = rJ[jj];
122 double s = rNow.begin()->second.first;
123 double e = rNow.rbegin()->second.second;
124
125 double dist1 = (s - End >= 0) ? (s - End) : 1e20;
126 double dist2 = (Start - e >= 0) ? (Start - e) : 1e20;
127 double dist = std::min(dist1, dist2);
128
129 if (dist < distMin) {
130 ipNow = calVec[j].pars.cnt.at(jj);
131 ipeNow = calVec[j].pars.cntUnc.at(jj);
132 sizeMatNow = calVec[j].pars.spreadMat;
133 distMin = dist;
134 }
135 }
136 }
137
138 //Store it to vectors
139 calVec[i].pars.cnt.resize(r.size());
140 calVec[i].pars.cntUnc.resize(r.size());
141 for (unsigned ii = 0; ii < r.size(); ++ii) {
142 calVec[i].pars.cnt.at(ii) = ipNow;
143 calVec[i].pars.cntUnc.at(ii) = ipeNow;
144 }
145 calVec[i].pars.spreadMat = sizeMatNow;
146 }
147
148 }

◆ filter()

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

Definition at line 38 of file Splitter.cc.

40 {
41 std::map<ExpRun, std::pair<double, double>> runsCopy;
42
43 for (auto r : runs) {
44 auto& I = r.second;
45 double d = I.second - I.first;
46 if (d > cut)
47 runsCopy[r.first] = r.second;
48 else
49 runsRemoved[r.first] = r.second;
50 }
51
52 return runsCopy;
53 }

◆ fitData()

std::pair< Pars, MatrixXd > fitData ( Pars  pars,
Limits  limits,
bool  UseCheb = true 
)

Fit the data with specified initial values of parameters and limits on them.

Definition at line 155 of file ChebFitter.cc.

156 {
157 m_useCheb = UseCheb;
158
159 ROOT::Math::Minimizer* minimum =
160 ROOT::Math::Factory::CreateMinimizer("Minuit2", "");
161
162 // set tolerance , etc...
163 minimum->SetMaxFunctionCalls(10000000); // for Minuit/Minuit2
164 minimum->SetMaxIterations(100000); // for GSL
165 minimum->SetTolerance(10.0);
166 //minimum->SetPrecision(1e-5);
167
168 //minimum->SetPrintLevel(3); //many outputs
169 minimum->SetPrintLevel(0); //few outputs
170 minimum->SetStrategy(2);
171 minimum->SetErrorDef(1);
172
173
174 // Set the free variables to be minimized !
175 m_parNames.clear();
176 int k = 0;
177 for (auto p : pars) {
178 std::string n = p.first;
179 double vCnt = p.second;
180 if (limits.count(n) == 1) {
181 double vMin = limits.at(n).first;
182 double vMax = limits.at(n).second;
183 double step = (vMax - vMin) / 100;
184 minimum->SetLimitedVariable(k, n, vCnt, step, vMin, vMax);
185 } else {
186 double step = 1;
187 minimum->SetVariable(k, n, vCnt, step);
188 }
189 m_parNames.push_back(n);
190 ++k;
191 }
192
193 // create function wrapper for minimizer
194 // a IMultiGenFunction type
195 ROOT::Math::Functor f(*this, pars.size());
196 minimum->SetFunction(f);
197
198
199 // do the minimization
200 minimum->Minimize();
201
202
203 Pars parsF;
204 for (unsigned i = 0; i < m_parNames.size(); ++i)
205 parsF[m_parNames[i]] = minimum->X()[i];
206
207 MatrixXd covMat(parsF.size(), parsF.size());
208 for (unsigned i = 0; i < parsF.size(); ++i)
209 for (unsigned j = 0; j < parsF.size(); ++j)
210 covMat(i, j) = minimum->CovMatrix(i, j);
211
212 // print pars
213 std::stringstream log;
214 log << "Minuit status : " << minimum->Status() << ", ";
215 for (auto p : parsF)
216 log << "\"" << p.first << "\" : " << p.second << ", ";
217
218 B2INFO(log.str());
219
220 delete minimum;
221
222 return std::make_pair(parsF, covMat);
223 }
bool m_useCheb
flag to use approximation based on Chebyshev polynomials
Definition: ChebFitter.h:87
std::vector< std::string > m_parNames
vector with names of the parameters
Definition: ChebFitter.h:86
std::map< std::string, double > Pars
values of parameters in ML fit
Definition: ChebFitter.h:25

◆ getCoefs()

MatrixXd getCoefs ( int  Size,
bool  isInverse = false 
)

Transformation matrix between Cheb nodes and coefficients of the Cheb polynomials Notice, that there are two alternative ways defining polynomial interpolation:

  • coefficients c_i in of the Cheb polynomials, i.e. f(x) = sum_i c_i P_i(x)
  • Values of the f(x) in the Cheb nodes, i.e. d_j = f(x_j), where x_j are the nodes The Chebyshev polynomials are defined for x between 0 and 1

Definition at line 127 of file nodes.cc.

128 {
129 const int N = Size - 1;
130 assert(N % 2 == 0);
131
132 MatrixXd coef(Size, Size);
133
134 double mul = 1;
135 double C = 1. / N;
136 if (isInverse == true) {C = 1. / 2; }
137
138 for (int k = 0; k <= N; ++k) {
139 if (!isInverse) {
140 coef(k, N) = C;
141 coef(k, 0) = C * (k % 2 == 1 ? -1 : 1);
142 } else {
143 mul = k % 2 == 1 ? -1 : 1;
144 coef(N - k, N) = C * mul;
145 coef(N - k, 0) = C ;
146 }
147
148 for (int n = 1; n <= N - 1; ++n) {
149 double el = cos(n * k * M_PI / N) * 2.*C * mul;
150 if (!isInverse) coef(k, N - n) = el;
151 else coef(N - k, N - n) = el;
152 }
153 }
154
155 return coef;
156 }

◆ getCoefsCheb()

MatrixXd getCoefsCheb ( int  oldSize)

Transformation matrix between Cheb nodes and Cheb coefficients with better normalization of the borders.

Definition at line 162 of file nodes.cc.

163 {
164 auto coef = getCoefs(oldSize);
165
166 coef.row(0) *= 0.5;
167 coef.row(coef.rows() - 1) *= 0.5;
168
169 return coef;
170 }

◆ getDataGrid()

VectorXd getDataGrid ( ) const
private

Calculate Chebyshev coefficients for the data set.

Definition at line 109 of file ChebFitter.cc.

110 {
111 double a = m_nodes[0];
112 double b = m_nodes[m_nodes.size() - 1];
113
114
115 VectorXd polSum = VectorXd::Zero(m_nodes.size());
116 for (double x : m_data) {
117 double xx = (x - a) / (b - a); //normalize between 0 and 1
118 polSum += getPols(m_nodes.size(), xx);
119 }
120
121
122 //transform to the basis of the cheb m_nodes
123 VectorXd gridVals = m_coefsMat * polSum;
124
125 return gridVals;
126 }
std::vector< double > m_data
vector with the data points to be fitted
Definition: ChebFitter.h:78
Eigen::VectorXd m_nodes
vector with cheb nodes
Definition: ChebFitter.h:83
Eigen::MatrixXd m_coefsMat
transformation matrix from chebPol to gridPoints
Definition: ChebFitter.h:81
Eigen::VectorXd getPols(int Size, double x)
Evaluate Chebyshev polynomials up to Size at point x It returns a vector of the P_i(x) for i=0....
Definition: nodes.cc:77

◆ getDataGridWithCov()

std::pair< VectorXd, MatrixXd > getDataGridWithCov ( ) const
private

Calculate Chebyshev coefficients with the covariance (useful for toy studies)

Definition at line 130 of file ChebFitter.cc.

131 {
132 double a = m_nodes[0];
133 double b = m_nodes[m_nodes.size() - 1];
134
135 MatrixXd polSum2 = MatrixXd::Zero(m_nodes.size(), m_nodes.size());
136 VectorXd polSum = VectorXd::Zero(m_nodes.size());
137 for (double x : m_data) {
138 double xx = (x - a) / (b - a); //normalize between 0 and 1
139 VectorXd pol = getPols(m_nodes.size(), xx);
140 polSum += pol;
141 polSum2 += pol * pol.transpose();
142 }
143
144
145 //transform to the basis of the cheb nodes
146 MatrixXd coefs = getCoefsCheb(polSum.size()).transpose();
147 VectorXd gridVals = coefs * polSum;
148 MatrixXd gridValsCov = coefs * polSum2 * coefs.transpose();
149
150 return std::make_pair(gridVals, gridValsCov);
151 }
Eigen::MatrixXd getCoefsCheb(int oldSize)
Transformation matrix between Cheb nodes and Cheb coefficients with better normalization of the borde...
Definition: nodes.cc:162

◆ getDedxFromParticle()

CDCDedxTrack const * getDedxFromParticle ( Particle const *  particle)

CDC dEdx value from particle.

Definition at line 37 of file DedxVariables.cc.

38 {
39 const Track* track = particle->getTrack();
40 if (!track) {
41 return nullptr;
42 }
43
44 const CDCDedxTrack* dedxTrack = track->getRelatedTo<CDCDedxTrack>();
45 if (!dedxTrack) {
46 return nullptr;
47 }
48
49 return dedxTrack;
50 }

◆ getFunctionFast()

double getFunctionFast ( const Pars pars,
double  x 
)
private

Evaluate the fitted function approximated with the Chebyshev polynomial, typically runs faster.

Definition at line 226 of file ChebFitter.cc.

227 {
228 static VectorXd funVals = getLogFunction(pars);
229
230
231 return exp(-0.5 * interpol(m_nodes, funVals, x));
232
233 }
Eigen::VectorXd getLogFunction(Pars pars) const
Get the -2*log(p(x)) on the Cheb nodes.
Definition: ChebFitter.cc:43
Eigen::VectorXd interpol(const Eigen::VectorXd &xi, double x)
Get Interpolation vector k_i for point x from the function values at points xi (polynomial interpolat...

◆ getID()

int getID ( const std::vector< double > &  breaks,
double  t 
)
inline

get id of the time point t

Definition at line 60 of file calibTools.h.

61 {
62 for (int i = 0; i < int(breaks.size()) + 1; ++i) {
63 double s = (i == 0) ? 0 : breaks[i - 1];
64 double e = (i == int(breaks.size())) ? 1e20 : breaks[i];
65 if (s <= t && t < e)
66 return i;
67 }
68 return -1;
69 }

◆ getLogFunction()

VectorXd getLogFunction ( Pars  pars) const
private

Get the -2*log(p(x)) on the Cheb nodes.

Definition at line 43 of file ChebFitter.cc.

44 {
45 VectorXd fVals(m_nodes.size());
46
47 //calc function values
48 fVals = m_nodes.unaryExpr([&](double x) { return m_myFun(x, pars); });
49
50 //normalize the function
51 double I = fVals.dot(m_weights);
52
53 fVals = -2 * log(fVals.array() / I); //normalize by integral
54
55 return fVals;
56
57 }
Eigen::VectorXd m_weights
vector with cheb weights for integration
Definition: ChebFitter.h:84
std::function< double(double, Pars)> m_myFun
function to fit
Definition: ChebFitter.h:89

◆ getLogLikelihoodFast()

double getLogLikelihoodFast ( const Pars pars) const
private

Calculate log likelihood using approximation based on Chebyshev polynomials (typically faster)

Definition at line 89 of file ChebFitter.cc.

90 {
91 VectorXd funVals = getLogFunction(pars);
92 double LL = funVals.dot(m_dataGrid);
93
94 return LL;
95 }
Eigen::VectorXd m_dataGrid
vector with the data points related to the cheb nodes (m_dataGrid.size = nodes.size)
Definition: ChebFitter.h:79

◆ getLogLikelihoodSlow()

double getLogLikelihoodSlow ( const Pars pars) const
private

Calculate log likelihood using exact formula.

Definition at line 76 of file ChebFitter.cc.

77 {
78
79 double L = 0;
80 for (double d : m_data) {
81 double v = m_myFun(d, pars);
82 L += -2 * log(v);
83 }
84
85 return L;
86 }

◆ getMinima()

std::pair< double, double > getMinima ( std::vector< std::vector< double > >  vals,
int  i0,
int  j0 
)
inline

Get minimum inside and outside of the smaller window defined by i0, j0.

Definition at line 23 of file minimizer.h.

24 {
25 int N = vals.size();
26 int Nzoom = (N + 1) / 2;
27
28
29 double minIn = std::numeric_limits<double>::max();
30 double minOut = std::numeric_limits<double>::max();
31
32 for (int i = 0; i < N; ++i)
33 for (int j = 0; j < N; ++j) {
34 if ((i0 <= i && i < i0 + Nzoom) &&
35 (j0 <= j && j < j0 + Nzoom))
36 minIn = std::min(minIn, vals[i][j]);
37 else
38 minOut = std::min(minOut, vals[i][j]);
39 }
40 return {minIn, minOut};
41 }

◆ getMinimum()

std::vector< double > getMinimum ( std::function< double(double, double)>  fun,
double  xMin,
double  xMax,
double  yMin,
double  yMax 
)
inline

Get minimum of 2D function in the rectangular domain defined by xMin,xMax & yMin,yMax.

Definition at line 44 of file minimizer.h.

45 {
46 const int N = 17; //the grid has size N x N
47 const int Nzoom = (N + 1) / 2; // the size of the zoomed rectangle where minimum is
48
49
50 const int kMax = 35; //number of iterations
51
52 std::vector<std::vector<double>> vals(N);
53 for (auto& v : vals) v.resize(N);
54
55 for (int k = 0; k < kMax; ++k) {
56
57 // get values of the function
58 for (int i = 0; i < N; ++i)
59 for (int j = 0; j < N; ++j) {
60 double x = xMin + i * (xMax - xMin) / (N - 1.);
61 double y = yMin + j * (yMax - yMin) / (N - 1.);
62 vals[i][j] = fun(x, y);
63 }
64
65 if (k == kMax - 1) break;
66
67 double mOutMax = - std::numeric_limits<double>::max();
68 int iOpt = -1, jOpt = -1;
69 //find optimal rectangle
70 for (int i = 0; i < N - Nzoom; ++i)
71 for (int j = 0; j < N - Nzoom; ++j) {
72 double mIn, mOut;
73 std::tie(mIn, mOut) = getMinima(vals, i, j);
74
75 if (mOut > mOutMax) {
76 mOutMax = mOut;
77 iOpt = i;
78 jOpt = j;
79 }
80 }
81
82 //Zoom to the optimal rectangle
83 // get values of the function
84
85 double xMinNow = xMin + iOpt * (xMax - xMin) / (N - 1.);
86 double xMaxNow = xMin + (iOpt + Nzoom - 1) * (xMax - xMin) / (N - 1.);
87
88 double yMinNow = yMin + jOpt * (yMax - yMin) / (N - 1.);
89 double yMaxNow = yMin + (jOpt + Nzoom - 1) * (yMax - yMin) / (N - 1.);
90
91 xMin = xMinNow;
92 xMax = xMaxNow;
93 yMin = yMinNow;
94 yMax = yMaxNow;
95
96 }
97
98
99 //get the overall minimum
100 double minTot = std::numeric_limits<double>::max();
101 int iOpt = -1, jOpt = -1;
102
103 for (int i = 0; i < N; ++i)
104 for (int j = 0; j < N; ++j) {
105 if (vals[i][j] < minTot) {
106 minTot = vals[i][j];
107 iOpt = i;
108 jOpt = j;
109 }
110 }
111
112 double xMinNow = xMin + iOpt * (xMax - xMin) / (N - 1.);
113 double yMinNow = yMin + jOpt * (yMax - yMin) / (N - 1.);
114
115 return {xMinNow, yMinNow};
116 }

◆ getMinLoss()

double getMinLoss ( const std::vector< Atom > &  vec,
int  e,
std::vector< int > &  breaks 
)
private

Recursive function to evaluate minimal sum of the lossFuctions for the optimal clustering.

It return the minimum of the lossFunction and optimal break points giving such value. It acts only on atoms with indexes between 0 and e

Parameters
vecVector of atoms, where each atom is an intervals in time
ethe index of the last atom included in the optimisation problem
[out]breaksOutput vector with indexes of optimal break points
Returns
: Minimal value of the summed loss function

Definition at line 204 of file Splitter.cc.

205 {
206 // If entry in cache (speed up)
207 if (cache[e].first >= 0) {
208 breaks = cache[e].second;
209 return cache[e].first;
210 }
211
212
213 std::vector<int> breaksOpt;
214 double minVal = 1e30;
215 int iMin = -10;
216 for (int i = -1; i <= e - 1; ++i) {
217 auto breaksNow = breaks;
218 double r1 = 0;
219 if (i != -1)
220 r1 = getMinLoss(vec, i, breaksNow);
221 double r2 = lossFunction(vec, i + 1, e);
222 double tot = r1 + r2;
223
224 if (tot < minVal) { //store minimum
225 minVal = tot;
226 iMin = i;
227 breaksOpt = breaksNow;
228 }
229 }
230
231 if (iMin != -1)
232 breaksOpt.push_back(iMin);
233
234
235 breaks = breaksOpt;
236 cache[e] = std::make_pair(minVal, breaks); //store solution to cache
237 return minVal;
238 }
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...
Definition: Splitter.cc:135

◆ getNodes()

VectorXd getNodes ( int  Size)

Get the vector of positions of the Chebyshev nodes The nodes are by definition between 0 and 1, there are Size nodes For the corresponding weights use getWeights.

Definition at line 65 of file nodes.cc.

66 {
67 assert((Size - 1) % 2 == 0);
68 VectorXd xi = VectorXd::Zero(Size);
69 for (int i = 0; i < Size; ++i) {
70 double Cos = cos(i / (Size - 1.) * M_PI);
71 xi[i] = (1 - Cos) / 2;
72 }
73 return xi;
74 }

◆ getPols()

VectorXd getPols ( int  Size,
double  x 
)

Evaluate Chebyshev polynomials up to Size at point x It returns a vector of the P_i(x) for i=0..Size-1 The polynomial is defined for x between 0 and 1.

Definition at line 77 of file nodes.cc.

78 {
79 VectorXd pol(Size);
80 double C = 2 * (2 * x - 1);
81
82 if (Size >= 1) pol[0] = 1;
83 if (Size >= 2) pol[1] = C / 2;
84
85 for (int i = 2; i < Size; ++i)
86 pol[i] = C * pol[i - 1] - pol[i - 2];
87 return pol;
88 }

◆ getPolsSum()

VectorXd getPolsSum ( int  Size,
VectorXd  x 
)

Calculate the Chebyshev polynomials of order i=0..Size-1 at points given in vector x_j and sum it over point index j It returns sum_j P_i(x_j) for i=0..Size-1 The Chebyshev polynomials are defined for x between 0 and 1.

Definition at line 96 of file nodes.cc.

97 {
98 assert(Size > 2);
99
100 VectorXd polSum(Size);
101
102 VectorXd pol0 = 0 * x.array() + 1;
103 VectorXd pol1 = 2 * x.array() - 1;
104 VectorXd C = 2 * pol1;
105
106 VectorXd pol2(x.size());
107 for (int i = 2; i < Size; ++i) {
108 polSum(i - 2) = pol0.sum();
109
110 pol2 = C.array() * pol1.array() - pol0.array();
111
112 pol0 = pol1;
113 pol1 = pol2;
114 }
115
116 polSum(Size - 2) = pol0.sum();
117 polSum(Size - 1) = pol1.sum();
118
119 return polSum;
120 }

◆ getPosition()

ExpRunEvt getPosition ( const std::vector< Evt > &  events,
double  tEdge 
)
inline

Get the exp-run-evt number from the event time [hours].

Parameters
eventsvector of events
tEdgethe event time of the event of interest [hours]
Returns
the position of the time point in the exp-run-evt format

Definition at line 341 of file Splitter.h.

342 {
343 ExpRunEvt evt(-1, -1, -1);
344 double tBreak = -1e10;
345 for (auto& e : events) {
346 if (e.t < tEdge) {
347 if (e.t > tBreak) {
348 tBreak = e.t;
349 evt = ExpRunEvt(e.exp, e.run, e.evtNo);
350 }
351 }
352 }
353 return evt;
354 }

◆ getRangeLin()

std::vector< double > getRangeLin ( int  nVals,
double  xMin,
double  xMax 
)
inline

Equidistant range between xMin and xMax for spline of the first order.

Definition at line 83 of file tools.h.

84 {
85 B2ASSERT("At least one value in the spline required", nVals >= 1);
86 if (nVals == 1) return {};
87 std::vector<double> v(nVals);
88 for (int i = 0; i < nVals; ++i)
89 v[i] = xMin + i * (xMax - xMin) / (nVals - 1);
90 return v;
91 }

◆ getRangeZero()

std::vector< double > getRangeZero ( int  nVals,
double  xMin,
double  xMax 
)
inline

Equidistant range between xMin and xMax for spline of the zero order.

Definition at line 94 of file tools.h.

95 {
96 B2ASSERT("At least one value in the spline required", nVals >= 1);
97 if (nVals == 1) return {};
98 std::vector<double> v(nVals - 1);
99 for (int i = 1; i < nVals; ++i)
100 v[i - 1] = xMin + i * (xMax - xMin) / (nVals);
101 return v;
102 }

◆ getRun()

static ExpRun getRun ( std::map< ExpRun, std::pair< double, double > >  runs,
double  t 
)
static

Get exp number + run number from time.

Parameters
runsmap, where key contain the exp-run number and value the start- and end-time of the run
ttime of interest [hours]
Returns
: the exp-run number at the input time t

Definition at line 262 of file Splitter.cc.

263 {
264 ExpRun rFound(-1, -1);
265 int nFound = 0;
266 for (auto r : runs) { //Linear search over runs
267 if (r.second.first <= t && t < r.second.second) {
268 rFound = r.first;
269 ++nFound;
270 }
271 }
272
273 B2ASSERT("Exactly one interval should be found", nFound == 1);
274 B2ASSERT("Assert that something was found", rFound != ExpRun(-1, -1));
275 return rFound;
276 }

◆ getRunInfo()

std::map< ExpRun, std::pair< double, double > > getRunInfo ( const std::vector< Evt > &  evts)
inline

Get the map of runs, where each run contains pair with start/end time [hours].

Parameters
evtsvector of events
Returns
a map where the key is exp-run and value start/end time of the particular run [hours]

Definition at line 312 of file Splitter.h.

313 {
314 std::map<ExpRun, std::pair<double, double>> runsInfo;
315
316 for (auto& evt : evts) {
317 int Exp = evt.exp;
318 int Run = evt.run;
319 double time = evt.t;
320 if (runsInfo.count(ExpRun(Exp, Run))) {
321 double tMin, tMax;
322 std::tie(tMin, tMax) = runsInfo.at(ExpRun(Exp, Run));
323 tMin = std::min(tMin, time);
324 tMax = std::max(tMax, time);
325 runsInfo.at(ExpRun(Exp, Run)) = {tMin, tMax};
326 } else {
327 runsInfo[ExpRun(Exp, Run)] = {time, time};
328 }
329
330 }
331 return runsInfo;
332 }

◆ getStartEndIndexes()

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 at line 83 of file Splitter.cc.

84 {
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); //There is one more interval than #breaks
88 int s = (indx == 0) ? 0 : breaks[indx - 1] + 1;
89 int e = (indx == int(breaks.size())) ? nIntervals - 1 : breaks[indx];
90 return {s, e};
91 }

◆ getSVDDedxFromParticle()

VXDDedxTrack const * getSVDDedxFromParticle ( Particle const *  particle)

SVD dEdx value from particle.

Definition at line 55 of file DedxVariables.cc.

56 {
57 const Track* track = particle->getTrack();
58 if (!track) {
59 return nullptr;
60 }
61
62 const VXDDedxTrack* dedxTrack = track->getRelatedTo<VXDDedxTrack>();
63 if (!dedxTrack) {
64 return nullptr;
65 }
66 return dedxTrack;
67 }

◆ getWeights()

VectorXd getWeights ( int  Size)

Get the vector of weights to calculate the integral over the Chebyshev nodes The nodes are by definition between 0 and 1, there are Size nodes To get their positions, use getNodes.

Definition at line 29 of file nodes.cc.

30 {
31 const int N = Size - 1;
32 assert(N % 2 == 0);
33
34 std::vector<std::vector<double>> coef(Size);
35 for (auto& el : coef) el.resize(Size);
36
37
38 for (int k = 0; k <= N / 2; ++k) {
39 coef[2 * k][N] = 1. / N;
40 coef[2 * k][0] = 1. / N ;
41
42 coef[2 * k][N / 2] = 2. / N * (2 * ((k + 1) % 2) - 1);
43
44 for (int n = 1; n <= N / 2 - 1; ++n)
45 coef[2 * k][n] = coef[2 * k][N - n] = 2. / N * cos(n * k * M_PI * 2 / N);
46 }
47
48 VectorXd wgt = VectorXd::Zero(Size);
49
50
51 for (int i = 0; i < Size; ++i) {
52 wgt[i] += coef[0][i];
53 wgt[i] += coef[N][i] / (1. - N * N);
54 for (int k = 1; k <= N / 2 - 1; ++k) {
55 double w = 2. / (1 - 4 * k * k);
56 wgt[i] += w * coef[2 * k][i];
57 }
58
59 wgt[i] *= 0.5; //for interval (0,1)
60 }
61 return wgt;
62 }

◆ init()

void init ( int  Size,
double  xMin,
double  xMax 
)

Initialize the fitter (the Chebyshev coefficients)

Definition at line 59 of file ChebFitter.cc.

60 {
61 // loading the Cheb nodes
62 m_nodes = (xMax - xMin) * getNodes(Size).array() + xMin;
63
64 // loading the weights for integration
65 m_weights = (xMax - xMin) * getWeights(Size);
66
67
68 // calculate the transformation matrix from pol coefs to grid points
69 m_coefsMat = getCoefsCheb(Size).transpose();
70
72
73 }
Eigen::VectorXd getWeights(int Size)
Get the vector of weights to calculate the integral over the Chebyshev nodes The nodes are by definit...
Definition: nodes.cc:29
Eigen::VectorXd getDataGrid() const
Calculate Chebyshev coefficients for the data set.
Definition: ChebFitter.cc:109
Eigen::VectorXd getNodes(int Size)
Get the vector of positions of the Chebyshev nodes The nodes are by definition between 0 and 1,...
Definition: nodes.cc:65

◆ interpol() [1/3]

VectorXd interpol ( const VectorXd &  xi,
double  x 
)

Get Interpolation vector k_i for point x from the function values at points xi (polynomial interpolation) In the second step, the function value at x can be evaluated as sum_i vals_i k_i.

Definition at line 189 of file nodes.cc.

190 {
191 double Norm = (xi[xi.size() - 1] - xi[0]) / 2;
192 VectorXd coefs(xi.size());
193 for (int i = 0; i < xi.size(); ++i) {
194 double num = 1, den = 1;
195 for (int j = 0; j < xi.size(); ++j)
196 if (j != i) {
197 num *= (x - xi(j)) / Norm;
198 den *= (xi(i) - xi(j)) / Norm;
199 }
200 coefs(i) = num / den;
201 }
202 return coefs;
203 }

◆ interpol() [2/3]

double interpol ( Eigen::VectorXd  xi,
Eigen::VectorXd  vals,
double  x 
)

Get interpolated function value at point x when function values vals at points xi are provided.

If the points xi are fixed and only vals are different between interpol calls, use interpol(xi, x) to speed up the evaluation.

◆ interpol() [3/3]

double interpol ( VectorXd  xi,
VectorXd  vals,
double  x 
)

Get interpolated function value at point x when function values vals at points xi are provided.

If the points xi are fixed and only vals are different between interpol calls, use interpol(xi, x) to speed up the evaluation.

Definition at line 210 of file nodes.cc.

211 {
212 VectorXd coefs = interpol(xi, x);
213 return coefs.dot(vals);
214 }

◆ lossFunction()

double lossFunction ( const std::vector< Atom > &  vec,
int  s,
int  e 
) const
private

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

the lossFunction formula (it can be modified according to the user's taste)

Parameters
vecVector of atoms, where each atom is an intervals in time
sFirst index of the calib. interval
eLast index of the calib. interval
Returns
: A value of the loss function

Definition at line 135 of file Splitter.cc.

136 {
137
138 //raw time
139 double rawTime = vec[e].t2 - vec[s].t1;
140
141 //max gap
142 double maxGap = 0;
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);
146 }
147
148 //net time
149 double netTime = 0;
150 for (int i = s; i <= e; ++i) {
151 netTime += vec[i].t2 - vec[i].t1;
152 }
153
154 // Number of events
155 double nEv = 0;
156 for (int i = s; i <= e; ++i) {
157 nEv += vec[i].nEv;
158 }
159 if (nEv == 0) nEv = 0.1;
160
161 //double loss = pow(rawTime - tBest, 2) + gapPenalty * pow(maxGap, 2);
162 //double loss = 1./nEv + timePenalty * pow(rawTime, 2);
163
164 lossFun->SetParameters(rawTime, netTime, maxGap, nEv);
165 double lossNew = lossFun->Eval(0);
166
167 return lossNew;
168 }
TF1 * lossFun
loss function used for clustering of Atoms
Definition: Splitter.h:297

◆ merge()

std::vector< std::vector< double > > merge ( std::vector< std::vector< std::vector< double > > >  toMerge)
inline

merge { vector<double> a, vector<double> b} into {a, b}

Definition at line 41 of file tools.h.

42 {
43 std::vector<std::vector<double>> allVecs;
44 for (const auto& v : toMerge)
45 allVecs.insert(allVecs.end(), v.begin(), v.end());
46 return allVecs;
47 }

◆ mergeIntervals()

std::map< ExpRun, std::pair< double, double > > mergeIntervals ( std::map< ExpRun, std::pair< double, double > >  I1,
std::map< ExpRun, std::pair< double, double > >  I2 
)
static

Merge two subintervals into one subinterval.

Parameters
I1First subinterval to merge
I2Second subinterval to merge
Returns
: The resulting subinterval

Definition at line 306 of file Splitter.cc.

308 {
309 std::map<ExpRun, std::pair<double, double>> I = I1;
310 for (auto r : I2) {
311 ExpRun run = r.first;
312 if (I.count(run) == 0)
313 I[run] = r.second;
314 else {
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));
316 }
317 }
318 return I;
319 }

◆ operator!=()

bool operator!= ( ExpRun  a,
ExpRun  b 
)
inline

Not equal for ExpRun.

Definition at line 71 of file Splitter.h.

71{ return (a.exp != b.exp || a.run != b.run); }

◆ operator()()

double operator() ( const double *  par) const

Evaluate the log likelihood.

Definition at line 97 of file ChebFitter.cc.

98 {
99 Pars pars;
100 for (unsigned i = 0; i < m_parNames.size(); ++i)
101 pars[m_parNames[i]] = par[i];
102
103 double LL = m_useCheb ? getLogLikelihoodFast(pars) : getLogLikelihoodSlow(pars);
104
105 return LL;
106 }
double getLogLikelihoodSlow(const Pars &pars) const
Calculate log likelihood using exact formula.
Definition: ChebFitter.cc:76
double getLogLikelihoodFast(const Pars &pars) const
Calculate log likelihood using approximation based on Chebyshev polynomials (typically faster)
Definition: ChebFitter.cc:89

◆ operator<()

bool operator< ( ExpRun  a,
ExpRun  b 
)
inline

less than for ExpRun

Definition at line 74 of file Splitter.h.

74{return ((a.exp < b.exp) || (a.exp == b.exp && a.run < b.run));}

◆ plotRuns()

void plotRuns ( std::vector< std::pair< double, double > >  runs)

plot runs on time axis

Definition at line 57 of file Splitter.cc.

58 {
59 TGraphErrors* gr = new TGraphErrors();
60
61
62 for (auto r : runs) {
63 double m = (r.first + r.second) / 2;
64 double e = (r.second - r.first) / 2;
65
66 gr->SetPoint(gr->GetN(), m, 1);
67 gr->SetPointError(gr->GetN() - 1, e, 0);
68 }
69
70 gStyle->SetEndErrorSize(6);
71
72 gr->SetLineWidth(1);
73 gr->SetMarkerSize(40);
74 gr->Draw("ape");
75 gr->GetYaxis()->SetRangeUser(-10, 10);
76 gr->GetXaxis()->SetRangeUser(0, 256);
77
78 gr->GetXaxis()->SetTitle("time [hours]");
79
80 }

◆ plotSRuns()

void plotSRuns ( std::vector< std::pair< double, double > >  runs,
std::vector< int >  breaks,
int  offset = 2 
)

plot clusters or runs on time axis

Definition at line 94 of file Splitter.cc.

95 {
96 TGraphErrors* gr = new TGraphErrors();
97
98 for (int i = 0; i < int(breaks.size()) + 1; ++i) {
99 int s, e;
100 std::tie(s, e) = getStartEndIndexes(runs.size(), breaks, i);
101 double a = runs[s].first;
102 double b = runs[e].second;
103
104 double m = (a + b) / 2;
105 double err = (b - a) / 2;
106
107 gr->SetPoint(gr->GetN(), m, offset);
108 gr->SetPointError(gr->GetN() - 1, err, 0);
109
110 }
111
112 gr->SetLineColor(kRed);
113 gr->SetMarkerColor(kRed);
114 gr->Draw("pe same");
115 }

◆ printBySize()

void printBySize ( std::vector< std::pair< double, double > >  runs)

print sorted lengths of the runs

Definition at line 119 of file Splitter.cc.

120 {
121 std::vector<double> dist;
122 for (auto r : runs) {
123 double d = r.second - r.first;
124 dist.push_back(d);
125 }
126
127 sort(dist.begin(), dist.end());
128
129 for (auto d : dist)
130 B2INFO(d);
131
132 }

◆ rn()

TString rn ( )
inline

Get random string.

Definition at line 38 of file tools.h.

38{return Form("%d", gRandom->Integer(1000000000)); }

◆ runAlgorithm()

CalibrationData runAlgorithm ( const std::vector< Evt > &  evts,
std::vector< std::map< ExpRun, std::pair< double, double > > >  range,
Fun  runCalibAnalysis 
)
inline

run calibration algorithm for single calibration interval

Definition at line 335 of file calibTools.h.

338 {
339 CalibrationData calD;
340 auto& r = range;
341 double rStart, rEnd;
342 std::tie(rStart, rEnd) = Splitter::getStartEnd(r);
343 B2INFO("Start of loop startTime endTime : " << rStart << " " << rEnd);
344
345 auto breaks = Splitter::getBreaks(r);
346
347 std::vector<Evt> evtsNow;
348
349 std::vector<int> Counts(breaks.size() + 1, 0);
350 // Select events belonging to the interval
351 for (const auto& ev : evts) {
352 if (rStart <= ev.t && ev.t < rEnd) {
353 evtsNow.push_back(ev);
354 ++Counts.at(getID(breaks, ev.t));
355 }
356 }
357
358 B2ASSERT("Number of intervals vs number of breakPoints", r.size() == breaks.size() + 1);
359
360 //Merge smallest interval if with low stat (try it 10times)
361 for (int k = 0; k < 10; ++k) {
362 int iMin = min_element(Counts.begin(), Counts.end()) - Counts.begin();
363 if (Counts.size() >= 2 && Counts[iMin] < 50) { //merge with neighbor if possible
364 auto iM = -1;
365 if (iMin == 0)
366 iM = iMin + 1;
367 else if (iMin == int(Counts.size()) - 1)
368 iM = iMin - 1;
369 else {
370 if (Counts[iMin + 1] < Counts[iMin - 1])
371 iM = iMin + 1;
372 else
373 iM = iMin - 1;
374 }
375 B2ASSERT("Number of intervals equal to size of counters", r.size() == Counts.size());
376
377 r.at(iM) = Splitter::mergeIntervals(r[iM], r[iMin]);
378 r.erase(r.begin() + iMin);
379 breaks = Splitter::getBreaks(r);
380 Counts[iM] += Counts[iMin];
381 Counts.erase(Counts.begin() + iMin);
382 }
383 }
384
385 B2INFO("#events " << " : " << evtsNow.size());
386 B2INFO("Breaks size " << " : " << breaks.size());
387
388 calD.breakPoints = convertSplitPoints(evtsNow, breaks);
389
390 calD.subIntervals = r;
391
392 if (breaks.size() > 0)
393 B2INFO("StartOfCalibInterval (run,evtNo,vtxIntervalsSize) " << calD.breakPoints.at(0).run << " " <<
394 calD.breakPoints.at(0).evt << " " << calD.breakPoints.size());
395
396
397 //If too few events, let have the output empty
398 //Will be filled with the closest neighbor at the next stage
399 if (evtsNow.size() < 50) {
400 return calD;
401 }
402
403 // Run the calibration
404 B2INFO("Start of running calibration over calibration interval");
405 tie(calD.pars.cnt, calD.pars.cntUnc, calD.pars.spreadMat) = runCalibAnalysis(evtsNow, breaks);
406 calD.pars.pulls.resize(calD.pars.cnt.size());
407 B2INFO("End of running analysis - SpreadMatX : " << sqrt(abs(calD.pars.spreadMat(0, 0))));
408 B2ASSERT("All subintervals have calibration of the mean value", calD.pars.cnt.size() == r.size());
409 B2ASSERT("All subintervals have calibration of the unc. of mean", calD.pars.cntUnc.size() == r.size());
410
411 calD.isCalibrated = true;
412
413 return calD;
414 }
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
std::vector< ExpRunEvt > convertSplitPoints(const std::vector< Evt > &events, std::vector< double > splitPoints)
Convert splitPoints [hours] to breakPoints in ExpRunEvt.
Definition: Splitter.h:363

◆ runCalibration()

CalibrationAlgorithm::EResult runCalibration ( TTree *  tracks,
const std::string &  calibName,
Fun1  GetEvents,
Fun2  calibAnalysis,
std::function< TObject *(Eigen::VectorXd, Eigen::MatrixXd, Eigen::MatrixXd)>  calibObjCreator,
TString  m_lossFunctionOuter,
TString  m_lossFunctionInner 
)

Run the the calibration over the whole event sample.

Parameters
tracksTTree object with mu-mu events
calibNamename of the calibration payload
GetEventsfunction that transforms TTree to std::vector
calibAnalysisfunction that performs the calibration on a single calibration interval
calibObjCreatorfunction that stores results to the payload class which inherits from TObject
m_lossFunctionOuterLost function for the calibration intervals of the spread parameters
m_lossFunctionInnerLost function for the calibration subintervals (for the mean value parameters)
Returns
State of the calibration run, i.e. EResult::c_OK if everything OK

Definition at line 428 of file calibTools.h.

431 {
432 // Check that there are at least some data
433 if (!tracks || tracks->GetEntries() < 15) {
434 if (tracks)
435 B2WARNING("Too few data : " << tracks->GetEntries());
436 return CalibrationAlgorithm::EResult::c_NotEnoughData;
437 }
438 B2INFO("Number of tracks: " << tracks->GetEntries());
439
440 // Tree to vector of Events
441 auto evts = GetEvents(tracks);
442
443 //Time range for each ExpRun
444 std::map<ExpRun, std::pair<double, double>> runsInfoOrg = getRunInfo(evts);
445 std::map<ExpRun, std::pair<double, double>> runsRemoved; //map with time intervals of very short runs
446 auto runsInfo = filter(runsInfoOrg, 2. / 60, runsRemoved); //include only runs longer than 2mins
447
448 // If nothing remains
449 if (runsInfo.size() == 0) {
450 B2WARNING("Too short run");
451 return CalibrationAlgorithm::EResult::c_NotEnoughData;
452 }
453
454 // Get intervals based on the input loss functions
455 Splitter splt;
456 auto splits = splt.getIntervals(runsInfo, evts, m_lossFunctionOuter, m_lossFunctionInner);
457
458 //Loop over all calibration intervals
459 std::vector<CalibrationData> calVec;
460 for (auto s : splits) {
461 CalibrationData calD = runAlgorithm(evts, s, calibAnalysis); // run the calibration over the interval s
462 calVec.push_back(calD);
463 }
464
465 // extrapolate results to the low-stat intervals
467
468 // Include removed short runs
469 for (auto shortRun : runsRemoved) {
470 addShortRun(calVec, shortRun);
471 }
472
473 // Store Payloads to files
474 storePayloads(evts, calVec, calibName, calibObjCreator);
475
476 return CalibrationAlgorithm::EResult::c_OK;
477 }
void addShortRun(std::vector< CalibrationData > &calVec, std::pair< ExpRun, std::pair< double, double > > shortRun)
Extrapolate calibration to the very short runs which were filtered before.
Definition: calibTools.h:151
void storePayloads(const std::vector< Evt > &evts, const std::vector< CalibrationData > &calVecConst, std::string objName, std::function< TObject *(Eigen::VectorXd, Eigen::MatrixXd, Eigen::MatrixXd) > getCalibObj)
Store payloads to files.
Definition: calibTools.h:225
CalibrationData runAlgorithm(const std::vector< Evt > &evts, std::vector< std::map< ExpRun, std::pair< double, double > > > range, Fun runCalibAnalysis)
run calibration algorithm for single calibration interval
Definition: calibTools.h:335
void extrapolateCalibration(std::vector< CalibrationData > &calVec)
Extrapolate calibration to intervals where it failed.
Definition: calibTools.h:102

◆ slice() [1/2]

std::vector< Atom > slice ( std::vector< Atom vec,
int  s,
int  e 
)
inline

Slice the vector to contain only elements with indexes s .. e (included)

Definition at line 85 of file Splitter.h.

86 {
87 return std::vector<Atom>(vec.begin() + s, vec.begin() + e + 1);
88 }

◆ slice() [2/2]

std::vector< double > slice ( std::vector< double >  v,
unsigned  ind,
unsigned  n 
)
inline

put slice of original vector v[ind:ind+n] into new one, n is number of elements

Definition at line 106 of file tools.h.

107 {
108 std::vector<double> vNew;
109 for (unsigned i = ind; i < ind + n && i < v.size(); ++i)
110 vNew.push_back(v[i]);
111 return vNew;
112 }

◆ splitToSmall()

std::vector< std::pair< double, double > > splitToSmall ( std::map< ExpRun, std::pair< double, double > >  runs,
double  intSize = 1. / 60 
)
staticprivate

Split the runs into small calibration intervals (atoms) of a specified size.

By definition each of these intervals spans only over single run. These will be clustered into larger intervals in the next steps

Parameters
runsRuns to split into the atoms
intSizeIntended size of the small intervals
Returns
: A vector with resulting time boundaries of the atoms

Definition at line 173 of file Splitter.cc.

174 {
175 // split into small intervals
176 std::vector<std::pair<double, double>> smallRuns;
177
178 for (auto r : runs) {
179 auto& I = r.second;
180 if (intSize < 0) {
181 smallRuns.push_back(I);
182 continue;
183 }
184
185 double runTime = I.second - I.first;
186 int nSplits = runTime / intSize; //1-m intervals
187 nSplits = std::max(1, nSplits); //at least 1 interval
188
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});
193 }
194 }
195 return smallRuns;
196 }

◆ storePayloads()

void storePayloads ( const std::vector< Evt > &  evts,
const std::vector< CalibrationData > &  calVecConst,
std::string  objName,
std::function< TObject *(Eigen::VectorXd, Eigen::MatrixXd, Eigen::MatrixXd) >  getCalibObj 
)
inline

Store payloads to files.

Definition at line 225 of file calibTools.h.

227 {
228 auto calVec = calVecConst;
229
230 // Loop to store payloads
231 ExpRun exprunLast(-1, -1); //last exprun
232 EventDependency* intraRun = nullptr;
233
234 // Loop over calibration intervals
235 for (unsigned i = 0; i < calVec.size(); ++i) {
236 const auto& r = calVec[i].subIntervals; // splits[i];
237 // Loop over calibration subintervals
238 for (int k = 0; k < int(r.size()); ++k) {
239
240 for (auto I : r[k]) { //interval required to be within single run
241 ExpRun exprun = I.first;
242
243 //Encode Start+End time in seconds of the payload
244 if (calVec[i].pars.cntUnc.at(k).rows() == 3) {
245 calVec[i].pars.cntUnc.at(k)(0, 1) = calVec[i].pars.cntUnc.at(k)(1, 0) = encodeNumber(calVec[i].pars.cntUnc.at(k)(0, 1),
246 round(I.second.first * 3600));
247 calVec[i].pars.cntUnc.at(k)(0, 2) = calVec[i].pars.cntUnc.at(k)(2, 0) = encodeNumber(calVec[i].pars.cntUnc.at(k)(0, 2),
248 round(I.second.second * 3600));
249 } else {
250 calVec[i].pars.cntUnc.at(k)(0, 0) = encodeNumber(calVec[i].pars.cntUnc.at(k)(0, 0), round(I.second.first * 3600));
251 calVec[i].pars.spreadMat(0, 0) = encodeNumber(calVec[i].pars.spreadMat(0, 0), round(I.second.second * 3600));
252 }
253
254 TObject* obj = getCalibObj(calVec[i].pars.cnt.at(k), calVec[i].pars.cntUnc.at(k), calVec[i].pars.spreadMat);
255 if (exprun != exprunLast) { //if new run
256 if (intraRun) { //if not first -> store
257 auto m_iov = IntervalOfValidity(exprunLast.exp, exprunLast.run, exprunLast.exp, exprunLast.run);
258 Database::Instance().storeData(objName, intraRun, m_iov);
259 }
260
261 intraRun = new EventDependency(obj);
262 } else {
263 int breakPoint;
264 if (k - 1 >= 0) {
265 breakPoint = calVec[i].breakPoints.at(k - 1).evt;
266 B2ASSERT("Payload saving consistency", calVec[i].breakPoints.at(k - 1).run == exprun.run);
267 } else {
268 B2ASSERT("Payload saving consistency", i != 0);
269 double rStart, rEnd;
270 std::tie(rStart, rEnd) = Splitter::getStartEnd(r);
271 auto pos = getPosition(evts, rStart);
272 breakPoint = pos.evt;
273 B2ASSERT("Payload saving consistency", pos.run == exprun.run);
274 }
275 intraRun->add(breakPoint, obj);
276 }
277 exprunLast = exprun;
278 }
279 } //end loop over calibration subintervals
280
281 } //end loop over calibration intervals
282
283 //Store the last entry
284 auto m_iov = IntervalOfValidity(exprunLast.exp, exprunLast.run, exprunLast.exp, exprunLast.run);
285 Database::Instance().storeData(objName, intraRun, m_iov);
286 }
ExpRunEvt getPosition(const std::vector< Evt > &events, double tEdge)
Get the exp-run-evt number from the event time [hours].
Definition: Splitter.h:341

◆ storePayloadsNoIntraRun()

void storePayloadsNoIntraRun ( const std::vector< CalibrationData > &  calVecConst,
std::string  objName,
std::function< TObject *(Eigen::VectorXd, Eigen::MatrixXd, Eigen::MatrixXd) >  getCalibObj 
)
inline

Store payloads to files, where calib data have no intra-run dependence.

Definition at line 290 of file calibTools.h.

292 {
293 auto calVec = calVecConst;
294
295 // Check that there is no intra-run dependence
296 std::set<ExpRun> existingRuns;
297 for (unsigned i = 0; i < calVec.size(); ++i) {
298 const auto& r = calVec[i].subIntervals;
299 // Loop over calibration subintervals
300 for (int k = 0; k < int(r.size()); ++k) {
301
302 for (auto I : r[k]) {
303 ExpRun exprun = I.first;
304 // make sure that the run isn't already in the list, to avoid duplicity
305 if (existingRuns.count(exprun) != 0)
306 B2FATAL("Intra-run dependence exists");
307 existingRuns.insert(exprun);
308 }
309 }
310 }
311
312
313 // Loop over calibration intervals
314 for (unsigned i = 0; i < calVec.size(); ++i) {
315 const auto& r = calVec[i].subIntervals; // splits[i];
316 // Loop over calibration subintervals
317 for (unsigned k = 0; k < r.size(); ++k) {
318
319 TObject* obj = getCalibObj(calVec[i].pars.cnt.at(k), calVec[i].pars.cntUnc.at(k), calVec[i].pars.spreadMat);
320
321 ExpRun start = (r[k].cbegin()->first);
322 ExpRun last = (r[k].crbegin()->first);
323
324 auto iov = IntervalOfValidity(start.exp, start.run, last.exp, last.run);
325 Database::Instance().storeData(objName, obj, iov);
326
327
328 } //end loop over calibration subintervals
329 } //end loop over calibration intervals
330
331 }

◆ toB2Vector3()

B2Vector3D toB2Vector3 ( Eigen::VectorXd  vIn)
inline

Function that converts Eigen vector to ROOT vector.

Definition at line 54 of file calibTools.h.

55 {
56 return B2Vector3D(vIn(0), vIn(1), vIn(2));
57 }

◆ toTMatrixDSym()

TMatrixDSym toTMatrixDSym ( Eigen::MatrixXd  mIn)
inline

Function that converts Eigen symmetric matrix to ROOT matrix.

Definition at line 44 of file calibTools.h.

45 {
46 TMatrixDSym mOut(mIn.rows());
47 for (int i = 0; i < mIn.rows(); ++i)
48 for (int j = 0; j < mIn.cols(); ++j)
49 mOut(i, j) = (mIn(i, j) + mIn(j, i)) / 2.;
50 return mOut;
51 }

◆ vec2vec() [1/2]

std::vector< double > vec2vec ( Eigen::VectorXd  v)
inline

ROOT vector -> std vector.

Definition at line 61 of file tools.h.

62 {
63 std::vector<double> vNew(v.rows());
64 for (int i = 0; i < v.rows(); ++i)
65 vNew[i] = v(i);
66 return vNew;
67 }

◆ vec2vec() [2/2]

Eigen::VectorXd vec2vec ( std::vector< double >  vec)
inline

std vector -> ROOT vector

Definition at line 51 of file tools.h.

52 {
53 Eigen::VectorXd v(vec.size());
54 for (unsigned i = 0; i < vec.size(); ++i) {
55 v[i] = vec[i];
56 }
57 return v;
58 }

◆ vecs2mat()

Eigen::MatrixXd vecs2mat ( std::vector< std::vector< double > >  vecs)
inline

merge columns (from std::vectors) into ROOT matrix

Definition at line 72 of file tools.h.

73 {
74 Eigen::MatrixXd m(vecs[0].size(), vecs.size());
75 for (unsigned i = 0; i < vecs[0].size(); ++i)
76 for (unsigned j = 0; j < vecs.size(); ++j) {
77 m(i, j) = vecs[j][i];
78 }
79 return m;
80 }