Belle II Software  release-08-01-10
KLMDQM2Module.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 header. */
10 #include <klm/modules/KLMDQM2/KLMDQM2Module.h>
11 
12 /* Basf2 headers. */
13 #include <mdst/dataobjects/Track.h>
14 
15 /* ROOT headers. */
16 #include <TDirectory.h>
17 
18 /* CLHEP headers. */
19 #include <CLHEP/Vector/ThreeVector.h>
20 
21 /* C++ headers. */
22 #include <string>
23 
24 using namespace Belle2;
25 
26 //-----------------------------------------------------------------
27 // Register the Module
28 //-----------------------------------------------------------------
29 REG_MODULE(KLMDQM2);
30 
31 //-----------------------------------------------------------------
32 // Implementation
33 //-----------------------------------------------------------------
34 
36  HistoModule(),
37  m_ElementNumbers(&(KLMElementNumbers::Instance())),
38  m_eklmElementNumbers(&(EKLMElementNumbers::Instance())),
39  m_PlaneArrayIndex(&(KLMPlaneArrayIndex::Instance())),
40  m_GeometryBKLM{nullptr},
41  m_MatchedHitsBKLM{nullptr},
42  m_AllExtHitsBKLM{nullptr},
43  m_MatchedHitsEKLM{nullptr},
44  m_AllExtHitsEKLM{nullptr},
45  m_MatchedHitsBKLMSector{nullptr},
46  m_AllExtHitsBKLMSector{nullptr},
47  m_MatchedHitsEKLMSector{nullptr},
48  m_AllExtHitsEKLMSector{nullptr}
49 {
50  // Set module properties
51  setDescription(R"DOC(Additional Module for KLMDQM plots after HLT filters
52 
53  An additional module developed to display plane efficiencies for the KLM during runs (i.e. for online analyses).
54  This module would be called after HLT filter in order to use mumu-tight skim to select reasonable events.
55  The output histograms would be plane efficiencies = MatchedDigits/AllExtits.
56  )DOC");
57 
58  // Parameter definitions
59  addParam("MuonListName", m_MuonListName, "Muon list name.",
60  std::string("mu+:all"));
61  addParam("AllowedDistance1D", m_AllowedDistance1D,
62  "Maximal distance in the units of strip number from ExtHit to "
63  "matching KLMDigit (not for multi-strip hits).", double(8));
64  addParam("MinimalMatchingDigits", m_MinimalMatchingDigits,
65  "Minimal number of matching digits.", 0);
66  addParam("MinimalMatchingDigitsOuterLayers",
68  "Minimal number of matching digits in outer layers.", 0);
69  addParam("MinimalMomentumNoOuterLayers", m_MinimalMomentumNoOuterLayers,
70  "Minimal momentum in case there are no hits in outer layers.", 0.0);
71  addParam("RemoveUnusedMuons", m_RemoveUnusedMuons,
72  "Whether to remove unused muons.", false);
73  addParam("IgnoreBackwardPropagation", m_IgnoreBackwardPropagation,
74  "Whether to ignore ExtHits with backward propagation.", false);
76  addParam("histogramDirectoryName", m_HistogramDirectoryName,
77  "Directory for KLM DQM histograms in ROOT file.",
78  std::string("KLMEfficiencyDQM"));
79  addParam("SoftwareTriggerName", m_SoftwareTriggerName,
80  "Software Trigger for event selection",
81  std::string("software_trigger_cut&skim&accept_mumutight"));
82 
83 }
84 
86 {
87 }
88 
90 {
91 
92  TDirectory* newDirectory{gDirectory->mkdir(m_HistogramDirectoryName.c_str())};
93  TDirectory::TContext context{gDirectory, newDirectory};
94 
95 
98 
99 
100  /* Number of hits per channel. */
101  /* KLM General Related. */
102  m_MatchedHitsBKLM = new TH1F("matched_hitsBKLM",
103  "Matched Hits in BKLM Plane",
104  m_PlaneNumBKLM, 0.5, 0.5 + m_PlaneNumBKLM);
105  m_MatchedHitsBKLM->GetXaxis()->SetTitle("Layer Number");
106 
107  m_AllExtHitsBKLM = new TH1F("all_ext_hitsBKLM",
108  "All ExtHits in BKLM Plane",
109  m_PlaneNumBKLM, 0.5, 0.5 + m_PlaneNumBKLM);
110  m_AllExtHitsBKLM->GetXaxis()->SetTitle("Layer number");
111 
112 
113 
114  m_MatchedHitsEKLM = new TH1F("matched_hitsEKLM",
115  "Matched Hits in EKLM Plane",
116  m_PlaneNumEKLM, 0.5, m_PlaneNumEKLM + 0.5);
117  m_MatchedHitsEKLM->GetXaxis()->SetTitle("Plane number");
118 
119  m_AllExtHitsEKLM = new TH1F("all_ext_hitsEKLM",
120  "All ExtHits in EKLM Plane",
121  m_PlaneNumEKLM, 0.5, m_PlaneNumEKLM + 0.5);
122  m_AllExtHitsEKLM->GetXaxis()->SetTitle("Plane number");
123 
124 
125 
126  /********************/
127  /********************/
128  /* binned by sector */
129  /********************/
130  /********************/
131 
132  m_MatchedHitsBKLMSector = new TH1F("matched_hitsBKLMSector",
133  "Matched Hits in BKLM Sector",
134  BKLMMaxSectors, 0.5, 0.5 + BKLMMaxSectors);
135  m_MatchedHitsBKLMSector->GetXaxis()->SetTitle("Sector Number");
136 
137  m_AllExtHitsBKLMSector = new TH1F("all_ext_hitsBKLMSector",
138  "All ExtHits in BKLM Sector",
139  BKLMMaxSectors, 0.5, 0.5 + BKLMMaxSectors);
140  m_AllExtHitsBKLMSector->GetXaxis()->SetTitle("Sector number");
141 
142 
143 
144  m_MatchedHitsEKLMSector = new TH1F("matched_hitsEKLMSector",
145  "Matched Hits in EKLM Sector",
146  EKLMMaxSectors, 0.5, EKLMMaxSectors + 0.5);
147  m_MatchedHitsEKLMSector->GetXaxis()->SetTitle("Sector number");
148 
149  m_AllExtHitsEKLMSector = new TH1F("all_ext_hitsEKLMSector",
150  "All ExtHits in EKLM Sector",
151  EKLMMaxSectors, 0.5, EKLMMaxSectors + 0.5);
152  m_AllExtHitsEKLMSector->GetXaxis()->SetTitle("Sector number");
153 
154 }//end of defineHisto
155 
157 {
158  REG_HISTOGRAM;
159  //inputs
160  m_softwareTriggerResult.isOptional();
161  m_MuonList.isRequired(m_MuonListName);
162  m_Digits.isOptional();
164 }
165 
167 {
168  //start by restarting histograms
169 
170  /* KLM General Related. */
171  m_MatchedHitsBKLM->Reset();
172  m_AllExtHitsBKLM->Reset();
173  m_MatchedHitsEKLM->Reset();
174  m_AllExtHitsEKLM->Reset();
175 
176  m_MatchedHitsBKLMSector->Reset();
177  m_AllExtHitsBKLMSector->Reset();
178  m_MatchedHitsEKLMSector->Reset();
179  m_AllExtHitsEKLMSector->Reset();
180 }
181 
183 {
184  if (triggerFlag() || m_SoftwareTriggerName == "") {
185  unsigned int nMuons = m_MuonList->getListSize();
186  for (unsigned int i = 0; i < nMuons; ++i) {
187  const Particle* muon = m_MuonList->getParticle(i);
192 
193  }
194 
195  }
196 }
197 
199 {
200 }
201 
203 {
204 }
205 
207 {
208 
209  bool passed = false;
211  try {
213  true : false;
214  } catch (const std::out_of_range&) {
215  passed = false;
216  }
217  }
218  return passed;
219 
220 }
221 
223  struct HitData* hitData)
224 {
225  for (const KLMDigit& digit : m_Digits) {
226  if (!(digit.getSubdetector() == hitData->subdetector &&
227  digit.getSection() == hitData->section &&
228  digit.getLayer() == hitData->layer &&
229  digit.getSector() == hitData->sector &&
230  digit.getPlane() == hitData->plane))
231  continue;
232 
233  // Defining quantities for distance cut
234  double stripPosition = digit.getStrip();
235  double allowedDistance1D = m_AllowedDistance1D;
236 
237  if (digit.isMultiStrip()) {
238  // Due to a firmware bug, we have to be wary with the allowed distance...
239  stripPosition = 0.5 * (digit.getLastStrip() + digit.getStrip());
240  allowedDistance1D *= (digit.getLastStrip() - digit.getStrip() + 1);
241  }
242  if (fabs(stripPosition - hitData->strip) < allowedDistance1D) {
243  hitData->digit = &digit;
244  return;
245  }
246  }
247 }
248 
250  std::map<KLMPlaneNumber, struct HitData>& hitMap,
251  KLMPlaneNumber planeGlobal, struct HitData* hitData)
252 {
253  std::map<KLMPlaneNumber, struct HitData>::iterator it;
254  it = hitMap.find(planeGlobal);
255  /*
256  * There may be more than one such hit e.g. if track crosses the edge
257  * of the strips or WLS fiber groove. Select only one hit per plane.
258  */
259  if (it == hitMap.end()) {
260  hitMap.insert(std::pair<KLMPlaneNumber, struct HitData>(
261  planeGlobal, *hitData));
262  }
263 }
264 
266  const Particle* muon, TH1F* matchedHitsBKLM, TH1F* allExtHitsBKLM,
267  TH1F* matchedHitsEKLM, TH1F* allExtHitsEKLM, TH1F* matchedHitsBKLMSec, TH1F* allExtHitsBKLMSec,
268  TH1F* matchedHitsEKLMSec, TH1F* allExtHitsEKLMSec)
269 {
270  const int nExtrapolationLayers =
272  const Track* track = muon->getTrack();
273  RelationVector<ExtHit> extHits = track->getRelationsTo<ExtHit>();
274  std::map<KLMPlaneNumber, struct HitData> selectedHits;
275  std::map<KLMPlaneNumber, struct HitData>::iterator it;
276  KLMChannelNumber channel;
278  struct HitData hitData, hitDataPrevious;
279  ROOT::Math::XYZVector extHitPosition;
280  CLHEP::Hep3Vector extHitPositionCLHEP, localPosition;
281  int layer;
282  int extHitLayer[nExtrapolationLayers] = {0};
283  int digitLayer[nExtrapolationLayers] = {0};
284  // initialize hitDataPrevious components
285  hitDataPrevious.subdetector = -1;
286  hitDataPrevious.section = -1;
287  hitDataPrevious.sector = -1;
288  hitDataPrevious.layer = -1;
289  for (const ExtHit& hit : extHits) {
290  /*
291  * Choose hits that exit the sensitive volume.
292  * It is not possible to use entry hits because of a bug in Geant4:
293  * the step does not always have a correct status (fGeomBoundary),
294  * and, consequently, ExtHits are not created.
295  */
296  if (hit.getStatus() != EXT_EXIT)
297  continue;
298  /*
299  * Ignore ExtHits with backward propagation. This affects cosmic events
300  * only. The removal of hits with backward propagation is normally
301  * not needed, however, it is added because of backward error propagation
302  * bug in Geant4 10.6.
303  */
305  if (hit.isBackwardPropagated())
306  continue;
307  }
308  KLMPlaneNumber planeGlobal = 0;
309  hitData.hit = &hit;
310  hitData.digit = nullptr;
311  if (hit.getDetectorID() == Const::EDetector::EKLM) {
312  int stripGlobal = hit.getCopyID();
315  stripGlobal, &hitData.section, &hitData.layer, &hitData.sector,
316  &hitData.plane, &hitData.strip);
318  hitData.section, hitData.sector, hitData.layer,
319  hitData.plane, hitData.strip);
320  status = m_ChannelStatus->getChannelStatus(channel);
321  if (status == KLMChannelStatus::c_Unknown)
322  B2FATAL("Incomplete KLM channel status data.");
323  if (status == KLMChannelStatus::c_Normal) {
325  hitData.subdetector, hitData.layer);
326  extHitLayer[layer - 1]++;
327  planeGlobal = m_ElementNumbers->planeNumberEKLM(
328  hitData.section, hitData.sector, hitData.layer,
329  hitData.plane);
330  addHit(selectedHits, planeGlobal, &hitData);
331  }
332  } else if (hit.getDetectorID() == Const::EDetector::BKLM) {
333  int moduleNumber = hit.getCopyID();
336  moduleNumber, &hitData.section, &hitData.sector, &hitData.layer);
338  /*
339  * For scintillators, the plane and strip numbers are recorded
340  * in the copy number.
341  */
343  moduleNumber, &hitData.section, &hitData.sector, &hitData.layer,
344  &hitData.plane, &hitData.strip);
346  hitData.section, hitData.sector, hitData.layer,
347  hitData.plane, hitData.strip);
348  status = m_ChannelStatus->getChannelStatus(channel);
349  if (status == KLMChannelStatus::c_Unknown)
350  B2FATAL("Incomplete KLM channel status data.");
351  if (status == KLMChannelStatus::c_Normal) {
353  hitData.subdetector, hitData.layer);
354  extHitLayer[layer - 1]++;
355  planeGlobal = m_ElementNumbers->planeNumberBKLM(
356  hitData.section, hitData.sector, hitData.layer,
357  hitData.plane);
358  addHit(selectedHits, planeGlobal, &hitData);
359  }
360  } else {
361  /* For RPCs, the sensitive volume corresponds to both readout planes. */
362  extHitPosition = hit.getPosition();
363  extHitPositionCLHEP.setX(extHitPosition.X());
364  extHitPositionCLHEP.setY(extHitPosition.Y());
365  extHitPositionCLHEP.setZ(extHitPosition.Z());
366  const bklm::Module* module =
367  m_GeometryBKLM->findModule(hitData.section, hitData.sector,
368  hitData.layer);
369  localPosition = module->globalToLocal(extHitPositionCLHEP);
371  hitData.strip = module->getZStrip(localPosition);
372  /*
373  * FIXME:
374  * There are 2 hits per module in RPCs, but the plane information is
375  * not available in ExtHit. For now, 2 entries are created (one for
376  * each plane) for the first hit, and the second one is removed.
377  */
378  if ((hitData.subdetector == hitDataPrevious.subdetector) &&
379  (hitData.section == hitDataPrevious.section) &&
380  (hitData.sector == hitDataPrevious.sector) &&
381  (hitData.layer == hitDataPrevious.layer))
382  continue;
383  std::memcpy(&hitDataPrevious, &hitData, sizeof(struct HitData));
384  /* The returned strip may be out of the valid range. */
386  hitData.section, hitData.sector, hitData.layer, hitData.plane,
387  hitData.strip, false)) {
389  hitData.section, hitData.sector, hitData.layer,
390  hitData.plane, hitData.strip);
391  status = m_ChannelStatus->getChannelStatus(channel);
392  if (status == KLMChannelStatus::c_Unknown)
393  B2FATAL("Incomplete KLM channel status data.");
394  if (status == KLMChannelStatus::c_Normal) {
396  hitData.subdetector, hitData.layer);
397  extHitLayer[layer - 1]++;
398  hitData.localPosition = localPosition.z();
399  planeGlobal = m_ElementNumbers->planeNumberBKLM(
400  hitData.section, hitData.sector, hitData.layer,
401  hitData.plane);
402  addHit(selectedHits, planeGlobal, &hitData);
403  }
404  }
406  hitData.strip = module->getPhiStrip(localPosition);
407  /* The returned strip may be out of the valid range. */
409  hitData.section, hitData.sector, hitData.layer, hitData.plane,
410  hitData.strip, false)) {
412  hitData.section, hitData.sector, hitData.layer,
413  hitData.plane, hitData.strip);
414  status = m_ChannelStatus->getChannelStatus(channel);
415  if (status == KLMChannelStatus::c_Unknown)
416  B2FATAL("Incomplete KLM channel status data.");
417  if (status == KLMChannelStatus::c_Normal) {
419  hitData.subdetector, hitData.layer);
420  extHitLayer[layer - 1]++;
421  hitData.localPosition = localPosition.y();
422  planeGlobal = m_ElementNumbers->planeNumberBKLM(
423  hitData.section, hitData.sector, hitData.layer,
424  hitData.plane);
425  addHit(selectedHits, planeGlobal, &hitData);
426  } //end of channel status check
427  } //end of channel number check
428  }//end of detector condition
429  } else
430  continue;
431  }
432  /* Find matching digits. */
433  int nDigits = 0;
434  for (it = selectedHits.begin(); it != selectedHits.end(); ++it) {
435  findMatchingDigit(&(it->second));
436  if (it->second.digit != nullptr) {
437  nDigits++;
439  it->second.subdetector, it->second.layer);
440  digitLayer[layer - 1]++;
441  }
442  }
443  if (nDigits < m_MinimalMatchingDigits)
444  return false;
445  /* Write efficiency histograms */
446  for (it = selectedHits.begin(); it != selectedHits.end(); ++it) {
447  int matchingDigits = 0;
448  int matchingDigitsOuterLayers = 0;
449  int extHitsOuterLayers = 0;
451  it->second.subdetector, it->second.layer) - 1;
452  for (int i = 0; i < nExtrapolationLayers; ++i) {
453  if (i != layer)
454  matchingDigits += digitLayer[i];
455  if (i > layer) {
456  matchingDigitsOuterLayers += digitLayer[i];
457  extHitsOuterLayers += extHitLayer[i];
458  }
459  }
460  /* Check the number of matching digits in other layers. */
461  if (matchingDigits < m_MinimalMatchingDigits)
462  continue;
463  /*
464  * Check the number of matching digits in outer layers relatively to
465  * this hit.
466  */
467  if (matchingDigitsOuterLayers < m_MinimalMatchingDigitsOuterLayers) {
473  if (extHitsOuterLayers >= m_MinimalMatchingDigitsOuterLayers)
474  continue;
475  /*
476  * If the number of ExtHits is insufficient, then check the momentum.
477  * The muons with sufficiently large momentum have a very small
478  * probability to get absorbed in the detector.
479  */
480  if (muon->getP() < m_MinimalMomentumNoOuterLayers)
481  continue;
482  }
483  //Filling AddExtHits and MatchedHits histograms
484  if (it->second.subdetector == KLMElementNumbers::c_EKLM) {
485  int planeNum = m_eklmElementNumbers->planeNumber(it->second.section, it->second.layer, it->second.sector, it->second.plane);
486  int sectorNum = (it->second.section - 1) * EKLMElementNumbers::getMaximalSectorNumber() + it->second.sector;
487  allExtHitsEKLM->Fill(planeNum);
488  allExtHitsEKLMSec->Fill(sectorNum);
489  if (it->second.digit) {
490  matchedHitsEKLM->Fill(planeNum);
491  matchedHitsEKLMSec->Fill(sectorNum);
492  }
493  }//end of if loop
494 
495 
496  else if (it->second.subdetector == KLMElementNumbers::c_BKLM) {
497  int layerGlobal = BKLMElementNumbers::layerGlobalNumber(
498  it->second.section, it->second.sector, it->second.layer);
499  int sectorGlobal = it->second.section * BKLMElementNumbers::getMaximalSectorNumber() + (it->second.sector);
500  allExtHitsBKLM->Fill(layerGlobal);
501  allExtHitsBKLMSec->Fill(sectorGlobal);
502  if (it->second.digit) {
503  matchedHitsBKLM->Fill(layerGlobal);
504  matchedHitsBKLMSec->Fill(sectorGlobal);
505  }
506  } else {
507  B2FATAL("Had a hit that didn't come from BKLM or EKLM?");
508  }
509 
510  } //end of selectedHits for loop
511  return true;
512 } //end of collectTrackData
static void channelNumberToElementNumbers(KLMChannelNumber channel, int *section, int *sector, int *layer, int *plane, int *strip)
Get element numbers by channel number.
static bool checkChannelNumber(int section, int sector, int layer, int plane, int strip, bool fatalError=true)
Check channel number.
@ c_FirstRPCLayer
First RPC layer.
static constexpr int getMaximalSectorNumber()
Get maximal sector number (1-based).
static void moduleNumberToElementNumbers(KLMModuleNumber module, int *section, int *sector, int *layer)
Get element numbers by module number.
static constexpr int getMaximalSectorGlobalNumber()
Get maximal sector global number.
static int layerGlobalNumber(int section, int sector, int layer)
Get layer global number.
EKLM element numbers.
void stripNumberToElementNumbers(int stripGlobal, int *section, int *layer, int *sector, int *plane, int *strip) const
Get element numbers by strip global number.
static constexpr int getMaximalSectorNumber()
Get maximal sector number.
int planeNumber(int section, int layer, int sector, int plane) const
Get plane number.
static constexpr int getMaximalSectorGlobalNumberKLMOrder()
Get maximal sector global number with KLM ordering (section, sector).
Store one Ext hit as a ROOT object.
Definition: ExtHit.h:32
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
ChannelStatus
Channel status.
@ c_Normal
Normally operating channel.
@ c_Unknown
Unknown status (no data).
TH1F * m_MatchedHitsBKLMSector
Matched hits in sector for BKLM.
TH1F * m_MatchedHitsEKLM
Matched hits in plane for EKLM.
const bklm::GeometryPar * m_GeometryBKLM
BKLM geometry.
void addHit(std::map< KLMPlaneNumber, struct HitData > &hitMap, KLMPlaneNumber planeGlobal, struct HitData *hitData)
Add hit to map.
StoreArray< KLMDigit > m_Digits
KLM digits.
DBObjPtr< KLMChannelStatus > m_ChannelStatus
Channel status.
TH1F * m_AllExtHitsBKLMSector
Extrapolated hits in sector for BKLM.
virtual void initialize() override
Register input and output data.
double m_MinimalMomentumNoOuterLayers
Minimal momentum in case there are no hits in outer layers.
bool m_RemoveUnusedMuons
Whether to remove unused muons.
virtual void event() override
Selection for mumu_tight_skim, then DQM plot filling
const KLMElementNumbers * m_ElementNumbers
KLM element numbers.
virtual void endRun() override
Called if the current run ends.
const EKLMElementNumbers * m_eklmElementNumbers
EKLM element numbers.
virtual void terminate() override
Called at the end of the event processing.
void findMatchingDigit(struct HitData *hitData)
Find matching digit.
std::string m_SoftwareTriggerName
Software Trigger Name.
TH1F * m_AllExtHitsEKLMSector
Extrapolated hits in sector for EKLM.
TH1F * m_MatchedHitsBKLM
Matched hits in plane for BKLM.
StoreObjPtr< ParticleList > m_MuonList
Muons.
virtual void beginRun() override
Called when entering a new run.
bool m_IgnoreBackwardPropagation
Whether to ignore ExtHits with backward propagation.
double m_AllowedDistance1D
Maximal distance in the units of strip number from ExtHit to matching KLMDigit.
KLMDQM2Module()
Constructor: Sets the description, the properties and the parameters of the module.
~KLMDQM2Module()
Destructor.
TH1F * m_AllExtHitsBKLM
Extrapolated hits in plane for BKLM.
TH1F * m_MatchedHitsEKLMSector
Matched hits in sector for EKLM.
const int m_PlaneNumBKLM
Number of layers/planes for BKLM.
int m_MinimalMatchingDigitsOuterLayers
Minimal number of matching digits in outer layers.
std::string m_MuonListName
Muon list name.
int m_MinimalMatchingDigits
Minimal number of matching digits.
StoreObjPtr< SoftwareTriggerResult > m_softwareTriggerResult
Trigger Information.
const int m_PlaneNumEKLM
Number of layers/planes for EKLM.
TH1F * m_AllExtHitsEKLM
Extrapolated hits in plane for EKLM.
bool collectDataTrack(const Particle *muon, TH1F *matchedHitsBKLM, TH1F *allExtHitsBKLM, TH1F *matchedHitsEKLM, TH1F *allExtHitsEKLM, TH1F *matchedHitsBKLMSec, TH1F *allExtHitsBKLMSec, TH1F *matchedHitsEKLMSec, TH1F *allExtHitsEKLMSec)
Collect the data for one muon.
bool triggerFlag()
Uses TrigResult along with desired software cut to determine whether histograms are filled or not for...
std::string m_HistogramDirectoryName
Directory for KLM DQM histograms in ROOT file.
void defineHisto() override
Definition of the histograms.
KLM digit (class representing a digitized hit in RPCs or scintillators).
Definition: KLMDigit.h:29
KLM element numbers.
KLMChannelNumber channelNumberBKLM(int section, int sector, int layer, int plane, int strip) const
Get channel number for BKLM.
KLMPlaneNumber planeNumberEKLM(int section, int sector, int layer, int plane) const
Get channel number for EKLM.
int getExtrapolationLayer(int subdetector, int layer) const
Get extrapolation layer number (BKLM - from 1 to 15, EKLM - from 16 to 29).
KLMPlaneNumber planeNumberBKLM(int section, int sector, int layer, int plane) const
Get plane number for BKLM.
static constexpr int getMaximalExtrapolationLayer()
Get maximal extrapolation layer.
KLMChannelNumber channelNumberEKLM(int section, int sector, int layer, int plane, int strip) const
Get channel number for EKLM.
KLM plane array index.
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
@ 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
Class to store reconstructed particles.
Definition: Particle.h:75
Class for type safe access to objects that are referred to in relations.
Class that bundles various TrackFitResults.
Definition: Track.h:25
const Module * findModule(int section, int sector, int layer) const
Get the pointer to the definition of a module.
Definition: GeometryPar.cc:721
static GeometryPar * instance(void)
Static method to get a reference to the singleton GeometryPar instance.
Definition: GeometryPar.cc:27
Define the geometry of a BKLM module Each sector [octant] contains Modules.
Definition: Module.h:76
REG_MODULE(arichBtest)
Register the Module.
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
uint16_t KLMChannelNumber
Channel number.
uint16_t KLMPlaneNumber
Plane number.
@ c_accept
Accept this event.
Abstract base class for different kinds of events.
double localPosition
Local coordinate.
Definition: KLMDQM2Module.h:76
const KLMDigit * digit
Digit.
Definition: KLMDQM2Module.h:82
const ExtHit * hit
Extrapolation hit.
Definition: KLMDQM2Module.h:79