Belle II Software development
ECLSplitterN1Module.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/eclSplitterN1/ECLSplitterN1Module.h>
11
12/* ECL header. */
13#include <ecl/dataobjects/ECLCalDigit.h>
14#include <ecl/dataobjects/ECLConnectedRegion.h>
15#include <ecl/dataobjects/ECLElementNumbers.h>
16#include <ecl/dataobjects/ECLLocalMaximum.h>
17#include <ecl/dataobjects/ECLShower.h>
18#include <ecl/dbobjects/ECLnOptimal.h>
19#include <ecl/geometry/ECLGeometryPar.h>
20#include <ecl/geometry/ECLNeighbours.h>
21#include <ecl/utility/Position.h>
22
23/* Basf2 headers. */
24#include <framework/geometry/B2Vector3.h>
25#include <framework/logging/Logger.h>
26#include <mdst/dataobjects/EventLevelClusteringInfo.h>
27
28/* C++ headers. */
29#include <algorithm>
30#include <string>
31#include <utility>
32
33// NAMESPACES
34using namespace Belle2;
35using namespace ECL;
36
37//-----------------------------------------------------------------
38// Register the Module(s)
39//-----------------------------------------------------------------
40REG_MODULE(ECLSplitterN1);
41REG_MODULE(ECLSplitterN1PureCsI);
42
43//-----------------------------------------------------------------
44// Implementation
45//-----------------------------------------------------------------
46
48 m_eclCalDigits(eclCalDigitArrayName()),
49 m_eclConnectedRegions(eclConnectedRegionArrayName()),
50 m_eclShowers(eclShowerArrayName()),
51 m_eclLocalMaximums(eclLocalMaximumArrayName()),
52 m_eventLevelClusteringInfo(eventLevelClusteringInfoName())
53{
54 // Set description.
55 setDescription("ECLSplitterN1Module: Baseline reconstruction splitter code for the n photon hypothesis.");
56 addParam("fullBkgdCount", m_fullBkgdCount,
57 "Number of background digits at full background (as provided by EventLevelClusteringInfo).",
58 182);
59
60 // Set module parameters.
61
62 // Splitter.
63 addParam("threshold", m_threshold, "Threshold energy after splitting.", 7.5 * Belle2::Unit::MeV);
64 addParam("expConstant", m_expConstant, "Constant a from exp(-a*dist/RM), typical: 1.5 to 3.5?", 2.5); // to be optimized!
65 addParam("maxIterations", m_maxIterations, "Maximum number of iterations for centroid shifts.", 100);
66 addParam("shiftTolerance", m_shiftTolerance, "Tolerance level for centroid shifts.", 1.0 * Belle2::Unit::mm);
67 addParam("minimumSharedEnergy", m_minimumSharedEnergy, "Minimum shared energy.", 25.0 * Belle2::Unit::keV);
68 addParam("maxSplits", m_maxSplits, "Maximum number of splits within one connected region.", 10);
69 addParam("cutDigitEnergyForEnergy", m_cutDigitEnergyForEnergy,
70 "Minimum digit energy to be included in the shower energy calculation. (NOT USED)", 0.5 * Belle2::Unit::MeV);
71 addParam("cutDigitTimeResidualForEnergy", m_cutDigitTimeResidualForEnergy,
72 "Maximum time residual to be included in the shower energy calculation. (NOT USED)", 5.0);
73
74 // Neighbour definitions
75 addParam("useOptimalNumberOfDigitsForEnergy", m_useOptimalNumberOfDigitsForEnergy,
76 "Optimize the number of digits for energy calculations.", 1);
77
78 // Position.
79 addParam("positionMethod", m_positionMethod, "Position determination method.", std::string("lilo"));
80 addParam("liloParameterA", m_liloParameterA, "Position determination linear-log. parameter A.", 4.0);
81 addParam("liloParameterB", m_liloParameterB, "Position determination linear-log. parameter B.", 0.0);
82 addParam("liloParameterC", m_liloParameterC, "Position determination linear-log. parameter C.", 0.0);
83
84 // Set parallel processing flag.
86}
87
89{
90 // do not delete objects here, do it in terminate()!
91}
92
94{
95
96 // Geometry instance.
98
99 // Check and format user input.
100 m_liloParameters.resize(3);
104
105 // ECL dataobjects.
109 m_eclShowers.registerInDataStore(eclShowerArrayName());
110
111 // mDST dataobjects.
113
114 // Register relations (we probably dont need all, but keep them for now for debugging).
115 m_eclShowers.registerRelationTo(m_eclConnectedRegions);
116 m_eclShowers.registerRelationTo(m_eclCalDigits);
117 m_eclShowers.registerRelationTo(m_eclLocalMaximums);
118 m_eclLocalMaximums.registerRelationTo(m_eclCalDigits);
119 m_eclConnectedRegions.requireRelationTo(m_eclLocalMaximums);
120 m_eclConnectedRegions.requireRelationTo(m_eclCalDigits);
121
122 // Initialize neighbour maps (we will optimize the endcaps later, there is more than just a certain energy containment to be considered)
123 m_NeighbourMap9 = new ECLNeighbours("N", 1); // N: 3x3 = 9
124 m_NeighbourMap21 = new ECLNeighbours("NC", 2); // NC: 5x5 excluding corners = 21
125
126 // initialize the vector that gives the relation between cellid and store array position
129
130}
131
133{
134 //..Read in nOptimal crystal payload from database
135 if (m_eclNOptimal.hasChanged()) {
136
137 //..Vectors of energy boundaries for each region
138 std::vector<float> eBoundariesFwd = m_eclNOptimal->getUpperBoundariesFwd();
139 std::vector<float> eBoundariesBrl = m_eclNOptimal->getUpperBoundariesBrl();
140 std::vector<float> eBoundariesBwd = m_eclNOptimal->getUpperBoundariesBwd();
141
142 //..Adjust the size of the vector of boundaries to match the number in the payload
143 m_nEnergyBins = eBoundariesBrl.size();
144 B2INFO("ECLSplitterN1 beginRun: number of nOptimal payload energies = " << m_nEnergyBins);
145 m_eBoundaries.resize(m_nLeakReg, std::vector<float>(m_nEnergyBins, 0.));
146
147 //..Copy values to m_eBoundaries
148 for (int ie = 0; ie < m_nEnergyBins; ie++) {
149 m_eBoundaries[0][ie] = eBoundariesFwd[ie];
150 m_eBoundaries[1][ie] = eBoundariesBrl[ie];
151 m_eBoundaries[2][ie] = eBoundariesBwd[ie];
152 B2INFO(" upper boundaries for energy point " << ie << " " << m_eBoundaries[0][ie] << " " << m_eBoundaries[1][ie] << " " <<
153 m_eBoundaries[2][ie]);
154 }
155
156 //..Group number of each crystal
157 m_groupNumber = m_eclNOptimal->getGroupNumber();
158
159 //..2D histogram of nOptimal for each group and energy point
160 m_nOptimal2D = m_eclNOptimal->getNOptimal();
161 }
162}
163
165{
166 B2DEBUG(175, "ECLCRSplitterModule::event()");
167
168 // Fill a vector that can be used to map cellid -> store array position for eclCalDigits.
169 memset(&m_StoreArrPosition[0], -1, m_StoreArrPosition.size() * sizeof m_StoreArrPosition[0]);
170 for (int i = 0; i < m_eclCalDigits.getEntries(); i++) {
171 m_StoreArrPosition[m_eclCalDigits[i]->getCellId()] = i;
172 }
173
174 // Fill a vector that can be used to map cellid -> store array position for eclLocalMaximums.
175 memset(&m_StoreArrPositionLM[0], -1, m_StoreArrPositionLM.size() * sizeof m_StoreArrPositionLM[0]);
176 for (int i = 0; i < m_eclLocalMaximums.getEntries(); i++) {
177 m_StoreArrPositionLM[m_eclLocalMaximums[i]->getCellId()] = i;
178 }
179
180 // Loop over all connected regions
181 for (auto& aCR : m_eclConnectedRegions) {
182 // list theat will hold all cellids in this connected region
183 m_cellIdInCR.clear();
184
185 const unsigned int entries = (aCR.getRelationsWith<ECLCalDigit>(eclCalDigitArrayName())).size();
186
187 m_cellIdInCR.resize(entries);
188
189 // Fill all calDigits ids in this CR into a vector to make them 'find'-able.
190 int i = 0;
191 for (const auto& caldigit : aCR.getRelationsWith<ECLCalDigit>(eclCalDigitArrayName())) {
192 m_cellIdInCR[i] = caldigit.getCellId();
193 ++i;
194 }
195
196 // Split and reconstruct the showers in this connected regions.
198
199 } // end auto& aCR
200
201}
202
203
205{
206// if (m_tg2OptimalNumberOfDigitsForEnergy) delete m_tg2OptimalNumberOfDigitsForEnergy;
207}
208
209
211{
214}
215
217{
218
219 // Get the event background level
220 const int bkgdcount = m_eventLevelClusteringInfo->getNECLCalDigitsOutOfTime();
221 double backgroundLevel = 0.0; // from out of time digit counting
222 if (m_fullBkgdCount > 0) {
223 backgroundLevel = static_cast<double>(bkgdcount) / static_cast<double>(m_fullBkgdCount);
224 }
225
226 // Get the number of LMs in this CR
227 const int nLocalMaximums = aCR.getRelationsWith<ECLLocalMaximum>(eclLocalMaximumArrayName()).size();
228
229 B2DEBUG(170, "ECLCRSplitterModule::splitConnectedRegion: nLocalMaximums = " << nLocalMaximums);
230
231 // Three cases:
232 // 1) There is no local maximum (can only happen if the CR seed is not a LM itself).
233 // 2) There is exactly one local maximum.
234 // 3) There are more than one, typically two or three, local maxima and we have to share energy between them.
235 // If there are more than m_maxSplits local maxima, the m_maxSplits highest energy local maxima will be used.
236
237 // ---------------------------------------------------------------------
238 if (nLocalMaximums == 1) {
239
240 // Create a shower.
241 const auto aECLShower = m_eclShowers.appendNew();
242
243 // Add relation to the CR.
244 aECLShower->addRelationTo(&aCR);
245
246 // Find the highest energetic crystal in this CR or use the LM.
247 double weightSum = 0.0;
248
249 // Add relation to the LM.
251 aECLShower->addRelationTo(locmaxvector[0]);
252
253 const int locmaxcellid = locmaxvector[0]->getCellId();
254 const int pos = m_StoreArrPosition[locmaxcellid];
255 double highestEnergyID = (m_eclCalDigits[pos])->getCellId();
256 double highestEnergy = (m_eclCalDigits[pos])->getEnergy();
257 double highestEnergyTime = (m_eclCalDigits[pos])->getTime();
258 double highestEnergyTimeResolution = (m_eclCalDigits[pos])->getTimeResolution();
259
260 // Get a first estimation of the energy using 3x3 neighbours.
261 const double energyEstimation = estimateEnergy(highestEnergyID);
262
263 // Check if 21 would be better in the present background conditions:
264 ECLNeighbours* neighbourMap; // FIXME pointer needed?
265 int nNeighbours = getNeighbourMap(energyEstimation, backgroundLevel);
266 if (nNeighbours == 9 and !m_useOptimalNumberOfDigitsForEnergy) neighbourMap = m_NeighbourMap9;
267 else neighbourMap = m_NeighbourMap21;
268
269 // Add neighbours and weights for the shower.
270 std::vector<ECLCalDigit> digits;
271 std::vector<double> weights;
272 for (auto& neighbourId : neighbourMap->getNeighbours(highestEnergyID)) {
273 const auto it = std::find(m_cellIdInCR.begin(), m_cellIdInCR.end(),
274 neighbourId); // check if the neighbour is in the list for this CR
275 if (it == m_cellIdInCR.end()) continue; // not in this CR
276
277 const int neighbourpos = m_StoreArrPosition[neighbourId];
278 digits.push_back(*m_eclCalDigits[neighbourpos]); // list of digits for position reconstruction
279 weights.push_back(1.0); // list of weights (all 1 in this case for now)
280 weightSum += 1.0;
281
282 aECLShower->addRelationTo(m_eclCalDigits[neighbourpos], 1.0); // add digits to this shower, weight = 1
283 }
284
285 // Get position.
286 const B2Vector3D& showerposition = Belle2::ECL::computePositionLiLo(digits, weights, m_liloParameters);
287 aECLShower->setTheta(showerposition.Theta());
288 aECLShower->setPhi(showerposition.Phi());
289 aECLShower->setR(showerposition.Mag());
290
291 // Get Energy, if requested, set some weights to zero for energy calculation.
292 double showerEnergy = 0.0;
294
295 // Get the optimal number of neighbours for this crystal and energy
296 std::vector<int> nOptimalVec = getOptimalNumberOfDigits(highestEnergyID, energyEstimation);
297 const unsigned int nOptimal = static_cast<unsigned int>(nOptimalVec[0]);
298 aECLShower->setNominalNumberOfCrystalsForEnergy(static_cast<double>(nOptimal));
299
300 // Store the indices used; will be needed later for energy corrections.
301 // Also store energy, used along with cellID to find leakage corrections.
302 aECLShower->setNOptimalGroupIndex(nOptimalVec[1]);
303 aECLShower->setNOptimalEnergyBin(nOptimalVec[2]);
304 aECLShower->setNOptimalEnergy(energyEstimation);
305
306 // Get the list of crystals used for the energy calculation
307 std::vector< std::pair<unsigned int, double>> listCrystalPairs; // cell id and weighted reconstructed energy
308 listCrystalPairs.resize(digits.size()); //resize to number of all crystals in cluster
309
310 std::vector < std::pair<double, double> > weighteddigits;
311 weighteddigits.resize(digits.size());
312 for (unsigned int i = 0; i < digits.size(); ++i) {
313 weighteddigits.at(i) = std::make_pair((digits.at(i)).getEnergy(), weights.at(i));
314 listCrystalPairs.at(i) = std::make_pair((digits.at(i)).getCellId(), weights.at(i) * (digits.at(i)).getEnergy());
315 }
316
317 // sort the listCrystals and keep the n highest in descending order
318 std::sort(listCrystalPairs.begin(), listCrystalPairs.end(), [](const auto & left, const auto & right) {
319 return left.second > right.second;
320 });
321 std::vector< unsigned int> listCrystals; //cell id
322
323 for (unsigned int i = 0; i < digits.size(); ++i) {
324 if (i < nOptimal) {
325 listCrystals.push_back(listCrystalPairs[i].first);
326 }
327 }
328
329 aECLShower->setNumberOfCrystalsForEnergy(static_cast<double>(listCrystals.size()));
330 aECLShower->setListOfCrystalsForEnergy(listCrystals);
331
332 showerEnergy = getEnergySum(weighteddigits, nOptimal);
333 B2DEBUG(175, "Shower Energy (1): " << showerEnergy);
334
335 } else {
336 showerEnergy = Belle2::ECL::computeEnergySum(digits, weights);
337 }
338
339 aECLShower->setEnergy(showerEnergy);
340 aECLShower->setEnergyRaw(showerEnergy);
341 aECLShower->setEnergyHighestCrystal(highestEnergy);
342 aECLShower->setTime(highestEnergyTime);
343 aECLShower->setDeltaTime99(highestEnergyTimeResolution);
344 aECLShower->setNumberOfCrystals(weightSum);
345 aECLShower->setCentralCellId(highestEnergyID);
346
347 B2DEBUG(175, "theta = " << showerposition.Theta());
348 B2DEBUG(175, "phi = " << showerposition.Phi());
349 B2DEBUG(175, "R = " << showerposition.Mag());
350 B2DEBUG(175, "energy = " << showerEnergy);
351 B2DEBUG(175, "time = " << highestEnergyTime);
352 B2DEBUG(175, "time resolution = " << highestEnergyTimeResolution);
353 B2DEBUG(175, "neighbours = " << nNeighbours);
354 B2DEBUG(175, "backgroundLevel = " << backgroundLevel);
355
356 // Fill shower Ids
357 aECLShower->setShowerId(1); // always one (only this single shower in the CR)
358 aECLShower->setHypothesisId(Belle2::ECLShower::c_nPhotons);
359 aECLShower->setConnectedRegionId(aCR.getCRId());
360
361 // Add relations of all CalDigits of the CR to the local maximum (here: all weights = 1).
362 const int posLM = m_StoreArrPositionLM[locmaxcellid];
363 for (const auto& aDigit : aCR.getRelationsWith<ECLCalDigit>()) {
364 const int posDigit = m_StoreArrPosition[aDigit.getCellId()];
365 m_eclLocalMaximums[posLM]->addRelationTo(m_eclCalDigits[posDigit], 1.0);
366 }
367
368 } // end case with one LM
369 else { // More than one LM, energy must be split. This algorithm is inspired by BaBar code.
370
371 // check if we have too many local maximums. if yes: limit to user set maximum
372 // create a vector with all local maximums and its crystal energies
373 std::vector<std::pair<ECLLocalMaximum, double>> lm_energy_vector;
374 for (auto& aLocalMaximum : aCR.getRelationsWith<ECLLocalMaximum>(eclLocalMaximumArrayName())) {
375 const int cellid = aLocalMaximum.getCellId();
376 const int pos = m_StoreArrPosition[cellid];
377 const double digitenergy = m_eclCalDigits[pos]->getEnergy();
378 lm_energy_vector.push_back(std::pair<ECLLocalMaximum, double>(aLocalMaximum, digitenergy));
379 };
380
381 // sort this vector in descending order and keep only up to m_maxSplits entries
382 if (lm_energy_vector.size() >= static_cast<size_t>(m_maxSplits)) {
383 std::sort(lm_energy_vector.begin(), lm_energy_vector.end(), [](const std::pair<ECLLocalMaximum, double>& x,
384 const std::pair<ECLLocalMaximum, double>& y) {
385 return x.second > y.second;
386 });
387
388 lm_energy_vector.resize(m_maxSplits);
389 }
390
391 std::vector<ECLCalDigit> digits;
392 std::vector<double> weights;
393 std::map<int, B2Vector3D> centroidList; // key = cellid, value = centroid position
394 std::map<int, double> centroidEnergyList; // key = cellid, value = centroid position
395 std::map<int, B2Vector3D> allPoints; // key = cellid, value = digit position
396 std::map<int, std::vector < double > > weightMap; // key = locmaxid, value = vector of weights
397 std::vector < ECLCalDigit > digitVector; // the order of weights in weightMap must be the same
398
399 // Fill the maxima positions in a map
400 std::map<int, B2Vector3D> localMaximumsPoints; // key = locmaxid, value = maximum position
401 std::map<int, B2Vector3D> centroidPoints; // key = locmaxid (as index), value = centroid position
402
403 for (auto& aLocalMaximum : lm_energy_vector) {
404
405 int cellid = aLocalMaximum.first.getCellId();
406
407 // Get the position of this crystal and fill it in two maps.
408 B2Vector3D vectorPosition = m_geom->GetCrystalPos(cellid - 1);
409 localMaximumsPoints.insert(std::map<int, B2Vector3D>::value_type(cellid, vectorPosition));
410 centroidPoints.insert(std::map<int, B2Vector3D>::value_type(cellid, vectorPosition));
411 }
412
413 // The following will be done iteratively. Empty clusters after splitting will be removed, and the procedure will be repeated.
414 bool iterateclusters = true;
415 do {
416 digits.clear();
417 weights.clear();
418 centroidList.clear();
419 centroidEnergyList.clear();
420 allPoints.clear();
421 weightMap.clear();
422 digitVector.clear();
423
424 // Fill all digits from this CR in a map
425 for (auto& aCalDigit : aCR.getRelationsWith<ECLCalDigit>(eclCalDigitArrayName())) {
426 const int cellid = aCalDigit.getCellId();
427 // get the position of this crystal and fill them in a map
428 B2Vector3D vectorPosition = m_geom->GetCrystalPos(cellid - 1);
429 allPoints.insert(std::map<int, B2Vector3D>::value_type(cellid, vectorPosition));
430 digits.push_back(aCalDigit);
431 }
432
433 for (const auto& digitpoint : allPoints) {
434 const int cellid = digitpoint.first;
435 const int pos = m_StoreArrPosition[cellid];
436 digitVector.push_back(*m_eclCalDigits[pos]);
437 }
438
439
440 // -----------------------------------------------------------------------------------------
441 // The 'heart' of the splitter
442 // Start with each local maximum and set it to the first centroid position. Then iterate over all ECLCalDigits
443 // in this CR and calculate the weighted distances to the local maximum
444 int nIterations = 0;
445 double centroidShiftAverage = 0.0;
446 // double lastcentroidShiftAverage = 0.0;
447
448 do {
449 B2DEBUG(175, "Iteration: #" << nIterations << " (of max. " << m_maxIterations << ")");
450
451 centroidShiftAverage = 0.0;
452 if (nIterations == 0) {
453 centroidList.clear();
454 centroidEnergyList.clear();
455 weightMap.clear();
456 }
457
458 // Loop over all local maximums points, each one will become a shower!
459 for (const auto& locmaxpoint : localMaximumsPoints) {
460
461 // cell id of this local maximum
462 const int locmaxcellid = locmaxpoint.first;
463
464 // clear weights vector
465 weights.clear();
466
467 // if this is the first iteration the shower energy is not know, take the local maximum energy * 1.5.
468 if (nIterations == 0) {
469 const int pos = m_StoreArrPosition[locmaxcellid];
470 const double locmaxenergy = m_eclCalDigits[pos]->getEnergy();
471 centroidEnergyList[locmaxcellid] = 1.5 * locmaxenergy;
472 }
473
474 B2DEBUG(175, "local maximum cellid: " << locmaxcellid);
475
476 //-------------------------------------------------------------------
477 // Loop over all digits. They get a weight using the distance to the respective centroid.
478 for (const auto& digitpoint : allPoints) {
479
480 // cellid and position of this digit
481 const int digitcellid = digitpoint.first;
482 B2Vector3D digitpos = digitpoint.second;
483
484 const int pos = m_StoreArrPosition[digitcellid];
485 const double digitenergy = m_eclCalDigits[pos]->getEnergy();
486
487 double weight = 0.0;
488 double energy = 0.0;
489 double distance = 0.0;
490 double distanceEnergySum = 0.0;
491
492 // Loop over all centroids
493 for (const auto& centroidpoint : centroidPoints) {
494
495 // cell id and position of this centroid
496 const int centroidcellid = centroidpoint.first;
497 B2Vector3D centroidpos = centroidpoint.second;
498
499 double thisdistance = 0.;
500
501 // in the first iteration, this distance is really zero, avoid floating point problems
502 if (nIterations == 0 and digitcellid == centroidcellid) {
503 thisdistance = 0.0;
504 } else {
505 B2Vector3D vectorDistance = ((centroidpos) - (digitpos));
506 thisdistance = vectorDistance.Mag();
507 }
508
509 // energy of the centroid aka locmax
510 const int thispos = m_StoreArrPosition[centroidcellid];
511 const double thisenergy = m_eclCalDigits[thispos]->getEnergy();
512
513 // Not the most efficienct way to get this information, but not worth the thinking yet:
514 if (locmaxcellid == centroidcellid) {
515 distance = thisdistance;
516 energy = thisenergy;
517 }
518
519 // Get the product of distance and energy.
520 const double expfactor = exp(-m_expConstant * thisdistance / c_molierRadius);
521 distanceEnergySum += (thisenergy * expfactor);
522
523 } // end centroidPoints
524
525 // Calculate the weight for this digits for this local maximum.
526 if (distanceEnergySum > 0.0) {
527 const double expfactor = exp(-m_expConstant * distance / c_molierRadius);
528 weight = energy * expfactor / distanceEnergySum;
529 } else {
530 weight = 0.0;
531 }
532
533 // Check if the weighted energy is above threshold
534 if ((digitenergy * weight) < m_minimumSharedEnergy) {
535 weight = 0.0;
536 }
537
538 // Check for rounding problems larger than unity
539 if (weight > 1.0) {
540 B2WARNING("ECLCRSplitterModule::splitConnectedRegion: Floating point glitch, weight for this digit " << weight <<
541 ", resetting it to 1.0.");
542 weight = 1.0;
543 }
544
545 // Fill the weight for this digits and this local maximum.
546 B2DEBUG(175, " cellid: " << digitcellid << ", energy: " << digitenergy << ", weight: " << weight << ", distance: " << distance);
547 weights.push_back(weight);
548
549 } // end allPoints
550
551 // Get the old centroid position.
552 B2Vector3D oldCentroidPos = (centroidPoints.find(locmaxcellid))->second;
553
554 // Calculate the new centroid position.
555 B2Vector3D newCentroidPos = Belle2::ECL::computePositionLiLo(digits, weights, m_liloParameters);
556
557 // Calculate new energy
558 const double newEnergy = Belle2::ECL::computeEnergySum(digits, weights);
559
560 // Calculate the shift of the centroid position for this local maximum.
561 const B2Vector3D centroidShift = (oldCentroidPos - newCentroidPos);
562
563 // Save the new centroid position (but dont update yet!), also save the weights and energy.
564 centroidList[locmaxcellid] = newCentroidPos;
565 weightMap[locmaxcellid] = weights;
566
567 B2DEBUG(175, "--> inserting new energy: " << newEnergy << " for local maximum " << locmaxcellid);
568 centroidEnergyList[locmaxcellid] = newEnergy;
569 double showerenergy = (*centroidEnergyList.find(locmaxcellid)).second / Belle2::Unit::MeV;
570 B2DEBUG(175, "--> new energy = " << showerenergy << " MeV");
571
572 // Add this to the average centroid shift.
573 centroidShiftAverage += centroidShift.Mag();
574
575 // Debugging output
576 B2DEBUG(175, " old centroid: " << oldCentroidPos.X() << " cm, " << oldCentroidPos.Y() << " cm, " << oldCentroidPos.Z() <<
577 "cm");
578 B2DEBUG(175, " new centroid: " << newCentroidPos.X() << " cm, " << newCentroidPos.Y() << " cm, " << newCentroidPos.Z() <<
579 "cm");
580 B2DEBUG(175, " centroid shift: " << centroidShift.Mag() << " cm");
581
582 } // end localMaximumsPoints
583
584 // Get the average centroid shift.
585 centroidShiftAverage /= static_cast<double>(nLocalMaximums);
586 // lastcentroidShiftAverage = centroidShiftAverage;
587 B2DEBUG(175, "--> average centroid shift: " << centroidShiftAverage << " cm (tolerance is " << m_shiftTolerance << " cm)");
588
589 // Update centroid positions for the next round
590 for (const auto& locmaxpoint : localMaximumsPoints) {
591 centroidPoints[locmaxpoint.first] = (centroidList.find(locmaxpoint.first))->second;
592 }
593
594 ++nIterations;
595
596 } while (nIterations < m_maxIterations and centroidShiftAverage > m_shiftTolerance);
597 // DONE!
598
599 // check that local maxima are still local maxima
600 std::vector<int> markfordeletion;
601 iterateclusters = false;
602 for (const auto& locmaxpoint : localMaximumsPoints) {
603
604 // Get locmax cellid
605 const int locmaxcellid = locmaxpoint.first;
606 const int pos = m_StoreArrPosition[locmaxcellid];
607
608 B2DEBUG(175, "locmaxcellid: " << locmaxcellid);
609 const double lmenergy = m_eclCalDigits[pos]->getEnergy();
610 B2DEBUG(175, "ok: ");
611
612 // Get the weight vector.
613 std::vector < double > myWeights = (*weightMap.find(locmaxcellid)).second;
614
615 for (unsigned int i = 0; i < digitVector.size(); ++i) {
616
617 const ECLCalDigit dig = digitVector[i];
618 const double weight = myWeights[i];
619 const int cellid = dig.getCellId();
620 const double energy = dig.getEnergy();
621
622 // two ways to fail:
623 // 1) another cell has more energy: energy*weight > lmenergy and cellid != locmaxcellid
624 // 2) local maximum has cell has less than threshold energy left: energy*weight < m_threshold and cellid == locmaxcellid
625 if ((energy * weight > lmenergy and cellid != locmaxcellid) or (energy * weight < m_threshold and cellid == locmaxcellid)) {
626 markfordeletion.push_back(locmaxcellid);
627 iterateclusters = true;
628 continue;
629 }
630 }
631 }
632
633 // delete LMs
634 for (const auto lmid : markfordeletion) {
635 localMaximumsPoints.erase(lmid);
636 centroidPoints.erase(lmid);
637 }
638
639 } while (iterateclusters);
640
641 // Create the ECLShower objects, one per LocalMaximumPoints
642 unsigned int iShower = 1;
643 for (const auto& locmaxpoint : localMaximumsPoints) {
644
645 const int locmaxcellid = locmaxpoint.first;
646 const int posLM = m_StoreArrPositionLM[locmaxcellid];
647
648 // Create a shower
649 const auto aECLShower = m_eclShowers.appendNew();
650
651 // Use the same method for the estimate (3x3).
652 const double energyEstimation = estimateEnergy(locmaxcellid);
653
654 // Get the neighbour list.
655 ECLNeighbours* neighbourMap; // FIXME need pointer?
656 int nNeighbours = getNeighbourMap(energyEstimation, backgroundLevel);
657 if (nNeighbours == 9 and !m_useOptimalNumberOfDigitsForEnergy) neighbourMap = m_NeighbourMap9;
658 else neighbourMap = m_NeighbourMap21;
659
660 // Get the neighbour list.
661 std::vector<short int> neighbourlist = neighbourMap->getNeighbours(locmaxcellid);
662
663 // Get the weight vector.
664 std::vector < double > myWeights = (*weightMap.find(locmaxcellid)).second;
665
666 // Loop over all digits.
667 std::vector<ECLCalDigit> newdigits;
668 std::vector<double> newweights;
669 double highestEnergy = 0.;
670 double highestEnergyTime = 0.;
671 double highestEnergyTimeResolution = 0.;
672 double weightSum = 0.0;
673
674 for (unsigned int i = 0; i < digitVector.size(); ++i) {
675
676 const ECLCalDigit dig = digitVector[i];
677 const double weight = myWeights[i];
678
679 const int cellid = dig.getCellId();
680 const int pos = m_StoreArrPosition[cellid];
681
682 // Add weighted relations of all CalDigits to the local maximum.
683 m_eclLocalMaximums[posLM]->addRelationTo(m_eclCalDigits[pos], weight);
684
685 // Positive weight and in allowed neighbour list?
686 if (weight > 0.0) {
687 if (std::find(neighbourlist.begin(), neighbourlist.end(), cellid) != neighbourlist.end()) {
688
689 aECLShower->addRelationTo(m_eclCalDigits[pos], weight);
690
691 newdigits.push_back(dig);
692 newweights.push_back(weight);
693
694 weightSum += weight;
695
696 const double energy = dig.getEnergy();
697
698 if (energy * weight > highestEnergy) {
699 highestEnergy = energy * weight;
700 highestEnergyTime = dig.getTime();
701 highestEnergyTimeResolution = dig.getTimeResolution();
702 }
703 }
704 }
705 }
706
707 // Old position:
708 B2Vector3D* oldshowerposition = new B2Vector3D((centroidList.find(locmaxcellid))->second);
709
710 B2DEBUG(175, "old theta: " << oldshowerposition->Theta());
711 B2DEBUG(175, "old phi: " << oldshowerposition->Phi());
712 B2DEBUG(175, "old R: " << oldshowerposition->Mag());
713 B2DEBUG(175, "old energy: " << energyEstimation);
714 delete oldshowerposition;
715
716 // New position (with reduced number of neighbours)
717 // There are some cases where high backgrounds fake local maxima and the splitted centroid position is far
718 // away from the original LM cell... this will throw a (non fatal) error, and create a cluster with zero energy now).
719 B2Vector3D* showerposition = new B2Vector3D(Belle2::ECL::computePositionLiLo(newdigits, newweights, m_liloParameters));
720 aECLShower->setTheta(showerposition->Theta());
721 aECLShower->setPhi(showerposition->Phi());
722 aECLShower->setR(showerposition->Mag());
723
724 B2DEBUG(175, "new theta: " << showerposition->Theta());
725 B2DEBUG(175, "new phi: " << showerposition->Phi());
726 B2DEBUG(175, "new R: " << showerposition->Mag());
727 delete showerposition;
728
729 // Get Energy, if requested, set weights to zero for energy calculation.
730 double showerEnergy = 0.0;
732
733
734 // Get the optimal number of neighbours for this crystal and energy
735 std::vector<int> nOptimalVec = getOptimalNumberOfDigits(locmaxcellid, energyEstimation);
736 const unsigned int nOptimal = static_cast<unsigned int>(nOptimalVec[0]);
737 aECLShower->setNominalNumberOfCrystalsForEnergy(static_cast<double>(nOptimal));
738
739 // Store the indices used; will be needed later for energy corrections.
740 // Also store energy, used along with cellID to find leakage corrections.
741 aECLShower->setNOptimalGroupIndex(nOptimalVec[1]);
742 aECLShower->setNOptimalEnergyBin(nOptimalVec[2]);
743 aECLShower->setNOptimalEnergy(energyEstimation);
744
745 // Get the list of crystals used for the energy calculation
746 std::vector< std::pair<unsigned int, double>> listCrystalPairs; // cell id and weighted reconstructed energy
747 listCrystalPairs.resize(newdigits.size()); //resize to number of all crystals in cluster
748
749 std::vector < std::pair<double, double> > weighteddigits;
750 weighteddigits.resize(newdigits.size());
751 for (unsigned int i = 0; i < newdigits.size(); ++i) {
752 weighteddigits.at(i) = std::make_pair((newdigits.at(i)).getEnergy(), newweights.at(i));
753 listCrystalPairs.at(i) = std::make_pair((newdigits.at(i)).getCellId(), newweights.at(i) * (newdigits.at(i)).getEnergy());
754 }
755
756 // sort the listCrystals and keep the n highest in descending order
757 std::sort(listCrystalPairs.begin(), listCrystalPairs.end(), [](const auto & left, const auto & right) {
758 return left.second > right.second;
759 });
760
761 std::vector< unsigned int> listCrystals; //cell id
762
763 for (unsigned int i = 0; i < newdigits.size(); ++i) {
764 if (i < nOptimal) {
765 listCrystals.push_back(listCrystalPairs[i].first);
766 }
767 }
768
769 aECLShower->setNumberOfCrystalsForEnergy(static_cast<double>(listCrystals.size()));
770 aECLShower->setListOfCrystalsForEnergy(listCrystals);
771
772 showerEnergy = getEnergySum(weighteddigits, nOptimal);
773 B2DEBUG(175, "Shower Energy (2): " << showerEnergy);
774
775 } else {
776 showerEnergy = Belle2::ECL::computeEnergySum(newdigits, newweights);
777 }
778
779 aECLShower->setEnergy(showerEnergy);
780 aECLShower->setEnergyRaw(showerEnergy);
781 aECLShower->setEnergyHighestCrystal(highestEnergy);
782 aECLShower->setTime(highestEnergyTime);
783 aECLShower->setDeltaTime99(highestEnergyTimeResolution);
784 aECLShower->setNumberOfCrystals(weightSum);
785 aECLShower->setCentralCellId(locmaxcellid);
786
787 B2DEBUG(175, "new energy: " << showerEnergy);
788
789 // Get unique ID
790 aECLShower->setShowerId(iShower);
791 ++iShower;
792 aECLShower->setHypothesisId(Belle2::ECLShower::c_nPhotons);
793 aECLShower->setConnectedRegionId(aCR.getCRId());
794
795 // Add relation to the CR.
796 aECLShower->addRelationTo(&aCR);
797
798 // Add relation to the LM.
799 aECLShower->addRelationTo(m_eclLocalMaximums[posLM]);
800 }
801 }
802}
803
804int ECLSplitterN1Module::getNeighbourMap(const double energy, const double background)
805{
806 if (background <= 0.1) return 21;
807 else {
808 if (energy > 0.06 + 0.4 * background) return 21; // based on preliminary study TF, valid in barrel only (TF).
809 else return 9;
810 }
811}
812
813std::vector<int> ECLSplitterN1Module::getOptimalNumberOfDigits(const int cellid, const double energy)
814{
815
816 //..nOptimal depends on the energy bin, which in turn depends on ECL region
817 int iRegion = 1; // barrel
818 if (ECLElementNumbers::isForward(cellid)) {iRegion = 0;}
819 if (ECLElementNumbers::isBackward(cellid)) {iRegion = 2;}
820
821 //..Find the energy bin
822 int iEnergy = 0;
823 while (energy > m_eBoundaries[iRegion][iEnergy] and iEnergy < m_nEnergyBins - 1) {iEnergy++;}
824
825 //..Group number just depends on the cellID
826 int iGroup = m_groupNumber[cellid - 1];
827
828 //..Optimal number of crystals from energy and group numbers
829 int nOptimalNeighbours = (int)(0.5 + m_nOptimal2D.GetBinContent(iGroup + 1, iEnergy + 1));
830
831 //..Store these in a vector to return
832 std::vector<int> nOptimalVector;
833 nOptimalVector.push_back(nOptimalNeighbours);
834 nOptimalVector.push_back(iGroup);
835 nOptimalVector.push_back(iEnergy);
836
837 B2DEBUG(175, "ECLSplitterN1Module::getOptimalNumberOfDigits: cellID: " << cellid << " energy: " << energy << " iG: " << iGroup <<
838 " iE: " << iEnergy << " nOpt: " << nOptimalNeighbours);
839
840 return nOptimalVector;
841
842}
843
844double ECLSplitterN1Module::getEnergySum(std::vector < std::pair<double, double> >& weighteddigits, const unsigned int n)
845{
846
847 double energysum = 0.;
848
849 std::sort(weighteddigits.begin(), weighteddigits.end(), [](const auto & left, const auto & right) {
850 return left.first * left.second > right.first * right.second;
851 });
852
853
854 unsigned int min = n;
855 if (weighteddigits.size() < n) min = weighteddigits.size();
856
857 for (unsigned int i = 0; i < min; ++i) {
858 B2DEBUG(175, "getEnergySum: " << weighteddigits.at(i).first << " " << weighteddigits.at(i).second);
859 energysum += (weighteddigits.at(i).first * weighteddigits.at(i).second);
860 }
861 B2DEBUG(175, "getEnergySum: energysum=" << energysum);
862
863 return energysum;
864}
865
866
867double ECLSplitterN1Module::estimateEnergy(const int centerid)
868{
869
870 double energyEstimation = 0.0;
871
872 for (auto& neighbourId : m_NeighbourMap9->getNeighbours(centerid)) {
873
874 // Check if this neighbour is in this CR
875 const auto it = std::find(m_cellIdInCR.begin(), m_cellIdInCR.end(),
876 neighbourId); // check if the neighbour is in the list for this CR
877 if (it == m_cellIdInCR.end()) continue; // not in this CR
878
879 const int pos = m_StoreArrPosition[neighbourId];
880 const double energyNeighbour = m_eclCalDigits[pos]->getEnergy();
881
882 energyEstimation += energyNeighbour;
883 }
884
885 return energyEstimation;
886}
DataType Phi() const
The azimuth angle.
Definition: B2Vector3.h:151
DataType Z() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:435
DataType Theta() const
The polar angle.
Definition: B2Vector3.h:153
DataType X() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:431
DataType Y() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:433
DataType Mag() const
The magnitude (rho in spherical coordinate system).
Definition: B2Vector3.h:159
Class to store calibrated ECLDigits: ECLCalDigits.
Definition: ECLCalDigit.h:23
int getCellId() const
Get Cell ID.
Definition: ECLCalDigit.h:118
double getEnergy() const
Get Calibrated Energy.
Definition: ECLCalDigit.h:123
double getTimeResolution() const
Get Calibrated Time Resolution.
Definition: ECLCalDigit.h:168
double getTime() const
Get Calibrated Time.
Definition: ECLCalDigit.h:163
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)
Definition: ECLShower.h:42
virtual const char * eclShowerArrayName() const
Default name ECLShowers.
ECL::ECLNeighbours * m_NeighbourMap9
Neighbour maps.
virtual const char * eventLevelClusteringInfoName() const
Name to be used for default option: EventLevelClusteringInfo.
double m_liloParameterB
lin-log parameter B
int m_maxIterations
Maximum number of iterations.
StoreArray< ECLShower > m_eclShowers
Store array: ECLShower.
ECL::ECLNeighbours * m_NeighbourMap21
5x5 neighbours excluding corners = 21
StoreArray< ECLConnectedRegion > m_eclConnectedRegions
Store array: ECLConnectedRegion.
double m_minimumSharedEnergy
Minimum shared energy.
const int m_nLeakReg
3 ECL regions: 0 = forward, 1 = barrel, 2 = backward
std::string m_positionMethod
Position calculation: lilo or linear.
void splitConnectedRegion(ECLConnectedRegion &aCR)
Split connected region into showers.
virtual void initialize() override
Initialize.
double m_cutDigitTimeResidualForEnergy
Maximum time residual to be included in the shower energy calculation.
StoreArray< ECLLocalMaximum > m_eclLocalMaximums
Store array: ECLLocalMaximum.
virtual void event() override
Event.
virtual const char * eclLocalMaximumArrayName() const
Default name ECLLocalMaximums.
double m_threshold
Local maximum threshold after splitting.
double estimateEnergy(const int centerid)
Estimate energy using 3x3 around central crystal.
virtual void endRun() override
End run.
virtual void terminate() override
Terminate.
std::vector< int > getOptimalNumberOfDigits(const int cellid, const double energy)
Get optimal number of digits (out of 21) based on first energy estimation and background level per ev...
ECL::ECLGeometryPar * m_geom
Geometry.
std::vector< int > m_cellIdInCR
list with all cellid of this connected region
double m_shiftTolerance
Tolerance level for centroid shifts.
DBObjPtr< ECLnOptimal > m_eclNOptimal
nOptimal payload
virtual void beginRun() override
Begin run.
int m_fullBkgdCount
Number of expected background digits at full background, FIXME: move to database.
int m_nEnergyBins
number of energies bins in nOptimal payload
const double c_molierRadius
Constant RM (Molier Radius) from exp(-a*dist/RM), http://pdg.lbl.gov/2009/AtomicNuclearProperties/HTM...
virtual const char * eclConnectedRegionArrayName() const
Default name ECLConnectedRegions.
int m_maxSplits
Maximum number of splits.
int m_useOptimalNumberOfDigitsForEnergy
Optimize the number of neighbours for energy calculations.
int getNeighbourMap(const double energy, const double background)
Get number of neighbours based on first energy estimation and background level per event.
double m_cutDigitEnergyForEnergy
Minimum digit energy to be included in the shower energy calculation.
virtual const char * eclCalDigitArrayName() const
Default name ECLCalDigits.
std::vector< double > m_liloParameters
lin-log parameters A, B, and C
std::vector< int > m_StoreArrPosition
vector (ECLElementNumbers::c_NCrystals + 1 entries) with cell id to store array positions
std::vector< int > m_StoreArrPositionLM
vector (ECLElementNumbers::c_NCrystals + 1 entries) with cell id to store array positions for LM
double m_liloParameterA
lin-log parameter A
StoreObjPtr< EventLevelClusteringInfo > m_eventLevelClusteringInfo
Store object pointer: EventLevelClusteringInfo.
std::vector< std::vector< float > > m_eBoundaries
energy boundaries each region
TH2F m_nOptimal2D
2D hist of nOptimal for Ebin vs groupID
std::vector< int > m_groupNumber
group number for each crystal
double getEnergySum(std::vector< std::pair< double, double > > &weighteddigits, const unsigned int n)
Get energy sum for weighted entries.
StoreArray< ECLCalDigit > m_eclCalDigits
Store array: ECLCalDigit.
double m_expConstant
Constant a from exp(-a*dist/RM), 1.5 to 2.5.
double m_liloParameterC
lin-log parameter C
static ECLGeometryPar * Instance()
Static method to get a reference to the ECLGeometryPar instance.
ROOT::Math::XYZVector GetCrystalPos(int cid)
The Position of crystal.
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.
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
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]
Definition: Unit.h:70
static const double keV
[kiloelectronvolt]
Definition: Unit.h:113
static const double MeV
[megaelectronvolt]
Definition: Unit.h:114
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
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.
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.