1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //_________________________________________________________________________
19 // Implementation version v1 of the PHOS particle identifier
20 // Particle identification based on the
21 // - RCPV: distance from CPV recpoint to EMCA recpoint.
23 // - PCA: Principal Components Analysis..
24 // The identified particle has an identification number corresponding
25 // to a 9 bits number:
26 // -Bit 0 to 2: bit set if RCPV > CpvEmcDistance (each bit corresponds
27 // to a different efficiency-purity point of the photon identification)
28 // -Bit 3 to 5: bit set if TOF < TimeGate (each bit corresponds
29 // to a different efficiency-purity point of the photon identification)
30 // -Bit 6 to 9: bit set if Principal Components are
31 // inside an ellipse defined by the parameters a, b, c, x0 and y0.
32 // (each bit corresponds to a different efficiency-purity point of the
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
36 // file is $ALICE_ROOT/PHOS/PCA8pa15_0.5-100.root. Analysis done for
37 // energies between 0.5 and 100 GeV.
38 // A calibrated energy is calculated. The energy of the reconstructed
39 // cluster is corrected with the formula A + B * E + C * E^2, whose
40 // parameters where obtained through the study of the reconstructed
41 // energy distribution of monoenergetic photons.
43 // All the parameters (RCPV(2 rows-3 columns),TOF(1r-3c),PCA(5r-4c)
44 // and calibration(1r-3c))are stored in a file called
45 // $ALICE_ROOT/PHOS/Parameters.dat. Each time that AliPHOSPIDv1 is
46 // initialized, this parameters are copied to a Matrix (9,4), a
50 // root [0] AliPHOSPIDv1 * p = new AliPHOSPIDv1("galice1.root")
51 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
52 // // reading headers from file galice1.root and create RecParticles
53 // TrackSegments and RecPoints are used
54 // // set file name for the branch RecParticles
55 // root [1] p->ExecuteTask("deb all time")
56 // // available options
57 // // "deb" - prints # of reconstructed particles
58 // // "deb all" - prints # and list of RecParticles
59 // // "time" - prints benchmarking results
61 // root [2] AliPHOSPIDv1 * p2 = new AliPHOSPIDv1("galice1.root","v1",kTRUE)
62 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
64 // root [3] p2->ExecuteTask()
68 //*-- Author: Yves Schutz (SUBATECH) & Gines Martinez (SUBATECH) &
69 // Gustavo Conesa April 2002
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
74 // of fixing 3 ellipses for different ranges of energy, it has been
75 // studied the dependency of the ellipses parameters with the
76 // energy, and they are implemented in the code as a funtion
81 // --- ROOT system ---
86 #include "TBenchmark.h"
88 #include "TPrincipal.h"
90 // --- Standard library ---
92 //#include <Riostream.h>
94 // --- AliRoot header files ---
96 #include "AliGenerator.h"
97 #include "AliPHOSPIDv1.h"
98 #include "AliPHOSTrackSegment.h"
99 #include "AliPHOSRecParticle.h"
100 #include "AliPHOSGeometry.h"
101 #include "AliPHOSGetter.h"
103 ClassImp( AliPHOSPIDv1)
105 //____________________________________________________________________________
106 AliPHOSPIDv1::AliPHOSPIDv1():AliPHOSPID()
111 fDefaultInit = kTRUE ;
115 //____________________________________________________________________________
116 AliPHOSPIDv1::AliPHOSPIDv1(AliPHOSPIDv1 & pid ):AliPHOSPID(pid)
121 fDefaultInit = kFALSE ;
125 //____________________________________________________________________________
126 AliPHOSPIDv1::AliPHOSPIDv1(const char * headerFile,const char * name, const Bool_t toSplit)
127 :AliPHOSPID(headerFile, name,toSplit)
129 //ctor with the indication on where to look for the track segments
134 fDefaultInit = kFALSE ;
138 //____________________________________________________________________________
139 AliPHOSPIDv1::~AliPHOSPIDv1()
142 // fDefaultInit = kTRUE if PID created by default ctor (to get just the parameters)
144 delete [] fX ; // Principal input
145 delete [] fP ; // Principal components
146 delete [] fPPi0 ; // Pi0 Principal components
149 // AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
150 // remove the task from the folder list
151 // gime->RemoveTask("P",GetName()) ;
152 // TString name(GetName()) ;
153 // name.ReplaceAll("pid", "clu") ;
154 // gime->RemoveTask("C",name) ;
156 // // remove the data from the folder list
157 // name = GetName() ;
158 // name.Remove(name.Index(":")) ;
159 // gime->RemoveObjects("RE", name) ; // EMCARecPoints
160 // gime->RemoveObjects("RC", name) ; // CPVRecPoints
161 // gime->RemoveObjects("T", name) ; // TrackSegments
162 // gime->RemoveObjects("P", name) ; // RecParticles
165 // gime->CloseFile() ;
171 //____________________________________________________________________________
172 const TString AliPHOSPIDv1::BranchName() const
174 // gives the name of the current branch
175 TString branchName(GetName() ) ;
176 branchName.Remove(branchName.Index(Version())-1) ;
180 //____________________________________________________________________________
181 void AliPHOSPIDv1::Init()
183 // Make all memory allocations that are not possible in default constructor
184 // Add the PID task to the list of PHOS tasks
186 if ( strcmp(GetTitle(), "") == 0 )
187 SetTitle("galice.root") ;
189 TString branchname(GetName()) ;
190 branchname.Remove(branchname.Index(Version())-1) ;
191 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(),branchname.Data(),fToSplit ) ;
193 // gime->SetRecParticlesTitle(BranchName()) ;
195 Error("Init", "Could not obtain the Getter object !" ) ;
201 //First - extract full path if necessary
202 TString fileName(GetTitle()) ;
203 Ssiz_t islash = fileName.Last('/') ;
204 if(islash<fileName.Length())
205 fileName.Remove(islash+1,fileName.Length()) ;
208 fileName+="PHOS.RecData." ;
209 if((strcmp(branchname.Data(),"Default")!=0)&&(strcmp(branchname.Data(),"")!=0)){
210 fileName+=branchname.Data() ;
214 fSplitFile = static_cast<TFile*>(gROOT->GetFile(fileName.Data()));
216 fSplitFile = TFile::Open(fileName.Data(),"update") ;
219 gime->PostPID(this) ;
220 // create a folder on the white board //YSAlice/WhiteBoard/RecParticles/PHOS/recparticlesName
221 gime->PostRecParticles(branchname) ;
225 //____________________________________________________________________________
226 void AliPHOSPIDv1::InitParameters()
229 // fHeaderFileName = GetTitle() ;
230 // TString name(GetName()) ;
231 // if (name.IsNull())
232 // name = "Default" ;
233 // fTrackSegmentsTitle = name ;
234 // fRecPointsTitle = name ;
235 // fRecParticlesTitle = name ;
236 // name.Append(":") ;
237 // name.Append(Version()) ;
239 fRecParticlesInRun = 0 ;
241 // fClusterizer = 0 ;
243 fRecParticlesInRun = 0 ;
244 TString pidName( GetName()) ;
245 if (pidName.IsNull() )
246 pidName = "Default" ;
247 pidName.Append(":") ;
248 pidName.Append(Version()) ;
250 fPi0Analysis = kFALSE ;
251 SetParameters() ; // fill the parameters matrix from parameters file
254 //____________________________________________________________________________
255 const Float_t AliPHOSPIDv1::GetCpvtoEmcDistanceCut(const Float_t e, const TString Axis) const
257 // Get CpvtoEmcDistance Cut depending on the cluster energy, axis and
258 // Purity-Efficiency point
261 if (Axis.Contains("X")) i = 1;
262 else if (Axis.Contains("Z")) i = 2;
264 Error("GetCpvtoEmcDistanceCut"," Invalid axis option ");
266 Float_t a = (*fParameters)(i,0) ;
267 Float_t b = (*fParameters)(i,1) ;
268 Float_t c = (*fParameters)(i,2) ;
270 Float_t sig = a + TMath::Exp(b-c*e);
273 //____________________________________________________________________________
275 const Double_t AliPHOSPIDv1::GetTimeGate(const Int_t effpur) const
277 // Get TimeGate parameter depending on Purity-Efficiency point
279 if(effpur>2 || effpur<0)
280 Error("GetTimeGate","Invalid Efficiency-Purity choice %d",effpur);
281 return (*fParameters)(3,effpur) ;
284 //_____________________________________________________________________________
285 const Float_t AliPHOSPIDv1::GetDistance(AliPHOSEmcRecPoint * emc,AliPHOSRecPoint * cpv, Option_t * Axis)const
287 // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
289 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
293 emc->GetLocalPosition(vecEmc) ;
294 cpv->GetLocalPosition(vecCpv) ;
295 if(emc->GetPHOSMod() == cpv->GetPHOSMod()){
296 // Correct to difference in CPV and EMC position due to different distance to center.
297 // we assume, that particle moves from center
298 Float_t dCPV = geom->GetIPtoOuterCoverDistance();
299 Float_t dEMC = geom->GetIPtoCrystalSurface() ;
301 vecCpv = dEMC * vecCpv - vecEmc ;
302 if (Axis == "X") return vecCpv.X();
303 if (Axis == "Y") return vecCpv.Y();
304 if (Axis == "Z") return vecCpv.Z();
305 if (Axis == "R") return vecCpv.Mag();
312 //____________________________________________________________________________
313 const Int_t AliPHOSPIDv1::GetCPVBit(AliPHOSEmcRecPoint * emc,AliPHOSRecPoint * cpv,const Int_t EffPur, const Float_t e) const
315 if(EffPur>2 || EffPur<0)
316 Error("GetCPVBit","Invalid Efficiency-Purity choice %d",EffPur);
318 Float_t sigX = GetCpvtoEmcDistanceCut(e,"X");
319 Float_t sigZ = GetCpvtoEmcDistanceCut(e,"Z");
321 Float_t deltaX = TMath::Abs(GetDistance(emc, cpv, "X"));
322 Float_t deltaZ = TMath::Abs(GetDistance(emc, cpv, "Z"));
324 if((deltaX>sigX*(EffPur+1)) || (deltaZ>sigZ*(EffPur+1)))
331 //____________________________________________________________________________
332 const Double_t AliPHOSPIDv1::GetCalibratedEnergy(const Float_t e) const
334 // It calibrates Energy depending on the recpoint energy.
335 // The energy of the reconstructed cluster is corrected with
336 // the formula A + B* E + C* E^2, whose parameters where obtained
337 // through the study of the reconstructed energy distribution of
338 // monoenergetic photons.
340 Double_t p[]={0.,0.,0.};
342 for(i=0;i<3;i++) p[i]= (*fParameters)(0,i);
343 Double_t enerec = p[0] + p[1]* e+ p[2] * e * e;
347 //____________________________________________________________________________
348 const Int_t AliPHOSPIDv1::GetPrincipalBit(const Double_t* p ,const Int_t effpur, const Float_t e)const
350 //Is the particle inside de PCA ellipse?
353 Double_t a = GetEllipseParameter("a", e);
354 Double_t b = GetEllipseParameter("b", e);
355 Double_t c = GetEllipseParameter("c", e);
356 Double_t xCenter = GetEllipseParameter("x0", e);
357 Double_t yCenter = GetEllipseParameter("y0", e);
359 Double_t r = TMath::Power((p[0] - xCenter)/a,2) +
360 TMath::Power((p[1] - yCenter)/b,2) +
361 c*(p[0] - xCenter)*(p[1] - yCenter)/(a*b) ;
362 //3 different ellipses defined
363 if((effpur==2)&&(r <1./2.)) prinbit= 1;
364 if((effpur==1)&&(r <2. )) prinbit= 1;
365 if((effpur==0)&&(r <9./2.)) prinbit= 1;
368 Error("GetPrincipalBit", "Negative square? R=%f \n",r) ;
373 //____________________________________________________________________________
374 const Int_t AliPHOSPIDv1::GetPrincipalPi0Bit(const Double_t* p, const Int_t effpur, const Float_t e)const
376 //Is the particle inside de Pi0 PCA ellipse?
379 Double_t a = GetEllipseParameterPi0("a", e);
380 Double_t b = GetEllipseParameterPi0("b", e);
381 Double_t c = GetEllipseParameterPi0("c", e);
382 Double_t xCenter = GetEllipseParameterPi0("x0", e);
383 Double_t yCenter = GetEllipseParameterPi0("y0", e);
385 Double_t r = TMath::Power((p[0] - xCenter)/a,2) +
386 TMath::Power((p[1] - yCenter)/b,2) +
387 c*(p[0] - xCenter)*(p[1] - yCenter)/(a*b) ;
388 //3 different ellipses defined
389 if((effpur==2)&&(r <1./2.)) prinbit= 1;
390 if((effpur==1)&&(r <2. )) prinbit= 1;
391 if((effpur==0)&&(r <9./2.)) prinbit= 1;
394 Error("GetPrincipalPi0Bit", "Negative square?") ;
399 //_____________________________________________________________________________
400 void AliPHOSPIDv1::SetCpvtoEmcDistanceCutParameters(Float_t e, Int_t effpur, TString Axis,Float_t cut)
402 // Set the parameters to calculate Cpvto EmcDistanceCut depending on the cluster energy and
403 // Purity-Efficiency point
405 if(effpur>2 || effpur<0)
406 Error("SetCpvtoEmcDistanceCutParameters","Invalid Efficiency-Purity choice %d",effpur);
409 if (Axis.Contains("X")) i = 1;
410 else if(Axis.Contains("Z")) i = 2;
412 Error("SetCpvtoEmcDistanceCutParameters"," Invalid axis option");
414 (*fParameters)(i,effpur) = cut ;
416 //_____________________________________________________________________________
417 void AliPHOSPIDv1::SetTimeGate(Int_t effpur, Float_t gate)
419 // Set the parameter TimeGate depending on the cluster energy and
420 // Purity-Efficiency point
421 if(effpur>2 || effpur<0)
422 Error("SetTimeGate","Invalid Efficiency-Purity choice %d",effpur);
424 (*fParameters)(3,effpur)= gate ;
426 //_____________________________________________________________________________
427 void AliPHOSPIDv1::SetParameters()
429 // PCA : To do the Principal Components Analysis it is necessary
430 // the Principal file, which is opened here
431 fX = new double[7]; // Data for the PCA
432 fP = new double[7]; // Eigenvalues of the PCA
433 fPPi0 = new double[7]; // Eigenvalues of the Pi0 PCA
435 // Read photon principals from the photon file
437 fFileName = "$ALICE_ROOT/PHOS/PCA8pa15_0.5-100.root" ;
438 TFile f( fFileName.Data(), "read" ) ;
439 fPrincipal = dynamic_cast<TPrincipal*> (f.Get("principal")) ;
442 // Read pi0 principals from the pi0 file
444 fFileNamePi0 = "$ALICE_ROOT/PHOS/PCA_pi0_40-120.root" ;
445 TFile fPi0( fFileNamePi0.Data(), "read" ) ;
446 fPrincipalPi0 = dynamic_cast<TPrincipal*> (fPi0.Get("principal")) ;
449 // Open parameters file and initialization of the Parameters matrix.
450 // In the File Parameters.dat are all the parameters. These are introduced
451 // in a matrix of 9x4
453 // All the parameters defined in this file are, in order of row:
454 // -CpvtoEmcDistanceCut (2 row (x and z) and 3 columns, each one depending
455 // on the parameter of the funtion that sets the cut in x or z.
456 // -TimeGate, 1 row and 3 columns (3 efficiency-purty cuts)
457 // -PCA, parameters of the functions that
458 // calculate the ellipse parameters, x0,y0,a, b, c. These 5 parameters
459 // (5 rows) depend on 4 parameters (columns).
460 // -Finally there is a row with the energy calibration parameters,
463 fFileNamePar = gSystem->ExpandPathName("$ALICE_ROOT/PHOS/Parameters.dat");
464 fParameters = new TMatrixD(14,4) ;
465 const Int_t kmaxLeng=255;
466 char string[kmaxLeng];
468 // Open a text file with PID parameters
469 FILE *fd = fopen(fFileNamePar.Data(),"r");
471 Fatal("SetParameter","File %s with a PID parameters cannot be opened\n",
472 fFileNamePar.Data());
475 // Read parameter file line-by-line and skip empty line and comments
476 while (fgets(string,kmaxLeng,fd) != NULL) {
477 if (string[0] == '\n' ) continue;
478 if (string[0] == '!' ) continue;
479 sscanf(string, "%lf %lf %lf %lf",
480 &(*fParameters)(i,0), &(*fParameters)(i,1),
481 &(*fParameters)(i,2), &(*fParameters)(i,3));
488 //________________________________________________________________________
489 void AliPHOSPIDv1::SetEllipseParameter(TString Param, Int_t i, Double_t par)
491 // Set the parameter "i" that is needed to calculate the ellipse
492 // parameter "Param".
495 if (Param.Contains("a")) p=4;
496 else if(Param.Contains("b")) p=5;
497 else if(Param.Contains("c")) p=6;
498 else if(Param.Contains("x0"))p=7;
499 else if(Param.Contains("y0"))p=8;
501 Error("SetEllipseParameter", "No parameter with index %d", i) ;
503 Error("SetEllipseParameter", "No parameter with name %s", Param.Data() ) ;
505 (*fParameters)(p,i) = par ;
507 //________________________________________________________________________
508 void AliPHOSPIDv1::SetEllipseParameterPi0(TString Param, Int_t i, Double_t par)
510 // Set the parameter "i" that is needed to calculate the ellipse
511 // parameter "Param".
512 if(!fPi0Analysis) Error("SetPi0EllipseParameter", "Pi 0 Analysis is off") ;
514 if (Param.Contains("a")) p=9;
515 else if(Param.Contains("b")) p=10;
516 else if(Param.Contains("c")) p=11;
517 else if(Param.Contains("x0"))p=12;
518 else if(Param.Contains("y0"))p=13;
520 Error("SetPi0EllipseParameter", "No parameter with index %d", i) ;
522 Error("SetPi0EllipseParameter", "No parameter with name %s", Param.Data() ) ;
524 (*fParameters)(p,i) = par ;
526 //________________________________________________________________________
527 const Double_t AliPHOSPIDv1::GetParameterToCalculateEllipse(const TString Param, const Int_t i) const
529 // Get the parameter "i" that is needed to calculate the ellipse
530 // parameter "Param".
535 if (Param.Contains("a")) p=4;
536 else if(Param.Contains("b")) p=5;
537 else if(Param.Contains("c")) p=6;
538 else if(Param.Contains("x0"))p=7;
539 else if(Param.Contains("y0"))p=8;
542 Error("GetParameterToCalculateEllipse", "No parameter with index", i) ;
544 Error("GetParameterToCalculateEllipse", "No parameter with name %s", Param.Data() ) ;
546 par = (*fParameters)(p,i) ;
551 //____________________________________________________________________________
552 const Double_t AliPHOSPIDv1::GetParameterToCalculatePi0Ellipse(const TString Param, const Int_t i) const
554 // Get the parameter "i" that is needed to calculate the ellipse
555 // parameter "Param".
557 if(!fPi0Analysis) Error("GetParameterToCalculatePi0Ellipse", "Pi 0 Analysis is off") ;
562 if(Param.Contains("a")) p=9;
563 if(Param.Contains("b")) p=10;
564 if(Param.Contains("c")) p=11;
565 if(Param.Contains("x0"))p=12;
566 if(Param.Contains("y0"))p=13;
569 Error("GetParameterToCalculatePi0Ellipse", "No parameter with index", i) ;
571 Error("GetParameterToCalculatePi0Ellipse", "No parameter with name %s", Param.Data() ) ;
573 par = (*fParameters)(p,i) ;
578 //____________________________________________________________________________
579 void AliPHOSPIDv1::SetCalibrationParameter(Int_t i,Double_t param) const
581 (*fParameters)(0,i) = param ;
583 //____________________________________________________________________________
584 const Double_t AliPHOSPIDv1::GetCalibrationParameter(const Int_t i) const
586 Float_t param = (*fParameters)(0,i);
589 //____________________________________________________________________________
590 const Double_t AliPHOSPIDv1::GetEllipseParameter(const TString Param,Float_t E) const
592 // Calculates the parameter Param of the ellipse
594 Double_t p[4]={0.,0.,0.,0.};
595 Double_t value = 0.0;
598 if(Param.Contains("a")){
599 for(i=0;i<4;i++)p[i]=(*fParameters)(4,i);
603 else if(Param.Contains("b")){
604 for(i=0;i<4;i++)p[i]=(*fParameters)(5,i);
608 else if(Param.Contains("c"))
609 for(i=0;i<4;i++)p[i]=(*fParameters)(6,i);
611 else if(Param.Contains("x0")){
612 for(i=0;i<4;i++)p[i]=(*fParameters)(7,i);
615 else if(Param.Contains("y0"))
616 for(i=0;i<4;i++)p[i]=(*fParameters)(8,i);
618 value = p[0]/TMath::Sqrt(E)+p[1]*E+p[2]*E*E+p[3];
622 //____________________________________________________________________________
623 // const Double_t AliPHOSPIDv1::GetEllipseParameter(const TString Param,Float_t E) const
625 // // Calculates the parameter Param of the pi0 ellipse
627 // Double_t p[3] = {0.,0.,0.};
628 // Double_t value = 0.0;
631 // if(Param.Contains("a"))
632 // for(i=0;i<3;i++)p[i]=(*fParameters)(4,i);
633 // else if(Param.Contains("b"))
634 // for(i=0;i<3;i++)p[i]=(*fParameters)(5,i);
635 // else if(Param.Contains("c"))
636 // for(i=0;i<3;i++)p[i]=(*fParameters)(6,i);
637 // else if(Param.Contains("x0"))
638 // for(i=0;i<3;i++)p[i]=(*fParameters)(7,i);
639 // else if(Param.Contains("y0"))
640 // for(i=0;i<3;i++)p[i]=(*fParameters)(8,i);
642 // value = p[0] + p[1]*E + p[2]*E*E;
645 //____________________________________________________________________________
646 const Double_t AliPHOSPIDv1::GetEllipseParameterPi0(const TString Param,Float_t E) const
648 // Calculates the parameter Param of the pi0 ellipse
650 Double_t p[3] = {0.,0.,0.};
651 Double_t value = 0.0;
654 if(Param.Contains("a"))
655 for(i=0;i<3;i++)p[i]=(*fParameters)(9,i);
656 else if(Param.Contains("b"))
657 for(i=0;i<3;i++)p[i]=(*fParameters)(10,i);
658 else if(Param.Contains("c"))
659 for(i=0;i<3;i++)p[i]=(*fParameters)(11,i);
660 else if(Param.Contains("x0"))
661 for(i=0;i<3;i++)p[i]=(*fParameters)(12,i);
662 else if(Param.Contains("y0"))
663 for(i=0;i<3;i++)p[i]=(*fParameters)(13,i);
665 value = p[0] + p[1]*E + p[2]*E*E;
668 //____________________________________________________________________________
670 void AliPHOSPIDv1::Exec(Option_t * option)
674 if( strcmp(GetName(), "")== 0 )
677 if(strstr(option,"tim"))
678 gBenchmark->Start("PHOSPID");
680 if(strstr(option,"print")) {
686 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
687 if(gime->BranchExists("RecParticles") )
689 Int_t nevents = gime->MaxEvent() ;
692 for(ievent = 0; ievent < nevents; ievent++){
693 gime->Event(ievent,"R") ;
695 if(gime->TrackSegments() && //Skip events, where no track segments made
696 gime->TrackSegments()->GetEntriesFast()) {
698 WriteRecParticles(ievent);
699 if(strstr(option,"deb"))
700 PrintRecParticles(option) ;
701 //increment the total number of rec particles per run
702 fRecParticlesInRun+=gime->RecParticles(BranchName())->GetEntriesFast() ;
706 if(strstr(option,"tim")){
707 gBenchmark->Stop("PHOSPID");
708 Info("Exec", "took %f seconds for PID %f seconds per event",
709 gBenchmark->GetCpuTime("PHOSPID"),
710 gBenchmark->GetCpuTime("PHOSPID")/nevents) ;
714 //____________________________________________________________________________
715 void AliPHOSPIDv1::MakeRecParticles(){
717 // Makes a RecParticle out of a TrackSegment
719 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
720 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
721 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
722 TClonesArray * trackSegments = gime->TrackSegments() ;
723 if ( !emcRecPoints || !cpvRecPoints || !trackSegments ) {
724 Fatal("MakeRecParticles", "RecPoints or TrackSegments not found !") ;
726 TClonesArray * recParticles = gime->RecParticles() ;
727 recParticles->Clear();
729 TIter next(trackSegments) ;
730 AliPHOSTrackSegment * ts ;
732 AliPHOSRecParticle * rp ;
733 while ( (ts = (AliPHOSTrackSegment *)next()) ) {
735 new( (*recParticles)[index] ) AliPHOSRecParticle() ;
736 rp = (AliPHOSRecParticle *)recParticles->At(index) ;
737 rp->SetTrackSegment(index) ;
738 rp->SetIndexInList(index) ;
740 AliPHOSEmcRecPoint * emc = 0 ;
741 if(ts->GetEmcIndex()>=0)
742 emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ;
744 AliPHOSRecPoint * cpv = 0 ;
745 if(ts->GetCpvIndex()>=0)
746 cpv = (AliPHOSRecPoint *) cpvRecPoints->At(ts->GetCpvIndex()) ;
748 // Now set type (reconstructed) of the particle
750 // Choose the cluster energy range
753 Fatal("MakeRecParticles", "-> emc(%d) = %d", ts->GetEmcIndex(), emc ) ;
756 Float_t e = emc->GetEnergy() ;
759 emc->GetElipsAxis(lambda) ;
761 if((lambda[0]>0.01) && (lambda[1]>0.01)){
762 // Looking PCA. Define and calculate the data (X),
763 // introduce in the function X2P that gives the components (P).
766 Float_t emaxdTotal = 0. ;
768 if((lambda[0]+lambda[1])!=0)
769 spher=fabs(lambda[0]-lambda[1])/(lambda[0]+lambda[1]);
771 emaxdTotal=emc->GetMaximalEnergy()/emc->GetEnergy();
775 fX[2] = emc->GetDispersion() ;
777 fX[4] = emc->GetMultiplicity() ;
779 fX[6] = emc->GetCoreEnergy() ;
781 fPrincipal->X2P(fX,fP);
783 fPrincipalPi0->X2P(fX,fPPi0);
787 fP[0]=-100.0; //We do not accept clusters with
788 fP[1]=-100.0; //one cell as a photon-like
795 Float_t time =emc->GetTime() ;
797 // Loop of Efficiency-Purity (the 3 points of purity or efficiency
798 // are taken into account to set the particle identification)
799 for(Int_t eff_pur = 0; eff_pur < 3 ; eff_pur++){
801 // Looking at the CPV detector. If RCPV greater than CpvEmcDistance,
802 // 1st,2nd or 3rd bit (depending on the efficiency-purity point )
804 if(GetCPVBit(emc, cpv, eff_pur,e) == 1 )
805 rp->SetPIDBit(eff_pur) ;
807 // Looking the TOF. If TOF smaller than gate, 4th, 5th or 6th
808 // bit (depending on the efficiency-purity point )is set to 1
809 if(time< (*fParameters)(2,eff_pur))
810 rp->SetPIDBit(eff_pur+3) ;
812 //If we are inside the ellipse, 7th, 8th or 9th
813 // bit (depending on the efficiency-purity point )is set to 1
814 if(GetPrincipalBit(fP,eff_pur,e) == 1)
815 rp->SetPIDBit(eff_pur+6) ;
818 //If we are inside the ellipse, 10th, 11th or 12th
819 // bit (depending on the efficiency-purity point )is set to 1
821 if(GetPrincipalPi0Bit(fPPi0,eff_pur,e) == 1)
822 rp->SetPIDBit(eff_pur+9) ;
827 //Set momentum, energy and other parameters
828 Float_t encal = GetCalibratedEnergy(e);
829 TVector3 dir = GetMomentumDirection(emc,cpv) ;
831 rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),encal) ;
833 rp->Name(); //If photon sets the particle pdg name to gamma
834 rp->SetProductionVertex(0,0,0,0);
835 rp->SetFirstMother(-1);
836 rp->SetLastMother(-1);
837 rp->SetFirstDaughter(-1);
838 rp->SetLastDaughter(-1);
839 rp->SetPolarisation(0,0,0);
845 //____________________________________________________________________________
846 void AliPHOSPIDv1::Print()
848 // Print the parameters used for the particle type identification
851 message = "\n=============== AliPHOSPID1 ================\n" ;
852 message += "Making PID\n";
853 message += " Pricipal analysis file from 0.5 to 100 %s\n" ;
854 message += " Name of parameters file %s\n" ;
855 message += " Matrix of Parameters: 14x4\n" ;
856 message += " Energy Calibration 1x3 [3 parametres to calibrate energy: A + B* E + C * E^2]\n" ;
857 message += " RCPV 2x3 rows x and z, columns function cut parameters\n" ;
858 message += " TOF 1x3 [High Eff-Low Pur,Medium Eff-Pur, Low Eff-High Pur]\n" ;
859 message += " PCA 5x4 [5 ellipse parametres and 4 parametres to calculate them: A/Sqrt(E) + B* E + C * E^2 + D]\n" ;
860 message += " Pi0 PCA 5x3 [5 ellipse parametres and 3 parametres to calculate them: A + B* E + C * E^2]\n" ;
861 Info("Print", message.Data(), fFileName.Data(), fFileNamePar.Data() ) ;
862 fParameters->Print() ;
865 //____________________________________________________________________________
866 void AliPHOSPIDv1::WriteRecParticles(Int_t event)
868 // writes the reconstructed particles to file
869 AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ;
871 TClonesArray * recParticles = gime->RecParticles() ;
872 recParticles->Expand(recParticles->GetEntriesFast() ) ;
880 sprintf(name,"%s%d", "TreeR",event) ;
881 treeR = dynamic_cast<TTree*>(fSplitFile->Get(name));
884 treeR = gAlice->TreeR();
888 gAlice->MakeTree("R", fSplitFile);
889 treeR = gAlice->TreeR() ;
893 Int_t bufferSize = 32000 ;
894 TBranch * rpBranch = treeR->Branch("PHOSRP",&recParticles,bufferSize);
895 rpBranch->SetTitle(BranchName());
899 Int_t splitlevel = 0 ;
900 AliPHOSPIDv1 * pid = this ;
901 TBranch * pidBranch = treeR->Branch("AliPHOSPID","AliPHOSPIDv1",&pid,bufferSize,splitlevel);
902 pidBranch->SetTitle(BranchName());
907 treeR->AutoSave() ; //Write(0,kOverwrite) ;
908 if(gAlice->TreeR()!=treeR){
913 //____________________________________________________________________________
914 TVector3 AliPHOSPIDv1::GetMomentumDirection(AliPHOSEmcRecPoint * emc, AliPHOSRecPoint * cpv)const
916 // Calculates the momentum direction:
917 // 1. if only a EMC RecPoint, direction is given by IP and this RecPoint
918 // 2. if a EMC RecPoint and CPV RecPoint, direction is given by the line through the 2 recpoints
919 // However because of the poor position resolution of PPSD the direction is always taken as if we were
922 TVector3 dir(0,0,0) ;
924 TVector3 emcglobalpos ;
927 emc->GetGlobalPosition(emcglobalpos, dummy) ;
930 dir.SetZ( -dir.Z() ) ; // why ?
932 //account correction to the position of IP
933 Float_t xo,yo,zo ; //Coordinates of the origin
934 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
935 TVector3 origin(xo,yo,zo);
940 //____________________________________________________________________________
941 void AliPHOSPIDv1::PrintRecParticles(Option_t * option)
943 // Print table of reconstructed particles
945 AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ;
947 TClonesArray * recParticles = gime->RecParticles(BranchName()) ;
950 message = "\nevent " ;
951 message += gAlice->GetEvNumber() ;
952 message += " found " ;
953 message += recParticles->GetEntriesFast();
954 message += " RecParticles\n" ;
956 if(strstr(option,"all")) { // printing found TS
957 message += "\n PARTICLE Index \n" ;
960 for (index = 0 ; index < recParticles->GetEntries() ; index++) {
961 AliPHOSRecParticle * rp = (AliPHOSRecParticle * ) recParticles->At(index) ;
963 message += rp->Name().Data() ;
965 message += rp->GetIndexInList() ;
967 message += rp->GetType() ;
970 Info("Print", message.Data() ) ;