Belle II Software  release-06-02-00
KLMReconstructorModule.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/KLMReconstructor/KLMReconstructorModule.h>
11 
12 /* KLM headers. */
13 #include <klm/bklm/geometry/Module.h>
14 #include <klm/dataobjects/eklm/EKLMElementNumbers.h>
15 
16 /* Belle 2 headers. */
17 #include <framework/gearbox/Const.h>
18 #include <framework/logging/Logger.h>
19 
20 /* CLHEP headers. */
21 #include <CLHEP/Vector/ThreeVector.h>
22 
23 /* C++ headers. */
24 #include <utility>
25 
26 using namespace Belle2;
27 using namespace Belle2::bklm;
28 
29 REG_MODULE(KLMReconstructor)
30 
31 static bool compareSector(KLMDigit* d1, KLMDigit* d2)
32 {
33  int s1, s2;
34  static const EKLMElementNumbers& elementNumbers =
36  s1 = elementNumbers.sectorNumber(d1->getSection(), d1->getLayer(),
37  d1->getSector());
38  s2 = elementNumbers.sectorNumber(d2->getSection(), d2->getLayer(),
39  d2->getSector());
40  return s1 < s2;
41 }
42 
43 static bool comparePlane(KLMDigit* d1, KLMDigit* d2)
44 {
45  return d1->getPlane() < d2->getPlane();
46 }
47 
48 static bool compareStrip(KLMDigit* d1, KLMDigit* d2)
49 {
50  return d1->getStrip() < d2->getStrip();
51 }
52 
53 static bool compareTime(KLMDigit* d1, KLMDigit* d2)
54 {
55  return d1->getTime() < d2->getTime();
56 }
57 
58 static bool sameSector(KLMDigit* d1, KLMDigit* d2)
59 {
60  return ((d1->getSection() == d2->getSection()) &&
61  (d1->getLayer() == d2->getLayer()) &&
62  (d1->getSector() == d2->getSector()));
63 }
64 
66  Module(),
67  m_CoincidenceWindow(0),
68  m_PromptTime(0),
69  m_PromptWindow(0),
70  m_EventT0Value(0.),
71  m_ElementNumbers(&(KLMElementNumbers::Instance())),
72  m_bklmGeoPar(nullptr),
73  m_eklmElementNumbers(&(EKLMElementNumbers::Instance())),
74  m_eklmGeoDat(nullptr),
75  m_eklmNStrip(0),
76  m_eklmTransformData{nullptr}
77 {
78  setDescription("Create BKLMHit1ds from KLMDigits and then create BKLMHit2ds from BKLMHit1ds; create EKLMHit2ds from KLMDigits.");
80  addParam("TimeCableDelayCorrection", m_TimeCableDelayCorrection,
81  "Perform cable delay time correction (true) or not (false).", true);
82  addParam("EventT0Correction", m_EventT0Correction,
83  "Perform EventT0 correction (true) or not (false)", true);
84  addParam("IfAlign", m_bklmIfAlign,
85  "Perform alignment correction (true) or not (false).",
86  bool(true));
87  addParam("IgnoreScintillators", m_bklmIgnoreScintillators,
88  "Ignore scintillators (to debug their electronics mapping).",
89  false);
90  addParam("CheckSegmentIntersection", m_eklmCheckSegmentIntersection,
91  "Check if segments intersect.", true);
92  addParam("IgnoreHotChannels", m_IgnoreHotChannels,
93  "Use only Normal and Dead (for debugging) channels during 2d hit reconstruction",
94  true);
95 }
96 
98 {
99 }
100 
102 {
103  m_Digits.isRequired();
105  m_EventT0.isRequired();
106  /* BKLM. */
107  m_bklmHit1ds.registerInDataStore();
108  m_bklmHit2ds.registerInDataStore();
109  m_bklmHit1ds.registerRelationTo(m_Digits);
110  m_bklmHit2ds.registerRelationTo(m_bklmHit1ds);
112  /* EKLM. */
113  m_eklmHit2ds.registerInDataStore();
114  m_eklmAlignmentHits.registerInDataStore();
115  m_eklmHit2ds.registerRelationTo(m_Digits);
116  m_eklmAlignmentHits.registerRelationTo(m_eklmHit2ds);
120  if (m_eklmGeoDat->getNPlanes() != 2)
121  B2FATAL("It is not possible to run EKLM reconstruction with 1 plane.");
123 }
124 
126 {
128  if (!m_TimeConstants.isValid())
129  B2FATAL("KLM time constants data are not available.");
130  if (!m_TimeCableDelay.isValid())
131  B2FATAL("KLM time cable decay data are not available.");
132  if (!m_TimeResolution.isValid())
133  B2ERROR("KLM time resolution data are not available. "
134  "The error is non-fatal because the data are only used to set "
135  "chi^2 of 2d hit, which is informational only now.");
136  }
137  if (!m_ChannelStatus.isValid())
138  B2FATAL("KLM channel status data are not available.");
139  if (!m_TimeWindow.isValid())
140  B2FATAL("KLM time window data are not available.");
141  m_CoincidenceWindow = m_TimeWindow->getCoincidenceWindow();
142  m_PromptTime = m_TimeWindow->getPromptTime();
143  m_PromptWindow = m_TimeWindow->getPromptWindow();
151  }
152 }
153 
155 {
156  m_EventT0Value = 0.;
157  if (m_EventT0.isValid())
158  if (m_EventT0->hasEventT0())
159  m_EventT0Value = m_EventT0->getEventT0();
162 }
163 
165 {
166  unsigned int cID = d->getUniqueChannelID();
167  ct -= m_TimeCableDelay->getTimeDelay(cID);
168 }
169 
170 /*
171 double KLMReconstructorModule::getTime(KLMDigit* d, double dist)
172 {
173  int strip;
174  strip = m_eklmElementNumbers->stripNumber(
175  d->getSection(), d->getLayer(), d->getSector(),
176  d->getPlane(), d->getStrip()) - 1;
177  return d->getTime() -
178  (dist / m_eklmTimeCalibration->getEffectiveLightSpeed() +
179  m_eklmTimeCalibrationData[strip]->getTimeShift());
180 
181  // TODO: Subtract time correction given by
182  // m_eklmTimeCalibration->getAmplitudeTimeConstant() / sqrt(d->getNPhotoelectrons()).
183  // It requires a new firmware version that will be able to extract amplitude.
184 
185 }
186 */
187 
189 {
190  int subdetector = digit->getSubdetector();
191  int section = digit->getSection();
192  int sector = digit->getSector();
193  int layer = digit->getLayer();
194  int plane = digit->getPlane();
195  int strip = digit->getStrip();
196  KLMChannelNumber channel = m_ElementNumbers->channelNumber(subdetector, section, sector, layer, plane, strip);
197  enum KLMChannelStatus::ChannelStatus status = m_ChannelStatus->getChannelStatus(channel);
198  if (status == KLMChannelStatus::c_Unknown)
199  B2FATAL("Incomplete KLM channel status data.");
200  if (status == KLMChannelStatus::c_Normal || status == KLMChannelStatus::c_Dead)
201  return true;
202  return false;
203 }
204 
206 {
207  /* Construct BKLMHit1Ds from KLMDigits. */
208  /* Sort KLMDigits by module and strip number. */
209  std::map<KLMChannelNumber, int> channelDigitMap;
210  for (int index = 0; index < m_Digits.getEntries(); ++index) {
211  const KLMDigit* digit = m_Digits[index];
213  continue;
214  if (digit->isMultiStrip())
215  continue;
216  if (m_bklmIgnoreScintillators && !digit->inRPC())
217  continue;
218  if (m_IgnoreHotChannels && !isNormal(digit))
219  continue;
220  if (digit->inRPC() || digit->isGood()) {
222  digit->getSection(), digit->getSector(),
223  digit->getLayer(), digit->getPlane(),
224  digit->getStrip());
225  channelDigitMap.insert(std::pair<KLMChannelNumber, int>(channel, index));
226  }
227  }
228  if (channelDigitMap.empty())
229  return;
230  std::vector<std::pair<const KLMDigit*, double>> digitCluster;
231  KLMChannelNumber previousChannel = channelDigitMap.begin()->first;
232  double averageTime = m_Digits[channelDigitMap.begin()->second]->getTime();
234  correctCableDelay(averageTime, m_Digits[channelDigitMap.begin()->second]);
235  for (std::map<KLMChannelNumber, int>::iterator it = channelDigitMap.begin(); it != channelDigitMap.end(); ++it) {
236  const KLMDigit* digit = m_Digits[it->second];
237  double digitTime = digit->getTime();
239  correctCableDelay(digitTime, digit);
240  if ((it->first > previousChannel + 1) || (std::fabs(digitTime - averageTime) > m_CoincidenceWindow)) {
241  m_bklmHit1ds.appendNew(digitCluster); // Also sets relation BKLMHit1d -> KLMDigit
242  digitCluster.clear();
243  }
244  previousChannel = it->first;
245  double n = (double)(digitCluster.size());
246  averageTime = (n * averageTime + digitTime) / (n + 1.0);
247  digitCluster.emplace_back(std::make_pair(digit, digitTime));
248  }
249  m_bklmHit1ds.appendNew(digitCluster); // Also sets relation BKLMHit1d -> KLMDigit
250 
251  /* Construct BKLMHit2Ds from orthogonal same-module BKLMHit1Ds. */
252  for (int i = 0; i < m_bklmHit1ds.getEntries(); ++i) {
253  int moduleID = m_bklmHit1ds[i]->getModuleID();
254  const bklm::Module* m = m_bklmGeoPar->findModule(m_bklmHit1ds[i]->getSection(), m_bklmHit1ds[i]->getSector(),
255  m_bklmHit1ds[i]->getLayer());
256  bool isPhiReadout = m_bklmHit1ds[i]->isPhiReadout();
257  for (int j = i + 1; j < m_bklmHit1ds.getEntries(); ++j) {
259  moduleID, m_bklmHit1ds[j]->getModuleID()))
260  continue;
261  if (isPhiReadout == m_bklmHit1ds[j]->isPhiReadout())
262  continue;
263  int phiIndex = isPhiReadout ? i : j;
264  int zIndex = isPhiReadout ? j : i;
265  const BKLMHit1d* phiHit = m_bklmHit1ds[phiIndex];
266  const BKLMHit1d* zHit = m_bklmHit1ds[zIndex];
267  CLHEP::Hep3Vector local = m->getLocalPosition(phiHit->getStripAve(), zHit->getStripAve());
268  CLHEP::Hep3Vector propagationDist;
269  if (m_bklmHit1ds[i]->getLayer() < BKLMElementNumbers::c_FirstRPCLayer) {
270  if (isPhiReadout) {
271  propagationDist = m->getPropagationDistance(
272  local, m_bklmHit1ds[j]->getStripMin(),
273  m_bklmHit1ds[i]->getStripMin());
274  } else {
275  propagationDist = m->getPropagationDistance(
276  local, m_bklmHit1ds[i]->getStripMin(),
277  m_bklmHit1ds[j]->getStripMin());
278  }
279  } else {
280  propagationDist = m->getPropagationTimes(local);
281  }
282  double delayPhi, delayZ;
283  if (phiHit->inRPC())
284  delayPhi = m_DelayRPCPhi;
285  else
286  delayPhi = m_DelayBKLMScintillators;
287  if (zHit->inRPC())
288  delayZ = m_DelayRPCZ;
289  else
290  delayZ = m_DelayBKLMScintillators;
291  double phiTime = phiHit->getTime() - propagationDist.y() * delayPhi;
292  double zTime = zHit->getTime() - propagationDist.z() * delayZ;
293  if (std::fabs(phiTime - zTime) > m_CoincidenceWindow)
294  continue;
295  // The second param in localToGlobal is whether do the alignment correction (true) or not (false)
296  CLHEP::Hep3Vector global = m->localToGlobal(local + m->getLocalReconstructionShift(), m_bklmIfAlign);
297  double time = 0.5 * (phiTime + zTime);
299  time -= m_EventT0Value;
300  BKLMHit2d* hit2d = m_bklmHit2ds.appendNew(phiHit, zHit, global, time); // Also sets relation BKLMHit2d -> BKLMHit1d
301  if (std::fabs(time - m_PromptTime) > m_PromptWindow)
302  hit2d->isOutOfTime(true);
303  }
304  }
305 }
306 
307 
309 {
310  int i, n;
311  double d1, d2, time, t1, t2, sd;
312  std::vector<KLMDigit*> digitVector;
313  std::vector<KLMDigit*>::iterator it1, it2, it3, it4, it5, it6, it7, it8, it9;
314  n = m_Digits.getEntries();
315  for (i = 0; i < n; i++) {
316  KLMDigit* digit = m_Digits[i];
318  continue;
319  if (digit->isMultiStrip())
320  continue;
321  if (m_IgnoreHotChannels && !isNormal(digit))
322  continue;
323  if (digit->isGood())
324  digitVector.push_back(digit);
325  }
326  /* Sort by sector. */
327  sort(digitVector.begin(), digitVector.end(), compareSector);
328  it1 = digitVector.begin();
329  while (it1 != digitVector.end()) {
330  it2 = it1;
331  while (1) {
332  ++it2;
333  if (it2 == digitVector.end())
334  break;
335  if (!sameSector(*it1, *it2))
336  break;
337  }
338  /* Now it1 .. it2 - hits in a sector. Sort by plane. */
339  sort(it1, it2, comparePlane);
340  /* If all hits are form the second plane, then continue. */
341  if ((*it1)->getPlane() != 1) {
342  it1 = it2;
343  continue;
344  }
345  it3 = it1;
346  while (1) {
347  ++it3;
348  if (it3 == it2)
349  break;
350  if ((*it3)->getPlane() != (*it1)->getPlane())
351  break;
352  }
353  /*
354  * Now it1 .. it3 - hits from the first plane, it3 .. it2 - hits from the
355  * second plane. If there are no hits from the second plane, then continue.
356  */
357  if (it3 == it2) {
358  it1 = it2;
359  continue;
360  }
361  /* Sort by strip. */
362  sort(it1, it3, compareStrip);
363  sort(it3, it2, compareStrip);
364  it4 = it1;
365  while (it4 != it2) {
366  it5 = it4;
367  while (1) {
368  ++it5;
369  if (it5 == it2)
370  break;
371  /* This loop is for both planes so it is necessary to compare planes. */
372  if ((*it5)->getStrip() != (*it4)->getStrip() ||
373  (*it5)->getPlane() != (*it4)->getPlane())
374  break;
375  }
376  /* Now it4 .. it5 - hits from the same strip. Sort by time. */
377  sort(it4, it5, compareTime);
378  it4 = it5;
379  }
380  /* Strip loop. */
381  it4 = it1;
382  while (it4 != it3) {
383  it5 = it4;
384  while (1) {
385  ++it5;
386  if (it5 == it3)
387  break;
388  if ((*it5)->getStrip() != (*it4)->getStrip())
389  break;
390  }
391  it6 = it3;
392  while (it6 != it2) {
393  it7 = it6;
394  while (1) {
395  ++it7;
396  if (it7 == it2)
397  break;
398  if ((*it7)->getStrip() != (*it6)->getStrip())
399  break;
400  }
401  /*
402  * Now it4 .. it5 - hits from a single fisrt-plane strip amd
403  * it6 .. it7 - hits from a single second-plane strip.
404  * If strips do not intersect, then continue.
405  */
406  HepGeom::Point3D<double> crossPoint(0, 0, 0);
407  if (!m_eklmTransformData->intersection(*it4, *it6, &crossPoint,
408  &d1, &d2, &sd,
410  it6 = it7;
411  continue;
412  }
413  for (it8 = it4; it8 != it5; ++it8) {
414  for (it9 = it6; it9 != it7; ++it9) {
415  t1 = (*it8)->getTime() - d1 * m_DelayEKLMScintillators
416  + 0.5 * sd / Const::speedOfLight;
417  t2 = (*it9)->getTime() - d2 * m_DelayEKLMScintillators
418  - 0.5 * sd / Const::speedOfLight;
420  correctCableDelay(t1, *it8);
421  correctCableDelay(t2, *it9);
422  }
423  if (std::fabs(t1 - t2) > m_CoincidenceWindow)
424  continue;
425  time = (t1 + t2) / 2;
427  time -= m_EventT0Value;
428  EKLMHit2d* hit2d = m_eklmHit2ds.appendNew(*it8);
429  hit2d->setEnergyDeposit((*it8)->getEnergyDeposit() +
430  (*it9)->getEnergyDeposit());
431  hit2d->setPosition(crossPoint.x(), crossPoint.y(), crossPoint.z());
432  double timeResolution = 1.0;
433  if (m_TimeResolution.isValid()) {
434  timeResolution *= m_TimeResolution->getTimeResolution(
435  (*it8)->getUniqueChannelID());
436  timeResolution *= m_TimeResolution->getTimeResolution(
437  (*it9)->getUniqueChannelID());
438  }
439  hit2d->setChiSq((t1 - t2) * (t1 - t2) / timeResolution);
440  hit2d->setTime(time);
441  hit2d->setMCTime(((*it8)->getMCTime() + (*it9)->getMCTime()) / 2);
442  hit2d->addRelationTo(*it8);
443  hit2d->addRelationTo(*it9);
444  for (i = 0; i < 2; i++) {
445  EKLMAlignmentHit* alignmentHit = m_eklmAlignmentHits.appendNew(i);
446  alignmentHit->addRelationTo(hit2d);
447  }
448  /* Exit the loop. Equivalent to selection of the earliest hit. */
449  break;
450  }
451  /* Exit the loop. Equivalent to selection of the earliest hit. */
452  break;
453  }
454  it6 = it7;
455  }
456  it4 = it5;
457  }
458  it1 = it2;
459  }
460 }
461 
463 {
464 }
465 
467 {
468  delete m_eklmTransformData;
469 }
static KLMChannelNumber channelNumber(int section, int sector, int layer, int plane, int strip)
Get channel number.
@ c_FirstRPCLayer
First RPC layer.
static bool hitsFromSameModule(int module1, int module2)
Check whether the hits are from the same module.
Store one reconstructed BKLM 1D hit as a ROOT object.
Definition: BKLMHit1d.h:30
bool inRPC() const
Determine whether this 1D hit is in RPC or scintillator.
Definition: BKLMHit1d.h:54
float getTime() const
Get reconstructed hit time.
Definition: BKLMHit1d.h:126
double getStripAve() const
Get average strip number.
Definition: BKLMHit1d.h:111
Store one BKLM strip hit as a ROOT object.
Definition: BKLMHit2d.h:32
bool isOutOfTime()
Determine whether this 2D hit is outside the trigger-coincidence window.
Definition: BKLMHit2d.h:137
static const double speedOfLight
[cm/ns]
Definition: Const.h:575
This dataobject is used only for EKLM alignment.
EKLM element numbers.
static const EKLMElementNumbers & Instance()
Instantiation.
int sectorNumber(int section, int layer, int sector) const
Get sector number.
static constexpr int getMaximalStripGlobalNumber()
Get maximal strip global number.
int getNPlanes() const
Get number of planes.
Class for 2d hits handling.
Definition: EKLMHit2d.h:30
void setChiSq(float chisq)
Set Chi^2 of the crossing point.
Definition: EKLMHit2d.h:63
void setEnergyDeposit(float eDep)
Set EnergyDeposit.
Definition: EKLMHitBase.h:109
void setTime(float time)
Set hit time.
Definition: EKLMHitBase.h:127
void setPosition(float x, float y, float z)
Set hit global position.
void setMCTime(float t)
Set MC time.
Definition: EKLMHitMCTime.h:45
static const GeometryData & Instance(enum DataSource dataSource=c_Database, const GearDir *gearDir=nullptr)
Instantiation.
Definition: GeometryData.cc:33
Transformation data.
Definition: TransformData.h:35
@ c_Alignment
Use alignment data (for everything else).
Definition: TransformData.h:45
bool intersection(KLMDigit *hit1, KLMDigit *hit2, HepGeom::Point3D< double > *cross, double *d1, double *d2, double *sd, bool segments=true) const
Check if strips intersect, and find intersection point if yes.
ChannelStatus
Channel status.
@ c_Normal
Normally operating channel.
@ c_Dead
Dead channel (no signal).
@ c_Unknown
Unknown status (no data).
KLM digit (class representing a digitized hit in RPCs or scintillators).
Definition: KLMDigit.h:30
bool inRPC() const
Determine whether the hit is in RPC or scintillator.
Definition: KLMDigit.h:213
int getSubdetector() const
Get subdetector number.
Definition: KLMDigit.h:79
int getLayer() const
Get layer number.
Definition: KLMDigit.h:133
bool isGood() const
Whether hit could be used late (if it passed discriminator threshold)
Definition: KLMDigit.h:355
float getTime() const
Get hit time.
Definition: KLMDigit.h:283
int getSection() const
Get section number.
Definition: KLMDigit.h:97
int getPlane() const
Get plane number.
Definition: KLMDigit.h:151
int getStrip() const
Get strip number.
Definition: KLMDigit.h:169
bool isMultiStrip() const
Determine whether this digit is a multi-strip one or not.
Definition: KLMDigit.h:204
int getSector() const
Get sector number.
Definition: KLMDigit.h:115
KLM element numbers.
KLMChannelNumber channelNumber(int subdetector, int section, int sector, int layer, int plane, int strip) const
Get channel number.
OptionalDBObjPtr< KLMTimeResolution > m_TimeResolution
KLM time resolution.
void reconstructEKLMHits()
Reconstruct EKLMHit2d.
bool m_TimeCableDelayCorrection
Perform cable delay time correction (true) or not (false).
StoreArray< KLMDigit > m_Digits
KLM digits.
DBObjPtr< KLMChannelStatus > m_ChannelStatus
Channel status.
EKLM::TransformData * m_eklmTransformData
Transformation data.
void initialize() override
Initializer.
void event() override
Called for each event.
const KLMElementNumbers * m_ElementNumbers
KLM element numbers.
bklm::GeometryPar * m_bklmGeoPar
BKLM GeometryPar singleton.
StoreArray< EKLMAlignmentHit > m_eklmAlignmentHits
Alignment Hits.
StoreArray< BKLMHit2d > m_bklmHit2ds
BKLM 2d hits.
double m_CoincidenceWindow
Half-width of the time coincidence window used to create a 2D hit from 1D digits/hits.
void endRun() override
Called if the current run ends.
const EKLMElementNumbers * m_eklmElementNumbers
EKLM element numbers.
double m_DelayBKLMScintillators
Delay (ns / cm) for BKLM scintillators.
void terminate() override
Called at the end of the event processing.
bool m_bklmIfAlign
Perform alignment correction (true) or not (false).
double m_EventT0Value
Value of the EventT0.
void reconstructBKLMHits()
Reconstruct BKLMHit1d and BKLMHit2d.
double m_PromptTime
Nominal time of prompt BKLMHit2ds.
double m_DelayRPCPhi
Delay (ns / cm) for RPC phi plane.
void beginRun() override
Called when entering a new run.
bool isNormal(const KLMDigit *digit) const
Check if channel is normal or dead.
OptionalDBObjPtr< KLMTimeCableDelay > m_TimeCableDelay
KLM time cable delay.
void correctCableDelay(double &td, const KLMDigit *digit)
Time correction by subtract cable delay.
double m_PromptWindow
Half-width of the time window relative to the prompt time for BKLMHit2ds.
StoreArray< BKLMHit1d > m_bklmHit1ds
BKLM 1d hits.
const EKLM::GeometryData * m_eklmGeoDat
Geometry data.
double m_DelayEKLMScintillators
Delay (ns / cm) for EKLM scintillators.
StoreObjPtr< EventT0 > m_EventT0
EventT0.
StoreArray< EKLMHit2d > m_eklmHit2ds
EKLM 2d hits.
bool m_EventT0Correction
Perform EventT0 correction (true) or not (false).
bool m_bklmIgnoreScintillators
Ignore scintillators (to debug their electronics mapping).
bool m_IgnoreHotChannels
Use only normal and dead (for debugging) channels during 2d hit reconstruction.
OptionalDBObjPtr< KLMTimeConstants > m_TimeConstants
KLM time constants.
double m_DelayRPCZ
Delay (ns / cm) for RPC Z plane.
DBObjPtr< KLMTimeWindow > m_TimeWindow
KLM time window.
bool m_eklmCheckSegmentIntersection
Check if segments intersect.
@ c_BKLM
BKLM scintillator.
@ c_EKLM
EKLM scintillator.
Base class for Modules.
Definition: Module.h:72
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
void addRelationTo(const RelationsInterface< BASE > *object, float weight=1.0, const std::string &namedRelation="") const
Add a relation from this object to another object (with caching).
const Module * findModule(int section, int sector, int layer) const
Get the pointer to the definition of a module.
Definition: GeometryPar.cc:727
static GeometryPar * instance(void)
Static method to get a reference to the singleton GeometryPar instance.
Definition: GeometryPar.cc:28
Define the geometry of a BKLM module Each sector [octant] contains Modules.
Definition: Module.h:76
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
uint16_t KLMChannelNumber
Channel number.
Abstract base class for different kinds of events.