Included the weight to hadrons as a function of their reconstructed energy and a...
[u/mrichter/AliRoot.git] / PHOS / AliPHOSPIDv1.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-commercial 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$ */
17
18 //_________________________________________________________________________
19 // Implementation version v1 of the PHOS particle identifier 
20 // Particle identification based on the 
21 //     - RCPV: distance from CPV recpoint to EMCA recpoint.
22 //     - TOF 
23 //     - PCA: Principal Components Analysis..
24 // The identified particle has an identification number corresponding 
25 // to a 9 bits number:
26 //     -Bit 0 to 2: bit set if RCPV > CpvEmcDistance (each bit corresponds
27 //      to a different efficiency-purity point of the photon identification) 
28 //     -Bit 3 to 5: bit set if TOF  < TimeGate (each bit corresponds
29 //      to a different efficiency-purity point of the photon identification) 
30 //     -Bit 6 to 9: bit set if Principal Components are 
31 //      inside an ellipse defined by the parameters a, b, c, x0 and y0.
32 //      (each bit corresponds to a different efficiency-purity point of the 
33 //      photon identification)
34 //      The PCA (Principal components analysis) needs a file that contains
35 //      a previous analysis of the correlations between the particles. This 
36 //      file is $ALICE_ROOT/PHOS/PCA8pa15_0.5-100.root. Analysis done for 
37 //      energies between 0.5 and 100 GeV.
38 //      A calibrated energy is calculated. The energy of the reconstructed
39 //      cluster is corrected with the formula A + B * E  + C * E^2, whose 
40 //      parameters where obtained through the study of the reconstructed 
41 //      energy distribution of monoenergetic photons. 
42 //
43 //      All the parameters (RCPV(2 rows-3 columns),TOF(1r-3c),PCA(5r-4c) 
44 //      and calibration(1r-3c))are stored in a file called 
45 //      $ALICE_ROOT/PHOS/Parameters.dat. Each time that AliPHOSPIDv1 is 
46 //      initialized, this parameters are copied to a Matrix (9,4), a 
47 //      TMatrixD object.  
48 //
49 // use case:
50 //  root [0] AliPHOSPIDv1 * p = new AliPHOSPIDv1("galice1.root")
51 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
52 //          // reading headers from file galice1.root and create  RecParticles 
53             // TrackSegments and RecPoints are used 
54 //          // set file name for the branch RecParticles
55 //  root [1] p->ExecuteTask("deb all time")
56 //          // available options
57 //          // "deb" - prints # of reconstructed particles
58 //          // "deb all" -  prints # and list of RecParticles
59 //          // "time" - prints benchmarking results
60 //                  
61 //  root [2] AliPHOSPIDv1 * p2 = new AliPHOSPIDv1("galice1.root","v1",kTRUE)
62 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
63 //                //Split mode.  
64 //  root [3] p2->ExecuteTask()
65 //
66
67
68 //*-- Author: Yves Schutz (SUBATECH)  & Gines Martinez (SUBATECH) & 
69 //            Gustavo Conesa April 2002
70 //            PCA redesigned by Gustavo Conesa October 2002:
71 //            The way of using the PCA has changed. Instead of 2
72 //            files with the PCA, each one with different energy ranges 
73 //            of application, we use the wide one (0.5-100 GeV), and instead
74 //            of fixing 3 ellipses for different ranges of energy, it has been
75 //            studied the dependency of the ellipses parameters with the 
76 //            energy, and they are implemented in the code as a funtion 
77 //            of the energy. 
78 //
79 //
80 //
81 // --- ROOT system ---
82
83
84 // --- Standard library ---
85 #include "TFormula.h"
86 #include "TBenchmark.h"
87 #include "TPrincipal.h"
88 #include "TFile.h" 
89 #include "TSystem.h"
90
91 // --- AliRoot header files ---
92               //#include "AliLog.h"
93 #include "AliGenerator.h"
94 #include "AliPHOS.h"
95 #include "AliPHOSPIDv1.h"
96 #include "AliPHOSGetter.h"
97
98 ClassImp( AliPHOSPIDv1) 
99
100 //____________________________________________________________________________
101 AliPHOSPIDv1::AliPHOSPIDv1():AliPHOSPID()
102
103   // default ctor
104  
105   InitParameters() ; 
106   fDefaultInit = kTRUE ; 
107 }
108
109 //____________________________________________________________________________
110 AliPHOSPIDv1::AliPHOSPIDv1(const AliPHOSPIDv1 & pid ):AliPHOSPID(pid)
111
112   // ctor
113   InitParameters() ; 
114   Init() ;
115
116 }
117
118 //____________________________________________________________________________
119 AliPHOSPIDv1::AliPHOSPIDv1(const TString alirunFileName, const TString eventFolderName):AliPHOSPID(alirunFileName, eventFolderName)
120
121   //ctor with the indication on where to look for the track segments
122  
123   InitParameters() ; 
124   Init() ;
125   fDefaultInit = kFALSE ; 
126 }
127
128 //____________________________________________________________________________
129 AliPHOSPIDv1::~AliPHOSPIDv1()
130
131   // dtor
132   fPrincipalPhoton = 0;
133   fPrincipalPi0 = 0;
134
135   delete [] fX ;       // Principal input 
136   delete [] fPPhoton ; // Photon Principal components
137   delete [] fPPi0 ;    // Pi0 Principal components
138
139   delete fParameters;
140   delete fTFphoton;
141   delete fTFpiong;
142   delete fTFkaong;
143   delete fTFkaonl;
144   delete fTFhhadrong;
145   delete fTFhhadronl;
146   delete fDFmuon;
147 }
148 //____________________________________________________________________________
149 const TString AliPHOSPIDv1::BranchName() const 
150 {  
151
152   return GetName() ;
153 }
154  
155 //____________________________________________________________________________
156 void AliPHOSPIDv1::Init()
157 {
158   // Make all memory allocations that are not possible in default constructor
159   // Add the PID task to the list of PHOS tasks
160
161   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
162   if(!gime)
163     gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data()) ; 
164
165   if ( !gime->PID() ) 
166     gime->PostPID(this) ;
167 }
168
169 //____________________________________________________________________________
170 void AliPHOSPIDv1::InitParameters()
171 {
172   // Initialize PID parameters
173   fWrite                   = kTRUE ;
174   fRecParticlesInRun = 0 ; 
175   fNEvent            = 0 ;            
176   fRecParticlesInRun = 0 ;
177   fBayesian          = kTRUE ;
178   SetParameters() ; // fill the parameters matrix from parameters file
179   SetEventRange(0,-1) ;
180
181   // initialisation of response function parameters
182   // Tof
183   // Photons
184   fTphoton[0] = 0.218    ;
185   //fTphoton[0] = 1.    ;
186   fTphoton[1] = 1.55E-8  ; 
187   fTphoton[2] = 5.05E-10 ;
188   fTFphoton = new TFormula("ToF response to photons" , "gaus") ; 
189   fTFphoton->SetParameters( fTphoton[0], fTphoton[1], fTphoton[2]) ; 
190 //   // Electrons
191 //   fTelectron[0] = 0.2      ;
192 //   fTelectron[1] = 1.55E-8  ; 
193 //   fTelectron[2] = 5.35E-10 ;
194 //   fTFelectron = new TFormula("ToF response to electrons" , "gaus") ; 
195 //   fTFelectron->SetParameters( fTelectron[0], fTelectron[1], fTelectron[2]) ; 
196 //   // Muons
197 //   fTmuon[0] = 0.2     ;
198 //   fTmuon[1] = 1.55E-8 ; 
199 //   fTmuon[2] = 5.1E-10 ;
200 //   fTFmuon = new TFormula("ToF response to muons" , "gaus") ; 
201 //   fTFmuon->SetParameters( fTmuon[0], fTmuon[1], fTmuon[2]) ; 
202
203   // Pions
204   //Gaus (0 to max probability)
205   fTpiong[0] = 0.0971    ; 
206   //fTpiong[0] = 1.       ;
207   fTpiong[1] = 1.58E-8  ; 
208   fTpiong[2] = 5.69E-10 ;
209   fTFpiong = new TFormula("ToF response to pions" , "gaus") ; 
210   fTFpiong->SetParameters( fTpiong[0], fTpiong[1], fTpiong[2]) ; 
211   // Landau (max probability to inf) 
212 //   fTpionl[0] = 0.05     ; 
213 //   //fTpionl[0] = 5.53     ; 
214 //   fTpionl[1] = 1.68E-8  ; 
215 //   fTpionl[2] = 5.38E-10 ;
216 //   fTFpionl = new TFormula("ToF response to pions" , "landau") ; 
217 //   fTFpionl->SetParameters( fTpionl[0], fTpionl[1], fTpionl[2]) ; 
218
219
220   // Kaons
221   //Gaus (0 to max probability)
222   fTkaong[0] = 0.0542  ; 
223   //fTkaong[0] = 1.      ;
224   fTkaong[1] = 1.64E-8 ; 
225   fTkaong[2] = 6.07-10 ;
226   fTFkaong = new TFormula("ToF response to kaon" , "gaus") ; 
227   fTFkaong->SetParameters( fTkaong[0], fTkaong[1], fTkaong[2]) ; 
228   //Landau (max probability to inf) 
229   fTkaonl[0] = 0.264   ;
230   //fTkaonl[0] = 5.53     ;
231   fTkaonl[1] = 1.68E-8  ; 
232   fTkaonl[2] = 4.10E-10 ;
233   fTFkaonl = new TFormula("ToF response to kaon" , "landau") ; 
234   fTFkaonl->SetParameters( fTkaonl[0], fTkaonl[1], fTkaonl[2]) ; 
235
236   //Heavy Hadrons
237   //Gaus (0 to max probability)
238   fThhadrong[0] = 0.0302   ;  
239   //fThhadrong[0] = 1.       ; 
240   fThhadrong[1] = 1.73E-8  ; 
241   fThhadrong[2] = 9.52E-10 ;
242   fTFhhadrong = new TFormula("ToF response to heavy hadrons" , "gaus") ; 
243   fTFhhadrong->SetParameters( fThhadrong[0], fThhadrong[1], fThhadrong[2]) ; 
244   //Landau (max probability to inf) 
245   fThhadronl[0] = 0.139    ;  
246   //fThhadronl[0] = 5.53      ;   
247   fThhadronl[1] = 1.745E-8  ; 
248   fThhadronl[2] = 1.00E-9  ;
249   fTFhhadronl = new TFormula("ToF response to heavy hadrons" , "landau") ; 
250   fTFhhadronl->SetParameters( fThhadronl[0], fThhadronl[1], fThhadronl[2]) ; 
251
252 /// /gaussian parametrization for pions
253 //   fTpion[0] = 3.93E-2 ;  fTpion[1] = 0.130   ; fTpion[2] =-6.37E-2 ;//constant
254 //   fTpion[3] = 1.65E-8 ;  fTpion[4] =-1.40E-9 ; fTpion[5] = 5.96E-10;//mean
255 //   fTpion[6] = 8.09E-10;  fTpion[7] =-4.65E-10; fTpion[8] = 1.50E-10;//sigma
256
257 // //landau parametrization for kaons
258 //   fTkaon[0] = 0.107   ;  fTkaon[1] = 0.166   ; fTkaon[2] = 0.243   ;//constant
259 //   fTkaon[3] = 1.80E-8 ;  fTkaon[4] =-2.96E-9 ; fTkaon[5] = 9.60E-10;//mean
260 //   fTkaon[6] = 1.37E-9 ;  fTkaon[7] =-1.80E-9 ; fTkaon[8] = 6.74E-10;//sigma
261
262 // //landau parametrization for nucleons
263 //   fThhadron[0] = 6.33E-2 ;  fThhadron[1] = 2.52E-2 ; fThhadron[2] = 2.16E-2 ;//constant
264 //   fThhadron[3] = 1.94E-8 ;  fThhadron[4] =-7.06E-10; fThhadron[5] =-4.69E-10;//mean
265 //   fThhadron[6] = 2.55E-9 ;  fThhadron[7] =-1.90E-9 ; fThhadron[8] = 5.41E-10;//sigma
266
267
268   // Shower shape: dispersion gaussian parameters
269   // Photons
270
271 //   fDphoton[0] = 3.84e-2;  fDphoton[1] = 4.46e-3 ; fDphoton[2] = -2.36e-2;//constant
272 //   //fDphoton[0] = 1.0    ;  fDphoton[1] = 0.      ; fDphoton[2] = 0.     ;//constant
273 //   fDphoton[3] = 1.55   ;  fDphoton[4] =-0.0863  ; fDphoton[5] = 0.287   ;//mean
274 //   fDphoton[6] = 0.0451 ;  fDphoton[7] =-0.0803  ; fDphoton[8] = 0.314   ;//sigma
275
276    fDphoton[0] = 4.62e-2;  fDphoton[1] = 1.39e-2 ; fDphoton[2] = -3.80e-2;//constant
277    //fDphoton[0] = 1.0    ;  fDphoton[1] = 0.      ; fDphoton[2] = 0.     ;//constant
278    fDphoton[3] = 1.53   ;  fDphoton[4] =-6.62e-2 ; fDphoton[5] = 0.339   ;//mean
279    fDphoton[6] = 6.89e-2;  fDphoton[7] =-6.59e-2 ; fDphoton[8] = 0.194   ;//sigma
280
281    fDpi0[0] = 0.0586  ;  fDpi0[1] = 1.06E-3 ; fDpi0[2] = 0.      ;//constant
282    //fDpi0[0] = 1.0     ;  fDpi0[1] = 0.0     ; fDpi0[2] = 0.      ;//constant
283   fDpi0[3] = 2.67    ;  fDpi0[4] =-2.00E-2 ; fDpi0[5] = 9.37E-5 ;//mean
284   fDpi0[6] = 0.153   ;  fDpi0[7] = 9.34E-4 ; fDpi0[8] =-1.49E-5 ;//sigma
285   //landau
286 //   fDhadron[0] = 0.007  ;  fDhadron[1] = 0.      ; fDhadron[2] = 0.    ;//constant
287 //   //fDhadron[0] = 5.53   ;  fDhadron[1] = 0.      ; fDhadron[2] = 0.    ;//constant
288 //   fDhadron[3] = 3.38   ;  fDhadron[4] = 0.0833  ; fDhadron[5] =-0.845 ;//mean
289 //   fDhadron[6] = 0.627  ;  fDhadron[7] = 0.012   ; fDhadron[8] =-0.170 ;//sigma
290
291   fDhadron[0] = 1.61E-2 ;  fDhadron[1] = 3.03E-3 ; fDhadron[2] = 1.01E-2 ;//constant
292   fDhadron[3] = 3.81    ;  fDhadron[4] = 0.232   ; fDhadron[5] =-1.25    ;//mean
293   fDhadron[6] = 0.897   ;  fDhadron[7] = 0.0987  ; fDhadron[8] =-0.534   ;//sigma
294   // Muons
295   fDmuon[0] = 0.0631     ;
296   fDmuon[1] = 1.4 ; 
297   fDmuon[2] = 0.0557 ;
298   fDFmuon = new TFormula("Shower shape response to muons" , "landau") ; 
299   fDFmuon->SetParameters( fDmuon[0], fDmuon[1], fDmuon[2]) ; 
300
301
302   // x(CPV-EMC) distance gaussian parameters
303   
304   fXelectron[0] = 8.06e-2 ;  fXelectron[1] = 1.00e-2; fXelectron[2] =-5.14e-2;//constant
305   //fXelectron[0] = 1.0     ;  fXelectron[1] = 0.     ; fXelectron[2] = 0.    ;//constant
306   fXelectron[3] = 0.202   ;  fXelectron[4] = 8.15e-3; fXelectron[5] = 4.55   ;//mean
307   fXelectron[6] = 0.334   ;  fXelectron[7] = 0.186  ; fXelectron[8] = 4.32e-2;//sigma
308   
309   //charged hadrons gaus
310   fXcharged[0] = 6.43e-3 ;  fXcharged[1] =-4.19e-5; fXcharged[2] = 1.42e-3;//constant
311   fXcharged[3] = 2.75    ;  fXcharged[4] =-0.40   ; fXcharged[5] = 1.68   ;//mean
312   fXcharged[6] = 3.135   ;  fXcharged[7] =-9.41e-2; fXcharged[8] = 1.31e-2;//sigma
313   
314   // z(CPV-EMC) distance gaussian parameters
315   
316   fZelectron[0] = 8.22e-2 ;  fZelectron[1] = 5.11e-3; fZelectron[2] =-3.05e-2;//constant
317   //fZelectron[0] = 1.0     ;  fZelectron[1] = 0.     ; fZelectron[2] = 0.    ;//constant
318   fZelectron[3] = 3.09e-2 ;  fZelectron[4] = 5.87e-2; fZelectron[5] =-9.49e-2;//mean
319   fZelectron[6] = 0.263   ;  fZelectron[7] =-9.02e-3; fZelectron[8] = 0.151 ;//sigma
320   
321   //charged hadrons gaus
322   
323   fZcharged[0] = 1.00e-2 ;  fZcharged[1] = 2.82E-4 ; fZcharged[2] = 2.87E-3 ;//constant
324   fZcharged[3] =-4.68e-2 ;  fZcharged[4] =-9.21e-3 ; fZcharged[5] = 4.91e-2 ;//mean
325   fZcharged[6] = 1.425   ;  fZcharged[7] =-5.90e-2 ; fZcharged[8] = 5.07e-2 ;//sigma
326   
327   //Threshold to differentiate between charged and neutral
328   fChargedNeutralThreshold = 1e-5;
329
330   //Weight to hadrons recontructed energy
331
332   fERecWeightPar[0] = 0.32 ; 
333   fERecWeightPar[1] = 3.8  ;
334   fERecWeightPar[2] = 5.4E-3 ; 
335   fERecWeightPar[3] = 5.6E-2 ;
336   fERecWeight = new TFormula("Weight for hadrons" , "[0]*exp(-x*[1])+[2]*exp(-x*[3])") ; 
337   fERecWeight ->SetParameters(fERecWeightPar[0],fERecWeightPar[1] ,fERecWeightPar[2] ,fERecWeightPar[3]) ; 
338
339
340   for (Int_t i =0; i<  AliPID::kSPECIESN ; i++)
341     fInitPID[i] = 1.;
342  
343 }
344
345 //________________________________________________________________________
346 void  AliPHOSPIDv1::Exec(Option_t *option)
347 {
348   // Steering method to perform particle reconstruction and identification
349   // for the event range from fFirstEvent to fLastEvent.
350   // This range is optionally set by SetEventRange().
351   // if fLastEvent=-1 (by default), then process events until the end.
352   
353   if(strstr(option,"tim"))
354     gBenchmark->Start("PHOSPID");
355   
356   if(strstr(option,"print")) {
357     Print() ; 
358     return ; 
359   }
360
361
362   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
363  
364   if (fLastEvent == -1) 
365     fLastEvent = gime->MaxEvent() - 1 ;
366   else 
367     fLastEvent = TMath::Min(fLastEvent,gime->MaxEvent());
368   Int_t nEvents   = fLastEvent - fFirstEvent + 1;
369
370   Int_t ievent ; 
371   for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
372     gime->Event(ievent,"TR") ;
373     if(gime->TrackSegments() && //Skip events, where no track segments made
374        gime->TrackSegments()->GetEntriesFast()) {
375
376       MakeRecParticles() ;
377
378       if(fBayesian)
379         MakePID() ; 
380       
381       WriteRecParticles();
382       if(strstr(option,"deb"))
383         PrintRecParticles(option) ;
384       //increment the total number of rec particles per run 
385       fRecParticlesInRun += gime->RecParticles()->GetEntriesFast() ; 
386     }
387   }
388   if(strstr(option,"deb"))
389       PrintRecParticles(option);
390   if(strstr(option,"tim")){
391     gBenchmark->Stop("PHOSPID");
392     AliInfo(Form("took %f seconds for PID %f seconds per event", 
393          gBenchmark->GetCpuTime("PHOSPID"),  
394          gBenchmark->GetCpuTime("PHOSPID")/nEvents)) ;
395   }
396   if(fWrite)
397     Unload();
398 }
399
400 //________________________________________________________________________
401 Double_t  AliPHOSPIDv1::GausF(Double_t  x, Double_t  y, Double_t * par)
402 {
403   //Given the energy x and the parameter y (tof, shower dispersion or cpv-emc distance), 
404   //this method returns a density probability of this parameter, given by a gaussian 
405   //function whose parameters depend with the energy  with a function: a/(x*x)+b/x+b
406   Double_t cnt    = par[1] / (x*x) + par[2] / x + par[0] ;
407   Double_t mean   = par[4] / (x*x) + par[5] / x + par[3] ;
408   Double_t sigma  = par[7] / (x*x) + par[8] / x + par[6] ;
409  
410   //  Double_t arg    = - (y-mean) * (y-mean) / (2*sigma*sigma) ;
411   //  return cnt * TMath::Exp(arg) ;
412   if(TMath::Abs(sigma) > 1.e-10){
413     return cnt*TMath::Gaus(y,mean,sigma);
414   }
415   else
416     return 0.;
417  
418 }
419 //________________________________________________________________________
420 Double_t  AliPHOSPIDv1::GausPol2(Double_t  x, Double_t y, Double_t * par)
421 {
422   //Given the energy x and the parameter y (tof, shower dispersion or cpv-emc distance), 
423   //this method returns a density probability of this parameter, given by a gaussian 
424   //function whose parameters depend with the energy like second order polinomial
425
426   Double_t cnt    = par[0] + par[1] * x + par[2] * x * x ;
427   Double_t mean   = par[3] + par[4] * x + par[5] * x * x ;
428   Double_t sigma  = par[6] + par[7] * x + par[8] * x * x ;
429
430   if(TMath::Abs(sigma) > 1.e-10){
431     return cnt*TMath::Gaus(y,mean,sigma);
432   }
433   else
434     return 0.;
435  
436
437
438 }
439
440 //____________________________________________________________________________
441 const TString AliPHOSPIDv1::GetFileNamePrincipal(TString particle) const
442 {
443   //Get file name that contains the PCA for a particle ("photon or pi0")
444   particle.ToLower();
445   TString name;
446   if      (particle=="photon") 
447     name = fFileNamePrincipalPhoton ;
448   else if (particle=="pi0"   ) 
449     name = fFileNamePrincipalPi0    ;
450   else    
451     AliError(Form("Wrong particle name: %s (choose from pi0/photon)\n",
452                   particle.Data()));
453   return name;
454 }
455
456 //____________________________________________________________________________
457 Float_t  AliPHOSPIDv1::GetParameterCalibration(Int_t i) const 
458 {
459   // Get the i-th parameter "Calibration"
460   Float_t param = 0.;
461   if (i>2 || i<0) { 
462     AliError(Form("Invalid parameter number: %d",i));
463   } else
464     param = (*fParameters)(0,i);
465   return param;
466 }
467
468 //____________________________________________________________________________
469 Float_t  AliPHOSPIDv1::GetCalibratedEnergy(Float_t e) const
470 {
471 //      It calibrates Energy depending on the recpoint energy.
472 //      The energy of the reconstructed cluster is corrected with 
473 //      the formula A + B* E  + C* E^2, whose parameters where obtained 
474 //      through the study of the reconstructed energy distribution of 
475 //      monoenergetic photons.
476  
477   Float_t p[]={0.,0.,0.};
478   for (Int_t i=0; i<3; i++) p[i] = GetParameterCalibration(i);
479   Float_t enerec = p[0] +  p[1]*e + p[2]*e*e;
480   return enerec ;
481
482 }
483
484 //____________________________________________________________________________
485 Float_t  AliPHOSPIDv1::GetParameterCpv2Emc(Int_t i, TString axis) const 
486 {
487   // Get the i-th parameter "CPV-EMC distance" for the specified axis
488   Float_t param = 0.;
489   if(i>2 || i<0) {
490     AliError(Form("Invalid parameter number: %d",i));
491   } else {
492     axis.ToLower();
493     if      (axis == "x") 
494       param = (*fParameters)(1,i);
495     else if (axis == "z") 
496       param = (*fParameters)(2,i);
497     else { 
498       AliError(Form("Invalid axis name: %s",axis.Data()));
499     }
500   }
501   return  param;
502 }
503
504 //____________________________________________________________________________
505 Float_t  AliPHOSPIDv1::GetCpv2EmcDistanceCut(TString axis, Float_t e) const
506 {
507   // Get CpvtoEmcDistance Cut depending on the cluster energy, axis and 
508   // Purity-Efficiency point 
509
510   axis.ToLower();
511   Float_t p[]={0.,0.,0.};
512   for (Int_t i=0; i<3; i++) p[i] = GetParameterCpv2Emc(i,axis);
513   Float_t sig = p[0] + TMath::Exp(p[1] - p[2]*e);
514   return sig;
515 }
516
517 //____________________________________________________________________________
518 Float_t  AliPHOSPIDv1::GetEllipseParameter(TString particle, TString param, Float_t e) const 
519 {
520   // Calculates the parameter param of the ellipse
521
522   particle.ToLower();
523   param.   ToLower();
524   Float_t p[4]={0.,0.,0.,0.};
525   Float_t value = 0.0;
526   for (Int_t i=0; i<4; i++) p[i] = GetParameterToCalculateEllipse(particle,param,i);
527   if (particle == "photon") {
528     if      (param.Contains("a"))  e = TMath::Min((Double_t)e,70.);
529     else if (param.Contains("b"))  e = TMath::Min((Double_t)e,70.);
530     else if (param.Contains("x0")) e = TMath::Max((Double_t)e,1.1);
531   }
532
533  if (particle == "photon")
534     value = p[0]/TMath::Sqrt(e) + p[1]*e + p[2]*e*e + p[3];
535   else if (particle == "pi0")
536     value = p[0] + p[1]*e + p[2]*e*e;
537
538   return value;
539 }
540
541 //_____________________________________________________________________________
542 Float_t  AliPHOSPIDv1::GetParameterPhotonBoundary (Int_t i) const
543
544   // Get the parameter "i" to calculate the boundary on the moment M2x
545   // for photons at high p_T
546   Float_t param = 0;
547   if (i>3 || i<0) {
548     AliError(Form("Wrong parameter number: %d\n",i));
549   } else
550     param = (*fParameters)(14,i) ;
551   return param;
552 }
553
554 //____________________________________________________________________________
555 Float_t  AliPHOSPIDv1::GetParameterPi0Boundary (Int_t i) const
556
557   // Get the parameter "i" to calculate the boundary on the moment M2x
558   // for pi0 at high p_T
559   Float_t param = 0;
560   if (i>2 || i<0) {
561     AliError(Form("Wrong parameter number: %d\n",i));
562   } else
563     param = (*fParameters)(15,i) ;
564   return param;
565 }
566
567 //____________________________________________________________________________
568 Float_t  AliPHOSPIDv1::GetParameterTimeGate(Int_t i) const
569 {
570   // Get TimeGate parameter depending on Purity-Efficiency i:
571   // i=0 - Low purity, i=1 - Medium purity, i=2 - High purity
572   Float_t param = 0.;
573   if(i>2 || i<0) {
574     AliError(Form("Invalid Efficiency-Purity choice %d",i));
575   } else
576     param = (*fParameters)(3,i) ; 
577   return param;
578 }
579
580 //_____________________________________________________________________________
581 Float_t  AliPHOSPIDv1::GetParameterToCalculateEllipse(TString particle, TString param, Int_t i) const
582
583   // Get the parameter "i" that is needed to calculate the ellipse 
584   // parameter "param" for the particle "particle" ("photon" or "pi0")
585
586   particle.ToLower();
587   param.   ToLower();
588   Int_t offset = -1;
589   if      (particle == "photon") 
590     offset=0;
591   else if (particle == "pi0")    
592     offset=5;
593   else
594     AliError(Form("Wrong particle name: %s (choose from pi0/photon)\n",
595                   particle.Data()));
596
597   Int_t p= -1;
598   Float_t par = 0;
599
600   if     (param.Contains("a")) p=4+offset; 
601   else if(param.Contains("b")) p=5+offset; 
602   else if(param.Contains("c")) p=6+offset; 
603   else if(param.Contains("x0"))p=7+offset; 
604   else if(param.Contains("y0"))p=8+offset;
605
606   if      (i>4 || i<0) {
607     AliError(Form("No parameter with index %d", i)) ; 
608   } else if (p==-1) {
609     AliError(Form("No parameter with name %s", param.Data() )) ; 
610   } else
611     par = (*fParameters)(p,i) ;
612   
613   return par;
614 }
615
616
617 //____________________________________________________________________________
618 Float_t  AliPHOSPIDv1::GetDistance(AliPHOSEmcRecPoint * emc,AliPHOSCpvRecPoint * cpv, Option_t *  axis)const
619 {
620   // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
621   
622   const AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ; 
623   TVector3 vecEmc ;
624   TVector3 vecCpv ;
625   if(cpv){
626     emc->GetLocalPosition(vecEmc) ;
627     cpv->GetLocalPosition(vecCpv) ; 
628     
629     if(emc->GetPHOSMod() == cpv->GetPHOSMod()){      
630       // Correct to difference in CPV and EMC position due to different distance to center.
631       // we assume, that particle moves from center
632       Float_t dCPV = geom->GetIPtoOuterCoverDistance();
633       Float_t dEMC = geom->GetIPtoCrystalSurface() ;
634       dEMC         = dEMC / dCPV ;
635       vecCpv = dEMC * vecCpv  - vecEmc ; 
636       if (axis == "X") return vecCpv.X();
637       if (axis == "Y") return vecCpv.Y();
638       if (axis == "Z") return vecCpv.Z();
639       if (axis == "R") return vecCpv.Mag();
640     }
641     return 100000000 ;
642   }
643   return 100000000 ;
644 }
645 //____________________________________________________________________________
646 Int_t  AliPHOSPIDv1::GetCPVBit(AliPHOSEmcRecPoint * emc,AliPHOSCpvRecPoint * cpv, Int_t effPur, Float_t e) const
647 {
648   //Calculates the pid bit for the CPV selection per each purity.
649   if(effPur>2 || effPur<0)
650     AliError(Form("Invalid Efficiency-Purity choice %d",effPur));
651   
652   Float_t sigX = GetCpv2EmcDistanceCut("X",e);
653   Float_t sigZ = GetCpv2EmcDistanceCut("Z",e);
654   
655   Float_t deltaX = TMath::Abs(GetDistance(emc, cpv,  "X"));
656   Float_t deltaZ = TMath::Abs(GetDistance(emc, cpv,  "Z"));
657   //Info("GetCPVBit"," xdist %f, sigx %f, zdist %f, sigz %f",deltaX, sigX, deltaZ,sigZ) ;
658  
659   //if(deltaX>sigX*(effPur+1))
660   //if((deltaX>sigX*(effPur+1)) || (deltaZ>sigZ*(effPur+1)))
661   if((deltaX>sigX*(effPur+1)) && (deltaZ>sigZ*(effPur+1)))
662     return 1;//Neutral
663   else
664     return 0;//Charged
665 }
666
667 //____________________________________________________________________________
668 Int_t  AliPHOSPIDv1::GetPrincipalBit(TString particle, const Double_t* p, Int_t effPur, Float_t e)const
669 {
670   //Is the particle inside de PCA ellipse?
671   
672   particle.ToLower();
673   Int_t    prinbit  = 0 ;
674   Float_t a  = GetEllipseParameter(particle,"a" , e); 
675   Float_t b  = GetEllipseParameter(particle,"b" , e);
676   Float_t c  = GetEllipseParameter(particle,"c" , e);
677   Float_t x0 = GetEllipseParameter(particle,"x0", e); 
678   Float_t y0 = GetEllipseParameter(particle,"y0", e);
679   
680   Float_t r = TMath::Power((p[0] - x0)/a,2) + 
681               TMath::Power((p[1] - y0)/b,2) +
682             c*(p[0] -  x0)*(p[1] - y0)/(a*b) ;
683   //3 different ellipses defined
684   if((effPur==2) && (r<1./2.)) prinbit= 1;
685   if((effPur==1) && (r<2.   )) prinbit= 1;
686   if((effPur==0) && (r<9./2.)) prinbit= 1;
687
688   if(r<0)
689     AliError("Negative square?") ;
690
691   return prinbit;
692
693 }
694 //____________________________________________________________________________
695 Int_t  AliPHOSPIDv1::GetHardPhotonBit(AliPHOSEmcRecPoint * emc) const
696 {
697   // Set bit for identified hard photons (E > 30 GeV)
698   // if the second moment M2x is below the boundary
699
700   Float_t e   = emc->GetEnergy();
701   if (e < 30.0) return 0;
702   Float_t m2x = emc->GetM2x();
703   Float_t m2xBoundary = GetParameterPhotonBoundary(0) *
704     TMath::Exp(-TMath::Power(e-GetParameterPhotonBoundary(1),2)/2.0/
705                 TMath::Power(GetParameterPhotonBoundary(2),2)) +
706     GetParameterPhotonBoundary(3);
707   AliDebug(1, Form("GetHardPhotonBit","E=%f, m2x=%f, boundary=%f",
708                        e,m2x,m2xBoundary));
709   if (m2x < m2xBoundary)
710     return 1;// A hard photon
711   else
712     return 0;// Not a hard photon
713 }
714
715 //____________________________________________________________________________
716 Int_t  AliPHOSPIDv1::GetHardPi0Bit(AliPHOSEmcRecPoint * emc) const
717 {
718   // Set bit for identified hard pi0  (E > 30 GeV)
719   // if the second moment M2x is above the boundary
720
721   Float_t e   = emc->GetEnergy();
722   if (e < 30.0) return 0;
723   Float_t m2x = emc->GetM2x();
724   Float_t m2xBoundary = GetParameterPi0Boundary(0) +
725                     e * GetParameterPi0Boundary(1);
726   AliDebug(1,Form("E=%f, m2x=%f, boundary=%f",e,m2x,m2xBoundary));
727   if (m2x > m2xBoundary)
728     return 1;// A hard pi0
729   else
730     return 0;// Not a hard pi0
731 }
732
733 //____________________________________________________________________________
734 TVector3 AliPHOSPIDv1::GetMomentumDirection(AliPHOSEmcRecPoint * emc, AliPHOSCpvRecPoint * )const 
735
736   // Calculates the momentum direction:
737   //   1. if only a EMC RecPoint, direction is given by IP and this RecPoint
738   //   2. if a EMC RecPoint and CPV RecPoint, direction is given by the line through the 2 recpoints 
739   //  However because of the poor position resolution of PPSD the direction is always taken as if we were 
740   //  in case 1.
741
742   TVector3 dir(0,0,0) ; 
743   TMatrix  dummy ;
744   
745   emc->GetGlobalPosition(dir, dummy) ;
746
747   //account correction to the position of IP
748   Float_t xo,yo,zo ; //Coordinates of the origin
749   if(gAlice && gAlice->GetMCApp() && gAlice->Generator()){
750     gAlice->Generator()->GetOrigin(xo,yo,zo) ;
751   }
752   else{
753     xo=yo=zo=0.;
754   }
755   TVector3 origin(xo,yo,zo);
756   dir = dir - origin ;
757   dir.SetMag(1.) ;
758
759   return dir ;  
760 }
761
762 //________________________________________________________________________
763 Double_t  AliPHOSPIDv1::LandauF(Double_t  x, Double_t y, Double_t * par)
764 {
765   //Given the energy x and the parameter y (tof, shower dispersion or cpv-emc distance), 
766   //this method returns a density probability of this parameter, given by a landau 
767   //function whose parameters depend with the energy  with a function: a/(x*x)+b/x+b
768
769   Double_t cnt    = par[1] / (x*x) + par[2] / x + par[0] ;
770   Double_t mean   = par[4] / (x*x) + par[5] / x + par[3] ;
771   Double_t sigma  = par[7] / (x*x) + par[8] / x + par[6] ;
772
773   if(TMath::Abs(sigma) > 1.e-10){
774     return cnt*TMath::Landau(y,mean,sigma);
775   }
776   else
777     return 0.;
778
779 }
780 //________________________________________________________________________
781 Double_t  AliPHOSPIDv1::LandauPol2(Double_t  x, Double_t y, Double_t * par)
782 {
783
784   //Given the energy x and the parameter y (tof, shower dispersion or cpv-emc distance), 
785   //this method returns a density probability of this parameter, given by a landau 
786   //function whose parameters depend with the energy like second order polinomial
787
788   Double_t cnt    = par[2] * (x*x) + par[1] * x + par[0] ;
789   Double_t mean   = par[5] * (x*x) + par[4] * x + par[3] ;
790   Double_t sigma  = par[8] * (x*x) + par[7] * x + par[6] ;
791
792    if(TMath::Abs(sigma) > 1.e-10){
793     return cnt*TMath::Landau(y,mean,sigma);
794   }
795   else
796     return 0.;
797
798
799 }
800 // //________________________________________________________________________
801 // Double_t  AliPHOSPIDv1::ChargedHadronDistProb(Double_t  x, Double_t y, Double_t * parg, Double_t * parl)
802 // {
803 //   Double_t cnt   = 0.0 ;
804 //   Double_t mean  = 0.0 ;
805 //   Double_t sigma = 0.0 ;
806 //   Double_t arg   = 0.0 ;
807 //   if (y < parl[4] / (x*x) + parl[5] / x + parl[3]){
808 //     cnt    = parg[1] / (x*x) + parg[2] / x + parg[0] ;
809 //     mean   = parg[4] / (x*x) + parg[5] / x + parg[3] ;
810 //     sigma  = parg[7] / (x*x) + parg[8] / x + parg[6] ;
811 //     TF1 * f = new TF1("gaus","gaus",0.,100.);
812 //     f->SetParameters(cnt,mean,sigma);
813 //     arg  = f->Eval(y) ;
814 //   }
815 //   else{
816 //     cnt    = parl[1] / (x*x) + parl[2] / x + parl[0] ;
817 //     mean   = parl[4] / (x*x) + parl[5] / x + parl[3] ;
818 //     sigma  = parl[7] / (x*x) + parl[8] / x + parl[6] ;
819 //     TF1 * f = new TF1("landau","landau",0.,100.);
820 //     f->SetParameters(cnt,mean,sigma);
821 //     arg  = f->Eval(y) ;
822 //   }
823 //   //  Double_t mean   = par[3] + par[4] * x + par[5] * x * x ;
824 //   //   Double_t sigma  = par[6] + par[7] * x + par[8] * x * x ;
825   
826 //   //Double_t arg    = -(y-mean)*(y-mean)/(2*sigma*sigma) ;
827 //   //return cnt * TMath::Exp(arg) ;
828   
829 //   return arg;
830   
831 // }
832 //____________________________________________________________________________
833 void  AliPHOSPIDv1::MakePID()
834 {
835   // construct the PID weight from a Bayesian Method
836   
837   const Int_t kSPECIES = AliPID::kSPECIESN ;
838  
839   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
840
841   Int_t nparticles = gime->RecParticles()->GetEntriesFast() ;
842
843   //   const Int_t kMAXPARTICLES = 2000 ; 
844   //   if (nparticles >= kMAXPARTICLES) 
845   //     Error("MakePID", "Change size of MAXPARTICLES") ; 
846   //   Double_t stof[kSPECIES][kMAXPARTICLES] ;
847   
848 //   const Int_t kMAXPARTICLES = 2000 ; 
849 //   if (nparticles >= kMAXPARTICLES) 
850 //     AliError("Change size of MAXPARTICLES") ; 
851 //   Double_t stof[kSPECIES][kMAXPARTICLES] ;
852
853
854   // make the normalized distribution of pid for this event 
855   // w(pid) in the Bayesian formulation
856 //   for(index = 0 ; index < nparticles ; index ++) {
857     
858 //     cout<<">>>>>>>>>>>>>>>Bayes Index "<<index<<endl;
859
860
861 //     AliPHOSEmcRecPoint * emc     = AliPHOSGetter::Instance()->EmcRecPoint(index) ;
862 //     AliPHOSCpvRecPoint * cpv     = AliPHOSGetter::Instance()->CpvRecPoint(index) ;
863
864   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
865   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
866   TClonesArray * trackSegments = gime->TrackSegments() ;
867   if ( !emcRecPoints || !cpvRecPoints || !trackSegments ) {
868     AliFatal("RecPoints or TrackSegments not found !") ;  
869   }
870   TIter next(trackSegments) ; 
871   AliPHOSTrackSegment * ts ; 
872   Int_t index = 0 ; 
873
874   Double_t * stof[kSPECIES] ;
875   Double_t * sdp [kSPECIES]  ;
876   Double_t * scpv[kSPECIES] ;
877   Double_t * sw  [kSPECIES] ;
878   //Info("MakePID","Begin MakePID"); 
879   
880   for (Int_t i =0; i< kSPECIES; i++){
881     stof[i] = new Double_t[nparticles] ;
882     sdp [i] = new Double_t[nparticles] ;
883     scpv[i] = new Double_t[nparticles] ;
884     sw  [i] = new Double_t[nparticles] ;
885   }
886   
887
888   while ( (ts = (AliPHOSTrackSegment *)next()) ) {
889     
890     //cout<<">>>>>> Bayesian Index "<<index<<endl;
891
892     AliPHOSEmcRecPoint * emc = 0 ;
893     if(ts->GetEmcIndex()>=0)
894       emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ;
895     
896     AliPHOSCpvRecPoint * cpv = 0 ;
897     if(ts->GetCpvIndex()>=0)
898       cpv = (AliPHOSCpvRecPoint *) cpvRecPoints->At(ts->GetCpvIndex()) ;
899     
900 //     Int_t track = 0 ; 
901 //     track = ts->GetTrackIndex() ; //TPC tracks ?
902     
903     if (!emc) {
904       AliFatal(Form("-> emc(%d) = %d", ts->GetEmcIndex(), emc )) ;
905     }
906
907     // ############Tof#############################
908
909     //    Info("MakePID", "TOF");
910     Float_t  en   = emc->GetEnergy();    
911     Double_t time = emc->GetTime() ;
912     //    cout<<">>>>>>>Energy "<<en<<"Time "<<time<<endl;
913    
914     // now get the signals probability
915     // s(pid) in the Bayesian formulation
916     
917     stof[AliPID::kPhoton][index]   = 1.; 
918     stof[AliPID::kElectron][index] = 1.;
919     stof[AliPID::kEleCon][index]   = 1.;
920     //We assing the same prob to charged hadrons, sum is 1
921     stof[AliPID::kPion][index]     = 1./3.; 
922     stof[AliPID::kKaon][index]     = 1./3.; 
923     stof[AliPID::kProton][index]   = 1./3.;
924     //We assing the same prob to neutral hadrons, sum is 1
925     stof[AliPID::kNeutron][index]  = 1./2.;
926     stof[AliPID::kKaon0][index]    = 1./2.;
927
928     stof[AliPID::kMuon][index]     = 1.; 
929
930     
931     if(en < 2.) {
932
933       Double_t pTofPion = fTFpiong ->Eval(time) ; //gaus distribution
934       Double_t pTofKaon = 0;
935
936       if(time < fTkaonl[1])
937         pTofKaon = fTFkaong  ->Eval(time) ; //gaus distribution
938       else 
939         pTofKaon = fTFkaonl  ->Eval(time) ; //landau distribution
940
941       Double_t pTofNucleon = 0;
942
943       if(time < fThhadronl[1])
944         pTofNucleon = fTFhhadrong   ->Eval(time) ; //gaus distribution
945       else
946         pTofNucleon = fTFhhadronl   ->Eval(time) ; //landau distribution
947       //We assing the same prob to neutral hadrons, sum is the average prob
948       Double_t pTofNeHadron =  (pTofKaon + pTofNucleon)/2. ;
949       //We assing the same prob to charged hadrons, sum is the average prob
950       Double_t pTofChHadron =  (pTofPion + pTofKaon + pTofNucleon)/3. ;
951
952       stof[AliPID::kPhoton][index]   = fTFphoton     ->Eval(time) ; 
953       //gaus distribution
954       stof[AliPID::kEleCon][index]   = stof[AliPID::kPhoton][index] ; 
955       //a conversion electron has the photon ToF
956       stof[AliPID::kMuon][index]     = stof[AliPID::kPhoton][index] ;
957  
958       stof[AliPID::kElectron][index] = pTofPion  ;                             
959
960       stof[AliPID::kPion][index]     =  pTofChHadron ; 
961       stof[AliPID::kKaon][index]     =  pTofChHadron ;
962       stof[AliPID::kProton][index]   =  pTofChHadron ;
963
964       stof[AliPID::kKaon0][index]    =  pTofNeHadron ;     
965       stof[AliPID::kNeutron][index]  =  pTofNeHadron ;            
966     } 
967     
968     //    Info("MakePID", "Dispersion");
969     
970     // ###########Shower shape: Dispersion####################
971     Float_t dispersion = emc->GetDispersion();
972     //dispersion is not well defined if the cluster is only in few crystals
973     
974     sdp[AliPID::kPhoton][index]   = 1. ;
975     sdp[AliPID::kElectron][index] = 1. ;
976     sdp[AliPID::kPion][index]     = 1. ; 
977     sdp[AliPID::kKaon][index]     = 1. ; 
978     sdp[AliPID::kProton][index]   = 1. ;
979     sdp[AliPID::kNeutron][index]  = 1. ;
980     sdp[AliPID::kEleCon][index]   = 1. ; 
981     sdp[AliPID::kKaon0][index]    = 1. ; 
982     sdp[AliPID::kMuon][index]     = 1. ; 
983     
984     if(en > 0.5 && emc->GetMultiplicity() > 3){ 
985       sdp[AliPID::kPhoton][index]   = GausF(en , dispersion, fDphoton) ;
986       sdp[AliPID::kElectron][index] = sdp[AliPID::kPhoton][index] ;
987       sdp[AliPID::kPion][index]     = LandauF(en , dispersion, fDhadron ) ; 
988       sdp[AliPID::kKaon][index]     = sdp[AliPID::kPion][index]  ; 
989       sdp[AliPID::kProton][index]   = sdp[AliPID::kPion][index]  ;
990       sdp[AliPID::kNeutron][index]  = sdp[AliPID::kPion][index]  ;
991       sdp[AliPID::kEleCon][index]   = sdp[AliPID::kPhoton][index]; 
992       sdp[AliPID::kKaon0][index]    = sdp[AliPID::kPion][index]  ; 
993       sdp[AliPID::kMuon][index]     = fDFmuon ->Eval(dispersion) ; 
994       //landau distribution
995     }
996     
997 //      Info("MakePID","multiplicity %d, dispersion %f", emc->GetMultiplicity(), dispersion);
998 //      Info("MakePID","ss: photon %f, hadron %f ",  sdp[AliPID::kPhoton][index],  sdp[AliPID::kPion][index]);
999 //       cout<<">>>>>multiplicity "<<emc->GetMultiplicity()<<", dispersion "<< dispersion<<endl ;
1000 //       cout<<"<<<<<ss: photon   "<<sdp[AliPID::kPhoton][index]<<", hadron    "<<sdp[AliPID::kPion][index]<<endl;
1001
1002     //########## CPV-EMC  Distance#######################
1003     //     Info("MakePID", "Distance");
1004
1005     Float_t x             = TMath::Abs(GetDistance(emc, cpv,  "X")) ;
1006     Float_t z             = GetDistance(emc, cpv,  "Z") ;
1007    
1008     Double_t pcpv         = 0 ;
1009     Double_t pcpvneutral  = 0. ;
1010    
1011     Double_t elprobx      = GausF(en , x, fXelectron) ;
1012     Double_t elprobz      = GausF(en , z, fZelectron) ;
1013     Double_t chprobx      = GausF(en , x, fXcharged)  ;
1014     Double_t chprobz      = GausF(en , z, fZcharged)  ;
1015     Double_t pcpvelectron = elprobx * elprobz;
1016     Double_t pcpvcharged  = chprobx * chprobz;
1017   
1018 //     cout<<">>>>energy "<<en<<endl;
1019 //     cout<<">>>>electron : x "<<x<<" xprob "<<elprobx<<" z "<<z<<" zprob "<<elprobz<<endl;
1020 //     cout<<">>>>hadron   : x "<<x<<" xprob "<<chprobx<<" z "<<z<<" zprob "<<chprobz<<endl;
1021 //     cout<<">>>>electron : px*pz "<<pcpvelectron <<" hadron: px*pz "<<pcpvcharged<<endl;  
1022
1023     // Is neutral or charged
1024     if(pcpvelectron >= pcpvcharged)  
1025       pcpv = pcpvelectron ;
1026     else
1027       pcpv = pcpvcharged ;
1028     
1029     if(pcpv < fChargedNeutralThreshold)
1030       {
1031         pcpvneutral  = 1. ;
1032         pcpvcharged  = 0. ;
1033         pcpvelectron = 0. ;
1034       }
1035     //    else
1036     //      cout<<">>>>>>>>>>>CHARGED>>>>>>>>>>>"<<endl;
1037     
1038     scpv[AliPID::kPion][index]     =  pcpvcharged  ; 
1039     scpv[AliPID::kKaon][index]     =  pcpvcharged  ; 
1040     scpv[AliPID::kProton][index]   =  pcpvcharged  ;
1041
1042     scpv[AliPID::kMuon][index]     =  pcpvelectron ; 
1043     scpv[AliPID::kElectron][index] =  pcpvelectron ;
1044     scpv[AliPID::kEleCon][index]   =  pcpvelectron ; 
1045
1046     scpv[AliPID::kPhoton][index]   =  pcpvneutral  ;
1047     scpv[AliPID::kNeutron][index]  =  pcpvneutral  ; 
1048     scpv[AliPID::kKaon0][index]    =  pcpvneutral  ; 
1049
1050     
1051     //   Info("MakePID", "CPV passed");
1052
1053     //############## Pi0 #############################
1054     stof[AliPID::kPi0][index]      = 0. ;  
1055     scpv[AliPID::kPi0][index]      = 0. ;
1056     sdp [AliPID::kPi0][index]      = 0. ;
1057
1058     if(en > 30.){
1059       // pi0 are detected via decay photon
1060       stof[AliPID::kPi0][index]  = fTFphoton  ->Eval(time) ;
1061       scpv[AliPID::kPi0][index]  = pcpvneutral  ;
1062       sdp [AliPID::kPi0][index]  = 1. ;
1063       if(emc->GetMultiplicity() > 3)
1064         sdp [AliPID::kPi0][index]  = GausPol2(en , dispersion, fDpi0) ;
1065     }
1066     
1067
1068     //############## muon #############################
1069
1070     if(en > 0.5){
1071       //Muons deposit few energy
1072       scpv[AliPID::kMuon][index]     =  0 ;
1073       stof[AliPID::kMuon][index]     =  0 ;
1074       sdp [AliPID::kMuon][index]     =  0 ;
1075     }
1076
1077     //Weight to apply to hadrons due to energy reconstruction
1078
1079     Float_t weight = fERecWeight ->Eval(en) ;
1080  
1081     sw[AliPID::kPhoton][index]   = 1. ;
1082     sw[AliPID::kElectron][index] = 1. ;
1083     sw[AliPID::kPion][index]     = weight ; 
1084     sw[AliPID::kKaon][index]     = weight ; 
1085     sw[AliPID::kProton][index]   = weight ;
1086     sw[AliPID::kNeutron][index]  = weight ;
1087     sw[AliPID::kEleCon][index]   = 1. ; 
1088     sw[AliPID::kKaon0][index]    = weight ; 
1089     sw[AliPID::kMuon][index]     = weight ; 
1090     sw[AliPID::kPi0][index]      = 1. ;
1091
1092 //     if(en > 0.5){
1093 //       cout<<"######################################################"<<endl;
1094 //       //cout<<"MakePID: energy "<<en<<", tof "<<time<<", distance "<<distance<<", dispersion "<<dispersion<<endl ;
1095 //       cout<<"MakePID: energy "<<en<<", tof "<<time<<", dispersion "<<dispersion<<", x "<<x<<", z "<<z<<endl ;
1096 //       cout<<">>>>>multiplicity "<<emc->GetMultiplicity()<<endl;
1097 //       cout<<">>>>electron : xprob "<<elprobx<<" zprob "<<elprobz<<endl;
1098 //       cout<<">>>>hadron   : xprob "<<chprobx<<" zprob "<<chprobz<<endl;
1099 //       cout<<">>>>electron : px*pz "<<pcpvelectron <<" hadron: px*pz "<<pcpvcharged<<endl;  
1100       
1101 //       cout<<"Photon   , pid "<< fInitPID[AliPID::kPhoton]<<" tof "<<stof[AliPID::kPhoton][index]
1102 //        <<", cpv "<<scpv[AliPID::kPhoton][index]<<", ss "<<sdp[AliPID::kPhoton][index]<<endl;
1103 //       cout<<"EleCon   , pid "<< fInitPID[AliPID::kEleCon]<<", tof "<<stof[AliPID::kEleCon][index]
1104 //        <<", cpv "<<scpv[AliPID::kEleCon][index]<<" ss "<<sdp[AliPID::kEleCon][index]<<endl;
1105 //       cout<<"Electron , pid "<< fInitPID[AliPID::kElectron]<<", tof "<<stof[AliPID::kElectron][index]
1106 //        <<", cpv "<<scpv[AliPID::kElectron][index]<<" ss "<<sdp[AliPID::kElectron][index]<<endl;
1107 //       cout<<"Muon     , pid "<< fInitPID[AliPID::kMuon]<<", tof "<<stof[AliPID::kMuon][index]
1108 //        <<", cpv "<<scpv[AliPID::kMuon][index]<<" ss "<<sdp[AliPID::kMuon][index]<<endl;
1109 //       cout<<"Pi0      , pid "<< fInitPID[AliPID::kPi0]<<", tof "<<stof[AliPID::kPi0][index]
1110 //        <<", cpv "<<scpv[AliPID::kPi0][index]<<" ss "<<sdp[AliPID::kPi0][index]<<endl;
1111 //       cout<<"Pion     , pid "<< fInitPID[AliPID::kPion]<<", tof "<<stof[AliPID::kPion][index]
1112 //        <<", cpv "<<scpv[AliPID::kPion][index]<<" ss "<<sdp[AliPID::kPion][index]<<endl;
1113 //       cout<<"Kaon0    , pid "<< fInitPID[AliPID::kKaon0]<<", tof "<<stof[AliPID::kKaon0][index]
1114 //        <<", cpv "<<scpv[AliPID::kKaon0][index]<<" ss "<<sdp[AliPID::kKaon0][index]<<endl;
1115 //       cout<<"Kaon     , pid "<< fInitPID[AliPID::kKaon]<<", tof "<<stof[AliPID::kKaon][index]
1116 //        <<", cpv "<<scpv[AliPID::kKaon][index]<<" ss "<<sdp[AliPID::kKaon][index]<<endl;
1117 //       cout<<"Neutron  , pid "<< fInitPID[AliPID::kNeutron]<<", tof "<<stof[AliPID::kNeutron][index]
1118 //        <<", cpv "<<scpv[AliPID::kNeutron][index]<<" ss "<<sdp[AliPID::kNeutron][index]<<endl;
1119 //       cout<<"Proton   , pid "<< fInitPID[AliPID::kProton]<<", tof "<<stof[AliPID::kProton][index]
1120 //        <<", cpv "<<scpv[AliPID::kProton][index]<<" ss "<<sdp[AliPID::kProton][index]<<endl;
1121 //       cout<<"######################################################"<<endl;
1122 //     }
1123       index++;
1124   }
1125   
1126   //for (index = 0 ; index < kSPECIES ; index++) 
1127   // pid[index] /= nparticles ; 
1128   
1129
1130   //  Info("MakePID", "Total Probability calculation");
1131   
1132   for(index = 0 ; index < nparticles ; index ++) {
1133     // calculates the Bayesian weight
1134     
1135     Int_t jndex ;
1136     Double_t wn = 0.0 ; 
1137     for (jndex = 0 ; jndex < kSPECIES ; jndex++) 
1138       wn += stof[jndex][index] * sdp[jndex][index]  * scpv[jndex][index] * sw[jndex][index] * fInitPID[jndex] ;
1139    
1140     //    cout<<"*************wn "<<wn<<endl;
1141     AliPHOSRecParticle * recpar = gime->RecParticle(index) ;  
1142     if (TMath::Abs(wn)>0)
1143       for (jndex = 0 ; jndex < kSPECIES ; jndex++) {
1144
1145 //      if(recpar->IsEleCon()){
1146 //       fInitPID[AliPID::kEleCon]   = 1. ;
1147 //       fInitPID[AliPID::kPhoton]   = 0. ;
1148 //      }
1149 //      else{
1150 //       fInitPID[AliPID::kEleCon]   = 0. ;
1151 //       fInitPID[AliPID::kPhoton]   = 1. ;
1152 //      }
1153         fInitPID[AliPID::kEleCon]   = 0. ;
1154         //cout<<"jndex "<<jndex<<" wn "<<wn<<" SetPID * wn"
1155         //<<stof[jndex][index] * sdp[jndex][index] * pid[jndex]  << endl;
1156         //cout<<" tof "<<stof[jndex][index] << " disp " <<sdp[jndex][index] << " pid "<< fInitPID[jndex] << endl;
1157 //      cout<<"Particle "<<jndex<<"  final prob * wn   "
1158 //          <<stof[jndex][index] * sdp[jndex][index] * scpv[jndex][index] * fInitPID[jndex] <<"  wn  "<< wn<<endl;
1159         recpar->SetPID(jndex, stof[jndex][index] * sdp[jndex][index] * sw[jndex][index] *
1160                        scpv[jndex][index] * fInitPID[jndex] / wn) ; 
1161 //      cout<<"final prob "<<stof[jndex][index] * sdp[jndex][index] * scpv[jndex][index] * fInitPID[jndex] / wn<<endl;
1162         //recpar->SetPID(jndex, stof[jndex][index] * fInitPID[jndex] / wn) ; 
1163         //cout<<"After SetPID"<<endl;
1164         //recpar->Print();
1165       }
1166   }
1167   //  Info("MakePID", "Delete");
1168   
1169     for (Int_t i =0; i< kSPECIES; i++){
1170       delete [] stof[i];
1171       delete [] sdp [i];
1172       delete [] scpv[i];
1173       delete [] sw  [i];
1174     }
1175   //  Info("MakePID","End MakePID"); 
1176 }
1177
1178 //____________________________________________________________________________
1179 void  AliPHOSPIDv1::MakeRecParticles()
1180 {
1181   // Makes a RecParticle out of a TrackSegment
1182   
1183   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
1184   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
1185   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
1186   TClonesArray * trackSegments = gime->TrackSegments() ; 
1187   if ( !emcRecPoints || !cpvRecPoints || !trackSegments ) {
1188     AliFatal("RecPoints or TrackSegments not found !") ;  
1189   }
1190   TClonesArray * recParticles  = gime->RecParticles() ; 
1191   recParticles->Clear();
1192
1193   TIter next(trackSegments) ; 
1194   AliPHOSTrackSegment * ts ; 
1195   Int_t index = 0 ; 
1196   AliPHOSRecParticle * rp ; 
1197   while ( (ts = (AliPHOSTrackSegment *)next()) ) {
1198     //  cout<<">>>>>>>>>>>>>>>PCA Index "<<index<<endl;
1199     new( (*recParticles)[index] ) AliPHOSRecParticle() ;
1200     rp = (AliPHOSRecParticle *)recParticles->At(index) ; 
1201     rp->SetTrackSegment(index) ;
1202     rp->SetIndexInList(index) ;
1203         
1204     AliPHOSEmcRecPoint * emc = 0 ;
1205     if(ts->GetEmcIndex()>=0)
1206       emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ;
1207     
1208     AliPHOSCpvRecPoint * cpv = 0 ;
1209     if(ts->GetCpvIndex()>=0)
1210       cpv = (AliPHOSCpvRecPoint *) cpvRecPoints->At(ts->GetCpvIndex()) ;
1211     
1212     Int_t track = 0 ; 
1213     track = ts->GetTrackIndex() ; 
1214       
1215     // Now set type (reconstructed) of the particle
1216
1217     // Choose the cluster energy range
1218     
1219     if (!emc) {
1220       AliFatal(Form("-> emc(%d) = %d", ts->GetEmcIndex(), emc )) ;
1221     }
1222
1223     Float_t e = emc->GetEnergy() ;   
1224     
1225     Float_t  lambda[2] ;
1226     emc->GetElipsAxis(lambda) ;
1227     
1228     if((lambda[0]>0.01) && (lambda[1]>0.01)){
1229       // Looking PCA. Define and calculate the data (X),
1230       // introduce in the function X2P that gives the components (P).  
1231
1232       Float_t  spher = 0. ;
1233       Float_t  emaxdtotal = 0. ; 
1234       
1235       if((lambda[0]+lambda[1])!=0) 
1236         spher=fabs(lambda[0]-lambda[1])/(lambda[0]+lambda[1]); 
1237       
1238       emaxdtotal=emc->GetMaximalEnergy()/emc->GetEnergy(); 
1239       
1240       fX[0] = lambda[0] ;  
1241       fX[1] = lambda[1] ; 
1242       fX[2] = emc->GetDispersion() ; 
1243       fX[3] = spher ; 
1244       fX[4] = emc->GetMultiplicity() ;  
1245       fX[5] = emaxdtotal ;  
1246       fX[6] = emc->GetCoreEnergy() ;  
1247       
1248       fPrincipalPhoton->X2P(fX,fPPhoton);
1249       fPrincipalPi0   ->X2P(fX,fPPi0);
1250
1251     }
1252     else{
1253       fPPhoton[0]=-100.0;  //We do not accept clusters with 
1254       fPPhoton[1]=-100.0;  //one cell as a photon-like
1255       fPPi0[0]   =-100.0;
1256       fPPi0[1]   =-100.0;
1257     }
1258     
1259     Float_t time = emc->GetTime() ;
1260     rp->SetTof(time) ; 
1261     
1262     // Loop of Efficiency-Purity (the 3 points of purity or efficiency 
1263     // are taken into account to set the particle identification)
1264     for(Int_t effPur = 0; effPur < 3 ; effPur++){
1265       
1266       // Looking at the CPV detector. If RCPV greater than CpvEmcDistance, 
1267       // 1st,2nd or 3rd bit (depending on the efficiency-purity point )
1268       // is set to 1
1269       if(GetCPVBit(emc, cpv, effPur,e) == 1 ){  
1270         rp->SetPIDBit(effPur) ;
1271         //cout<<"CPV bit "<<effPur<<endl;
1272       }
1273       // Looking the TOF. If TOF smaller than gate,  4th, 5th or 6th 
1274       // bit (depending on the efficiency-purity point )is set to 1             
1275       if(time< (*fParameters)(3,effPur)) 
1276         rp->SetPIDBit(effPur+3) ;                   
1277   
1278       //Photon PCA
1279       //If we are inside the ellipse, 7th, 8th or 9th 
1280       // bit (depending on the efficiency-purity point )is set to 1 
1281       if(GetPrincipalBit("photon",fPPhoton,effPur,e) == 1) 
1282         rp->SetPIDBit(effPur+6) ;
1283
1284       //Pi0 PCA
1285       //If we are inside the ellipse, 10th, 11th or 12th 
1286       // bit (depending on the efficiency-purity point )is set to 1 
1287       if(GetPrincipalBit("pi0"   ,fPPi0   ,effPur,e) == 1) 
1288         rp->SetPIDBit(effPur+9) ;
1289     }
1290     if(GetHardPhotonBit(emc))
1291       rp->SetPIDBit(12) ;
1292     if(GetHardPi0Bit   (emc))
1293       rp->SetPIDBit(13) ;
1294     
1295     if(track >= 0) 
1296       rp->SetPIDBit(14) ; 
1297
1298     //Set momentum, energy and other parameters 
1299     Float_t  encal = GetCalibratedEnergy(e);
1300     TVector3 dir   = GetMomentumDirection(emc,cpv) ; 
1301     dir.SetMag(encal) ;
1302     rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),encal) ;
1303     rp->SetCalcMass(0);
1304     rp->Name(); //If photon sets the particle pdg name to gamma
1305     rp->SetProductionVertex(0,0,0,0);
1306     rp->SetFirstMother(-1);
1307     rp->SetLastMother(-1);
1308     rp->SetFirstDaughter(-1);
1309     rp->SetLastDaughter(-1);
1310     rp->SetPolarisation(0,0,0);
1311     //Set the position in global coordinate system from the RecPoint
1312     AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
1313     AliPHOSTrackSegment * ts  = gime->TrackSegment(rp->GetPHOSTSIndex()) ; 
1314     AliPHOSEmcRecPoint  * erp = gime->EmcRecPoint(ts->GetEmcIndex()) ; 
1315     TVector3 pos ; 
1316     geom->GetGlobal(erp, pos) ; 
1317     rp->SetPos(pos);
1318     index++ ; 
1319   }
1320 }
1321   
1322 //____________________________________________________________________________
1323 void  AliPHOSPIDv1::Print() const
1324 {
1325   // Print the parameters used for the particle type identification
1326
1327     AliInfo("=============== AliPHOSPIDv1 ================") ;
1328     printf("Making PID\n") ;
1329     printf("    Pricipal analysis file from 0.5 to 100 %s\n", fFileNamePrincipalPhoton.Data() )   ; 
1330     printf("    Name of parameters file     %s\n", fFileNameParameters.Data() )  ;
1331     printf("    Matrix of Parameters: 14x4\n") ;
1332     printf("        Energy Calibration  1x3 [3 parametres to calibrate energy: A + B* E + C * E^2]\n") ;
1333     printf("        RCPV 2x3 rows x and z, columns function cut parameters\n") ;
1334     printf("        TOF  1x3 [High Eff-Low Pur,Medium Eff-Pur, Low Eff-High Pur]\n") ;
1335     printf("        PCA  5x4 [5 ellipse parametres and 4 parametres to calculate them: A/Sqrt(E) + B* E + C * E^2 + D]\n") ;
1336     Printf("    Pi0 PCA  5x3 [5 ellipse parametres and 3 parametres to calculate them: A + B* E + C * E^2]\n") ;
1337     fParameters->Print() ;
1338 }
1339
1340
1341
1342 //____________________________________________________________________________
1343 void AliPHOSPIDv1::PrintRecParticles(Option_t * option)
1344 {
1345   // Print table of reconstructed particles
1346
1347   AliPHOSGetter *gime = AliPHOSGetter::Instance() ; 
1348
1349   TClonesArray * recParticles = gime->RecParticles() ; 
1350
1351   TString message ; 
1352   message  = "\nevent " ;
1353   message += gAlice->GetEvNumber() ; 
1354   message += "       found " ; 
1355   message += recParticles->GetEntriesFast(); 
1356   message += " RecParticles\n" ; 
1357
1358   if(strstr(option,"all")) {  // printing found TS
1359     message += "\n  PARTICLE         Index    \n" ; 
1360     
1361     Int_t index ;
1362     for (index = 0 ; index < recParticles->GetEntries() ; index++) {
1363       AliPHOSRecParticle * rp = (AliPHOSRecParticle * ) recParticles->At(index) ;       
1364       message += "\n" ;
1365       message += rp->Name().Data() ;  
1366       message += " " ;
1367       message += rp->GetIndexInList() ;  
1368       message += " " ;
1369       message += rp->GetType()  ;
1370     }
1371   }
1372   AliInfo(message.Data() ) ; 
1373 }
1374
1375 //____________________________________________________________________________
1376 void  AliPHOSPIDv1::SetParameters() 
1377 {
1378   // PCA : To do the Principal Components Analysis it is necessary 
1379   // the Principal file, which is opened here
1380   fX       = new double[7]; // Data for the PCA 
1381   fPPhoton = new double[7]; // Eigenvalues of the PCA
1382   fPPi0    = new double[7]; // Eigenvalues of the Pi0 PCA
1383
1384   // Read photon principals from the photon file
1385   
1386   fFileNamePrincipalPhoton = "$ALICE_ROOT/PHOS/PCA8pa15_0.5-100.root" ; 
1387   TFile f( fFileNamePrincipalPhoton.Data(), "read" ) ;
1388   fPrincipalPhoton = dynamic_cast<TPrincipal*> (f.Get("principal")) ; 
1389   f.Close() ; 
1390
1391   // Read pi0 principals from the pi0 file
1392
1393   fFileNamePrincipalPi0    = "$ALICE_ROOT/PHOS/PCA_pi0_40-120.root" ;
1394   TFile fPi0( fFileNamePrincipalPi0.Data(), "read" ) ;
1395   fPrincipalPi0    = dynamic_cast<TPrincipal*> (fPi0.Get("principal")) ; 
1396   fPi0.Close() ;
1397
1398   // Open parameters file and initialization of the Parameters matrix. 
1399   // In the File Parameters.dat are all the parameters. These are introduced 
1400   // in a matrix of 16x4  
1401   // 
1402   // All the parameters defined in this file are, in order of row: 
1403   // line   0   : calibration 
1404   // lines  1,2 : CPV rectangular cat for X and Z
1405   // line   3   : TOF cut
1406   // lines  4-8 : parameters to calculate photon PCA ellipse
1407   // lines  9-13: parameters to calculate pi0 PCA ellipse
1408   // lines 14-15: parameters to calculate border for high-pt photons and pi0
1409
1410   fFileNameParameters = gSystem->ExpandPathName("$ALICE_ROOT/PHOS/Parameters.dat");
1411   fParameters = new TMatrix(16,4) ;
1412   const Int_t kMaxLeng=255;
1413   char string[kMaxLeng];
1414
1415   // Open a text file with PID parameters
1416   FILE *fd = fopen(fFileNameParameters.Data(),"r");
1417   if (!fd)
1418     AliFatal(Form("File %s with a PID parameters cannot be opened\n",
1419           fFileNameParameters.Data()));
1420
1421   Int_t i=0;
1422   // Read parameter file line-by-line and skip empty line and comments
1423   while (fgets(string,kMaxLeng,fd) != NULL) {
1424     if (string[0] == '\n' ) continue;
1425     if (string[0] == '!'  ) continue;
1426     sscanf(string, "%f %f %f %f",
1427            &(*fParameters)(i,0), &(*fParameters)(i,1), 
1428            &(*fParameters)(i,2), &(*fParameters)(i,3));
1429     i++;
1430     AliDebug(1, Form("SetParameters", "line %d: %s",i,string));
1431   }
1432   fclose(fd);
1433 }
1434
1435 //____________________________________________________________________________
1436 void  AliPHOSPIDv1::SetParameterCalibration(Int_t i,Float_t param) 
1437 {
1438   // Set parameter "Calibration" i to a value param
1439   if(i>2 || i<0) {
1440     AliError(Form("Invalid parameter number: %d",i));
1441   } else
1442     (*fParameters)(0,i) = param ;
1443 }
1444
1445 //____________________________________________________________________________
1446 void  AliPHOSPIDv1::SetParameterCpv2Emc(Int_t i, TString axis, Float_t cut) 
1447 {
1448   // Set the parameters to calculate Cpv-to-Emc Distance Cut depending on 
1449   // Purity-Efficiency point i
1450
1451   if(i>2 || i<0) {
1452     AliError(Form("Invalid parameter number: %d",i));
1453   } else {
1454     axis.ToLower();
1455     if      (axis == "x") (*fParameters)(1,i) = cut;
1456     else if (axis == "z") (*fParameters)(2,i) = cut;
1457     else { 
1458       AliError(Form("Invalid axis name: %s",axis.Data()));
1459     }
1460   }
1461 }
1462
1463 //____________________________________________________________________________
1464 void  AliPHOSPIDv1::SetParameterPhotonBoundary(Int_t i,Float_t param) 
1465 {
1466   // Set parameter "Hard photon boundary" i to a value param
1467   if(i>4 || i<0) {
1468     AliError(Form("Invalid parameter number: %d",i));
1469   } else
1470     (*fParameters)(14,i) = param ;
1471 }
1472
1473 //____________________________________________________________________________
1474 void  AliPHOSPIDv1::SetParameterPi0Boundary(Int_t i,Float_t param) 
1475 {
1476   // Set parameter "Hard pi0 boundary" i to a value param
1477   if(i>1 || i<0) {
1478     AliError(Form("Invalid parameter number: %d",i));
1479   } else
1480     (*fParameters)(15,i) = param ;
1481 }
1482
1483 //_____________________________________________________________________________
1484 void  AliPHOSPIDv1::SetParameterTimeGate(Int_t i, Float_t gate) 
1485 {
1486   // Set the parameter TimeGate depending on Purity-Efficiency point i 
1487   if (i>2 || i<0) {
1488     AliError(Form("Invalid Efficiency-Purity choice %d",i));
1489   } else
1490     (*fParameters)(3,i)= gate ; 
1491
1492
1493 //_____________________________________________________________________________
1494 void  AliPHOSPIDv1::SetParameterToCalculateEllipse(TString particle, TString param, Int_t i, Float_t par) 
1495 {  
1496   // Set the parameter "i" that is needed to calculate the ellipse 
1497   // parameter "param" for a particle "particle"
1498   
1499   particle.ToLower();
1500   param.   ToLower();
1501   Int_t p= -1;
1502   Int_t offset=0;
1503
1504   if      (particle == "photon") offset=0;
1505   else if (particle == "pi0")    offset=5;
1506   else
1507     AliError(Form("Wrong particle name: %s (choose from pi0/photon)\n",
1508                   particle.Data()));
1509
1510   if     (param.Contains("a")) p=4+offset; 
1511   else if(param.Contains("b")) p=5+offset; 
1512   else if(param.Contains("c")) p=6+offset; 
1513   else if(param.Contains("x0"))p=7+offset; 
1514   else if(param.Contains("y0"))p=8+offset;
1515   if((i>4)||(i<0)) {
1516     AliError(Form("No parameter with index %d", i)) ; 
1517   } else if(p==-1) {
1518     AliError(Form("No parameter with name %s", param.Data() )) ; 
1519   } else
1520     (*fParameters)(p,i) = par ;
1521
1522
1523 //____________________________________________________________________________
1524 void AliPHOSPIDv1::Unload() 
1525 {
1526   //Unloads RecPoints, Tracks and RecParticles
1527   AliPHOSGetter * gime = AliPHOSGetter::Instance() ;  
1528   gime->PhosLoader()->UnloadRecPoints() ;
1529   gime->PhosLoader()->UnloadTracks() ;
1530   gime->PhosLoader()->UnloadRecParticles() ;
1531 }
1532
1533 //____________________________________________________________________________
1534 void  AliPHOSPIDv1::WriteRecParticles()
1535 {
1536   //It writes reconstructed particles and pid to file
1537
1538   AliPHOSGetter *gime = AliPHOSGetter::Instance() ; 
1539
1540   TClonesArray * recParticles = gime->RecParticles() ; 
1541   recParticles->Expand(recParticles->GetEntriesFast() ) ;
1542   if(fWrite){
1543     TTree * treeP =  gime->TreeP();
1544     
1545     //First rp
1546     Int_t bufferSize = 32000 ;
1547     TBranch * rpBranch = treeP->Branch("PHOSRP",&recParticles,bufferSize);
1548     rpBranch->SetTitle(BranchName());
1549     
1550     rpBranch->Fill() ;
1551     
1552     gime->WriteRecParticles("OVERWRITE");
1553     gime->WritePID("OVERWRITE");
1554   }
1555 }
1556
1557
1558 //_______________________________________________________________________
1559 void AliPHOSPIDv1::SetInitPID(const Double_t *p) {
1560   // Sets values for the initial population of each particle type 
1561   for (Int_t i=0; i<AliPID::kSPECIESN; i++) fInitPID[i] = p[i];
1562 }
1563 //_______________________________________________________________________
1564 void AliPHOSPIDv1::GetInitPID(Double_t *p) const {
1565   // Gets values for the initial population of each particle type 
1566   for (Int_t i=0; i<AliPID::kSPECIESN; i++) p[i] = fInitPID[i];
1567 }