R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandPrimaryInteractionFinder.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
15#include "FairLogger.h"
16#include "TFile.h"
17#include <FairRootManager.h>
18#include <TH1D.h>
19#include <TH2D.h>
20#include <vector>
21
22namespace
23{
24 constexpr auto NEUTRON_PID = 2112;
25 constexpr auto MIN_TOF_DEFAULT = 1e99;
26
27 auto IsPrimaryTrack(const R3BMCTrack& track) -> bool
28 {
29 // TODO: The test can be modified to rely on some other information,
30 // e.g. if the neutrons were created in a reaction in the target simulated with Geant
31 return track.GetMotherId() == -1 && track.GetPdgCode() == NEUTRON_PID;
32 }
33
34 // Map a hit to every point if possible (else -1)
35 auto MapPointsToHits(const std::vector<R3BNeulandPoint>& points, const std::vector<R3BNeulandHit>& hits)
36 -> std::map<const R3BNeulandPoint*, const R3BNeulandHit*>
37 {
38 auto point_to_hit_map = std::map<const R3BNeulandPoint*, const R3BNeulandHit*>{};
39 for (const auto& point : points)
40 {
41 point_to_hit_map[&point] = nullptr;
42 for (const auto& hit : hits)
43 {
44 if (point.GetPaddle() == hit.GetPaddle())
45 {
46 point_to_hit_map[&point] = &hit;
47 }
48 }
49 }
50 return point_to_hit_map;
51 };
52
53 // Map a primary track to every point if possible (else -1)
54 auto MapPointsToPrimaryTracks(const std::vector<R3BNeulandPoint>& points, const std::vector<R3BMCTrack>& tracks)
55 -> std::map<const R3BNeulandPoint*, const R3BMCTrack*>
56 {
57 std::map<const R3BNeulandPoint*, const R3BMCTrack*> point_to_track_map;
58 for (const auto& point : points)
59 {
60 point_to_track_map[&point] = nullptr;
61 auto iTrack = point.GetTrackID();
62 while (iTrack > -1)
63 {
64 const auto& track = tracks.at(iTrack);
65 if (IsPrimaryTrack(track))
66 {
67 point_to_track_map[&point] = &track;
68 break;
69 }
70 // Else, start tracing back:
71 iTrack = track.GetMotherId();
72 }
73 }
74 return point_to_track_map;
75 }
76
77 auto FindFirstPoint(const R3BMCTrack& track,
78 const std::map<const R3BNeulandPoint*, const R3BMCTrack*>& point_to_track_map)
79 -> const R3BNeulandPoint*
80 {
81 const R3BNeulandPoint* minPoint{ nullptr };
82 // Search for min ToF:
83 auto minToF = MIN_TOF_DEFAULT;
84
85 for (auto [map_point, map_track] : point_to_track_map)
86 {
87 // Only look at points traced back to this track:
88 if (map_track == &track)
89 {
90 const auto ToF = map_point->GetTime();
91 if (ToF < minToF)
92 {
93 minToF = ToF;
94 minPoint = map_point;
95 }
96 }
97 }
98 return minPoint;
99 }
100
101 auto FindFirstHit(const R3BMCTrack& track,
102 const std::map<const R3BNeulandPoint*, const R3BMCTrack*>& point_to_track_map,
103 const std::map<const R3BNeulandPoint*, const R3BNeulandHit*>& point_to_hit_map)
104 -> const R3BNeulandHit*
105 {
106 const R3BNeulandPoint* minHitPoint = nullptr;
107 // Search for min ToF:
108 double minHitToF = MIN_TOF_DEFAULT;
109
110 for (const auto [map_point, map_track] : point_to_track_map)
111 {
112 // Only look at points traced back to this track:
113 // AND where a hit is registered
114 if (map_track == &track && point_to_hit_map.find(map_point) != point_to_hit_map.end())
115 {
116 const auto ToF = map_point->GetTime();
117 if (ToF < minHitToF)
118 {
119 minHitToF = ToF;
120 minHitPoint = map_point;
121 }
122 }
123 }
124 return minHitPoint != nullptr ? point_to_hit_map.at(minHitPoint) : nullptr;
125 }
126
127} // namespace
128
129// NOLINTNEXTLINE: bugprone-easily-swappable-parameters
131 std::string_view hitsIn,
132 std::string_view pointsOut,
133 std::string_view hitsOut,
134 std::string_view tracksOut)
135 : FairTask("R3BNeulandPrimaryInteractionFinder")
136 , fTracksIn("MCTrack")
137 , fPointsIn(pointsIn)
138 , fHitsIn(hitsIn)
139 , fTracksOut(tracksOut)
140 , fPointsOut(pointsOut)
141 , fHitsOut(hitsOut)
142 , fhDistance(new TH1D("fhDistance", "Distance firstPoint to firstHit", 10000, 0, 1000))
143 , fhPointsVsHits(new TH2D("fhPointsVsHits", "Number of Primary Points vs Number of Primary Hits", 6, 0, 6, 6, 0, 6))
144 , fhPointVsHitPaddle(
145 new TH2D("fhPointVsHitPaddle", "First Point Paddle vs First Hit Paddle", 3001, -1, 3000, 3001, -1, 3000))
146{
147 fhDistance->GetXaxis()->SetTitle("Distance [cm]");
148 fhPointsVsHits->GetXaxis()->SetTitle("# First Points");
149 fhPointsVsHits->GetYaxis()->SetTitle("# First Hits");
150 fhPointVsHitPaddle->GetXaxis()->SetTitle("First Point Paddle");
151 fhPointVsHitPaddle->GetYaxis()->SetTitle("First Hit Paddle");
152}
153
155{
156 fTracksIn.init();
157 fPointsIn.init();
158 fHitsIn.init();
159 fTracksOut.init();
160 fPointsOut.init();
161 fHitsOut.init();
162 return kSUCCESS;
163}
164
166{
167 const auto& tracks = fTracksIn.read();
168 const auto& points = fPointsIn.get();
169 const auto& hits = fHitsIn.get();
170 fTracksOut.clear();
171 fPointsOut.clear();
172 fHitsOut.clear();
173
174 const auto point_to_hit_map = MapPointsToHits(points, hits);
175 const auto point_to_track_map = MapPointsToPrimaryTracks(points, tracks);
176
177 if (FairLogger::GetLogger()->IsLogNeeded(fair::Severity::debug))
178 {
179 LOG(debug) << "R3BNeulandPrimaryInteractionFinder: Points without Hit in: ";
180 for (const auto [map_point, map_hit] : point_to_hit_map)
181 {
182 if (map_hit == nullptr)
183 {
184 const auto& point = map_point;
185 LOG(debug) << point->GetDetectorID() << ":" << tracks.at(point->GetTrackID()).GetPdgCode() << ":"
186 << point->GetLightYield() << ":" << point->GetEnergyLoss() << "\t";
187 }
188 }
189 }
190
191 for (const auto& track : tracks)
192 {
193 if (IsPrimaryTrack(track))
194 {
195 fTracksOut.get().push_back(track);
196
197 const auto* firstPoint = FindFirstPoint(track, point_to_track_map);
198 const auto* firstHit = FindFirstHit(track, point_to_track_map, point_to_hit_map);
199
200 if (firstPoint != nullptr)
201 {
202 fPointsOut.get().push_back(*firstPoint);
203 }
204
205 if (firstHit != nullptr)
206 {
207 // NOTE: This may add the same hit multiple times.
208 fHitsOut.get().push_back(*firstHit);
209 }
210
211 if ((firstHit != nullptr) && (firstPoint != nullptr))
212 {
213 fhDistance->Fill((firstPoint->GetPosition() - firstHit->GetPosition()).Mag());
214 }
215
216 fhPointVsHitPaddle->Fill((firstPoint != nullptr) ? firstPoint->GetPaddle() : -1,
217 (firstHit != nullptr) ? firstHit->GetPaddle() : -1);
218 }
219 }
220
221 fhPointsVsHits->Fill(static_cast<double>(fPointsOut.get().size()), static_cast<double>(fHitsOut.get().size()));
222}
223
225{
226 TDirectory* tmp = gDirectory;
227 FairRootManager::Instance()->GetOutFile()->cd();
228 gDirectory->mkdir("R3BNeulandPrimaryInteractionFinder");
229 gDirectory->cd("R3BNeulandPrimaryInteractionFinder");
230
231 fhDistance->Write();
232 fhPointsVsHits->Write();
233 fhPointVsHitPaddle->Write();
234
235 gDirectory = tmp;
236}
R3BNeulandPrimaryInteractionFinder(std::string_view pointsIn="NeulandPoints", std::string_view hitsIn="NeulandHits", std::string_view pointsOut="NeulandPrimaryPoints", std::string_view hitsOut="NeulandPrimaryHits", std::string_view tracksOut="NeulandPrimaryTracks")