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