]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSPIDv1.cxx
Methods made virtual to allow On-Flight reconstruction with Light version of getter
[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 #include "TROOT.h"
83 #include "TTree.h"
84 #include "TFile.h"
85 #include "TF2.h"
86 #include "TFormula.h"
87 #include "TCanvas.h"
88 #include "TFolder.h"
89 #include "TSystem.h"
90 #include "TBenchmark.h"
91 #include "TMatrixD.h"
92 #include "TPrincipal.h"
93 #include "TSystem.h"
94
95 // --- Standard library ---
96
97
98 // --- AliRoot header files ---
99
100 #include "AliGenerator.h"
101 #include "AliPHOS.h"
102 #include "AliPHOSPIDv1.h"
103 #include "AliPHOSClusterizerv1.h"
104 #include "AliPHOSEmcRecPoint.h"
105 #include "AliPHOSTrackSegment.h"
106 #include "AliPHOSTrackSegmentMakerv1.h"
107 #include "AliPHOSRecParticle.h"
108 #include "AliPHOSGeometry.h"
109 #include "AliPHOSGetter.h"
110
111 ClassImp( AliPHOSPIDv1) 
112
113 //____________________________________________________________________________
114 AliPHOSPIDv1::AliPHOSPIDv1():AliPHOSPID()
115
116   // default ctor
117  
118   InitParameters() ; 
119   fDefaultInit = kTRUE ; 
120 }
121
122 //____________________________________________________________________________
123 AliPHOSPIDv1::AliPHOSPIDv1(const AliPHOSPIDv1 & pid ):AliPHOSPID(pid)
124
125   // ctor
126   InitParameters() ; 
127   Init() ;
128
129 }
130
131 //____________________________________________________________________________
132 AliPHOSPIDv1::AliPHOSPIDv1(const TString alirunFileName, const TString eventFolderName):AliPHOSPID(alirunFileName, eventFolderName)
133
134   //ctor with the indication on where to look for the track segments
135  
136   InitParameters() ; 
137   Init() ;
138   fDefaultInit = kFALSE ; 
139 }
140
141 //____________________________________________________________________________
142 AliPHOSPIDv1::~AliPHOSPIDv1()
143
144   // dtor
145
146   delete [] fX ;       // Principal input 
147   delete [] fPPhoton ; // Photon Principal components
148   delete [] fPPi0 ;    // Pi0 Principal components
149 }
150 //____________________________________________________________________________
151 const TString AliPHOSPIDv1::BranchName() const 
152 {  
153
154   return GetName() ;
155 }
156  
157 //____________________________________________________________________________
158 void AliPHOSPIDv1::Init()
159 {
160   // Make all memory allocations that are not possible in default constructor
161   // Add the PID task to the list of PHOS tasks
162
163   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
164   if(!gime)
165     gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data()) ; 
166
167   if ( !gime->PID() ) 
168     gime->PostPID(this) ;
169 }
170
171 //____________________________________________________________________________
172 void AliPHOSPIDv1::InitParameters()
173 {
174   // Initialize PID parameters
175   fWrite                   = kTRUE ;
176   fRecParticlesInRun = 0 ; 
177   fNEvent            = 0 ;            
178   fRecParticlesInRun = 0 ;
179   SetParameters() ; // fill the parameters matrix from parameters file
180   SetEventRange(0,-1) ;
181   // initialisation of response function parameters
182   // Tof
183   // Photons
184   fTphoton[0] = 2.97E-1 ; 
185   fTphoton[1] = 1.55E-8 ; 
186   fTphoton[2] = 5.40E-10 ;
187   fTFphoton = new TFormula("ToF response to photon" , "gaus") ; 
188   fTFphoton->SetParameters( fTphoton[0], fTphoton[1], fTphoton[2]) ; 
189   // Electrons
190   fTelectron[0] = 2.73E-1 ; 
191   fTelectron[1] = 1.56E-8 ; 
192   fTelectron[2] = 5.70E-10 ;
193   fTFelectron = new TFormula("ToF response to electron" , "gaus") ; 
194   fTFelectron->SetParameters( fTelectron[0], fTelectron[1], fTelectron[2]) ; 
195   //Charged Hadrons
196   fTchargedhadron[0] = 1.58E-1 ; 
197   fTchargedhadron[1] = 1.64E-8 ; 
198   fTchargedhadron[2] = 3.59E-10 ;
199   fTFchargedhadron = new TFormula("ToF response to charged hadron" , "landau") ; 
200   fTFchargedhadron->SetParameters( fTchargedhadron[0], fTchargedhadron[1], fTchargedhadron[2]) ; 
201   //Neutral Hadrons
202   fTneutralhadron[0] = 9.62E-1 ; 
203   fTneutralhadron[1] = 1.65E-8 ; 
204   fTneutralhadron[2] = 6.46E-10 ;
205   fTFneutralhadron = new TFormula("ToF response to neutral hadron" , "landau") ; 
206   fTFneutralhadron->SetParameters( fTneutralhadron[0], fTneutralhadron[1], fTneutralhadron[2]) ;   
207 }
208
209 //________________________________________________________________________
210 void  AliPHOSPIDv1::Exec(Option_t *option)
211 {
212   // Steering method to perform particle reconstruction and identification
213   // for the event range from fFirstEvent to fLastEvent.
214   // This range is optionally set by SetEventRange().
215   // if fLastEvent=-1 (by default), then process events until the end.
216   
217   if(strstr(option,"tim"))
218     gBenchmark->Start("PHOSPID");
219   
220   if(strstr(option,"print")) {
221     Print() ; 
222     return ; 
223   }
224
225
226   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
227  
228   if (fLastEvent == -1) 
229     fLastEvent = gime->MaxEvent() - 1 ;
230   else 
231     fLastEvent = TMath::Min(fLastEvent,gime->MaxEvent());
232   Int_t nEvents   = fLastEvent - fFirstEvent + 1;
233
234   Int_t ievent ; 
235   for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
236     gime->Event(ievent,"TR") ;
237     if(gime->TrackSegments() && //Skip events, where no track segments made
238        gime->TrackSegments()->GetEntriesFast()) {
239       MakeRecParticles() ;
240       MakePID() ; 
241       WriteRecParticles();
242       if(strstr(option,"deb"))
243         PrintRecParticles(option) ;
244       //increment the total number of rec particles per run 
245       fRecParticlesInRun += gime->RecParticles()->GetEntriesFast() ; 
246     }
247   }
248   if(strstr(option,"deb"))
249       PrintRecParticles(option);
250   if(strstr(option,"tim")){
251     gBenchmark->Stop("PHOSPID");
252     Info("Exec", "took %f seconds for PID %f seconds per event", 
253          gBenchmark->GetCpuTime("PHOSPID"),  
254          gBenchmark->GetCpuTime("PHOSPID")/nEvents) ;
255   }
256   if(fWrite)
257     Unload();
258 }
259
260 //____________________________________________________________________________
261 const TString AliPHOSPIDv1::GetFileNamePrincipal(TString particle) const
262 {
263   //Get file name that contains the PCA for a particle ("photon or pi0")
264   particle.ToLower();
265   TString name;
266   if      (particle=="photon") name = fFileNamePrincipalPhoton ;
267   else if (particle=="pi0"   ) name = fFileNamePrincipalPi0    ;
268   else    Error("GetFileNamePrincipal","Wrong particle name: %s (choose from pi0/photon)\n",particle.Data());
269   return name;
270 }
271
272 //____________________________________________________________________________
273 Float_t  AliPHOSPIDv1::GetParameterCalibration(Int_t i) const 
274 {
275   // Get the i-th parameter "Calibration"
276   Float_t param = 0.;
277   if (i>2 || i<0)
278     Error("GetParameterCalibration","Invalid parameter number: %d",i);
279   else
280     param = (*fParameters)(0,i);
281   return param;
282 }
283
284 //____________________________________________________________________________
285 Float_t  AliPHOSPIDv1::GetCalibratedEnergy(Float_t e) const
286 {
287 //      It calibrates Energy depending on the recpoint energy.
288 //      The energy of the reconstructed cluster is corrected with 
289 //      the formula A + B* E  + C* E^2, whose parameters where obtained 
290 //      through the study of the reconstructed energy distribution of 
291 //      monoenergetic photons.
292  
293   Float_t p[]={0.,0.,0.};
294   for (Int_t i=0; i<3; i++) p[i] = GetParameterCalibration(i);
295   Float_t enerec = p[0] +  p[1]*e + p[2]*e*e;
296   return enerec ;
297
298 }
299
300 //____________________________________________________________________________
301 Float_t  AliPHOSPIDv1::GetParameterCpv2Emc(Int_t i, TString axis) const 
302 {
303   // Get the i-th parameter "CPV-EMC distance" for the specified axis
304   Float_t param = 0.;
305   if(i>2 || i<0)
306     Error("GetParameterCpv2Emc","Invalid parameter number: %d",i);
307   else {
308     axis.ToLower();
309     if      (axis == "x") param = (*fParameters)(1,i);
310     else if (axis == "z") param = (*fParameters)(2,i);
311     else Error("GetParameterCpv2Emc","Invalid axis name: %s",axis.Data());
312   }
313   return  param;
314 }
315
316 //____________________________________________________________________________
317 Float_t  AliPHOSPIDv1::GetCpv2EmcDistanceCut(TString axis, Float_t e) const
318 {
319   // Get CpvtoEmcDistance Cut depending on the cluster energy, axis and 
320   // Purity-Efficiency point 
321
322   axis.ToLower();
323   Float_t p[]={0.,0.,0.};
324   for (Int_t i=0; i<3; i++) p[i] = GetParameterCpv2Emc(i,axis);
325   Float_t sig = p[0] + TMath::Exp(p[1] - p[2]*e);
326   return sig;
327 }
328
329 //____________________________________________________________________________
330 Float_t  AliPHOSPIDv1::GetEllipseParameter(TString particle, TString param, Float_t e) const 
331 {
332   // Calculates the parameter param of the ellipse
333
334   particle.ToLower();
335   param.   ToLower();
336   Float_t p[4]={0.,0.,0.,0.};
337   Float_t value = 0.0;
338   for (Int_t i=0; i<4; i++) p[i] = GetParameterToCalculateEllipse(particle,param,i);
339   if (particle == "photon") {
340     if      (param.Contains("a"))  e = TMath::Min((Double_t)e,70.);
341     else if (param.Contains("b"))  e = TMath::Min((Double_t)e,70.);
342     else if (param.Contains("x0")) e = TMath::Max((Double_t)e,1.1);
343   }
344
345  if (particle == "photon")
346     value = p[0]/TMath::Sqrt(e) + p[1]*e + p[2]*e*e + p[3];
347   else if (particle == "pi0")
348     value = p[0] + p[1]*e + p[2]*e*e;
349
350   return value;
351 }
352
353 //_____________________________________________________________________________
354 Float_t  AliPHOSPIDv1::GetParameterPhotonBoundary (Int_t i) const
355
356   // Get the parameter "i" to calculate the boundary on the moment M2x
357   // for photons at high p_T
358   Float_t param = 0;
359   if (i>3 || i<0)
360     Error("GetParameterPhotonBoundary","Wrong parameter number: %d\n",i);
361   else
362     param = (*fParameters)(14,i) ;
363   return param;
364 }
365
366 //____________________________________________________________________________
367 Float_t  AliPHOSPIDv1::GetParameterPi0Boundary (Int_t i) const
368
369   // Get the parameter "i" to calculate the boundary on the moment M2x
370   // for pi0 at high p_T
371   Float_t param = 0;
372   if (i>2 || i<0)
373     Error("GetParameterPi0Boundary","Wrong parameter number: %d\n",i);
374   else
375     param = (*fParameters)(15,i) ;
376   return param;
377 }
378
379 //____________________________________________________________________________
380 Float_t  AliPHOSPIDv1::GetParameterTimeGate(Int_t i) const
381 {
382   // Get TimeGate parameter depending on Purity-Efficiency i:
383   // i=0 - Low purity, i=1 - Medium purity, i=2 - High purity
384   Float_t param = 0.;
385   if(i>2 || i<0)
386     Error("GetParameterTimeGate","Invalid Efficiency-Purity choice %d",i);
387   else
388     param = (*fParameters)(3,i) ; 
389   return param;
390 }
391
392 //_____________________________________________________________________________
393 Float_t  AliPHOSPIDv1::GetParameterToCalculateEllipse(TString particle, TString param, Int_t i) const
394
395   // Get the parameter "i" that is needed to calculate the ellipse 
396   // parameter "param" for the particle "particle" ("photon" or "pi0")
397
398   particle.ToLower();
399   param.   ToLower();
400   Int_t offset = -1;
401   if      (particle == "photon") offset=0;
402   else if (particle == "pi0")    offset=5;
403   else
404     Error("GetParameterToCalculateEllipse","Wrong particle name: %s (choose from pi0/photon)\n",particle.Data());
405
406   Int_t p= -1;
407   Float_t par = 0;
408
409   if     (param.Contains("a")) p=4+offset; 
410   else if(param.Contains("b")) p=5+offset; 
411   else if(param.Contains("c")) p=6+offset; 
412   else if(param.Contains("x0"))p=7+offset; 
413   else if(param.Contains("y0"))p=8+offset;
414
415   if      (i>4 || i<0)
416     Error("GetParameterToCalculateEllipse", "No parameter with index", i) ; 
417   else if (p==-1)
418     Error("GetParameterToCalculateEllipse", "No parameter with name %s", param.Data() ) ; 
419   else
420     par = (*fParameters)(p,i) ;
421   
422   return par;
423 }
424
425
426 //____________________________________________________________________________
427 Float_t  AliPHOSPIDv1::GetDistance(AliPHOSEmcRecPoint * emc,AliPHOSCpvRecPoint * cpv, Option_t *  axis)const
428 {
429   // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
430   
431   const AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ; 
432   TVector3 vecEmc ;
433   TVector3 vecCpv ;
434   if(cpv){
435     emc->GetLocalPosition(vecEmc) ;
436     cpv->GetLocalPosition(vecCpv) ; 
437     
438     if(emc->GetPHOSMod() == cpv->GetPHOSMod()){      
439       // Correct to difference in CPV and EMC position due to different distance to center.
440       // we assume, that particle moves from center
441       Float_t dCPV = geom->GetIPtoOuterCoverDistance();
442       Float_t dEMC = geom->GetIPtoCrystalSurface() ;
443       dEMC         = dEMC / dCPV ;
444       vecCpv = dEMC * vecCpv  - vecEmc ; 
445       if (axis == "X") return vecCpv.X();
446       if (axis == "Y") return vecCpv.Y();
447       if (axis == "Z") return vecCpv.Z();
448       if (axis == "R") return vecCpv.Mag();
449     }
450     return 100000000 ;
451   }
452   return 100000000 ;
453 }
454 //____________________________________________________________________________
455 Int_t  AliPHOSPIDv1::GetCPVBit(AliPHOSEmcRecPoint * emc,AliPHOSCpvRecPoint * cpv, Int_t effPur, Float_t e) const
456 {
457   if(effPur>2 || effPur<0)
458     Error("GetCPVBit","Invalid Efficiency-Purity choice %d",effPur);
459   
460   Float_t sigX = GetCpv2EmcDistanceCut("X",e);
461   Float_t sigZ = GetCpv2EmcDistanceCut("Z",e);
462   
463   Float_t deltaX = TMath::Abs(GetDistance(emc, cpv,  "X"));
464   Float_t deltaZ = TMath::Abs(GetDistance(emc, cpv,  "Z"));
465
466   if((deltaX>sigX*(effPur+1))|(deltaZ>sigZ*(effPur+1)))
467     return 1;//Neutral
468   else
469     return 0;//Charged
470 }
471
472 //____________________________________________________________________________
473 Int_t  AliPHOSPIDv1::GetPrincipalBit(TString particle, const Double_t* p, Int_t effPur, Float_t e)const
474 {
475   //Is the particle inside de PCA ellipse?
476   
477   particle.ToLower();
478   Int_t    prinbit  = 0 ;
479   Float_t a        = GetEllipseParameter(particle,"a" , e); 
480   Float_t b        = GetEllipseParameter(particle,"b" , e);
481   Float_t c        = GetEllipseParameter(particle,"c" , e);
482   Float_t x0 = GetEllipseParameter(particle,"x0", e); 
483   Float_t y0 = GetEllipseParameter(particle,"y0", e);
484   
485   Float_t r = TMath::Power((p[0] - x0)/a,2) + 
486               TMath::Power((p[1] - y0)/b,2) +
487             c*(p[0] -  x0)*(p[1] - y0)/(a*b) ;
488   //3 different ellipses defined
489   if((effPur==2) && (r<1./2.)) prinbit= 1;
490   if((effPur==1) && (r<2.   )) prinbit= 1;
491   if((effPur==0) && (r<9./2.)) prinbit= 1;
492
493   if(r<0)
494     Error("GetPrincipalBit", "Negative square?") ;
495
496   return prinbit;
497
498 }
499 //____________________________________________________________________________
500 Int_t  AliPHOSPIDv1::GetHardPhotonBit(AliPHOSEmcRecPoint * emc) const
501 {
502   // Set bit for identified hard photons (E > 30 GeV)
503   // if the second moment M2x is below the boundary
504
505   Float_t e   = emc->GetEnergy();
506   if (e < 30.0) return 0;
507   Float_t m2x = emc->GetM2x();
508   Float_t m2xBoundary = GetParameterPhotonBoundary(0) *
509     TMath::Exp(-TMath::Power(e-GetParameterPhotonBoundary(1),2)/2.0/
510                 TMath::Power(GetParameterPhotonBoundary(2),2)) +
511     GetParameterPhotonBoundary(3);
512   //Info("GetHardPhotonBit","E=%f, m2x=%f, boundary=%f",e,m2x,m2xBoundary);
513   if (m2x < m2xBoundary)
514     return 1;// A hard photon
515   else
516     return 0;// Not a hard photon
517 }
518
519 //____________________________________________________________________________
520 Int_t  AliPHOSPIDv1::GetHardPi0Bit(AliPHOSEmcRecPoint * emc) const
521 {
522   // Set bit for identified hard pi0  (E > 30 GeV)
523   // if the second moment M2x is above the boundary
524
525   Float_t e   = emc->GetEnergy();
526   if (e < 30.0) return 0;
527   Float_t m2x = emc->GetM2x();
528   Float_t m2xBoundary = GetParameterPi0Boundary(0) +
529                     e * GetParameterPi0Boundary(1);
530   //Info("GetHardPi0Bit","E=%f, m2x=%f, boundary=%f",e,m2x,m2xBoundary);
531   if (m2x > m2xBoundary)
532     return 1;// A hard pi0
533   else
534     return 0;// Not a hard pi0
535 }
536
537 //____________________________________________________________________________
538 TVector3 AliPHOSPIDv1::GetMomentumDirection(AliPHOSEmcRecPoint * emc, AliPHOSCpvRecPoint * )const 
539
540   // Calculates the momentum direction:
541   //   1. if only a EMC RecPoint, direction is given by IP and this RecPoint
542   //   2. if a EMC RecPoint and CPV RecPoint, direction is given by the line through the 2 recpoints 
543   //  However because of the poor position resolution of PPSD the direction is always taken as if we were 
544   //  in case 1.
545
546   TVector3 dir(0,0,0) ; 
547   
548   TVector3 emcglobalpos ;
549   TMatrix  dummy ;
550   
551   emc->GetGlobalPosition(emcglobalpos, dummy) ;
552   
553
554   dir = emcglobalpos ;  
555   dir.SetMag(1.) ;
556
557   //account correction to the position of IP
558   Float_t xo,yo,zo ; //Coordinates of the origin
559   if(gAlice && gAlice->GetMCApp() && gAlice->Generator())
560     gAlice->Generator()->GetOrigin(xo,yo,zo) ;
561   else{
562     xo=yo=zo=0.;
563   }
564   TVector3 origin(xo,yo,zo);
565   dir = dir - origin ;
566
567   return dir ;  
568 }
569
570 //____________________________________________________________________________
571 void  AliPHOSPIDv1::MakePID()
572 {
573   // construct the PID weight from a Bayesian Method
574
575   Int_t index ;
576   const Int_t kSPECIES = AliESDtrack::kSPECIESN ;
577   Double_t pid[kSPECIES] = {0., 0., 0., 0., 0., 0.} ;  
578   Int_t nparticles = AliPHOSGetter::Instance()->RecParticles()->GetEntriesFast() ;
579 //   const Int_t kMAXPARTICLES = 2000 ; 
580 //   if (nparticles >= kMAXPARTICLES) 
581 //     Error("MakePID", "Change size of MAXPARTICLES") ; 
582 //   Double_t stof[kSPECIES][kMAXPARTICLES] ;
583   Double_t * stof[kSPECIES];
584   for (Int_t i =0; i< kSPECIES; i++)
585     stof[i] = new Double_t[nparticles];
586
587   // make the normalized distribution of pid for this event 
588   // w(pid) in the Bayesian formulation
589   for(index = 0 ; index < nparticles ; index ++) {
590     AliPHOSRecParticle * recpar = AliPHOSGetter::Instance()->RecParticle(index) ;  
591     
592     pid[AliESDtrack::kKaon0] = 0  ; 
593     
594     if (recpar->IsPhoton() || recpar->IsHardPhoton())  
595       pid[AliESDtrack::kPhoton]++ ; 
596     
597     else if (recpar->IsPi0() || recpar->IsHardPi0())   
598       pid[AliESDtrack::kPi0]++ ; 
599     
600     else if (recpar->IsElectron()) {
601       pid[AliESDtrack::kElectron]++ ; 
602       pid[AliESDtrack::kMuon]++ ; 
603     }
604     
605     else if (recpar->IsChargedHadron()){
606       pid[AliESDtrack::kPion]++ ; 
607       pid[AliESDtrack::kKaon]++ ; 
608       pid[AliESDtrack::kProton]++ ;
609     } 
610     
611     else if (recpar->IsNeutralHadron())
612       pid[AliESDtrack::kNeutron]++ ; 
613     
614     else if (recpar->IsEleCon())   
615       pid[AliESDtrack::kEleCon]++ ;
616     
617     // now get the signals probability
618     // s(pid) in the Bayesian formulation
619     // Tof
620
621     stof[AliESDtrack::kPhoton][index]   = fTFphoton->Eval(recpar->ToF()) ; 
622     stof[AliESDtrack::kPi0][index]      = fTFphoton->Eval(recpar->ToF()) ; // pi0 are detected via decay photon 
623     stof[AliESDtrack::kElectron][index] = fTFelectron->Eval(recpar->ToF()) ; 
624     stof[AliESDtrack::kPion][index]     = fTFchargedhadron->Eval(recpar->ToF()) ; 
625     stof[AliESDtrack::kKaon][index]     = fTFchargedhadron->Eval(recpar->ToF()) ; 
626     stof[AliESDtrack::kProton][index]   = fTFchargedhadron->Eval(recpar->ToF()) ; 
627     stof[AliESDtrack::kNeutron][index]  = fTFneutralhadron->Eval(recpar->ToF()) ; 
628     stof[AliESDtrack::kEleCon][index]   = fTFphoton->Eval(recpar->ToF()) ; // a conversion electron has the photon ToF
629     stof[AliESDtrack::kKaon0][index]    = 0 ; // do not know yet what to to with K0
630     stof[AliESDtrack::kMuon][index]     = 0 ; // do not know yet what to do with muon
631     
632   }
633   for (index = 0 ; index < kSPECIES ; index++) 
634     pid[index] /= nparticles ; 
635  
636
637   for(index = 0 ; index < nparticles ; index ++) {
638     // calculates the Bayesian weight
639     Int_t jndex ;
640     Double_t wn = 0.0 ; 
641     for (jndex = 0 ; jndex < kSPECIES ; jndex++) 
642       wn += stof[jndex][index] * pid[jndex] ;
643     AliPHOSRecParticle * recpar = AliPHOSGetter::Instance()->RecParticle(index) ;  
644     if (TMath::Abs(wn)>0)
645       for (jndex = 0 ; jndex < kSPECIES ; jndex++) {
646         recpar->SetPID(jndex, stof[jndex][index] * pid[jndex] / wn) ; 
647       }
648   }
649   for (Int_t i =0; i< kSPECIES; i++)
650     delete [] stof[i];
651 }
652
653 //____________________________________________________________________________
654 void  AliPHOSPIDv1::MakeRecParticles()
655 {
656   // Makes a RecParticle out of a TrackSegment
657   
658   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
659   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
660   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
661   TClonesArray * trackSegments = gime->TrackSegments() ; 
662   if ( !emcRecPoints || !cpvRecPoints || !trackSegments ) {
663     Fatal("MakeRecParticles", "RecPoints or TrackSegments not found !") ;  
664   }
665   TClonesArray * recParticles  = gime->RecParticles() ; 
666   recParticles->Clear();
667
668   TIter next(trackSegments) ; 
669   AliPHOSTrackSegment * ts ; 
670   Int_t index = 0 ; 
671   AliPHOSRecParticle * rp ; 
672   while ( (ts = (AliPHOSTrackSegment *)next()) ) {
673     
674     new( (*recParticles)[index] ) AliPHOSRecParticle() ;
675     rp = (AliPHOSRecParticle *)recParticles->At(index) ; 
676     rp->SetTrackSegment(index) ;
677     rp->SetIndexInList(index) ;
678         
679     AliPHOSEmcRecPoint * emc = 0 ;
680     if(ts->GetEmcIndex()>=0)
681       emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ;
682     
683     AliPHOSCpvRecPoint * cpv = 0 ;
684     if(ts->GetCpvIndex()>=0)
685       cpv = (AliPHOSCpvRecPoint *) cpvRecPoints->At(ts->GetCpvIndex()) ;
686     
687     Int_t track = 0 ; 
688     track = ts->GetTrackIndex() ; 
689       
690     // Now set type (reconstructed) of the particle
691
692     // Choose the cluster energy range
693     
694     if (!emc) {
695       Fatal("MakeRecParticles", "-> emc(%d) = %d", ts->GetEmcIndex(), emc ) ;
696     }
697
698     Float_t e = emc->GetEnergy() ;   
699     
700     Float_t  lambda[2] ;
701     emc->GetElipsAxis(lambda) ;
702     
703     if((lambda[0]>0.01) && (lambda[1]>0.01)){
704       // Looking PCA. Define and calculate the data (X),
705       // introduce in the function X2P that gives the components (P).  
706
707       Float_t  Spher = 0. ;
708       Float_t  Emaxdtotal = 0. ; 
709       
710       if((lambda[0]+lambda[1])!=0) 
711         Spher=fabs(lambda[0]-lambda[1])/(lambda[0]+lambda[1]); 
712       
713       Emaxdtotal=emc->GetMaximalEnergy()/emc->GetEnergy(); 
714       
715       fX[0] = lambda[0] ;  
716       fX[1] = lambda[1] ; 
717       fX[2] = emc->GetDispersion() ; 
718       fX[3] = Spher ; 
719       fX[4] = emc->GetMultiplicity() ;  
720       fX[5] = Emaxdtotal ;  
721       fX[6] = emc->GetCoreEnergy() ;  
722       
723       fPrincipalPhoton->X2P(fX,fPPhoton);
724       fPrincipalPi0   ->X2P(fX,fPPi0);
725
726     }
727     else{
728       fPPhoton[0]=-100.0;  //We do not accept clusters with 
729       fPPhoton[1]=-100.0;  //one cell as a photon-like
730       fPPi0[0]   =-100.0;
731       fPPi0[1]   =-100.0;
732     }
733     
734     Float_t time = emc->GetTime() ;
735     rp->SetTof(time) ; 
736     
737     // Loop of Efficiency-Purity (the 3 points of purity or efficiency 
738     // are taken into account to set the particle identification)
739     for(Int_t effPur = 0; effPur < 3 ; effPur++){
740       
741       // Looking at the CPV detector. If RCPV greater than CpvEmcDistance, 
742       // 1st,2nd or 3rd bit (depending on the efficiency-purity point )
743       // is set to 1
744       if(GetCPVBit(emc, cpv, effPur,e) == 1 )  
745         rp->SetPIDBit(effPur) ;
746       
747       // Looking the TOF. If TOF smaller than gate,  4th, 5th or 6th 
748       // bit (depending on the efficiency-purity point )is set to 1             
749       if(time< (*fParameters)(3,effPur)) 
750         rp->SetPIDBit(effPur+3) ;                   
751   
752       //Photon PCA
753       //If we are inside the ellipse, 7th, 8th or 9th 
754       // bit (depending on the efficiency-purity point )is set to 1 
755       if(GetPrincipalBit("photon",fPPhoton,effPur,e) == 1) 
756         rp->SetPIDBit(effPur+6) ;
757
758       //Pi0 PCA
759       //If we are inside the ellipse, 10th, 11th or 12th 
760       // bit (depending on the efficiency-purity point )is set to 1 
761       if(GetPrincipalBit("pi0"   ,fPPi0   ,effPur,e) == 1) 
762         rp->SetPIDBit(effPur+9) ;
763     }
764     if(GetHardPhotonBit(emc))
765       rp->SetPIDBit(12) ;
766     if(GetHardPi0Bit   (emc))
767       rp->SetPIDBit(13) ;
768     
769     if(track >= 0) 
770       rp->SetPIDBit(14) ; 
771
772     //Set momentum, energy and other parameters 
773     Float_t  encal = GetCalibratedEnergy(e);
774     TVector3 dir   = GetMomentumDirection(emc,cpv) ; 
775     dir.SetMag(encal) ;
776     rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),encal) ;
777     rp->SetCalcMass(0);
778     rp->Name(); //If photon sets the particle pdg name to gamma
779     rp->SetProductionVertex(0,0,0,0);
780     rp->SetFirstMother(-1);
781     rp->SetLastMother(-1);
782     rp->SetFirstDaughter(-1);
783     rp->SetLastDaughter(-1);
784     rp->SetPolarisation(0,0,0);
785     //Set the position in global coordinate system from the RecPoint
786     AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
787     AliPHOSTrackSegment * ts  = gime->TrackSegment(rp->GetPHOSTSIndex()) ; 
788     AliPHOSEmcRecPoint  * erp = gime->EmcRecPoint(ts->GetEmcIndex()) ; 
789     TVector3 pos ; 
790     geom->GetGlobal(erp, pos) ; 
791     rp->SetPos(pos);
792     index++ ; 
793   }
794 }
795   
796 //____________________________________________________________________________
797 void  AliPHOSPIDv1::Print() const
798 {
799   // Print the parameters used for the particle type identification
800
801     Info("Print", "=============== AliPHOSPIDv1 ================") ;
802     printf("Making PID\n") ;
803     printf("    Pricipal analysis file from 0.5 to 100 %s\n", fFileNamePrincipalPhoton.Data() )   ; 
804     printf("    Name of parameters file     %s\n", fFileNameParameters.Data() )  ;
805     printf("    Matrix of Parameters: 14x4\n") ;
806     printf("        Energy Calibration  1x3 [3 parametres to calibrate energy: A + B* E + C * E^2]\n") ;
807     printf("        RCPV 2x3 rows x and z, columns function cut parameters\n") ;
808     printf("        TOF  1x3 [High Eff-Low Pur,Medium Eff-Pur, Low Eff-High Pur]\n") ;
809     printf("        PCA  5x4 [5 ellipse parametres and 4 parametres to calculate them: A/Sqrt(E) + B* E + C * E^2 + D]\n") ;
810     Printf("    Pi0 PCA  5x3 [5 ellipse parametres and 3 parametres to calculate them: A + B* E + C * E^2]\n") ;
811     fParameters->Print() ;
812 }
813
814
815
816 //____________________________________________________________________________
817 void AliPHOSPIDv1::PrintRecParticles(Option_t * option)
818 {
819   // Print table of reconstructed particles
820
821   AliPHOSGetter *gime = AliPHOSGetter::Instance() ; 
822
823   TClonesArray * recParticles = gime->RecParticles() ; 
824
825   TString message ; 
826   message  = "\nevent " ;
827   message += gAlice->GetEvNumber() ; 
828   message += "       found " ; 
829   message += recParticles->GetEntriesFast(); 
830   message += " RecParticles\n" ; 
831
832   if(strstr(option,"all")) {  // printing found TS
833     message += "\n  PARTICLE         Index    \n" ; 
834     
835     Int_t index ;
836     for (index = 0 ; index < recParticles->GetEntries() ; index++) {
837       AliPHOSRecParticle * rp = (AliPHOSRecParticle * ) recParticles->At(index) ;       
838       message += "\n" ;
839       message += rp->Name().Data() ;  
840       message += " " ;
841       message += rp->GetIndexInList() ;  
842       message += " " ;
843       message += rp->GetType()  ;
844     }
845   }
846   Info("Print", message.Data() ) ; 
847 }
848
849 //____________________________________________________________________________
850 void  AliPHOSPIDv1::SetParameters() 
851 {
852   // PCA : To do the Principal Components Analysis it is necessary 
853   // the Principal file, which is opened here
854   fX       = new double[7]; // Data for the PCA 
855   fPPhoton = new double[7]; // Eigenvalues of the PCA
856   fPPi0    = new double[7]; // Eigenvalues of the Pi0 PCA
857
858   // Read photon principals from the photon file
859   
860   fFileNamePrincipalPhoton = "$ALICE_ROOT/PHOS/PCA8pa15_0.5-100.root" ; 
861   TFile f( fFileNamePrincipalPhoton.Data(), "read" ) ;
862   fPrincipalPhoton = dynamic_cast<TPrincipal*> (f.Get("principal")) ; 
863   f.Close() ; 
864
865   // Read pi0 principals from the pi0 file
866
867   fFileNamePrincipalPi0    = "$ALICE_ROOT/PHOS/PCA_pi0_40-120.root" ;
868   TFile fPi0( fFileNamePrincipalPi0.Data(), "read" ) ;
869   fPrincipalPi0    = dynamic_cast<TPrincipal*> (fPi0.Get("principal")) ; 
870   fPi0.Close() ;
871
872   // Open parameters file and initialization of the Parameters matrix. 
873   // In the File Parameters.dat are all the parameters. These are introduced 
874   // in a matrix of 16x4  
875   // 
876   // All the parameters defined in this file are, in order of row: 
877   // line   0   : calibration 
878   // lines  1,2 : CPV rectangular cat for X and Z
879   // line   3   : TOF cut
880   // lines  4-8 : parameters to calculate photon PCA ellipse
881   // lines  9-13: parameters to calculate pi0 PCA ellipse
882   // lines 14-15: parameters to calculate border for high-pt photons and pi0
883
884   fFileNameParameters = gSystem->ExpandPathName("$ALICE_ROOT/PHOS/Parameters.dat");
885   fParameters = new TMatrix(16,4) ;
886   const Int_t maxLeng=255;
887   char string[maxLeng];
888
889   // Open a text file with PID parameters
890   FILE *fd = fopen(fFileNameParameters.Data(),"r");
891   if (!fd)
892     Fatal("SetParameter","File %s with a PID parameters cannot be opened\n",
893           fFileNameParameters.Data());
894
895   Int_t i=0;
896   // Read parameter file line-by-line and skip empty line and comments
897   while (fgets(string,maxLeng,fd) != NULL) {
898     if (string[0] == '\n' ) continue;
899     if (string[0] == '!'  ) continue;
900     sscanf(string, "%f %f %f %f",
901            &(*fParameters)(i,0), &(*fParameters)(i,1), 
902            &(*fParameters)(i,2), &(*fParameters)(i,3));
903     i++;
904     //Info("SetParameters", "line %d: %s",i,string);
905   }
906   fclose(fd);
907 }
908
909 //____________________________________________________________________________
910 void  AliPHOSPIDv1::SetParameterCalibration(Int_t i,Float_t param) 
911 {
912   // Set parameter "Calibration" i to a value param
913   if(i>2 || i<0)
914     Error("SetParameterCalibration","Invalid parameter number: %d",i);
915   else
916     (*fParameters)(0,i) = param ;
917 }
918
919 //____________________________________________________________________________
920 void  AliPHOSPIDv1::SetParameterCpv2Emc(Int_t i, TString axis, Float_t cut) 
921 {
922   // Set the parameters to calculate Cpv-to-Emc Distance Cut depending on 
923   // Purity-Efficiency point i
924
925   if(i>2 || i<0)
926     Error("SetParameterCpv2Emc","Invalid parameter number: %d",i);
927   else {
928     axis.ToLower();
929     if      (axis == "x") (*fParameters)(1,i) = cut;
930     else if (axis == "z") (*fParameters)(2,i) = cut;
931     else Error("SetParameterCpv2Emc","Invalid axis name: %s",axis.Data());
932   }
933 }
934
935 //____________________________________________________________________________
936 void  AliPHOSPIDv1::SetParameterPhotonBoundary(Int_t i,Float_t param) 
937 {
938   // Set parameter "Hard photon boundary" i to a value param
939   if(i>4 || i<0)
940     Error("SetParameterPhotonBoundary","Invalid parameter number: %d",i);
941   else
942     (*fParameters)(14,i) = param ;
943 }
944
945 //____________________________________________________________________________
946 void  AliPHOSPIDv1::SetParameterPi0Boundary(Int_t i,Float_t param) 
947 {
948   // Set parameter "Hard pi0 boundary" i to a value param
949   if(i>1 || i<0)
950     Error("SetParameterPi0Boundary","Invalid parameter number: %d",i);
951   else
952     (*fParameters)(15,i) = param ;
953 }
954
955 //_____________________________________________________________________________
956 void  AliPHOSPIDv1::SetParameterTimeGate(Int_t i, Float_t gate) 
957 {
958   // Set the parameter TimeGate depending on Purity-Efficiency point i 
959   if (i>2 || i<0)
960     Error("SetParameterTimeGate","Invalid Efficiency-Purity choice %d",i);
961   else
962     (*fParameters)(3,i)= gate ; 
963
964
965 //_____________________________________________________________________________
966 void  AliPHOSPIDv1::SetParameterToCalculateEllipse(TString particle, TString param, Int_t i, Float_t par) 
967 {  
968   // Set the parameter "i" that is needed to calculate the ellipse 
969   // parameter "param" for a particle "particle"
970   
971   particle.ToLower();
972   param.   ToLower();
973   Int_t p= -1;
974   Int_t offset=0;
975
976   if      (particle == "photon") offset=0;
977   else if (particle == "pi0")    offset=5;
978   else
979     Error("SetParameterToCalculateEllipse","Wrong particle name: %s (choose from pi0/photon)\n",particle.Data());
980
981   if     (param.Contains("a")) p=4+offset; 
982   else if(param.Contains("b")) p=5+offset; 
983   else if(param.Contains("c")) p=6+offset; 
984   else if(param.Contains("x0"))p=7+offset; 
985   else if(param.Contains("y0"))p=8+offset;
986   if((i>4)||(i<0))
987     Error("SetEllipseParameter", "No parameter with index %d", i) ; 
988   else if(p==-1)
989     Error("SetEllipseParameter", "No parameter with name %s", param.Data() ) ; 
990   else
991     (*fParameters)(p,i) = par ;
992
993
994 //____________________________________________________________________________
995 void AliPHOSPIDv1::Unload() 
996 {
997   AliPHOSGetter * gime = AliPHOSGetter::Instance() ;  
998   gime->PhosLoader()->UnloadRecPoints() ;
999   gime->PhosLoader()->UnloadTracks() ;
1000   gime->PhosLoader()->UnloadRecParticles() ;
1001 }
1002
1003 //____________________________________________________________________________
1004 void  AliPHOSPIDv1::WriteRecParticles()
1005 {
1006  
1007   AliPHOSGetter *gime = AliPHOSGetter::Instance() ; 
1008
1009   TClonesArray * recParticles = gime->RecParticles() ; 
1010   recParticles->Expand(recParticles->GetEntriesFast() ) ;
1011   if(fWrite){
1012     TTree * treeP =  gime->TreeP();
1013     
1014     //First rp
1015     Int_t bufferSize = 32000 ;
1016     TBranch * rpBranch = treeP->Branch("PHOSRP",&recParticles,bufferSize);
1017     rpBranch->SetTitle(BranchName());
1018     
1019     rpBranch->Fill() ;
1020     
1021     gime->WriteRecParticles("OVERWRITE");
1022     gime->WritePID("OVERWRITE");
1023   }
1024 }
1025