Belle II Software  release-05-02-19
ECLBackgroundModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2010 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Sam de Jong *
7  * Alexandre Beaulieu *
8  * Poyuan Chen *
9  * *
10  * This software is provided "as is" without any warranty. *
11  **************************************************************************/
12 //This module
13 #include <ecl/modules/eclBackgroundStudy/ECLBackgroundModule.h>
14 
15 //Root
16 #include <TVector3.h>
17 #include <TH1F.h>
18 #include <TH2F.h>
19 
20 //Framework
21 #include <framework/logging/Logger.h>
22 
23 //ECL
24 #include <ecl/modules/eclBackgroundStudy/ECLCrystalData.h>
25 #include <ecl/dataobjects/ECLShower.h>
26 #include <ecl/dataobjects/ECLSimHit.h>
27 
28 #ifdef DOARICH
29 #include <arich/geometry/ARICHGeometryPar.h>
30 #endif
31 
32 //Simulation
33 #include <simulation/dataobjects/BeamBackHit.h>
34 
35 //MDST
36 #include <mdst/dataobjects/MCParticle.h>
37 
38 #define PI 3.14159265358979323846
39 
40 using namespace std;
41 using namespace Belle2;
42 
43 //-----------------------------------------------------------------
44 // Register the Module
45 //-----------------------------------------------------------------
46 REG_MODULE(ECLBackground)
47 
48 //-----------------------------------------------------------------
49 // Implementation
50 //-----------------------------------------------------------------
51 
52 
54 {
55  //Set module properties
56  setDescription("Processes background campaigns and produces histograms. Requires HistoManager");
57 
58  std::vector<int> empty;
59  addParam("sampleTime", m_sampleTime, "Length of sample, in us", 1000);
60  addParam("doARICH", m_doARICH, "If true, some ARICH plots (for shielding studies) will be produced", false);
61  addParam("crystalsOfInterest", m_CryInt, "Cell ID of crystals of interest. Dose will be printed at end of run", empty);
62 
63 }
64 
65 ECLBackgroundModule::~ECLBackgroundModule()
66 {
67 }
68 
69 //If a histogram is initalized here, it will be saved.
70 void ECLBackgroundModule::defineHisto()
71 {
72  std::ostringstream s;
73  s << m_sampleTime;
74 
75  //initalize histograms
76  h_nECLSimHits = new TH1F("ECL_Sim_Hits", "ECL Sim Hits", 100, 0, 100);
77 
78 
79  //Radiation dose
80  h_CrystalRadDose = new TH1F("Crystal_Rad_Dose", "Crystal Radiation Dose vs #theta_{ID};#theta_{ID};Gy/yr", 69, -0.5, 68.5);
81  h_CrystalRadDoseTheta = new TH1F("Crystal_Rad_Dose_Theta", "Crystal Radiation Dose vs #theta;#theta (deg);Gy/yr", 100, 12, 152);
82  h_CrystalThetaID2 = new TH1F("Crystal_Dose_ThetaID_2", "Crystal Radiation Dose vs #phi_{ID}, #theta_{ID}=2; #phi_{ID};Gy/yr", 64,
83  -0.5, 63.5);
84  h_CrystalThetaID67 = new TH1F("Crystal_Dose_ThetaID_67", "Crystal Radiation Dose vs #phi_{ID}, #theta_{ID}=67; #phi_{ID};Gy/yr", 64,
85  -0.5, 63.5);
86  h_BarrelDose = new TH1F("Crystal_Dose_Barrel", "Crystal Radiation Dose in Barrel, 12<#theta_{ID}<59; #phi_{ID}; Gy/yr", 144, -0.5,
87  143.5);
88  h_DiodeRadDose = new TH1F("Diode_Rad_Dose", "Diode Radiation Dose vs #theta_{ID};#theta_{ID};Gy/yr", 69, -0.5, 68.5);
89 
90  //hit locations
91  h_ProdVert = new TH1F("MCProd_Vert", "Production Vertex;z (cm)", 125, -200, 300);
92  h_HitLocations = new TH2F("Hit_Locations", "Hit locations;z (cm); r (cm)", 250, -200, 300, 80, 0, 160);
93  h_ProdVertvsThetaId = new TH2F("MCProd_Vert_vs_ThetaID", "Production Vertex vs #theta_{ID};#theta_{ID};z (cm)", 69, -0.5, 68.5, 125,
94  -200, 300);
95  hEdepPerRing = new TH1F("hEdepPerRing", "Energy deposited per #theta_{ID};#theta_{ID}; GeV", 69, -0.5, 68.5);
96  hNevtPerRing = new TH1F("hNevtPerRing", "Number of events #theta_{ID} (for pile-up);#theta_{ID};N_{event}", 69, -0.5, 68.5);
97 
98 
99 
100  //Neutrons
101  h_NeutronFluxThetaID2 = new TH1F("Neutron_Flux_ThetaID_2", "Diode Neutron Flux, #theta_{ID}=2;#phi_{ID}; yr^{-1}/cm^{-2}", 64, -0.5,
102  63.5);
103  h_NeutronFluxThetaID67 = new TH1F("Neutron_Flux_ThetaID_67", "Diode Neutron Flux, #theta_{ID}=67;#phi_{ID}; yr^{-1}/cm^{-2}", 64,
104  -0.5, 63.5);
105  h_NeutronFlux = new TH1F("Neutron_Flux", "Diode Neutron Flux vs #theta_{ID};#theta_{ID}; yr^{-1}/cm^{-2}", 69, -0.5, 68.5);
106  h_NeutronE = new TH1F("Neutron_Energy", "Neutron Energy; Energy (MeV)", 200, 0, 0.5);
107  h_NeutronEThetaID0 = new TH1F("Neutron_Energy_ThetaID0", "Neutron Energy, First Crystal; Energy (MeV)", 50, 0, 0.5);
108 
109  h_PhotonE = new TH1F("Photon_Energy", "Energy of photons creating hits in ECL; Energy (MeV)", 200, 0, 10);
110 
111  //showers
112  TString stime = s.str();
113  h_Shower = new TH1F("Shower_E_Dist", "Shower Energy distribution " + stime + " #mu s;GeV;# of showers", 100, 0, 0.5);
114  h_ShowerVsTheta = new TH2F("Shower_E_Dist_vs_theta", "Shower Energy distribution " + stime + " #mu s;GeV;#theta (deg)", 100, 0, 0.5,
115  180, 0, 180);
116 
117 
118 
119  //
120  // Below are for the ECL shields studies
121  //
123 
124  //Doses
125  hEMDose = new TH1F("hEMDose", "Crystal Radiation Dose; Cell ID ; Gy/yr", 8736, 0, 8736);
126  hEnergyPerCrystal = new TH1F("hEnergyPerCrystal", "Energy per crystal; Cell ID; GeV", 8736, 0, 8736);
127 
128  //Diodes
129  hDiodeFlux = new TH1F("hDiodeFlux", "Diode Neutron Flux ; Cell ID ; 1MeV-equiv / cm^{2} yr", 8736, 0, 8736);
130 
131  //Radiation spectra
132  hEgamma = new TH1F("hEgamma", "Log Spectrum of the photons hitting the crystals / 1 MeV; log_{10}(E_{#gamma}/1MeV) ", 500, -4, 3);
133  hEneu = new TH1F("hEneu", "Log Spectrum of the neutrons hitting the diodes / 1 MeV; log_{10}(E_{n}/1MeV)", 500, -10, 2);
134 
135  //ARICH plots
136  if (m_doARICH) {
137  hARICHDoseBB = new TH1F("hARICHDoseBB", "Radiation dose in ARICH boards (cBB); Ring-ID; Gy/yr", 7, -0.5, 6.5);
138  hHAPDFlux = new TH1F("hARICHnFlux", "1-MeV equivalent neutron flux in ARICH diodes (BB) ; Ring-ID ; 1-MeV-equiv / cm^{2} yr", 7,
139  -0.5, 6.5);
140  }
141 
142  hEMDoseECF = new TH2F();
143  hEMDoseECB = new TH2F();
144  hEMDoseBAR = new TH2F();
145  hEMDoseWideTID = new TH1F();
146 
147  hDiodeFluxECF = new TH2F();
148  hDiodeFluxECB = new TH2F();
149  hDiodeFluxBAR = new TH2F();
150  hDiodeFluxWideTID = new TH1F();
151 
152  hEnergyPerCrystalECF = new TH2F();
153  hEnergyPerCrystalECB = new TH2F();
154  hEnergyPerCrystalBAR = new TH2F();
155  hEnergyPerCrystalWideTID = new TH1F();
156 
157 
158 
159 }
160 
161 void ECLBackgroundModule::initialize()
162 {
163 
164  REG_HISTOGRAM
165 
166  if (m_doARICH) B2INFO("ECLBackgroundModule: ARICH plots are being produced");
167 
168  // Initialize variables
169 #ifdef DOARICH
170  if (m_doARICH) m_arichgp = ARICHGeometryPar::Instance();
171 #endif
172 
173  m_nEvent = 0;
174  BuildECL();
175 
176 }
177 
178 void ECLBackgroundModule::beginRun()
179 {
180 }
181 
182 void ECLBackgroundModule::event()
183 {
184 
185 
186 
187  //some variables that will be used many times
188  int m_cellID, m_thetaID, m_phiID, pid, NperRing;
189  double edep, theta, Energy, diodeDose, weightedFlux;
190  TVector3 rHit;
191 
192  //ignore events with huge number of SimHits (usually a glitchy event)
193  if (m_eclArray.getEntries() > 4000) {
194  B2INFO("ECLBackgroundModule: Skipping event #" << m_nEvent << " due to large number of ECLSimHits");
195  m_nEvent++;
196  return;
197  }
198 
199  bool isE = false;
200  //bool EinTheta[nECLThetaID] = {false};
201 
202  double edepSum = 0;
203  //double edepSumTheta[nECLThetaID] = {0};
204  //double E_tot[nECLCrystalTot] = {0};
205 
206 
207  auto edepSumTheta = new double[nECLThetaID]();
208  auto E_tot = new double[nECLCrystalTot]();
209 
210  auto EinTheta = new bool[nECLThetaID]();
211  std::fill_n(EinTheta, nECLThetaID, false);
212 
213 
214  h_nECLSimHits->Fill(m_eclArray.getEntries()); //number of ECL hits in an event
215 
216  //MC ID of photon hits
217  vector<int> MCPhotonIDs;
218 
219  int hitNum = m_eclArray.getEntries();
220  for (int i = 0; i < hitNum; i++) { //loop over ECLSimHits
221  ECLSimHit* aECLHit = m_eclArray[i];
222  m_cellID = aECLHit->getCellId() - 1; //cell ID
223  edep = aECLHit->getEnergyDep(); //energy deposited
224  G4ThreeVector hitPosn = aECLHit->getPosition(); //position of hit
225  pid = aECLHit->getPDGCode();
226  float Mass = Crystal[m_cellID]->GetMass();
227  m_thetaID = Crystal[m_cellID]->GetThetaID();
228  m_phiID = Crystal[m_cellID]->GetPhiID();
229  NperRing = Crystal[m_cellID]->GetNperThetaID(); //number of crystals in this theta ring
230  theta = Crystal[m_cellID]->GetTheta();
231 
232 
233  //get Track ID of photons which create the SimHits
234  if (pid == 22) MCPhotonIDs.push_back(aECLHit->getTrackId());
235 
236  edepSum = edepSum + edep;
237  E_tot[m_cellID] = edep + E_tot[m_cellID]; //sum energy deposited in this crystal
238  edepSumTheta[m_thetaID] = edepSumTheta[m_thetaID] + edep; //sum of energy for this thetaID value
239  EinTheta[m_thetaID] = true; //there is an energy deposit in this theta ring. used later
240  isE = true;
241 
242 
243  //fill histograms
244  //radiation dose for this SimHit
245  double dose = edep * GeVtoJ * usInYr / (m_sampleTime * Mass);
246 
247  h_CrystalRadDoseTheta->Fill(theta, dose / NperRing);
248  h_CrystalRadDose->AddBinContent(m_thetaID + 1, dose / NperRing);
249  hEMDose->AddBinContent(m_cellID + 1, dose);
250  hEnergyPerCrystal->AddBinContent(m_cellID + 1, edep);
251 
252  //2nd thetaID ring
253  if (m_thetaID == 2) {
254  h_CrystalThetaID2->AddBinContent(m_phiID + 1, dose);
255  }
256  //67th thetaID ring
257  if (m_thetaID == 67) {
258  h_CrystalThetaID67->AddBinContent(m_phiID + 1, dose);
259  }
260  //Barrel
261  if (m_thetaID < 59 && m_thetaID > 12) {
262  h_BarrelDose->AddBinContent(m_phiID + 1, dose / 46);
263  }
264 
265  //location of the hit
266  h_HitLocations->Fill(hitPosn.z(), hitPosn.perp());
267 
268  }
269 
270 
271  //for pileup noise estimation. To properly produce pileup noise plot, see comment at EOF.
272  for (int iECLCell = 0; iECLCell < nECLCrystalTot; iECLCell++) {
273  edep = E_tot[iECLCell];
274  m_thetaID = Crystal[iECLCell]->GetThetaID();
275  NperRing = Crystal[iECLCell]->GetNperThetaID();
276  if (edep > 0.000000000001) {
277  hNevtPerRing->Fill(m_thetaID, 1.0 / NperRing);
278  hEdepPerRing->Fill(m_thetaID, edep / NperRing);
279  }
280 
281  }
282 
283 
284  //One track can create several ECLSimHits. Remove the duplicates
285  sort(MCPhotonIDs.begin(), MCPhotonIDs.end());
286  vector<int>::iterator it;
287  it = std::unique(MCPhotonIDs.begin(), MCPhotonIDs.end());
288  MCPhotonIDs.resize(std::distance(MCPhotonIDs.begin() , it));
289 
290 
291  //loop over MCParticles to find the photons that caused the simhits
292  for (int i = 0; i < (int)MCPhotonIDs.size(); i++) {
293  for (int j = 0; j < m_mcParticles.getEntries(); j++) {
294  if (m_mcParticles[j]->getIndex() == MCPhotonIDs[i]) {
295  h_PhotonE->Fill(m_mcParticles[j]->getEnergy() * 1000);
296  hEgamma->Fill(log10(m_mcParticles[j]->getEnergy() * 1000));
297  break; //once the correct MCParticle is found, stop looping over MCParticles
298  }
299  }
300  }
301 
302  //*****************end of crystal analysis
303 
304  //start of diode analysis
305  int neuHits = m_BeamBackArray.getEntries();
306  for (int iHits = 0; iHits < neuHits; iHits++) { //loop over m_BeamBackArray
307  BeamBackHit* aBeamBackSimHit = m_BeamBackArray[iHits];
308 
309  //get relevant values
310  m_cellID = aBeamBackSimHit->getIdentifier();
311  double damage = aBeamBackSimHit->getNeutronWeight();
312  edep = aBeamBackSimHit->getEnergyDeposit();
313  pid = aBeamBackSimHit->getPDG();
314  int SubDet = aBeamBackSimHit->getSubDet();
315  Energy = aBeamBackSimHit->getEnergy();
316  rHit = aBeamBackSimHit->getPosition();
317 
318 
319  if (SubDet == 6) { //ECL
320 
321  m_thetaID = Crystal[m_cellID]->GetThetaID();
322  m_phiID = Crystal[m_cellID]->GetPhiID();
323  NperRing = Crystal[m_cellID]->GetNperThetaID();
324  diodeDose = edep * GeVtoJ * usInYr / (m_sampleTime * DiodeMass);
325  h_DiodeRadDose->AddBinContent(m_thetaID + 1, diodeDose / NperRing); //diode radiation dose plot
326 
327  if (pid == 2112) {
328 
329  weightedFlux = damage * usInYr / (m_sampleTime * DiodeArea) ; //neutrons per cm^2 per year
330 
331  //neutron plots
332  if (m_thetaID == 0) h_NeutronEThetaID0->Fill(Energy * 1000);
333  h_NeutronE->Fill(Energy * 1000);
334  hEneu->Fill(log10(Energy * 1000));
335 
336  h_NeutronFlux->AddBinContent(m_thetaID + 1, weightedFlux / NperRing);
337  hDiodeFlux->AddBinContent(m_cellID + 1, weightedFlux);
338 
339  if (m_thetaID == 2) h_NeutronFluxThetaID2->AddBinContent(m_phiID + 1 , weightedFlux);
340  if (m_thetaID == 67) h_NeutronFluxThetaID67->AddBinContent(m_phiID + 1 , weightedFlux);
341 
342 
343  }
344 
345  } else if (SubDet == 4 && m_doARICH) { //ARICH
346  FillARICHBeamBack(aBeamBackSimHit);
347  }
348  }
349 
350  int nShower = m_eclShowerArray.getEntries();
351  for (int i = 0; i < nShower; i++) {
352  ECLShower* aShower = m_eclShowerArray[i];
353 
354  Energy = aShower->getEnergy();
355  theta = aShower->getTheta();
356 
357  //get number of background showers with energy above 20MeV
358  if (Energy > 0.02) {
359  h_Shower->Fill(Energy);
360  h_ShowerVsTheta->Fill(Energy, theta * 180 / PI);
361 
362  }
363  }
364 
365 
366  for (int i = 0; i < nECLThetaID; i++) {
367  if (EinTheta[i]) {
368  //0th McParticle in an event is the origin of all particles
369  h_ProdVertvsThetaId->Fill(i, m_mcParticles[0]->getProductionVertex().z(), edepSumTheta[i]);
370  }
371  }
372 
373 
374  if (isE) {
375  h_ProdVert->Fill(m_mcParticles[0]->getProductionVertex().z(), edepSum);
376  }
377 
378  if (m_nEvent % ((int)m_sampleTime * 100) == 0) B2INFO("ECLBackgroundModule: At Event #" << m_nEvent);
379  m_nEvent++;
380 
381  delete[] edepSumTheta;
382  delete[] E_tot;
383  delete[] EinTheta;
384 
385 }
386 
387 void ECLBackgroundModule::endRun()
388 {
389  B2INFO("ECLBackgroundModule: Total Number of events: " << m_nEvent);
390 
391  //print doses of crystals of interest
392  for (int i = 0; i < (int)m_CryInt.size(); i++) {
393  if (m_CryInt[i] > 8736) {
394  B2WARNING("ECLBackgroundModule: Invalid cell ID. must be less than 8736");
395  continue;
396  }
397  double dose = hEMDose->GetBinContent(m_CryInt[i] + 1); //add 1 since bin #1 corrosponds to cell ID #0
398  int thetaID = Crystal[m_CryInt[i]]->GetThetaID();
399  int phiID = Crystal[m_CryInt[i]]->GetPhiID();
400  B2RESULT("Dose in Crystal " << m_CryInt[i] << ": " << dose << " ThetaID=" << thetaID << ", PhiID=" << phiID);
401  }
402 
403 
404  hEnergyPerCrystalECF = BuildPosHisto(hEnergyPerCrystal, "forward");
405  hEnergyPerCrystalECB = BuildPosHisto(hEnergyPerCrystal, "backward");
406  hEnergyPerCrystalBAR = BuildPosHisto(hEnergyPerCrystal, "barrel");
407  hEnergyPerCrystalWideTID = BuildThetaIDWideHisto(hEnergyPerCrystal);
408 
409 
410  hEMDoseECF = BuildPosHisto(hEMDose, "forward");
411  hEMDoseECB = BuildPosHisto(hEMDose, "backward");
412  hEMDoseBAR = BuildPosHisto(hEMDose, "barrel");
413  hEMDoseWideTID = BuildThetaIDWideHisto(hEMDose);
414 
415  hDiodeFluxECF = BuildPosHisto(hDiodeFlux, "forward");
416  hDiodeFluxECB = BuildPosHisto(hDiodeFlux, "backward");
417  hDiodeFluxBAR = BuildPosHisto(hDiodeFlux, "barrel");
418  hDiodeFluxWideTID = BuildThetaIDWideHisto(hDiodeFlux);
419 
420  hEMDose->SetTitle("Crystal Radiation Dose vs Cell ID");
421  hDiodeFlux->SetTitle("Diode Neutron Flux vs Cell ID");
422 
423 }
424 
425 void ECLBackgroundModule::terminate()
426 {
427 }
428 
429 
430 //
431 // Methods to study performance of ECL shields
432 // and potential impact on ARICH doses
434 #ifdef DOARICH
435 int ECLBackgroundModule::FillARICHBeamBack(BeamBackHit* aBBHit)
436 {
437 
438  double _damage = aBBHit->getNeutronWeight();
439  double _eDep = aBBHit->getEnergyDeposit();
440  float _trlen = aBBHit->getTrackLength();
441  int _pid = aBBHit->getPDG();
442  TVector3 _posHit = aBBHit->getPosition();
443 
444  int _moduleID = m_arichgp->getCopyNo(_posHit);
445 
446  double r = _posHit.Perp();
447  int _ring = 0;
448 
449  _ring = ARICHmod2row(_moduleID);
450 
451  B2DEBUG(200, "Filling ARICH BeamBackHit");
452  B2DEBUG(200, " PDG = " << _pid);
453  B2DEBUG(200, " Edep = " << _eDep);
454  B2DEBUG(200, " Ring = " << _ring);
455  B2DEBUG(200, " Radius = " << r);
456  B2DEBUG(200, " Module = " << _moduleID);
457 
458  if (2112 == _pid) {
459  hHAPDFlux->Fill(_ring, _damage * _trlen / HAPDthickness * usInYr / (m_sampleTime * HAPDarea * nHAPDperRing[_ring]));
460  } else {
461  hARICHDoseBB->Fill(_ring, _eDep / (HAPDmass * nHAPDperRing[_ring]) * GeVtoJ * usInYr / m_sampleTime);
462  }
463 
464  return 1;
465 }
466 
467 #else
468 int ECLBackgroundModule::FillARICHBeamBack(BeamBackHit* aBBHit) { return 1;}
469 #endif
470 
471 int ECLBackgroundModule::BuildECL()
472 {
473  for (int i = 0; i < nECLCrystalTot; i++) {
474  Crystal[i] = new ECLCrystalData(i);
475  }
476  return 1;
477 }
478 
479 //Method used for debugging.
480 int ECLBackgroundModule::SetPosHistos(TH1F* h, TH2F* hFWD, TH2F* hBAR, TH2F* hBWD)
481 {
482  char FWDtitle[100];
483  char BWDtitle[100];
484  char BARtitle[100];
485 
486  char FWDname[16];
487  char BWDname[16];
488  char BARname[16];
489 
490  sprintf(FWDtitle, "%s -- Forward Endcap" , h->GetTitle());
491  sprintf(BWDtitle, "%s -- Backward Endcap", h->GetTitle());
492  sprintf(BARtitle, "%s -- Barrel", h->GetTitle());
493 
494  sprintf(FWDname, "%sFWD", h->GetName());
495  sprintf(BWDname, "%sBWD", h->GetName());
496  sprintf(BARname, "%sBAR", h->GetName());
497 
498  // Fill 2D histograms with the values in the 1D histogram
499  for (int i = 0; i < nECLCrystalTot; i++) {
500  float value = h->GetBinContent(i + 1);
501 
502  if (i < nECLCrystalECF) {
503  hFWD->Fill(floor(Crystal[i]->GetX()), floor(Crystal[i]->GetY()), value);
504 
505  } else if (i >= (nECLCrystalBAR + nECLCrystalECF)) {
506  hBWD->Fill(floor(Crystal[i]->GetX()), floor(Crystal[i]->GetY()), value);
507 
508  } else
509  hBAR->Fill(floor(Crystal[i]->GetZ()), floor(Crystal[i]->GetR() * (Crystal[i]->GetPhi() - 180) * PI / 180) , value);
510  }
511 
512  return 1;
513 }
514 
515 
516 TH2F* ECLBackgroundModule::BuildPosHisto(TH1F* h, const char* sub)
517 {
518 
519  // Initialize variables
520  char _title[100];
521  char _name[16];
522  TH2F* h_out;
523 
524  double value = 0;
525 
526 // Forward endcap value vs (x,y)
527  if (!strcmp(sub, "forward")) {
528  sprintf(_name, "%sFWD", h->GetName());
529  sprintf(_title, "%s -- Forward Endcap;x(cm);y(cm)" , h->GetTitle());
530  h_out = new TH2F(_name, _title, 90, -150, 150, 90, -150, 150); //position in cm
531  h_out->Sumw2();
532  for (int i = 0; i < nECLCrystalECF; i++) {
533  value = h->GetBinContent(i + 1);
534  h_out->Fill(floor(Crystal[i]->GetX()),
535  floor(Crystal[i]->GetY()),
536  value);
537  }
538 
539  // Backward endcap value vs (x,y)
540  } else if (!strcmp(sub, "backward")) {
541  sprintf(_name, "%sBWD", h->GetName());
542  sprintf(_title, "%s -- Backward Endcap;x(cm);y(cm)", h->GetTitle());
543  h_out = new TH2F(_name, _title, 90, -150, 150, 90, -150, 150); //position in cm
544  h_out->Sumw2();
545  for (int i = (nECLCrystalBAR + nECLCrystalECF); i < nECLCrystalTot; i++) {
546  value = h->GetBinContent(i + 1);
547  h_out->Fill(floor(Crystal[i]->GetX()),
548  floor(Crystal[i]->GetY()),
549  value);
550  }
551 
552 
553  // The rest: barrel value vs (theta_ID, phi_ID)
554  } else if (!strcmp(sub, "barrel")) {
555  sprintf(_name, "%sBAR", h->GetName());
556  sprintf(_title, "%s -- Barrel;#theta_{ID};#phi_{ID}", h->GetTitle());
557  h_out = new TH2F(_name, _title, 47, 12, 59, 144, 0, 144); //position in cm (along z and along r*phi)
558  h_out->Sumw2();
559  for (int i = nECLCrystalECF; i < (nECLCrystalBAR + nECLCrystalECF); i++) {
560  value = h->GetBinContent(i + 1);
561  h_out->Fill(Crystal[i]->GetThetaID(), Crystal[i]->GetPhiID(), value);
562  }
563 
564  } else {
565  B2WARNING("ECLBackgroundModule: Unable to BuildPosHisto. Check Arguments.");
566  h_out = new TH2F("(empty)", "(empty)", 1, 0, 1, 1, 0, 1);
567  }
568 
569  return h_out;
570 }
571 
572 
573 TH1F* ECLBackgroundModule::BuildThetaIDWideHisto(TH1F* h_cry)
574 {
575 
576  char _title[100];
577  char _name[64];
578 
579  //Define the boundaries of the bins
580  static const int _nbins = 21;
581  static const double _xbins[] = { -0.5, 0.5, 4.5, 8.5, 11.5, 12.5,
582  16.5, 20.5, 24.5, 28.5, 32.5,
583  36.5, 40.5, 44.5, 48.5, 52.5,
584  56.5, 58.5, 59.5, 63.5, 67.5, 68.5
585  };
586 
587 
588  sprintf(_title, "%s vs #theta_{ID} -- averages" , h_cry->GetTitle());
589  sprintf(_name, "%svsTheWide", h_cry->GetName());
590 
591  //New pointer to the returned histogram ...
592  TH1F* h_out = new TH1F(_name, _title, 1, 0, 1);
593  // ... but only temp variables to the temporary ones
594  TH1F h_mass("h_mass", "Total Mass per Theta-ID", 1, 0, 1);
595  TH1F h_N("h_N", "Entries (unweighted) per Theta-ID bin", 1, 0, 1);
596 
597  //Apply all the same binning
598  h_out->SetBins(_nbins, _xbins);
599  h_mass.SetBins(_nbins, _xbins);
600  h_N.SetBins(_nbins, _xbins);
601 
602  h_out->SetTitle(_title);
603  h_out->Sumw2();
604 
605  //Make histo for total mass, then divide!
606  for (int i = 0; i < nECLCrystalTot; i++) {
607  h_out->Fill(Crystal[i]->GetThetaID(), h_cry->GetBinContent(i + 1) * Crystal[i]->GetMass());
608  h_mass.Fill(Crystal[i]->GetThetaID(), Crystal[i]->GetMass());
609  h_mass.SetBinError(Crystal[i]->GetThetaID(), 0);
610  h_N.Fill(Crystal[i]->GetThetaID());
611  }
612  h_out->SetXTitle("#theta_{ID}");
613  h_out->SetYTitle(h_cry->GetYaxis()->GetTitle());
614  h_out->Divide(&h_mass);
615 
616  return h_out;
617 }
618 
619 
620 int ECLBackgroundModule::ARICHmod2row(int modID)
621 {
622  if (modID <= 42) return 0;
623  else if (modID <= 90) return 1;
624  else if (modID <= 144) return 2;
625  else if (modID <= 204) return 3;
626  else if (modID <= 270) return 4;
627  else if (modID <= 342) return 5;
628  else if (modID <= 420) return 6;
629 
630  B2WARNING("ECLBackgroundModule: ARICHmod2row: modID out of bound; can't get ring index");
631  return -1;
632 }
633 
634 
635 
636 
637 /*
638 // In order to produce the pileup noise estimate, use this function (paste into another file):
639 void PileUpNoise(){
640 
641  const int numOfTypes = 8;
642  TString Filetypes[] = {"Touschek_HER", "Touschek_LER", "Coulomb_HER", "Coulomb_LER", "RBB_HER", "RBB_LER", "BHWide_HER", "BHWide_LER"};
643  TFile *f;
644 
645  double hEdepPerRing[69] = {};
646  double hNevtPerRing[69] = {};
647  double sampletime=1000;
648 
649  TH1F *h_Pileup = new TH1F("Pile_up", "Estimated Pile up Noise vs #theta_{ID}; #theta_{ID}; MeV", 69, -0.5, 68.5);
650 
651  for(int i=0; i<numOfTypes; i++){
652  f = new TFile(Filetypes[i]+".root"); //loads the sample files (eg, Touschek_HER.root, RBB_LER.root, etc)
653  TH1F *hNevt = (TH1F*)gROOT->FindObject("hNevtPerRing");
654  TH1F *hEdep = (TH1F*)gROOT->FindObject("hEdepPerRing");
655  for(int j=1; j<70; j++){
656  hEdepPerRing[j-1] = hEdepPerRing[j-1] + hEdep->GetBinContent(j);
657  hNevtPerRing[j-1] = hNevtPerRing[j-1] + hNevt->GetBinContent(j);
658  }
659  }
660 
661  for(int i=1; i<70; i++){
662  double Eavg = hEdepPerRing[i-1] / hNevtPerRing[i-1];
663  double pileup = sqrt( hNevtPerRing[i-1] / sampletime ) * Eavg * 1000;
664  h_Pileup->SetBinContent(i, pileup);
665  }
666 
667 }
668 */
Belle2::ECLCrystalData
Class for obtaining crystal details for a given crystal cell An evolved look-up table.
Definition: ECLCrystalData.h:31
Belle2::ECLBackgroundModule
A module to study background campaigns and produce histograms.
Definition: ECLBackgroundModule.h:46
Belle2::ECLSimHit::getPDGCode
int getPDGCode() const
Get Particle PDG (can be one of secondaries)
Definition: ECLSimHit.h:109
Belle2::ECLSimHit
ClassECLSimHit - Geant4 simulated hit for the ECL.
Definition: ECLSimHit.h:42
Belle2::ECLShower::getTheta
double getTheta() const
Get Theta.
Definition: ECLShower.h:296
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::BeamBackHit::getEnergyDeposit
double getEnergyDeposit() const
Get particle energy deposit in sensitive volume.
Definition: BeamBackHit.h:116
Belle2::ECLSimHit::getPosition
G4ThreeVector getPosition() const
Get Position.
Definition: ECLSimHit.h:139
Belle2::BeamBackHit::getEnergy
double getEnergy() const
Get energy of the particle.
Definition: BeamBackHit.h:110
Belle2::BeamBackHit::getPDG
int getPDG() const
Get the lund code of the particle that hit the sensitive area.
Definition: BeamBackHit.h:95
Belle2::ECLSimHit::getCellId
int getCellId() const
Get Cell ID.
Definition: ECLSimHit.h:99
Belle2::BeamBackHit::getTrackLength
double getTrackLength() const
the length of the track in the volume
Definition: BeamBackHit.h:119
Belle2::BeamBackHit::getPosition
const TVector3 & getPosition() const
Get global position of the particle hit.
Definition: BeamBackHit.h:101
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::ECLShower::getEnergy
double getEnergy() const
Get Energy.
Definition: ECLShower.h:286
Belle2::ECLSimHit::getEnergyDep
double getEnergyDep() const
Get Deposit energy.
Definition: ECLSimHit.h:119
Belle2::BeamBackHit::getNeutronWeight
double getNeutronWeight() const
get the effective neutron weigth
Definition: BeamBackHit.h:122
Belle2::ECLSimHit::getTrackId
int getTrackId() const
Get Track ID.
Definition: ECLSimHit.h:104
Belle2::BeamBackHit::getIdentifier
int getIdentifier() const
Get the identifier of subdetector component in which hit occured.
Definition: BeamBackHit.h:89
Belle2::BeamBackHit::getSubDet
int getSubDet() const
Det the index of subdetector in which hit occured.
Definition: BeamBackHit.h:92
Belle2::BeamBackHit
Class BeamBackHit - Stores hits from beam backgound simulation.
Definition: BeamBackHit.h:39
Belle2::HistoModule
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
Belle2::ECLShower
Class to store ECL Showers.
Definition: ECLShower.h:42