Belle II Software  release-05-02-19
V0ObjectsDQMModule.cc
1 //+
2 // File : V0ObjectsDQMModule.cc
3 // Description : Module to monitor displaced vertices on HLT
4 //
5 // Author : Bryan Fulsom, PNNL
6 // Date : 2019-01-17
7 //-
8 
9 #include <dqm/modules/V0ObjectsDQM/V0ObjectsDQMModule.h>
10 #include <framework/datastore/StoreObjPtr.h>
11 #include <analysis/dataobjects/ParticleList.h>
12 #include <TDirectory.h>
13 
14 using namespace Belle2;
15 
16 //-----------------------------------------------------------------
17 // Register the Module
18 //-----------------------------------------------------------------
19 REG_MODULE(V0ObjectsDQM)
20 
21 //-----------------------------------------------------------------
22 // Implementation
23 //-----------------------------------------------------------------
24 
25 V0ObjectsDQMModule::V0ObjectsDQMModule() : HistoModule()
26 {
27  //Set module properties
28 
29  setDescription("Monitor displaced vertices");
31 
32  addParam("V0PListName", m_V0PListName, "Name of the vertexed particle list", std::string("K_S0:V0DQM"));
33 }
34 
36 {
37  TDirectory* oldDir = gDirectory;
38  oldDir->mkdir("V0Objects");
39  oldDir->cd("V0Objects");
40 
41  for (int j = 0; j < 32; j++) {
42  m_h_xvsy[j] = new TH2F(Form("xvsy[%i]", j), Form("xvsy[%i]", j), 200, -10, 10, 200, -10, 10);
43  m_h_xvsy[j]->SetXTitle("x [cm]");
44  m_h_xvsy[j]->SetYTitle("y [cm]");
45  m_h_xvsy[j]->SetStats(kFALSE);
46  }
47  m_h_xvsz = new TH2F("xvsz", "xvsz", 1500, -75, 75, 400, -10, 10);
48  m_h_xvsz->SetXTitle("z [cm]");
49  m_h_xvsz->SetYTitle("x [cm]");
50  m_h_xvsz->SetStats(kFALSE);
51 
52  oldDir->cd();
53 }
54 
55 
57 {
58  REG_HISTOGRAM
59 
60 }
61 
62 
64 {
65  for (int j = 0; j < 32; j++) {
66  m_h_xvsy[j]->Reset();
67  }
68  m_h_xvsz->Reset();
69 }
70 
71 
73 {
74 
76 
77  if (V0Particles.isValid()) {
78  for (unsigned int i = 0; i < V0Particles->getListSize(); i++) {
79  Particle* V0 = V0Particles->getParticle(i);
80  //Get the vertex position, fill accordingly
81  float vtxx = V0->getX();
82  float vtxy = V0->getY();
83  float vtxz = V0->getZ();
84  if (fabs(vtxz) < 75 && fabs(vtxx) < 10 && fabs(vtxy) < 10) {
85  m_h_xvsy[int(floor((vtxz + 75) / 5))]->Fill(vtxx, vtxy);
86  if (vtxz <= -5. || vtxz >= 8.) m_h_xvsz->Fill(vtxz, vtxx);
87  }
88  }
89  }
90 
91 }
Belle2::V0ObjectsDQMModule::initialize
void initialize() override final
Function for dynamic initialization of module.
Definition: V0ObjectsDQMModule.cc:56
Belle2::Module::setDescription
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:216
Belle2::V0
Object holding information for V0s.
Definition: V0.h:40
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::V0ObjectsDQMModule::beginRun
void beginRun() override final
Function to process begin_run record.
Definition: V0ObjectsDQMModule.cc:63
Belle2::Module::setPropertyFlags
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:210
Belle2::V0ObjectsDQMModule::m_V0PListName
std::string m_V0PListName
Name of the V0 particle list.
Definition: V0ObjectsDQMModule.h:43
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::V0ObjectsDQMModule::event
void event() override final
Function to process event record.
Definition: V0ObjectsDQMModule.cc:72
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::V0ObjectsDQMModule::defineHisto
void defineHisto() override final
Function to define histograms.
Definition: V0ObjectsDQMModule.cc:35
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