9 #include <ecl/modules/eclBackgroundStudy/ECLBackgroundModule.h>
17 #include <framework/logging/Logger.h>
20 #include <ecl/modules/eclBackgroundStudy/ECLCrystalData.h>
21 #include <ecl/dataobjects/ECLShower.h>
22 #include <ecl/dataobjects/ECLSimHit.h>
25 #include <arich/geometry/ARICHGeometryPar.h>
29 #include <simulation/dataobjects/BeamBackHit.h>
32 #include <mdst/dataobjects/MCParticle.h>
34 #define PI 3.14159265358979323846
52 setDescription(
"Processes background campaigns and produces histograms. Requires HistoManager");
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);
61 ECLBackgroundModule::~ECLBackgroundModule()
66 void ECLBackgroundModule::defineHisto()
72 h_nECLSimHits =
new TH1F(
"ECL_Sim_Hits",
"ECL Sim Hits", 100, 0, 100);
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,
80 h_CrystalThetaID67 =
new TH1F(
"Crystal_Dose_ThetaID_67",
"Crystal Radiation Dose vs #phi_{ID}, #theta_{ID}=67; #phi_{ID};Gy/yr", 64,
82 h_BarrelDose =
new TH1F(
"Crystal_Dose_Barrel",
"Crystal Radiation Dose in Barrel, 12<#theta_{ID}<59; #phi_{ID}; Gy/yr", 144, -0.5,
84 h_DiodeRadDose =
new TH1F(
"Diode_Rad_Dose",
"Diode Radiation Dose vs #theta_{ID};#theta_{ID};Gy/yr", 69, -0.5, 68.5);
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,
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);
97 h_NeutronFluxThetaID2 =
new TH1F(
"Neutron_Flux_ThetaID_2",
"Diode Neutron Flux, #theta_{ID}=2;#phi_{ID}; yr^{-1}/cm^{-2}", 64, -0.5,
99 h_NeutronFluxThetaID67 =
new TH1F(
"Neutron_Flux_ThetaID_67",
"Diode Neutron Flux, #theta_{ID}=67;#phi_{ID}; yr^{-1}/cm^{-2}", 64,
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);
105 h_PhotonE =
new TH1F(
"Photon_Energy",
"Energy of photons creating hits in ECL; Energy (MeV)", 200, 0, 10);
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,
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);
125 hDiodeFlux =
new TH1F(
"hDiodeFlux",
"Diode Neutron Flux ; Cell ID ; 1MeV-equiv / cm^{2} yr", 8736, 0, 8736);
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);
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,
138 hEMDoseECF =
new TH2F();
139 hEMDoseECB =
new TH2F();
140 hEMDoseBAR =
new TH2F();
141 hEMDoseWideTID =
new TH1F();
143 hDiodeFluxECF =
new TH2F();
144 hDiodeFluxECB =
new TH2F();
145 hDiodeFluxBAR =
new TH2F();
146 hDiodeFluxWideTID =
new TH1F();
148 hEnergyPerCrystalECF =
new TH2F();
149 hEnergyPerCrystalECB =
new TH2F();
150 hEnergyPerCrystalBAR =
new TH2F();
151 hEnergyPerCrystalWideTID =
new TH1F();
157 void ECLBackgroundModule::initialize()
163 B2INFO(
"ECLBackgroundModule: ARICH plots are being produced");
166 m_arichgp = ARICHGeometryPar::Instance();
175 void ECLBackgroundModule::beginRun()
179 void ECLBackgroundModule::event()
185 int m_cellID, m_thetaID, m_phiID, pid, NperRing;
186 double edep, theta, Energy, diodeDose, weightedFlux;
189 if (m_eclArray.getEntries() > 4000) {
190 B2INFO(
"ECLBackgroundModule: Skipping event #" << m_nEvent <<
" due to large number of ECLSimHits");
203 auto edepSumTheta =
new double[nECLThetaID]();
204 auto E_tot =
new double[nECLCrystalTot]();
206 auto EinTheta =
new bool[nECLThetaID]();
207 std::fill_n(EinTheta, nECLThetaID,
false);
210 h_nECLSimHits->Fill(m_eclArray.getEntries());
213 vector<int> MCPhotonIDs;
215 int hitNum = m_eclArray.getEntries();
216 for (
int i = 0; i < hitNum; i++) {
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();
226 theta = Crystal[m_cellID]->GetTheta();
230 if (pid == 22) MCPhotonIDs.push_back(aECLHit->
getTrackId());
232 edepSum = edepSum + edep;
233 E_tot[m_cellID] = edep + E_tot[m_cellID];
234 edepSumTheta[m_thetaID] = edepSumTheta[m_thetaID] + edep;
235 EinTheta[m_thetaID] =
true;
241 double dose = edep * GeVtoJ * usInYr / (m_sampleTime * Mass);
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);
249 if (m_thetaID == 2) {
250 h_CrystalThetaID2->AddBinContent(m_phiID + 1, dose);
253 if (m_thetaID == 67) {
254 h_CrystalThetaID67->AddBinContent(m_phiID + 1, dose);
257 if (m_thetaID < 59 && m_thetaID > 12) {
258 h_BarrelDose->AddBinContent(m_phiID + 1, dose / 46);
262 h_HitLocations->Fill(hitPosn.z(), hitPosn.perp());
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);
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));
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));
301 int neuHits = m_BeamBackArray.getEntries();
302 for (
int iHits = 0; iHits < neuHits; iHits++) {
303 BeamBackHit* aBeamBackSimHit = m_BeamBackArray[iHits];
309 pid = aBeamBackSimHit->
getPDG();
310 int SubDet = aBeamBackSimHit->
getSubDet();
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);
325 weightedFlux = damage * usInYr / (m_sampleTime * DiodeArea) ;
328 if (m_thetaID == 0) h_NeutronEThetaID0->Fill(Energy * 1000);
329 h_NeutronE->Fill(Energy * 1000);
330 hEneu->Fill(log10(Energy * 1000));
332 h_NeutronFlux->AddBinContent(m_thetaID + 1, weightedFlux / NperRing);
333 hDiodeFlux->AddBinContent(m_cellID + 1, weightedFlux);
335 if (m_thetaID == 2) h_NeutronFluxThetaID2->AddBinContent(m_phiID + 1 , weightedFlux);
336 if (m_thetaID == 67) h_NeutronFluxThetaID67->AddBinContent(m_phiID + 1 , weightedFlux);
341 }
else if (SubDet == 4 && m_doARICH) {
342 FillARICHBeamBack(aBeamBackSimHit);
346 int nShower = m_eclShowerArray.getEntries();
347 for (
int i = 0; i < nShower; i++) {
348 ECLShower* aShower = m_eclShowerArray[i];
355 h_Shower->Fill(Energy);
356 h_ShowerVsTheta->Fill(Energy, theta * 180 / PI);
362 for (
int i = 0; i < nECLThetaID; i++) {
365 h_ProdVertvsThetaId->Fill(i, m_mcParticles[0]->getProductionVertex().z(), edepSumTheta[i]);
371 h_ProdVert->Fill(m_mcParticles[0]->getProductionVertex().z(), edepSum);
374 if (m_nEvent % ((
int)m_sampleTime * 100) == 0) B2INFO(
"ECLBackgroundModule: At Event #" << m_nEvent);
377 delete[] edepSumTheta;
383 void ECLBackgroundModule::endRun()
385 B2INFO(
"ECLBackgroundModule: Total Number of events: " << m_nEvent);
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");
393 double dose = hEMDose->GetBinContent(m_CryInt[i] + 1);
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);
400 hEnergyPerCrystalECF = BuildPosHisto(hEnergyPerCrystal,
"forward");
401 hEnergyPerCrystalECB = BuildPosHisto(hEnergyPerCrystal,
"backward");
402 hEnergyPerCrystalBAR = BuildPosHisto(hEnergyPerCrystal,
"barrel");
403 hEnergyPerCrystalWideTID = BuildThetaIDWideHisto(hEnergyPerCrystal);
406 hEMDoseECF = BuildPosHisto(hEMDose,
"forward");
407 hEMDoseECB = BuildPosHisto(hEMDose,
"backward");
408 hEMDoseBAR = BuildPosHisto(hEMDose,
"barrel");
409 hEMDoseWideTID = BuildThetaIDWideHisto(hEMDose);
411 hDiodeFluxECF = BuildPosHisto(hDiodeFlux,
"forward");
412 hDiodeFluxECB = BuildPosHisto(hDiodeFlux,
"backward");
413 hDiodeFluxBAR = BuildPosHisto(hDiodeFlux,
"barrel");
414 hDiodeFluxWideTID = BuildThetaIDWideHisto(hDiodeFlux);
416 hEMDose->SetTitle(
"Crystal Radiation Dose vs Cell ID");
417 hDiodeFlux->SetTitle(
"Diode Neutron Flux vs Cell ID");
421 void ECLBackgroundModule::terminate()
431 int ECLBackgroundModule::FillARICHBeamBack(
BeamBackHit* aBBHit)
437 int _pid = aBBHit->
getPDG();
440 int _moduleID = m_arichgp->getCopyNo(_posHit);
442 double r = _posHit.Perp();
445 _ring = ARICHmod2row(_moduleID);
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);
455 hHAPDFlux->Fill(_ring, _damage * _trlen / HAPDthickness * usInYr / (m_sampleTime * HAPDarea * nHAPDperRing[_ring]));
457 hARICHDoseBB->Fill(_ring, _eDep / (HAPDmass * nHAPDperRing[_ring]) * GeVtoJ * usInYr / m_sampleTime);
464 int ECLBackgroundModule::FillARICHBeamBack(
BeamBackHit* aBBHit) {
return 1;}
467 int ECLBackgroundModule::BuildECL()
469 for (
int i = 0; i < nECLCrystalTot; i++) {
476 int ECLBackgroundModule::SetPosHistos(TH1F* h, TH2F* hFWD, TH2F* hBAR, TH2F* hBWD)
487 for (
int i = 0; i < nECLCrystalTot; i++) {
488 float value = h->GetBinContent(i + 1);
490 if (i < nECLCrystalECF) {
491 hFWD->Fill(floor(Crystal[i]->GetX()), floor(Crystal[i]->GetY()), value);
493 }
else if (i >= (nECLCrystalBAR + nECLCrystalECF)) {
494 hBWD->Fill(floor(Crystal[i]->GetX()), floor(Crystal[i]->GetY()), value);
497 hBAR->Fill(floor(Crystal[i]->GetZ()), floor(Crystal[i]->GetR() * (Crystal[i]->GetPhi() - 180) * PI / 180) , value);
504 TH2F* ECLBackgroundModule::BuildPosHisto(TH1F* h,
const char* sub)
508 TH2F* h_out =
nullptr;
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);
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()),
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);
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()),
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);
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);
549 B2WARNING(
"ECLBackgroundModule: Unable to BuildPosHisto. Check Arguments.");
550 h_out =
new TH2F(
"(empty)",
"(empty)", 1, 0, 1, 1, 0, 1);
557 TH1F* ECLBackgroundModule::BuildThetaIDWideHisto(TH1F* h_cry)
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
569 std::string _title = h_cry->GetTitle() + std::string(
" vs #theta_{ID} -- averages");
570 std::string _name = h_cry->GetName() + std::string(
"vsTheWide");
573 TH1F* h_out =
new TH1F(_name.c_str(), _title.c_str(), 1, 0, 1);
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);
579 h_out->SetBins(_nbins, _xbins);
580 h_mass.SetBins(_nbins, _xbins);
581 h_N.SetBins(_nbins, _xbins);
583 h_out->SetTitle(_title.c_str());
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());
593 h_out->SetXTitle(
"#theta_{ID}");
594 h_out->SetYTitle(h_cry->GetYaxis()->GetTitle());
595 h_out->Divide(&h_mass);
601 int ECLBackgroundModule::ARICHmod2row(
int modID)
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;
611 B2WARNING(
"ECLBackgroundModule: ARICHmod2row: modID out of bound; can't get ring index");
Class BeamBackHit - Stores hits from beam backgound simulation.
double getNeutronWeight() const
get the effective neutron weigth
const TVector3 & getPosition() const
Get global position of the particle hit.
double getEnergy() const
Get energy of the particle.
double getEnergyDeposit() const
Get particle energy deposit in sensitive volume.
int getPDG() const
Get the lund code of the particle that hit the sensitive area.
double getTrackLength() const
the length of the track in the volume
int getIdentifier() const
Get the identifier of subdetector component in which hit occured.
int getSubDet() const
Det the index of subdetector in which hit occured.
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.
double getEnergy() const
Get Energy.
double getTheta() const
Get Theta.
ClassECLSimHit - Geant4 simulated hit for the ECL.
int getPDGCode() const
Get Particle PDG (can be one of secondaries)
int getTrackId() const
Get Track ID.
int getCellId() const
Get Cell ID.
double getEnergyDep() const
Get Deposit energy.
G4ThreeVector getPosition() const
Get Position.
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Abstract base class for different kinds of events.