]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/STEERBase/AliPIDResponse.cxx
Revert wrong commit
[u/mrichter/AliRoot.git] / STEER / STEERBase / AliPIDResponse.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id: AliPIDResponse.cxx 46193 2010-12-21 09:00:14Z wiechula $ */
17
18 //-----------------------------------------------------------------
19 //        Base class for handling the pid response               //
20 //        functions of all detectors                             //
21 //        and give access to the nsigmas                         //
22 //                                                               //
23 //   Origin: Jens Wiechula, Uni Tuebingen, jens.wiechula@cern.ch //
24 //-----------------------------------------------------------------
25
26 #include <TList.h>
27 #include <TObjArray.h>
28 #include <TPRegexp.h>
29 #include <TF1.h>
30 #include <TH2D.h>
31 #include <TSpline.h>
32 #include <TFile.h>
33 #include <TArrayI.h>
34 #include <TArrayF.h>
35 #include <TLinearFitter.h>
36
37 #include <AliVEvent.h>
38 #include <AliVTrack.h>
39 #include <AliLog.h>
40 #include <AliPID.h>
41 #include <AliOADBContainer.h>
42 #include <AliTRDPIDResponseObject.h>
43 #include <AliTOFPIDParams.h>
44 #include <AliHMPIDPIDParams.h>
45
46 #include "AliPIDResponse.h"
47 #include "AliDetectorPID.h"
48
49 #include "AliCentrality.h"
50
51 ClassImp(AliPIDResponse);
52
53 AliPIDResponse::AliPIDResponse(Bool_t isMC/*=kFALSE*/) :
54 TNamed("PIDResponse","PIDResponse"),
55 fITSResponse(isMC),
56 fTPCResponse(),
57 fTRDResponse(),
58 fTOFResponse(),
59 fHMPIDResponse(),
60 fEMCALResponse(),
61 fRange(5.),
62 fITSPIDmethod(kITSTruncMean),
63 fTuneMConData(kFALSE),
64 fTuneMConDataMask(kDetTOF|kDetTPC),
65 fIsMC(isMC),
66 fCachePID(kTRUE),
67 fOADBPath(),
68 fCustomTPCpidResponse(),
69 fBeamType("PP"),
70 fLHCperiod(),
71 fMCperiodTPC(),
72 fMCperiodUser(),
73 fCurrentFile(),
74 fCurrentAliRootRev(-1),
75 fRecoPass(0),
76 fRecoPassUser(-1),
77 fRun(-1),
78 fOldRun(-1),
79 fResT0A(75.),
80 fResT0C(65.),
81 fResT0AC(55.),
82 fArrPidResponseMaster(NULL),
83 fResolutionCorrection(NULL),
84 fOADBvoltageMaps(NULL),
85 fUseTPCEtaCorrection(kFALSE),
86 fUseTPCMultiplicityCorrection(kFALSE),
87 fTRDPIDResponseObject(NULL),
88 fTOFtail(0.9),
89 fTOFPIDParams(NULL),
90 fHMPIDPIDParams(NULL),
91 fEMCALPIDParams(NULL),
92 fCurrentEvent(NULL),
93 fCurrCentrality(0.0)
94 {
95   //
96   // default ctor
97   //
98   AliLog::SetClassDebugLevel("AliPIDResponse",0);
99   AliLog::SetClassDebugLevel("AliESDpid",0);
100   AliLog::SetClassDebugLevel("AliAODpidUtil",0);
101
102 }
103
104 //______________________________________________________________________________
105 AliPIDResponse::~AliPIDResponse()
106 {
107   //
108   // dtor
109   //
110   delete fArrPidResponseMaster;
111   delete fTRDPIDResponseObject;
112   delete fTOFPIDParams;
113 }
114
115 //______________________________________________________________________________
116 AliPIDResponse::AliPIDResponse(const AliPIDResponse &other) :
117 TNamed(other),
118 fITSResponse(other.fITSResponse),
119 fTPCResponse(other.fTPCResponse),
120 fTRDResponse(other.fTRDResponse),
121 fTOFResponse(other.fTOFResponse),
122 fHMPIDResponse(other.fHMPIDResponse),
123 fEMCALResponse(other.fEMCALResponse),
124 fRange(other.fRange),
125 fITSPIDmethod(other.fITSPIDmethod),
126 fTuneMConData(other.fTuneMConData),
127 fTuneMConDataMask(other.fTuneMConDataMask),
128 fIsMC(other.fIsMC),
129 fCachePID(other.fCachePID),
130 fOADBPath(other.fOADBPath),
131 fCustomTPCpidResponse(other.fCustomTPCpidResponse),
132 fBeamType("PP"),
133 fLHCperiod(),
134 fMCperiodTPC(),
135 fMCperiodUser(other.fMCperiodUser),
136 fCurrentFile(),
137 fCurrentAliRootRev(other.fCurrentAliRootRev),
138 fRecoPass(0),
139 fRecoPassUser(other.fRecoPassUser),
140 fRun(-1),
141 fOldRun(-1),
142 fResT0A(75.),
143 fResT0C(65.),
144 fResT0AC(55.),
145 fArrPidResponseMaster(NULL),
146 fResolutionCorrection(NULL),
147 fOADBvoltageMaps(NULL),
148 fUseTPCEtaCorrection(other.fUseTPCEtaCorrection),
149 fUseTPCMultiplicityCorrection(other.fUseTPCMultiplicityCorrection),
150 fTRDPIDResponseObject(NULL),
151 fTOFtail(0.9),
152 fTOFPIDParams(NULL),
153 fHMPIDPIDParams(NULL),
154 fEMCALPIDParams(NULL),
155 fCurrentEvent(NULL),
156 fCurrCentrality(0.0)
157 {
158   //
159   // copy ctor
160   //
161 }
162
163 //______________________________________________________________________________
164 AliPIDResponse& AliPIDResponse::operator=(const AliPIDResponse &other)
165 {
166   //
167   // copy ctor
168   //
169   if(this!=&other) {
170     delete fArrPidResponseMaster;
171     TNamed::operator=(other);
172     fITSResponse=other.fITSResponse;
173     fTPCResponse=other.fTPCResponse;
174     fTRDResponse=other.fTRDResponse;
175     fTOFResponse=other.fTOFResponse;
176     fHMPIDResponse=other.fHMPIDResponse;
177     fEMCALResponse=other.fEMCALResponse;
178     fRange=other.fRange;
179     fITSPIDmethod=other.fITSPIDmethod;
180     fOADBPath=other.fOADBPath;
181     fCustomTPCpidResponse=other.fCustomTPCpidResponse;
182     fTuneMConData=other.fTuneMConData;
183     fTuneMConDataMask=other.fTuneMConDataMask;
184     fIsMC=other.fIsMC;
185     fCachePID=other.fCachePID;
186     fBeamType="PP";
187     fLHCperiod="";
188     fMCperiodTPC="";
189     fMCperiodUser=other.fMCperiodUser;
190     fCurrentFile="";
191     fCurrentAliRootRev=other.fCurrentAliRootRev;
192     fRecoPass=0;
193     fRecoPassUser=other.fRecoPassUser;
194     fRun=-1;
195     fOldRun=-1;
196     fResT0A=75.;
197     fResT0C=65.;
198     fResT0AC=55.;
199     fArrPidResponseMaster=NULL;
200     fResolutionCorrection=NULL;
201     fOADBvoltageMaps=NULL;
202     fUseTPCEtaCorrection=other.fUseTPCEtaCorrection;
203     fUseTPCMultiplicityCorrection=other.fUseTPCMultiplicityCorrection;
204     fTRDPIDResponseObject=NULL;
205     fEMCALPIDParams=NULL;
206     fTOFtail=0.9;
207     fTOFPIDParams=NULL;
208     fHMPIDPIDParams=NULL;
209     fCurrentEvent=other.fCurrentEvent;
210
211   }
212   return *this;
213 }
214
215 //______________________________________________________________________________
216 Float_t AliPIDResponse::NumberOfSigmas(EDetector detector, const AliVParticle *vtrack, AliPID::EParticleType type) const
217 {
218   //
219   // NumberOfSigmas for 'detCode'
220   //
221   
222   const AliVTrack *track=static_cast<const AliVTrack*>(vtrack);
223   // look for cached value first
224   const AliDetectorPID *detPID=track->GetDetectorPID();
225   
226   if ( detPID && detPID->HasNumberOfSigmas(detector)){
227     return detPID->GetNumberOfSigmas(detector, type);
228   } else if (fCachePID) {
229     FillTrackDetectorPID(track, detector);
230     detPID=track->GetDetectorPID();
231     return detPID->GetNumberOfSigmas(detector, type);
232   }
233   
234   return GetNumberOfSigmas(detector, track, type);
235 }
236
237 //______________________________________________________________________________
238 AliPIDResponse::EDetPidStatus AliPIDResponse::NumberOfSigmas(EDetector detCode, const AliVParticle *track,
239                                                              AliPID::EParticleType type, Double_t &val) const
240 {
241   //
242   // NumberOfSigmas with detector status as return value
243   //
244   
245   val=NumberOfSigmas(detCode, track, type);
246   return CheckPIDStatus(detCode, (AliVTrack*)track);
247 }
248
249 //______________________________________________________________________________
250 // public buffered versions of the PID calculation
251 //
252
253 //______________________________________________________________________________
254 Float_t AliPIDResponse::NumberOfSigmasITS(const AliVParticle *vtrack, AliPID::EParticleType type) const
255 {
256   //
257   // Calculate the number of sigmas in the ITS
258   //
259   
260   return NumberOfSigmas(kITS, vtrack, type);
261 }
262
263 //______________________________________________________________________________
264 Float_t AliPIDResponse::NumberOfSigmasTPC(const AliVParticle *vtrack, AliPID::EParticleType type) const
265 {
266   //
267   // Calculate the number of sigmas in the TPC
268   //
269   
270   return NumberOfSigmas(kTPC, vtrack, type);
271 }
272
273 //______________________________________________________________________________
274 Float_t AliPIDResponse::NumberOfSigmasTPC( const AliVParticle *vtrack, 
275                                            AliPID::EParticleType type,
276                                            AliTPCPIDResponse::ETPCdEdxSource dedxSource) const
277 {
278   //get number of sigmas according the selected TPC gain configuration scenario
279   const AliVTrack *track=static_cast<const AliVTrack*>(vtrack);
280
281   Float_t nSigma=fTPCResponse.GetNumberOfSigmas(track, type, dedxSource, fUseTPCEtaCorrection, fUseTPCMultiplicityCorrection);
282
283   return nSigma;
284 }
285
286 //______________________________________________________________________________
287 Float_t AliPIDResponse::NumberOfSigmasTOF(const AliVParticle *vtrack, AliPID::EParticleType type) const
288 {
289   //
290   // Calculate the number of sigmas in the TOF
291   //
292   
293   return NumberOfSigmas(kTOF, vtrack, type);
294 }
295
296 //______________________________________________________________________________
297 Float_t AliPIDResponse::NumberOfSigmasHMPID(const AliVParticle *vtrack, AliPID::EParticleType type) const
298 {
299   //
300   // Calculate the number of sigmas in the EMCAL
301   //
302   
303   return NumberOfSigmas(kHMPID, vtrack, type);
304 }
305
306 //______________________________________________________________________________
307 Float_t AliPIDResponse::NumberOfSigmasEMCAL(const AliVParticle *vtrack, AliPID::EParticleType type) const
308 {
309   //
310   // Calculate the number of sigmas in the EMCAL
311   //
312   
313   return NumberOfSigmas(kEMCAL, vtrack, type);
314 }
315
316 //______________________________________________________________________________
317 Float_t  AliPIDResponse::NumberOfSigmasEMCAL(const AliVParticle *vtrack, AliPID::EParticleType type, Double_t &eop, Double_t showershape[4])  const
318 {
319   //
320   // emcal nsigma with eop and showershape
321   //
322   AliVTrack *track=(AliVTrack*)vtrack;
323   
324   AliVCluster *matchedClus = NULL;
325
326   Double_t mom     = -1.; 
327   Double_t pt      = -1.; 
328   Double_t EovP    = -1.;
329   Double_t fClsE   = -1.;
330
331   // initialize eop and shower shape parameters
332   eop = -1.;
333   for(Int_t i = 0; i < 4; i++){
334     showershape[i] = -1.;
335   }
336   
337   Int_t nMatchClus = -1;
338   Int_t charge     = 0;
339   
340   // Track matching
341   nMatchClus = track->GetEMCALcluster();
342   if(nMatchClus > -1){
343
344     mom    = track->P();
345     pt     = track->Pt();
346     charge = track->Charge();
347     
348     matchedClus = (AliVCluster*)fCurrentEvent->GetCaloCluster(nMatchClus);
349     
350     if(matchedClus){
351       
352       // matched cluster is EMCAL
353       if(matchedClus->IsEMCAL()){
354         
355         fClsE       = matchedClus->E();
356         EovP        = fClsE/mom;
357         
358         // fill used EMCAL variables here
359         eop            = EovP; // E/p
360         showershape[0] = matchedClus->GetNCells(); // number of cells in cluster
361         showershape[1] = matchedClus->GetM02(); // long axis
362         showershape[2] = matchedClus->GetM20(); // short axis
363         showershape[3] = matchedClus->GetDispersion(); // dispersion
364
365         // look for cached value first
366         const AliDetectorPID *detPID=track->GetDetectorPID();
367         const EDetector detector=kEMCAL;
368         
369         if ( detPID && detPID->HasNumberOfSigmas(detector)){
370           return detPID->GetNumberOfSigmas(detector, type);
371         } else if (fCachePID) {
372           FillTrackDetectorPID(track, detector);
373           detPID=track->GetDetectorPID();
374           return detPID->GetNumberOfSigmas(detector, type);
375         }
376         
377         // NSigma value really meaningful only for electrons!
378         return fEMCALResponse.GetNumberOfSigmas(pt,EovP,type,charge);
379       }
380     }
381   }
382   return -999;
383 }
384
385 //______________________________________________________________________________
386 AliPIDResponse::EDetPidStatus AliPIDResponse::GetSignalDelta(EDetector detector, const AliVParticle *track, AliPID::EParticleType type, Double_t &val, Bool_t ratio/*=kFALSE*/) const
387 {
388   //
389   //
390   //
391   val=-9999.;
392   switch (detector){
393     case kITS:   return GetSignalDeltaITS(track,type,val,ratio); break;
394     case kTPC:   return GetSignalDeltaTPC(track,type,val,ratio); break;
395     case kTOF:   return GetSignalDeltaTOF(track,type,val,ratio); break;
396     case kHMPID: return GetSignalDeltaHMPID(track,type,val,ratio); break;
397     default: return kDetNoSignal;
398   }
399   return kDetNoSignal;
400 }
401
402 //______________________________________________________________________________
403 Double_t AliPIDResponse::GetSignalDelta(EDetector detCode, const AliVParticle *track, AliPID::EParticleType type, Bool_t ratio/*=kFALSE*/) const
404 {
405   //
406   //
407   //
408   Double_t val=-9999.;
409   EDetPidStatus stat=GetSignalDelta(detCode, track, type, val, ratio);
410   if ( stat==kDetNoSignal ) val=-9999.;
411   return val;
412 }
413
414 //______________________________________________________________________________
415 AliPIDResponse::EDetPidStatus AliPIDResponse::ComputePIDProbability  (EDetCode  detCode, const AliVTrack *track, Int_t nSpecies, Double_t p[]) const
416 {
417   // Compute PID response of 'detCode'
418   
419   // find detector code from detector bit mask
420   Int_t detector=-1;
421   for (Int_t idet=0; idet<kNdetectors; ++idet) if ( (detCode&(1<<idet)) ) { detector=idet; break; }
422   if (detector==-1) return kDetNoSignal;
423
424   return ComputePIDProbability((EDetector)detector, track, nSpecies, p);
425 }
426
427 //______________________________________________________________________________
428 AliPIDResponse::EDetPidStatus AliPIDResponse::ComputePIDProbability  (EDetector detector,  const AliVTrack *track, Int_t nSpecies, Double_t p[]) const
429 {
430   //
431   // Compute PID response of 'detector'
432   //
433
434   const AliDetectorPID *detPID=track->GetDetectorPID();
435
436   if ( detPID && detPID->HasRawProbability(detector)){
437     return detPID->GetRawProbability(detector, p, nSpecies);
438   } else if (fCachePID) {
439     FillTrackDetectorPID(track, detector);
440     detPID=track->GetDetectorPID();
441     return detPID->GetRawProbability(detector, p, nSpecies);
442   }
443   
444   //if no caching return values calculated from scratch
445   return GetComputePIDProbability(detector, track, nSpecies, p);
446 }
447
448 //______________________________________________________________________________
449 AliPIDResponse::EDetPidStatus AliPIDResponse::ComputeITSProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[]) const
450 {
451   // Compute PID response for the ITS
452   return ComputePIDProbability(kITS, track, nSpecies, p);
453 }
454
455 //______________________________________________________________________________
456 AliPIDResponse::EDetPidStatus AliPIDResponse::ComputeTPCProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[]) const
457 {
458   // Compute PID response for the TPC
459   return ComputePIDProbability(kTPC, track, nSpecies, p);
460 }
461
462 //______________________________________________________________________________
463 AliPIDResponse::EDetPidStatus AliPIDResponse::ComputeTOFProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[]) const
464 {
465   // Compute PID response for the
466   return ComputePIDProbability(kTOF, track, nSpecies, p);
467 }
468
469 //______________________________________________________________________________
470 AliPIDResponse::EDetPidStatus AliPIDResponse::ComputeTRDProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[]) const
471 {
472   // Compute PID response for the
473   return ComputePIDProbability(kTRD, track, nSpecies, p);
474 }
475
476 //______________________________________________________________________________
477 AliPIDResponse::EDetPidStatus AliPIDResponse::ComputeEMCALProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[]) const
478 {
479   // Compute PID response for the EMCAL
480   return ComputePIDProbability(kEMCAL, track, nSpecies, p);
481 }
482 //______________________________________________________________________________
483 AliPIDResponse::EDetPidStatus AliPIDResponse::ComputePHOSProbability (const AliVTrack */*track*/, Int_t nSpecies, Double_t p[]) const
484 {
485   // Compute PID response for the PHOS
486   
487   // set flat distribution (no decision)
488   for (Int_t j=0; j<nSpecies; j++) p[j]=1./nSpecies;
489   return kDetNoSignal;
490 }
491
492 //______________________________________________________________________________
493 AliPIDResponse::EDetPidStatus AliPIDResponse::ComputeHMPIDProbability(const AliVTrack *track, Int_t nSpecies, Double_t p[]) const
494 {
495   // Compute PID response for the HMPID
496   return ComputePIDProbability(kHMPID, track, nSpecies, p);
497 }
498
499 //______________________________________________________________________________
500 AliPIDResponse::EDetPidStatus AliPIDResponse::ComputeTRDProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[],AliTRDPIDResponse::ETRDPIDMethod PIDmethod) const
501 {
502   // Compute PID response for the
503   return GetComputeTRDProbability(track, nSpecies, p, PIDmethod);
504 }
505
506 //______________________________________________________________________________
507 AliPIDResponse::EDetPidStatus AliPIDResponse::CheckPIDStatus(EDetector detector, const AliVTrack *track) const
508 {
509   // calculate detector pid status
510   
511   const Int_t iDetCode=(Int_t)detector;
512   if (iDetCode<0||iDetCode>=kNdetectors) return kDetNoSignal;
513   const AliDetectorPID *detPID=track->GetDetectorPID();
514   
515   if ( detPID ){
516     return detPID->GetPIDStatus(detector);
517   } else if (fCachePID) {
518     FillTrackDetectorPID(track, detector);
519     detPID=track->GetDetectorPID();
520     return detPID->GetPIDStatus(detector);
521   }
522   
523   // if not buffered and no buffering is requested
524   return GetPIDStatus(detector, track);
525 }
526
527 //______________________________________________________________________________
528 void AliPIDResponse::InitialiseEvent(AliVEvent *event, Int_t pass, Int_t run)
529 {
530   //
531   // Apply settings for the current event
532   //
533   fRecoPass=pass;
534   
535
536   fCurrentEvent=NULL;
537   if (!event) return;
538   fCurrentEvent=event;
539   if (run>0) fRun=run;
540   else fRun=event->GetRunNumber();
541   
542   if (fRun!=fOldRun){
543     ExecNewRun();
544     fOldRun=fRun;
545   }
546   
547   //TPC resolution parametrisation PbPb
548   if ( fResolutionCorrection ){
549     Double_t corrSigma=fResolutionCorrection->Eval(GetTPCMultiplicityBin(event));
550     fTPCResponse.SetSigma(3.79301e-03*corrSigma, 2.21280e+04);
551   }
552   
553   // Set up TPC multiplicity for PbPb
554   //TODO Will NOT give the desired number for AODs -> Needs new variable/function in future.
555   // Fatal, if AOD event and correction enabled
556   //printf("DETECTED class: %s (%d)\n\n\n\n", event->IsA()->GetName(), fUseTPCMultiplicityCorrection);//TODO
557   if (fUseTPCMultiplicityCorrection && strcmp(event->IsA()->GetName(), "AliESDEvent") != 0) {
558     AliFatal("TPC multiplicity correction is enabled, but will NOT work for AOD events, only for ESD => Disabled multiplicity correction!");
559     fUseTPCMultiplicityCorrection = kFALSE;
560   }
561   
562   if (fUseTPCMultiplicityCorrection)
563     fTPCResponse.SetCurrentEventMultiplicity(event->GetNumberOfTracks());
564   else
565     fTPCResponse.SetCurrentEventMultiplicity(0);
566   
567   //TOF resolution
568   SetTOFResponse(event, (AliPIDResponse::EStartTimeType_t)fTOFPIDParams->GetStartTimeMethod());
569
570
571   // Get and set centrality
572   AliCentrality *centrality = event->GetCentrality();
573   if(centrality){
574     fCurrCentrality = centrality->GetCentralityPercentile("V0M");
575   }
576   else{
577     fCurrCentrality = -1;
578   }
579
580   // Set centrality percentile for EMCAL
581   fEMCALResponse.SetCentrality(fCurrCentrality);
582
583 }
584
585 //______________________________________________________________________________
586 void AliPIDResponse::ExecNewRun()
587 {
588   //
589   // Things to Execute upon a new run
590   //
591   SetRecoInfo();
592   
593   SetITSParametrisation();
594   
595   SetTPCPidResponseMaster();
596   SetTPCParametrisation();
597   SetTPCEtaMaps();
598
599   SetTRDPidResponseMaster(); 
600   InitializeTRDResponse();
601
602   SetEMCALPidResponseMaster(); 
603   InitializeEMCALResponse();
604   
605   SetTOFPidResponseMaster();
606   InitializeTOFResponse();
607
608   SetHMPIDPidResponseMaster();
609   InitializeHMPIDResponse();
610
611   if (fCurrentEvent) fTPCResponse.SetMagField(fCurrentEvent->GetMagneticField());
612 }
613
614 //______________________________________________________________________________
615 Double_t AliPIDResponse::GetTPCMultiplicityBin(const AliVEvent * const event)
616 {
617   //
618   // Get TPC multiplicity in bins of 150
619   //
620   
621   const AliVVertex* vertexTPC = event->GetPrimaryVertex();
622   Double_t tpcMulti=0.;
623   if(vertexTPC){
624     Double_t vertexContribTPC=vertexTPC->GetNContributors();
625     tpcMulti=vertexContribTPC/150.;
626     if (tpcMulti>20.) tpcMulti=20.;
627   }
628   
629   return tpcMulti;
630 }
631
632 //______________________________________________________________________________
633 void AliPIDResponse::SetRecoInfo()
634 {
635   //
636   // Set reconstruction information
637   //
638   
639   //reset information
640   fLHCperiod="";
641   fMCperiodTPC="";
642   
643   fBeamType="";
644     
645   fBeamType="PP";
646
647   Bool_t hasProdInfo=(fCurrentFile.BeginsWith("LHC"));
648   
649   TPRegexp reg(".*(LHC1[1-3][a-z]+[0-9]+[a-z_]*)/.*");
650   if (hasProdInfo) reg=TPRegexp("LHC1[1-2][a-z]+[0-9]+[a-z_]*");
651   TPRegexp reg12a17("LHC1[2-3][a-z]");
652
653   //find the period by run number (UGLY, but not stored in ESD and AOD... )
654   if (fRun>=114737&&fRun<=117223)      { fLHCperiod="LHC10B"; fMCperiodTPC="LHC10D1";  }
655   else if (fRun>=118503&&fRun<=121040) { fLHCperiod="LHC10C"; fMCperiodTPC="LHC10D1";  }
656   else if (fRun>=122195&&fRun<=126437) { fLHCperiod="LHC10D"; fMCperiodTPC="LHC10F6A"; }
657   else if (fRun>=127710&&fRun<=130850) { fLHCperiod="LHC10E"; fMCperiodTPC="LHC10F6A"; }
658   else if (fRun>=133004&&fRun<=135029) { fLHCperiod="LHC10F"; fMCperiodTPC="LHC10F6A"; }
659   else if (fRun>=135654&&fRun<=136377) { fLHCperiod="LHC10G"; fMCperiodTPC="LHC10F6A"; }
660   else if (fRun>=136851&&fRun<=139846) {
661     fLHCperiod="LHC10H";
662     fMCperiodTPC="LHC10H8";
663     if (reg.MatchB(fCurrentFile)) fMCperiodTPC="LHC11A10";
664     fBeamType="PBPB";
665   }
666   else if (fRun>=139847&&fRun<=146974) { fLHCperiod="LHC11A"; fMCperiodTPC="LHC10F6A"; }
667   //TODO: periods 11B (146975-150721), 11C (150722-155837) are not yet treated assume 11d for the moment
668   else if (fRun>=146975&&fRun<=155837) { fLHCperiod="LHC11D"; fMCperiodTPC="LHC10F6A"; }
669   else if (fRun>=155838&&fRun<=159649) { fLHCperiod="LHC11D"; fMCperiodTPC="LHC10F6A"; }
670   // also for 11e (159650-162750),f(162751-165771) use 11d
671   else if (fRun>=159650&&fRun<=162750) { fLHCperiod="LHC11D"; fMCperiodTPC="LHC10F6A"; }
672   else if (fRun>=162751&&fRun<=165771) { fLHCperiod="LHC11D"; fMCperiodTPC="LHC10F6A"; }
673   
674   else if (fRun>=165772 && fRun<=170718) {
675     fLHCperiod="LHC11H";
676     fMCperiodTPC="LHC11A10";
677     fBeamType="PBPB";
678     if (reg12a17.MatchB(fCurrentFile)) fMCperiodTPC="LHC12A17";
679   }
680   if (fRun>=170719 && fRun<=177311) { fLHCperiod="LHC12A"; fBeamType="PP"; /*fMCperiodTPC="";*/ }
681   // for the moment use LHC12b parameters up to LHC12e
682   if (fRun>=177312 /*&& fRun<=179356*/) { fLHCperiod="LHC12B"; fBeamType="PP"; /*fMCperiodTPC="";*/ }
683 //   if (fRun>=179357 && fRun<=183173) { fLHCperiod="LHC12C"; fBeamType="PP"; /*fMCperiodTPC="";*/ }
684 //   if (fRun>=183174 && fRun<=186345) { fLHCperiod="LHC12D"; fBeamType="PP"; /*fMCperiodTPC="";*/ }
685 //   if (fRun>=186346 && fRun<=186635) { fLHCperiod="LHC12E"; fBeamType="PP"; /*fMCperiodTPC="";*/ }
686
687 //   if (fRun>=186636 && fRun<=188166) { fLHCperiod="LHC12F"; fBeamType="PP"; /*fMCperiodTPC="";*/ }
688 //   if (fRun >= 188167 && fRun <= 188355 ) { fLHCperiod="LHC12G"; fBeamType="PP"; /*fMCperiodTPC="";*/ }
689 //   if (fRun >= 188356 && fRun <= 188503 ) { fLHCperiod="LHC12G"; fBeamType="PPB"; /*fMCperiodTPC="";*/ }
690 // for the moment use 12g parametrisation for all full gain runs (LHC12f+)
691   if (fRun >= 186636 && fRun < 194480) { fLHCperiod="LHC12G"; fBeamType="PPB"; fMCperiodTPC="LHC12G"; }
692
693   // New parametrisation for 2013 pPb runs
694   if (fRun >= 194480) { 
695     fLHCperiod="LHC13B"; 
696     fBeamType="PPB";
697     fMCperiodTPC="LHC12G";
698   
699     if (fCurrentAliRootRev >= 61605)
700       fMCperiodTPC="LHC13B2_FIX";
701     if (fCurrentAliRootRev >= 62714)
702       fMCperiodTPC="LHC13B2_FIXn1";
703   }
704
705   //exception new pp MC productions from 2011
706   if (fBeamType=="PP" && reg.MatchB(fCurrentFile)) { fMCperiodTPC="LHC11B2"; fBeamType="PP"; }
707   // exception for 11f1
708   if (fCurrentFile.Contains("LHC11f1")) fMCperiodTPC="LHC11F1";
709   // exception for 12f1a, 12f1b and 12i3
710   if (fCurrentFile.Contains("LHC12f1") || fCurrentFile.Contains("LHC12i3")) fMCperiodTPC="LHC12F1";
711   // exception for 12c4
712   if (fCurrentFile.Contains("LHC12c4")) fMCperiodTPC="LHC12C4";
713 }
714
715 //______________________________________________________________________________
716 void AliPIDResponse::SetITSParametrisation()
717 {
718   //
719   // Set the ITS parametrisation
720   //
721 }
722
723  
724 //______________________________________________________________________________
725 void AliPIDResponse::AddPointToHyperplane(TH2D* h, TLinearFitter* linExtrapolation, Int_t binX, Int_t binY)
726 {
727   if (h->GetBinContent(binX, binY) <= 1e-4)
728     return; // Reject bins without content (within some numerical precision) or with strange content
729     
730   Double_t coord[2] = {0, 0};
731   coord[0] = h->GetXaxis()->GetBinCenter(binX);
732   coord[1] = h->GetYaxis()->GetBinCenter(binY);
733   Double_t binError = h->GetBinError(binX, binY);
734   if (binError <= 0) {
735     binError = 1000; // Should not happen because bins without content are rejected for the map (TH2D* h)
736     printf("ERROR: This should never happen: Trying to add bin in addPointToHyperplane with error not set....\n");
737   }
738   linExtrapolation->AddPoint(coord, h->GetBinContent(binX, binY, binError));
739 }
740
741
742 //______________________________________________________________________________
743 TH2D* AliPIDResponse::RefineHistoViaLinearInterpolation(TH2D* h, Double_t refineFactorX, Double_t refineFactorY)
744 {
745   if (!h)
746     return 0x0;
747   
748   // Interpolate to finer map
749   TLinearFitter* linExtrapolation = new TLinearFitter(2, "hyp2", "");
750   
751   Double_t upperMapBoundY = h->GetYaxis()->GetBinUpEdge(h->GetYaxis()->GetNbins());
752   Double_t lowerMapBoundY = h->GetYaxis()->GetBinLowEdge(1);
753   Int_t nBinsX = 30;
754   // Binning was find to yield good results, if 40 bins are chosen for the range 0.0016 to 0.02. For the new variable range,
755   // scale the number of bins correspondingly
756   Int_t nBinsY = TMath::Nint((upperMapBoundY - lowerMapBoundY) / (0.02 - 0.0016) * 40);
757   Int_t nBinsXrefined = nBinsX * refineFactorX;
758   Int_t nBinsYrefined = nBinsY * refineFactorY; 
759   
760   TH2D* hRefined = new TH2D(Form("%s_refined", h->GetName()),  Form("%s (refined)", h->GetTitle()),
761                             nBinsXrefined, h->GetXaxis()->GetBinLowEdge(1), h->GetXaxis()->GetBinUpEdge(h->GetXaxis()->GetNbins()),
762                             nBinsYrefined, lowerMapBoundY, upperMapBoundY);
763   
764   for (Int_t binX = 1; binX <= nBinsXrefined; binX++)  {
765     for (Int_t binY = 1; binY <= nBinsYrefined; binY++)  {
766       
767       hRefined->SetBinContent(binX, binY, 1); // Default value is 1
768       
769       Double_t centerX = hRefined->GetXaxis()->GetBinCenter(binX);
770       Double_t centerY = hRefined->GetYaxis()->GetBinCenter(binY);
771       
772       /*OLD
773       linExtrapolation->ClearPoints();
774       
775       // For interpolation: Just take the corresponding bin from the old histo.
776       // For extrapolation: take the last available bin from the old histo.
777       // If the boundaries are to be skipped, also skip the corresponding bins
778       Int_t oldBinX = h->GetXaxis()->FindBin(centerX);
779       if (oldBinX < 1)  
780         oldBinX = 1;
781       if (oldBinX > nBinsX)
782         oldBinX = nBinsX;
783       
784       Int_t oldBinY = h->GetYaxis()->FindBin(centerY);
785       if (oldBinY < 1)  
786         oldBinY = 1;
787       if (oldBinY > nBinsY)
788         oldBinY = nBinsY;
789       
790       // Neighbours left column
791       if (oldBinX >= 2) {
792         if (oldBinY >= 2) {
793           AddPointToHyperplane(h, linExtrapolation, oldBinX - 1, oldBinY - 1);
794         }
795         
796         AddPointToHyperplane(h, linExtrapolation, oldBinX - 1, oldBinY);
797         
798         if (oldBinY < nBinsY) {
799           AddPointToHyperplane(h, linExtrapolation, oldBinX - 1, oldBinY + 1);
800         }
801       }
802       
803       // Neighbours (and point itself) same column
804       if (oldBinY >= 2) {
805         AddPointToHyperplane(h, linExtrapolation, oldBinX, oldBinY - 1);
806       }
807         
808       AddPointToHyperplane(h, linExtrapolation, oldBinX, oldBinY);
809         
810       if (oldBinY < nBinsY) {
811         AddPointToHyperplane(h, linExtrapolation, oldBinX, oldBinY + 1);
812       }
813       
814       // Neighbours right column
815       if (oldBinX < nBinsX) {
816         if (oldBinY >= 2) {
817           AddPointToHyperplane(h, linExtrapolation, oldBinX + 1, oldBinY - 1);
818         }
819         
820         AddPointToHyperplane(h, linExtrapolation, oldBinX + 1, oldBinY);
821         
822         if (oldBinY < nBinsY) {
823           AddPointToHyperplane(h, linExtrapolation, oldBinX + 1, oldBinY + 1);
824         }
825       }
826       
827       
828       // Fit 2D-hyperplane
829       if (linExtrapolation->GetNpoints() <= 0)
830         continue;
831         
832       if (linExtrapolation->Eval() != 0)// EvalRobust -> Takes much, much, [...], much more time (~hours instead of seconds)
833         continue;
834       
835       // Fill the bin of the refined histogram with the extrapolated value
836       Double_t interpolatedValue = linExtrapolation->GetParameter(0) + linExtrapolation->GetParameter(1) * centerX
837                                  + linExtrapolation->GetParameter(2) * centerY;
838       */
839       Double_t interpolatedValue = h->Interpolate(centerX, centerY) ;
840       hRefined->SetBinContent(binX, binY, interpolatedValue);      
841     }
842   } 
843   
844   
845   // Problem: Interpolation does not work before/beyond center of first/last bin (as the name suggests).
846   // Therefore, for each row in dEdx: Take last bin from old map and interpolate values from center and edge.
847   // Assume line through these points and extropolate to last bin of refined map
848   const Double_t firstOldXbinUpEdge = h->GetXaxis()->GetBinUpEdge(1);
849   const Double_t firstOldXbinCenter = h->GetXaxis()->GetBinCenter(1);
850   
851   const Double_t oldXbinHalfWidth = firstOldXbinUpEdge - firstOldXbinCenter;
852   
853   const Double_t lastOldXbinLowEdge = h->GetXaxis()->GetBinLowEdge(h->GetNbinsX());
854   const Double_t lastOldXbinCenter = h->GetXaxis()->GetBinCenter(h->GetNbinsX());
855   
856   for (Int_t binY = 1; binY <= nBinsYrefined; binY++)  {
857     Double_t centerY = hRefined->GetYaxis()->GetBinCenter(binY);
858     
859     const Double_t interpolatedCenterFirstXbin = h->Interpolate(firstOldXbinCenter, centerY);
860     const Double_t interpolatedUpEdgeFirstXbin = h->Interpolate(firstOldXbinUpEdge, centerY);
861     
862     const Double_t extrapolationSlopeFirstXbin = (interpolatedUpEdgeFirstXbin - interpolatedCenterFirstXbin) / oldXbinHalfWidth;
863     const Double_t extrapolationOffsetFirstXbin = interpolatedCenterFirstXbin;
864     
865     
866     const Double_t interpolatedCenterLastXbin = h->Interpolate(lastOldXbinCenter, centerY);
867     const Double_t interpolatedLowEdgeLastXbin = h->Interpolate(lastOldXbinLowEdge, centerY);
868     
869     const Double_t extrapolationSlopeLastXbin = (interpolatedCenterLastXbin - interpolatedLowEdgeLastXbin) / oldXbinHalfWidth;
870     const Double_t extrapolationOffsetLastXbin = interpolatedCenterLastXbin;
871
872     for (Int_t binX = 1; binX <= nBinsXrefined; binX++)  {
873       Double_t centerX = hRefined->GetXaxis()->GetBinCenter(binX);
874      
875       if (centerX < firstOldXbinCenter) {
876         Double_t extrapolatedValue = extrapolationOffsetFirstXbin + (centerX - firstOldXbinCenter) * extrapolationSlopeFirstXbin;
877         hRefined->SetBinContent(binX, binY, extrapolatedValue);      
878       }
879       else if (centerX <= lastOldXbinCenter) {
880         continue;
881       }
882       else {
883         Double_t extrapolatedValue = extrapolationOffsetLastXbin + (centerX - lastOldXbinCenter) * extrapolationSlopeLastXbin;
884         hRefined->SetBinContent(binX, binY, extrapolatedValue);     
885       }
886     }
887   } 
888   
889   delete linExtrapolation;
890   
891   return hRefined;
892 }
893
894 //______________________________________________________________________________
895 void AliPIDResponse::SetTPCEtaMaps(Double_t refineFactorMapX, Double_t refineFactorMapY,
896                                    Double_t refineFactorSigmaMapX, Double_t refineFactorSigmaMapY)
897 {
898   //
899   // Load the TPC eta correction maps from the OADB
900   //
901   
902   if (fUseTPCEtaCorrection == kFALSE) {
903     // Disable eta correction via setting no maps
904     if (!fTPCResponse.SetEtaCorrMap(0x0))
905       AliInfo("Request to disable TPC eta correction -> Eta correction has been disabled"); 
906     else
907       AliError("Request to disable TPC eta correction -> Some error occured when unloading the correction maps");
908     
909     if (!fTPCResponse.SetSigmaParams(0x0, 0))
910       AliInfo("Request to disable TPC eta correction -> Using old parametrisation for sigma"); 
911     else
912       AliError("Request to disable TPC eta correction -> Some error occured when unloading the sigma maps");
913     
914     return;
915   }
916   
917   TString dataType = "DATA";
918   TString period = fLHCperiod.IsNull() ? "No period information" : fLHCperiod;
919   
920   if (fIsMC)  {
921     if (!(fTuneMConData && ((fTuneMConDataMask & kDetTPC) == kDetTPC))) {
922       period=fMCperiodTPC;
923       dataType="MC";
924     }
925     fRecoPass = 1;
926     
927     if (!(fTuneMConData && ((fTuneMConDataMask & kDetTPC) == kDetTPC)) && fMCperiodTPC.IsNull()) {
928       AliFatal("MC detected, but no MC period set -> Not changing eta maps!");
929       return;
930     }
931   }
932
933   Int_t recopass = fRecoPass;
934   if (fTuneMConData && ((fTuneMConDataMask & kDetTPC) == kDetTPC) )
935     recopass = fRecoPassUser;
936   
937   TString defaultObj = Form("Default_%s_pass%d", dataType.Data(), recopass);
938   
939   AliInfo(Form("Current period and reco pass: %s.pass%d", period.Data(), recopass));
940   
941   // Invalidate old maps
942   fTPCResponse.SetEtaCorrMap(0x0);
943   fTPCResponse.SetSigmaParams(0x0, 0);
944   
945   // Load the eta correction maps
946   AliOADBContainer etaMapsCont(Form("TPCetaMaps_%s_pass%d", dataType.Data(), recopass)); 
947   
948   Int_t statusCont = etaMapsCont.InitFromFile(Form("%s/COMMON/PID/data/TPCetaMaps.root", fOADBPath.Data()),
949                                               Form("TPCetaMaps_%s_pass%d", dataType.Data(), recopass));
950   if (statusCont) {
951     AliError("Failed initializing TPC eta correction maps from OADB -> Disabled eta correction");
952     fUseTPCEtaCorrection = kFALSE;
953   }
954   else {
955     AliInfo(Form("Loading TPC eta correction map from %s/COMMON/PID/data/TPCetaMaps.root", fOADBPath.Data()));
956     
957     TH2D* etaMap = 0x0;
958     
959     if (fIsMC && !(fTuneMConData && ((fTuneMConDataMask & kDetTPC) == kDetTPC))) {
960       TString searchMap = Form("TPCetaMaps_%s_%s_pass%d", dataType.Data(), period.Data(), recopass);
961       etaMap = dynamic_cast<TH2D *>(etaMapsCont.GetDefaultObject(searchMap.Data()));
962       if (!etaMap) {
963         // Try default object
964         etaMap = dynamic_cast<TH2D *>(etaMapsCont.GetDefaultObject(defaultObj.Data()));
965       }
966     }
967     else {
968       etaMap = dynamic_cast<TH2D *>(etaMapsCont.GetObject(fRun, defaultObj.Data()));
969     }
970     
971         
972     if (!etaMap) {
973       AliError(Form("TPC eta correction map not found for run %d and also no default map found -> Disabled eta correction!!!", fRun));
974       fUseTPCEtaCorrection = kFALSE;
975     }
976     else {
977       TH2D* etaMapRefined = RefineHistoViaLinearInterpolation(etaMap, refineFactorMapX, refineFactorMapY);
978       
979       if (etaMapRefined) {
980         if (!fTPCResponse.SetEtaCorrMap(etaMapRefined)) {
981           AliError(Form("Failed to set TPC eta correction map for run %d -> Disabled eta correction!!!", fRun));
982           fTPCResponse.SetEtaCorrMap(0x0);
983           fUseTPCEtaCorrection = kFALSE;
984         }
985         else {
986           AliInfo(Form("Loaded TPC eta correction map (refine factors %.2f/%.2f) from %s/COMMON/PID/data/TPCetaMaps.root: %s", 
987                        refineFactorMapX, refineFactorMapY, fOADBPath.Data(), fTPCResponse.GetEtaCorrMap()->GetTitle()));
988         }
989         
990         delete etaMapRefined;
991       }
992       else {
993         AliError(Form("Failed to set TPC eta correction map for run %d (map was loaded, but couldn't be refined) -> Disabled eta correction!!!", fRun));
994         fUseTPCEtaCorrection = kFALSE;
995       }
996     }
997   }
998   
999   // If there was some problem loading the eta maps, it makes no sense to load the sigma maps (that require eta corrected data)
1000   if (fUseTPCEtaCorrection == kFALSE) {
1001     AliError("Failed to load TPC eta correction map required by sigma maps -> Using old parametrisation for sigma"); 
1002     return;
1003   }
1004   
1005   // Load the sigma parametrisation (1/dEdx vs tanTheta_local (~eta))
1006   AliOADBContainer etaSigmaMapsCont(Form("TPCetaSigmaMaps_%s_pass%d", dataType.Data(), recopass)); 
1007   
1008   statusCont = etaSigmaMapsCont.InitFromFile(Form("%s/COMMON/PID/data/TPCetaMaps.root", fOADBPath.Data()),
1009                                              Form("TPCetaSigmaMaps_%s_pass%d", dataType.Data(), recopass));
1010   if (statusCont) {
1011     AliError("Failed initializing TPC eta sigma maps from OADB -> Using old sigma parametrisation");
1012   }
1013   else {
1014     AliInfo(Form("Loading TPC eta sigma map from %s/COMMON/PID/data/TPCetaMaps.root", fOADBPath.Data()));
1015     
1016     TObjArray* etaSigmaPars = 0x0;
1017     
1018     if (fIsMC && !(fTuneMConData && ((fTuneMConDataMask & kDetTPC) == kDetTPC))) {
1019       TString searchMap = Form("TPCetaSigmaMaps_%s_%s_pass%d", dataType.Data(), period.Data(), recopass);
1020       etaSigmaPars = dynamic_cast<TObjArray *>(etaSigmaMapsCont.GetDefaultObject(searchMap.Data()));
1021       if (!etaSigmaPars) {
1022         // Try default object
1023         etaSigmaPars = dynamic_cast<TObjArray *>(etaSigmaMapsCont.GetDefaultObject(defaultObj.Data()));
1024       }
1025     }
1026     else {
1027       etaSigmaPars = dynamic_cast<TObjArray *>(etaSigmaMapsCont.GetObject(fRun, defaultObj.Data()));
1028     }
1029     
1030     if (!etaSigmaPars) {
1031       AliError(Form("TPC eta sigma parametrisation not found for run %d -> Using old sigma parametrisation!!!", fRun));
1032     }
1033     else {
1034       TH2D* etaSigmaPar1Map = dynamic_cast<TH2D *>(etaSigmaPars->FindObject("sigmaPar1Map"));
1035       TNamed* sigmaPar0Info = dynamic_cast<TNamed *>(etaSigmaPars->FindObject("sigmaPar0"));
1036       Double_t sigmaPar0 = 0.0;
1037       
1038       if (sigmaPar0Info) {
1039         TString sigmaPar0String = sigmaPar0Info->GetTitle();
1040         sigmaPar0 = sigmaPar0String.Atof();
1041       }
1042       else {
1043         // Something is weired because the object for parameter 0 could not be loaded -> New sigma parametrisation can not be used!
1044         etaSigmaPar1Map = 0x0;
1045       }
1046       
1047       TH2D* etaSigmaPar1MapRefined = RefineHistoViaLinearInterpolation(etaSigmaPar1Map, refineFactorSigmaMapX, refineFactorSigmaMapY);
1048       
1049       
1050       if (etaSigmaPar1MapRefined) {
1051         if (!fTPCResponse.SetSigmaParams(etaSigmaPar1MapRefined, sigmaPar0)) {
1052           AliError(Form("Failed to set TPC eta sigma map for run %d -> Using old sigma parametrisation!!!", fRun));
1053           fTPCResponse.SetSigmaParams(0x0, 0);
1054         }
1055         else {
1056           AliInfo(Form("Loaded TPC sigma correction map (refine factors %.2f/%.2f) from %s/COMMON/PID/data/TPCetaMaps.root: %s", 
1057                        refineFactorSigmaMapX, refineFactorSigmaMapY, fOADBPath.Data(), fTPCResponse.GetSigmaPar1Map()->GetTitle()));
1058         }
1059         
1060         delete etaSigmaPar1MapRefined;
1061       }
1062       else {
1063         AliError(Form("Failed to set TPC eta sigma map for run %d (map was loaded, but couldn't be refined) -> Using old sigma parametrisation!!!",
1064                       fRun));
1065       }
1066     }
1067   }
1068 }
1069
1070 //______________________________________________________________________________
1071 void AliPIDResponse::SetTPCPidResponseMaster()
1072 {
1073   //
1074   // Load the TPC pid response functions from the OADB
1075   // Load the TPC voltage maps from OADB
1076   //
1077   //don't load twice for the moment
1078    if (fArrPidResponseMaster) return;
1079  
1080
1081   //reset the PID response functions
1082   delete fArrPidResponseMaster;
1083   fArrPidResponseMaster=NULL;
1084   
1085   TString fileName(Form("%s/COMMON/PID/data/TPCPIDResponse.root", fOADBPath.Data()));
1086   TFile *f=NULL;
1087   if (!fCustomTPCpidResponse.IsNull()) fileName=fCustomTPCpidResponse;
1088   
1089   TString fileNamePIDresponse(Form("%s/COMMON/PID/data/TPCPIDResponse.root", fOADBPath.Data()));
1090   f=TFile::Open(fileNamePIDresponse.Data());
1091   if (f && f->IsOpen() && !f->IsZombie()){
1092     fArrPidResponseMaster=dynamic_cast<TObjArray*>(f->Get("TPCPIDResponse"));
1093   }
1094   delete f;
1095
1096   TString fileNameVoltageMaps(Form("%s/COMMON/PID/data/TPCvoltageSettings.root", fOADBPath.Data()));
1097   f=TFile::Open(fileNameVoltageMaps.Data());
1098   if (f && f->IsOpen() && !f->IsZombie()){
1099     fOADBvoltageMaps=dynamic_cast<AliOADBContainer*>(f->Get("TPCvoltageSettings"));
1100   }
1101   delete f;
1102   
1103   if (!fArrPidResponseMaster){
1104     AliFatal(Form("Could not retrieve the TPC pid response from: %s",fileNamePIDresponse.Data()));
1105     return;
1106   }
1107   fArrPidResponseMaster->SetOwner();
1108
1109   if (!fOADBvoltageMaps)
1110   {
1111     AliFatal(Form("Could not retrieve the TPC voltage maps from: %s",fileNameVoltageMaps.Data()));
1112   }
1113   fArrPidResponseMaster->SetOwner();
1114 }
1115
1116 //______________________________________________________________________________
1117 void AliPIDResponse::SetTPCParametrisation()
1118 {
1119   //
1120   // Change BB parametrisation for current run
1121   //
1122   
1123   //
1124   //reset old splines
1125   //
1126   fTPCResponse.ResetSplines();
1127   
1128   if (fLHCperiod.IsNull()) {
1129     AliError("No period set, not changing parametrisation");
1130     return;
1131   }
1132   
1133   //
1134   // Set default parametrisations for data and MC
1135   //
1136   
1137   //data type
1138   TString datatype="DATA";
1139   //in case of mc fRecoPass is per default 1
1140   if (fIsMC) {
1141       if(!(fTuneMConData && ((fTuneMConDataMask & kDetTPC) == kDetTPC))) datatype="MC";
1142       fRecoPass=1;
1143   }
1144
1145   // period
1146   TString period=fLHCperiod;
1147   if (fIsMC && !(fTuneMConData && ((fTuneMConDataMask & kDetTPC) == kDetTPC))) period=fMCperiodTPC;
1148
1149   Int_t recopass = fRecoPass;
1150   if(fTuneMConData && ((fTuneMConDataMask & kDetTPC) == kDetTPC)) recopass = fRecoPassUser;
1151     
1152   AliInfo(Form("Searching splines for: %s %s PASS%d %s",datatype.Data(),period.Data(),recopass,fBeamType.Data()));
1153   Bool_t found=kFALSE;
1154   //
1155   //set the new PID splines
1156   //
1157   if (fArrPidResponseMaster){
1158     //for MC don't use period information
1159     //if (fIsMC) period="[A-Z0-9]*";
1160     //for MC use MC period information
1161     //pattern for the default entry (valid for all particles)
1162     TPRegexp reg(Form("TSPLINE3_%s_([A-Z]*)_%s_PASS%d_%s_MEAN(_*)([A-Z1-9]*)",datatype.Data(),period.Data(),recopass,fBeamType.Data()));
1163
1164     //find particle id and gain scenario
1165     for (Int_t igainScenario=0; igainScenario<AliTPCPIDResponse::fgkNumberOfGainScenarios; igainScenario++)
1166     {
1167       TObject *grAll=NULL;
1168       TString gainScenario = AliTPCPIDResponse::GainScenarioName(igainScenario);
1169       gainScenario.ToUpper();
1170       //loop over entries and filter them
1171       for (Int_t iresp=0; iresp<fArrPidResponseMaster->GetEntriesFast();++iresp)
1172       {
1173         TObject *responseFunction=fArrPidResponseMaster->At(iresp);
1174         if (responseFunction==NULL) continue;
1175         TString responseName=responseFunction->GetName();
1176          
1177         if (!reg.MatchB(responseName)) continue;
1178
1179         TObjArray *arr=reg.MatchS(responseName); if (!arr) continue;
1180         TObject* tmp=NULL;
1181         tmp=arr->At(1); if (!tmp) continue;
1182         TString particleName=tmp->GetName();
1183         tmp=arr->At(3); if (!tmp) continue;
1184         TString gainScenarioName=tmp->GetName();
1185         delete arr;
1186         if (particleName.IsNull()) continue;
1187         if (!grAll && particleName=="ALL" && gainScenarioName==gainScenario) grAll=responseFunction;
1188         else 
1189         {
1190           for (Int_t ispec=0; ispec<(AliTPCPIDResponse::fgkNumberOfParticleSpecies); ++ispec)
1191           {
1192             TString particle=AliPID::ParticleName(ispec);
1193             particle.ToUpper();
1194             //std::cout<<responseName<<" "<<particle<<" "<<particleName<<" "<<gainScenario<<" "<<gainScenarioName<<std::endl;
1195             if ( particle == particleName && gainScenario == gainScenarioName )
1196             {
1197               fTPCResponse.SetResponseFunction( responseFunction,
1198                                                 (AliPID::EParticleType)ispec,
1199                                                 (AliTPCPIDResponse::ETPCgainScenario)igainScenario );
1200               fTPCResponse.SetUseDatabase(kTRUE);
1201               AliInfo(Form("Adding graph: %d %d - %s",ispec,igainScenario,responseFunction->GetName()));
1202               found=kTRUE;
1203               break;
1204             }
1205           }
1206         }
1207       }
1208       
1209       // Retrieve responsefunction for pions - will (if available) be used for muons if there are no dedicated muon splines.
1210       // For light nuclei, try to set the proton spline, if no dedicated splines are available.
1211       // In both cases: Use default splines, if no dedicated splines and no pion/proton splines are available.
1212       TObject* responseFunctionPion = fTPCResponse.GetResponseFunction( (AliPID::EParticleType)AliPID::kPion,                             
1213                                                                         (AliTPCPIDResponse::ETPCgainScenario)igainScenario);
1214       TObject* responseFunctionProton = fTPCResponse.GetResponseFunction( (AliPID::EParticleType)AliPID::kProton,                             
1215                                                                           (AliTPCPIDResponse::ETPCgainScenario)igainScenario);
1216       
1217       for (Int_t ispec=0; ispec<(AliTPCPIDResponse::fgkNumberOfParticleSpecies); ++ispec)
1218       {
1219         if (!fTPCResponse.GetResponseFunction( (AliPID::EParticleType)ispec,
1220           (AliTPCPIDResponse::ETPCgainScenario)igainScenario))
1221         {
1222           if (ispec == AliPID::kMuon) { // Muons
1223             if (responseFunctionPion) {
1224               fTPCResponse.SetResponseFunction( responseFunctionPion,
1225                                                 (AliPID::EParticleType)ispec,
1226                                                 (AliTPCPIDResponse::ETPCgainScenario)igainScenario );
1227               fTPCResponse.SetUseDatabase(kTRUE);
1228               AliInfo(Form("Adding graph: %d %d - %s",ispec,igainScenario,responseFunctionPion->GetName()));
1229               found=kTRUE;  
1230             }
1231             else if (grAll) {
1232               fTPCResponse.SetResponseFunction( grAll,
1233                                                 (AliPID::EParticleType)ispec,
1234                                                 (AliTPCPIDResponse::ETPCgainScenario)igainScenario );
1235               fTPCResponse.SetUseDatabase(kTRUE);
1236               AliInfo(Form("Adding graph: %d %d - %s",ispec,igainScenario,grAll->GetName()));
1237               found=kTRUE;
1238             }
1239             //else
1240             //  AliError(Form("No splines found for muons (also no pion splines and no default splines) for gain scenario %d!", igainScenario));
1241           }
1242           else if (ispec >= AliPID::kSPECIES) { // Light nuclei
1243             if (responseFunctionProton) {
1244               fTPCResponse.SetResponseFunction( responseFunctionProton,
1245                                                 (AliPID::EParticleType)ispec,
1246                                                 (AliTPCPIDResponse::ETPCgainScenario)igainScenario );
1247               fTPCResponse.SetUseDatabase(kTRUE);
1248               AliInfo(Form("Adding graph: %d %d - %s",ispec,igainScenario,responseFunctionProton->GetName()));
1249               found=kTRUE;  
1250             }
1251             else if (grAll) {
1252               fTPCResponse.SetResponseFunction( grAll,
1253                                                 (AliPID::EParticleType)ispec,
1254                                                 (AliTPCPIDResponse::ETPCgainScenario)igainScenario );
1255               fTPCResponse.SetUseDatabase(kTRUE);
1256               AliInfo(Form("Adding graph: %d %d - %s",ispec,igainScenario,grAll->GetName()));
1257               found=kTRUE;
1258             }
1259             //else
1260             //  AliError(Form("No splines found for species %d (also no proton splines and no default splines) for gain scenario %d!",
1261             //                ispec, igainScenario));
1262           }
1263         }
1264       }
1265     }
1266   }
1267   else AliInfo("no fArrPidResponseMaster");
1268
1269   if (!found){
1270     AliError(Form("No splines found for: %s %s PASS%d %s",datatype.Data(),period.Data(),recopass,fBeamType.Data()));
1271   }
1272
1273
1274   //
1275   // Setup multiplicity correction
1276   //
1277   if (fUseTPCMultiplicityCorrection && (fBeamType.CompareTo("PBPB") == 0 || fBeamType.CompareTo("AA") == 0)) {
1278     AliInfo("Multiplicity correction enabled!");
1279     
1280     //TODO After testing, load parameters from outside       
1281     if (period.Contains("LHC11A10"))  {//LHC11A10A
1282       AliInfo("Using multiplicity correction parameters for 11a10!");
1283       fTPCResponse.SetParameterMultiplicityCorrection(0, 6.90133e-06);
1284       fTPCResponse.SetParameterMultiplicityCorrection(1, -1.22123e-03);
1285       fTPCResponse.SetParameterMultiplicityCorrection(2, 1.80220e-02);
1286       fTPCResponse.SetParameterMultiplicityCorrection(3, 0.1);
1287       fTPCResponse.SetParameterMultiplicityCorrection(4, 6.45306e-03);
1288       
1289       fTPCResponse.SetParameterMultiplicityCorrectionTanTheta(0, -2.85505e-07);
1290       fTPCResponse.SetParameterMultiplicityCorrectionTanTheta(1, -1.31911e-06);
1291       fTPCResponse.SetParameterMultiplicityCorrectionTanTheta(2, -0.5);
1292
1293       fTPCResponse.SetParameterMultiplicitySigmaCorrection(0, -4.29665e-05);
1294       fTPCResponse.SetParameterMultiplicitySigmaCorrection(1, 1.37023e-02);
1295       fTPCResponse.SetParameterMultiplicitySigmaCorrection(2, -6.36337e-01);
1296       fTPCResponse.SetParameterMultiplicitySigmaCorrection(3, 1.13479e-02);
1297     }
1298     else if (period.Contains("LHC10H") && recopass == 2) {    
1299       AliInfo("Using multiplicity correction parameters for 10h.pass2!");
1300       fTPCResponse.SetParameterMultiplicityCorrection(0, 3.21636e-07);
1301       fTPCResponse.SetParameterMultiplicityCorrection(1, -6.65876e-04);
1302       fTPCResponse.SetParameterMultiplicityCorrection(2, 1.28786e-03);
1303       fTPCResponse.SetParameterMultiplicityCorrection(3, 1.47677e-02);
1304       fTPCResponse.SetParameterMultiplicityCorrection(4, 0);
1305       
1306       fTPCResponse.SetParameterMultiplicityCorrectionTanTheta(0, 7.23591e-08);
1307       fTPCResponse.SetParameterMultiplicityCorrectionTanTheta(1, 2.7469e-06);
1308       fTPCResponse.SetParameterMultiplicityCorrectionTanTheta(2, -0.5);
1309       
1310       fTPCResponse.SetParameterMultiplicitySigmaCorrection(0, -1.22590e-05);
1311       fTPCResponse.SetParameterMultiplicitySigmaCorrection(1, 6.88888e-03);
1312       fTPCResponse.SetParameterMultiplicitySigmaCorrection(2, -3.20788e-01);
1313       fTPCResponse.SetParameterMultiplicitySigmaCorrection(3, 1.07345e-02);
1314     }
1315     else {
1316       AliError(Form("Multiplicity correction is enabled, but no multiplicity correction parameters have been found for period %s.pass%d -> Mulitplicity correction DISABLED!", period.Data(), recopass));
1317       fUseTPCMultiplicityCorrection = kFALSE;
1318       fTPCResponse.ResetMultiplicityCorrectionFunctions();
1319     }
1320   }
1321   else {
1322     // Just set parameters such that overall correction factor is 1, i.e. no correction.
1323     // This is just a reasonable choice for the parameters for safety reasons. Disabling
1324     // the multiplicity correction will anyhow skip the calculation of the corresponding
1325     // correction factor inside THIS class. Nevertheless, experts can access the TPCPIDResponse
1326     // directly and use it for calculations - which should still give valid results, even if
1327     // the multiplicity correction is explicitely enabled in such expert calls.
1328     
1329     AliInfo(Form("Multiplicity correction %sdisabled (%s)!", fUseTPCMultiplicityCorrection ? "automatically " : "",
1330                  fUseTPCMultiplicityCorrection ? "no PbPb or AA" : "requested by user"));
1331     
1332     fUseTPCMultiplicityCorrection = kFALSE;
1333     fTPCResponse.ResetMultiplicityCorrectionFunctions();
1334   }
1335   
1336   /*
1337   //TODO NOW start
1338   for (Int_t i = 0; i <= 4 + 1; i++) {
1339     printf("parMultCorr: %d, %e\n", i, fTPCResponse.GetMultiplicityCorrectionFunction()->GetParameter(i));
1340   }
1341   for (Int_t j = 0; j <= 2 + 1; j++) {
1342     printf("parMultCorrTanTheta: %d, %e\n", j, fTPCResponse.GetMultiplicityCorrectionFunctionTanTheta()->GetParameter(j));
1343   }
1344   for (Int_t j = 0; j <= 3 + 1; j++) {
1345     printf("parMultSigmaCorr: %d, %e\n", j, fTPCResponse.GetMultiplicitySigmaCorrectionFunction()->GetParameter(j));
1346   }
1347   
1348   //TODO NOW end
1349   */
1350   
1351   //
1352   // Setup old resolution parametrisation
1353   //
1354   
1355   //default
1356   fTPCResponse.SetSigma(3.79301e-03, 2.21280e+04);
1357   
1358   if (fRun>=122195){
1359     fTPCResponse.SetSigma(2.30176e-02, 5.60422e+02);
1360   }
1361
1362   if (fRun>=186636){
1363 //   if (fRun>=188356){
1364     fTPCResponse.SetSigma(8.62022e-04, 9.08156e+05);
1365   }
1366   
1367   if (fArrPidResponseMaster)
1368   fResolutionCorrection=(TF1*)fArrPidResponseMaster->FindObject(Form("TF1_%s_ALL_%s_PASS%d_%s_SIGMA",datatype.Data(),period.Data(),recopass,fBeamType.Data()));
1369   
1370   if (fResolutionCorrection) AliInfo(Form("Setting multiplicity correction function: %s",fResolutionCorrection->GetName()));
1371
1372   //read in the voltage map
1373   TVectorF* gsm = 0x0;
1374   if (fOADBvoltageMaps) gsm=dynamic_cast<TVectorF*>(fOADBvoltageMaps->GetObject(fRun));
1375   if (gsm) 
1376   {
1377     fTPCResponse.SetVoltageMap(*gsm);
1378     TString vals;
1379     AliInfo(Form("Reading the voltage map for run %d\n",fRun));
1380     vals="IROC A: "; for (Int_t i=0; i<18; i++){vals+=Form("%.2f ",(*gsm)[i]);}
1381     AliInfo(vals.Data());
1382     vals="IROC C: "; for (Int_t i=18; i<36; i++){vals+=Form("%.2f ",(*gsm)[i]);}
1383     AliInfo(vals.Data());
1384     vals="OROC A: "; for (Int_t i=36; i<54; i++){vals+=Form("%.2f ",(*gsm)[i]);}
1385     AliInfo(vals.Data());
1386     vals="OROC C: "; for (Int_t i=54; i<72; i++){vals+=Form("%.2f ",(*gsm)[i]);}
1387     AliInfo(vals.Data());
1388   }
1389   else AliInfo("no voltage map, ideal default assumed");
1390 }
1391
1392 //______________________________________________________________________________
1393 void AliPIDResponse::SetTRDPidResponseMaster()
1394 {
1395   //
1396   // Load the TRD pid params and references from the OADB
1397   //
1398   if(fTRDPIDResponseObject) return;
1399   AliOADBContainer contParams("contParams"); 
1400
1401   Int_t statusResponse = contParams.InitFromFile(Form("%s/COMMON/PID/data/TRDPIDResponse.root", fOADBPath.Data()), "AliTRDPIDResponseObject");
1402   if(statusResponse){
1403     AliError("Failed initializing PID Response Object from OADB");
1404   } else {
1405     AliInfo(Form("Loading TRD Response from %s/COMMON/PID/data/TRDPIDResponse.root", fOADBPath.Data()));
1406     fTRDPIDResponseObject = dynamic_cast<AliTRDPIDResponseObject *>(contParams.GetObject(fRun));
1407     if(!fTRDPIDResponseObject){
1408       AliError(Form("TRD Response not found in run %d", fRun));
1409     }
1410   }
1411 }
1412
1413 //______________________________________________________________________________
1414 void AliPIDResponse::InitializeTRDResponse(){
1415   //
1416   // Set PID Params and references to the TRD PID response
1417   // 
1418     fTRDResponse.SetPIDResponseObject(fTRDPIDResponseObject);
1419 }
1420
1421 //______________________________________________________________________________
1422 void AliPIDResponse::SetTRDSlices(UInt_t TRDslicesForPID[2],AliTRDPIDResponse::ETRDPIDMethod method) const{
1423
1424     if(fLHCperiod.Contains("LHC10D") || fLHCperiod.Contains("LHC10E")){
1425         // backward compatibility for setting with 8 slices
1426         TRDslicesForPID[0] = 0;
1427         TRDslicesForPID[1] = 7;
1428     }
1429     else{
1430         if(method==AliTRDPIDResponse::kLQ1D){
1431             TRDslicesForPID[0] = 0; // first Slice contains normalized dEdx
1432             TRDslicesForPID[1] = 0;
1433         }
1434         if(method==AliTRDPIDResponse::kLQ2D){
1435             TRDslicesForPID[0] = 1;
1436             TRDslicesForPID[1] = 7;
1437         }
1438     }
1439     AliDebug(1,Form("Slice Range set to %d - %d",TRDslicesForPID[0],TRDslicesForPID[1]));
1440 }
1441
1442 //______________________________________________________________________________
1443 void AliPIDResponse::SetTOFPidResponseMaster()
1444 {
1445   //
1446   // Load the TOF pid params from the OADB
1447   //
1448
1449   if (fTOFPIDParams) delete fTOFPIDParams;
1450   fTOFPIDParams=NULL;
1451
1452   TFile *oadbf = new TFile(Form("%s/COMMON/PID/data/TOFPIDParams.root",fOADBPath.Data()));
1453   if (oadbf && oadbf->IsOpen()) {
1454     AliInfo(Form("Loading TOF Params from %s/COMMON/PID/data/TOFPIDParams.root", fOADBPath.Data()));
1455     AliOADBContainer *oadbc = (AliOADBContainer *)oadbf->Get("TOFoadb");
1456     if (oadbc) fTOFPIDParams = dynamic_cast<AliTOFPIDParams *>(oadbc->GetObject(fRun,"TOFparams"));
1457     oadbf->Close();
1458     delete oadbc;
1459   }
1460   delete oadbf;
1461
1462   if (!fTOFPIDParams) AliFatal("TOFPIDParams could not be retrieved");
1463 }
1464
1465 //______________________________________________________________________________
1466 void AliPIDResponse::InitializeTOFResponse(){
1467   //
1468   // Set PID Params to the TOF PID response
1469   //
1470
1471   AliInfo("TOF PID Params loaded from OADB");
1472   AliInfo(Form("  TOF resolution %5.2f [ps]",fTOFPIDParams->GetTOFresolution()));
1473   AliInfo(Form("  StartTime method %d",fTOFPIDParams->GetStartTimeMethod()));
1474   AliInfo(Form("  TOF res. mom. params: %5.2f %5.2f %5.2f %5.2f",
1475                fTOFPIDParams->GetSigParams(0),fTOFPIDParams->GetSigParams(1),fTOFPIDParams->GetSigParams(2),fTOFPIDParams->GetSigParams(3)));
1476   
1477   for (Int_t i=0;i<4;i++) {
1478     fTOFResponse.SetTrackParameter(i,fTOFPIDParams->GetSigParams(i));
1479   }
1480   fTOFResponse.SetTimeResolution(fTOFPIDParams->GetTOFresolution());
1481
1482   AliInfo("TZERO resolution loaded from ESDrun/AODheader");
1483   Float_t t0Spread[4];
1484   for (Int_t i=0;i<4;i++) t0Spread[i]=fCurrentEvent->GetT0spread(i);
1485   AliInfo(Form("  TZERO spreads from data: (A+C)/2 %f A %f C %f (A'-C')/2: %f",t0Spread[0],t0Spread[1],t0Spread[2],t0Spread[3]));
1486   Float_t a = t0Spread[1]*t0Spread[1]-t0Spread[0]*t0Spread[0]+t0Spread[3]*t0Spread[3];
1487   Float_t c = t0Spread[2]*t0Spread[2]-t0Spread[0]*t0Spread[0]+t0Spread[3]*t0Spread[3];
1488   if ( (t0Spread[0] > 50. && t0Spread[0] < 400.) && (a > 0.) && (c>0.)) {
1489     fResT0AC=t0Spread[3];
1490     fResT0A=TMath::Sqrt(a);
1491     fResT0C=TMath::Sqrt(c);
1492   } else {
1493     AliInfo("  TZERO spreads not present or inconsistent, loading default");
1494     fResT0A=75.;
1495     fResT0C=65.;
1496     fResT0AC=55.;
1497   }
1498   AliInfo(Form("  TZERO resolution set to: T0A: %f [ps] T0C: %f [ps] T0AC %f [ps]",fResT0A,fResT0C,fResT0AC));
1499
1500 }
1501
1502 //______________________________________________________________________________
1503 void AliPIDResponse::SetHMPIDPidResponseMaster()
1504 {
1505   //
1506   // Load the HMPID pid params from the OADB
1507   //
1508
1509   if (fHMPIDPIDParams) delete fHMPIDPIDParams;
1510   fHMPIDPIDParams=NULL;
1511
1512   TFile *oadbf = new TFile(Form("%s/COMMON/PID/data/HMPIDPIDParams.root",fOADBPath.Data()));
1513   if (oadbf && oadbf->IsOpen()) {
1514     AliInfo(Form("Loading HMPID Params from %s/COMMON/PID/data/HMPIDPIDParams.root", fOADBPath.Data()));
1515     AliOADBContainer *oadbc = (AliOADBContainer *)oadbf->Get("HMPoadb");
1516     if (oadbc) fHMPIDPIDParams = dynamic_cast<AliHMPIDPIDParams *>(oadbc->GetObject(fRun,"HMPparams"));
1517     oadbf->Close();
1518     delete oadbc;
1519   }
1520   delete oadbf;
1521
1522   if (!fHMPIDPIDParams) AliFatal("HMPIDPIDParams could not be retrieved");
1523 }
1524
1525 //______________________________________________________________________________
1526 void AliPIDResponse::InitializeHMPIDResponse(){
1527   //
1528   // Set PID Params to the HMPID PID response
1529   //
1530
1531   fHMPIDResponse.SetRefIndexArray(fHMPIDPIDParams->GetHMPIDrefIndex());
1532 }
1533
1534 //______________________________________________________________________________
1535 Bool_t AliPIDResponse::IdentifiedAsElectronTRD(const AliVTrack *vtrack, Double_t efficiencyLevel,Double_t centrality,AliTRDPIDResponse::ETRDPIDMethod PIDmethod) const {
1536   //
1537   // Check whether track is identified as electron under a given electron efficiency hypothesis
1538     //
1539
1540   Double_t probs[AliPID::kSPECIES];
1541   ComputeTRDProbability(vtrack, AliPID::kSPECIES, probs,PIDmethod);
1542
1543   Int_t ntracklets = vtrack->GetTRDntrackletsPID();
1544   // Take mean of the TRD momenta in the given tracklets
1545   Float_t p = 0, trdmomenta[AliVTrack::kTRDnPlanes];
1546   Int_t nmomenta = 0;
1547   for(Int_t iPl=0;iPl<AliVTrack::kTRDnPlanes;iPl++){
1548     if(vtrack->GetTRDmomentum(iPl) > 0.){
1549       trdmomenta[nmomenta++] = vtrack->GetTRDmomentum(iPl); 
1550     }
1551   }
1552   p = TMath::Mean(nmomenta, trdmomenta);
1553
1554   return fTRDResponse.IdentifiedAsElectron(ntracklets, probs, p, efficiencyLevel,centrality,PIDmethod);
1555 }
1556
1557 //______________________________________________________________________________
1558 void AliPIDResponse::SetEMCALPidResponseMaster()
1559 {
1560   //
1561   // Load the EMCAL pid response functions from the OADB
1562   //
1563   TObjArray* fEMCALPIDParamsRun      = NULL;
1564   TObjArray* fEMCALPIDParamsPass     = NULL;
1565
1566   if(fEMCALPIDParams) return;
1567   AliOADBContainer contParams("contParams"); 
1568
1569   Int_t statusPars = contParams.InitFromFile(Form("%s/COMMON/PID/data/EMCALPIDParams.root", fOADBPath.Data()), "AliEMCALPIDParams");
1570   if(statusPars){
1571     AliError("Failed initializing PID Params from OADB");
1572   } 
1573   else {
1574     AliInfo(Form("Loading EMCAL Params from %s/COMMON/PID/data/EMCALPIDParams.root", fOADBPath.Data()));
1575
1576     fEMCALPIDParamsRun = dynamic_cast<TObjArray *>(contParams.GetObject(fRun));
1577     if(fEMCALPIDParamsRun)  fEMCALPIDParamsPass = dynamic_cast<TObjArray *>(fEMCALPIDParamsRun->FindObject(Form("pass%d",fRecoPass)));
1578     if(fEMCALPIDParamsPass) fEMCALPIDParams     = dynamic_cast<TObjArray *>(fEMCALPIDParamsPass->FindObject(Form("EMCALPIDParams_Particles")));
1579
1580     if(!fEMCALPIDParams){
1581       AliInfo(Form("EMCAL Params not found in run %d pass %d", fRun, fRecoPass));
1582       AliInfo("Will take the standard LHC11d instead ...");
1583
1584       fEMCALPIDParamsRun = dynamic_cast<TObjArray *>(contParams.GetObject(156477));
1585       if(fEMCALPIDParamsRun)  fEMCALPIDParamsPass = dynamic_cast<TObjArray *>(fEMCALPIDParamsRun->FindObject(Form("pass%d",1)));
1586       if(fEMCALPIDParamsPass) fEMCALPIDParams     = dynamic_cast<TObjArray *>(fEMCALPIDParamsPass->FindObject(Form("EMCALPIDParams_Particles")));
1587
1588       if(!fEMCALPIDParams){
1589         AliError(Form("DEFAULT EMCAL Params (LHC11d) not found in file %s/COMMON/PID/data/EMCALPIDParams.root", fOADBPath.Data()));     
1590       }
1591     }
1592   }
1593 }
1594
1595 //______________________________________________________________________________
1596 void AliPIDResponse::InitializeEMCALResponse(){
1597   //
1598   // Set PID Params to the EMCAL PID response
1599   // 
1600   fEMCALResponse.SetPIDParams(fEMCALPIDParams);
1601
1602 }
1603
1604 //______________________________________________________________________________
1605 void AliPIDResponse::FillTrackDetectorPID(const AliVTrack *track, EDetector detector) const
1606 {
1607   //
1608   // create detector PID information and setup the transient pointer in the track
1609   //
1610   
1611   // check if detector number is inside accepted range
1612   if (detector == kNdetectors) return;
1613   
1614   // get detector pid
1615   AliDetectorPID *detPID=const_cast<AliDetectorPID*>(track->GetDetectorPID());
1616   if (!detPID) {
1617     detPID=new AliDetectorPID;
1618     (const_cast<AliVTrack*>(track))->SetDetectorPID(detPID);
1619   }
1620   
1621   //check if values exist
1622   if (detPID->HasRawProbability(detector) && detPID->HasNumberOfSigmas(detector)) return;
1623   
1624   //TODO: which particles to include? See also the loops below...
1625   Double_t values[AliPID::kSPECIESC]={0};
1626
1627   //probabilities
1628   EDetPidStatus status=GetComputePIDProbability(detector,track,AliPID::kSPECIESC,values);
1629   detPID->SetRawProbability(detector, values, (Int_t)AliPID::kSPECIESC, status);
1630   
1631   //nsigmas
1632   for (Int_t ipart=0; ipart<AliPID::kSPECIESC; ++ipart)
1633     values[ipart]=GetNumberOfSigmas(detector,track,(AliPID::EParticleType)ipart);
1634   // the pid status is the same for probabilities and nSigmas, so it is
1635   // fine to use the one from the probabilities also here
1636   detPID->SetNumberOfSigmas(detector, values, (Int_t)AliPID::kSPECIESC, status);
1637   
1638 }
1639
1640 //______________________________________________________________________________
1641 void AliPIDResponse::FillTrackDetectorPID()
1642 {
1643   //
1644   // create detector PID information and setup the transient pointer in the track
1645   //
1646
1647   if (!fCurrentEvent) return;
1648   
1649   for (Int_t itrack=0; itrack<fCurrentEvent->GetNumberOfTracks(); ++itrack){
1650     AliVTrack *track=dynamic_cast<AliVTrack*>(fCurrentEvent->GetTrack(itrack));
1651     if (!track) continue;
1652
1653     for (Int_t idet=0; idet<kNdetectors; ++idet){
1654       FillTrackDetectorPID(track, (EDetector)idet);
1655     }
1656   }
1657 }
1658
1659 //______________________________________________________________________________
1660 void AliPIDResponse::SetTOFResponse(AliVEvent *vevent,EStartTimeType_t option){
1661   //
1662   // Set TOF response function
1663   // Input option for event_time used
1664   //
1665   
1666     Float_t t0spread = 0.; //vevent->GetEventTimeSpread();
1667     if(t0spread < 10) t0spread = 80;
1668
1669     // T0 from TOF algorithm
1670
1671     Bool_t flagT0TOF=kFALSE;
1672     Bool_t flagT0T0=kFALSE;
1673     Float_t *startTime = new Float_t[fTOFResponse.GetNmomBins()];
1674     Float_t *startTimeRes = new Float_t[fTOFResponse.GetNmomBins()];
1675     Int_t *startTimeMask = new Int_t[fTOFResponse.GetNmomBins()];
1676
1677     // T0-TOF arrays
1678     Float_t *estimatedT0event = new Float_t[fTOFResponse.GetNmomBins()];
1679     Float_t *estimatedT0resolution = new Float_t[fTOFResponse.GetNmomBins()];
1680     for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
1681       estimatedT0event[i]=0.0;
1682       estimatedT0resolution[i]=0.0;
1683       startTimeMask[i] = 0;
1684     }
1685
1686     Float_t resT0A=fResT0A;
1687     Float_t resT0C=fResT0C;
1688     Float_t resT0AC=fResT0AC;
1689     if(vevent->GetT0TOF()){ // check if T0 detector information is available
1690         flagT0T0=kTRUE;
1691     }
1692
1693
1694     AliTOFHeader *tofHeader = (AliTOFHeader*)vevent->GetTOFHeader();
1695
1696     if (tofHeader) { // read global info and T0-TOF
1697       fTOFResponse.SetTimeResolution(tofHeader->GetTOFResolution());
1698       t0spread = tofHeader->GetT0spread(); // read t0 sprad
1699       if(t0spread < 10) t0spread = 80;
1700
1701       flagT0TOF=kTRUE;
1702       for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){ // read T0-TOF default value
1703         startTime[i]=tofHeader->GetDefaultEventTimeVal();
1704         startTimeRes[i]=tofHeader->GetDefaultEventTimeRes();
1705         if(startTimeRes[i] < 1.e-5) startTimeRes[i] = t0spread;
1706       }
1707
1708       TArrayI *ibin=(TArrayI*)tofHeader->GetNvalues();
1709       TArrayF *t0Bin=(TArrayF*)tofHeader->GetEventTimeValues();
1710       TArrayF *t0ResBin=(TArrayF*)tofHeader->GetEventTimeRes();
1711       for(Int_t j=0;j < tofHeader->GetNbins();j++){ // fill T0-TOF in p-bins
1712         Int_t icurrent = (Int_t)ibin->GetAt(j);
1713         startTime[icurrent]=t0Bin->GetAt(j);
1714         startTimeRes[icurrent]=t0ResBin->GetAt(j);
1715         if(startTimeRes[icurrent] < 1.e-5) startTimeRes[icurrent] = t0spread;
1716       }
1717     }
1718
1719     // for cut of 3 sigma on t0 spread
1720     Float_t t0cut = 3 * t0spread;
1721     if(t0cut < 500) t0cut = 500;
1722
1723     if(option == kFILL_T0){ // T0-FILL is used
1724         for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
1725           estimatedT0event[i]=0.0;
1726           estimatedT0resolution[i]=t0spread;
1727         }
1728         fTOFResponse.SetT0event(estimatedT0event);
1729         fTOFResponse.SetT0resolution(estimatedT0resolution);
1730     }
1731
1732     if(option == kTOF_T0){ // T0-TOF is used when available (T0-FILL otherwise) from ESD
1733         if(flagT0TOF){
1734             fTOFResponse.SetT0event(startTime);
1735             fTOFResponse.SetT0resolution(startTimeRes);
1736             for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
1737               if(startTimeRes[i]<t0spread) startTimeMask[i]=1;
1738               fTOFResponse.SetT0binMask(i,startTimeMask[i]);
1739             }
1740         }
1741         else{
1742             for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
1743               estimatedT0event[i]=0.0;
1744               estimatedT0resolution[i]=t0spread;
1745               fTOFResponse.SetT0binMask(i,startTimeMask[i]);
1746             }
1747             fTOFResponse.SetT0event(estimatedT0event);
1748             fTOFResponse.SetT0resolution(estimatedT0resolution);
1749         }
1750     }
1751     else if(option == kBest_T0){ // T0-T0 or T0-TOF are used when available (T0-FILL otherwise) from ESD
1752         Float_t t0AC=-10000;
1753         Float_t t0A=-10000;
1754         Float_t t0C=-10000;
1755         if(flagT0T0){
1756             t0A= vevent->GetT0TOF()[1];
1757             t0C= vevent->GetT0TOF()[2];
1758         //      t0AC= vevent->GetT0TOF()[0];
1759         t0AC= t0A/resT0A/resT0A + t0C/resT0C/resT0C;
1760         resT0AC= TMath::Sqrt(1./resT0A/resT0A + 1./resT0C/resT0C);
1761         t0AC /= resT0AC*resT0AC;
1762         }
1763
1764         Float_t t0t0Best = 0;
1765         Float_t t0t0BestRes = 9999;
1766         Int_t t0used=0;
1767         if(TMath::Abs(t0A) < t0cut && TMath::Abs(t0C) < t0cut && TMath::Abs(t0C-t0A) < 500){
1768             t0t0Best = t0AC;
1769             t0t0BestRes = resT0AC;
1770             t0used=6;
1771         }
1772         else if(TMath::Abs(t0C) < t0cut){
1773             t0t0Best = t0C;
1774             t0t0BestRes = resT0C;
1775             t0used=4;
1776         }
1777         else if(TMath::Abs(t0A) < t0cut){
1778             t0t0Best = t0A;
1779             t0t0BestRes = resT0A;
1780             t0used=2;
1781         }
1782
1783         if(flagT0TOF){ // if T0-TOF info is available
1784             for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
1785                 if(t0t0BestRes < 999){
1786                   if(startTimeRes[i] < t0spread){
1787                     Double_t wtot = 1./startTimeRes[i]/startTimeRes[i] + 1./t0t0BestRes/t0t0BestRes;
1788                     Double_t t0best = startTime[i]/startTimeRes[i]/startTimeRes[i] + t0t0Best/t0t0BestRes/t0t0BestRes;
1789                     estimatedT0event[i]=t0best / wtot;
1790                     estimatedT0resolution[i]=1./TMath::Sqrt(wtot);
1791                     startTimeMask[i] = t0used+1;
1792                   }
1793                   else {
1794                     estimatedT0event[i]=t0t0Best;
1795                     estimatedT0resolution[i]=t0t0BestRes;
1796                     startTimeMask[i] = t0used;
1797                   }
1798                 }
1799                 else{
1800                   estimatedT0event[i]=startTime[i];
1801                   estimatedT0resolution[i]=startTimeRes[i];
1802                   if(startTimeRes[i]<t0spread) startTimeMask[i]=1;
1803                 }
1804                 fTOFResponse.SetT0binMask(i,startTimeMask[i]);
1805             }
1806             fTOFResponse.SetT0event(estimatedT0event);
1807             fTOFResponse.SetT0resolution(estimatedT0resolution);
1808         }
1809         else{ // if no T0-TOF info is available
1810             for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
1811               fTOFResponse.SetT0binMask(i,t0used);
1812               if(t0t0BestRes < 999){
1813                 estimatedT0event[i]=t0t0Best;
1814                 estimatedT0resolution[i]=t0t0BestRes;
1815               }
1816               else{
1817                 estimatedT0event[i]=0.0;
1818                 estimatedT0resolution[i]=t0spread;
1819               }
1820             }
1821             fTOFResponse.SetT0event(estimatedT0event);
1822             fTOFResponse.SetT0resolution(estimatedT0resolution);
1823         }
1824     }
1825
1826     else if(option == kT0_T0){ // T0-T0 is used when available (T0-FILL otherwise)
1827         Float_t t0AC=-10000;
1828         Float_t t0A=-10000;
1829         Float_t t0C=-10000;
1830         if(flagT0T0){
1831             t0A= vevent->GetT0TOF()[1];
1832             t0C= vevent->GetT0TOF()[2];
1833         //      t0AC= vevent->GetT0TOF()[0];
1834         t0AC= t0A/resT0A/resT0A + t0C/resT0C/resT0C;
1835         resT0AC= TMath::Sqrt(1./resT0A/resT0A + 1./resT0C/resT0C);
1836         t0AC /= resT0AC*resT0AC;
1837         }
1838
1839         if(TMath::Abs(t0A) < t0cut && TMath::Abs(t0C) < t0cut && TMath::Abs(t0C-t0A) < 500){
1840             for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
1841               estimatedT0event[i]=t0AC;
1842               estimatedT0resolution[i]=resT0AC;
1843               fTOFResponse.SetT0binMask(i,6);
1844             }
1845         }
1846         else if(TMath::Abs(t0C) < t0cut){
1847             for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
1848               estimatedT0event[i]=t0C;
1849               estimatedT0resolution[i]=resT0C;
1850               fTOFResponse.SetT0binMask(i,4);
1851             }
1852         }
1853         else if(TMath::Abs(t0A) < t0cut){
1854             for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
1855               estimatedT0event[i]=t0A;
1856               estimatedT0resolution[i]=resT0A;
1857               fTOFResponse.SetT0binMask(i,2);
1858             }
1859         }
1860         else{
1861             for(Int_t i=0;i<fTOFResponse.GetNmomBins();i++){
1862               estimatedT0event[i]=0.0;
1863               estimatedT0resolution[i]=t0spread;
1864               fTOFResponse.SetT0binMask(i,0);
1865             }
1866         }
1867         fTOFResponse.SetT0event(estimatedT0event);
1868         fTOFResponse.SetT0resolution(estimatedT0resolution);
1869     }
1870     delete [] startTime;
1871     delete [] startTimeRes;
1872     delete [] startTimeMask;
1873     delete [] estimatedT0event;
1874     delete [] estimatedT0resolution;
1875 }
1876
1877 //______________________________________________________________________________
1878 // private non cached versions of the PID calculation
1879 //
1880
1881
1882 //______________________________________________________________________________
1883 Float_t AliPIDResponse::GetNumberOfSigmas(EDetector detector, const AliVParticle *vtrack, AliPID::EParticleType type) const
1884 {
1885   //
1886   // NumberOfSigmas for 'detCode'
1887   //
1888
1889   const AliVTrack *track=static_cast<const AliVTrack*>(vtrack);
1890   
1891   switch (detector){
1892     case kITS:   return GetNumberOfSigmasITS(track, type);   break;
1893     case kTPC:   return GetNumberOfSigmasTPC(track, type);   break;
1894     case kTOF:   return GetNumberOfSigmasTOF(track, type);   break;
1895     case kHMPID: return GetNumberOfSigmasHMPID(track, type); break;
1896     case kEMCAL: return GetNumberOfSigmasEMCAL(track, type); break;
1897     default: return -999.;
1898   }
1899
1900   return -999.;
1901 }
1902
1903 //______________________________________________________________________________
1904 Float_t AliPIDResponse::GetNumberOfSigmasITS(const AliVParticle *vtrack, AliPID::EParticleType type) const
1905 {
1906   //
1907   // Calculate the number of sigmas in the ITS
1908   //
1909   
1910   AliVTrack *track=(AliVTrack*)vtrack;
1911
1912   const EDetPidStatus pidStatus=GetITSPIDStatus(track);
1913   if (pidStatus!=kDetPidOk) return -999.;
1914
1915   return fITSResponse.GetNumberOfSigmas(track,type);
1916 }
1917
1918 //______________________________________________________________________________
1919 Float_t AliPIDResponse::GetNumberOfSigmasTPC(const AliVParticle *vtrack, AliPID::EParticleType type) const
1920 {
1921   //
1922   // Calculate the number of sigmas in the TPC
1923   //
1924   
1925   AliVTrack *track=(AliVTrack*)vtrack;
1926
1927   const EDetPidStatus pidStatus=GetTPCPIDStatus(track);
1928   if (pidStatus!=kDetPidOk) return -999.;
1929
1930   // the following call is needed in order to fill the transient data member
1931   // fTPCsignalTuned which is used in the TPCPIDResponse to judge
1932   // if using tuned on data
1933   if (fTuneMConData && ((fTuneMConDataMask & kDetTPC) == kDetTPC))
1934     this->GetTPCsignalTunedOnData(track);
1935   
1936   return fTPCResponse.GetNumberOfSigmas(track, type, AliTPCPIDResponse::kdEdxDefault, fUseTPCEtaCorrection, fUseTPCMultiplicityCorrection);
1937 }
1938
1939 //______________________________________________________________________________
1940 Float_t AliPIDResponse::GetNumberOfSigmasTOF(const AliVParticle *vtrack, AliPID::EParticleType type) const
1941 {
1942   //
1943   // Calculate the number of sigmas in the TOF
1944   //
1945   
1946   AliVTrack *track=(AliVTrack*)vtrack;
1947
1948   const EDetPidStatus pidStatus=GetTOFPIDStatus(track);
1949   if (pidStatus!=kDetPidOk) return -999.;
1950   
1951   return GetNumberOfSigmasTOFold(vtrack, type);
1952 }
1953 //______________________________________________________________________________
1954
1955 Float_t AliPIDResponse::GetNumberOfSigmasHMPID(const AliVParticle *vtrack, AliPID::EParticleType type) const
1956 {
1957   //
1958   // Calculate the number of sigmas in the HMPID
1959   //  
1960   AliVTrack *track=(AliVTrack*)vtrack;
1961     
1962   const EDetPidStatus pidStatus=GetHMPIDPIDStatus(track);
1963   if (pidStatus!=kDetPidOk) return -999.; 
1964   
1965   return fHMPIDResponse.GetNumberOfSigmas(track, type);
1966 }
1967
1968 //______________________________________________________________________________
1969 Float_t AliPIDResponse::GetNumberOfSigmasEMCAL(const AliVParticle *vtrack, AliPID::EParticleType type) const
1970 {
1971   //
1972   // Calculate the number of sigmas in the EMCAL
1973   //
1974   
1975   AliVTrack *track=(AliVTrack*)vtrack;
1976
1977   const EDetPidStatus pidStatus=GetEMCALPIDStatus(track);
1978   if (pidStatus!=kDetPidOk) return -999.;
1979
1980   const Int_t nMatchClus = track->GetEMCALcluster();
1981   AliVCluster *matchedClus = (AliVCluster*)fCurrentEvent->GetCaloCluster(nMatchClus);
1982   
1983   const Double_t mom    = track->P();
1984   const Double_t pt     = track->Pt();
1985   const Int_t    charge = track->Charge();
1986   const Double_t fClsE  = matchedClus->E();
1987   const Double_t EovP   = fClsE/mom;
1988   
1989   return fEMCALResponse.GetNumberOfSigmas(pt,EovP,type,charge);
1990 }
1991
1992 //______________________________________________________________________________
1993 AliPIDResponse::EDetPidStatus AliPIDResponse::GetSignalDeltaITS(const AliVParticle *vtrack, AliPID::EParticleType type, Double_t &val, Bool_t ratio/*=kFALSE*/) const
1994 {
1995   //
1996   // Signal minus expected Signal for ITS
1997   //
1998   AliVTrack *track=(AliVTrack*)vtrack;
1999   val=fITSResponse.GetSignalDelta(track,type,ratio);
2000   
2001   return GetITSPIDStatus(track);
2002 }
2003
2004 //______________________________________________________________________________
2005 AliPIDResponse::EDetPidStatus AliPIDResponse::GetSignalDeltaTPC(const AliVParticle *vtrack, AliPID::EParticleType type, Double_t &val, Bool_t ratio/*=kFALSE*/) const
2006 {
2007   //
2008   // Signal minus expected Signal for TPC
2009   //
2010   AliVTrack *track=(AliVTrack*)vtrack;
2011   
2012   // the following call is needed in order to fill the transient data member
2013   // fTPCsignalTuned which is used in the TPCPIDResponse to judge
2014   // if using tuned on data
2015   if (fTuneMConData && ((fTuneMConDataMask & kDetTPC) == kDetTPC))
2016     this->GetTPCsignalTunedOnData(track);
2017   
2018   val=fTPCResponse.GetSignalDelta(track, type, AliTPCPIDResponse::kdEdxDefault, fUseTPCEtaCorrection, fUseTPCMultiplicityCorrection, ratio);
2019   
2020   return GetTPCPIDStatus(track);
2021 }
2022
2023 //______________________________________________________________________________
2024 AliPIDResponse::EDetPidStatus AliPIDResponse::GetSignalDeltaTOF(const AliVParticle *vtrack, AliPID::EParticleType type, Double_t &val, Bool_t ratio/*=kFALSE*/) const
2025 {
2026   //
2027   // Signal minus expected Signal for TOF
2028   //
2029   AliVTrack *track=(AliVTrack*)vtrack;
2030   val=GetSignalDeltaTOFold(track, type, ratio);
2031   
2032   return GetTOFPIDStatus(track);
2033 }
2034
2035 //______________________________________________________________________________
2036 AliPIDResponse::EDetPidStatus AliPIDResponse::GetSignalDeltaHMPID(const AliVParticle *vtrack, AliPID::EParticleType type, Double_t &val, Bool_t ratio/*=kFALSE*/) const
2037 {
2038   //
2039   // Signal minus expected Signal for HMPID
2040   //
2041   AliVTrack *track=(AliVTrack*)vtrack;
2042   val=fHMPIDResponse.GetSignalDelta(track, type, ratio);
2043   
2044   return GetHMPIDPIDStatus(track);
2045 }
2046
2047 //______________________________________________________________________________
2048 AliPIDResponse::EDetPidStatus AliPIDResponse::GetComputePIDProbability  (EDetector detCode,  const AliVTrack *track, Int_t nSpecies, Double_t p[]) const
2049 {
2050   //
2051   // Compute PID response of 'detCode'
2052   //
2053
2054   switch (detCode){
2055     case kITS: return GetComputeITSProbability(track, nSpecies, p); break;
2056     case kTPC: return GetComputeTPCProbability(track, nSpecies, p); break;
2057     case kTRD: return GetComputeTRDProbability(track, nSpecies, p); break;
2058     case kTOF: return GetComputeTOFProbability(track, nSpecies, p); break;
2059     case kPHOS: return GetComputePHOSProbability(track, nSpecies, p); break;
2060     case kEMCAL: return GetComputeEMCALProbability(track, nSpecies, p); break;
2061     case kHMPID: return GetComputeHMPIDProbability(track, nSpecies, p); break;
2062     default: return kDetNoSignal;
2063   }
2064 }
2065
2066 //______________________________________________________________________________
2067 AliPIDResponse::EDetPidStatus AliPIDResponse::GetComputeITSProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[]) const
2068 {
2069   //
2070   // Compute PID response for the ITS
2071   //
2072   
2073   // set flat distribution (no decision)
2074   for (Int_t j=0; j<nSpecies; j++) p[j]=1./nSpecies;
2075   
2076   const EDetPidStatus pidStatus=GetITSPIDStatus(track);
2077   if (pidStatus!=kDetPidOk) return pidStatus;
2078   
2079   if (track->GetDetectorPID()){
2080     return track->GetDetectorPID()->GetRawProbability(kITS, p, nSpecies);
2081   }
2082   
2083   //check for ITS standalone tracks
2084   Bool_t isSA=kTRUE;
2085   if( track->GetStatus() & AliVTrack::kTPCin ) isSA=kFALSE;
2086
2087   Double_t mom=track->P();
2088   Double_t dedx=track->GetITSsignal();
2089   Double_t momITS=mom;
2090   UChar_t clumap=track->GetITSClusterMap();
2091   Int_t nPointsForPid=0;
2092   for(Int_t i=2; i<6; i++){
2093     if(clumap&(1<<i)) ++nPointsForPid;
2094   }
2095
2096   Bool_t mismatch=kTRUE/*, heavy=kTRUE*/;
2097   for (Int_t j=0; j<nSpecies; j++) {
2098     Double_t mass=AliPID::ParticleMassZ(j);//GeV/c^2
2099     const Double_t chargeFactor = TMath::Power(AliPID::ParticleCharge(j),2.);
2100     Double_t bethe=fITSResponse.Bethe(momITS,mass)*chargeFactor;
2101     //TODO: in case of the electron, use the SA parametrisation,
2102     //      this needs to be changed if ITS provides a parametrisation
2103     //      for electrons also for ITS+TPC tracks
2104     Double_t sigma=fITSResponse.GetResolution(bethe,nPointsForPid,isSA || (j==(Int_t)AliPID::kElectron));
2105     if (TMath::Abs(dedx-bethe) > fRange*sigma) {
2106       p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma;
2107     } else {
2108       p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma;
2109       mismatch=kFALSE;
2110     }
2111   }
2112
2113   if (mismatch){
2114     for (Int_t j=0; j<nSpecies; j++) p[j]=1./nSpecies;
2115   }
2116
2117   return kDetPidOk;
2118 }
2119 //______________________________________________________________________________
2120 AliPIDResponse::EDetPidStatus AliPIDResponse::GetComputeTPCProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[]) const
2121 {
2122   //
2123   // Compute PID response for the TPC
2124   //
2125   
2126   // set flat distribution (no decision)
2127   for (Int_t j=0; j<nSpecies; j++) p[j]=1./nSpecies;
2128   
2129   const EDetPidStatus pidStatus=GetTPCPIDStatus(track);
2130   if (pidStatus!=kDetPidOk) return pidStatus;
2131   
2132   Double_t dedx=track->GetTPCsignal();
2133   Bool_t mismatch=kTRUE/*, heavy=kTRUE*/;
2134   
2135   if (fTuneMConData && ((fTuneMConDataMask & kDetTPC) == kDetTPC)) dedx = this->GetTPCsignalTunedOnData(track);
2136   
2137   Double_t bethe = 0.;
2138   Double_t sigma = 0.;
2139   
2140   for (Int_t j=0; j<nSpecies; j++) {
2141     AliPID::EParticleType type=AliPID::EParticleType(j);
2142     
2143     bethe=fTPCResponse.GetExpectedSignal(track, type, AliTPCPIDResponse::kdEdxDefault, fUseTPCEtaCorrection, fUseTPCMultiplicityCorrection);
2144     sigma=fTPCResponse.GetExpectedSigma(track, type, AliTPCPIDResponse::kdEdxDefault, fUseTPCEtaCorrection, fUseTPCMultiplicityCorrection);
2145     
2146     if (TMath::Abs(dedx-bethe) > fRange*sigma) {
2147       p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma;
2148     } else {
2149       p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma;
2150       mismatch=kFALSE;
2151     }
2152   }
2153   
2154   if (mismatch){
2155     for (Int_t j=0; j<nSpecies; j++) p[j]=1./nSpecies;
2156   }
2157   
2158   return kDetPidOk;
2159 }
2160 //______________________________________________________________________________
2161 AliPIDResponse::EDetPidStatus AliPIDResponse::GetComputeTOFProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[]) const
2162 {
2163   //
2164   // Compute PID probabilities for TOF
2165   //
2166   
2167   // set flat distribution (no decision)
2168   for (Int_t j=0; j<nSpecies; j++) p[j]=1./nSpecies;
2169   
2170   const EDetPidStatus pidStatus=GetTOFPIDStatus(track);
2171   if (pidStatus!=kDetPidOk) return pidStatus;
2172
2173   const Double_t meanCorrFactor = 0.07/fTOFtail; // Correction factor on the mean because of the tail (should be ~ 0.1 with tail = 1.1)
2174   
2175   for (Int_t j=0; j<nSpecies; j++) {
2176     AliPID::EParticleType type=AliPID::EParticleType(j);
2177     const Double_t nsigmas=GetNumberOfSigmasTOFold(track,type) + meanCorrFactor;
2178     
2179     const Double_t expTime = fTOFResponse.GetExpectedSignal(track,type);
2180     const Double_t sig     = fTOFResponse.GetExpectedSigma(track->P(),expTime,AliPID::ParticleMassZ(type));
2181     if (TMath::Abs(nsigmas) > (fRange+2)) {
2182       if(nsigmas < fTOFtail)
2183         p[j] = TMath::Exp(-0.5*(fRange+2)*(fRange+2))/sig;
2184       else
2185         p[j] = TMath::Exp(-(fRange+2 - fTOFtail*0.5)*fTOFtail)/sig;
2186     } else{
2187       if(nsigmas < fTOFtail)
2188         p[j] = TMath::Exp(-0.5*nsigmas*nsigmas)/sig;
2189       else
2190         p[j] = TMath::Exp(-(nsigmas - fTOFtail*0.5)*fTOFtail)/sig;
2191     }    
2192   }
2193   
2194   return kDetPidOk;
2195 }
2196 //______________________________________________________________________________
2197 AliPIDResponse::EDetPidStatus AliPIDResponse::GetComputeTRDProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[],AliTRDPIDResponse::ETRDPIDMethod PIDmethod/*=AliTRDPIDResponse::kLQ1D*/) const
2198 {
2199   //
2200   // Compute PID probabilities for the TRD
2201   //
2202   
2203   // set flat distribution (no decision)
2204   for (Int_t j=0; j<nSpecies; j++) p[j]=1./nSpecies;
2205   
2206   const EDetPidStatus pidStatus=GetTRDPIDStatus(track);
2207   if (pidStatus!=kDetPidOk) return pidStatus;
2208
2209   UInt_t TRDslicesForPID[2];
2210   SetTRDSlices(TRDslicesForPID,PIDmethod);
2211   
2212   Float_t mom[6]={0.};
2213   Double_t dedx[48]={0.};  // Allocate space for the maximum number of TRD slices
2214   Int_t nslices = TRDslicesForPID[1] - TRDslicesForPID[0] + 1;
2215   AliDebug(1, Form("First Slice: %d, Last Slice: %d, Number of slices: %d",  TRDslicesForPID[0], TRDslicesForPID[1], nslices));
2216   for(UInt_t ilayer = 0; ilayer < 6; ilayer++){
2217     mom[ilayer] = track->GetTRDmomentum(ilayer);
2218     for(UInt_t islice = TRDslicesForPID[0]; islice <= TRDslicesForPID[1]; islice++){
2219       dedx[ilayer*nslices+islice-TRDslicesForPID[0]] = track->GetTRDslice(ilayer, islice);
2220     }
2221   }
2222   
2223   fTRDResponse.GetResponse(nslices, dedx, mom, p,PIDmethod);
2224   return kDetPidOk;
2225 }
2226
2227 //______________________________________________________________________________
2228 AliPIDResponse::EDetPidStatus AliPIDResponse::GetComputeEMCALProbability  (const AliVTrack *track, Int_t nSpecies, Double_t p[]) const
2229 {
2230   //
2231   // Compute PID response for the EMCAL
2232   //
2233   
2234   for (Int_t j=0; j<nSpecies; j++) p[j]=1./nSpecies;
2235
2236   const EDetPidStatus pidStatus=GetEMCALPIDStatus(track);
2237   if (pidStatus!=kDetPidOk) return pidStatus;
2238
2239   const Int_t nMatchClus = track->GetEMCALcluster();
2240   AliVCluster *matchedClus = (AliVCluster*)fCurrentEvent->GetCaloCluster(nMatchClus);
2241   
2242   const Double_t mom    = track->P();
2243   const Double_t pt     = track->Pt();
2244   const Int_t    charge = track->Charge();
2245   const Double_t fClsE  = matchedClus->E();
2246   const Double_t EovP   = fClsE/mom;
2247   
2248   // compute the probabilities
2249   fEMCALResponse.ComputeEMCALProbability(nSpecies,pt,EovP,charge,p);
2250   return kDetPidOk;
2251 }
2252
2253 //______________________________________________________________________________
2254 AliPIDResponse::EDetPidStatus AliPIDResponse::GetComputePHOSProbability (const AliVTrack */*track*/, Int_t nSpecies, Double_t p[]) const
2255 {
2256   //
2257   // Compute PID response for the PHOS
2258   //
2259   
2260   // set flat distribution (no decision)
2261   for (Int_t j=0; j<nSpecies; j++) p[j]=1./nSpecies;
2262   return kDetNoSignal;
2263 }
2264
2265 //______________________________________________________________________________
2266 AliPIDResponse::EDetPidStatus AliPIDResponse::GetComputeHMPIDProbability(const AliVTrack *track, Int_t nSpecies, Double_t p[]) const
2267 {
2268   //
2269   // Compute PID response for the HMPID
2270   //
2271   
2272   // set flat distribution (no decision)
2273   for (Int_t j=0; j<nSpecies; j++) p[j]=1./nSpecies;
2274   
2275   const EDetPidStatus pidStatus=GetHMPIDPIDStatus(track);
2276   if (pidStatus!=kDetPidOk) return pidStatus;
2277   
2278   fHMPIDResponse.GetProbability(track,nSpecies,p);
2279     
2280   return kDetPidOk;
2281 }
2282
2283 //______________________________________________________________________________
2284 AliPIDResponse::EDetPidStatus AliPIDResponse::GetITSPIDStatus(const AliVTrack *track) const
2285 {
2286   // compute ITS pid status
2287
2288   // check status bits
2289   if ((track->GetStatus()&AliVTrack::kITSin)==0 &&
2290     (track->GetStatus()&AliVTrack::kITSout)==0) return kDetNoSignal;
2291
2292   const Float_t dEdx=track->GetITSsignal();
2293   if (dEdx<=0) return kDetNoSignal;
2294   
2295   // requite at least 3 pid clusters
2296   const UChar_t clumap=track->GetITSClusterMap();
2297   Int_t nPointsForPid=0;
2298   for(Int_t i=2; i<6; i++){
2299     if(clumap&(1<<i)) ++nPointsForPid;
2300   }
2301   
2302   if(nPointsForPid<3) { 
2303     return kDetNoSignal;
2304   }
2305   
2306   return kDetPidOk;
2307 }
2308
2309 //______________________________________________________________________________
2310 AliPIDResponse::EDetPidStatus AliPIDResponse:: GetTPCPIDStatus(const AliVTrack *track) const
2311 {
2312   // compute TPC pid status
2313   
2314   // check quality of the track
2315   if ( (track->GetStatus()&AliVTrack::kTPCin )==0 && (track->GetStatus()&AliVTrack::kTPCout)==0 ) return kDetNoSignal;
2316
2317   // check pid values
2318   const Double_t dedx=track->GetTPCsignal();
2319   const UShort_t signalN=track->GetTPCsignalN();
2320   if (signalN<10 || dedx<10) return kDetNoSignal;
2321
2322   if (!(fArrPidResponseMaster && fArrPidResponseMaster->At(AliPID::kPion))) return kDetNoParams;
2323   
2324   return kDetPidOk;
2325 }
2326
2327 //______________________________________________________________________________
2328 AliPIDResponse::EDetPidStatus AliPIDResponse::GetTRDPIDStatus(const AliVTrack *track) const
2329 {
2330   // compute TRD pid status
2331
2332   if((track->GetStatus()&AliVTrack::kTRDout)==0) return kDetNoSignal;
2333   return kDetPidOk;
2334 }
2335
2336 //______________________________________________________________________________
2337 AliPIDResponse::EDetPidStatus AliPIDResponse::GetTOFPIDStatus(const AliVTrack *track) const
2338 {
2339   // compute TOF pid status
2340
2341   if ((track->GetStatus()&AliVTrack::kTOFout)==0) return kDetNoSignal;
2342   if ((track->GetStatus()&AliVTrack::kTIME)==0) return kDetNoSignal;
2343
2344   return kDetPidOk;
2345 }
2346
2347 //______________________________________________________________________________
2348 Float_t AliPIDResponse::GetTOFMismatchProbability(const AliVTrack *track) const
2349 {
2350   // compute mismatch probability cross-checking at 5 sigmas with TPC
2351   // currently just implemented as a 5 sigma compatibility cut
2352
2353   // check pid status
2354   const EDetPidStatus tofStatus=GetTOFPIDStatus(track);
2355   if (tofStatus!=kDetPidOk) return 0.;
2356
2357   //mismatch
2358   const EDetPidStatus tpcStatus=GetTPCPIDStatus(track);
2359   if (tpcStatus!=kDetPidOk) return 0.;
2360   
2361   const Double_t meanCorrFactor = 0.11/fTOFtail; // Correction factor on the mean because of the tail (should be ~ 0.1 with tail = 1.1)
2362   Bool_t mismatch = kTRUE/*, heavy = kTRUE*/;
2363   for (Int_t j=0; j<AliPID::kSPECIESC; j++) {
2364     AliPID::EParticleType type=AliPID::EParticleType(j);
2365     const Double_t nsigmas=GetNumberOfSigmasTOFold(track,type) + meanCorrFactor;
2366     
2367     if (TMath::Abs(nsigmas)<5.){
2368       const Double_t nsigmasTPC=GetNumberOfSigmasTPC(track,type);
2369       if (TMath::Abs(nsigmasTPC)<5.) mismatch=kFALSE;
2370     }
2371   }
2372   
2373   if (mismatch){
2374     return 1.;
2375   }
2376   
2377   return 0.;
2378 }
2379
2380 //______________________________________________________________________________
2381 AliPIDResponse::EDetPidStatus AliPIDResponse:: GetHMPIDPIDStatus(const AliVTrack *track) const
2382 {
2383   // compute HMPID pid status
2384   
2385   Int_t ch = track->GetHMPIDcluIdx()/1000000;
2386   Double_t HMPIDsignal = track->GetHMPIDsignal(); 
2387   
2388   if((track->GetStatus()&AliVTrack::kHMPIDpid)==0 || ch<0 || ch>6 || HMPIDsignal<0) return kDetNoSignal;
2389   
2390   return kDetPidOk;
2391 }
2392
2393 //______________________________________________________________________________
2394 AliPIDResponse::EDetPidStatus AliPIDResponse:: GetPHOSPIDStatus(const AliVTrack */*track*/) const
2395 {
2396   // compute PHOS pid status
2397   return kDetNoSignal;  
2398 }
2399
2400 //______________________________________________________________________________
2401 AliPIDResponse::EDetPidStatus AliPIDResponse:: GetEMCALPIDStatus(const AliVTrack *track) const
2402 {
2403   // compute EMCAL pid status
2404
2405
2406   // Track matching
2407   const Int_t nMatchClus = track->GetEMCALcluster();
2408   if (nMatchClus<0) return kDetNoSignal;
2409
2410   AliVCluster *matchedClus = (AliVCluster*)fCurrentEvent->GetCaloCluster(nMatchClus);
2411
2412   if (!(matchedClus && matchedClus->IsEMCAL())) return kDetNoSignal;
2413
2414   const Int_t charge = track->Charge();
2415   if (TMath::Abs(charge)!=1) return kDetNoSignal;
2416
2417   if (!(fEMCALPIDParams && fEMCALPIDParams->At(AliPID::kElectron))) return kDetNoParams;
2418   
2419   return kDetPidOk;
2420
2421 }
2422
2423 //______________________________________________________________________________
2424 AliPIDResponse::EDetPidStatus AliPIDResponse::GetPIDStatus(EDetector detector, const AliVTrack *track) const
2425 {
2426   //
2427   // check pid status for a track
2428   //
2429
2430   switch (detector){
2431     case kITS:   return GetITSPIDStatus(track);   break;
2432     case kTPC:   return GetTPCPIDStatus(track);   break;
2433     case kTRD:   return GetTRDPIDStatus(track);   break;
2434     case kTOF:   return GetTOFPIDStatus(track);   break;
2435     case kPHOS:  return GetPHOSPIDStatus(track);  break;
2436     case kEMCAL: return GetEMCALPIDStatus(track); break;
2437     case kHMPID: return GetHMPIDPIDStatus(track); break;
2438     default: return kDetNoSignal;
2439   }
2440   return kDetNoSignal;
2441   
2442 }