Belle II Software  release-08-01-05
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 // CDC
13 
14 // Dataobject classes
15 #include <framework/database/DBObjPtr.h>
16 
17 #include <TF1.h>
18 #include <TVector3.h>
19 #include <TDirectory.h>
20 
21 #include <fstream>
22 #include <math.h>
23 
24 #include <cdc/dataobjects/WireID.h>
25 #include <cdc/geometry/CDCGeometryPar.h>
26 
27 using namespace std;
28 using namespace Belle2;
29 using namespace CDC;
30 
31 //-----------------------------------------------------------------
32 // Register module
33 //-----------------------------------------------------------------
34 
35 REG_MODULE(CDCDQM);
36 
37 CDCDQMModule::CDCDQMModule() : HistoModule()
38 {
39  // set module description (e.g. insert text)
40  setDescription("Make summary of data quality.");
41  addParam("MinHits", m_minHits, "Include only events with more than MinHits hits in ARICH", 0);
43 }
44 
46 {
47 }
48 
50 {
51 
52  TDirectory* oldDir = gDirectory;
53  oldDir->mkdir("CDC");
54  oldDir->cd("CDC");
55  m_hNEvents = new TH1F("hNEvents", "hNEvents", 10, 0, 10);
56  m_hNEvents->GetXaxis()->SetBinLabel(1, "number of events");
57  m_hBit = new TH2F("hBit", "m_hBit", 7, 0, 7.0, 48, 0, 48.0);
58  m_hBit->SetTitle("CDC:Removed Data Bit;CDCRawIndex;Channell Index");
59  m_hOcc = new TH1F("hOcc", "hOccupancy", 150, 0, 1.5);
60  m_hPhi = new TH1F("hPhi", "cdc:#phi distribution", 360, -180.0, 180.0);
61  m_hPhi->SetTitle("CDC-track-phi;track-phi;entries");
62  m_hADC = new TH2F("hADC", "hADC", 300, 0, 300, 200, 0, 1000);
63  m_hADC->SetTitle("ADC vs CDC-Boards;Board index;ADC");
64  m_hTDC = new TH2F("hTDC", "hTDC", 300, 0, 300, 1000, 4200, 5200);
65  m_hTDC->SetTitle("TDC vs CDC-Boards;Board index;TDC");
66  m_hHit = new TH2F("hHit", "hHit", 56, 0, 56, 400, 0, 400);
67  m_hHit->SetTitle("CDC-hits;layer index;nhits");
68  m_h2HitPhi = new TH2F("h2HitPhi", "h2HitPhi", 90, -180.0, 180.0, 56, 0, 56);
69  m_h2HitPhi->SetTitle("CDC-hits-map (#phi vs layer);Track-#phi;Layer index");
70 
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_hPhi->Reset();
98  m_hADC->Reset();
99  m_hTDC->Reset();
100  m_hHit->Reset();
101  m_h2HitPhi->Reset();
102 }
103 
105 {
106  static CDCGeometryPar& cdcgeo = CDCGeometryPar::Instance();
107  const int nWires = 14336;
108  setReturnValue(1);
109  if (!m_trgSummary.isValid() || (m_trgSummary->getTimType() == Belle2::TRGSummary::TTYP_RAND)) {
110  setReturnValue(0);
111  return;
112  }
113 
114  if (m_cdcHits.getEntries() < m_minHits) {
115  setReturnValue(0);
116  return;
117  }
118 
119  m_nEvents += 1;
120  m_hOcc->Fill(static_cast<float>(m_cdcHits.getEntries()) / nWires);
121 
122  for (const auto& hit : m_cdcHits) {
123  int lay = hit.getICLayer();
124  int wire = hit.getIWire();
125  m_hHit->Fill(lay, wire);
126  }
127 
128  // to record removed databits
129  const int nEntries = m_rawCDCs.getEntries();
130  B2DEBUG(99, "nEntries of RawCDCs : " << nEntries);
131  for (int i = 0; i < nEntries; ++i) {
132  const int nEntriesRawCDC = m_rawCDCs[i]->GetNumEntries();
133  B2DEBUG(99, LogVar("nEntries of rawCDC[i]", nEntriesRawCDC));
134  for (int j = 0; j < nEntriesRawCDC; ++j) {
135  int MaxNumOfCh = m_rawCDCs[i]->GetMaxNumOfCh(j);
136  if (MaxNumOfCh != 4 && MaxNumOfCh != 48) {
137  B2ERROR("CDCDQM: Invalid value of GetMaxNumOfCh");
138  } else if (MaxNumOfCh == 48) {
139  for (int k = 0; k < MaxNumOfCh; ++k) {
140  if (m_rawCDCs[i]->CheckOnlineRemovedDataBit(j, k) == true)m_hBit->SetBinContent(i + 1, k + 1, -0.5);
141  else m_hBit->SetBinContent(i + 1, k + 1, 0.5);
142  }
143  }
144  }
145  }
146 
147  // ADC vs layer 2D histogram with only good track related hits
148  for (const auto& b2track : m_Tracks) {
149 
150  const Belle2::TrackFitResult* fitresult = b2track.getTrackFitResultWithClosestMass(Const::pion);
151  if (!fitresult) {
152  B2WARNING("No track fit result found.");
153  continue;
154  }
155 
156  Belle2::RecoTrack* track = b2track.getRelatedTo<Belle2::RecoTrack>(m_recoTrackArrayName);
157  if (!track) {
158  B2WARNING("Can not access RecoTrack of this Belle2::Track");
159  continue;
160  }
161 
162  const genfit::FitStatus* fs = track->getTrackFitStatus();
163  if (!fs) {
164  B2WARNING("Can not access FitStatus of this Track");
165  continue;
166  }
167  // require high NDF track
168  int ndf = fs->getNdf();
169  if (ndf < 20) continue;
170 
171  m_hPhi->Fill(fitresult->getPhi() / Unit::deg);
172 
173  // Fill histograms of ADC/TDC if hits are associated with track
174  for (const RecoHitInformation::UsedCDCHit* hit : track->getCDCHitList()) {
175 
176  const genfit::TrackPoint* tp = track->getCreatedTrackPoint(track->getRecoHitInformation(hit));
177  if (!tp) {
178  B2WARNING("Can not access TrackPoint of this hit");
179  continue;
180  }
181 
182  UChar_t lay = hit->getICLayer();
183  UShort_t IWire = hit->getIWire();
184  UShort_t adc = hit->getADCCount();
185  unsigned short tdc = hit->getTDCCount();
186  WireID wireid(lay, IWire);
187  unsigned short bid = cdcgeo.getBoardID(wireid);
188 
189  m_hADC->Fill(bid, adc);
190  m_hTDC->Fill(bid, tdc);
191  m_h2HitPhi->Fill(fitresult->getPhi() / Unit::deg, lay);
192  }
193  }
194 }
195 
197 {
198  m_hNEvents->SetBinContent(1, m_nEvents);
199 }
200 
202 {
203 }
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:111
void initialize() override
Initialize the Module.
Definition: CDCDQMModule.cc:74
int m_minHits
Minimum hits for processing.
Definition: CDCDQMModule.h:102
void event() override
Event processor.
StoreArray< CDCRawHit > m_cdcRawHits
CDC raw hits.
Definition: CDCDQMModule.h:95
void endRun() override
End-of-run action.
StoreObjPtr< TRGSummary > m_trgSummary
Trigger summary.
Definition: CDCDQMModule.h:97
void terminate() override
Termination action.
TH1F * m_hNEvents
Histogram of num.
Definition: CDCDQMModule.h:104
TH2F * m_hTDC
Histogram of TDC with track associated hits for all boards (0-299)
Definition: CDCDQMModule.h:108
void beginRun() override
Called when entering a new run.
Definition: CDCDQMModule.cc:87
TH1F * m_hPhi
Histogram of cdc track phi.
Definition: CDCDQMModule.h:106
StoreArray< CDCHit > m_cdcHits
CDC hits.
Definition: CDCDQMModule.h:94
virtual ~CDCDQMModule()
Destructor.
Definition: CDCDQMModule.cc:45
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:105
TH2F * m_h2HitPhi
Histogram of track associated hits in phi vs layer
Definition: CDCDQMModule.h:110
StoreArray< RecoTrack > m_RecoTracks
RecoTracks.
Definition: CDCDQMModule.h:99
StoreArray< Track > m_Tracks
Tracks.
Definition: CDCDQMModule.h:98
TH2F * m_hHit
Histogram of Hits for all layers (0-55)
Definition: CDCDQMModule.h:109
void defineHisto() override
Histogram definitions.
Definition: CDCDQMModule.cc:49
Long64_t m_nEvents
Number of events processed.
Definition: CDCDQMModule.h:103
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:652
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
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.
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.
Class where important numbers and properties of a fit can be stored.
Definition: FitStatus.h:80
double getNdf() const
Get the degrees of freedom of the fit.
Definition: FitStatus.h:122
Object containing AbsMeasurement and AbsFitterInfo objects.
Definition: TrackPoint.h:46
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
Abstract base class for different kinds of events.