Belle II Software  release-08-01-10
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 #include <ecl/dataobjects/ECLElementNumbers.h>
20 
21 using namespace Belle2;
22 using namespace ECL;
23 
24 //-----------------------------------------------------------------
25 // Register the Module
26 //-----------------------------------------------------------------
27 REG_MODULE(ECLShowerCorrector);
28 REG_MODULE(ECLShowerCorrectorPureCsI);
29 
30 //-----------------------------------------------------------------
31 // Implementation
32 //-----------------------------------------------------------------
33 
35  m_eclShowers(eclShowerArrayName()),
36  m_eclLeakageCorrections("ECLLeakageCorrections")
37 {
38 
39  // Set description
40  setDescription("ECLShowerCorrectorModule: Corrects energy of ECLShowers and highest energy crystal for shower leakage, beam backgrounds, and clustering");
42 
43 }
44 
46 {
47  if (m_leakagePosition != nullptr)
48  delete m_leakagePosition;
49 }
50 
52 {
53  B2DEBUG(28, "ECLShowerCorrectorModule::initialize()");
54 
55  //..Register in datastore
56  m_eclShowers.registerInDataStore(eclShowerArrayName());
57  m_eventLevelClusteringInfo.registerInDataStore();
58 
59  //..Class to find cellID and position within crystal from theta and phi
61 
62 }
63 
65 {
66  //-----------------------------------------------------------------
67  //..Read in leakage corrections from database
68  if (m_eclLeakageCorrections.hasChanged()) {
69 
70  //..Vectors of log(E) for each region
71  std::vector<float> logEnergiesFwd = m_eclLeakageCorrections->getlogEnergiesFwd();
72  std::vector<float> logEnergiesBrl = m_eclLeakageCorrections->getlogEnergiesBrl();
73  std::vector<float> logEnergiesBwd = m_eclLeakageCorrections->getlogEnergiesBwd();
74 
75  //..Adjust the size of the vector of log energies to match the number in the payload
76  nEnergies = logEnergiesBrl.size();
77  B2INFO("ECLShowerCorrector beginRun: Number of energies = " << nEnergies);
78  leakLogE.resize(nLeakReg, std::vector<float>(nEnergies, 0.));
79 
80  //..Copy values to leakLogE
81  for (int ie = 0; ie < nEnergies; ie++) {
82  leakLogE[0][ie] = logEnergiesFwd[ie];
83  leakLogE[1][ie] = logEnergiesBrl[ie];
84  leakLogE[2][ie] = logEnergiesBwd[ie];
85  }
86 
87  //..Position dependent corrections
88  thetaCorrection = m_eclLeakageCorrections->getThetaCorrections();
89  phiCorrection = m_eclLeakageCorrections->getPhiCorrections();
90 
91  //..Relevant parameters
92  nPositionBins = thetaCorrection.GetNbinsY();
94  } else if (!m_eclLeakageCorrections.isValid()) {
95  B2FATAL("ECLShowerCorrectorModule: missing eclLeakageCorrections payload");
96  }
97 
98  //-----------------------------------------------------------------
99  //..Get correction histograms related to the nOptimal number of crystals.
100  // All three have energy bin as y, crystals group number as x.
101  // Energy bin and group number for the shower are found in
102  // eclSplitterN1 and are stored in the ECLShower dataobject.
103 
104  if (m_eclNOptimal.hasChanged()) {
105 
106  //..Bias is the difference between the peak energy in nOptimal crystals
107  // before bias correction and the mc true deposited energy.
108  m_bias = m_eclNOptimal->getBias();
109 
110  //..Log of the peak energy contained in nOptimal crystals after bias correction
111  m_logPeakEnergy = m_eclNOptimal->getLogPeakEnergy();
112 
113  //..peakFracEnergy is the peak energy after subtracting the beam bias
114  // divided by the generated photon energy.
115  m_peakFracEnergy = m_eclNOptimal->getPeakFracEnergy();
116  } else if (!m_eclNOptimal.isValid()) {
117  B2FATAL("ECLShowerCorrectorModule: missing eclNOptimal payload");
118  }
119 
120  //-----------------------------------------------------------------
121  //..Get correction histograms related to the nOptimal number of crystals.
122  // All three have energy bin as y, crystals group number as x.
123  // Energy bin and group number for the shower are found in
124  // eclSplitterN1 and are stored in the ECLShower dataobject.
125 
126  if (m_eclNOptimal.hasChanged()) {
127 
128  //..Bias is the difference between the peak energy in nOptimal crystals
129  // before bias correction and the mc true deposited energy.
130  m_bias = m_eclNOptimal->getBias();
131 
132  //..Log of the peak energy contained in nOptimal crystals after bias correction
133  m_logPeakEnergy = m_eclNOptimal->getLogPeakEnergy();
134 
135  //..peakFracEnergy is the peak energy after subtracting the beam bias
136  // divided by the generated photon energy.
137  m_peakFracEnergy = m_eclNOptimal->getPeakFracEnergy();
138  }
139 }
140 
142 {
143 
144  //-----------------------------------------------------------------
145  //..Loop over all ECLShowers.
146  for (auto& eclShower : m_eclShowers) {
147 
148  //..Only want to correct EM showers
149  if (eclShower.getHypothesisId() != ECLShower::c_nPhotons) {continue;}
150 
151  //..Will correct both raw cluster energy and energy of the center crystal
152  const double energyRaw = eclShower.getEnergy();
153  const double energyRawHighest = eclShower.getEnergyHighestCrystal();
154 
155  //-----------------------------------------------------------------
156  //..Correct for bias and peak value corresponding to nOptimal crystals
157  const int iGroup = eclShower.getNOptimalGroupIndex();
158  const int iEnergy = eclShower.getNOptimalEnergyBin();
159  const double e3x3 = eclShower.getNOptimalEnergy(); // only for debugging
160 
161  //..The optimal number of crystals for generated energy iEnergy and
162  // group of crystals iGroup is nOptimal. There are three corresponding bins for
163  // each of the 2D correction histograms m_logPeakEnergy, m_bias, and
164  // m_peakFracEnergy. For example, the bias for nOptimal crystals in group iGroup is
165  // m_bias(3*iGroup+1, iEnergy+1) for generated energy point iEnergy,
166  // m_bias(3*iGroup+2, iEnergy+1) for generated energy point iEnergy-1, and
167  // m_bias(3*iGroup+3, iEnergy+1) for generated energy point iEnergy+1
168  // This allows the correction to be interpolated to an arbitrary observed energy.
169  const int iy = iEnergy + 1;
170 
171  const int ixNom = 3 * iGroup + 1;
172  const int ixLowerE = ixNom + 1;
173  const int ixHigherE = ixNom + 2;
174 
175  const double logENom = m_logPeakEnergy.GetBinContent(ixNom, iy);
176  const double logELower = m_logPeakEnergy.GetBinContent(ixLowerE, iy);
177  const double logEHigher = m_logPeakEnergy.GetBinContent(ixHigherE, iy);
178 
179  const double biasNom = m_bias.GetBinContent(ixNom, iy);
180  const double biasLower = m_bias.GetBinContent(ixLowerE, iy);
181  const double biasHigher = m_bias.GetBinContent(ixHigherE, iy);
182 
183  const double peakNom = m_peakFracEnergy.GetBinContent(ixNom, iy);
184  const double peakLower = m_peakFracEnergy.GetBinContent(ixLowerE, iy);
185  const double peakHigher = m_peakFracEnergy.GetBinContent(ixHigherE, iy);
186 
187  //..Interpolate in log of raw energy
188  const double logESumN = log(energyRaw);
189 
190  double logEOther = logELower;
191  double biasOther = biasLower;
192  double peakOther = peakLower;
193  if (logESumN > logENom) {
194  logEOther = logEHigher;
195  biasOther = biasHigher;
196  peakOther = peakHigher;
197  }
198 
199  //..The nominal and "other" energies may be identical if this is the first or last energy
200  double bias = biasNom;
201  double peak = peakNom;
202  if (std::abs(logEOther - logENom) > 0.0001) {
203  bias = biasNom + (biasOther - biasNom) * (logESumN - logENom) / (logEOther - logENom);
204  peak = peakNom + (peakOther - peakNom) * (logESumN - logENom) / (logEOther - logENom);
205  }
206 
207  //..Correct raw energy for bias and peak
208  const double ePartialCorr = (energyRaw - bias) / peak;
209 
210  //-----------------------------------------------------------------
211  //..Shower quantities needed to find the correction.
212 
213  //..Location starting point
214  const int icellIDMaxE = eclShower.getCentralCellId();
215  const float thetaLocation = eclShower.getTheta();
216  const float phiLocation = eclShower.getPhi();
217 
218  //-----------------------------------------------------------------
219  //..Location of shower. cellID = positionVector[0] is for debugging only
220  std::vector<int> positionVector = m_leakagePosition->getLeakagePosition(icellIDMaxE, thetaLocation, phiLocation, nPositionBins);
221  const int iThetaID = positionVector[1];
222  const int iRegion = positionVector[2];
223  const int iThetaBin = positionVector[3];
224  const int iPhiBin = positionVector[4];
225  const int iPhiMech = positionVector[5];
226 
227  //-----------------------------------------------------------------
228  //..Energy points that bracket this value
229  float logEnergy = log(ePartialCorr);
230  int ie0 = 0;
231  if (logEnergy < leakLogE[iRegion][0]) {
232  ie0 = 0;
233  } else if (logEnergy > leakLogE[iRegion][nEnergies - 1]) {
234  ie0 = nEnergies - 2;
235  } else {
236  while (logEnergy > leakLogE[iRegion][ie0 + 1]) {ie0++;}
237  }
238 
239  //..Correction from lower energy point.
240  // The following include +1 because first histogram bin is 1 not 0.
241  int iXBin = iThetaID + ie0 * nThetaID + 1; // thetaID / energy bin
242  double thetaCor = thetaCorrection.GetBinContent(iXBin, iThetaBin + 1);
243  double phiCor = phiCorrection.GetBinContent(iXBin + iPhiMech * nXBins, iPhiBin + 1);
244  const double cor0 = thetaCor * phiCor;
245 
246  //..Correction from upper energy point
247  iXBin = iThetaID + (ie0 + 1) * nThetaID + 1;
248  thetaCor = thetaCorrection.GetBinContent(iXBin, iThetaBin + 1);
249  phiCor = phiCorrection.GetBinContent(iXBin + iPhiMech * nXBins, iPhiBin + 1);
250  const double cor1 = thetaCor * phiCor;
251 
252  //..Interpolate (in logE)
253  const double positionCor = cor0 + (cor1 - cor0) * (logEnergy - leakLogE[iRegion][ie0]) / (leakLogE[iRegion][ie0 + 1] -
254  leakLogE[iRegion][ie0]);
255 
256  //-----------------------------------------------------------------
257  //..Apply correction. Assume the same correction for the maximum energy crystal.
258  const double correctedEnergy = ePartialCorr / positionCor;
259  const double overallCorrection = energyRaw / correctedEnergy;
260  const double correctedEnergyHighest = energyRawHighest / overallCorrection;
261 
262  //..Set the corrected energies
263  eclShower.setEnergy(correctedEnergy);
264  eclShower.setEnergyHighestCrystal(correctedEnergyHighest);
265 
266  B2DEBUG(28, "ECLShowerCorrectorModule: cellID: " << positionVector[0] << " iG: " << iGroup << " iE: " << iEnergy << " Eraw: " <<
267  energyRaw << " E3x3: " << e3x3 << " Ecor: " << correctedEnergy);
268  B2DEBUG(28, " peakNom: " << peakNom << " biasNom: " << biasNom << " positionCor: " << positionCor << " overallCor: " <<
269  overallCorrection);
270 
271  } // end loop over showers
272 
273  //-----------------------------------------------------------------
274  //..Count number of showers in each region for EventLevelClusteringInfo
275  uint16_t nShowersPerRegion[nLeakReg] = {};
276  for (auto& eclShower : m_eclShowers) {
277  if (eclShower.getHypothesisId() != ECLShower::c_nPhotons) {continue;}
278  const int iCellId = eclShower.getCentralCellId();
279  if (ECLElementNumbers::isForward(iCellId)) {nShowersPerRegion[0]++;}
280  if (ECLElementNumbers::isBarrel(iCellId)) {nShowersPerRegion[1]++;}
281  if (ECLElementNumbers::isBackward(iCellId)) {nShowersPerRegion[2]++;}
282  }
283 
284  //..Store
285  m_eventLevelClusteringInfo->setNECLShowersFWD(nShowersPerRegion[0]);
286  m_eventLevelClusteringInfo->setNECLShowersBarrel(nShowersPerRegion[1]);
287  m_eventLevelClusteringInfo->setNECLShowersBWD(nShowersPerRegion[2]);
288 
289 }
290 
292 {
293  ;
294 }
295 
297 {
298  ;
299 }
static constexpr unsigned int nLeakReg
0 = forward, 1 = barrel, 2 = backward
TH2F m_logPeakEnergy
log of peak energy (GeV) contained in nOptimal crystals
TH2F phiCorrection
histogram of phi-dependent corrections
StoreArray< ECLShower > m_eclShowers
Store array: ECLShower.
virtual void initialize() override
Initialize.
virtual void event() override
Event.
const unsigned int nThetaID
69 thetaIDs
virtual const char * eclShowerArrayName() const
We need names for the data objects to differentiate between PureCsI and default.
virtual void endRun() override
End run.
int nEnergies
number of energies for which there are leakage corrections
TH2F m_bias
2D histogram of bias = sum of ECLCalDigit energy minus true (GeV)
virtual void terminate() override
Terminate.
TH2F thetaCorrection
histogram of theta-dependent corrections
std::vector< std::vector< float > > leakLogE
log(E) values for each region
int nXBins
number of thetaID x energy bins
TH2F m_peakFracEnergy
2D histogram of peak fractional contained energy
DBObjPtr< ECLnOptimal > m_eclNOptimal
nOptimal payload: bias from beam background, and correction from number of crystals
virtual void beginRun() override
Begin run.
DBObjPtr< ECLLeakageCorrections > m_eclLeakageCorrections
Leakage correction from database: leakage correction as a function of position.
ECL::ECLLeakagePosition * m_leakagePosition
location of cluster; cellID and position within the crystal
int nPositionBins
number of locations across crystal
StoreObjPtr< EventLevelClusteringInfo > m_eventLevelClusteringInfo
EventLevelClusteringInfo.
@ c_nPhotons
CR is split into n photons (N1)
Definition: ECLShower.h:42
Class to get position information for a cluster for leakage corrections.
std::vector< int > getLeakagePosition(const int cellIDFromEnergy, const float theta, const float phi, const int nPositions)
Return postion.
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
REG_MODULE(arichBtest)
Register the Module.
bool isBackward(int cellId)
Check whether the crystal is in backward ECL.
bool isForward(int cellId)
Check whether the crystal is in forward ECL.
bool isBarrel(int cellId)
Check whether the crystal is in barrel ECL.
Abstract base class for different kinds of events.