Belle II Software development
CDCDQMModule.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
9// Own include
10#include <cdc/modules/cdcDQM/CDCDQMModule.h>
11
12// Dataobject classes
13#include <framework/database/DBObjPtr.h>
14
15#include <TF1.h>
16#include <TDirectory.h>
17
18#include <fstream>
19#include <math.h>
20
21#include <cdc/dataobjects/WireID.h>
22#include <cdc/geometry/CDCGeometryPar.h>
23
24using namespace std;
25using namespace Belle2;
26using namespace CDC;
27
28//-----------------------------------------------------------------
29// Register module
30//-----------------------------------------------------------------
31
32REG_MODULE(CDCDQM);
33
35{
36 // set module description (e.g. insert text)
37 setDescription("Make summary of data quality.");
38 addParam("MinHits", m_minHits, "Include only events with more than MinHits hits in CDC", 0);
40}
41
43{
44}
45
47{
48
49 TDirectory* oldDir = gDirectory;
50 oldDir->mkdir("CDC");
51 oldDir->cd("CDC");
52 m_hNEvents = new TH1F("hNEvents", "hNEvents", 10, 0, 10);
53 m_hNEvents->GetXaxis()->SetBinLabel(1, "number of events");
54 m_hBit = new TH2F("hBit", "m_hBit", 7, 0, 7.0, 48, 0, 48.0);
55 m_hBit->SetTitle("CDC:Removed Data Bit;CDCRawIndex;Channell Index");
56 m_hOcc = new TH1F("hOcc", "hOccupancy", 150, 0, 1.5);
57 m_hADC = new TH2F("hADC", "hADC", 300, 0, 300, 200, 0, 1000);
58 m_hADC->SetTitle("ADC vs CDC-Boards;Board index;ADC");
59 m_hTDC = new TH2F("hTDC", "hTDC", 300, 0, 300, 1000, 4200, 5200);
60 m_hTDC->SetTitle("TDC vs CDC-Boards;Board index;TDC");
61 m_hHit = new TH2F("hHit", "hHit", 56, 0, 56, 400, 0, 400);
62 m_hHit->SetTitle("CDC-hits;layer index;nhits");
63 m_hPhi = new TH1F("hPhi", "", 360, -180.0, 180.0);
64 m_hPhi->SetTitle("CDC-track-#phi;cdctrack #phi (IP tracks + all events);entries");
65 m_hPhiIndex = new TH2F("hPhiIndex", "", 360, -180.0, 180.0, 8, 0, 8.0);
66 m_hPhiIndex->SetTitle("CDC-track-#phi;cdctrack #phi vs skims;selection-index");
67 m_hPhiEff = new TH2F("hPhiEff", "", 360, -180.0, 180.0, 100, 0, 100.0);
68 m_hPhiEff->SetTitle("CDC-track-#phi;cdctrack #phi vs cdchits;ncdchits");
69 m_hPhiHit = new TH2F("h2HitPhi", "h2HitPhi", 90, -180.0, 180.0, 56, 0, 56);
70 m_hPhiHit->SetTitle("CDC-hits-map (#phi vs layer);Track-#phi;Layer index");
71 oldDir->cd();
72}
73
75{
76 REG_HISTOGRAM
77 m_cdcHits.isOptional();
78 m_cdcRawHits.isOptional();
79 m_trgSummary.isOptional();
80
81 if (!m_Tracks.isOptional()) {
82 B2WARNING("Missing Tracks array");
83 return;
84 }
85}
86
88{
89 if (!m_RecoTracks.isOptional()) {
90 B2DEBUG(22, "Missing recoTracks array in beginRun() ");
91 return;
92 }
93
94 m_hNEvents->Reset();
95 m_hBit->Reset();
96 m_hOcc->Reset();
97 m_hADC->Reset();
98 m_hTDC->Reset();
99 m_hHit->Reset();
100 m_hPhi->Reset();
101 m_hPhiIndex->Reset();
102 m_hPhiEff->Reset();
103 m_hPhiHit->Reset();
104}
105
107{
108 static CDCGeometryPar& cdcgeo = CDCGeometryPar::Instance();
109 const int nWires = 14336;
111 if (!m_trgSummary.isValid() || (m_trgSummary->getTimType() == Belle2::TRGSummary::TTYP_RAND)) {
113 return;
114 }
115
116 if (!m_TrgResult.isValid()) {
117 B2WARNING("SoftwareTriggerResult object not available but require to select bhabha/mumu/hadron events skim");
118 return;
119 }
120
121 const std::map<std::string, int>& fresults = m_TrgResult->getResults();
122 if ((fresults.find("software_trigger_cut&skim&accept_bhabha") == fresults.end()) ||
123 (fresults.find("software_trigger_cut&skim&accept_mumu_tight_or_highm") == fresults.end()) ||
124 (fresults.find("software_trigger_cut&skim&accept_hadron") == fresults.end())) {
125 B2WARNING("CDCDQMModule: Can't find required bhabha or mumu or hadron trigger identifier");
126 return;
127 }
128
129 const bool IsBhabha = (m_TrgResult->getResult("software_trigger_cut&skim&accept_bhabha") ==
131 const bool IsHadron = (m_TrgResult->getResult("software_trigger_cut&skim&accept_hadron") ==
133 const bool IsMumu = (m_TrgResult->getResult("software_trigger_cut&skim&accept_mumu_tight_or_highm") ==
135
136 if (m_cdcHits.getEntries() < m_minHits) {
138 return;
139 }
140
141 m_nEvents += 1;
142 m_hOcc->Fill(static_cast<float>(m_cdcHits.getEntries()) / nWires);
143
144 for (const auto& hit : m_cdcHits) {
145 int lay = hit.getICLayer();
146 int wire = hit.getIWire();
147 m_hHit->Fill(lay, wire);
148 }
149
150 // to record removed databits
151 const int nEntries = m_rawCDCs.getEntries();
152 B2DEBUG(99, "nEntries of RawCDCs : " << nEntries);
153 for (int i = 0; i < nEntries; ++i) {
154 const int nEntriesRawCDC = m_rawCDCs[i]->GetNumEntries();
155 B2DEBUG(99, LogVar("nEntries of rawCDC[i]", nEntriesRawCDC));
156 for (int j = 0; j < nEntriesRawCDC; ++j) {
157 int MaxNumOfCh = m_rawCDCs[i]->GetMaxNumOfCh(j);
158 if (MaxNumOfCh != 4 && MaxNumOfCh != 48) {
159 B2ERROR("CDCDQM: Invalid value of GetMaxNumOfCh");
160 } else if (MaxNumOfCh == 48) {
161 for (int k = 0; k < MaxNumOfCh; ++k) {
162 if (m_rawCDCs[i]->CheckOnlineRemovedDataBit(j, k) == true)m_hBit->SetBinContent(i + 1, k + 1, -0.5);
163 else m_hBit->SetBinContent(i + 1, k + 1, 0.5);
164 }
165 }
166 }
167 }
168
169 // ADC vs layer 2D histogram with only good track related hits
170 double iselect = -1.0;
171
172 for (const auto& b2track : m_Tracks) {
173
174 const Belle2::TrackFitResult* fitresult = b2track.getTrackFitResultWithClosestMass(Const::pion);
175 if (!fitresult) {
176 B2WARNING("No track fit result found.");
177 continue;
178 }
179
180 Belle2::RecoTrack* track = b2track.getRelatedTo<Belle2::RecoTrack>(m_recoTrackArrayName);
181 if (!track) {
182 B2WARNING("Can not access RecoTrack of this Belle2::Track");
183 continue;
184 }
185
186 const genfit::FitStatus* fs = track->getTrackFitStatus();
187 if (!fs) {
188 B2WARNING("Can not access FitStatus of this Track");
189 continue;
190 }
191 // require high NDF track
192 int ndf = fs->getNdf();
193 if (ndf < 20) continue;
194
195 double phiDegree = fitresult->getPhi() / Unit::deg;
196
197 if (fabs(fitresult->getD0()) > 1.0 || fabs(fitresult->getZ0()) > 1.0) {
198 //Off IP tracks
199 m_hPhiIndex->Fill(phiDegree, 0.5); //all skims
200 if (IsBhabha) iselect = 1.5;
201 if (IsHadron) iselect = 2.5;
202 if (IsMumu) iselect = 3.5;
203 m_hPhiIndex->Fill(phiDegree, iselect);
204 } else {
205 //IP tracks
206 m_hPhi->Fill(phiDegree);
207 m_hPhiIndex->Fill(phiDegree, 4.5); //all skims
208 if (IsBhabha) iselect = 5.5;
209 if (IsHadron) iselect = 6.5;
210 if (IsMumu) iselect = 7.5;
211 m_hPhiIndex->Fill(phiDegree, iselect);
212
213 //for tracking efficiency part
214 double nsvdhits = fitresult->getHitPatternVXD().getNSVDHits();
215 double ncdchits = fitresult->getHitPatternCDC().getNHits();
216 if (nsvdhits > 6) {
217 if (ncdchits >= 100)ncdchits = 99.5; //push to last bin
218 m_hPhiEff->Fill(phiDegree, ncdchits);
219 }
220 }
221
222 // Fill histograms of ADC/TDC if hits are associated with track
223 for (const RecoHitInformation::UsedCDCHit* hit : track->getCDCHitList()) {
224
225 const genfit::TrackPoint* tp = track->getCreatedTrackPoint(track->getRecoHitInformation(hit));
226 if (!tp) {
227 B2WARNING("Can not access TrackPoint of this hit");
228 continue;
229 }
230
231 UChar_t lay = hit->getICLayer();
232 UShort_t IWire = hit->getIWire();
233 UShort_t adc = hit->getADCCount();
234 unsigned short tdc = hit->getTDCCount();
235 WireID wireid(lay, IWire);
236 unsigned short bid = cdcgeo.getBoardID(wireid);
237
238 m_hADC->Fill(bid, adc);
239 m_hTDC->Fill(bid, tdc);
240 m_hPhiHit->Fill(fitresult->getPhi() / Unit::deg, lay);
241
242 }
243 }
244}
245
247{
248 m_hNEvents->SetBinContent(1, m_nEvents);
249}
250
252{
253}
std::string m_recoTrackArrayName
Belle2::RecoTrack StoreArray name.
Definition: CDCDQMModule.h:100
StoreArray< RawCDC > m_rawCDCs
Input array for CDC Raw.
Definition: CDCDQMModule.h:96
TH2F * m_hBit
Histogram of online databit removed.
Definition: CDCDQMModule.h:110
void initialize() override
Initialize the Module.
Definition: CDCDQMModule.cc:74
int m_minHits
Minimum hits for processing.
Definition: CDCDQMModule.h:103
void event() override
Event processor.
CDCDQMModule()
Constructor.
Definition: CDCDQMModule.cc:34
StoreArray< CDCRawHit > m_cdcRawHits
CDC raw hits.
Definition: CDCDQMModule.h:95
void endRun() override
End-of-run action.
TH2F * m_hPhiIndex
Histogram of cdc phi of different IP + skims tracks.
Definition: CDCDQMModule.h:112
StoreObjPtr< TRGSummary > m_trgSummary
Trigger summary.
Definition: CDCDQMModule.h:97
void terminate() override
Termination action.
TH1F * m_hNEvents
Histogram of num.
Definition: CDCDQMModule.h:105
TH2F * m_hTDC
Histogram of TDC with track associated hits for all boards (0-299)
Definition: CDCDQMModule.h:108
TH2F * m_hPhiHit
Histogram of track associated hits in phi vs layer
Definition: CDCDQMModule.h:114
void beginRun() override
Called when entering a new run.
Definition: CDCDQMModule.cc:87
TH2F * m_hPhiEff
Histogram of cdc phi of tracking eff.
Definition: CDCDQMModule.h:113
TH1F * m_hPhi
Histogram of cdc phi of IP tracks.
Definition: CDCDQMModule.h:111
StoreArray< CDCHit > m_cdcHits
CDC hits.
Definition: CDCDQMModule.h:94
virtual ~CDCDQMModule()
Destructor.
Definition: CDCDQMModule.cc:42
TH2F * m_hADC
Histogram of ADC with track associated hits for all boards (0-299)
Definition: CDCDQMModule.h:107
TH1F * m_hOcc
Histogram of occupancy.
Definition: CDCDQMModule.h:106
StoreArray< RecoTrack > m_RecoTracks
RecoTracks.
Definition: CDCDQMModule.h:99
StoreArray< Track > m_Tracks
Tracks.
Definition: CDCDQMModule.h:98
StoreObjPtr< SoftwareTriggerResult > m_TrgResult
Store array for Trigger selection.
Definition: CDCDQMModule.h:101
TH2F * m_hHit
Histogram of Hits for all layers (0-55)
Definition: CDCDQMModule.h:109
void defineHisto() override
Histogram definitions.
Definition: CDCDQMModule.cc:46
Long64_t m_nEvents
Number of events processed.
Definition: CDCDQMModule.h:104
Class containing the result of the unpacker in raw data and the result of the digitizer in simulation...
Definition: CDCHit.h:40
The Class for CDC Geometry Parameters.
unsigned short getBoardID(const WireID &wID) const
Returns frontend board id. corresponding to the wire id.
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
static const ChargedStable pion
charged pion particle
Definition: Const.h:661
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
unsigned short getNHits() const
Get the total Number of CDC hits in the fit.
unsigned short getNSVDHits() const
Get total number of hits in the SVD.
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
void setReturnValue(int value)
Sets the return value for this module as integer.
Definition: Module.cc:220
@ 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
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
@ TTYP_RAND
random trigger events
Definition: TRGSummary.h:67
Values of the result of a track fit with a given particle hypothesis.
double getPhi() const
Getter for phi0 with CDF naming convention.
double getD0() const
Getter for d0.
double getZ0() const
Getter for z0.
HitPatternCDC getHitPatternCDC() const
Getter for the hit pattern in the CDC;.
HitPatternVXD getHitPatternVXD() const
Getter for the hit pattern in the VXD;.
static const double deg
degree to radians
Definition: Unit.h:109
Class to identify a wire inside the CDC.
Definition: WireID.h:34
Class to store variables with their name which were sent to the logging service.
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
@ c_accept
Accept this event.
Abstract base class for different kinds of events.
STL namespace.