Belle II Software  release-08-00-10
eclLeakageCollectorModule.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/eclLeakageCollector/eclLeakageCollectorModule.h>
11 
12 /* ECL headers. */
13 #include <ecl/dbobjects/ECLCrystalCalib.h>
14 #include <ecl/dataobjects/ECLShower.h>
15 #include <ecl/dataobjects/ECLCalDigit.h>
16 #include <ecl/geometry/ECLLeakagePosition.h>
17 #include <ecl/dataobjects/ECLElementNumbers.h>
18 
19 /* Basf2 headers. */
20 #include <framework/gearbox/Const.h>
21 #include <framework/dataobjects/EventMetaData.h>
22 #include <framework/geometry/VectorUtil.h>
23 #include <mdst/dataobjects/MCParticle.h>
24 
25 /* Root headers. */
26 #include <Math/Vector3D.h>
27 #include <Math/VectorUtil.h>
28 #include <TMath.h>
29 #include <TTree.h>
30 
31 /* C++ headers. */
32 #include <iostream>
33 
34 using namespace Belle2;
35 using namespace ECL;
36 
37 //-----------------------------------------------------------------
38 // Register the Module
39 //-----------------------------------------------------------------
40 REG_MODULE(eclLeakageCollector);
41 
42 //-----------------------------------------------------------------
43 // Implementation
44 //-----------------------------------------------------------------
45 
46 //-----------------------------------------------------------------
48  m_mcParticleArray("MCParticles"),
49  m_evtMetaData("EventMetaData")
50 {
52  setDescription("Store quantities from single photon MC used to calculated ECL energy leakage corrections");
54  addParam("position_bins", m_position_bins, "number of crystal subdivisions in theta and phi", 29);
55  addParam("number_energies", m_number_energies, "number of generated energies", 8);
56  addParam("energies_forward", m_energies_forward, "generated photon energies, forward", std::vector<double> {0.030, 0.050, 0.100, 0.200, 0.483, 1.166, 2.816, 6.800});
57  addParam("energies_barrel", m_energies_barrel, "generated photon energies, barrel", std::vector<double> {0.030, 0.050, 0.100, 0.200, 0.458, 1.049, 2.402, 5.500});
58  addParam("energies_backward", m_energies_backward, "generated photon energies, backward", std::vector<double> {0.030, 0.050, 0.100, 0.200, 0.428, 0.917, 1.962, 4.200});
59  addParam("showerArrayName", m_showerArrayName, "name of ECLShower data object", std::string("ECLTrimmedShowers"));
60 }
61 
62 
63 //-----------------------------------------------------------------
65 {
66  //-----------------------------------------------------------------
67  //..Parameters and other basic info
68  B2INFO("eclLeakageCollector: Experiment = " << m_evtMetaData->getExperiment() << " run = " << m_evtMetaData->getRun());
69 
70  //..Check that input parameters are consistent
71  const int n_e_forward = m_energies_forward.size();
72  const int n_e_barrel = m_energies_barrel.size();
73  const int n_e_backward = m_energies_backward.size();
74  if (n_e_forward != m_number_energies or n_e_barrel != m_number_energies or n_e_backward != m_number_energies) {
75  B2FATAL("eclLeakageCollector: length of energy vectors inconsistent with parameter number_energies: " << n_e_forward << " " <<
76  n_e_barrel << " " << n_e_backward << " " << m_number_energies);
77  }
78 
79  //..Store generated energies as integers in MeV
80  i_energies.resize(nLeakReg, std::vector<int>(m_number_energies, 0));
81  for (int ie = 0; ie < m_number_energies; ie++) {
82  i_energies[0][ie] = (int)(1000.*m_energies_forward[ie] + 0.5);
83  i_energies[1][ie] = (int)(1000.*m_energies_barrel[ie] + 0.5);
84  i_energies[2][ie] = (int)(1000.*m_energies_backward[ie] + 0.5);
85  }
86 
87  //..Require all energies are different, and that there are at least two
88  if (m_number_energies < 2) {
89  B2FATAL("eclLeakageCollector: require at least two energy points. m_number_energies = " << m_number_energies);
90  }
91  for (int ie = 0; ie < m_number_energies - 1; ie++) {
92  for (int ireg = 0; ireg < nLeakReg; ireg++) {
93  if (i_energies[ireg][ie] == i_energies[ireg][ie + 1]) {
94  B2FATAL("eclLeakageCollector: identical energies, ireg = " << ireg << ", " << i_energies[ireg][ie] << " MeV");
95  }
96  }
97  }
98 
99  //-----------------------------------------------------------------
100  //..Write out the input parameters
101  B2INFO("eclLeakageCollector parameters: ");
102  B2INFO("position_bins " << m_position_bins);
103  B2INFO("number_energies " << m_number_energies);
104  std::cout << "energies_forward ";
105  for (int ie = 0; ie < m_number_energies; ie++) {std::cout << m_energies_forward[ie] << " ";}
106  std::cout << std::endl;
107  std::cout << "energies_barrel ";
108  for (int ie = 0; ie < m_number_energies; ie++) {std::cout << m_energies_barrel[ie] << " ";}
109  std::cout << std::endl;
110  std::cout << "energies_backward ";
111  for (int ie = 0; ie < m_number_energies; ie++) {std::cout << m_energies_backward[ie] << " ";}
112  std::cout << std::endl;
113  B2INFO("showerArrayName " << m_showerArrayName);
114 
115  //-----------------------------------------------------------------
116  //..Define histogram to store parameters
117  const int nBinX = 3 + nLeakReg * m_number_energies;
118  auto inputParameters = new TH1F("inputParameters", "eclLeakageCollector job parameters", nBinX, 0, nBinX);
119  registerObject<TH1F>("inputParameters", inputParameters);
120 
121  //..TTree stores required quantities for each photon
122  auto tree = new TTree("single_photon_leakage", "");
123  tree->Branch("cellID", &t_cellID, "cellID/I");
124  tree->Branch("thetaID", &t_thetaID, "thetaID/I");
125  tree->Branch("region", &t_region, "region/I");
126  tree->Branch("thetaBin", &t_thetaBin, "thetaBin/I");
127  tree->Branch("phiBin", &t_phiBin, "phiBin/I");
128  tree->Branch("phiMech", &t_phiMech, "phiMech/I");
129  tree->Branch("energyBin", &t_energyBin, "energyBin/I");
130  tree->Branch("nCrys", &t_nCrys, "nCrys/I");
131  tree->Branch("energyFrac", &t_energyFrac, "energyFrac/F");
132  tree->Branch("origEnergyFrac", &t_origEnergyFrac, "origEnergyFrac/F");
133  tree->Branch("locationError", &t_locationError, "locationError/F");
134  registerObject<TTree>("tree", tree);
135 
136 
137  //-----------------------------------------------------------------
138  //..Class to find cellID and position within crystal from theta and phi
139  std::cout << "creating leakagePosition object " << std::endl;
141 
142  //-----------------------------------------------------------------
143  //..Required arrays
146 
147  //-----------------------------------------------------------------
148  //..nOptimal payload. Get optimal number of crystals, and
149  // corresponding correction factors.
150  m_nOptimal2D = m_eclNOptimal->getNOptimal();
151  m_peakFracEnergy = m_eclNOptimal->getPeakFracEnergy();
152  m_bias = m_eclNOptimal->getBias();
153  m_logPeakEnergy = m_eclNOptimal->getLogPeakEnergy();
154  m_groupNumber = m_eclNOptimal->getGroupNumber();
155 
156  //..Vectors of energy boundaries for each region
157  std::vector<float> eBoundariesFwd = m_eclNOptimal->getUpperBoundariesFwd();
158  std::vector<float> eBoundariesBrl = m_eclNOptimal->getUpperBoundariesBrl();
159  std::vector<float> eBoundariesBwd = m_eclNOptimal->getUpperBoundariesBwd();
160  m_nEnergyBins = eBoundariesBrl.size();
161 
162  //..Copy values to m_eBoundaries
163  m_eBoundaries.resize(m_nLeakReg, std::vector<float>(m_nEnergyBins, 0.));
164  for (int ie = 0; ie < m_nEnergyBins; ie++) {
165  m_eBoundaries[0][ie] = eBoundariesFwd[ie];
166  m_eBoundaries[1][ie] = eBoundariesBrl[ie];
167  m_eBoundaries[2][ie] = eBoundariesBwd[ie];
168  B2INFO(" nOptimal upper boundaries for energy point " << ie << " " << m_eBoundaries[0][ie] << " " << m_eBoundaries[1][ie] << " " <<
169  m_eBoundaries[2][ie]);
170  }
171 
172 }
173 
174 //-----------------------------------------------------------------
176 {
177 
178  //-----------------------------------------------------------------
179  //..First time, store the job parameters
180  if (storeCalib) {
181  getObjectPtr<TH1F>("inputParameters")->Fill(0.01, m_position_bins);
182  getObjectPtr<TH1F>("inputParameters")->Fill(1.01, m_number_energies);
183  double firstBin = 2.01;
184  for (int ie = 0; ie < m_number_energies; ie++) {
185  getObjectPtr<TH1F>("inputParameters")->Fill(firstBin + ie, m_energies_forward[ie]);
186  }
187  firstBin += m_number_energies;
188  for (int ie = 0; ie < m_number_energies; ie++) {
189  getObjectPtr<TH1F>("inputParameters")->Fill(firstBin + ie, m_energies_barrel[ie]);
190  }
191  firstBin += m_number_energies;
192  for (int ie = 0; ie < m_number_energies; ie++) {
193  getObjectPtr<TH1F>("inputParameters")->Fill(firstBin + ie, m_energies_backward[ie]);
194  }
195 
196  //..Keep track of how many times inputParameters was filled
197  int lastBin = getObjectPtr<TH1F>("inputParameters")->GetNbinsX();
198  getObjectPtr<TH1F>("inputParameters")->SetBinContent(lastBin, 1.);
199  storeCalib = false;
200  }
201 
202  //-----------------------------------------------------------------
203  //..Generated MC particle. Should only be one entry, but check.
204  const int nMC = m_mcParticleArray.getEntries();
205  if (nMC != 1) {return;}
206  double mcLabE = m_mcParticleArray[0]->getEnergy();
207  ROOT::Math::XYZVector mcp3 = m_mcParticleArray[0]->getMomentum();
208  ROOT::Math::XYZVector vertex = m_mcParticleArray[0]->getProductionVertex();
209 
210  //-----------------------------------------------------------------
211  //..Find the ECLShower (should only be one when using trimmed data object)
212  const int nShower = m_eclShowerArray.getEntries();
213  double minAngle = 999.;
214  int minShower = -1;
215  for (int is = 0; is < nShower; is++) {
216 
217  //..Only interested in photon hypothesis showers
218  if (m_eclShowerArray[is]->getHypothesisId() == ECLShower::c_nPhotons) {
219 
220  //..Make a position vector from theta, phi, and R
221  ROOT::Math::XYZVector position;
222  VectorUtil::setMagThetaPhi(
223  position, m_eclShowerArray[is]->getR(),
224  m_eclShowerArray[is]->getTheta(), m_eclShowerArray[is]->getPhi());
225 
226  //..Direction is difference between position and vertex
227  ROOT::Math::XYZVector direction = ROOT::Math::XYZVector(position) - vertex;
228 
229  double angle = ROOT::Math::VectorUtil::Angle(direction, mcp3);
230  if (angle < minAngle) {
231  minAngle = angle;
232  minShower = is;
233  }
234  }
235  }
236  if (minShower == -1) {return;}
237 
238  //-----------------------------------------------------------------
239  //..Location of selected shower.
240  const int maxECellId = m_eclShowerArray[minShower]->getCentralCellId();
241  const float thetaLocation = m_eclShowerArray[minShower]->getTheta();
242  const float phiLocation = m_eclShowerArray[minShower]->getPhi();
243  std::vector<int> positionVector = leakagePosition->getLeakagePosition(maxECellId, thetaLocation, phiLocation, m_position_bins);
244 
245  //..TTree items
246  t_cellID = positionVector[0];
247  t_thetaID = positionVector[1];
248  t_region = positionVector[2];
249  t_thetaBin = positionVector[3];
250  t_phiBin = positionVector[4];
251  t_phiMech = positionVector[5];
252 
253  //-----------------------------------------------------------------
254  //..Generated and reconstructed energy quantities
255 
256  //..Find the generated energy bin by requiring it be close enough to expected value
257  const float genEnergyMeV = 1000.*mcLabE;
258  const float tolerance = std::max(0.002 * genEnergyMeV, 1.0);
259  t_energyBin = -1;
260  for (int ie = 0; ie < m_number_energies; ie++) {
261  if (std::abs(genEnergyMeV - i_energies[t_region][ie]) < tolerance) {
262  t_energyBin = ie;
263  break;
264  }
265  }
266 
267  //..Quit if the true energy is not equal to a generated one.
268  // This can happen if the cluster is reconstructed in the wrong region.
269  if (t_energyBin == -1) {return;}
270 
271  //..Reconstructed energy after existing leakage correction, normalized to generated
272  t_origEnergyFrac = m_eclShowerArray[minShower]->getEnergy() / mcLabE;
273 
274  //..Sum of nOptimal crystals (without leakage correction), when events were generated
275  const double eRaw = m_eclShowerArray[minShower]->getEnergyRaw();
276 
277  //-----------------------------------------------------------------
278  //..Find nOptimal from cellID and 3x3 energy found by ECLSplitterN1Module.
279  // Note that payload may have been updated since events were generated, so
280  // we need to redo the sum of energies of the nOptimal crystals.
281  const int ig = m_groupNumber[maxECellId - 1];
282  const double e3x3Orig = m_eclShowerArray[minShower]->getNOptimalEnergy();
283  double e3x3 = e3x3Orig;
284 
285  //..Alternate e3x3, for samples produced before nOptimal was introduced
286  const double eHighestCorr = m_eclShowerArray[minShower]->getEnergyHighestCrystal();
287  const double correction = m_eclShowerArray[minShower]->getEnergy() / m_eclShowerArray[minShower]->getEnergyRaw();
288  const double eHighestRaw = eHighestCorr / correction;
289  const double e3x3Alt = eHighestRaw / m_eclShowerArray[minShower]->getE1oE9();
290  if (e3x3 < 0.001) {e3x3 = e3x3Alt;}
291 
292  //..Need the detector region to get the energy bin boundaries
293  int iRegion = 1; // barrel
294  if (ECLElementNumbers::isForward(maxECellId)) {iRegion = 0;}
295  if (ECLElementNumbers::isBackward(maxECellId)) {iRegion = 2;}
296 
297  //..nOptimal energy bin for this energy.
298  int iEnergy = 0;
299  while (e3x3 > m_eBoundaries[iRegion][iEnergy] and iEnergy < m_nEnergyBins - 1) {iEnergy++;}
300 
301  //..nOptimal
302  t_nCrys = m_nOptimal2D.GetBinContent(ig + 1, iEnergy + 1);
303 
304  //-----------------------------------------------------------------
305  //..Find the sum of the nOptimal most-energetic digits
306 
307  //..Get the ECLCalDigits and rank them by energy
308  std::vector<double> digitEnergies;
309  const auto showerDigitRelations = m_eclShowerArray[minShower]->getRelationsWith<ECLCalDigit>("ECLTrimmedDigits");
310  unsigned int nRelatedDigits = showerDigitRelations.size();
311  for (unsigned int ir = 0; ir < nRelatedDigits; ir++) {
312  const auto calDigit = showerDigitRelations.object(ir);
313  const auto weight = showerDigitRelations.weight(ir);
314  digitEnergies.push_back(calDigit->getEnergy() * weight);
315  }
316 
317  //..Rank from lowest to highest
318  std::sort(digitEnergies.begin(), digitEnergies.end());
319 
320  //..Sum the highest nOptimal digit energies (if there are that many)
321  double eSumOfN = 1.e-10; // cpp does not like 0
322  for (int isum = 0; isum < t_nCrys; isum++) {
323  const int i = (int)nRelatedDigits - 1 - isum;
324  if (i >= 0) {eSumOfN += digitEnergies[i];}
325  }
326 
327  //-----------------------------------------------------------------
328  //..Correct the sum of nOptimal crystals for bias and nCrystal peak energy
329  // We need to do this because the nOptimal payload may have changed
330  // since the events were generated.
331 
332  //..To allow for energy interpolation, there are three bins per group and energy:
333  // iy = iEnergy + 1 in all three cases
334  // ix = 3*ig + 1 = value for nominal energy, i.e logPeakEnergy(3*ig+1, iEnergy+1)
335  // ix = 3*ig + 2 = value for lower energy, i.e logPeakEnergy(3*ig+2, iEnergy+1)
336  // ix = 3*ig + 3 = value for higher energy, i.e logPeakEnergy(3*ig+3, iEnergy+1)
337  const int iy = iEnergy + 1;
338 
339  const int ixNom = 3 * ig + 1;
340  const int ixLowerE = ixNom + 1;
341  const int ixHigherE = ixNom + 2;
342 
343  const double logENom = m_logPeakEnergy.GetBinContent(ixNom, iy);
344  const double logELower = m_logPeakEnergy.GetBinContent(ixLowerE, iy);
345  const double logEHigher = m_logPeakEnergy.GetBinContent(ixHigherE, iy);
346 
347  const double biasNom = m_bias.GetBinContent(ixNom, iy);
348  const double biasLower = m_bias.GetBinContent(ixLowerE, iy);
349  const double biasHigher = m_bias.GetBinContent(ixHigherE, iy);
350 
351  const double peakNom = m_peakFracEnergy.GetBinContent(ixNom, iy);
352  const double peakLower = m_peakFracEnergy.GetBinContent(ixLowerE, iy);
353  const double peakHigher = m_peakFracEnergy.GetBinContent(ixHigherE, iy);
354 
355  //..Interpolate in log energy
356  const double logESumN = log(eSumOfN);
357 
358  double logEOther = logELower;
359  double biasOther = biasLower;
360  double peakOther = peakLower;
361  if (logESumN > logENom) {
362  logEOther = logEHigher;
363  biasOther = biasHigher;
364  peakOther = peakHigher;
365  }
366 
367  //..The nominal and "other" energies may be identical if this is the first or last energy
368  double bias = biasNom;
369  double peak = peakNom;
370  if (abs(logEOther - logENom) > 0.0001) {
371  bias = biasNom + (biasOther - biasNom) * (logESumN - logENom) / (logEOther - logENom);
372  peak = peakNom + (peakOther - peakNom) * (logESumN - logENom) / (logEOther - logENom);
373  }
374 
375  //..Normalized reconstructed energy after bias and nCrystal correction
376  t_energyFrac = (eSumOfN - bias) / peak / mcLabE;
377 
378  //-----------------------------------------------------------------
379  //..Distance between generated and reconstructed positions
380  const double radius = m_eclShowerArray[minShower]->getR();
381  ROOT::Math::XYZVector measuredLocation;
382  VectorUtil::setMagThetaPhi(
383  measuredLocation, radius, thetaLocation, phiLocation);
384  ROOT::Math::XYZVector measuredDirection = ROOT::Math::XYZVector(measuredLocation) - vertex;
385  t_locationError = radius * ROOT::Math::VectorUtil::Angle(measuredDirection, mcp3);
386 
387 
388  //-----------------------------------------------------------------
389  //..Debug: dump some events
390  if (m_nDump < 100) {
391  m_nDump++;
392  B2INFO(" ");
393  B2INFO("cellID " << t_cellID << " ig " << ig << " iEnergy " << iEnergy << " ie " << t_energyBin << " nOpt " << t_nCrys);
394  B2INFO(" e3x3Orig " << e3x3Orig << " e3x3Alt " << e3x3Alt << " Eraw " << eRaw << " ESum " << eSumOfN << " eFrac " <<
395  t_energyFrac);
396  B2INFO(" 3 log E " << logELower << " " << logENom << " " << logEHigher << " log E " << logESumN);
397  B2INFO(" 3 biases " << biasLower << " " << biasNom << " " << biasHigher << " bias " << bias);
398  B2INFO(" 3 peaks " << peakLower << " " << peakNom << " " << peakHigher << " peak " << peak);
399  }
400 
401  //-----------------------------------------------------------------
402  //..Done
403  getObjectPtr<TTree>("tree")->Fill();
404 }
405 
Calibration collector module base class.
Class to store calibrated ECLDigits: ECLCalDigits.
Definition: ECLCalDigit.h:23
@ 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.
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
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
StoreArray< MCParticle > m_mcParticleArray
Required input array of MCParticles.
eclLeakageCollectorModule()
Constructor: Sets the description, the properties and the parameters of the module.
StoreArray< ECLShower > m_eclShowerArray
Required input array of ECLShowers.
TH2F m_logPeakEnergy
log of peak contained energy in GeV
float t_energyFrac
measured energy without leakage correction divided by generated
const int m_nLeakReg
3 ECL regions: 0 = forward, 1 = barrel, 2 = backward
float t_locationError
reconstructed minus generated position (cm)
std::vector< double > m_energies_forward
generated photon energies, forward
const int nLeakReg
Some other useful quantities.
int m_number_energies
number of generated energies (8)
std::vector< std::vector< int > > i_energies
Generated energies in MeV in each region.
bool storeCalib
store parameters first event
TH2F m_bias
2D histogram of bias = sum of ECLCalDigit energy minus true (GeV)
TH2F m_peakFracEnergy
2D histogram of peak fractional contained energy
virtual void collect() override
Accumulate TTree.
std::string m_showerArrayName
Required arrays.
StoreObjPtr< EventMetaData > m_evtMetaData
dataStore EventMetaData
DBObjPtr< ECLnOptimal > m_eclNOptimal
nOptimal payload
std::vector< double > m_energies_backward
generated photon energies, backward
ECL::ECLLeakagePosition * leakagePosition
location of position of cluster
virtual void prepare() override
Define histograms and read payloads from DB.
int t_phiBin
binned location in phi relative to crystal edge.
int t_phiMech
0: mechanical structure next to phi edge; 1: no mech structure
int t_region
region of photon 0=forward 1=barrel 2=backward
int m_nEnergyBins
number of energies bins in nOptimal payload
int t_thetaBin
binned location in theta relative to crystal edge
int t_nCrys
number of crystals used to calculate energy
std::vector< double > m_energies_barrel
generated photon energies, barrel
int m_position_bins
Parameters to pass along to the algorithm.
float t_origEnergyFrac
original leakage-corrected energy / generated
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
int m_nDump
Number of events with diagnostic printouts.
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
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.