AliPhysics  1c9c77b (1c9c77b)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliCaloPID.cxx
Go to the documentation of this file.
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 // --- ROOT system ---
17 #include <TMath.h>
18 #include <TString.h>
19 #include <TList.h>
20 
21 // ---- ANALYSIS system ----
22 #include "AliCaloPID.h"
23 #include "AliAODCaloCluster.h"
24 #include "AliESDCaloCluster.h"
25 #include "AliVCaloCells.h"
26 #include "AliVTrack.h"
27 #include "AliAODPWG4Particle.h"
28 #include "AliCalorimeterUtils.h"
29 #include "AliFiducialCut.h" // detector enum definition
30 #include "AliVEvent.h"
31 #include "AliLog.h"
32 
33 // ---- Detector ----
34 #include "AliEMCALPIDUtils.h"
35 
39 
40 //________________________
43 //________________________
45 TObject(), fDebug(-1), fParticleFlux(kLow),
46 //Bayesian
47 fEMCALPIDUtils(), fUseBayesianWeights(kFALSE), fRecalculateBayesian(kFALSE),
48 fEMCALPhotonWeight(0.), fEMCALPi0Weight(0.),
49 fEMCALElectronWeight(0.), fEMCALChargeWeight(0.), fEMCALNeutralWeight(0.),
50 fPHOSPhotonWeight(0.), fPHOSPi0Weight(0.),
51 fPHOSElectronWeight(0.), fPHOSChargeWeight(0.) , fPHOSNeutralWeight(0.),
52 fPHOSWeightFormula(0), fPHOSPhotonWeightFormula(0), fPHOSPi0WeightFormula(0),
53 fPHOSPhotonWeightFormulaExpression(""),
54 fPHOSPi0WeightFormulaExpression(""),
55 //PID calculation
56 fEMCALL0CutMax(100.), fEMCALL0CutMin(0),
57 fEMCALDEtaCut(2000.), fEMCALDPhiCut(2000.),
58 fTOFCut(0.),
59 fPHOSDispersionCut(1000), fPHOSRCut(1000),
60 //Split
61 fUseSimpleMassCut(kFALSE),
62 fUseSimpleM02Cut(kFALSE),
63 fUseSplitAsyCut(kFALSE),
64 fUseSplitSSCut(kTRUE),
65 fSplitM02MaxCut(0), fSplitM02MinCut(0), fSplitMinNCells(0),
66 fMassEtaMin(0), fMassEtaMax(0),
67 fMassPi0Min(0), fMassPi0Max(0),
68 fMassPhoMin(0), fMassPhoMax(0),
69 fM02MaxParamShiftNLMN(0),
70 fSplitWidthSigma(0), fMassShiftHighECell(0)
71 {
73 }
74 
75 //________________________________________
80 //________________________________________
82 TObject(), fDebug(-1), fParticleFlux(flux),
83 //Bayesian
84 fEMCALPIDUtils(), fUseBayesianWeights(kFALSE), fRecalculateBayesian(kFALSE),
85 fEMCALPhotonWeight(0.), fEMCALPi0Weight(0.),
86 fEMCALElectronWeight(0.), fEMCALChargeWeight(0.), fEMCALNeutralWeight(0.),
87 fPHOSPhotonWeight(0.), fPHOSPi0Weight(0.),
88 fPHOSElectronWeight(0.), fPHOSChargeWeight(0.) , fPHOSNeutralWeight(0.),
89 fPHOSWeightFormula(0), fPHOSPhotonWeightFormula(0), fPHOSPi0WeightFormula(0),
90 fPHOSPhotonWeightFormulaExpression(""),
91 fPHOSPi0WeightFormulaExpression(""),
92 //PID calculation
93 fEMCALL0CutMax(100.), fEMCALL0CutMin(0),
94 fEMCALDEtaCut(2000.), fEMCALDPhiCut(2000.),
95 fTOFCut(0.),
96 fPHOSDispersionCut(1000), fPHOSRCut(1000),
97 //Split
98 fUseSimpleMassCut(kFALSE),
99 fUseSimpleM02Cut(kFALSE),
100 fUseSplitAsyCut(kFALSE),
101 fUseSplitSSCut(kTRUE),
102 fSplitM02MaxCut(0), fSplitM02MinCut(0), fSplitMinNCells(0),
103 fMassEtaMin(0), fMassEtaMax(0),
104 fMassPi0Min(0), fMassPi0Max(0),
105 fMassPhoMin(0), fMassPhoMax(0),
106 fM02MaxParamShiftNLMN(0),
107 fSplitWidthSigma(0), fMassShiftHighECell(0)
108 {
109  InitParameters();
110 }
111 
112 //_______________________________________________
117 //_______________________________________________
118 AliCaloPID::AliCaloPID(const TNamed * emcalpid) :
119 TObject(), fDebug(-1), fParticleFlux(kLow),
120 //Bayesian
121 fEMCALPIDUtils((AliEMCALPIDUtils*)emcalpid),
122 fUseBayesianWeights(kFALSE), fRecalculateBayesian(kFALSE),
123 fEMCALPhotonWeight(0.), fEMCALPi0Weight(0.),
124 fEMCALElectronWeight(0.), fEMCALChargeWeight(0.), fEMCALNeutralWeight(0.),
125 fPHOSPhotonWeight(0.), fPHOSPi0Weight(0.),
126 fPHOSElectronWeight(0.), fPHOSChargeWeight(0.) , fPHOSNeutralWeight(0.),
127 fPHOSWeightFormula(0), fPHOSPhotonWeightFormula(0), fPHOSPi0WeightFormula(0),
128 fPHOSPhotonWeightFormulaExpression(""),
129 fPHOSPi0WeightFormulaExpression(""),
130 //PID calculation
131 fEMCALL0CutMax(100.), fEMCALL0CutMin(0),
132 fEMCALDEtaCut(2000.), fEMCALDPhiCut(2000.),
133 fTOFCut(0.),
134 fPHOSDispersionCut(1000), fPHOSRCut(1000),
135 //Split
136 fUseSimpleMassCut(kFALSE),
137 fUseSimpleM02Cut(kFALSE),
138 fUseSplitAsyCut(kFALSE),
139 fUseSplitSSCut(kTRUE),
140 fSplitM02MaxCut(0), fSplitM02MinCut(0), fSplitMinNCells(0),
141 fMassEtaMin(0), fMassEtaMax(0),
142 fMassPi0Min(0), fMassPi0Max(0),
143 fMassPhoMin(0), fMassPhoMax(0),
144 fM02MaxParamShiftNLMN(0),
145 fSplitWidthSigma(0), fMassShiftHighECell(0)
146 
147 {
148  InitParameters();
149 }
150 
151 //_______________________
152 // Destructor.
153 //_______________________
155 {
156  delete fPHOSPhotonWeightFormula ;
157  delete fPHOSPi0WeightFormula ;
158  delete fEMCALPIDUtils ;
159 }
160 
161 //_______________________________
162 // Initialize the parameters of the PID.
163 //_______________________________
165 {
166  // Bayesian
167  fEMCALPhotonWeight = 0.6 ;
168  fEMCALPi0Weight = 0.6 ;
169  fEMCALElectronWeight = 0.6 ;
170  fEMCALChargeWeight = 0.6 ;
171  fEMCALNeutralWeight = 0.6 ;
172 
173  fPHOSPhotonWeight = 0.6 ;
174  fPHOSPi0Weight = 0.6 ;
175  fPHOSElectronWeight = 0.6 ;
176  fPHOSChargeWeight = 0.6 ;
177  fPHOSNeutralWeight = 0.6 ;
178 
179  //Formula to set the PID weight threshold for photon or pi0
180  fPHOSWeightFormula = kFALSE;
181  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))";
182  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))" ;
183 
185  {
186  if(fParticleFlux == kLow)
187  {
188  AliInfo("SetLOWFluxParam");
189  fEMCALPIDUtils->SetLowFluxParam() ;
190  }
191  else if (fParticleFlux == kHigh)
192  {
193  AliInfo("SetHighFluxParam");
194  fEMCALPIDUtils->SetHighFluxParam() ;
195  }
196  }
197 
198  //PID recalculation, not bayesian
199 
200  //EMCAL
201  fEMCALL0CutMax = 0.3 ;
202  fEMCALL0CutMin = 0.01;
203 
204  fEMCALDPhiCut = 0.05; // Same cut as in AliEMCALRecoUtils
205  fEMCALDEtaCut = 0.025;// Same cut as in AliEMCALRecoUtils
206 
207  // PHOS / EMCAL, not used
208  fTOFCut = 1.e-6;
209 
210  //PHOS
211  fPHOSRCut = 2. ;
212  fPHOSDispersionCut = 2.5;
213 
214  // Cluster splitting
215 
216  fSplitM02MinCut = 0.3 ;
217  fSplitM02MaxCut = 5 ;
218  fSplitMinNCells = 4 ;
219 
220  fMassEtaMin = 0.4;
221  fMassEtaMax = 0.6;
222 
223  fMassPi0Min = 0.11;
224  fMassPi0Max = 0.18;
225 
226  fMassPhoMin = 0.0;
227  fMassPhoMax = 0.08;
228 
229  fMassPi0Param[0][0] = 0 ; // Constant term on mass dependence
230  fMassPi0Param[0][1] = 0 ; // slope term on mass dependence
231  fMassPi0Param[0][2] = 0 ; // E function change
232  fMassPi0Param[0][3] = 0.044 ; // constant term on mass dependence
233  fMassPi0Param[0][4] = 0.0049; // slope term on mass dependence
234  fMassPi0Param[0][5] = 0.070 ; // Absolute low mass cut
235 
236  fMassPi0Param[1][0] = 0.115 ; // Constant term below 21 GeV
237  fMassPi0Param[1][1] = 0.00096; // slope term below 21 GeV
238  fMassPi0Param[1][2] = 21 ; // E function change
239  fMassPi0Param[1][3] = 0.10 ; // constant term on mass dependence
240  fMassPi0Param[1][4] = 0.0017; // slope term on mass dependence
241  fMassPi0Param[1][5] = 0.070 ; // Absolute low mass cut
242 
243  fWidthPi0Param[0][0] = 0.012 ; // Constant term on width dependence
244  fWidthPi0Param[0][1] = 0.0 ; // Slope term on width dependence
245  fWidthPi0Param[0][2] = 19 ; // E function change
246  fWidthPi0Param[0][3] = 0.0012; // Constant term on width dependence
247  fWidthPi0Param[0][4] = 0.0006; // Slope term on width dependence
248  fWidthPi0Param[0][5] = 0.0 ; // xx term
249 
250  fWidthPi0Param[1][0] = 0.009 ; // Constant term on width dependence
251  fWidthPi0Param[1][1] = 0.000 ; // Slope term on width dependence
252  fWidthPi0Param[1][2] = 10 ; // E function change
253  fWidthPi0Param[1][3] = 0.0023 ; // Constant term on width dependence
254  fWidthPi0Param[1][4] = 0.00067; // Slope term on width dependence
255  fWidthPi0Param[1][5] = 0.000 ;// xx term
256 
257  fMassShiftHighECell = 0; // Shift of cuts in case of higher energy threshold in cells, 5 MeV when Ecell>150 MeV
258 
259  //TF1 *lM02MinNLM1 = new TF1("M02MinNLM1","exp(2.135-0.245*x)",6,13.6);
260  fM02MinParam[0][0] = 2.135 ;
261  fM02MinParam[0][1] =-0.245 ;
262  fM02MinParam[0][2] = 0.0 ;
263  fM02MinParam[0][3] = 0.0 ;
264  fM02MinParam[0][4] = 0.0 ;
265 
266  // Same as NLM=1 for NLM=2
267  fM02MinParam[1][0] = 2.135 ;
268  fM02MinParam[1][1] =-0.245 ;
269  fM02MinParam[1][2] = 0.0 ;
270  fM02MinParam[1][3] = 0.0 ;
271  fM02MinParam[1][4] = 0.0 ;
272 
273  //TF1 *lM02MaxNLM1 = new TF1("M02MaxNLM1","exp(0.0662-0.0201*x)-0.0955+0.00186*x[0]+9.91/x[0]",6,100);
274  fM02MaxParam[0][0] = 0.0662 ;
275  fM02MaxParam[0][1] =-0.0201 ;
276  fM02MaxParam[0][2] =-0.0955 ;
277  fM02MaxParam[0][3] = 0.00186;
278  fM02MaxParam[0][4] = 9.91 ;
279 
280  //TF1 *lM02MaxNLM2 = new TF1("M02MaxNLM2","exp(0.353-0.0264*x)-0.524+0.00559*x[0]+21.9/x[0]",6,100);
281  fM02MaxParam[1][0] = 0.353 ;
282  fM02MaxParam[1][1] =-0.0264 ;
283  fM02MaxParam[1][2] =-0.524 ;
284  fM02MaxParam[1][3] = 0.00559;
285  fM02MaxParam[1][4] = 21.9 ;
286 
287  fM02MaxParamShiftNLMN = 0.75;
288 
289  //TF1 *lAsyNLM1 = new TF1("lAsyNLM1","0.96-879/(x*x*x)",5,100);
290  fAsyMinParam[0][0] = 0.96 ;
291  fAsyMinParam[0][1] = 0 ;
292  fAsyMinParam[0][2] =-879 ;
293  fAsyMinParam[0][3] = 0.96 ; // Absolute max
294 
295  //TF1 *lAsyNLM2 = new TF1("lAsyNLM2","0.95+0.0015*x-233/(x*x*x)",5,100);
296  fAsyMinParam[1][0] = 0.95 ;
297  fAsyMinParam[1][1] = 0.0015;
298  fAsyMinParam[1][2] =-233 ;
299  fAsyMinParam[1][3] = 1.0 ; // Absolute max
300 
301  fSplitEFracMin[0] = 0.0 ; // 0.96
302  fSplitEFracMin[1] = 0.0 ; // 0.96
303  fSplitEFracMin[2] = 0.0 ; // 0.7
304 
305  fSubClusterEMin[0] = 0.0; // 3 GeV
306  fSubClusterEMin[1] = 0.0; // 1 GeV
307  fSubClusterEMin[2] = 0.0; // 1 GeV
308 
309 
310  fSplitWidthSigma = 3. ;
311 }
312 
313 
314 //_________________________________________________________________________________________
320 //_________________________________________________________________________________________
322 {
323  if(!fUseSplitAsyCut) return kTRUE ;
324 
325  Float_t abasy = TMath::Abs(asy);
326 
327  Int_t inlm = nlm-1;
328  if(nlm > 2) inlm=1; // only 2 cases defined nlm=1 and nlm>=2
329 
330  // Get the parametrized min cut of asymmetry for NLM=2 up to 11 GeV
331 
332  Float_t cut = fAsyMinParam[inlm][0] + fAsyMinParam[inlm][1]*energy + fAsyMinParam[inlm][2]/energy/energy/energy ;
333 
334  // In any case and beyond validity energy range of the function,
335  // the parameter cannot be smaller than 1
336  if( cut > fAsyMinParam[inlm][3] ) cut = fAsyMinParam[inlm][3];
337 
338  //printf("energy %2.2f - nlm: %d (%d)- p0 %f, p1 %f, p2 %f, p3 %f ; cut: %2.2f\n",energy,nlm,inlm,
339  // fAsyMinParam[inlm][0],fAsyMinParam[inlm][1],fAsyMinParam[inlm][2],fAsyMinParam[inlm][3],cut);
340 
341  if(abasy < cut) return kTRUE;
342  else return kFALSE;
343 }
344 
345 //______________________________________________________________________________________
351 //______________________________________________________________________________________
353 {
355  {
356  if(mass < fMassPi0Max && mass > fMassPi0Min) return kTRUE;
357  else return kFALSE;
358  }
359 
360  // Get the selected mean value as reference for the mass
361  Int_t inlm = nlm-1;
362  if(nlm > 2) inlm=1; // only 2 cases defined nlm=1 and nlm>=2
363 
364  Float_t meanMass = energy * fMassPi0Param[inlm][1] + fMassPi0Param[inlm][0];
365  if(energy > fMassPi0Param[inlm][2]) meanMass = energy * fMassPi0Param[inlm][4] + fMassPi0Param[inlm][3];
366 
367  // In case of higher energy cell cut than 50 MeV, smaller mean mass 0-5 MeV, not really necessary
368  meanMass -= fMassShiftHighECell;
369 
370  // Get the parametrized width of the mass
371  Float_t width = 0.009;
372  if (energy > 8 && energy < fWidthPi0Param[inlm][2])
373  width = energy * fWidthPi0Param[inlm][1] + fWidthPi0Param[inlm][0];
374  else if( energy > fWidthPi0Param[inlm][2])
375  width = energy * energy * fWidthPi0Param[inlm][5] + energy * fWidthPi0Param[inlm][4] + fWidthPi0Param[inlm][3];
376 
377  // Calculate the 2 sigma cut
378  Float_t minMass = meanMass-fSplitWidthSigma*width;
379  Float_t maxMass = meanMass+fSplitWidthSigma*width;
380 
381  // In case of low energy, hard cut to avoid conversions
382  if(energy < 10 && minMass < fMassPi0Param[inlm][5] ) minMass = fMassPi0Param[inlm][5];
383 
384  //printf("E %2.2f, mass %1.1f, nlm %d: sigma %1.1f width %3.1f, mean Mass %3.0f, minMass %3.0f, maxMass %3.0f\n ",
385  // energy,mass *1000, inlm, fSplitWidthSigma, width*1000, meanMass*1000,minMass*1000,maxMass*1000);
386 
387  if(mass < maxMass && mass > minMass) return kTRUE;
388  else return kFALSE;
389 }
390 
391 //________________________________________________
395 //________________________________________________
397 {
398  Float_t minCut = fSplitM02MinCut;
399  Float_t maxCut = fSplitM02MaxCut;
400 
401  if(m02 < maxCut && m02 > minCut) return kTRUE;
402  else return kFALSE;
403 }
404 
405 //_______________________________________________________________________________
411 //_______________________________________________________________________________
413 {
414  if(!fUseSplitSSCut) return kTRUE ;
415 
416  //First check the absolute minimum and maximum
417  if(!IsInM02Range(m02)) return kFALSE ;
418 
419  //If requested, check the E dependent cuts
420  else if(!fUseSimpleM02Cut)
421  {
422  Int_t inlm = nlm-1;
423  if(nlm > 2) inlm=1; // only 2 cases defined nlm=1 and nlm>=2
424 
425  Float_t minCut = fSplitM02MinCut;
426  Float_t maxCut = fSplitM02MaxCut;
427 
428  //e^{a+bx} + c + dx + e/x
429  if(energy > 1) minCut = TMath::Exp( fM02MinParam[inlm][0] + fM02MinParam[inlm][1]*energy ) +
430  fM02MinParam[inlm][2] + fM02MinParam[inlm][3]*energy + fM02MinParam[inlm][4]/energy;
431 
432  if(energy > 1) maxCut = TMath::Exp( fM02MaxParam[inlm][0] + fM02MaxParam[inlm][1]*energy ) +
433  fM02MaxParam[inlm][2] + fM02MaxParam[inlm][3]*energy + fM02MaxParam[inlm][4]/energy;
434 
435  // In any case and beyond validity energy range of the function,
436  // the parameter cannot be smaller than 0.3 or larger than 4-5
437  if( minCut < fSplitM02MinCut) minCut = fSplitM02MinCut;
438  if( maxCut > fSplitM02MaxCut) maxCut = fSplitM02MaxCut;
439  if( nlm > 2 ) maxCut+=fM02MaxParamShiftNLMN;
440 
441  //if(energy > 7) printf("\t \t E %2.2f, nlm %d, m02 %2.2f, minM02 %2.2f, maxM02 %2.2f\n",energy, nlm, m02,minCut,maxCut);
442 
443  if(m02 < maxCut && m02 > minCut) return kTRUE;
444  else return kFALSE;
445 
446  }
447 
448  else return kTRUE;
449 }
450 
451 
452 //______________________________________________________________________________
453 // Select the appropriate shower shape range in splitting method to select eta's
454 // Use same parametrization as pi0, just shift the distributions (to be tuned)
459 //______________________________________________________________________________
461 {
462  if(!fUseSplitSSCut) return kTRUE ;
463 
464  //First check the absolute minimum and maximum
465  if(!IsInM02Range(m02)) return kFALSE ;
466 
467  //DO NOT USE, study parametrization
468 
469  //If requested, check the E dependent cuts
470  else if(!fUseSimpleM02Cut)
471  {
472  Int_t inlm = nlm-1;
473  if(nlm > 2) inlm=1; // only 2 cases defined nlm=1 and nlm>=2
474 
475  Float_t minCut = fSplitM02MinCut;
476  Float_t maxCut = fSplitM02MaxCut;
477 
478  Float_t shiftE = energy-20; // to be tuned
479  if(nlm==1) shiftE=energy-28;
480 
481  //e^{a+bx} + c + dx + e/x
482  if(shiftE > 1) minCut = TMath::Exp( fM02MinParam[inlm][0] + fM02MinParam[inlm][1]*shiftE ) +
483  fM02MinParam[inlm][2] + fM02MinParam[inlm][3]*shiftE + fM02MinParam[inlm][4]/shiftE;
484 
485  // In any case the parameter cannot be smaller than 0.3
486  if( minCut < fSplitM02MinCut) minCut = fSplitM02MinCut;
487 
488  shiftE = energy+20; // to be tuned
489 
490  if(shiftE > 1) maxCut = 1 + TMath::Exp( fM02MaxParam[inlm][0] + fM02MaxParam[inlm][1]*shiftE ) +
491  fM02MaxParam[inlm][2] + fM02MaxParam[inlm][3]*shiftE + fM02MaxParam[inlm][4]/shiftE;
492 
493  // In any case the parameter cannot be smaller than 4-5
494  if( maxCut > fSplitM02MaxCut) maxCut = fSplitM02MaxCut;
495  if( nlm > 2 ) maxCut+=fM02MaxParamShiftNLMN;
496 
497  //if(energy>6)printf("\t \t E %2.2f, nlm %d, m02 %2.2f, minM02 %2.2f, maxM02 %2.2f\n",energy, nlm, m02,minCut,maxCut);
498 
499  if(m02 < maxCut && m02 > minCut) return kTRUE;
500  else return kFALSE;
501 
502  }
503 
504  else return kTRUE;
505 }
506 
507 //______________________________________________________________________________
514 //______________________________________________________________________________
516 {
517  if(!fUseSplitSSCut) return kTRUE ;
518 
519  Float_t minCut = 0.1;
520  Float_t maxCut = fSplitM02MinCut;
521 
522  if(!fUseSimpleM02Cut)
523  {
524  Int_t inlm = nlm-1;
525  if(nlm > 2) inlm=1; // only 2 cases defined nlm=1 and nlm>=2
526 
527  //e^{a+bx} + c + dx + e/x
528  if(energy > 1) maxCut = TMath::Exp( fM02MinParam[inlm][0] + fM02MinParam[inlm][1]*energy ) +
529  fM02MinParam[inlm][2] + fM02MinParam[inlm][3]*energy + fM02MinParam[inlm][4]/energy;
530 
531  if( maxCut < fSplitM02MinCut) maxCut = fSplitM02MinCut;
532  }
533 
534  if(m02 < maxCut && m02 > minCut) return kTRUE;
535  else return kFALSE;
536 
537 }
538 
539 //______________________________________________
541 //______________________________________________
542 AliEMCALPIDUtils *AliCaloPID::GetEMCALPIDUtils()
543 {
544  if(!fEMCALPIDUtils) fEMCALPIDUtils = new AliEMCALPIDUtils ;
545 
546  return fEMCALPIDUtils ;
547 }
548 
549 
550 //________________________________________________________________
553 //________________________________________________________________
555 {
556  Float_t energy = cluster->E();
557  Float_t lambda0 = cluster->GetM02();
558  Float_t lambda1 = cluster->GetM20();
559 
560  // ---------------------
561  // Use bayesian approach
562  // ---------------------
563 
565  {
566  Double_t weights[AliPID::kSPECIESCN];
567 
568  if(cluster->IsEMCAL() && fRecalculateBayesian)
569  {
570  fEMCALPIDUtils->ComputePID(energy, lambda0);
571  for(Int_t i = 0; i < AliPID::kSPECIESCN; i++) weights[i] = fEMCALPIDUtils->GetPIDFinal(i);
572  }
573  else
574  {
575  for(Int_t i = 0; i < AliPID::kSPECIESCN; i++) weights[i] = cluster->GetPID()[i];
576  }
577 
578  if(fDebug > 0) PrintClusterPIDWeights(weights);
579 
580  return GetIdentifiedParticleTypeFromBayesWeights(cluster->IsEMCAL(), weights, energy);
581  }
582 
583  // -------------------------------------------------------
584  // Calculate PID SS from data, do not use bayesian weights
585  // -------------------------------------------------------
586 
587  AliDebug(1,Form("EMCAL %d?, E %3.2f, l0 %3.2f, l1 %3.2f, disp %3.2f, tof %1.11f, distCPV %3.2f, distToBC %1.1f, NMax %d",
588  cluster->IsEMCAL(),energy,lambda0,cluster->GetM20(),cluster->GetDispersion(),cluster->GetTOF(),
589  cluster->GetEmcCpvDistance(), cluster->GetDistanceToBadChannel(),cluster->GetNExMax()));
590 
591  if(cluster->IsEMCAL())
592  {
593  AliDebug(1,Form("EMCAL SS %f <%f < %f?",fEMCALL0CutMin, lambda0, fEMCALL0CutMax));
594 
595  if(lambda0 < fEMCALL0CutMax && lambda0 > fEMCALL0CutMin) return kPhoton ;
596  else return kNeutralUnknown ;
597  } // EMCAL
598  else // PHOS
599  {
600  if(TestPHOSDispersion(energy,lambda0,lambda1) < fPHOSDispersionCut) return kPhoton;
601  else return kNeutralUnknown;
602  }
603 }
604 
605 //_________________________________________________________________________________________________________
610 //_________________________________________________________________________________________________________
612 {
613  if(!pid)
614  {
615  AliFatal("pid pointer not initialized!!!");
616  return kNeutralUnknown; // not needed, added to content coverity
617  }
618 
619  Float_t wPh = fPHOSPhotonWeight ;
620  Float_t wPi0 = fPHOSPi0Weight ;
622  Float_t wCh = fPHOSChargeWeight ;
624 
625  if(!isEMCAL && fPHOSWeightFormula){
626  wPh = GetPHOSPhotonWeightFormula()->Eval(energy) ;
627  wPi0 = GetPHOSPi0WeightFormula() ->Eval(energy);
628  }
629  else
630  {
631  wPh = fEMCALPhotonWeight ;
632  wPi0 = fEMCALPi0Weight ;
633  wE = fEMCALElectronWeight ;
634  wCh = fEMCALChargeWeight ;
635  wNe = fEMCALNeutralWeight ;
636  }
637 
638  if(fDebug > 0) PrintClusterPIDWeights(pid);
639 
641  Float_t chargedHadronWeight = pid[AliVCluster::kProton]+pid[AliVCluster::kKaon]+
642  pid[AliVCluster::kPion]+pid[AliVCluster::kMuon];
643  Float_t neutralHadronWeight = pid[AliVCluster::kNeutron]+pid[AliVCluster::kKaon0];
644  Float_t allChargedWeight = pid[AliVCluster::kElectron]+pid[AliVCluster::kEleCon]+ chargedHadronWeight;
645  Float_t allNeutralWeight = pid[AliVCluster::kPhoton]+pid[AliVCluster::kPi0]+ neutralHadronWeight;
646 
647  //Select most probable ID
648  if(!isEMCAL) // PHOS
649  {
650  if(pid[AliVCluster::kPhoton] > wPh) pdg = kPhoton ;
651  else if(pid[AliVCluster::kPi0] > wPi0) pdg = kPi0 ;
652  else if(pid[AliVCluster::kElectron] > wE) pdg = kElectron ;
653  else if(pid[AliVCluster::kEleCon] > wE) pdg = kEleCon ;
654  else if(chargedHadronWeight > wCh) pdg = kChargedHadron ;
655  else if(neutralHadronWeight > wNe) pdg = kNeutralHadron ;
656  else if(allChargedWeight > allNeutralWeight)
657  pdg = kChargedUnknown ;
658  else
659  pdg = kNeutralUnknown ;
660  }
661  else //EMCAL
662  {
663  if(pid[AliVCluster::kPhoton] > wPh) pdg = kPhoton ;
664  else if(pid[AliVCluster::kElectron] > wE) pdg = kElectron ;
665  else if(pid[AliVCluster::kPhoton]+pid[AliVCluster::kElectron] > wPh) pdg = kPhoton ; //temporal sollution until track matching for electrons is considered
666  else if(pid[AliVCluster::kPi0] > wPi0) pdg = kPi0 ;
667  else if(chargedHadronWeight + neutralHadronWeight > wCh) pdg = kChargedHadron ;
668  else if(neutralHadronWeight + chargedHadronWeight > wNe) pdg = kNeutralHadron ;
669  else pdg = kNeutralUnknown ;
670  }
671 
672  AliDebug(1,Form("Final Pdg: %d, cluster energy %2.2f", pdg,energy));
673 
674  return pdg ;
675 
676 }
677 
678 //____________________________________________________________________________________________________________
697 //____________________________________________________________________________________________________________
699  AliVCaloCells* cells,
700  AliCalorimeterUtils * caloutils,
701  Double_t vertex[3],
702  Int_t & nMax,
703  Double_t & mass, Double_t & angle,
704  TLorentzVector & l1, TLorentzVector & l2,
705  Int_t & absId1, Int_t & absId2,
706  Float_t & distbad1, Float_t & distbad2,
707  Bool_t & fidcut1, Bool_t & fidcut2 ) const
708 {
709  Float_t eClus = cluster->E();
710  Float_t m02 = cluster->GetM02();
711  const Int_t nc = cluster->GetNCells();
712  Int_t absIdList[nc];
713  Float_t maxEList [nc];
714 
715  mass = -1.;
716  angle = -1.;
717 
718  //If too low number of cells, skip it
719  if ( nc < fSplitMinNCells) return kNeutralUnknown ;
720 
721  AliDebug(2,"\t pass nCells cut");
722 
723  // Get Number of local maxima
724  nMax = caloutils->GetNumberOfLocalMaxima(cluster, cells, absIdList, maxEList) ;
725 
726  AliDebug(1,Form("Cluster : E %1.1f, M02 %1.2f, NLM %d, N Cells %d",eClus,m02,nMax,nc));
727 
728  //---------------------------------------------------------------------
729  // Get the 2 max indeces and do inv mass
730  //---------------------------------------------------------------------
731 
733  if(cluster->IsPHOS()) calorimeter = AliFiducialCut::kPHOS;
734 
735  if ( nMax == 2 )
736  {
737  absId1 = absIdList[0];
738  absId2 = absIdList[1];
739 
740  //Order in energy
741  Float_t en1 = cells->GetCellAmplitude(absId1);
742  caloutils->RecalibrateCellAmplitude(en1,calorimeter,absId1);
743  Float_t en2 = cells->GetCellAmplitude(absId2);
744  caloutils->RecalibrateCellAmplitude(en2,calorimeter,absId2);
745  if(en1 < en2)
746  {
747  absId2 = absIdList[0];
748  absId1 = absIdList[1];
749  }
750  }
751  else if( nMax == 1 )
752  {
753 
754  absId1 = absIdList[0];
755 
756  //Find second highest energy cell
757 
758  Float_t enmax = 0 ;
759  for(Int_t iDigit = 0 ; iDigit < cluster->GetNCells() ; iDigit++)
760  {
761  Int_t absId = cluster->GetCellsAbsId()[iDigit];
762  if( absId == absId1 ) continue ;
763  Float_t endig = cells->GetCellAmplitude(absId);
764  caloutils->RecalibrateCellAmplitude(endig,calorimeter,absId);
765  if(endig > enmax)
766  {
767  enmax = endig ;
768  absId2 = absId ;
769  }
770  }// cell loop
771  }// 1 maxima
772  else
773  { // n max > 2
774  // loop on maxima, find 2 highest
775 
776  // First max
777  Float_t enmax = 0 ;
778  for(Int_t iDigit = 0 ; iDigit < nMax ; iDigit++)
779  {
780  Float_t endig = maxEList[iDigit];
781  if(endig > enmax)
782  {
783  enmax = endig ;
784  absId1 = absIdList[iDigit];
785  }
786  }// first maxima loop
787 
788  // Second max
789  Float_t enmax2 = 0;
790  for(Int_t iDigit = 0 ; iDigit < nMax ; iDigit++)
791  {
792  if(absIdList[iDigit]==absId1) continue;
793  Float_t endig = maxEList[iDigit];
794  if(endig > enmax2)
795  {
796  enmax2 = endig ;
797  absId2 = absIdList[iDigit];
798  }
799  }// second maxima loop
800 
801  } // n local maxima > 2
802 
803  if(absId2<0 || absId1<0)
804  {
805  AliDebug(1,Form("Bad index for local maxima : N max %d, i1 %d, i2 %d, cluster E %2.2f, ncells %d, m02 %2.2f",
806  nMax,absId1,absId2,eClus,nc,m02));
807  return kNeutralUnknown ;
808  }
809 
810  //---------------------------------------------------------------------
811  // Split the cluster energy in 2, around the highest 2 local maxima
812  //---------------------------------------------------------------------
813 
814  AliAODCaloCluster cluster1(0, 0,NULL,0.,NULL,NULL,1,0);
815  AliAODCaloCluster cluster2(1, 0,NULL,0.,NULL,NULL,1,0);
816 
817  caloutils->SplitEnergy(absId1,absId2,cluster, cells, &cluster1, &cluster2,nMax); /*absIdList, maxEList,*/
818 
819  fidcut1 = caloutils->GetEMCALRecoUtils()->CheckCellFiducialRegion(caloutils->GetEMCALGeometry(), &cluster1,cells);
820  fidcut2 = caloutils->GetEMCALRecoUtils()->CheckCellFiducialRegion(caloutils->GetEMCALGeometry(), &cluster2,cells);
821 
822  caloutils->GetEMCALRecoUtils()->RecalculateClusterDistanceToBadChannel(caloutils->GetEMCALGeometry(),cells,&cluster1);
823  caloutils->GetEMCALRecoUtils()->RecalculateClusterDistanceToBadChannel(caloutils->GetEMCALGeometry(),cells,&cluster2);
824 
825  distbad1 = cluster1.GetDistanceToBadChannel();
826  distbad2 = cluster2.GetDistanceToBadChannel();
827 // if(!fidcut2 || !fidcut1 || distbad1 < 2 || distbad2 < 2)
828 // printf("*** Dist to bad channel cl %f, cl1 %f, cl2 %f; fid cut cl %d, cl1 %d, cl2 %d \n",
829 // cluster->GetDistanceToBadChannel(),distbad1,distbad2,
830 // caloutils->GetEMCALRecoUtils()->CheckCellFiducialRegion(caloutils->GetEMCALGeometry(), cluster,cells),fidcut1,fidcut2);
831 
832  cluster1.GetMomentum(l1,vertex);
833  cluster2.GetMomentum(l2,vertex);
834 
835  mass = (l1+l2).M();
836  angle = l2.Angle(l1.Vect());
837  Float_t e1 = cluster1.E();
838  Float_t e2 = cluster2.E();
839 
840  // Consider clusters with splitted energy not too different to original cluster energy
841  Float_t splitFracCut = 0;
842  if(nMax < 3) splitFracCut = fSplitEFracMin[nMax-1];
843  else splitFracCut = fSplitEFracMin[2];
844  if((e1+e2)/eClus < splitFracCut) return kNeutralUnknown ;
845 
846  AliDebug(2,"\t pass Split E frac cut");
847 
848  // Consider sub-clusters with minimum energy
849  Float_t minECut = fSubClusterEMin[2];
850  if (nMax == 2) minECut = fSubClusterEMin[1];
851  else if(nMax == 1) minECut = fSubClusterEMin[0];
852  if(e1 < minECut || e2 < minECut)
853  {
854  //printf("Reject: e1 %2.1f, e2 %2.1f, cut %2.1f\n",e1,e2,minECut);
855  return kNeutralUnknown ;
856  }
857 
858  AliDebug(2,"\t pass min sub-cluster E cut");
859 
860  // Asymmetry of cluster
861  Float_t asy =-10;
862  if(e1+e2 > 0) asy = (e1-e2) / (e1+e2);
863 
864  if( !IsInPi0SplitAsymmetryRange(eClus,asy,nMax) ) return kNeutralUnknown ;
865 
866 
867  AliDebug(2,"\t pass asymmetry cut");
868 
869  Bool_t pi0OK = kFALSE;
870  Bool_t etaOK = kFALSE;
871  Bool_t conOK = kFALSE;
872 
873  //If too small or big M02, skip it
874  if (IsInPi0M02Range(eClus,m02,nMax)) pi0OK = kTRUE;
875  else if(IsInEtaM02Range(eClus,m02,nMax)) etaOK = kTRUE;
876  else if(IsInConM02Range(eClus,m02,nMax)) conOK = kTRUE;
877 
878  Float_t energy = eClus;
879  if(nMax > 2) energy = e1+e2; // In case of NLM>2 use mass cut for NLM=2 but for the split sum not the cluster energy that is not the pi0 E.
880 
881  // Check the mass, and set an ID to the splitted cluster
882  if ( conOK && mass < fMassPhoMax && mass > fMassPhoMin ) { AliDebug(2,"\t Split Conv"); return kPhoton ; }
883  else if( etaOK && mass < fMassEtaMax && mass > fMassEtaMin ) { AliDebug(2,"\t Split Eta" ); return kEta ; }
884  else if( pi0OK && IsInPi0SplitMassRange(energy,mass,nMax) ) { AliDebug(2,"\t Split Pi0" ); return kPi0 ; }
885  else return kNeutralUnknown ;
886 
887 }
888 
889 //_________________________________________
891 //_________________________________________
893 {
894  TString parList ; //this will be list of parameters used for this analysis.
895  const Int_t buffersize = 255;
896  char onePar[buffersize] ;
897  snprintf(onePar,buffersize,"--- AliCaloPID ---") ;
898  parList+=onePar ;
900  {
901  snprintf(onePar,buffersize,"fEMCALPhotonWeight =%2.2f (EMCAL bayesian weight for photons)",fEMCALPhotonWeight) ;
902  parList+=onePar ;
903  snprintf(onePar,buffersize,"fEMCALPi0Weight =%2.2f (EMCAL bayesian weight for pi0)",fEMCALPi0Weight) ;
904  parList+=onePar ;
905  snprintf(onePar,buffersize,"fEMCALElectronWeight =%2.2f(EMCAL bayesian weight for electrons)",fEMCALElectronWeight) ;
906  parList+=onePar ;
907  snprintf(onePar,buffersize,"fEMCALChargeWeight =%2.2f (EMCAL bayesian weight for charged hadrons)",fEMCALChargeWeight) ;
908  parList+=onePar ;
909  snprintf(onePar,buffersize,"fEMCALNeutralWeight =%2.2f (EMCAL bayesian weight for neutral hadrons)",fEMCALNeutralWeight) ;
910  parList+=onePar ;
911  snprintf(onePar,buffersize,"fPHOSPhotonWeight =%2.2f (PHOS bayesian weight for photons)",fPHOSPhotonWeight) ;
912  parList+=onePar ;
913  snprintf(onePar,buffersize,"fPHOSPi0Weight =%2.2f (PHOS bayesian weight for pi0)",fPHOSPi0Weight) ;
914  parList+=onePar ;
915  snprintf(onePar,buffersize,"fPHOSElectronWeight =%2.2f(PHOS bayesian weight for electrons)",fPHOSElectronWeight) ;
916  parList+=onePar ;
917  snprintf(onePar,buffersize,"fPHOSChargeWeight =%2.2f (PHOS bayesian weight for charged hadrons)",fPHOSChargeWeight) ;
918  parList+=onePar ;
919  snprintf(onePar,buffersize,"fPHOSNeutralWeight =%2.2f (PHOS bayesian weight for neutral hadrons)",fPHOSNeutralWeight) ;
920  parList+=onePar ;
921 
923  {
924  snprintf(onePar,buffersize,"PHOS Photon Weight Formula: %s",fPHOSPhotonWeightFormulaExpression.Data() ) ;
925  parList+=onePar;
926  snprintf(onePar,buffersize,"PHOS Pi0 Weight Formula: %s",fPHOSPi0WeightFormulaExpression.Data() ) ;
927  parList+=onePar;
928  }
929  }
930  else
931  {
932  snprintf(onePar,buffersize,"EMCAL: fEMCALL0CutMin =%2.2f, fEMCALL0CutMax =%2.2f (Cut on Shower Shape)",fEMCALL0CutMin, fEMCALL0CutMax) ;
933  parList+=onePar ;
934  snprintf(onePar,buffersize,"EMCAL: fEMCALDEtaCut =%2.2f, fEMCALDPhiCut =%2.2f (Cut on track matching)",fEMCALDEtaCut, fEMCALDPhiCut) ;
935  parList+=onePar ;
936  snprintf(onePar,buffersize,"fTOFCut =%e (Cut on TOF, used in PID evaluation)",fTOFCut) ;
937  parList+=onePar ;
938  snprintf(onePar,buffersize,"fPHOSRCut =%2.2f, fPHOSDispersionCut =%2.2f (Cut on Shower Shape and CPV)",fPHOSRCut,fPHOSDispersionCut) ;
939  parList+=onePar ;
940 
941  }
942 
943  if(fUseSimpleM02Cut)
944  {
945  snprintf(onePar,buffersize,"%2.2f< M02 < %2.2f", fSplitM02MinCut, fSplitM02MaxCut) ;
946  parList+=onePar ;
947  }
948  snprintf(onePar,buffersize,"fMinNCells =%d", fSplitMinNCells) ;
949  parList+=onePar ;
951  {
952  snprintf(onePar,buffersize,"pi0 : %2.1f < m <%2.1f", fMassPi0Min,fMassPi0Max);
953  parList+=onePar ;
954  }
955  snprintf(onePar,buffersize,"eta : %2.1f < m <%2.1f", fMassEtaMin,fMassEtaMax);
956  parList+=onePar ;
957  snprintf(onePar,buffersize,"conv: %2.1f < m <%2.1f", fMassPhoMin,fMassPhoMax);
958  parList+=onePar ;
959 
960 
961  return parList;
962 }
963 
964 //________________________________________________
966 //________________________________________________
967 void AliCaloPID::Print(const Option_t * opt) const
968 {
969  if(! opt)
970  return;
971 
972  printf("***** Print: %s %s ******\n", GetName(), GetTitle() ) ;
973 
975  {
976  printf("PHOS PID weight , photon %0.2f, pi0 %0.2f, e %0.2f, charge %0.2f, neutral %0.2f \n",
979  printf("EMCAL PID weight, photon %0.2f, pi0 %0.2f, e %0.2f, charge %0.2f, neutral %0.2f\n",
982 
983  printf("PHOS Parametrized weight on? = %d\n", fPHOSWeightFormula) ;
985  {
986  printf("Photon weight formula = %s\n", fPHOSPhotonWeightFormulaExpression.Data());
987  printf("Pi0 weight formula = %s\n", fPHOSPi0WeightFormulaExpression .Data());
988  }
989  if(fRecalculateBayesian) printf(" Recalculate bayesian with Particle Flux? = %d\n",fParticleFlux);
990  }
991  else
992  {
993  printf("TOF cut = %e\n", fTOFCut);
994  printf("EMCAL Lambda0 cut min = %2.2f; max = %2.2f\n", fEMCALL0CutMin,fEMCALL0CutMax);
995  printf("EMCAL cluster-track dEta < %2.3f; dPhi < %2.3f\n", fEMCALDEtaCut, fEMCALDPhiCut);
996  printf("PHOS Treac matching cut =%2.2f, Dispersion Cut =%2.2f \n",fPHOSRCut, fPHOSDispersionCut) ;
997 
998  }
999 
1000  printf("Min. N Cells =%d \n", fSplitMinNCells) ;
1001  if(fUseSimpleM02Cut) printf("%2.2f < lambda_0^2 <%2.2f \n",fSplitM02MinCut,fSplitM02MaxCut);
1002  if(fUseSimpleMassCut)printf("pi0 : %2.2f<m<%2.2f \n", fMassPi0Min,fMassPi0Max);
1003  printf("eta : %2.2f<m<%2.2f \n", fMassEtaMin,fMassEtaMax);
1004  printf("phot: %2.2f<m<%2.2f \n", fMassPhoMin,fMassPhoMax);
1005 
1006  printf(" \n");
1007 }
1008 
1009 //_________________________________________________________________
1010 // Print PID of cluster, (AliVCluster*)cluster->GetPID()
1011 //_________________________________________________________________
1013 {
1014  printf("AliCaloPID::PrintClusterPIDWeights() \n \t ph %0.2f, pi0 %0.2f, el %0.2f, conv el %0.2f, \n \t \
1015  pion %0.2f, kaon %0.2f, proton %0.2f , neutron %0.2f, kaon %0.2f \n",
1016  pid[AliVCluster::kPhoton], pid[AliVCluster::kPi0],
1017  pid[AliVCluster::kElectron], pid[AliVCluster::kEleCon],
1018  pid[AliVCluster::kPion], pid[AliVCluster::kKaon],
1019  pid[AliVCluster::kProton],
1020  pid[AliVCluster::kNeutron], pid[AliVCluster::kKaon0]);
1021 }
1022 
1023 //___________________________________________________________________________
1025 //___________________________________________________________________________
1026 void AliCaloPID::SetPIDBits(AliVCluster * cluster,
1027  AliAODPWG4Particle * ph, AliCalorimeterUtils* cu,
1028  AliVEvent* event)
1029 {
1030  // Dispersion/lambdas
1031  //Double_t disp= cluster->GetDispersion() ;
1032  Double_t l1 = cluster->GetM20() ;
1033  Double_t l0 = cluster->GetM02() ;
1034  Bool_t isDispOK = kTRUE ;
1035  if(cluster->IsPHOS()){
1036  if(TestPHOSDispersion(ph->Pt(),l0,l1) < fPHOSDispersionCut) isDispOK = kTRUE;
1037  else isDispOK = kFALSE;
1038  }
1039  else{//EMCAL
1040 
1041  if(l0 > fEMCALL0CutMin && l0 < fEMCALL0CutMax) isDispOK = kTRUE;
1042 
1043  }
1044 
1045  ph->SetDispBit(isDispOK) ;
1046 
1047  //TOF
1048  Double_t tof=cluster->GetTOF() ;
1049  ph->SetTOFBit(TMath::Abs(tof)<fTOFCut) ;
1050 
1051  //Charged
1052  Bool_t isNeutral = IsTrackMatched(cluster,cu,event);
1053 
1054  ph->SetChargedBit(isNeutral);
1055 
1056  //Set PID pdg
1057  ph->SetIdentifiedParticleType(GetIdentifiedParticleType(cluster));
1058 
1059  AliDebug(1,Form("TOF %e, Lambda0 %2.2f, Lambda1 %2.2f",tof , l0, l1));
1060  AliDebug(1,Form("pdg %d, bits: TOF %d, Dispersion %d, Charge %d",
1061  ph->GetIdentifiedParticleType(), ph->GetTOFBit() , ph->GetDispBit() , ph->GetChargedBit()));
1062 }
1063 
1064 //_________________________________________________________
1070 //_________________________________________________________
1071 Bool_t AliCaloPID::IsTrackMatched(AliVCluster* cluster,
1072  AliCalorimeterUtils * cu,
1073  AliVEvent* event) const
1074 {
1075  Int_t nMatches = cluster->GetNTracksMatched();
1076  AliVTrack * track = 0;
1077 
1078  // At least one match
1079  //
1080  if(nMatches <= 0) return kFALSE;
1081 
1082  // Select the track, depending on ESD or AODs
1083  //
1084  //In case of ESDs,
1085  //by default without match one entry with negative index, no match, reject.
1086  //
1087  if(!strcmp("AliESDCaloCluster",Form("%s",cluster->ClassName())))
1088  {
1089  Int_t iESDtrack = ((AliESDCaloCluster*)cluster)->GetTracksMatched()->At(0); //cluster->GetTrackMatchedIndex();
1090 
1091  if(iESDtrack >= 0) track = dynamic_cast<AliVTrack*> (event->GetTrack(iESDtrack));
1092  else return kFALSE;
1093 
1094  if (!track)
1095  {
1096  AliWarning(Form("Null matched track in ESD for index %d",iESDtrack));
1097  return kFALSE;
1098  }
1099  } // ESDs
1100  else
1101  { // AODs
1102  track = dynamic_cast<AliVTrack*> (cluster->GetTrackMatched(0));
1103  if (!track)
1104  {
1105  AliWarning("Null matched track in AOD!");
1106  return kFALSE;
1107  }
1108  } // AODs
1109 
1110  Float_t dZ = cluster->GetTrackDz();
1111  Float_t dR = cluster->GetTrackDx();
1112 
1113  // Comment out, new value already set in AliCalorimeterUtils::RecalculateClusterTrackMatching()
1114  // when executed in the reader.
1115  // // if track matching was recalculated
1116  // if(cluster->IsEMCAL() && cu && cu->IsRecalculationOfClusterTrackMatchingOn())
1117  // {
1118  // dR = 2000., dZ = 2000.;
1119  // cu->GetEMCALRecoUtils()->GetMatchedResiduals(cluster->GetID(),dZ,dR);
1120  // //AliDebug(2,"Residuals, (Old, New): z (%2.4f,%2.4f), x (%2.4f,%2.4f)\n", cluster->GetTrackDz(),dZ,cluster->GetTrackDx(),dR));
1121  // }
1122 
1123  if(cluster->IsPHOS())
1124  {
1125  Int_t charge = track->Charge();
1126  Double_t mf = event->GetMagneticField();
1127  if(TestPHOSChargedVeto(dR, dZ, track->Pt(), charge, mf ) < fPHOSRCut) return kTRUE;
1128  else return kFALSE;
1129 
1130  } // PHOS
1131  else // EMCAL
1132  {
1133  AliDebug(1,Form("EMCAL dR %f < %f, dZ %f < %f ",dR, fEMCALDPhiCut, dZ, fEMCALDEtaCut));
1134 
1135  if(TMath::Abs(dR) < fEMCALDPhiCut &&
1136  TMath::Abs(dZ) < fEMCALDEtaCut) return kTRUE;
1137  else return kFALSE;
1138  }// EMCAL cluster
1139 }
1140 
1141 //___________________________________________________________________________________________________
1149 //___________________________________________________________________________________________________
1151 {
1152  Double_t l2Mean = 1.53126+9.50835e+06/(1.+1.08728e+07*pt+1.73420e+06*pt*pt) ;
1153  Double_t l1Mean = 1.12365+0.123770*TMath::Exp(-pt*0.246551)+5.30000e-03*pt ;
1154  Double_t l2Sigma = 6.48260e-02+7.60261e+10/(1.+1.53012e+11*pt+5.01265e+05*pt*pt)+9.00000e-03*pt;
1155  Double_t l1Sigma = 4.44719e-04+6.99839e-01/(1.+1.22497e+00*pt+6.78604e-07*pt*pt)+9.00000e-03*pt;
1156  Double_t c =-0.35-0.550*TMath::Exp(-0.390730*pt) ;
1157  Double_t r2 = 0.5* (l1-l1Mean)*(l1-l1Mean)/l1Sigma/l1Sigma +
1158  0.5* (l2-l2Mean)*(l2-l2Mean)/l2Sigma/l2Sigma +
1159  0.5*c*(l1-l1Mean)*(l2-l2Mean)/l1Sigma/l2Sigma ;
1160 
1161  AliDebug(1,Form("PHOS SS R %f < %f?", TMath::Sqrt(r2), fPHOSDispersionCut));
1162 
1163  return TMath::Sqrt(r2) ;
1164 }
1165 
1166 //_______________________________________________________________________________________________
1181 //_______________________________________________________________________________________________
1183  Int_t charge, Double_t mf) const
1184 {
1185  Double_t meanX = 0.;
1186  Double_t meanZ = 0.;
1187  Double_t sx = TMath::Min(5.4,2.59719e+02*TMath::Exp(-pt/1.02053e-01)+
1188  6.58365e-01*5.91917e-01*5.91917e-01/((pt-9.61306e-01)*(pt-9.61306e-01)+5.91917e-01*5.91917e-01)+
1189  1.59219);
1190  Double_t sz = TMath::Min(2.75,4.90341e+02*1.91456e-02*1.91456e-02/(pt*pt+1.91456e-02*1.91456e-02)+
1191  1.60) ;
1192 
1193  if(mf<0.){ //field --
1194  meanZ = -0.468318 ;
1195  if(charge>0)
1196  meanX = TMath::Min(7.3, 3.89994*1.20679 *1.20679 /(pt*pt+1.20679*1.20679)+
1197  0.249029+2.49088e+07*TMath::Exp(-pt*3.33650e+01)) ;
1198  else
1199  meanX =-TMath::Min(7.7, 3.86040*0.912499*0.912499/(pt*pt+0.912499*0.912499)+
1200  1.23114 +4.48277e+05*TMath::Exp(-pt*2.57070e+01)) ;
1201  }
1202  else{ //Field ++
1203  meanZ = -0.468318;
1204  if(charge>0)
1205  meanX =-TMath::Min(8.0,3.86040*1.31357*1.31357/(pt*pt+1.31357*1.31357)+
1206  0.880579+7.56199e+06*TMath::Exp(-pt*3.08451e+01)) ;
1207  else
1208  meanX = TMath::Min(6.85, 3.89994*1.16240*1.16240/(pt*pt+1.16240*1.16240)-
1209  0.120787+2.20275e+05*TMath::Exp(-pt*2.40913e+01)) ;
1210  }
1211 
1212  Double_t rz = (dz-meanZ)/sz ;
1213  Double_t rx = (dx-meanX)/sx ;
1214 
1215  AliDebug(1,Form("PHOS Matching R %f < %f",TMath::Sqrt(rx*rx+rz*rz), fPHOSRCut));
1216 
1217  return TMath::Sqrt(rx*rx+rz*rz) ;
1218 }
1219 
Int_t charge
Int_t pdg
TString fPHOSPhotonWeightFormulaExpression
Photon weight formula in string.
Definition: AliCaloPID.h:298
Bool_t IsInEtaM02Range(Float_t energy, Float_t m02, Int_t nlm) const
Definition: AliCaloPID.cxx:460
Float_t fTOFCut
Cut on TOF, used in PID evaluation.
Definition: AliCaloPID.h:308
Float_t fMassPhoMin
Min Photon mass.
Definition: AliCaloPID.h:326
TFormula * GetPHOSPi0WeightFormula()
Definition: AliCaloPID.h:164
double Double_t
Definition: External.C:58
Float_t fPHOSNeutralWeight
Bayesian PID weight for neutral hadrons in PHOS.
Definition: AliCaloPID.h:293
Bool_t fUseBayesianWeights
Select clusters based on weights calculated in reconstruction.
Definition: AliCaloPID.h:281
Bool_t IsInPi0SplitAsymmetryRange(Float_t energy, Float_t asy, Int_t nlm) const
Definition: AliCaloPID.cxx:321
Float_t fMassPhoMax
Min Photon mass.
Definition: AliCaloPID.h:327
Float_t fMassPi0Max
Min Pi0 mass, simple cut case.
Definition: AliCaloPID.h:325
Float_t fEMCALElectronWeight
Bayesian PID weight for electrons in EMCAL.
Definition: AliCaloPID.h:286
virtual ~AliCaloPID()
Definition: AliCaloPID.cxx:154
AliEMCALRecoUtils * GetEMCALRecoUtils() const
Float_t fPHOSPhotonWeight
Bayesian PID weight for photons in PHOS.
Definition: AliCaloPID.h:289
Float_t fSplitEFracMin[3]
Definition: AliCaloPID.h:334
Float_t fM02MaxParamShiftNLMN
shift of max M02 for NLM>2.
Definition: AliCaloPID.h:332
Float_t fWidthPi0Param[2][6]
Width param, 2 regions in energy.
Definition: AliCaloPID.h:329
Double_t mass
Bool_t fUseSplitAsyCut
Remove splitted clusters with too large asymmetry.
Definition: AliCaloPID.h:317
Float_t fEMCALPi0Weight
Bayesian PID weight for pi0 in EMCAL.
Definition: AliCaloPID.h:285
Float_t fEMCALPhotonWeight
Bayesian PID weight for photons in EMCAL.
Definition: AliCaloPID.h:284
Float_t fSubClusterEMin[3]
Do not use sub-clusters with too low energy depeding on NLM.
Definition: AliCaloPID.h:336
Float_t fMassEtaMax
Max Eta mass.
Definition: AliCaloPID.h:323
TCanvas * c
Definition: TestFitELoss.C:172
TString fPHOSPi0WeightFormulaExpression
Pi0 weight formula in string.
Definition: AliCaloPID.h:299
TFormula * GetPHOSPhotonWeightFormula()
Definition: AliCaloPID.h:158
Float_t fEMCALDEtaCut
Track matching cut on Dz.
Definition: AliCaloPID.h:305
const TString calorimeter
Definition: anaM.C:35
Float_t fMassPi0Min
Min Pi0 mass, simple cut case.
Definition: AliCaloPID.h:324
Float_t fPHOSChargeWeight
Bayesian PID weight for charged hadrons in PHOS.
Definition: AliCaloPID.h:292
Float_t fMassPi0Param[2][6]
Mean mass param, 2 regions in energy.
Definition: AliCaloPID.h:328
TString GetPIDParametersList()
Put data member values in string to keep in output container.
Definition: AliCaloPID.cxx:892
void Print(const Option_t *opt) const
Print some relevant parameters set for the analysis.
Definition: AliCaloPID.cxx:967
Float_t fM02MaxParam[2][5]
5 param for expo + pol fit on M02 maximum for pi0 selection.
Definition: AliCaloPID.h:331
Bool_t fRecalculateBayesian
Recalculate PID bayesian or use simple PID?
Definition: AliCaloPID.h:282
AliEMCALGeometry * GetEMCALGeometry() const
int Int_t
Definition: External.C:63
Float_t fSplitWidthSigma
Cut on mass+-width*fSplitWidthSigma.
Definition: AliCaloPID.h:337
Int_t fDebug
Debug level.
Definition: AliCaloPID.h:275
float Float_t
Definition: External.C:68
Bool_t IsInPi0M02Range(Float_t energy, Float_t m02, Int_t nlm) const
Definition: AliCaloPID.cxx:412
Float_t fSplitM02MaxCut
Study clusters with l0 smaller than cut.
Definition: AliCaloPID.h:319
Float_t TestPHOSDispersion(Double_t pt, Double_t m20, Double_t m02) const
Bool_t fPHOSWeightFormula
Use parametrized weight threshold, function of energy.
Definition: AliCaloPID.h:295
Float_t TestPHOSChargedVeto(Double_t dx, Double_t dz, Double_t ptTrack, Int_t chargeTrack, Double_t mf) const
void SetPIDBits(AliVCluster *cluster, AliAODPWG4Particle *aodph, AliCalorimeterUtils *cu, AliVEvent *event)
Set Bits for PID selection.
Float_t fEMCALNeutralWeight
Bayesian PID weight for neutral hadrons in EMCAL.
Definition: AliCaloPID.h:288
Int_t GetNumberOfLocalMaxima(AliVCluster *cluster, AliVCaloCells *cells)
Find the number of local maxima in cluster.
Float_t fPHOSRCut
Track-Cluster distance cut for track matching in PHOS.
Definition: AliCaloPID.h:311
Int_t GetIdentifiedParticleType(AliVCluster *cluster)
Definition: AliCaloPID.cxx:554
Bool_t Data(TH1F *h, Double_t *rangefit, Bool_t writefit, Double_t &sgn, Double_t &errsgn, Double_t &bkg, Double_t &errbkg, Double_t &sgnf, Double_t &errsgnf, Double_t &sigmafit, Int_t &status)
Bool_t fUseSplitSSCut
Remove splitted clusters out of shower shape band.
Definition: AliCaloPID.h:318
energy
Float_t fMassShiftHighECell
Shift cuts 5 MeV for Ecell > 150 MeV, default Ecell > 50 MeV.
Definition: AliCaloPID.h:338
Bool_t IsInPi0SplitMassRange(Float_t energy, Float_t mass, Int_t nlm) const
Definition: AliCaloPID.cxx:352
Float_t fPHOSElectronWeight
Bayesian PID weight for electrons in PHOS.
Definition: AliCaloPID.h:291
Float_t fEMCALDPhiCut
Track matching cut on Dx.
Definition: AliCaloPID.h:306
Float_t fEMCALL0CutMax
Max Cut on shower shape lambda0, used in PID evaluation, only EMCAL.
Definition: AliCaloPID.h:303
ClassImp(AliAnalysisTaskCRC) AliAnalysisTaskCRC
Int_t GetIdentifiedParticleTypeFromClusterSplitting(AliVCluster *cluster, AliVCaloCells *cells, AliCalorimeterUtils *caloutils, Double_t vertex[3], Int_t &nLocMax, Double_t &mass, Double_t &angle, TLorentzVector &l1, TLorentzVector &l2, Int_t &absId1, Int_t &absId2, Float_t &distbad1, Float_t &distbad2, Bool_t &fidcut1, Bool_t &fidcut2) const
Definition: AliCaloPID.cxx:698
TFormula * fPHOSPhotonWeightFormula
Formula for photon weight.
Definition: AliCaloPID.h:296
Double_t minMass
Bool_t IsInConM02Range(Float_t energy, Float_t m02, Int_t nlm) const
Definition: AliCaloPID.cxx:515
Float_t fPHOSDispersionCut
Shower shape elipse radious cut.
Definition: AliCaloPID.h:310
Int_t GetIdentifiedParticleTypeFromBayesWeights(Bool_t isEMCAL, Double_t *pid, Float_t energy)
Definition: AliCaloPID.cxx:611
void RecalibrateCellAmplitude(Float_t &amp, Int_t calo, Int_t absId) const
Recalculate cell energy if recalibration factor.
const char Option_t
Definition: External.C:48
Bool_t fUseSimpleM02Cut
Use simple min-max M02 cut.
Definition: AliCaloPID.h:316
Class for PID selection with calorimeters.
Definition: AliCaloPID.h:51
Float_t fEMCALL0CutMin
Min Cut on shower shape lambda0, used in PID evaluation, only EMCAL.
Definition: AliCaloPID.h:304
AliEMCALPIDUtils * GetEMCALPIDUtils()
Definition: AliCaloPID.cxx:542
Double_t maxMass
Float_t fMassEtaMin
Min Eta mass.
Definition: AliCaloPID.h:322
bool Bool_t
Definition: External.C:53
Class with utils specific to calorimeter clusters/cells.
Float_t fEMCALChargeWeight
Bayesian PID weight for charged hadrons in EMCAL.
Definition: AliCaloPID.h:287
AliEMCALPIDUtils * fEMCALPIDUtils
Pointer to EMCALPID to redo the PID Bayesian calculation.
Definition: AliCaloPID.h:280
void InitParameters()
Definition: AliCaloPID.cxx:164
Int_t fParticleFlux
Particle flux for setting PID parameters.
Definition: AliCaloPID.h:276
Float_t fM02MinParam[2][5]
5 param for expo + pol fit on M02 minimum for pi0 selection (maximum for conversions).
Definition: AliCaloPID.h:330
void PrintClusterPIDWeights(const Double_t *pid) const
Bool_t fUseSimpleMassCut
Use simple min-max pi0 mass cut.
Definition: AliCaloPID.h:315
Float_t fPHOSPi0Weight
Bayesian PID weight for pi0 in PHOS.
Definition: AliCaloPID.h:290
TFormula * fPHOSPi0WeightFormula
Formula for pi0 weight.
Definition: AliCaloPID.h:297
Bool_t IsInM02Range(Float_t m02) const
Definition: AliCaloPID.cxx:396
Int_t fSplitMinNCells
Study clusters with ncells larger than cut.
Definition: AliCaloPID.h:321
void SplitEnergy(Int_t absId1, Int_t absId2, AliVCluster *cluster, AliVCaloCells *cells, AliAODCaloCluster *cluster1, AliAODCaloCluster *cluster2, Int_t nMax, Int_t eventNumber=0)
Bool_t IsTrackMatched(AliVCluster *cluster, AliCalorimeterUtils *cu, AliVEvent *event) const
Float_t fSplitM02MinCut
Study clusters with l0 larger than cut, simple case.
Definition: AliCaloPID.h:320
Float_t fAsyMinParam[2][4]
4 param for fit on asymmetry minimum, for 2 cases, NLM=1 and NLM>=2.
Definition: AliCaloPID.h:333