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