10 #include <top/modules/TOPBackground/TOPBackgroundModule.h>
12 #include <top/dataobjects/TOPBarHit.h>
13 #include <top/dataobjects/TOPDigit.h>
14 #include <top/dataobjects/TOPSimHit.h>
15 #include <top/geometry/TOPGeometryPar.h>
17 #include <mdst/dataobjects/MCParticle.h>
18 #include <simulation/dataobjects/BeamBackHit.h>
21 #include <framework/datastore/DataStore.h>
22 #include <framework/datastore/StoreArray.h>
25 #include <framework/logging/Logger.h>
26 #include <framework/gearbox/Const.h>
35 #include <geometry/Materials.h>
36 #include <G4Material.hh>
39 using namespace boost;
98 setDescription(
"A module to analyze beam background simulations regarding TOP");
101 addParam(
"Type", m_BkgType,
"Backgound type" ,
string(
"Backgound"));
102 addParam(
"Output", m_OutputFileName,
"Name of the output file",
103 string(
"Backgound.root"));
104 addParam(
"TimeOfSimulation", m_TimeOfSimulation,
105 "Real time in micro seconds that corresponds to simulated data", 5.);
108 TOPBackgroundModule::~TOPBackgroundModule()
113 void TOPBackgroundModule::initialize()
119 m_rootFile = TFile::Open(m_OutputFileName.c_str(),
"RECREATE");
120 origingamma =
new TTree(
"origingamma",
"tree");
121 originpe =
new TTree(
"originpe",
"tree2");
123 origingamma->Branch(
"x", &origingamma_x);
124 origingamma->Branch(
"y", &origingamma_y);
125 origingamma->Branch(
"z", &origingamma_z);
126 originpe->Branch(
"x", &originpe_x);
127 originpe->Branch(
"y", &originpe_y);
128 originpe->Branch(
"z", &originpe_z);
130 peflux =
new TH1F(
"Photoelectron flux",
"Photoelectron flux", 33, -11.25, 360);
131 nflux =
new TH1F(
"Neutron flux",
"Neutron flux", 33, -11.25, 360);
132 rdose =
new TH1F(
"Radiation dose",
"Radiation dose", 33, -11.25, 360);
133 zdist =
new TH1F(
"Photoelectron origin",
"Photoelectron origin", 200, -400, 400);
134 genergy =
new TH1F(
"Energy distribution of photons",
"Energy distribution of photons", 150, 0, 5);
135 genergy2 =
new TH1F(
"Energy distribution of gammas",
"Energy distribution of gammas", 500, 0, 5);
137 zdistg =
new TH1F(
"Photoelectron flux z",
"Photoelectron flux Z projection", 800, -400, 400);
138 originpt =
new TH1F(
"P_t of origin electron",
"P_t of origin electron", 300, 0.06, 0.14);
140 nflux_bar =
new TH2F(
"Neutron flux on bar",
"Neutron flux on bar", 32, -114.8, 211.5, 16, -0, 360);
141 gflux_bar =
new TH2F(
"Gamma flux on bar",
"Gamma flux on bar MHz/cm^{2}", 32, -114.8, 211.5, 16, -0, 360);
142 cflux_bar =
new TH2F(
"Charged flux on bar",
"Charged flux on bar MHz/cm^{2}", 32, -114.8, 211.5, 16, -0, 360);
144 norigin =
new TH1F(
"neutron(BAR) origin",
"neutron(BAR) origin", 200, -400, 400);
145 gorigin =
new TH1F(
"gamma(BAR) origin",
"gamma(BAR) origin", 200, -400, 400);
146 corigin =
new TH1F(
"charged(BAR) origin",
"charged(BAR) origin", 200, -400, 400);
148 nprim =
new TH1F(
"neutron(BAR) primary",
"neutron(BAR) primary", 200, -400, 400);
149 gprim =
new TH1F(
"gamma(BAR) primary",
"gamma(BAR) primary", 200, -400, 400);
150 cprim =
new TH1F(
"charged(BAR) primary",
"charged(BAR) primary", 200, -400, 400);
152 origin_zx =
new TGraph();
153 origin_zy =
new TGraph();
155 prim_zx =
new TGraph();
156 prim_zy =
new TGraph();
157 module_occupancy =
new TGraph();
159 origin_zy->SetName(
"originZY");
160 origin_zx->SetName(
"originZX");
161 module_occupancy->SetName(
"occupancy");
163 const auto& geo = TOP::TOPGeometryPar::Instance()->getGeometry()->getFrontEnd();
164 double S1 = geo.getFrontBoardWidth() * geo.getFrontBoardHeight();
165 double S2 = geo.getHVBoardWidth() * geo.getHVBoardLength();
166 double V1 = S1 * geo.getFrontBoardThickness();
167 double V2 = S2 * geo.getHVBoardThickness();
168 G4Material* material1 = geometry::Materials::get(geo.getFrontBoardMaterial());
169 if (!material1) B2FATAL(
"Material '" << geo.getFrontBoardMaterial() <<
"' not found");
170 G4Material* material2 = geometry::Materials::get(geo.getHVBoardMaterial());
171 if (!material2) B2FATAL(
"Material '" << geo.getHVBoardMaterial() <<
"' not found");
172 double density1 = material1->GetDensity() / CLHEP::kg * CLHEP::cm * CLHEP::cm * CLHEP::cm;
173 double density2 = material2->GetDensity() / CLHEP::kg * CLHEP::cm * CLHEP::cm * CLHEP::cm;
176 PCBmass = V1 * density1 + V2 * density2;
179 evtoJ = 1.60217653 * 1e-10;
180 mtoc = 1.97530864197531;
192 void TOPBackgroundModule::beginRun()
195 B2INFO(
"TOPBackground: Processing:");
199 void TOPBackgroundModule::event()
207 for (
int i = 0; i < nHits; i++) {
211 peflux->AddBinContent(barID * 2, 1. / m_TimeOfSimulation / 32.0);
213 const TOPSimHit* simHit = DataStore::getRelated<TOPSimHit>(aDigit);
214 if (!simHit)
continue;
217 module_occupancy->SetPoint(count_occ, PMTID , barID);
222 const MCParticle* particle = DataStore::getRelated<MCParticle>(simHit);
231 const MCParticle* pommother = mother->getMother();
234 zdist->Fill(mother->getVertex().Z());
235 zdistg->Fill(mother->getVertex().Z(), 1. / m_TimeOfSimulation / 32.0 / 16.0);
236 originpe_x = mother->getVertex().X();
237 originpe_y = mother->getVertex().Y();
238 originpe_z = mother->getVertex().Z();
240 TVector3 momentum = mother->getMomentum();
241 if (m_BkgType.at(m_BkgType.size() - 3) ==
'L') momentum.RotateY(0.0415);
242 else if (m_BkgType.at(m_BkgType.size() - 3) ==
'H') momentum.RotateY(-0.0415);
243 double px = momentum.X();
244 double py = momentum.Y();
245 double pt = sqrt(px * px + py * py);
258 for (
int iHit = 0; iHit < nHits; ++iHit) {
261 if (subdet != 5)
continue;
264 double phi = pos.XYvector().Phi_0_2pi(pos.XYvector().Phi()) / 3.14159265358979 * 180.;
265 int barID = int (phi / 22.5 + 0.5);
272 if (tophit->
getPDG() == Const::neutron.getPDGCode()) {
276 nflux->AddBinContent(barID * 2, w / m_TimeOfSimulation / PCBarea * yearns * tlen / 0.2);
280 rdose->AddBinContent(barID * 2, edep / m_TimeOfSimulation * yearns / PCBmass * evtoJ);
285 for (
int iHit = 0; iHit < nHits; ++iHit) {
288 int PDG = toptrk->
getPDG();
291 if (PDG == Const::neutron.getPDGCode()) {
292 nflux_bar->Fill(toptrk->
getPosition().Z(), (barID - 1) * 22.5,
293 1. / 917.65 / m_TimeOfSimulation * yearns * 2.);
296 if (PDG == Const::photon.getPDGCode()) {
297 gflux_bar->Fill(toptrk->
getPosition().Z(), (barID - 1) * 22.5,
298 1. / 917.65 / m_TimeOfSimulation * 2.);
300 genergy2->Fill(toptrk->
getMomentum().Mag() * 1000);
312 cflux_bar->Fill(toptrk->
getPosition().Z(), (barID - 1) * 22.5,
313 1. / 917.65 / m_TimeOfSimulation * 2.);
322 void TOPBackgroundModule::myprint(TH1F* histo,
const char* path,
const char* xtit =
"",
const char* ytit =
"",
double tresh = 0)
326 gStyle->SetOptStat(
"");
327 gStyle->SetOptFit(1111);
329 gStyle->SetCanvasColor(-1);
330 gStyle->SetPadColor(-1);
331 gStyle->SetFrameFillColor(-1);
332 gStyle->SetHistFillColor(-1);
333 gStyle->SetTitleFillColor(-1);
334 gStyle->SetFillColor(-1);
335 gStyle->SetFillStyle(4000);
336 gStyle->SetStatStyle(0);
337 gStyle->SetTitleStyle(0);
338 gStyle->SetCanvasBorderSize(0);
339 gStyle->SetCanvasBorderMode(0);
340 gStyle->SetPadBorderMode(0);
342 gStyle->SetFrameBorderSize(0);
343 gStyle->SetLegendBorderSize(0);
344 gStyle->SetStatBorderSize(0);
345 gStyle->SetTitleBorderSize(0);
350 TCanvas* c1 =
new TCanvas(
"c1",
"", 1920, 1200);
352 double x1 = histo->GetBinLowEdge(1);
353 double nb = histo->GetNbinsX();
354 double bin = histo->GetBinWidth(1);
355 double x2 = x1 + bin * nb;
357 double max = histo->GetBinContent(histo->GetMaximumBin());
360 histo->GetYaxis()->SetRangeUser(0, tresh * 1.1);
363 TLine* line =
new TLine(x1, tresh, x2, tresh);
364 line->SetLineColor(1);
365 line->SetLineWidth(3);
366 line->SetLineStyle(2);
369 histo->SetFillColor(2);
370 histo->SetLineColor(1);
372 gPad->SetTopMargin(0.08);
373 gPad->SetBottomMargin(0.15);
376 histo->GetXaxis()->SetLabelSize(0.06);
377 histo->GetYaxis()->SetLabelSize(0.06);
378 histo->GetXaxis()->SetTitleSize(0.06);
379 histo->GetYaxis()->SetTitleSize(0.06);
380 histo->GetXaxis()->SetTitle(xtit);
381 histo->GetYaxis()->SetTitle(ytit);
382 histo->GetXaxis()->SetTitleOffset(0.9);
383 histo->GetYaxis()->SetTitleOffset(0.7);
387 TLegend* leg =
new TLegend(0.75, 0.95, 0.90, 1.00);
388 leg->AddEntry(histo, m_BkgType.c_str(),
"pf");
399 void TOPBackgroundModule::endRun()
401 B2INFO(
"TOPBackground: Finished:");
404 void TOPBackgroundModule::terminate()
414 origingamma->Write();
426 module_occupancy->Write();
434 B2INFO(
"TOPBackground finished");
439 void TOPBackgroundModule::printModuleParams()
const
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 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 getSubDet() const
Det the index of subdetector in which hit occured.
A Class to store the Monte Carlo particle information.
Accessor to arrays stored in the data store.
int getEntries() const
Get the number of objects in the array.
Class to store track parameters of incoming MC particles relation to MCParticle filled in top/simulat...
TVector3 getMomentum() const
Returns impact momentum.
int getModuleID() const
Returns module ID.
int getPDG() const
Returns PDG code of particle.
TVector3 getPosition() const
Returns impact point.
TVector3 getProductionPoint() const
Returns production point.
Class to store TOP digitized hits (output of TOPDigitizer or raw data unpacker) relations to TOPSimHi...
int getModuleID() const
Returns module ID.
Class to store simulated hits of Cherenkov photons on PMT's input for digitization module (TOPDigitiz...
double getEnergy() const
Returns photon energy.
int getPmtID() const
Returns PMT ID.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
MCParticle * getMother() const
Returns a pointer to the mother particle.
Abstract base class for different kinds of events.