R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandMultiplicityCalorimetricTrain.cxx
Go to the documentation of this file.
2#include "FairLogger.h"
3#include "FairRootManager.h"
4#include "FairRuntimeDb.h"
5#include "Math/Factory.h"
6#include "Math/Functor.h"
7#include "Math/Minimizer.h"
8#include "R3BNeulandCluster.h"
10#include "TDirectory.h"
11#include <FairTask.h>
12#include <RtypesCore.h>
13#include <TCutG.h>
14#include <TH2.h>
15#include <TString.h>
16#include <cmath>
17#include <fairlogger/Logger.h>
18#include <iostream>
19#include <numeric>
20#include <string>
21#include <string_view>
22#include <utility>
23
24/*
25 * ^
26 * |
27 * n |
28 * c 3
29 * l X X
30 * u X X
31 * s X X
32 * t X X
33 * e 0 X
34 * r | X X
35 * |---1 X X X X 2----->
36 * edep
37 */
38
40 std::string_view tracks,
41 std::string_view phits)
42 : FairTask("R3BNeulandMultiplicityCalorimetricTrain")
43 , fClusters(clusters)
44 , fTracks(tracks)
45 , fPHits(std::move(phits))
46 , fPar(nullptr)
47 , fUseHits(false)
48 , fEdepOpt({ 200, 25, 50, 1500 })
49 , fEdepOffOpt({ 5, 1, 0, 250 })
50 , fNclusterOpt({ 10, 5, 5, 50 })
51 , fNclusterOffOpt({ 2, 1, 0, 10 })
52 , fWeight(0)
53{
54}
55
57{
58 for (auto& nc : fCuts)
59 {
60 delete nc.second;
61 }
62}
63
65{
66 // Input
67 fClusters.init();
68 fTracks.init();
69 fPHits.init(true);
70
71 // Output Parameter Container
72 auto* ioman = FairRootManager::Instance();
73 if (ioman == nullptr)
74 {
75 LOG(fatal) << "R3BNeulandMultiplicityCalorimetricTrain:Init: No FairRootManager";
76 return kFATAL;
77 }
78
79 auto* rtdb = FairRuntimeDb::instance();
80 if (rtdb == nullptr)
81 {
82 LOG(fatal) << "R3BNeulandMultiplicityCalorimetricTrain::Init: No FairRuntimeDb!";
83 return kFATAL;
84 }
86 rtdb->getContainer("R3BNeulandMultiplicityCalorimetricPar"));
87 if (fPar == nullptr)
88 {
89 LOG(fatal) << "R3BNeulandMultiplicityCalorimetricTrain::Init: No R3BNeulandMultiplicityCalorimetricPar!";
90 return kFATAL;
91 }
92 // FIXME: FairRuntimeDB needs to be forced to load the Data from the second file with Run Id 1
93 // rtdb->initContainers(rtdb->getCurrentRun()->getRunId(), 1);
94
95 return kSUCCESS;
96}
97
99{
100 const auto nPN = fUseHits ? fPHits.get().size() : fTracks.get().size();
101
102 const auto& clusters = fClusters.get();
103 const auto nClusters = clusters.size();
104
105 if (nClusters == 0)
106 {
107 return;
108 }
109
110 const auto Edep =
111 std::accumulate(clusters.cbegin(),
112 clusters.cend(),
113 double{},
114 [](double sum, const R3BNeulandCluster& cluster) { return sum + cluster.GetE(); });
115
116 GetOrBuildHist(nPN)->Fill(static_cast<int>(std::ceil(Edep)), static_cast<double>(nClusters));
117}
118
120{
121 Optimize();
122
123 if (FairLogger::GetLogger()->IsLogNeeded(fair::Severity::debug))
124 {
125 Print();
126 }
127
128 auto* rtdb = FairRuntimeDb::instance();
129 rtdb->addRun(1);
130 fPar->SetNeutronCuts(fCuts);
131 fPar->setChanged();
132 rtdb->writeContainer(fPar, rtdb->getRun(1), nullptr);
133
134 TDirectory* tmp = gDirectory;
135 FairRootManager::Instance()->GetOutFile()->cd();
136
137 gDirectory->mkdir("NeulandMultiplicityCalorimetricTrain");
138 gDirectory->cd("NeulandMultiplicityCalorimetricTrain");
139
140 for (auto& nh : fHists)
141 {
142 nh.second->Write();
143 }
144
145 gDirectory = tmp;
146}
147
149{
150 ROOT::Math::Minimizer* min = ROOT::Math::Factory::CreateMinimizer("Genetic");
151
152 const ROOT::Math::Functor minimizor_functor([&](const double* cut) { return WastedEfficiency(cut); }, 4);
153 min->SetFunction(minimizor_functor);
154
155 min->SetLimitedVariable(0, "edep", fEdepOpt.at(0), fEdepOpt.at(1), fEdepOpt.at(2), fEdepOpt.at(3));
156 min->SetLimitedVariable(1, "edepoff", fEdepOffOpt.at(0), fEdepOffOpt.at(1), fEdepOffOpt.at(2), fEdepOffOpt.at(3));
157 min->SetLimitedVariable(
158 2, "ncluster", fNclusterOpt.at(0), fNclusterOpt.at(1), fNclusterOpt.at(2), fNclusterOpt.at(3));
159 min->SetLimitedVariable(
160 3, "nclusteroff", fNclusterOffOpt.at(0), fNclusterOffOpt.at(1), fNclusterOffOpt.at(2), fNclusterOffOpt.at(3));
161
162 min->Minimize();
163
164 LOG(info) << "R3BNeulandMultiplicityCalorimetricTrain::Optimize done!";
165}
166
167TCutG* R3BNeulandMultiplicityCalorimetricTrain::GetCut(const unsigned int nNeutrons,
168 const double edep,
169 const double edepoff,
170 const double ncluster,
171 const double nclusteroff)
172{
173 if (!fCuts[nNeutrons])
174 {
175 fCuts[nNeutrons] = new TCutG(TString(std::to_string(nNeutrons)), 4);
176 fCuts.at(nNeutrons)->SetVarX("Total Energy [MeV]");
177 fCuts.at(nNeutrons)->SetVarY("Number of Clusters");
178 }
179 auto cut = fCuts[nNeutrons];
180
181 // nmin: p0 = p1 = (-1,-1)
182 if (nNeutrons == fHists.begin()->first)
183 {
184 cut->SetPoint(0, -1, -1);
185 cut->SetPoint(1, -1, -1);
186 }
187 else
188 {
189 cut->SetPoint(0, -1, ncluster * (nNeutrons - 1) + nclusteroff);
190 cut->SetPoint(1, edep * (nNeutrons - 1) + edepoff, -1);
191 }
192
193 // nmax: p2 (inf, -1) p3 (-1, inf)
194 if (nNeutrons == fHists.rbegin()->first)
195 {
196 cut->SetPoint(2, 100000, -1);
197 cut->SetPoint(3, -1, 100000);
198 }
199 else
200 {
201 cut->SetPoint(2, edep * nNeutrons + edepoff, -1);
202 cut->SetPoint(3, -1, ncluster * nNeutrons + nclusteroff);
203 }
204
205 return cut;
206}
207
209{
210 double edep = d[0];
211 double edepoff = d[1];
212 double ncluster = d[2];
213 double nclusteroff = d[3];
214
215 double wasted_efficiency = 0;
216 for (auto& nh : fHists)
217 {
218 const unsigned int nNeutrons = nh.first;
219 const double efficiency =
220 ((double)GetCut(nNeutrons, edep, edepoff, ncluster, nclusteroff)->IntegralHist(nh.second) /
221 (double)nh.second->GetEntries());
222 wasted_efficiency += (1. - efficiency) * (1. + fWeight * nNeutrons);
223 }
224 return wasted_efficiency;
225}
226
228{
229 std::cout << "\t";
230 for (const auto& nh : fHists)
231 {
232 std::cout << nh.first << "n\t";
233 }
234 std::cout << std::endl;
235
236 for (const auto& nc : fCuts)
237 {
238 const unsigned int nOut = nc.first;
239 const TCutG* cut = nc.second;
240
241 std::cout << nOut << "n:\t";
242 for (const auto& nh : fHists)
243 {
244 std::cout << ((double)cut->IntegralHist(nh.second) / (double)nh.second->GetEntries()) << "\t";
245 }
246 std::cout << std::endl;
247 }
248 std::cout << std::endl;
249}
250
252{
253 if (fHists.find(i) == fHists.end())
254 {
255 const TString name = "hnPN" + std::to_string(i);
256 fHists[i] = new TH2D(name, name, 150, 0, 5000, 50, 0, 100);
257 fHists.at(i)->GetXaxis()->SetTitle("Total Energy [MeV]");
258 fHists.at(i)->GetYaxis()->SetTitle("Number of Clusters");
259 }
260 return fHists.at(i);
261}
NeuLAND number of clusters / energy - neutron multiplicity parameter storage.
R3BNeulandMultiplicityCalorimetricTrain(std::string_view clusters="NeulandClusters", std::string_view tracks="NeulandPrimaryTracks", std::string_view phits="NeulandPrimaryHits")
auto GetCut(unsigned int nNeutrons, double edep, double edepoff, double ncluster, double nclusteroff) -> TCutG *
R3B::InputVectorConnector< R3BNeulandCluster > fClusters