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