Belle II Software  release-05-02-19
PhysicsObjectsDQMModule.cc
1 //+
2 // File : PhysicsObjectsDQMModule.cc
3 // Description : Module to monitor physics objects on HLT
4 //
5 // Author : Boqun Wang, University of Cincinnati
6 // Date : May - 2018
7 //-
8 
9 #include <dqm/modules/PhysicsObjectsDQM/PhysicsObjectsDQMModule.h>
10 #include <analysis/dataobjects/ParticleList.h>
11 #include <analysis/variables/EventShapeVariables.h>
12 #include <framework/datastore/StoreObjPtr.h>
13 #include <mdst/dataobjects/SoftwareTriggerResult.h>
14 #include <TDirectory.h>
15 #include <map>
16 
17 using namespace Belle2;
18 
19 //-----------------------------------------------------------------
20 // Register the Module
21 //-----------------------------------------------------------------
22 REG_MODULE(PhysicsObjectsDQM)
23 
24 //-----------------------------------------------------------------
25 // Implementation
26 //-----------------------------------------------------------------
27 
28 PhysicsObjectsDQMModule::PhysicsObjectsDQMModule() : HistoModule()
29 {
30  //Set module properties
31 
32  setDescription("Monitor Physics Objects Quality");
34 
35  addParam("TriggerIdentifier", m_triggerIdentifier,
36  "Trigger identifier string used to select events for the histograms", std::string("software_trigger_cut&skim&accept_hadron"));
37  addParam("PI0PListName", m_pi0PListName, "Name of the pi0 particle list", std::string("pi0:physDQM"));
38  addParam("KS0PListName", m_ks0PListName, "Name of the KS0 particle list", std::string("K_S0:physDQM"));
39 }
40 
42 {
43  TDirectory* oldDir = gDirectory;
44  oldDir->mkdir("PhysicsObjects")->cd();
45 
46  m_h_mKS0 = new TH1F("mKS0", "KS0 Invariant Mass", 20, 0.48, 0.52);
47  m_h_mKS0->SetXTitle("M(K_{S}^{0}) [GeV]");
48 
49  m_h_mPI0 = new TH1F("mPI0", "pi0 Invariant Mass", 25, 0.10, 0.15);
50  m_h_mPI0->SetXTitle("M(#pi^{0}) [GeV]");
51 
52  m_h_R2 = new TH1F("R2", "Event Level R2", 36, 0, 1.2);
53  m_h_R2->SetXTitle("R2");
54 
55  oldDir->cd();
56 }
57 
58 
60 {
61  REG_HISTOGRAM
62 
64  result.isOptional();
65 }
66 
67 
69 {
70  m_h_mKS0->Reset();
71  m_h_mPI0->Reset();
72  m_h_R2->Reset();
73 }
74 
75 
77 {
78 }
79 
80 
82 {
83 }
84 
85 
87 {
89  if (!result.isValid()) {
90  B2WARNING("SoftwareTriggerResult object not available but needed to select events for the histograms.");
91  return;
92  }
93 
94  const std::map<std::string, int>& results = result->getResults();
95  if (results.find(m_triggerIdentifier) == results.end()) {
96  B2WARNING("PhysicsObjectsDQM: Can't find trigger identifier: " << m_triggerIdentifier);
97  return;
98  }
99 
100  const bool accepted = (result->getResult(m_triggerIdentifier) == SoftwareTriggerCutResult::c_accept);
101 
102  if (accepted == false) return;
103 
106 
107  double R2 = Belle2::Variable::foxWolframR2(nullptr);
108  m_h_R2->Fill(R2);
109 
110  if (pi0Particles.isValid() && abs(pi0Particles->getPDGCode()) == 111) {
111  for (unsigned int i = 0; i < pi0Particles->getListSize(); i++) {
112  Particle* pi0 = pi0Particles->getParticle(i);
113  m_h_mPI0->Fill(pi0->getMass());
114  }
115  }
116  if (ks0Particles.isValid() && abs(ks0Particles->getPDGCode()) == 310) {
117  for (unsigned int i = 0; i < ks0Particles->getListSize(); i++) {
118  Particle* ks0 = ks0Particles->getParticle(i);
119  m_h_mKS0->Fill(ks0->getMass());
120  }
121  }
122 }
Belle2::Module::setDescription
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:216
Belle2::PhysicsObjectsDQMModule::defineHisto
void defineHisto() override
Function to define histograms.
Definition: PhysicsObjectsDQMModule.cc:41
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::Module::c_ParallelProcessingCertified
@ 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:82
Belle2::PhysicsObjectsDQMModule::beginRun
void beginRun() override
Function to process begin_run record.
Definition: PhysicsObjectsDQMModule.cc:68
Belle2::PhysicsObjectsDQMModule::m_h_mPI0
TH1F * m_h_mPI0
PI0 invariant mass.
Definition: PhysicsObjectsDQMModule.h:42
Belle2::SoftwareTriggerCutResult::c_accept
@ c_accept
Accept this event.
Belle2::PhysicsObjectsDQMModule::m_ks0PListName
std::string m_ks0PListName
Name of the KS0 particle list.
Definition: PhysicsObjectsDQMModule.h:54
Belle2::PhysicsObjectsDQMModule::m_h_mKS0
TH1F * m_h_mKS0
KS0 invariant mass.
Definition: PhysicsObjectsDQMModule.h:39
Belle2::PhysicsObjectsDQMModule::m_pi0PListName
std::string m_pi0PListName
Name of the pi0 particle list.
Definition: PhysicsObjectsDQMModule.h:51
Belle2::Module::setPropertyFlags
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:210
Belle2::Particle::getMass
float getMass() const
Returns invariant mass (= nominal for FS particles)
Definition: Particle.h:433
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::StoreObjPtr
Type-safe access to single objects in the data store.
Definition: ParticleList.h:33
Belle2::PhysicsObjectsDQMModule::event
void event() override
Function to process event record.
Definition: PhysicsObjectsDQMModule.cc:86
Belle2::Particle
Class to store reconstructed particles.
Definition: Particle.h:77
Belle2::Module::addParam
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:562
Belle2::PhysicsObjectsDQMModule::m_h_R2
TH1F * m_h_R2
R2.
Definition: PhysicsObjectsDQMModule.h:45
Belle2::PhysicsObjectsDQMModule::initialize
void initialize() override
Function for dynamic initialization of module.
Definition: PhysicsObjectsDQMModule.cc:59
Belle2::PhysicsObjectsDQMModule::m_triggerIdentifier
std::string m_triggerIdentifier
Trigger identifier string used to select events for the histograms.
Definition: PhysicsObjectsDQMModule.h:48
Belle2::PhysicsObjectsDQMModule::terminate
void terminate() override
Function to terminate module.
Definition: PhysicsObjectsDQMModule.cc:81
Belle2::PhysicsObjectsDQMModule::endRun
void endRun() override
Function to process end_run record.
Definition: PhysicsObjectsDQMModule.cc:76
Belle2::HistoModule
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
Belle2::StoreObjPtr::isValid
bool isValid() const
Check whether the object was created.
Definition: StoreObjPtr.h:120