]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSPIDv0.cxx
updating
[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 "TROOT.h"
61 #include "TTree.h"
62 #include "TFile.h"
63 #include "TF2.h"
64 #include "TFormula.h"
65 #include "TCanvas.h"
66 #include "TFolder.h"
67 #include "TSystem.h"
68 #include "TBenchmark.h"
69 // --- Standard library ---
70
71 // --- AliRoot header files ---
72
73 #include "AliRun.h"
74 #include "AliGenerator.h"
75 #include "AliPHOS.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"
84
85 ClassImp( AliPHOSPIDv0) 
86
87 //____________________________________________________________________________
88 AliPHOSPIDv0::AliPHOSPIDv0():AliPHOSPID()
89
90   // default ctor
91   fFormula           = 0 ;
92   fDispersion        = 0. ; 
93   fCpvEmcDistance    = 0 ; 
94   fTimeGate          = 2.e-9 ;
95   fEventFolderName    = "" ; 
96   fTrackSegmentsTitle= "" ; 
97   fRecPointsTitle    = "" ; 
98   fRecParticlesTitle = "" ; 
99   fIDOptions         = "dis time" ; 
100   fRecParticlesInRun = 0 ;
101   fClusterizer = 0;
102   fTSMaker = 0;
103 }
104
105 //____________________________________________________________________________
106 AliPHOSPIDv0::AliPHOSPIDv0(const char * evFolderName,const char * name) : AliPHOSPID(evFolderName, name)
107
108   //ctor with the indication on where to look for the track segments
109
110   fFormula        = new TFormula("LambdaCuts","(x>1)*(x<2.5)*(y>0)*(y<x)") ;   
111   fDispersion     = 2.0 ; 
112   fCpvEmcDistance = 3.0 ;
113   fTimeGate          = 2.e-9 ;
114  
115   fEventFolderName     = GetTitle() ; 
116   fTrackSegmentsTitle = GetName();
117   fRecPointsTitle     = GetName();
118   fRecParticlesTitle  = GetName();
119   fIDOptions          = "dis time" ;
120     
121   fRecParticlesInRun = 0 ; 
122
123   Init() ;
124
125 }
126
127 //____________________________________________________________________________
128 AliPHOSPIDv0::~AliPHOSPIDv0()
129
130 }
131
132 //____________________________________________________________________________
133 Float_t  AliPHOSPIDv0::GetDistance(AliPHOSEmcRecPoint * emc,AliPHOSRecPoint * cpv, Option_t *  Axis)const
134 {
135   // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
136  
137   const AliPHOSGeometry * geom = AliPHOSLoader::GetPHOSGeometry() ; 
138   TVector3 vecEmc ;
139   TVector3 vecCpv ;
140   
141   emc->GetLocalPosition(vecEmc) ;
142   cpv->GetLocalPosition(vecCpv) ; 
143   if(emc->GetPHOSMod() == cpv->GetPHOSMod()){ 
144     
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() ;
149     dEMC         = dEMC / dCPV ;
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();
155   } 
156  
157   return 100000000 ;
158 }
159
160 //____________________________________________________________________________
161 void  AliPHOSPIDv0::Exec(Option_t * option) 
162 {
163   //Steering method
164   
165   if( strcmp(GetName(), "")== 0 ) 
166     Init() ;
167   
168   if(strstr(option,"tim"))
169     gBenchmark->Start("PHOSPID");
170   
171   if(strstr(option,"print")) {
172     Print() ; 
173     return ; 
174   }
175
176   AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
177   if(runget == 0x0) 
178    {
179      Error("Exec","Can not find run getter in event folder \"%s\"",GetTitle());
180      return;
181    }
182   
183   AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
184   if ( gime == 0 ) 
185    {
186      Error("Exec","Could not obtain the Loader object !"); 
187      return ;
188    } 
189
190   if(gime->BranchExists("RecParticles") )
191     return ;
192
193   
194   Int_t nevents = runget->GetNumberOfEvents() ;       //(Int_t) gAlice->TreeE()->GetEntries() ;
195
196   Int_t ievent ;
197   
198   for(ievent = 0; ievent < nevents; ievent++){
199     runget->GetEvent(ievent);
200     Info("Exec", "event %d %d %d", ievent, gime->EmcRecPoints(), gime->TrackSegments()) ;
201     MakeRecParticles() ;
202     
203     WriteRecParticles();
204     
205     if(strstr(option,"deb"))
206       PrintRecParticles(option) ;
207
208     //increment the total number of rec particles per run 
209     fRecParticlesInRun += gime->RecParticles()->GetEntriesFast() ; 
210
211   }
212   
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) ; 
217   }
218   
219 }
220 //____________________________________________________________________________
221 void AliPHOSPIDv0::Init()
222 {
223   // Make all memory allocations that are not possible in default constructor
224   // Add the PID task to the list of PHOS tasks
225   
226   AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
227   if(runget == 0x0) 
228    {
229      Error("Exec","Can not find run getter in event folder \"%s\"",GetTitle());
230      return;
231    }
232   
233   AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
234   if ( gime == 0 ) 
235    {
236      Error("Exec","Could not obtain the Loader object !"); 
237      return ;
238    } 
239    
240   gime->PostPID(this);
241   gime->LoadRecParticles("UPDATE");
242   
243 }
244
245 //____________________________________________________________________________
246 void  AliPHOSPIDv0::MakeRecParticles(){
247
248   TString taskName(GetName()) ; 
249   taskName.Remove(taskName.Index(Version())-1) ;
250
251   AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
252   if(runget == 0x0) 
253    {
254      Error("Exec","Can not find run getter in event folder \"%s\"",GetTitle());
255      return;
256    }
257   
258   AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
259   if ( gime == 0 ) 
260    {
261      Error("Exec","Could not obtain the Loader object !"); 
262      return ;
263    } 
264
265   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
266   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
267   TClonesArray * trackSegments = gime->TrackSegments() ; 
268   TClonesArray * recParticles  = gime->RecParticles() ; 
269
270   recParticles->Clear();
271   
272   TIter next(trackSegments) ; 
273   AliPHOSTrackSegment * ts ; 
274   Int_t index = 0 ; 
275   AliPHOSRecParticle * rp ; 
276   
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 ) ;
280   
281   while ( (ts = (AliPHOSTrackSegment *)next()) ) {
282     
283     new( (*recParticles)[index] ) AliPHOSRecParticle() ;
284     rp = (AliPHOSRecParticle *)recParticles->At(index) ; 
285     rp->SetTrackSegment(index) ;
286     rp->SetIndexInList(index) ;
287     
288     AliPHOSEmcRecPoint * emc = 0 ;
289     if(ts->GetEmcIndex()>=0)
290       emc = (AliPHOSEmcRecPoint *) emcRecPoints->At(ts->GetEmcIndex()) ;
291     
292     AliPHOSRecPoint    * cpv = 0 ;
293     if(ts->GetCpvIndex()>=0)
294       cpv = (AliPHOSRecPoint *)   cpvRecPoints->At(ts->GetCpvIndex()) ;
295     
296     //set momentum and energy first
297     Float_t    e = emc->GetEnergy() ;
298     TVector3 dir = GetMomentumDirection(emc,cpv) ; 
299     dir.SetMag(e) ;
300
301     rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),e) ;
302     rp->SetCalcMass(0);
303
304     //now set type (reconstructed) of the particle    
305     Int_t showerprofile = 0;  // 0 narrow and 1 wide
306     
307     if(ellips){
308       Float_t lambda[2] ;
309       emc->GetElipsAxis(lambda) ;
310       if(fFormula->Eval(lambda[0],lambda[1]) <= 0 )
311         showerprofile = 1 ;  // not narrow
312     }
313     
314     if(disp)
315       if(emc->GetDispersion() > fDispersion )
316         showerprofile = 1 ;  // not narrow
317     
318     Int_t slow = 0 ;
319     if(time)
320       if(emc->GetTime() > fTimeGate )
321         slow = 0 ; 
322         
323     // Looking at the CPV detector
324     Int_t cpvdetector= 0 ;  //1 hit and 0 no hit     
325     if(cpv)
326       if(GetDistance(emc, cpv,  "R") < fCpvEmcDistance) 
327         cpvdetector = 1 ;  
328     
329     Int_t type = showerprofile + 2 * slow  + 4 * cpvdetector ;
330     rp->SetType(type) ; 
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);
337     index++ ; 
338   }
339   
340 }
341
342 //____________________________________________________________________________
343 void  AliPHOSPIDv0:: Print() const
344 {
345   // Print the parameters used for the particle type identification
346   TString message ; 
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(),  
356        GetTitle(), 
357        fRecPointsTitle.Data(), 
358        fTrackSegmentsTitle.Data(), 
359        fRecParticlesTitle.Data(), 
360        fCpvEmcDistance );
361
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) ;
368 }
369
370 //____________________________________________________________________________
371 void  AliPHOSPIDv0::SetShowerProfileCut(const char * formula)
372 {
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.
375   if(fFormula) 
376     delete fFormula; 
377   fFormula = new TFormula("Lambda Cut",formula) ;
378 }
379 //____________________________________________________________________________
380 void  AliPHOSPIDv0::WriteRecParticles()
381 {
382  
383   AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
384   if(runget == 0x0) 
385    {
386      Error("Exec","Can not find run getter in event folder \"%s\"",GetTitle());
387      return;
388    }
389   
390   AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
391   if ( gime == 0 ) 
392    {
393      Error("Exec","Could not obtain the Loader object !"); 
394      return ;
395    } 
396
397   TClonesArray * recParticles = gime->RecParticles() ; 
398   recParticles->Expand(recParticles->GetEntriesFast() ) ;
399
400   TTree * treeR = gime->TreeR();
401   
402   if(!treeR){
403     gime->MakeTree("R");
404     treeR = gime->TreeR() ;
405   }
406   
407   //First rp
408   Int_t bufferSize = 32000 ;    
409   TBranch * rpBranch = treeR->Branch("PHOSRP",&recParticles,bufferSize);
410   rpBranch->SetTitle(fRecParticlesTitle);
411   
412   //second, pid
413   Int_t splitlevel = 0 ; 
414   AliPHOSPIDv0 * pid = this ;
415   TBranch * pidBranch = treeR->Branch("AliPHOSPID","AliPHOSPIDv0",&pid,bufferSize,splitlevel);
416   pidBranch->SetTitle(fRecParticlesTitle.Data());
417   
418   rpBranch->Fill() ;
419   pidBranch->Fill() ;
420
421   gime->WriteRecParticles("OVERWRITE");
422   gime->WritePID("OVERWRITE");
423
424 }
425
426 //____________________________________________________________________________
427 void  AliPHOSPIDv0::PlotDispersionCuts()const
428 {
429   // produces a plot of the dispersion cut
430   TCanvas*  lambdas = new TCanvas("lambdas","Cuts on the ellipse axis",200,10,700,500);
431  
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) ;
438     ell->Draw() ;
439   }
440   
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) ;
449     dsp->SetNpx(200) ;
450     dsp->SetNpy(200) ;
451     if(fIDOptions.Contains("ell",TString::kIgnoreCase ) )
452       dsp->Draw("same") ;
453     else
454       dsp->Draw() ;
455   }
456   lambdas->Update();
457 }
458
459 //____________________________________________________________________________
460 TVector3 AliPHOSPIDv0::GetMomentumDirection(AliPHOSEmcRecPoint * emc, AliPHOSRecPoint * )const 
461
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 
466   //  in case 1.
467
468   TVector3 dir(0,0,0) ; 
469   
470   TVector3 emcglobalpos ;
471   TMatrix  dummy ;
472   
473   emc->GetGlobalPosition(emcglobalpos, dummy) ;
474   
475  
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 ; 
486   //      }
487   //   }
488   //   else { // looks like a neutral
489   //    dir = emcglobalpos ;  
490   //  }
491   
492   dir = emcglobalpos ;  
493   dir.SetZ( -dir.Z() ) ;   // why ?  
494   dir.SetMag(1.) ;
495
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);
500   dir = dir - origin ;
501
502   return dir ;  
503 }
504 //____________________________________________________________________________
505 void AliPHOSPIDv0::PrintRecParticles(Option_t * option)
506 {
507   // Print table of reconstructed particles
508
509   AliRunLoader* runget = AliRunLoader::GetRunLoader(GetTitle());
510   if(runget == 0x0) 
511    {
512      Error("WriteRecParticles","Can not find run getter in event folder \"%s\"",GetTitle());
513      return;
514    }
515   
516   AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(runget->GetLoader("PHOSLoader"));
517   if ( gime == 0 ) 
518    {
519      Error("WriteRecParticles","Could not obtain the Loader object !"); 
520      return ;
521    } 
522
523   TString taskName(GetName()) ; 
524   taskName.Remove(taskName.Index(Version())-1) ;
525   TClonesArray * recParticles = gime->RecParticles() ; 
526   
527   TString message ; 
528   message  = "event %d\n" ; 
529   message += "       found %d RecParticles\n" ; 
530   Info("PrintRecParticles", message.Data(), gAlice->GetEvNumber(), recParticles->GetEntriesFast() ) ;   
531
532   if(strstr(option,"all")) {  // printing found TS
533     Info("PrintRecParticles","  PARTICLE   Index    \n"  ) ; 
534    
535     Int_t index ;
536     for (index = 0 ; index < recParticles->GetEntries() ; index++) {
537       AliPHOSRecParticle * rp = (AliPHOSRecParticle * ) recParticles->At(index) ;       
538       
539       Text_t particle[11];
540       switch(rp->GetType()) {
541       case  AliPHOSFastRecParticle::kNEUTRALEMFAST:
542         strcpy( particle, "NEUTRAL EM FAST");
543         break;
544       case  AliPHOSFastRecParticle::kNEUTRALHAFAST:
545         strcpy(particle, "NEUTRAL HA FAST");
546         break;
547       case  AliPHOSFastRecParticle::kNEUTRALEMSLOW:
548         strcpy(particle, "NEUTRAL EM SLOW");
549         break ;
550       case  AliPHOSFastRecParticle::kNEUTRALHASLOW: 
551         strcpy(particle, "NEUTRAL HA SLOW");
552         break ;
553       case  AliPHOSFastRecParticle::kCHARGEDEMFAST:
554         strcpy(particle, "CHARGED EM FAST") ;
555         break ;
556       case  AliPHOSFastRecParticle::kCHARGEDHAFAST:
557         strcpy(particle, "CHARGED HA FAST") ;
558         break ; 
559       case  AliPHOSFastRecParticle::kCHARGEDEMSLOW:
560         strcpy(particle, "CHARGEDEMSLOW") ;
561         break ;
562       case  AliPHOSFastRecParticle::kCHARGEDHASLOW:
563         strcpy(particle, "CHARGED HA SLOW") ;
564         break ; 
565       }
566       
567       //    Int_t * primaries; 
568       //    Int_t nprimaries;
569       //    primaries = rp->GetPrimaries(nprimaries);
570       
571       Info("PrintRecParticles", "          %s     %d\n",  particle, rp->GetIndexInList()) ;
572     }
573   }  
574 }
575
576
577