]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/PartCorrBase/AliCaloPID.cxx
Centrality update (Alberica)
[u/mrichter/AliRoot.git] / PWG4 / PartCorrBase / AliCaloPID.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 /* $Id: AliCaloPID.cxx 21839 2007-10-29 13:49:42Z gustavo $ */
16
17 //_________________________________________________________________________
18 // Class for PID selection with calorimeters
19 // The Output of the 2 main methods GetIdentifiedParticleType is a PDG number identifying the cluster, 
20 // being kPhoton, kElectron, kPi0 ... as defined in the header file
21 //   - GetIdentifiedParticleType(const TString calo, const Double_t * pid, const Float_t energy)
22 //      Reads the PID weights array of the ESDs and depending on its magnitude identifies the particle
23 //   - GetIdentifiedParticleType(const TString calo,const TLorentzVector mom, const AliVCluster * cluster)
24 //      Recalcultes PID, the bayesian or any new one to be implemented in the future
25 //      Right now only the possibility to recalculate EMCAL with bayesian and simple PID.
26 //      In order to recalculate Bayesian, it is necessary to load the EMCALUtils library
27 //      and do SwitchOnBayesianRecalculation().
28 //      To change the PID parameters from Low to High like the ones by default, use the constructor 
29 //      AliCaloPID(flux)
30 //      where flux is AliCaloPID::kLow or AliCaloPID::kHigh
31 //      If it is necessary to change the parameters use the constructor 
32 //      AliCaloPID(AliEMCALPIDUtils *utils) and set the parameters before.
33 //   - SetPIDBits: Simple PID, depending on the thresholds fL0Cut fTOFCut and even the
34 //     result of the PID bayesian a different PID bit is set. 
35 //
36 //  All these methods can be called in the analysis you are interested.
37 //
38 //*-- Author: Gustavo Conesa (LNF-INFN) 
39 //////////////////////////////////////////////////////////////////////////////
40
41
42 // --- ROOT system ---
43 #include <TMath.h>
44 #include <TString.h>
45 #include <TList.h>
46
47 //---- ANALYSIS system ----
48 #include "AliCaloPID.h"
49 #include "AliVCluster.h"
50 #include "AliVTrack.h"
51 #include "AliAODPWG4Particle.h"
52 #include "AliCalorimeterUtils.h"
53
54 ClassImp(AliCaloPID)
55
56
57 //________________________________________________
58 AliCaloPID::AliCaloPID() : 
59 TObject(), 
60 fEMCALPhotonWeight(0.),   fEMCALPi0Weight(0.),  
61 fEMCALElectronWeight(0.), fEMCALChargeWeight(0.),      fEMCALNeutralWeight(0.),
62 fPHOSPhotonWeight(0.),    fPHOSPi0Weight(0.),  
63 fPHOSElectronWeight(0.),  fPHOSChargeWeight(0.) ,      fPHOSNeutralWeight(0.), 
64 fPHOSWeightFormula(0),    fPHOSPhotonWeightFormula(0), fPHOSPi0WeightFormula(0),
65 fPHOSPhotonWeightFormulaExpression(""), fPHOSPi0WeightFormulaExpression(""),
66 fL0CutMax(100.),          fL0CutMin(0),                fTOFCut(0.), 
67 fRcutPHOS(0),             fDebug(-1), 
68 fRecalculateBayesian(kFALSE), fParticleFlux(kLow),     fEMCALPIDUtils(),
69 fHistoNEBins(100),        fHistoEMax(100.),            fHistoEMin(0.),
70 fHistoNDEtaBins(100),     fHistoDEtaMax(0.15),         fHistoDEtaMin(-0.15),
71 fHistoNDPhiBins(100),     fHistoDPhiMax(0.15),         fHistoDPhiMin(-0.15),
72 fhTrackMatchedDEta(0x0),  fhTrackMatchedDPhi(0x0),     fhTrackMatchedDEtaDPhi(0x0)
73 {
74   //Ctor
75   
76   //Initialize parameters
77   InitParameters();
78 }
79
80 //________________________________________________
81 AliCaloPID::AliCaloPID(const Int_t flux) : 
82 TObject(), 
83 fEMCALPhotonWeight(0.),   fEMCALPi0Weight(0.),  
84 fEMCALElectronWeight(0.), fEMCALChargeWeight(0.),      fEMCALNeutralWeight(0.),
85 fPHOSPhotonWeight(0.),    fPHOSPi0Weight(0.),  
86 fPHOSElectronWeight(0.),  fPHOSChargeWeight(0.) ,      fPHOSNeutralWeight(0.), 
87 fPHOSWeightFormula(0),    fPHOSPhotonWeightFormula(0), fPHOSPi0WeightFormula(0),
88 fPHOSPhotonWeightFormulaExpression(""), fPHOSPi0WeightFormulaExpression(""),
89 fL0CutMax(100.),          fL0CutMin(0),                fTOFCut(0.), 
90 fRcutPHOS(0),             fDebug(-1), 
91 fRecalculateBayesian(kFALSE), fParticleFlux(flux),     fEMCALPIDUtils(),
92 fHistoNEBins(100),        fHistoEMax(100.),            fHistoEMin(0.),
93 fHistoNDEtaBins(100),     fHistoDEtaMax(0.15),         fHistoDEtaMin(-0.15),
94 fHistoNDPhiBins(100),     fHistoDPhiMax(0.15),         fHistoDPhiMin(-0.15),
95 fhTrackMatchedDEta(0x0),  fhTrackMatchedDPhi(0x0),     fhTrackMatchedDEtaDPhi(0x0)
96 {
97   //Ctor
98         
99   //Initialize parameters
100   InitParameters();
101
102 }
103
104 //________________________________________________
105 AliCaloPID::AliCaloPID(const TNamed * emcalpid) : 
106 TObject(), 
107 fEMCALPhotonWeight(0.),   fEMCALPi0Weight(0.),  
108 fEMCALElectronWeight(0.), fEMCALChargeWeight(0.),      fEMCALNeutralWeight(0.),
109 fPHOSPhotonWeight(0.),    fPHOSPi0Weight(0.),  
110 fPHOSElectronWeight(0.),  fPHOSChargeWeight(0.) ,      fPHOSNeutralWeight(0.), 
111 fPHOSWeightFormula(0),    fPHOSPhotonWeightFormula(0), fPHOSPi0WeightFormula(0),
112 fPHOSPhotonWeightFormulaExpression(""), fPHOSPi0WeightFormulaExpression(""),
113 fL0CutMax(100.),          fL0CutMin(0),                fTOFCut(0.), 
114 fRcutPHOS(0),             fDebug(-1), 
115 fRecalculateBayesian(kFALSE), fParticleFlux(kLow),     fEMCALPIDUtils((AliEMCALPIDUtils*) emcalpid),
116 fHistoNEBins(100),        fHistoEMax(100.),            fHistoEMin(0.),
117 fHistoNDEtaBins(100),     fHistoDEtaMax(0.15),         fHistoDEtaMin(-0.15),
118 fHistoNDPhiBins(100),     fHistoDPhiMax(0.15),         fHistoDPhiMin(-0.15),
119 fhTrackMatchedDEta(0x0),  fhTrackMatchedDPhi(0x0),     fhTrackMatchedDEtaDPhi(0x0)
120 {
121   //Ctor
122         
123   //Initialize parameters
124   InitParameters();
125 }
126
127 //_________________________________
128 AliCaloPID::~AliCaloPID() {
129   //Dtor
130   
131   delete fPHOSPhotonWeightFormula ;
132   delete fPHOSPi0WeightFormula ;
133   delete fEMCALPIDUtils ;
134
135 }
136
137 //________________________________________________________________________
138 TList *  AliCaloPID::GetCreateOutputObjects()
139 {  
140   // Create histograms to be saved in output file and 
141   // store them in outputContainer of the analysis class that calls this class.
142   
143   TList * outputContainer = new TList() ; 
144   outputContainer->SetName("CaloPIDHistos") ; 
145   
146   outputContainer->SetOwner(kFALSE);
147   
148   fhTrackMatchedDEta  = new TH2F
149   ("TrackMatchedDEta",
150    "d#eta of cluster-track vs cluster energy",
151    fHistoNEBins,fHistoEMin,fHistoEMax,fHistoNDEtaBins,fHistoDEtaMin,fHistoDEtaMax); 
152   fhTrackMatchedDEta->SetYTitle("d#eta");
153   fhTrackMatchedDEta->SetXTitle("E_{cluster} (GeV)");
154   
155   fhTrackMatchedDPhi  = new TH2F
156   ("TrackMatchedDPhi",
157    "d#phi of cluster-track vs cluster energy"
158    ,fHistoNEBins,fHistoEMin,fHistoEMax,fHistoNDPhiBins,fHistoDPhiMin,fHistoDPhiMax); 
159   fhTrackMatchedDPhi->SetYTitle("d#phi (rad)");
160   fhTrackMatchedDPhi->SetXTitle("E_{cluster} (GeV)");
161   
162   fhTrackMatchedDEtaDPhi  = new TH2F
163   ("TrackMatchedDEtaDPhi",
164    "d#eta vs d#phi of cluster-track vs cluster energy"
165    ,fHistoNDEtaBins,fHistoDEtaMin,fHistoDEtaMax,fHistoNDPhiBins,fHistoDPhiMin,fHistoDPhiMax); 
166   fhTrackMatchedDEtaDPhi->SetYTitle("d#phi (rad)");
167   fhTrackMatchedDEtaDPhi->SetXTitle("d#eta");   
168   
169   outputContainer->Add(fhTrackMatchedDEta) ; 
170   outputContainer->Add(fhTrackMatchedDPhi) ;
171   outputContainer->Add(fhTrackMatchedDEtaDPhi) ; 
172   
173   return outputContainer;
174 }
175
176
177
178 //_______________________________________________________________
179 void AliCaloPID::InitParameters()
180 {
181   //Initialize the parameters of the PID.
182   
183   fEMCALPhotonWeight   = 0.6 ;
184   fEMCALPi0Weight      = 0.6 ;
185   fEMCALElectronWeight = 0.6 ;
186   fEMCALChargeWeight   = 0.6 ;
187   fEMCALNeutralWeight  = 0.6 ;
188   
189   fPHOSPhotonWeight    = 0.6 ;
190   fPHOSPi0Weight       = 0.6 ;
191   fPHOSElectronWeight  = 0.6 ;
192   fPHOSChargeWeight    = 0.6 ;
193   fPHOSNeutralWeight   = 0.6 ;
194   
195   //Formula to set the PID weight threshold for photon or pi0
196   fPHOSWeightFormula       = kFALSE;
197   fPHOSPhotonWeightFormulaExpression = "0.98*(x<40)+ 0.68*(x>=100)+(x>=40 && x<100)*(0.98+x*(6e-3)-x*x*(2e-04)+x*x*x*(1.1e-06))";
198   fPHOSPi0WeightFormulaExpression    = "0.98*(x<65)+ 0.915*(x>=100)+(x>=65 && x-x*(1.95e-3)-x*x*(4.31e-05)+x*x*x*(3.61e-07))"   ;
199   
200   fL0CutMax = 0.3;
201   fL0CutMin = 0.01;
202   fTOFCut   = 5.e-9;
203   fRcutPHOS = 10.;
204   
205   fDebug    =-1;
206         
207   if(fRecalculateBayesian){
208         if(fParticleFlux == kLow){
209                 printf("AliCaloPID::Init() - SetLOWFluxParam\n");
210                 fEMCALPIDUtils->SetLowFluxParam() ;
211         }
212         else if (fParticleFlux == kHigh){
213                 printf("AliCaloPID::Init() - SetHIGHFluxParam\n");
214                 fEMCALPIDUtils->SetHighFluxParam() ;
215         }
216   }
217 }
218
219 //_______________________________________________________________
220 Int_t AliCaloPID::GetIdentifiedParticleType(const TString calo, const Double_t * pid, const Float_t energy) {
221   //Return most probable identity of the particle.
222   
223   if(!pid){ 
224     printf("AliCaloPID::GetIdentifiedParticleType() - pid pointer not initialized!!!\n");
225     abort();
226   }
227   
228   Float_t wPh  =  fPHOSPhotonWeight ;
229   Float_t wPi0 =  fPHOSPi0Weight ;
230   Float_t wE   =  fPHOSElectronWeight ;
231   Float_t wCh  =  fPHOSChargeWeight ;
232   Float_t wNe  =  fPHOSNeutralWeight ;
233     
234   if(calo == "PHOS" && fPHOSWeightFormula){
235     wPh  = GetPHOSPhotonWeightFormula()->Eval(energy) ;
236     wPi0 = GetPHOSPi0WeightFormula()   ->Eval(energy);
237   }
238   
239   if(calo == "EMCAL"){
240     
241     wPh  =  fEMCALPhotonWeight ;
242     wPi0 =  fEMCALPi0Weight ;
243     wE   =  fEMCALElectronWeight ;
244     wCh  =  fEMCALChargeWeight ;
245     wNe  =  fEMCALNeutralWeight ;
246     
247   }
248   
249   if(fDebug > 0)  printf("AliCaloPID::GetIdentifiedParticleType: calo %s, ph %0.2f, pi0 %0.2f, el %0.2f, conv el %0.2f, hadrons: pion %0.2f, kaon %0.2f, proton %0.2f , neutron %0.2f, kaon %0.2f \n",
250                          calo.Data(),pid[AliVCluster::kPhoton], pid[AliVCluster::kPi0],
251                          pid[AliVCluster::kElectron], pid[AliVCluster::kEleCon],
252                          pid[AliVCluster::kPion], pid[AliVCluster::kKaon], pid[AliVCluster::kProton],
253                          pid[AliVCluster::kNeutron], pid[AliVCluster::kKaon0]);
254   
255   Int_t pdg = kNeutralUnknown ;
256   Float_t chargedHadronWeight = pid[AliVCluster::kProton]+pid[AliVCluster::kKaon]+
257     pid[AliVCluster::kPion]+pid[AliVCluster::kMuon];
258   Float_t neutralHadronWeight = pid[AliVCluster::kNeutron]+pid[AliVCluster::kKaon0];
259   Float_t allChargedWeight    = pid[AliVCluster::kElectron]+pid[AliVCluster::kEleCon]+ chargedHadronWeight;
260   Float_t allNeutralWeight    = pid[AliVCluster::kPhoton]+pid[AliVCluster::kPi0]+ neutralHadronWeight;
261   
262   //Select most probable ID
263   if(calo=="PHOS"){
264     if(pid[AliVCluster::kPhoton] > wPh)        pdg = kPhoton ;
265     else if(pid[AliVCluster::kPi0] > wPi0)     pdg = kPi0 ; 
266     else if(pid[AliVCluster::kElectron] > wE)  pdg = kElectron ;
267     else if(pid[AliVCluster::kEleCon] >  wE)   pdg = kEleCon ;
268     else if(chargedHadronWeight > wCh)         pdg = kChargedHadron ;  
269     else if(neutralHadronWeight > wNe)         pdg = kNeutralHadron ; 
270     else if(allChargedWeight >  allNeutralWeight)
271       pdg = kChargedUnknown ; 
272     else 
273       pdg = kNeutralUnknown ;
274   }
275   else{//EMCAL
276     
277     if(pid[AliVCluster::kPhoton]  > wPh)                     pdg = kPhoton ;
278     else if(pid[AliVCluster::kElectron]  > wE)               pdg = kElectron ;
279     else if(pid[AliVCluster::kPhoton]+pid[AliVCluster::kElectron]  > wPh) pdg = kPhoton ; //temporal sollution until track matching for electrons is considered
280     else if(pid[AliVCluster::kPi0] > wPi0)                   pdg = kPi0 ; 
281     else if(chargedHadronWeight + neutralHadronWeight > wCh) pdg = kChargedHadron ;  
282     else if(neutralHadronWeight + chargedHadronWeight > wNe) pdg = kNeutralHadron ; 
283     else                                                     pdg = kNeutralUnknown ;
284   }
285   
286   if(fDebug > 0)printf("AliCaloPID::GetIdentifiedParticleType:Final Pdg: %d, cluster energy %2.2f \n", pdg,energy);
287    //printf("PDG1\n");
288   return pdg ;
289   
290 }
291
292 //_______________________________________________________________
293 Int_t AliCaloPID::GetIdentifiedParticleType(const TString calo,const TLorentzVector mom, const AliVCluster * cluster) {
294   //Recalculated PID with all parameters
295   
296   Float_t lambda0 = cluster->GetM02();
297   Float_t lambda1 = cluster->GetM20();
298   Float_t energy  = mom.E();    
299   
300   if(fDebug > 0) printf("AliCaloPID::GetIdentifiedParticleType: Calorimeter %s, E %3.2f, l0 %3.2f, l1 %3.2f, disp %3.2f, tof %1.11f, distCPV %3.2f, distToBC %1.1f, NMax %d\n",
301                         calo.Data(),energy,lambda0,cluster->GetM20(),cluster->GetDispersion(),cluster->GetTOF(), 
302                         cluster->GetEmcCpvDistance(), cluster->GetDistanceToBadChannel(),cluster->GetNExMax());
303   
304   if(calo == "EMCAL") {
305           //Recalculate Bayesian
306           if(fRecalculateBayesian){       
307                   if(fDebug > 0)  {
308                           const Double_t  *pid0 = cluster->GetPID();
309                           printf("AliCaloPID::GetIdentifiedParticleType: BEFORE calo %s, ph %0.2f, pi0 %0.2f, el %0.2f, conv el %0.2f, hadrons: pion %0.2f, kaon %0.2f, proton %0.2f , neutron %0.2f, kaon %0.2f \n",
310                calo.Data(),pid0[AliVCluster::kPhoton], pid0[AliVCluster::kPi0],
311                pid0[AliVCluster::kElectron], pid0[AliVCluster::kEleCon],
312                pid0[AliVCluster::kPion], pid0[AliVCluster::kKaon], pid0[AliVCluster::kProton],
313                pid0[AliVCluster::kNeutron], pid0[AliVCluster::kKaon0]);
314                   }
315                   
316       fEMCALPIDUtils->ComputePID(energy, lambda0);
317       Double_t pid[AliPID::kSPECIESN];
318       for(Int_t i = 0; i < AliPID::kSPECIESN; i++) pid[i] = fEMCALPIDUtils->GetPIDFinal(i);
319       return GetIdentifiedParticleType(calo, pid, energy);
320                   
321     }
322     
323     // If no use of bayesian, simplest PID  
324     if(lambda0 < fL0CutMax && lambda0 > fL0CutMin) return kPhoton ;
325     else return  kPi0 ; // Wild guess, it can be hadron but not photon
326   
327   }//EMCAL
328   else {//PHOS
329     
330     // Do not use bayesian, cut based on shower ellipse
331     return IsPHOSPhoton(lambda0,lambda1);
332   
333   }
334     
335 }
336
337 //__________________________________________________
338 TString  AliCaloPID::GetPIDParametersList()  {
339   //Put data member values in string to keep in output container
340   
341   TString parList ; //this will be list of parameters used for this analysis.
342   const Int_t buffersize = 255;
343   char onePar[buffersize] ;
344   snprintf(onePar,buffersize,"--- AliCaloPID ---\n") ;
345   parList+=onePar ;     
346   snprintf(onePar,buffersize,"fL0CutMin =%2.2f, fL0CutMax =%2.2f  (Cut on Shower Shape, used in PID evaluation) \n",fL0CutMin, fL0CutMax) ;
347   parList+=onePar ;
348   snprintf(onePar,buffersize,"fTOFCut  =%e (Cut on TOF, used in PID evaluation) \n",fTOFCut) ;
349   parList+=onePar ;
350   snprintf(onePar,buffersize,"fEMCALPhotonWeight =%2.2f (EMCAL bayesian weight for photons)\n",fEMCALPhotonWeight) ;
351   parList+=onePar ;
352   snprintf(onePar,buffersize,"fEMCALPi0Weight =%2.2f (EMCAL bayesian weight for pi0)\n",fEMCALPi0Weight) ;
353   parList+=onePar ;
354   snprintf(onePar,buffersize,"fEMCALElectronWeight =%2.2f(EMCAL bayesian weight for electrons)\n",fEMCALElectronWeight) ;
355   parList+=onePar ;
356   snprintf(onePar,buffersize,"fEMCALChargeWeight =%2.2f (EMCAL bayesian weight for charged hadrons)\n",fEMCALChargeWeight) ;
357   parList+=onePar ;
358   snprintf(onePar,buffersize,"fEMCALNeutralWeight =%2.2f (EMCAL bayesian weight for neutral hadrons)\n",fEMCALNeutralWeight) ;
359   parList+=onePar ;
360   snprintf(onePar,buffersize,"fPHOSPhotonWeight =%2.2f (PHOS bayesian weight for photons)\n",fPHOSPhotonWeight) ;
361   parList+=onePar ;
362   snprintf(onePar,buffersize,"fPHOSPi0Weight =%2.2f (PHOS bayesian weight for pi0)\n",fPHOSPi0Weight) ;
363   parList+=onePar ;
364   snprintf(onePar,buffersize,"fPHOSElectronWeight =%2.2f(PHOS bayesian weight for electrons)\n",fPHOSElectronWeight) ;
365   parList+=onePar ;
366   snprintf(onePar,buffersize,"fPHOSChargeWeight =%2.2f (PHOS bayesian weight for charged hadrons)\n",fPHOSChargeWeight) ;
367   parList+=onePar ;
368   snprintf(onePar,buffersize,"fPHOSNeutralWeight =%2.2f (PHOS bayesian weight for neutral hadrons)\n",fPHOSNeutralWeight) ;
369   parList+=onePar ;
370   
371   if(fPHOSWeightFormula){
372     snprintf(onePar,buffersize,"PHOS Photon Weight Formula: %s\n",fPHOSPhotonWeightFormulaExpression.Data() ) ;
373     parList+=onePar;
374     snprintf(onePar,buffersize,"PHOS Pi0    Weight Formula: %s\n",fPHOSPi0WeightFormulaExpression.Data()    ) ;
375     parList+=onePar;      
376   }
377   
378   return parList; 
379   
380 }
381
382 //________________________________________________________________
383 void AliCaloPID::Print(const Option_t * opt) const
384 {
385   
386   //Print some relevant parameters set for the analysis
387   if(! opt)
388     return;
389   
390   printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
391   
392   printf("PHOS PID weight , photon %0.2f, pi0 %0.2f, e %0.2f, charge %0.2f, neutral %0.2f \n",  
393          fPHOSPhotonWeight,  fPHOSPi0Weight, 
394          fPHOSElectronWeight,  fPHOSChargeWeight,   fPHOSNeutralWeight) ; 
395   printf("EMCAL PID weight, photon %0.2f, pi0 %0.2f, e %0.2f, charge %0.2f, neutral %0.2f\n",   
396          fEMCALPhotonWeight,  fEMCALPi0Weight, 
397          fEMCALElectronWeight,  fEMCALChargeWeight,  fEMCALNeutralWeight) ; 
398   
399   printf("PHOS Parametrized weight on?  =     %d\n",  fPHOSWeightFormula) ; 
400   if(fPHOSWeightFormula){
401     printf("Photon weight formula = %s\n", fPHOSPhotonWeightFormulaExpression.Data());
402     printf("Pi0    weight formula = %s\n", fPHOSPi0WeightFormulaExpression   .Data());
403   }
404   
405   printf("TOF cut        = %e\n",fTOFCut);
406   printf("Lambda0 cut min = %2.2f; max = %2.2f\n",fL0CutMin, fL0CutMax);
407   printf("Debug level    = %d\n",fDebug);
408   printf("Recalculate Bayesian?    = %d\n",fRecalculateBayesian);
409  if(fRecalculateBayesian) printf("Particle Flux?    = %d\n",fParticleFlux);
410   printf(" \n");
411   
412
413
414 //_______________________________________________________________
415 void AliCaloPID::SetPIDBits(const TString calo, const AliVCluster * cluster, AliAODPWG4Particle * ph, const AliCalorimeterUtils* cu) {
416   //Set Bits for PID selection
417   
418   //Dispersion/lambdas
419   //Double_t disp= cluster->GetDispersion()  ;
420   Double_t l1  = cluster->GetM20() ;
421   Double_t l0  = cluster->GetM02() ; 
422   Bool_t isDispOK = kTRUE ;
423   if(cluster->IsPHOS()){ 
424     
425    isDispOK =  IsPHOSPhoton(l0,l1);
426     
427   }
428   else{//EMCAL
429     
430     if(l0 > fL0CutMin && l0 < fL0CutMax) isDispOK = kTRUE;
431     
432   }
433   
434   ph->SetDispBit(isDispOK) ;
435   
436   //TOF
437   Double_t tof=cluster->GetTOF()  ;
438   ph->SetTOFBit(TMath::Abs(tof)<fTOFCut) ; 
439   
440   //Charged veto  
441   //Bool_t isNeutral = kTRUE ;
442   //if(cluster->IsPHOS())  isNeutral = cluster->GetEmcCpvDistance() > 5. ;
443   //else 
444   Bool_t isNeutral = IsTrackMatched(cluster,cu);
445   
446   ph->SetChargedBit(isNeutral);
447   
448   //Set PID pdg
449   ph->SetIdentifiedParticleType(GetIdentifiedParticleType(calo,cluster->GetPID(),ph->E()));
450   
451   if(fDebug > 0){ 
452     printf("AliCaloPID::SetPIDBits: TOF %e, Lambda0 %2.2f, Lambda1 %2.2f\n",tof , l0, l1);      
453     printf("AliCaloPID::SetPIDBits: pdg %d, bits: TOF %d, Dispersion %d, Charge %d\n",
454              ph->GetIdentifiedParticleType(), ph->GetTOFBit() , ph->GetDispBit() , ph->GetChargedBit()); 
455   }
456 }
457
458 //__________________________________________________________________________
459 Bool_t AliCaloPID::IsTrackMatched(const AliVCluster* cluster, const AliCalorimeterUtils * cu) const {
460   //Check if there is any track attached to this cluster
461   
462   Int_t nMatches = cluster->GetNTracksMatched();
463   
464   //  if(nMatches>0){
465   //    printf("N matches %d, first match (ESD) %d or (AOD) %d\n",nMatches,cluster->GetTrackMatchedIndex(), cluster->GetTrackMatched(0));
466   //    if     (cluster->GetTrackMatched(0))        printf("\t matched track id %d\n",((AliVTrack*)cluster->GetTrackMatched(0)) ->GetID() ) ;
467   //  }
468   //  else {
469   //    printf("Not Matched");
470   //  }
471   
472   //If EMCAL track matching needs to be recalculated
473   if(cluster->IsEMCAL() && cu && cu->IsRecalculationOfClusterTrackMatchingOn()){
474     Float_t dR = 999., dZ = 999.;
475     cu->GetEMCALRecoUtils()->GetMatchedResiduals(cluster->GetID(),dR,dZ);
476     
477     if(dR < 999) {     
478       
479       if(fhTrackMatchedDEta){
480         fhTrackMatchedDEta->Fill(cluster->E(),dZ);
481         fhTrackMatchedDPhi->Fill(cluster->E(),dR);
482         if(cluster->E() > 0.5)fhTrackMatchedDEtaDPhi->Fill(dZ,dR);
483       }
484       //printf("dR %f, dZ %f \n",dR,dZ);
485       return kTRUE;
486     }
487     else         
488       return kFALSE;
489   }//EMCAL cluster and recalculation of matching on
490   
491   if(fhTrackMatchedDEta){
492     fhTrackMatchedDEta->Fill(cluster->GetTrackDz(),cluster->E());
493     fhTrackMatchedDPhi->Fill(cluster->GetTrackDx(),cluster->E());
494     if(cluster->E() > 0.5)fhTrackMatchedDEtaDPhi->Fill(cluster->GetTrackDz(),cluster->GetTrackDx());
495   }
496   
497   if(!strcmp("AliESDCaloCluster",Form("%s",cluster->ClassName())))
498   {    
499     if (nMatches > 0) {
500       if (nMatches == 1 ) {
501         Int_t iESDtrack = cluster->GetTrackMatchedIndex();
502         //printf("Track Matched index %d\n",iESDtrack);
503         if(iESDtrack==-1) return kFALSE ;// Default value of array, there is no match
504         else  {
505           if(cluster->IsPHOS()) {
506             Double_t r = TMath::Sqrt(cluster->GetTrackDx()*cluster->GetTrackDx()+cluster->GetTrackDz()*cluster->GetTrackDz());
507             if( r < fRcutPHOS) return kTRUE ;// Default value of array, there is no match
508             else               return kFALSE;
509           }//PHOS
510           else return kTRUE; // EMCAL and not default
511         }
512       }//Just one, check
513       else return kTRUE ;//More than one, there is a match.
514     }// > 0
515     else return kFALSE; //It does not happen, but in case
516     
517   }//ESDs
518   else
519   {    
520     if(nMatches > 0) return kTRUE; //There is at least one match.
521     else             return kFALSE;
522     
523   }//AODs or MC (copy into AOD)
524   
525   return kFALSE;
526   
527 }
528
529 //__________________________________________________________________________
530 Bool_t AliCaloPID::IsPHOSPhoton(const Double_t l0, const Double_t l1) {
531   // Check the shape of the PHOS cluster
532   // Return true if photon like, from Dmitri P.
533   // TO DO, move parameters to data members
534   
535   Double_t l0Mean  = 1.22 ;
536   Double_t l1Mean  = 2.0  ;
537   Double_t l0Sigma = 0.42 ;
538   Double_t l1Sigma = 0.71 ;
539   Double_t c       =-0.59 ;
540   
541   Double_t R2=(l1-l1Mean)*(l0-l0Mean)/l0Sigma/l0Sigma+(l1-l1Mean)*(l1-l1Mean)/l1Sigma/l1Sigma-c*(l0-l0Mean)*(l1-l1Mean)/l0Sigma/l1Sigma ;
542   
543   return (R2<9.) ;
544   
545 }
546