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