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 v0 of the PHOS particle identifier
20 // Particle identification based on the
22 // - Preshower information (in MIXT or GPS2 geometries)
25 // CPV or Preshower clusters should be closer in PHOS plane than fCpvEmcDistance (in cm).
26 // This parameter can be set by method SetCpvtoEmcDistanceCut(Float_t cut)
28 // One can set desirable ID method by the function SetIdentificationMethod(option).
29 // Presently the following options can be used together or separately :
30 // - "disp": use dispersion cut on shower width
31 // (width can be set by method SetDispersionCut(Float_t cut)
32 // - "ell" : use cut on the axis of the ellipse, drawn around shower
33 // (this cut can be changed by SetShowerProfileCut(char* formula),
34 // where formula - any function of two variables f(lambda[0],lambda[1]).
35 // Shower is considered as EM if f() > 0 )
36 // One can visualize current cuts calling method PlotDispersionCuts().
39 // root [0] AliPHOSPIDv0 * p1 = new AliPHOSPIDv0("galice.root")
40 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
41 // root [1] p1->SetIdentificationMethod("disp ellipse")
42 // root [2] p1->ExecuteTask()
43 // root [3] AliPHOSPIDv0 * p2 = new AliPHOSPIDv0("galice1.root","ts1")
44 // Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
45 // // reading headers from file galice1.root and TrackSegments
46 // // with title "ts1"
47 // root [4] p2->SetRecParticlesBranch("rp1")
48 // // set file name for the branch RecParticles
49 // root [5] p2->ExecuteTask("deb all time")
50 // // available options
51 // // "deb" - prints # of reconstructed particles
52 // // "deb all" - prints # and list of RecParticles
53 // // "time" - prints benchmarking results
55 //*-- Author: Yves Schutz (SUBATECH) & Gines Martinez (SUBATECH) &
56 // Dmitri Peressounko (SUBATECH & Kurchatov Institute)
57 // Completely redesined by Dmitri Peressounko, March 2001
59 // --- ROOT system ---
68 #include "TBenchmark.h"
69 // --- Standard library ---
71 // --- AliRoot header files ---
74 #include "AliGenerator.h"
76 #include "AliPHOSPIDv0.h"
77 #include "AliPHOSEmcRecPoint.h"
78 #include "AliPHOSClusterizerv1.h"
79 #include "AliPHOSTrackSegment.h"
80 #include "AliPHOSTrackSegmentMakerv1.h"
81 #include "AliPHOSRecParticle.h"
82 #include "AliPHOSGeometry.h"
83 #include "AliPHOSLoader.h"
85 ClassImp( AliPHOSPIDv0)
87 //____________________________________________________________________________
88 AliPHOSPIDv0::AliPHOSPIDv0():AliPHOSPID()
95 fEventFolderName = "" ;
96 fTrackSegmentsTitle= "" ;
97 fRecPointsTitle = "" ;
98 fRecParticlesTitle = "" ;
99 fIDOptions = "dis time" ;
100 fRecParticlesInRun = 0 ;
105 //____________________________________________________________________________
106 AliPHOSPIDv0::AliPHOSPIDv0(const char * evFolderName,const char * name) : AliPHOSPID(evFolderName, name)
108 //ctor with the indication on where to look for the track segments
110 fFormula = new TFormula("LambdaCuts","(x>1)*(x<2.5)*(y>0)*(y<x)") ;
112 fCpvEmcDistance = 3.0 ;
115 fEventFolderName = GetTitle() ;
116 fTrackSegmentsTitle = GetName();
117 fRecPointsTitle = GetName();
118 fRecParticlesTitle = GetName();
119 fIDOptions = "dis time" ;
121 fRecParticlesInRun = 0 ;
127 //____________________________________________________________________________
128 AliPHOSPIDv0::~AliPHOSPIDv0()
132 //____________________________________________________________________________
133 Float_t AliPHOSPIDv0::GetDistance(AliPHOSEmcRecPoint * emc,AliPHOSRecPoint * cpv, Option_t * Axis)const
135 // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
137 const AliPHOSGeometry * geom = AliPHOSLoader::GetPHOSGeometry() ;
141 emc->GetLocalPosition(vecEmc) ;
142 cpv->GetLocalPosition(vecCpv) ;
143 if(emc->GetPHOSMod() == cpv->GetPHOSMod()){
145 // Correct to difference in CPV and EMC position due to different distance to center.
146 // we assume, that particle moves from center
147 Float_t dCPV = geom->GetIPtoOuterCoverDistance();
148 Float_t dEMC = geom->GetIPtoCrystalSurface() ;
150 vecCpv = dEMC * vecCpv - vecEmc ;
151 if (Axis == "X") return vecCpv.X();
152 if (Axis == "Y") return vecCpv.Y();
153 if (Axis == "Z") return vecCpv.Z();
154 if (Axis == "R") return vecCpv.Mag();
160 //____________________________________________________________________________
161 void AliPHOSPIDv0::Exec(Option_t * option)
165 if( strcmp(GetName(), "")== 0 )
168 if(strstr(option,"tim"))
169 gBenchmark->Start("PHOSPID");
171 if(strstr(option,"print")) {
176 AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
179 Error("Exec","Can not find run getter in event folder \"%s\"",GetTitle());
183 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
186 Error("Exec","Could not obtain the Loader object !");
190 if(gime->BranchExists("RecParticles") )
194 Int_t nevents = runget->GetNumberOfEvents() ; //(Int_t) gAlice->TreeE()->GetEntries() ;
198 for(ievent = 0; ievent < nevents; ievent++){
199 runget->GetEvent(ievent);
200 Info("Exec", "event %d %d %d", ievent, gime->EmcRecPoints(), gime->TrackSegments()) ;
205 if(strstr(option,"deb"))
206 PrintRecParticles(option) ;
208 //increment the total number of rec particles per run
209 fRecParticlesInRun += gime->RecParticles()->GetEntriesFast() ;
213 if(strstr(option,"tim")){
214 gBenchmark->Stop("PHOSPID");
215 Info("Exec", "took %f seconds for PID %f seconds per event",
216 gBenchmark->GetCpuTime("PHOSPID"), gBenchmark->GetCpuTime("PHOSPID")/nevents) ;
220 //____________________________________________________________________________
221 void AliPHOSPIDv0::Init()
223 // Make all memory allocations that are not possible in default constructor
224 // Add the PID task to the list of PHOS tasks
226 AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
229 Error("Exec","Can not find run getter in event folder \"%s\"",GetTitle());
233 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
236 Error("Exec","Could not obtain the Loader object !");
241 gime->LoadRecParticles("UPDATE");
245 //____________________________________________________________________________
246 void AliPHOSPIDv0::MakeRecParticles(){
248 TString taskName(GetName()) ;
249 taskName.Remove(taskName.Index(Version())-1) ;
251 AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
254 Error("Exec","Can not find run getter in event folder \"%s\"",GetTitle());
258 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
261 Error("Exec","Could not obtain the Loader object !");
265 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
266 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
267 TClonesArray * trackSegments = gime->TrackSegments() ;
268 TClonesArray * recParticles = gime->RecParticles() ;
270 recParticles->Clear();
272 TIter next(trackSegments) ;
273 AliPHOSTrackSegment * ts ;
275 AliPHOSRecParticle * rp ;
277 Bool_t ellips = fIDOptions.Contains("ell",TString::kIgnoreCase ) ;
278 Bool_t disp = fIDOptions.Contains("dis",TString::kIgnoreCase ) ;
279 Bool_t time = fIDOptions.Contains("tim",TString::kIgnoreCase ) ;
281 while ( (ts = (AliPHOSTrackSegment *)next()) ) {
283 new( (*recParticles)[index] ) AliPHOSRecParticle() ;
284 rp = (AliPHOSRecParticle *)recParticles->At(index) ;
285 rp->SetTrackSegment(index) ;
286 rp->SetIndexInList(index) ;
288 AliPHOSEmcRecPoint * emc = 0 ;
289 if(ts->GetEmcIndex()>=0)
290 emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ;
292 AliPHOSRecPoint * cpv = 0 ;
293 if(ts->GetCpvIndex()>=0)
294 cpv = (AliPHOSRecPoint *) cpvRecPoints->At(ts->GetCpvIndex()) ;
296 //set momentum and energy first
297 Float_t e = emc->GetEnergy() ;
298 TVector3 dir = GetMomentumDirection(emc,cpv) ;
301 rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),e) ;
304 //now set type (reconstructed) of the particle
305 Int_t showerprofile = 0; // 0 narrow and 1 wide
309 emc->GetElipsAxis(lambda) ;
310 if(fFormula->Eval(lambda[0],lambda[1]) <= 0 )
311 showerprofile = 1 ; // not narrow
315 if(emc->GetDispersion() > fDispersion )
316 showerprofile = 1 ; // not narrow
320 if(emc->GetTime() > fTimeGate )
323 // Looking at the CPV detector
324 Int_t cpvdetector= 0 ; //1 hit and 0 no hit
326 if(GetDistance(emc, cpv, "R") < fCpvEmcDistance)
329 Int_t type = showerprofile + 2 * slow + 4 * cpvdetector ;
331 rp->SetProductionVertex(0,0,0,0);
332 rp->SetFirstMother(-1);
333 rp->SetLastMother(-1);
334 rp->SetFirstDaughter(-1);
335 rp->SetLastDaughter(-1);
336 rp->SetPolarisation(0,0,0);
342 //____________________________________________________________________________
343 void AliPHOSPIDv0:: Print() const
345 // Print the parameters used for the particle type identification
347 message = "=============== AliPHOSPIDv0 ================\n" ;
348 message += "Making PID\n" ;
349 message += " Headers file: %s\n" ;
350 message += " RecPoints branch title: %s\n" ;
351 message += " TrackSegments Branch title: %s\n" ;
352 message += " RecParticles Branch title %s\n" ;
353 message += "with parameters:\n" ;
354 message += " Maximal EMC - CPV distance (cm) %f\n" ;
355 Info("Print", message.Data(),
357 fRecPointsTitle.Data(),
358 fTrackSegmentsTitle.Data(),
359 fRecParticlesTitle.Data(),
362 if(fIDOptions.Contains("dis",TString::kIgnoreCase ))
363 Info("Print", " dispersion cut %f", fDispersion ) ;
364 if(fIDOptions.Contains("ell",TString::kIgnoreCase ))
365 Info("Print", " Eliptic cuts function: %s", fFormula->GetTitle() ) ;
366 if(fIDOptions.Contains("tim",TString::kIgnoreCase ))
367 Info("Print", " Time Gate used: %f", fTimeGate) ;
370 //____________________________________________________________________________
371 void AliPHOSPIDv0::SetShowerProfileCut(const char * formula)
373 //set shape of the cut on the axis of ellipce, drown around shouer
374 //shower considered "narrow" if Formula(lambda[0],lambda[1]) > 0.
377 fFormula = new TFormula("Lambda Cut",formula) ;
379 //____________________________________________________________________________
380 void AliPHOSPIDv0::WriteRecParticles()
383 AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
386 Error("Exec","Can not find run getter in event folder \"%s\"",GetTitle());
390 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
393 Error("Exec","Could not obtain the Loader object !");
397 TClonesArray * recParticles = gime->RecParticles() ;
398 recParticles->Expand(recParticles->GetEntriesFast() ) ;
400 TTree * treeR = gime->TreeR();
404 treeR = gime->TreeR() ;
408 Int_t bufferSize = 32000 ;
409 TBranch * rpBranch = treeR->Branch("PHOSRP",&recParticles,bufferSize);
410 rpBranch->SetTitle(fRecParticlesTitle);
413 Int_t splitlevel = 0 ;
414 AliPHOSPIDv0 * pid = this ;
415 TBranch * pidBranch = treeR->Branch("AliPHOSPID","AliPHOSPIDv0",&pid,bufferSize,splitlevel);
416 pidBranch->SetTitle(fRecParticlesTitle.Data());
421 gime->WriteRecParticles("OVERWRITE");
422 gime->WritePID("OVERWRITE");
426 //____________________________________________________________________________
427 void AliPHOSPIDv0::PlotDispersionCuts()const
429 // produces a plot of the dispersion cut
430 TCanvas* lambdas = new TCanvas("lambdas","Cuts on the ellipse axis",200,10,700,500);
432 if(fIDOptions.Contains("ell",TString::kIgnoreCase ) ){
433 TF2 * ell = new TF2("Elliptic Cuts",fFormula->GetName(),0,3,0,3) ;
434 ell->SetMinimum(0.0000001) ;
435 ell->SetMaximum(0.001) ;
436 ell->SetLineStyle(1) ;
437 ell->SetLineWidth(2) ;
441 if( fIDOptions.Contains("dis",TString::kIgnoreCase ) ){
442 TF2 * dsp = new TF2("dispersion","(y<x)*(x*x+y*y < [0]*[0])",0,3,0,3) ;
443 dsp->SetParameter(0,fDispersion) ;
444 dsp->SetMinimum(0.0000001) ;
445 dsp->SetMaximum(0.001) ;
446 dsp->SetLineStyle(1) ;
447 dsp->SetLineColor(2) ;
448 dsp->SetLineWidth(2) ;
451 if(fIDOptions.Contains("ell",TString::kIgnoreCase ) )
459 //____________________________________________________________________________
460 TVector3 AliPHOSPIDv0::GetMomentumDirection(AliPHOSEmcRecPoint * emc, AliPHOSRecPoint * )const
462 // Calculates the momentum direction:
463 // 1. if only a EMC RecPoint, direction is given by IP and this RecPoint
464 // 2. if a EMC RecPoint and CPV RecPoint, direction is given by the line through the 2 recpoints
465 // However because of the poor position resolution of PPSD the direction is always taken as if we were
468 TVector3 dir(0,0,0) ;
470 TVector3 emcglobalpos ;
473 emc->GetGlobalPosition(emcglobalpos, dummy) ;
476 // The following commented code becomes valid once the PPSD provides
477 // a reasonable position resolution, at least as good as EMC !
478 // TVector3 ppsdlglobalpos ;
479 // TVector3 ppsduglobalpos ;
480 // if( fPpsdLowRecPoint ){ // certainly a photon that has concerted
481 // fPpsdLowRecPoint->GetGlobalPosition(ppsdlglobalpos, mdummy) ;
482 // dir = emcglobalpos - ppsdlglobalpos ;
483 // if( fPpsdUpRecPoint ){ // not looks like a charged
484 // fPpsdUpRecPoint->GetGlobalPosition(ppsduglobalpos, mdummy) ;
485 // dir = ( dir + emcglobalpos - ppsduglobalpos ) * 0.5 ;
488 // else { // looks like a neutral
489 // dir = emcglobalpos ;
493 dir.SetZ( -dir.Z() ) ; // why ?
496 //account correction to the position of IP
497 Float_t xo,yo,zo ; //Coordinates of the origin
498 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
499 TVector3 origin(xo,yo,zo);
504 //____________________________________________________________________________
505 void AliPHOSPIDv0::PrintRecParticles(Option_t * option)
507 // Print table of reconstructed particles
509 AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
512 Error("WriteRecParticles","Can not find run getter in event folder \"%s\"",GetTitle());
516 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
519 Error("WriteRecParticles","Could not obtain the Loader object !");
523 TString taskName(GetName()) ;
524 taskName.Remove(taskName.Index(Version())-1) ;
525 TClonesArray * recParticles = gime->RecParticles() ;
528 message = "event %d\n" ;
529 message += " found %d RecParticles\n" ;
530 Info("PrintRecParticles", message.Data(), gAlice->GetEvNumber(), recParticles->GetEntriesFast() ) ;
532 if(strstr(option,"all")) { // printing found TS
533 Info("PrintRecParticles"," PARTICLE Index \n" ) ;
536 for (index = 0 ; index < recParticles->GetEntries() ; index++) {
537 AliPHOSRecParticle * rp = (AliPHOSRecParticle * ) recParticles->At(index) ;
540 switch(rp->GetType()) {
541 case AliPHOSFastRecParticle::kNEUTRALEMFAST:
542 strcpy( particle, "NEUTRAL EM FAST");
544 case AliPHOSFastRecParticle::kNEUTRALHAFAST:
545 strcpy(particle, "NEUTRAL HA FAST");
547 case AliPHOSFastRecParticle::kNEUTRALEMSLOW:
548 strcpy(particle, "NEUTRAL EM SLOW");
550 case AliPHOSFastRecParticle::kNEUTRALHASLOW:
551 strcpy(particle, "NEUTRAL HA SLOW");
553 case AliPHOSFastRecParticle::kCHARGEDEMFAST:
554 strcpy(particle, "CHARGED EM FAST") ;
556 case AliPHOSFastRecParticle::kCHARGEDHAFAST:
557 strcpy(particle, "CHARGED HA FAST") ;
559 case AliPHOSFastRecParticle::kCHARGEDEMSLOW:
560 strcpy(particle, "CHARGEDEMSLOW") ;
562 case AliPHOSFastRecParticle::kCHARGEDHASLOW:
563 strcpy(particle, "CHARGED HA SLOW") ;
567 // Int_t * primaries;
569 // primaries = rp->GetPrimaries(nprimaries);
571 Info("PrintRecParticles", " %s %d\n", particle, rp->GetIndexInList()) ;