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