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 <TMath.h>
17#include <TDirectory.h>
18
19#include <fstream>
20#include <math.h>
21#include <set>
22
23#include <cdc/dataobjects/WireID.h>
24#include <cdc/geometry/CDCGeometryPar.h>
25
26using namespace std;
27using namespace Belle2;
28using namespace CDC;
29
30//-----------------------------------------------------------------
31// Register module
32//-----------------------------------------------------------------
33
34REG_MODULE(CDCDQM);
35
37{
38 // set module description (e.g. insert text)
39 setDescription("Make summary of data quality.");
40 addParam("MinHits", m_minHits, "Include only events with more than MinHits hits in CDC", 0);
41 addParam("MinNDF", m_minNdf, "Include only tracks with more than this NDF", 20);
42 addParam("AdjustWireShift", m_adjustWireShift,
43 "If true, gets the corrected phi view of the wires in stereo layers", m_adjustWireShift);
45}
46
48{
49}
50
52{
53 TDirectory* oldDir = gDirectory;
54 oldDir->mkdir("CDC");
55 oldDir->cd("CDC");
56 m_hNEvents = new TH1F("hNEvents", "hNEvents", 10, 0, 10);
57 m_hNEvents->GetXaxis()->SetBinLabel(1, "number of events");
58 m_hBit = new TH2F("hBit", "m_hBit", 7, 0, 7.0, 48, 0, 48.0);
59 m_hBit->SetTitle("CDC:Removed Data Bit;CDCRawIndex;Channell Index");
60 m_hOcc = new TH1F("hOcc", "hOccupancy", 150, 0, 1.5);
61 m_hADC = new TH2F("hADC", "hADC", 300, 0, 300, 200, 0, 1000);
62 m_hADC->SetTitle("ADC vs CDC-Boards;Board index;ADC");
63 m_hTDC = new TH2F("hTDC", "hTDC", 300, 0, 300, 1000, 4200, 5200);
64 m_hTDC->SetTitle("TDC vs CDC-Boards;Board index;TDC");
65 m_hHit = new TH2F("hHit", "hHit", 56, 0, 56, 400, 0, 400);
66 m_hHit->SetTitle("CDC-hits;layer index;nhits");
67 m_hPhi = new TH1F("hPhi", "", 360, -180.0, 180.0);
68 m_hPhi->SetTitle("CDC-track-#phi;cdctrack #phi (IP tracks + all events);entries");
69 m_hPhiIndex = new TH2F("hPhiIndex", "", 360, -180.0, 180.0, 8, 0, 8.0);
70 m_hPhiIndex->SetTitle("CDC-track-#phi;cdctrack #phi vs skims;selection-index");
71 m_hPhiEff = new TH2F("hPhiEff", "", 360, -180.0, 180.0, 100, 0, 100.0);
72 m_hPhiEff->SetTitle("CDC-track-#phi;cdctrack #phi vs cdchits;ncdchits");
73 m_hPhiHit = new TH2F("h2HitPhi", "h2HitPhi", 90, -180.0, 180.0, 56, 0, 56);
74 m_hPhiHit->SetTitle("CDC-hits-map (#phi vs layer);Track-#phi;Layer index");
75 m_hPhiNCDC = new TH2F("hPhiNCDC", "hPhiNCDC", 45, -180.0, 180.0, 61, -0.5, 60.5);
76 m_hPhiNCDC->SetTitle("nCDCHits vs #phi;Track-#phi;nCDCHits;Track / bin");
77 m_hTrackingWireEff = new TH2F("hTrackingWireEff", "title", 400, 0.5, 400 + 0.5, 56 * 2, -0.5, 56 * 2 - 0.5);
78 m_hTrackingWireEff->SetTitle("Attached vs Expected wires for all layers (backplate view);wire;layer;Track / bin");
79 oldDir->cd();
80}
81
83{
84 REG_HISTOGRAM
85 m_cdcHits.isOptional();
86 m_cdcRawHits.isOptional();
87 m_trgSummary.isOptional();
88
89 if (!m_Tracks.isOptional()) {
90 B2WARNING("Missing Tracks array");
91 return;
92 }
93}
94
96{
97 if (!m_RecoTracks.isOptional()) {
98 B2DEBUG(22, "Missing recoTracks array in beginRun() ");
99 return;
100 }
101
102 m_hNEvents->Reset();
103 m_hBit->Reset();
104 m_hOcc->Reset();
105 m_hADC->Reset();
106 m_hTDC->Reset();
107 m_hHit->Reset();
108 m_hPhi->Reset();
109 m_hPhiIndex->Reset();
110 m_hPhiEff->Reset();
111 m_hPhiHit->Reset();
112 m_hPhiNCDC->Reset();
113 m_hTrackingWireEff->Reset();
114}
115
117{
118 static CDCGeometryPar& cdcgeo = CDCGeometryPar::Instance();
119 const int nWires = 14336;
121 if (!m_trgSummary.isValid() || (m_trgSummary->getTimType() == Belle2::TRGSummary::TTYP_RAND)) {
123 return;
124 }
125
126 if (!m_TrgResult.isValid()) {
127 B2WARNING("SoftwareTriggerResult object not available but require to select bhabha/mumu/hadron events skim");
128 return;
129 }
130
131 const std::map<std::string, int>& fresults = m_TrgResult->getResults();
132 if ((fresults.find("software_trigger_cut&skim&accept_bhabha") == fresults.end()) ||
133 (fresults.find("software_trigger_cut&skim&accept_mumu_tight_or_highm") == fresults.end()) ||
134 (fresults.find("software_trigger_cut&skim&accept_hadron") == fresults.end())) {
135 B2WARNING("CDCDQMModule: Can't find required bhabha or mumu or hadron trigger identifier");
136 return;
137 }
138
139 const bool IsBhabha = (m_TrgResult->getResult("software_trigger_cut&skim&accept_bhabha") ==
141 const bool IsHadron = (m_TrgResult->getResult("software_trigger_cut&skim&accept_hadron") ==
143 const bool IsMumu = (m_TrgResult->getResult("software_trigger_cut&skim&accept_mumu_tight_or_highm") ==
145
146 if (m_cdcHits.getEntries() < m_minHits) {
148 return;
149 }
150
151 m_nEvents += 1;
152 m_hOcc->Fill(static_cast<float>(m_cdcHits.getEntries()) / nWires);
153
154 for (const auto& hit : m_cdcHits) {
155 int lay = hit.getICLayer();
156 int wire = hit.getIWire();
157 m_hHit->Fill(lay, wire);
158 }
159
160 // to record removed databits
161 const int nEntries = m_rawCDCs.getEntries();
162 B2DEBUG(99, "nEntries of RawCDCs : " << nEntries);
163 for (int i = 0; i < nEntries; ++i) {
164 const int nEntriesRawCDC = m_rawCDCs[i]->GetNumEntries();
165 B2DEBUG(99, LogVar("nEntries of rawCDC[i]", nEntriesRawCDC));
166 for (int j = 0; j < nEntriesRawCDC; ++j) {
167 int MaxNumOfCh = m_rawCDCs[i]->GetMaxNumOfCh(j);
168 if (MaxNumOfCh != 4 && MaxNumOfCh != 48) {
169 B2ERROR("CDCDQM: Invalid value of GetMaxNumOfCh");
170 } else if (MaxNumOfCh == 48) {
171 for (int k = 0; k < MaxNumOfCh; ++k) {
172 if (m_rawCDCs[i]->CheckOnlineRemovedDataBit(j, k) == true)m_hBit->SetBinContent(i + 1, k + 1, -0.5);
173 else m_hBit->SetBinContent(i + 1, k + 1, 0.5);
174 }
175 }
176 }
177 }
178
179 // ADC vs layer 2D histogram with only good track related hits
180 double iselect = -1.0;
181
182 for (const auto& b2track : m_Tracks) {
183
184 const Belle2::TrackFitResult* fitresult = b2track.getTrackFitResultWithClosestMass(Const::pion);
185 if (!fitresult) {
186 B2WARNING("No track fit result found.");
187 continue;
188 }
189
190 Belle2::RecoTrack* track = b2track.getRelatedTo<Belle2::RecoTrack>(m_recoTrackArrayName);
191 if (!track) {
192 B2WARNING("Can not access RecoTrack of this Belle2::Track");
193 continue;
194 }
195
196 const genfit::FitStatus* fs = track->getTrackFitStatus();
197 if (!fs) {
198 B2WARNING("Can not access FitStatus of this Track");
199 continue;
200 }
201
202 // use tracks from IP for layer effi calculation
203 if (std::fabs(fitresult->getD0()) < 1.0 && std::fabs(fitresult->getZ0()) < 1.0) {
204 std::set<int> hitInSLayer; // list of contributing layers for this track
205 for (const RecoHitInformation::UsedCDCHit* hit : track->getCDCHitList()) {
206 const genfit::TrackPoint* tp = track->getCreatedTrackPoint(track->getRecoHitInformation(hit));
207 if (!tp) continue;
208 hitInSLayer.insert(hit->getICLayer());
209 }
210 if (hitInSLayer.empty()) continue;
211 auto helix = fitresult->getHelix();
212 int nSLayers = cdcgeo.getNumberOfSenseLayers();
213 for (int lay = 0; lay < nSLayers; lay++) {
214 double layerR = cdcgeo.senseWireR(lay);
215 // extrapolate to layer
216 double arcLength = helix.getArcLength2DAtCylindricalR(layerR);
217 if (std::isnan(arcLength)) continue;
218 const auto& result = helix.getPositionAtArcLength2D(arcLength);
219 if (result.Z() > cdcgeo.senseWireFZ(lay) || result.Z() < cdcgeo.senseWireBZ(lay)) continue;
220 // get phi of extrapolation and fill
221 double phi = getShiftedPhi(result, lay);
222 int fillXbin = findPhiBin(phi, lay);
223 m_hTrackingWireEff->Fill(fillXbin, lay);
224 if (hitInSLayer.count(lay)) // if hit is attached in this layer
225 m_hTrackingWireEff->Fill(fillXbin, lay + nSLayers);
226 }
227 }
228
229 double phiDegree = fitresult->getPhi() / Unit::deg;
230 m_hPhiNCDC->Fill(phiDegree, TMath::Min(int(track->getNumberOfCDCHits()), 60));
231
232 // require high NDF track
233 int ndf = fs->getNdf();
234 if (ndf < m_minNdf) continue;
235
236 if (fabs(fitresult->getD0()) > 1.0 || fabs(fitresult->getZ0()) > 1.0) {
237 //Off IP tracks
238 m_hPhiIndex->Fill(phiDegree, 0.5); //all skims
239 if (IsBhabha) iselect = 1.5;
240 if (IsHadron) iselect = 2.5;
241 if (IsMumu) iselect = 3.5;
242 m_hPhiIndex->Fill(phiDegree, iselect);
243 } else {
244 //IP tracks
245 m_hPhi->Fill(phiDegree);
246 m_hPhiIndex->Fill(phiDegree, 4.5); //all skims
247 if (IsBhabha) iselect = 5.5;
248 if (IsHadron) iselect = 6.5;
249 if (IsMumu) iselect = 7.5;
250 m_hPhiIndex->Fill(phiDegree, iselect);
251
252 //for tracking efficiency part
253 double nsvdhits = fitresult->getHitPatternVXD().getNSVDHits();
254 double ncdchits = fitresult->getHitPatternCDC().getNHits();
255 if (nsvdhits > 6) {
256 if (ncdchits >= 100)ncdchits = 99.5; //push to last bin
257 m_hPhiEff->Fill(phiDegree, ncdchits);
258 }
259 }
260
261 // Fill histograms of ADC/TDC if hits are associated with track
262 for (const RecoHitInformation::UsedCDCHit* hit : track->getCDCHitList()) {
263
264 const genfit::TrackPoint* tp = track->getCreatedTrackPoint(track->getRecoHitInformation(hit));
265 if (!tp) {
266 B2WARNING("Can not access TrackPoint of this hit");
267 continue;
268 }
269
270 UChar_t lay = hit->getICLayer();
271 UShort_t IWire = hit->getIWire();
272 UShort_t adc = hit->getADCCount();
273 unsigned short tdc = hit->getTDCCount();
274 WireID wireid(lay, IWire);
275 unsigned short bid = cdcgeo.getBoardID(wireid);
276
277 m_hADC->Fill(bid, adc);
278 m_hTDC->Fill(bid, tdc);
279 m_hPhiHit->Fill(phiDegree, lay);
280 }
281 }
282}
283
285{
286 m_hNEvents->SetBinContent(1, m_nEvents);
287}
288
290{
291}
292
293int CDCDQMModule::findPhiBin(double phi, const int& lay)
294{
295 static CDCGeometryPar& cdcgeo = CDCGeometryPar::Instance();
296 int nWires = cdcgeo.nWiresInLayer(lay);
297 double offset = cdcgeo.offset(lay);
298 double shift = (offset - 0.5) * 2 * TMath::Pi() / nWires;
299 double binLow = shift;
300 double binHigh = 2 * TMath::Pi() + shift;
301 // wrap around
302 if (phi < binLow) phi += 2 * TMath::Pi(); // Wrap around for values less than binLow
303 else if (phi >= binHigh) phi -= 2 * TMath::Pi(); // Wrap around for values >= binHigh
304 double binWidth = (binHigh - binLow) / nWires;
305 int binIndex = (phi - binLow) / binWidth;
306 // Return the bin index (1-based for ROOT histograms)
307 return binIndex + 1;
308}
309
310double CDCDQMModule::getShiftedPhi(const ROOT::Math::XYZVector& position, const int& lay)
311{
312 double phi = TMath::ATan2(position.Y(), position.X());
313 if (m_adjustWireShift) {
314 static CDCGeometryPar& cdcgeo = CDCGeometryPar::Instance();
315 int nShifts = cdcgeo.nShifts(lay);
316 if (nShifts) {
317 int nWires = cdcgeo.nWiresInLayer(lay);
318 double fZ = cdcgeo.senseWireFZ(lay);
319 double bZ = cdcgeo.senseWireBZ(lay);
320 double R = cdcgeo.senseWireR(lay);
321 double phiSize = 2 * TMath::Pi() / nWires;
322 double phiF = phiSize * 0.5 * nShifts;
323 B2Vector3D f(R * TMath::Cos(phiF), R * TMath::Sin(phiF), fZ);
324 B2Vector3D b(R, 0, bZ);
325 B2Vector3D v = f - b;
326 B2Vector3D u = v.Unit();
327 double beta = (position.Z() - b.Z()) / u.Z();
328 B2Vector3D p = b + beta * u;
329 phi -= TMath::ATan2(p.Y(), p.X());
330 }
331 }
332 while (phi < 0) phi += (2 * TMath::Pi());
333 while (phi >= 2 * TMath::Pi()) phi -= (2 * TMath::Pi());
334 return phi;
335}
double R
typedef autogenerated by FFTW
B2Vector3< DataType > Unit() const
Unit vector parallel to this.
Definition: B2Vector3.h:269
std::string m_recoTrackArrayName
Belle2::RecoTrack StoreArray name.
Definition: CDCDQMModule.h:111
StoreArray< RawCDC > m_rawCDCs
Input array for CDC Raw.
Definition: CDCDQMModule.h:107
TH2F * m_hPhiNCDC
Histogram of track associated nCDCHits vs phi.
Definition: CDCDQMModule.h:128
TH2F * m_hBit
Histogram of online databit removed.
Definition: CDCDQMModule.h:123
void initialize() override
Initialize the Module.
Definition: CDCDQMModule.cc:82
int m_minHits
Minimum hits for processing.
Definition: CDCDQMModule.h:114
int m_minNdf
Minimum track NDF for cdc dqm plots.
Definition: CDCDQMModule.h:115
void event() override
Event processor.
CDCDQMModule()
Constructor.
Definition: CDCDQMModule.cc:36
StoreArray< CDCRawHit > m_cdcRawHits
CDC raw hits.
Definition: CDCDQMModule.h:106
void endRun() override
End-of-run action.
TH2F * m_hPhiIndex
Histogram of cdc phi of different IP + skims tracks.
Definition: CDCDQMModule.h:125
StoreObjPtr< TRGSummary > m_trgSummary
Trigger summary.
Definition: CDCDQMModule.h:108
void terminate() override
Termination action.
double getShiftedPhi(const ROOT::Math::XYZVector &position, const int &lay)
Compute and shift phi if it is stereo layer return [0,2pi].
int findPhiBin(double phi, const int &lay)
Find bin corresponds to a specific phi in a layer phi must be in [0,2pi].
TH1F * m_hNEvents
Histogram of num.
Definition: CDCDQMModule.h:118
TH2F * m_hTDC
Histogram of TDC with track associated hits for all boards (0-299)
Definition: CDCDQMModule.h:121
TH2F * m_hPhiHit
Histogram of track associated hits in phi vs layer
Definition: CDCDQMModule.h:127
void beginRun() override
Called when entering a new run.
Definition: CDCDQMModule.cc:95
TH2F * m_hPhiEff
Histogram of cdc phi of tracking eff.
Definition: CDCDQMModule.h:126
bool m_adjustWireShift
If true, gets the correct phi view of the boards.
Definition: CDCDQMModule.h:116
TH1F * m_hPhi
Histogram of cdc phi of IP tracks.
Definition: CDCDQMModule.h:124
StoreArray< CDCHit > m_cdcHits
CDC hits.
Definition: CDCDQMModule.h:105
virtual ~CDCDQMModule()
Destructor.
Definition: CDCDQMModule.cc:47
TH2F * m_hADC
Histogram of ADC with track associated hits for all boards (0-299)
Definition: CDCDQMModule.h:120
TH1F * m_hOcc
Histogram of occupancy.
Definition: CDCDQMModule.h:119
StoreArray< RecoTrack > m_RecoTracks
RecoTracks.
Definition: CDCDQMModule.h:110
StoreArray< Track > m_Tracks
Tracks.
Definition: CDCDQMModule.h:109
StoreObjPtr< SoftwareTriggerResult > m_TrgResult
Store array for Trigger selection.
Definition: CDCDQMModule.h:112
TH2F * m_hHit
Histogram of Hits for all layers (0-55)
Definition: CDCDQMModule.h:122
void defineHisto() override
Histogram definitions.
Definition: CDCDQMModule.cc:51
Long64_t m_nEvents
Number of events processed.
Definition: CDCDQMModule.h:117
TH2F * m_hTrackingWireEff
Histogram of attached and expected CDC wires in layer.
Definition: CDCDQMModule.h:129
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.
int nShifts(int layerId) const
Returns number shift.
double offset(int layerID) const
Return wire offset in phi direction at endplate.
double senseWireBZ(int layerId) const
Returns backward z position of sense wire in each layer.
unsigned nWiresInLayer(int layerId) const
Returns wire numbers in a layer.
double senseWireFZ(int layerId) const
Returns forward z position of sense wire in each layer.
ushort getNumberOfSenseLayers() const
Get the number of sense layers.
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.
double senseWireR(int layerId) const
Returns radius of sense wire in each layer.
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.
Helix getHelix() const
Conversion to framework Helix (without covariance).
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.