Belle II Software release-09-00-00
PhysicsObjectsMiraBelleEcmsBBModule.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 <dqm/modules/PhysicsObjectsMiraBelle/PhysicsObjectsMiraBelleEcmsBBModule.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#include <mdst/dataobjects/SoftwareTriggerResult.h>
16#include <framework/datastore/StoreObjPtr.h>
17#include <framework/datastore/StoreArray.h>
18#include <analysis/dataobjects/ParticleList.h>
19
20#include <TTree.h>
21#include <TH1D.h>
22
23
24using namespace Belle2;
25
26
27
28//-----------------------------------------------------------------
29// Register the Module
30//-----------------------------------------------------------------
31REG_MODULE(PhysicsObjectsMiraBelleEcmsBB);
32
33//-----------------------------------------------------------------
34// Implementation
35//-----------------------------------------------------------------
36
38{
39 //Set module properties
40 setDescription("Monitor of the CMS collision energy based on hadronic B decays");
42
43 addParam("TriggerIdentifier", m_triggerIdentifier,
44 "Trigger identifier string used to select events for the histograms", std::string("software_trigger_cut&skim&accept_btocharm"));
45 addParam("BmListName", m_BmListName, "Name of the B- particle list", std::string("B-:combined"));
46 addParam("B0ListName", m_B0ListName, "Name of the B0 particle list", std::string("B0:combined"));
47}
48
49
51{
52 REG_HISTOGRAM
53}
54
56{
57 TDirectory* oldDir = gDirectory;
58 oldDir->mkdir("PhysicsObjectsMiraBelleEcmsBB");
59 oldDir->cd("PhysicsObjectsMiraBelleEcmsBB");
60
61 const double cMBp = EvtGenDatabasePDG::Instance()->GetParticle("B+")->Mass();
62
63 // In the original ML-based analysis there are 40 bins for visualisation
64 // Since binned fit is used here, 80 bins are used instead.
65 // For simplicity, the same histogram lower bound (B+ mass) is used
66 m_hB0 = new TH1D("hB0", "", 80, cMBp, 5.37);
67 m_hBp = new TH1D("hBp", "", 80, cMBp, 5.37);
68 oldDir->cd();
69}
70
72{
73 if (m_hBp) m_hBp->Reset();
74 if (m_hB0) m_hB0->Reset();
75}
76
77
78
79
80
82{
84 if (!result.isValid()) {
85 B2WARNING("SoftwareTriggerResult object not available but needed to select events for the histograms.");
86 return;
87 }
88
89 const bool accepted = (result->getResult(m_triggerIdentifier) == SoftwareTriggerCutResult::c_accept);
90 if (accepted == false) return;
91
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
110 if (Bparts.size() == 0) return;
111
112 for (unsigned i = 0; i < Bparts.size(); ++i) {
113 const Particle* Bpart = Bparts[i];
114
115 const Particle* D = nullptr;
116 double dmDstar = std::numeric_limits<double>::quiet_NaN();
117
118 static const int c_PdgD0 = abs(EvtGenDatabasePDG::Instance()->GetParticle("D0")->PdgCode());
119 static const int c_PdgDplus = abs(EvtGenDatabasePDG::Instance()->GetParticle("D+")->PdgCode());
120
121 static const int c_PdgDstar0 = abs(EvtGenDatabasePDG::Instance()->GetParticle("D*0")->PdgCode());
122 static const int c_PdgDstarPlus = abs(EvtGenDatabasePDG::Instance()->GetParticle("D*+")->PdgCode());
123
124 //if D0 or D+ meson
125 if (abs(Bpart->getDaughter(0)->getPDGCode()) == c_PdgD0 || abs(Bpart->getDaughter(0)->getPDGCode()) == c_PdgDplus) {
126 D = Bpart->getDaughter(0);
127 } else if (abs(Bpart->getDaughter(0)->getPDGCode()) == c_PdgDstar0 || abs(Bpart->getDaughter(0)->getPDGCode()) == c_PdgDstarPlus) {
128 const Particle* Dstar = Bpart->getDaughter(0);
129 D = Dstar->getDaughter(0);
130 dmDstar = Dstar->getMass() - D->getMass();
131 } else {
132 B2INFO("No D meson found");
133 }
134 double mD = D ? D->getMass() : std::numeric_limits<double>::quiet_NaN();
135
136 //Convert mBC and deltaE to the Y4S reference
137 double pBcms = PCmsLabTransform::labToCms(Bpart->get4Vector()).P();
138 double mInv = Bpart->getMass();
139 double pdg = Bpart->getPDGCode();
140 double R2 = Variable::R2(Bpart);
141
142 //get mass of B+- or B0
143 double mB = EvtGenDatabasePDG::Instance()->GetParticle(abs(pdg))->Mass();
144
145
146 // Filling the histograms
147 if (c_mDmin < mD && mD < c_mDmax)
148 if (abs(mInv - mB) < c_mBwindow)
149 if (R2 < c_R2max)
150 if (std::isnan(dmDstar) || (c_dmDstarMin < dmDstar && dmDstar < c_dmDstarMax)) {
151 double eBC = sqrt(pBcms * pBcms + mB * mB); // beam constrained energy
152 if (abs(pdg) == 511) {
153 m_hB0->Fill(eBC);
154 } else {
155 m_hBp->Fill(eBC);
156 }
157 }
158 }
159}
static EvtGenDatabasePDG * Instance()
Instance method that loads the EvtGen table.
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
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.
Class to store reconstructed particles.
Definition: Particle.h:75
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 getMass() const
Returns invariant mass (= nominal for FS particles)
Definition: Particle.h:507
void initialize() override final
Register the histograms.
std::string m_B0ListName
List name for neutral B candidates.
static constexpr double c_mDmax
Maximal value of the D meson inv mass.
void defineHisto() override final
Initialize the histograms.
std::string m_triggerIdentifier
Trigger identifier string used to select events for the histograms.
void event() override final
Event processor Filling of the histograms.
static constexpr double c_mDmin
Minimal value of the D meson inv mass.
static constexpr double c_dmDstarMin
Minimal value of the m(D*)-(mD)
std::string m_BmListName
List name for charged B candidates.
static constexpr double c_dmDstarMax
Maximal value of the m(D*)-(mD)
void beginRun() override final
Reset the histograms.
static constexpr double c_mBwindow
Maximal deviation of B meson inv mass from PDG value.
static constexpr double c_R2max
Maximal allowed R2 value (to suppress continuum)
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
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
@ c_accept
Accept this event.
Abstract base class for different kinds of events.