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