Correct: wrong matching between CPV and EMC, initial pi0 population was set to 0
[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   for (Int_t i =0; i<  AliPID::kSPECIESN ; i++)
328     fInitPID[i] = 1.;
329   
330 }
331
332 //________________________________________________________________________
333 void  AliPHOSPIDv1::Exec(Option_t *option)
334 {
335   // Steering method to perform particle reconstruction and identification
336   // for the event range from fFirstEvent to fLastEvent.
337   // This range is optionally set by SetEventRange().
338   // if fLastEvent=-1 (by default), then process events until the end.
339   
340   if(strstr(option,"tim"))
341     gBenchmark->Start("PHOSPID");
342   
343   if(strstr(option,"print")) {
344     Print() ; 
345     return ; 
346   }
347
348
349   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
350  
351   if (fLastEvent == -1) 
352     fLastEvent = gime->MaxEvent() - 1 ;
353   else 
354     fLastEvent = TMath::Min(fLastEvent,gime->MaxEvent());
355   Int_t nEvents   = fLastEvent - fFirstEvent + 1;
356
357   Int_t ievent ; 
358   for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
359     gime->Event(ievent,"TR") ;
360     if(gime->TrackSegments() && //Skip events, where no track segments made
361        gime->TrackSegments()->GetEntriesFast()) {
362
363       MakeRecParticles() ;
364
365       if(fBayesian)
366         MakePID() ; 
367       
368       WriteRecParticles();
369       if(strstr(option,"deb"))
370         PrintRecParticles(option) ;
371       //increment the total number of rec particles per run 
372       fRecParticlesInRun += gime->RecParticles()->GetEntriesFast() ; 
373     }
374   }
375   if(strstr(option,"deb"))
376       PrintRecParticles(option);
377   if(strstr(option,"tim")){
378     gBenchmark->Stop("PHOSPID");
379     AliInfo(Form("took %f seconds for PID %f seconds per event", 
380          gBenchmark->GetCpuTime("PHOSPID"),  
381          gBenchmark->GetCpuTime("PHOSPID")/nEvents)) ;
382   }
383   if(fWrite)
384     Unload();
385 }
386
387 //________________________________________________________________________
388 Double_t  AliPHOSPIDv1::GausF(Double_t  x, Double_t  y, Double_t * par)
389 {
390   //Given the energy x and the parameter y (tof, shower dispersion or cpv-emc distance), 
391   //this method returns a density probability of this parameter, given by a gaussian 
392   //function whose parameters depend with the energy  with a function: a/(x*x)+b/x+b
393   Double_t cnt    = par[1] / (x*x) + par[2] / x + par[0] ;
394   Double_t mean   = par[4] / (x*x) + par[5] / x + par[3] ;
395   Double_t sigma  = par[7] / (x*x) + par[8] / x + par[6] ;
396  
397   //  Double_t arg    = - (y-mean) * (y-mean) / (2*sigma*sigma) ;
398   //  return cnt * TMath::Exp(arg) ;
399   if(TMath::Abs(sigma) > 1.e-10){
400     return cnt*TMath::Gaus(y,mean,sigma);
401   }
402   else
403     return 0.;
404  
405 }
406 //________________________________________________________________________
407 Double_t  AliPHOSPIDv1::GausPol2(Double_t  x, Double_t y, Double_t * par)
408 {
409   //Given the energy x and the parameter y (tof, shower dispersion or cpv-emc distance), 
410   //this method returns a density probability of this parameter, given by a gaussian 
411   //function whose parameters depend with the energy like second order polinomial
412
413   Double_t cnt    = par[0] + par[1] * x + par[2] * x * x ;
414   Double_t mean   = par[3] + par[4] * x + par[5] * x * x ;
415   Double_t sigma  = par[6] + par[7] * x + par[8] * x * x ;
416
417   if(TMath::Abs(sigma) > 1.e-10){
418     return cnt*TMath::Gaus(y,mean,sigma);
419   }
420   else
421     return 0.;
422  
423
424
425 }
426
427 //____________________________________________________________________________
428 const TString AliPHOSPIDv1::GetFileNamePrincipal(TString particle) const
429 {
430   //Get file name that contains the PCA for a particle ("photon or pi0")
431   particle.ToLower();
432   TString name;
433   if      (particle=="photon") 
434     name = fFileNamePrincipalPhoton ;
435   else if (particle=="pi0"   ) 
436     name = fFileNamePrincipalPi0    ;
437   else    
438     AliError(Form("Wrong particle name: %s (choose from pi0/photon)\n",
439                   particle.Data()));
440   return name;
441 }
442
443 //____________________________________________________________________________
444 Float_t  AliPHOSPIDv1::GetParameterCalibration(Int_t i) const 
445 {
446   // Get the i-th parameter "Calibration"
447   Float_t param = 0.;
448   if (i>2 || i<0) { 
449     AliError(Form("Invalid parameter number: %d",i));
450   } else
451     param = (*fParameters)(0,i);
452   return param;
453 }
454
455 //____________________________________________________________________________
456 Float_t  AliPHOSPIDv1::GetCalibratedEnergy(Float_t e) const
457 {
458 //      It calibrates Energy depending on the recpoint energy.
459 //      The energy of the reconstructed cluster is corrected with 
460 //      the formula A + B* E  + C* E^2, whose parameters where obtained 
461 //      through the study of the reconstructed energy distribution of 
462 //      monoenergetic photons.
463  
464   Float_t p[]={0.,0.,0.};
465   for (Int_t i=0; i<3; i++) p[i] = GetParameterCalibration(i);
466   Float_t enerec = p[0] +  p[1]*e + p[2]*e*e;
467   return enerec ;
468
469 }
470
471 //____________________________________________________________________________
472 Float_t  AliPHOSPIDv1::GetParameterCpv2Emc(Int_t i, TString axis) const 
473 {
474   // Get the i-th parameter "CPV-EMC distance" for the specified axis
475   Float_t param = 0.;
476   if(i>2 || i<0) {
477     AliError(Form("Invalid parameter number: %d",i));
478   } else {
479     axis.ToLower();
480     if      (axis == "x") 
481       param = (*fParameters)(1,i);
482     else if (axis == "z") 
483       param = (*fParameters)(2,i);
484     else { 
485       AliError(Form("Invalid axis name: %s",axis.Data()));
486     }
487   }
488   return  param;
489 }
490
491 //____________________________________________________________________________
492 Float_t  AliPHOSPIDv1::GetCpv2EmcDistanceCut(TString axis, Float_t e) const
493 {
494   // Get CpvtoEmcDistance Cut depending on the cluster energy, axis and 
495   // Purity-Efficiency point 
496
497   axis.ToLower();
498   Float_t p[]={0.,0.,0.};
499   for (Int_t i=0; i<3; i++) p[i] = GetParameterCpv2Emc(i,axis);
500   Float_t sig = p[0] + TMath::Exp(p[1] - p[2]*e);
501   return sig;
502 }
503
504 //____________________________________________________________________________
505 Float_t  AliPHOSPIDv1::GetEllipseParameter(TString particle, TString param, Float_t e) const 
506 {
507   // Calculates the parameter param of the ellipse
508
509   particle.ToLower();
510   param.   ToLower();
511   Float_t p[4]={0.,0.,0.,0.};
512   Float_t value = 0.0;
513   for (Int_t i=0; i<4; i++) p[i] = GetParameterToCalculateEllipse(particle,param,i);
514   if (particle == "photon") {
515     if      (param.Contains("a"))  e = TMath::Min((Double_t)e,70.);
516     else if (param.Contains("b"))  e = TMath::Min((Double_t)e,70.);
517     else if (param.Contains("x0")) e = TMath::Max((Double_t)e,1.1);
518   }
519
520  if (particle == "photon")
521     value = p[0]/TMath::Sqrt(e) + p[1]*e + p[2]*e*e + p[3];
522   else if (particle == "pi0")
523     value = p[0] + p[1]*e + p[2]*e*e;
524
525   return value;
526 }
527
528 //_____________________________________________________________________________
529 Float_t  AliPHOSPIDv1::GetParameterPhotonBoundary (Int_t i) const
530
531   // Get the parameter "i" to calculate the boundary on the moment M2x
532   // for photons at high p_T
533   Float_t param = 0;
534   if (i>3 || i<0) {
535     AliError(Form("Wrong parameter number: %d\n",i));
536   } else
537     param = (*fParameters)(14,i) ;
538   return param;
539 }
540
541 //____________________________________________________________________________
542 Float_t  AliPHOSPIDv1::GetParameterPi0Boundary (Int_t i) const
543
544   // Get the parameter "i" to calculate the boundary on the moment M2x
545   // for pi0 at high p_T
546   Float_t param = 0;
547   if (i>2 || i<0) {
548     AliError(Form("Wrong parameter number: %d\n",i));
549   } else
550     param = (*fParameters)(15,i) ;
551   return param;
552 }
553
554 //____________________________________________________________________________
555 Float_t  AliPHOSPIDv1::GetParameterTimeGate(Int_t i) const
556 {
557   // Get TimeGate parameter depending on Purity-Efficiency i:
558   // i=0 - Low purity, i=1 - Medium purity, i=2 - High purity
559   Float_t param = 0.;
560   if(i>2 || i<0) {
561     AliError(Form("Invalid Efficiency-Purity choice %d",i));
562   } else
563     param = (*fParameters)(3,i) ; 
564   return param;
565 }
566
567 //_____________________________________________________________________________
568 Float_t  AliPHOSPIDv1::GetParameterToCalculateEllipse(TString particle, TString param, Int_t i) const
569
570   // Get the parameter "i" that is needed to calculate the ellipse 
571   // parameter "param" for the particle "particle" ("photon" or "pi0")
572
573   particle.ToLower();
574   param.   ToLower();
575   Int_t offset = -1;
576   if      (particle == "photon") 
577     offset=0;
578   else if (particle == "pi0")    
579     offset=5;
580   else
581     AliError(Form("Wrong particle name: %s (choose from pi0/photon)\n",
582                   particle.Data()));
583
584   Int_t p= -1;
585   Float_t par = 0;
586
587   if     (param.Contains("a")) p=4+offset; 
588   else if(param.Contains("b")) p=5+offset; 
589   else if(param.Contains("c")) p=6+offset; 
590   else if(param.Contains("x0"))p=7+offset; 
591   else if(param.Contains("y0"))p=8+offset;
592
593   if      (i>4 || i<0) {
594     AliError(Form("No parameter with index %d", i)) ; 
595   } else if (p==-1) {
596     AliError(Form("No parameter with name %s", param.Data() )) ; 
597   } else
598     par = (*fParameters)(p,i) ;
599   
600   return par;
601 }
602
603
604 //____________________________________________________________________________
605 Float_t  AliPHOSPIDv1::GetDistance(AliPHOSEmcRecPoint * emc,AliPHOSCpvRecPoint * cpv, Option_t *  axis)const
606 {
607   // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
608   
609   const AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ; 
610   TVector3 vecEmc ;
611   TVector3 vecCpv ;
612   if(cpv){
613     emc->GetLocalPosition(vecEmc) ;
614     cpv->GetLocalPosition(vecCpv) ; 
615     
616     if(emc->GetPHOSMod() == cpv->GetPHOSMod()){      
617       // Correct to difference in CPV and EMC position due to different distance to center.
618       // we assume, that particle moves from center
619       Float_t dCPV = geom->GetIPtoOuterCoverDistance();
620       Float_t dEMC = geom->GetIPtoCrystalSurface() ;
621       dEMC         = dEMC / dCPV ;
622       vecCpv = dEMC * vecCpv  - vecEmc ; 
623       if (axis == "X") return vecCpv.X();
624       if (axis == "Y") return vecCpv.Y();
625       if (axis == "Z") return vecCpv.Z();
626       if (axis == "R") return vecCpv.Mag();
627     }
628     return 100000000 ;
629   }
630   return 100000000 ;
631 }
632 //____________________________________________________________________________
633 Int_t  AliPHOSPIDv1::GetCPVBit(AliPHOSEmcRecPoint * emc,AliPHOSCpvRecPoint * cpv, Int_t effPur, Float_t e) const
634 {
635   //Calculates the pid bit for the CPV selection per each purity.
636   if(effPur>2 || effPur<0)
637     AliError(Form("Invalid Efficiency-Purity choice %d",effPur));
638   
639   Float_t sigX = GetCpv2EmcDistanceCut("X",e);
640   Float_t sigZ = GetCpv2EmcDistanceCut("Z",e);
641   
642   Float_t deltaX = TMath::Abs(GetDistance(emc, cpv,  "X"));
643   Float_t deltaZ = TMath::Abs(GetDistance(emc, cpv,  "Z"));
644   //Info("GetCPVBit"," xdist %f, sigx %f, zdist %f, sigz %f",deltaX, sigX, deltaZ,sigZ) ;
645  
646   //if(deltaX>sigX*(effPur+1))
647   //if((deltaX>sigX*(effPur+1)) || (deltaZ>sigZ*(effPur+1)))
648   if((deltaX>sigX*(effPur+1)) && (deltaZ>sigZ*(effPur+1)))
649     return 1;//Neutral
650   else
651     return 0;//Charged
652 }
653
654 //____________________________________________________________________________
655 Int_t  AliPHOSPIDv1::GetPrincipalBit(TString particle, const Double_t* p, Int_t effPur, Float_t e)const
656 {
657   //Is the particle inside de PCA ellipse?
658   
659   particle.ToLower();
660   Int_t    prinbit  = 0 ;
661   Float_t a  = GetEllipseParameter(particle,"a" , e); 
662   Float_t b  = GetEllipseParameter(particle,"b" , e);
663   Float_t c  = GetEllipseParameter(particle,"c" , e);
664   Float_t x0 = GetEllipseParameter(particle,"x0", e); 
665   Float_t y0 = GetEllipseParameter(particle,"y0", e);
666   
667   Float_t r = TMath::Power((p[0] - x0)/a,2) + 
668               TMath::Power((p[1] - y0)/b,2) +
669             c*(p[0] -  x0)*(p[1] - y0)/(a*b) ;
670   //3 different ellipses defined
671   if((effPur==2) && (r<1./2.)) prinbit= 1;
672   if((effPur==1) && (r<2.   )) prinbit= 1;
673   if((effPur==0) && (r<9./2.)) prinbit= 1;
674
675   if(r<0)
676     AliError("Negative square?") ;
677
678   return prinbit;
679
680 }
681 //____________________________________________________________________________
682 Int_t  AliPHOSPIDv1::GetHardPhotonBit(AliPHOSEmcRecPoint * emc) const
683 {
684   // Set bit for identified hard photons (E > 30 GeV)
685   // if the second moment M2x is below the boundary
686
687   Float_t e   = emc->GetEnergy();
688   if (e < 30.0) return 0;
689   Float_t m2x = emc->GetM2x();
690   Float_t m2xBoundary = GetParameterPhotonBoundary(0) *
691     TMath::Exp(-TMath::Power(e-GetParameterPhotonBoundary(1),2)/2.0/
692                 TMath::Power(GetParameterPhotonBoundary(2),2)) +
693     GetParameterPhotonBoundary(3);
694   AliDebug(1, Form("GetHardPhotonBit","E=%f, m2x=%f, boundary=%f",
695                        e,m2x,m2xBoundary));
696   if (m2x < m2xBoundary)
697     return 1;// A hard photon
698   else
699     return 0;// Not a hard photon
700 }
701
702 //____________________________________________________________________________
703 Int_t  AliPHOSPIDv1::GetHardPi0Bit(AliPHOSEmcRecPoint * emc) const
704 {
705   // Set bit for identified hard pi0  (E > 30 GeV)
706   // if the second moment M2x is above the boundary
707
708   Float_t e   = emc->GetEnergy();
709   if (e < 30.0) return 0;
710   Float_t m2x = emc->GetM2x();
711   Float_t m2xBoundary = GetParameterPi0Boundary(0) +
712                     e * GetParameterPi0Boundary(1);
713   AliDebug(1,Form("E=%f, m2x=%f, boundary=%f",e,m2x,m2xBoundary));
714   if (m2x > m2xBoundary)
715     return 1;// A hard pi0
716   else
717     return 0;// Not a hard pi0
718 }
719
720 //____________________________________________________________________________
721 TVector3 AliPHOSPIDv1::GetMomentumDirection(AliPHOSEmcRecPoint * emc, AliPHOSCpvRecPoint * )const 
722
723   // Calculates the momentum direction:
724   //   1. if only a EMC RecPoint, direction is given by IP and this RecPoint
725   //   2. if a EMC RecPoint and CPV RecPoint, direction is given by the line through the 2 recpoints 
726   //  However because of the poor position resolution of PPSD the direction is always taken as if we were 
727   //  in case 1.
728
729   TVector3 dir(0,0,0) ; 
730   
731   TVector3 emcglobalpos ;
732   TMatrix  dummy ;
733   
734   emc->GetGlobalPosition(emcglobalpos, dummy) ;
735
736   dir = emcglobalpos ;  
737
738   //account correction to the position of IP
739   Float_t xo,yo,zo ; //Coordinates of the origin
740   if(gAlice && gAlice->GetMCApp() && gAlice->Generator()){
741     gAlice->Generator()->GetOrigin(xo,yo,zo) ;
742   }
743   else{
744     xo=yo=zo=0.;
745   }
746   TVector3 origin(xo,yo,zo);
747   dir = dir - origin ;
748   dir.SetMag(1.) ;
749
750   return dir ;  
751 }
752
753 //________________________________________________________________________
754 Double_t  AliPHOSPIDv1::LandauF(Double_t  x, Double_t y, Double_t * par)
755 {
756   //Given the energy x and the parameter y (tof, shower dispersion or cpv-emc distance), 
757   //this method returns a density probability of this parameter, given by a landau 
758   //function whose parameters depend with the energy  with a function: a/(x*x)+b/x+b
759
760   Double_t cnt    = par[1] / (x*x) + par[2] / x + par[0] ;
761   Double_t mean   = par[4] / (x*x) + par[5] / x + par[3] ;
762   Double_t sigma  = par[7] / (x*x) + par[8] / x + par[6] ;
763
764   if(TMath::Abs(sigma) > 1.e-10){
765     return cnt*TMath::Landau(y,mean,sigma);
766   }
767   else
768     return 0.;
769
770 }
771 //________________________________________________________________________
772 Double_t  AliPHOSPIDv1::LandauPol2(Double_t  x, Double_t y, Double_t * par)
773 {
774
775   //Given the energy x and the parameter y (tof, shower dispersion or cpv-emc distance), 
776   //this method returns a density probability of this parameter, given by a landau 
777   //function whose parameters depend with the energy like second order polinomial
778
779   Double_t cnt    = par[2] * (x*x) + par[1] * x + par[0] ;
780   Double_t mean   = par[5] * (x*x) + par[4] * x + par[3] ;
781   Double_t sigma  = par[8] * (x*x) + par[7] * x + par[6] ;
782
783    if(TMath::Abs(sigma) > 1.e-10){
784     return cnt*TMath::Landau(y,mean,sigma);
785   }
786   else
787     return 0.;
788
789
790 }
791 // //________________________________________________________________________
792 // Double_t  AliPHOSPIDv1::ChargedHadronDistProb(Double_t  x, Double_t y, Double_t * parg, Double_t * parl)
793 // {
794 //   Double_t cnt   = 0.0 ;
795 //   Double_t mean  = 0.0 ;
796 //   Double_t sigma = 0.0 ;
797 //   Double_t arg   = 0.0 ;
798 //   if (y < parl[4] / (x*x) + parl[5] / x + parl[3]){
799 //     cnt    = parg[1] / (x*x) + parg[2] / x + parg[0] ;
800 //     mean   = parg[4] / (x*x) + parg[5] / x + parg[3] ;
801 //     sigma  = parg[7] / (x*x) + parg[8] / x + parg[6] ;
802 //     TF1 * f = new TF1("gaus","gaus",0.,100.);
803 //     f->SetParameters(cnt,mean,sigma);
804 //     arg  = f->Eval(y) ;
805 //   }
806 //   else{
807 //     cnt    = parl[1] / (x*x) + parl[2] / x + parl[0] ;
808 //     mean   = parl[4] / (x*x) + parl[5] / x + parl[3] ;
809 //     sigma  = parl[7] / (x*x) + parl[8] / x + parl[6] ;
810 //     TF1 * f = new TF1("landau","landau",0.,100.);
811 //     f->SetParameters(cnt,mean,sigma);
812 //     arg  = f->Eval(y) ;
813 //   }
814 //   //  Double_t mean   = par[3] + par[4] * x + par[5] * x * x ;
815 //   //   Double_t sigma  = par[6] + par[7] * x + par[8] * x * x ;
816   
817 //   //Double_t arg    = -(y-mean)*(y-mean)/(2*sigma*sigma) ;
818 //   //return cnt * TMath::Exp(arg) ;
819   
820 //   return arg;
821   
822 // }
823 //____________________________________________________________________________
824 void  AliPHOSPIDv1::MakePID()
825 {
826   // construct the PID weight from a Bayesian Method
827   
828   const Int_t kSPECIES = AliPID::kSPECIESN ;
829  
830   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
831
832   Int_t nparticles = gime->RecParticles()->GetEntriesFast() ;
833
834   //   const Int_t kMAXPARTICLES = 2000 ; 
835   //   if (nparticles >= kMAXPARTICLES) 
836   //     Error("MakePID", "Change size of MAXPARTICLES") ; 
837   //   Double_t stof[kSPECIES][kMAXPARTICLES] ;
838   
839 //   const Int_t kMAXPARTICLES = 2000 ; 
840 //   if (nparticles >= kMAXPARTICLES) 
841 //     AliError("Change size of MAXPARTICLES") ; 
842 //   Double_t stof[kSPECIES][kMAXPARTICLES] ;
843
844
845   // make the normalized distribution of pid for this event 
846   // w(pid) in the Bayesian formulation
847 //   for(index = 0 ; index < nparticles ; index ++) {
848     
849 //     cout<<">>>>>>>>>>>>>>>Bayes Index "<<index<<endl;
850
851
852 //     AliPHOSEmcRecPoint * emc     = AliPHOSGetter::Instance()->EmcRecPoint(index) ;
853 //     AliPHOSCpvRecPoint * cpv     = AliPHOSGetter::Instance()->CpvRecPoint(index) ;
854
855   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
856   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
857   TClonesArray * trackSegments = gime->TrackSegments() ;
858   if ( !emcRecPoints || !cpvRecPoints || !trackSegments ) {
859     AliFatal("RecPoints or TrackSegments not found !") ;  
860   }
861   TIter next(trackSegments) ; 
862   AliPHOSTrackSegment * ts ; 
863   Int_t index = 0 ; 
864
865   Double_t * stof[kSPECIES] ;
866   Double_t * sdp [kSPECIES]  ;
867   Double_t * scpv[kSPECIES] ;
868   
869   //Info("MakePID","Begin MakePID"); 
870   
871   for (Int_t i =0; i< kSPECIES; i++){
872     stof[i] = new Double_t[nparticles] ;
873     sdp [i] = new Double_t[nparticles] ;
874     scpv[i] = new Double_t[nparticles] ;
875   }
876   
877
878   while ( (ts = (AliPHOSTrackSegment *)next()) ) {
879     
880     //cout<<">>>>>> Bayesian Index "<<index<<endl;
881
882     AliPHOSEmcRecPoint * emc = 0 ;
883     if(ts->GetEmcIndex()>=0)
884       emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ;
885     
886     AliPHOSCpvRecPoint * cpv = 0 ;
887     if(ts->GetCpvIndex()>=0)
888       cpv = (AliPHOSCpvRecPoint *) cpvRecPoints->At(ts->GetCpvIndex()) ;
889     
890 //     Int_t track = 0 ; 
891 //     track = ts->GetTrackIndex() ; //TPC tracks ?
892     
893     if (!emc) {
894       AliFatal(Form("-> emc(%d) = %d", ts->GetEmcIndex(), emc )) ;
895     }
896
897     // ############Tof#############################
898
899     //    Info("MakePID", "TOF");
900     Float_t en   = emc->GetEnergy();    
901     Double_t time = emc->GetTime() ;
902     //    cout<<">>>>>>>Energy "<<en<<"Time "<<time<<endl;
903     //Conversion Electrons initial population. TO BE REMOVED
904     fInitPID[AliPID::kEleCon]   = 0. ;
905
906     
907     // now get the signals probability
908     // s(pid) in the Bayesian formulation
909     
910     stof[AliPID::kPhoton][index]   = 1.; 
911     stof[AliPID::kElectron][index] = 1.;
912     stof[AliPID::kEleCon][index]   = 1.;
913     //We assing the same prob to charged hadrons, sum is 1
914     stof[AliPID::kPion][index]     = 1./3.; 
915     stof[AliPID::kKaon][index]     = 1./3.; 
916     stof[AliPID::kProton][index]   = 1./3.;
917     //We assing the same prob to neutral hadrons, sum is 1
918     stof[AliPID::kNeutron][index]  = 1./2.;
919     stof[AliPID::kKaon0][index]    = 1./2.;
920
921     stof[AliPID::kMuon][index]     = 1.; 
922
923     
924     if(en < 2.) {
925
926       Double_t pTofPion = fTFpiong ->Eval(time) ; //gaus distribution
927       Double_t pTofKaon = 0;
928
929       if(time < fTkaonl[1])
930         pTofKaon = fTFkaong      ->Eval(time) ; //gaus distribution
931       else 
932         pTofKaon = fTFkaonl      ->Eval(time) ;  //landau distribution
933
934       Double_t pTofNucleon = 0;
935
936       if(time < fThhadronl[1])
937         pTofNucleon = fTFhhadrong   ->Eval(time) ; //gaus distribution
938       else
939         pTofNucleon = fTFhhadronl   ->Eval(time) ; //landau distribution
940       //We assing the same prob to charged hadrons, sum is the average prob
941       Double_t pTofNeHadron =  (pTofPion + pTofKaon + pTofNucleon)/2. ;
942       //We assing the same prob to neutral hadrons, sum is the average prob
943       Double_t pTofChHadron =  (pTofKaon + pTofNucleon)/3. ;
944
945       stof[AliPID::kPhoton][index]   = fTFphoton     ->Eval(time) ; //gaus distribution
946       stof[AliPID::kEleCon][index]   = stof[AliPID::kPhoton][index] ; // a conversion electron has the photon ToF
947       stof[AliPID::kMuon][index]     = stof[AliPID::kPhoton][index] ;
948  
949       stof[AliPID::kElectron][index] = pTofPion  ;                               
950
951       stof[AliPID::kPion][index]     =  pTofChHadron ; 
952       stof[AliPID::kKaon][index]     =  pTofChHadron ;
953       stof[AliPID::kProton][index]   =  pTofChHadron ;
954
955       stof[AliPID::kKaon0][index]    =  pTofNeHadron ;     
956       stof[AliPID::kNeutron][index]  =  pTofNeHadron ;            
957     } 
958     
959     //    Info("MakePID", "Dispersion");
960     
961     // ###########Shower shape: Dispersion####################
962     Float_t dispersion = emc->GetDispersion();
963     //dispersion is not well defined if the cluster is only in few crystals
964     
965     sdp[AliPID::kPhoton][index]   = 1. ;
966     sdp[AliPID::kElectron][index] = 1. ;
967     sdp[AliPID::kPion][index]     = 1. ; 
968     sdp[AliPID::kKaon][index]     = 1. ; 
969     sdp[AliPID::kProton][index]   = 1. ;
970     sdp[AliPID::kNeutron][index]  = 1. ;
971     sdp[AliPID::kEleCon][index]   = 1. ; 
972     sdp[AliPID::kKaon0][index]    = 1. ; 
973     sdp[AliPID::kMuon][index]     = 1. ; 
974     
975     if(en > 0.5 && emc->GetMultiplicity() > 3){ 
976       sdp[AliPID::kPhoton][index]   = GausF(en , dispersion, fDphoton) ;
977       sdp[AliPID::kElectron][index] = sdp[AliPID::kPhoton][index] ;
978       sdp[AliPID::kPion][index]     = LandauF(en , dispersion, fDhadron ) ; 
979       sdp[AliPID::kKaon][index]     = sdp[AliPID::kPion][index]  ; 
980       sdp[AliPID::kProton][index]   = sdp[AliPID::kPion][index]  ;
981       sdp[AliPID::kNeutron][index]  = sdp[AliPID::kPion][index]  ;
982       sdp[AliPID::kEleCon][index]   = sdp[AliPID::kPhoton][index]; 
983       sdp[AliPID::kKaon0][index]    = sdp[AliPID::kPion][index]  ; 
984       sdp[AliPID::kMuon][index]     = fDFmuon ->Eval(dispersion) ; //landau distribution
985     }
986     
987 //      Info("MakePID","multiplicity %d, dispersion %f", emc->GetMultiplicity(), dispersion);
988 //      Info("MakePID","ss: photon %f, hadron %f ",  sdp[AliPID::kPhoton][index],  sdp[AliPID::kPion][index]);
989 //       cout<<">>>>>multiplicity "<<emc->GetMultiplicity()<<", dispersion "<< dispersion<<endl ;
990 //       cout<<"<<<<<ss: photon   "<<sdp[AliPID::kPhoton][index]<<", hadron    "<<sdp[AliPID::kPion][index]<<endl;
991
992     //########## CPV-EMC  Distance#######################
993     //     Info("MakePID", "Distance");
994     //    Float_t distance = GetDistance(emc, cpv,  "R") ;
995     Float_t x             = TMath::Abs(GetDistance(emc, cpv,  "X")) ;
996     Float_t z             = GetDistance(emc, cpv,  "Z") ;
997     //    Info("MakePID", "Distance %f", distance);
998     Double_t pcpv         = 0 ;
999     Double_t pcpvneutral  = 0. ;
1000    
1001     Double_t elprobx      = GausF(en , x, fXelectron) ;
1002     Double_t elprobz      = GausF(en , z, fZelectron) ;
1003     Double_t chprobx      = GausF(en , x, fXcharged)  ;
1004     Double_t chprobz      = GausF(en , z, fZcharged)  ;
1005     Double_t pcpvelectron = elprobx * elprobz;
1006     Double_t pcpvcharged  = chprobx * chprobz;
1007   
1008 //     cout<<">>>>energy "<<en<<endl;
1009 //     cout<<">>>>electron : x "<<x<<" xprob "<<elprobx<<" z "<<z<<" zprob "<<elprobz<<endl;
1010 //     cout<<">>>>hadron   : x "<<x<<" xprob "<<chprobx<<" z "<<z<<" zprob "<<chprobz<<endl;
1011 //     cout<<">>>>electron : px*pz "<<pcpvelectron <<" hadron: px*pz "<<pcpvcharged<<endl;  
1012
1013     // Is neutral or charged
1014     if(pcpvelectron >= pcpvcharged)  
1015       pcpv = pcpvelectron ;
1016     else
1017       pcpv = pcpvcharged ;
1018     
1019     if(pcpv < 1e-7)
1020       {
1021         pcpvneutral  = 1. ;
1022         pcpvcharged  = 0. ;
1023         pcpvelectron = 0. ;
1024       }
1025     //    else
1026     //      cout<<">>>>>>>>>>>CHARGED>>>>>>>>>>>"<<endl;
1027     
1028     scpv[AliPID::kPion][index]     =  pcpvcharged  ; 
1029     scpv[AliPID::kKaon][index]     =  pcpvcharged  ; 
1030     scpv[AliPID::kProton][index]   =  pcpvcharged  ;
1031
1032     scpv[AliPID::kMuon][index]     =  pcpvelectron ; 
1033     scpv[AliPID::kElectron][index] =  pcpvelectron ;
1034     scpv[AliPID::kEleCon][index]   =  pcpvelectron ; 
1035
1036     scpv[AliPID::kPhoton][index]   =  pcpvneutral  ;
1037     scpv[AliPID::kNeutron][index]  =  pcpvneutral  ; 
1038     scpv[AliPID::kKaon0][index]    =  pcpvneutral  ; 
1039
1040     
1041     //   Info("MakePID", "CPV passed");
1042
1043     //############## Pi0 #############################
1044     stof[AliPID::kPi0][index]      = 0. ;  
1045     scpv[AliPID::kPi0][index]      = 0. ;
1046     sdp [AliPID::kPi0][index]      = 0. ;
1047
1048     if(en > 30.){
1049       // pi0 are detected via decay photon
1050       stof[AliPID::kPi0][index]  = fTFphoton  ->Eval(time) ;
1051       scpv[AliPID::kPi0][index]  = pcpvneutral  ;
1052       sdp [AliPID::kPi0][index]  = 1. ;
1053       if(emc->GetMultiplicity() > 3)
1054         sdp [AliPID::kPi0][index]  = GausPol2(en , dispersion, fDpi0) ;
1055     }
1056     
1057
1058     //############## muon #############################
1059
1060     if(en > 0.5){
1061       //Muons deposit few energy
1062       scpv[AliPID::kMuon][index]     =  0 ;
1063       stof[AliPID::kMuon][index]     =  0 ;
1064       sdp [AliPID::kMuon][index]     =  0 ;
1065     }
1066
1067 //     if(en > 0.5){
1068 //       cout<<"######################################################"<<endl;
1069 //       //cout<<"MakePID: energy "<<en<<", tof "<<time<<", distance "<<distance<<", dispersion "<<dispersion<<endl ;
1070 //       cout<<"MakePID: energy "<<en<<", tof "<<time<<", dispersion "<<dispersion<<", x "<<x<<", z "<<z<<endl ;
1071 //       cout<<">>>>>multiplicity "<<emc->GetMultiplicity()<<endl;
1072 //       cout<<">>>>electron : xprob "<<elprobx<<" zprob "<<elprobz<<endl;
1073 //       cout<<">>>>hadron   : xprob "<<chprobx<<" zprob "<<chprobz<<endl;
1074 //       cout<<">>>>electron : px*pz "<<pcpvelectron <<" hadron: px*pz "<<pcpvcharged<<endl;  
1075       
1076 //       cout<<"Photon   , pid "<< fInitPID[AliPID::kPhoton]<<" tof "<<stof[AliPID::kPhoton][index]
1077 //        <<", cpv "<<scpv[AliPID::kPhoton][index]<<", ss "<<sdp[AliPID::kPhoton][index]<<endl;
1078 //       cout<<"EleCon   , pid "<< fInitPID[AliPID::kEleCon]<<", tof "<<stof[AliPID::kEleCon][index]
1079 //        <<", cpv "<<scpv[AliPID::kEleCon][index]<<" ss "<<sdp[AliPID::kEleCon][index]<<endl;
1080 //       cout<<"Electron , pid "<< fInitPID[AliPID::kElectron]<<", tof "<<stof[AliPID::kElectron][index]
1081 //        <<", cpv "<<scpv[AliPID::kElectron][index]<<" ss "<<sdp[AliPID::kElectron][index]<<endl;
1082 //       cout<<"Muon     , pid "<< fInitPID[AliPID::kMuon]<<", tof "<<stof[AliPID::kMuon][index]
1083 //        <<", cpv "<<scpv[AliPID::kMuon][index]<<" ss "<<sdp[AliPID::kMuon][index]<<endl;
1084 //       cout<<"Pi0      , pid "<< fInitPID[AliPID::kPi0]<<", tof "<<stof[AliPID::kPi0][index]
1085 //        <<", cpv "<<scpv[AliPID::kPi0][index]<<" ss "<<sdp[AliPID::kPi0][index]<<endl;
1086 //       cout<<"Pion     , pid "<< fInitPID[AliPID::kPion]<<", tof "<<stof[AliPID::kPion][index]
1087 //        <<", cpv "<<scpv[AliPID::kPion][index]<<" ss "<<sdp[AliPID::kPion][index]<<endl;
1088 //       cout<<"Kaon0    , pid "<< fInitPID[AliPID::kKaon0]<<", tof "<<stof[AliPID::kKaon0][index]
1089 //        <<", cpv "<<scpv[AliPID::kKaon0][index]<<" ss "<<sdp[AliPID::kKaon0][index]<<endl;
1090 //       cout<<"Kaon     , pid "<< fInitPID[AliPID::kKaon]<<", tof "<<stof[AliPID::kKaon][index]
1091 //        <<", cpv "<<scpv[AliPID::kKaon][index]<<" ss "<<sdp[AliPID::kKaon][index]<<endl;
1092 //       cout<<"Neutron  , pid "<< fInitPID[AliPID::kNeutron]<<", tof "<<stof[AliPID::kNeutron][index]
1093 //        <<", cpv "<<scpv[AliPID::kNeutron][index]<<" ss "<<sdp[AliPID::kNeutron][index]<<endl;
1094 //       cout<<"Proton   , pid "<< fInitPID[AliPID::kProton]<<", tof "<<stof[AliPID::kProton][index]
1095 //        <<", cpv "<<scpv[AliPID::kProton][index]<<" ss "<<sdp[AliPID::kProton][index]<<endl;
1096 //       cout<<"######################################################"<<endl;
1097 //     }
1098       index++;
1099   }
1100   
1101   //for (index = 0 ; index < kSPECIES ; index++) 
1102   // pid[index] /= nparticles ; 
1103   
1104
1105   //  Info("MakePID", "Total Probability calculation");
1106   
1107   for(index = 0 ; index < nparticles ; index ++) {
1108     // calculates the Bayesian weight
1109     
1110     Int_t jndex ;
1111     Double_t wn = 0.0 ; 
1112     for (jndex = 0 ; jndex < kSPECIES ; jndex++) 
1113       wn += stof[jndex][index] * sdp[jndex][index]  * scpv[jndex][index] * fInitPID[jndex] ;
1114    
1115     //    cout<<"*************wn "<<wn<<endl;
1116     AliPHOSRecParticle * recpar = gime->RecParticle(index) ;  
1117     if (TMath::Abs(wn)>0)
1118       for (jndex = 0 ; jndex < kSPECIES ; jndex++) {
1119         //cout<<"jndex "<<jndex<<" wn "<<wn<<" SetPID * wn"
1120         //<<stof[jndex][index] * sdp[jndex][index] * pid[jndex]  << endl;
1121         //cout<<" tof "<<stof[jndex][index] << " disp " <<sdp[jndex][index] << " pid "<< fInitPID[jndex] << endl;
1122 //      cout<<"Particle "<<jndex<<"  final prob * wn   "
1123 //          <<stof[jndex][index] * sdp[jndex][index] * scpv[jndex][index] * fInitPID[jndex] <<"  wn  "<< wn<<endl;
1124         recpar->SetPID(jndex, stof[jndex][index] * sdp[jndex][index] * 
1125                        scpv[jndex][index] * fInitPID[jndex] / wn) ; 
1126 //      cout<<"final prob "<<stof[jndex][index] * sdp[jndex][index] * scpv[jndex][index] * fInitPID[jndex] / wn<<endl;
1127         //recpar->SetPID(jndex, stof[jndex][index] * fInitPID[jndex] / wn) ; 
1128         //cout<<"After SetPID"<<endl;
1129         //recpar->Print();
1130       }
1131   }
1132   //  Info("MakePID", "Delete");
1133   
1134     for (Int_t i =0; i< kSPECIES; i++){
1135       delete [] stof[i];
1136       delete [] sdp[i];
1137       delete [] scpv[i];
1138     }
1139   //  Info("MakePID","End MakePID"); 
1140 }
1141
1142 //____________________________________________________________________________
1143 void  AliPHOSPIDv1::MakeRecParticles()
1144 {
1145   // Makes a RecParticle out of a TrackSegment
1146   
1147   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
1148   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
1149   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
1150   TClonesArray * trackSegments = gime->TrackSegments() ; 
1151   if ( !emcRecPoints || !cpvRecPoints || !trackSegments ) {
1152     AliFatal("RecPoints or TrackSegments not found !") ;  
1153   }
1154   TClonesArray * recParticles  = gime->RecParticles() ; 
1155   recParticles->Clear();
1156
1157   TIter next(trackSegments) ; 
1158   AliPHOSTrackSegment * ts ; 
1159   Int_t index = 0 ; 
1160   AliPHOSRecParticle * rp ; 
1161   while ( (ts = (AliPHOSTrackSegment *)next()) ) {
1162     //  cout<<">>>>>>>>>>>>>>>PCA Index "<<index<<endl;
1163     new( (*recParticles)[index] ) AliPHOSRecParticle() ;
1164     rp = (AliPHOSRecParticle *)recParticles->At(index) ; 
1165     rp->SetTrackSegment(index) ;
1166     rp->SetIndexInList(index) ;
1167         
1168     AliPHOSEmcRecPoint * emc = 0 ;
1169     if(ts->GetEmcIndex()>=0)
1170       emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ;
1171     
1172     AliPHOSCpvRecPoint * cpv = 0 ;
1173     if(ts->GetCpvIndex()>=0)
1174       cpv = (AliPHOSCpvRecPoint *) cpvRecPoints->At(ts->GetCpvIndex()) ;
1175     
1176     Int_t track = 0 ; 
1177     track = ts->GetTrackIndex() ; 
1178       
1179     // Now set type (reconstructed) of the particle
1180
1181     // Choose the cluster energy range
1182     
1183     if (!emc) {
1184       AliFatal(Form("-> emc(%d) = %d", ts->GetEmcIndex(), emc )) ;
1185     }
1186
1187     Float_t e = emc->GetEnergy() ;   
1188     
1189     Float_t  lambda[2] ;
1190     emc->GetElipsAxis(lambda) ;
1191     
1192     if((lambda[0]>0.01) && (lambda[1]>0.01)){
1193       // Looking PCA. Define and calculate the data (X),
1194       // introduce in the function X2P that gives the components (P).  
1195
1196       Float_t  spher = 0. ;
1197       Float_t  emaxdtotal = 0. ; 
1198       
1199       if((lambda[0]+lambda[1])!=0) 
1200         spher=fabs(lambda[0]-lambda[1])/(lambda[0]+lambda[1]); 
1201       
1202       emaxdtotal=emc->GetMaximalEnergy()/emc->GetEnergy(); 
1203       
1204       fX[0] = lambda[0] ;  
1205       fX[1] = lambda[1] ; 
1206       fX[2] = emc->GetDispersion() ; 
1207       fX[3] = spher ; 
1208       fX[4] = emc->GetMultiplicity() ;  
1209       fX[5] = emaxdtotal ;  
1210       fX[6] = emc->GetCoreEnergy() ;  
1211       
1212       fPrincipalPhoton->X2P(fX,fPPhoton);
1213       fPrincipalPi0   ->X2P(fX,fPPi0);
1214
1215     }
1216     else{
1217       fPPhoton[0]=-100.0;  //We do not accept clusters with 
1218       fPPhoton[1]=-100.0;  //one cell as a photon-like
1219       fPPi0[0]   =-100.0;
1220       fPPi0[1]   =-100.0;
1221     }
1222     
1223     Float_t time = emc->GetTime() ;
1224     rp->SetTof(time) ; 
1225     
1226     // Loop of Efficiency-Purity (the 3 points of purity or efficiency 
1227     // are taken into account to set the particle identification)
1228     for(Int_t effPur = 0; effPur < 3 ; effPur++){
1229       
1230       // Looking at the CPV detector. If RCPV greater than CpvEmcDistance, 
1231       // 1st,2nd or 3rd bit (depending on the efficiency-purity point )
1232       // is set to 1
1233       if(GetCPVBit(emc, cpv, effPur,e) == 1 ){  
1234         rp->SetPIDBit(effPur) ;
1235         //cout<<"CPV bit "<<effPur<<endl;
1236       }
1237       // Looking the TOF. If TOF smaller than gate,  4th, 5th or 6th 
1238       // bit (depending on the efficiency-purity point )is set to 1             
1239       if(time< (*fParameters)(3,effPur)) 
1240         rp->SetPIDBit(effPur+3) ;                   
1241   
1242       //Photon PCA
1243       //If we are inside the ellipse, 7th, 8th or 9th 
1244       // bit (depending on the efficiency-purity point )is set to 1 
1245       if(GetPrincipalBit("photon",fPPhoton,effPur,e) == 1) 
1246         rp->SetPIDBit(effPur+6) ;
1247
1248       //Pi0 PCA
1249       //If we are inside the ellipse, 10th, 11th or 12th 
1250       // bit (depending on the efficiency-purity point )is set to 1 
1251       if(GetPrincipalBit("pi0"   ,fPPi0   ,effPur,e) == 1) 
1252         rp->SetPIDBit(effPur+9) ;
1253     }
1254     if(GetHardPhotonBit(emc))
1255       rp->SetPIDBit(12) ;
1256     if(GetHardPi0Bit   (emc))
1257       rp->SetPIDBit(13) ;
1258     
1259     if(track >= 0) 
1260       rp->SetPIDBit(14) ; 
1261
1262     //Set momentum, energy and other parameters 
1263     Float_t  encal = GetCalibratedEnergy(e);
1264     TVector3 dir   = GetMomentumDirection(emc,cpv) ; 
1265     dir.SetMag(encal) ;
1266     rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),encal) ;
1267     rp->SetCalcMass(0);
1268     rp->Name(); //If photon sets the particle pdg name to gamma
1269     rp->SetProductionVertex(0,0,0,0);
1270     rp->SetFirstMother(-1);
1271     rp->SetLastMother(-1);
1272     rp->SetFirstDaughter(-1);
1273     rp->SetLastDaughter(-1);
1274     rp->SetPolarisation(0,0,0);
1275     //Set the position in global coordinate system from the RecPoint
1276     AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
1277     AliPHOSTrackSegment * ts  = gime->TrackSegment(rp->GetPHOSTSIndex()) ; 
1278     AliPHOSEmcRecPoint  * erp = gime->EmcRecPoint(ts->GetEmcIndex()) ; 
1279     TVector3 pos ; 
1280     geom->GetGlobal(erp, pos) ; 
1281     rp->SetPos(pos);
1282     index++ ; 
1283   }
1284 }
1285   
1286 //____________________________________________________________________________
1287 void  AliPHOSPIDv1::Print() const
1288 {
1289   // Print the parameters used for the particle type identification
1290
1291     AliInfo("=============== AliPHOSPIDv1 ================") ;
1292     printf("Making PID\n") ;
1293     printf("    Pricipal analysis file from 0.5 to 100 %s\n", fFileNamePrincipalPhoton.Data() )   ; 
1294     printf("    Name of parameters file     %s\n", fFileNameParameters.Data() )  ;
1295     printf("    Matrix of Parameters: 14x4\n") ;
1296     printf("        Energy Calibration  1x3 [3 parametres to calibrate energy: A + B* E + C * E^2]\n") ;
1297     printf("        RCPV 2x3 rows x and z, columns function cut parameters\n") ;
1298     printf("        TOF  1x3 [High Eff-Low Pur,Medium Eff-Pur, Low Eff-High Pur]\n") ;
1299     printf("        PCA  5x4 [5 ellipse parametres and 4 parametres to calculate them: A/Sqrt(E) + B* E + C * E^2 + D]\n") ;
1300     Printf("    Pi0 PCA  5x3 [5 ellipse parametres and 3 parametres to calculate them: A + B* E + C * E^2]\n") ;
1301     fParameters->Print() ;
1302 }
1303
1304
1305
1306 //____________________________________________________________________________
1307 void AliPHOSPIDv1::PrintRecParticles(Option_t * option)
1308 {
1309   // Print table of reconstructed particles
1310
1311   AliPHOSGetter *gime = AliPHOSGetter::Instance() ; 
1312
1313   TClonesArray * recParticles = gime->RecParticles() ; 
1314
1315   TString message ; 
1316   message  = "\nevent " ;
1317   message += gAlice->GetEvNumber() ; 
1318   message += "       found " ; 
1319   message += recParticles->GetEntriesFast(); 
1320   message += " RecParticles\n" ; 
1321
1322   if(strstr(option,"all")) {  // printing found TS
1323     message += "\n  PARTICLE         Index    \n" ; 
1324     
1325     Int_t index ;
1326     for (index = 0 ; index < recParticles->GetEntries() ; index++) {
1327       AliPHOSRecParticle * rp = (AliPHOSRecParticle * ) recParticles->At(index) ;       
1328       message += "\n" ;
1329       message += rp->Name().Data() ;  
1330       message += " " ;
1331       message += rp->GetIndexInList() ;  
1332       message += " " ;
1333       message += rp->GetType()  ;
1334     }
1335   }
1336   AliInfo(message.Data() ) ; 
1337 }
1338
1339 //____________________________________________________________________________
1340 void  AliPHOSPIDv1::SetParameters() 
1341 {
1342   // PCA : To do the Principal Components Analysis it is necessary 
1343   // the Principal file, which is opened here
1344   fX       = new double[7]; // Data for the PCA 
1345   fPPhoton = new double[7]; // Eigenvalues of the PCA
1346   fPPi0    = new double[7]; // Eigenvalues of the Pi0 PCA
1347
1348   // Read photon principals from the photon file
1349   
1350   fFileNamePrincipalPhoton = "$ALICE_ROOT/PHOS/PCA8pa15_0.5-100.root" ; 
1351   TFile f( fFileNamePrincipalPhoton.Data(), "read" ) ;
1352   fPrincipalPhoton = dynamic_cast<TPrincipal*> (f.Get("principal")) ; 
1353   f.Close() ; 
1354
1355   // Read pi0 principals from the pi0 file
1356
1357   fFileNamePrincipalPi0    = "$ALICE_ROOT/PHOS/PCA_pi0_40-120.root" ;
1358   TFile fPi0( fFileNamePrincipalPi0.Data(), "read" ) ;
1359   fPrincipalPi0    = dynamic_cast<TPrincipal*> (fPi0.Get("principal")) ; 
1360   fPi0.Close() ;
1361
1362   // Open parameters file and initialization of the Parameters matrix. 
1363   // In the File Parameters.dat are all the parameters. These are introduced 
1364   // in a matrix of 16x4  
1365   // 
1366   // All the parameters defined in this file are, in order of row: 
1367   // line   0   : calibration 
1368   // lines  1,2 : CPV rectangular cat for X and Z
1369   // line   3   : TOF cut
1370   // lines  4-8 : parameters to calculate photon PCA ellipse
1371   // lines  9-13: parameters to calculate pi0 PCA ellipse
1372   // lines 14-15: parameters to calculate border for high-pt photons and pi0
1373
1374   fFileNameParameters = gSystem->ExpandPathName("$ALICE_ROOT/PHOS/Parameters.dat");
1375   fParameters = new TMatrix(16,4) ;
1376   const Int_t kMaxLeng=255;
1377   char string[kMaxLeng];
1378
1379   // Open a text file with PID parameters
1380   FILE *fd = fopen(fFileNameParameters.Data(),"r");
1381   if (!fd)
1382     AliFatal(Form("File %s with a PID parameters cannot be opened\n",
1383           fFileNameParameters.Data()));
1384
1385   Int_t i=0;
1386   // Read parameter file line-by-line and skip empty line and comments
1387   while (fgets(string,kMaxLeng,fd) != NULL) {
1388     if (string[0] == '\n' ) continue;
1389     if (string[0] == '!'  ) continue;
1390     sscanf(string, "%f %f %f %f",
1391            &(*fParameters)(i,0), &(*fParameters)(i,1), 
1392            &(*fParameters)(i,2), &(*fParameters)(i,3));
1393     i++;
1394     AliDebug(1, Form("SetParameters", "line %d: %s",i,string));
1395   }
1396   fclose(fd);
1397 }
1398
1399 //____________________________________________________________________________
1400 void  AliPHOSPIDv1::SetParameterCalibration(Int_t i,Float_t param) 
1401 {
1402   // Set parameter "Calibration" i to a value param
1403   if(i>2 || i<0) {
1404     AliError(Form("Invalid parameter number: %d",i));
1405   } else
1406     (*fParameters)(0,i) = param ;
1407 }
1408
1409 //____________________________________________________________________________
1410 void  AliPHOSPIDv1::SetParameterCpv2Emc(Int_t i, TString axis, Float_t cut) 
1411 {
1412   // Set the parameters to calculate Cpv-to-Emc Distance Cut depending on 
1413   // Purity-Efficiency point i
1414
1415   if(i>2 || i<0) {
1416     AliError(Form("Invalid parameter number: %d",i));
1417   } else {
1418     axis.ToLower();
1419     if      (axis == "x") (*fParameters)(1,i) = cut;
1420     else if (axis == "z") (*fParameters)(2,i) = cut;
1421     else { 
1422       AliError(Form("Invalid axis name: %s",axis.Data()));
1423     }
1424   }
1425 }
1426
1427 //____________________________________________________________________________
1428 void  AliPHOSPIDv1::SetParameterPhotonBoundary(Int_t i,Float_t param) 
1429 {
1430   // Set parameter "Hard photon boundary" i to a value param
1431   if(i>4 || i<0) {
1432     AliError(Form("Invalid parameter number: %d",i));
1433   } else
1434     (*fParameters)(14,i) = param ;
1435 }
1436
1437 //____________________________________________________________________________
1438 void  AliPHOSPIDv1::SetParameterPi0Boundary(Int_t i,Float_t param) 
1439 {
1440   // Set parameter "Hard pi0 boundary" i to a value param
1441   if(i>1 || i<0) {
1442     AliError(Form("Invalid parameter number: %d",i));
1443   } else
1444     (*fParameters)(15,i) = param ;
1445 }
1446
1447 //_____________________________________________________________________________
1448 void  AliPHOSPIDv1::SetParameterTimeGate(Int_t i, Float_t gate) 
1449 {
1450   // Set the parameter TimeGate depending on Purity-Efficiency point i 
1451   if (i>2 || i<0) {
1452     AliError(Form("Invalid Efficiency-Purity choice %d",i));
1453   } else
1454     (*fParameters)(3,i)= gate ; 
1455
1456
1457 //_____________________________________________________________________________
1458 void  AliPHOSPIDv1::SetParameterToCalculateEllipse(TString particle, TString param, Int_t i, Float_t par) 
1459 {  
1460   // Set the parameter "i" that is needed to calculate the ellipse 
1461   // parameter "param" for a particle "particle"
1462   
1463   particle.ToLower();
1464   param.   ToLower();
1465   Int_t p= -1;
1466   Int_t offset=0;
1467
1468   if      (particle == "photon") offset=0;
1469   else if (particle == "pi0")    offset=5;
1470   else
1471     AliError(Form("Wrong particle name: %s (choose from pi0/photon)\n",
1472                   particle.Data()));
1473
1474   if     (param.Contains("a")) p=4+offset; 
1475   else if(param.Contains("b")) p=5+offset; 
1476   else if(param.Contains("c")) p=6+offset; 
1477   else if(param.Contains("x0"))p=7+offset; 
1478   else if(param.Contains("y0"))p=8+offset;
1479   if((i>4)||(i<0)) {
1480     AliError(Form("No parameter with index %d", i)) ; 
1481   } else if(p==-1) {
1482     AliError(Form("No parameter with name %s", param.Data() )) ; 
1483   } else
1484     (*fParameters)(p,i) = par ;
1485
1486
1487 //____________________________________________________________________________
1488 void AliPHOSPIDv1::Unload() 
1489 {
1490   //Unloads RecPoints, Tracks and RecParticles
1491   AliPHOSGetter * gime = AliPHOSGetter::Instance() ;  
1492   gime->PhosLoader()->UnloadRecPoints() ;
1493   gime->PhosLoader()->UnloadTracks() ;
1494   gime->PhosLoader()->UnloadRecParticles() ;
1495 }
1496
1497 //____________________________________________________________________________
1498 void  AliPHOSPIDv1::WriteRecParticles()
1499 {
1500   //It writes reconstructed particles and pid to file
1501
1502   AliPHOSGetter *gime = AliPHOSGetter::Instance() ; 
1503
1504   TClonesArray * recParticles = gime->RecParticles() ; 
1505   recParticles->Expand(recParticles->GetEntriesFast() ) ;
1506   if(fWrite){
1507     TTree * treeP =  gime->TreeP();
1508     
1509     //First rp
1510     Int_t bufferSize = 32000 ;
1511     TBranch * rpBranch = treeP->Branch("PHOSRP",&recParticles,bufferSize);
1512     rpBranch->SetTitle(BranchName());
1513     
1514     rpBranch->Fill() ;
1515     
1516     gime->WriteRecParticles("OVERWRITE");
1517     gime->WritePID("OVERWRITE");
1518   }
1519 }
1520
1521
1522 //_______________________________________________________________________
1523 void AliPHOSPIDv1::SetInitPID(const Double_t *p) {
1524   // Sets values for the initial population of each particle type 
1525   for (Int_t i=0; i<AliPID::kSPECIESN; i++) fInitPID[i] = p[i];
1526 }
1527 //_______________________________________________________________________
1528 void AliPHOSPIDv1::GetInitPID(Double_t *p) const {
1529   // Gets values for the initial population of each particle type 
1530   for (Int_t i=0; i<AliPID::kSPECIESN; i++) p[i] = fInitPID[i];
1531 }