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 "AliPHOSClusterizerv1.h"
78 #include "AliPHOSTrackSegment.h"
79 #include "AliPHOSTrackSegmentMakerv1.h"
80 #include "AliPHOSRecParticle.h"
81 #include "AliPHOSGeometry.h"
82 #include "AliPHOSLoader.h"
84 ClassImp( AliPHOSPIDv0)
86 //____________________________________________________________________________
87 AliPHOSPIDv0::AliPHOSPIDv0():AliPHOSPID()
94 fEventFolderName = "" ;
95 fTrackSegmentsTitle= "" ;
96 fRecPointsTitle = "" ;
97 fRecParticlesTitle = "" ;
98 fIDOptions = "dis time" ;
99 fRecParticlesInRun = 0 ;
104 //____________________________________________________________________________
105 AliPHOSPIDv0::AliPHOSPIDv0(const char * evFolderName,const char * name) : AliPHOSPID(evFolderName, name)
107 //ctor with the indication on where to look for the track segments
109 fFormula = new TFormula("LambdaCuts","(x>1)*(x<2.5)*(y>0)*(y<x)") ;
111 fCpvEmcDistance = 3.0 ;
114 fEventFolderName = GetTitle() ;
115 fTrackSegmentsTitle = GetName();
116 fRecPointsTitle = GetName();
117 fRecParticlesTitle = GetName();
118 fIDOptions = "dis time" ;
120 fRecParticlesInRun = 0 ;
126 //____________________________________________________________________________
127 AliPHOSPIDv0::~AliPHOSPIDv0()
131 //____________________________________________________________________________
132 Float_t AliPHOSPIDv0::GetDistance(AliPHOSEmcRecPoint * emc,AliPHOSRecPoint * cpv, Option_t * Axis)const
134 // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
136 const AliPHOSGeometry * geom = AliPHOSLoader::GetPHOSGeometry() ;
140 emc->GetLocalPosition(vecEmc) ;
141 cpv->GetLocalPosition(vecCpv) ;
142 if(emc->GetPHOSMod() == cpv->GetPHOSMod()){
144 // Correct to difference in CPV and EMC position due to different distance to center.
145 // we assume, that particle moves from center
146 Float_t dCPV = geom->GetIPtoOuterCoverDistance();
147 Float_t dEMC = geom->GetIPtoCrystalSurface() ;
149 vecCpv = dEMC * vecCpv - vecEmc ;
150 if (Axis == "X") return vecCpv.X();
151 if (Axis == "Y") return vecCpv.Y();
152 if (Axis == "Z") return vecCpv.Z();
153 if (Axis == "R") return vecCpv.Mag();
159 //____________________________________________________________________________
160 void AliPHOSPIDv0::Exec(Option_t * option)
164 if( strcmp(GetName(), "")== 0 )
167 if(strstr(option,"tim"))
168 gBenchmark->Start("PHOSPID");
170 if(strstr(option,"print")) {
175 AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
178 Error("Exec","Can not find run getter in event folder \"%s\"",GetTitle());
182 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
185 Error("Exec","Could not obtain the Loader object !");
189 if(gime->BranchExists("RecParticles") )
193 Int_t nevents = runget->GetNumberOfEvents() ; //(Int_t) gAlice->TreeE()->GetEntries() ;
197 for(ievent = 0; ievent < nevents; ievent++){
198 runget->GetEvent(ievent);
199 Info("Exec", "event %d %d %d", ievent, gime->EmcRecPoints(), gime->TrackSegments()) ;
204 if(strstr(option,"deb"))
205 PrintRecParticles(option) ;
207 //increment the total number of rec particles per run
208 fRecParticlesInRun += gime->RecParticles()->GetEntriesFast() ;
212 if(strstr(option,"tim")){
213 gBenchmark->Stop("PHOSPID");
214 Info("Exec", "took %f seconds for PID %f seconds per event",
215 gBenchmark->GetCpuTime("PHOSPID"), gBenchmark->GetCpuTime("PHOSPID")/nevents) ;
219 //____________________________________________________________________________
220 void AliPHOSPIDv0::Init()
222 // Make all memory allocations that are not possible in default constructor
223 // Add the PID task to the list of PHOS tasks
225 AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
228 Error("Exec","Can not find run getter in event folder \"%s\"",GetTitle());
232 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
235 Error("Exec","Could not obtain the Loader object !");
240 gime->LoadRecParticles("UPDATE");
244 //____________________________________________________________________________
245 void AliPHOSPIDv0::MakeRecParticles(){
247 TString taskName(GetName()) ;
248 taskName.Remove(taskName.Index(Version())-1) ;
250 AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
253 Error("Exec","Can not find run getter in event folder \"%s\"",GetTitle());
257 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
260 Error("Exec","Could not obtain the Loader object !");
264 TObjArray * emcRecPoints = gime->EmcRecPoints() ;
265 TObjArray * cpvRecPoints = gime->CpvRecPoints() ;
266 TClonesArray * trackSegments = gime->TrackSegments() ;
267 TClonesArray * recParticles = gime->RecParticles() ;
269 recParticles->Clear();
271 TIter next(trackSegments) ;
272 AliPHOSTrackSegment * ts ;
274 AliPHOSRecParticle * rp ;
276 Bool_t ellips = fIDOptions.Contains("ell",TString::kIgnoreCase ) ;
277 Bool_t disp = fIDOptions.Contains("dis",TString::kIgnoreCase ) ;
278 Bool_t time = fIDOptions.Contains("tim",TString::kIgnoreCase ) ;
280 while ( (ts = (AliPHOSTrackSegment *)next()) ) {
282 new( (*recParticles)[index] ) AliPHOSRecParticle() ;
283 rp = (AliPHOSRecParticle *)recParticles->At(index) ;
284 rp->SetTrackSegment(index) ;
285 rp->SetIndexInList(index) ;
287 AliPHOSEmcRecPoint * emc = 0 ;
288 if(ts->GetEmcIndex()>=0)
289 emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ;
291 AliPHOSRecPoint * cpv = 0 ;
292 if(ts->GetCpvIndex()>=0)
293 cpv = (AliPHOSRecPoint *) cpvRecPoints->At(ts->GetCpvIndex()) ;
295 //set momentum and energy first
296 Float_t e = emc->GetEnergy() ;
297 TVector3 dir = GetMomentumDirection(emc,cpv) ;
300 rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),e) ;
303 //now set type (reconstructed) of the particle
304 Int_t showerprofile = 0; // 0 narrow and 1 wide
308 emc->GetElipsAxis(lambda) ;
309 if(fFormula->Eval(lambda[0],lambda[1]) <= 0 )
310 showerprofile = 1 ; // not narrow
314 if(emc->GetDispersion() > fDispersion )
315 showerprofile = 1 ; // not narrow
319 if(emc->GetTime() > fTimeGate )
322 // Looking at the CPV detector
323 Int_t cpvdetector= 0 ; //1 hit and 0 no hit
325 if(GetDistance(emc, cpv, "R") < fCpvEmcDistance)
328 Int_t type = showerprofile + 2 * slow + 4 * cpvdetector ;
330 rp->SetProductionVertex(0,0,0,0);
331 rp->SetFirstMother(-1);
332 rp->SetLastMother(-1);
333 rp->SetFirstDaughter(-1);
334 rp->SetLastDaughter(-1);
335 rp->SetPolarisation(0,0,0);
341 //____________________________________________________________________________
342 void AliPHOSPIDv0:: Print() const
344 // Print the parameters used for the particle type identification
346 message = "=============== AliPHOSPIDv0 ================\n" ;
347 message += "Making PID\n" ;
348 message += " Headers file: %s\n" ;
349 message += " RecPoints branch title: %s\n" ;
350 message += " TrackSegments Branch title: %s\n" ;
351 message += " RecParticles Branch title %s\n" ;
352 message += "with parameters:\n" ;
353 message += " Maximal EMC - CPV distance (cm) %f\n" ;
354 Info("Print", message.Data(),
356 fRecPointsTitle.Data(),
357 fTrackSegmentsTitle.Data(),
358 fRecParticlesTitle.Data(),
361 if(fIDOptions.Contains("dis",TString::kIgnoreCase ))
362 Info("Print", " dispersion cut %f", fDispersion ) ;
363 if(fIDOptions.Contains("ell",TString::kIgnoreCase ))
364 Info("Print", " Eliptic cuts function: %s", fFormula->GetTitle() ) ;
365 if(fIDOptions.Contains("tim",TString::kIgnoreCase ))
366 Info("Print", " Time Gate used: %f", fTimeGate) ;
369 //____________________________________________________________________________
370 void AliPHOSPIDv0::SetShowerProfileCut(char * formula)
372 //set shape of the cut on the axis of ellipce, drown around shouer
373 //shower considered "narrow" if Formula(lambda[0],lambda[1]) > 0.
376 fFormula = new TFormula("Lambda Cut",formula) ;
378 //____________________________________________________________________________
379 void AliPHOSPIDv0::WriteRecParticles()
382 AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
385 Error("Exec","Can not find run getter in event folder \"%s\"",GetTitle());
389 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
392 Error("Exec","Could not obtain the Loader object !");
396 TClonesArray * recParticles = gime->RecParticles() ;
397 recParticles->Expand(recParticles->GetEntriesFast() ) ;
399 TTree * treeR = gime->TreeR();
403 treeR = gime->TreeR() ;
407 Int_t bufferSize = 32000 ;
408 TBranch * rpBranch = treeR->Branch("PHOSRP",&recParticles,bufferSize);
409 rpBranch->SetTitle(fRecParticlesTitle);
412 Int_t splitlevel = 0 ;
413 AliPHOSPIDv0 * pid = this ;
414 TBranch * pidBranch = treeR->Branch("AliPHOSPID","AliPHOSPIDv0",&pid,bufferSize,splitlevel);
415 pidBranch->SetTitle(fRecParticlesTitle.Data());
420 gime->WriteRecParticles("OVERWRITE");
421 gime->WritePID("OVERWRITE");
425 //____________________________________________________________________________
426 void AliPHOSPIDv0::PlotDispersionCuts()const
428 // produces a plot of the dispersion cut
429 TCanvas* lambdas = new TCanvas("lambdas","Cuts on the ellipse axis",200,10,700,500);
431 if(fIDOptions.Contains("ell",TString::kIgnoreCase ) ){
432 TF2 * ell = new TF2("Elliptic Cuts",fFormula->GetName(),0,3,0,3) ;
433 ell->SetMinimum(0.0000001) ;
434 ell->SetMaximum(0.001) ;
435 ell->SetLineStyle(1) ;
436 ell->SetLineWidth(2) ;
440 if( fIDOptions.Contains("dis",TString::kIgnoreCase ) ){
441 TF2 * dsp = new TF2("dispersion","(y<x)*(x*x+y*y < [0]*[0])",0,3,0,3) ;
442 dsp->SetParameter(0,fDispersion) ;
443 dsp->SetMinimum(0.0000001) ;
444 dsp->SetMaximum(0.001) ;
445 dsp->SetLineStyle(1) ;
446 dsp->SetLineColor(2) ;
447 dsp->SetLineWidth(2) ;
450 if(fIDOptions.Contains("ell",TString::kIgnoreCase ) )
458 //____________________________________________________________________________
459 TVector3 AliPHOSPIDv0::GetMomentumDirection(AliPHOSEmcRecPoint * emc, AliPHOSRecPoint * cpv)const
461 // Calculates the momentum direction:
462 // 1. if only a EMC RecPoint, direction is given by IP and this RecPoint
463 // 2. if a EMC RecPoint and CPV RecPoint, direction is given by the line through the 2 recpoints
464 // However because of the poor position resolution of PPSD the direction is always taken as if we were
467 TVector3 dir(0,0,0) ;
469 TVector3 emcglobalpos ;
472 emc->GetGlobalPosition(emcglobalpos, dummy) ;
475 // The following commented code becomes valid once the PPSD provides
476 // a reasonable position resolution, at least as good as EMC !
477 // TVector3 ppsdlglobalpos ;
478 // TVector3 ppsduglobalpos ;
479 // if( fPpsdLowRecPoint ){ // certainly a photon that has concerted
480 // fPpsdLowRecPoint->GetGlobalPosition(ppsdlglobalpos, mdummy) ;
481 // dir = emcglobalpos - ppsdlglobalpos ;
482 // if( fPpsdUpRecPoint ){ // not looks like a charged
483 // fPpsdUpRecPoint->GetGlobalPosition(ppsduglobalpos, mdummy) ;
484 // dir = ( dir + emcglobalpos - ppsduglobalpos ) * 0.5 ;
487 // else { // looks like a neutral
488 // dir = emcglobalpos ;
492 dir.SetZ( -dir.Z() ) ; // why ?
495 //account correction to the position of IP
496 Float_t xo,yo,zo ; //Coordinates of the origin
497 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
498 TVector3 origin(xo,yo,zo);
503 //____________________________________________________________________________
504 void AliPHOSPIDv0::PrintRecParticles(Option_t * option)
506 // Print table of reconstructed particles
508 AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
511 Error("WriteRecParticles","Can not find run getter in event folder \"%s\"",GetTitle());
515 AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
518 Error("WriteRecParticles","Could not obtain the Loader object !");
522 TString taskName(GetName()) ;
523 taskName.Remove(taskName.Index(Version())-1) ;
524 TClonesArray * recParticles = gime->RecParticles() ;
527 message = "event %d\n" ;
528 message += " found %d RecParticles\n" ;
529 Info("PrintRecParticles", message.Data(), gAlice->GetEvNumber(), recParticles->GetEntriesFast() ) ;
531 if(strstr(option,"all")) { // printing found TS
532 Info("PrintRecParticles"," PARTICLE Index \n" ) ;
535 for (index = 0 ; index < recParticles->GetEntries() ; index++) {
536 AliPHOSRecParticle * rp = (AliPHOSRecParticle * ) recParticles->At(index) ;
539 switch(rp->GetType()) {
540 case AliPHOSFastRecParticle::kNEUTRALEMFAST:
541 strcpy( particle, "NEUTRAL EM FAST");
543 case AliPHOSFastRecParticle::kNEUTRALHAFAST:
544 strcpy(particle, "NEUTRAL HA FAST");
546 case AliPHOSFastRecParticle::kNEUTRALEMSLOW:
547 strcpy(particle, "NEUTRAL EM SLOW");
549 case AliPHOSFastRecParticle::kNEUTRALHASLOW:
550 strcpy(particle, "NEUTRAL HA SLOW");
552 case AliPHOSFastRecParticle::kCHARGEDEMFAST:
553 strcpy(particle, "CHARGED EM FAST") ;
555 case AliPHOSFastRecParticle::kCHARGEDHAFAST:
556 strcpy(particle, "CHARGED HA FAST") ;
558 case AliPHOSFastRecParticle::kCHARGEDEMSLOW:
559 strcpy(particle, "CHARGEDEMSLOW") ;
561 case AliPHOSFastRecParticle::kCHARGEDHASLOW:
562 strcpy(particle, "CHARGED HA SLOW") ;
566 // Int_t * primaries;
568 // primaries = rp->GetPrimaries(nprimaries);
570 Info("PrintRecParticles", " %s %d\n", particle, rp->GetIndexInList()) ;