Belle II Software  release-08-02-06
ECLLocalMaximumFinderModule.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 <ecl/modules/eclLocalMaximumFinder/ECLLocalMaximumFinderModule.h>
11 
12 /* ECL headers. */
13 #include <ecl/dataobjects/ECLCalDigit.h>
14 #include <ecl/dataobjects/ECLConnectedRegion.h>
15 #include <ecl/dataobjects/ECLDigit.h>
16 #include <ecl/dataobjects/ECLElementNumbers.h>
17 #include <ecl/dataobjects/ECLHit.h>
18 #include <ecl/dataobjects/ECLLocalMaximum.h>
19 #include <ecl/geometry/ECLGeometryPar.h>
20 #include <ecl/geometry/ECLNeighbours.h>
21 
22 /* ROOT headers. */
23 #include <TFile.h>
24 #include <TTree.h>
25 
26 /* Basf2 headers. */
27 #include <framework/gearbox/Const.h>
28 #include <framework/logging/Logger.h>
29 #include <mdst/dataobjects/MCParticle.h>
30 
31 // NAMESPACE(S)
32 using namespace Belle2;
33 using namespace ECL;
34 
35 //-----------------------------------------------------------------
36 // Register the Module
37 //-----------------------------------------------------------------
38 REG_MODULE(ECLLocalMaximumFinder);
39 REG_MODULE(ECLLocalMaximumFinderPureCsI);
40 
41 //-----------------------------------------------------------------
42 // Implementation
43 //-----------------------------------------------------------------
45  m_mcParticles(mcParticleArrayName()),
46  m_eclHits(eclHitArrayName()),
47  m_eclDigits(eclDigitArrayName()),
48  m_eclCalDigits(eclCalDigitArrayName()),
49  m_eclConnectedRegions(eclConnectedRegionArrayName()),
50  m_eclLocalMaximums(eclLocalMaximumArrayName())
51 {
52  // Set description.
53  setDescription("ECLLocalMaximumFinderModule");
54 
55  // Parallel processing certification.
57 
58  // Add module parameters.
59  addParam("energyCut", m_energyCut, "Seed energy cut [MeV], minimum is 5.0 MeV.", 10.0 * Belle2::Unit::MeV);
60  addParam("isTrainingMode", m_isTrainingMode,
61  "Run in training mode (i.e. fill file with MVA input variables and determine MC truth of LM.).", 0);
62  addParam("outfileName", m_outfileName, "Output file name for training file.", std::string("ECLLocalMaximumFinderOutput.root"));
63  addParam("method", m_method, "Method to determine the LM (cut, none, fastbdt).", std::string("none"));
64  addParam("truthFraction", m_truthFraction, "Minimum matched energy fraction truth/rec for the LM.", 0.51);
65  addParam("cutOffset", m_cutOffset, "Cut method specific: Offset. (BaBar: 2.5, high eff: 1.40)", 2.5);
66  addParam("cutSlope", m_cutSlope, "Cut method specific: Slope. (BaBar: 0.5, high eff: 3.0)", 0.5);
67  addParam("cutRatioCorrection", m_cutRatioCorrection, "Cut method specific: Ratio correction.", 0.0);
68 
69 }
70 
72 {
73  ;
74 }
75 
77 {
78  B2DEBUG(200, "ECLLocalMaximumFinderModule::initialize()");
79 
80  m_eclLocalMaximums.registerInDataStore(eclLocalMaximumArrayName());
81  m_eclCalDigits.registerInDataStore(eclCalDigitArrayName());
83  m_eclConnectedRegions.registerRelationTo(m_eclLocalMaximums);
84  m_eclCalDigits.registerRelationTo(m_eclLocalMaximums);
85  m_eventLevelClusteringInfo.registerInDataStore();
86 
87  // Check user input.
89  B2WARNING("ECLLocalMaximumFinderModule::initialize: Energy threshold too small, resetting to " << c_minEnergyCut << " GeV");
91  }
92 
93  if (m_truthFraction < 1e-6) {
94  B2WARNING("ECLLocalMaximumFinderModule::initialize: Matching fraction must be above 1e-6, input: " << m_truthFraction <<
95  " will be reset to 1e-6.");
96  m_truthFraction = 1e-6;
97  }
98 
99  if (m_truthFraction > 1.0) {
100  B2WARNING("ECLLocalMaximumFinderModule::initialize: Matching fraction must be below 1.0, input: " << m_truthFraction <<
101  " will be reset to 1.0.");
102  m_truthFraction = 1.0;
103  }
104 
105  // Geometry instance.
107 
108  // Initialize neighbour map.
109  m_neighbourMap = new ECLNeighbours("N", 1);
110 
111  // Reset all variables.
113 
114  // Open output files and declare branches if in training mode. Each file will hold a flat ntuple of training data.
115  m_outfile = nullptr;
116  m_tree = nullptr;
117  if (m_isTrainingMode > 0) {
118  const int cBufferLength = 500;
119  char tmpBuffer[cBufferLength];
120  int tmpBufferLength = sprintf(tmpBuffer, "%s", m_outfileName.c_str());
121 
122  if (tmpBufferLength >= cBufferLength) {
123  B2FATAL("ECLLocalMaximumFinderModule::initialize(): Output training filename length too long!");
124  }
125 
126  m_outfile = new TFile(tmpBuffer, "RECREATE");
127  m_tree = new TTree("locmax", "locmax");
128  m_tree->Branch("energyRatioNeighbour0", &m_energyRatioNeighbour[0], "energyRatioNeighbour0/F");
129  m_tree->Branch("energyRatioNeighbour1", &m_energyRatioNeighbour[1], "energyRatioNeighbour1/F");
130  m_tree->Branch("energyRatioNeighbour2", &m_energyRatioNeighbour[2], "energyRatioNeighbour2/F");
131  m_tree->Branch("energyRatioNeighbour3", &m_energyRatioNeighbour[3], "energyRatioNeighbour3/F");
132  m_tree->Branch("energyRatioNeighbour4", &m_energyRatioNeighbour[4], "energyRatioNeighbour4/F");
133  m_tree->Branch("energyRatioNeighbour5", &m_energyRatioNeighbour[5], "energyRatioNeighbour5/F");
134  m_tree->Branch("energyRatioNeighbour6", &m_energyRatioNeighbour[6], "energyRatioNeighbour6/F");
135  m_tree->Branch("energyRatioNeighbour7", &m_energyRatioNeighbour[7], "energyRatioNeighbour7/F");
136  m_tree->Branch("energyRatioNeighbour8", &m_energyRatioNeighbour[8], "energyRatioNeighbour8/F");
137  m_tree->Branch("energyRatioNeighbour9", &m_energyRatioNeighbour[9], "energyRatioNeighbour9/F");
138  m_tree->Branch("energyRatioNeighbour10", &m_energyRatioNeighbour[10], "energyRatioNeighbour10/F");
139  m_tree->Branch("energyRatioNeighbour11", &m_energyRatioNeighbour[11], "energyRatioNeighbour11/F");
140  m_tree->Branch("energy", &m_energy, "energy/F");
141  m_tree->Branch("target", &m_target, "target/F");
142  m_tree->Branch("targetindex", &m_targetindex, "targetindex/F");
143  m_tree->Branch("targetpi0index", &m_targetpi0index, "targetpi0index/F");
144  m_tree->Branch("thetaId", &m_thetaId, "thetaId/F");
145  m_tree->Branch("phiId", &m_phiId, "phiId/F");
146  m_tree->Branch("cellId", &m_cellId, "cellId/F");
147  m_tree->Branch("maxNeighbourEnergy", &m_maxNeighbourEnergy, "maxNeighbourEnergy/F");
148  m_tree->Branch("nNeighbours10", &m_nNeighbours10, "nNeighbours10/F");
149  m_tree->Branch("time", &m_time, "time/F");
150  m_tree->Branch("timeResolution", &m_timeResolution, "timeResolution/F");
151  m_tree->Branch("timeFitFailed", &m_timeFitFailed, "timeFitFailed/F");
152  m_tree->Branch("CRId", &m_CRId, "CRId/F");
153  m_tree->Branch("LMId", &m_LMId, "LMId/F");
154  }
155 
156  // initialize the vector that gives the relation between cellid and store array position
158 
159 }
160 
162 {
163  ;
164 }
165 
167 {
168  B2DEBUG(200, "ECLLocalMaximumFinderModule::event()");
169 
170  // Fill a vector that can be used to map cellid -> store array position
171  std::fill_n(m_StoreArrPosition.begin(), m_StoreArrPosition.size(), -1);
172  for (int i = 0; i < m_eclCalDigits.getEntries(); i++) {
173  m_StoreArrPosition[m_eclCalDigits[i]->getCellId()] = i;
174  }
175 
176  // Vector with neighbour ids.
177  std::vector< double > vNeighourEnergies;
178  vNeighourEnergies.resize(c_nMaxNeighbours);
179 
180  // Loop over connected regions.
181  for (const ECLConnectedRegion& aCR : m_eclConnectedRegions) {
182  const int crId = aCR.getCRId();
183  int iLM = 1;
184 
185  // Loop over all entries in this CR.
186  for (const ECLCalDigit& aECLCalDigit : aCR.getRelationsTo<ECLCalDigit>()) {
187 
188  // Check seed energy cut.
189  if (aECLCalDigit.getEnergy() >= m_energyCut) {
190 
191  // Clean up for this candiate (MVA is trained per LM, regardless of CR)
192  std::fill_n(vNeighourEnergies.begin(), vNeighourEnergies.size(),
193  -999); // -999 means later: this digit is just not available in this neighbour definition.
196 
197  // Check neighbours: Must be a local energy maximum.
198  bool isLocMax = 1;
199  int neighbourCount = 0;
200  for (auto& neighbourId : m_neighbourMap->getNeighbours(aECLCalDigit.getCellId())) {
201  if (neighbourId == aECLCalDigit.getCellId()) continue; // Skip the center cell to avoid possible floating point issues.
202 
203  const int pos = m_StoreArrPosition[neighbourId]; // Get position in the store array for this digit.
204 
205  double energyNeighbour = 0.0;
206  if (pos >= 0) {
207  energyNeighbour = m_eclCalDigits[pos]->getEnergy(); // Get the energy directly from the store array.
208  vNeighourEnergies[neighbourCount] = energyNeighbour;
209  } else {
210  // Digit does not belong to this CR
211  vNeighourEnergies[neighbourCount] = 0.0;
212  }
213  ++neighbourCount;
214 
215  if (energyNeighbour > aECLCalDigit.getEnergy()) {
216  isLocMax = 0;
217  break;
218  }
219  }
220 
221  // It is a local maximum. Get all variables needed for classification.
222  if (isLocMax) {
223 
224  for (unsigned int npos = 0; npos < vNeighourEnergies.size(); ++npos) {
225  if (vNeighourEnergies[npos] >= 0) {
226  m_energyRatioNeighbour[npos] = static_cast < float >(vNeighourEnergies[npos] / aECLCalDigit.getEnergy());
227  if (vNeighourEnergies[npos] > m_maxNeighbourEnergy) m_maxNeighbourEnergy = vNeighourEnergies[npos];
228  if (vNeighourEnergies[npos] > 1.0 * Belle2::Unit::MeV) ++m_nNeighbours10;
229  } else m_energyRatioNeighbour[npos] = 0.0;
230  }
231  m_time = aECLCalDigit.getTime();
232  m_maxEnergyRatio = m_maxNeighbourEnergy / aECLCalDigit.getEnergy();
233 
234  // Fill training monitoring variables and MC truth information.
235  if (m_isTrainingMode > 0) {
236 
237  m_energy = static_cast < float >(aECLCalDigit.getEnergy());
238  m_cellId = static_cast < float >(aECLCalDigit.getCellId());
239  m_timeResolution = static_cast < float >(aECLCalDigit.getTimeResolution());
240  m_timeFitFailed = static_cast < float >(aECLCalDigit.isFailedFit());
241  m_CRId = static_cast < float >(crId);
242  m_LMId = static_cast < float >(iLM);
243 
244  m_geom->Mapping(m_cellId - 1);
245  m_thetaId = static_cast < float >(m_geom->GetThetaID());
246  m_phiId = static_cast < float >(m_geom->GetPhiID());
247 
248  // This requires MC matching before this stage!
249  int pi0index = -1;
250  int maxtype = 0;
251  int maxpos = 0;
252 
253  auto relatedParticlePairs = aECLCalDigit.getRelationsWith<MCParticle>();
254 
255  for (unsigned int irel = 0; irel < relatedParticlePairs.size(); irel++) {
256  const auto particle = relatedParticlePairs.object(irel);
257  const double weight = relatedParticlePairs.weight(irel);
258 
259  int motherpdg = -1;
260  int motherindex = -1;
261  pi0index = -1;
262  getEnteringMother(*particle, motherpdg, motherindex, pi0index);
263  addToSignalEnergy(motherpdg, motherindex, pi0index, weight);
264  }
265 
266  B2DEBUG(175, " -> digt energy: " << aECLCalDigit.getEnergy());
267  B2DEBUG(175, "photon: " << m_signalEnergy[0][0] << " " << m_signalEnergy[0][1]);
268  B2DEBUG(175, "pi0: " << m_signalEnergy[1][0] << " " << m_signalEnergy[1][1]);
269  B2DEBUG(175, "electron: " << m_signalEnergy[2][0] << " " << m_signalEnergy[2][1]);
270  B2DEBUG(175, "muon: " << m_signalEnergy[3][0] << " " << m_signalEnergy[3][1]);
271  B2DEBUG(175, "neutral hadron: " << m_signalEnergy[4][0] << " " << m_signalEnergy[4][1]);
272  B2DEBUG(175, "charged hadron: " << m_signalEnergy[5][0] << " " << m_signalEnergy[5][1]);
273  B2DEBUG(175, "other: " << m_signalEnergy[6][0]);
274 
275  maxtype = 0;
276  maxpos = 0;
277  getMax(maxtype, maxpos);
278 
279  if (maxtype >= 0) {
280  if (m_signalEnergy[maxtype][maxpos] >= m_truthFraction * aECLCalDigit.getEnergy()) {
281  m_target = maxtype;
282  m_targetindex = m_signalId[maxtype][maxpos];
283  m_targetpi0index = pi0index;
284  } else {
285  m_target = 7;
286  }
287  } else {
288  m_target = 7;
289  }
290 
291  m_tree->Fill();
292 
293  } // end training
294 
295  if (m_method == "cut") {
296 
297  B2DEBUG(200, "m_cutSlope: " << m_cutSlope << ", m_nNeighbours10: " << m_nNeighbours10 << ", m_cutOffset: " << m_cutOffset <<
298  ", m_maxNeighbourEnergy: " << m_maxNeighbourEnergy << ", m_cutRatioCorrection: " << m_cutRatioCorrection <<
299  ", aECLCalDigit.getEnergy(): " << aECLCalDigit.getEnergy());
300  B2DEBUG(200, "m_cutSlope * (m_nNeighbours10 - m_cutOffset): " << m_cutSlope * (m_nNeighbours10 - m_cutOffset));
301  B2DEBUG(200, "(m_maxNeighbourEnergy - m_cutRatioCorrection) / (aECLCalDigit.getEnergy() - m_cutRatioCorrection)" <<
302  (m_maxNeighbourEnergy - m_cutRatioCorrection) / (aECLCalDigit.getEnergy() - m_cutRatioCorrection) << "\n");
303 
305  (aECLCalDigit.getEnergy() - m_cutRatioCorrection)) {
306  makeLocalMaximum(aCR, aECLCalDigit.getCellId(), iLM);
307  ++iLM;
308  }
309  } else if (m_method == "none") { // All energy local maximums will become local maximums.
310  makeLocalMaximum(aCR, aECLCalDigit.getCellId(), iLM);
311  ++iLM;
312  }
313 
314  }
315 
316  } // end check energy
317 
318  } // end CalDigit loop
319 
320  // Check if there is at least one local maximum in the CR. If not, make the highest energetic crystal one.
321  if (iLM == 1) {
322 
323  int highestEnergyCellId = -1;
324  double highestEnergy = 0.0;
325 
326  // Loop over all entries in this CR.
327  for (const ECLCalDigit& aECLCalDigit : aCR.getRelationsTo<ECLCalDigit>(eclCalDigitArrayName())) {
328  if (aECLCalDigit.getEnergy() > highestEnergy) {
329  highestEnergyCellId = aECLCalDigit.getCellId();
330  highestEnergy = aECLCalDigit.getEnergy();
331  }
332  } // end CalDigit loop
333 
334  makeLocalMaximum(aCR, highestEnergyCellId, 0);
335 
336  }
337 
338  } // end CR loop
339 
340  // Find the number of local maximums in each ECL region in mdst
341  uint16_t nLMPerRegion[3] = {};
342  for (const ECLLocalMaximum& aLM : m_eclLocalMaximums) {
343  const int iCellId = aLM.getCellId();
344  if (ECLElementNumbers::isForward(iCellId)) {nLMPerRegion[0]++;}
345  if (ECLElementNumbers::isBarrel(iCellId)) {nLMPerRegion[1]++;}
346  if (ECLElementNumbers::isBackward(iCellId)) {nLMPerRegion[2]++;}
347  }
348 
349  // Store numbers in EventLevelClusteringInfo mdst object
351  m_eventLevelClusteringInfo->setNECLLocalMaximumsFWD(nLMPerRegion[0]);
352  m_eventLevelClusteringInfo->setNECLLocalMaximumsBarrel(nLMPerRegion[1]);
353  m_eventLevelClusteringInfo->setNECLLocalMaximumsBWD(nLMPerRegion[2]);
354 
355 }
356 
357 
359 {
360  B2DEBUG(200, "ECLLocalMaximumFinderModule::endRun()");
361 }
362 
363 
365 {
366  B2DEBUG(200, "ECLLocalMaximumFinderModule::terminate()");
367 
368  // If run in trainingmode, write the training file.
369  if (m_outfile and m_tree) {
370  m_outfile->cd();
371  m_tree->Write();
372  m_outfile->Write();
373  m_outfile->Close();
374  delete m_outfile;
375  }
376 
377  if (m_neighbourMap) delete m_neighbourMap;
378 
379 }
380 
381 void ECLLocalMaximumFinderModule::makeLocalMaximum(const ECLConnectedRegion& aCR, const int cellId, const int lmId)
382 {
383  // Set the local maximum dataobject.
384  const auto aLocalMaximum = m_eclLocalMaximums.appendNew();
385 
386  B2DEBUG(175, "ECLLocalMaximumFinderModule::makeLocalMaximum(): local maximum cellid: " << cellId);
387 
388  // Set the id of this local maximum.
389  aLocalMaximum->setLMId(lmId);
390 
391  // Set the cell Id of the digit.
392  aLocalMaximum->setCellId(cellId);
393 
394  // Add relations to ECLConnectedRegion.
395  aCR.addRelationTo(aLocalMaximum);
396 
397 }
398 
400 {
401  m_target = -1;
402  m_targetindex = -1;
403  m_targetpi0index = -1;
404 
405  m_totalSignalEnergy = 0.;
406 
407  for (unsigned int i = 0; i < 10; ++i) {
408  for (unsigned int j = 0; j < 5; ++j) {
409  m_signalEnergy[i][j] = 0.;
410  m_signalId[i][j] = -1;
411  }
412  }
413 
414 }
415 
417 {
418 
419  for (unsigned i = 0; i < c_nMaxNeighbours; ++i) {
420  m_energyRatioNeighbour[i] = 0.;
421  }
422 
423  m_energy = 0.;
424  m_cellId = -1;
425  m_thetaId = -1;
426  m_phiId = -1;
427  m_nNeighbours10 = 0;
428  m_maxNeighbourEnergy = 0.0;
429  m_time = -9999;
430  m_timeResolution = -9999.;
431  m_timeFitFailed = -9999.;
432  m_CRId = -1;
433  m_LMId = -1;
434  m_maxEnergyRatio = 0.0;
435 }
436 
437 
438 void ECLLocalMaximumFinderModule::getEnteringMother(const MCParticle& particle, int& pdg, int& arrayindex, int& pi0arrayindex)
439 {
440 
441  int index = particle.getArrayIndex();
442  int pi0index = -1;
443 
444  while (!isEnteringECL(m_mcParticles[index]->getProductionVertex())) {
445  if (m_mcParticles[index]->getMother()) index = m_mcParticles[index]->getMother()->getArrayIndex();
446  else index = -1;
447  };
448 
449  // For photon mother: are they from a pi0? This can be used to improved overlap/merged pi0 reconstruction.
450  if (m_mcParticles[index]->getPDG() == Const::photon.getPDGCode()) {
451  if (m_mcParticles[index]->getMother()->getPDG() == Const::pi0.getPDGCode()) {
452  pi0index = m_mcParticles[index]->getMother()->getArrayIndex();
453  }
454  }
455 
456  // Dont include mother if its energy is too low or if its a photon from a neutron interaction.
457  if ((m_mcParticles[index]->getEnergy() < 5.0 * Belle2::Unit::MeV)
458  or (m_mcParticles[index]->getPDG() == Const::photon.getPDGCode()
459  and abs(m_mcParticles[index]->getMother()->getPDG()) == Const::proton.getPDGCode())) {
460  pdg = -1;
461  arrayindex = -1;
462  pi0arrayindex = -1;
463  } else {
464  pdg = m_mcParticles[index]->getPDG();
465  arrayindex = index;
466  pi0arrayindex = pi0index;
467  }
468 }
469 
471 {
472 
473  const double theta = vertex.Theta();
474 
475  if (theta > 0.555015 and theta < 2.26369) { //barrel
476  double radius = vertex.Perp();
477  if (radius < 125 * Belle2::Unit::cm) return true;
478  } else if (theta <= 0.555015) { //fwd
479  if (vertex.Z() < 196.16 * Belle2::Unit::cm) return true;
480  } else if (theta >= 2.26369) { //bwd
481  if (vertex.Z() > -102.16 * Belle2::Unit::cm) return true;
482  }
483 
484  return false;
485 }
486 
487 int ECLLocalMaximumFinderModule::getIdPosition(const int type, const int id)
488 {
489 
490  for (int i = 0; i < 5; i++) {
491  if (m_signalId[type][i] == id) return i;
492  if (m_signalId[type][i] == -1) {
493  m_signalId[type][i] = id;
494  return i; //next free one
495  }
496  }
497 
498  return -1;
499 }
500 
501 // of all entries, get the maximum!
502 void ECLLocalMaximumFinderModule::getMax(int& type, int& id)
503 {
504  double maxe = 0.;
505  int maxtype = -1;
506  int maxid = -1;
507 
508  for (unsigned int i = 0; i < 10; ++i) {
509  for (unsigned int j = 0; j < 5; ++j) {
510  if (m_signalEnergy[i][j] > maxe) {
511  maxe = m_signalEnergy[i][j];
512  maxtype = i;
513  maxid = j;
514  }
515  }
516  }
517 
518  type = maxtype;
519  id = maxid;
520 
521 }
522 
523 void ECLLocalMaximumFinderModule::addToSignalEnergy(int motherpdg, int motherindex, int pi0index, double weight)
524 {
525 
526  // for the LM training and CR/LM debugging
527  if (motherpdg == Const::photon.getPDGCode()) {
528  if (pi0index >= 0) { // photon from pi0
529  int idpos = getIdPosition(1, motherindex);
530  m_signalEnergy[1][idpos] += weight;
531  } else {
532  int idpos = getIdPosition(0, motherindex);
533  m_signalEnergy[0][idpos] += weight; // photon from another source
534  }
535  } else if (abs(motherpdg) == Const::electron.getPDGCode()) { // electron
536  int idpos = getIdPosition(2, motherindex);
537  m_signalEnergy[2][idpos] += weight;
538  } else if (abs(motherpdg) == Const::muon.getPDGCode()) { // muon
539  int idpos = getIdPosition(3, motherindex);
540  m_signalEnergy[3][idpos] += weight;
541  } else if (abs(motherpdg) == Const::Klong.getPDGCode() or abs(motherpdg) == Const::neutron.getPDGCode()) { // neutral hadron
542  int idpos = getIdPosition(4, motherindex);
543  m_signalEnergy[4][idpos] += weight;
544  } else if (abs(motherpdg) == Const::pion.getPDGCode() or abs(motherpdg) == Const::kaon.getPDGCode()
545  or abs(motherpdg) == Const::proton.getPDGCode()) { // charged hadron
546  int idpos = getIdPosition(5, motherindex);
547  m_signalEnergy[5][idpos] += weight;
548  } else { // everything else
549  int idpos = getIdPosition(6, motherindex);
550  m_signalEnergy[6][idpos] += weight;
551  }
552 }
static const ParticleType neutron
neutron particle
Definition: Const.h:666
static const ParticleType pi0
neutral pion particle
Definition: Const.h:665
static const ChargedStable muon
muon particle
Definition: Const.h:651
static const ChargedStable pion
charged pion particle
Definition: Const.h:652
static const ParticleType Klong
K^0_L particle.
Definition: Const.h:669
static const ChargedStable proton
proton particle
Definition: Const.h:654
static const ChargedStable kaon
charged kaon particle
Definition: Const.h:653
static const ParticleType photon
photon particle
Definition: Const.h:664
static const ChargedStable electron
electron particle
Definition: Const.h:650
Class to store calibrated ECLDigits: ECLCalDigits.
Definition: ECLCalDigit.h:23
Class to store connected regions (CRs)
TTree * m_tree
tree that contain information for MVA training
double m_cutRatioCorrection
correction for nominator and denominator of the ratio.
double m_signalEnergy[10][5]
total energy per MC matching type of this digit
void getEnteringMother(const MCParticle &particle, int &pdg, int &arrayindex, int &pi0arrayindex)
Get enterging mother of this particle.
StoreArray< ECLConnectedRegion > m_eclConnectedRegions
Store array: ECLConnectedRegion.
float m_energyRatioNeighbour[c_nMaxNeighbours]
energy ratio of neighbour 0..9 to center
double m_totalSignalEnergy
total energy of this digit
float m_phiId
local maximum center theta Id
virtual void initialize() override
Initialize.
StoreArray< ECLLocalMaximum > m_eclLocalMaximums
Store array: ECLLocalMaximum.
float m_cellId
local maximum center cell Id
int getIdPosition(const int type, const int id)
Get Id position in the vector.
float m_energy
Variables to monitor the MVA training.
std::string m_method
Method to find the local maximum.
virtual void terminate() override
Terminate (close ROOT files here if you have opened any).
std::string m_outfileName
file name prefix of the training output file
ECL::ECLNeighbours * m_neighbourMap
Neighbour maps.
void makeLocalMaximum(const ECLConnectedRegion &aCR, const int cellId, const int lmId)
Make local maximum for a given connected region.
void resetClassifierVariables()
Reset Classifier Variables.
float m_thetaId
local maximum center theta Id
float m_maxNeighbourEnergy
highest energy of all neighbours
bool isEnteringECL(const B2Vector3D &vec)
Check if particle is produced outside of the ECL.
TFile * m_outfile
Output training files and trees.
void addToSignalEnergy(int motherpdg, int motherindex, int pi0index, double weight)
Add energy to vector.
static const unsigned c_nMaxNeighbours
Variables (possibly) used for MVA classification.
StoreArray< MCParticle > m_mcParticles
Store array: MCParticle.
float m_maxEnergyRatio
Highest energetic neighbour energy divided by LM energy.
int m_signalId[10][5]
total energy per MC matching type of this digit
const double c_minEnergyCut
Minimum LM energy.
virtual const char * eclCalDigitArrayName() const
Name to be used for default or PureCsI option: ECLCalDigits.
std::vector< int > m_StoreArrPosition
vector (ECLElementNumbers::c_NCrystals + 1 entries) with cell id to store array positions
StoreObjPtr< EventLevelClusteringInfo > m_eventLevelClusteringInfo
EventLevelClusteringInfo.
void resetTrainingVariables()
Reset Debug Variables.
void getMax(int &maxtype, int &maxpos)
Get the highest energy deposition particle type.
StoreArray< ECLCalDigit > m_eclCalDigits
Store array: ECLCalDigit.
int m_isTrainingMode
training mode for MVA methods (i.e.
float m_nNeighbours10
Variables (possibly) used for cut classification.
virtual const char * eclConnectedRegionArrayName() const
Name to be used for default option: ECLConnectedRegions.
virtual const char * eclLocalMaximumArrayName() const
Name to be used for default option: ECLLocalMaximums.
Class to store local maxima (LM)
static ECLGeometryPar * Instance()
Static method to get a reference to the ECLGeometryPar instance.
void Mapping(int cid)
Mapping theta, phi Id.
int GetThetaID()
Get Theta Id.
Class to get the neighbours for a given cell id.
Definition: ECLNeighbours.h:25
const std::vector< short int > & getNeighbours(short int cid) const
Return the neighbours for a given cell ID.
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
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).
RelationVector< TO > getRelationsTo(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from this object to another store array.
static const double MeV
[megaelectronvolt]
Definition: Unit.h:114
static const double cm
Standard units with the value = 1.
Definition: Unit.h:47
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
const int c_NCrystals
Number of crystals.
bool isBackward(int cellId)
Check whether the crystal is in backward ECL.
bool isForward(int cellId)
Check whether the crystal is in forward ECL.
bool isBarrel(int cellId)
Check whether the crystal is in barrel ECL.
Abstract base class for different kinds of events.