R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandNeutronReconstructionMon.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
16#include <FairTask.h>
17#include <Rtypes.h>
18#include <RtypesCore.h>
19#include <TH1.h>
20#include <TH2.h>
21#include <TLorentzVector.h>
22#include <TString.h>
23#include <TVector3.h>
24#include <algorithm>
25// #include <cmath>
26#include <cstddef>
27#include <fairlogger/Logger.h>
28#include <functional>
29// #include <iostream>
30// #include <numeric>
31#include "TClonesArray.h"
32#include "TDirectory.h"
33#include <TFile.h>
34#include <utility>
35#include <vector>
36// #include "FairMCPoint.h"
37#include "FairRootManager.h"
38#include "R3BMCTrack.h"
39#include "R3BNeulandNeutron.h"
40#include "R3BNeulandPoint.h"
41
42static const Double_t c2 = 898.75517873681758374; // cm²/ns²
43static const Double_t massNeutron = 939.565379; // MeV/c²
44
45namespace
46{
47 // auto Distance(const R3BNeulandNeutron& neutron, const FairMCPoint& mc_point) -> Double_t
48 // {
49 // auto pos_vec = ROOT::Math::XYZVector{ mc_point.GetX(), mc_point.GetY(), mc_point.GetZ() };
50 // pos_vec -= neutron.GetPosition();
51 // return std::sqrt(pos_vec.Dot(pos_vec));
52 // }
53
54 // auto Score(const std::vector<std::pair<R3BNeulandNeutron, FairMCPoint>>& combination) -> Double_t
55 // {
56 // return std::accumulate(combination.begin(),
57 // combination.end(),
58 // 0.,
59 // [](const Double_t sum, const std::pair<R3BNeulandNeutron, FairMCPoint>& pair)
60 // { return sum + Distance(pair.first, pair.second); });
61 // }
62 template <typename T, typename U>
63 auto GetAllCombinations(std::vector<T> one_series /* Copy! */,
64 std::vector<U> other_series /* Copy! */,
65 std::function<bool(const U&, const U&)> comparator)
66 -> std::vector<std::vector<std::pair<T, U>>>
67 {
68 std::vector<std::vector<std::pair<T, U>>> out;
69
70 // Bring both inputs up to the same length
71 if (one_series.size() < other_series.size())
72 {
73 one_series.resize(other_series.size());
74 }
75 if (other_series.size() < one_series.size())
76 {
77 other_series.resize(one_series.size());
78 }
79
80 std::sort(other_series.begin(), other_series.end(), comparator);
81 // NOLINTBEGIN(cppcoreguidelines-avoid-do-while)
82 do
83 {
84 std::vector<std::pair<T, U>> tmp;
85 for (size_t i = 0; i < one_series.size(); i++)
86 {
87 tmp.push_back({ one_series.at(i), other_series.at(i) });
88 }
89 out.push_back(std::move(tmp));
90 } while (std::next_permutation(other_series.begin(), other_series.end(), comparator));
91 // NOLINTEND(cppcoreguidelines-avoid-do-while)
92
93 return out;
94 }
95} // namespace
96
98 : FairTask("R3B Neuland Neutron Reconstruction Evaluation")
99 , fInput(input)
100 , fOutput(output)
101{
102 LOG(info) << "Using R3B Neuland Neutron Reconstruction Evaluation";
103}
104
106
108{
109 FairRootManager* ioman = FairRootManager::Instance();
110 if (!ioman)
111 {
112 LOG(fatal) << "R3BNeulandNeutronReconstructionMon::Init: No FairRootManager";
113 return kFATAL;
114 }
115
116 // Set Input: TClonesArray of R3BNeulandNeutrons
117 if (dynamic_cast<TClonesArray*>(ioman->GetObject(fInput)) == nullptr)
118 {
119 LOG(fatal) << "R3BNeulandNeutronReconstructionMon::Init No " << fInput << "!";
120 return kFATAL;
121 }
122 if (!TString((dynamic_cast<TClonesArray*>(ioman->GetObject(fInput)))->GetClass()->GetName())
123 .EqualTo("R3BNeulandNeutron"))
124 {
125 LOG(fatal) << "R3BNeulandNeutronReconstructionMon::Init Branch " << fInput
126 << " does not contain "
127 "R3BNeulandNeutrons!";
128 return kFATAL;
129 }
130 fReconstructedNeutrons = dynamic_cast<TClonesArray*>(ioman->GetObject(fInput));
131
132 // Set Input: TClonesArray of FairMCPoints
133 if (dynamic_cast<TClonesArray*>(ioman->GetObject("NeulandPrimaryPoints")) == nullptr)
134 {
135 LOG(fatal) << "R3BNeulandNeutronReconstructionMon::Init No NeulandPrimaryPoints!";
136 return kFATAL;
137 }
138 if (!TString((dynamic_cast<TClonesArray*>(ioman->GetObject("NeulandPrimaryPoints")))->GetClass()->GetName())
139 .EqualTo("R3BNeulandPoint"))
140 {
141 LOG(fatal) << "R3BNeulandNeutronReconstructionMon::Init Branch NeulandPrimaryPoints "
142 "does not contain R3BNeulandPoint!";
143 return kFATAL;
144 }
145 fPrimaryNeutronInteractionPoints = dynamic_cast<TClonesArray*>(ioman->GetObject("NeulandPrimaryPoints"));
146
147 // Set Input: TClonesArray of R3BMCTrack
148 if (dynamic_cast<TClonesArray*>(ioman->GetObject("MCTrack")) == nullptr)
149 {
150 LOG(fatal) << "R3BNeulandNeutronReconstructionMon::Init No MCTrack!";
151 return kFATAL;
152 }
153 if (!TString((dynamic_cast<TClonesArray*>(ioman->GetObject("MCTrack")))->GetClass()->GetName())
154 .EqualTo("R3BMCTrack"))
155 {
156 LOG(fatal) << "R3BNeulandNeutronReconstructionMon::Init Branch MCTrack "
157 "does not contain FairMCPoints!";
158 return kFATAL;
159 }
160 fMCTracks = dynamic_cast<TClonesArray*>(ioman->GetObject("MCTrack"));
161
162 TH1::AddDirectory(kFALSE);
163
164 fhScore = new TH1D("fhScore", "Neuland Neutron Reconstruction Score (lower is better)", 5000, 0, 5000);
165 fhCountN = new TH1D("fhCountN", "Number of reconstructed Neutrons", 11, -0.5, 10.5);
167 new TH1D("fhCountNdiff", "Number of reacted primary Neutrons - Number of reconstructed Neutrons", 11, -5, 5);
168 fhEdiff = new TH1D("fhEdiff", "Energy of primary Neutron - Energy of reconstructed Neutron", 2001, -1000, 1000);
169
170 fhErel = new TH1D("fhErel", "fhErel", 5000, 0, 5000);
171 fhErelMC = new TH1D("fhErelMC", "fhErelMC", 5000, 0, 5000);
172 fhErelVSnNreco = new TH2D("fhErelVSnNreco", "fhErelVSnNreco", 5000, 0, 5000, 11, -0.5, 10.5);
173
174 fhErelVSnNrecoNPNIPs = new TH2D("fhErelVSnNrecoNPNIPs", "fhErelVSnNrecoNPNIPs", 5000, 0, 5000, 11, -0.5, 10.5);
175
176 fhNreacNreco = new TH2D("fhNreacNreco", "fhNreacNreco", 11, -0.5, 10.5, 11, -0.5, 10.5);
177 fhNreacNreco->GetXaxis()->SetTitle("nReac");
178 fhNreacNreco->GetYaxis()->SetTitle("nReco");
179
180 return kSUCCESS;
181}
182
184{
185 const Int_t nReconstructedNeutrons = fReconstructedNeutrons->GetEntries();
186 std::vector<R3BNeulandNeutron> neutrons;
187 neutrons.reserve(nReconstructedNeutrons);
188 for (Int_t i = 0; i < nReconstructedNeutrons; i++)
189 {
190 neutrons.push_back(*(dynamic_cast<R3BNeulandNeutron*>(fReconstructedNeutrons->At(i))));
191 }
192
193 const Int_t nNPNIPS = fPrimaryNeutronInteractionPoints->GetEntries();
194 std::vector<R3BNeulandPoint> npnips;
195 npnips.reserve(nReconstructedNeutrons);
196 for (Int_t i = 0; i < nNPNIPS; i++)
197 {
198 npnips.push_back(*(dynamic_cast<R3BNeulandPoint*>(fPrimaryNeutronInteractionPoints->At(i))));
199 }
200 fhCountN->Fill(nReconstructedNeutrons);
201 fhCountNdiff->Fill((Int_t)nNPNIPS - (Int_t)nReconstructedNeutrons);
202 fhNreacNreco->Fill(nNPNIPS, nReconstructedNeutrons);
203
204 // Find the lowes score by brute force creating all combinations of MC-known and reconstructed neutron hits.
205 // Score is the sum of distances. Note that detecting a wrong amount of neutrons will create a huge score
206 // WiP. Slow because of all the copying
207 /*const auto allCombinations =
208 GetAllCombinations<R3BNeulandNeutron, FairMCPoint>(neutrons,
209 npnips,
210 [](const FairMCPoint& a, const FairMCPoint& b)
211 {
212 return a.GetTime() < b.GetTime();
213 });
214
215 const auto bestCombination = *std::min_element(allCombinations.begin(),
216 allCombinations.end(),
217 [](const std::vector<std::pair<R3BNeulandNeutron, FairMCPoint>>& a,
218 const std::vector<std::pair<R3BNeulandNeutron, FairMCPoint>>& b)
219 {
220 return Score(a) < Score(b);
221 });
222
223 const Double_t minScore = Score(bestCombination);
224 fhScore->Fill(minScore);
225
226 for (const auto& pair : bestCombination)
227 {
228 if (pair.second.GetEnergyLoss() > 0 || pair.first.GetEkin() > 0)
229 {
230 // TODO: Only for massNeutron
231 if (pair.second.GetEnergyLoss() > 0)
232 {
233 fhEdiff->Fill(pair.second.GetEnergyLoss() * 1000. - massNeutron - pair.first.GetEkin());
234 }
235 else
236 {
237 fhEdiff->Fill(-pair.first.GetEkin());
238 }
239 }
240 }*/
241
242 // Erel
243 // if (nReconstructedNeutrons == 0)
244 {
245 TLorentzVector p4_reco;
246 Double_t m0_reco = 0;
247
248 TLorentzVector p4_mc;
249 Double_t m0_mc = 0;
250
251 TLorentzVector p4_npnips;
252 Double_t m0_npnips = 0;
253
254 for (const auto& neutron : neutrons)
255 {
256 auto pos = neutron.GetP();
257 p4_reco += TLorentzVector(TVector3{ pos.X(), pos.Y(), pos.Z() }, neutron.GetEtot());
258 m0_reco += massNeutron;
259 }
260
261 for (const auto& npnip : npnips)
262 {
263 // Workaround mom = 0
264 R3BMCTrack* MCTrack = dynamic_cast<R3BMCTrack*>(fMCTracks->At(npnip.GetTrackID()));
265 p4_npnips += TLorentzVector(MCTrack->GetPx() * 1000.,
266 MCTrack->GetPy() * 1000.,
267 MCTrack->GetPz() * 1000.,
268 MCTrack->GetEnergy() * 1000.);
269 m0_npnips += MCTrack->GetMass() * 1000.;
270 }
271
272 R3BMCTrack* MCTrack;
273 const Int_t nMCTracks = fMCTracks->GetEntries();
274 for (Int_t i = 0; i < nMCTracks; i++)
275 {
276 MCTrack = dynamic_cast<R3BMCTrack*>(fMCTracks->At(i));
277 // We are only interested in primary particles, these should come first,
278 // thus we can stop once the Mother Id != -1
279 if (MCTrack->GetMotherId() != -1)
280 {
281 // break;
282 continue;
283 }
284
285 // For E_rel with reconstructed neutrons, only sum up non-neutron primary particles
286 // This information has to come from the other detectors later.
287 if (MCTrack->GetPdgCode() != 2112)
288 {
289 p4_reco += TLorentzVector(MCTrack->GetPx() * 1000.,
290 MCTrack->GetPy() * 1000.,
291 MCTrack->GetPz() * 1000.,
292 MCTrack->GetEnergy() * 1000.);
293 m0_reco += MCTrack->GetMass() * 1000.;
294
295 p4_npnips += TLorentzVector(MCTrack->GetPx() * 1000.,
296 MCTrack->GetPy() * 1000.,
297 MCTrack->GetPz() * 1000.,
298 MCTrack->GetEnergy() * 1000.);
299 m0_npnips += MCTrack->GetMass() * 1000.;
300 }
301
302 // For comparison, calculate E_rel for all primary particles from MC
303 p4_mc += TLorentzVector(MCTrack->GetPx() * 1000.,
304 MCTrack->GetPy() * 1000.,
305 MCTrack->GetPz() * 1000.,
306 MCTrack->GetEnergy() * 1000.);
307 m0_mc += MCTrack->GetMass() * 1000.;
308 }
309
310 const Double_t minv = p4_reco.Mag() - m0_reco;
311 fhErel->Fill(minv * 1000.);
312 fhErelVSnNreco->Fill(minv * 1000., nReconstructedNeutrons);
313
314 if (fhmErelnReco[nReconstructedNeutrons] == nullptr)
315 {
316 fhmErelnReco[nReconstructedNeutrons] = new TH1D("fhErel" + TString::Itoa(nReconstructedNeutrons, 10),
317 "fhErel" + TString::Itoa(nReconstructedNeutrons, 10),
318 5000,
319 0,
320 5000);
321 }
322 fhmErelnReco.at(nReconstructedNeutrons)->Fill(minv * 1000.);
323
324 const Double_t minv_mc = p4_mc.Mag() - m0_mc;
325 fhErelMC->Fill(minv_mc * 1000.);
326
327 const Double_t minv_npnips = p4_npnips.Mag() - m0_npnips;
328 fhErelVSnNrecoNPNIPs->Fill(minv_npnips * 1000., nNPNIPS);
329 }
330}
331
333{
334 TDirectory* tmp = gDirectory;
335 FairRootManager::Instance()->GetOutFile()->cd();
336
337 gDirectory->mkdir(fOutput);
338 gDirectory->cd(fOutput);
339
340 fhCountN->Write();
341 fhCountNdiff->Write();
342 fhScore->Write();
343 fhEdiff->Write();
344 fhErel->Write();
345 fhErelMC->Write();
346 fhErelVSnNreco->Write();
347 fhErelVSnNrecoNPNIPs->Write();
348 fhNreacNreco->Write();
349
350 for (auto& nh : fhmErelnReco)
351 {
352 nh.second->Write();
353 }
354
355 gDirectory = tmp;
356}
357
ClassImp(R3B::Neuland::Cal2HitPar)
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 GetPy() const
Definition R3BMCTrack.h:58
double GetMass() const
Definition R3BMCTrack.h:60
double GetPx() const
Definition R3BMCTrack.h:57
int GetMotherId() const
Definition R3BMCTrack.h:49
R3BNeulandNeutronReconstructionMon(const TString input="NeulandNeutrons", const TString output="NeulandNeutronReconstructionMon")