]>
Commit | Line | Data |
---|---|---|
b8bfee30 | 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 | AliPIDCombined::AliPIDCombined(): | |
36 | TNamed("CombinedPID","CombinedPID"), | |
37 | fDetectorMask(0), | |
38 | fRejectMismatchMask(0x7F), | |
39 | fEnablePriors(kTRUE), | |
40 | fSelectedSpecies(AliPID::kSPECIES) | |
41 | { | |
42 | // | |
43 | // default constructor | |
44 | // | |
45 | for (Int_t i=0;i<AliPID::kSPECIES+AliPID::kSPECIESLN;i++) fPriorsDistributions[i]=NULL; | |
46 | AliLog::SetClassDebugLevel("AliPIDCombined",10); | |
47 | } | |
48 | ||
49 | //------------------------------------------------------------------------------------------------- | |
50 | AliPIDCombined::AliPIDCombined(const TString& name,const TString& title): | |
51 | TNamed(name,title), | |
52 | fDetectorMask(0), | |
53 | fRejectMismatchMask(0x7F), | |
54 | fEnablePriors(kTRUE), | |
55 | fSelectedSpecies(AliPID::kSPECIES) | |
56 | { | |
57 | // | |
58 | // default constructor with name and title | |
59 | // | |
60 | for (Int_t i=0;i<AliPID::kSPECIES+AliPID::kSPECIESLN;i++) fPriorsDistributions[i]=NULL; | |
61 | AliLog::SetClassDebugLevel("AliPIDCombined",10); | |
62 | ||
63 | } | |
64 | ||
65 | //------------------------------------------------------------------------------------------------- | |
66 | AliPIDCombined::~AliPIDCombined() { | |
67 | ||
68 | for(Int_t i=0;i < AliPID::kSPECIES+AliPID::kSPECIESLN;i++){ | |
69 | if(fPriorsDistributions[i]) | |
70 | delete fPriorsDistributions[i]; | |
71 | } | |
72 | } | |
73 | ||
74 | //------------------------------------------------------------------------------------------------- | |
75 | void AliPIDCombined::SetPriorDistribution(AliPID::EParticleType type,TH1F *prior) { | |
76 | if ( (type < 0) || ( type >= (AliPID::kSPECIESN+AliPID::kSPECIESLN) ) ){ | |
77 | AliError(Form("Invalid EParticleType setting prior (offending type: %d)",type)); | |
78 | return; | |
79 | } | |
80 | if(prior) { | |
81 | Int_t i = type; | |
82 | if (type >= (AliPID::EParticleType)AliPID::kSPECIES ) { | |
83 | if (type < AliPID::kDeuteron) { | |
84 | AliError(Form("Invalid EParticleType setting prior. Type: %d (neutral) not supported",type)); | |
85 | return; | |
86 | } else i = (Int_t)type - (AliPID::kSPECIESN-AliPID::kSPECIES); | |
87 | } | |
88 | if (fPriorsDistributions[i]) { | |
89 | delete fPriorsDistributions[i]; | |
90 | } | |
91 | fPriorsDistributions[i]=new TH1F(*prior); | |
92 | } | |
93 | } | |
94 | ||
95 | //------------------------------------------------------------------------------------------------- | |
96 | UInt_t AliPIDCombined::ComputeProbabilities(const AliVTrack *track, const AliPIDResponse *response, Double_t* bayesProbabilities) const { | |
97 | Double_t pITS[fSelectedSpecies]; | |
98 | Double_t pTPC[fSelectedSpecies]; | |
99 | Double_t pTRD[fSelectedSpecies]; | |
100 | Double_t pTOF[fSelectedSpecies]; | |
101 | Double_t pHMPID[fSelectedSpecies]; | |
102 | Double_t pEMCAL[fSelectedSpecies]; | |
103 | Double_t pPHOS[fSelectedSpecies]; | |
104 | for (Int_t i=0;i<fSelectedSpecies;i++) { | |
105 | pITS[i]=1.; | |
106 | pTPC[i]=1.; | |
107 | pTRD[i]=1.; | |
108 | pTOF[i]=1.; | |
109 | pHMPID[i]=1.; | |
110 | pEMCAL[i]=1.; | |
111 | pPHOS[i]=1.; | |
112 | } | |
113 | UInt_t usedMask=0; | |
114 | AliPIDResponse::EDetPidStatus status=AliPIDResponse::kDetNoSignal; | |
115 | Double_t p[fSelectedSpecies]; | |
116 | ||
117 | // getting probability distributions for selected detectors only | |
118 | if (fDetectorMask & AliPIDResponse::kDetITS) { | |
119 | status = response->ComputeITSProbability(track, fSelectedSpecies, pITS); | |
120 | SetCombinedStatus(status,&usedMask,AliPIDResponse::kDetITS,pITS); | |
121 | } | |
122 | ||
123 | if (fDetectorMask & AliPIDResponse::kDetTPC) { | |
124 | status = response->ComputeTPCProbability(track, fSelectedSpecies, pTPC); | |
125 | SetCombinedStatus(status,&usedMask,AliPIDResponse::kDetTPC,pTPC); | |
126 | } | |
127 | ||
128 | ||
129 | if (fDetectorMask & AliPIDResponse::kDetTRD) { | |
130 | status = response->ComputeTRDProbability(track, fSelectedSpecies, pTRD); | |
131 | SetCombinedStatus(status,&usedMask,AliPIDResponse::kDetTRD,pTRD); | |
132 | } | |
133 | ||
134 | if (fDetectorMask & AliPIDResponse::kDetTOF) { | |
135 | status = response->ComputeTOFProbability(track, fSelectedSpecies, pTOF); | |
136 | SetCombinedStatus(status,&usedMask,AliPIDResponse::kDetTOF,pTOF); | |
137 | } | |
138 | ||
139 | if (fDetectorMask & AliPIDResponse::kDetHMPID) { | |
140 | status = response->ComputeHMPIDProbability(track, fSelectedSpecies, pHMPID); | |
141 | SetCombinedStatus(status,&usedMask,AliPIDResponse::kDetHMPID,pHMPID); | |
142 | } | |
143 | ||
144 | ||
145 | if (fDetectorMask & AliPIDResponse::kDetEMCAL) { | |
146 | status = response->ComputeEMCALProbability(track, fSelectedSpecies, pEMCAL); | |
147 | SetCombinedStatus(status,&usedMask,AliPIDResponse::kDetEMCAL,pEMCAL); | |
148 | } | |
149 | ||
150 | ||
151 | if (fDetectorMask & AliPIDResponse::kDetPHOS) { | |
152 | status = response->ComputePHOSProbability(track, fSelectedSpecies, pPHOS); | |
153 | SetCombinedStatus(status,&usedMask,AliPIDResponse::kDetPHOS,pPHOS); | |
154 | } | |
155 | ||
156 | ||
157 | for (Int_t i =0; i<fSelectedSpecies; i++){ | |
158 | p[i] = pITS[i]*pTPC[i]*pTRD[i]*pTOF[i]*pHMPID[i]*pEMCAL[i]*pPHOS[i]; | |
159 | } | |
160 | Double_t priors[fSelectedSpecies]; | |
161 | if (fEnablePriors) GetPriors(track,priors); | |
162 | else { for (Int_t i=0;i<fSelectedSpecies;i++) priors[i]=1.;} | |
163 | ComputeBayesProbabilities(bayesProbabilities,p,priors); | |
164 | return usedMask; | |
165 | } | |
166 | ||
167 | ||
168 | //------------------------------------------------------------------------------------------------- | |
169 | void AliPIDCombined::GetPriors(const AliVTrack *track, Double_t* p) const { | |
170 | ||
171 | // | |
172 | // get priors from distributions | |
173 | // | |
174 | ||
175 | Double_t pt=TMath::Abs(track->Pt()); | |
176 | Double_t sumPriors = 0; | |
177 | for (Int_t i=0;i<fSelectedSpecies;++i){ | |
178 | if (fPriorsDistributions[i]){ | |
179 | p[i]=fPriorsDistributions[i]->Interpolate(pt); | |
180 | } | |
181 | else { | |
182 | p[i]=0.; | |
183 | } | |
184 | sumPriors+=p[i]; | |
185 | } | |
186 | ||
187 | // normalizing priors | |
188 | if (sumPriors == 0) return; | |
189 | for (Int_t i=0;i<fSelectedSpecies;++i){ | |
190 | p[i]=p[i]/sumPriors; | |
191 | } | |
192 | return; | |
193 | } | |
194 | ||
195 | //------------------------------------------------------------------------------------------------- | |
196 | void AliPIDCombined::ComputeBayesProbabilities(Double_t* probabilities, const Double_t* probDensity, const Double_t* prior) const { | |
197 | ||
198 | ||
199 | // | |
200 | // calculate Bayesian probabilities | |
201 | // | |
202 | Double_t sum = 0.; | |
203 | for (Int_t i = 0; i < fSelectedSpecies; i++) { | |
204 | sum += probDensity[i] * prior[i]; | |
205 | } | |
206 | if (sum <= 0) { | |
207 | AliError("Invalid probability densities or priors"); | |
208 | for (Int_t i = 0; i < fSelectedSpecies; i++) probabilities[i] = -1; | |
209 | return; | |
210 | } | |
211 | for (Int_t i = 0; i < fSelectedSpecies; i++) { | |
212 | probabilities[i] = probDensity[i] * prior[i] / sum; | |
213 | } | |
214 | ||
215 | ||
216 | } | |
217 | ||
218 | //---------------------------------------------------------------------------------------- | |
219 | void AliPIDCombined::SetCombinedStatus(AliPIDResponse::EDetPidStatus status, UInt_t *maskDetIn, AliPIDResponse::EDetCode bit, Double_t* p) const { | |
220 | switch (status) { | |
221 | case AliPIDResponse::kDetNoSignal: | |
222 | break; | |
223 | case AliPIDResponse::kDetPidOk: | |
224 | *maskDetIn = *maskDetIn | bit; | |
225 | break; | |
226 | case AliPIDResponse::kDetMismatch: | |
227 | if ( fRejectMismatchMask & bit) for (Int_t j=0;j<fSelectedSpecies;j++) p[j]=1./fSelectedSpecies; | |
228 | else *maskDetIn = *maskDetIn | bit; | |
229 | break; | |
230 | } | |
231 | ||
232 | } | |
233 | ||
234 | ClassImp(AliPIDCombined); | |
235 |