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