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>
36static TObject* getInvariantMassObj(VectorXd vMass, MatrixXd vMassUnc, MatrixXd vMassSpread)
40 double mass = vMass(0);
41 double unc = vMassUnc(0, 0);
42 double spread = vMassSpread(0, 0);
44 payload->setMass(mass, unc, spread);
45 TObject* obj =
static_cast<TObject*
>(payload);
58std::vector<CalibrationData> runMuMuCalibration(
const std::vector<Belle2::InvariantMassMuMuCalib::Event>& evts, Fun calibAnalysis,
59 TString lossFunctionOuter, TString lossFunctionInner)
62 std::map<ExpRun, std::pair<double, double>> runsInfoOrg =
getRunInfo(evts);
63 std::map<ExpRun, std::pair<double, double>> runsRemoved;
64 auto runsInfo =
filter(runsInfoOrg, 2. / 60, runsRemoved);
67 if (runsInfo.size() == 0) {
68 B2WARNING(
"Too short run");
74 auto splits = splt.
getIntervals(runsInfo, evts, lossFunctionOuter, lossFunctionInner, 100 );
77 std::vector<CalibrationData> calVec;
78 for (
auto s : splits) {
80 calVec.push_back(calD);
87 for (
auto shortRun : runsRemoved) {
95std::vector<CalibrationData> addSpreadAndOffset(
const std::vector<CalibrationData>& mumuCalResults,
double spread,
double spreadUnc,
96 double offset,
double offsetUnc)
98 std::vector<CalibrationData> mumuCalResultsNew = mumuCalResults;
100 for (
auto& calibNew : mumuCalResultsNew) {
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;
107 calibNew.pars.shift = offset;
108 calibNew.pars.shiftUnc = offsetUnc;
109 calibNew.pars.pulls[j] = 0;
112 return mumuCalResultsNew;
118static std::vector<TTree*> getTrees(TString tag, TFile* f)
123 f->cd(tag +
"/events");
124 TList* l = f->CurrentDirectory().load()->GetListOfKeys();
127 std::vector<TTree*> treeVec;
129 for (
const auto&& obj : *l) {
130 TString n = obj->GetName();
132 TString trDir = tag +
"/events/" + TString(obj->GetName()) +
"/events_1";
134 treeVec.push_back((TTree*) f->Get(trDir));
141static std::vector<InvariantMassMuMuCalib::Event> readMuMuFiles(std::vector<std::string> files,
bool is4S)
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);
148 for (
auto tr : trees) {
150 std::vector<InvariantMassMuMuCalib::Event> eventsNow = InvariantMassMuMuCalib::getEvents(tr, is4S);
151 events.insert(events.end(), eventsNow.begin(), eventsNow.end());
162static std::vector<InvariantMassBhadCalib::Event> readBhadFiles(std::vector<std::string> files)
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);
169 for (
auto tr : trees) {
171 std::vector<InvariantMassBhadCalib::Event> eventsNow = InvariantMassBhadCalib::getEvents(tr);
172 events.insert(events.end(), eventsNow.begin(), eventsNow.end());
190 if (files.size() == 0)
199 TTree* eventsTr = getObjectPtr<TTree>(
"events").get();
202 if (!eventsTr || eventsTr->GetEntries() < 15) {
205 B2WARNING(
"Small number of events in the 4S mumu sample, only " << eventsTr->GetEntries());
207 B2WARNING(
"Small number of events in the off-resonance mumu sample, only " << eventsTr->GetEntries());
209 B2WARNING(
"Pointer to the \"events\" object not defined for the mumu sample");
213 B2INFO(
"Number of mumu events: " << eventsTr->GetEntries());
216 return InvariantMassMuMuCalib::getEvents(eventsTr, is4S);
222 if (files.size() == 0)
230 TTree* eventsTr = getObjectPtr<TTree>(
"events").get();
233 if (!eventsTr || eventsTr->GetEntries() < 15) {
235 B2WARNING(
"Too few data, only " << eventsTr->GetEntries() <<
" events found in hadronic B decay sample");
237 B2WARNING(
"Pointer to the \"events\" object not defined for hadronic B decay sample");
241 B2INFO(
"Number of mumu events: " << eventsTr->GetEntries());
244 return InvariantMassBhadCalib::getEvents(eventsTr);
250double weightAvg(
double x,
double xe,
double y,
double ye)
252 double xe2 = xe * xe;
253 double ye2 = ye * ye;
254 double res = (ye2 * x + xe2 * y) / (xe2 + ye2);
260void printToFile(
const std::vector<CalibrationData>& CalResults, TString outFileName)
263 std::ofstream finalOut(outFileName);
265 "t1 t2 exp1 run1 exp2 run2 state Ecms EcmsUnc pull shift shiftUnc spread spreadUnc" <<
268 for (
auto cal : CalResults) {
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;
276 double shift = cal.pars.shift;
277 double shiftUnc = cal.pars.shiftUnc;
278 double pull = cal.pars.pulls[j];
280 double s = cal.subIntervals[j].begin()->second.first;
281 double e = cal.subIntervals[j].rbegin()->second.second;
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;
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;
304 const std::vector<std::vector<InvariantMassMuMuCalib::Event>>& evtsMuMuBlocks)
307 for (
int b = 0; b < int(CalResultsBlocks.size()); ++b) {
309 if (evtsMuMuBlocks[b][0].is4S)
continue;
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);
326 double spreadUnc = 0.2e-3;
328 double shiftUnc = 0.3e-3;
331 if (shifts.size() == 1) {
333 spreadUnc = spreadsUnc[0];
335 shiftUnc = shiftsUnc[0];
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]);
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);
346 CalResultsBlocks[b] = addSpreadAndOffset(CalResultsBlocks[b], spread, spreadUnc, shift, shiftUnc);
348 return CalResultsBlocks;
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);
372 B2FATAL(
"Unrecognised data type");
377 for (
auto r : filesMuMu4S)
378 B2INFO(
"MuMu4SFile name " << r);
379 for (
auto r : filesMuMuOff) {
380 B2INFO(
"MuMuOffFile name " << r);
382 for (
auto r : filesHad4S)
383 B2INFO(
"Had4SFile name " << r);
388 const auto evtsMuMuOff = readMuMuFiles(filesMuMuOff,
false );
389 const auto evtsMuMu4S = readMuMuFiles(filesMuMu4S,
true );
391 if (evtsMuMuOff.size() + evtsMuMu4S.size() < 50) {
392 B2WARNING(
"Not enough data, there are only " << evtsMuMuOff.size() + evtsMuMu4S.size() <<
" mumu events");
398 const auto evtsHad = readBhadFiles(filesHad4S);
401 std::vector<InvariantMassMuMuCalib::Event> evtsMuMuTemp = evtsMuMu4S;
402 evtsMuMuTemp.insert(evtsMuMuTemp.end(), evtsMuMuOff.begin(), evtsMuMuOff.end());
404 if (evtsMuMuTemp.size() < 50) {
405 B2WARNING(
"Not enough mumu data, there are only " << evtsMuMuTemp.size() <<
" mumu events");
415 std::vector<std::vector<InvariantMassMuMuCalib::Event>> evtsMuMuBlocks;
416 bool is4Sold = evtsMuMuTemp[0].is4S;
417 evtsMuMuBlocks.push_back({});
418 for (
const auto& ev : evtsMuMuTemp) {
419 if (is4Sold != ev.is4S) {
420 evtsMuMuBlocks.push_back({});
423 evtsMuMuBlocks.back().push_back(ev);
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());
430 B2INFO(
"Number of main calibration blocks " << evtsMuMuBlocks.size());
432 std::ofstream BonlyOut(
"BonlyEcmsCalib.txt");
433 BonlyOut <<
"t1 t2 exp1 run1 exp2 run2 Ecms EcmsUnc spread spreadUnc" << std::endl;
435 std::vector<std::vector<CalibrationData>> CalResultsBlocks;
437 for (
const auto& evtsMuMu : evtsMuMuBlocks) {
440 const std::vector<CalibrationData> mumuCalResults = runMuMuCalibration(evtsMuMu,
441 InvariantMassMuMuCalib::runMuMuInvariantMassAnalysis,
446 if (!
m_runHadB || evtsMuMu[0].is4S ==
false) {
447 CalResultsBlocks.push_back(mumuCalResults);
453 auto combCalResults = mumuCalResults;
455 for (
unsigned i = 0; i < mumuCalResults.size(); ++i) {
456 const auto& calibMuMu = mumuCalResults[i];
457 auto& calibComb = combCalResults[i];
459 std::vector<std::pair<double, double>> limits, mumuVals;
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});
466 double val = calibMuMu.pars.cnt[j](0);
467 double unc = calibMuMu.pars.cntUnc[j](0, 0);
468 mumuVals.push_back({val, unc});
475 auto resB = InvariantMassBhadCalib::doBhadOnlyFit(evtsHad, limits);
478 double s = calibMuMu.subIntervals.front().begin()->second.first;
479 double e = calibMuMu.subIntervals.back().rbegin()->second.second;
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;
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;
494 auto res = InvariantMassBhadCalib::doBhadFit(evtsHad, limits, mumuVals, resB);
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];
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];
511 CalResultsBlocks.push_back(combCalResults);
515 assert(CalResultsBlocks.size() == evtsMuMuBlocks.size());
520 std::vector<CalibrationData> CalResults;
521 for (
auto cb : CalResultsBlocks)
523 CalResults.push_back(c);
526 printToFile(CalResults,
"finalEcmsCalib.txt");
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.
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.
double sqrt(double a)
sqrt for double
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.
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.
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::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].
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
void extrapolateCalibration(std::vector< CalibrationData > &calVec)
Extrapolate calibration to intervals where it failed.
Abstract base class for different kinds of events.
Parameters and data relevant for single calibration interval.
structure containing variables relevant for the hadronic B decays
Event containing two tracks.