Belle II Software development
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
17using namespace Belle2;
18
19
20
21//-----------------------------------------------------------------
22// Register the Module
23//-----------------------------------------------------------------
24REG_MODULE(EcmsCollector);
25
26//-----------------------------------------------------------------
27// Implementation
28//-----------------------------------------------------------------
29
31{
32 //Set module properties
33
34 setDescription("Collect data for eCMS calibration algorithm using the momenta of the hadronic events");
36}
37
39{
40 B2INFO("Init of the trees");
41 TString objectName = "events";
42 //Data object creation --------------------------------------------------
43 TTree* tree = new TTree(objectName, "");
44
45 tree->Branch<int>("event", &m_evt);
46 tree->Branch<int>("exp", &m_exp);
47 tree->Branch<int>("run", &m_run);
48 tree->Branch<double>("time", &m_time);
49
50 tree->Branch<std::vector<double>>("pBcms", &m_pBcms);
51 tree->Branch<std::vector<double>>("mB", &m_mB);
52
53 tree->Branch<std::vector<int>>("pdg", &m_pdg);
54 tree->Branch<std::vector<int>>("mode", &m_mode);
55 tree->Branch<std::vector<double>>("Kpid", &m_Kpid);
56 tree->Branch<std::vector<double>>("R2", &m_R2);
57 tree->Branch<std::vector<double>>("mD", &m_mD);
58 tree->Branch<std::vector<double>>("dmDstar", &m_dmDstar);
59
60
61 // We register the objects so that our framework knows about them.
62 // Don't try and hold onto the pointers or fill these objects directly
63 // Use the getObjectPtr functions to access collector objects
64 registerObject<TTree>(objectName.Data(), tree);
65}
66
67
69{
70 m_pBcms.resize(n);
71 m_mB.resize(n);
72 m_pdg.resize(n);
73 m_mode.resize(n);
74 m_Kpid.resize(n);
75 m_R2.resize(n);
76 m_mD.resize(n);
77 m_dmDstar.resize(n);
78}
79
80
81
83{
84 // store event info
85 m_evt = m_emd->getEvent();
86 m_run = m_emd->getRun();
87 m_exp = m_emd->getExperiment();
88 m_time = m_emd->getTime() / 1e9 / 3600.; //from ns to hours
89
90
91 StoreObjPtr<ParticleList> B0("B0:merged");
92 StoreObjPtr<ParticleList> Bm("B-:merged");
93
94
95 //put all the B candidates into the vector
96 std::vector<const Particle*> Bparts;
97
98 if (B0.isValid()) {
99 for (unsigned i = 0; i < B0->getListSize(); ++i)
100 if (B0->getParticle(i))
101 Bparts.push_back(B0->getParticle(i));
102 }
103 if (Bm.isValid()) {
104 for (unsigned i = 0; i < Bm->getListSize(); ++i)
105 if (Bm->getParticle(i))
106 Bparts.push_back(Bm->getParticle(i));
107 }
108
109 if (Bparts.size() == 0) return;
110
111
112 resize(Bparts.size());
113
114 for (unsigned i = 0; i < Bparts.size(); ++i) {
115 const Particle* Bpart = Bparts[i];
116
117 //Convert mBC and deltaE to the Y4S reference
119 m_mB[i] = Bpart->getMass();
120 m_pdg[i] = Bpart->getPDGCode();
121 m_mode[i] = Bpart->getExtraInfo("decayModeID");
122 m_R2[i] = Variable::R2(Bpart);
123
124
125 const Particle* D = nullptr;
126 m_dmDstar[i] = -99;
127
128
129 static const int c_PdgD0 = abs(EvtGenDatabasePDG::Instance()->GetParticle("D0")->PdgCode());
130 static const int c_PdgDplus = abs(EvtGenDatabasePDG::Instance()->GetParticle("D+")->PdgCode());
131
132 static const int c_PdgDstar0 = abs(EvtGenDatabasePDG::Instance()->GetParticle("D*0")->PdgCode());
133 static const int c_PdgDstarPlus = abs(EvtGenDatabasePDG::Instance()->GetParticle("D*+")->PdgCode());
134
135 //if D0 or D+ meson
136 if (abs(Bpart->getDaughter(0)->getPDGCode()) == c_PdgD0 || abs(Bpart->getDaughter(0)->getPDGCode()) == c_PdgDplus) {
137 D = Bpart->getDaughter(0);
138 } else if (abs(Bpart->getDaughter(0)->getPDGCode()) == c_PdgDstar0 || abs(Bpart->getDaughter(0)->getPDGCode()) == c_PdgDstarPlus) {
139 const Particle* Dstar = Bpart->getDaughter(0);
140 D = Dstar->getDaughter(0);
141 m_dmDstar[i] = Dstar->getMass() - D->getMass();
142 } else {
143 B2INFO("No D meson found");
144 }
145 if (D != nullptr) {
146 m_mD[i] = D->getMass();
147 const Particle* Kaon = D->getDaughter(0);
148
149 m_Kpid[i] = -99;
150 if (Kaon && Kaon->getPIDLikelihood()) {
152 }
153 } else {
154 m_mD[i] = -99;
155 m_Kpid[i] = -99;
156 }
157 }
158
159
160 getObjectPtr<TTree>("events")->Fill();
161
162}
Calibration collector module base class.
StoreObjPtr< EventMetaData > m_emd
Current EventMetaData.
static const ChargedStable pion
charged pion particle
Definition: Const.h:661
static const ChargedStable kaon
charged kaon particle
Definition: Const.h:662
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.
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
static ROOT::Math::PxPyPzMVector labToCms(const ROOT::Math::PxPyPzMVector &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...
Class to store reconstructed particles.
Definition: Particle.h:75
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:871
int getPDGCode(void) const
Returns PDG code.
Definition: Particle.h:454
ROOT::Math::PxPyPzEVector get4Vector() const
Returns Lorentz vector.
Definition: Particle.h:547
const Particle * getDaughter(unsigned i) const
Returns a pointer to the i-th daughter particle.
Definition: Particle.cc:631
double getExtraInfo(const std::string &name) const
Return given value if set.
Definition: Particle.cc:1289
double getMass() const
Returns invariant mass (= nominal for FS particles)
Definition: Particle.h:507
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
bool isValid() const
Check whether the object was created.
Definition: StoreObjPtr.h:111
#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.