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