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 ---
90 #include "TBenchmark.h"
92 #include "TPrincipal.h"
95 // --- Standard library ---
98 // --- AliRoot header files ---
100 #include "AliGenerator.h"
102 #include "AliPHOSPIDv1.h"
103 #include "AliPHOSClusterizerv1.h"
104 #include "AliPHOSEmcRecPoint.h"
105 #include "AliPHOSTrackSegment.h"
106 #include "AliPHOSTrackSegmentMakerv1.h"
107 #include "AliPHOSRecParticle.h"
108 #include "AliPHOSGeometry.h"
109 #include "AliPHOSGetter.h"
111 ClassImp( AliPHOSPIDv1)
113 //____________________________________________________________________________
114 AliPHOSPIDv1::AliPHOSPIDv1():AliPHOSPID()
119 fDefaultInit = kTRUE ;
122 //____________________________________________________________________________
123 AliPHOSPIDv1::AliPHOSPIDv1(const AliPHOSPIDv1 & pid ):AliPHOSPID(pid)
131 //____________________________________________________________________________
132 AliPHOSPIDv1::AliPHOSPIDv1(const TString alirunFileName, const TString eventFolderName):AliPHOSPID(alirunFileName, eventFolderName)
134 //ctor with the indication on where to look for the track segments
138 fDefaultInit = kFALSE ;
141 //____________________________________________________________________________
142 AliPHOSPIDv1::~AliPHOSPIDv1()
146 delete [] fX ; // Principal input
147 delete [] fPPhoton ; // Photon Principal components
148 delete [] fPPi0 ; // Pi0 Principal components
150 //____________________________________________________________________________
151 const TString AliPHOSPIDv1::BranchName() const
157 //____________________________________________________________________________
158 void AliPHOSPIDv1::Init()
160 // Make all memory allocations that are not possible in default constructor
161 // Add the PID task to the list of PHOS tasks
163 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data()) ;
166 gime->PostPID(this) ;
169 //____________________________________________________________________________
170 void AliPHOSPIDv1::InitParameters()
172 // Initialize PID parameters
173 fRecParticlesInRun = 0 ;
175 fRecParticlesInRun = 0 ;
176 SetParameters() ; // fill the parameters matrix from parameters file
177 SetEventRange(0,-1) ;
180 //________________________________________________________________________
181 void AliPHOSPIDv1::Exec(Option_t *option)
183 // Steering method to perform particle reconstruction and identification
184 // for the event range from fFirstEvent to fLastEvent.
185 // This range is optionally set by SetEventRange().
186 // if fLastEvent=-1 (by default), then process events until the end.
188 if(strstr(option,"tim"))
189 gBenchmark->Start("PHOSPID");
191 if(strstr(option,"print")) {
197 AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;
199 if (fLastEvent == -1)
200 fLastEvent = gime->MaxEvent() - 1 ;
202 fLastEvent = TMath::Min(fLastEvent,gime->MaxEvent());
203 Int_t nEvents = fLastEvent - fFirstEvent + 1;
206 for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
207 gime->Event(ievent,"TR") ;
208 if(gime->TrackSegments() && //Skip events, where no track segments made
209 gime->TrackSegments()->GetEntriesFast()) {
212 if(strstr(option,"deb"))
213 PrintRecParticles(option) ;
214 //increment the total number of rec particles per run
215 fRecParticlesInRun += gime->RecParticles()->GetEntriesFast() ;
218 if(strstr(option,"tim")){
219 gBenchmark->Stop("PHOSPID");
220 Info("Exec", "took %f seconds for PID %f seconds per event",
221 gBenchmark->GetCpuTime("PHOSPID"),
222 gBenchmark->GetCpuTime("PHOSPID")/nEvents) ;
228 //____________________________________________________________________________
229 const TString AliPHOSPIDv1::GetFileNamePrincipal(TString particle) const
231 //Get file name that contains the PCA for a particle ("photon or pi0")
234 if (particle=="photon") name = fFileNamePrincipalPhoton ;
235 else if (particle=="pi0" ) name = fFileNamePrincipalPi0 ;
236 else Error("GetFileNamePrincipal","Wrong particle name: %s (choose from pi0/photon)\n",particle.Data());
240 //____________________________________________________________________________
241 Float_t AliPHOSPIDv1::GetParameterCalibration(Int_t i) const
243 // Get the i-th parameter "Calibration"
246 Error("GetParameterCalibration","Invalid parameter number: %d",i);
248 param = (*fParameters)(0,i);
252 //____________________________________________________________________________
253 Float_t AliPHOSPIDv1::GetCalibratedEnergy(Float_t e) const
255 // It calibrates Energy depending on the recpoint energy.
256 // The energy of the reconstructed cluster is corrected with
257 // the formula A + B* E + C* E^2, whose parameters where obtained
258 // through the study of the reconstructed energy distribution of
259 // monoenergetic photons.
261 Float_t p[]={0.,0.,0.};
262 for (Int_t i=0; i<3; i++) p[i] = GetParameterCalibration(i);
263 Float_t enerec = p[0] + p[1]*e + p[2]*e*e;
268 //____________________________________________________________________________
269 Float_t AliPHOSPIDv1::GetParameterCpv2Emc(Int_t i, TString axis) const
271 // Get the i-th parameter "CPV-EMC distance" for the specified axis
274 Error("GetParameterCpv2Emc","Invalid parameter number: %d",i);
277 if (axis == "x") param = (*fParameters)(1,i);
278 else if (axis == "z") param = (*fParameters)(2,i);
279 else Error("GetParameterCpv2Emc","Invalid axis name: %s",axis.Data());
284 //____________________________________________________________________________
285 Float_t AliPHOSPIDv1::GetCpv2EmcDistanceCut(TString axis, Float_t e) const
287 // Get CpvtoEmcDistance Cut depending on the cluster energy, axis and
288 // Purity-Efficiency point
291 Float_t p[]={0.,0.,0.};
292 for (Int_t i=0; i<3; i++) p[i] = GetParameterCpv2Emc(i,axis);
293 Float_t sig = p[0] + TMath::Exp(p[1] - p[2]*e);
297 //____________________________________________________________________________
298 Float_t AliPHOSPIDv1::GetEllipseParameter(TString particle, TString param, Float_t e) const
300 // Calculates the parameter param of the ellipse
304 Float_t p[4]={0.,0.,0.,0.};
306 for (Int_t i=0; i<4; i++) p[i] = GetParameterToCalculateEllipse(particle,param,i);
307 if (particle == "photon") {
308 if (param.Contains("a")) e = TMath::Min((Double_t)e,70.);
309 else if (param.Contains("b")) e = TMath::Min((Double_t)e,70.);
310 else if (param.Contains("x0")) e = TMath::Max((Double_t)e,1.1);
313 value = p[0]/TMath::Sqrt(e) + p[1]*e + p[2]*e*e + p[3];
317 //_____________________________________________________________________________
318 Float_t AliPHOSPIDv1::GetParameterPhotonBoundary (Int_t i) const
320 // Get the parameter "i" to calculate the boundary on the moment M2x
321 // for photons at high p_T
324 Error("GetParameterPhotonBoundary","Wrong parameter number: %d\n",i);
326 param = (*fParameters)(14,i) ;
330 //____________________________________________________________________________
331 Float_t AliPHOSPIDv1::GetParameterPi0Boundary (Int_t i) const
333 // Get the parameter "i" to calculate the boundary on the moment M2x
334 // for pi0 at high p_T
337 Error("GetParameterPi0Boundary","Wrong parameter number: %d\n",i);
339 param = (*fParameters)(15,i) ;
343 //____________________________________________________________________________
344 Float_t AliPHOSPIDv1::GetParameterTimeGate(Int_t i) const
346 // Get TimeGate parameter depending on Purity-Efficiency i:
347 // i=0 - Low purity, i=1 - Medium purity, i=2 - High purity
350 Error("GetParameterTimeGate","Invalid Efficiency-Purity choice %d",i);
352 param = (*fParameters)(3,i) ;
356 //_____________________________________________________________________________
357 Float_t AliPHOSPIDv1::GetParameterToCalculateEllipse(TString particle, TString param, Int_t i) const
359 // Get the parameter "i" that is needed to calculate the ellipse
360 // parameter "param" for the particle "particle" ("photon" or "pi0")
365 if (particle == "photon") offset=0;
366 else if (particle == "pi0") offset=5;
368 Error("GetParameterToCalculateEllipse","Wrong particle name: %s (choose from pi0/photon)\n",particle.Data());
373 if (param.Contains("a")) p=4+offset;
374 else if(param.Contains("b")) p=5+offset;
375 else if(param.Contains("c")) p=6+offset;
376 else if(param.Contains("x0"))p=7+offset;
377 else if(param.Contains("y0"))p=8+offset;
380 Error("GetParameterToCalculateEllipse", "No parameter with index", i) ;
382 Error("GetParameterToCalculateEllipse", "No parameter with name %s", param.Data() ) ;
384 par = (*fParameters)(p,i) ;
390 //____________________________________________________________________________
391 Float_t AliPHOSPIDv1::GetDistance(AliPHOSEmcRecPoint * emc,AliPHOSRecPoint * cpv, Option_t * axis)const
393 // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
395 const AliPHOSGeometry * geom = AliPHOSGetter::Instance()->PHOSGeometry() ;
399 emc->GetLocalPosition(vecEmc) ;
400 cpv->GetLocalPosition(vecCpv) ;
401 if(emc->GetPHOSMod() == cpv->GetPHOSMod()){
402 // Correct to difference in CPV and EMC position due to different distance to center.
403 // we assume, that particle moves from center
404 Float_t dCPV = geom->GetIPtoOuterCoverDistance();
405 Float_t dEMC = geom->GetIPtoCrystalSurface() ;
407 vecCpv = dEMC * vecCpv - vecEmc ;
408 if (axis == "X") return vecCpv.X();
409 if (axis == "Y") return vecCpv.Y();
410 if (axis == "Z") return vecCpv.Z();
411 if (axis == "R") return vecCpv.Mag();
417 //____________________________________________________________________________
418 Int_t AliPHOSPIDv1::GetCPVBit(AliPHOSEmcRecPoint * emc,AliPHOSRecPoint * cpv, Int_t effPur, Float_t e) const
420 if(effPur>2 || effPur<0)
421 Error("GetCPVBit","Invalid Efficiency-Purity choice %d",effPur);
423 Float_t sigX = GetCpv2EmcDistanceCut("X",e);
424 Float_t sigZ = GetCpv2EmcDistanceCut("Z",e);
426 Float_t deltaX = TMath::Abs(GetDistance(emc, cpv, "X"));
427 Float_t deltaZ = TMath::Abs(GetDistance(emc, cpv, "Z"));
429 if((deltaX>sigX*(effPur+1))|(deltaZ>sigZ*(effPur+1)))
435 //____________________________________________________________________________
436 Int_t AliPHOSPIDv1::GetPrincipalBit(TString particle, const Double_t* p, Int_t effPur, Float_t e)const
438 //Is the particle inside de PCA ellipse?
442 Float_t a = GetEllipseParameter(particle,"a" , e);
443 Float_t b = GetEllipseParameter(particle,"b" , e);
444 Float_t c = GetEllipseParameter(particle,"c" , e);
445 Float_t x0 = GetEllipseParameter(particle,"x0", e);
446 Float_t y0 = GetEllipseParameter(particle,"y0", e);
448 Float_t r = TMath::Power((p[0] - x0)/a,2) +
449 TMath::Power((p[1] - y0)/b,2) +
450 c*(p[0] - x0)*(p[1] - y0)/(a*b) ;
451 //3 different ellipses defined
452 if((effPur==2) && (r<1./2.)) prinbit= 1;
453 if((effPur==1) && (r<2. )) prinbit= 1;
454 if((effPur==0) && (r<9./2.)) prinbit= 1;
457 Error("GetPrincipalBit", "Negative square?") ;
462 //____________________________________________________________________________
463 Int_t AliPHOSPIDv1::GetHardPhotonBit(AliPHOSEmcRecPoint * emc) const
465 // Set bit for identified hard photons (E > 30 GeV)
466 // if the second moment M2x is below the boundary
468 Float_t e = emc->GetEnergy();
469 if (e < 30.0) return 0;
470 Float_t m2x = emc->GetM2x();
471 Float_t m2xBoundary = GetParameterPhotonBoundary(0) *
472 TMath::Exp(-TMath::Power(e-GetParameterPhotonBoundary(1),2)/2.0/
473 TMath::Power(GetParameterPhotonBoundary(2),2)) +
474 GetParameterPhotonBoundary(3);
475 Info("GetHardPhotonBit","E=%f, m2x=%f, boundary=%f",e,m2x,m2xBoundary);
476 if (m2x < m2xBoundary)
477 return 1;// A hard photon
479 return 0;// Not a hard photon
482 //____________________________________________________________________________
483 Int_t AliPHOSPIDv1::GetHardPi0Bit(AliPHOSEmcRecPoint * emc) const
485 // Set bit for identified hard pi0 (E > 30 GeV)
486 // if the second moment M2x is above the boundary
488 Float_t e = emc->GetEnergy();
489 if (e < 30.0) return 0;
490 Float_t m2x = emc->GetM2x();
491 Float_t m2xBoundary = GetParameterPi0Boundary(0) +
492 e * GetParameterPi0Boundary(1);
493 Info("GetHardPi0Bit","E=%f, m2x=%f, boundary=%f",e,m2x,m2xBoundary);
494 if (m2x > m2xBoundary)
495 return 1;// A hard pi0
497 return 0;// Not a hard pi0
500 //____________________________________________________________________________
501 TVector3 AliPHOSPIDv1::GetMomentumDirection(AliPHOSEmcRecPoint * emc, AliPHOSRecPoint * )const
503 // Calculates the momentum direction:
504 // 1. if only a EMC RecPoint, direction is given by IP and this RecPoint
505 // 2. if a EMC RecPoint and CPV RecPoint, direction is given by the line through the 2 recpoints
506 // However because of the poor position resolution of PPSD the direction is always taken as if we were
509 TVector3 dir(0,0,0) ;
511 TVector3 emcglobalpos ;
514 emc->GetGlobalPosition(emcglobalpos, dummy) ;
520 //account correction to the position of IP
521 Float_t xo,yo,zo ; //Coordinates of the origin
522 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
523 TVector3 origin(xo,yo,zo);
529 //____________________________________________________________________________
530 void AliPHOSPIDv1::MakeRecParticles()
532 // Makes a RecParticle out of a TrackSegment
534 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
535 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
536 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
537 TClonesArray * trackSegments = gime->TrackSegments() ;
538 if ( !emcRecPoints || !cpvRecPoints || !trackSegments ) {
539 Fatal("MakeRecParticles", "RecPoints or TrackSegments not found !") ;
541 TClonesArray * recParticles = gime->RecParticles() ;
542 recParticles->Clear();
544 TIter next(trackSegments) ;
545 AliPHOSTrackSegment * ts ;
547 AliPHOSRecParticle * rp ;
548 while ( (ts = (AliPHOSTrackSegment *)next()) ) {
550 new( (*recParticles)[index] ) AliPHOSRecParticle() ;
551 rp = (AliPHOSRecParticle *)recParticles->At(index) ;
552 rp->SetTrackSegment(index) ;
553 rp->SetIndexInList(index) ;
555 AliPHOSEmcRecPoint * emc = 0 ;
556 if(ts->GetEmcIndex()>=0)
557 emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ;
559 AliPHOSRecPoint * cpv = 0 ;
560 if(ts->GetCpvIndex()>=0)
561 cpv = (AliPHOSRecPoint *) cpvRecPoints->At(ts->GetCpvIndex()) ;
563 // Now set type (reconstructed) of the particle
565 // Choose the cluster energy range
568 Fatal("MakeRecParticles", "-> emc(%d) = %d", ts->GetEmcIndex(), emc ) ;
571 Float_t e = emc->GetEnergy() ;
574 emc->GetElipsAxis(lambda) ;
576 if((lambda[0]>0.01) && (lambda[1]>0.01)){
577 // Looking PCA. Define and calculate the data (X),
578 // introduce in the function X2P that gives the components (P).
581 Float_t Emaxdtotal = 0. ;
583 if((lambda[0]+lambda[1])!=0)
584 Spher=fabs(lambda[0]-lambda[1])/(lambda[0]+lambda[1]);
586 Emaxdtotal=emc->GetMaximalEnergy()/emc->GetEnergy();
590 fX[2] = emc->GetDispersion() ;
592 fX[4] = emc->GetMultiplicity() ;
594 fX[6] = emc->GetCoreEnergy() ;
596 fPrincipalPhoton->X2P(fX,fPPhoton);
597 fPrincipalPi0 ->X2P(fX,fPPi0);
601 fPPhoton[0]=-100.0; //We do not accept clusters with
602 fPPhoton[1]=-100.0; //one cell as a photon-like
607 Float_t time =emc->GetTime() ;
609 // Loop of Efficiency-Purity (the 3 points of purity or efficiency
610 // are taken into account to set the particle identification)
611 for(Int_t effPur = 0; effPur < 3 ; effPur++){
613 // Looking at the CPV detector. If RCPV greater than CpvEmcDistance,
614 // 1st,2nd or 3rd bit (depending on the efficiency-purity point )
616 if(GetCPVBit(emc, cpv, effPur,e) == 1 )
617 rp->SetPIDBit(effPur) ;
619 // Looking the TOF. If TOF smaller than gate, 4th, 5th or 6th
620 // bit (depending on the efficiency-purity point )is set to 1
621 if(time< (*fParameters)(2,effPur))
622 rp->SetPIDBit(effPur+3) ;
625 //If we are inside the ellipse, 7th, 8th or 9th
626 // bit (depending on the efficiency-purity point )is set to 1
627 if(GetPrincipalBit("photon",fPPhoton,effPur,e) == 1)
628 rp->SetPIDBit(effPur+6) ;
631 //If we are inside the ellipse, 10th, 11th or 12th
632 // bit (depending on the efficiency-purity point )is set to 1
633 if(GetPrincipalBit("pi0" ,fPPi0 ,effPur,e) == 1)
634 rp->SetPIDBit(effPur+9) ;
636 if(GetHardPhotonBit(emc))
638 if(GetHardPi0Bit (emc))
641 //Set momentum, energy and other parameters
642 Float_t encal = GetCalibratedEnergy(e);
643 TVector3 dir = GetMomentumDirection(emc,cpv) ;
645 rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),encal) ;
647 rp->Name(); //If photon sets the particle pdg name to gamma
648 rp->SetProductionVertex(0,0,0,0);
649 rp->SetFirstMother(-1);
650 rp->SetLastMother(-1);
651 rp->SetFirstDaughter(-1);
652 rp->SetLastDaughter(-1);
653 rp->SetPolarisation(0,0,0);
654 //Set the position in global coordinate system from the RecPoint
655 AliPHOSGeometry * geom = gime->PHOSGeometry() ;
656 AliPHOSTrackSegment * ts = gime->TrackSegment(rp->GetPHOSTSIndex()) ;
657 AliPHOSEmcRecPoint * erp = gime->EmcRecPoint(ts->GetEmcIndex()) ;
659 geom->GetGlobal(erp, pos) ;
666 //____________________________________________________________________________
667 void AliPHOSPIDv1::Print() const
669 // Print the parameters used for the particle type identification
671 Info("Print", "=============== AliPHOSPIDv1 ================") ;
672 printf("Making PID\n") ;
673 printf(" Pricipal analysis file from 0.5 to 100 %s\n", fFileNamePrincipalPhoton.Data() ) ;
674 printf(" Name of parameters file %s\n", fFileNameParameters.Data() ) ;
675 printf(" Matrix of Parameters: 14x4\n") ;
676 printf(" Energy Calibration 1x3 [3 parametres to calibrate energy: A + B* E + C * E^2]\n") ;
677 printf(" RCPV 2x3 rows x and z, columns function cut parameters\n") ;
678 printf(" TOF 1x3 [High Eff-Low Pur,Medium Eff-Pur, Low Eff-High Pur]\n") ;
679 printf(" PCA 5x4 [5 ellipse parametres and 4 parametres to calculate them: A/Sqrt(E) + B* E + C * E^2 + D]\n") ;
680 Printf(" Pi0 PCA 5x3 [5 ellipse parametres and 3 parametres to calculate them: A + B* E + C * E^2]\n") ;
681 fParameters->Print() ;
686 //____________________________________________________________________________
687 void AliPHOSPIDv1::PrintRecParticles(Option_t * option)
689 // Print table of reconstructed particles
691 AliPHOSGetter *gime = AliPHOSGetter::Instance() ;
693 TClonesArray * recParticles = gime->RecParticles() ;
696 message = "\nevent " ;
697 message += gAlice->GetEvNumber() ;
698 message += " found " ;
699 message += recParticles->GetEntriesFast();
700 message += " RecParticles\n" ;
702 if(strstr(option,"all")) { // printing found TS
703 message += "\n PARTICLE Index \n" ;
706 for (index = 0 ; index < recParticles->GetEntries() ; index++) {
707 AliPHOSRecParticle * rp = (AliPHOSRecParticle * ) recParticles->At(index) ;
709 message += rp->Name().Data() ;
711 message += rp->GetIndexInList() ;
713 message += rp->GetType() ;
716 Info("Print", message.Data() ) ;
719 //____________________________________________________________________________
720 void AliPHOSPIDv1::SetParameters()
722 // PCA : To do the Principal Components Analysis it is necessary
723 // the Principal file, which is opened here
724 fX = new double[7]; // Data for the PCA
725 fPPhoton = new double[7]; // Eigenvalues of the PCA
726 fPPi0 = new double[7]; // Eigenvalues of the Pi0 PCA
728 // Read photon principals from the photon file
730 fFileNamePrincipalPhoton = "$ALICE_ROOT/PHOS/PCA8pa15_0.5-100.root" ;
731 TFile f( fFileNamePrincipalPhoton.Data(), "read" ) ;
732 fPrincipalPhoton = dynamic_cast<TPrincipal*> (f.Get("principal")) ;
735 // Read pi0 principals from the pi0 file
737 fFileNamePrincipalPi0 = "$ALICE_ROOT/PHOS/PCA_pi0_40-120.root" ;
738 TFile fPi0( fFileNamePrincipalPi0.Data(), "read" ) ;
739 fPrincipalPi0 = dynamic_cast<TPrincipal*> (fPi0.Get("principal")) ;
742 // Open parameters file and initialization of the Parameters matrix.
743 // In the File Parameters.dat are all the parameters. These are introduced
744 // in a matrix of 16x4
746 // All the parameters defined in this file are, in order of row:
747 // line 0 : calibration
748 // lines 1,2 : CPV rectangular cat for X and Z
750 // lines 4-8 : parameters to calculate photon PCA ellipse
751 // lines 9-13: parameters to calculate pi0 PCA ellipse
752 // lines 14-15: parameters to calculate border for high-pt photons and pi0
754 fFileNameParameters = gSystem->ExpandPathName("$ALICE_ROOT/PHOS/Parameters.dat");
755 fParameters = new TMatrix(16,4) ;
756 const Int_t maxLeng=255;
757 char string[maxLeng];
759 // Open a text file with PID parameters
760 FILE *fd = fopen(fFileNameParameters.Data(),"r");
762 Fatal("SetParameter","File %s with a PID parameters cannot be opened\n",
763 fFileNameParameters.Data());
766 // Read parameter file line-by-line and skip empty line and comments
767 while (fgets(string,maxLeng,fd) != NULL) {
768 if (string[0] == '\n' ) continue;
769 if (string[0] == '!' ) continue;
770 sscanf(string, "%f %f %f %f",
771 &(*fParameters)(i,0), &(*fParameters)(i,1),
772 &(*fParameters)(i,2), &(*fParameters)(i,3));
774 //printf("line %d: %s",i,string);
779 //____________________________________________________________________________
780 void AliPHOSPIDv1::SetParameterCalibration(Int_t i,Float_t param)
782 // Set parameter "Calibration" i to a value param
784 Error("SetParameterCalibration","Invalid parameter number: %d",i);
786 (*fParameters)(0,i) = param ;
789 //____________________________________________________________________________
790 void AliPHOSPIDv1::SetParameterCpv2Emc(Int_t i, TString axis, Float_t cut)
792 // Set the parameters to calculate Cpv-to-Emc Distance Cut depending on
793 // Purity-Efficiency point i
796 Error("SetParameterCpv2Emc","Invalid parameter number: %d",i);
799 if (axis == "x") (*fParameters)(1,i) = cut;
800 else if (axis == "z") (*fParameters)(2,i) = cut;
801 else Error("SetParameterCpv2Emc","Invalid axis name: %s",axis.Data());
805 //____________________________________________________________________________
806 void AliPHOSPIDv1::SetParameterPhotonBoundary(Int_t i,Float_t param)
808 // Set parameter "Hard photon boundary" i to a value param
810 Error("SetParameterPhotonBoundary","Invalid parameter number: %d",i);
812 (*fParameters)(14,i) = param ;
815 //____________________________________________________________________________
816 void AliPHOSPIDv1::SetParameterPi0Boundary(Int_t i,Float_t param)
818 // Set parameter "Hard pi0 boundary" i to a value param
820 Error("SetParameterPi0Boundary","Invalid parameter number: %d",i);
822 (*fParameters)(15,i) = param ;
825 //_____________________________________________________________________________
826 void AliPHOSPIDv1::SetParameterTimeGate(Int_t i, Float_t gate)
828 // Set the parameter TimeGate depending on Purity-Efficiency point i
830 Error("SetParameterTimeGate","Invalid Efficiency-Purity choice %d",i);
832 (*fParameters)(3,i)= gate ;
835 //_____________________________________________________________________________
836 void AliPHOSPIDv1::SetParameterToCalculateEllipse(TString particle, TString param, Int_t i, Float_t par)
838 // Set the parameter "i" that is needed to calculate the ellipse
839 // parameter "param" for a particle "particle"
846 if (particle == "photon") offset=0;
847 else if (particle == "pi0") offset=5;
849 Error("SetParameterToCalculateEllipse","Wrong particle name: %s (choose from pi0/photon)\n",particle.Data());
851 if (param.Contains("a")) p=4+offset;
852 else if(param.Contains("b")) p=5+offset;
853 else if(param.Contains("c")) p=6+offset;
854 else if(param.Contains("x0"))p=7+offset;
855 else if(param.Contains("y0"))p=8+offset;
857 Error("SetEllipseParameter", "No parameter with index %d", i) ;
859 Error("SetEllipseParameter", "No parameter with name %s", param.Data() ) ;
861 (*fParameters)(p,i) = par ;
864 //____________________________________________________________________________
865 void AliPHOSPIDv1::Unload()
867 AliPHOSGetter * gime = AliPHOSGetter::Instance() ;
868 gime->PhosLoader()->UnloadRecPoints() ;
869 gime->PhosLoader()->UnloadTracks() ;
870 gime->PhosLoader()->UnloadRecParticles() ;
873 //____________________________________________________________________________
874 void AliPHOSPIDv1::WriteRecParticles()
877 AliPHOSGetter *gime = AliPHOSGetter::Instance() ;
879 TClonesArray * recParticles = gime->RecParticles() ;
880 recParticles->Expand(recParticles->GetEntriesFast() ) ;
881 TTree * treeP = gime->TreeP();
884 Int_t bufferSize = 32000 ;
885 TBranch * rpBranch = treeP->Branch("PHOSRP",&recParticles,bufferSize);
886 rpBranch->SetTitle(BranchName());
890 gime->WriteRecParticles("OVERWRITE");
891 gime->WritePID("OVERWRITE");