Belle II Software  release-06-02-00
KLMStripEfficiencyCollectorModule.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/KLMStripEfficiencyCollector/KLMStripEfficiencyCollectorModule.h>
11 #include <klm/dataobjects/KLMChannelIndex.h>
12 
13 /* CLHEP headers. */
14 #include <CLHEP/Vector/ThreeVector.h>
15 
16 using namespace Belle2;
17 
18 REG_MODULE(KLMStripEfficiencyCollector)
19 
22  m_ElementNumbers(&(KLMElementNumbers::Instance())),
23  m_eklmElementNumbers(&(EKLMElementNumbers::Instance())),
24  m_GeometryBKLM(nullptr),
25  m_PlaneArrayIndex(&(KLMPlaneArrayIndex::Instance())),
26  m_MatchingFile(nullptr),
27  m_MatchingTree(nullptr),
28  m_MatchingHitData( {0, 0, 0, 0, 0, 0, 0., nullptr, nullptr}),
29  m_MatchedStrip(0)
30 {
31  setDescription("Module for KLM strip efficiency data collection.");
32  addParam("MuonListName", m_MuonListName, "Muon list name.",
33  std::string("mu+:all"));
34  addParam("AllowedDistance1D", m_AllowedDistance1D,
35  "Maximal distance in the units of strip number from ExtHit to "
36  "matching KLMDigit.", double(8));
37  addParam("MinimalMatchingDigits", m_MinimalMatchingDigits,
38  "Minimal number of matching digits.", 0);
39  addParam("MinimalMatchingDigitsOuterLayers",
40  m_MinimalMatchingDigitsOuterLayers,
41  "Minimal number of matching digits in outer layers.", 0);
42  addParam("MinimalMomentumNoOuterLayers", m_MinimalMomentumNoOuterLayers,
43  "Minimal momentum in case there are no hits in outer layers.", 0.0);
44  addParam("RemoveUnusedMuons", m_RemoveUnusedMuons,
45  "Whether to remove unused muons.", false);
46  addParam("IgnoreBackwardPropagation", m_IgnoreBackwardPropagation,
47  "Whether to ignore ExtHits with backward propagation.", false);
48  addParam("Debug", m_Debug, "Debug mode.", false);
49  addParam("DebugFileName", m_MatchingFileName, "Debug file name.", std::string("matching.root"));
50  setPropertyFlags(c_ParallelProcessingCertified);
51 }
52 
54 {
55 }
56 
58 {
59  m_Digits.isRequired();
60  m_tracks.isRequired();
61  m_extHits.isRequired();
62  m_MuonList.isRequired(m_MuonListName);
63  int nPlanes = m_PlaneArrayIndex->getNElements();
64  TH1F* matchedDigitsInPlane = new TH1F(
65  "matchedDigitsInPlane", "Number of matching KLMDigits",
66  nPlanes, -0.5, double(nPlanes) - 0.5);
67  TH1F* allExtHitsInPlane = new TH1F(
68  "allExtHitsInPlane", "Number of ExtHits",
69  nPlanes, -0.5, double(nPlanes) - 0.5);
70  registerObject<TH1F>("matchedDigitsInPlane", matchedDigitsInPlane);
71  registerObject<TH1F>("allExtHitsInPlane", allExtHitsInPlane);
73  if (m_Debug) {
74  m_MatchingFile = new TFile(m_MatchingFileName.c_str(), "recreate");
75  m_MatchingTree = new TTree("t_matching", "");
76  m_MatchingTree->Branch("subdetector", &m_MatchingHitData.subdetector,
77  "subdetector/I");
78  m_MatchingTree->Branch("section", &m_MatchingHitData.section, "section/I");
79  m_MatchingTree->Branch("sector", &m_MatchingHitData.sector, "sector/I");
80  m_MatchingTree->Branch("layer", &m_MatchingHitData.layer, "layer/I");
81  m_MatchingTree->Branch("plane", &m_MatchingHitData.plane, "plane/I");
82  m_MatchingTree->Branch("strip", &m_MatchingHitData.strip, "strip/I");
83  m_MatchingTree->Branch("matchedStrip", &m_MatchedStrip, "matchedStrip/I");
84  }
85 }
86 
88 {
89  if (m_Debug) {
90  m_MatchingFile->cd();
91  m_MatchingTree->Write();
92  delete m_MatchingTree;
93  delete m_MatchingFile;
94  }
95 }
96 
98 {
99  if (!m_ChannelStatus.isValid())
100  B2FATAL("KLM channel status data are not available.");
101  int minimalActivePlanes = -1;
103  for (KLMChannelIndex& klmSector : klmSectors) {
104  KLMChannelIndex klmNextSector(klmSector);
105  ++klmNextSector;
106  int activePlanes = 0;
107  KLMChannelIndex klmPlane(klmSector);
109  for (; klmPlane != klmNextSector; ++klmPlane) {
110  KLMChannelIndex klmNextPlane(klmPlane);
111  ++klmNextPlane;
112  bool isActive = false;
113  KLMChannelIndex klmChannel(klmPlane);
115  for (; klmChannel != klmNextPlane; ++ klmChannel) {
116  KLMChannelNumber channel = klmChannel.getKLMChannelNumber();
117  enum KLMChannelStatus::ChannelStatus status =
118  m_ChannelStatus->getChannelStatus(channel);
119  if (status == KLMChannelStatus::c_Unknown)
120  B2FATAL("Incomplete KLM channel status data.");
121  if (status == KLMChannelStatus::c_Normal) {
122  isActive = true;
123  break;
124  }
125  }
126  if (isActive)
127  ++activePlanes;
128  }
129  /*
130  * If a sector is completely off, it is not necessary to reduce
131  * the minimal number of digits, because the efficiencies cannot be
132  * measured for the entire sector anyway.
133  */
134  if (activePlanes == 0)
135  continue;
136  if (minimalActivePlanes < 0)
137  minimalActivePlanes = activePlanes;
138  else if (minimalActivePlanes < activePlanes)
139  minimalActivePlanes = activePlanes;
140  }
141  if ((minimalActivePlanes >= 0) &&
142  (minimalActivePlanes < m_MinimalMatchingDigits)) {
143  B2WARNING("The minimal number of active planes (" << minimalActivePlanes <<
144  ") is less than the minimal number of matching digits (" <<
145  m_MinimalMatchingDigits << "). The minimal number of "
146  "matching digits is reduced to make the calibration possible.");
147  m_MinimalMatchingDigits = minimalActivePlanes;
148  }
149 }
150 
152 {
153 }
154 
156 {
157  std::vector<unsigned int> toRemove;
158  unsigned int nMuons = m_MuonList->getListSize();
159  /*
160  * The getter functions create the object for particular experiment and run.
161  * It is not guaranteed that collectDataTrack() is called if there are
162  * no tracks (e.g. for too short or bad runs). Thus, they are called here
163  * to guarantee the creation of histograms.
164  */
165  TH1F* matchedDigitsInPlane = getObjectPtr<TH1F>("matchedDigitsInPlane");
166  TH1F* allExtHitsInPlane = getObjectPtr<TH1F>("allExtHitsInPlane");
167  for (unsigned int i = 0; i < nMuons; ++i) {
168  const Particle* muon = m_MuonList->getParticle(i);
169  bool trackUsed = collectDataTrack(muon, matchedDigitsInPlane,
170  allExtHitsInPlane);
171  if (m_RemoveUnusedMuons && !trackUsed)
172  toRemove.push_back(muon->getArrayIndex());
173  }
175  m_MuonList->removeParticles(toRemove);
176 }
177 
179  std::map<KLMPlaneNumber, struct HitData>& hitMap,
180  KLMPlaneNumber planeGlobal, struct HitData* hitData)
181 {
182  std::map<KLMPlaneNumber, struct HitData>::iterator it;
183  it = hitMap.find(planeGlobal);
184  /*
185  * There may be more than one such hit e.g. if track crosses the edge
186  * of the strips or WLS fiber groove. Select only one hit per plane.
187  */
188  if (it == hitMap.end()) {
189  hitMap.insert(std::pair<KLMPlaneNumber, struct HitData>(
190  planeGlobal, *hitData));
191  }
192 }
193 
195  struct HitData* hitData)
196 {
197  for (const KLMDigit& digit : m_Digits) {
198  /*
199  * TODO: multi-strip digits are ugnored for now.
200  * It is necessary to take them into account.
201  */
202  if (digit.isMultiStrip())
203  continue;
204  if (!(digit.getSubdetector() == hitData->subdetector &&
205  digit.getSection() == hitData->section &&
206  digit.getLayer() == hitData->layer &&
207  digit.getSector() == hitData->sector &&
208  digit.getPlane() == hitData->plane))
209  continue;
210  if (fabs(digit.getStrip() - hitData->strip) < m_AllowedDistance1D) {
211  hitData->digit = &digit;
212  return;
213  }
214  }
215 }
216 
218  const Particle* muon, TH1F* matchedDigitsInPlane, TH1F* allExtHitsInPlane)
219 {
220  const int nExtrapolationLayers =
222  const Track* track = muon->getTrack();
223  RelationVector<ExtHit> extHits = track->getRelationsTo<ExtHit>();
224  std::map<KLMPlaneNumber, struct HitData> selectedHits;
225  std::map<KLMPlaneNumber, struct HitData>::iterator it;
226  KLMChannelNumber channel;
228  struct HitData hitData, hitDataPrevious;
229  TVector3 extHitPosition;
230  CLHEP::Hep3Vector extHitPositionCLHEP, localPosition;
231  int layer;
232  int extHitLayer[nExtrapolationLayers] = {0};
233  int digitLayer[nExtrapolationLayers] = {0};
234  hitDataPrevious.subdetector = -1;
235  for (const ExtHit& hit : extHits) {
236  /*
237  * Choose hits that exit the sensitive volume.
238  * It is not possible to use entry hits because of a bug in Geant4:
239  * the step does not always have a correct status (fGeomBoundary),
240  * and, consequently, ExtHits are not created.
241  */
242  if (hit.getStatus() != EXT_EXIT)
243  continue;
244  /*
245  * Ignore ExtHits with backward propagation. This affects cosmic events
246  * only. The removal of hits with backward propagation is normally
247  * not needed, however, it is added because of backward error propagation
248  * bug in Geant4 10.6.
249  */
251  if (hit.isBackwardPropagated())
252  continue;
253  }
254  KLMPlaneNumber planeGlobal = 0;
255  hitData.hit = &hit;
256  hitData.digit = nullptr;
257  if (hit.getDetectorID() == Const::EDetector::EKLM) {
258  int stripGlobal = hit.getCopyID();
261  stripGlobal, &hitData.section, &hitData.layer, &hitData.sector,
262  &hitData.plane, &hitData.strip);
264  hitData.section, hitData.sector, hitData.layer,
265  hitData.plane, hitData.strip);
266  status = m_ChannelStatus->getChannelStatus(channel);
267  if (status == KLMChannelStatus::c_Unknown)
268  B2FATAL("Incomplete KLM channel status data.");
269  if (status == KLMChannelStatus::c_Normal) {
271  hitData.subdetector, hitData.layer);
272  extHitLayer[layer - 1]++;
273  planeGlobal = m_ElementNumbers->planeNumberEKLM(
274  hitData.section, hitData.sector, hitData.layer,
275  hitData.plane);
276  addHit(selectedHits, planeGlobal, &hitData);
277  }
278  } else if (hit.getDetectorID() == Const::EDetector::BKLM) {
279  int moduleNumber = hit.getCopyID();
282  moduleNumber, &hitData.section, &hitData.sector, &hitData.layer);
284  /*
285  * For scintillators, the plane and strip numbers are recorded
286  * in the copy number.
287  */
289  moduleNumber, &hitData.section, &hitData.sector, &hitData.layer,
290  &hitData.plane, &hitData.strip);
292  hitData.section, hitData.sector, hitData.layer,
293  hitData.plane, hitData.strip);
294  status = m_ChannelStatus->getChannelStatus(channel);
295  if (status == KLMChannelStatus::c_Unknown)
296  B2FATAL("Incomplete KLM channel status data.");
297  if (status == KLMChannelStatus::c_Normal) {
299  hitData.subdetector, hitData.layer);
300  extHitLayer[layer - 1]++;
301  planeGlobal = m_ElementNumbers->planeNumberBKLM(
302  hitData.section, hitData.sector, hitData.layer,
303  hitData.plane);
304  addHit(selectedHits, planeGlobal, &hitData);
305  }
306  } else {
307  /* For RPCs, the sensitive volume corresponds to both readout planes. */
308  extHitPosition = hit.getPosition();
309  extHitPositionCLHEP.setX(extHitPosition.X());
310  extHitPositionCLHEP.setY(extHitPosition.Y());
311  extHitPositionCLHEP.setZ(extHitPosition.Z());
312  const bklm::Module* module =
313  m_GeometryBKLM->findModule(hitData.section, hitData.sector,
314  hitData.layer);
315  localPosition = module->globalToLocal(extHitPositionCLHEP);
317  hitData.strip = module->getZStrip(localPosition);
318  /*
319  * FIXME
320  * There are 2 hits per module in RPCs, but the plane information is
321  * not available in ExtHit. For now, 2 entries are created (one for
322  * each plane) for the first hit, and the second one is removed.
323  */
324  if (hitData.subdetector == hitDataPrevious.subdetector &&
325  hitData.section == hitDataPrevious.section &&
326  hitData.sector == hitDataPrevious.sector &&
327  hitData.layer == hitDataPrevious.layer)
328  continue;
329  std::memcpy(&hitDataPrevious, &hitData, sizeof(struct HitData));
330  /* The returned strip may be out of the valid range. */
332  hitData.section, hitData.sector, hitData.layer, hitData.plane,
333  hitData.strip, false)) {
335  hitData.section, hitData.sector, hitData.layer,
336  hitData.plane, hitData.strip);
337  status = m_ChannelStatus->getChannelStatus(channel);
338  if (status == KLMChannelStatus::c_Unknown)
339  B2FATAL("Incomplete KLM channel status data.");
340  if (status == KLMChannelStatus::c_Normal) {
342  hitData.subdetector, hitData.layer);
343  extHitLayer[layer - 1]++;
344  hitData.localPosition = localPosition.z();
345  planeGlobal = m_ElementNumbers->planeNumberBKLM(
346  hitData.section, hitData.sector, hitData.layer,
347  hitData.plane);
348  addHit(selectedHits, planeGlobal, &hitData);
349  }
350  }
352  hitData.strip = module->getPhiStrip(localPosition);
353  /* The returned strip may be out of the valid range. */
355  hitData.section, hitData.sector, hitData.layer, hitData.plane,
356  hitData.strip, false)) {
358  hitData.section, hitData.sector, hitData.layer,
359  hitData.plane, hitData.strip);
360  status = m_ChannelStatus->getChannelStatus(channel);
361  if (status == KLMChannelStatus::c_Unknown)
362  B2FATAL("Incomplete KLM channel status data.");
363  if (status == KLMChannelStatus::c_Normal) {
365  hitData.subdetector, hitData.layer);
366  extHitLayer[layer - 1]++;
367  hitData.localPosition = localPosition.y();
368  planeGlobal = m_ElementNumbers->planeNumberBKLM(
369  hitData.section, hitData.sector, hitData.layer,
370  hitData.plane);
371  addHit(selectedHits, planeGlobal, &hitData);
372  }
373  }
374  }
375  } else
376  continue;
377  }
378  /* Find matching digits. */
379  int nDigits = 0;
380  std::map<int, int>::iterator it2;
381  for (it = selectedHits.begin(); it != selectedHits.end(); ++it) {
382  findMatchingDigit(&(it->second));
383  if (it->second.digit != nullptr) {
384  nDigits++;
386  it->second.subdetector, it->second.layer);
387  digitLayer[layer - 1]++;
388  }
389  }
390  if (nDigits < m_MinimalMatchingDigits)
391  return false;
392  /* Write efficiency histograms */
393  for (it = selectedHits.begin(); it != selectedHits.end(); ++it) {
394  int matchingDigits = 0;
395  int matchingDigitsOuterLayers = 0;
396  int extHitsOuterLayers = 0;
398  it->second.subdetector, it->second.layer) - 1;
399  for (int i = 0; i < nExtrapolationLayers; ++i) {
400  if (i != layer)
401  matchingDigits += digitLayer[i];
402  if (i > layer) {
403  matchingDigitsOuterLayers += digitLayer[i];
404  extHitsOuterLayers += extHitLayer[i];
405  }
406  }
407  /* Check the number of matching digits in other layers. */
408  if (matchingDigits < m_MinimalMatchingDigits)
409  continue;
410  /*
411  * Check the number of matching digits in outer layers relatively to
412  * this hit.
413  */
414  if (matchingDigitsOuterLayers < m_MinimalMatchingDigitsOuterLayers) {
420  if (extHitsOuterLayers >= m_MinimalMatchingDigitsOuterLayers)
421  continue;
422  /*
423  * If the number of ExtHits is insufficient, then check the momentum.
424  * The muons with sufficiently large momentum have a very small
425  * probability to get absorbed in the detector.
426  */
427  if (muon->getMomentum().Mag() < m_MinimalMomentumNoOuterLayers)
428  continue;
429  }
430  allExtHitsInPlane->Fill(m_PlaneArrayIndex->getIndex(it->first));
431  if (it->second.digit) {
432  matchedDigitsInPlane->Fill(m_PlaneArrayIndex->getIndex(it->first));
433  if (m_Debug) {
434  std::memcpy(&m_MatchingHitData, &(it->second), sizeof(struct HitData));
435  m_MatchedStrip = it->second.digit->getStrip();
436  m_MatchingTree->Fill();
437  }
438  }
439  }
440  return true;
441 }
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 void moduleNumberToElementNumbers(KLMModuleNumber module, int *section, int *sector, int *layer)
Get element numbers by module number.
Calibration collector module base class.
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.
Store one Ext hit as a ROOT object.
Definition: ExtHit.h:30
KLM channel index.
KLMChannelNumber getKLMChannelNumber() const
Get KLM channel number.
void setIndexLevel(enum IndexLevel indexLevel)
Set index level.
ChannelStatus
Channel status.
@ c_Normal
Normally operating channel.
@ c_Unknown
Unknown status (no data).
KLM digit (class representing a digitized hit in RPCs or scintillators).
Definition: KLMDigit.h:30
uint16_t getNElements() const
Get number of elements.
uint16_t getIndex(uint16_t number) const
Get element index.
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.
Module KLMStripEfficiencyCollectorModule.
const bklm::GeometryPar * m_GeometryBKLM
BKLM geometry.
void addHit(std::map< KLMPlaneNumber, struct HitData > &hitMap, KLMPlaneNumber planeGlobal, struct HitData *hitData)
Add hit to map.
DBObjPtr< KLMChannelStatus > m_ChannelStatus
Channel status.
void startRun() override
This method is called at the beginning of the run.
double m_MinimalMomentumNoOuterLayers
Minimal momentum in case there are no hits in outer layers.
const KLMElementNumbers * m_ElementNumbers
Element numbers.
const KLMPlaneArrayIndex * m_PlaneArrayIndex
Plane array index.
const EKLMElementNumbers * m_eklmElementNumbers
EKLM element numbers.
void findMatchingDigit(struct HitData *hitData)
Find matching digit.
void collect() override
This method is called for each event.
bool collectDataTrack(const Particle *muon, TH1F *matchedDigitsInPlane, TH1F *allExtHitsInPlane)
Collect the data for one muon.
void closeRun() override
This method is called at the end of 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.
int m_MinimalMatchingDigitsOuterLayers
Minimal number of matching digits in outer layers.
int m_MinimalMatchingDigits
Minimal number of matching digits.
Class to store reconstructed particles.
Definition: Particle.h:74
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: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
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
uint16_t KLMChannelNumber
Channel number.
uint16_t KLMPlaneNumber
Plane number.
Abstract base class for different kinds of events.