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