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