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