1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 $ */
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
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.
36 // All these methods can be called in the analysis you are interested.
38 //*-- Author: Gustavo Conesa (LNF-INFN)
39 //////////////////////////////////////////////////////////////////////////////
42 // --- ROOT system ---
46 //---- ANALYSIS system ----
47 #include "AliCaloPID.h"
48 #include "AliVCluster.h"
49 #include "AliVTrack.h"
50 #include "AliAODPWG4Particle.h"
51 #include "AliCalorimeterUtils.h"
56 //________________________________________________
57 AliCaloPID::AliCaloPID() :
59 fEMCALPhotonWeight(0.), fEMCALPi0Weight(0.),
60 fEMCALElectronWeight(0.), fEMCALChargeWeight(0.), fEMCALNeutralWeight(0.),
61 fPHOSPhotonWeight(0.), fPHOSPi0Weight(0.),
62 fPHOSElectronWeight(0.), fPHOSChargeWeight(0.) , fPHOSNeutralWeight(0.),
63 fPHOSWeightFormula(0), fPHOSPhotonWeightFormula(0), fPHOSPi0WeightFormula(0),
64 fPHOSPhotonWeightFormulaExpression(""), fPHOSPi0WeightFormulaExpression(""),
65 fL0CutMax(100.), fL0CutMin(0), fTOFCut(0.),
66 fRcutPHOS(0), fDebug(-1),
67 fRecalculateBayesian(kFALSE), fParticleFlux(kLow), fEMCALPIDUtils(),
68 fHistoNEBins(100), fHistoEMax(100.), fHistoEMin(0.),
69 fHistoNDEtaBins(100), fHistoDEtaMax(0.15), fHistoDEtaMin(-0.15),
70 fHistoNDPhiBins(100), fHistoDPhiMax(0.15), fHistoDPhiMin(-0.15),
71 fhTrackMatchedDEta(0x0), fhTrackMatchedDPhi(0x0), fhTrackMatchedDEtaDPhi(0x0)
75 //Initialize parameters
79 //________________________________________________
80 AliCaloPID::AliCaloPID(const Int_t flux) :
82 fEMCALPhotonWeight(0.), fEMCALPi0Weight(0.),
83 fEMCALElectronWeight(0.), fEMCALChargeWeight(0.), fEMCALNeutralWeight(0.),
84 fPHOSPhotonWeight(0.), fPHOSPi0Weight(0.),
85 fPHOSElectronWeight(0.), fPHOSChargeWeight(0.) , fPHOSNeutralWeight(0.),
86 fPHOSWeightFormula(0), fPHOSPhotonWeightFormula(0), fPHOSPi0WeightFormula(0),
87 fPHOSPhotonWeightFormulaExpression(""), fPHOSPi0WeightFormulaExpression(""),
88 fL0CutMax(100.), fL0CutMin(0), fTOFCut(0.),
89 fRcutPHOS(0), fDebug(-1),
90 fRecalculateBayesian(kFALSE), fParticleFlux(flux), fEMCALPIDUtils(),
91 fHistoNEBins(100), fHistoEMax(100.), fHistoEMin(0.),
92 fHistoNDEtaBins(100), fHistoDEtaMax(0.15), fHistoDEtaMin(-0.15),
93 fHistoNDPhiBins(100), fHistoDPhiMax(0.15), fHistoDPhiMin(-0.15),
94 fhTrackMatchedDEta(0x0), fhTrackMatchedDPhi(0x0), fhTrackMatchedDEtaDPhi(0x0)
98 //Initialize parameters
103 //________________________________________________
104 AliCaloPID::AliCaloPID(const TTask * emcalpid) :
106 fEMCALPhotonWeight(0.), fEMCALPi0Weight(0.),
107 fEMCALElectronWeight(0.), fEMCALChargeWeight(0.), fEMCALNeutralWeight(0.),
108 fPHOSPhotonWeight(0.), fPHOSPi0Weight(0.),
109 fPHOSElectronWeight(0.), fPHOSChargeWeight(0.) , fPHOSNeutralWeight(0.),
110 fPHOSWeightFormula(0), fPHOSPhotonWeightFormula(0), fPHOSPi0WeightFormula(0),
111 fPHOSPhotonWeightFormulaExpression(""), fPHOSPi0WeightFormulaExpression(""),
112 fL0CutMax(100.), fL0CutMin(0), fTOFCut(0.),
113 fRcutPHOS(0), fDebug(-1),
114 fRecalculateBayesian(kFALSE), fParticleFlux(kLow), fEMCALPIDUtils((AliEMCALPIDUtils*) emcalpid),
115 fHistoNEBins(100), fHistoEMax(100.), fHistoEMin(0.),
116 fHistoNDEtaBins(100), fHistoDEtaMax(0.15), fHistoDEtaMin(-0.15),
117 fHistoNDPhiBins(100), fHistoDPhiMax(0.15), fHistoDPhiMin(-0.15),
118 fhTrackMatchedDEta(0x0), fhTrackMatchedDPhi(0x0), fhTrackMatchedDEtaDPhi(0x0)
122 //Initialize parameters
126 //_________________________________
127 AliCaloPID::~AliCaloPID() {
130 delete fPHOSPhotonWeightFormula ;
131 delete fPHOSPi0WeightFormula ;
132 delete fEMCALPIDUtils ;
136 //________________________________________________________________________
137 TList * AliCaloPID::GetCreateOutputObjects()
139 // Create histograms to be saved in output file and
140 // store them in outputContainer of the analysis class that calls this class.
142 TList * outputContainer = new TList() ;
143 outputContainer->SetName("CaloPIDHistos") ;
145 outputContainer->SetOwner(kFALSE);
147 fhTrackMatchedDEta = new TH2F
149 "d#eta of cluster-track vs cluster energy",
150 fHistoNEBins,fHistoEMin,fHistoEMax,fHistoNDEtaBins,fHistoDEtaMin,fHistoDEtaMax);
151 fhTrackMatchedDEta->SetYTitle("d#eta");
152 fhTrackMatchedDEta->SetXTitle("E_{cluster} (GeV)");
154 fhTrackMatchedDPhi = new TH2F
156 "d#phi of cluster-track vs cluster energy"
157 ,fHistoNEBins,fHistoEMin,fHistoEMax,fHistoNDPhiBins,fHistoDPhiMin,fHistoDPhiMax);
158 fhTrackMatchedDPhi->SetYTitle("d#phi (rad)");
159 fhTrackMatchedDPhi->SetXTitle("E_{cluster} (GeV)");
161 fhTrackMatchedDEtaDPhi = new TH2F
162 ("TrackMatchedDEtaDPhi",
163 "d#eta vs d#phi of cluster-track vs cluster energy"
164 ,fHistoNDEtaBins,fHistoDEtaMin,fHistoDEtaMax,fHistoNDPhiBins,fHistoDPhiMin,fHistoDPhiMax);
165 fhTrackMatchedDEtaDPhi->SetYTitle("d#phi (rad)");
166 fhTrackMatchedDEtaDPhi->SetXTitle("d#eta");
168 outputContainer->Add(fhTrackMatchedDEta) ;
169 outputContainer->Add(fhTrackMatchedDPhi) ;
170 outputContainer->Add(fhTrackMatchedDEtaDPhi) ;
172 return outputContainer;
177 //_______________________________________________________________
178 void AliCaloPID::InitParameters()
180 //Initialize the parameters of the PID.
182 fEMCALPhotonWeight = 0.6 ;
183 fEMCALPi0Weight = 0.6 ;
184 fEMCALElectronWeight = 0.6 ;
185 fEMCALChargeWeight = 0.6 ;
186 fEMCALNeutralWeight = 0.6 ;
188 fPHOSPhotonWeight = 0.6 ;
189 fPHOSPi0Weight = 0.6 ;
190 fPHOSElectronWeight = 0.6 ;
191 fPHOSChargeWeight = 0.6 ;
192 fPHOSNeutralWeight = 0.6 ;
194 //Formula to set the PID weight threshold for photon or pi0
195 fPHOSWeightFormula = kFALSE;
196 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))";
197 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))" ;
206 if(fRecalculateBayesian){
207 if(fParticleFlux == kLow){
208 printf("AliCaloPID::Init() - SetLOWFluxParam\n");
209 fEMCALPIDUtils->SetLowFluxParam() ;
211 else if (fParticleFlux == kHigh){
212 printf("AliCaloPID::Init() - SetHIGHFluxParam\n");
213 fEMCALPIDUtils->SetHighFluxParam() ;
218 //_______________________________________________________________
219 Int_t AliCaloPID::GetIdentifiedParticleType(const TString calo, const Double_t * pid, const Float_t energy) {
220 //Return most probable identity of the particle.
223 printf("AliCaloPID::GetIdentifiedParticleType() - pid pointer not initialized!!!\n");
227 Float_t wPh = fPHOSPhotonWeight ;
228 Float_t wPi0 = fPHOSPi0Weight ;
229 Float_t wE = fPHOSElectronWeight ;
230 Float_t wCh = fPHOSChargeWeight ;
231 Float_t wNe = fPHOSNeutralWeight ;
233 if(calo == "PHOS" && fPHOSWeightFormula){
234 wPh = GetPHOSPhotonWeightFormula()->Eval(energy) ;
235 wPi0 = GetPHOSPi0WeightFormula() ->Eval(energy);
240 wPh = fEMCALPhotonWeight ;
241 wPi0 = fEMCALPi0Weight ;
242 wE = fEMCALElectronWeight ;
243 wCh = fEMCALChargeWeight ;
244 wNe = fEMCALNeutralWeight ;
248 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",
249 calo.Data(),pid[AliVCluster::kPhoton], pid[AliVCluster::kPi0],
250 pid[AliVCluster::kElectron], pid[AliVCluster::kEleCon],
251 pid[AliVCluster::kPion], pid[AliVCluster::kKaon], pid[AliVCluster::kProton],
252 pid[AliVCluster::kNeutron], pid[AliVCluster::kKaon0]);
254 Int_t pdg = kNeutralUnknown ;
255 Float_t chargedHadronWeight = pid[AliVCluster::kProton]+pid[AliVCluster::kKaon]+
256 pid[AliVCluster::kPion]+pid[AliVCluster::kMuon];
257 Float_t neutralHadronWeight = pid[AliVCluster::kNeutron]+pid[AliVCluster::kKaon0];
258 Float_t allChargedWeight = pid[AliVCluster::kElectron]+pid[AliVCluster::kEleCon]+ chargedHadronWeight;
259 Float_t allNeutralWeight = pid[AliVCluster::kPhoton]+pid[AliVCluster::kPi0]+ neutralHadronWeight;
261 //Select most probable ID
263 if(pid[AliVCluster::kPhoton] > wPh) pdg = kPhoton ;
264 else if(pid[AliVCluster::kPi0] > wPi0) pdg = kPi0 ;
265 else if(pid[AliVCluster::kElectron] > wE) pdg = kElectron ;
266 else if(pid[AliVCluster::kEleCon] > wE) pdg = kEleCon ;
267 else if(chargedHadronWeight > wCh) pdg = kChargedHadron ;
268 else if(neutralHadronWeight > wNe) pdg = kNeutralHadron ;
269 else if(allChargedWeight > allNeutralWeight)
270 pdg = kChargedUnknown ;
272 pdg = kNeutralUnknown ;
276 if(pid[AliVCluster::kPhoton] > wPh) pdg = kPhoton ;
277 else if(pid[AliVCluster::kElectron] > wE) pdg = kElectron ;
278 else if(pid[AliVCluster::kPhoton]+pid[AliVCluster::kElectron] > wPh) pdg = kPhoton ; //temporal sollution until track matching for electrons is considered
279 else if(pid[AliVCluster::kPi0] > wPi0) pdg = kPi0 ;
280 else if(chargedHadronWeight + neutralHadronWeight > wCh) pdg = kChargedHadron ;
281 else if(neutralHadronWeight + chargedHadronWeight > wNe) pdg = kNeutralHadron ;
282 else pdg = kNeutralUnknown ;
285 if(fDebug > 0)printf("AliCaloPID::GetIdentifiedParticleType:Final Pdg: %d, cluster energy %2.2f \n", pdg,energy);
291 //_______________________________________________________________
292 Int_t AliCaloPID::GetIdentifiedParticleType(const TString calo,const TLorentzVector mom, const AliVCluster * cluster) {
293 //Recalculated PID with all parameters
295 Float_t lambda0 = cluster->GetM02();
296 Float_t lambda1 = cluster->GetM20();
297 Float_t energy = mom.E();
299 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",
300 calo.Data(),energy,lambda0,cluster->GetM20(),cluster->GetDispersion(),cluster->GetTOF(),
301 cluster->GetEmcCpvDistance(), cluster->GetDistanceToBadChannel(),cluster->GetNExMax());
303 if(calo == "EMCAL") {
304 //Recalculate Bayesian
305 if(fRecalculateBayesian){
307 const Double_t *pid0 = cluster->GetPID();
308 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",
309 calo.Data(),pid0[AliVCluster::kPhoton], pid0[AliVCluster::kPi0],
310 pid0[AliVCluster::kElectron], pid0[AliVCluster::kEleCon],
311 pid0[AliVCluster::kPion], pid0[AliVCluster::kKaon], pid0[AliVCluster::kProton],
312 pid0[AliVCluster::kNeutron], pid0[AliVCluster::kKaon0]);
315 fEMCALPIDUtils->ComputePID(energy, lambda0);
316 Double_t pid[AliPID::kSPECIESN];
317 for(Int_t i = 0; i < AliPID::kSPECIESN; i++) pid[i] = fEMCALPIDUtils->GetPIDFinal(i);
318 return GetIdentifiedParticleType(calo, pid, energy);
322 // If no use of bayesian, simplest PID
323 if(lambda0 < fL0CutMax && lambda0 > fL0CutMin) return kPhoton ;
324 else return kPi0 ; // Wild guess, it can be hadron but not photon
329 // Do not use bayesian, cut based on shower ellipse
330 return IsPHOSPhoton(lambda0,lambda1);
336 //__________________________________________________
337 TString AliCaloPID::GetPIDParametersList() {
338 //Put data member values in string to keep in output container
340 TString parList ; //this will be list of parameters used for this analysis.
341 const Int_t buffersize = 255;
342 char onePar[buffersize] ;
343 snprintf(onePar,buffersize,"--- AliCaloPID ---\n") ;
345 snprintf(onePar,buffersize,"fL0CutMin =%2.2f, fL0CutMax =%2.2f (Cut on Shower Shape, used in PID evaluation) \n",fL0CutMin, fL0CutMax) ;
347 snprintf(onePar,buffersize,"fTOFCut =%e (Cut on TOF, used in PID evaluation) \n",fTOFCut) ;
349 snprintf(onePar,buffersize,"fEMCALPhotonWeight =%2.2f (EMCAL bayesian weight for photons)\n",fEMCALPhotonWeight) ;
351 snprintf(onePar,buffersize,"fEMCALPi0Weight =%2.2f (EMCAL bayesian weight for pi0)\n",fEMCALPi0Weight) ;
353 snprintf(onePar,buffersize,"fEMCALElectronWeight =%2.2f(EMCAL bayesian weight for electrons)\n",fEMCALElectronWeight) ;
355 snprintf(onePar,buffersize,"fEMCALChargeWeight =%2.2f (EMCAL bayesian weight for charged hadrons)\n",fEMCALChargeWeight) ;
357 snprintf(onePar,buffersize,"fEMCALNeutralWeight =%2.2f (EMCAL bayesian weight for neutral hadrons)\n",fEMCALNeutralWeight) ;
359 snprintf(onePar,buffersize,"fPHOSPhotonWeight =%2.2f (PHOS bayesian weight for photons)\n",fPHOSPhotonWeight) ;
361 snprintf(onePar,buffersize,"fPHOSPi0Weight =%2.2f (PHOS bayesian weight for pi0)\n",fPHOSPi0Weight) ;
363 snprintf(onePar,buffersize,"fPHOSElectronWeight =%2.2f(PHOS bayesian weight for electrons)\n",fPHOSElectronWeight) ;
365 snprintf(onePar,buffersize,"fPHOSChargeWeight =%2.2f (PHOS bayesian weight for charged hadrons)\n",fPHOSChargeWeight) ;
367 snprintf(onePar,buffersize,"fPHOSNeutralWeight =%2.2f (PHOS bayesian weight for neutral hadrons)\n",fPHOSNeutralWeight) ;
370 if(fPHOSWeightFormula){
371 snprintf(onePar,buffersize,"PHOS Photon Weight Formula: %s\n",fPHOSPhotonWeightFormulaExpression.Data() ) ;
373 snprintf(onePar,buffersize,"PHOS Pi0 Weight Formula: %s\n",fPHOSPi0WeightFormulaExpression.Data() ) ;
381 //________________________________________________________________
382 void AliCaloPID::Print(const Option_t * opt) const
385 //Print some relevant parameters set for the analysis
389 printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
391 printf("PHOS PID weight , photon %0.2f, pi0 %0.2f, e %0.2f, charge %0.2f, neutral %0.2f \n",
392 fPHOSPhotonWeight, fPHOSPi0Weight,
393 fPHOSElectronWeight, fPHOSChargeWeight, fPHOSNeutralWeight) ;
394 printf("EMCAL PID weight, photon %0.2f, pi0 %0.2f, e %0.2f, charge %0.2f, neutral %0.2f\n",
395 fEMCALPhotonWeight, fEMCALPi0Weight,
396 fEMCALElectronWeight, fEMCALChargeWeight, fEMCALNeutralWeight) ;
398 printf("PHOS Parametrized weight on? = %d\n", fPHOSWeightFormula) ;
399 if(fPHOSWeightFormula){
400 printf("Photon weight formula = %s\n", fPHOSPhotonWeightFormulaExpression.Data());
401 printf("Pi0 weight formula = %s\n", fPHOSPi0WeightFormulaExpression .Data());
404 printf("TOF cut = %e\n",fTOFCut);
405 printf("Lambda0 cut min = %2.2f; max = %2.2f\n",fL0CutMin, fL0CutMax);
406 printf("Debug level = %d\n",fDebug);
407 printf("Recalculate Bayesian? = %d\n",fRecalculateBayesian);
408 if(fRecalculateBayesian) printf("Particle Flux? = %d\n",fParticleFlux);
413 //_______________________________________________________________
414 void AliCaloPID::SetPIDBits(const TString calo, const AliVCluster * cluster, AliAODPWG4Particle * ph, const AliCalorimeterUtils* cu) {
415 //Set Bits for PID selection
418 //Double_t disp= cluster->GetDispersion() ;
419 Double_t l1 = cluster->GetM20() ;
420 Double_t l0 = cluster->GetM02() ;
421 Bool_t isDispOK = kTRUE ;
422 if(cluster->IsPHOS()){
424 isDispOK = IsPHOSPhoton(l0,l1);
429 if(l0 > fL0CutMin && l0 < fL0CutMax) isDispOK = kTRUE;
433 ph->SetDispBit(isDispOK) ;
436 Double_t tof=cluster->GetTOF() ;
437 ph->SetTOFBit(TMath::Abs(tof)<fTOFCut) ;
440 //Bool_t isNeutral = kTRUE ;
441 //if(cluster->IsPHOS()) isNeutral = cluster->GetEmcCpvDistance() > 5. ;
443 Bool_t isNeutral = IsTrackMatched(cluster,cu);
445 ph->SetChargedBit(isNeutral);
448 ph->SetIdentifiedParticleType(GetIdentifiedParticleType(calo,cluster->GetPID(),ph->E()));
451 printf("AliCaloPID::SetPIDBits: TOF %e, Lambda0 %2.2f, Lambda1 %2.2f\n",tof , l0, l1);
452 printf("AliCaloPID::SetPIDBits: pdg %d, bits: TOF %d, Dispersion %d, Charge %d\n",
453 ph->GetIdentifiedParticleType(), ph->GetTOFBit() , ph->GetDispBit() , ph->GetChargedBit());
457 //__________________________________________________________________________
458 Bool_t AliCaloPID::IsTrackMatched(const AliVCluster* cluster, const AliCalorimeterUtils * cu) const {
459 //Check if there is any track attached to this cluster
461 Int_t nMatches = cluster->GetNTracksMatched();
464 // printf("N matches %d, first match (ESD) %d or (AOD) %d\n",nMatches,cluster->GetTrackMatchedIndex(), cluster->GetTrackMatched(0));
465 // if (cluster->GetTrackMatched(0)) printf("\t matched track id %d\n",((AliVTrack*)cluster->GetTrackMatched(0)) ->GetID() ) ;
468 // printf("Not Matched");
471 //If EMCAL track matching needs to be recalculated
472 if(cluster->IsEMCAL() && cu && cu->IsRecalculationOfClusterTrackMatchingOn()){
473 Float_t dR = 999., dZ = 999.;
474 cu->GetEMCALRecoUtils()->GetMatchedResiduals(cluster->GetID(),dR,dZ);
478 if(fhTrackMatchedDEta){
479 fhTrackMatchedDEta->Fill(cluster->E(),dZ);
480 fhTrackMatchedDPhi->Fill(cluster->E(),dR);
481 if(cluster->E() > 0.5)fhTrackMatchedDEtaDPhi->Fill(dZ,dR);
483 //printf("dR %f, dZ %f \n",dR,dZ);
488 }//EMCAL cluster and recalculation of matching on
490 if(fhTrackMatchedDEta){
491 fhTrackMatchedDEta->Fill(cluster->GetTrackDz(),cluster->E());
492 fhTrackMatchedDPhi->Fill(cluster->GetTrackDx(),cluster->E());
493 if(cluster->E() > 0.5)fhTrackMatchedDEtaDPhi->Fill(cluster->GetTrackDz(),cluster->GetTrackDx());
496 if(!strcmp("AliESDCaloCluster",Form("%s",cluster->ClassName())))
499 if (nMatches == 1 ) {
500 Int_t iESDtrack = cluster->GetTrackMatchedIndex();
501 //printf("Track Matched index %d\n",iESDtrack);
502 if(iESDtrack==-1) return kFALSE ;// Default value of array, there is no match
504 if(cluster->IsPHOS()) {
505 Double_t r = TMath::Sqrt(cluster->GetTrackDx()*cluster->GetTrackDx()+cluster->GetTrackDz()*cluster->GetTrackDz());
506 if( r < fRcutPHOS) return kTRUE ;// Default value of array, there is no match
509 else return kTRUE; // EMCAL and not default
512 else return kTRUE ;//More than one, there is a match.
514 else return kFALSE; //It does not happen, but in case
519 if(nMatches > 0) return kTRUE; //There is at least one match.
522 }//AODs or MC (copy into AOD)
528 //__________________________________________________________________________
529 Bool_t AliCaloPID::IsPHOSPhoton(const Double_t l0, const Double_t l1) {
530 // Check the shape of the PHOS cluster
531 // Return true if photon like, from Dmitri P.
532 // TO DO, move parameters to data members
534 Double_t l0Mean = 1.22 ;
535 Double_t l1Mean = 2.0 ;
536 Double_t l0Sigma = 0.42 ;
537 Double_t l1Sigma = 0.71 ;
540 Double_t R2=(l1-l1Mean)*(l0-l0Mean)/l0Sigma/l0Sigma+(l1-l1Mean)*(l1-l1Mean)/l1Sigma/l1Sigma-c*(l0-l0Mean)*(l1-l1Mean)/l0Sigma/l1Sigma ;