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