Belle II Software  release-06-02-00
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 //..This module`
10 #include <ecl/modules/eclLeakageCollector/eclLeakageCollectorModule.h>
11 
12 //..Framework
13 #include <framework/gearbox/Const.h>
14 #include <framework/dataobjects/EventMetaData.h>
15 
16 //..Root
17 #include <TH2F.h>
18 #include <TVector3.h>
19 #include <TMath.h>
20 #include <TTree.h>
21 
22 //..mdst
23 #include <mdst/dataobjects/MCParticle.h>
24 
25 //..ECL
26 #include <ecl/dbobjects/ECLCrystalCalib.h>
27 #include <ecl/dataobjects/ECLShower.h>
28 #include <ecl/geometry/ECLLeakagePosition.h>
29 
30 
31 //..Other
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");
53  setPropertyFlags(c_ParallelProcessingCertified);
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("ECLShowers"));
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.001);
83  i_energies[1][ie] = (int)(1000.*m_energies_barrel[ie] + 0.001);
84  i_energies[2][ie] = (int)(1000.*m_energies_backward[ie] + 0.001);
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;
140  leakagePosition = new ECLLeakagePosition();
141 
142  //-----------------------------------------------------------------
143  //..Required arrays
144  m_eclShowerArray.isRequired(m_showerArrayName);
145  m_mcParticleArray.isRequired();
146 }
147 
148 //-----------------------------------------------------------------
150 {
151 
152  //-----------------------------------------------------------------
153  //..First time, store the job parameters
154  if (storeCalib) {
155  getObjectPtr<TH1F>("inputParameters")->Fill(0.01, m_position_bins);
156  getObjectPtr<TH1F>("inputParameters")->Fill(1.01, m_number_energies);
157  double firstBin = 2.01;
158  for (int ie = 0; ie < m_number_energies; ie++) {
159  getObjectPtr<TH1F>("inputParameters")->Fill(firstBin + ie, m_energies_forward[ie]);
160  }
161  firstBin += m_number_energies;
162  for (int ie = 0; ie < m_number_energies; ie++) {
163  getObjectPtr<TH1F>("inputParameters")->Fill(firstBin + ie, m_energies_barrel[ie]);
164  }
165  firstBin += m_number_energies;
166  for (int ie = 0; ie < m_number_energies; ie++) {
167  getObjectPtr<TH1F>("inputParameters")->Fill(firstBin + ie, m_energies_backward[ie]);
168  }
169 
170  //..Keep track of how many times inputParameters was filled
171  int lastBin = getObjectPtr<TH1F>("inputParameters")->GetNbinsX();
172  getObjectPtr<TH1F>("inputParameters")->SetBinContent(lastBin, 1.);
173  storeCalib = false;
174  }
175 
176  //-----------------------------------------------------------------
177  //..Generated MC particle (always the first entry in the list)
178  double mcLabE = m_mcParticleArray[0]->getEnergy();
179  TVector3 mcp3 = m_mcParticleArray[0]->getMomentum();
180  TVector3 vertex = m_mcParticleArray[0]->getProductionVertex();
181 
182  //-----------------------------------------------------------------
183  //..Find the shower closest to the MC true angle
184  const int nShower = m_eclShowerArray.getEntries();
185  double minAngle = 999.;
186  int minShower = -1;
187  for (int is = 0; is < nShower; is++) {
188 
189  //..Only interested in photon hypothesis showers
190  if (m_eclShowerArray[is]->getHypothesisId() == ECLShower::c_nPhotons) {
191 
192  //..Make a position vector from theta, phi, and R
193  TVector3 position(0., 0., 1.);
194  position.SetTheta(m_eclShowerArray[is]->getTheta());
195  position.SetPhi(m_eclShowerArray[is]->getPhi());
196  position.SetMag(m_eclShowerArray[is]->getR());
197 
198  //..Direction is difference between position and vertex
199  TVector3 direction = position - vertex;
200 
201  double angle = direction.Angle(mcp3);
202  if (angle < minAngle) {
203  minAngle = angle;
204  minShower = is;
205  }
206  }
207  }
208  if (minShower == -1) {return;}
209 
210  //-----------------------------------------------------------------
211  //..Location of selected shower.
212  const int maxECellId = m_eclShowerArray[minShower]->getCentralCellId();
213  const float thetaLocation = m_eclShowerArray[minShower]->getTheta();
214  const float phiLocation = m_eclShowerArray[minShower]->getPhi();
215  std::vector<int> positionVector = leakagePosition->getLeakagePosition(maxECellId, thetaLocation, phiLocation, m_position_bins);
216 
217  //..TTree items
218  t_cellID = positionVector[0];
219  t_thetaID = positionVector[1];
220  t_region = positionVector[2];
221  t_thetaBin = positionVector[3];
222  t_phiBin = positionVector[4];
223  t_phiMech = positionVector[5];
224 
225  //-----------------------------------------------------------------
226  //..Generated and reconstructed energy quantities
227 
228  //..Find the generated energy bin
229  const int iGenEnergyMeV = (int)(1000.*mcLabE + 0.001);
230  t_energyBin = -1;
231  for (int ie = 0; ie < m_number_energies; ie++) {
232  if (iGenEnergyMeV == i_energies[t_region][ie]) {
233  t_energyBin = ie;
234  break;
235  }
236  }
237 
238  //..Crystals used to calculate energy
239  int nForEnergy = (int)(m_eclShowerArray[minShower]->getNumberOfCrystalsForEnergy() + 0.001);
240  t_nCrys = std::min(nCrysMax, nForEnergy);
241 
242  //..Reconstructed (without leakage corrections), normalized to generated
243  t_energyFrac = m_eclShowerArray[minShower]->getEnergyRaw() / mcLabE;
244 
245  //..Reconstructed energy after existing leakage correction, normalized to generated
246  t_origEnergyFrac = m_eclShowerArray[minShower]->getEnergy() / mcLabE;
247 
248  //-----------------------------------------------------------------
249  //..Distance between generated and reconstructed positions
250  const double radius = m_eclShowerArray[minShower]->getR();
251  TVector3 measuredLocation(0., 0., 1.);
252  measuredLocation.SetTheta(thetaLocation);
253  measuredLocation.SetPhi(phiLocation);
254  measuredLocation.SetMag(radius);
255  TVector3 measuredDirection = measuredLocation - vertex;
256  t_locationError = radius * measuredDirection.Angle(mcp3);
257 
258  //-----------------------------------------------------------------
259  //..Done
260  getObjectPtr<TTree>("tree")->Fill();
261 }
Calibration collector module base class.
@ c_nPhotons
CR is split into n photons (N1)
Definition: ECLShower.h:37
Class to get position information for a cluster for leakage corrections.
Store event, run, and experiment numbers.
Definition: EventMetaData.h:33
Store information needed to calculate ECL energy leakage corrections.
virtual void collect() override
Accumulate TTree.
virtual void prepare() override
Define histograms and read payloads from DB.
#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.