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