R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandMCMon.cxx
Go to the documentation of this file.
1/******************************************************************************
2 * Copyright (C) 2019 GSI Helmholtzzentrum für Schwerionenforschung GmbH *
3 * Copyright (C) 2019-2025 Members of R3B Collaboration *
4 * *
5 * This software is distributed under the terms of the *
6 * GNU General Public Licence (GPL) version 3, *
7 * copied verbatim in the file "LICENSE". *
8 * *
9 * In applying this license GSI does not waive the privileges and immunities *
10 * granted to it by virtue of its status as an Intergovernmental Organization *
11 * or submit itself to any jurisdiction. *
12 ******************************************************************************/
13
14#include "R3BNeulandMCMon.h"
15#include "FairLogger.h"
16#include "FairRootManager.h"
17#include "FairRun.h"
18#include "Math/Vector4D.h"
19#include "TDirectory.h"
20#include "TFile.h"
21#include "TH1D.h"
22#include "TH2D.h"
23#include "TH3D.h"
24#include <algorithm>
25#include <iostream>
26#include <numeric>
27#include <string>
28
29inline Bool_t IsPrimaryNeutron(const R3BMCTrack* mcTrack)
30{
31 return (mcTrack->GetPdgCode() == 2112 && mcTrack->GetMotherId() == -1);
32}
33
34inline Bool_t IsMotherPrimaryNeutron(const R3BMCTrack* mcTrack, const std::vector<R3BMCTrack*>& tracks)
35{
36 if (mcTrack->GetMotherId() < 0 || size_t(mcTrack->GetMotherId()) >= tracks.size())
37 {
38 return false;
39 }
40 return IsPrimaryNeutron(tracks.at(mcTrack->GetMotherId()));
41}
42
43inline Double_t GetKineticEnergy(const R3BMCTrack* mcTrack)
44{
45 return (mcTrack->GetEnergy() - mcTrack->GetMass()) * 1000.;
46}
47
48inline Double_t GetTheta(const R3BMCTrack* mcTrack) { return std::acos(mcTrack->GetPz() / mcTrack->GetP()); }
49
50static std::map<int, std::string> lookup_table = {
51 { -211, "pion" }, { 22, "gamma" }, { 111, "pion" }, { 211, "pion" }, { 2112, "n" },
52 { 2212, "p" }, { 1000010020, "d, t" }, { 1000010030, "d, t" }, { 1000020030, "alpha" }, { 1000020040, "alpha" },
53};
54
55std::string lookup(int pid)
56{
57 const auto l = lookup_table.find(pid);
58 if (l == lookup_table.end())
59 {
60 return std::string("heavy");
61 }
62 else
63 {
64 return l->second;
65 }
66}
67
68template <typename T>
69void writeout(std::map<T, TH1D*>& map, const std::string& what, const int nEvents = 0)
70{
71 if (!map.empty())
72 {
73 std::ostream* out = &std::cout;
74 std::ofstream of;
75 if (FairRun::Instance() != nullptr)
76 {
77 const std::string f =
78 std::string{ FairRun::Instance()->GetUserOutputFileName().Data() } + std::string(".") + what + ".dat";
79 std::cout << "Writing to file " << f << std::endl;
80 of.open(f);
81 out = &of;
82 }
83 gDirectory->mkdir(what.c_str());
84 gDirectory->cd(what.c_str());
85 for (const auto& kv : map)
86 {
87 auto h = kv.second;
88 const auto c = h->GetEntries();
89 if (nEvents > 0 && nEvents > c)
90 {
91 h->Fill(0., nEvents - c);
92 }
93 *out << kv.first << "\t" << c << "\t" << h->GetMean() << "\t" << h->GetMeanError() << std::endl;
94 kv.second->Write();
95 }
96 gDirectory->cd("..");
97 }
98}
99
100R3BNeulandMCMon::R3BNeulandMCMon(const Option_t* option)
101 : FairTask("R3B NeuLAND Neuland Monte Carlo Monitor")
102 , fIs3DTrackEnabled(false)
103 , fIsFullSimAnaEnabled(false)
104 , fPrimaryNeutronInteractionPoints("NeulandPrimaryPoints")
105 , fMCTracks("MCTrack")
106 , fNeulandPoints("NeulandPoints")
107 , nEvents(0)
108{
109 LOG(info) << "Using R3B NeuLAND Neuland Monte Carlo Monitor";
110
111 TString opt = option;
112 opt.ToUpper();
113 if (opt.Contains("3DTRACK"))
114 {
115 fIs3DTrackEnabled = true;
116 LOG(info) << "... with 3D track visualization";
117 }
118 if (opt.Contains("FULLSIMANA"))
119 {
120 fIsFullSimAnaEnabled = true;
121 LOG(info) << "... with full simulation neutron reaction product analysis";
122 }
123}
124
126{
127 fPrimaryNeutronInteractionPoints.init();
128 fMCTracks.Init();
129 fNeulandPoints.init();
130
131 fhNPNIPsEToFVSTime = new TH2D("NPNIPEToFVSTime", "NPNIP E_{ToF} vs. NPNIP Time", 100, 0, 1000, 500, 0, 500);
132 fhNPNIPsEToFVSTime->GetXaxis()->SetTitle("NPNIP E_{ToF} [MeV]");
133 fhNPNIPsEToFVSTime->GetYaxis()->SetTitle("NPNIP t [ns]");
134
135 fhPDG = new TH1D("hPDG",
136 "Number of particles by PDG code created by primary neutron interaction resulting in a point",
137 8000,
138 -4000,
139 4000);
140 fhEPrimarys = new TH1D("hE_primarys", "Energy of primary particles", 10000, 0, 10000);
141 fhEPrimarys2 = new TH1D("hE_primarys2", "Tot. Energy of primary particles", 100000, 0, 100000);
142 fhEPrimaryNeutrons = new TH1D("hE_primary_neutrons", "Energy of primary Neutrons", 10000, 0, 10000);
143 fhErelMC = new TH1D("hErelMC", "Erel", 5000, 0, 5000);
144 fhEtotPrim = new TH1D("hE_tot_prim",
145 "Total Light Yield of non-neutron LandPoints created by primary neutron interaction(s)",
146 6000,
147 0,
148 6000);
149 fhEtot = new TH1D("hE_tot", "Total Light Yield of non-neutron LandPoints", 6000, 0, 6000);
150 fhESecondaryNeutrons = new TH1D(
151 "hE_secondary_neutrons", "Energy of neutron tracks created by primary neutron interaction", 6000, 0, 6000);
152 fhMotherIDs = new TH1D("hmotherIDs", "MotherIDs", 6010, -10, 6000);
153 fhPrimaryDaughterIDs = new TH1D("hprimary_daughter_IDs", "IDs of tracks with a primary mother", 6001, -1, 6000);
154 fhMCToF = new TH1D("fhMCToF", "Energy of primary Neutron - ToF Energy from PNIPS", 2001, -1000, 1000);
155
156 fhNPNIPSrvsz = new TH2D("fhNPNIPSrvsz", "NPNIPS R = #sqrt(x**2+y**2) vs z", 600, 1000, 4000, 40, 0, 200);
157 fhNPNIPSrvsz->GetXaxis()->SetTitle("NPNIP Z [cm]");
158 fhNPNIPSrvsz->GetYaxis()->SetTitle("NPNIP R [cm]");
159
160 fhNPNIPSxy = new TH2D("fhNPNIPSxy", "NPNIPS xy", 101, -200, 200, 101, -200, 200);
161 fhNPNIPSxy->GetXaxis()->SetTitle("NPNIP X [cm]");
162 fhNPNIPSxy->GetYaxis()->SetTitle("NPNIP Y [cm]");
163
164 fhnNPNIPs = new TH1D("fhnNPNIPs", "Number of NPNIPs per event", 10, -0.5, 9.5);
165
166 fhnProducts = new TH1D("nProducts", "Number of products created by primary neutron interaction", 20, -0.5, 19.5);
167 fhnProductsCharged =
168 new TH1D("nProductsCharged",
169 "Number of products created by primary neutron interaction w/o gammas, neutrons, pion0",
170 20,
171 -0.5,
172 19.5);
173 fhSumProductEnergy = new TH1D("SumProductEnergy", "Sum product energy", 1000, 0, 1000);
174 fhnSecondaryNeutrons = new TH1D("NumberSecondaryNeutrons", "Number of Secondary Neutrons", 10, -0.5, 9.5);
175 fhnSecondaryProtons = new TH1D("NumberSecondaryProtons", "Number of Secondary Protons", 10, -0.5, 9.5);
176
177 // bins = 10**numpy.arange(-2,2.42,0.03)
178 std::array<double, 148> bins = {
179 1.00000000e-02, 1.07151931e-02, 1.14815362e-02, 1.23026877e-02, 1.31825674e-02, 1.41253754e-02, 1.51356125e-02,
180 1.62181010e-02, 1.73780083e-02, 1.86208714e-02, 1.99526231e-02, 2.13796209e-02, 2.29086765e-02, 2.45470892e-02,
181 2.63026799e-02, 2.81838293e-02, 3.01995172e-02, 3.23593657e-02, 3.46736850e-02, 3.71535229e-02, 3.98107171e-02,
182 4.26579519e-02, 4.57088190e-02, 4.89778819e-02, 5.24807460e-02, 5.62341325e-02, 6.02559586e-02, 6.45654229e-02,
183 6.91830971e-02, 7.41310241e-02, 7.94328235e-02, 8.51138038e-02, 9.12010839e-02, 9.77237221e-02, 1.04712855e-01,
184 1.12201845e-01, 1.20226443e-01, 1.28824955e-01, 1.38038426e-01, 1.47910839e-01, 1.58489319e-01, 1.69824365e-01,
185 1.81970086e-01, 1.94984460e-01, 2.08929613e-01, 2.23872114e-01, 2.39883292e-01, 2.57039578e-01, 2.75422870e-01,
186 2.95120923e-01, 3.16227766e-01, 3.38844156e-01, 3.63078055e-01, 3.89045145e-01, 4.16869383e-01, 4.46683592e-01,
187 4.78630092e-01, 5.12861384e-01, 5.49540874e-01, 5.88843655e-01, 6.30957344e-01, 6.76082975e-01, 7.24435960e-01,
188 7.76247117e-01, 8.31763771e-01, 8.91250938e-01, 9.54992586e-01, 1.02329299e+00, 1.09647820e+00, 1.17489755e+00,
189 1.25892541e+00, 1.34896288e+00, 1.44543977e+00, 1.54881662e+00, 1.65958691e+00, 1.77827941e+00, 1.90546072e+00,
190 2.04173794e+00, 2.18776162e+00, 2.34422882e+00, 2.51188643e+00, 2.69153480e+00, 2.88403150e+00, 3.09029543e+00,
191 3.31131121e+00, 3.54813389e+00, 3.80189396e+00, 4.07380278e+00, 4.36515832e+00, 4.67735141e+00, 5.01187234e+00,
192 5.37031796e+00, 5.75439937e+00, 6.16595002e+00, 6.60693448e+00, 7.07945784e+00, 7.58577575e+00, 8.12830516e+00,
193 8.70963590e+00, 9.33254301e+00, 1.00000000e+01, 1.07151931e+01, 1.14815362e+01, 1.23026877e+01, 1.31825674e+01,
194 1.41253754e+01, 1.51356125e+01, 1.62181010e+01, 1.73780083e+01, 1.86208714e+01, 1.99526231e+01, 2.13796209e+01,
195 2.29086765e+01, 2.45470892e+01, 2.63026799e+01, 2.81838293e+01, 3.01995172e+01, 3.23593657e+01, 3.46736850e+01,
196 3.71535229e+01, 3.98107171e+01, 4.26579519e+01, 4.57088190e+01, 4.89778819e+01, 5.24807460e+01, 5.62341325e+01,
197 6.02559586e+01, 6.45654229e+01, 6.91830971e+01, 7.41310241e+01, 7.94328235e+01, 8.51138038e+01, 9.12010839e+01,
198 9.77237221e+01, 1.04712855e+02, 1.12201845e+02, 1.20226443e+02, 1.28824955e+02, 1.38038426e+02, 1.47910839e+02,
199 1.58489319e+02, 1.69824365e+02, 1.81970086e+02, 1.94984460e+02, 2.08929613e+02, 2.23872114e+02, 2.39883292e+02,
200 2.57039578e+02,
201 };
202
203 fhElossVSLight = new TH2D("fhElossVSLight",
204 "Deposited Energy vs Generated Light",
205 bins.size() - 1,
206 bins.data(),
207 bins.size() - 1,
208 bins.data());
209 fhElossVSLightLog = new TH2D("fhElossVSLightLog", "Deposited Energy vs Generated Light", 600, -3, 3, 600, -3, 3);
210
211 fhThetaLight = new TH2D("fhThetaLight", "fhThetaLight", 200, -100, 100, 400, 0, 400);
212
213 if (fIs3DTrackEnabled)
214 {
215 // XYZ -> ZXY (side view)
216
217 fh3 = new TH3D("hMCTracks", "hMCTracks", 60, 1400, 1700, 50, -125, 125, 50, -125, 125);
218 fh3->SetTitle("NeuLAND MCTracks");
219 fh3->GetXaxis()->SetTitle("Z");
220 fh3->GetYaxis()->SetTitle("X");
221 fh3->GetZaxis()->SetTitle("Y");
222 FairRootManager::Instance()->Register("NeulandMCMon", "MC Tracks in NeuLAND", fh3, kTRUE);
223
224 fh3PNIP = new TH3D("h3DNPNIPS", "h3DNPNIPS", 60, 1400, 1700, 50, -125, 125, 50, -125, 125);
225 fh3PNIP->SetTitle("NeuLAND 3D NPNIPS");
226 fh3PNIP->GetXaxis()->SetTitle("Z");
227 fh3PNIP->GetYaxis()->SetTitle("X");
228 fh3PNIP->GetZaxis()->SetTitle("Y");
229 FairRootManager::Instance()->Register("Neuland3DPNIP", "First interactions in NeuLAND", fh3PNIP, kTRUE);
230 }
231
232 return kSUCCESS;
233}
234
236{
237 nEvents++;
238
239 const auto& npnips = fPrimaryNeutronInteractionPoints.get();
240 const auto& mcTracks = fMCTracks.Retrieve();
241 const auto& points = fNeulandPoints.get();
242
243 for (const auto& mcTrack : mcTracks)
244 {
245 // Distribution of MC Track mother id's
246 fhMotherIDs->Fill(mcTrack->GetMotherId());
247
248 // Energy of primary Particles
249 if (mcTrack->GetMotherId() == -1)
250 {
251 fhEPrimarys->Fill(GetKineticEnergy(mcTrack));
252 fhEPrimarys2->Fill(mcTrack->GetEnergy() * 1000.);
253 }
254
255 // Energy of primary Neutrons
256 if (IsPrimaryNeutron(mcTrack))
257 {
258 fhEPrimaryNeutrons->Fill(GetKineticEnergy(mcTrack));
259 }
260
261 // Remaining energy of neutrons after first interaction
262 if (IsMotherPrimaryNeutron(mcTrack, mcTracks) && mcTrack->GetPdgCode() == 2112)
263 {
264 fhESecondaryNeutrons->Fill(GetKineticEnergy(mcTrack));
265 }
266 }
267
268 fhnNPNIPs->Fill(npnips.size());
269 for (const auto& npnip : npnips)
270 {
271 // WIP: ToF Calculation -> Should respect other origin than 0,0,0,0.
272 const Double_t s2 = std::pow(npnip.GetX(), 2) + std::pow(npnip.GetY(), 2) + std::pow(npnip.GetZ(), 2); // cm²
273 const Double_t v2 = s2 / std::pow(npnip.GetTime(), 2); // ns²
274
275 const Double_t c2 = 898.75517873681758374898; // cm²/ns²
276 const Double_t massNeutron = 939.565379; // MeV/c²
277 const Double_t ETimeOfFlight = massNeutron * ((1. / std::sqrt(1 - (v2 / c2))) - 1);
278
279 auto mcTrack = mcTracks.at(npnip.GetTrackID());
280 fhNPNIPsEToFVSTime->Fill(ETimeOfFlight, npnip.GetTime());
281 fhMCToF->Fill(GetKineticEnergy(mcTrack) - ETimeOfFlight);
282 fhNPNIPSrvsz->Fill(npnip.GetZ(), std::sqrt(std::pow(npnip.GetX(), 2) + std::pow(npnip.GetY(), 2)));
283 fhNPNIPSxy->Fill(npnip.GetX(), npnip.GetY());
284 }
285
286 {
287 std::map<Int_t, Double_t> EtotPDG;
288 Double_t Etot = 0.;
289 Double_t EtotPrim = 0.;
290
291 for (const auto& point : points)
292 {
293 const R3BMCTrack* mcTrack = mcTracks.at(point.GetTrackID());
294
295 Etot += point.GetLightYield() * 1000.;
296
297 // Select tracks with a primary neutron mother
298 if (IsMotherPrimaryNeutron(mcTrack, mcTracks))
299 {
300
301 // Total energy of non-neutron secondary particles where mother is a primary neutron
302 if (mcTrack->GetPdgCode() != 2112)
303 {
304 EtotPrim += point.GetLightYield() * 1000.;
305 }
306
307 // Distribution of secondary particles
308 fhPDG->Fill(mcTrack->GetPdgCode());
309 fhPrimaryDaughterIDs->Fill(point.GetTrackID());
310
311 // Buld Histograms for each particle PDG if it doesn't exist
312 if (!fhmEPdg[mcTrack->GetPdgCode()])
313 {
314 TString name = TString("Light Yield of PID ") + TString::Itoa(mcTrack->GetPdgCode(), 10) +
315 TString(" with a primary neutron mother");
316 fhmEPdg[mcTrack->GetPdgCode()] =
317 new TH1D("hE_PDG_" + TString::Itoa(mcTrack->GetPdgCode(), 10), name, 3000, 0, 3000);
318 }
319 // Get Energy py particle where the mother is a primary neutron
320 fhmEPdg[mcTrack->GetPdgCode()]->Fill(point.GetLightYield() * 1000.); // point.GetEnergyLoss()*1000.);
321 } // end primary neutron mother
322
323 // Sum energy per particle type per event
324 if (!EtotPDG[mcTrack->GetPdgCode()])
325 {
326 EtotPDG[mcTrack->GetPdgCode()] = 0.;
327 }
328 EtotPDG[mcTrack->GetPdgCode()] += point.GetLightYield() * 1000.; // point.GetEnergyLoss()*1000.;
329
330 fhElossVSLight->Fill(point.GetEnergyLoss() * 1000., point.GetLightYield() * 1000.);
331 fhElossVSLightLog->Fill(std::log10(point.GetEnergyLoss() * 1000.),
332 std::log10(point.GetLightYield() * 1000.));
333
334 if (!fhmElossVSLightLogPdg[mcTrack->GetPdgCode()])
335 {
336 fhmElossVSLightLogPdg[mcTrack->GetPdgCode()] =
337 new TH2D("fhElossVSLightLog_" + TString::Itoa(mcTrack->GetPdgCode(), 10),
338 "Deposited Energy vs Generated Light",
339 600,
340 -3,
341 3,
342 600,
343 -3,
344 3);
345 }
346 fhmElossVSLightLogPdg[mcTrack->GetPdgCode()]->Fill(std::log10(point.GetEnergyLoss() * 1000.),
347 std::log10(point.GetLightYield() * 1000.));
348
349 fhThetaLight->Fill(GetTheta(mcTrack), point.GetLightYield() * 1000.);
350 }
351
352 fhEtot->Fill(Etot);
353 fhEtotPrim->Fill(EtotPrim);
354
355 for (const auto& kv : EtotPDG)
356 {
357 if (!fhmEtotPdg[kv.first])
358 {
359 TString name = TString("Sum Light Yield of PID ") + TString::Itoa(kv.first, 10);
360 fhmEtotPdg[kv.first] = new TH1D("fhEtotPDG_" + TString::Itoa(kv.first, 10), name, 3000, 0, 3000);
361 }
362 fhmEtotPdg[kv.first]->Fill(kv.second);
363
364 if (!fhmEtotPdgRel[kv.first])
365 {
366 TString name = TString("Percent Light Yield of PID ") + TString::Itoa(kv.first, 10) +
367 TString(" to total Light Yield.");
368 fhmEtotPdgRel[kv.first] = new TH1D("fhEtotPDGRel_" + TString::Itoa(kv.first, 10), name, 110, 0, 110);
369 }
370 if (Etot == 0)
371 {
372 fhmEtotPdgRel[kv.first]->Fill(0);
373 }
374 else
375 {
376 fhmEtotPdgRel[kv.first]->Fill(kv.second / Etot * 100.);
377 }
378 }
379 }
380
381 {
382 std::vector<int> primaryNeutronIDs;
383 std::map<int, std::vector<R3BMCTrack*>> tracksByPrimaryNeutronID;
384
385 const Int_t nTracks = mcTracks.size();
386 for (Int_t j = 0; j < nTracks; j++)
387 {
388 R3BMCTrack* track = mcTracks.at(j);
389 if (track->GetMotherId() == -1 && track->GetPdgCode() == 2112)
390 {
391 primaryNeutronIDs.push_back(j);
392 }
393 }
394
395 for (Int_t j = 0; j < nTracks; j++)
396 {
397 R3BMCTrack* track = mcTracks.at(j);
398 for (const Int_t primaryNeutronID : primaryNeutronIDs)
399 {
400 if (track->GetMotherId() == primaryNeutronID)
401 {
402 tracksByPrimaryNeutronID[primaryNeutronID].push_back(track);
403 }
404 }
405 }
406
407 for (const Int_t primaryNeutronID : primaryNeutronIDs)
408 {
409 auto& tracks = tracksByPrimaryNeutronID[primaryNeutronID];
410 fhnProducts->Fill(tracks.size());
411
412 fhnProductsCharged->Fill(std::accumulate(tracks.begin(),
413 tracks.end(),
414 0,
415 [](const int a, const R3BMCTrack* b)
416 {
417 if (b->GetPdgCode() == 2112 || b->GetPdgCode() == 22 ||
418 b->GetPdgCode() == 111)
419 {
420 return a;
421 }
422 else
423 {
424 return a + 1;
425 }
426 }));
427
428 fhnSecondaryNeutrons->Fill(std::accumulate(tracks.begin(),
429 tracks.end(),
430 0,
431 [](const int c, const R3BMCTrack* t)
432 {
433 if (t->GetPdgCode() == 2112)
434 {
435 return c + 1;
436 }
437 return c;
438 }));
439
440 fhnSecondaryProtons->Fill(std::accumulate(tracks.begin(),
441 tracks.end(),
442 0,
443 [](const int c, const R3BMCTrack* t)
444 {
445 if (t->GetPdgCode() == 2212)
446 {
447 return c + 1;
448 }
449 return c;
450 }));
451
452 const Double_t sumProductEnergy = std::accumulate(tracks.begin(),
453 tracks.end(),
454 0.,
455 [](const Double_t a, const R3BMCTrack* b)
456 { return a + 1000. * (b->GetEnergy() - b->GetMass()); });
457
458 fhSumProductEnergy->Fill(sumProductEnergy);
459
460 if (fIsFullSimAnaEnabled)
461 {
462 std::map<int, int> mCountByProductPdg;
463
464 for (const auto& track : tracks)
465 {
466 if (mCountByProductPdg.find(track->GetPdgCode()) == mCountByProductPdg.end())
467 {
468 mCountByProductPdg[track->GetPdgCode()] = 1;
469 }
470 else
471 {
472 mCountByProductPdg[track->GetPdgCode()]++;
473 }
474
475 auto& hEnergy = fhmEnergyByProductPdg[track->GetPdgCode()];
476 if (hEnergy == nullptr)
477 {
478 hEnergy = new TH1D(TString::Itoa(track->GetPdgCode(), 10),
479 TString::Itoa(track->GetPdgCode(), 10),
480 1000,
481 0,
482 1000);
483 }
484 hEnergy->Fill(1000. * (track->GetEnergy() - track->GetMass()));
485
486 const auto l = lookup(track->GetPdgCode());
487 auto& hEnergyRed = fhmEnergyByProductPdgReduced[l];
488 if (hEnergyRed == nullptr)
489 {
490 hEnergyRed = new TH1D(l.c_str(), l.c_str(), 1000, 0, 1000);
491 }
492 hEnergyRed->Fill(1000. * (track->GetEnergy() - track->GetMass()));
493 }
494
495 for (const auto& kv : mCountByProductPdg)
496 {
497 const int pdg = kv.first;
498 auto& hCount = fhmCountByProductPdg[pdg];
499 if (hCount == nullptr)
500 {
501 hCount = new TH1D(TString::Itoa(pdg, 10), TString::Itoa(pdg, 10), 10, 0, 10);
502 }
503 hCount->Fill(kv.second);
504
505 const auto l = lookup(kv.first);
506 auto& hCountRed = fhmCountByProductPdgReduced[l];
507 if (hCountRed == nullptr)
508 {
509 hCountRed = new TH1D(l.c_str(), l.c_str(), 10, 0, 10);
510 }
511 hCountRed->Fill(kv.second);
512 }
513
514 std::sort(tracks.begin(),
515 tracks.end(),
516 [](const R3BMCTrack* a, const R3BMCTrack* b) { return a->GetPdgCode() < b->GetPdgCode(); });
517
518 // std::vector<int> filteredPdgCodes{ 22, -211, 211, 111 };
519 std::vector<int> filteredPdgCodes{ 22, 111 };
520 const TString reaction = std::accumulate(
521 tracks.begin(),
522 tracks.end(),
523 TString(),
524 [&](TString s, const R3BMCTrack* b) -> TString
525 {
526 if (std::find(filteredPdgCodes.begin(), filteredPdgCodes.end(), b->GetPdgCode()) !=
527 filteredPdgCodes.end())
528 {
529 return s;
530 }
531 return s += "_" + TString::Itoa(b->GetPdgCode(), 10);
532 });
533
534 if (fhmEnergyByReaction[reaction] == nullptr)
535 {
536 fhmEnergyByReaction[reaction] = new TH1D(reaction, reaction, 1000, 0, 1000);
537 }
538 fhmEnergyByReaction[reaction]->Fill(sumProductEnergy);
539
540 for (const auto& track : tracks)
541 {
542 auto& hEnergy = fhmmEnergyByReactionByProductPdg[reaction][track->GetPdgCode()];
543 if (hEnergy == nullptr)
544 {
545 hEnergy = new TH1D(TString::Itoa(track->GetPdgCode(), 10) + reaction,
546 TString::Itoa(track->GetPdgCode(), 10) + reaction,
547 1000,
548 0,
549 1000);
550 }
551 hEnergy->Fill(1000. * (track->GetEnergy() - track->GetMass()));
552 }
553 }
554 }
555 }
556
557 {
558 ROOT::Math::PxPyPzEVector p4;
559 double m0 = 0;
560 for (const auto track : mcTracks)
561 {
562 if (track->GetMotherId() != -1)
563 {
564 break;
565 }
566 p4 += track->GetFourMomentum();
567 m0 += track->GetMass();
568 }
569 fhErelMC->Fill((p4.M() - m0) * 1e6);
570 }
571
572 if (fIs3DTrackEnabled)
573 {
574 fh3->Reset("ICES");
575 for (const auto& point : points)
576 {
577 if (point.GetLightYield() > 0)
578 {
579 const auto pos = point.GetPosition();
580 fh3->Fill(pos.Z(), pos.X(), pos.Y(), point.GetEnergyLoss() * 1000.);
581 }
582 }
583
584 fh3PNIP->Reset("ICES");
585 for (const auto& npnip : npnips)
586 {
587 fh3PNIP->Fill(npnip.GetZ(), npnip.GetX(), npnip.GetY(), npnip.GetTime());
588 }
589 }
590}
591
593{
594 TDirectory* tmp = gDirectory;
595 FairRootManager::Instance()->GetOutFile()->cd();
596
597 gDirectory->mkdir("NeulandMCMon");
598 gDirectory->cd("NeulandMCMon");
599
600 fhPDG->Write();
601 fhEPrimarys->Write();
602 fhEPrimarys2->Write();
603 fhEPrimaryNeutrons->Write();
604 fhErelMC->Write();
605 fhEtot->Write();
606 fhEtotPrim->Write();
607 fhESecondaryNeutrons->Write();
608 fhMotherIDs->Write();
609 fhPrimaryDaughterIDs->Write();
610 fhMCToF->Write();
611 fhNPNIPsEToFVSTime->Write();
612 fhNPNIPSrvsz->Write();
613 fhNPNIPSxy->Write();
614 fhnNPNIPs->Write();
615 fhThetaLight->Write();
616 fhElossVSLight->Write();
617 fhElossVSLightLog->Write();
618
619 gDirectory = tmp;
620 gDirectory->cd("NeulandMCMon");
621 gDirectory->mkdir("LightYield");
622 gDirectory->cd("LightYield");
623
624 gDirectory->mkdir("LightYieldByProductPdg");
625 gDirectory->cd("LightYieldByProductPdg");
626 for (const auto& kv : fhmEPdg)
627 {
628 kv.second->Write();
629 }
630 for (const auto& kv : fhmElossVSLightLogPdg)
631 {
632 kv.second->Write();
633 }
634 gDirectory->cd("..");
635
636 gDirectory->mkdir("SumLightYieldByProductPdg");
637 gDirectory->cd("SumLightYieldByProductPdg");
638 for (const auto& kv : fhmEtotPdg)
639 {
640 kv.second->Write();
641 }
642 gDirectory->cd("..");
643
644 gDirectory->mkdir("RelSumLightYieldByProductPdg");
645 gDirectory->cd("RelSumLightYieldByProductPdg");
646 for (const auto& kv : fhmEtotPdgRel)
647 {
648 kv.second->Write();
649 }
650 gDirectory->cd("..");
651
652 gDirectory = tmp;
653 gDirectory->cd("NeulandMCMon");
654 gDirectory->mkdir("Products");
655 gDirectory->cd("Products");
656
657 fhnProducts->Write();
658 fhnProductsCharged->Write();
659 fhSumProductEnergy->Write();
660 fhnSecondaryProtons->Write();
661 fhnSecondaryNeutrons->Write();
662
663 writeout(fhmEnergyByProductPdg, "EnergyByProduct");
664 writeout(fhmEnergyByProductPdgReduced, "EnergyByProductReduced");
665 writeout(fhmCountByProductPdg, "CountByProduct", nEvents);
666 writeout(fhmCountByProductPdgReduced, "CountByProductReduced", nEvents);
667 writeout(fhmEnergyByReaction, "EnergyByReaction");
668
669 if (!fhmmEnergyByReactionByProductPdg.empty())
670 {
671 gDirectory->mkdir("EnergyByReactionByProductPdg");
672 gDirectory->cd("EnergyByReactionByProductPdg");
673 for (const auto& kv : fhmmEnergyByReactionByProductPdg)
674 {
675 gDirectory->mkdir(kv.first);
676 gDirectory->cd(kv.first);
677 for (const auto& kv2 : kv.second)
678 {
679 kv2.second->Write();
680 }
681 gDirectory->cd("..");
682 }
683 gDirectory->cd("..");
684 }
685
686 gDirectory = tmp;
687}
688
ClassImp(R3B::Neuland::Cal2HitPar)
void writeout(std::map< T, TH1D * > &map, const std::string &what, const int nEvents=0)
Bool_t IsMotherPrimaryNeutron(const R3BMCTrack *mcTrack, const std::vector< R3BMCTrack * > &tracks)
Double_t GetKineticEnergy(const R3BMCTrack *mcTrack)
Double_t GetTheta(const R3BMCTrack *mcTrack)
Bool_t IsPrimaryNeutron(const R3BMCTrack *mcTrack)
std::string lookup(int pid)
static std::map< int, std::string > lookup_table
static const Double_t c
static const Double_t c2
static const Double_t massNeutron
int GetPdgCode() const
Accessors.
Definition R3BMCTrack.h:48
double GetEnergy() const
Definition R3BMCTrack.h:63
double GetPz() const
Definition R3BMCTrack.h:59
double GetP() const
Definition R3BMCTrack.h:65
double GetMass() const
Definition R3BMCTrack.h:60
int GetMotherId() const
Definition R3BMCTrack.h:49
R3BNeulandMCMon(const Option_t *option="")
void Finish() override
void Exec(Option_t *) override
InitStatus Init() override
STL class.
constexpr UInt_t nEvents