]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/STEERBase/AliPIDCombined.cxx
0ac1c114c28b722ae00a992daa1e84fb94167be4
[u/mrichter/AliRoot.git] / STEER / STEERBase / AliPIDCombined.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
17 //-----------------------------------------------------------------
18 //        Base class for combining PID of different detectors    //
19 //        (user selected) and compute Bayesian probabilities     //
20 //                                                               //
21 //                                                               //
22 //   Origin: Pietro Antonioli, INFN-BO Pietro.Antonioli@cern.ch  //
23 //                                                               //
24 //-----------------------------------------------------------------
25
26 #include <TH1.h>
27
28 #include <AliVTrack.h>
29 #include <AliLog.h>
30 #include <AliPID.h>
31 #include <AliPIDResponse.h>
32
33 #include "AliPIDCombined.h"
34
35 #include "TMath.h"
36 #include "TFile.h"
37
38 #include "AliOADBContainer.h"
39
40 // initialize static members
41 TH2F* AliPIDCombined::fDefaultPriorsTPC[]={0x0};
42 TH2F* AliPIDCombined::fDefaultPriorsTPCpPb[]={0x0};
43 Float_t AliPIDCombined::fTOFmismProb = 0;
44
45 ClassImp(AliPIDCombined);
46
47 AliPIDCombined::AliPIDCombined():
48         TNamed("CombinedPID","CombinedPID"),
49         fDetectorMask(0),
50         fRejectMismatchMask(0x7F),
51         fEnablePriors(kTRUE),
52         fSelectedSpecies(AliPID::kSPECIES),
53         fUseDefaultTPCPriors(kFALSE)
54 {       
55   //
56   // default constructor
57   //    
58   for (Int_t i=0;i<AliPID::kSPECIESC;i++) fPriorsDistributions[i]=NULL;
59   AliLog::SetClassDebugLevel("AliPIDCombined",10);
60 }
61
62 //-------------------------------------------------------------------------------------------------     
63 AliPIDCombined::AliPIDCombined(const TString& name,const TString& title):
64         TNamed(name,title),
65         fDetectorMask(0),
66         fRejectMismatchMask(0x7F),
67         fEnablePriors(kTRUE),
68         fSelectedSpecies(AliPID::kSPECIES),
69         fUseDefaultTPCPriors(kFALSE)
70 {
71   //
72   // default constructor with name and title
73   //
74   for (Int_t i=0;i<AliPID::kSPECIESC;i++) fPriorsDistributions[i]=NULL;
75   AliLog::SetClassDebugLevel("AliPIDCombined",10);
76
77 }
78
79 //-------------------------------------------------------------------------------------------------     
80 AliPIDCombined::~AliPIDCombined() {
81
82   for(Int_t i=0;i < AliPID::kSPECIESC;i++){
83     if(fPriorsDistributions[i])
84       delete fPriorsDistributions[i];
85   }
86 }
87
88 //-------------------------------------------------------------------------------------------------     
89 void AliPIDCombined::SetPriorDistribution(AliPID::EParticleType type,TH1F *prior) {
90   if ( (type < 0) || ( type >= ((AliPID::EParticleType)AliPID::kSPECIESC) ) ){
91     AliError(Form("Invalid EParticleType setting prior (offending type: %d)",type));
92     return;
93   }
94   if(prior) {
95     Int_t i = (Int_t)type;
96     if (fPriorsDistributions[i]) {
97       delete fPriorsDistributions[i]; 
98     }
99     fPriorsDistributions[i]=new TH1F(*prior);
100   }
101 }
102
103 //-------------------------------------------------------------------------------------------------     
104 UInt_t AliPIDCombined::ComputeProbabilities(const AliVTrack *track, const AliPIDResponse *response, Double_t* bayesProbabilities,Double_t *priorsOwn) const {
105   //
106   // (1) Get raw probabilities of requested detectors and combine
107   // (2) priors passed as argument
108   // (3) Compute Bayes probabilities
109   //
110
111
112   // (1) Get raw probabilities of selected detectors and combine
113   UInt_t usedMask=0;             // detectors actually used for track
114   fTOFmismProb = 0; // reset TOF mismatch weights
115
116   AliPIDResponse::EDetPidStatus status=AliPIDResponse::kDetNoSignal;
117   Double_t p[AliPID::kSPECIESC];  // combined probabilities of selected detectors
118   Double_t pMismTOF[AliPID::kSPECIESC];  // combined TOF mismatch probabilities using selected detectors
119   for (Int_t i=0;i<fSelectedSpecies;i++){ p[i]=1.;pMismTOF[i]=1.;} // no decision
120   for (Int_t ibit = 0; ibit < 7 ; ibit++) {
121     AliPIDResponse::EDetCode detBit = (AliPIDResponse::EDetCode)(1<<ibit);
122     if (fDetectorMask & detBit) {           // getting probabilities for requested detectors only
123       Double_t detProb[AliPID::kSPECIESC];
124       status = response->ComputePIDProbability(detBit,track,fSelectedSpecies,detProb);
125       if (status == AliPIDResponse::kDetPidOk) {
126         if (fRejectMismatchMask & detBit) {         // mismatch check (currently just for TOF)
127           if (detBit == AliPIDResponse::kDetTOF) {
128             fTOFmismProb = response->GetTOFMismatchProbability(); // mismatch weights computed with TOF probs (no arguments)
129             //Float_t probMis = response->GetTOFMismatchProbability(track); // mismatch compatibility TPC-TOF cut
130             SetCombinedStatus(status,&usedMask,detBit,detProb,fTOFmismProb);
131           } else {
132             SetCombinedStatus(status,&usedMask,detBit);
133           }
134         } else {
135           SetCombinedStatus(status,&usedMask,detBit);
136         }
137         for (Int_t i=0;i<fSelectedSpecies;i++){
138           p[i] *= detProb[i];
139           if(detBit == AliPIDResponse::kDetTOF){
140             Float_t pt = track->Pt();
141             Float_t mismPropagationFactor[] = {1.,1.,1.,1 + TMath::Exp(1 - 1.12*pt),1 + 1./(4.71114 - 5.72372*pt + 2.94715*pt*pt),1.,1.,1.,1.,1.}; // correction for kaons and protons which has to be alligned with the one in AliPIDResponse
142             pMismTOF[i] *= fTOFmismProb*mismPropagationFactor[i];
143           }
144           else pMismTOF[i] *= detProb[i];
145         }
146       }
147     }
148   }
149   // if no detectors available there is no point to go further
150   if (usedMask == 0) return usedMask;
151
152   // (2) Get priors and propagate depending on detectors used
153   Double_t priors[AliPID::kSPECIESC];
154   memset(priors,0,fSelectedSpecies*sizeof(Double_t));
155
156   if(priorsOwn){
157     for(Int_t ipr=0;ipr < fSelectedSpecies;ipr++)
158       priors[ipr] = priorsOwn[ipr];
159   }
160
161   else if (fEnablePriors){
162     Bool_t isPPB = (response->GetBeamType() == AliPIDResponse::kPPB);
163     GetPriors(track,priors,response->GetCurrentCentrality(),isPPB);
164
165     // We apply the propagation factors of the more external detector
166     // 
167     // TPC+HMPID    --> apply HMPID propagation factors (TRD and TOF may be present)
168     // TPC+EMCAL    --> apply EMCAL propagation factors (TRD and TOF may be present)
169     // TPC+TOF      --> apply TOF propagation factors (TRD may be present, HMPID and EMCAL not (if requested))
170     // TPC+TRD      --> apply TRD propagation factors (TOF, HMPID and EMCAL not present (if requested) )
171     // 
172     if(fUseDefaultTPCPriors) {
173
174       Double_t pt=TMath::Abs(track->Pt());
175       if ( ( (usedMask & AliPIDResponse::kDetEMCAL)==AliPIDResponse::kDetEMCAL) || ( (usedMask & AliPIDResponse::kDetHMPID)==AliPIDResponse::kDetHMPID) ) {
176       // we assume EMCAL and HMPID cannot be simultaneously present
177       if ( (usedMask & AliPIDResponse::kDetEMCAL)==AliPIDResponse::kDetEMCAL ) {
178         // EMCal case (for the moment only in combination with TPC)
179         // propagation factors determined from LHC11d MC (LHC12a15f)
180         // v2 clusterizer, dEta < 0.015, dPhi < 0.03, NonLinearityFunction = 6
181                
182         Float_t electronEMCALfactor=0.1;
183         Float_t kaonEMCALfactor = 0.1;
184         Float_t protonEMCALfactor = 0.1;
185                
186         if(track->Charge() > 0){
187           // positiv charge (start parametrization at 0.75 GeV/c and stop at 20 GeV/c for the moment)
188           if(pt > 0.75 && pt < 20.0){
189             electronEMCALfactor =  ( 0.214159 * ( 1 - TMath::Exp(-TMath::Power(pt,0.484512)/0.700499)*TMath::Power(pt,-0.669644)) ) /  ( 0.210436 * ( 1 - TMath::Exp(-TMath::Power(pt,-0.219228)/0.947432)*TMath::Power(pt,-0.700792)) );
190             kaonEMCALfactor =  ( 0.208686 * ( 1 - TMath::Exp(-TMath::Power(pt,-3.98149e-05)/1.28447)*TMath::Power(pt,-0.629191)) ) /  ( 0.210436 * ( 1 - TMath::Exp(-TMath::Power(pt,-0.219228)/0.947432)*TMath::Power(pt,-0.700792)) );
191             protonEMCALfactor =  ( 0.27555 * ( 1 - TMath::Exp(-TMath::Power(pt,-1.97226e-05)/1.52719)*TMath::Power(pt,-0.209068)) ) /  ( 0.210436 * ( 1 - TMath::Exp(-TMath::Power(pt,-0.219228)/0.947432)*TMath::Power(pt,-0.700792)) );
192                    
193           }
194         }
195                
196         if(track->Charge() < 0){
197           // negative charge  (start parametrization at 0.75 GeV/c and stop at 20 GeV/c for the moment)
198           if(pt > 0.75 && pt < 20.0){              
199             electronEMCALfactor =  ( 0.216895 * ( 1 - TMath::Exp(-TMath::Power(pt,0.000105924)/0.865938)*TMath::Power(pt,-1.32787)) ) /  ( 0.210385 * ( 1 - TMath::Exp(-TMath::Power(pt,4.41206e-07)/1.08984)*TMath::Power(pt,-0.544375)) );
200             kaonEMCALfactor =  ( 0.204117 * ( 1 - TMath::Exp(-TMath::Power(pt,-1.6853e-05)/1.61765)*TMath::Power(pt,-0.738355)) ) /  ( 0.210385 * ( 1 - TMath::Exp(-TMath::Power(pt,4.41206e-07)/1.08984)*TMath::Power(pt,-0.544375)) );
201             protonEMCALfactor =  ( 0.215679 * ( 1 - TMath::Exp(-TMath::Power(pt,-4.10015e-05)/1.40921)*TMath::Power(pt,-0.533752)) ) /  ( 0.210385 * ( 1 - TMath::Exp(-TMath::Power(pt,4.41206e-07)/1.08984)*TMath::Power(pt,-0.544375)) );
202           }
203         }
204         priors[0] *= electronEMCALfactor;
205         priors[3] *= kaonEMCALfactor;
206         priors[4] *= protonEMCALfactor;       
207       } // end of EMCAL case
208
209       else if ( (usedMask & AliPIDResponse::kDetHMPID)==AliPIDResponse::kDetHMPID ) {  // HMPID case
210
211         Float_t kaonHMPIDfactor   = 0.;
212         Float_t protonHMPIDfactor = 0.;                     
213         if(pt>1. && pt<6.)  kaonHMPIDfactor   = (-0.0729337 + pt*0.0999531 - pt*pt*0.0371803 + pt*pt*pt*0.00706436 - pt*pt*pt*pt*0.000643619 + pt*pt*pt*pt*pt*2.21853e-05)/(-0.00896231+ pt*0.0330702 - pt*pt*0.0109562+ pt*pt*pt*0.00232895 - pt*pt*pt*pt*0.000246143 + pt*pt*pt*pt*pt*9.59812e-06); 
214         if(pt>1.4 && pt<6.) protonHMPIDfactor = (-0.0444188 + pt*0.0681506 - pt*pt*0.0231819 + pt*pt*pt*0.00400771 - pt*pt*pt*pt*0.000339315 + pt*pt*pt*pt*pt*1.12616e-05)/(-0.00896231+ pt*0.0330702 - pt*pt*0.0109562+ pt*pt*pt*0.00232895 - pt*pt*pt*pt*0.000246143 + pt*pt*pt*pt*pt*9.59812e-06);        
215         if(pt>6. && pt<8.)  kaonHMPIDfactor   = (-0.0729337 + pt*0.0999531 - pt*pt*0.0371803 + pt*pt*pt*0.00706436 - pt*pt*pt*pt*0.000643619 + pt*pt*pt*pt*pt*2.21853e-05)/0.0530456;        
216         if(pt>8.)           kaonHMPIDfactor   = 0.0550432/0.0530456; 
217         if(pt>6. && pt<8.5) protonHMPIDfactor = (-0.0444188 + pt*0.0681506 - pt*pt*0.0231819 + pt*pt*pt*0.00400771 - pt*pt*pt*pt*0.000339315 + pt*pt*pt*pt*pt*1.12616e-05)/0.0530456;      
218         if(pt>8.5)          protonHMPIDfactor = 0.0530071/0.0530456;       
219                                  
220         if(track->Charge() < 0){ 
221          if(pt>0.4 && pt<6.) protonHMPIDfactor = (-0.0351485 + pt*0.0473821 - pt*pt*0.0147947 + pt*pt*pt*0.00254811- pt*pt*pt*pt*0.000224724 + pt*pt*pt*pt*pt*7.9303e-06)/(-0.00896231+ pt*0.0330702 - pt*pt*0.0109562+ pt*pt*pt*0.00232895 - pt*pt*pt*pt*0.000246143 + pt*pt*pt*pt*pt*9.59812e-06);
222          if(pt>6. && pt<8.5) protonHMPIDfactor = (-0.0351485 + pt*0.0473821 - pt*pt*0.0147947 + pt*pt*pt*0.00254811- pt*pt*pt*pt*0.000224724 + pt*pt*pt*pt*pt*7.9303e-06)/0.0530456; 
223          if(pt>8.5)          protonHMPIDfactor = 0.0457756/0.0530456; 
224         } 
225       
226         priors[3] *= kaonHMPIDfactor;
227         priors[4] *= protonHMPIDfactor;
228                
229       }
230
231     } // end of outer cases: EMCAL/HMPID
232       else if ( (usedMask & AliPIDResponse::kDetTOF) == AliPIDResponse::kDetTOF ){
233         Float_t kaonTOFfactor = 0.1;
234         if(pt > 0.29){
235           kaonTOFfactor = 1 - TMath::Exp(-TMath::Power(pt,4.19618E-07)/4.96502e-01)*TMath::Power(pt,-1.50705);
236           if(pt < 0.4) kaonTOFfactor *= 1-1.29854e-04/TMath::Power(pt,7.5);
237         }
238
239         Float_t protonTOFfactor = 0.1;
240         if(pt > 0.4) protonTOFfactor = 1 - TMath::Exp(-TMath::Power(pt,3.30978)/8.57396E-02)*TMath::Power(pt,-4.42661E-01);
241
242         if(track->Charge() < 0){ // for negative tracks
243           kaonTOFfactor *= 1 - TMath::Exp(-TMath::Power(pt,4.87912E-07)/3.26431E-01)*TMath::Power(pt,-1.22893);
244           protonTOFfactor *= 1 - TMath::Exp(-TMath::Power(pt,2.00575E-07)/4.95605E-01)*TMath::Power(pt,-6.71305E-01);
245         }
246         // TODO: we may need an electron factor for TOF as well, especially if TRD is present!
247         priors[3] *= kaonTOFfactor;
248         priors[4] *= protonTOFfactor;
249       }  // end of TOF case
250       else if ( (usedMask & AliPIDResponse::kDetTRD)==AliPIDResponse::kDetTRD ) {
251         Float_t electronTRDfactor=0.1;
252         Float_t kaonTRDfactor = 0.1;
253         Float_t protonTRDfactor = 0.1;
254                  
255         if(track->Charge() > 0){
256           // positiv charge
257           if(pt > 0.25) electronTRDfactor =  1 - TMath::Exp(-TMath::Power(pt,5.13315e-03)/2.11145e-01)*TMath::Power(pt,-2.97659e+00);
258           if(pt > 0.35) kaonTRDfactor = 1 - TMath::Exp(-TMath::Power(pt,-4.29549e-02)/4.87989e-01)*TMath::Power(pt,-1.54270e+00);
259           if(pt > 0.35) protonTRDfactor = 1 - TMath::Exp(-TMath::Power(pt,2.81238e+00)/7.57082e-02)*TMath::Power(pt,-8.12595e-01);
260         }
261                  
262         if(track->Charge() < 0){
263           // negative charge
264           if(pt > 0.25) electronTRDfactor =  1 - TMath::Exp(-TMath::Power(pt,2.45537e-02)/1.90397e-01)*TMath::Power(pt,-3.33121e+00);
265           if(pt > 0.35) kaonTRDfactor = 1 - TMath::Exp(-TMath::Power(pt, -3.42831e-03)/5.57013e-01)*TMath::Power(pt,-1.39202e+00);
266           if(pt > 0.35) protonTRDfactor = 1 - TMath::Exp(-TMath::Power(pt,3.36631e+00)/7.18819e-02)*TMath::Power(pt,-2.00577e-01);
267         }
268         // what about electrons
269         priors[0] *= electronTRDfactor;
270         priors[3] *= kaonTRDfactor;
271         priors[4] *= protonTRDfactor;         
272       } // end of TRD case
273     } // end of fUseDefaultTPCPriors
274   }   // end of use priors
275   else { for (Int_t i=0;i<fSelectedSpecies;i++) priors[i]=1.;}
276
277   //  Compute Bayes probabilities
278   ComputeBayesProbabilities(bayesProbabilities,p,priors,pMismTOF);
279
280   // compute TOF probability contribution from mismatch
281   fTOFmismProb = 0; 
282   for (Int_t i=0;i<fSelectedSpecies;i++) fTOFmismProb += pMismTOF[i];
283
284   return usedMask;
285 }
286
287
288 //-------------------------------------------------------------------------------------------------
289 void AliPIDCombined::GetPriors(const AliVTrack *track, Double_t* p,Float_t centrality,Bool_t isPPB) const {
290         
291         //
292         // get priors from distributions
293         //
294         
295
296         Double_t pt=TMath::Abs(track->Pt());
297
298         if(fUseDefaultTPCPriors){ // use default priors if requested
299           Float_t usedCentr = centrality+5*(centrality>0);
300           if(usedCentr < -0.99) usedCentr = -0.99;
301           else if(usedCentr > 99) usedCentr = 99;
302           if(pt > 9.99) pt = 9.99;
303           else if(pt < 0.01)  pt = 0.01;
304
305           if(! isPPB)
306             for(Int_t i=0;i<AliPID::kSPECIESC;i++) p[i] = fDefaultPriorsTPC[i]->Interpolate(usedCentr,pt);
307           else
308             for(Int_t i=0;i<AliPID::kSPECIESC;i++) p[i] = fDefaultPriorsTPCpPb[i]->Interpolate(usedCentr,pt);
309
310           return;
311         }
312         
313         Double_t sumPriors = 0;
314         for (Int_t i=0;i<fSelectedSpecies;++i){
315                 if (fPriorsDistributions[i]){
316                         p[i]=fPriorsDistributions[i]->Interpolate(pt);
317                 }
318                 else {
319                         p[i]=0.;
320                 }
321                 sumPriors+=p[i];                
322         }
323
324         // normalizing priors
325         if (sumPriors == 0) return;
326         for (Int_t i=0;i<fSelectedSpecies;++i){
327            p[i]=p[i]/sumPriors;
328         }
329         return;
330 }
331 //-------------------------------------------------------------------------------------------------
332 void AliPIDCombined::GetPriors(const AliVTrack *track,Double_t* p,const AliPIDResponse *response,UInt_t detUsed) const{
333         
334         //
335         // get priors from distributions
336         //
337         Double_t pt=TMath::Abs(track->Pt());
338         
339         if(fUseDefaultTPCPriors){ // use default priors if requested
340           Float_t usedCentr = response->GetCurrentCentrality()+5*(response->GetCurrentCentrality()>0);
341           if(usedCentr < -0.99) usedCentr = -0.99;
342           else if(usedCentr > 99) usedCentr = 99;
343           if(pt > 9.99) pt = 9.99;
344           else if(pt < 0.01)  pt = 0.01;
345           
346           for(Int_t i=0;i<AliPID::kSPECIESC;i++) p[i] = fDefaultPriorsTPC[i]->Interpolate(usedCentr,pt);
347
348           // Extra factor if TOF matching was required
349           if(detUsed & AliPIDResponse::kDetTOF){
350             Float_t kaonTOFfactor = 0.1;
351             if(pt > 0.29){
352               kaonTOFfactor = 1 - TMath::Exp(-TMath::Power(pt,4.19618E-07)/4.96502e-01)*TMath::Power(pt,-1.50705);
353               if(pt < 0.4) kaonTOFfactor *= 1-1.29854e-04/TMath::Power(pt,7.5);
354             }
355             Float_t protonTOFfactor = 0.1;
356             if(pt > 0.4) protonTOFfactor = 1 - TMath::Exp(-TMath::Power(pt,3.30978)/8.57396E-02)*TMath::Power(pt,-4.42661E-01);
357             
358             if(track->Charge() < 0){ // for negative tracks
359               kaonTOFfactor *= 1 - TMath::Exp(-TMath::Power(pt,4.87912E-07)/3.26431E-01)*TMath::Power(pt,-1.22893);
360               protonTOFfactor *= 1 - TMath::Exp(-TMath::Power(pt,2.00575E-07)/4.95605E-01)*TMath::Power(pt,-6.71305E-01);
361             }
362             
363             p[3] *= kaonTOFfactor;
364             p[4] *= protonTOFfactor;
365           }
366
367           return;
368         }
369
370
371         Double_t sumPriors = 0;
372         for (Int_t i=0;i<fSelectedSpecies;++i){
373                 if (fPriorsDistributions[i]){
374                         p[i]=fPriorsDistributions[i]->Interpolate(pt);
375                 }
376                 else {
377                         p[i]=0.;
378                 }
379                 sumPriors+=p[i];                
380         }
381
382         // normalizing priors
383         if (sumPriors == 0) return;
384         for (Int_t i=0;i<fSelectedSpecies;++i){
385            p[i]=p[i]/sumPriors;
386         }
387         return;
388 }
389 //-------------------------------------------------------------------------------------------------     
390 void AliPIDCombined::ComputeBayesProbabilities(Double_t* probabilities, const Double_t* probDensity, const Double_t* prior, Double_t* probDensityMism) const {
391
392
393   //
394   // calculate Bayesian probabilities
395   //
396   Double_t sum = 0.;
397   for (Int_t i = 0; i < fSelectedSpecies; i++) {
398     sum += probDensity[i] * prior[i];
399   }
400   if (sum <= 0) {
401
402     AliError("Invalid probability densities or priors");
403     for (Int_t i = 0; i < fSelectedSpecies; i++) probabilities[i] = -1;
404     return;
405   }
406   for (Int_t i = 0; i < fSelectedSpecies; i++) {
407     probabilities[i] = probDensity[i] * prior[i] / sum;
408     if(probDensityMism) probDensityMism[i] *= prior[i] / sum;
409   }
410
411
412 }
413
414
415 //----------------------------------------------------------------------------------------
416 void AliPIDCombined::SetCombinedStatus(AliPIDResponse::EDetPidStatus status, UInt_t *maskDetIn, AliPIDResponse::EDetCode bit) const {
417   switch (status) {
418   case AliPIDResponse::kDetNoParams:
419   case AliPIDResponse::kDetNoSignal:
420   case AliPIDResponse::kDetMismatch: // for backward compatibilty, we need then to remove kDetMismatch from AliPIDResponse
421     break;
422   case AliPIDResponse::kDetPidOk:
423     *maskDetIn = *maskDetIn | bit;
424     break;
425   }
426
427 }
428
429 //----------------------------------------------------------------------------------------
430 void AliPIDCombined::SetCombinedStatus(AliPIDResponse::EDetPidStatus status, UInt_t *maskDetIn, AliPIDResponse::EDetCode bit, Double_t* /*p*/, const Float_t /*probMis*/) const {
431   switch (status) {
432   case AliPIDResponse::kDetNoParams:
433   case AliPIDResponse::kDetNoSignal:
434   case AliPIDResponse::kDetMismatch: // for backward compatibility, we need then to remove kDeteMismatch from AliPIDResponse
435     break;
436   case AliPIDResponse::kDetPidOk:
437     //      if (probMis > 0.5) for (Int_t j=0;j<fSelectedSpecies;j++) p[j]=1./fSelectedSpecies; // no longer used because mismatch is in the framework now
438     //else 
439     *maskDetIn = *maskDetIn | bit;
440     break;
441   }
442
443 }
444
445
446
447 //----------------------------------------------------------------------------------------
448 void AliPIDCombined::SetDefaultTPCPriors(){
449   fEnablePriors=kTRUE;
450   fUseDefaultTPCPriors = kTRUE;
451
452   //check if priors are already initialized
453   if (fDefaultPriorsTPC[0]) return;
454   
455   TString oadbfilename("$ALICE_ROOT/OADB/COMMON/PID/data/");
456   oadbfilename += "/PIDdefaultPriors.root";
457   TFile * foadb = TFile::Open(oadbfilename.Data());
458   if(!foadb || !foadb->IsOpen()) {
459     delete foadb;
460     AliFatal(Form("Cannot open OADB file %s", oadbfilename.Data()));
461     return;
462   }
463   
464   AliOADBContainer * priorsContainer = (AliOADBContainer*) foadb->Get("priorsTPC");
465   if (!priorsContainer) AliFatal("Cannot fetch OADB container for Priors");
466   
467   const char *namespecies[AliPID::kSPECIESC] = {"El","Mu","Pi","Ka","Pr","De","Tr","He3","He3"};
468   char name[100];
469
470   for(Int_t i=0;i < AliPID::kSPECIESC;i++){
471     snprintf(name,99,"hDefault%sPriors",namespecies[i]);
472     fDefaultPriorsTPC[i] = (TH2F*) priorsContainer->GetDefaultObject(name);
473     if (!fDefaultPriorsTPC[i]) AliFatal(Form("Cannot find priors for %s", namespecies[i]));
474     snprintf(name,99,"hDefault%sPriorsPPb",namespecies[i]);
475     fDefaultPriorsTPCpPb[i] = (TH2F*) priorsContainer->GetDefaultObject(name);
476     if (!fDefaultPriorsTPCpPb[i]) {
477         fDefaultPriorsTPCpPb[i] = fDefaultPriorsTPC[i]; // use PbPb ones
478     }
479   }
480
481   delete foadb;
482 }