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