]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSPIDv0.cxx
Important comment added.
[u/mrichter/AliRoot.git] / PHOS / AliPHOSPIDv0.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id$ */
17
18 //_________________________________________________________________________
19 // Implementation version v0 of the PHOS particle identifier 
20 // Particle identification based on the 
21 //     - CPV information, 
22 //     - Preshower information (in MIXT or GPS2 geometries)
23 //     - shower width.
24 //
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)  
27 //
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().    
37 //
38 // use case:
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
54 //                  
55 //*-- Author: Yves Schutz (SUBATECH)  & Gines Martinez (SUBATECH) & 
56 //            Dmitri Peressounko (SUBATECH & Kurchatov Institute)
57 //            Completely redesined by Dmitri Peressounko, March 2001
58
59 // --- ROOT system ---
60 #include "TTree.h"
61 #include "TF2.h"
62 #include "TFormula.h"
63 #include "TCanvas.h"
64
65 #include "TBenchmark.h"
66 // --- Standard library ---
67
68 // --- AliRoot header files ---
69
70 #include "AliRun.h"
71 #include "AliGenerator.h"
72 #include "AliPHOSPIDv0.h"
73 #include "AliPHOSEmcRecPoint.h"
74 #include "AliPHOSTrackSegment.h"
75 #include "AliPHOSRecParticle.h"
76 #include "AliPHOSGeometry.h"
77 #include "AliPHOSLoader.h"
78
79 ClassImp( AliPHOSPIDv0) 
80
81 //____________________________________________________________________________
82 AliPHOSPIDv0::AliPHOSPIDv0():AliPHOSPID()
83
84   // default ctor
85   fFormula           = 0 ;
86   fDispersion        = 0. ; 
87   fCpvEmcDistance    = 0 ; 
88   fTimeGate          = 2.e-9 ;
89   fEventFolderName    = "" ; 
90   fTrackSegmentsTitle= "" ; 
91   fRecPointsTitle    = "" ; 
92   fRecParticlesTitle = "" ; 
93   fIDOptions         = "dis time" ; 
94   fRecParticlesInRun = 0 ;
95   fClusterizer = 0;
96   fTSMaker = 0;
97 }
98
99 //____________________________________________________________________________
100 AliPHOSPIDv0::AliPHOSPIDv0(const char * evFolderName,const char * name) : AliPHOSPID(evFolderName, name)
101
102   //ctor with the indication on where to look for the track segments
103
104   fFormula        = new TFormula("LambdaCuts","(x>1)*(x<2.5)*(y>0)*(y<x)") ;   
105   fDispersion     = 2.0 ; 
106   fCpvEmcDistance = 3.0 ;
107   fTimeGate          = 2.e-9 ;
108  
109   fEventFolderName     = GetTitle() ; 
110   fTrackSegmentsTitle = GetName();
111   fRecPointsTitle     = GetName();
112   fRecParticlesTitle  = GetName();
113   fIDOptions          = "dis time" ;
114     
115   fRecParticlesInRun = 0 ; 
116
117   Init() ;
118
119 }
120
121 //____________________________________________________________________________
122 AliPHOSPIDv0::~AliPHOSPIDv0()
123
124 }
125
126 //____________________________________________________________________________
127 Float_t  AliPHOSPIDv0::GetDistance(AliPHOSEmcRecPoint * emc,AliPHOSRecPoint * cpv, Option_t *  Axis)const
128 {
129   // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
130  
131   const AliPHOSGeometry * geom = AliPHOSLoader::GetPHOSGeometry() ; 
132   TVector3 vecEmc ;
133   TVector3 vecCpv ;
134   
135   emc->GetLocalPosition(vecEmc) ;
136   cpv->GetLocalPosition(vecCpv) ; 
137   if(emc->GetPHOSMod() == cpv->GetPHOSMod()){ 
138     
139     // Correct to difference in CPV and EMC position due to different distance to center.
140     // we assume, that particle moves from center
141     Float_t dCPV = geom->GetIPtoOuterCoverDistance();
142     Float_t dEMC = geom->GetIPtoCrystalSurface() ;
143     dEMC         = dEMC / dCPV ;
144     vecCpv = dEMC * vecCpv  - vecEmc ; 
145     if (Axis == "X") return vecCpv.X();
146     if (Axis == "Y") return vecCpv.Y();
147     if (Axis == "Z") return vecCpv.Z();
148     if (Axis == "R") return vecCpv.Mag();
149   } 
150  
151   return 100000000 ;
152 }
153
154 //____________________________________________________________________________
155 void  AliPHOSPIDv0::Exec(Option_t * option) 
156 {
157   //Steering method
158   
159   if( strcmp(GetName(), "")== 0 ) 
160     Init() ;
161   
162   if(strstr(option,"tim"))
163     gBenchmark->Start("PHOSPID");
164   
165   if(strstr(option,"print")) {
166     Print() ; 
167     return ; 
168   }
169
170   AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
171   if(runget == 0x0) 
172    {
173      Error("Exec","Can not find run getter in event folder \"%s\"",GetTitle());
174      return;
175    }
176   
177   AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
178   if ( gime == 0 ) 
179    {
180      Error("Exec","Could not obtain the Loader object !"); 
181      return ;
182    } 
183
184   if(gime->BranchExists("RecParticles") )
185     return ;
186
187   
188   Int_t nevents = runget->GetNumberOfEvents() ;       //(Int_t) gAlice->TreeE()->GetEntries() ;
189
190   Int_t ievent ;
191   
192   for(ievent = 0; ievent < nevents; ievent++){
193     runget->GetEvent(ievent);
194     Info("Exec", "event %d %d %d", ievent, gime->EmcRecPoints(), gime->TrackSegments()) ;
195     MakeRecParticles() ;
196     
197     WriteRecParticles();
198     
199     if(strstr(option,"deb"))
200       PrintRecParticles(option) ;
201
202     //increment the total number of rec particles per run 
203     fRecParticlesInRun += gime->RecParticles()->GetEntriesFast() ; 
204
205   }
206   
207   if(strstr(option,"tim")){
208     gBenchmark->Stop("PHOSPID");
209     Info("Exec", "took %f seconds for PID %f seconds per event", 
210          gBenchmark->GetCpuTime("PHOSPID"), gBenchmark->GetCpuTime("PHOSPID")/nevents) ; 
211   }
212   
213 }
214 //____________________________________________________________________________
215 void AliPHOSPIDv0::Init()
216 {
217   // Make all memory allocations that are not possible in default constructor
218   // Add the PID task to the list of PHOS tasks
219   
220   AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
221   if(runget == 0x0) 
222    {
223      Error("Exec","Can not find run getter in event folder \"%s\"",GetTitle());
224      return;
225    }
226   
227   AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
228   if ( gime == 0 ) 
229    {
230      Error("Exec","Could not obtain the Loader object !"); 
231      return ;
232    } 
233    
234   gime->PostPID(this);
235   gime->LoadRecParticles("UPDATE");
236   
237 }
238
239 //____________________________________________________________________________
240 void  AliPHOSPIDv0::MakeRecParticles()
241 {
242   // Reconstructs the particles from the tracksegments
243
244   TString taskName(GetName()) ; 
245   taskName.Remove(taskName.Index(Version())-1) ;
246
247   AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
248   if(runget == 0x0) 
249    {
250      Error("Exec","Can not find run getter in event folder \"%s\"",GetTitle());
251      return;
252    }
253   
254   AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
255   if ( gime == 0 ) 
256    {
257      Error("Exec","Could not obtain the Loader object !"); 
258      return ;
259    } 
260
261   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
262   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
263   TClonesArray * trackSegments = gime->TrackSegments() ; 
264   TClonesArray * recParticles  = gime->RecParticles() ; 
265
266   recParticles->Clear();
267   
268   TIter next(trackSegments) ; 
269   AliPHOSTrackSegment * ts ; 
270   Int_t index = 0 ; 
271   AliPHOSRecParticle * rp ; 
272   
273   Bool_t ellips = fIDOptions.Contains("ell",TString::kIgnoreCase ) ;
274   Bool_t disp   = fIDOptions.Contains("dis",TString::kIgnoreCase ) ;
275   Bool_t time   = fIDOptions.Contains("tim",TString::kIgnoreCase ) ;
276   
277   while ( (ts = (AliPHOSTrackSegment *)next()) ) {
278     
279     new( (*recParticles)[index] ) AliPHOSRecParticle() ;
280     rp = (AliPHOSRecParticle *)recParticles->At(index) ; 
281     rp->SetTrackSegment(index) ;
282     rp->SetIndexInList(index) ;
283     
284     AliPHOSEmcRecPoint * emc = 0 ;
285     if(ts->GetEmcIndex()>=0)
286       emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ;
287     
288     AliPHOSRecPoint    * cpv = 0 ;
289     if(ts->GetCpvIndex()>=0)
290       cpv = (AliPHOSRecPoint *)   cpvRecPoints->At(ts->GetCpvIndex()) ;
291     
292     //set momentum and energy first
293     Float_t    e = emc->GetEnergy() ;
294     TVector3 dir = GetMomentumDirection(emc,cpv) ; 
295     dir.SetMag(e) ;
296
297     rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),e) ;
298     rp->SetCalcMass(0);
299
300     //now set type (reconstructed) of the particle    
301     Int_t showerprofile = 0;  // 0 narrow and 1 wide
302     
303     if(ellips){
304       Float_t lambda[2] ;
305       emc->GetElipsAxis(lambda) ;
306       if(fFormula->Eval(lambda[0],lambda[1]) <= 0 )
307         showerprofile = 1 ;  // not narrow
308     }
309     
310     if(disp)
311       if(emc->GetDispersion() > fDispersion )
312         showerprofile = 1 ;  // not narrow
313     
314     Int_t slow = 0 ;
315     if(time)
316       if(emc->GetTime() > fTimeGate )
317         slow = 0 ; 
318         
319     // Looking at the CPV detector
320     Int_t cpvdetector= 0 ;  //1 hit and 0 no hit     
321     if(cpv)
322       if(GetDistance(emc, cpv,  "R") < fCpvEmcDistance) 
323         cpvdetector = 1 ;  
324     
325     Int_t type = showerprofile + 2 * slow  + 4 * cpvdetector ;
326     rp->SetType(type) ; 
327     rp->SetProductionVertex(0,0,0,0);
328     rp->SetFirstMother(-1);
329     rp->SetLastMother(-1);
330     rp->SetFirstDaughter(-1);
331     rp->SetLastDaughter(-1);
332     rp->SetPolarisation(0,0,0);
333     index++ ; 
334   }
335   
336 }
337
338 //____________________________________________________________________________
339 void  AliPHOSPIDv0:: Print() const
340 {
341   // Print the parameters used for the particle type identification
342   TString message ; 
343   message  = "=============== AliPHOSPIDv0 ================\n" ;
344   message += "Making PID\n" ;
345   message += "    Headers file:               %s\n" ; 
346   message += "    RecPoints branch title:     %s\n" ;
347   message += "    TrackSegments Branch title: %s\n" ; 
348   message += "    RecParticles Branch title   %s\n" ;  
349   message += "with parameters:\n"  ;
350   message += "    Maximal EMC - CPV  distance (cm) %f\n" ;
351   Info("Print", message.Data(),  
352        GetTitle(), 
353        fRecPointsTitle.Data(), 
354        fTrackSegmentsTitle.Data(), 
355        fRecParticlesTitle.Data(), 
356        fCpvEmcDistance );
357
358   if(fIDOptions.Contains("dis",TString::kIgnoreCase ))
359     Info("Print", "                    dispersion cut %f",  fDispersion ) ;
360   if(fIDOptions.Contains("ell",TString::kIgnoreCase ))
361     Info("Print", "             Eliptic cuts function: %s",  fFormula->GetTitle() ) ;
362   if(fIDOptions.Contains("tim",TString::kIgnoreCase ))
363     Info("Print",  "             Time Gate used: %f",  fTimeGate) ;
364 }
365
366 //____________________________________________________________________________
367 void  AliPHOSPIDv0::SetShowerProfileCut(const char * formula)
368 {
369   //set shape of the cut on the axis of ellipce, drown around shouer
370   //shower considered "narrow" if Formula(lambda[0],lambda[1]) > 0.
371   if(fFormula) 
372     delete fFormula; 
373   fFormula = new TFormula("Lambda Cut",formula) ;
374 }
375 //____________________________________________________________________________
376 void  AliPHOSPIDv0::WriteRecParticles()
377 {
378   // Saves the reconstructed particles too a file
379  
380   AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
381   if(runget == 0x0) 
382    {
383      Error("Exec","Can not find run getter in event folder \"%s\"",GetTitle());
384      return;
385    }
386   
387   AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
388   if ( gime == 0 ) 
389    {
390      Error("Exec","Could not obtain the Loader object !"); 
391      return ;
392    } 
393
394   TClonesArray * recParticles = gime->RecParticles() ; 
395   recParticles->Expand(recParticles->GetEntriesFast() ) ;
396
397   TTree * treeR = gime->TreeR();
398   
399   if(!treeR){
400     gime->MakeTree("R");
401     treeR = gime->TreeR() ;
402   }
403   
404   //First rp
405   Int_t bufferSize = 32000 ;    
406   TBranch * rpBranch = treeR->Branch("PHOSRP",&recParticles,bufferSize);
407   rpBranch->SetTitle(fRecParticlesTitle);
408   
409   //second, pid
410   Int_t splitlevel = 0 ; 
411   AliPHOSPIDv0 * pid = this ;
412   TBranch * pidBranch = treeR->Branch("AliPHOSPID","AliPHOSPIDv0",&pid,bufferSize,splitlevel);
413   pidBranch->SetTitle(fRecParticlesTitle.Data());
414   
415   rpBranch->Fill() ;
416   pidBranch->Fill() ;
417
418   gime->WriteRecParticles("OVERWRITE");
419   gime->WritePID("OVERWRITE");
420
421 }
422
423 //____________________________________________________________________________
424 void  AliPHOSPIDv0::PlotDispersionCuts()const
425 {
426   // produces a plot of the dispersion cut
427   TCanvas*  lambdas = new TCanvas("lambdas","Cuts on the ellipse axis",200,10,700,500);
428  
429   if(fIDOptions.Contains("ell",TString::kIgnoreCase ) ){
430     TF2 * ell = new TF2("Elliptic Cuts",fFormula->GetName(),0,3,0,3) ;
431     ell->SetMinimum(0.0000001) ;
432     ell->SetMaximum(0.001) ;
433     ell->SetLineStyle(1) ;
434     ell->SetLineWidth(2) ;
435     ell->Draw() ;
436   }
437   
438   if( fIDOptions.Contains("dis",TString::kIgnoreCase ) ){
439     TF2 * dsp = new TF2("dispersion","(y<x)*(x*x+y*y < [0]*[0])",0,3,0,3) ;
440     dsp->SetParameter(0,fDispersion) ;
441     dsp->SetMinimum(0.0000001) ;
442     dsp->SetMaximum(0.001) ;
443     dsp->SetLineStyle(1) ;
444     dsp->SetLineColor(2) ;
445     dsp->SetLineWidth(2) ;
446     dsp->SetNpx(200) ;
447     dsp->SetNpy(200) ;
448     if(fIDOptions.Contains("ell",TString::kIgnoreCase ) )
449       dsp->Draw("same") ;
450     else
451       dsp->Draw() ;
452   }
453   lambdas->Update();
454 }
455
456 //____________________________________________________________________________
457 TVector3 AliPHOSPIDv0::GetMomentumDirection(AliPHOSEmcRecPoint * emc, AliPHOSRecPoint * )const 
458
459   // Calculates the momentum direction:
460   //   1. if only a EMC RecPoint, direction is given by IP and this RecPoint
461   //   2. if a EMC RecPoint and CPV RecPoint, direction is given by the line through the 2 recpoints 
462   //  However because of the poor position resolution of PPSD the direction is always taken as if we were 
463   //  in case 1.
464
465   TVector3 dir(0,0,0) ; 
466   
467   TVector3 emcglobalpos ;
468   TMatrix  dummy ;
469   
470   emc->GetGlobalPosition(emcglobalpos, dummy) ;
471   
472  
473   // The following commented code becomes valid once the PPSD provides 
474   // a reasonable position resolution, at least as good as EMC ! 
475   //   TVector3 ppsdlglobalpos ;
476   //   TVector3 ppsduglobalpos ;
477   //   if( fPpsdLowRecPoint ){ // certainly a photon that has concerted
478   //     fPpsdLowRecPoint->GetGlobalPosition(ppsdlglobalpos, mdummy) ; 
479   //     dir = emcglobalpos -  ppsdlglobalpos ; 
480   //     if( fPpsdUpRecPoint ){ // not looks like a charged       
481   //        fPpsdUpRecPoint->GetGlobalPosition(ppsduglobalpos, mdummy) ; 
482   //        dir = ( dir +  emcglobalpos -  ppsduglobalpos ) * 0.5 ; 
483   //      }
484   //   }
485   //   else { // looks like a neutral
486   //    dir = emcglobalpos ;  
487   //  }
488   
489   dir = emcglobalpos ;  
490   dir.SetZ( -dir.Z() ) ;   // why ?  
491   dir.SetMag(1.) ;
492
493   //account correction to the position of IP
494   Float_t xo,yo,zo ; //Coordinates of the origin
495   gAlice->Generator()->GetOrigin(xo,yo,zo) ;
496   TVector3 origin(xo,yo,zo);
497   dir = dir - origin ;
498
499   return dir ;  
500 }
501 //____________________________________________________________________________
502 void AliPHOSPIDv0::PrintRecParticles(Option_t * option)
503 {
504   // Print table of reconstructed particles
505
506   AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
507   if(runget == 0x0) 
508    {
509      Error("WriteRecParticles","Can not find run getter in event folder \"%s\"",GetTitle());
510      return;
511    }
512   
513   AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
514   if ( gime == 0 ) 
515    {
516      Error("WriteRecParticles","Could not obtain the Loader object !"); 
517      return ;
518    } 
519
520   TString taskName(GetName()) ; 
521   taskName.Remove(taskName.Index(Version())-1) ;
522   TClonesArray * recParticles = gime->RecParticles() ; 
523   
524   TString message ; 
525   message  = "event %d\n" ; 
526   message += "       found %d RecParticles\n" ; 
527   Info("PrintRecParticles", message.Data(), gAlice->GetEvNumber(), recParticles->GetEntriesFast() ) ;   
528
529   if(strstr(option,"all")) {  // printing found TS
530     Info("PrintRecParticles","  PARTICLE   Index    \n"  ) ; 
531    
532     Int_t index ;
533     for (index = 0 ; index < recParticles->GetEntries() ; index++) {
534       AliPHOSRecParticle * rp = (AliPHOSRecParticle * ) recParticles->At(index) ;       
535       
536       Text_t particle[11];
537       switch(rp->GetType()) {
538       case  AliPHOSFastRecParticle::kNEUTRALEMFAST:
539         strcpy( particle, "NEUTRAL EM FAST");
540         break;
541       case  AliPHOSFastRecParticle::kNEUTRALHAFAST:
542         strcpy(particle, "NEUTRAL HA FAST");
543         break;
544       case  AliPHOSFastRecParticle::kNEUTRALEMSLOW:
545         strcpy(particle, "NEUTRAL EM SLOW");
546         break ;
547       case  AliPHOSFastRecParticle::kNEUTRALHASLOW: 
548         strcpy(particle, "NEUTRAL HA SLOW");
549         break ;
550       case  AliPHOSFastRecParticle::kCHARGEDEMFAST:
551         strcpy(particle, "CHARGED EM FAST") ;
552         break ;
553       case  AliPHOSFastRecParticle::kCHARGEDHAFAST:
554         strcpy(particle, "CHARGED HA FAST") ;
555         break ; 
556       case  AliPHOSFastRecParticle::kCHARGEDEMSLOW:
557         strcpy(particle, "CHARGEDEMSLOW") ;
558         break ;
559       case  AliPHOSFastRecParticle::kCHARGEDHASLOW:
560         strcpy(particle, "CHARGED HA SLOW") ;
561         break ; 
562       }
563       
564       //    Int_t * primaries; 
565       //    Int_t nprimaries;
566       //    primaries = rp->GetPrimaries(nprimaries);
567       
568       Info("PrintRecParticles", "          %s     %d\n",  particle, rp->GetIndexInList()) ;
569     }
570   }  
571 }
572
573
574