13 #include <ecl/modules/eclBackgroundStudy/ECLBackgroundModule.h>
21 #include <framework/logging/Logger.h>
24 #include <ecl/modules/eclBackgroundStudy/ECLCrystalData.h>
25 #include <ecl/dataobjects/ECLShower.h>
26 #include <ecl/dataobjects/ECLSimHit.h>
29 #include <arich/geometry/ARICHGeometryPar.h>
33 #include <simulation/dataobjects/BeamBackHit.h>
36 #include <mdst/dataobjects/MCParticle.h>
38 #define PI 3.14159265358979323846
56 setDescription(
"Processes background campaigns and produces histograms. Requires HistoManager");
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);
65 ECLBackgroundModule::~ECLBackgroundModule()
70 void ECLBackgroundModule::defineHisto()
76 h_nECLSimHits =
new TH1F(
"ECL_Sim_Hits",
"ECL Sim Hits", 100, 0, 100);
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,
84 h_CrystalThetaID67 =
new TH1F(
"Crystal_Dose_ThetaID_67",
"Crystal Radiation Dose vs #phi_{ID}, #theta_{ID}=67; #phi_{ID};Gy/yr", 64,
86 h_BarrelDose =
new TH1F(
"Crystal_Dose_Barrel",
"Crystal Radiation Dose in Barrel, 12<#theta_{ID}<59; #phi_{ID}; Gy/yr", 144, -0.5,
88 h_DiodeRadDose =
new TH1F(
"Diode_Rad_Dose",
"Diode Radiation Dose vs #theta_{ID};#theta_{ID};Gy/yr", 69, -0.5, 68.5);
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,
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);
101 h_NeutronFluxThetaID2 =
new TH1F(
"Neutron_Flux_ThetaID_2",
"Diode Neutron Flux, #theta_{ID}=2;#phi_{ID}; yr^{-1}/cm^{-2}", 64, -0.5,
103 h_NeutronFluxThetaID67 =
new TH1F(
"Neutron_Flux_ThetaID_67",
"Diode Neutron Flux, #theta_{ID}=67;#phi_{ID}; yr^{-1}/cm^{-2}", 64,
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);
109 h_PhotonE =
new TH1F(
"Photon_Energy",
"Energy of photons creating hits in ECL; Energy (MeV)", 200, 0, 10);
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,
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);
129 hDiodeFlux =
new TH1F(
"hDiodeFlux",
"Diode Neutron Flux ; Cell ID ; 1MeV-equiv / cm^{2} yr", 8736, 0, 8736);
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);
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,
142 hEMDoseECF =
new TH2F();
143 hEMDoseECB =
new TH2F();
144 hEMDoseBAR =
new TH2F();
145 hEMDoseWideTID =
new TH1F();
147 hDiodeFluxECF =
new TH2F();
148 hDiodeFluxECB =
new TH2F();
149 hDiodeFluxBAR =
new TH2F();
150 hDiodeFluxWideTID =
new TH1F();
152 hEnergyPerCrystalECF =
new TH2F();
153 hEnergyPerCrystalECB =
new TH2F();
154 hEnergyPerCrystalBAR =
new TH2F();
155 hEnergyPerCrystalWideTID =
new TH1F();
161 void ECLBackgroundModule::initialize()
166 if (m_doARICH) B2INFO(
"ECLBackgroundModule: ARICH plots are being produced");
170 if (m_doARICH) m_arichgp = ARICHGeometryPar::Instance();
178 void ECLBackgroundModule::beginRun()
182 void ECLBackgroundModule::event()
188 int m_cellID, m_thetaID, m_phiID, pid, NperRing;
189 double edep, theta, Energy, diodeDose, weightedFlux;
193 if (m_eclArray.getEntries() > 4000) {
194 B2INFO(
"ECLBackgroundModule: Skipping event #" << m_nEvent <<
" due to large number of ECLSimHits");
207 auto edepSumTheta =
new double[nECLThetaID]();
208 auto E_tot =
new double[nECLCrystalTot]();
210 auto EinTheta =
new bool[nECLThetaID]();
211 std::fill_n(EinTheta, nECLThetaID,
false);
214 h_nECLSimHits->Fill(m_eclArray.getEntries());
217 vector<int> MCPhotonIDs;
219 int hitNum = m_eclArray.getEntries();
220 for (
int i = 0; i < hitNum; i++) {
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();
230 theta = Crystal[m_cellID]->GetTheta();
234 if (pid == 22) MCPhotonIDs.push_back(aECLHit->
getTrackId());
236 edepSum = edepSum + edep;
237 E_tot[m_cellID] = edep + E_tot[m_cellID];
238 edepSumTheta[m_thetaID] = edepSumTheta[m_thetaID] + edep;
239 EinTheta[m_thetaID] =
true;
245 double dose = edep * GeVtoJ * usInYr / (m_sampleTime * Mass);
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);
253 if (m_thetaID == 2) {
254 h_CrystalThetaID2->AddBinContent(m_phiID + 1, dose);
257 if (m_thetaID == 67) {
258 h_CrystalThetaID67->AddBinContent(m_phiID + 1, dose);
261 if (m_thetaID < 59 && m_thetaID > 12) {
262 h_BarrelDose->AddBinContent(m_phiID + 1, dose / 46);
266 h_HitLocations->Fill(hitPosn.z(), hitPosn.perp());
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);
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));
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));
305 int neuHits = m_BeamBackArray.getEntries();
306 for (
int iHits = 0; iHits < neuHits; iHits++) {
307 BeamBackHit* aBeamBackSimHit = m_BeamBackArray[iHits];
313 pid = aBeamBackSimHit->
getPDG();
314 int SubDet = aBeamBackSimHit->
getSubDet();
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);
329 weightedFlux = damage * usInYr / (m_sampleTime * DiodeArea) ;
332 if (m_thetaID == 0) h_NeutronEThetaID0->Fill(Energy * 1000);
333 h_NeutronE->Fill(Energy * 1000);
334 hEneu->Fill(log10(Energy * 1000));
336 h_NeutronFlux->AddBinContent(m_thetaID + 1, weightedFlux / NperRing);
337 hDiodeFlux->AddBinContent(m_cellID + 1, weightedFlux);
339 if (m_thetaID == 2) h_NeutronFluxThetaID2->AddBinContent(m_phiID + 1 , weightedFlux);
340 if (m_thetaID == 67) h_NeutronFluxThetaID67->AddBinContent(m_phiID + 1 , weightedFlux);
345 }
else if (SubDet == 4 && m_doARICH) {
346 FillARICHBeamBack(aBeamBackSimHit);
350 int nShower = m_eclShowerArray.getEntries();
351 for (
int i = 0; i < nShower; i++) {
352 ECLShower* aShower = m_eclShowerArray[i];
359 h_Shower->Fill(Energy);
360 h_ShowerVsTheta->Fill(Energy, theta * 180 / PI);
366 for (
int i = 0; i < nECLThetaID; i++) {
369 h_ProdVertvsThetaId->Fill(i, m_mcParticles[0]->getProductionVertex().z(), edepSumTheta[i]);
375 h_ProdVert->Fill(m_mcParticles[0]->getProductionVertex().z(), edepSum);
378 if (m_nEvent % ((
int)m_sampleTime * 100) == 0) B2INFO(
"ECLBackgroundModule: At Event #" << m_nEvent);
381 delete[] edepSumTheta;
387 void ECLBackgroundModule::endRun()
389 B2INFO(
"ECLBackgroundModule: Total Number of events: " << m_nEvent);
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");
397 double dose = hEMDose->GetBinContent(m_CryInt[i] + 1);
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);
404 hEnergyPerCrystalECF = BuildPosHisto(hEnergyPerCrystal,
"forward");
405 hEnergyPerCrystalECB = BuildPosHisto(hEnergyPerCrystal,
"backward");
406 hEnergyPerCrystalBAR = BuildPosHisto(hEnergyPerCrystal,
"barrel");
407 hEnergyPerCrystalWideTID = BuildThetaIDWideHisto(hEnergyPerCrystal);
410 hEMDoseECF = BuildPosHisto(hEMDose,
"forward");
411 hEMDoseECB = BuildPosHisto(hEMDose,
"backward");
412 hEMDoseBAR = BuildPosHisto(hEMDose,
"barrel");
413 hEMDoseWideTID = BuildThetaIDWideHisto(hEMDose);
415 hDiodeFluxECF = BuildPosHisto(hDiodeFlux,
"forward");
416 hDiodeFluxECB = BuildPosHisto(hDiodeFlux,
"backward");
417 hDiodeFluxBAR = BuildPosHisto(hDiodeFlux,
"barrel");
418 hDiodeFluxWideTID = BuildThetaIDWideHisto(hDiodeFlux);
420 hEMDose->SetTitle(
"Crystal Radiation Dose vs Cell ID");
421 hDiodeFlux->SetTitle(
"Diode Neutron Flux vs Cell ID");
425 void ECLBackgroundModule::terminate()
435 int ECLBackgroundModule::FillARICHBeamBack(
BeamBackHit* aBBHit)
441 int _pid = aBBHit->
getPDG();
444 int _moduleID = m_arichgp->getCopyNo(_posHit);
446 double r = _posHit.Perp();
449 _ring = ARICHmod2row(_moduleID);
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);
459 hHAPDFlux->Fill(_ring, _damage * _trlen / HAPDthickness * usInYr / (m_sampleTime * HAPDarea * nHAPDperRing[_ring]));
461 hARICHDoseBB->Fill(_ring, _eDep / (HAPDmass * nHAPDperRing[_ring]) * GeVtoJ * usInYr / m_sampleTime);
468 int ECLBackgroundModule::FillARICHBeamBack(
BeamBackHit* aBBHit) {
return 1;}
471 int ECLBackgroundModule::BuildECL()
473 for (
int i = 0; i < nECLCrystalTot; i++) {
480 int ECLBackgroundModule::SetPosHistos(TH1F* h, TH2F* hFWD, TH2F* hBAR, TH2F* hBWD)
490 sprintf(FWDtitle,
"%s -- Forward Endcap" , h->GetTitle());
491 sprintf(BWDtitle,
"%s -- Backward Endcap", h->GetTitle());
492 sprintf(BARtitle,
"%s -- Barrel", h->GetTitle());
494 sprintf(FWDname,
"%sFWD", h->GetName());
495 sprintf(BWDname,
"%sBWD", h->GetName());
496 sprintf(BARname,
"%sBAR", h->GetName());
499 for (
int i = 0; i < nECLCrystalTot; i++) {
500 float value = h->GetBinContent(i + 1);
502 if (i < nECLCrystalECF) {
503 hFWD->Fill(floor(Crystal[i]->GetX()), floor(Crystal[i]->GetY()), value);
505 }
else if (i >= (nECLCrystalBAR + nECLCrystalECF)) {
506 hBWD->Fill(floor(Crystal[i]->GetX()), floor(Crystal[i]->GetY()), value);
509 hBAR->Fill(floor(Crystal[i]->GetZ()), floor(Crystal[i]->GetR() * (Crystal[i]->GetPhi() - 180) * PI / 180) , value);
516 TH2F* ECLBackgroundModule::BuildPosHisto(TH1F* h,
const char* sub)
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);
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()),
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);
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()),
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);
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);
565 B2WARNING(
"ECLBackgroundModule: Unable to BuildPosHisto. Check Arguments.");
566 h_out =
new TH2F(
"(empty)",
"(empty)", 1, 0, 1, 1, 0, 1);
573 TH1F* ECLBackgroundModule::BuildThetaIDWideHisto(TH1F* h_cry)
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
588 sprintf(_title,
"%s vs #theta_{ID} -- averages" , h_cry->GetTitle());
589 sprintf(_name,
"%svsTheWide", h_cry->GetName());
592 TH1F* h_out =
new TH1F(_name, _title, 1, 0, 1);
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);
598 h_out->SetBins(_nbins, _xbins);
599 h_mass.SetBins(_nbins, _xbins);
600 h_N.SetBins(_nbins, _xbins);
602 h_out->SetTitle(_title);
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());
612 h_out->SetXTitle(
"#theta_{ID}");
613 h_out->SetYTitle(h_cry->GetYaxis()->GetTitle());
614 h_out->Divide(&h_mass);
620 int ECLBackgroundModule::ARICHmod2row(
int modID)
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;
630 B2WARNING(
"ECLBackgroundModule: ARICHmod2row: modID out of bound; can't get ring index");