Belle II Software  release-06-02-00
EcmsCollectorModule.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 #include <tracking/modules/EcmsCollector/EcmsCollectorModule.h>
9 
10 #include <analysis/dataobjects/ParticleList.h>
11 #include <analysis/variables/ContinuumSuppressionVariables.h>
12 #include <analysis/utility/PCmsLabTransform.h>
13 #include <mdst/dataobjects/PIDLikelihood.h>
14 #include <framework/particledb/EvtGenDatabasePDG.h>
15 
16 
17 using namespace Belle2;
18 using std::vector;
19 
20 
21 
22 //-----------------------------------------------------------------
23 // Register the Module
24 //-----------------------------------------------------------------
25 REG_MODULE(EcmsCollector)
26 
27 //-----------------------------------------------------------------
28 // Implementation
29 //-----------------------------------------------------------------
30 
32 {
33  //Set module properties
34 
35  setDescription("Collect data for eCMS calibration algorithm using the momenta of the hadronic events");
36  setPropertyFlags(c_ParallelProcessingCertified);
37 }
38 
40 {
41  B2INFO("Init of the trees");
42  TString objectName = "events";
43  //Data object creation --------------------------------------------------
44  TTree* tree = new TTree(objectName, "");
45 
46  tree->Branch<int>("event", &m_evt);
47  tree->Branch<int>("exp", &m_exp);
48  tree->Branch<int>("run", &m_run);
49  tree->Branch<double>("time", &m_time);
50 
51  tree->Branch<vector<double>>("pBcms", &m_pBcms);
52  tree->Branch<vector<double>>("mB", &m_mB);
53 
54  tree->Branch<vector<int>>("pdg", &m_pdg);
55  tree->Branch<vector<int>>("mode", &m_mode);
56  tree->Branch<vector<double>>("Kpid", &m_Kpid);
57  tree->Branch<vector<double>>("R2", &m_R2);
58  tree->Branch<vector<double>>("mD", &m_mD);
59  tree->Branch<vector<double>>("dmDstar", &m_dmDstar);
60 
61 
62  // We register the objects so that our framework knows about them.
63  // Don't try and hold onto the pointers or fill these objects directly
64  // Use the getObjectPtr functions to access collector objects
65  registerObject<TTree>(objectName.Data(), tree);
66 }
67 
68 
70 {
71  m_pBcms.resize(n);
72  m_mB.resize(n);
73  m_pdg.resize(n);
74  m_mode.resize(n);
75  m_Kpid.resize(n);
76  m_R2.resize(n);
77  m_mD.resize(n);
78  m_dmDstar.resize(n);
79 }
80 
81 
82 
84 {
85  // store event info
86  m_evt = m_emd->getEvent();
87  m_run = m_emd->getRun();
88  m_exp = m_emd->getExperiment();
89  m_time = m_emd->getTime() / 1e9 / 3600.; //from ns to hours
90 
91 
92  StoreObjPtr<ParticleList> B0("B0:merged");
93  StoreObjPtr<ParticleList> Bm("B-:merged");
94 
95 
96  //put all the B candidates into the vector
97  vector<const Particle*> Bparts;
98 
99  if (B0.isValid()) {
100  for (unsigned i = 0; i < B0->getListSize(); ++i)
101  if (B0->getParticle(i))
102  Bparts.push_back(B0->getParticle(i));
103  }
104  if (Bm.isValid()) {
105  for (unsigned i = 0; i < Bm->getListSize(); ++i)
106  if (Bm->getParticle(i))
107  Bparts.push_back(Bm->getParticle(i));
108  }
109 
110  if (Bparts.size() == 0) return;
111 
112 
113  resize(Bparts.size());
114 
115  for (unsigned i = 0; i < Bparts.size(); ++i) {
116  const Particle* Bpart = Bparts[i];
117 
118  //Convert mBC and deltaE to the Y4S reference
119  m_pBcms[i] = PCmsLabTransform::labToCms(Bpart->get4Vector()).P();
120  m_mB[i] = Bpart->getMass();
121  m_pdg[i] = Bpart->getPDGCode();
122  m_mode[i] = Bpart->getExtraInfo("decayModeID");
123  m_R2[i] = Variable::R2(Bpart);
124 
125 
126  const Particle* D = nullptr;
127  m_dmDstar[i] = -99;
128 
129 
130  static const int c_PdgD0 = abs(EvtGenDatabasePDG::Instance()->GetParticle("D0")->PdgCode());
131  static const int c_PdgDplus = abs(EvtGenDatabasePDG::Instance()->GetParticle("D+")->PdgCode());
132 
133  static const int c_PdgDstar0 = abs(EvtGenDatabasePDG::Instance()->GetParticle("D*0")->PdgCode());
134  static const int c_PdgDstarPlus = abs(EvtGenDatabasePDG::Instance()->GetParticle("D*+")->PdgCode());
135 
136  //if D0 or D+ meson
137  if (abs(Bpart->getDaughter(0)->getPDGCode()) == c_PdgD0 || abs(Bpart->getDaughter(0)->getPDGCode()) == c_PdgDplus) {
138  D = Bpart->getDaughter(0);
139  } else if (abs(Bpart->getDaughter(0)->getPDGCode()) == c_PdgDstar0 || abs(Bpart->getDaughter(0)->getPDGCode()) == c_PdgDstarPlus) {
140  const Particle* Dstar = Bpart->getDaughter(0);
141  D = Dstar->getDaughter(0);
142  m_dmDstar[i] = Dstar->getMass() - D->getMass();
143  } else {
144  B2INFO("No D meson found");
145  }
146  m_mD[i] = D->getMass();
147  const Particle* Kaon = D->getDaughter(0);
148 
149 
150 
151  m_Kpid[i] = -99;
152  if (Kaon && Kaon->getPIDLikelihood()) {
154  }
155  }
156 
157 
158  getObjectPtr<TTree>("events")->Fill();
159 
160 }
Calibration collector module base class.
StoreObjPtr< EventMetaData > m_emd
Current EventMetaData.
static const ChargedStable pion
charged pion particle
Definition: Const.h:542
static const ChargedStable kaon
charged kaon particle
Definition: Const.h:543
Collector for the collision energy calibration based on the hadronic modes.
void prepare() override final
Initialize the module.
std::vector< double > m_dmDstar
D*-D0 mass.
std::vector< double > m_Kpid
Kaon PID.
std::vector< double > m_pBcms
B mesons CMS momentum.
std::vector< int > m_mode
decay mode ID
double m_time
event time [hours]
void resize(int n)
resize the event members to n candidates
std::vector< double > m_mB
B mesons mass.
void collect() override final
Event processor The filling of the tree.
std::vector< double > m_mD
D meson mass.
std::vector< double > m_R2
the R2 variable used for the continuum suppression
std::vector< int > m_pdg
B meson PDG code (can neutral or charged)
static EvtGenDatabasePDG * Instance()
Instance method that loads the EvtGen table.
static TLorentzVector labToCms(const TLorentzVector &vec)
Transforms Lorentz vector into CM System.
double getProbability(const Const::ChargedStable &p1, const Const::ChargedStable &p2, Const::PIDDetectorSet set=Const::PIDDetectorSet::set()) const
Return combined likelihood probability for a particle being p1 and not p2, assuming equal prior proba...
Definition: PIDLikelihood.h:83
Class to store reconstructed particles.
Definition: Particle.h:74
float getExtraInfo(const std::string &name) const
Return given value if set.
Definition: Particle.cc:1242
float getMass() const
Returns invariant mass (= nominal for FS particles)
Definition: Particle.h:445
const PIDLikelihood * getPIDLikelihood() const
Returns the pointer to the PIDLikelihood object that is related to the Track, which was used to creat...
Definition: Particle.cc:840
int getPDGCode(void) const
Returns PDG code.
Definition: Particle.h:392
TLorentzVector get4Vector() const
Returns Lorentz vector.
Definition: Particle.h:477
const Particle * getDaughter(unsigned i) const
Returns a pointer to the i-th daughter particle.
Definition: Particle.cc:639
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:95
bool isValid() const
Check whether the object was created.
Definition: StoreObjPtr.h:110
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.