comments added
[u/mrichter/AliRoot.git] / PHOS / AliPHOSPIDv1.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 v1 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] AliPHOSPIDv1 * p1 = new AliPHOSPIDv1("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] AliPHOSPIDv1 * p2 = new AliPHOSPIDv1("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 //            Complitely 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 #include <iostream.h>
72 #include <iomanip.h>
73
74 // --- AliRoot header files ---
75
76 #include "AliRun.h"
77 #include "AliGenerator.h"
78 #include "AliPHOS.h"
79 #include "AliPHOSPIDv1.h"
80 #include "AliPHOSClusterizerv1.h"
81 #include "AliPHOSTrackSegment.h"
82 #include "AliPHOSTrackSegmentMakerv1.h"
83 #include "AliPHOSRecParticle.h"
84
85 ClassImp( AliPHOSPIDv1) 
86
87 //____________________________________________________________________________
88 AliPHOSPIDv1::AliPHOSPIDv1():AliPHOSPID()
89
90   // default ctor
91   fIsInitialized = kFALSE ;
92 }
93
94 //____________________________________________________________________________
95 AliPHOSPIDv1::AliPHOSPIDv1(const char * headeFile,const char * tsBranchTitle):AliPHOSPID()
96
97   //ctor with the indication on where to look for the track segments
98   fHeaderFileName = headeFile ;
99
100   fTSTitle = tsBranchTitle ;
101
102   SetName("AliPHOSPID") ;
103   SetTitle("version1") ;
104
105   TFile * file = (TFile*) gROOT->GetFile(fHeaderFileName.Data() ) ;
106
107   if(file == 0){
108     if(fHeaderFileName.Contains("rfio")) // if we read file using HPSS
109       file = TFile::Open(fHeaderFileName.Data(),"update") ;
110     else
111       file = new TFile(fHeaderFileName.Data(),"update") ;
112     gAlice = (AliRun *) file->Get("gAlice") ;
113   }
114   
115   AliPHOS * phos = (AliPHOS *) gAlice->GetDetector("PHOS") ;    
116   fGeom  = AliPHOSGeometry::GetInstance(phos->GetGeometry()->GetName(),phos->GetGeometry()->GetTitle() );
117   
118   fTrackSegments = new TClonesArray("AliPHOSTrackSegment",1) ;
119   fTSMaker       = 0 ;
120   fEmcRecPoints  = new TObjArray(1) ;
121   fCpvRecPoints  = new TObjArray(1) ;
122   fClusterizer   = 0 ;
123   fRecParticles  = new TClonesArray("AliPHOSRecParticle",100) ;
124
125   fFormula = new TFormula("LambdaCuts","(x>1)*(x<3)*(y>0)*(y<x)") ;
126   
127   // add Task to //root/Tasks folder
128   TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
129   roottasks->Add(this) ; 
130
131   fDispersion = 2.0; 
132   fCpvEmcDistance = 3.0 ;
133   fIsInitialized = kTRUE ;
134
135 }
136 //____________________________________________________________________________
137 AliPHOSPIDv1::~AliPHOSPIDv1()
138
139   //dtor 
140 }
141 //____________________________________________________________________________
142 void AliPHOSPIDv1::Init()
143 {
144   // Make all memory allocations that are not possible in default constructor
145   if(!fIsInitialized){
146     if(fHeaderFileName.IsNull())
147       fHeaderFileName = "galice.root" ;
148     
149     TFile * file = (TFile*) gROOT->GetFile(fHeaderFileName.Data() ) ;
150
151     if(file == 0){
152       file = new TFile(fHeaderFileName.Data(),"update") ;
153       gAlice = (AliRun *) file->Get("gAlice") ;
154     }
155
156     AliPHOS * phos = (AliPHOS *) gAlice->GetDetector("PHOS") ;    
157     fGeom  = AliPHOSGeometry::GetInstance(phos->GetGeometry()->GetName(),phos->GetGeometry()->GetTitle() );
158
159     fTrackSegments = new TClonesArray("AliPHOSTrackSegment",1) ;
160     fTSMaker       = new AliPHOSTrackSegmentMakerv1() ;
161     fEmcRecPoints  = new TObjArray(1) ;
162     fCpvRecPoints  = new TObjArray(1) ;
163     fClusterizer   = new AliPHOSClusterizerv1() ;
164     fRecParticles  = new TClonesArray("AliPHOSRecParticle",100) ;
165     
166     fFormula = new TFormula("LambdaCuts","(x>1)*(x<2.5)*(y>0)*(y<x)") ;
167     
168     // add Task to //root/Tasks folder
169     TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ; 
170     roottasks->Add(this) ; 
171
172     fDispersion = 2.0; 
173     fCpvEmcDistance = 3.0 ;
174
175     fIsInitialized = kTRUE ;
176   }
177 }
178 //____________________________________________________________________________
179 Bool_t AliPHOSPIDv1::ReadTrackSegments()
180 {
181   // Reads TrackSegments an extracts the title of the RecPoints 
182   // branch from which TS were made of.
183   // Then reads both TrackSegments and RecPoints.
184
185   //Fist read Track Segment Branch and extract RecPointsBranch from fTSMaker
186   fTrackSegments->Clear() ; 
187   fEmcRecPoints->Clear() ;
188   fCpvRecPoints->Clear() ;
189   fRecParticles->Clear() ;
190
191   gAlice->GetEvent(fNEvent) ;
192
193   TTree * treeR = gAlice->TreeR()  ; 
194
195   if(treeR==0){
196     char treeName[20]; 
197     sprintf(treeName,"TreeR%d",fNEvent);
198     cout << "Error in AliPHOSClusterizerv1 : no "<<treeName << endl  ;
199     cout << "   Do nothing " << endl ;
200     return kFALSE ;
201   }
202
203   //first read TSMaker branch and extract information about RecPoints Branches
204   TBranch * tsMakerBranch = 0;
205   TBranch * tsBranch = 0;
206
207   TObjArray * branches = treeR->GetListOfBranches() ;
208   Int_t ibranch;
209   Bool_t tsMakerNotFound = kTRUE ;
210   Bool_t tsNotFound = kTRUE ;
211   
212   for(ibranch = 0;(ibranch <branches->GetEntries())&&(tsMakerNotFound||tsNotFound);ibranch++){
213     if(tsMakerNotFound){
214       tsMakerBranch=(TBranch *) branches->At(ibranch) ;
215       if( fTSTitle.CompareTo(tsMakerBranch->GetTitle())==0 )
216         if( strcmp(tsMakerBranch->GetName(),"AliPHOSTrackSegmentMaker") == 0) 
217           tsMakerNotFound = kFALSE ;
218     }
219     if(tsNotFound){
220       tsBranch=(TBranch *) branches->At(ibranch) ;
221       if( fTSTitle.CompareTo(tsBranch->GetTitle())==0 )
222         if( strcmp(tsBranch->GetName(),"PHOSTS") == 0) 
223           tsNotFound = kFALSE ;
224     }
225   }
226   
227   if(tsMakerNotFound ||tsNotFound ){
228     cout << "Can't find Branch with TrackSegmentMaker and TrackSegments " ;
229     cout << "Do nothing" <<endl  ;
230     return kFALSE ;
231   }
232
233   tsMakerBranch->SetAddress(&fTSMaker) ;
234   tsBranch->SetAddress(&fTrackSegments) ;
235
236   tsMakerBranch->GetEntry(0) ;
237   tsBranch->GetEntry(0) ;
238
239   fRecPointsTitle = fTSMaker->GetRecPointsBranch() ;
240
241   //reading now recponts branches
242   TBranch * emcBranch = 0;
243   TBranch * cpvBranch = 0;
244   TBranch * cluBranch = 0;
245
246   Bool_t emcNotFound = kTRUE ;
247   Bool_t cpvNotFound = kTRUE ;
248   Bool_t cluNotFound = kTRUE ;
249  
250   for(ibranch = 0;(ibranch <branches->GetEntries())&&(emcNotFound||cpvNotFound||cluNotFound);ibranch++){
251     if(emcNotFound){
252       emcBranch=(TBranch *) branches->At(ibranch) ;
253       if( fRecPointsTitle.CompareTo(emcBranch->GetTitle())==0 )
254         if( strcmp(emcBranch->GetName(),"PHOSEmcRP") == 0) 
255           emcNotFound = kFALSE ;
256     }
257     if(cpvNotFound){
258       cpvBranch=(TBranch *) branches->At(ibranch) ;
259       if( fRecPointsTitle.CompareTo(cpvBranch->GetTitle())==0 )
260         if( strcmp(cpvBranch->GetName(),"PHOSCpvRP") == 0) 
261           cpvNotFound = kFALSE ;
262     }
263     if(cluNotFound){
264       cluBranch=(TBranch *) branches->At(ibranch) ;
265       if( fRecPointsTitle.CompareTo(cluBranch->GetTitle())==0 )
266         if( strcmp(cluBranch->GetName(),"AliPHOSClusterizer") == 0) 
267           cluNotFound = kFALSE ;
268     }
269   }
270   
271   if(emcNotFound ||cpvNotFound ||cluNotFound ){
272     cout << "Can't find Branch with RecPoints or AliPHOSClusterizer " ;
273     cout << "Do nothing" <<endl  ;
274     return kFALSE ;
275   }
276   
277   emcBranch->SetAddress(&fEmcRecPoints) ;
278   cpvBranch->SetAddress(&fCpvRecPoints) ;
279   cluBranch->SetAddress(&fClusterizer) ;
280
281   emcBranch->GetEntry(0) ;
282   cpvBranch->GetEntry(0) ;
283   cluBranch->GetEntry(0) ;
284
285   return kTRUE ;
286
287
288
289 }
290 //____________________________________________________________________________
291 Float_t  AliPHOSPIDv1::GetDistance(AliPHOSEmcRecPoint * emc,AliPHOSRecPoint * cpv, Option_t *  Axis)const
292 {
293   // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
294  
295   TVector3 vecEmc ;
296   TVector3 vecCpv ;
297   
298   emc->GetLocalPosition(vecEmc) ;
299   cpv->GetLocalPosition(vecCpv) ; 
300   if(emc->GetPHOSMod() == cpv->GetPHOSMod()){ 
301     
302     // Correct to difference in CPV and EMC position due to different distance to center.
303     // we assume, that particle moves from center
304     Float_t dCPV = fGeom->GetIPtoOuterCoverDistance();
305     Float_t dEMC = fGeom->GetIPtoCrystalSurface() ;
306     dEMC         = dEMC / dCPV ;
307     vecCpv = dEMC * vecCpv  - vecEmc ; 
308     if (Axis == "X") return vecCpv.X();
309     if (Axis == "Y") return vecCpv.Y();
310     if (Axis == "Z") return vecCpv.Z();
311     if (Axis == "R") return vecCpv.Mag();
312   } 
313  
314   return 100000000 ;
315 }
316
317 //____________________________________________________________________________
318 void  AliPHOSPIDv1::Exec(Option_t * option) 
319 {
320   //Steering method
321
322   if(!fIsInitialized) 
323     Init() ;
324
325   if(strstr(option,"tim"))
326     gBenchmark->Start("PHOSPID");
327
328
329   Int_t nEvents = (Int_t) gAlice->TreeE()->GetEntries() ;
330
331   for(fNEvent = 0 ;fNEvent <nEvents; fNEvent++){
332     if(!ReadTrackSegments())
333       return ;
334     MakeRecParticles() ;
335     WriteRecParticles();
336     if(strstr(option,"deb"))
337       PrintRecParticles(option) ;
338   }
339
340   if(strstr(option,"tim")){
341     gBenchmark->Stop("PHOSPID");
342     cout << "AliPHOSPID:" << endl ;
343     cout << "  took " << gBenchmark->GetCpuTime("PHOSPID") << " seconds for PID " 
344          <<  gBenchmark->GetCpuTime("PHOSPID")/nEvents << " seconds per event " << endl ;
345     cout << endl ;
346   }
347
348 }
349 //____________________________________________________________________________
350 void  AliPHOSPIDv1::MakeRecParticles(){
351
352   // Makes a RecParticle out of a TrackSegment
353
354   TIter next(fTrackSegments) ; 
355   AliPHOSTrackSegment * ts ; 
356   Int_t index = 0 ; 
357   AliPHOSRecParticle * rp ; 
358   
359   Bool_t ellips = fIDOptions.Contains("ell",TString::kIgnoreCase ) ;
360   Bool_t disp   = fIDOptions.Contains("dis",TString::kIgnoreCase ) ;
361   
362   while ( (ts = (AliPHOSTrackSegment *)next()) ) {
363     
364     new( (*fRecParticles)[index] ) AliPHOSRecParticle() ;
365     rp = (AliPHOSRecParticle *)fRecParticles->At(index) ; 
366     rp->SetTraskSegment(index) ;
367     
368     AliPHOSEmcRecPoint * emc = 0 ;
369     if(ts->GetEmcIndex()>=0)
370       emc = (AliPHOSEmcRecPoint *) fEmcRecPoints->At(ts->GetEmcIndex()) ;
371     
372     AliPHOSRecPoint    * cpv = 0 ;
373     if(ts->GetCpvIndex()>=0)
374       cpv = (AliPHOSRecPoint *)   fCpvRecPoints->At(ts->GetCpvIndex()) ;
375     
376     AliPHOSRecPoint    * ppsd = 0 ;
377     if(ts->GetPpsdIndex()>=0)
378       ppsd= (AliPHOSRecPoint *)   fCpvRecPoints->At(ts->GetPpsdIndex()) ;
379
380     //set momentum and energy first
381     Float_t    e = emc->GetEnergy() ;
382     TVector3 dir = GetMomentumDirection(emc,cpv,ppsd) ; 
383     dir.SetMag(e) ;
384
385     rp->SetMomentum(dir.X(),dir.Y(),dir.Z(),e) ;
386     rp->SetCalcMass(0);
387
388     //now set type (reconstructed) of the particle    
389     Int_t showerprofile = 0;  // 0 narrow and 1 wide
390     
391     if(ellips){
392       Float_t lambda[2] ;
393       emc->GetElipsAxis(lambda) ;
394       if(fFormula->Eval(lambda[0],lambda[1]) <= 0 )
395         showerprofile = 1 ;  // not narrow
396     }
397     
398     if(disp)
399       if(emc->GetDispersion() > fDispersion )
400         showerprofile = 1 ;  // not narrow
401     
402     
403     // Looking at the photon conversion detector
404     Int_t pcdetector= 0 ;  //1 hit and 0 no hit
405     if(ppsd)
406       if(GetDistance(emc, ppsd, "R") < fCpvEmcDistance) 
407         pcdetector = 1 ;  
408     
409     // Looking at the CPV detector
410     Int_t cpvdetector= 0 ;  //1 hit and 0 no hit     
411     if(cpv)
412       if(GetDistance(emc, cpv,  "R") < fCpvEmcDistance) 
413         cpvdetector = 1 ;  
414     
415     Int_t type = showerprofile + 2 * pcdetector + 4 * cpvdetector ;
416     rp->SetType(type) ; 
417     index++ ; 
418   }
419   
420 }
421
422 //____________________________________________________________________________
423 void  AliPHOSPIDv1:: Print(Option_t * option) const
424 {
425   // Print the parameters used for the particle type identification
426     cout <<  "=============== AliPHOSPID1 ================" << endl ;
427     cout <<  "Making PID "<< endl ;
428     cout <<  "    Headers file:               " << fHeaderFileName.Data() << endl ;
429     cout <<  "    RecPoints branch title:     " << fRecPointsTitle.Data() << endl ;
430     cout <<  "    TrackSegments Branch title: " << fTSTitle.Data() << endl ;
431     cout <<  "    RecParticles Branch title   " << fRecparticlesTitle.Data() << endl;
432     cout <<  "with parameters: " << endl ;
433     cout <<  "    Maximal EMC - CPV (PPSD) distance (cm) " << fCpvEmcDistance << endl ;
434     if(fIDOptions.Contains("dis",TString::kIgnoreCase ))
435       cout <<  "                            dispersion cut " << fDispersion << endl ;
436     if(fIDOptions.Contains("ell",TString::kIgnoreCase )){
437       cout << "Eliptic cuts function: " << endl ;
438       cout << fFormula->GetTitle() << endl ;
439     }
440     cout <<  "============================================" << endl ;
441 }
442
443 //____________________________________________________________________________
444 void  AliPHOSPIDv1::SetShowerProfileCut(char * formula)
445 {
446   //set shape of the cut on the axis of ellipce, drown around shouer
447   //shower considered "narrow" if Formula(lambda[0],lambda[1]) > 0.
448   if(fFormula) 
449     delete fFormula; 
450   fFormula = new TFormula("Lambda Cut",formula) ;
451 }
452 //____________________________________________________________________________
453 void  AliPHOSPIDv1::WriteRecParticles()
454 {
455   //check, if these branches already exist  
456   TBranch * pidBranch = 0;
457   TBranch * rpBranch = 0;
458
459   TObjArray * branches = gAlice->TreeR()->GetListOfBranches() ;
460   Int_t ibranch;
461   Bool_t pidNotFound = kTRUE ;
462   Bool_t rpNotFound = kTRUE ;
463   
464   for(ibranch = 0;(ibranch <branches->GetEntries())&& pidNotFound && rpNotFound;ibranch++){
465     if(pidNotFound){
466       pidBranch=(TBranch *) branches->At(ibranch) ;
467       if( (strcmp(pidBranch->GetName(),"PHOSPID") == 0) &&
468           (fRecparticlesTitle.CompareTo(pidBranch->GetTitle()) == 0) )
469         pidNotFound = kFALSE ;
470     }
471     if(rpNotFound){
472       rpBranch=(TBranch *) branches->At(ibranch) ;
473       if( (strcmp(rpBranch->GetName(),"PHOSRP") == 0) &&
474           (fRecparticlesTitle.CompareTo(rpBranch->GetTitle())==0 ))
475         rpNotFound = kFALSE ;
476     }
477   }
478   
479   if(!pidNotFound || !rpNotFound) {
480     cout << "AliPHOSPIDv1 error: " << endl ;
481     cout << "       Branch PHOSRP and PHOSPID with title '"<<fRecparticlesTitle.Data()<<"' already exist "<< endl ;
482     cout << "       can not overwrite " << endl ;
483     return ;
484   }
485
486   //Make branch in TreeR for TrackSegments 
487   char * filename = 0;
488   if(gSystem->Getenv("CONFIG_SPLIT_FILE")!=0){   //generating file name
489     filename = new char[strlen(gAlice->GetBaseFile())+20] ;
490     sprintf(filename,"%s/PHOS.Reco.root",gAlice->GetBaseFile()) ; 
491   }
492
493   TDirectory *cwd = gDirectory;
494   
495   //First rp
496   Int_t bufferSize = 32000 ;    
497   rpBranch = gAlice->TreeR()->Branch("PHOSRP",&fRecParticles,bufferSize);
498   rpBranch->SetTitle(fRecparticlesTitle.Data());
499   if (filename) {
500     rpBranch->SetFile(filename);
501     TIter next( rpBranch->GetListOfBranches());
502     TBranch * sb ;
503     while ((sb=(TBranch*)next())) {
504       sb->SetFile(filename);
505     }   
506     cwd->cd();
507   }
508
509   //second, pid
510   Int_t splitlevel = 0 ; 
511   AliPHOSPIDv1 * pid = this ;
512   pidBranch = gAlice->TreeR()->Branch("AliPHOSPID","AliPHOSPIDv1",&pid,bufferSize,splitlevel);
513   pidBranch->SetTitle(fRecparticlesTitle.Data());
514   if (filename) {
515     pidBranch->SetFile(filename);
516     TIter next( pidBranch->GetListOfBranches());
517     TBranch * sb ;
518     while ((sb=(TBranch*)next())) {
519       sb->SetFile(filename);
520     }   
521     cwd->cd();
522   }    
523   
524   rpBranch->Fill() ;
525   pidBranch->Fill() ;
526   gAlice->TreeR()->Fill() ;  
527   gAlice->TreeR()->Write(0,kOverwrite) ;  
528   
529 }
530 //____________________________________________________________________________
531 void  AliPHOSPIDv1::PlotDispersionCuts()const
532 {
533   // produces a plot of the dispersion cut
534   TCanvas*  lambdas = new TCanvas("lambdas","Cuts on the ellipse axis",200,10,700,500);
535   
536   if(fIDOptions.Contains("ell",TString::kIgnoreCase ) ){
537     TF2 * ell = new TF2("Elliptic Cuts",fFormula->GetName(),0,3,0,3) ;
538     ell->SetMinimum(0.0000001) ;
539     ell->SetMaximum(0.001) ;
540     ell->SetLineStyle(1) ;
541     ell->SetLineWidth(2) ;
542     ell->Draw() ;
543   }
544   
545   if( fIDOptions.Contains("dis",TString::kIgnoreCase ) ){
546     TF2 * dsp = new TF2("dispersion","(y<x)*(x*x+y*y < [0]*[0])",0,3,0,3) ;
547     dsp->SetParameter(0,fDispersion) ;
548     dsp->SetMinimum(0.0000001) ;
549     dsp->SetMaximum(0.001) ;
550     dsp->SetLineStyle(1) ;
551     dsp->SetLineColor(2) ;
552     dsp->SetLineWidth(2) ;
553     dsp->SetNpx(200) ;
554     dsp->SetNpy(200) ;
555     if(fIDOptions.Contains("ell",TString::kIgnoreCase ) )
556       dsp->Draw("same") ;
557     else
558       dsp->Draw() ;
559   }
560   lambdas->Update();
561 }
562
563 //____________________________________________________________________________
564 TVector3 AliPHOSPIDv1::GetMomentumDirection(AliPHOSEmcRecPoint * emc, AliPHOSRecPoint * cpv,AliPHOSRecPoint * ppsd)const 
565
566   // Calculates the momentum direction:
567   //   1. if only a EMC RecPoint, direction is given by IP and this RecPoint
568   //   2. if a EMC RecPoint and one PPSD RecPoint, direction is given by the line through the 2 recpoints 
569   //   3. if a EMC RecPoint and two PPSD RecPoints, dirrection is given by the average line through 
570   //      the 2 pairs of recpoints  
571   // However because of the poor position resolution of PPSD the direction is always taken as if we were 
572   //  in case 1.
573
574   TVector3 dir(0,0,0) ; 
575   
576   TVector3 emcglobalpos ;
577   TMatrix  dummy ;
578   
579   emc->GetGlobalPosition(emcglobalpos, dummy) ;
580   
581  
582   // The following commented code becomes valid once the PPSD provides 
583   // a reasonable position resolution, at least as good as EMC ! 
584   //   TVector3 ppsdlglobalpos ;
585   //   TVector3 ppsduglobalpos ;
586   //   if( fPpsdLowRecPoint ){ // certainly a photon that has concerted
587   //     fPpsdLowRecPoint->GetGlobalPosition(ppsdlglobalpos, mdummy) ; 
588   //     dir = emcglobalpos -  ppsdlglobalpos ; 
589   //     if( fPpsdUpRecPoint ){ // not looks like a charged       
590   //        fPpsdUpRecPoint->GetGlobalPosition(ppsduglobalpos, mdummy) ; 
591   //        dir = ( dir +  emcglobalpos -  ppsduglobalpos ) * 0.5 ; 
592   //      }
593   //   }
594   //   else { // looks like a neutral
595   //    dir = emcglobalpos ;  
596   //  }
597   
598   dir = emcglobalpos ;  
599   dir.SetZ( -dir.Z() ) ;   // why ?  
600   dir.SetMag(1.) ;
601
602   //account correction to the position of IP
603   Float_t xo,yo,zo ; //Coordinates of the origin
604   gAlice->Generator()->GetOrigin(xo,yo,zo) ;
605   TVector3 origin(xo,yo,zo);
606   dir = dir - origin ;
607
608   return dir ;  
609 }
610 //____________________________________________________________________________
611 void AliPHOSPIDv1::PrintRecParticles(Option_t * option)
612 {
613   // Prints the list of reconstructed particles
614   cout << "AliPHOSPIDv1: " << endl ;
615   cout << "       found " << fRecParticles->GetEntriesFast() << " RecParticles " << endl ;
616
617   if(strstr(option,"all")) {  // printing found TS
618     
619     cout << "  PARTICLE "   
620          << "  Index    "  << endl ;
621       //         << "  X        "     
622       //         << "  Y        " 
623       //         << "  Z        "    
624       //         << " # of primaries "          
625       //         << " Primaries list "    <<  endl;      
626     
627     Int_t index ;
628     for (index = 0 ; index < fRecParticles->GetEntries() ; index++) {
629       AliPHOSRecParticle * rp = (AliPHOSRecParticle * ) fRecParticles->At(index) ;       
630       
631       Text_t particle[11];
632       switch(rp->GetType()) {
633       case  AliPHOSFastRecParticle::kNEUTRALEM:
634         strcpy( particle, "NEUTRAL_EM");
635         break;
636       case  AliPHOSFastRecParticle::kNEUTRALHA:
637         strcpy(particle, "NEUTRAL_HA");
638         break;
639       case  AliPHOSFastRecParticle::kGAMMA:
640         strcpy(particle, "GAMMA");
641         break ;
642       case  AliPHOSFastRecParticle::kGAMMAHA: 
643         strcpy(particle, "GAMMA_H");
644         break ;
645       case  AliPHOSFastRecParticle::kABSURDEM:
646         strcpy(particle, "ABSURD_EM") ;
647         break ;
648       case  AliPHOSFastRecParticle::kABSURDHA:
649         strcpy(particle, "ABSURD_HA") ;
650         break ; 
651       case  AliPHOSFastRecParticle::kELECTRON:
652         strcpy(particle, "ELECTRON") ;
653         break ;
654       case  AliPHOSFastRecParticle::kCHARGEDHA:
655         strcpy(particle, "CHARGED_HA") ;
656         break ; 
657       }
658       
659       //    Int_t * primaries; 
660       //    Int_t nprimaries;
661       //    primaries = rp->GetPrimaries(nprimaries);
662       
663       cout << setw(15) << particle << "  "
664            << setw(3) <<  rp->GetIndexInList() << " "  ;
665         //         << setw(4) <<  nprimaries << "  ";
666         //      for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
667         //      cout << setw(4)  <<  primaries[iprimary] << " ";
668       cout << endl;      
669     }
670     cout << "-------------------------------------------" << endl ;
671   }
672   
673 }
674
675
676