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