Belle II Software  release-05-02-19
InvariantMassAlgorithm.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2021 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Radek Zlebcik *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <mdst/dbobjects/CollisionInvariantMass.h>
12 #include <tracking/calibration/InvariantMassAlgorithm.h>
13 #include <tracking/calibration/InvariantMassStandAlone.h>
14 #include <tracking/calibration/calibTools.h>
15 
16 #include <Eigen/Dense>
17 
18 using Eigen::VectorXd;
19 using Eigen::MatrixXd;
20 
21 using namespace Belle2;
22 
23 using Belle2::InvariantMassCalib::getEvents;
24 using Belle2::InvariantMassCalib::runInvariantMassAnalysis;
25 
26 //Using boostVector collector for the input
28 {
29  setDescription("Collision invariant mass calibration algorithm");
30 }
31 
32 
34 static TObject* getInvariantMassObj(VectorXd vMass, MatrixXd vMassUnc, MatrixXd vMassSpread)
35 {
36  auto payload = new CollisionInvariantMass();
37 
38  double mass = vMass(0);
39  double unc2 = vMassUnc(0, 0);
40  double spread2 = vMassSpread(0, 0);
41 
42  payload->setMass(mass, sqrt(unc2), sqrt(spread2));
43  TObject* obj = static_cast<TObject*>(payload);
44  return obj;
45 }
46 
47 
48 
49 /* Main calibration method calling dedicated functions */
51 {
52  TTree* tracks = getObjectPtr<TTree>("events").get();
53  return runCalibration(tracks, "CollisionInvariantMass", getEvents,
54  runInvariantMassAnalysis, getInvariantMassObj,
56 }
Belle2::InvariantMassAlgorithm::m_lossFunctionInner
TString m_lossFunctionInner
Inner loss function (for calibration subintervals with constant InvariantMass)
Definition: InvariantMassAlgorithm.h:56
Belle2::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.
Definition: calibTools.h:386
Belle2::InvariantMassAlgorithm::calibrate
virtual EResult calibrate() override
Run algo on data.
Definition: InvariantMassAlgorithm.cc:50
Belle2::InvariantMassAlgorithm::InvariantMassAlgorithm
InvariantMassAlgorithm()
Constructor set the prefix to BoostVectorCollector.
Definition: InvariantMassAlgorithm.cc:27
Belle2::InvariantMassAlgorithm::m_lossFunctionOuter
TString m_lossFunctionOuter
Outer loss function (for calibration intervals with constant InvarinatMass spread)
Definition: InvariantMassAlgorithm.h:53
Belle2::CalibrationAlgorithm::setDescription
void setDescription(const std::string &description)
Set algorithm description (in constructor)
Definition: CalibrationAlgorithm.h:331
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::CollisionInvariantMass
This class contains the measured average center-of-mass energy, which is equal to the invariant mass ...
Definition: CollisionInvariantMass.h:31
Belle2::CalibrationAlgorithm::EResult
EResult
The result of calibration.
Definition: CalibrationAlgorithm.h:50
Belle2::CalibrationAlgorithm
Base class for calibration algorithms.
Definition: CalibrationAlgorithm.h:47