CbmRoot
Loading...
Searching...
No Matches
CbmRichTrainAnnSelect.cxx
Go to the documentation of this file.
1/* Copyright (C) 2005-2017 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Semen Lebedev [committer] */
4
13
14#include "CbmDrawHist.h"
15#include "CbmMCTrack.h"
16#include "CbmRichConverter.h"
17#include "CbmRichHit.h"
18#include "CbmRichRing.h"
21#include "CbmTrackMatchNew.h"
22#include "FairRootManager.h"
23#include "TCanvas.h"
24#include "TClonesArray.h"
25#include "TH1D.h"
26#include "TMultiLayerPerceptron.h"
27
28#include <boost/assign/list_of.hpp>
29
30#include <iostream>
31#include <sstream>
32#include <vector>
33
34using boost::assign::list_of;
35using std::cout;
36using std::endl;
37using std::string;
38using std::vector;
39
41 : fRichRings(NULL)
42 , fMcTracks(NULL)
43 , fRichRingMatches(NULL)
44 ,
45
46 fEventNumber(0)
47 , fQuota(0.6)
48 , fMaxNofTrainSamples(5000)
49 , fNofFakeLikeTrue(0)
50 , fNofTrueLikeFake(0)
51 , fAnnCut(-0.5)
52 ,
53
54 fhNofHits()
55 , fhAngle()
56 , fhNofHitsOnRing()
57 , fhChi2()
58 , fhRadPos()
59 , fhRadius()
60 ,
61
62 fhAnnOutput()
63 , fhCumProb()
64 ,
65
66 fRSParams()
67 ,
68
69 fFitCOP(NULL)
70 , fSelectImpl(NULL)
71 ,
72
73 fHists()
74{
75 fhNofHits.resize(2);
76 fhAngle.resize(2);
77 fhNofHitsOnRing.resize(2);
78 fhChi2.resize(2);
79 fhRadPos.resize(2);
80 fhRadius.resize(2);
81 fhAnnOutput.resize(2);
82 fhCumProb.resize(2);
83 string ss;
84 for (int i = 0; i < 2; i++) {
85 if (i == 0) ss = "True";
86 if (i == 1) ss = "Fake";
87 // Difference Fake and True rings histograms BEGIN
88 fhNofHits[i] =
89 new TH1D(string("fhNofHits" + ss).c_str(), "Number of hits in ring;Nof hits in ring;Counter", 50, 0, 50);
90 fHists.push_back(fhNofHits[i]);
91 fhAngle[i] = new TH1D(string("fhAngle" + ss).c_str(), "Biggest angle in ring;Angle [rad];Counter", 50, 0, 6.5);
92 fHists.push_back(fhAngle[i]);
94 new TH1D(string("fhNofHitsOnRing" + ss).c_str(), "Number of hits on ring;Nof hits on ring;Counter", 50, 0, 50);
95 fHists.push_back(fhNofHitsOnRing[i]);
96 fhChi2[i] = new TH1D(string("fhFakeChi2" + ss).c_str(), "Chi2;Chi2;Counter", 100, 0., 1.0);
97 fHists.push_back(fhChi2[i]);
98 fhRadPos[i] =
99 new TH1D(string("fhRadPos" + ss).c_str(), "Radial position;Radial position [cm];Counter", 150, 0, 150);
100 fHists.push_back(fhRadPos[i]);
101 fhRadius[i] = new TH1D(string("fhRadius" + ss).c_str(), "Radius;Radius [cm];Counter", 80, 0., 9.);
102 fHists.push_back(fhRadius[i]);
103 // ANN outputs
104 fhAnnOutput[i] = new TH1D(string("fhAnnOutput" + ss).c_str(), "ANN output;ANN output;Counter", 100, -1.2, 1.2);
105 fHists.push_back(fhAnnOutput[i]);
106 fhCumProb[i] =
107 new TH1D(string("fhCumProb" + ss).c_str(), "ANN output;ANN output;Cumulative probability", 100, -1.2, 1.2);
108 fHists.push_back(fhCumProb[i]);
109 }
110
111 fRSParams.resize(2);
112}
113
115
117{
118 cout << "InitStatus CbmRichTrainAnnSelect::Init()" << endl;
119 FairRootManager* ioman = FairRootManager::Instance();
120 if (NULL == ioman) {
121 LOG(fatal) << GetName() << "::Init: CbmRootManager is not instantiated";
122 }
123
124 fRichRings = (TClonesArray*) ioman->GetObject("RichRing");
125 if (NULL == fRichRings) {
126 LOG(fatal) << GetName() << "::Init: No RichRing array!";
127 }
128
129 fMcTracks = (TClonesArray*) ioman->GetObject("MCTrack");
130 if (NULL == fMcTracks) {
131 LOG(fatal) << GetName() << "::Init: No MCTrack array!";
132 }
133
134 fRichRingMatches = (TClonesArray*) ioman->GetObject("RichRingMatch");
135 if (NULL == fRichRingMatches) {
136 LOG(fatal) << GetName() << "::Init: No RichRingMatch array!";
137 }
138
142
143 return kSUCCESS;
144}
145
146void CbmRichTrainAnnSelect::Exec(Option_t* /*option*/)
147{
148 fEventNumber++;
149 cout << "CbmRichRingQa Event No. " << fEventNumber << endl;
150
151 SetRecFlag();
153
154 cout << "nof trues = " << fRSParams[0].size() << endl;
155 cout << "nof fakes = " << fRSParams[1].size() << endl;
156}
157
159{
160 Int_t nMatches = fRichRingMatches->GetEntriesFast();
161 vector<Int_t> clone;
162 clone.clear();
163
164 for (Int_t iMatches = 0; iMatches < nMatches; iMatches++) {
165 CbmTrackMatchNew* match = (CbmTrackMatchNew*) fRichRingMatches->At(iMatches);
166 if (NULL == match) continue;
167
168 CbmRichRing* ring = (CbmRichRing*) fRichRings->At(iMatches);
169 if (NULL == ring) continue;
170
171 Double_t lPercTrue = match->GetTrueOverAllHitsRatio();
172
173 Int_t trackID = match->GetMatchedLink().GetIndex();
174 if (trackID > fMcTracks->GetEntriesFast() || trackID < 0) continue;
175 CbmMCTrack* track = (CbmMCTrack*) fMcTracks->At(trackID);
176 if (NULL == track) continue;
177
178 Int_t gcode = TMath::Abs(track->GetPdgCode());
179 //Double_t momentum = track->GetP();
180 Int_t motherId = track->GetMotherId();
181
182 ring->SetRecFlag(-1);
183
184 if (lPercTrue < fQuota) { // fake ring
185 ring->SetRecFlag(1);
186 }
187 else { // true rings
188 if (gcode == 11 && motherId == -1) { // select primary electrons
189 ring->SetRecFlag(3);
190 }
191 }
192 } // iMatches
193}
194
196{
197 Int_t nMatches = fRichRingMatches->GetEntriesFast();
198
199 for (Int_t iMatches = 0; iMatches < nMatches; iMatches++) {
200 CbmTrackMatchNew* match = (CbmTrackMatchNew*) fRichRingMatches->At(iMatches);
201 if (NULL == match) continue;
202
203 CbmRichRing* ring = (CbmRichRing*) fRichRings->At(iMatches);
204 if (NULL == ring) continue;
205
206 CbmRichRingLight ringLight;
208 fFitCOP->DoFit(&ringLight);
209
210 Int_t recFlag = ring->GetRecFlag();
211 Double_t angle = fSelectImpl->GetAngle(&ringLight);
212 Int_t hitsOnRing = fSelectImpl->GetNofHitsOnRingCircle(&ringLight);
213 Double_t chi2 = ringLight.GetChi2() / ringLight.GetNofHits();
214 Int_t nHits = ringLight.GetNofHits();
215 Double_t radPos = ringLight.GetRadialPosition();
216 Double_t radius = ringLight.GetRadius();
217
219 p.fAngle = angle;
220 p.fChi2 = chi2;
221 p.fHitsOnRing = hitsOnRing;
222 p.fNofHits = nHits;
223 p.fRadPos = radPos;
224 p.fRadius = radius;
225
226 if (recFlag == 1) {
227 fhNofHits[1]->Fill(nHits);
228 fhAngle[1]->Fill(angle);
229 fhNofHitsOnRing[1]->Fill(hitsOnRing);
230 fhRadPos[1]->Fill(radPos);
231 fhRadius[1]->Fill(radius);
232 fhChi2[1]->Fill(chi2);
233 fRSParams[1].push_back(p);
234 }
235
236 if (recFlag == 3) {
237 fhNofHits[0]->Fill(nHits);
238 fhAngle[0]->Fill(angle);
239 fhNofHitsOnRing[0]->Fill(hitsOnRing);
240 fhRadPos[0]->Fill(radPos);
241 fhRadius[0]->Fill(radius);
242 fhChi2[0]->Fill(chi2);
243 fRSParams[0].push_back(p);
244 }
245 } // iMatches
246}
247
249{
250 TTree* simu = new TTree("MonteCarlo", "MontecarloData");
251 Double_t x[6];
252 Double_t xOut;
253
254 simu->Branch("x0", &x[0], "x0/D");
255 simu->Branch("x1", &x[1], "x1/D");
256 simu->Branch("x2", &x[2], "x2/D");
257 simu->Branch("x3", &x[3], "x3/D");
258 simu->Branch("x4", &x[4], "x4/D");
259 simu->Branch("x5", &x[5], "x5/D");
260 simu->Branch("xOut", &xOut, "xOut/D");
261
262 for (int j = 0; j < 2; j++) {
263 for (unsigned int i = 0; i < fRSParams[j].size(); i++) {
264 x[0] = fRSParams[j][i].fNofHits / 45.;
265 x[1] = fRSParams[j][i].fAngle / 6.28;
266 x[2] = fRSParams[j][i].fHitsOnRing / 45.;
267 x[3] = fRSParams[j][i].fRadPos / 110.;
268 x[4] = fRSParams[j][i].fRadius / 10.;
269 x[5] = fRSParams[j][i].fChi2 / 0.4;
270
271 for (int k = 0; k < 6; k++) {
272 if (x[k] < 0.) x[k] = 0.;
273 if (x[k] > 1.) x[k] = 1.;
274 }
275
276 if (j == 0) xOut = 1.;
277 if (j == 1) xOut = -1.;
278 simu->Fill();
279 if (i >= fMaxNofTrainSamples) break;
280 }
281 }
282
283 TMultiLayerPerceptron network("x0,x1,x2,x3,x4,x5:10:xOut", simu, "Entry$+1");
284 //network.LoadWeights("/u/slebedev/JUL09/trunk/parameters/rich/NeuralNet_RingSelection_Weights_Compact.txt");
285 network.Train(300, "text,update=10");
286 network.DumpWeights("rich_v17a_select_ann_weights.txt");
287 //network.Export();
288
289 Double_t params[6];
290
293
294 for (int j = 0; j < 2; j++) {
295 for (unsigned int i = 0; i < fRSParams[j].size(); i++) {
296 params[0] = fRSParams[j][i].fNofHits / 45.;
297 params[1] = fRSParams[j][i].fAngle / 6.28;
298 params[2] = fRSParams[j][i].fHitsOnRing / 45.;
299 params[3] = fRSParams[j][i].fRadPos / 110.;
300 params[4] = fRSParams[j][i].fRadius / 10.;
301 params[5] = fRSParams[j][i].fChi2 / 0.4;
302
303 for (int k = 0; k < 6; k++) {
304 if (params[k] < 0.) params[k] = 0.;
305 if (params[k] > 1.) params[k] = 1.;
306 }
307
308 Double_t netEval = network.Evaluate(0, params);
309
310 //if (netEval > maxEval) netEval = maxEval - 0.01;
311 // if (netEval < minEval) netEval = minEval + 0.01;
312
313 fhAnnOutput[j]->Fill(netEval);
314 if (netEval >= fAnnCut && j == 1) fNofFakeLikeTrue++;
315 if (netEval < fAnnCut && j == 0) fNofTrueLikeFake++;
316 }
317 }
318}
319
321{
322 cout << "nof Trues = " << fRSParams[0].size() << endl;
323 cout << "nof Fakes = " << fRSParams[1].size() << endl;
324 cout << "Fake like True = " << fNofFakeLikeTrue
325 << ", fake remain = " << (double) fNofFakeLikeTrue / fRSParams[1].size() << endl;
326 cout << "True like Fake = " << fNofTrueLikeFake << ", true loss = " << (double) fNofTrueLikeFake / fRSParams[0].size()
327 << endl;
328 cout << "ANN cut = " << fAnnCut << endl;
329
330 Double_t cumProbFake = 0.;
331 Double_t cumProbTrue = 0.;
332 Int_t nofFake = (Int_t) fhAnnOutput[1]->GetEntries();
333 Int_t nofTrue = (Int_t) fhAnnOutput[0]->GetEntries();
334 for (Int_t i = 1; i <= fhAnnOutput[1]->GetNbinsX(); i++) {
335 cumProbFake += fhAnnOutput[1]->GetBinContent(i);
336 fhCumProb[1]->SetBinContent(i, (Double_t) cumProbFake / nofFake);
337
338 cumProbTrue += fhAnnOutput[0]->GetBinContent(i);
339 fhCumProb[0]->SetBinContent(i, 1. - (Double_t) cumProbTrue / nofTrue);
340 }
341
342 // TCanvas* c1 = new TCanvas("ann_select_ann_output", "ann_select_ann_output", 500, 500);
343 new TCanvas("ann_select_ann_output", "ann_select_ann_output", 500, 500);
344 fhAnnOutput[0]->Scale(1. / fhAnnOutput[0]->Integral());
345 fhAnnOutput[1]->Scale(1. / fhAnnOutput[1]->Integral());
346 DrawH1(list_of(fhAnnOutput[0])(fhAnnOutput[1]), list_of("True")("Fake"), kLinear, kLog, true, 0.8, 0.8, 0.99, 0.99);
347
348 // TCanvas* c2 = new TCanvas("ann_select_cum_prob", "ann_select_cum_prob", 500, 500);
349 new TCanvas("ann_select_cum_prob", "ann_select_cum_prob", 500, 500);
350 DrawH1(list_of(fhCumProb[0])(fhCumProb[1]), list_of("True")("Fake"), kLinear, kLinear, true, 0.8, 0.8, 0.99, 0.99);
351
352 TCanvas* c3 = new TCanvas("ann_select_params", "ann_select_params", 900, 600);
353 c3->Divide(3, 2);
354 c3->cd(1);
355 DrawH1(list_of(fhNofHits[0])(fhNofHits[1]), list_of("True")("Fake"), kLinear, kLog, true, 0.8, 0.8, 0.99, 0.99);
356 c3->cd(2);
357 DrawH1(list_of(fhAngle[0])(fhAngle[1]), list_of("True")("Fake"), kLinear, kLog, true, 0.8, 0.8, 0.99, 0.99);
358 c3->cd(3);
359 DrawH1(list_of(fhNofHitsOnRing[0])(fhNofHitsOnRing[1]), list_of("True")("Fake"), kLinear, kLog, true, 0.8, 0.8, 0.99,
360 0.99);
361 c3->cd(4);
362 DrawH1(list_of(fhRadPos[0])(fhRadPos[1]), list_of("True")("Fake"), kLinear, kLog, true, 0.8, 0.8, 0.99, 0.99);
363 c3->cd(5);
364 DrawH1(list_of(fhChi2[0])(fhChi2[1]), list_of("True")("Fake"), kLinear, kLog, true, 0.8, 0.8, 0.99, 0.99);
365 c3->cd(6);
366 DrawH1(list_of(fhRadius[0])(fhRadius[1]), list_of("True")("Fake"), kLinear, kLog, true, 0.8, 0.8, 0.99, 0.99);
367}
368
370{
372 Draw();
373 for (unsigned int i = 0; i < fHists.size(); i++) {
374 fHists[i]->Scale(1. / fHists[i]->Integral());
375 }
376
377 TDirectory* current = gDirectory;
378 TDirectory* rich = current->mkdir("CbmRichTrainAnnSelect");
379 rich->cd();
380
381 for (unsigned int i = 0; i < fHists.size(); i++) {
382 fHists[i]->Write();
383 }
384
385 rich->cd();
386 current->cd();
387
388 delete fFitCOP;
389 delete fSelectImpl;
390}
391
ClassImp(CbmConverterManager)
void DrawH1(TH1 *hist, HistScale logx, HistScale logy, const string &drawOpt, Int_t color, Int_t lineWidth, Int_t lineStyle, Int_t markerSize, Int_t markerStyle)
Helper functions for drawing 1D and 2D histograms and graphs.
@ kLinear
Definition CbmDrawHist.h:69
@ kLog
Definition CbmDrawHist.h:68
Convert internal data classes to cbmroot common data classes.
Here the ring is fitted with the COP algorithm from A. Ayriyan/G. Ososkov.
Train ANN for fake rejection.
static constexpr size_t size()
Definition KfSimdPseudo.h:2
int32_t GetMotherId() const
Definition CbmMCTrack.h:69
int32_t GetPdgCode() const
Definition CbmMCTrack.h:68
const CbmLink & GetMatchedLink() const
Definition CbmMatch.h:41
static void Init()
Initialize array of RICH hits.
static void CopyHitsToRingLight(const CbmRichRing *ring1, CbmRichRingLight *ring2)
Copy hits from CbmRichRing to CbmRichRingLight.
Here the ring is fitted with the COP algorithm from A. Ayriyan/G. Ososkov.
virtual void DoFit(CbmRichRingLight *ring)
Inherited from CbmRichRingFitterBase.
float GetRadius() const
int GetNofHits() const
Return number of hits in ring.
float GetChi2() const
float GetRadialPosition() const
Return radial position of the ring.
int GetNofHitsOnRingCircle(CbmRichRingLight *ring)
float GetAngle(CbmRichRingLight *ring)
int32_t GetRecFlag() const
Definition CbmRichRing.h:98
void SetRecFlag(int32_t recflag)
Definition CbmRichRing.h:62
Train ANN for fake rejection.
virtual void FinishTask()
Inherited from FairTask.
void TrainAndTestAnn()
Train and test ANN.
CbmRichRingFitterCOP * fFitCOP
vector< vector< RingSelectParam > > fRSParams
CbmRichTrainAnnSelect()
Default constructor.
CbmRichRingSelectImpl * fSelectImpl
void Draw(Option_t *="")
Draw results.
void SetRecFlag()
Set recFlag weather ring was found correctly or not.
virtual void Exec(Option_t *option)
Inherited from FairTask.
void DiffFakeTrueCircle()
Fill ring selection parameters in array and histograms.
virtual InitStatus Init()
Inherited from FairTask.
virtual ~CbmRichTrainAnnSelect()
Destructor.
double GetTrueOverAllHitsRatio() const
Input Parameters for ANN.