Belle II Software development
InvariantMassAlgorithm.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 <mdst/dbobjects/CollisionInvariantMass.h>
11#include <reconstruction/calibration/BeamSpotBoostInvMass/InvariantMassAlgorithm.h>
12#include <reconstruction/calibration/BeamSpotBoostInvMass/InvariantMassMuMuStandAlone.h>
13#include <reconstruction/calibration/BeamSpotBoostInvMass/InvariantMassBhadStandAlone.h>
14#include <reconstruction/calibration/BeamSpotBoostInvMass/calibTools.h>
15#include <TFile.h>
16
17#include <Eigen/Dense>
18#include <iostream>
19#include <fstream>
20#include <iomanip>
21
22using Eigen::VectorXd;
23using Eigen::MatrixXd;
24
25using namespace Belle2;
26
27
28//Using boostVector collector for the input
30{
31 setDescription("Collision invariant mass calibration algorithm");
32}
33
34
36static TObject* getInvariantMassObj(VectorXd vMass, MatrixXd vMassUnc, MatrixXd vMassSpread)
37{
38 auto payload = new CollisionInvariantMass();
39
40 double mass = vMass(0);
41 double unc = vMassUnc(0, 0);
42 double spread = vMassSpread(0, 0);
43
44 payload->setMass(mass, unc, spread);
45 TObject* obj = static_cast<TObject*>(payload);
46 return obj;
47}
48
49
57template<typename Fun>
58std::vector<CalibrationData> runMuMuCalibration(const std::vector<Belle2::InvariantMassMuMuCalib::Event>& evts, Fun calibAnalysis,
59 TString lossFunctionOuter, TString lossFunctionInner)
60{
61 //Time range for each ExpRun
62 std::map<ExpRun, std::pair<double, double>> runsInfoOrg = getRunInfo(evts);
63 std::map<ExpRun, std::pair<double, double>> runsRemoved; //map with time intervals of very short runs
64 auto runsInfo = filter(runsInfoOrg, 2. / 60, runsRemoved); //include only runs longer than 2mins
65
66 // If nothing remains
67 if (runsInfo.size() == 0) {
68 B2WARNING("Too short run");
69 return {};
70 }
71
72 // Get intervals based on the input loss functions
73 Splitter splt;
74 auto splits = splt.getIntervals(runsInfo, evts, lossFunctionOuter, lossFunctionInner, 100 /*maximal length of the run */);
75
76 //Loop over all calibration intervals
77 std::vector<CalibrationData> calVec;
78 for (auto s : splits) {
79 CalibrationData calD = runAlgorithm(evts, s, calibAnalysis); // run the calibration over the interval s
80 calVec.push_back(calD);
81 }
82
83 // extrapolate results to the low-stat intervals
85
86 // Include removed short runs
87 for (auto shortRun : runsRemoved) {
88 addShortRun(calVec, shortRun);
89 }
90
91 return calVec;
92}
93
94
95std::vector<CalibrationData> addSpreadAndOffset(const std::vector<CalibrationData>& mumuCalResults, double spread, double spreadUnc,
96 double offset, double offsetUnc)
97{
98 std::vector<CalibrationData> mumuCalResultsNew = mumuCalResults;
99
100 for (auto& calibNew : mumuCalResultsNew) {
101
102 for (unsigned j = 0; j < calibNew.subIntervals.size(); ++j) {
103 calibNew.pars.cnt[j](0) += offset;
104 calibNew.pars.spreadMat(0, 0) = spread;
105 calibNew.pars.spreadUnc = spreadUnc;
106
107 calibNew.pars.shift = offset;
108 calibNew.pars.shiftUnc = offsetUnc;
109 calibNew.pars.pulls[j] = 0;
110 }
111 }
112 return mumuCalResultsNew;
113}
114
115
116
117
118static std::vector<TTree*> getTrees(TString tag, TFile* f)
119{
120 if (f == nullptr)
121 return {};
122
123 f->cd(tag + "/events");
124 TList* l = f->CurrentDirectory().load()->GetListOfKeys();
125
126
127 std::vector<TTree*> treeVec;
128
129 for (const auto&& obj : *l) {
130 TString n = obj->GetName();
131
132 TString trDir = tag + "/events/" + TString(obj->GetName()) + "/events_1";
133
134 treeVec.push_back((TTree*) f->Get(trDir));
135 }
136
137 return treeVec;
138}
139
140
141static std::vector<InvariantMassMuMuCalib::Event> readMuMuFiles(std::vector<std::string> files, bool is4S)
142{
143 std::vector<InvariantMassMuMuCalib::Event> events;
144 for (auto fName : files) {
145 TFile* f = TFile::Open(fName.c_str(), "READ");
146 auto trees = getTrees("BoostVectorCollector", f);
147
148 for (auto tr : trees) {
149 if (!tr) continue;
150 std::vector<InvariantMassMuMuCalib::Event> eventsNow = InvariantMassMuMuCalib::getEvents(tr, is4S);
151 events.insert(events.end(), eventsNow.begin(), eventsNow.end());
152 }
153
154 f->Close();
155 }
156
157 sort(events.begin(), events.end(), [](InvariantMassMuMuCalib::Event e1, InvariantMassMuMuCalib::Event e2) {return e1.t < e2.t;});
158 return events;
159}
160
161
162static std::vector<InvariantMassBhadCalib::Event> readBhadFiles(std::vector<std::string> files)
163{
164 std::vector<InvariantMassBhadCalib::Event> events;
165 for (auto fName : files) {
166 TFile* f = TFile::Open(fName.c_str(), "READ");
167 auto trees = getTrees("EcmsCollector", f);
168
169 for (auto tr : trees) {
170 if (!tr) continue;
171 std::vector<InvariantMassBhadCalib::Event> eventsNow = InvariantMassBhadCalib::getEvents(tr);
172 events.insert(events.end(), eventsNow.begin(), eventsNow.end());
173 }
174
175 f->Close();
176 }
177
178 sort(events.begin(), events.end(), [](InvariantMassBhadCalib::Event e1, InvariantMassBhadCalib::Event e2) {return e1.t < e2.t;});
179 return events;
180}
181
182
183
184
185
186
187
188std::vector<InvariantMassMuMuCalib::Event> InvariantMassAlgorithm::getDataMuMu(const std::vector<std::string>& files, bool is4S)
189{
190 if (files.size() == 0)
191 return {};
192
194 setPrefix("BoostVectorCollector");
195
196 setInputFileNames(files);
198
199 TTree* eventsTr = getObjectPtr<TTree>("events").get();
200
201 // Check that there are at least some mumu data
202 if (!eventsTr || eventsTr->GetEntries() < 15) {
203 if (eventsTr) {
204 if (is4S)
205 B2WARNING("Small number of events in the 4S mumu sample, only " << eventsTr->GetEntries());
206 else
207 B2WARNING("Small number of events in the off-resonance mumu sample, only " << eventsTr->GetEntries());
208 } else
209 B2WARNING("Pointer to the \"events\" object not defined for the mumu sample");
210
211 return {};
212 }
213 B2INFO("Number of mumu events: " << eventsTr->GetEntries());
214
215
216 return InvariantMassMuMuCalib::getEvents(eventsTr, is4S);
217
218}
219
220std::vector<InvariantMassBhadCalib::Event> InvariantMassAlgorithm::getDataHadB(const std::vector<std::string>& files)
221{
222 if (files.size() == 0)
223 return {};
224
226 setPrefix("EcmsCollector");
227
228 setInputFileNames(files);
230 TTree* eventsTr = getObjectPtr<TTree>("events").get();
231
232 // Check that there are at least some mumu data
233 if (!eventsTr || eventsTr->GetEntries() < 15) {
234 if (eventsTr)
235 B2WARNING("Too few data, only " << eventsTr->GetEntries() << " events found in hadronic B decay sample");
236 else
237 B2WARNING("Pointer to the \"events\" object not defined for hadronic B decay sample");
238
239 return {};
240 }
241 B2INFO("Number of mumu events: " << eventsTr->GetEntries());
242
243
244 return InvariantMassBhadCalib::getEvents(eventsTr);
245
246}
247
248
250double weightAvg(double x, double xe, double y, double ye)
251{
252 double xe2 = xe * xe;
253 double ye2 = ye * ye;
254 double res = (ye2 * x + xe2 * y) / (xe2 + ye2);
255 return res;
256}
257
258
260void printToFile(const std::vector<CalibrationData>& CalResults, TString outFileName)
261{
262 //Store info to the text file
263 std::ofstream finalOut(outFileName);
264 finalOut <<
265 "t1 t2 exp1 run1 exp2 run2 state Ecms EcmsUnc pull shift shiftUnc spread spreadUnc" <<
266 std::endl;
267
268 for (auto cal : CalResults) {
269
270 for (unsigned j = 0; j < cal.subIntervals.size(); ++j) {
271 double cnt = cal.pars.cnt[j](0);
272 double cntUnc = cal.pars.cntUnc[j](0, 0);
273 double spread = cal.pars.spreadMat(0, 0);
274 double spreadUnc = cal.pars.spreadUnc;
275
276 double shift = cal.pars.shift;
277 double shiftUnc = cal.pars.shiftUnc;
278 double pull = cal.pars.pulls[j];
279
280 double s = cal.subIntervals[j].begin()->second.first;
281 double e = cal.subIntervals[j].rbegin()->second.second;
282
283 double exp1 = cal.subIntervals[j].begin()->first.exp;
284 double exp2 = cal.subIntervals[j].rbegin()->first.exp;
285 double run1 = cal.subIntervals[j].begin()->first.run;
286 double run2 = cal.subIntervals[j].rbegin()->first.run;
287
288
289 finalOut << std::setprecision(8) << s << " " << e << " " << exp1 << " " << run1 << " " << exp2 << " " << run2 << " " << j << " "
290 << 1e3 * cnt << " " << 1e3 * cntUnc << " " << pull <<
291 " " << 1e3 * shift << " " << 1e3 *
292 shiftUnc << " " << 1e3 * spread << " " << 1e3 * spreadUnc << std::endl;
293
294 }
295
296 }
297 finalOut.close();
298}
299
300
301
302std::vector<std::vector<CalibrationData>> InvariantMassAlgorithm::adjustOffResonanceEnergy(std::vector<std::vector<CalibrationData>>
303 CalResultsBlocks,
304 const std::vector<std::vector<InvariantMassMuMuCalib::Event>>& evtsMuMuBlocks)
305{
306 // Adjust the energy in the off-resoncance runs
307 for (int b = 0; b < int(CalResultsBlocks.size()); ++b) {
308 // only deal with off-res runs
309 if (evtsMuMuBlocks[b][0].is4S) continue;
310
311 std::vector<CalibPars> parsEdges;
312 if (b > 0 && evtsMuMuBlocks[b - 1][0].is4S)
313 parsEdges.push_back(CalResultsBlocks[b - 1].back().pars);
314 if (b < int(CalResultsBlocks.size()) - 1 && evtsMuMuBlocks[b + 1][0].is4S)
315 parsEdges.push_back(CalResultsBlocks[b + 1].front().pars);
316 std::vector<double> shifts, shiftsUnc, spreads, spreadsUnc;
317 for (auto p : parsEdges) {
318 shifts.push_back(p.shift);
319 shiftsUnc.push_back(p.shiftUnc);
320 spreads.push_back(p.spreadMat(0, 0));
321 spreadsUnc.push_back(p.spreadUnc);
322 }
323
324 // by default take the offset from the steering
325 double spread = m_eCMSmumuSpread;
326 double spreadUnc = 0.2e-3; // default spread for off-resonance [GeV]
327 double shift = m_eCMSmumuShift;
328 double shiftUnc = 0.3e-3; // default shift unc [GeV]
329
330 // in case of single neighboring 4S block
331 if (shifts.size() == 1) {
332 spread = spreads[0];
333 spreadUnc = spreadsUnc[0];
334 shift = shifts[0];
335 shiftUnc = shiftsUnc[0];
336 }
337 // if neighbors are from both sides take the average
338 else if (shifts.size() == 2) {
339 spread = weightAvg(spreads[0], spreadsUnc[0], spreads[1], spreadsUnc[1]);
340 shift = weightAvg(shifts[0], shiftsUnc[0], shifts[1], shiftsUnc[1]);
341
342 spreadUnc = sqrt(pow(spreads[0] - spreads[1], 2) + (pow(spreadsUnc[0], 2) + pow(spreadsUnc[1], 2)) / 2);
343 shiftUnc = sqrt((pow(shift - shifts[0], 2) + pow(shift - shifts[1], 2)) / 2 + (pow(shiftsUnc[0], 2) + pow(shiftsUnc[1], 2)) / 2);
344 }
345
346 CalResultsBlocks[b] = addSpreadAndOffset(CalResultsBlocks[b], spread, spreadUnc, shift, shiftUnc);
347 }
348 return CalResultsBlocks;
349}
350
351
352
353
354
355
356
357
358/* Main calibration method calling dedicated functions */
360{
361 std::vector<std::string> files = getVecInputFileNames();
362
363 std::vector<std::string> filesHad4S, filesMuMu4S, filesMuMuOff;
364 for (auto f : files) {
365 if (f.find("/dimuon_4S/") != std::string::npos)
366 filesMuMu4S.push_back(f);
367 else if (f.find("/dimuon_Off/") != std::string::npos)
368 filesMuMuOff.push_back(f);
369 else if (f.find("/hadB_4S/") != std::string::npos)
370 filesHad4S.push_back(f);
371 else {
372 B2FATAL("Unrecognised data type");
373 }
374 }
375
376
377 for (auto r : filesMuMu4S)
378 B2INFO("MuMu4SFile name " << r);
379 for (auto r : filesMuMuOff) {
380 B2INFO("MuMuOffFile name " << r);
381 }
382 for (auto r : filesHad4S)
383 B2INFO("Had4SFile name " << r);
384
385
386
387 // load the mumu data
388 const auto evtsMuMuOff = readMuMuFiles(filesMuMuOff, false /*is4S*/);
389 const auto evtsMuMu4S = readMuMuFiles(filesMuMu4S, true /*is4S*/);
390
391 if (evtsMuMuOff.size() + evtsMuMu4S.size() < 50) {
392 B2WARNING("Not enough data, there are only " << evtsMuMuOff.size() + evtsMuMu4S.size() << " mumu events");
394 }
395
396
397 // load hadB data
398 const auto evtsHad = readBhadFiles(filesHad4S);
399
400 // merge mumu data
401 std::vector<InvariantMassMuMuCalib::Event> evtsMuMuTemp = evtsMuMu4S;
402 evtsMuMuTemp.insert(evtsMuMuTemp.end(), evtsMuMuOff.begin(), evtsMuMuOff.end());
403
404 if (evtsMuMuTemp.size() < 50) {
405 B2WARNING("Not enough mumu data, there are only " << evtsMuMuTemp.size() << " mumu events");
407 }
408
409
410 // sort by time
411 std::sort(evtsMuMuTemp.begin(), evtsMuMuTemp.end(), [](const InvariantMassMuMuCalib::Event & evt1,
412 const InvariantMassMuMuCalib::Event & evt2) { return evt1.t < evt2.t; });
413
414 //split the mumu events into the blocks with identical run type
415 std::vector<std::vector<InvariantMassMuMuCalib::Event>> evtsMuMuBlocks;
416 bool is4Sold = evtsMuMuTemp[0].is4S;
417 evtsMuMuBlocks.push_back({}); // empty vec to start
418 for (const auto& ev : evtsMuMuTemp) {
419 if (is4Sold != ev.is4S) {
420 evtsMuMuBlocks.push_back({}); // adding new entry
421 is4Sold = ev.is4S;
422 }
423 evtsMuMuBlocks.back().push_back(ev);
424 }
425 B2INFO("Number of mumu 4S events " << evtsMuMu4S.size());
426 B2INFO("Number of mumu off-res events " << evtsMuMuOff.size());
427 B2INFO("Total number of mumu events " << evtsMuMuTemp.size());
428 B2INFO("Number of hadronic B events " << evtsHad.size());
429
430 B2INFO("Number of main calibration blocks " << evtsMuMuBlocks.size());
431
432 std::ofstream BonlyOut("BonlyEcmsCalib.txt");
433 BonlyOut << "t1 t2 exp1 run1 exp2 run2 Ecms EcmsUnc spread spreadUnc" << std::endl;
434
435 std::vector<std::vector<CalibrationData>> CalResultsBlocks;
436 //calibrate each run-type block separately
437 for (const auto& evtsMuMu : evtsMuMuBlocks) {
438
439 // Run the mumuCalibration
440 const std::vector<CalibrationData> mumuCalResults = runMuMuCalibration(evtsMuMu,
441 InvariantMassMuMuCalib::runMuMuInvariantMassAnalysis,
443
444
445 //if off-res block or hadB is turned off
446 if (!m_runHadB || evtsMuMu[0].is4S == false) {
447 CalResultsBlocks.push_back(mumuCalResults);
448 continue;
449 }
450
451
452 // If its 4S type and hadB allowed, run the combined Calib
453 auto combCalResults = mumuCalResults;
454 //Update the calibration with Bhad data
455 for (unsigned i = 0; i < mumuCalResults.size(); ++i) {
456 const auto& calibMuMu = mumuCalResults[i];
457 auto& calibComb = combCalResults[i];
458
459 std::vector<std::pair<double, double>> limits, mumuVals;
460
461 for (unsigned j = 0; j < calibMuMu.subIntervals.size(); ++j) {
462 double s = calibMuMu.subIntervals[j].begin()->second.first;
463 double e = calibMuMu.subIntervals[j].rbegin()->second.second;
464 limits.push_back({s, e});
465
466 double val = calibMuMu.pars.cnt[j](0);
467 double unc = calibMuMu.pars.cntUnc[j](0, 0);
468 mumuVals.push_back({val, unc});
469 }
470
471
473 // run B mesons stand-alone calibration
475 auto resB = InvariantMassBhadCalib::doBhadOnlyFit(evtsHad, limits);
476
477 // store it to file
478 double s = calibMuMu.subIntervals.front().begin()->second.first; // start time in hours
479 double e = calibMuMu.subIntervals.back().rbegin()->second.second;// end time in hours
480
481 double exp1 = calibMuMu.subIntervals.front().begin()->first.exp;
482 double exp2 = calibMuMu.subIntervals.back().rbegin()->first.exp;
483 double run1 = calibMuMu.subIntervals.front().begin()->first.run;
484 double run2 = calibMuMu.subIntervals.back().rbegin()->first.run;
485
486 BonlyOut << std::setprecision(8) << s << " " << e << " " << exp1 << " " << run1 << " " << exp2 << " " << run2 << " " <<
487 1e3 * resB[0] << " " << 1e3 * resB[1] << " " << 1e3 * resB[2] << " " << 1e3 * resB[3] << std::endl;
488
489
490
492 //run the combined calibration
494 auto res = InvariantMassBhadCalib::doBhadFit(evtsHad, limits, mumuVals, resB);
495
496 // update calibration data
497 for (unsigned j = 0; j < calibMuMu.subIntervals.size(); ++j) {
498 calibComb.pars.cnt[j](0) = res[j][0];
499 calibComb.pars.cntUnc[j](0, 0) = res[j][1];
500 calibComb.pars.spreadMat(0, 0) = res[j][2];
501
502 calibComb.pars.spreadUnc = res[j][3];
503 calibComb.pars.shift = res[j][4];
504 calibComb.pars.shiftUnc = res[j][5];
505 calibComb.pars.pulls[j] = res[j][6];
506 }
507
508
509 }
510
511 CalResultsBlocks.push_back(combCalResults);
512
513 }
514
515 assert(CalResultsBlocks.size() == evtsMuMuBlocks.size());
516
517 CalResultsBlocks = adjustOffResonanceEnergy(CalResultsBlocks, evtsMuMuBlocks);
518
519
520 std::vector<CalibrationData> CalResults;
521 for (auto cb : CalResultsBlocks)
522 for (auto c : cb)
523 CalResults.push_back(c);
524
525 //Store info to the text file
526 printToFile(CalResults, "finalEcmsCalib.txt");
527
528
529 storePayloadsNoIntraRun(CalResults, "CollisionInvariantMass", getInvariantMassObj);
530
531
533}
Base class for calibration algorithms.
void setDescription(const std::string &description)
Set algorithm description (in constructor)
EResult
The result of calibration.
@ c_OK
Finished successfully =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
void setInputFileNames(PyObject *inputFileNames)
Set the input file names used for this algorithm from a Python list.
void setPrefix(const std::string &prefix)
Set the prefix used to identify datastore objects.
void fillRunToInputFilesMap()
Fill the mapping of ExpRun -> Files.
std::vector< std::string > getVecInputFileNames() const
Get the input file names used for this algorithm as a STL vector.
void clearCalibrationData()
Clear calibration data.
This class contains the measured average center-of-mass energy, which is equal to the invariant mass ...
bool m_runHadB
Run the calibration from had-B decays.
TString m_lossFunctionOuter
Outer loss function (for calibration intervals with constant InvarinatMass spread)
InvariantMassAlgorithm()
Constructor set the prefix to BoostVectorCollector.
std::vector< std::vector< CalibrationData > > adjustOffResonanceEnergy(std::vector< std::vector< CalibrationData > > CalResultsBlocks, const std::vector< std::vector< InvariantMassMuMuCalib::Event > > &evtsMuMuBlocks)
Adjust the energy of the off-resonance runs based on the energy offset between mumu and hadB method i...
double m_eCMSmumuShift
Shift between the energy from the mumu events and the real value.
std::vector< InvariantMassMuMuCalib::Event > getDataMuMu(const std::vector< std::string > &files, bool is4S)
Load the mumu data from files.
std::vector< InvariantMassBhadCalib::Event > getDataHadB(const std::vector< std::string > &files)
Load the hadB data from files.
virtual EResult calibrate() override
Run algo on data.
TString m_lossFunctionInner
Inner loss function (for calibration subintervals with constant InvariantMass)
double m_eCMSmumuSpread
Energy spread for mumu only run (m_runHadB == false)
Class that allows to split runs into the intervals of intended properties given by the lossFunction.
Definition: Splitter.h:108
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.
Definition: Splitter.h:163
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
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 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.
Definition: calibTools.h:290
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: Splitter.cc:38
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].
Definition: Splitter.h:312
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
Abstract base class for different kinds of events.
Parameters and data relevant for single calibration interval.
Definition: calibTools.h:87
structure containing variables relevant for the hadronic B decays