Belle II Software  release-05-01-25
CDCDedxDQM.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2012 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Jitendra Kumar, Jake Bennett
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <reconstruction/modules/CDCDedxDQM/CDCDedxDQM.h>
12 #include <framework/core/HistoModule.h>
13 
14 #include <TDirectory.h>
15 
16 using namespace Belle2;
17 
18 REG_MODULE(CDCDedxDQM)
19 
20 //---------------------------------
22 {
23  setPropertyFlags(c_ParallelProcessingCertified); // parallel processing
24  setDescription("Make data quality monitoring plots for dE/dx: means and resolutions for bhabha skim, dedx band plots for bhabha/hadron skim.");
25 }
26 
27 //---------------------------------
29 
30 
31 //---------------------------------
33 {
34 
35  TDirectory* oldDir = gDirectory;
36  oldDir->mkdir("CDCDedx")->cd();
37 
38  StoreObjPtr<EventMetaData> eventMetaDataPtr;
39  fCurrentEventNum = eventMetaDataPtr->getRun();
40 
41  temp1D = new TH1D("hdEdx_PerRun", Form("hdEdx_PerRun%d", fCurrentEventNum), 300, 0., 3.);
42  temp1D->GetXaxis()->SetTitle("dEdx trucMean of bhabha tracks");
43  temp1D->GetYaxis()->SetTitle("Entries");
44 
45  temp2D = new TH2D("hdEdxVsP_PerRun", Form("hdEdxVsP_PerRun%d", fCurrentEventNum), 400, 0.02, 8.0, 500, 0.10, 15.0);
46  temp2D->GetXaxis()->SetTitle("p (GeV) of hadron tracks");
47  temp2D->GetYaxis()->SetTitle("dEdx");
48 
49  oldDir->cd();
50 
51 }
52 
53 
54 //---------------------------------
56 {
57 
58  if (!m_cdcDedxTracks.isOptional()) {
59  B2WARNING("Missing CDCDedxTracks array, CDCDedxDQM is skipped.");
60  return;
61  }
62 
63  m_TrgResult.isOptional();
64  m_cdcDedxTracks.isRequired();
65  REG_HISTOGRAM
66 }
67 
68 
69 
70 //---------------------------------
72 {
73  if (!m_cdcDedxTracks.isOptional()) {
74  B2WARNING("Missing CDCDedxTracks array, CDCDedxDQM is skipped.");
75  return;
76  }
77 
78  temp1D->Reset();
79  temp2D->Reset();
80 
81 }
82 
83 
84 //---------------------------------
86 {
87 
88  if (!m_TrgResult.isValid()) {
89  B2WARNING("SoftwareTriggerResult object not available but require to select skim events for this module: going back");
90  return;
91  }
92 
93  const std::map<std::string, int>& fresults = m_TrgResult->getResults();
94  if (fresults.find("software_trigger_cut&skim&accept_bhabha") == fresults.end()
95  and fresults.find("software_trigger_cut&skim&accept_hadron") == fresults.end()) {
96  //B2WARNING("CDCDedxDQMModule: Can't find required trigger identifiers: going back ");
97  return;
98  }
99 
100  const bool IsBhabhaEvtAccepted = (m_TrgResult->getResult("software_trigger_cut&skim&accept_bhabha") ==
102  const bool IsHadronEvtAccepted = (m_TrgResult->getResult("software_trigger_cut&skim&accept_hadron") ==
104  if (!IsBhabhaEvtAccepted and !IsHadronEvtAccepted) {
105  //B2WARNING("CDCDedxDQMModule: not an bhabha or hadron event: going back");
106  return;
107  }
108  // fill histograms for bhabha-events only
109  if (!m_cdcDedxTracks.isOptional()) {
110  B2WARNING("Missing CDCDedxTracks array, CDCDedxDQM is skipped.");
111  return;
112  }
113 
114  for (Int_t idedx = 0; idedx < m_cdcDedxTracks.getEntries(); idedx++) {
115 
116  CDCDedxTrack* dedxTrack = m_cdcDedxTracks[idedx];
117  //per run
118  if (IsBhabhaEvtAccepted)temp1D->Fill(float(dedxTrack->getDedxNoSat()));
119  if (IsHadronEvtAccepted)temp2D->Fill(float(dedxTrack->getMomentum()), float(dedxTrack->getDedxNoSat()));
120  }
121 }
122 
123 
124 //---------------------------------
126 {
127 
128 }
129 
130 
131 //---------------------------------
133 {
134 
135 }
Belle2::CDCDedxDQMModule::fCurrentEventNum
Int_t fCurrentEventNum
variable to get run number
Definition: CDCDedxDQM.h:81
Belle2::CDCDedxTrack::getDedxNoSat
double getDedxNoSat() const
Get dE/dx truncated mean without the saturation correction for this track.
Definition: CDCDedxTrack.h:117
Belle2::CDCDedxDQMModule::endRun
virtual void endRun() override
This method is called at the end of each run.
Definition: CDCDedxDQM.cc:125
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::CDCDedxTrack::getMomentum
double getMomentum() const
Return the track momentum valid in the CDC.
Definition: CDCDedxTrack.h:129
Belle2::SoftwareTriggerCutResult::c_accept
@ c_accept
Accept this event.
Belle2::CDCDedxDQMModule::initialize
virtual void initialize() override
Initialize the module.
Definition: CDCDedxDQM.cc:55
Belle2::CDCDedxDQMModule::temp2D
TH2D * temp2D
Dedx vs P histogram per run.
Definition: CDCDedxDQM.h:84
Belle2::CDCDedxDQMModule::terminate
virtual void terminate() override
End of the event processing.
Definition: CDCDedxDQM.cc:132
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::CDCDedxDQMModule::event
virtual void event() override
This method is called for each event.
Definition: CDCDedxDQM.cc:85
Belle2::CDCDedxDQMModule::m_cdcDedxTracks
StoreArray< CDCDedxTrack > m_cdcDedxTracks
Store array for CDCDedxTrack.
Definition: CDCDedxDQM.h:79
Belle2::CDCDedxDQMModule::~CDCDedxDQMModule
virtual ~CDCDedxDQMModule()
Destructor.
Definition: CDCDedxDQM.cc:28
Belle2::CDCDedxDQMModule::beginRun
virtual void beginRun() override
This method is called for each run.
Definition: CDCDedxDQM.cc:71
Belle2::CDCDedxDQMModule::defineHisto
virtual void defineHisto() override
Defination of histograms.
Definition: CDCDedxDQM.cc:32
Belle2::CDCDedxDQMModule
This module to design collect CDC dEdx monitoring for DQM and only minimal information are stored.
Definition: CDCDedxDQM.h:45
Belle2::CDCDedxDQMModule::temp1D
TH1D * temp1D
Dedx histogram per run.
Definition: CDCDedxDQM.h:83
Belle2::CDCDedxDQMModule::m_TrgResult
StoreObjPtr< SoftwareTriggerResult > m_TrgResult
Store array for Trigger selection.
Definition: CDCDedxDQM.h:78
Belle2::CDCDedxTrack
Debug output for CDCDedxPID module.
Definition: CDCDedxTrack.h:36
Belle2::HistoModule
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29