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