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