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