Belle II Software  release-08-01-10
eclNOptimalCollectorModule.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/eclNOptimalCollector/eclNOptimalCollectorModule.h>
11 
12 /* ECL headers. */
13 #include <ecl/dataobjects/ECLCalDigit.h>
14 #include <ecl/dataobjects/ECLShower.h>
15 #include <ecl/geometry/ECLNeighbours.h>
16 
17 /* Basf2 headers. */
18 #include <mdst/dataobjects/MCParticle.h>
19 #include <mdst/dataobjects/ECLCluster.h>
20 
21 /* ROOT headers. */
22 #include <TH2F.h>
23 
24 /* C++ headers. */
25 #include <iostream>
26 
27 using namespace std;
28 using namespace Belle2;
29 using namespace ECL;
30 
31 //-----------------------------------------------------------------
32 // Register the Module
33 //-----------------------------------------------------------------
34 REG_MODULE(eclNOptimalCollector);
35 
36 //-----------------------------------------------------------------
37 // Implementation
38 //-----------------------------------------------------------------
39 
42 eclNOptimalCollectorModule::eclNOptimalCollectorModule() : CalibrationCollectorModule()
43 {
44  // Set module properties
45  setDescription("Calibration collector module to find optimal number of crystal for cluster energies");
47  addParam("numberEnergies", m_numberEnergies, "number of generated energies", 8);
48  addParam("energiesForward", m_energiesForward, "generated photon energies, forward", std::vector<double> {0.030, 0.050, 0.100, 0.200, 0.483, 1.166, 2.816, 6.800});
49  addParam("energiesBarrel", m_energiesBarrel, "generated photon energies, barrel", std::vector<double> {0.030, 0.050, 0.100, 0.200, 0.458, 1.049, 2.402, 5.500});
50  addParam("energiesBackward", m_energiesBackward, "generated photon energies, backward", std::vector<double> {0.030, 0.050, 0.100, 0.200, 0.428, 0.917, 1.962, 4.200});
51  addParam("digitArrayName", m_digitArrayName, "name of ECLCalDigit data object", std::string("ECLTrimmedDigits"));
52  addParam("showerArrayName", m_showerArrayName, "name of ECLShower data object", std::string("ECLTrimmedShowers"));
53  addParam("nGroupPerThetaID", m_nGroupPerThetaID, "groups per standard thetaID", 8);
54 }
55 
59 {
60 
61  //--------------------------------------------------------
62  //..Check that input parameters are consistent
63  const int nForward = m_energiesForward.size();
64  const int nBarrel = m_energiesBarrel.size();
65  const int nBackward = m_energiesBackward.size();
66  if (nForward != m_numberEnergies or nBarrel != m_numberEnergies or nBackward != m_numberEnergies) {
67  B2FATAL("eclNOptimalCollector: length of energy vectors inconsistent with parameter numberEnergies: " << nForward << " " <<
68  nBarrel << " " << nBackward << " " << m_numberEnergies);
69  }
70 
71  //..Store generated energies as integers in MeV
72  iEnergies.resize(nLeakReg, std::vector<int>(m_numberEnergies, 0));
73  for (int ie = 0; ie < m_numberEnergies; ie++) {
74  iEnergies[0][ie] = (int)(1000.*m_energiesForward[ie] + 0.5);
75  iEnergies[1][ie] = (int)(1000.*m_energiesBarrel[ie] + 0.5);
76  iEnergies[2][ie] = (int)(1000.*m_energiesBackward[ie] + 0.5);
77  }
78 
79  //--------------------------------------------------------
80  //..Required data objects
82  m_eclClusterArray.isRequired("ECLClusters");
85 
86  //--------------------------------------------------------
87  //..Sort the crystals into groups of similar performance
88 
89  //..Record the thetaID of each cellID
90  neighbours = new ECLNeighbours("N", 1);
91  std::vector<int> nCrysPerRing;
92  for (int thID = 0; thID < 69; thID++) {
93  const int nCrys = neighbours->getCrystalsPerRing(thID);
94  nCrysPerRing.push_back(nCrys);
95  for (int phID = 0; phID < nCrys; phID++) {
96  thetaIDofCrysID.push_back(thID);
97  }
98  }
99 
100  //..First crystalID of each thetaID
101  std::vector<int> firstCrysIdPerRing;
102  firstCrysIdPerRing.push_back(0);
103  for (int thID = 1; thID < 69; thID++) {
104  const int iCrysID = firstCrysIdPerRing[thID - 1] + nCrysPerRing[thID - 1];
105  firstCrysIdPerRing.push_back(iCrysID);
106  }
107 
108  //..Range of thetaID from first and last cellID
109  const int firstThetaId = thetaIDofCrysID[iFirstCellId - 1];
110  const int lastThetaId = thetaIDofCrysID[iLastCellId - 1];
111 
112  //..Nominal groups per thetaID, but some thetaID have double that
113  const int specialThetaID[4] = {5, 11, 60, 65};
114  int iSp = 0;
115  int firstGroupOfThetaID = -m_nGroupPerThetaID;
116  for (int thID = firstThetaId; thID <= lastThetaId; thID++) {
117  firstGroupOfThetaID += m_nGroupPerThetaID;
118  const int nCrysPerGroup = nCrysPerRing[thID] / m_nGroupPerThetaID;
119  for (int phID = 0; phID < nCrysPerRing[thID]; phID++) {
120  const int iCrysID = firstCrysIdPerRing[thID] + phID;
121  const int iLocalGroup = phID / nCrysPerGroup;
122  iGroupOfCrystal[iCrysID] = firstGroupOfThetaID + iLocalGroup;
123  }
124 
125  //..Special thetaID with double the number of groups
126  if (thID == specialThetaID[iSp]) {
127  iSp++;
128  firstGroupOfThetaID += m_nGroupPerThetaID;
129  for (int phID = 1; phID < nCrysPerRing[thID]; phID += 3) {
130  const int iCrysID = firstCrysIdPerRing[thID] + phID;
131  const int iLocalGroup = phID / nCrysPerGroup;
132  iGroupOfCrystal[iCrysID] = firstGroupOfThetaID + iLocalGroup;
133  }
134 
135  }
136  }
137 
138  //..Assign crystals outside of the useful range to the first or last group
139  for (int ic = 1; ic < iFirstCellId; ic++) {
140  const int iCrysID = ic - 1;
142  }
143 
144  for (int ic = iLastCellId + 1; ic <= ECLElementNumbers::c_NCrystals; ic++) {
145  const int iCrysID = ic - 1;
147  }
148 
149  //..Number of groups
151 
152  //-----------------------------------------------------------------
153  //..Define histogram to store parameters
154  const int nBinX = 3 + nLeakReg * m_numberEnergies;
155  auto inputParameters = new TH1F("inputParameters", "eclNOptimalCollector job parameters", nBinX, 0, nBinX);
156  registerObject<TH1F>("inputParameters", inputParameters);
157 
158  //..Store group number for each crystal
159  auto groupNumberOfEachCellID = new TH1F("groupNumberOfEachCellID", "group number of each cellID;cellID",
161  registerObject<TH1F>("groupNumberOfEachCellID", groupNumberOfEachCellID);
162 
163  //-----------------------------------------------------------------
164  //..Histograms to check MC sample properties
165  auto entriesPerThetaIdEnergy = new TH2F("entriesPerThetaIdEnergy", "Entries per thetaID/energy point;thetaID;energy point", 69, 0,
167  registerObject<TH2F>("entriesPerThetaIdEnergy", entriesPerThetaIdEnergy);
168 
169  auto mcEnergyDiff = new TH2F("mcEnergyDiff", "mc E minus nominal (MeV) per thetaID/energy point;thetaID;energy point", 69, 0, 69,
171  mcEnergyDiff->Sumw2();
172  registerObject<TH2F>("mcEnergyDiff", mcEnergyDiff);
173 
174 
175  //-----------------------------------------------------------------
176  //..Histograms to store the sum of n crystals. One per group / energy.
177  // Include an extra nCrystal bin at the end to store the current raw energy
178  const int nHist = nCrystalGroups * m_numberEnergies;
179  vector<TH2F*> eSum(nHist);
180  vector<TH2F*> biasSum(nHist);
181  int iHist = -1;
182  for (int ig = 0; ig < nCrystalGroups; ig++) {
183  for (int ie = 0; ie < m_numberEnergies; ie++) {
184  iHist++;
185  std::string name = "eSum_" + std::to_string(ig) + "_" + std::to_string(ie);
186  TString hname = name;
187  TString title = "energy summed over nCrys divided by eMC, group " + std::to_string(ig) + ", E point " + std::to_string(
188  ie) + ";nCrys;energy sum / Etrue";
189  eSum[iHist] = new TH2F(hname, title, nCrysMax + 1, 1., nCrysMax + 2., 2000, 0., 2.);
190  registerObject<TH2F>(name, eSum[iHist]);
191 
192  name = "biasSum_" + std::to_string(ig) + "_" + std::to_string(ie);
193  hname = name;
194  title = "energy minus mc true summing over nCrys, group " + std::to_string(ig) + ", E point " + std::to_string(
195  ie) + ";nCrys;bias = energy minus mc truth (GeV)";
196  biasSum[iHist] = new TH2F(hname, title, nCrysMax + 1, 1., nCrysMax + 2., 1000, -0.1, 0.1);
197  registerObject<TH2F>(name, biasSum[iHist]);
198  }
199  }
200 }
201 
205 {
206 
207  //-----------------------------------------------------------------
208  //..First time, store the job parameters
209  if (storeParameters) {
210  getObjectPtr<TH1F>("inputParameters")->Fill(0.01, nCrystalGroups);
211  getObjectPtr<TH1F>("inputParameters")->Fill(1.01, m_numberEnergies);
212  double firstBin = 2.01;
213  for (int ie = 0; ie < m_numberEnergies; ie++) {
214  getObjectPtr<TH1F>("inputParameters")->Fill(firstBin + ie, m_energiesForward[ie]);
215  }
216  firstBin += m_numberEnergies;
217  for (int ie = 0; ie < m_numberEnergies; ie++) {
218  getObjectPtr<TH1F>("inputParameters")->Fill(firstBin + ie, m_energiesBarrel[ie]);
219  }
220  firstBin += m_numberEnergies;
221  for (int ie = 0; ie < m_numberEnergies; ie++) {
222  getObjectPtr<TH1F>("inputParameters")->Fill(firstBin + ie, m_energiesBackward[ie]);
223  }
224 
225  //..Keep track of how many times inputParameters was filled
226  int lastBin = getObjectPtr<TH1F>("inputParameters")->GetNbinsX();
227  getObjectPtr<TH1F>("inputParameters")->SetBinContent(lastBin, 1.);
228 
229  //..Store the group number for every cellID
230  for (int ic = 1; ic < 8737; ic++) {
231  getObjectPtr<TH1F>("groupNumberOfEachCellID")->Fill(ic + 0.01, iGroupOfCrystal[ic - 1]);
232  }
233 
234  //..Call each eSum hist once to ensure they exist (with non-real nCrys)
235  for (int ig = 0; ig < nCrystalGroups; ig++) {
236  for (int ie = 0; ie < m_numberEnergies; ie++) {
237  std::string histName = "eSum_" + std::to_string(ig) + "_" + std::to_string(ie);
238  getObjectPtr<TH2F>(histName)->Fill(0.01, 0.96);
239 
240  std::string histNameBias = "biasSum_" + std::to_string(ig) + "_" + std::to_string(ie);
241  getObjectPtr<TH2F>(histNameBias)->Fill(0.01, 0.96);
242  }
243  }
244 
245  storeParameters = false;
246  }
247 
248  //--------------------------------------------------------
249  //..Find the ECLShower (should only be one when using trimmed data object)
250  int iMax = -1;
251  double showerMaxE = 0.; // shower energy before leakage correction
252  const int nShower = m_eclShowerArray.getEntries();
253  for (int i = 0; i < nShower; i++) {
254  if (m_eclShowerArray[i]->getHypothesisId() == ECLShower::c_nPhotons) {
255  const double nominalE = m_eclShowerArray[i]->getEnergyRaw();
256  if (nominalE > showerMaxE) {
257  showerMaxE = nominalE;
258  iMax = i;
259  }
260  }
261  }
262 
263  //..Quit now if no shower
264  if (iMax == -1) {return;}
265 
266  //--------------------------------------------------------
267  //..Get cellID from related cluster
268  int iCellId = -1;
269  const auto showerClusterRelations = m_eclShowerArray[iMax]->getRelationsWith<ECLCluster>();
270  const unsigned int nRelatedClusters = showerClusterRelations.size();
271  if (nRelatedClusters > 0) {
272  const auto cluster = showerClusterRelations.object(0);
273  iCellId = cluster->getMaxECellId();
274  }
275 
276  //..Quit if cellID is not in the range being calibrated
277  if (iCellId<iFirstCellId or iCellId>iLastCellId) {return;}
278 
279  //..ECL region forward/barrel/backward
280  int iRegion = 1;
281  if (ECLElementNumbers::isForward(iCellId)) {
282  iRegion = 0;
283  } else if (ECLElementNumbers::isBackward(iCellId)) {
284  iRegion = 2;
285  }
286 
287  //--------------------------------------------------------
288  //..Get the generated energy from the MCParticles block
289 
290  //..Should only be one entry. Quit if this is not the case.
291  const int nMC = m_mcParticleArray.getEntries();
292  if (nMC != 1) {return;}
293 
294  //..Energy. Convert to MeV to get corresponding bin
295  const double eTrue = m_mcParticleArray[0]->getEnergy();
296  const float genEnergyMeV = 1000.*eTrue;
297  const float tolerance = std::max(0.002 * genEnergyMeV, 1.0);
298 
299  //..Look for generated energy to be close enough to the expected one
300  int iEnergy = -1;
301  double minDiff = 9999.;
302  for (int ie = 0; ie < m_numberEnergies; ie++) {
303  double diff = std::abs(genEnergyMeV - iEnergies[iRegion][ie]);
304  if (diff < std::abs(minDiff)) {
305  minDiff = genEnergyMeV - iEnergies[iRegion][ie];
306  }
307  if (diff < tolerance) {
308  iEnergy = ie;
309  break;
310  }
311  }
312 
313  //..Quit if the true energy is not equal to a generated one.
314  // This happens if the cluster is reconstructed in the wrong region.
315  if (iEnergy == -1) {return;}
316 
317  //..Some statistics for validation
318  getObjectPtr<TH2F>("entriesPerThetaIdEnergy")->Fill(thetaIDofCrysID[iCellId - 1] + 0.001, iEnergy + 0.001);
319  getObjectPtr<TH2F>("mcEnergyDiff")->Fill(thetaIDofCrysID[iCellId - 1] + 0.001, iEnergy + 0.001, minDiff);
320 
321  //--------------------------------------------------------
322  //..Get the ECLCalDigits and weights associated with the cluster,
323  // plus MC true energy
324  std::vector<double> digitEnergies;
325  std::vector < std::pair<double, double> > energies;
326 
327  const auto showerDigitRelations = m_eclShowerArray[iMax]->getRelationsWith<ECLCalDigit>(m_digitArrayName);
328  unsigned int nRelatedDigits = showerDigitRelations.size();
329  for (unsigned int ir = 0; ir < nRelatedDigits; ir++) {
330  const auto calDigit = showerDigitRelations.object(ir);
331  const auto weight = showerDigitRelations.weight(ir);
332  digitEnergies.push_back(calDigit->getEnergy() * weight);
333  const double eCalDigit = weight * calDigit->getEnergy();
334 
335  //..MC energy via relation to ECLCalDigit
336  auto digitMCRelations = calDigit->getRelationsTo<MCParticle>();
337  double eMC = 0.;
338  for (unsigned int i = 0; i < digitMCRelations.size(); i++) {
339  eMC += digitMCRelations.weight(i);
340  }
341  std::pair<double, double> pTemp = std::make_pair(eCalDigit, eMC);
342  energies.push_back(pTemp);
343  }
344 
345  //..Sort digit and mc energies from lowest to highest
346  std::sort(energies.begin(), energies.end());
347 
348  //--------------------------------------------------------
349  //..Store the energy sum of 1, 2, ..., n crystals (max 21), and corresponding bias
350  std::string histName = "eSum_" + std::to_string(iGroupOfCrystal[iCellId - 1]) + "_" + std::to_string(iEnergy);
351  std::string histNameBias = "biasSum_" + std::to_string(iGroupOfCrystal[iCellId - 1]) + "_" + std::to_string(iEnergy);
352 
353  double eSumOfN = 0.;
354  double biasSumOfN = 0.;
355  for (int isum = 0; isum < nCrysMax; isum++) {
356  int i = (int)nRelatedDigits - 1 - isum;
357  if (i >= 0) {
358  eSumOfN += energies[i].first;
359  biasSumOfN += (energies[i].first - energies[i].second);
360  }
361  getObjectPtr<TH2F>(histName)->Fill(isum + 1.01, eSumOfN / eTrue);
362  getObjectPtr<TH2F>(histNameBias)->Fill(isum + 1.01, biasSumOfN);
363  }
364 
365  //..Also store the current raw energy for monitoring purposes
366  getObjectPtr<TH2F>(histName)->Fill(nCrysMax + 1.01, showerMaxE / eTrue);
367 
368 }
Calibration collector module base class.
Class to store calibrated ECLDigits: ECLCalDigits.
Definition: ECLCalDigit.h:23
ECL cluster data.
Definition: ECLCluster.h:27
@ c_nPhotons
CR is split into n photons (N1)
Definition: ECLShower.h:42
Class to get the neighbours for a given cell id.
Definition: ECLNeighbours.h:25
short int getCrystalsPerRing(const short int thetaid) const
return number of crystals in a given theta ring
Definition: ECLNeighbours.h:39
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
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
Array of MCParticles.
StoreArray< ECLShower > m_eclShowerArray
Required arrays.
const int iFirstCellId
first useful cellID (first of thetaID 3)
int m_nGroupPerThetaID
number of groups per standard thetaID
std::vector< std::vector< int > > iEnergies
Some other useful quantities.
std::vector< double > m_energiesBarrel
generated photon energies, barrel
const int iLastCellId
first useful cellID (last of thetaID 66)
const int nCrysMax
max number of crystals used to calculate energy
int nCrystalGroups
sort the crystals into this many groups
const int nLeakReg
0 = forward, 1 = barrel, 2 = backward
StoreArray< ECLCluster > m_eclClusterArray
Array of ECLClusters.
int iGroupOfCrystal[ECLElementNumbers::c_NCrystals]
group number of each crystal
ECL::ECLNeighbours * neighbours
neighbours to crystal
std::vector< int > thetaIDofCrysID
thetaID of each crystal
void collect() override
Select events and crystals and accumulate histograms.
std::string m_showerArrayName
Name of ECLShower StoreArray.
std::string m_digitArrayName
Name of ECLCalDigit StoreArray.
StoreArray< ECLCalDigit > m_eclCalDigitArray
Array of ECLCalDigits.
void prepare() override
Define histograms.
std::vector< double > m_energiesForward
generated photon energies, forward
bool storeParameters
store parameters first event
std::vector< double > m_energiesBackward
generated photon energies, backward
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
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
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.