6ad0bfa0 |
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 | |
b2a60966 |
16 | /* $Id$ */ |
17 | |
6ad0bfa0 |
18 | //_________________________________________________________________________ |
b2a60966 |
19 | // Implementation version v1 of the PHOS particle identifier |
7acf6008 |
20 | // Particle identification based on the |
148b2bba |
21 | // - RCPV: distance from CPV recpoint to EMCA recpoint. |
22 | // - TOF |
23 | // - PCA: Principal Components Analysis.. |
24 | // The identified particle has an identification number corresponding |
25 | // to a 9 bits number: |
bc0c084c |
26 | // -Bit 0 to 2: bit set if RCPV > CpvEmcDistance (each bit corresponds |
148b2bba |
27 | // to a different efficiency-purity point of the photon identification) |
bc0c084c |
28 | // -Bit 3 to 5: bit set if TOF < TimeGate (each bit corresponds |
148b2bba |
29 | // to a different efficiency-purity point of the photon identification) |
30 | // -Bit 6 to 9: bit set if Principal Components are |
50739f15 |
31 | // inside an ellipse defined by the parameters a, b, c, x0 and y0. |
148b2bba |
32 | // (each bit corresponds to a different efficiency-purity point of the |
50739f15 |
33 | // photon identification) |
34 | // The PCA (Principal components analysis) needs a file that contains |
35 | // a previous analysis of the correlations between the particles. This |
bc0c084c |
36 | // file is $ALICE_ROOT/PHOS/PCA8pa15_0.5-100.root. Analysis done for |
50739f15 |
37 | // energies between 0.5 and 100 GeV. |
9fa5f1d0 |
38 | // A calibrated energy is calculated. The energy of the reconstructed |
50739f15 |
39 | // cluster is corrected with the formula A + B * E + C * E^2, whose |
bc0c084c |
40 | // parameters where obtained through the study of the reconstructed |
50739f15 |
41 | // energy distribution of monoenergetic photons. |
a4e98857 |
42 | // |
bc0c084c |
43 | // All the parameters (RCPV(2 rows-3 columns),TOF(1r-3c),PCA(5r-4c) |
50739f15 |
44 | // and calibration(1r-3c))are stored in a file called |
45 | // $ALICE_ROOT/PHOS/Parameters.dat. Each time that AliPHOSPIDv1 is |
bc0c084c |
46 | // initialized, this parameters are copied to a Matrix (9,4), a |
50739f15 |
47 | // TMatrixD object. |
7acf6008 |
48 | // |
a4e98857 |
49 | // use case: |
50739f15 |
50 | // root [0] AliPHOSPIDv1 * p = new AliPHOSPIDv1("galice1.root") |
a4e98857 |
51 | // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated |
50739f15 |
52 | // // reading headers from file galice1.root and create RecParticles |
53 | // TrackSegments and RecPoints are used |
54 | // // set file name for the branch RecParticles |
f0a4c9e9 |
55 | // root [1] p->ExecuteTask("deb all time") |
50739f15 |
56 | // // available options |
57 | // // "deb" - prints # of reconstructed particles |
58 | // // "deb all" - prints # and list of RecParticles |
59 | // // "time" - prints benchmarking results |
7acf6008 |
60 | // |
50739f15 |
61 | // root [2] AliPHOSPIDv1 * p2 = new AliPHOSPIDv1("galice1.root","v1",kTRUE) |
148b2bba |
62 | // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated |
50739f15 |
63 | // //Split mode. |
f0a4c9e9 |
64 | // root [3] p2->ExecuteTask() |
65 | // |
50739f15 |
66 | |
f0a4c9e9 |
67 | |
7acf6008 |
68 | //*-- Author: Yves Schutz (SUBATECH) & Gines Martinez (SUBATECH) & |
148b2bba |
69 | // Gustavo Conesa April 2002 |
50739f15 |
70 | // PCA redesigned by Gustavo Conesa October 2002: |
71 | // The way of using the PCA has changed. Instead of 2 |
72 | // files with the PCA, each one with different energy ranges |
73 | // of application, we use the wide one (0.5-100 GeV), and instead |
bc0c084c |
74 | // of fixing 3 ellipses for different ranges of energy, it has been |
50739f15 |
75 | // studied the dependency of the ellipses parameters with the |
76 | // energy, and they are implemented in the code as a funtion |
77 | // of the energy. |
78 | // |
79 | // |
80 | // |
6ad0bfa0 |
81 | // --- ROOT system --- |
c947e71a |
82 | |
83 | |
84 | // --- Standard library --- |
acb5beb7 |
85 | #include "TFormula.h" |
7acf6008 |
86 | #include "TBenchmark.h" |
148b2bba |
87 | #include "TPrincipal.h" |
c947e71a |
88 | #include "TFile.h" |
e3817e5f |
89 | #include "TSystem.h" |
148b2bba |
90 | |
6ad0bfa0 |
91 | // --- AliRoot header files --- |
c947e71a |
92 | //#include "AliLog.h" |
7acf6008 |
93 | #include "AliGenerator.h" |
e3817e5f |
94 | #include "AliPHOS.h" |
26d4b141 |
95 | #include "AliPHOSPIDv1.h" |
7b7c1533 |
96 | #include "AliPHOSGetter.h" |
6ad0bfa0 |
97 | |
26d4b141 |
98 | ClassImp( AliPHOSPIDv1) |
6ad0bfa0 |
99 | |
1cb7c1ee |
100 | //____________________________________________________________________________ |
101 | AliPHOSPIDv1::AliPHOSPIDv1():AliPHOSPID() |
102 | { |
a4e98857 |
103 | // default ctor |
148b2bba |
104 | |
8d0f3f77 |
105 | InitParameters() ; |
92f521a9 |
106 | fDefaultInit = kTRUE ; |
7acf6008 |
107 | } |
108 | |
581354c5 |
109 | //____________________________________________________________________________ |
88cb7938 |
110 | AliPHOSPIDv1::AliPHOSPIDv1(const AliPHOSPIDv1 & pid ):AliPHOSPID(pid) |
581354c5 |
111 | { |
386aef34 |
112 | // ctor |
581354c5 |
113 | InitParameters() ; |
581354c5 |
114 | Init() ; |
581354c5 |
115 | |
116 | } |
117 | |
7acf6008 |
118 | //____________________________________________________________________________ |
88cb7938 |
119 | AliPHOSPIDv1::AliPHOSPIDv1(const TString alirunFileName, const TString eventFolderName):AliPHOSPID(alirunFileName, eventFolderName) |
7acf6008 |
120 | { |
a4e98857 |
121 | //ctor with the indication on where to look for the track segments |
7b7c1533 |
122 | |
8d0f3f77 |
123 | InitParameters() ; |
2bd5457f |
124 | Init() ; |
92f521a9 |
125 | fDefaultInit = kFALSE ; |
7acf6008 |
126 | } |
7b7c1533 |
127 | |
7acf6008 |
128 | //____________________________________________________________________________ |
129 | AliPHOSPIDv1::~AliPHOSPIDv1() |
130 | { |
79bb1b62 |
131 | // dtor |
acb5beb7 |
132 | fPrincipalPhoton = 0; |
133 | fPrincipalPi0 = 0; |
9fa5f1d0 |
134 | |
e3817e5f |
135 | delete [] fX ; // Principal input |
136 | delete [] fPPhoton ; // Photon Principal components |
137 | delete [] fPPi0 ; // Pi0 Principal components |
acb5beb7 |
138 | |
139 | delete fParameters; |
140 | delete fTFphoton; |
141 | delete fTFpiong; |
142 | delete fTFkaong; |
143 | delete fTFkaonl; |
144 | delete fTFhhadrong; |
145 | delete fTFhhadronl; |
146 | delete fDFmuon; |
7acf6008 |
147 | } |
a496c46c |
148 | //____________________________________________________________________________ |
149 | const TString AliPHOSPIDv1::BranchName() const |
150 | { |
88cb7938 |
151 | |
152 | return GetName() ; |
a496c46c |
153 | } |
154 | |
148b2bba |
155 | //____________________________________________________________________________ |
156 | void AliPHOSPIDv1::Init() |
157 | { |
158 | // Make all memory allocations that are not possible in default constructor |
159 | // Add the PID task to the list of PHOS tasks |
a496c46c |
160 | |
adcca1e6 |
161 | AliPHOSGetter * gime = AliPHOSGetter::Instance() ; |
162 | if(!gime) |
163 | gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data()) ; |
88cb7938 |
164 | |
165 | if ( !gime->PID() ) |
166 | gime->PostPID(this) ; |
148b2bba |
167 | } |
8d0f3f77 |
168 | |
169 | //____________________________________________________________________________ |
170 | void AliPHOSPIDv1::InitParameters() |
171 | { |
e3817e5f |
172 | // Initialize PID parameters |
adcca1e6 |
173 | fWrite = kTRUE ; |
8d0f3f77 |
174 | fRecParticlesInRun = 0 ; |
8d0f3f77 |
175 | fNEvent = 0 ; |
8d0f3f77 |
176 | fRecParticlesInRun = 0 ; |
35adb638 |
177 | fBayesian = kTRUE ; |
9fa5f1d0 |
178 | SetParameters() ; // fill the parameters matrix from parameters file |
eabde521 |
179 | SetEventRange(0,-1) ; |
35adb638 |
180 | |
cc1fe362 |
181 | // initialisation of response function parameters |
182 | // Tof |
183 | // Photons |
35adb638 |
184 | fTphoton[0] = 0.218 ; |
185 | //fTphoton[0] = 1. ; |
186 | fTphoton[1] = 1.55E-8 ; |
187 | fTphoton[2] = 5.05E-10 ; |
188 | fTFphoton = new TFormula("ToF response to photons" , "gaus") ; |
cc1fe362 |
189 | fTFphoton->SetParameters( fTphoton[0], fTphoton[1], fTphoton[2]) ; |
35adb638 |
190 | // // Electrons |
191 | // fTelectron[0] = 0.2 ; |
192 | // fTelectron[1] = 1.55E-8 ; |
193 | // fTelectron[2] = 5.35E-10 ; |
194 | // fTFelectron = new TFormula("ToF response to electrons" , "gaus") ; |
195 | // fTFelectron->SetParameters( fTelectron[0], fTelectron[1], fTelectron[2]) ; |
196 | // // Muons |
197 | // fTmuon[0] = 0.2 ; |
198 | // fTmuon[1] = 1.55E-8 ; |
199 | // fTmuon[2] = 5.1E-10 ; |
200 | // fTFmuon = new TFormula("ToF response to muons" , "gaus") ; |
201 | // fTFmuon->SetParameters( fTmuon[0], fTmuon[1], fTmuon[2]) ; |
202 | |
203 | // Pions |
204 | //Gaus (0 to max probability) |
205 | fTpiong[0] = 0.0971 ; |
206 | //fTpiong[0] = 1. ; |
207 | fTpiong[1] = 1.58E-8 ; |
208 | fTpiong[2] = 5.69E-10 ; |
209 | fTFpiong = new TFormula("ToF response to pions" , "gaus") ; |
210 | fTFpiong->SetParameters( fTpiong[0], fTpiong[1], fTpiong[2]) ; |
211 | // Landau (max probability to inf) |
212 | // fTpionl[0] = 0.05 ; |
213 | // //fTpionl[0] = 5.53 ; |
214 | // fTpionl[1] = 1.68E-8 ; |
215 | // fTpionl[2] = 5.38E-10 ; |
216 | // fTFpionl = new TFormula("ToF response to pions" , "landau") ; |
217 | // fTFpionl->SetParameters( fTpionl[0], fTpionl[1], fTpionl[2]) ; |
218 | |
219 | |
220 | // Kaons |
221 | //Gaus (0 to max probability) |
222 | fTkaong[0] = 0.0542 ; |
223 | //fTkaong[0] = 1. ; |
224 | fTkaong[1] = 1.64E-8 ; |
225 | fTkaong[2] = 6.07-10 ; |
226 | fTFkaong = new TFormula("ToF response to kaon" , "gaus") ; |
227 | fTFkaong->SetParameters( fTkaong[0], fTkaong[1], fTkaong[2]) ; |
228 | //Landau (max probability to inf) |
229 | fTkaonl[0] = 0.264 ; |
230 | //fTkaonl[0] = 5.53 ; |
231 | fTkaonl[1] = 1.68E-8 ; |
232 | fTkaonl[2] = 4.10E-10 ; |
233 | fTFkaonl = new TFormula("ToF response to kaon" , "landau") ; |
234 | fTFkaonl->SetParameters( fTkaonl[0], fTkaonl[1], fTkaonl[2]) ; |
235 | |
236 | //Heavy Hadrons |
237 | //Gaus (0 to max probability) |
238 | fThhadrong[0] = 0.0302 ; |
239 | //fThhadrong[0] = 1. ; |
240 | fThhadrong[1] = 1.73E-8 ; |
241 | fThhadrong[2] = 9.52E-10 ; |
242 | fTFhhadrong = new TFormula("ToF response to heavy hadrons" , "gaus") ; |
243 | fTFhhadrong->SetParameters( fThhadrong[0], fThhadrong[1], fThhadrong[2]) ; |
244 | //Landau (max probability to inf) |
245 | fThhadronl[0] = 0.139 ; |
246 | //fThhadronl[0] = 5.53 ; |
247 | fThhadronl[1] = 1.745E-8 ; |
248 | fThhadronl[2] = 1.00E-9 ; |
249 | fTFhhadronl = new TFormula("ToF response to heavy hadrons" , "landau") ; |
250 | fTFhhadronl->SetParameters( fThhadronl[0], fThhadronl[1], fThhadronl[2]) ; |
251 | |
252 | /// /gaussian parametrization for pions |
253 | // fTpion[0] = 3.93E-2 ; fTpion[1] = 0.130 ; fTpion[2] =-6.37E-2 ;//constant |
254 | // fTpion[3] = 1.65E-8 ; fTpion[4] =-1.40E-9 ; fTpion[5] = 5.96E-10;//mean |
255 | // fTpion[6] = 8.09E-10; fTpion[7] =-4.65E-10; fTpion[8] = 1.50E-10;//sigma |
256 | |
257 | // //landau parametrization for kaons |
258 | // fTkaon[0] = 0.107 ; fTkaon[1] = 0.166 ; fTkaon[2] = 0.243 ;//constant |
259 | // fTkaon[3] = 1.80E-8 ; fTkaon[4] =-2.96E-9 ; fTkaon[5] = 9.60E-10;//mean |
260 | // fTkaon[6] = 1.37E-9 ; fTkaon[7] =-1.80E-9 ; fTkaon[8] = 6.74E-10;//sigma |
261 | |
262 | // //landau parametrization for nucleons |
263 | // fThhadron[0] = 6.33E-2 ; fThhadron[1] = 2.52E-2 ; fThhadron[2] = 2.16E-2 ;//constant |
264 | // fThhadron[3] = 1.94E-8 ; fThhadron[4] =-7.06E-10; fThhadron[5] =-4.69E-10;//mean |
265 | // fThhadron[6] = 2.55E-9 ; fThhadron[7] =-1.90E-9 ; fThhadron[8] = 5.41E-10;//sigma |
266 | |
267 | |
268 | // Shower shape: dispersion gaussian parameters |
269 | // Photons |
270 | |
c947e71a |
271 | // fDphoton[0] = 3.84e-2; fDphoton[1] = 4.46e-3 ; fDphoton[2] = -2.36e-2;//constant |
272 | // //fDphoton[0] = 1.0 ; fDphoton[1] = 0. ; fDphoton[2] = 0. ;//constant |
273 | // fDphoton[3] = 1.55 ; fDphoton[4] =-0.0863 ; fDphoton[5] = 0.287 ;//mean |
274 | // fDphoton[6] = 0.0451 ; fDphoton[7] =-0.0803 ; fDphoton[8] = 0.314 ;//sigma |
275 | |
276 | fDphoton[0] = 4.62e-2; fDphoton[1] = 1.39e-2 ; fDphoton[2] = -3.80e-2;//constant |
277 | //fDphoton[0] = 1.0 ; fDphoton[1] = 0. ; fDphoton[2] = 0. ;//constant |
278 | fDphoton[3] = 1.53 ; fDphoton[4] =-6.62e-2 ; fDphoton[5] = 0.339 ;//mean |
279 | fDphoton[6] = 6.89e-2; fDphoton[7] =-6.59e-2 ; fDphoton[8] = 0.194 ;//sigma |
35adb638 |
280 | |
281 | fDpi0[0] = 0.0586 ; fDpi0[1] = 1.06E-3 ; fDpi0[2] = 0. ;//constant |
282 | //fDpi0[0] = 1.0 ; fDpi0[1] = 0.0 ; fDpi0[2] = 0. ;//constant |
283 | fDpi0[3] = 2.67 ; fDpi0[4] =-2.00E-2 ; fDpi0[5] = 9.37E-5 ;//mean |
284 | fDpi0[6] = 0.153 ; fDpi0[7] = 9.34E-4 ; fDpi0[8] =-1.49E-5 ;//sigma |
285 | //landau |
286 | // fDhadron[0] = 0.007 ; fDhadron[1] = 0. ; fDhadron[2] = 0. ;//constant |
287 | // //fDhadron[0] = 5.53 ; fDhadron[1] = 0. ; fDhadron[2] = 0. ;//constant |
288 | // fDhadron[3] = 3.38 ; fDhadron[4] = 0.0833 ; fDhadron[5] =-0.845 ;//mean |
289 | // fDhadron[6] = 0.627 ; fDhadron[7] = 0.012 ; fDhadron[8] =-0.170 ;//sigma |
290 | |
c947e71a |
291 | fDhadron[0] = 1.61E-2 ; fDhadron[1] = 3.03E-3 ; fDhadron[2] = 1.01E-2 ;//constant |
292 | fDhadron[3] = 3.81 ; fDhadron[4] = 0.232 ; fDhadron[5] =-1.25 ;//mean |
293 | fDhadron[6] = 0.897 ; fDhadron[7] = 0.0987 ; fDhadron[8] =-0.534 ;//sigma |
35adb638 |
294 | // Muons |
295 | fDmuon[0] = 0.0631 ; |
296 | fDmuon[1] = 1.4 ; |
297 | fDmuon[2] = 0.0557 ; |
298 | fDFmuon = new TFormula("Shower shape response to muons" , "landau") ; |
299 | fDFmuon->SetParameters( fDmuon[0], fDmuon[1], fDmuon[2]) ; |
300 | |
35adb638 |
301 | |
c947e71a |
302 | // x(CPV-EMC) distance gaussian parameters |
303 | |
304 | fXelectron[0] = 8.06e-2 ; fXelectron[1] = 1.00e-2; fXelectron[2] =-5.14e-2;//constant |
305 | //fXelectron[0] = 1.0 ; fXelectron[1] = 0. ; fXelectron[2] = 0. ;//constant |
306 | fXelectron[3] = 0.202 ; fXelectron[4] = 8.15e-3; fXelectron[5] = 4.55 ;//mean |
307 | fXelectron[6] = 0.334 ; fXelectron[7] = 0.186 ; fXelectron[8] = 4.32e-2;//sigma |
308 | |
309 | //charged hadrons gaus |
310 | fXcharged[0] = 6.43e-3 ; fXcharged[1] =-4.19e-5; fXcharged[2] = 1.42e-3;//constant |
311 | fXcharged[3] = 2.75 ; fXcharged[4] =-0.40 ; fXcharged[5] = 1.68 ;//mean |
312 | fXcharged[6] = 3.135 ; fXcharged[7] =-9.41e-2; fXcharged[8] = 1.31e-2;//sigma |
313 | |
314 | // z(CPV-EMC) distance gaussian parameters |
315 | |
316 | fZelectron[0] = 8.22e-2 ; fZelectron[1] = 5.11e-3; fZelectron[2] =-3.05e-2;//constant |
317 | //fZelectron[0] = 1.0 ; fZelectron[1] = 0. ; fZelectron[2] = 0. ;//constant |
318 | fZelectron[3] = 3.09e-2 ; fZelectron[4] = 5.87e-2; fZelectron[5] =-9.49e-2;//mean |
319 | fZelectron[6] = 0.263 ; fZelectron[7] =-9.02e-3; fZelectron[8] = 0.151 ;//sigma |
320 | |
321 | //charged hadrons gaus |
322 | |
323 | fZcharged[0] = 1.00e-2 ; fZcharged[1] = 2.82E-4 ; fZcharged[2] = 2.87E-3 ;//constant |
324 | fZcharged[3] =-4.68e-2 ; fZcharged[4] =-9.21e-3 ; fZcharged[5] = 4.91e-2 ;//mean |
325 | fZcharged[6] = 1.425 ; fZcharged[7] =-5.90e-2 ; fZcharged[8] = 5.07e-2 ;//sigma |
326 | |
fb7b51ad |
327 | //Threshold to differentiate between charged and neutral |
328 | fChargedNeutralThreshold = 1e-5; |
329 | |
330 | //Weight to hadrons recontructed energy |
331 | |
332 | fERecWeightPar[0] = 0.32 ; |
333 | fERecWeightPar[1] = 3.8 ; |
334 | fERecWeightPar[2] = 5.4E-3 ; |
335 | fERecWeightPar[3] = 5.6E-2 ; |
336 | fERecWeight = new TFormula("Weight for hadrons" , "[0]*exp(-x*[1])+[2]*exp(-x*[3])") ; |
337 | fERecWeight ->SetParameters(fERecWeightPar[0],fERecWeightPar[1] ,fERecWeightPar[2] ,fERecWeightPar[3]) ; |
338 | |
339 | |
304864ab |
340 | for (Int_t i =0; i< AliPID::kSPECIESN ; i++) |
35adb638 |
341 | fInitPID[i] = 1.; |
fb7b51ad |
342 | |
8d0f3f77 |
343 | } |
344 | |
88cb7938 |
345 | //________________________________________________________________________ |
eabde521 |
346 | void AliPHOSPIDv1::Exec(Option_t *option) |
88cb7938 |
347 | { |
eabde521 |
348 | // Steering method to perform particle reconstruction and identification |
349 | // for the event range from fFirstEvent to fLastEvent. |
350 | // This range is optionally set by SetEventRange(). |
351 | // if fLastEvent=-1 (by default), then process events until the end. |
61d3d6aa |
352 | |
88cb7938 |
353 | if(strstr(option,"tim")) |
354 | gBenchmark->Start("PHOSPID"); |
355 | |
356 | if(strstr(option,"print")) { |
357 | Print() ; |
358 | return ; |
359 | } |
360 | |
361 | |
adcca1e6 |
362 | AliPHOSGetter * gime = AliPHOSGetter::Instance() ; |
88cb7938 |
363 | |
71cee46d |
364 | if (fLastEvent == -1) |
365 | fLastEvent = gime->MaxEvent() - 1 ; |
366 | else |
367 | fLastEvent = TMath::Min(fLastEvent,gime->MaxEvent()); |
eabde521 |
368 | Int_t nEvents = fLastEvent - fFirstEvent + 1; |
88cb7938 |
369 | |
71cee46d |
370 | Int_t ievent ; |
371 | for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) { |
88cb7938 |
372 | gime->Event(ievent,"TR") ; |
373 | if(gime->TrackSegments() && //Skip events, where no track segments made |
374 | gime->TrackSegments()->GetEntriesFast()) { |
7fb46731 |
375 | |
88cb7938 |
376 | MakeRecParticles() ; |
35adb638 |
377 | |
378 | if(fBayesian) |
379 | MakePID() ; |
380 | |
90cceaf6 |
381 | WriteRecParticles(); |
88cb7938 |
382 | if(strstr(option,"deb")) |
383 | PrintRecParticles(option) ; |
384 | //increment the total number of rec particles per run |
385 | fRecParticlesInRun += gime->RecParticles()->GetEntriesFast() ; |
386 | } |
387 | } |
ff417097 |
388 | if(strstr(option,"deb")) |
389 | PrintRecParticles(option); |
88cb7938 |
390 | if(strstr(option,"tim")){ |
391 | gBenchmark->Stop("PHOSPID"); |
351dd634 |
392 | AliInfo(Form("took %f seconds for PID %f seconds per event", |
88cb7938 |
393 | gBenchmark->GetCpuTime("PHOSPID"), |
351dd634 |
394 | gBenchmark->GetCpuTime("PHOSPID")/nEvents)) ; |
88cb7938 |
395 | } |
adcca1e6 |
396 | if(fWrite) |
397 | Unload(); |
88cb7938 |
398 | } |
399 | |
35adb638 |
400 | //________________________________________________________________________ |
17323043 |
401 | Double_t AliPHOSPIDv1::GausF(Double_t x, Double_t y, Double_t * par) |
35adb638 |
402 | { |
c947e71a |
403 | //Given the energy x and the parameter y (tof, shower dispersion or cpv-emc distance), |
404 | //this method returns a density probability of this parameter, given by a gaussian |
405 | //function whose parameters depend with the energy with a function: a/(x*x)+b/x+b |
406 | Double_t cnt = par[1] / (x*x) + par[2] / x + par[0] ; |
35adb638 |
407 | Double_t mean = par[4] / (x*x) + par[5] / x + par[3] ; |
408 | Double_t sigma = par[7] / (x*x) + par[8] / x + par[6] ; |
c947e71a |
409 | |
35adb638 |
410 | // Double_t arg = - (y-mean) * (y-mean) / (2*sigma*sigma) ; |
411 | // return cnt * TMath::Exp(arg) ; |
c947e71a |
412 | if(TMath::Abs(sigma) > 1.e-10){ |
acb5beb7 |
413 | return cnt*TMath::Gaus(y,mean,sigma); |
35adb638 |
414 | } |
415 | else |
416 | return 0.; |
c947e71a |
417 | |
35adb638 |
418 | } |
419 | //________________________________________________________________________ |
17323043 |
420 | Double_t AliPHOSPIDv1::GausPol2(Double_t x, Double_t y, Double_t * par) |
35adb638 |
421 | { |
c947e71a |
422 | //Given the energy x and the parameter y (tof, shower dispersion or cpv-emc distance), |
423 | //this method returns a density probability of this parameter, given by a gaussian |
424 | //function whose parameters depend with the energy like second order polinomial |
425 | |
35adb638 |
426 | Double_t cnt = par[0] + par[1] * x + par[2] * x * x ; |
427 | Double_t mean = par[3] + par[4] * x + par[5] * x * x ; |
428 | Double_t sigma = par[6] + par[7] * x + par[8] * x * x ; |
429 | |
c947e71a |
430 | if(TMath::Abs(sigma) > 1.e-10){ |
acb5beb7 |
431 | return cnt*TMath::Gaus(y,mean,sigma); |
35adb638 |
432 | } |
433 | else |
434 | return 0.; |
c947e71a |
435 | |
436 | |
437 | |
35adb638 |
438 | } |
439 | |
69183710 |
440 | //____________________________________________________________________________ |
e3817e5f |
441 | const TString AliPHOSPIDv1::GetFileNamePrincipal(TString particle) const |
148b2bba |
442 | { |
e3817e5f |
443 | //Get file name that contains the PCA for a particle ("photon or pi0") |
444 | particle.ToLower(); |
445 | TString name; |
351dd634 |
446 | if (particle=="photon") |
447 | name = fFileNamePrincipalPhoton ; |
448 | else if (particle=="pi0" ) |
449 | name = fFileNamePrincipalPi0 ; |
450 | else |
451 | AliError(Form("Wrong particle name: %s (choose from pi0/photon)\n", |
452 | particle.Data())); |
e3817e5f |
453 | return name; |
454 | } |
bc0c084c |
455 | |
e3817e5f |
456 | //____________________________________________________________________________ |
fc7e2f43 |
457 | Float_t AliPHOSPIDv1::GetParameterCalibration(Int_t i) const |
e3817e5f |
458 | { |
459 | // Get the i-th parameter "Calibration" |
460 | Float_t param = 0.; |
351dd634 |
461 | if (i>2 || i<0) { |
462 | AliError(Form("Invalid parameter number: %d",i)); |
463 | } else |
e3817e5f |
464 | param = (*fParameters)(0,i); |
465 | return param; |
466 | } |
bc0c084c |
467 | |
88cb7938 |
468 | //____________________________________________________________________________ |
fc7e2f43 |
469 | Float_t AliPHOSPIDv1::GetCalibratedEnergy(Float_t e) const |
88cb7938 |
470 | { |
471 | // It calibrates Energy depending on the recpoint energy. |
472 | // The energy of the reconstructed cluster is corrected with |
473 | // the formula A + B* E + C* E^2, whose parameters where obtained |
474 | // through the study of the reconstructed energy distribution of |
475 | // monoenergetic photons. |
476 | |
477 | Float_t p[]={0.,0.,0.}; |
478 | for (Int_t i=0; i<3; i++) p[i] = GetParameterCalibration(i); |
479 | Float_t enerec = p[0] + p[1]*e + p[2]*e*e; |
480 | return enerec ; |
481 | |
482 | } |
483 | |
e3817e5f |
484 | //____________________________________________________________________________ |
fc7e2f43 |
485 | Float_t AliPHOSPIDv1::GetParameterCpv2Emc(Int_t i, TString axis) const |
e3817e5f |
486 | { |
487 | // Get the i-th parameter "CPV-EMC distance" for the specified axis |
488 | Float_t param = 0.; |
351dd634 |
489 | if(i>2 || i<0) { |
490 | AliError(Form("Invalid parameter number: %d",i)); |
491 | } else { |
e3817e5f |
492 | axis.ToLower(); |
351dd634 |
493 | if (axis == "x") |
494 | param = (*fParameters)(1,i); |
495 | else if (axis == "z") |
496 | param = (*fParameters)(2,i); |
497 | else { |
498 | AliError(Form("Invalid axis name: %s",axis.Data())); |
499 | } |
e3817e5f |
500 | } |
501 | return param; |
502 | } |
503 | |
504 | //____________________________________________________________________________ |
fc7e2f43 |
505 | Float_t AliPHOSPIDv1::GetCpv2EmcDistanceCut(TString axis, Float_t e) const |
e3817e5f |
506 | { |
88cb7938 |
507 | // Get CpvtoEmcDistance Cut depending on the cluster energy, axis and |
508 | // Purity-Efficiency point |
509 | |
510 | axis.ToLower(); |
511 | Float_t p[]={0.,0.,0.}; |
512 | for (Int_t i=0; i<3; i++) p[i] = GetParameterCpv2Emc(i,axis); |
513 | Float_t sig = p[0] + TMath::Exp(p[1] - p[2]*e); |
514 | return sig; |
e3817e5f |
515 | } |
516 | |
88cb7938 |
517 | //____________________________________________________________________________ |
fc7e2f43 |
518 | Float_t AliPHOSPIDv1::GetEllipseParameter(TString particle, TString param, Float_t e) const |
88cb7938 |
519 | { |
520 | // Calculates the parameter param of the ellipse |
e3817e5f |
521 | |
522 | particle.ToLower(); |
523 | param. ToLower(); |
88cb7938 |
524 | Float_t p[4]={0.,0.,0.,0.}; |
525 | Float_t value = 0.0; |
526 | for (Int_t i=0; i<4; i++) p[i] = GetParameterToCalculateEllipse(particle,param,i); |
527 | if (particle == "photon") { |
528 | if (param.Contains("a")) e = TMath::Min((Double_t)e,70.); |
529 | else if (param.Contains("b")) e = TMath::Min((Double_t)e,70.); |
530 | else if (param.Contains("x0")) e = TMath::Max((Double_t)e,1.1); |
531 | } |
e3817e5f |
532 | |
443caba9 |
533 | if (particle == "photon") |
534 | value = p[0]/TMath::Sqrt(e) + p[1]*e + p[2]*e*e + p[3]; |
535 | else if (particle == "pi0") |
536 | value = p[0] + p[1]*e + p[2]*e*e; |
537 | |
88cb7938 |
538 | return value; |
e3817e5f |
539 | } |
540 | |
541 | //_____________________________________________________________________________ |
fc7e2f43 |
542 | Float_t AliPHOSPIDv1::GetParameterPhotonBoundary (Int_t i) const |
e3817e5f |
543 | { |
544 | // Get the parameter "i" to calculate the boundary on the moment M2x |
545 | // for photons at high p_T |
546 | Float_t param = 0; |
351dd634 |
547 | if (i>3 || i<0) { |
548 | AliError(Form("Wrong parameter number: %d\n",i)); |
549 | } else |
e3817e5f |
550 | param = (*fParameters)(14,i) ; |
551 | return param; |
148b2bba |
552 | } |
e3817e5f |
553 | |
148b2bba |
554 | //____________________________________________________________________________ |
fc7e2f43 |
555 | Float_t AliPHOSPIDv1::GetParameterPi0Boundary (Int_t i) const |
e3817e5f |
556 | { |
557 | // Get the parameter "i" to calculate the boundary on the moment M2x |
558 | // for pi0 at high p_T |
559 | Float_t param = 0; |
351dd634 |
560 | if (i>2 || i<0) { |
561 | AliError(Form("Wrong parameter number: %d\n",i)); |
562 | } else |
e3817e5f |
563 | param = (*fParameters)(15,i) ; |
564 | return param; |
565 | } |
148b2bba |
566 | |
e3817e5f |
567 | //____________________________________________________________________________ |
fc7e2f43 |
568 | Float_t AliPHOSPIDv1::GetParameterTimeGate(Int_t i) const |
e3817e5f |
569 | { |
88cb7938 |
570 | // Get TimeGate parameter depending on Purity-Efficiency i: |
571 | // i=0 - Low purity, i=1 - Medium purity, i=2 - High purity |
572 | Float_t param = 0.; |
351dd634 |
573 | if(i>2 || i<0) { |
574 | AliError(Form("Invalid Efficiency-Purity choice %d",i)); |
575 | } else |
88cb7938 |
576 | param = (*fParameters)(3,i) ; |
577 | return param; |
e3817e5f |
578 | } |
579 | |
148b2bba |
580 | //_____________________________________________________________________________ |
fc7e2f43 |
581 | Float_t AliPHOSPIDv1::GetParameterToCalculateEllipse(TString particle, TString param, Int_t i) const |
88cb7938 |
582 | { |
583 | // Get the parameter "i" that is needed to calculate the ellipse |
584 | // parameter "param" for the particle "particle" ("photon" or "pi0") |
585 | |
e3817e5f |
586 | particle.ToLower(); |
587 | param. ToLower(); |
88cb7938 |
588 | Int_t offset = -1; |
351dd634 |
589 | if (particle == "photon") |
590 | offset=0; |
591 | else if (particle == "pi0") |
592 | offset=5; |
e3817e5f |
593 | else |
351dd634 |
594 | AliError(Form("Wrong particle name: %s (choose from pi0/photon)\n", |
595 | particle.Data())); |
88cb7938 |
596 | |
597 | Int_t p= -1; |
598 | Float_t par = 0; |
e3817e5f |
599 | |
600 | if (param.Contains("a")) p=4+offset; |
601 | else if(param.Contains("b")) p=5+offset; |
602 | else if(param.Contains("c")) p=6+offset; |
603 | else if(param.Contains("x0"))p=7+offset; |
604 | else if(param.Contains("y0"))p=8+offset; |
12022e83 |
605 | |
351dd634 |
606 | if (i>4 || i<0) { |
607 | AliError(Form("No parameter with index %d", i)) ; |
608 | } else if (p==-1) { |
609 | AliError(Form("No parameter with name %s", param.Data() )) ; |
610 | } else |
88cb7938 |
611 | par = (*fParameters)(p,i) ; |
612 | |
613 | return par; |
12022e83 |
614 | } |
615 | |
12022e83 |
616 | |
617 | //____________________________________________________________________________ |
8d4608b5 |
618 | Float_t AliPHOSPIDv1::GetDistance(AliPHOSEmcRecPoint * emc,AliPHOSCpvRecPoint * cpv, Option_t * axis)const |
69183710 |
619 | { |
620 | // Calculates the distance between the EMC RecPoint and the PPSD RecPoint |
148b2bba |
621 | |
88cb7938 |
622 | const AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ; |
69183710 |
623 | TVector3 vecEmc ; |
7acf6008 |
624 | TVector3 vecCpv ; |
148b2bba |
625 | if(cpv){ |
626 | emc->GetLocalPosition(vecEmc) ; |
627 | cpv->GetLocalPosition(vecCpv) ; |
2bb500e5 |
628 | |
148b2bba |
629 | if(emc->GetPHOSMod() == cpv->GetPHOSMod()){ |
630 | // Correct to difference in CPV and EMC position due to different distance to center. |
631 | // we assume, that particle moves from center |
632 | Float_t dCPV = geom->GetIPtoOuterCoverDistance(); |
633 | Float_t dEMC = geom->GetIPtoCrystalSurface() ; |
634 | dEMC = dEMC / dCPV ; |
635 | vecCpv = dEMC * vecCpv - vecEmc ; |
e3817e5f |
636 | if (axis == "X") return vecCpv.X(); |
637 | if (axis == "Y") return vecCpv.Y(); |
638 | if (axis == "Z") return vecCpv.Z(); |
639 | if (axis == "R") return vecCpv.Mag(); |
640 | } |
148b2bba |
641 | return 100000000 ; |
642 | } |
7acf6008 |
643 | return 100000000 ; |
69183710 |
644 | } |
bc0c084c |
645 | //____________________________________________________________________________ |
8d4608b5 |
646 | Int_t AliPHOSPIDv1::GetCPVBit(AliPHOSEmcRecPoint * emc,AliPHOSCpvRecPoint * cpv, Int_t effPur, Float_t e) const |
bc0c084c |
647 | { |
c947e71a |
648 | //Calculates the pid bit for the CPV selection per each purity. |
e3817e5f |
649 | if(effPur>2 || effPur<0) |
351dd634 |
650 | AliError(Form("Invalid Efficiency-Purity choice %d",effPur)); |
bc0c084c |
651 | |
e3817e5f |
652 | Float_t sigX = GetCpv2EmcDistanceCut("X",e); |
653 | Float_t sigZ = GetCpv2EmcDistanceCut("Z",e); |
bc0c084c |
654 | |
655 | Float_t deltaX = TMath::Abs(GetDistance(emc, cpv, "X")); |
656 | Float_t deltaZ = TMath::Abs(GetDistance(emc, cpv, "Z")); |
7fb46731 |
657 | //Info("GetCPVBit"," xdist %f, sigx %f, zdist %f, sigz %f",deltaX, sigX, deltaZ,sigZ) ; |
658 | |
659 | //if(deltaX>sigX*(effPur+1)) |
660 | //if((deltaX>sigX*(effPur+1)) || (deltaZ>sigZ*(effPur+1))) |
661 | if((deltaX>sigX*(effPur+1)) && (deltaZ>sigZ*(effPur+1))) |
bc0c084c |
662 | return 1;//Neutral |
663 | else |
664 | return 0;//Charged |
bc0c084c |
665 | } |
69183710 |
666 | |
6ad0bfa0 |
667 | //____________________________________________________________________________ |
fc7e2f43 |
668 | Int_t AliPHOSPIDv1::GetPrincipalBit(TString particle, const Double_t* p, Int_t effPur, Float_t e)const |
148b2bba |
669 | { |
50739f15 |
670 | //Is the particle inside de PCA ellipse? |
581354c5 |
671 | |
e3817e5f |
672 | particle.ToLower(); |
673 | Int_t prinbit = 0 ; |
7fb46731 |
674 | Float_t a = GetEllipseParameter(particle,"a" , e); |
675 | Float_t b = GetEllipseParameter(particle,"b" , e); |
676 | Float_t c = GetEllipseParameter(particle,"c" , e); |
e3817e5f |
677 | Float_t x0 = GetEllipseParameter(particle,"x0", e); |
678 | Float_t y0 = GetEllipseParameter(particle,"y0", e); |
679 | |
680 | Float_t r = TMath::Power((p[0] - x0)/a,2) + |
681 | TMath::Power((p[1] - y0)/b,2) + |
682 | c*(p[0] - x0)*(p[1] - y0)/(a*b) ; |
50739f15 |
683 | //3 different ellipses defined |
e3817e5f |
684 | if((effPur==2) && (r<1./2.)) prinbit= 1; |
685 | if((effPur==1) && (r<2. )) prinbit= 1; |
686 | if((effPur==0) && (r<9./2.)) prinbit= 1; |
50739f15 |
687 | |
581354c5 |
688 | if(r<0) |
351dd634 |
689 | AliError("Negative square?") ; |
1f0e7ccd |
690 | |
691 | return prinbit; |
148b2bba |
692 | |
148b2bba |
693 | } |
1f0e7ccd |
694 | //____________________________________________________________________________ |
fc7e2f43 |
695 | Int_t AliPHOSPIDv1::GetHardPhotonBit(AliPHOSEmcRecPoint * emc) const |
1f0e7ccd |
696 | { |
e3817e5f |
697 | // Set bit for identified hard photons (E > 30 GeV) |
698 | // if the second moment M2x is below the boundary |
699 | |
700 | Float_t e = emc->GetEnergy(); |
701 | if (e < 30.0) return 0; |
702 | Float_t m2x = emc->GetM2x(); |
703 | Float_t m2xBoundary = GetParameterPhotonBoundary(0) * |
704 | TMath::Exp(-TMath::Power(e-GetParameterPhotonBoundary(1),2)/2.0/ |
705 | TMath::Power(GetParameterPhotonBoundary(2),2)) + |
706 | GetParameterPhotonBoundary(3); |
351dd634 |
707 | AliDebug(1, Form("GetHardPhotonBit","E=%f, m2x=%f, boundary=%f", |
708 | e,m2x,m2xBoundary)); |
e3817e5f |
709 | if (m2x < m2xBoundary) |
710 | return 1;// A hard photon |
711 | else |
712 | return 0;// Not a hard photon |
1f0e7ccd |
713 | } |
92f521a9 |
714 | |
e3817e5f |
715 | //____________________________________________________________________________ |
fc7e2f43 |
716 | Int_t AliPHOSPIDv1::GetHardPi0Bit(AliPHOSEmcRecPoint * emc) const |
e3817e5f |
717 | { |
718 | // Set bit for identified hard pi0 (E > 30 GeV) |
719 | // if the second moment M2x is above the boundary |
720 | |
721 | Float_t e = emc->GetEnergy(); |
722 | if (e < 30.0) return 0; |
723 | Float_t m2x = emc->GetM2x(); |
724 | Float_t m2xBoundary = GetParameterPi0Boundary(0) + |
725 | e * GetParameterPi0Boundary(1); |
351dd634 |
726 | AliDebug(1,Form("E=%f, m2x=%f, boundary=%f",e,m2x,m2xBoundary)); |
e3817e5f |
727 | if (m2x > m2xBoundary) |
728 | return 1;// A hard pi0 |
bc0c084c |
729 | else |
e3817e5f |
730 | return 0;// Not a hard pi0 |
f0a4c9e9 |
731 | } |
e3817e5f |
732 | |
733 | //____________________________________________________________________________ |
8d4608b5 |
734 | TVector3 AliPHOSPIDv1::GetMomentumDirection(AliPHOSEmcRecPoint * emc, AliPHOSCpvRecPoint * )const |
88cb7938 |
735 | { |
736 | // Calculates the momentum direction: |
737 | // 1. if only a EMC RecPoint, direction is given by IP and this RecPoint |
738 | // 2. if a EMC RecPoint and CPV RecPoint, direction is given by the line through the 2 recpoints |
739 | // However because of the poor position resolution of PPSD the direction is always taken as if we were |
740 | // in case 1. |
f0a4c9e9 |
741 | |
88cb7938 |
742 | TVector3 dir(0,0,0) ; |
88cb7938 |
743 | TMatrix dummy ; |
bf8f1fbd |
744 | |
fb7b51ad |
745 | emc->GetGlobalPosition(dir, dummy) ; |
e3817e5f |
746 | |
88cb7938 |
747 | //account correction to the position of IP |
748 | Float_t xo,yo,zo ; //Coordinates of the origin |
7fb46731 |
749 | if(gAlice && gAlice->GetMCApp() && gAlice->Generator()){ |
adcca1e6 |
750 | gAlice->Generator()->GetOrigin(xo,yo,zo) ; |
7fb46731 |
751 | } |
adcca1e6 |
752 | else{ |
753 | xo=yo=zo=0.; |
754 | } |
88cb7938 |
755 | TVector3 origin(xo,yo,zo); |
756 | dir = dir - origin ; |
7fb46731 |
757 | dir.SetMag(1.) ; |
e3817e5f |
758 | |
88cb7938 |
759 | return dir ; |
7b7c1533 |
760 | } |
761 | |
35adb638 |
762 | //________________________________________________________________________ |
17323043 |
763 | Double_t AliPHOSPIDv1::LandauF(Double_t x, Double_t y, Double_t * par) |
35adb638 |
764 | { |
c947e71a |
765 | //Given the energy x and the parameter y (tof, shower dispersion or cpv-emc distance), |
766 | //this method returns a density probability of this parameter, given by a landau |
767 | //function whose parameters depend with the energy with a function: a/(x*x)+b/x+b |
768 | |
769 | Double_t cnt = par[1] / (x*x) + par[2] / x + par[0] ; |
35adb638 |
770 | Double_t mean = par[4] / (x*x) + par[5] / x + par[3] ; |
771 | Double_t sigma = par[7] / (x*x) + par[8] / x + par[6] ; |
c947e71a |
772 | |
773 | if(TMath::Abs(sigma) > 1.e-10){ |
acb5beb7 |
774 | return cnt*TMath::Landau(y,mean,sigma); |
35adb638 |
775 | } |
776 | else |
777 | return 0.; |
778 | |
779 | } |
780 | //________________________________________________________________________ |
17323043 |
781 | Double_t AliPHOSPIDv1::LandauPol2(Double_t x, Double_t y, Double_t * par) |
35adb638 |
782 | { |
c947e71a |
783 | |
784 | //Given the energy x and the parameter y (tof, shower dispersion or cpv-emc distance), |
785 | //this method returns a density probability of this parameter, given by a landau |
786 | //function whose parameters depend with the energy like second order polinomial |
787 | |
35adb638 |
788 | Double_t cnt = par[2] * (x*x) + par[1] * x + par[0] ; |
c947e71a |
789 | Double_t mean = par[5] * (x*x) + par[4] * x + par[3] ; |
790 | Double_t sigma = par[8] * (x*x) + par[7] * x + par[6] ; |
35adb638 |
791 | |
c947e71a |
792 | if(TMath::Abs(sigma) > 1.e-10){ |
acb5beb7 |
793 | return cnt*TMath::Landau(y,mean,sigma); |
35adb638 |
794 | } |
795 | else |
796 | return 0.; |
c947e71a |
797 | |
798 | |
35adb638 |
799 | } |
800 | // //________________________________________________________________________ |
801 | // Double_t AliPHOSPIDv1::ChargedHadronDistProb(Double_t x, Double_t y, Double_t * parg, Double_t * parl) |
802 | // { |
803 | // Double_t cnt = 0.0 ; |
804 | // Double_t mean = 0.0 ; |
805 | // Double_t sigma = 0.0 ; |
806 | // Double_t arg = 0.0 ; |
807 | // if (y < parl[4] / (x*x) + parl[5] / x + parl[3]){ |
808 | // cnt = parg[1] / (x*x) + parg[2] / x + parg[0] ; |
809 | // mean = parg[4] / (x*x) + parg[5] / x + parg[3] ; |
810 | // sigma = parg[7] / (x*x) + parg[8] / x + parg[6] ; |
811 | // TF1 * f = new TF1("gaus","gaus",0.,100.); |
812 | // f->SetParameters(cnt,mean,sigma); |
813 | // arg = f->Eval(y) ; |
814 | // } |
815 | // else{ |
816 | // cnt = parl[1] / (x*x) + parl[2] / x + parl[0] ; |
817 | // mean = parl[4] / (x*x) + parl[5] / x + parl[3] ; |
818 | // sigma = parl[7] / (x*x) + parl[8] / x + parl[6] ; |
819 | // TF1 * f = new TF1("landau","landau",0.,100.); |
820 | // f->SetParameters(cnt,mean,sigma); |
821 | // arg = f->Eval(y) ; |
822 | // } |
823 | // // Double_t mean = par[3] + par[4] * x + par[5] * x * x ; |
824 | // // Double_t sigma = par[6] + par[7] * x + par[8] * x * x ; |
825 | |
826 | // //Double_t arg = -(y-mean)*(y-mean)/(2*sigma*sigma) ; |
827 | // //return cnt * TMath::Exp(arg) ; |
828 | |
829 | // return arg; |
830 | |
831 | // } |
2cc71c1e |
832 | //____________________________________________________________________________ |
833 | void AliPHOSPIDv1::MakePID() |
834 | { |
835 | // construct the PID weight from a Bayesian Method |
c947e71a |
836 | |
304864ab |
837 | const Int_t kSPECIES = AliPID::kSPECIESN ; |
7fb46731 |
838 | |
839 | AliPHOSGetter * gime = AliPHOSGetter::Instance() ; |
840 | |
841 | Int_t nparticles = gime->RecParticles()->GetEntriesFast() ; |
842 | |
c947e71a |
843 | // const Int_t kMAXPARTICLES = 2000 ; |
844 | // if (nparticles >= kMAXPARTICLES) |
845 | // Error("MakePID", "Change size of MAXPARTICLES") ; |
846 | // Double_t stof[kSPECIES][kMAXPARTICLES] ; |
847 | |
53ab54c8 |
848 | // const Int_t kMAXPARTICLES = 2000 ; |
849 | // if (nparticles >= kMAXPARTICLES) |
351dd634 |
850 | // AliError("Change size of MAXPARTICLES") ; |
53ab54c8 |
851 | // Double_t stof[kSPECIES][kMAXPARTICLES] ; |
53ab54c8 |
852 | |
7fb46731 |
853 | |
854 | // make the normalized distribution of pid for this event |
855 | // w(pid) in the Bayesian formulation |
856 | // for(index = 0 ; index < nparticles ; index ++) { |
857 | |
858 | // cout<<">>>>>>>>>>>>>>>Bayes Index "<<index<<endl; |
859 | |
860 | |
861 | // AliPHOSEmcRecPoint * emc = AliPHOSGetter::Instance()->EmcRecPoint(index) ; |
862 | // AliPHOSCpvRecPoint * cpv = AliPHOSGetter::Instance()->CpvRecPoint(index) ; |
863 | |
864 | TObjArray * emcRecPoints = gime->EmcRecPoints() ; |
865 | TObjArray * cpvRecPoints = gime->CpvRecPoints() ; |
866 | TClonesArray * trackSegments = gime->TrackSegments() ; |
867 | if ( !emcRecPoints || !cpvRecPoints || !trackSegments ) { |
868 | AliFatal("RecPoints or TrackSegments not found !") ; |
869 | } |
870 | TIter next(trackSegments) ; |
871 | AliPHOSTrackSegment * ts ; |
872 | Int_t index = 0 ; |
873 | |
35adb638 |
874 | Double_t * stof[kSPECIES] ; |
875 | Double_t * sdp [kSPECIES] ; |
876 | Double_t * scpv[kSPECIES] ; |
fb7b51ad |
877 | Double_t * sw [kSPECIES] ; |
35adb638 |
878 | //Info("MakePID","Begin MakePID"); |
879 | |
880 | for (Int_t i =0; i< kSPECIES; i++){ |
881 | stof[i] = new Double_t[nparticles] ; |
882 | sdp [i] = new Double_t[nparticles] ; |
883 | scpv[i] = new Double_t[nparticles] ; |
fb7b51ad |
884 | sw [i] = new Double_t[nparticles] ; |
35adb638 |
885 | } |
886 | |
7fb46731 |
887 | |
888 | while ( (ts = (AliPHOSTrackSegment *)next()) ) { |
889 | |
890 | //cout<<">>>>>> Bayesian Index "<<index<<endl; |
891 | |
892 | AliPHOSEmcRecPoint * emc = 0 ; |
893 | if(ts->GetEmcIndex()>=0) |
894 | emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ; |
895 | |
896 | AliPHOSCpvRecPoint * cpv = 0 ; |
897 | if(ts->GetCpvIndex()>=0) |
898 | cpv = (AliPHOSCpvRecPoint *) cpvRecPoints->At(ts->GetCpvIndex()) ; |
cc1fe362 |
899 | |
7fb46731 |
900 | // Int_t track = 0 ; |
901 | // track = ts->GetTrackIndex() ; //TPC tracks ? |
cc1fe362 |
902 | |
7fb46731 |
903 | if (!emc) { |
904 | AliFatal(Form("-> emc(%d) = %d", ts->GetEmcIndex(), emc )) ; |
905 | } |
c947e71a |
906 | |
7fb46731 |
907 | // ############Tof############################# |
c947e71a |
908 | |
7fb46731 |
909 | // Info("MakePID", "TOF"); |
fb7b51ad |
910 | Float_t en = emc->GetEnergy(); |
7fb46731 |
911 | Double_t time = emc->GetTime() ; |
912 | // cout<<">>>>>>>Energy "<<en<<"Time "<<time<<endl; |
fb7b51ad |
913 | |
35adb638 |
914 | // now get the signals probability |
915 | // s(pid) in the Bayesian formulation |
cc1fe362 |
916 | |
304864ab |
917 | stof[AliPID::kPhoton][index] = 1.; |
918 | stof[AliPID::kElectron][index] = 1.; |
304864ab |
919 | stof[AliPID::kEleCon][index] = 1.; |
7fb46731 |
920 | //We assing the same prob to charged hadrons, sum is 1 |
921 | stof[AliPID::kPion][index] = 1./3.; |
922 | stof[AliPID::kKaon][index] = 1./3.; |
923 | stof[AliPID::kProton][index] = 1./3.; |
924 | //We assing the same prob to neutral hadrons, sum is 1 |
925 | stof[AliPID::kNeutron][index] = 1./2.; |
926 | stof[AliPID::kKaon0][index] = 1./2.; |
927 | |
304864ab |
928 | stof[AliPID::kMuon][index] = 1.; |
7fb46731 |
929 | |
c947e71a |
930 | |
35adb638 |
931 | if(en < 2.) { |
7fb46731 |
932 | |
933 | Double_t pTofPion = fTFpiong ->Eval(time) ; //gaus distribution |
934 | Double_t pTofKaon = 0; |
935 | |
35adb638 |
936 | if(time < fTkaonl[1]) |
fb7b51ad |
937 | pTofKaon = fTFkaong ->Eval(time) ; //gaus distribution |
35adb638 |
938 | else |
fb7b51ad |
939 | pTofKaon = fTFkaonl ->Eval(time) ; //landau distribution |
7fb46731 |
940 | |
941 | Double_t pTofNucleon = 0; |
942 | |
35adb638 |
943 | if(time < fThhadronl[1]) |
7fb46731 |
944 | pTofNucleon = fTFhhadrong ->Eval(time) ; //gaus distribution |
35adb638 |
945 | else |
7fb46731 |
946 | pTofNucleon = fTFhhadronl ->Eval(time) ; //landau distribution |
7fb46731 |
947 | //We assing the same prob to neutral hadrons, sum is the average prob |
fb7b51ad |
948 | Double_t pTofNeHadron = (pTofKaon + pTofNucleon)/2. ; |
949 | //We assing the same prob to charged hadrons, sum is the average prob |
950 | Double_t pTofChHadron = (pTofPion + pTofKaon + pTofNucleon)/3. ; |
7fb46731 |
951 | |
fb7b51ad |
952 | stof[AliPID::kPhoton][index] = fTFphoton ->Eval(time) ; |
953 | //gaus distribution |
954 | stof[AliPID::kEleCon][index] = stof[AliPID::kPhoton][index] ; |
955 | //a conversion electron has the photon ToF |
304864ab |
956 | stof[AliPID::kMuon][index] = stof[AliPID::kPhoton][index] ; |
7fb46731 |
957 | |
fb7b51ad |
958 | stof[AliPID::kElectron][index] = pTofPion ; |
7fb46731 |
959 | |
960 | stof[AliPID::kPion][index] = pTofChHadron ; |
961 | stof[AliPID::kKaon][index] = pTofChHadron ; |
962 | stof[AliPID::kProton][index] = pTofChHadron ; |
963 | |
964 | stof[AliPID::kKaon0][index] = pTofNeHadron ; |
965 | stof[AliPID::kNeutron][index] = pTofNeHadron ; |
cc1fe362 |
966 | } |
c947e71a |
967 | |
968 | // Info("MakePID", "Dispersion"); |
cc1fe362 |
969 | |
7fb46731 |
970 | // ###########Shower shape: Dispersion#################### |
35adb638 |
971 | Float_t dispersion = emc->GetDispersion(); |
972 | //dispersion is not well defined if the cluster is only in few crystals |
cc1fe362 |
973 | |
304864ab |
974 | sdp[AliPID::kPhoton][index] = 1. ; |
975 | sdp[AliPID::kElectron][index] = 1. ; |
976 | sdp[AliPID::kPion][index] = 1. ; |
977 | sdp[AliPID::kKaon][index] = 1. ; |
978 | sdp[AliPID::kProton][index] = 1. ; |
979 | sdp[AliPID::kNeutron][index] = 1. ; |
980 | sdp[AliPID::kEleCon][index] = 1. ; |
981 | sdp[AliPID::kKaon0][index] = 1. ; |
982 | sdp[AliPID::kMuon][index] = 1. ; |
c947e71a |
983 | |
984 | if(en > 0.5 && emc->GetMultiplicity() > 3){ |
7fb46731 |
985 | sdp[AliPID::kPhoton][index] = GausF(en , dispersion, fDphoton) ; |
304864ab |
986 | sdp[AliPID::kElectron][index] = sdp[AliPID::kPhoton][index] ; |
987 | sdp[AliPID::kPion][index] = LandauF(en , dispersion, fDhadron ) ; |
988 | sdp[AliPID::kKaon][index] = sdp[AliPID::kPion][index] ; |
989 | sdp[AliPID::kProton][index] = sdp[AliPID::kPion][index] ; |
990 | sdp[AliPID::kNeutron][index] = sdp[AliPID::kPion][index] ; |
991 | sdp[AliPID::kEleCon][index] = sdp[AliPID::kPhoton][index]; |
992 | sdp[AliPID::kKaon0][index] = sdp[AliPID::kPion][index] ; |
fb7b51ad |
993 | sdp[AliPID::kMuon][index] = fDFmuon ->Eval(dispersion) ; |
994 | //landau distribution |
35adb638 |
995 | } |
cc1fe362 |
996 | |
7fb46731 |
997 | // Info("MakePID","multiplicity %d, dispersion %f", emc->GetMultiplicity(), dispersion); |
998 | // Info("MakePID","ss: photon %f, hadron %f ", sdp[AliPID::kPhoton][index], sdp[AliPID::kPion][index]); |
c947e71a |
999 | // cout<<">>>>>multiplicity "<<emc->GetMultiplicity()<<", dispersion "<< dispersion<<endl ; |
304864ab |
1000 | // cout<<"<<<<<ss: photon "<<sdp[AliPID::kPhoton][index]<<", hadron "<<sdp[AliPID::kPion][index]<<endl; |
c947e71a |
1001 | |
7fb46731 |
1002 | //########## CPV-EMC Distance####################### |
1003 | // Info("MakePID", "Distance"); |
fb7b51ad |
1004 | |
7fb46731 |
1005 | Float_t x = TMath::Abs(GetDistance(emc, cpv, "X")) ; |
1006 | Float_t z = GetDistance(emc, cpv, "Z") ; |
fb7b51ad |
1007 | |
7fb46731 |
1008 | Double_t pcpv = 0 ; |
c947e71a |
1009 | Double_t pcpvneutral = 0. ; |
7fb46731 |
1010 | |
1011 | Double_t elprobx = GausF(en , x, fXelectron) ; |
1012 | Double_t elprobz = GausF(en , z, fZelectron) ; |
1013 | Double_t chprobx = GausF(en , x, fXcharged) ; |
1014 | Double_t chprobz = GausF(en , z, fZcharged) ; |
c947e71a |
1015 | Double_t pcpvelectron = elprobx * elprobz; |
1016 | Double_t pcpvcharged = chprobx * chprobz; |
7fb46731 |
1017 | |
1018 | // cout<<">>>>energy "<<en<<endl; |
c947e71a |
1019 | // cout<<">>>>electron : x "<<x<<" xprob "<<elprobx<<" z "<<z<<" zprob "<<elprobz<<endl; |
1020 | // cout<<">>>>hadron : x "<<x<<" xprob "<<chprobx<<" z "<<z<<" zprob "<<chprobz<<endl; |
1021 | // cout<<">>>>electron : px*pz "<<pcpvelectron <<" hadron: px*pz "<<pcpvcharged<<endl; |
1022 | |
7fb46731 |
1023 | // Is neutral or charged |
35adb638 |
1024 | if(pcpvelectron >= pcpvcharged) |
1025 | pcpv = pcpvelectron ; |
1026 | else |
1027 | pcpv = pcpvcharged ; |
1028 | |
fb7b51ad |
1029 | if(pcpv < fChargedNeutralThreshold) |
35adb638 |
1030 | { |
1031 | pcpvneutral = 1. ; |
1032 | pcpvcharged = 0. ; |
1033 | pcpvelectron = 0. ; |
1034 | } |
c947e71a |
1035 | // else |
1036 | // cout<<">>>>>>>>>>>CHARGED>>>>>>>>>>>"<<endl; |
35adb638 |
1037 | |
304864ab |
1038 | scpv[AliPID::kPion][index] = pcpvcharged ; |
1039 | scpv[AliPID::kKaon][index] = pcpvcharged ; |
1040 | scpv[AliPID::kProton][index] = pcpvcharged ; |
7fb46731 |
1041 | |
1042 | scpv[AliPID::kMuon][index] = pcpvelectron ; |
304864ab |
1043 | scpv[AliPID::kElectron][index] = pcpvelectron ; |
304864ab |
1044 | scpv[AliPID::kEleCon][index] = pcpvelectron ; |
7fb46731 |
1045 | |
1046 | scpv[AliPID::kPhoton][index] = pcpvneutral ; |
1047 | scpv[AliPID::kNeutron][index] = pcpvneutral ; |
304864ab |
1048 | scpv[AliPID::kKaon0][index] = pcpvneutral ; |
7fb46731 |
1049 | |
35adb638 |
1050 | |
1051 | // Info("MakePID", "CPV passed"); |
c947e71a |
1052 | |
7fb46731 |
1053 | //############## Pi0 ############################# |
304864ab |
1054 | stof[AliPID::kPi0][index] = 0. ; |
1055 | scpv[AliPID::kPi0][index] = 0. ; |
1056 | sdp [AliPID::kPi0][index] = 0. ; |
c947e71a |
1057 | |
35adb638 |
1058 | if(en > 30.){ |
1059 | // pi0 are detected via decay photon |
304864ab |
1060 | stof[AliPID::kPi0][index] = fTFphoton ->Eval(time) ; |
1061 | scpv[AliPID::kPi0][index] = pcpvneutral ; |
1062 | sdp [AliPID::kPi0][index] = 1. ; |
c947e71a |
1063 | if(emc->GetMultiplicity() > 3) |
304864ab |
1064 | sdp [AliPID::kPi0][index] = GausPol2(en , dispersion, fDpi0) ; |
35adb638 |
1065 | } |
1066 | |
7fb46731 |
1067 | |
1068 | //############## muon ############################# |
1069 | |
35adb638 |
1070 | if(en > 0.5){ |
1071 | //Muons deposit few energy |
304864ab |
1072 | scpv[AliPID::kMuon][index] = 0 ; |
1073 | stof[AliPID::kMuon][index] = 0 ; |
1074 | sdp [AliPID::kMuon][index] = 0 ; |
c947e71a |
1075 | } |
1076 | |
fb7b51ad |
1077 | //Weight to apply to hadrons due to energy reconstruction |
1078 | |
1079 | Float_t weight = fERecWeight ->Eval(en) ; |
1080 | |
1081 | sw[AliPID::kPhoton][index] = 1. ; |
1082 | sw[AliPID::kElectron][index] = 1. ; |
1083 | sw[AliPID::kPion][index] = weight ; |
1084 | sw[AliPID::kKaon][index] = weight ; |
1085 | sw[AliPID::kProton][index] = weight ; |
1086 | sw[AliPID::kNeutron][index] = weight ; |
1087 | sw[AliPID::kEleCon][index] = 1. ; |
1088 | sw[AliPID::kKaon0][index] = weight ; |
1089 | sw[AliPID::kMuon][index] = weight ; |
1090 | sw[AliPID::kPi0][index] = 1. ; |
1091 | |
7fb46731 |
1092 | // if(en > 0.5){ |
1093 | // cout<<"######################################################"<<endl; |
1094 | // //cout<<"MakePID: energy "<<en<<", tof "<<time<<", distance "<<distance<<", dispersion "<<dispersion<<endl ; |
1095 | // cout<<"MakePID: energy "<<en<<", tof "<<time<<", dispersion "<<dispersion<<", x "<<x<<", z "<<z<<endl ; |
1096 | // cout<<">>>>>multiplicity "<<emc->GetMultiplicity()<<endl; |
1097 | // cout<<">>>>electron : xprob "<<elprobx<<" zprob "<<elprobz<<endl; |
1098 | // cout<<">>>>hadron : xprob "<<chprobx<<" zprob "<<chprobz<<endl; |
1099 | // cout<<">>>>electron : px*pz "<<pcpvelectron <<" hadron: px*pz "<<pcpvcharged<<endl; |
c947e71a |
1100 | |
7fb46731 |
1101 | // cout<<"Photon , pid "<< fInitPID[AliPID::kPhoton]<<" tof "<<stof[AliPID::kPhoton][index] |
1102 | // <<", cpv "<<scpv[AliPID::kPhoton][index]<<", ss "<<sdp[AliPID::kPhoton][index]<<endl; |
1103 | // cout<<"EleCon , pid "<< fInitPID[AliPID::kEleCon]<<", tof "<<stof[AliPID::kEleCon][index] |
1104 | // <<", cpv "<<scpv[AliPID::kEleCon][index]<<" ss "<<sdp[AliPID::kEleCon][index]<<endl; |
1105 | // cout<<"Electron , pid "<< fInitPID[AliPID::kElectron]<<", tof "<<stof[AliPID::kElectron][index] |
1106 | // <<", cpv "<<scpv[AliPID::kElectron][index]<<" ss "<<sdp[AliPID::kElectron][index]<<endl; |
1107 | // cout<<"Muon , pid "<< fInitPID[AliPID::kMuon]<<", tof "<<stof[AliPID::kMuon][index] |
1108 | // <<", cpv "<<scpv[AliPID::kMuon][index]<<" ss "<<sdp[AliPID::kMuon][index]<<endl; |
1109 | // cout<<"Pi0 , pid "<< fInitPID[AliPID::kPi0]<<", tof "<<stof[AliPID::kPi0][index] |
1110 | // <<", cpv "<<scpv[AliPID::kPi0][index]<<" ss "<<sdp[AliPID::kPi0][index]<<endl; |
1111 | // cout<<"Pion , pid "<< fInitPID[AliPID::kPion]<<", tof "<<stof[AliPID::kPion][index] |
1112 | // <<", cpv "<<scpv[AliPID::kPion][index]<<" ss "<<sdp[AliPID::kPion][index]<<endl; |
1113 | // cout<<"Kaon0 , pid "<< fInitPID[AliPID::kKaon0]<<", tof "<<stof[AliPID::kKaon0][index] |
1114 | // <<", cpv "<<scpv[AliPID::kKaon0][index]<<" ss "<<sdp[AliPID::kKaon0][index]<<endl; |
1115 | // cout<<"Kaon , pid "<< fInitPID[AliPID::kKaon]<<", tof "<<stof[AliPID::kKaon][index] |
1116 | // <<", cpv "<<scpv[AliPID::kKaon][index]<<" ss "<<sdp[AliPID::kKaon][index]<<endl; |
1117 | // cout<<"Neutron , pid "<< fInitPID[AliPID::kNeutron]<<", tof "<<stof[AliPID::kNeutron][index] |
1118 | // <<", cpv "<<scpv[AliPID::kNeutron][index]<<" ss "<<sdp[AliPID::kNeutron][index]<<endl; |
1119 | // cout<<"Proton , pid "<< fInitPID[AliPID::kProton]<<", tof "<<stof[AliPID::kProton][index] |
1120 | // <<", cpv "<<scpv[AliPID::kProton][index]<<" ss "<<sdp[AliPID::kProton][index]<<endl; |
1121 | // cout<<"######################################################"<<endl; |
1122 | // } |
1123 | index++; |
cc1fe362 |
1124 | } |
35adb638 |
1125 | |
1126 | //for (index = 0 ; index < kSPECIES ; index++) |
1127 | // pid[index] /= nparticles ; |
1128 | |
7fb46731 |
1129 | |
35adb638 |
1130 | // Info("MakePID", "Total Probability calculation"); |
1131 | |
cc1fe362 |
1132 | for(index = 0 ; index < nparticles ; index ++) { |
1133 | // calculates the Bayesian weight |
7fb46731 |
1134 | |
cc1fe362 |
1135 | Int_t jndex ; |
1136 | Double_t wn = 0.0 ; |
1137 | for (jndex = 0 ; jndex < kSPECIES ; jndex++) |
fb7b51ad |
1138 | wn += stof[jndex][index] * sdp[jndex][index] * scpv[jndex][index] * sw[jndex][index] * fInitPID[jndex] ; |
7fb46731 |
1139 | |
1140 | // cout<<"*************wn "<<wn<<endl; |
1141 | AliPHOSRecParticle * recpar = gime->RecParticle(index) ; |
e74ea0e9 |
1142 | if (TMath::Abs(wn)>0) |
1143 | for (jndex = 0 ; jndex < kSPECIES ; jndex++) { |
fb7b51ad |
1144 | |
1145 | // if(recpar->IsEleCon()){ |
1146 | // fInitPID[AliPID::kEleCon] = 1. ; |
1147 | // fInitPID[AliPID::kPhoton] = 0. ; |
1148 | // } |
1149 | // else{ |
1150 | // fInitPID[AliPID::kEleCon] = 0. ; |
1151 | // fInitPID[AliPID::kPhoton] = 1. ; |
1152 | // } |
1153 | fInitPID[AliPID::kEleCon] = 0. ; |
35adb638 |
1154 | //cout<<"jndex "<<jndex<<" wn "<<wn<<" SetPID * wn" |
1155 | //<<stof[jndex][index] * sdp[jndex][index] * pid[jndex] << endl; |
1156 | //cout<<" tof "<<stof[jndex][index] << " disp " <<sdp[jndex][index] << " pid "<< fInitPID[jndex] << endl; |
7fb46731 |
1157 | // cout<<"Particle "<<jndex<<" final prob * wn " |
1158 | // <<stof[jndex][index] * sdp[jndex][index] * scpv[jndex][index] * fInitPID[jndex] <<" wn "<< wn<<endl; |
fb7b51ad |
1159 | recpar->SetPID(jndex, stof[jndex][index] * sdp[jndex][index] * sw[jndex][index] * |
35adb638 |
1160 | scpv[jndex][index] * fInitPID[jndex] / wn) ; |
1161 | // cout<<"final prob "<<stof[jndex][index] * sdp[jndex][index] * scpv[jndex][index] * fInitPID[jndex] / wn<<endl; |
1162 | //recpar->SetPID(jndex, stof[jndex][index] * fInitPID[jndex] / wn) ; |
1163 | //cout<<"After SetPID"<<endl; |
1164 | //recpar->Print(); |
e74ea0e9 |
1165 | } |
2cc71c1e |
1166 | } |
35adb638 |
1167 | // Info("MakePID", "Delete"); |
1168 | |
acb5beb7 |
1169 | for (Int_t i =0; i< kSPECIES; i++){ |
1170 | delete [] stof[i]; |
fb7b51ad |
1171 | delete [] sdp [i]; |
acb5beb7 |
1172 | delete [] scpv[i]; |
fb7b51ad |
1173 | delete [] sw [i]; |
acb5beb7 |
1174 | } |
35adb638 |
1175 | // Info("MakePID","End MakePID"); |
2cc71c1e |
1176 | } |
1177 | |
7acf6008 |
1178 | //____________________________________________________________________________ |
e3817e5f |
1179 | void AliPHOSPIDv1::MakeRecParticles() |
1180 | { |
b2a60966 |
1181 | // Makes a RecParticle out of a TrackSegment |
148b2bba |
1182 | |
88cb7938 |
1183 | AliPHOSGetter * gime = AliPHOSGetter::Instance() ; |
fbf811ec |
1184 | TObjArray * emcRecPoints = gime->EmcRecPoints() ; |
1185 | TObjArray * cpvRecPoints = gime->CpvRecPoints() ; |
1186 | TClonesArray * trackSegments = gime->TrackSegments() ; |
148b2bba |
1187 | if ( !emcRecPoints || !cpvRecPoints || !trackSegments ) { |
351dd634 |
1188 | AliFatal("RecPoints or TrackSegments not found !") ; |
148b2bba |
1189 | } |
fbf811ec |
1190 | TClonesArray * recParticles = gime->RecParticles() ; |
01a599c9 |
1191 | recParticles->Clear(); |
148b2bba |
1192 | |
7b7c1533 |
1193 | TIter next(trackSegments) ; |
7acf6008 |
1194 | AliPHOSTrackSegment * ts ; |
6ad0bfa0 |
1195 | Int_t index = 0 ; |
09fc14a0 |
1196 | AliPHOSRecParticle * rp ; |
7acf6008 |
1197 | while ( (ts = (AliPHOSTrackSegment *)next()) ) { |
7fb46731 |
1198 | // cout<<">>>>>>>>>>>>>>>PCA Index "<<index<<endl; |
7b7c1533 |
1199 | new( (*recParticles)[index] ) AliPHOSRecParticle() ; |
1200 | rp = (AliPHOSRecParticle *)recParticles->At(index) ; |
f0a4c9e9 |
1201 | rp->SetTrackSegment(index) ; |
9688c1dd |
1202 | rp->SetIndexInList(index) ; |
148b2bba |
1203 | |
7acf6008 |
1204 | AliPHOSEmcRecPoint * emc = 0 ; |
1205 | if(ts->GetEmcIndex()>=0) |
7b7c1533 |
1206 | emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ; |
fad3e5b9 |
1207 | |
8d4608b5 |
1208 | AliPHOSCpvRecPoint * cpv = 0 ; |
7acf6008 |
1209 | if(ts->GetCpvIndex()>=0) |
8d4608b5 |
1210 | cpv = (AliPHOSCpvRecPoint *) cpvRecPoints->At(ts->GetCpvIndex()) ; |
fad3e5b9 |
1211 | |
bd76890a |
1212 | Int_t track = 0 ; |
1213 | track = ts->GetTrackIndex() ; |
1214 | |
148b2bba |
1215 | // Now set type (reconstructed) of the particle |
1216 | |
1217 | // Choose the cluster energy range |
9fa5f1d0 |
1218 | |
fbf811ec |
1219 | if (!emc) { |
351dd634 |
1220 | AliFatal(Form("-> emc(%d) = %d", ts->GetEmcIndex(), emc )) ; |
fbf811ec |
1221 | } |
50739f15 |
1222 | |
e3817e5f |
1223 | Float_t e = emc->GetEnergy() ; |
bc0c084c |
1224 | |
6f969528 |
1225 | Float_t lambda[2] ; |
1226 | emc->GetElipsAxis(lambda) ; |
50739f15 |
1227 | |
1228 | if((lambda[0]>0.01) && (lambda[1]>0.01)){ |
1229 | // Looking PCA. Define and calculate the data (X), |
bc0c084c |
1230 | // introduce in the function X2P that gives the components (P). |
1231 | |
c947e71a |
1232 | Float_t spher = 0. ; |
1233 | Float_t emaxdtotal = 0. ; |
50739f15 |
1234 | |
bc0c084c |
1235 | if((lambda[0]+lambda[1])!=0) |
c947e71a |
1236 | spher=fabs(lambda[0]-lambda[1])/(lambda[0]+lambda[1]); |
50739f15 |
1237 | |
c947e71a |
1238 | emaxdtotal=emc->GetMaximalEnergy()/emc->GetEnergy(); |
50739f15 |
1239 | |
1240 | fX[0] = lambda[0] ; |
1241 | fX[1] = lambda[1] ; |
1242 | fX[2] = emc->GetDispersion() ; |
c947e71a |
1243 | fX[3] = spher ; |
50739f15 |
1244 | fX[4] = emc->GetMultiplicity() ; |
c947e71a |
1245 | fX[5] = emaxdtotal ; |
50739f15 |
1246 | fX[6] = emc->GetCoreEnergy() ; |
1247 | |
e3817e5f |
1248 | fPrincipalPhoton->X2P(fX,fPPhoton); |
1249 | fPrincipalPi0 ->X2P(fX,fPPi0); |
1f0e7ccd |
1250 | |
50739f15 |
1251 | } |
1252 | else{ |
e3817e5f |
1253 | fPPhoton[0]=-100.0; //We do not accept clusters with |
1254 | fPPhoton[1]=-100.0; //one cell as a photon-like |
1255 | fPPi0[0] =-100.0; |
1256 | fPPi0[1] =-100.0; |
50739f15 |
1257 | } |
1258 | |
2cc71c1e |
1259 | Float_t time = emc->GetTime() ; |
1260 | rp->SetTof(time) ; |
9fa5f1d0 |
1261 | |
bc0c084c |
1262 | // Loop of Efficiency-Purity (the 3 points of purity or efficiency |
1263 | // are taken into account to set the particle identification) |
e3817e5f |
1264 | for(Int_t effPur = 0; effPur < 3 ; effPur++){ |
50739f15 |
1265 | |
bc0c084c |
1266 | // Looking at the CPV detector. If RCPV greater than CpvEmcDistance, |
1267 | // 1st,2nd or 3rd bit (depending on the efficiency-purity point ) |
1268 | // is set to 1 |
35adb638 |
1269 | if(GetCPVBit(emc, cpv, effPur,e) == 1 ){ |
e3817e5f |
1270 | rp->SetPIDBit(effPur) ; |
35adb638 |
1271 | //cout<<"CPV bit "<<effPur<<endl; |
1272 | } |
50739f15 |
1273 | // Looking the TOF. If TOF smaller than gate, 4th, 5th or 6th |
1274 | // bit (depending on the efficiency-purity point )is set to 1 |
2cc71c1e |
1275 | if(time< (*fParameters)(3,effPur)) |
e3817e5f |
1276 | rp->SetPIDBit(effPur+3) ; |
2cc71c1e |
1277 | |
e3817e5f |
1278 | //Photon PCA |
50739f15 |
1279 | //If we are inside the ellipse, 7th, 8th or 9th |
1280 | // bit (depending on the efficiency-purity point )is set to 1 |
e3817e5f |
1281 | if(GetPrincipalBit("photon",fPPhoton,effPur,e) == 1) |
1282 | rp->SetPIDBit(effPur+6) ; |
1f0e7ccd |
1283 | |
e3817e5f |
1284 | //Pi0 PCA |
1f0e7ccd |
1285 | //If we are inside the ellipse, 10th, 11th or 12th |
1286 | // bit (depending on the efficiency-purity point )is set to 1 |
e3817e5f |
1287 | if(GetPrincipalBit("pi0" ,fPPi0 ,effPur,e) == 1) |
1288 | rp->SetPIDBit(effPur+9) ; |
f0a4c9e9 |
1289 | } |
e3817e5f |
1290 | if(GetHardPhotonBit(emc)) |
1291 | rp->SetPIDBit(12) ; |
1292 | if(GetHardPi0Bit (emc)) |
1293 | rp->SetPIDBit(13) ; |
1f0e7ccd |
1294 | |
bd76890a |
1295 | if(track >= 0) |
1296 | rp->SetPIDBit(14) ; |
1297 | |
9fa5f1d0 |
1298 | //Set momentum, energy and other parameters |
50739f15 |
1299 | Float_t encal = GetCalibratedEnergy(e); |
9fa5f1d0 |
1300 | TVector3 dir = GetMomentumDirection(emc,cpv) ; |
1301 | dir.SetMag(encal) ; |
1302 | rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),encal) ; |
1303 | rp->SetCalcMass(0); |
e0ed2e49 |
1304 | rp->Name(); //If photon sets the particle pdg name to gamma |
e747b8da |
1305 | rp->SetProductionVertex(0,0,0,0); |
1306 | rp->SetFirstMother(-1); |
1307 | rp->SetLastMother(-1); |
1308 | rp->SetFirstDaughter(-1); |
1309 | rp->SetLastDaughter(-1); |
1310 | rp->SetPolarisation(0,0,0); |
d956e9b7 |
1311 | //Set the position in global coordinate system from the RecPoint |
1312 | AliPHOSGeometry * geom = gime->PHOSGeometry() ; |
1313 | AliPHOSTrackSegment * ts = gime->TrackSegment(rp->GetPHOSTSIndex()) ; |
1314 | AliPHOSEmcRecPoint * erp = gime->EmcRecPoint(ts->GetEmcIndex()) ; |
1315 | TVector3 pos ; |
1316 | geom->GetGlobal(erp, pos) ; |
1317 | rp->SetPos(pos); |
6ad0bfa0 |
1318 | index++ ; |
1319 | } |
6ad0bfa0 |
1320 | } |
e3817e5f |
1321 | |
09fc14a0 |
1322 | //____________________________________________________________________________ |
88cb7938 |
1323 | void AliPHOSPIDv1::Print() const |
09fc14a0 |
1324 | { |
b2a60966 |
1325 | // Print the parameters used for the particle type identification |
bc0c084c |
1326 | |
351dd634 |
1327 | AliInfo("=============== AliPHOSPIDv1 ================") ; |
88cb7938 |
1328 | printf("Making PID\n") ; |
1329 | printf(" Pricipal analysis file from 0.5 to 100 %s\n", fFileNamePrincipalPhoton.Data() ) ; |
1330 | printf(" Name of parameters file %s\n", fFileNameParameters.Data() ) ; |
1331 | printf(" Matrix of Parameters: 14x4\n") ; |
1332 | printf(" Energy Calibration 1x3 [3 parametres to calibrate energy: A + B* E + C * E^2]\n") ; |
1333 | printf(" RCPV 2x3 rows x and z, columns function cut parameters\n") ; |
1334 | printf(" TOF 1x3 [High Eff-Low Pur,Medium Eff-Pur, Low Eff-High Pur]\n") ; |
1335 | printf(" PCA 5x4 [5 ellipse parametres and 4 parametres to calculate them: A/Sqrt(E) + B* E + C * E^2 + D]\n") ; |
1336 | Printf(" Pi0 PCA 5x3 [5 ellipse parametres and 3 parametres to calculate them: A + B* E + C * E^2]\n") ; |
50739f15 |
1337 | fParameters->Print() ; |
09fc14a0 |
1338 | } |
1339 | |
a496c46c |
1340 | |
69183710 |
1341 | |
7acf6008 |
1342 | //____________________________________________________________________________ |
a4e98857 |
1343 | void AliPHOSPIDv1::PrintRecParticles(Option_t * option) |
1344 | { |
dd5c4038 |
1345 | // Print table of reconstructed particles |
1346 | |
88cb7938 |
1347 | AliPHOSGetter *gime = AliPHOSGetter::Instance() ; |
bf8f1fbd |
1348 | |
88cb7938 |
1349 | TClonesArray * recParticles = gime->RecParticles() ; |
21cd0c07 |
1350 | |
1351 | TString message ; |
3bf72d32 |
1352 | message = "\nevent " ; |
1353 | message += gAlice->GetEvNumber() ; |
1354 | message += " found " ; |
1355 | message += recParticles->GetEntriesFast(); |
1356 | message += " RecParticles\n" ; |
1357 | |
7acf6008 |
1358 | if(strstr(option,"all")) { // printing found TS |
3bf72d32 |
1359 | message += "\n PARTICLE Index \n" ; |
7acf6008 |
1360 | |
1361 | Int_t index ; |
7b7c1533 |
1362 | for (index = 0 ; index < recParticles->GetEntries() ; index++) { |
21cd0c07 |
1363 | AliPHOSRecParticle * rp = (AliPHOSRecParticle * ) recParticles->At(index) ; |
3bf72d32 |
1364 | message += "\n" ; |
1365 | message += rp->Name().Data() ; |
1366 | message += " " ; |
1367 | message += rp->GetIndexInList() ; |
1368 | message += " " ; |
1369 | message += rp->GetType() ; |
7acf6008 |
1370 | } |
3bf72d32 |
1371 | } |
351dd634 |
1372 | AliInfo(message.Data() ) ; |
69183710 |
1373 | } |
88cb7938 |
1374 | |
1375 | //____________________________________________________________________________ |
1376 | void AliPHOSPIDv1::SetParameters() |
1377 | { |
1378 | // PCA : To do the Principal Components Analysis it is necessary |
1379 | // the Principal file, which is opened here |
1380 | fX = new double[7]; // Data for the PCA |
1381 | fPPhoton = new double[7]; // Eigenvalues of the PCA |
1382 | fPPi0 = new double[7]; // Eigenvalues of the Pi0 PCA |
1383 | |
1384 | // Read photon principals from the photon file |
1385 | |
1386 | fFileNamePrincipalPhoton = "$ALICE_ROOT/PHOS/PCA8pa15_0.5-100.root" ; |
1387 | TFile f( fFileNamePrincipalPhoton.Data(), "read" ) ; |
1388 | fPrincipalPhoton = dynamic_cast<TPrincipal*> (f.Get("principal")) ; |
1389 | f.Close() ; |
1390 | |
1391 | // Read pi0 principals from the pi0 file |
1392 | |
1393 | fFileNamePrincipalPi0 = "$ALICE_ROOT/PHOS/PCA_pi0_40-120.root" ; |
1394 | TFile fPi0( fFileNamePrincipalPi0.Data(), "read" ) ; |
1395 | fPrincipalPi0 = dynamic_cast<TPrincipal*> (fPi0.Get("principal")) ; |
1396 | fPi0.Close() ; |
1397 | |
1398 | // Open parameters file and initialization of the Parameters matrix. |
1399 | // In the File Parameters.dat are all the parameters. These are introduced |
1400 | // in a matrix of 16x4 |
1401 | // |
1402 | // All the parameters defined in this file are, in order of row: |
1403 | // line 0 : calibration |
1404 | // lines 1,2 : CPV rectangular cat for X and Z |
1405 | // line 3 : TOF cut |
1406 | // lines 4-8 : parameters to calculate photon PCA ellipse |
1407 | // lines 9-13: parameters to calculate pi0 PCA ellipse |
1408 | // lines 14-15: parameters to calculate border for high-pt photons and pi0 |
1409 | |
1410 | fFileNameParameters = gSystem->ExpandPathName("$ALICE_ROOT/PHOS/Parameters.dat"); |
1411 | fParameters = new TMatrix(16,4) ; |
c947e71a |
1412 | const Int_t kMaxLeng=255; |
1413 | char string[kMaxLeng]; |
88cb7938 |
1414 | |
1415 | // Open a text file with PID parameters |
1416 | FILE *fd = fopen(fFileNameParameters.Data(),"r"); |
1417 | if (!fd) |
351dd634 |
1418 | AliFatal(Form("File %s with a PID parameters cannot be opened\n", |
1419 | fFileNameParameters.Data())); |
88cb7938 |
1420 | |
1421 | Int_t i=0; |
1422 | // Read parameter file line-by-line and skip empty line and comments |
c947e71a |
1423 | while (fgets(string,kMaxLeng,fd) != NULL) { |
88cb7938 |
1424 | if (string[0] == '\n' ) continue; |
1425 | if (string[0] == '!' ) continue; |
1426 | sscanf(string, "%f %f %f %f", |
1427 | &(*fParameters)(i,0), &(*fParameters)(i,1), |
1428 | &(*fParameters)(i,2), &(*fParameters)(i,3)); |
1429 | i++; |
351dd634 |
1430 | AliDebug(1, Form("SetParameters", "line %d: %s",i,string)); |
88cb7938 |
1431 | } |
1432 | fclose(fd); |
1433 | } |
1434 | |
1435 | //____________________________________________________________________________ |
1436 | void AliPHOSPIDv1::SetParameterCalibration(Int_t i,Float_t param) |
1437 | { |
1438 | // Set parameter "Calibration" i to a value param |
351dd634 |
1439 | if(i>2 || i<0) { |
1440 | AliError(Form("Invalid parameter number: %d",i)); |
1441 | } else |
88cb7938 |
1442 | (*fParameters)(0,i) = param ; |
1443 | } |
1444 | |
1445 | //____________________________________________________________________________ |
1446 | void AliPHOSPIDv1::SetParameterCpv2Emc(Int_t i, TString axis, Float_t cut) |
1447 | { |
1448 | // Set the parameters to calculate Cpv-to-Emc Distance Cut depending on |
1449 | // Purity-Efficiency point i |
1450 | |
351dd634 |
1451 | if(i>2 || i<0) { |
1452 | AliError(Form("Invalid parameter number: %d",i)); |
1453 | } else { |
88cb7938 |
1454 | axis.ToLower(); |
1455 | if (axis == "x") (*fParameters)(1,i) = cut; |
1456 | else if (axis == "z") (*fParameters)(2,i) = cut; |
351dd634 |
1457 | else { |
1458 | AliError(Form("Invalid axis name: %s",axis.Data())); |
1459 | } |
88cb7938 |
1460 | } |
1461 | } |
1462 | |
1463 | //____________________________________________________________________________ |
1464 | void AliPHOSPIDv1::SetParameterPhotonBoundary(Int_t i,Float_t param) |
1465 | { |
1466 | // Set parameter "Hard photon boundary" i to a value param |
351dd634 |
1467 | if(i>4 || i<0) { |
1468 | AliError(Form("Invalid parameter number: %d",i)); |
1469 | } else |
88cb7938 |
1470 | (*fParameters)(14,i) = param ; |
1471 | } |
1472 | |
1473 | //____________________________________________________________________________ |
1474 | void AliPHOSPIDv1::SetParameterPi0Boundary(Int_t i,Float_t param) |
1475 | { |
1476 | // Set parameter "Hard pi0 boundary" i to a value param |
351dd634 |
1477 | if(i>1 || i<0) { |
1478 | AliError(Form("Invalid parameter number: %d",i)); |
1479 | } else |
88cb7938 |
1480 | (*fParameters)(15,i) = param ; |
1481 | } |
1482 | |
1483 | //_____________________________________________________________________________ |
1484 | void AliPHOSPIDv1::SetParameterTimeGate(Int_t i, Float_t gate) |
1485 | { |
1486 | // Set the parameter TimeGate depending on Purity-Efficiency point i |
351dd634 |
1487 | if (i>2 || i<0) { |
1488 | AliError(Form("Invalid Efficiency-Purity choice %d",i)); |
1489 | } else |
88cb7938 |
1490 | (*fParameters)(3,i)= gate ; |
1491 | } |
1492 | |
1493 | //_____________________________________________________________________________ |
1494 | void AliPHOSPIDv1::SetParameterToCalculateEllipse(TString particle, TString param, Int_t i, Float_t par) |
1495 | { |
1496 | // Set the parameter "i" that is needed to calculate the ellipse |
1497 | // parameter "param" for a particle "particle" |
1498 | |
1499 | particle.ToLower(); |
1500 | param. ToLower(); |
1501 | Int_t p= -1; |
1502 | Int_t offset=0; |
1503 | |
1504 | if (particle == "photon") offset=0; |
1505 | else if (particle == "pi0") offset=5; |
1506 | else |
351dd634 |
1507 | AliError(Form("Wrong particle name: %s (choose from pi0/photon)\n", |
1508 | particle.Data())); |
88cb7938 |
1509 | |
1510 | if (param.Contains("a")) p=4+offset; |
1511 | else if(param.Contains("b")) p=5+offset; |
1512 | else if(param.Contains("c")) p=6+offset; |
1513 | else if(param.Contains("x0"))p=7+offset; |
1514 | else if(param.Contains("y0"))p=8+offset; |
351dd634 |
1515 | if((i>4)||(i<0)) { |
1516 | AliError(Form("No parameter with index %d", i)) ; |
1517 | } else if(p==-1) { |
1518 | AliError(Form("No parameter with name %s", param.Data() )) ; |
1519 | } else |
88cb7938 |
1520 | (*fParameters)(p,i) = par ; |
1521 | } |
1522 | |
1523 | //____________________________________________________________________________ |
1524 | void AliPHOSPIDv1::Unload() |
1525 | { |
c947e71a |
1526 | //Unloads RecPoints, Tracks and RecParticles |
88cb7938 |
1527 | AliPHOSGetter * gime = AliPHOSGetter::Instance() ; |
1528 | gime->PhosLoader()->UnloadRecPoints() ; |
1529 | gime->PhosLoader()->UnloadTracks() ; |
1530 | gime->PhosLoader()->UnloadRecParticles() ; |
1531 | } |
1532 | |
1533 | //____________________________________________________________________________ |
90cceaf6 |
1534 | void AliPHOSPIDv1::WriteRecParticles() |
88cb7938 |
1535 | { |
c947e71a |
1536 | //It writes reconstructed particles and pid to file |
1537 | |
88cb7938 |
1538 | AliPHOSGetter *gime = AliPHOSGetter::Instance() ; |
1539 | |
1540 | TClonesArray * recParticles = gime->RecParticles() ; |
1541 | recParticles->Expand(recParticles->GetEntriesFast() ) ; |
adcca1e6 |
1542 | if(fWrite){ |
1543 | TTree * treeP = gime->TreeP(); |
1544 | |
1545 | //First rp |
1546 | Int_t bufferSize = 32000 ; |
1547 | TBranch * rpBranch = treeP->Branch("PHOSRP",&recParticles,bufferSize); |
1548 | rpBranch->SetTitle(BranchName()); |
1549 | |
1550 | rpBranch->Fill() ; |
1551 | |
1552 | gime->WriteRecParticles("OVERWRITE"); |
1553 | gime->WritePID("OVERWRITE"); |
1554 | } |
88cb7938 |
1555 | } |
1556 | |
35adb638 |
1557 | |
1558 | //_______________________________________________________________________ |
1559 | void AliPHOSPIDv1::SetInitPID(const Double_t *p) { |
1560 | // Sets values for the initial population of each particle type |
304864ab |
1561 | for (Int_t i=0; i<AliPID::kSPECIESN; i++) fInitPID[i] = p[i]; |
35adb638 |
1562 | } |
1563 | //_______________________________________________________________________ |
1564 | void AliPHOSPIDv1::GetInitPID(Double_t *p) const { |
1565 | // Gets values for the initial population of each particle type |
304864ab |
1566 | for (Int_t i=0; i<AliPID::kSPECIESN; i++) p[i] = fInitPID[i]; |
35adb638 |
1567 | } |