CbmRoot
Loading...
Searching...
No Matches
CbmTrdElectronsTrainAnn.cxx
Go to the documentation of this file.
1/* Copyright (C) 2009-2021 UGiessen, JINR-LIT
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Semen Lebedev [committer] */
4
6
7#include "CbmDrawHist.h"
8#include "CbmMCTrack.h"
9#include "CbmTrackMatchNew.h"
10#include "CbmTrdHit.h"
11#include "CbmTrdPoint.h"
12#include "CbmTrdTrack.h"
13#include "TCanvas.h"
14#include "TClonesArray.h"
15#include "TCut.h"
16#include "TFile.h"
17#include "TGraph.h"
18#include "TLine.h"
19#include "TMVA/Config.h"
20#include "TMVA/PDF.h"
21#include "TMath.h"
22#include "TRandom.h"
23#include "TSystem.h"
24
25#include <FairRootManager.h>
26#include <Logger.h>
27
28#include <boost/assign/list_of.hpp>
29
30#include <iostream>
31#include <sstream>
32#include <string>
33#include <vector>
34
35using boost::assign::list_of;
36using std::cout;
37using std::endl;
38using std::string;
39using std::stringstream;
40using std::vector;
41
43 : FairTask()
44 , fMCTracks(NULL)
45 , fTrdPoints(NULL)
46 , fTrdTracks(NULL)
47 , fTrdTrackMatches(NULL)
48 , fTrdHits(NULL)
49 , fEloss()
50 , fHists()
51 , fhResults()
52 , fhMeanEloss()
53 , fhEloss()
54 , fhElossSort()
55 , fEventNum(0)
56 , fOutputDir("results/")
57 , fSigmaError(0.0)
58 , fIsDoTrain(true)
59 , fTransformType(1)
60 , fBeamDataFile("")
61 , fBeamDataPiHist("")
62 , fBeamDataElHist("")
63 , fAnnInput()
64 , fXOut(-1.)
65 , fNofTrdLayers(nofTrdLayers)
66 , fMaxEval(1.3)
67 , fMinEval(-1.3)
68 , fNN(NULL)
69 , fReader(NULL)
70 , fIdMethod(kANN)
71 , fNofAnnEpochs(50)
72 , fNofTrainSamples(2500)
73 , fElIdEfficiency(0.9)
74 , fhOutput()
75 , fhCumProbOutput()
76 , fhInput()
77{
78 fEloss.resize(2);
79 fhMeanEloss.resize(2);
80 fhEloss.resize(2);
81 string s;
82
83 TH1::SetDefaultBufferSize(30000);
84
85 fhResults = new TH1D("fhResults", "fhResults", 2, 0, 2);
86 fHists.push_back(fhResults);
87
88 for (int i = 0; i < 2; i++) {
89 if (i == 0) s = "El";
90 if (i == 1) s = "Pi";
91 fhMeanEloss[i] = new TH1D(("fhMeanEloss" + s).c_str(), "fhMeanEloss;Mean energy loss [a.u.];Yield", 100, 0., 50e-6);
92 fHists.push_back(fhMeanEloss[i]);
93 fhEloss[i] = new TH1D(("fhEloss" + s).c_str(), "fhEloss;Energy loss [a.u.];Yield", 100, 0., 50e-6);
94 fHists.push_back(fhEloss[i]);
95 }
96
97 Int_t nofSortBins = 200;
98 fhElossSort.resize(2);
99 for (int j = 0; j < 2; j++) {
100 fhElossSort[j].resize(fNofTrdLayers);
101 for (Int_t i = 0; i < fNofTrdLayers; i++) {
102 std::stringstream ss1;
103 if (j == 0) ss1 << "fhElossSortEl" << i;
104 if (j == 1) ss1 << "fhElossSortPi" << i;
105 fhElossSort[j][i] =
106 new TH1D(ss1.str().c_str(), (ss1.str() + ";Energy loss [a.u.];Counters").c_str(), nofSortBins, 0., 0.);
107 }
108 }
109
110 // Histograms for algorithm testing
111 Int_t nofBins = 10000;
112 fhOutput.resize(2);
113 fhCumProbOutput.resize(2);
114 fhInput.resize(2);
115 for (Int_t i = 0; i < 2; i++) {
116 s = (i == 0) ? "El" : "Pi";
117 fhOutput[i] = new TH1D(("fhOutput" + s).c_str(), "fhOutput;Output;Counter", nofBins, fMinEval, fMaxEval);
118 fhCumProbOutput[i] = new TH1D(("fhCumProbOutput" + s).c_str(), "fhCumProbOutput;Output;Cumulative probability",
119 nofBins, fMinEval, fMaxEval);
120
121 fhInput[i].resize(fNofTrdLayers);
122 for (Int_t j = 0; j < fNofTrdLayers; j++) {
123 stringstream ss;
124 ss << s << j;
125 fhInput[i][j] = new TH1D(("fhInput" + ss.str()).c_str(), "fhInput", 100, 0.0, 0.0);
126 }
127 }
128}
129
131
133{
134 FairRootManager* ioman = FairRootManager::Instance();
135 if (NULL == ioman) {
136 LOG(fatal) << GetName() << "::Init: RootManager not instantised";
137 }
138
139 fMCTracks = (TClonesArray*) ioman->GetObject("MCTrack");
140 if (NULL == fMCTracks) {
141 LOG(fatal) << GetName() << "::Init: No MCTrack array";
142 }
143
144 fTrdPoints = (TClonesArray*) ioman->GetObject("TrdPoint");
145 if (NULL == fTrdPoints) {
146 LOG(fatal) << GetName() << "::Init: No TRDPoint array";
147 }
148
149 fTrdTracks = (TClonesArray*) ioman->GetObject("TrdTrack");
150 if (NULL == fTrdTracks) {
151 LOG(fatal) << GetName() << "::Init: No TrdTrack array";
152 }
153
154 fTrdTrackMatches = (TClonesArray*) ioman->GetObject("TrdTrackMatch");
155 if (NULL == fTrdTrackMatches) {
156 LOG(fatal) << GetName() << "::Init: No TrdTrackMatch array";
157 }
158
159 fTrdHits = (TClonesArray*) ioman->GetObject("TrdHit");
160 if (NULL == fTrdHits) {
161 LOG(fatal) << GetName() << "::Init: No TrdHit array";
162 }
163
164 return kSUCCESS;
165}
166
168{
169 fEventNum++;
170 cout << "CbmTrdElectronsTrainAnn, event " << fEventNum << endl;
171
175
176 cout << "Nof electrons = " << fEloss[0].size() << endl;
177 cout << "Nof pions = " << fEloss[1].size() << endl;
178
179 // Finish();
180}
181
183{
184 Run();
185 Draw();
186
187 if (fOutputDir != "") {
188 gSystem->mkdir(fOutputDir.c_str(), true);
189 }
190
192 TFile* oldFile = gFile;
193 TDirectory* oldDir = gDirectory;
194
195 TFile* f = new TFile(string(fOutputDir + "/trd_elid_hist.root").c_str(), "RECREATE");
196 for (unsigned int i = 0; i < fHists.size(); i++) {
197 fHists[i]->Write();
198 }
199
201 gFile = oldFile;
202 gDirectory = oldDir;
203
204 f->Close();
205}
206
208{
212
213 cout << "Nof electrons = " << fEloss[0].size() << endl;
214 cout << "Nof pions = " << fEloss[1].size() << endl;
215
216 Finish();
217}
218
220{
221 if (fBeamDataFile == "" || fBeamDataPiHist == "" || fBeamDataElHist == "") {
222 LOG(fatal) << GetName() << "::FillElossVectorReal: Set input file for beam data and histogram names!";
223 }
224
226 TFile* oldFile = gFile;
227 TDirectory* oldDir = gDirectory;
228
229 TFile* file = new TFile(fBeamDataFile.c_str(), "READ");
230 LOG_IF(fatal, !file) << "Could not open file " << fBeamDataFile;
231
232 TH1F* hPion = static_cast<TH1F*>(file->Get<TH1F>(fBeamDataPiHist.c_str())->Clone());
233 LOG_IF(fatal, !hPion) << "Histogram " << fBeamDataPiHist << " not found in file " << fBeamDataFile;
234
235 TH1F* hElectron = static_cast<TH1F*>(file->Get<TH1F>(fBeamDataElHist.c_str())->Clone());
236 LOG_IF(fatal, !hElectron) << "Histogram " << fBeamDataElHist << " not found in file " << fBeamDataFile;
237
238 double scaleX =
239 fhEloss[0]->GetXaxis()->GetBinUpEdge(fhEloss[0]->GetNbinsX()) / hPion->GetXaxis()->GetBinUpEdge(hPion->GetNbinsX());
240
241 int nofSimulatedParticles = 1000000;
242 fEloss[0].resize(nofSimulatedParticles);
243 fEloss[1].resize(nofSimulatedParticles);
244
245 for (int iP = 0; iP < 2; iP++) {
246 TH1F* hist = hElectron;
247 if (iP == 1) hist = hPion;
248 for (int iT = 0; iT < nofSimulatedParticles; iT++) {
249 fEloss[iP][iT].resize(fNofTrdLayers);
250 for (int iH = 0; iH < fNofTrdLayers; iH++) {
251 double eloss = scaleX * hist->GetRandom();
252 TrdEloss e(eloss);
253 fEloss[iP][iT][iH] = e;
254 }
255 }
256 }
257
259 gFile = oldFile;
260 gDirectory = oldDir;
261}
262
264{
265 Int_t nofTrdTracks = fTrdTracks->GetEntriesFast();
266 for (Int_t iTrdTrack = 0; iTrdTrack < nofTrdTracks; iTrdTrack++) {
267 CbmTrdTrack* trdtrack = (CbmTrdTrack*) fTrdTracks->At(iTrdTrack);
268 Int_t nHits = trdtrack->GetNofHits();
269
270
271 CbmTrackMatchNew* trdmatch = (CbmTrackMatchNew*) fTrdTrackMatches->At(iTrdTrack);
272 Int_t iMC = trdmatch->GetMatchedLink().GetIndex();
273 if (iMC < 0 || iMC > fMCTracks->GetEntriesFast()) continue;
274
275 CbmMCTrack* mctrack = (CbmMCTrack*) fMCTracks->At(iMC);
276 Int_t pdg = TMath::Abs(mctrack->GetPdgCode());
277 // Double_t momMC = mctrack->GetP();
278 Double_t motherId = mctrack->GetMotherId();
279
281 if (motherId != -1) continue;
282
283 int iP = 0; //[0] = electron
284 if (pdg == 211) iP = 1; //[1] = pion
285 vector<TrdEloss> v;
286
287 if (nHits < fNofTrdLayers) continue;
288
289 for (Int_t iHit = 0; iHit < nHits; iHit++) {
290 Int_t hitIndex = trdtrack->GetHitIndex(iHit);
291 CbmTrdHit* trdhit = (CbmTrdHit*) fTrdHits->At(hitIndex);
292 TrdEloss e(trdhit->GetELoss(), trdhit->GetELoss(), 0 /*trdhit->GetELossdEdX(), trdhit->GetELossTR()*/);
293 if (static_cast<UInt_t>(fNofTrdLayers) == v.size()) break;
294 v.push_back(e);
295 } //iHit
296 fEloss[iP].push_back(v);
297 } //iTrdTrack
298}
299
301{
302 // [0]=electron, [1]=pion
303 for (int i = 0; i < 2; i++) {
304 for (unsigned int iT = 0; iT < fEloss[i].size(); iT++) {
305 Int_t nHits = fEloss[i][iT].size();
306 double sumEloss = 0.;
307 for (int iH = 0; iH < nHits; iH++) {
308 sumEloss += fEloss[i][iT][iH].fEloss;
309 fhEloss[i]->Fill(fEloss[i][iT][iH].fEloss);
310 } // iH
311 fhMeanEloss[i]->Fill(sumEloss / nHits);
312 } //iT
313 }
314}
315
317{
318 vector<double> v;
319 v.resize(fNofTrdLayers);
320 for (int iP = 0; iP < 2; iP++) {
321 for (unsigned int iT = 0; iT < fEloss[iP].size(); iT++) {
322 for (int iH = 0; iH < fNofTrdLayers; iH++) {
323 v[iH] = fEloss[iP][iT][iH].fEloss;
324 }
325 std::sort(v.begin(), v.end());
326 for (int iH = 0; iH < fNofTrdLayers; iH++) {
327 fhElossSort[iP][iH]->Fill(v[iH]);
328 }
329 } //iT
330
331 for (int iH = 0; iH < fNofTrdLayers; iH++) {
332 fhElossSort[iP][iH]->Scale(1. / fhElossSort[iP][iH]->Integral());
333 }
334 }
335}
336
338{
339 if (fIdMethod == kBDT || fIdMethod == kANN) {
340 if (fIsDoTrain) DoTrain();
341 DoTest();
342 }
343 else if (fIdMethod == kMEDIAN || fIdMethod == kLIKELIHOOD || fIdMethod == kMeanCut) {
344 DoTest();
345 }
346}
347
349{
350 if (fIdMethod == kBDT || fIdMethod == kANN) {
351 if (fTransformType == 1) {
352 Transform1();
353 }
354 else if (fTransformType == 2) {
355 Transform2();
356 }
357 }
358}
359
360
362{
363 Double_t ANNCoef1 = 1.06;
364 Double_t ANNCoef2 = 0.57;
365 //Double_t ANNCoef1[] = {1.04,1.105,1.154,1.277,1.333,1.394,1.47,1.50,1.54,1.58};
366 //Double_t ANNCoef2[] = {0.548,0.567,0.585,0.63,0.645,0.664,0.69,0.705,0.716,0.723};
367 for (UInt_t i = 0; i < fAnnInput.size(); i++) {
368 fAnnInput[i] = fAnnInput[i] * 1e6;
369 fAnnInput[i] = (fAnnInput[i] - ANNCoef1) / ANNCoef2 - 0.225;
370 }
371 sort(fAnnInput.begin(), fAnnInput.end());
372 for (UInt_t i = 0; i < fAnnInput.size(); i++) {
373 fAnnInput[i] = TMath::LandauI(fAnnInput[i]);
374 }
375}
376
378{
379 sort(fAnnInput.begin(), fAnnInput.end());
380
381 Int_t size = fAnnInput.size();
382 for (Int_t j = 0; j < size; j++) {
383 Int_t binNumEl = fhElossSort[0][j]->FindBin(fAnnInput[j]);
384 if (binNumEl > fhElossSort[0][j]->GetNbinsX()) binNumEl = fhElossSort[0][j]->GetNbinsX();
385 Double_t probEl = fhElossSort[0][j]->GetBinContent(binNumEl);
386
387 Int_t binNumPi = fhElossSort[1][j]->FindBin(fAnnInput[j]);
388 if (binNumPi > fhElossSort[1][j]->GetNbinsX()) binNumPi = fhElossSort[1][j]->GetNbinsX();
389 Double_t probPi = fhElossSort[1][j]->GetBinContent(binNumPi);
390
391 if (TMath::IsNaN(probPi / (probEl + probPi))) {
392 fAnnInput[j] = 0.;
393 }
394 else {
395 fAnnInput[j] = probPi / (probEl + probPi);
396 }
397 }
398}
399
401{
402 Double_t lPi = 1.;
403 Double_t lEl = 1.;
404 // sort(fAnnInput.begin(), fAnnInput.end());
405 for (UInt_t j = 0; j < fAnnInput.size(); j++) {
406 Int_t binNumEl = fhEloss[0]->FindBin(fAnnInput[j]);
407 if (binNumEl > fhEloss[0]->GetNbinsX()) binNumEl = fhEloss[0]->GetNbinsX();
408 Double_t probEl = fhEloss[0]->GetBinContent(binNumEl);
409 lEl = lEl * probEl;
410
411 Int_t binNumPi = fhEloss[1]->FindBin(fAnnInput[j]);
412 if (binNumPi > fhEloss[1]->GetNbinsX()) binNumPi = fhEloss[1]->GetNbinsX();
413 Double_t probPi = fhEloss[1]->GetBinContent(binNumPi);
414 lPi = lPi * probPi;
415 }
416 return lEl / (lEl + lPi);
417}
418
420{
421 sort(fAnnInput.begin(), fAnnInput.end());
422 double eval = 0.;
423 if (fNofTrdLayers % 2 == 0) {
424 eval = (fAnnInput[fNofTrdLayers / 2] + fAnnInput[fNofTrdLayers / 2 + 1]) / 2.;
425 }
426 else {
427 eval = fAnnInput[fNofTrdLayers / 2 + 1];
428 }
429 double scaleX = 1. / fhEloss[0]->GetXaxis()->GetBinUpEdge(fhEloss[0]->GetNbinsX());
430 return eval * scaleX;
431}
432
434{
435 double eval = 0.;
436 for (int i = 0; i < fNofTrdLayers; i++) {
437 eval += fAnnInput[i];
438 }
439 eval = eval / fNofTrdLayers;
440 double scaleX = 1. / fhEloss[0]->GetXaxis()->GetBinUpEdge(fhEloss[0]->GetNbinsX());
441 return eval * scaleX;
442}
443
445{
446 if (fIdMethod == kBDT) {
447 return fReader->EvaluateMVA("BDT");
448 }
449 else if (fIdMethod == kANN) {
450 Double_t par[fNofTrdLayers];
451 for (UInt_t i = 0; i < fAnnInput.size(); i++)
452 par[i] = fAnnInput[i];
453 return fNN->Evaluate(0, par);
454 }
455 else if (fIdMethod == kMEDIAN) {
456 return Median();
457 }
458 else if (fIdMethod == kLIKELIHOOD) {
459 return Likelihood();
460 }
461 else if (fIdMethod == kMeanCut) {
462 return MeanCut();
463 }
464 return -1.;
465}
466
467
469{
470 fAnnInput.clear();
471 fAnnInput.resize(fNofTrdLayers);
472
473 TTree* simu = CreateTree();
474 for (Int_t i = 0; i < 2; i++) {
475 for (UInt_t iT = 0; iT < fEloss[i].size() - 2; iT++) {
476 for (int iH = 0; iH < fNofTrdLayers; iH++) {
477 fAnnInput[iH] = fEloss[i][iT][iH].fEloss;
478 }
479 fXOut = (i == 0) ? 1. : -1.;
480 Transform();
481 simu->Fill();
482 if (Int_t(iT) >= fNofTrainSamples) break;
483 } //iT
484 } //i
485
486 if (fIdMethod == kANN) {
487 if (NULL != fNN) delete fNN;
488 TString mlpString = CreateAnnString();
489 cout << "-I- Create ANN: " << mlpString << endl;
490 cout << "-I- Number of training epochs = " << fNofAnnEpochs << endl;
491 fNN = new TMultiLayerPerceptron(mlpString, simu, "(Entry$+1)");
492 fNN->Train(fNofAnnEpochs, "+text,update=10");
493 stringstream ss;
494 ss << fOutputDir + "/ann_weights_" << fNofTrdLayers << ".txt";
495 fNN->DumpWeights(ss.str().c_str());
496 }
497 else if (fIdMethod == kBDT) {
498 CreateFactory(simu);
499 (TMVA::gConfig().GetIONames()).fWeightFileDir = fOutputDir;
500 TCut mycuts = "";
501 TCut mycutb = "";
502 // factory->PrepareTrainingAndTestTree(mycuts, mycutb,"SplitMode=Random:NormMode=NumEvents:!V");
503 //factory->BookMethod(TMVA::Types::kTMlpANN, "TMlpANN","!H:!V:NCycles=50:HiddenLayers=N+1");
504 stringstream ss;
505 ss << "nTrain_Signal=" << fNofTrainSamples - 500 << ":nTrain_Background=" << fNofTrainSamples - 500
506 << ":nTest_Signal=0:nTest_Background=0";
507 //TMVA_API
508 // factory->PrepareTrainingAndTestTree("", ss.str().c_str());
509 // factory->BookMethod(TMVA::Types::kBDT, "BDT", "!H:!V:NTrees=400:BoostType=AdaBoost:SeparationType=GiniIndex:nCuts=20:PruneMethod=CostComplexity:PruneStrength=4.5");
510 // factory->TrainAllMethods();
511 }
512}
513
515{
516 cout << "-I- Start pretesting " << endl;
517 if (fIdMethod == kBDT) {
518 cout << "-I- IdMethod = kBDT " << endl;
519 }
520 else if (fIdMethod == kANN) {
521 cout << "-I- IdMethod = kANN " << endl;
522 }
523 else if (fIdMethod == kMEDIAN) {
524 cout << "-I- IdMethod = kMEDIANA " << endl;
525 }
526 else if (fIdMethod == kLIKELIHOOD) {
527 cout << "-I- IdMethod = kLIKELIHOOD " << endl;
528 }
529
530 fAnnInput.clear();
531 fAnnInput.resize(fNofTrdLayers);
532
533 if (fIdMethod == kBDT) {
535 fReader->BookMVA("BDT", fOutputDir + "/TMVAnalysis_BDT.weights.xml");
536 }
537 else if (fIdMethod == kANN) {
538 if (fNN != NULL) delete fNN;
539 TTree* simu = CreateTree();
540 TString mlpString = CreateAnnString();
541 fNN = new TMultiLayerPerceptron(mlpString, simu, "(Entry$+1)");
542 stringstream ss;
543 ss << fOutputDir + "/ann_weights_" << fNofTrdLayers << ".txt";
544 fNN->LoadWeights(ss.str().c_str());
545 }
546
547 for (Int_t i = 0; i < 2; i++) {
548 for (UInt_t iT = 0; iT < fEloss[i].size(); iT++) {
549 if (Int_t(iT) < fNofTrainSamples) continue; // exclude training samples
550 if (iT % 100000 == 0) cout << "-I- read number: " << iT << endl;
551
552 //change electron hit to pion hit and vice versa
553 //int randPi = rand()%fElossPi.size();
554 //for (int i1 = 0; i1 < 5; i1++){
555 // fElossEl[iE][i1] = fElossPi[randPi][i1];
556 // }
557 for (int iH = 0; iH < fNofTrdLayers; iH++) {
558 if (fSigmaError != 0.) fEloss[i][iT][iH].fEloss += gRandom->Gaus(0., fSigmaError);
559 fAnnInput[iH] = fEloss[i][iT][iH].fEloss;
560 }
561
562 Transform();
563 Bool_t isEl = (i == 0) ? true : false;
564 FillAnnInputHist(isEl);
565 Double_t nnEval = Eval(isEl);
566 if (nnEval > fMaxEval) nnEval = fMaxEval - 0.01;
567 if (nnEval < fMinEval) nnEval = fMinEval + 0.01;
568 fhOutput[i]->Fill(nnEval);
569 }
570 }
572}
573
575{
576 fAnnInput.clear();
577 fAnnInput.resize(fNofTrdLayers);
578
579 DoPreTest();
580 double cut = FindOptimalCut();
581 cout << "-I- Optimal cut = " << cut << " for " << 100 * fElIdEfficiency << "% electron efficiency" << endl;
582 cout << "-I- Start testing" << endl;
583 fAnnInput.clear();
584 fAnnInput.resize(fNofTrdLayers);
585
586 if (fIdMethod == kBDT) {
588 //reader->BookMVA("TMlpANN", "weights/TMVAnalysis_TMlpANN.weights.txt");
589 fReader->BookMVA("BDT", fOutputDir + "/TMVAnalysis_BDT.weights.xml");
590 }
591 else if (fIdMethod == kANN) {
592 if (fNN != NULL) delete fNN;
593 TTree* simu = CreateTree();
594 TString mlpString = CreateAnnString();
595 fNN = new TMultiLayerPerceptron(mlpString, simu, "(Entry$+1)");
596 stringstream ss;
597 ss << fOutputDir + "/ann_weights_" << fNofTrdLayers << ".txt";
598 fNN->LoadWeights(ss.str().c_str());
599 }
600
601 int nofPiLikeEl = 0;
602 int nofElLikePi = 0;
603 int nofElTest = 0;
604 int nofPiTest = 0;
605
606 for (Int_t i = 0; i < 2; i++) {
607 for (UInt_t iT = 0; iT < fEloss[i].size(); iT++) {
608 if (Int_t(iT) < fNofTrainSamples) continue; //exclude training samples
609 if (iT % 100000 == 0) cout << "-I- Read number: " << iT << endl;
610
611 for (Int_t iH = 0; iH < fNofTrdLayers; iH++) {
612 fAnnInput[iH] = fEloss[i][iT][iH].fEloss;
613 }
614 Transform();
615
616 Bool_t isEl = (i == 0) ? true : false;
617 Double_t nnEval = Eval(isEl);
618 if (nnEval > fMaxEval) nnEval = fMaxEval - 0.01;
619 if (nnEval < fMinEval) nnEval = fMinEval + 0.01;
620 if (i == 0) {
621 nofElTest++;
622 if (nnEval < cut) nofElLikePi++;
623 }
624 else if (i == 1) {
625 nofPiTest++;
626 if (nnEval > cut) nofPiLikeEl++;
627 }
628 }
629 }
630 double piSupp = -1;
631 if (nofPiLikeEl != 0) piSupp = (double) nofPiTest / nofPiLikeEl;
632 double elEff = (double) nofElLikePi / nofElTest * 100.;
633
634 cout << "Testing results:" << endl;
635 cout << "nof TRD layers " << fNofTrdLayers << endl;
636 cout << "cut = " << cut << endl;
637 cout << "nof El = " << nofElTest << endl;
638 cout << "nof Pi = " << nofPiTest << endl;
639 cout << "nof Pi identified as El = " << nofPiLikeEl << endl;
640 cout << "nof El identified as Pi = " << nofElLikePi << endl;
641 cout << "pion suppression = " << nofPiTest << "/" << nofPiLikeEl << " = " << piSupp << endl;
642 cout << "electron efficiency loss in % = " << nofElLikePi << "/" << nofElTest << " = " << elEff << endl;
643
644 // write results to histogramm
645 fhResults->SetBinContent(1, piSupp);
646 fhResults->SetBinContent(2, elEff);
647}
648
650{
651 for (Int_t i = 0; i < 2; i++) {
652 Double_t cumProb = 0.;
653 Double_t nofTr = fhOutput[i]->GetEntries();
654 for (Int_t iH = 1; iH <= fhOutput[i]->GetNbinsX(); iH++) {
655 cumProb += fhOutput[i]->GetBinContent(iH);
656 Double_t pr = (i == 0) ? 1. - cumProb / nofTr : cumProb / nofTr;
657 fhCumProbOutput[i]->SetBinContent(iH, pr);
658 }
659 }
660}
661
663{
664 Int_t nBins = fhCumProbOutput[0]->GetNbinsX();
665 Double_t x[nBins], y[nBins];
666
667 for (Int_t i = 1; i <= nBins; i++) {
668 Double_t cpi = 100. * fhCumProbOutput[1]->GetBinContent(i);
669 if (cpi == 100.) cpi = 100. - 0.0001;
670 y[i - 1] = 100. / (100. - cpi);
671 x[i - 1] = 100. * fhCumProbOutput[0]->GetBinContent(i);
672 }
673 TGraph* rocGraph = new TGraph(nBins, x, y);
674 rocGraph->SetTitle("ROC diagram");
675 rocGraph->GetXaxis()->SetTitle("Efficiency [%]");
676 rocGraph->GetYaxis()->SetTitle("Pion suppression");
677 rocGraph->SetLineWidth(3);
678 return rocGraph;
679}
680
682{
683 Double_t optimalCut = -1;
684 for (Int_t i = 1; i <= fhCumProbOutput[0]->GetNbinsX(); i++) {
685 if (fhCumProbOutput[0]->GetBinContent(i) <= fElIdEfficiency) {
686 optimalCut = fhCumProbOutput[0]->GetBinCenter(i);
687 return optimalCut;
688 }
689 }
690 return optimalCut;
691}
692
694{
695 fAnnInput.clear();
696 fAnnInput.resize(fNofTrdLayers);
697 TTree* simu = new TTree("MonteCarlo", "MontecarloData");
698 size_t buf_size = 100;
699 char txt1[buf_size], txt2[buf_size];
700 for (Int_t i = 0; i < fNofTrdLayers; i++) {
701 snprintf(txt1, buf_size - 1, "x%d", i);
702 snprintf(txt2, buf_size - 1, "x%d/F", i);
703 simu->Branch(txt1, &fAnnInput[i], txt2);
704 }
705 simu->Branch("xOut", &fXOut, "xOut/F");
706
707 return simu;
708}
709
711{
712 string st = "";
713 size_t buf_size = 50;
714 char txt[buf_size];
715 for (Int_t i = 0; i < fNofTrdLayers - 1; i++) {
716 snprintf(txt, buf_size - 1, "x%d", i);
717 st = st + txt + ",";
718 }
719 snprintf(txt, buf_size - 1, "x%d", fNofTrdLayers - 1);
720 st = st + txt + "";
721 snprintf(txt, buf_size - 1, "%d", 2 * fNofTrdLayers);
722 st = st + ":" + txt + ":xOut";
723 return st;
724}
725
727{
728 if (fOutputDir != "") gSystem->mkdir(fOutputDir.c_str(), true);
729 TFile* outputFile = TFile::Open(string(fOutputDir + "/tmva_output.root").c_str(), "RECREATE");
730
731 TMVA::Factory* factory = new TMVA::Factory("TMVAnalysis", outputFile);
732
733 TCut piCut = "xOut<0";
734 TCut elCut = "xOut>0";
735 //TMVA_API
736 //factory->SetInputTrees(simu, elCut, piCut);
737 size_t buf_size = 100;
738 char txt1[buf_size];
739 for (Int_t i = 0; i < fNofTrdLayers; i++) {
740 snprintf(txt1, buf_size - 1, "x%d", i);
742 //factory->AddVariable(txt1, 'F');
743 }
744
745 return factory;
746}
747
749{
750 if (!fReader) delete fReader;
751
752 fAnnInput.clear();
753 fAnnInput.resize(fNofTrdLayers);
754
755 TMVA::Reader* reader = new TMVA::Reader();
756
757 size_t buf_size = 100;
758 char txt1[buf_size];
759 for (Int_t i = 0; i < fNofTrdLayers; i++) {
760 snprintf(txt1, buf_size - 1, "x%d", i);
761 reader->AddVariable(txt1, &fAnnInput[i]);
762 }
763 return reader;
764}
765
767{
768 for (UInt_t j = 0; j < fAnnInput.size(); j++) {
769 if (isEl) {
770 fhInput[0][j]->Fill(fAnnInput[j]);
771 }
772 else {
773 fhInput[1][j]->Fill(fAnnInput[j]);
774 }
775 }
776}
777
779{
781 TCanvas* cEloss = new TCanvas("trd_elid_eloss", "trd_elid_eloss", 1200, 600);
782 cEloss->Divide(2, 1);
783 cEloss->cd(1);
784 DrawH1(fhMeanEloss, list_of("e^{#pm}")("#pi^{#pm}"), kLinear, kLog, true, 0.8, 0.8, 0.99, 0.99);
785 cEloss->cd(2);
786 DrawH1(fhEloss, list_of("e^{#pm}")("#pi^{#pm}"), kLinear, kLog, true, 0.8, 0.8, 0.99, 0.99);
787
788
789 TCanvas* cElossSort = new TCanvas("trd_elid_eloss_sort", "trd_elid_eloss_sort", 1200, 900);
790 cElossSort->Divide(4, 3);
791 for (int iL = 0; iL < fNofTrdLayers; iL++) {
792 cElossSort->cd(iL + 1);
793 DrawH1(list_of(fhElossSort[0][iL])(fhElossSort[1][iL]), list_of("e^{#pm}")("#pi^{#pm}"), kLinear, kLog, true, 0.8,
794 0.8, 0.99, 0.99);
795 }
796
797 // TCanvas* cClassifierOutput = new TCanvas("trd_elid_classifier_output","trd_elid_classifier_output", 500, 500);
798 TH1D* out0 = (TH1D*) fhOutput[0]->Clone();
799 TH1D* out1 = (TH1D*) fhOutput[1]->Clone();
800 out0->Rebin(50);
801 out1->Rebin(50);
802 DrawH1(list_of(out0)(out1), list_of("e^{#pm}")("#pi^{#pm}"), kLinear, kLog, true, 0.8, 0.8, 0.99, 0.99);
803 out0->Scale(1. / out0->Integral());
804 out1->Scale(1. / out1->Integral());
805
806 // TCanvas* cCumProbOutput = new TCanvas("trd_elid_cum_prob_output","trd_elid_cum_prob_output", 500,500);
807 DrawH1(fhCumProbOutput, list_of("e^{#pm}")("#pi^{#pm}"), kLinear, kLinear, true, 0.8, 0.8, 0.99, 0.99);
808
809 // TCanvas* cRoc = new TCanvas("trd_elid_roc","trd_elid_roc", 500, 500);
810 TGraph* rocGraph = CreateRocDiagramm();
811 DrawGraph(rocGraph);
812 gPad->SetLogy();
813
814 TCanvas* cInput = new TCanvas("trd_elid_input", "trd_elid_ann_input", 10, 10, 800, 800);
815 cInput->Divide(4, 3);
816 for (int i = 0; i < fNofTrdLayers; i++) {
817 cInput->cd(i + 1);
818 DrawH1(list_of(fhInput[0][i])(fhInput[1][i]), list_of("e^{#pm}")("#pi^{#pm}"), kLinear, kLog, true, 0.8, 0.8, 0.99,
819 0.99);
820 }
821}
822
ClassImp(CbmConverterManager)
void SetDefaultDrawStyle()
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)
void DrawGraph(TGraph *graph, 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
Class for hits in TRD detector.
static constexpr size_t size()
Definition KfSimdPseudo.h:2
fscal v[fmask::Size]
Definition KfSimdPseudo.h:4
int32_t GetMotherId() const
Definition CbmMCTrack.h:69
int32_t GetPdgCode() const
Definition CbmMCTrack.h:68
const CbmLink & GetMatchedLink() const
Definition CbmMatch.h:41
virtual int32_t GetNofHits() const
Definition CbmTrack.h:58
int32_t GetHitIndex(int32_t iHit) const
Definition CbmTrack.h:59
virtual ~CbmTrdElectronsTrainAnn()
Destructor.
void FillElossHist()
Fill histograms with energy loss information.
void Draw(Option_t *="")
Draw results.
virtual InitStatus Init()
Inherited from FairTask.
virtual void Finish()
Inherited from FairTask.
CbmTrdElectronsTrainAnn(int nofTrdLayers)
Default constructor.
void FillElossVectorSim()
Fill vector with energy loss information for simulated data.
void SortElossAndFillHist()
Sort energy losses and fill histograms.
std::vector< TH1 * > fhCumProbOutput
void FillElossVectorReal()
Fill vector with energy loss information simulated from real data spectra.
std::vector< std::vector< TH1 * > > fhInput
std::vector< std::vector< std::vector< TrdEloss > > > fEloss
TMVA::Factory * CreateFactory(TTree *simu)
std::vector< Float_t > fAnnInput
virtual void Exec(Option_t *opt)
Inherited from FairTask.
TMultiLayerPerceptron * fNN
std::vector< std::vector< TH1 * > > fhElossSort
data class for a reconstructed Energy-4D measurement in the TRD
Definition CbmTrdHit.h:40
double GetELoss() const
Definition CbmTrdHit.h:79
Represents information about energy losses in one layer.