Belle II Software  release-06-00-14
ECLShowerCorrectorModule.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/eclShowerCorrection/ECLShowerCorrectorModule.h>
11 
12 //..Framework
13 #include <framework/logging/Logger.h>
14 
15 //..ECL
16 #include <ecl/dbobjects/ECLLeakageCorrections.h>
17 #include <ecl/dataobjects/ECLShower.h>
18 #include <ecl/geometry/ECLLeakagePosition.h>
19 
20 using namespace Belle2;
21 using namespace ECL;
22 
23 //-----------------------------------------------------------------
24 // Register the Module
25 //-----------------------------------------------------------------
26 REG_MODULE(ECLShowerCorrector)
27 REG_MODULE(ECLShowerCorrectorPureCsI)
28 
29 //-----------------------------------------------------------------
30 // Implementation
31 //-----------------------------------------------------------------
32 
34  m_eclShowers(eclShowerArrayName()),
35  m_eclLeakageCorrections("ECLLeakageCorrections")
36 {
37 
38  // Set description
39  setDescription("ECLShowerCorrectorModule: Corrects energy of ECLShowers and highest energy crystal for shower leakage, beam backgrounds, and clustering");
40  setPropertyFlags(c_ParallelProcessingCertified);
41 
42 }
43 
45 {
46  ;
47 }
48 
50 {
51  B2DEBUG(28, "ECLShowerCorrectorModule::initialize()");
52 
53  //..Register in datastore
54  m_eclShowers.registerInDataStore(eclShowerArrayName());
55 
56  //..Class to find cellID and position within crystal from theta and phi
57  leakagePosition = new ECLLeakagePosition();
58 
59 }
60 
62 {
63  //-----------------------------------------------------------------
64  //..Read in leakage corrections from database
65  if (m_eclLeakageCorrections.hasChanged()) {
66 
67  //..Vectors of log(E) for each region
68  std::vector<float> logEnergiesFwd = m_eclLeakageCorrections->getlogEnergiesFwd();
69  std::vector<float> logEnergiesBrl = m_eclLeakageCorrections->getlogEnergiesBrl();
70  std::vector<float> logEnergiesBwd = m_eclLeakageCorrections->getlogEnergiesBwd();
71 
72  //..Adjust the size of the vector of log energies to match the number in the payload
73  nEnergies = logEnergiesBrl.size();
74  B2INFO("ECLShowerCorrector beginRun: Number of energies = " << nEnergies);
75  leakLogE.resize(nLeakReg, std::vector<float>(nEnergies, 0.));
76 
77  //..Copy values to leakLogE
78  for (int ie = 0; ie < nEnergies; ie++) {
79  leakLogE[0][ie] = logEnergiesFwd[ie];
80  leakLogE[1][ie] = logEnergiesBrl[ie];
81  leakLogE[2][ie] = logEnergiesBwd[ie];
82  }
83 
84  //..Position and nCrys dependent corrections
85  thetaCorrection = m_eclLeakageCorrections->getThetaCorrections();
86  phiCorrection = m_eclLeakageCorrections->getPhiCorrections();
87  crysCorrection = m_eclLeakageCorrections->getnCrystalCorrections();
88 
89  //..Relevant parameters
90  nPositionBins = thetaCorrection.GetNbinsY();
91  nCrysMax = crysCorrection.GetNbinsY();
92  nXBins = nThetaID * nEnergies;
93  }
94 }
95 
97 {
98 
99  //-----------------------------------------------------------------
100  //..Loop over all ECLShowers.
101  for (auto& eclShower : m_eclShowers) {
102 
103  //..Only want to correct EM showers
104  if (eclShower.getHypothesisId() != ECLShower::c_nPhotons) {break;}
105 
106  //-----------------------------------------------------------------
107  //..Shower quantities needed to find the correction.
108  const double energyRaw = eclShower.getEnergy();
109  const double energyRawHighest = eclShower.getEnergyHighestCrystal();
110 
111  //..Crystals used to calculate energy
112  int nForEnergy = (int)(eclShower.getNumberOfCrystalsForEnergy() + 0.001);
113  const int nCrys = std::min(nCrysMax, nForEnergy);
114 
115  //..Location starting point
116  const int icellIDMaxE = eclShower.getCentralCellId();
117  const float thetaLocation = eclShower.getTheta();
118  const float phiLocation = eclShower.getPhi();
119 
120  //-----------------------------------------------------------------
121  //..Location of shower. cellID = positionVector[0] is for debugging only
122  std::vector<int> positionVector = leakagePosition->getLeakagePosition(icellIDMaxE, thetaLocation, phiLocation, nPositionBins);
123  const int iThetaID = positionVector[1];
124  const int iRegion = positionVector[2];
125  const int iThetaBin = positionVector[3];
126  const int iPhiBin = positionVector[4];
127  const int iPhiMech = positionVector[5];
128 
129  //-----------------------------------------------------------------
130  //..The correction is a function of corrected energy, so will need to iterate
131  double correction = 0.96; // typical correction as starting point
132  for (int iter = 0; iter < 2; iter++) {
133 
134  //..Energy points that bracket this value
135  float logEnergy = log(energyRaw / correction);
136  int ie0 = 0;
137  if (logEnergy < leakLogE[iRegion][0]) {
138  ie0 = 0;
139  } else if (logEnergy > leakLogE[iRegion][nEnergies - 1]) {
140  ie0 = nEnergies - 2;
141  } else {
142  while (logEnergy > leakLogE[iRegion][ie0 + 1]) {ie0++;}
143  }
144 
145  //..Correction from lower energy point.
146  // The following include +1 because first histogram bin is 1 not 0.
147  int iXBin = iThetaID + ie0 * nThetaID + 1; // thetaID / energy bin
148  double thetaCor = thetaCorrection.GetBinContent(iXBin, iThetaBin + 1);
149  double phiCor = phiCorrection.GetBinContent(iXBin + iPhiMech * nXBins, iPhiBin + 1);
150  double crysCor = crysCorrection.GetBinContent(iXBin, nCrys + 1);
151  const double cor0 = thetaCor * phiCor * crysCor;
152 
153  //..Correction from upper energy point
154  iXBin = iThetaID + (ie0 + 1) * nThetaID + 1;
155  thetaCor = thetaCorrection.GetBinContent(iXBin, iThetaBin + 1);
156  phiCor = phiCorrection.GetBinContent(iXBin + iPhiMech * nXBins, iPhiBin + 1);
157  crysCor = crysCorrection.GetBinContent(iXBin, nCrys + 1);
158  const double cor1 = thetaCor * phiCor * crysCor;
159 
160  //..Interpolate (in logE)
161  correction = cor0 + (cor1 - cor0) * (logEnergy - leakLogE[iRegion][ie0]) / (leakLogE[iRegion][ie0 + 1] - leakLogE[iRegion][ie0]);
162  }
163 
164  //-----------------------------------------------------------------
165  //..Apply correction
166  const double correctedEnergy = energyRaw / correction;
167  const double correctedEnergyHighest = energyRawHighest / correction;
168  B2DEBUG(28, "Correction factor=" << correction << ", correctedEnergy=" << correctedEnergy << ", correctedEnergyHighest=" <<
169  correctedEnergyHighest);
170 
171  //..Set the corrected energies
172  eclShower.setEnergy(correctedEnergy);
173  eclShower.setEnergyHighestCrystal(correctedEnergyHighest);
174 
175  } // end loop over showers
176 
177 }
178 
180 {
181  ;
182 }
183 
185 {
186  ;
187 }
DB object to store leakage corrections, including nCrys dependence
Class to perform the shower correction.
virtual void initialize() override
Initialize.
virtual void event() override
Event.
virtual void endRun() override
End run.
virtual void terminate() override
Terminate.
virtual void beginRun() override
Begin run.
@ c_nPhotons
CR is split into n photons (N1)
Definition: ECLShower.h:37
Class to get position information for a cluster for leakage corrections.
Base class for Modules.
Definition: Module.h:72
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.