]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaRec/AliTRDpidRefMakerNN.cxx
New directory for the VMC tests (Ivana, Eva)
[u/mrichter/AliRoot.git] / TRD / qaRec / AliTRDpidRefMakerNN.cxx
1 /*************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
7 * Permission to use, copy, modify and distribute this software and its   *
8 * documentation strictly for non-commercialf purposes is hereby granted  *
9 * without fee, provided that the above copyright notice appears in all   *
10 * copies and that both the copyright notice and this permission notice   *
11 * appear in the supporting documentation. The authors make no claims     *
12 * about the suitability of this software for any purpose. It is          *
13 * provided "as is" without express or implied warranty.                  *
14 **************************************************************************/
15
16 /* $Id: AliTRDpidRefMakerNN.cxx 27496 2008-07-22 08:35:45Z cblume $ */
17
18 ////////////////////////////////////////////////////////////////////////////
19 //                                                                        //
20 //  Builds the reference tree for the training of neural networks         //
21 //                                                                        //
22 ////////////////////////////////////////////////////////////////////////////
23
24 #include "TSystem.h"
25 #include "TDatime.h"
26 #include "TPDGCode.h"
27 #include "TH1F.h"
28 #include "TH2F.h"
29 #include "TFile.h"
30 #include "TGraphErrors.h"
31 #include "TTree.h"
32 #include "TEventList.h"
33 #include "TMultiLayerPerceptron.h"
34
35 #include "AliPID.h"
36 #include "AliESDtrack.h"
37 #include "AliTrackReference.h"
38
39 #include "AliTRDtrackV1.h"
40 #include "AliTRDReconstructor.h"
41 #include "AliTRDpidUtil.h"
42 #include "AliTRDpidRefMakerNN.h"
43 #include "AliTRDpidUtil.h"
44
45 #include "../Cal/AliTRDCalPID.h"
46 #include "../Cal/AliTRDCalPIDNN.h"
47 #include "info/AliTRDtrackInfo.h"
48 #include "info/AliTRDv0Info.h"
49
50 ClassImp(AliTRDpidRefMakerNN)
51
52 //________________________________________________________________________
53 AliTRDpidRefMakerNN::AliTRDpidRefMakerNN() 
54   :AliTRDpidRefMaker("PidRefMakerNN", "PID(NN) Reference Maker")
55 //   :AliTRDrecoTask("PidRefMakerNN", "PID(NN) Reference Maker")
56   ,fTrainMomBin(kAll)
57   ,fEpochs(1000)
58   ,fMinTrain(100)
59   ,fDate(0)
60   ,fDoTraining(0)
61   ,fContinueTraining(0)
62   ,fTrainPath(0x0)
63   ,fScale(0)
64 {
65   //
66   // Default constructor
67   //
68
69   memset(fTrain, 0, AliTRDCalPID::kNMom*AliTRDgeometry::kNlayer*sizeof(TEventList*));
70   memset(fTest, 0, AliTRDCalPID::kNMom*AliTRDgeometry::kNlayer*sizeof(TEventList*));
71   memset(fNet, 0, AliTRDgeometry::kNlayer*sizeof(TMultiLayerPerceptron*));
72
73   SetAbundance(.67);
74   SetScaledEdx(Float_t(AliTRDCalPIDNN::kMLPscale));
75   TDatime datime;
76   fDate = datime.GetDate();
77
78   DefineInput(1, TObjArray::Class());
79   DefineOutput(1, TTree::Class());
80 }
81
82
83 //________________________________________________________________________
84 AliTRDpidRefMakerNN::~AliTRDpidRefMakerNN() 
85 {
86 }
87
88
89 //________________________________________________________________________
90 void AliTRDpidRefMakerNN::CreateOutputObjects()
91 {
92   // Create histograms
93   // Called once
94
95   AliTRDpidRefMaker::CreateOutputObjects();
96   TGraphErrors *gEffisTrain = new TGraphErrors(kMoniTrain);
97   gEffisTrain -> SetLineColor(4);
98   gEffisTrain -> SetMarkerColor(4);
99   gEffisTrain -> SetMarkerStyle(29);
100   gEffisTrain -> SetMarkerSize(1);
101
102   TGraphErrors *gEffisTest = new TGraphErrors(kMoniTrain);
103   gEffisTest -> SetLineColor(2);
104   gEffisTest -> SetMarkerColor(2);
105   gEffisTest -> SetMarkerStyle(29);
106   gEffisTest -> SetMarkerSize(1);
107
108   fContainer -> AddAt(gEffisTrain,kGraphTrain);
109   fContainer -> AddAt(gEffisTest,kGraphTest);
110 }
111
112
113
114 //________________________________________________________________________
115 Bool_t AliTRDpidRefMakerNN::PostProcess()
116 {
117   // Draw result to the screen
118   // Called once at the end of the query
119
120   // build the training andthe test list for the neural networks
121   MakeTrainingLists();        
122   if(!fDoTraining) return kTRUE;
123
124   // train the neural networks and build the refrence histos for 2-dim LQ
125   gSystem->Exec(Form("mkdir ./Networks_%d/",fDate));
126   AliDebug(2, Form("TrainMomBin [%d] [%d]", fTrainMomBin, kAll));
127
128   // train single network for a single momentum (recommended)
129   if(!(fTrainMomBin == kAll)){
130     if(fTrain[fTrainMomBin][0] -> GetN() < fMinTrain){
131       AliDebug(2, Form("Warning in AliTRDpidRefMakerNN::PostProcess : Not enough events for training available! Please check Data sample!"));
132       return kFALSE;
133     }
134     MakeRefs(fTrainMomBin);
135 //     TrainNetworks(fTrainMomBin);
136     MonitorTraining(fTrainMomBin);
137   }
138   // train all momenta
139   else{
140     for(Int_t iMomBin = 0; iMomBin < AliTRDCalPID::kNMom; iMomBin++){
141       if(fTrain[iMomBin][0] -> GetN() < fMinTrain){
142         AliDebug(2, Form("Warning in AliTRDpidRefMakerNN::PostProcess : Not enough events for training available for momentum bin [%d]! Please check Data sample!", iMomBin));
143   continue;
144       }
145       MakeRefs(fTrainMomBin);
146 //       TrainNetworks(iMomBin);
147       MonitorTraining(iMomBin);
148     }
149   }
150
151   return kTRUE; // testing protection
152 }
153
154
155 //________________________________________________________________________
156 void AliTRDpidRefMakerNN::MakeTrainingLists() 
157 {
158   //
159   // build the training lists for the neural networks
160   //
161
162   if (!fData) {
163     LoadFile("TRD.CalibPidRefMakerNN.root");
164   }
165
166   if (!fData) {
167     Printf("ERROR tree for training list not available");
168     return;
169   }
170
171   AliDebug(2, "  Making training lists! \n");
172
173   Int_t nPart[AliPID::kSPECIES][AliTRDCalPID::kNMom];
174   memset(nPart, 0, AliPID::kSPECIES*AliTRDCalPID::kNMom*sizeof(Int_t));
175
176   // set needed branches
177   LinkPIDdata();
178
179   // start first loop to check total number of each particle type
180   for(Int_t iEv=0; iEv < fData -> GetEntries(); iEv++){
181     fData -> GetEntry(iEv);
182
183     // use only events with goes through 6 layers TRD
184     if(fPIDdataArray->fNtracklets != AliTRDgeometry::kNlayer) continue;
185
186     for(Int_t ily=AliTRDgeometry::kNlayer; ily--;) nPart[fPIDbin][fPIDdataArray->fData[ily].fPLbin & 0xf]++;
187   }
188
189   AliDebug(2, "Particle multiplicities:");
190   for(Int_t iMomBin = 0; iMomBin <AliTRDCalPID::kNMom; iMomBin++)
191     AliDebug(2, Form("Momentum[%d]  Elecs[%d] Muons[%d] Pions[%d] Kaons[%d] Protons[%d]", iMomBin, nPart[AliPID::kElectron][iMomBin], nPart[AliPID::kMuon][iMomBin], nPart[AliPID::kPion][iMomBin], nPart[AliPID::kKaon][iMomBin], nPart[AliPID::kProton][iMomBin]));
192
193   // implement counter of training and test sample size
194   Int_t iTrain[AliTRDCalPID::kNMom], iTest[AliTRDCalPID::kNMom];
195   memset(iTrain, 0, AliTRDCalPID::kNMom*sizeof(Int_t));
196   memset(iTest, 0, AliTRDCalPID::kNMom*sizeof(Int_t));
197
198   // set training sample size per momentum interval to 2/3 
199   // of smallest particle counter and test sample to 1/3
200   for(Int_t iMomBin = 0; iMomBin < AliTRDCalPID::kNMom; iMomBin++){
201     iTrain[iMomBin] = nPart[0][iMomBin];
202     for(Int_t iPart = 1; iPart < AliPID::kSPECIES; iPart++){
203       if(iTrain[iMomBin] > nPart[iPart][iMomBin])
204   iTrain[iMomBin] = nPart[iPart][iMomBin];
205     } 
206     iTrain[iMomBin] = Int_t(iTrain[iMomBin] * .66);
207     iTest[iMomBin] = Int_t( iTrain[iMomBin] * .5);
208     AliDebug(2, Form("Momentum[%d]  Train[%d] Test[%d]", iMomBin, iTrain[iMomBin], iTest[iMomBin]));
209   }
210
211
212   // reset couters
213   memset(nPart, 0, AliPID::kSPECIES*AliTRDCalPID::kNMom*sizeof(Int_t));
214
215   // start second loop to set the event lists
216   for(Int_t iEv = 0; iEv < fData -> GetEntries(); iEv++){
217     fData -> GetEntry(iEv);
218
219     // use only events with goes through 6 layers TRD
220     if(fPIDdataArray->fNtracklets != AliTRDgeometry::kNlayer) continue;
221
222     for(Int_t ily=AliTRDgeometry::kNlayer; ily--;){ 
223       Int_t iMomBin = fPIDdataArray->fData[ily].fPLbin & 0xf;
224     
225       // set event list
226       if(nPart[fPIDbin][iMomBin] < iTrain[iMomBin]){
227         fTrain[iMomBin][ily] -> Enter(iEv + ily);
228         nPart[fPIDbin][iMomBin]++;
229       } else if(nPart[fPIDbin][iMomBin] < iTest[iMomBin]+iTrain[iMomBin]){
230         fTest[iMomBin][ily] -> Enter(iEv + ily);
231         nPart[fPIDbin][iMomBin]++;
232       } else continue;
233     }
234   }
235   
236   AliDebug(2, "Particle multiplicities in both lists:");
237   for(Int_t iMomBin = 0; iMomBin <AliTRDCalPID::kNMom; iMomBin++)
238     AliDebug(2, Form("Momentum[%d]  Elecs[%d] Muons[%d] Pions[%d] Kaons[%d] Protons[%d]", iMomBin, nPart[AliPID::kElectron][iMomBin], nPart[AliPID::kMuon][iMomBin], nPart[AliPID::kPion][iMomBin], nPart[AliPID::kKaon][iMomBin], nPart[AliPID::kProton][iMomBin]));
239 }
240
241
242 //________________________________________________________________________
243 void AliTRDpidRefMakerNN::MakeRefs(Int_t mombin) 
244 {
245   //
246   // train the neural networks
247   //
248   
249   
250   if (!fData) LoadFile(Form("TRD.Calib%s.root", GetName()));
251
252   if (!fData) {
253     AliError("Tree for training list not available");
254     return;
255   }
256
257   TDatime datime;
258   fDate = datime.GetDate();
259
260   AliDebug(2, Form("Training momentum bin %d", mombin));
261
262   // set variable to monitor the training and to save the development of the networks
263   Int_t nEpochs = fEpochs/kMoniTrain;       
264   AliDebug(2, Form("Training %d times %d epochs", kMoniTrain, nEpochs));
265
266   // make directories to save the networks 
267   gSystem->Exec(Form("rm -r ./Networks_%d/MomBin_%d",fDate, mombin));
268   gSystem->Exec(Form("mkdir ./Networks_%d/MomBin_%d",fDate, mombin));
269
270   // variable to check if network can load weights from previous training
271   Bool_t bFirstLoop[AliTRDgeometry::kNlayer];
272   memset(bFirstLoop, kTRUE, AliTRDgeometry::kNlayer*sizeof(Bool_t));
273
274   // train networks over several loops and save them after each loop
275   for(Int_t iLoop = 0; iLoop < kMoniTrain; iLoop++){
276     // loop over chambers
277     for(Int_t iChamb = 0; iChamb < AliTRDgeometry::kNlayer; iChamb++){
278       // set the event lists
279       fData -> SetEventList(fTrain[mombin][iChamb]);
280       fData -> SetEventList(fTest[mombin][iChamb]);
281       
282       AliDebug(2, Form("Trainingloop[%d] Chamber[%d]", iLoop, iChamb));
283       
284       // check if network is already implemented
285       if(bFirstLoop[iChamb] == kTRUE){
286   fNet[iChamb] = new TMultiLayerPerceptron("fdEdx[0],fdEdx[1],fdEdx[2],fdEdx[3],fdEdx[4],fdEdx[5],fdEdx[6],fdEdx[7]:15:7:fPID[0],fPID[1],fPID[2],fPID[3],fPID[4]!",fData,fTrain[mombin][iChamb],fTest[mombin][iChamb]);
287   fNet[iChamb] -> SetLearningMethod(TMultiLayerPerceptron::kStochastic);       // set learning method
288   fNet[iChamb] -> TMultiLayerPerceptron::SetEta(0.001);                        // set learning speed
289   if(!fContinueTraining){
290     if(DebugLevel()>=2) fNet[iChamb] -> Train(nEpochs,"text update=10, graph");
291     else fNet[iChamb] -> Train(nEpochs,"");
292   }
293   else{
294     fNet[iChamb] -> LoadWeights(Form("./Networks_%d/MomBin_%d/Net%d_%d",fTrainPath, mombin, iChamb, kMoniTrain - 1));
295     if(DebugLevel()>=2) fNet[iChamb] -> Train(nEpochs,"text update=10, graph+");      
296     else fNet[iChamb] -> Train(nEpochs,"+");                   
297   }
298   bFirstLoop[iChamb] = kFALSE;
299       }
300       else{    
301   if(DebugLevel()>=2) fNet[iChamb] -> Train(nEpochs,"text update=10, graph+");      
302   else fNet[iChamb] -> Train(nEpochs,"+");                   
303       }
304       
305       // save weights for monitoring of the training
306       fNet[iChamb] -> DumpWeights(Form("./Networks_%d/MomBin_%d/Net%d_%d",fDate, mombin, iChamb, iLoop));
307     } // end chamber loop
308   }   // end training loop
309 }
310
311
312
313 //________________________________________________________________________
314 void AliTRDpidRefMakerNN::MonitorTraining(Int_t mombin) 
315 {
316   //
317   // train the neural networks
318   //
319   
320   if (!fData) LoadFile(Form("TRD.Calib%s.root", GetName()));
321   if (!fData) {
322     AliError("Tree for training list not available");
323     return;
324   }
325
326   // init networks and set event list
327   for(Int_t iChamb = 0; iChamb < AliTRDgeometry::kNlayer; iChamb++){
328   fNet[iChamb] = new TMultiLayerPerceptron("fdEdx[0],fdEdx[1],fdEdx[2],fdEdx[3],fdEdx[4],fdEdx[5],fdEdx[6],fdEdx[7]:15:7:fPID[0],fPID[1],fPID[2],fPID[3],fPID[4]!",fData,fTrain[mombin][iChamb],fTest[mombin][iChamb]);   
329   fData -> SetEventList(fTrain[mombin][iChamb]);
330   fData -> SetEventList(fTest[mombin][iChamb]);
331   }
332
333   // implement variables for likelihoods
334   Float_t like[AliPID::kSPECIES][AliTRDgeometry::kNlayer];
335   memset(like, 0, AliPID::kSPECIES*AliTRDgeometry::kNlayer*sizeof(Float_t));
336   Float_t likeAll[AliPID::kSPECIES], totProb;
337
338   Double_t pionEffiTrain[kMoniTrain], pionEffiErrTrain[kMoniTrain];
339   Double_t pionEffiTest[kMoniTrain], pionEffiErrTest[kMoniTrain];
340   memset(pionEffiTrain, 0, kMoniTrain*sizeof(Double_t));
341   memset(pionEffiErrTrain, 0, kMoniTrain*sizeof(Double_t));
342   memset(pionEffiTest, 0, kMoniTrain*sizeof(Double_t));
343   memset(pionEffiErrTest, 0, kMoniTrain*sizeof(Double_t));
344
345   // init histos
346   const Float_t epsilon = 1/(2*(AliTRDpidUtil::kBins-1));     // get nice histos with bin center at 0 and 1
347   TH1F *hElecs, *hPions;
348   hElecs = new TH1F("hElecs","Likelihood for electrons", AliTRDpidUtil::kBins, 0.-epsilon, 1.+epsilon);
349   hPions = new TH1F("hPions","Likelihood for pions", AliTRDpidUtil::kBins, 0.-epsilon, 1.+epsilon);
350
351   TGraphErrors *gEffisTrain=0x0, *gEffisTest=0x0;
352   gEffisTrain = (TGraphErrors*)fContainer->At(kGraphTrain);
353   gEffisTest = (TGraphErrors*)fContainer->At(kGraphTest);
354
355   AliTRDpidUtil *util = new AliTRDpidUtil();
356   
357   // monitor training progress
358   for(Int_t iLoop = 0; iLoop < kMoniTrain; iLoop++){
359
360     // load weights
361     for(Int_t iChamb = 0; iChamb < AliTRDgeometry::kNlayer; iChamb++){
362       fNet[iChamb] -> LoadWeights(Form("./Networks_%d/MomBin_%d/Net%d_%d",fDate, mombin, iChamb, iLoop));
363     }
364
365     // event loop training list
366     for(Int_t iEvent = 0; iEvent < fTrain[mombin][0] -> GetN(); iEvent++ ){
367
368       // reset particle probabilities
369       for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
370   likeAll[iPart] = 1./AliPID::kSPECIES;
371       }
372       totProb = 0.;
373
374       fData -> GetEntry(fTrain[mombin][0] -> GetEntry(iEvent));
375       // use event only if it is electron or pion
376       if(!((fPID[AliPID::kElectron] == 1.0) || (fPID[AliPID::kPion] == 1.0))) continue;
377
378       // get the probabilities for each particle type in each chamber
379       for(Int_t iChamb = 0; iChamb < AliTRDgeometry::kNlayer; iChamb++){
380   for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
381     like[iPart][iChamb] = fNet[iChamb] -> Result(fTrain[mombin][iChamb] -> GetEntry(iEvent), iPart);
382     likeAll[iPart] *=  like[iPart][iChamb];
383   }
384       }
385
386       // get total probability and normalize it
387       for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
388   totProb += likeAll[iPart];
389       }
390       for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
391   likeAll[iPart] /= totProb;
392       }
393
394       // fill likelihood distributions
395       if(fPID[AliPID::kElectron] == 1)      
396   hElecs -> Fill(likeAll[AliPID::kElectron]);
397       if(fPID[AliPID::kPion] == 1)      
398   hPions -> Fill(likeAll[AliPID::kElectron]);
399     } // end event loop
400
401
402     // calculate the pion efficiency and fill the graph
403     util -> CalculatePionEffi(hElecs, hPions);
404     pionEffiTrain[iLoop] = util -> GetPionEfficiency();
405     pionEffiErrTrain[iLoop] = util -> GetError();
406
407     gEffisTrain -> SetPoint(iLoop, iLoop+1, pionEffiTrain[iLoop]);
408     gEffisTrain -> SetPointError(iLoop, 0, pionEffiErrTrain[iLoop]);
409     hElecs -> Reset();
410     hPions -> Reset();
411     AliDebug(2, Form("TrainingLoop[%d] PionEfficiency[%f +/- %f]", iLoop, pionEffiTrain[iLoop], pionEffiErrTrain[iLoop]));
412     // end training loop
413     
414
415
416     // event loop test list
417     for(Int_t iEvent = 0; iEvent < fTest[mombin][0] -> GetN(); iEvent++ ){
418
419       // reset particle probabilities
420       for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
421   likeAll[iPart] = 1./AliTRDgeometry::kNlayer;
422       }
423       totProb = 0.;
424
425       fData -> GetEntry(fTest[mombin][0] -> GetEntry(iEvent));
426       // use event only if it is electron or pion
427       if(!((fPID[AliPID::kElectron] == 1.0) || (fPID[AliPID::kPion] == 1.0))) continue;
428
429       // get the probabilities for each particle type in each chamber
430       for(Int_t iChamb = 0; iChamb < AliTRDgeometry::kNlayer; iChamb++){
431   for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
432     like[iPart][iChamb] = fNet[iChamb] -> Result(fTest[mombin][iChamb] -> GetEntry(iEvent), iPart);
433     likeAll[iPart] *=  like[iPart][iChamb];
434   }
435       }
436
437       // get total probability and normalize it
438       for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
439   totProb += likeAll[iPart];
440       }
441       for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
442   likeAll[iPart] /= totProb;
443       }
444
445       // fill likelihood distributions
446       if(fPID[AliPID::kElectron] == 1)      
447   hElecs -> Fill(likeAll[AliPID::kElectron]);
448       if(fPID[AliPID::kPion] == 1)      
449   hPions -> Fill(likeAll[AliPID::kElectron]);
450     } // end event loop
451
452     // calculate the pion efficiency and fill the graph
453     util -> CalculatePionEffi(hElecs, hPions);
454     pionEffiTest[iLoop] = util -> GetPionEfficiency();
455     pionEffiErrTest[iLoop] = util -> GetError();
456
457     gEffisTest -> SetPoint(iLoop, iLoop+1, pionEffiTest[iLoop]);
458     gEffisTest -> SetPointError(iLoop, 0, pionEffiErrTest[iLoop]);
459     hElecs -> Reset();
460     hPions -> Reset();
461     AliDebug(2, Form("TestLoop[%d] PionEfficiency[%f +/- %f] \n", iLoop, pionEffiTest[iLoop], pionEffiErrTest[iLoop]));
462     
463   } //   end training loop
464
465   util -> Delete();
466
467   gEffisTest -> Draw("PAL");
468   gEffisTrain -> Draw("PL");
469
470 }
471
472
473 //________________________________________________________________________
474 void AliTRDpidRefMakerNN::LoadFile(const Char_t *InFileNN) 
475 {
476   //
477   // Loads the files and sets the event list
478   // for neural network training and 
479   // building of the 2-dim reference histograms.
480   // Useable for training outside of the makeResults.C macro
481   //
482
483   TFile *fInFileNN;
484   fInFileNN = new TFile(InFileNN, "READ");
485   fData = (TTree*)fInFileNN -> Get("NN");
486
487   for(Int_t iMom = 0; iMom < AliTRDCalPID::kNMom; iMom++){
488     for(Int_t ily = 0; ily < AliTRDgeometry::kNlayer; ily++){
489       fTrain[iMom][ily] = new TEventList(Form("fTrainMom%d_%d", iMom, ily), Form("Training list for momentum intervall %d and plane %d", iMom, ily));
490       fTest[iMom][ily] = new TEventList(Form("fTestMom%d_%d", iMom, ily), Form("Test list for momentum intervall %d and plane %d", iMom, ily));
491     }
492   }
493 }
494
495
496 // //________________________________________________________________________
497 // void AliTRDpidRefMakerNN::LoadContainer(const Char_t *InFileCont) 
498 // {
499
500 //   //
501 //   // Loads the container if no container is there.
502 //   // Useable for training outside of the makeResults.C macro
503 //   //
504
505 //   TFile *fInFileCont;
506 //   fInFileCont = new TFile(InFileCont, "READ");
507 //   fContainer = (TObjArray*)fInFileCont -> Get("PidRefMaker");
508
509 // }
510
511
512 // //________________________________________________________________________
513 // void AliTRDpidRefMakerNN::CreateGraphs()
514 // {
515 //   // Create histograms
516 //   // Called once
517
518 //   OpenFile(0, "RECREATE");
519 //   fContainer = new TObjArray();
520 //   fContainer->AddAt(new TH1F("hPDG","hPDG",AliPID::kSPECIES,-0.5,5.5),0);
521
522 //   TGraphErrors *gEffisTrain = new TGraphErrors(kMoniTrain);
523 //   gEffisTrain -> SetLineColor(4);
524 //   gEffisTrain -> SetMarkerColor(4);
525 //   gEffisTrain -> SetMarkerStyle(29);
526 //   gEffisTrain -> SetMarkerSize(2);
527
528 //   TGraphErrors *gEffisTest = new TGraphErrors(kMoniTrain);
529 //   gEffisTest -> SetLineColor(2);
530 //   gEffisTest -> SetMarkerColor(2);
531 //   gEffisTest -> SetMarkerSize(2);
532
533 //   fContainer -> AddAt(gEffisTrain,kGraphTrain);
534 //   fContainer -> AddAt(gEffisTest,kGraphTest);
535 // }
536