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 > fCpvEmcDistance (each bit corresponds
27 // to a different efficiency-purity point of the photon identification)
28 // -Bit 3 to 5: bit set if TOF < fTimeGate (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 don 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 thourgh the study of the reconstructed
41 // energy distribution of monoenergetic photons.
43 // All the parameters (RCPV(6 rows-3 columns),TOF(6r-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 (18,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 elipses 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 ---
97 // --- AliRoot header files ---
100 #include "AliGenerator.h"
102 #include "AliPHOSPIDv1.h"
103 #include "AliPHOSClusterizerv1.h"
104 #include "AliPHOSTrackSegment.h"
105 #include "AliPHOSTrackSegmentMakerv1.h"
106 #include "AliPHOSRecParticle.h"
107 #include "AliPHOSGeometry.h"
108 #include "AliPHOSGetter.h"
110 ClassImp( AliPHOSPIDv1)
112 //____________________________________________________________________________
113 AliPHOSPIDv1::AliPHOSPIDv1():AliPHOSPID()
118 fDefaultInit = kTRUE ;
122 //____________________________________________________________________________
123 AliPHOSPIDv1::AliPHOSPIDv1(const char * headerFile,const char * name, const Bool_t toSplit)
124 :AliPHOSPID(headerFile, name,toSplit)
126 //ctor with the indication on where to look for the track segments
131 fDefaultInit = kFALSE ;
135 //____________________________________________________________________________
136 AliPHOSPIDv1::~AliPHOSPIDv1()
139 // fDefaultInit = kTRUE if PID created by default ctor (to get just the parameters)
141 delete [] fX ; // Principal input
142 delete [] fP ; // Principal components
143 // delete fParameters ; // Matrix of Parameters
148 // AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
149 // remove the task from the folder list
150 // gime->RemoveTask("P",GetName()) ;
151 // TString name(GetName()) ;
152 // name.ReplaceAll("pid", "clu") ;
153 // gime->RemoveTask("C",name) ;
155 // // remove the data from the folder list
156 // name = GetName() ;
157 // name.Remove(name.Index(":")) ;
158 // gime->RemoveObjects("RE", name) ; // EMCARecPoints
159 // gime->RemoveObjects("RC", name) ; // CPVRecPoints
160 // gime->RemoveObjects("T", name) ; // TrackSegments
161 // gime->RemoveObjects("P", name) ; // RecParticles
164 // gime->CloseFile() ;
170 //____________________________________________________________________________
171 const TString AliPHOSPIDv1::BranchName() const
173 TString branchName(GetName() ) ;
174 branchName.Remove(branchName.Index(Version())-1) ;
178 //____________________________________________________________________________
179 void AliPHOSPIDv1::Init()
181 // Make all memory allocations that are not possible in default constructor
182 // Add the PID task to the list of PHOS tasks
184 if ( strcmp(GetTitle(), "") == 0 )
185 SetTitle("galice.root") ;
187 TString branchname(GetName()) ;
188 branchname.Remove(branchname.Index(Version())-1) ;
189 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(),branchname.Data(),fToSplit ) ;
191 // gime->SetRecParticlesTitle(BranchName()) ;
193 Error("Init", "Could not obtain the Getter object !" ) ;
199 //First - extract full path if necessary
200 TString fileName(GetTitle()) ;
201 Ssiz_t islash = fileName.Last('/') ;
202 if(islash<fileName.Length())
203 fileName.Remove(islash+1,fileName.Length()) ;
206 fileName+="PHOS.RecData." ;
207 if((strcmp(branchname.Data(),"Default")!=0)&&(strcmp(branchname.Data(),"")!=0)){
208 fileName+=branchname.Data() ;
212 fSplitFile = static_cast<TFile*>(gROOT->GetFile(fileName.Data()));
214 fSplitFile = TFile::Open(fileName.Data(),"update") ;
217 gime->PostPID(this) ;
218 // create a folder on the white board //YSAlice/WhiteBoard/RecParticles/PHOS/recparticlesName
219 gime->PostRecParticles(branchname) ;
223 //____________________________________________________________________________
224 void AliPHOSPIDv1::InitParameters()
227 // fHeaderFileName = GetTitle() ;
228 // TString name(GetName()) ;
229 // if (name.IsNull())
230 // name = "Default" ;
231 // fTrackSegmentsTitle = name ;
232 // fRecPointsTitle = name ;
233 // fRecParticlesTitle = name ;
234 // name.Append(":") ;
235 // name.Append(Version()) ;
237 fRecParticlesInRun = 0 ;
239 // fClusterizer = 0 ;
241 fRecParticlesInRun = 0 ;
242 TString pidName( GetName()) ;
243 if (pidName.IsNull() )
244 pidName = "Default" ;
245 pidName.Append(":") ;
246 pidName.Append(Version()) ;
248 SetParameters() ; // fill the parameters matrix from parameters file
251 //____________________________________________________________________________
252 const Double_t AliPHOSPIDv1::GetCpvtoEmcDistanceCut(const Float_t Cluster_En, const TString Eff_Pur) const
254 // Get CpvtoEmcDistanceCut parameter depending on the cluster energy and
255 // Purity-Efficiency point (possible options "HIGH EFFICIENCY"
256 // "MEDIUM EFFICIENCY" "LOW EFFICIENCY" and 3 more options changing
257 // EFFICIENCY by PURITY)
259 Int_t eff_pur = GetEffPurOption(Eff_Pur);
260 Int_t cluster = GetClusterOption(Cluster_En) ;
261 if((cluster!= -1)&&(eff_pur != -1))
262 return (*fParameters)(cluster,eff_pur) ;
266 //____________________________________________________________________________
268 const Double_t AliPHOSPIDv1::GetTimeGate(const Float_t Cluster_En, const TString Eff_Pur) const
270 // Get TimeGate parameter depending on the cluster energy and
271 // Purity-Efficiency point (possible options "HIGH EFFICIENCY"
272 // "MEDIUM EFFICIENCY" "LOW EFFICIENCY" and 3 more options changing
273 // EFFICIENCY by PURITY)
275 Int_t eff_pur = GetEffPurOption(Eff_Pur);
276 Int_t cluster = GetClusterOption(Cluster_En) ;
277 if((cluster!= -1)&&(eff_pur != -1))
278 return (*fParameters)(cluster+6,eff_pur) ;
283 //_____________________________________________________________________________
284 const Float_t AliPHOSPIDv1::GetDistance(AliPHOSEmcRecPoint * emc,AliPHOSRecPoint * cpv, Option_t * Axis)const
286 // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
288 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
292 emc->GetLocalPosition(vecEmc) ;
293 cpv->GetLocalPosition(vecCpv) ;
294 if(emc->GetPHOSMod() == cpv->GetPHOSMod()){
295 // Correct to difference in CPV and EMC position due to different distance to center.
296 // we assume, that particle moves from center
297 Float_t dCPV = geom->GetIPtoOuterCoverDistance();
298 Float_t dEMC = geom->GetIPtoCrystalSurface() ;
300 vecCpv = dEMC * vecCpv - vecEmc ;
301 if (Axis == "X") return vecCpv.X();
302 if (Axis == "Y") return vecCpv.Y();
303 if (Axis == "Z") return vecCpv.Z();
304 if (Axis == "R") return vecCpv.Mag();
312 //____________________________________________________________________________
313 const Double_t AliPHOSPIDv1::GetCalibratedEnergy(const Float_t e) const
315 // It calibrates Energy depending on the recpoint energy.
316 // The energy of the reconstructed cluster is corrected with
317 // the formula A + B* E + C* E^2, whose parameters where obtained
318 // through the study of the reconstructed energy distribution of
319 // monoenergetic photons.
321 Double_t p[]={0.,0.,0.};
323 for(i=0;i<3;i++) p[i]= (*fParameters)(17,i);
324 Double_t enerec = p[0] + p[1]* e+ p[2] * e * e;
328 //____________________________________________________________________________
329 const Int_t AliPHOSPIDv1::GetPrincipalSign(const Double_t* P,const Int_t eff_pur, const Float_t E)const
331 //Is the particle inside de PCA ellipse?
334 Double_t A = GetEllipseParameter("a", E);
335 Double_t B = GetEllipseParameter("b", E);
336 Double_t C = GetEllipseParameter("c", E);
337 Double_t X_center = GetEllipseParameter("x0", E);
338 Double_t Y_center = GetEllipseParameter("y0", E);
340 Double_t R = TMath::Power((P[0] - X_center)/A,2) +
341 TMath::Power((P[1] - Y_center)/B,2) +
342 C*(P[0] - X_center)*(P[1] - Y_center)/(A*B) ;
343 //3 different ellipses defined
344 if((eff_pur==2)&&(R <1./2.)) prinsign= 1;
345 if((eff_pur==1)&&(R <2. )) prinsign= 1;
346 if((eff_pur==0)&&(R <9./2.)) prinsign= 1;
349 Error("GetPrincipalSign", "Negative square?") ;
354 //_____________________________________________________________________________
355 void AliPHOSPIDv1::SetCpvtoEmcDistanceCut(Float_t Cluster_En, TString Eff_Pur, Float_t cut)
358 // Set the parameter Cpvto EmcDistanceCut depending on the cluster energy and
359 // Purity-Efficiency point (possible options "HIGH EFFICIENCY"
360 // "MEDIUM EFFICIENCY" "LOW EFFICIENCY" and 3 more options changing
361 // EFFICIENCY by PURITY)
364 Int_t eff_pur = GetEffPurOption(Eff_Pur);
365 Int_t cluster = GetClusterOption(Cluster_En) ;
366 if((cluster!= -1)&&(eff_pur != -1))
367 (*fParameters)(cluster,eff_pur) = cut ;
369 //_____________________________________________________________________________
370 void AliPHOSPIDv1::SetTimeGate(Float_t Cluster_En, TString Eff_Pur, Float_t gate)
373 // Set the parameter TimeGate depending on the cluster energy and
374 // Purity-Efficiency point (possible options "HIGH EFFICIENCY"
375 // "MEDIUM EFFICIENCY" "LOW EFFICIENCY" and 3 more options changing
376 // EFFICIENCY by PURITY)
378 Int_t eff_pur = GetEffPurOption(Eff_Pur);
379 Int_t cluster = GetClusterOption(Cluster_En) ;
380 if((cluster!= -1)&&(eff_pur != -1))
381 (*fParameters)(cluster+6,eff_pur) = gate ;
383 //_____________________________________________________________________________
384 void AliPHOSPIDv1::SetParameters()
386 // PCA : To do the Principal Components Analysis it is necessary
387 // the Principal file, which is opened here
388 fX = new double[7]; // Data for the PCA
389 fP = new double[7]; // Eigenvalues of the PCA
391 // Open principal and parameters files to be used
393 fFileName = "$ALICE_ROOT/PHOS/PCA8pa15_0.5-100.root" ;
394 fFileNamePar = gSystem->ExpandPathName("$ALICE_ROOT/PHOS/Parameters.dat");
395 TFile f( fFileName.Data(), "read" ) ;
396 fPrincipal = dynamic_cast<TPrincipal*> (f.Get("principal")) ;
399 // Initialization of the Parameters matrix. In the File Parameters.dat
400 // are all the parameters. These are introduced in a matrix of 18x4
402 // All the parameters defined in this file are, in order of row:
403 // CpvtoEmcDistanceCut (6 rows, each one depends on the energy range of the
404 // particle, and 3 columns, each one depending on the efficiency-purity
405 // point), TimeGate (the same) and the parameters of the functions that
406 // calculate the ellipse parameters, x0,y0,a, b, c. These 5 parameters
407 // (5 rows) depend on 4 parameters (columns). Finally there is a row with
408 // the energy calibration parameters, 3 parameters.
410 fParameters = new TMatrixD(18,4) ;
412 ifstream paramFile(fFileNamePar) ;
414 for(h = 0; h< 18; h++){
415 for(n = 0; n< 4; n++){
416 paramFile >> (*fParameters)(h,n) ;
421 //_____________________________________________________________________________
422 const Int_t AliPHOSPIDv1::GetClusterOption(const Float_t Cluster_En) const
424 // Gives the cluster energy range, for each range there is associated a TOF or RCPV
427 if((Cluster_En > 0.0 )&&(Cluster_En <= 2.0 )) cluster = 0 ;
428 if((Cluster_En > 2.0 )&&(Cluster_En <= 5.0 )) cluster = 1 ;
429 if((Cluster_En > 5.0 )&&(Cluster_En <= 10.0)) cluster = 2 ;
430 if((Cluster_En > 10.0)&&(Cluster_En <= 20.0)) cluster = 3 ;
431 if((Cluster_En > 20.0)&&(Cluster_En <= 30.0)) cluster = 4 ;
432 if( Cluster_En > 30.0) cluster = 5 ;
436 //____________________________________________________________________________
437 const Int_t AliPHOSPIDv1::GetEffPurOption(const TString Eff_Pur) const
440 // Looks for the Purity-Efficiency point (possible options "HIGH EFFICIENCY"
441 // "MEDIUM EFFICIENCY" "LOW EFFICIENCY" and 3 more options changing
442 // EFFICIENCY by PURITY)
446 if(Eff_Pur.Contains("HIGH EFFICIENCY") ||Eff_Pur.Contains("LOW PURITY") )
448 else if(Eff_Pur.Contains("MEDIUM EFFICIENCY") ||Eff_Pur.Contains("MEDIUM PURITY") )
450 else if(Eff_Pur.Contains("LOW EFFICIENCY")||Eff_Pur.Contains("HIGH PURITY") )
455 message = "Invalid Efficiency-Purity option\n";
456 message += "Possible options: HIGH EFFICIENCY = LOW PURITY\n" ;
457 message += " MEDIUM EFFICIENCY = MEDIUM PURITY\n" ;
458 message += " LOW EFFICIENCY = HIGH PURITY\n" ;
463 //________________________________________________________________________
464 void AliPHOSPIDv1::SetEllipseParameter(TString Param, Int_t i, Double_t par)
466 // Set the parameter "i" that is needed to calculate the ellipse
467 // parameter "Param".
471 if(Param.Contains("a"))p=12;
472 if(Param.Contains("b"))p=13;
473 if(Param.Contains("c"))p=14;
474 if(Param.Contains("x0"))p=15;
475 if(Param.Contains("y0"))p=16;
477 Error("SetEllipseParameter", "No parameter with index %d", i) ;
479 Error("SetEllipseParameter", "No parameter with name %s", Param.Data() ) ;
481 (*fParameters)(p,i) = par ;
483 //________________________________________________________________________
484 const Double_t AliPHOSPIDv1::GetParameterToCalculateEllipse(const TString Param, const Int_t i) const
486 // Get the parameter "i" that is needed to calculate the ellipse
487 // parameter "Param".
492 if(Param.Contains("a"))p=12;
493 if(Param.Contains("b"))p=13;
494 if(Param.Contains("c"))p=14;
495 if(Param.Contains("x0"))p=15;
496 if(Param.Contains("y0"))p=16;
499 Error("GetParameterToCalculateEllipse", "No parameter with index", i) ;
501 Error("GetParameterToCalculateEllipse", "No parameter with name %s", Param.Data() ) ;
503 par = (*fParameters)(p,i) ;
508 //____________________________________________________________________________
509 void AliPHOSPIDv1::SetCalibrationParameter(Int_t i,Double_t param)
511 (*fParameters)(17,i) = param ;
513 //____________________________________________________________________________
514 const Double_t AliPHOSPIDv1::GetCalibrationParameter(const Int_t i) const
516 Float_t param = (*fParameters)(17,i);
519 //____________________________________________________________________________
520 const Double_t AliPHOSPIDv1::GetEllipseParameter(const TString Param,Float_t E) const
522 Double_t p[4]={0.,0.,0.,0.};
523 Double_t value = 0.0;
526 if(Param.Contains("a")){
527 for(i=0;i<4;i++)p[i]=(*fParameters)(12,i);
531 if(Param.Contains("b")){
532 for(i=0;i<4;i++)p[i]=(*fParameters)(13,i);
536 if(Param.Contains("c"))
537 for(i=0;i<4;i++)p[i]=(*fParameters)(14,i);
539 if(Param.Contains("x0")){
540 for(i=0;i<4;i++)p[i]=(*fParameters)(15,i);
543 if(Param.Contains("y0"))
544 for(i=0;i<4;i++)p[i]=(*fParameters)(16,i);
546 value = p[0]/TMath::Sqrt(E)+p[1]*E+p[2]*E*E+p[3];
549 //____________________________________________________________________________
551 void AliPHOSPIDv1::Exec(Option_t * option)
555 if( strcmp(GetName(), "")== 0 )
558 if(strstr(option,"tim"))
559 gBenchmark->Start("PHOSPID");
561 if(strstr(option,"print")) {
567 // gAlice->GetEvent(0) ;
569 // //check, if the branch with name of this" already exits?
570 // if (gAlice->TreeR()) {
571 // TObjArray * lob = (TObjArray*)gAlice->TreeR()->GetListOfBranches() ;
573 // TBranch * branch = 0 ;
574 // Bool_t phospidfound = kFALSE, pidfound = kFALSE ;
576 // TString taskName(GetName()) ;
577 // taskName.Remove(taskName.Index(Version())-1) ;
579 // while ( (branch = (TBranch*)next()) && (!phospidfound || !pidfound) ) {
580 // if ( (strcmp(branch->GetName(), "PHOSPID")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) )
581 // phospidfound = kTRUE ;
583 // else if ( (strcmp(branch->GetName(), "AliPHOSPID")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) )
584 // pidfound = kTRUE ;
587 // if ( phospidfound || pidfound ) {
588 // Error("Exec", "RecParticles and/or PIDtMaker branch with name %s already exists", taskName.Data() ) ;
593 // Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
595 // AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
597 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
598 if(gime->BranchExists("RecParticles") )
600 Int_t nevents = gime->MaxEvent() ; //(Int_t) gAlice->TreeE()->GetEntries() ;
604 for(ievent = 0; ievent < nevents; ievent++){
605 gime->Event(ievent,"R") ;
609 WriteRecParticles(ievent);
611 if(strstr(option,"deb"))
612 PrintRecParticles(option) ;
614 //increment the total number of rec particles per run
615 fRecParticlesInRun += gime->RecParticles(BranchName())->GetEntriesFast() ;
619 if(strstr(option,"tim")){
620 gBenchmark->Stop("PHOSPID");
621 Info("Exec", "took %f seconds for PID %f seconds per event",
622 gBenchmark->GetCpuTime("PHOSPID"),
623 gBenchmark->GetCpuTime("PHOSPID")/nevents) ;
627 //____________________________________________________________________________
628 void AliPHOSPIDv1::MakeRecParticles(){
630 // Makes a RecParticle out of a TrackSegment
632 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
633 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
634 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
635 TClonesArray * trackSegments = gime->TrackSegments() ;
636 if ( !emcRecPoints || !cpvRecPoints || !trackSegments ) {
637 Error("MakeRecParticles", "RecPoints or TrackSegments not found !") ;
640 TClonesArray * recParticles = gime->RecParticles() ;
641 recParticles->Clear();
644 TIter next(trackSegments) ;
645 AliPHOSTrackSegment * ts ;
647 AliPHOSRecParticle * rp ;
648 while ( (ts = (AliPHOSTrackSegment *)next()) ) {
650 new( (*recParticles)[index] ) AliPHOSRecParticle() ;
651 rp = (AliPHOSRecParticle *)recParticles->At(index) ;
652 rp->SetTrackSegment(index) ;
653 rp->SetIndexInList(index) ;
655 AliPHOSEmcRecPoint * emc = 0 ;
656 if(ts->GetEmcIndex()>=0)
657 emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ;
659 AliPHOSRecPoint * cpv = 0 ;
660 if(ts->GetCpvIndex()>=0)
661 cpv = (AliPHOSRecPoint *) cpvRecPoints->At(ts->GetCpvIndex()) ;
663 // Now set type (reconstructed) of the particle
665 // Choose the cluster energy range
667 // YK: check if (emc != 0) !!!
669 Error("MakeRecParticles", "-> emc(%d) = %d", ts->GetEmcIndex(), emc ) ;
672 Float_t e = emc->GetEnergy() ;
673 Int_t cluster = GetClusterOption(e) ;// Gives value to cluster that defines the energy range parameter to be used in de RCPV, TOF and used in the PCA.
674 if(cluster== -1) continue ;
677 emc->GetElipsAxis(lambda) ;
679 if((lambda[0]>0.01) && (lambda[1]>0.01)){
680 // Looking PCA. Define and calculate the data (X),
681 // introduce in the function
682 // X2P that gives the components (P).
684 Float_t Emaxdtotal = 0. ;
686 if((lambda[0]+lambda[1])!=0) Spher=fabs(lambda[0]-lambda[1])/(lambda[0]+lambda[1]);
688 Emaxdtotal=emc->GetMaximalEnergy()/emc->GetEnergy();
692 fX[2] = emc->GetDispersion() ;
694 fX[4] = emc->GetMultiplicity() ;
696 fX[6] = emc->GetCoreEnergy() ;
698 fPrincipal->X2P(fX,fP);
701 fP[0]=-100.0; //We do not accept clusters with
702 fP[1]=-100.0; //one cell as a photon-like
705 Float_t time =emc->GetTime() ;
707 // Loop of Efficiency-Purity (the 3 points of purity or efficiency are taken
708 // into account to set the particle identification)
709 for(Int_t eff_pur = 0; eff_pur < 3 ; eff_pur++){
711 // Looking at the CPV detector. If RCPV greater than CpvEmcDistance, 1st,
712 // 2nd or 3rd bit (depending on the efficiency-purity point )is set to 1 .
714 if(GetDistance(emc, cpv, "R") > (*fParameters)(cluster,eff_pur) )
715 rp->SetPIDBit(eff_pur) ;
717 // Looking the TOF. If TOF smaller than gate, 4th, 5th or 6th
718 // bit (depending on the efficiency-purity point )is set to 1
719 if(time< (*fParameters)(cluster+6,eff_pur)) {
720 rp->SetPIDBit(eff_pur+3) ;
723 //If we are inside the ellipse, 7th, 8th or 9th
724 // bit (depending on the efficiency-purity point )is set to 1
725 if(GetPrincipalSign(fP,eff_pur,e) == 1)
726 rp->SetPIDBit(eff_pur+6) ;
729 //Set momentum, energy and other parameters
730 Float_t encal = GetCalibratedEnergy(e);
731 TVector3 dir = GetMomentumDirection(emc,cpv) ;
733 rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),encal) ;
735 rp->Name(); //If photon sets the particle pdg name to gamma
736 rp->SetProductionVertex(0,0,0,0);
737 rp->SetFirstMother(-1);
738 rp->SetLastMother(-1);
739 rp->SetFirstDaughter(-1);
740 rp->SetLastDaughter(-1);
741 rp->SetPolarisation(0,0,0);
747 //____________________________________________________________________________
748 void AliPHOSPIDv1:: Print()
750 // Print the parameters used for the particle type identification
752 message = "\n=============== AliPHOSPID1 ================\n" ;
753 message += "Making PID\n";
754 message += " Pricipal analysis file from 0.5 to 100 %s\n" ;
755 message += " Name of parameters file %s\n" ;
756 message += " Matrix of Parameters: 18x4\n" ;
757 message += " RCPV 6x3 [High Eff-Low Pur,Medium Eff-Pur, Low Eff-High Pur]\n" ;
758 message += " TOF 6x3 [High Eff-Low Pur,Medium Eff-Pur, Low Eff-High Pur]\n" ;
759 message += " PCA 5x4 [5 ellipse parametres and 4 parametres to calculate them: A/Sqrt(E) + B* E + C * E^2 + D]\n" ;
760 message += " Energy Calibration 1x3 [3 parametres to calibrate energy: A + B* E + C * E^2]\n" ;
761 Info("Print", message.Data(), fFileName.Data(), fFileNamePar.Data() ) ;
762 fParameters->Print() ;
765 //____________________________________________________________________________
766 void AliPHOSPIDv1::WriteRecParticles(Int_t event)
769 AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ;
771 TClonesArray * recParticles = gime->RecParticles() ;
772 recParticles->Expand(recParticles->GetEntriesFast() ) ;
780 sprintf(name,"%s%d", "TreeR",event) ;
781 treeR = dynamic_cast<TTree*>(fSplitFile->Get(name));
784 treeR = gAlice->TreeR();
788 gAlice->MakeTree("R", fSplitFile);
789 treeR = gAlice->TreeR() ;
793 Int_t bufferSize = 32000 ;
794 TBranch * rpBranch = treeR->Branch("PHOSRP",&recParticles,bufferSize);
795 rpBranch->SetTitle(BranchName());
799 Int_t splitlevel = 0 ;
800 AliPHOSPIDv1 * pid = this ;
801 TBranch * pidBranch = treeR->Branch("AliPHOSPID","AliPHOSPIDv1",&pid,bufferSize,splitlevel);
802 pidBranch->SetTitle(BranchName());
807 treeR->AutoSave() ; //Write(0,kOverwrite) ;
808 if(gAlice->TreeR()!=treeR){
813 //____________________________________________________________________________
814 TVector3 AliPHOSPIDv1::GetMomentumDirection(AliPHOSEmcRecPoint * emc, AliPHOSRecPoint * cpv)const
816 // Calculates the momentum direction:
817 // 1. if only a EMC RecPoint, direction is given by IP and this RecPoint
818 // 2. if a EMC RecPoint and CPV RecPoint, direction is given by the line through the 2 recpoints
819 // However because of the poor position resolution of PPSD the direction is always taken as if we were
822 TVector3 dir(0,0,0) ;
824 TVector3 emcglobalpos ;
827 emc->GetGlobalPosition(emcglobalpos, dummy) ;
831 dir.SetZ( -dir.Z() ) ; // why ?
834 //account correction to the position of IP
835 Float_t xo,yo,zo ; //Coordinates of the origin
836 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
837 TVector3 origin(xo,yo,zo);
842 //____________________________________________________________________________
843 void AliPHOSPIDv1::PrintRecParticles(Option_t * option)
845 // Print table of reconstructed particles
847 AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ;
849 TClonesArray * recParticles = gime->RecParticles(BranchName()) ;
852 message = "\nevent " ;
853 message += gAlice->GetEvNumber() ;
854 message += " found " ;
855 message += recParticles->GetEntriesFast();
856 message += " RecParticles\n" ;
858 if(strstr(option,"all")) { // printing found TS
859 message += "\n PARTICLE Index \n" ;
862 for (index = 0 ; index < recParticles->GetEntries() ; index++) {
863 AliPHOSRecParticle * rp = (AliPHOSRecParticle * ) recParticles->At(index) ;
865 message += rp->Name().Data() ;
867 message += rp->GetIndexInList() ;
869 message += rp->GetType() ;
872 Info("Print", message.Data() ) ;