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 ---
74 // --- AliRoot header files ---
77 #include "AliGenerator.h"
79 #include "AliPHOSPIDv0.h"
80 #include "AliPHOSClusterizerv1.h"
81 #include "AliPHOSTrackSegment.h"
82 #include "AliPHOSTrackSegmentMakerv1.h"
83 #include "AliPHOSRecParticle.h"
84 #include "AliPHOSGeometry.h"
85 #include "AliPHOSGetter.h"
87 ClassImp( AliPHOSPIDv0)
89 //____________________________________________________________________________
90 AliPHOSPIDv0::AliPHOSPIDv0():AliPHOSPID()
97 fHeaderFileName = "" ;
98 fTrackSegmentsTitle= "" ;
99 fRecPointsTitle = "" ;
100 fRecParticlesTitle = "" ;
101 fIDOptions = "dis time" ;
102 fRecParticlesInRun = 0 ;
107 //____________________________________________________________________________
108 AliPHOSPIDv0::AliPHOSPIDv0(const char * headerFile,const char * name, const Bool_t toSplit) : AliPHOSPID(headerFile, name,toSplit)
110 //ctor with the indication on where to look for the track segments
112 fFormula = new TFormula("LambdaCuts","(x>1)*(x<2.5)*(y>0)*(y<x)") ;
114 fCpvEmcDistance = 3.0 ;
117 fHeaderFileName = GetTitle() ;
118 fTrackSegmentsTitle = GetName() ;
119 fRecPointsTitle = GetName() ;
120 fRecParticlesTitle = GetName() ;
121 fIDOptions = "dis time" ;
123 TString tempo(GetName()) ;
125 tempo.Append(Version()) ;
127 fRecParticlesInRun = 0 ;
133 //____________________________________________________________________________
134 AliPHOSPIDv0::~AliPHOSPIDv0()
138 //____________________________________________________________________________
139 Float_t AliPHOSPIDv0::GetDistance(AliPHOSEmcRecPoint * emc,AliPHOSRecPoint * cpv, Option_t * Axis)const
141 // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
143 const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
147 emc->GetLocalPosition(vecEmc) ;
148 cpv->GetLocalPosition(vecCpv) ;
149 if(emc->GetPHOSMod() == cpv->GetPHOSMod()){
151 // Correct to difference in CPV and EMC position due to different distance to center.
152 // we assume, that particle moves from center
153 Float_t dCPV = geom->GetIPtoOuterCoverDistance();
154 Float_t dEMC = geom->GetIPtoCrystalSurface() ;
156 vecCpv = dEMC * vecCpv - vecEmc ;
157 if (Axis == "X") return vecCpv.X();
158 if (Axis == "Y") return vecCpv.Y();
159 if (Axis == "Z") return vecCpv.Z();
160 if (Axis == "R") return vecCpv.Mag();
166 //____________________________________________________________________________
167 void AliPHOSPIDv0::Exec(Option_t * option)
171 if( strcmp(GetName(), "")== 0 )
174 if(strstr(option,"tim"))
175 gBenchmark->Start("PHOSPID");
177 if(strstr(option,"print")) {
182 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
183 if(gime->BranchExists("RecParticles") )
186 // gAlice->GetEvent(0) ;
187 // //check, if the branch with name of this" already exits?
188 // TObjArray * lob = (TObjArray*)gAlice->TreeR()->GetListOfBranches() ;
190 // TBranch * branch = 0 ;
191 // Bool_t phospidfound = kFALSE, pidfound = kFALSE ;
193 // TString taskName(GetName()) ;
194 // taskName.Remove(taskName.Index(Version())-1) ;
196 // while ( (branch = (TBranch*)next()) && (!phospidfound || !pidfound) ) {
197 // if ( (strcmp(branch->GetName(), "PHOSPID")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) )
198 // phospidfound = kTRUE ;
200 // else if ( (strcmp(branch->GetName(), "AliPHOSPID")==0) && (strcmp(branch->GetTitle(), taskName.Data())==0) )
201 // pidfound = kTRUE ;
204 // if ( phospidfound || pidfound ) {
205 // cerr << "WARNING: AliPHOSPIDv0::Exec -> RecParticles and/or PIDtMaker branch with name "
206 // << taskName.Data() << " already exits" << endl ;
210 Int_t nevents = gime->MaxEvent() ; //(Int_t) gAlice->TreeE()->GetEntries() ;
213 for(ievent = 0; ievent < nevents; ievent++){
214 gime->Event(ievent,"R") ;
215 cout << "event " << ievent << " " << gime->EmcRecPoints() << " " << gime->TrackSegments() << endl ;
218 WriteRecParticles(ievent);
220 if(strstr(option,"deb"))
221 PrintRecParticles(option) ;
223 //increment the total number of rec particles per run
224 fRecParticlesInRun += gime->RecParticles()->GetEntriesFast() ;
228 if(strstr(option,"tim")){
229 gBenchmark->Stop("PHOSPID");
230 cout << "AliPHOSPID:" << endl ;
231 cout << " took " << gBenchmark->GetCpuTime("PHOSPID") << " seconds for PID "
232 << gBenchmark->GetCpuTime("PHOSPID")/nevents << " seconds per event " << endl ;
237 //____________________________________________________________________________
238 void AliPHOSPIDv0::Init()
240 // Make all memory allocations that are not possible in default constructor
241 // Add the PID task to the list of PHOS tasks
243 if ( strcmp(GetTitle(), "") == 0 )
244 SetTitle("galice.root") ;
246 TString taskName(GetName()) ;
247 taskName.Remove(taskName.Index(Version())-1) ;
249 AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), taskName.Data(),fToSplit) ;
251 cerr << "ERROR: AliPHOSPIDv0::Init -> Could not obtain the Getter object !" << endl ;
257 //First - extract full path if necessary
258 TString fileName(GetTitle()) ;
259 Ssiz_t islash = fileName.Last('/') ;
260 if(islash<fileName.Length())
261 fileName.Remove(islash+1,fileName.Length()) ;
264 fileName+="PHOS.RecData." ;
265 if((strcmp(taskName.Data(),"Default")!=0)&&(strcmp(taskName.Data(),"")!=0)){
270 fSplitFile = static_cast<TFile*>(gROOT->GetFile(fileName.Data()));
272 fSplitFile = TFile::Open(fileName.Data(),"update") ;
276 gime->PostPID(this) ;
277 // create a folder on the white board //YSAlice/WhiteBoard/RecParticles/PHOS/recparticlesName
278 gime->PostRecParticles(taskName.Data() ) ;
282 //____________________________________________________________________________
283 void AliPHOSPIDv0::MakeRecParticles(){
285 // Makes a RecParticle out of a TrackSegment
286 TString taskName(GetName()) ;
287 taskName.Remove(taskName.Index(Version())-1) ;
289 AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
290 TObjArray * emcRecPoints = gime->EmcRecPoints(taskName) ;
291 TObjArray * cpvRecPoints = gime->CpvRecPoints(taskName) ;
292 TClonesArray * trackSegments = gime->TrackSegments(taskName) ;
293 TClonesArray * recParticles = gime->RecParticles(taskName) ;
294 recParticles->Clear();
296 TIter next(trackSegments) ;
297 AliPHOSTrackSegment * ts ;
299 AliPHOSRecParticle * rp ;
301 Bool_t ellips = fIDOptions.Contains("ell",TString::kIgnoreCase ) ;
302 Bool_t disp = fIDOptions.Contains("dis",TString::kIgnoreCase ) ;
303 Bool_t time = fIDOptions.Contains("tim",TString::kIgnoreCase ) ;
305 while ( (ts = (AliPHOSTrackSegment *)next()) ) {
307 new( (*recParticles)[index] ) AliPHOSRecParticle() ;
308 rp = (AliPHOSRecParticle *)recParticles->At(index) ;
309 rp->SetTrackSegment(index) ;
310 rp->SetIndexInList(index) ;
312 AliPHOSEmcRecPoint * emc = 0 ;
313 if(ts->GetEmcIndex()>=0)
314 emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ;
316 AliPHOSRecPoint * cpv = 0 ;
317 if(ts->GetCpvIndex()>=0)
318 cpv = (AliPHOSRecPoint *) cpvRecPoints->At(ts->GetCpvIndex()) ;
320 //set momentum and energy first
321 Float_t e = emc->GetEnergy() ;
322 TVector3 dir = GetMomentumDirection(emc,cpv) ;
325 rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),e) ;
328 //now set type (reconstructed) of the particle
329 Int_t showerprofile = 0; // 0 narrow and 1 wide
333 emc->GetElipsAxis(lambda) ;
334 if(fFormula->Eval(lambda[0],lambda[1]) <= 0 )
335 showerprofile = 1 ; // not narrow
339 if(emc->GetDispersion() > fDispersion )
340 showerprofile = 1 ; // not narrow
344 if(emc->GetTime() > fTimeGate )
347 // Looking at the CPV detector
348 Int_t cpvdetector= 0 ; //1 hit and 0 no hit
350 if(GetDistance(emc, cpv, "R") < fCpvEmcDistance)
353 Int_t type = showerprofile + 2 * slow + 4 * cpvdetector ;
355 rp->SetProductionVertex(0,0,0,0);
356 rp->SetFirstMother(-1);
357 rp->SetLastMother(-1);
358 rp->SetFirstDaughter(-1);
359 rp->SetLastDaughter(-1);
360 rp->SetPolarisation(0,0,0);
366 //____________________________________________________________________________
367 void AliPHOSPIDv0:: Print(Option_t * option) const
369 // Print the parameters used for the particle type identification
370 cout << "=============== AliPHOSPID1 ================" << endl ;
371 cout << "Making PID "<< endl ;
372 cout << " Headers file: " << fHeaderFileName.Data() << endl ;
373 cout << " RecPoints branch title: " << fRecPointsTitle.Data() << endl ;
374 cout << " TrackSegments Branch title: " << fTrackSegmentsTitle.Data() << endl ;
375 cout << " RecParticles Branch title " << fRecParticlesTitle.Data() << endl;
376 cout << "with parameters: " << endl ;
377 cout << " Maximal EMC - CPV distance (cm) " << fCpvEmcDistance << endl ;
378 if(fIDOptions.Contains("dis",TString::kIgnoreCase ))
379 cout << " dispersion cut " << fDispersion << endl ;
380 if(fIDOptions.Contains("ell",TString::kIgnoreCase )){
381 cout << " Eliptic cuts function: " << endl ;
382 cout << fFormula->GetTitle() << endl ;
384 if(fIDOptions.Contains("tim",TString::kIgnoreCase ))
385 cout << " Time Gate uzed: " << fTimeGate << endl ;
386 cout << "============================================" << endl ;
389 //____________________________________________________________________________
390 void AliPHOSPIDv0::SetShowerProfileCut(char * formula)
392 //set shape of the cut on the axis of ellipce, drown around shouer
393 //shower considered "narrow" if Formula(lambda[0],lambda[1]) > 0.
396 fFormula = new TFormula("Lambda Cut",formula) ;
398 //____________________________________________________________________________
399 void AliPHOSPIDv0::WriteRecParticles(Int_t event)
402 AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ;
403 TString taskName(GetName()) ;
404 taskName.Remove(taskName.Index(Version())-1) ;
405 TClonesArray * recParticles = gime->RecParticles(taskName) ;
406 recParticles->Expand(recParticles->GetEntriesFast() ) ;
415 sprintf(name,"%s%d", "TreeR",event) ;
416 treeR = dynamic_cast<TTree*>(fSplitFile->Get(name));
419 treeR = gAlice->TreeR();
423 gAlice->MakeTree("R", fSplitFile);
424 treeR = gAlice->TreeR() ;
427 // //Make branch in TreeR for RecParticles
428 // char * filename = 0;
429 // if(gSystem->Getenv("CONFIG_SPLIT_FILE")!=0){ //generating file name
430 // filename = new char[strlen(gAlice->GetBaseFile())+20] ;
431 // sprintf(filename,"%s/PHOS.Reco.root",gAlice->GetBaseFile()) ;
434 // TDirectory *cwd = gDirectory;
437 Int_t bufferSize = 32000 ;
438 TBranch * rpBranch = treeR->Branch("PHOSRP",&recParticles,bufferSize);
439 rpBranch->SetTitle(fRecParticlesTitle);
441 // rpBranch->SetFile(filename);
442 // TIter next( rpBranch->GetListOfBranches());
444 // while ((sb=(TBranch*)next())) {
445 // sb->SetFile(filename);
451 Int_t splitlevel = 0 ;
452 AliPHOSPIDv0 * pid = this ;
453 TBranch * pidBranch = treeR->Branch("AliPHOSPID","AliPHOSPIDv0",&pid,bufferSize,splitlevel);
454 pidBranch->SetTitle(fRecParticlesTitle.Data());
456 // pidBranch->SetFile(filename);
457 // TIter next( pidBranch->GetListOfBranches());
459 // while ((sb=(TBranch*)next())) {
460 // sb->SetFile(filename);
468 treeR->AutoSave() ; //Write(0,kOverwrite) ;
472 //____________________________________________________________________________
473 void AliPHOSPIDv0::PlotDispersionCuts()const
475 // produces a plot of the dispersion cut
476 TCanvas* lambdas = new TCanvas("lambdas","Cuts on the ellipse axis",200,10,700,500);
478 if(fIDOptions.Contains("ell",TString::kIgnoreCase ) ){
479 TF2 * ell = new TF2("Elliptic Cuts",fFormula->GetName(),0,3,0,3) ;
480 ell->SetMinimum(0.0000001) ;
481 ell->SetMaximum(0.001) ;
482 ell->SetLineStyle(1) ;
483 ell->SetLineWidth(2) ;
487 if( fIDOptions.Contains("dis",TString::kIgnoreCase ) ){
488 TF2 * dsp = new TF2("dispersion","(y<x)*(x*x+y*y < [0]*[0])",0,3,0,3) ;
489 dsp->SetParameter(0,fDispersion) ;
490 dsp->SetMinimum(0.0000001) ;
491 dsp->SetMaximum(0.001) ;
492 dsp->SetLineStyle(1) ;
493 dsp->SetLineColor(2) ;
494 dsp->SetLineWidth(2) ;
497 if(fIDOptions.Contains("ell",TString::kIgnoreCase ) )
505 //____________________________________________________________________________
506 TVector3 AliPHOSPIDv0::GetMomentumDirection(AliPHOSEmcRecPoint * emc, AliPHOSRecPoint * cpv)const
508 // Calculates the momentum direction:
509 // 1. if only a EMC RecPoint, direction is given by IP and this RecPoint
510 // 2. if a EMC RecPoint and CPV RecPoint, direction is given by the line through the 2 recpoints
511 // However because of the poor position resolution of PPSD the direction is always taken as if we were
514 TVector3 dir(0,0,0) ;
516 TVector3 emcglobalpos ;
519 emc->GetGlobalPosition(emcglobalpos, dummy) ;
522 // The following commented code becomes valid once the PPSD provides
523 // a reasonable position resolution, at least as good as EMC !
524 // TVector3 ppsdlglobalpos ;
525 // TVector3 ppsduglobalpos ;
526 // if( fPpsdLowRecPoint ){ // certainly a photon that has concerted
527 // fPpsdLowRecPoint->GetGlobalPosition(ppsdlglobalpos, mdummy) ;
528 // dir = emcglobalpos - ppsdlglobalpos ;
529 // if( fPpsdUpRecPoint ){ // not looks like a charged
530 // fPpsdUpRecPoint->GetGlobalPosition(ppsduglobalpos, mdummy) ;
531 // dir = ( dir + emcglobalpos - ppsduglobalpos ) * 0.5 ;
534 // else { // looks like a neutral
535 // dir = emcglobalpos ;
539 dir.SetZ( -dir.Z() ) ; // why ?
542 //account correction to the position of IP
543 Float_t xo,yo,zo ; //Coordinates of the origin
544 gAlice->Generator()->GetOrigin(xo,yo,zo) ;
545 TVector3 origin(xo,yo,zo);
550 //____________________________________________________________________________
551 void AliPHOSPIDv0::PrintRecParticles(Option_t * option)
553 // Print table of reconstructed particles
555 AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ;
557 TString taskName(GetName()) ;
558 taskName.Remove(taskName.Index(Version())-1) ;
559 TClonesArray * recParticles = gime->RecParticles(taskName) ;
561 cout << "AliPHOSPIDv0: event "<<gAlice->GetEvNumber() << endl ;
562 cout << " found " << recParticles->GetEntriesFast() << " RecParticles " << endl ;
564 if(strstr(option,"all")) { // printing found TS
567 << " Index " << endl ;
571 // << " # of primaries "
572 // << " Primaries list " << endl;
575 for (index = 0 ; index < recParticles->GetEntries() ; index++) {
576 AliPHOSRecParticle * rp = (AliPHOSRecParticle * ) recParticles->At(index) ;
579 switch(rp->GetType()) {
580 case AliPHOSFastRecParticle::kNEUTRALEMFAST:
581 strcpy( particle, "NEUTRAL EM FAST");
583 case AliPHOSFastRecParticle::kNEUTRALHAFAST:
584 strcpy(particle, "NEUTRAL HA FAST");
586 case AliPHOSFastRecParticle::kNEUTRALEMSLOW:
587 strcpy(particle, "NEUTRAL EM SLOW");
589 case AliPHOSFastRecParticle::kNEUTRALHASLOW:
590 strcpy(particle, "NEUTRAL HA SLOW");
592 case AliPHOSFastRecParticle::kCHARGEDEMFAST:
593 strcpy(particle, "CHARGED EM FAST") ;
595 case AliPHOSFastRecParticle::kCHARGEDHAFAST:
596 strcpy(particle, "CHARGED HA FAST") ;
598 case AliPHOSFastRecParticle::kCHARGEDEMSLOW:
599 strcpy(particle, "CHARGEDEMSLOW") ;
601 case AliPHOSFastRecParticle::kCHARGEDHASLOW:
602 strcpy(particle, "CHARGED HA SLOW") ;
606 // Int_t * primaries;
608 // primaries = rp->GetPrimaries(nprimaries);
610 cout << setw(10) << particle << " "
611 << setw(5) << rp->GetIndexInList() << " " ;
612 // << setw(4) << nprimaries << " ";
613 // for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
614 // cout << setw(4) << primaries[iprimary] << " ";
617 cout << "-------------------------------------------" << endl ;