10 #include <ecl/modules/eclSplitterN1/ECLSplitterN1Module.h>
23 #include <framework/logging/Logger.h>
24 #include <framework/utilities/FileSystem.h>
25 #include <framework/geometry/B2Vector3.h>
28 #include <ecl/utility/Position.h>
29 #include <ecl/dataobjects/ECLCalDigit.h>
30 #include <ecl/dataobjects/ECLConnectedRegion.h>
31 #include <ecl/dataobjects/ECLLocalMaximum.h>
32 #include <ecl/dataobjects/ECLShower.h>
33 #include <ecl/geometry/ECLNeighbours.h>
34 #include <ecl/geometry/ECLGeometryPar.h>
37 #include <mdst/dataobjects/EventLevelClusteringInfo.h>
54 m_eclCalDigits(eclCalDigitArrayName()),
55 m_eclConnectedRegions(eclConnectedRegionArrayName()),
56 m_eclShowers(eclShowerArrayName()),
57 m_eclLocalMaximums(eclLocalMaximumArrayName()),
58 m_eventLevelClusteringInfo(eventLevelClusteringInfoName())
61 setDescription(
"ECLSplitterN1Module: Baseline reconstruction splitter code for the n photon hypothesis.");
62 addParam(
"fullBkgdCount", m_fullBkgdCount,
63 "Number of background digits at full background (as provided by EventLevelClusteringInfo).",
69 addParam(
"threshold", m_threshold,
"Threshold energy after splitting.", 7.5 *
Belle2::Unit::MeV);
70 addParam(
"expConstant", m_expConstant,
"Constant a from exp(-a*dist/RM), typical: 1.5 to 3.5?", 2.5);
71 addParam(
"maxIterations", m_maxIterations,
"Maximum number of iterations for centroid shifts.", 100);
72 addParam(
"shiftTolerance", m_shiftTolerance,
"Tolerance level for centroid shifts.", 1.0 *
Belle2::Unit::mm);
73 addParam(
"minimumSharedEnergy", m_minimumSharedEnergy,
"Minimum shared energy.", 25.0 *
Belle2::Unit::keV);
74 addParam(
"maxSplits", m_maxSplits,
"Maximum number of splits within one connected region.", 10);
75 addParam(
"cutDigitEnergyForEnergy", m_cutDigitEnergyForEnergy,
76 "Minimum digit energy to be included in the shower energy calculation. (NOT USED)", 0.5 *
Belle2::Unit::MeV);
77 addParam(
"cutDigitTimeResidualForEnergy", m_cutDigitTimeResidualForEnergy,
78 "Maximum time residual to be included in the shower energy calculation. (NOT USED)", 5.0);
81 addParam(
"useOptimalNumberOfDigitsForEnergy", m_useOptimalNumberOfDigitsForEnergy,
82 "Optimize the number of digits for energy calculations.", 1);
83 addParam(
"fileBackgroundNormName", m_fileBackgroundNormName,
"Background filename.",
85 addParam(
"fileNOptimalFWDName", m_fileNOptimalFWDName,
"FWD number of optimal neighbours filename.",
87 addParam(
"fileNOptimalBarrelName", m_fileNOptimalBarrelName,
"Barrel number of optimal neighbours filename.",
89 addParam(
"fileNOptimalBWDName", m_fileNOptimalBWDName,
"BWD number of optimal neighbours filename.",
93 addParam(
"positionMethod", m_positionMethod,
"Position determination method.", std::string(
"lilo"));
94 addParam(
"liloParameterA", m_liloParameterA,
"Position determination linear-log. parameter A.", 4.0);
95 addParam(
"liloParameterB", m_liloParameterB,
"Position determination linear-log. parameter B.", 0.0);
96 addParam(
"liloParameterC", m_liloParameterC,
"Position determination linear-log. parameter C.", 0.0);
99 setPropertyFlags(c_ParallelProcessingCertified);
114 m_liloParameters.resize(3);
115 m_liloParameters.at(0) = m_liloParameterA;
116 m_liloParameters.at(1) = m_liloParameterB;
117 m_liloParameters.at(2) = m_liloParameterC;
120 m_eclCalDigits.registerInDataStore(eclCalDigitArrayName());
121 m_eclConnectedRegions.registerInDataStore(eclConnectedRegionArrayName());
122 m_eclShowers.registerInDataStore(eclShowerArrayName());
123 m_eclLocalMaximums.registerInDataStore(eclLocalMaximumArrayName());
124 m_eventLevelClusteringInfo.registerInDataStore(eventLevelClusteringInfoName());
127 m_eclShowers.registerRelationTo(m_eclConnectedRegions);
128 m_eclShowers.registerRelationTo(m_eclCalDigits);
129 m_eclShowers.registerRelationTo(m_eclLocalMaximums);
130 m_eclLocalMaximums.registerRelationTo(m_eclCalDigits);
137 m_StoreArrPosition.resize(8736 + 1);
138 m_StoreArrPositionLM.resize(8736 + 1);
141 m_fileBackgroundNorm =
new TFile(m_fileBackgroundNormName.c_str(),
"READ");
142 if (!m_fileBackgroundNorm) B2FATAL(
"Could not find file: " << m_fileBackgroundNormName);
143 m_th1fBackgroundNorm =
dynamic_cast<TH1F*
>(m_fileBackgroundNorm->Get(
"background_norm"));
144 if (!m_th1fBackgroundNorm) B2FATAL(
"Could not find m_th1fBackgroundNorm");
147 m_fileNOptimalFWD =
new TFile(m_fileNOptimalFWDName.c_str(),
"READ");
148 if (!m_fileNOptimalFWD) B2FATAL(
"Could not find file: " << m_fileNOptimalFWDName);
149 for (
unsigned t = 0; t < 13; ++t) {
150 for (
unsigned s = 0; s < c_nSectorCellIdFWD[t]; ++s) {
151 m_tg2dNOptimalFWD[t][s] =
dynamic_cast<TGraph2D*
>(m_fileNOptimalFWD->Get(Form(
"thetaid-%i_sectorcellid-%i", t, s)));
152 if (!m_tg2dNOptimalFWD[t][s]) B2FATAL(
"Could not find TGraph2D m_tg2dNOptimalFWD!");
156 m_fileNOptimalBarrel =
new TFile(m_fileNOptimalBarrelName.c_str(),
"READ");
157 if (!m_fileNOptimalBarrel) B2FATAL(
"Could not find file: " << m_fileNOptimalBarrelName);
159 m_tg2dNOptimalBarrel =
dynamic_cast<TGraph2D*
>(m_fileNOptimalBarrel->Get(
"thetaid-50_sectorcellid-8"));
160 if (!m_tg2dNOptimalBarrel) B2FATAL(
"Could not find TGraph2D m_tg2dNOptimalBarrel!");
162 m_fileNOptimalBWD =
new TFile(m_fileNOptimalBWDName.c_str(),
"READ");
163 if (!m_fileNOptimalBWD) B2FATAL(
"Could not find file: " << m_fileNOptimalBWDName);
164 for (
unsigned t = 0; t < 10; ++t) {
165 for (
unsigned s = 0; s < c_nSectorCellIdBWD[t]; ++s) {
166 m_tg2dNOptimalBWD[t][s] =
dynamic_cast<TGraph2D*
>(m_fileNOptimalBWD->Get(Form(
"thetaid-%i_sectorcellid-%i", t + 59, s)));
167 if (!m_tg2dNOptimalBWD[t][s]) B2FATAL(
"Could not find TGraph2D m_tg2dNOptimalBWD!");
180 B2DEBUG(175,
"ECLCRSplitterModule::event()");
183 memset(&m_StoreArrPosition[0], -1, m_StoreArrPosition.size() *
sizeof m_StoreArrPosition[0]);
184 for (
int i = 0; i < m_eclCalDigits.getEntries(); i++) {
185 m_StoreArrPosition[m_eclCalDigits[i]->getCellId()] = i;
189 memset(&m_StoreArrPositionLM[0], -1, m_StoreArrPositionLM.size() *
sizeof m_StoreArrPositionLM[0]);
190 for (
int i = 0; i < m_eclLocalMaximums.getEntries(); i++) {
191 m_StoreArrPositionLM[m_eclLocalMaximums[i]->getCellId()] = i;
195 for (
auto& aCR : m_eclConnectedRegions) {
197 m_cellIdInCR.clear();
199 const unsigned int entries = (aCR.getRelationsWith<
ECLCalDigit>(eclCalDigitArrayName())).size();
201 m_cellIdInCR.resize(entries);
205 for (
const auto& caldigit : aCR.getRelationsWith<
ECLCalDigit>(eclCalDigitArrayName())) {
206 m_cellIdInCR[i] = caldigit.getCellId();
211 splitConnectedRegion(aCR);
227 if (m_fileBackgroundNorm)
delete m_fileBackgroundNorm;
228 if (m_fileNOptimalFWD)
delete m_fileNOptimalFWD;
229 if (m_fileNOptimalBarrel)
delete m_fileNOptimalBarrel;
230 if (m_fileNOptimalBWD)
delete m_fileNOptimalBWD;
232 if (m_NeighbourMap9)
delete m_NeighbourMap9;
233 if (m_NeighbourMap21)
delete m_NeighbourMap21;
240 const int bkgdcount = m_eventLevelClusteringInfo->getNECLCalDigitsOutOfTime();
241 double backgroundLevel = 0.0;
242 if (m_fullBkgdCount > 0) {
243 backgroundLevel =
static_cast<double>(bkgdcount) /
static_cast<double>(m_fullBkgdCount);
249 B2DEBUG(175,
"ECLCRSplitterModule::splitConnectedRegion: nLocalMaximums = " << nLocalMaximums);
258 if (nLocalMaximums == 1 or nLocalMaximums >= m_maxSplits) {
261 const auto aECLShower = m_eclShowers.appendNew();
264 aECLShower->addRelationTo(&aCR);
267 double weightSum = 0.0;
273 const int locmaxcellid = locmaxvector[0]->getCellId();
274 const int pos = m_StoreArrPosition[locmaxcellid];
275 double highestEnergyID = (m_eclCalDigits[pos])->getCellId();
276 double highestEnergy = (m_eclCalDigits[pos])->getEnergy();
277 double highestEnergyTime = (m_eclCalDigits[pos])->getTime();
278 double highestEnergyTimeResolution = (m_eclCalDigits[pos])->getTimeResolution();
281 const double energyEstimation = estimateEnergy(highestEnergyID);
285 int nNeighbours = getNeighbourMap(energyEstimation, backgroundLevel);
286 if (nNeighbours == 9 and !m_useOptimalNumberOfDigitsForEnergy) neighbourMap = m_NeighbourMap9;
287 else neighbourMap = m_NeighbourMap21;
290 std::vector<ECLCalDigit> digits;
291 std::vector<double> weights;
292 for (
auto& neighbourId : neighbourMap->
getNeighbours(highestEnergyID)) {
293 const auto it = std::find(m_cellIdInCR.begin(), m_cellIdInCR.end(),
295 if (it == m_cellIdInCR.end())
continue;
297 const int neighbourpos = m_StoreArrPosition[neighbourId];
298 digits.push_back(*m_eclCalDigits[neighbourpos]);
299 weights.push_back(1.0);
302 aECLShower->addRelationTo(m_eclCalDigits[neighbourpos], 1.0);
306 const B2Vector3D& showerposition = Belle2::ECL::computePositionLiLo(digits, weights, m_liloParameters);
307 aECLShower->setTheta(showerposition.
Theta());
308 aECLShower->setPhi(showerposition.
Phi());
309 aECLShower->setR(showerposition.
Mag());
312 double showerEnergy = 0.0;
313 if (m_useOptimalNumberOfDigitsForEnergy) {
316 const unsigned int nOptimal = getOptimalNumberOfDigits(highestEnergyID, energyEstimation, backgroundLevel);
317 aECLShower->setNominalNumberOfCrystalsForEnergy(
static_cast<double>(nOptimal));
320 std::vector< std::pair<unsigned int, double>> listCrystalPairs;
321 listCrystalPairs.resize(digits.size());
323 std::vector < std::pair<double, double> > weighteddigits;
324 weighteddigits.resize(digits.size());
325 for (
unsigned int i = 0; i < digits.size(); ++i) {
326 weighteddigits.at(i) = std::make_pair((digits.at(i)).getEnergy(), weights.at(i));
327 listCrystalPairs.at(i) = std::make_pair((digits.at(i)).getCellId(), weights.at(i) * (digits.at(i)).getEnergy());
331 std::sort(listCrystalPairs.begin(), listCrystalPairs.end(), [](
const auto & left,
const auto & right) {
332 return left.second < right.second;
334 std::vector< unsigned int> listCrystals;
336 for (
unsigned int i = 0; i < digits.size(); ++i) {
338 listCrystals.push_back(listCrystalPairs[i].first);
342 aECLShower->setNumberOfCrystalsForEnergy(
static_cast<double>(listCrystals.size()));
343 aECLShower->setListOfCrystalsForEnergy(listCrystals);
345 showerEnergy = getEnergySum(weighteddigits, nOptimal);
346 B2DEBUG(150,
"Shower Energy (1): " << showerEnergy);
349 showerEnergy = Belle2::ECL::computeEnergySum(digits, weights);
352 aECLShower->setEnergy(showerEnergy);
353 aECLShower->setEnergyRaw(showerEnergy);
354 aECLShower->setEnergyHighestCrystal(highestEnergy);
355 aECLShower->setTime(highestEnergyTime);
356 aECLShower->setDeltaTime99(highestEnergyTimeResolution);
357 aECLShower->setNumberOfCrystals(weightSum);
358 aECLShower->setCentralCellId(highestEnergyID);
360 B2DEBUG(175,
"theta = " << showerposition.
Theta());
361 B2DEBUG(175,
"phi = " << showerposition.
Phi());
362 B2DEBUG(175,
"R = " << showerposition.
Mag());
363 B2DEBUG(175,
"energy = " << showerEnergy);
364 B2DEBUG(175,
"time = " << highestEnergyTime);
365 B2DEBUG(175,
"time resolution = " << highestEnergyTimeResolution);
366 B2DEBUG(175,
"neighbours = " << nNeighbours);
367 B2DEBUG(175,
"backgroundLevel = " << backgroundLevel);
370 aECLShower->setShowerId(1);
372 aECLShower->setConnectedRegionId(aCR.
getCRId());
375 const int posLM = m_StoreArrPositionLM[locmaxcellid];
377 const int posDigit = m_StoreArrPosition[aDigit.getCellId()];
378 m_eclLocalMaximums[posLM]->addRelationTo(m_eclCalDigits[posDigit], 1.0);
384 std::vector<ECLCalDigit> digits;
385 std::vector<double> weights;
386 std::map<int, B2Vector3D> centroidList;
387 std::map<int, double> centroidEnergyList;
388 std::map<int, B2Vector3D> allPoints;
389 std::map<int, std::vector < double > > weightMap;
390 std::vector < ECLCalDigit > digitVector;
393 std::map<int, B2Vector3D> localMaximumsPoints;
394 std::map<int, B2Vector3D> centroidPoints;
397 int cellid = aLocalMaximum.getCellId();
400 B2Vector3D vectorPosition = m_geom->GetCrystalPos(cellid - 1);
401 localMaximumsPoints.insert(std::map<int, B2Vector3D>::value_type(cellid, vectorPosition));
402 centroidPoints.insert(std::map<int, B2Vector3D>::value_type(cellid, vectorPosition));
406 bool iterateclusters =
true;
410 centroidList.clear();
411 centroidEnergyList.clear();
418 const int cellid = aCalDigit.getCellId();
420 B2Vector3D vectorPosition = m_geom->GetCrystalPos(cellid - 1);
421 allPoints.insert(std::map<int, B2Vector3D>::value_type(cellid, vectorPosition));
422 digits.push_back(aCalDigit);
425 for (
const auto& digitpoint : allPoints) {
426 const int cellid = digitpoint.first;
427 const int pos = m_StoreArrPosition[cellid];
428 digitVector.push_back(*m_eclCalDigits[pos]);
437 double centroidShiftAverage = 0.0;
441 B2DEBUG(175,
"Iteration: #" << nIterations <<
" (of max. " << m_maxIterations <<
")");
443 centroidShiftAverage = 0.0;
444 if (nIterations == 0) {
445 centroidList.clear();
446 centroidEnergyList.clear();
451 for (
const auto& locmaxpoint : localMaximumsPoints) {
454 const int locmaxcellid = locmaxpoint.first;
460 if (nIterations == 0) {
461 const int pos = m_StoreArrPosition[locmaxcellid];
462 const double locmaxenergy = m_eclCalDigits[pos]->getEnergy();
463 centroidEnergyList[locmaxcellid] = 1.5 * locmaxenergy;
466 B2DEBUG(175,
"local maximum cellid: " << locmaxcellid);
471 for (
const auto& digitpoint : allPoints) {
474 const int digitcellid = digitpoint.first;
477 const int pos = m_StoreArrPosition[digitcellid];
478 const double digitenergy = m_eclCalDigits[pos]->getEnergy();
482 double distance = 0.0;
483 double distanceEnergySum = 0.0;
486 for (
const auto& centroidpoint : centroidPoints) {
489 const int centroidcellid = centroidpoint.first;
490 B2Vector3D centroidpos = centroidpoint.second;
492 double thisdistance = 0.;
495 if (nIterations == 0 and digitcellid == centroidcellid) {
498 B2Vector3D vectorDistance = ((centroidpos) - (digitpos));
499 thisdistance = vectorDistance.
Mag();
503 const int thispos = m_StoreArrPosition[centroidcellid];
504 const double thisenergy = m_eclCalDigits[thispos]->getEnergy();
507 if (locmaxcellid == centroidcellid) {
508 distance = thisdistance;
513 const double expfactor = exp(-m_expConstant * thisdistance / c_molierRadius);
514 distanceEnergySum += (thisenergy * expfactor);
519 if (distanceEnergySum > 0.0) {
520 const double expfactor = exp(-m_expConstant * distance / c_molierRadius);
521 weight = energy * expfactor / distanceEnergySum;
527 if ((digitenergy * weight) < m_minimumSharedEnergy) {
533 B2WARNING(
"ECLCRSplitterModule::splitConnectedRegion: Floating point glitch, weight for this digit " << weight <<
534 ", resetting it to 1.0.");
539 B2DEBUG(175,
" cellid: " << digitcellid <<
", energy: " << digitenergy <<
", weight: " << weight <<
", distance: " << distance);
540 weights.push_back(weight);
546 B2Vector3D oldCentroidPos = (centroidPoints.find(locmaxcellid))->second;
549 B2Vector3D newCentroidPos = Belle2::ECL::computePositionLiLo(digits, weights, m_liloParameters);
552 const double newEnergy = Belle2::ECL::computeEnergySum(digits, weights);
555 const B2Vector3D centroidShift = (oldCentroidPos - newCentroidPos);
558 centroidList[locmaxcellid] = newCentroidPos;
559 weightMap[locmaxcellid] = weights;
561 B2DEBUG(175,
"--> inserting new energy: " << newEnergy <<
" for local maximum " << locmaxcellid);
562 centroidEnergyList[locmaxcellid] = newEnergy;
563 double showerenergy = (*centroidEnergyList.find(locmaxcellid)).second /
Belle2::Unit::MeV;
564 B2DEBUG(175,
"--> new energy = " << showerenergy <<
" MeV");
567 centroidShiftAverage += centroidShift.
Mag();
570 B2DEBUG(175,
" old centroid: " << oldCentroidPos.
x() <<
" cm, " << oldCentroidPos.
y() <<
" cm, " << oldCentroidPos.
z() <<
572 B2DEBUG(175,
" new centroid: " << newCentroidPos.
x() <<
" cm, " << newCentroidPos.
y() <<
" cm, " << newCentroidPos.
z() <<
574 B2DEBUG(175,
" centroid shift: " << centroidShift.
Mag() <<
" cm");
579 centroidShiftAverage /=
static_cast<double>(nLocalMaximums);
581 B2DEBUG(175,
"--> average centroid shift: " << centroidShiftAverage <<
" cm (tolerance is " << m_shiftTolerance <<
" cm)");
584 for (
const auto& locmaxpoint : localMaximumsPoints) {
585 centroidPoints[locmaxpoint.first] = (centroidList.find(locmaxpoint.first))->second;
590 }
while (nIterations < m_maxIterations and centroidShiftAverage > m_shiftTolerance);
594 std::vector<int> markfordeletion;
595 iterateclusters =
false;
596 for (
const auto& locmaxpoint : localMaximumsPoints) {
599 const int locmaxcellid = locmaxpoint.first;
600 const int pos = m_StoreArrPosition[locmaxcellid];
602 B2DEBUG(175,
"locmaxcellid: " << locmaxcellid);
603 const double lmenergy = m_eclCalDigits[pos]->getEnergy();
604 B2DEBUG(175,
"ok: ");
607 std::vector < double > myWeights = (*weightMap.find(locmaxcellid)).second;
609 for (
unsigned int i = 0; i < digitVector.size(); ++i) {
612 const double weight = myWeights[i];
619 if ((energy * weight > lmenergy and cellid != locmaxcellid) or (energy * weight < m_threshold and cellid == locmaxcellid)) {
620 markfordeletion.push_back(locmaxcellid);
621 iterateclusters =
true;
628 for (
const auto lmid : markfordeletion) {
629 localMaximumsPoints.erase(lmid);
630 centroidPoints.erase(lmid);
633 }
while (iterateclusters);
636 unsigned int iShower = 1;
637 for (
const auto& locmaxpoint : localMaximumsPoints) {
639 const int locmaxcellid = locmaxpoint.first;
640 const int posLM = m_StoreArrPositionLM[locmaxcellid];
643 const auto aECLShower = m_eclShowers.appendNew();
646 const double energyEstimation = estimateEnergy(locmaxcellid);
650 int nNeighbours = getNeighbourMap(energyEstimation, backgroundLevel);
651 if (nNeighbours == 9 and !m_useOptimalNumberOfDigitsForEnergy) neighbourMap = m_NeighbourMap9;
652 else neighbourMap = m_NeighbourMap21;
655 std::vector<short int> neighbourlist = neighbourMap->
getNeighbours(locmaxcellid);
658 std::vector < double > myWeights = (*weightMap.find(locmaxcellid)).second;
661 std::vector<ECLCalDigit> newdigits;
662 std::vector<double> newweights;
663 double highestEnergy = 0.;
664 double highestEnergyTime = 0.;
665 double highestEnergyTimeResolution = 0.;
666 double weightSum = 0.0;
668 for (
unsigned int i = 0; i < digitVector.size(); ++i) {
671 const double weight = myWeights[i];
674 const int pos = m_StoreArrPosition[cellid];
677 m_eclLocalMaximums[posLM]->addRelationTo(m_eclCalDigits[pos], weight);
681 if (std::find(neighbourlist.begin(), neighbourlist.end(), cellid) != neighbourlist.end()) {
683 aECLShower->addRelationTo(m_eclCalDigits[pos], weight);
685 newdigits.push_back(dig);
686 newweights.push_back(weight);
692 if (energy * weight > highestEnergy) {
693 highestEnergy = energy * weight;
694 highestEnergyTime = dig.
getTime();
704 B2DEBUG(175,
"old theta: " << oldshowerposition->
Theta());
705 B2DEBUG(175,
"old phi: " << oldshowerposition->
Phi());
706 B2DEBUG(175,
"old R: " << oldshowerposition->
Mag());
707 B2DEBUG(175,
"old energy: " << energyEstimation);
708 delete oldshowerposition;
713 B2Vector3D* showerposition =
new B2Vector3D(Belle2::ECL::computePositionLiLo(newdigits, newweights, m_liloParameters));
714 aECLShower->setTheta(showerposition->
Theta());
715 aECLShower->setPhi(showerposition->
Phi());
716 aECLShower->setR(showerposition->
Mag());
718 B2DEBUG(175,
"new theta: " << showerposition->
Theta());
719 B2DEBUG(175,
"new phi: " << showerposition->
Phi());
720 B2DEBUG(175,
"new R: " << showerposition->
Mag());
721 delete showerposition;
724 double showerEnergy = 0.0;
725 if (m_useOptimalNumberOfDigitsForEnergy) {
728 const unsigned int nOptimal = getOptimalNumberOfDigits(locmaxcellid, energyEstimation, backgroundLevel);
729 aECLShower->setNominalNumberOfCrystalsForEnergy(
static_cast<double>(nOptimal));
732 std::vector< std::pair<unsigned int, double>> listCrystalPairs;
733 listCrystalPairs.resize(newdigits.size());
735 std::vector < std::pair<double, double> > weighteddigits;
736 weighteddigits.resize(newdigits.size());
737 for (
unsigned int i = 0; i < newdigits.size(); ++i) {
738 weighteddigits.at(i) = std::make_pair((newdigits.at(i)).getEnergy(), newweights.at(i));
739 listCrystalPairs.at(i) = std::make_pair((newdigits.at(i)).getCellId(), newweights.at(i) * (newdigits.at(i)).getEnergy());
743 std::sort(listCrystalPairs.begin(), listCrystalPairs.end(), [](
const auto & left,
const auto & right) {
744 return left.second < right.second;
747 std::vector< unsigned int> listCrystals;
749 for (
unsigned int i = 0; i < newdigits.size(); ++i) {
751 listCrystals.push_back(listCrystalPairs[i].first);
755 aECLShower->setNumberOfCrystalsForEnergy(
static_cast<double>(listCrystals.size()));
756 aECLShower->setListOfCrystalsForEnergy(listCrystals);
758 showerEnergy = getEnergySum(weighteddigits, nOptimal);
759 B2DEBUG(150,
"Shower Energy (2): " << showerEnergy);
762 showerEnergy = Belle2::ECL::computeEnergySum(newdigits, newweights);
765 aECLShower->setEnergy(showerEnergy);
766 aECLShower->setEnergyRaw(showerEnergy);
767 aECLShower->setEnergyHighestCrystal(highestEnergy);
768 aECLShower->setTime(highestEnergyTime);
769 aECLShower->setDeltaTime99(highestEnergyTimeResolution);
770 aECLShower->setNumberOfCrystals(weightSum);
771 aECLShower->setCentralCellId(locmaxcellid);
773 B2DEBUG(175,
"new energy: " << showerEnergy);
776 aECLShower->setShowerId(iShower);
779 aECLShower->setConnectedRegionId(aCR.
getCRId());
782 aECLShower->addRelationTo(&aCR);
785 aECLShower->addRelationTo(m_eclLocalMaximums[posLM]);
792 if (background <= 0.1)
return 21;
794 if (energy > 0.06 + 0.4 * background)
return 21;
801 int nOptimalNeighbours = 21;
804 const double bgCorrected = bg * m_th1fBackgroundNorm->GetBinContent(cellid);
807 if (bgCorrected > 0.025) {
809 double energyChecked = energy;
810 double bgChecked = bgCorrected;
811 if (bgCorrected > 1.0) bgChecked = 1.0;
816 m_geom->Mapping(cellid - 1);
817 const int thetaId = m_geom->GetThetaID();
818 const int crystalsPerSector = c_crystalsPerRing[thetaId] / 16;
819 const int phiIdInSector = m_geom->GetPhiID() % crystalsPerSector;
822 nOptimalNeighbours =
static_cast<unsigned int>(m_tg2dNOptimalFWD[thetaId][phiIdInSector]->Interpolate(energyChecked, bgChecked));
823 }
else if (thetaId < 59) {
824 nOptimalNeighbours =
static_cast<unsigned int>(m_tg2dNOptimalBarrel->Interpolate(energyChecked, bgChecked));
826 nOptimalNeighbours =
static_cast<unsigned int>(m_tg2dNOptimalBWD[thetaId - 59][phiIdInSector]->Interpolate(energyChecked,
829 B2DEBUG(175,
"ECLSplitterN1Module::getOptimalNumberOfDigits: theta ID: " << thetaId <<
", bg: " << bg <<
" (after corr.: " <<
830 bgCorrected <<
"), energy: " << energy <<
", n: " << nOptimalNeighbours);
831 }
else B2DEBUG(175,
"ECLSplitterN1Module::getOptimalNumberOfDigits: small bg: " << bg <<
" (after corr.: " << bgCorrected <<
832 "), energy: " << energy <<
", n: " << nOptimalNeighbours);
834 return nOptimalNeighbours;
841 double energysum = 0.;
843 std::sort(weighteddigits.begin(), weighteddigits.end(), std::greater<std::pair<double, double>>());
845 unsigned int min = n;
846 if (weighteddigits.size() < n) min = weighteddigits.size();
848 for (
unsigned int i = 0; i < min; ++i) {
849 B2DEBUG(150,
"getEnergySum: " << weighteddigits.at(i).first <<
" " << weighteddigits.at(i).second);
850 energysum += (weighteddigits.at(i).first * weighteddigits.at(i).second);
852 B2DEBUG(150,
"getEnergySum: energysum=" << energysum);
861 double energyEstimation = 0.0;
863 for (
auto& neighbourId : m_NeighbourMap9->getNeighbours(centerid)) {
866 const auto it = std::find(m_cellIdInCR.begin(), m_cellIdInCR.end(),
868 if (it == m_cellIdInCR.end())
continue;
870 const int pos = m_StoreArrPosition[neighbourId];
871 const double energyNeighbour = m_eclCalDigits[pos]->getEnergy();
873 energyEstimation += energyNeighbour;
876 return energyEstimation;
DataType Phi() const
The azimuth angle.
DataType y() const
access variable Y (= .at(1) without boundary check)
DataType Theta() const
The polar angle.
DataType z() const
access variable Z (= .at(2) without boundary check)
DataType Mag() const
The magnitude (rho in spherical coordinate system).
DataType x() const
access variable X (= .at(0) without boundary check)
Class to store calibrated ECLDigits: ECLCalDigits.
int getCellId() const
Get Cell ID.
double getEnergy() const
Get Calibrated Energy.
double getTimeResolution() const
Get Calibrated Time Resolution.
double getTime() const
Get Calibrated Time.
Class to store connected regions (CRs)
int getCRId() const
Get CR identifieer.
Class to store local maxima (LM)
@ c_nPhotons
CR is split into n photons (N1)
Class to perform the shower correction.
~ECLSplitterN1Module()
Destructor.
void splitConnectedRegion(ECLConnectedRegion &aCR)
Split connected region into showers.
virtual void initialize() override
Initialize.
virtual void event() override
Event.
double estimateEnergy(const int centerid)
Estimate energy using 3x3 around central crystal.
virtual void endRun() override
End run.
virtual void terminate() override
Terminate.
virtual void beginRun() override
Begin run.
int getNeighbourMap(const double energy, const double background)
Get number of neighbours based on first energy estimation and background level per event.
unsigned int getOptimalNumberOfDigits(const int cellid, const double energy, const double background)
Get optimal number of digits (out of 21) based on first energy estimation and background level per ev...
double getEnergySum(std::vector< std::pair< double, double > > &weighteddigits, const unsigned int n)
Get energy sum for weighted entries.
static ECLGeometryPar * Instance()
Static method to get a reference to the ECLGeometryPar instance.
Class to get the neighbours for a given cell id.
const std::vector< short int > & getNeighbours(short int cid) const
Return the neighbours for a given cell ID.
static std::string findFile(const std::string &path, bool silent=false)
Search for given file or directory in local or central release directory, and return absolute path if...
Class for type safe access to objects that are referred to in relations.
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< T > getRelationsWith(const std::string &name="", const std::string &namedRelation="") const
Get the relations between this object and another store array.
static const double mm
[millimeters]
static const double keV
[kiloelectronvolt]
static const double MeV
[megaelectronvolt]
static const double GeV
Standard of [energy, momentum, mass].
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
B2Vector3< double > B2Vector3D
typedef for common usage with double
Abstract base class for different kinds of events.