]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaRec/AliTRDpidRefMakerNN.cxx
new strategy for TRD-PID ref maker (AlexW & AlexB)
[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   if(fDebugLevel>=2) Printf("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       if(fDebugLevel>=2) Printf("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   if(fDebugLevel>=2) Printf("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   if(fDebugLevel>=2) Printf("\n 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   if(fDebugLevel>=2){ 
190     Printf("Particle multiplicities:");
191     for(Int_t iMomBin = 0; iMomBin <AliTRDCalPID::kNMom; iMomBin++)
192       Printf("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]);
193     Printf("\n");
194   }
195
196   // implement counter of training and test sample size
197   Int_t iTrain[AliTRDCalPID::kNMom], iTest[AliTRDCalPID::kNMom];
198   memset(iTrain, 0, AliTRDCalPID::kNMom*sizeof(Int_t));
199   memset(iTest, 0, AliTRDCalPID::kNMom*sizeof(Int_t));
200
201   // set training sample size per momentum interval to 2/3 
202   // of smallest particle counter and test sample to 1/3
203   for(Int_t iMomBin = 0; iMomBin < AliTRDCalPID::kNMom; iMomBin++){
204     iTrain[iMomBin] = nPart[0][iMomBin];
205     for(Int_t iPart = 1; iPart < AliPID::kSPECIES; iPart++){
206       if(iTrain[iMomBin] > nPart[iPart][iMomBin])
207   iTrain[iMomBin] = nPart[iPart][iMomBin];
208     } 
209     iTrain[iMomBin] = Int_t(iTrain[iMomBin] * .66);
210     iTest[iMomBin] = Int_t( iTrain[iMomBin] * .5);
211     if(fDebugLevel>=2) Printf("Momentum[%d]  Train[%d] Test[%d]", iMomBin, iTrain[iMomBin], iTest[iMomBin]);
212   }
213   if(fDebugLevel>=2) Printf("\n");
214
215
216   // reset couters
217   memset(nPart, 0, AliPID::kSPECIES*AliTRDCalPID::kNMom*sizeof(Int_t));
218
219   // start second loop to set the event lists
220   for(Int_t iEv = 0; iEv < fData -> GetEntries(); iEv++){
221     fData -> GetEntry(iEv);
222
223     // use only events with goes through 6 layers TRD
224     if(fPIDdataArray->fNtracklets != AliTRDgeometry::kNlayer) continue;
225
226     for(Int_t ily=AliTRDgeometry::kNlayer; ily--;){ 
227       Int_t iMomBin = fPIDdataArray->fData[ily].fPLbin & 0xf;
228     
229       // set event list
230       if(nPart[fPIDbin][iMomBin] < iTrain[iMomBin]){
231         fTrain[iMomBin][ily] -> Enter(iEv + ily);
232         nPart[fPIDbin][iMomBin]++;
233       } else if(nPart[fPIDbin][iMomBin] < iTest[iMomBin]+iTrain[iMomBin]){
234         fTest[iMomBin][ily] -> Enter(iEv + ily);
235         nPart[fPIDbin][iMomBin]++;
236       } else continue;
237     }
238   }
239   
240   if(fDebugLevel>=2){ 
241     Printf("Particle multiplicities in both lists:");
242     for(Int_t iMomBin = 0; iMomBin <AliTRDCalPID::kNMom; iMomBin++)
243       Printf("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]);
244     Printf("\n");
245   }
246 }
247
248
249 //________________________________________________________________________
250 void AliTRDpidRefMakerNN::MakeRefs(Int_t mombin) 
251 {
252   //
253   // train the neural networks
254   //
255   
256   
257   if (!fData) LoadFile(Form("TRD.Calib%s.root", GetName()));
258
259   if (!fData) {
260     AliError("Tree for training list not available");
261     return;
262   }
263
264   TDatime datime;
265   fDate = datime.GetDate();
266
267   AliDebug(2, Form("Training momentum bin %d", mombin));
268
269   // set variable to monitor the training and to save the development of the networks
270   Int_t nEpochs = fEpochs/kMoniTrain;       
271   AliDebug(2, Form("Training %d times %d epochs", kMoniTrain, nEpochs));
272
273   // make directories to save the networks 
274   gSystem->Exec(Form("rm -r ./Networks_%d/MomBin_%d",fDate, mombin));
275   gSystem->Exec(Form("mkdir ./Networks_%d/MomBin_%d",fDate, mombin));
276
277   // variable to check if network can load weights from previous training
278   Bool_t bFirstLoop[AliTRDgeometry::kNlayer];
279   memset(bFirstLoop, kTRUE, AliTRDgeometry::kNlayer*sizeof(Bool_t));
280
281   // train networks over several loops and save them after each loop
282   for(Int_t iLoop = 0; iLoop < kMoniTrain; iLoop++){
283     // loop over chambers
284     for(Int_t iChamb = 0; iChamb < AliTRDgeometry::kNlayer; iChamb++){
285       // set the event lists
286       fData -> SetEventList(fTrain[mombin][iChamb]);
287       fData -> SetEventList(fTest[mombin][iChamb]);
288       
289       AliDebug(2, Form("Trainingloop[%d] Chamber[%d]", iLoop, iChamb));
290       
291       // check if network is already implemented
292       if(bFirstLoop[iChamb] == kTRUE){
293   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]);
294   fNet[iChamb] -> SetLearningMethod(TMultiLayerPerceptron::kStochastic);       // set learning method
295   fNet[iChamb] -> TMultiLayerPerceptron::SetEta(0.001);                        // set learning speed
296   if(!fContinueTraining){
297     if(fDebugLevel>=2) fNet[iChamb] -> Train(nEpochs,"text update=10, graph");
298     else fNet[iChamb] -> Train(nEpochs,"");
299   }
300   else{
301     fNet[iChamb] -> LoadWeights(Form("./Networks_%d/MomBin_%d/Net%d_%d",fTrainPath, mombin, iChamb, kMoniTrain - 1));
302     if(fDebugLevel>=2) fNet[iChamb] -> Train(nEpochs,"text update=10, graph+");      
303     else fNet[iChamb] -> Train(nEpochs,"+");                   
304   }
305   bFirstLoop[iChamb] = kFALSE;
306       }
307       else{    
308   if(fDebugLevel>=2) fNet[iChamb] -> Train(nEpochs,"text update=10, graph+");      
309   else fNet[iChamb] -> Train(nEpochs,"+");                   
310       }
311       
312       // save weights for monitoring of the training
313       fNet[iChamb] -> DumpWeights(Form("./Networks_%d/MomBin_%d/Net%d_%d",fDate, mombin, iChamb, iLoop));
314     } // end chamber loop
315   }   // end training loop
316 }
317
318
319
320 //________________________________________________________________________
321 void AliTRDpidRefMakerNN::MonitorTraining(Int_t mombin) 
322 {
323   //
324   // train the neural networks
325   //
326   
327   if (!fData) LoadFile(Form("TRD.Calib%s.root", GetName()));
328   if (!fData) {
329     AliError("Tree for training list not available");
330     return;
331   }
332
333   // init networks and set event list
334   for(Int_t iChamb = 0; iChamb < AliTRDgeometry::kNlayer; iChamb++){
335   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]);   
336   fData -> SetEventList(fTrain[mombin][iChamb]);
337   fData -> SetEventList(fTest[mombin][iChamb]);
338   }
339
340   // implement variables for likelihoods
341   Float_t like[AliPID::kSPECIES][AliTRDgeometry::kNlayer];
342   memset(like, 0, AliPID::kSPECIES*AliTRDgeometry::kNlayer*sizeof(Float_t));
343   Float_t likeAll[AliPID::kSPECIES], totProb;
344
345   Double_t pionEffiTrain[kMoniTrain], pionEffiErrTrain[kMoniTrain];
346   Double_t pionEffiTest[kMoniTrain], pionEffiErrTest[kMoniTrain];
347   memset(pionEffiTrain, 0, kMoniTrain*sizeof(Double_t));
348   memset(pionEffiErrTrain, 0, kMoniTrain*sizeof(Double_t));
349   memset(pionEffiTest, 0, kMoniTrain*sizeof(Double_t));
350   memset(pionEffiErrTest, 0, kMoniTrain*sizeof(Double_t));
351
352   // init histos
353   const Float_t epsilon = 1/(2*(AliTRDpidUtil::kBins-1));     // get nice histos with bin center at 0 and 1
354   TH1F *hElecs, *hPions;
355   hElecs = new TH1F("hElecs","Likelihood for electrons", AliTRDpidUtil::kBins, 0.-epsilon, 1.+epsilon);
356   hPions = new TH1F("hPions","Likelihood for pions", AliTRDpidUtil::kBins, 0.-epsilon, 1.+epsilon);
357
358   TGraphErrors *gEffisTrain=0x0, *gEffisTest=0x0;
359   gEffisTrain = (TGraphErrors*)fContainer->At(kGraphTrain);
360   gEffisTest = (TGraphErrors*)fContainer->At(kGraphTest);
361
362   AliTRDpidUtil *util = new AliTRDpidUtil();
363   
364   // monitor training progress
365   for(Int_t iLoop = 0; iLoop < kMoniTrain; iLoop++){
366
367     // load weights
368     for(Int_t iChamb = 0; iChamb < AliTRDgeometry::kNlayer; iChamb++){
369       fNet[iChamb] -> LoadWeights(Form("./Networks_%d/MomBin_%d/Net%d_%d",fDate, mombin, iChamb, iLoop));
370     }
371
372     // event loop training list
373     for(Int_t iEvent = 0; iEvent < fTrain[mombin][0] -> GetN(); iEvent++ ){
374
375       // reset particle probabilities
376       for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
377   likeAll[iPart] = 1./AliPID::kSPECIES;
378       }
379       totProb = 0.;
380
381       fData -> GetEntry(fTrain[mombin][0] -> GetEntry(iEvent));
382       // use event only if it is electron or pion
383       if(!((fPID[AliPID::kElectron] == 1.0) || (fPID[AliPID::kPion] == 1.0))) continue;
384
385       // get the probabilities for each particle type in each chamber
386       for(Int_t iChamb = 0; iChamb < AliTRDgeometry::kNlayer; iChamb++){
387   for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
388     like[iPart][iChamb] = fNet[iChamb] -> Result(fTrain[mombin][iChamb] -> GetEntry(iEvent), iPart);
389     likeAll[iPart] *=  like[iPart][iChamb];
390   }
391       }
392
393       // get total probability and normalize it
394       for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
395   totProb += likeAll[iPart];
396       }
397       for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
398   likeAll[iPart] /= totProb;
399       }
400
401       // fill likelihood distributions
402       if(fPID[AliPID::kElectron] == 1)      
403   hElecs -> Fill(likeAll[AliPID::kElectron]);
404       if(fPID[AliPID::kPion] == 1)      
405   hPions -> Fill(likeAll[AliPID::kElectron]);
406     } // end event loop
407
408
409     // calculate the pion efficiency and fill the graph
410     util -> CalculatePionEffi(hElecs, hPions);
411     pionEffiTrain[iLoop] = util -> GetPionEfficiency();
412     pionEffiErrTrain[iLoop] = util -> GetError();
413
414     gEffisTrain -> SetPoint(iLoop, iLoop+1, pionEffiTrain[iLoop]);
415     gEffisTrain -> SetPointError(iLoop, 0, pionEffiErrTrain[iLoop]);
416     hElecs -> Reset();
417     hPions -> Reset();
418     if(fDebugLevel>=2) Printf("TrainingLoop[%d] PionEfficiency[%f +/- %f]", iLoop, pionEffiTrain[iLoop], pionEffiErrTrain[iLoop]);
419     // end training loop
420     
421
422
423     // event loop test list
424     for(Int_t iEvent = 0; iEvent < fTest[mombin][0] -> GetN(); iEvent++ ){
425
426       // reset particle probabilities
427       for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
428   likeAll[iPart] = 1./AliTRDgeometry::kNlayer;
429       }
430       totProb = 0.;
431
432       fData -> GetEntry(fTest[mombin][0] -> GetEntry(iEvent));
433       // use event only if it is electron or pion
434       if(!((fPID[AliPID::kElectron] == 1.0) || (fPID[AliPID::kPion] == 1.0))) continue;
435
436       // get the probabilities for each particle type in each chamber
437       for(Int_t iChamb = 0; iChamb < AliTRDgeometry::kNlayer; iChamb++){
438   for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
439     like[iPart][iChamb] = fNet[iChamb] -> Result(fTest[mombin][iChamb] -> GetEntry(iEvent), iPart);
440     likeAll[iPart] *=  like[iPart][iChamb];
441   }
442       }
443
444       // get total probability and normalize it
445       for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
446   totProb += likeAll[iPart];
447       }
448       for(Int_t iPart = 0; iPart < AliPID::kSPECIES; iPart++){
449   likeAll[iPart] /= totProb;
450       }
451
452       // fill likelihood distributions
453       if(fPID[AliPID::kElectron] == 1)      
454   hElecs -> Fill(likeAll[AliPID::kElectron]);
455       if(fPID[AliPID::kPion] == 1)      
456   hPions -> Fill(likeAll[AliPID::kElectron]);
457     } // end event loop
458
459     // calculate the pion efficiency and fill the graph
460     util -> CalculatePionEffi(hElecs, hPions);
461     pionEffiTest[iLoop] = util -> GetPionEfficiency();
462     pionEffiErrTest[iLoop] = util -> GetError();
463
464     gEffisTest -> SetPoint(iLoop, iLoop+1, pionEffiTest[iLoop]);
465     gEffisTest -> SetPointError(iLoop, 0, pionEffiErrTest[iLoop]);
466     hElecs -> Reset();
467     hPions -> Reset();
468     if(fDebugLevel>=2) Printf("TestLoop[%d] PionEfficiency[%f +/- %f] \n", iLoop, pionEffiTest[iLoop], pionEffiErrTest[iLoop]);
469     
470   } //   end training loop
471
472   util -> Delete();
473
474   gEffisTest -> Draw("PAL");
475   gEffisTrain -> Draw("PL");
476
477 }
478
479
480 //________________________________________________________________________
481 void AliTRDpidRefMakerNN::LoadFile(const Char_t *InFileNN) 
482 {
483   //
484   // Loads the files and sets the event list
485   // for neural network training and 
486   // building of the 2-dim reference histograms.
487   // Useable for training outside of the makeResults.C macro
488   //
489
490   TFile *fInFileNN;
491   fInFileNN = new TFile(InFileNN, "READ");
492   fData = (TTree*)fInFileNN -> Get("NN");
493
494   for(Int_t iMom = 0; iMom < AliTRDCalPID::kNMom; iMom++){
495     for(Int_t ily = 0; ily < AliTRDgeometry::kNlayer; ily++){
496       fTrain[iMom][ily] = new TEventList(Form("fTrainMom%d_%d", iMom, ily), Form("Training list for momentum intervall %d and plane %d", iMom, ily));
497       fTest[iMom][ily] = new TEventList(Form("fTestMom%d_%d", iMom, ily), Form("Test list for momentum intervall %d and plane %d", iMom, ily));
498     }
499   }
500 }
501
502
503 // //________________________________________________________________________
504 // void AliTRDpidRefMakerNN::LoadContainer(const Char_t *InFileCont) 
505 // {
506
507 //   //
508 //   // Loads the container if no container is there.
509 //   // Useable for training outside of the makeResults.C macro
510 //   //
511
512 //   TFile *fInFileCont;
513 //   fInFileCont = new TFile(InFileCont, "READ");
514 //   fContainer = (TObjArray*)fInFileCont -> Get("PidRefMaker");
515
516 // }
517
518
519 // //________________________________________________________________________
520 // void AliTRDpidRefMakerNN::CreateGraphs()
521 // {
522 //   // Create histograms
523 //   // Called once
524
525 //   OpenFile(0, "RECREATE");
526 //   fContainer = new TObjArray();
527 //   fContainer->AddAt(new TH1F("hPDG","hPDG",AliPID::kSPECIES,-0.5,5.5),0);
528
529 //   TGraphErrors *gEffisTrain = new TGraphErrors(kMoniTrain);
530 //   gEffisTrain -> SetLineColor(4);
531 //   gEffisTrain -> SetMarkerColor(4);
532 //   gEffisTrain -> SetMarkerStyle(29);
533 //   gEffisTrain -> SetMarkerSize(2);
534
535 //   TGraphErrors *gEffisTest = new TGraphErrors(kMoniTrain);
536 //   gEffisTest -> SetLineColor(2);
537 //   gEffisTest -> SetMarkerColor(2);
538 //   gEffisTest -> SetMarkerSize(2);
539
540 //   fContainer -> AddAt(gEffisTrain,kGraphTrain);
541 //   fContainer -> AddAt(gEffisTest,kGraphTest);
542 // }
543