]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSTrackSegmentMakerv2.cxx
Runloader is updated when moving to next file (quick fix).
[u/mrichter/AliRoot.git] / PHOS / AliPHOSTrackSegmentMakerv2.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 /* $Id$ */
16
17 /* History of cvs commits:
18  *
19  * $Log$
20  *
21  */
22
23 //_________________________________________________________________________
24 // Implementation version 2 of algorithm class to construct PHOS track segments
25 // Track segment for PHOS is list of 
26 //        EMC RecPoint + (possibly) projection of TPC track
27 // To find TrackSegments we do the following: 
28 //  for each EMC RecPoints we look at
29 //   TPC projections radius fRtpc. 
30 //  If there is such a track
31 //   we make "Link" it is just indexes of EMC and TPC track and distance
32 //   between them in the PHOS plane. 
33 //  Then we sort "Links" and starting from the 
34 //   least "Link" pointing to the unassined EMC and TPC assing them to 
35 //   new TrackSegment. 
36 // If there is no TPC track we make TrackSegment 
37 // consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
38 //// In principle this class should be called from AliPHOSReconstructor, but 
39 // one can use it as well in standalone mode.
40 // Use  case:
41 //  root [0] AliPHOSTrackSegmentMakerv2 * t = new AliPHOSTrackSegmentMaker("galice.root", "tracksegmentsname", "recpointsname")
42 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
43 //               // reads gAlice from header file "galice.root", uses recpoints stored in the branch names "recpointsname" (default = "Default")
44 //               // and saves recpoints in branch named "tracksegmentsname" (default = "recpointsname")                       
45 //  root [1] t->ExecuteTask()
46 //  root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
47 //  root [4] t->ExecuteTask("deb all time") 
48 //                 
49 //*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH) & Yves Schutz (SUBATECH) 
50 //
51
52 // --- ROOT system ---
53 #include "TFile.h"
54 #include "TTree.h"
55 #include "TBenchmark.h"
56
57 // --- Standard library ---
58 #include "Riostream.h"
59 // --- AliRoot header files ---
60 #include "AliPHOSGeometry.h"
61 #include "AliPHOSTrackSegmentMakerv2.h"
62 #include "AliPHOSTrackSegment.h"
63 #include "AliPHOSLink.h"
64 #include "AliPHOSGetter.h"
65 #include "AliESD.h"
66 #include "AliESDtrack.h"
67
68 ClassImp( AliPHOSTrackSegmentMakerv2) 
69
70
71 //____________________________________________________________________________
72 AliPHOSTrackSegmentMakerv2::AliPHOSTrackSegmentMakerv2() :
73   AliPHOSTrackSegmentMaker(),
74   fDefaultInit(kTRUE),
75   fWrite(kFALSE),
76   fNTrackSegments(0),
77   fRtpc(0.f),
78   fLinkUpArray(0),
79   fEmcFirst(0),
80   fEmcLast(0),
81   fModule(0),
82   fTrackSegmentsInRun(0)
83   
84 {
85   // default ctor (to be used mainly by Streamer)
86   InitParameters() ; 
87   for(Int_t i =0 ; i<5; i++)fTPC[i]=0 ;
88 }
89
90 //____________________________________________________________________________
91 AliPHOSTrackSegmentMakerv2::AliPHOSTrackSegmentMakerv2(const TString & alirunFileName, const TString & eventFolderName) :
92   AliPHOSTrackSegmentMaker(alirunFileName, eventFolderName),
93   fDefaultInit(kFALSE),
94   fWrite(kFALSE),
95   fNTrackSegments(0),
96   fRtpc(0.f),
97   fLinkUpArray(0),
98   fEmcFirst(0),
99   fEmcLast(0),
100   fModule(0),
101   fTrackSegmentsInRun(0)
102 {
103   // ctor
104   InitParameters() ; 
105   Init() ;
106   fESD = 0;
107 }
108
109 //____________________________________________________________________________
110 AliPHOSTrackSegmentMakerv2::AliPHOSTrackSegmentMakerv2(const AliPHOSTrackSegmentMakerv2 & tsm) :
111   AliPHOSTrackSegmentMaker(tsm),
112   fDefaultInit(kFALSE),
113   fWrite(kFALSE),
114   fNTrackSegments(0),
115   fRtpc(0.f),
116   fLinkUpArray(0),
117   fEmcFirst(0),
118   fEmcLast(0),
119   fModule(0),
120   fTrackSegmentsInRun(0)
121 {
122   // cpy ctor: no implementation yet
123   // requested by the Coding Convention
124   Fatal("cpy ctor", "not implemented") ;
125 }
126
127
128 //____________________________________________________________________________
129  AliPHOSTrackSegmentMakerv2::~AliPHOSTrackSegmentMakerv2()
130
131   // dtor
132   // fDefaultInit = kTRUE if TrackSegmentMaker created by default ctor (to get just the parameters)
133   if (!fDefaultInit)  
134     delete fLinkUpArray ;
135   for(Int_t imod=0; imod<5; imod++){
136     if(fTPC[imod]) delete fTPC[imod] ;
137   }
138 }
139
140 //____________________________________________________________________________
141 const TString AliPHOSTrackSegmentMakerv2::BranchName() const 
142 {  
143  
144   return GetName() ;
145 }
146
147 //____________________________________________________________________________
148 void  AliPHOSTrackSegmentMakerv2::FillOneModule()
149 {
150   // Finds first and last indexes between which 
151   // clusters from one PHOS module are
152
153   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
154   const AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
155   
156   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
157  
158   //First EMC clusters
159   Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
160   for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&  
161         ((dynamic_cast<AliPHOSRecPoint *>(emcRecPoints->At(fEmcLast)))->GetPHOSMod() == fModule ); 
162       fEmcLast ++)  ;
163   
164   //Now TPC tracks
165   if(fESD){
166     //Do it ones, only first time
167     if(fModule==1){
168       Int_t nTracks = fESD->GetNumberOfTracks();
169       TVector3 inPHOS ;
170    
171       //In this particular case we use fixed vertex position at zero
172       Double_t vtx[3]={0.,0.,0.} ;
173       AliESDtrack *track;
174       Double_t xyz[3] ;
175       Int_t nPHOSmod = geom->GetNModules() ;
176       for(Int_t imod=0 ; imod< nPHOSmod; imod++){
177         fTPC[imod]->Clear() ;
178       }
179       Double_t rEMC = geom->GetIPtoCrystalSurface() ; //Use here ideal geometry 
180       for (Int_t iTrack=0; iTrack<nTracks; iTrack++) {
181         track = fESD->GetTrack(iTrack);
182         if (!track->GetXYZAt(rEMC, fESD->GetMagneticField(), xyz))
183           continue; //track coord on the cylinder of PHOS radius
184         if ((TMath::Abs(xyz[0])+TMath::Abs(xyz[1])+TMath::Abs(xyz[2]))<=0)
185           continue;
186         //Check if this track hits PHOS
187         inPHOS.SetXYZ(xyz[0],xyz[1],xyz[2]);
188         Int_t modNum ; 
189         Double_t x,z ;
190         geom->ImpactOnEmc(vtx, inPHOS.Theta(), inPHOS.Phi(), modNum, x,z) ;
191         if(modNum>0 && modNum<=nPHOSmod){
192           //Mark this track as one belonging to module
193           fTPC[modNum-1]->AddLast(track) ;
194         }
195       }
196     }
197   } 
198
199 }
200
201 //____________________________________________________________________________
202 void  AliPHOSTrackSegmentMakerv2::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,
203                                                          AliESDtrack *track,
204                                                          Float_t &dx, Float_t &dz) const
205 {
206   // Calculates the distance between the EMC RecPoint and the CPV RecPoint
207   // Clusters are sorted in "rows" and "columns" of width 1 cm
208
209   const AliPHOSGeometry * geom = AliPHOSGeometry::GetInstance() ;
210   TVector3 emcGlobal; // Global position of the CPV recpoint
211   geom->GetGlobal((AliRecPoint*)emcClu,emcGlobal);
212   Double_t rEMC = emcGlobal.Pt() ;// Radius from IP to current point 
213
214 //  printf("EMC: x=%f, y=%f, z=%f \n",emcGlobal.X(),emcGlobal.Y(),emcGlobal.Z()) ;
215
216   Double_t xyz[3] ;
217   if (track->GetXYZAt(rEMC, fESD->GetMagneticField(), xyz)){ //calculate distance
218 // printf("xyz: x=%f, y=%f, z=%f \n",xyz[0],xyz[1],xyz[2]) ;
219     dx=TMath::Sqrt((emcGlobal.X()-xyz[0])*(emcGlobal.X()-xyz[0])+(emcGlobal.Y()-xyz[1])*(emcGlobal.Y()-xyz[1])) ;
220     dx=TMath::Sign(dx,(Float_t)(emcGlobal.X()-xyz[0])) ; //set direction
221     dz=emcGlobal.Z()-xyz[2] ;
222 // printf("  dx=%f, dz=%f \n",dx,dz) ;
223   }
224   else{
225     dx=999. ; 
226     dz=999. ;
227   }
228
229   return ;
230 }
231 //____________________________________________________________________________
232 void  AliPHOSTrackSegmentMakerv2::Init()
233 {
234   // Make all memory allocations that are not possible in default constructor
235   
236   AliPHOSGetter* gime = AliPHOSGetter::Instance();
237   if(!gime)
238     gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
239   
240   fLinkUpArray  = new TClonesArray("AliPHOSLink", 1000); 
241   if ( !gime->TrackSegmentMaker() ) {
242     gime->PostTrackSegmentMaker(this);
243   }
244   AliPHOSGeometry * geom = gime->PHOSGeometry() ;
245   Int_t nMod = geom->GetNModules() ;
246   for(Int_t imod=0; imod<nMod ; imod++){
247    fTPC[imod]=new TList() ;
248   }
249 }
250
251 //____________________________________________________________________________
252 void  AliPHOSTrackSegmentMakerv2::InitParameters()
253 {
254   //Initializes parameters
255   fRtpc      = 4. ;
256   fEmcFirst  = 0 ;    
257   fEmcLast   = 0 ;   
258   fLinkUpArray = 0 ;
259   fWrite                   = kTRUE ;
260   fTrackSegmentsInRun       = 0 ; 
261   SetEventRange(0,-1) ;
262 }
263
264
265 //____________________________________________________________________________
266 void  AliPHOSTrackSegmentMakerv2::MakeLinks()const
267
268   // Finds distances (links) between all EMC and CPV clusters, 
269   // which are not further apart from each other than fRcpv 
270   // and sort them in accordance with this distance
271   
272   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
273   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
274
275   fLinkUpArray->Clear() ;    
276
277   AliPHOSEmcRecPoint * emcclu ;
278
279   Int_t iLinkUp  = 0 ;
280   
281   Int_t iEmcRP;
282   for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
283     emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
284     Int_t mod=emcclu->GetPHOSMod() ;
285     TIter next(fTPC[mod-1]) ; 
286     AliESDtrack *track ;
287     Int_t itrack=0 ;
288     while((track= static_cast<AliESDtrack *>(next()))){
289       itrack = track->GetID() ;
290       Float_t dx,dz ;
291       GetDistanceInPHOSPlane(emcclu, track, dx,dz) ;     
292       if(TMath::Sqrt(dx*dx+dz*dz) < fRtpc ){ 
293         new ((*fLinkUpArray)[iLinkUp++])  AliPHOSLink(dx, dz, iEmcRP, itrack, -1) ;
294       }      
295     }
296   } 
297   
298   fLinkUpArray->Sort() ;  //first links with smallest distances
299 }
300
301 //____________________________________________________________________________
302 void  AliPHOSTrackSegmentMakerv2::MakePairs()
303
304   // Using the previously made list of "links", we found the smallest link - i.e. 
305   // link with the least distance between EMC and CPV and pointing to still 
306   // unassigned RecParticles. We assign these RecPoints to TrackSegment and 
307   // remove them from the list of "unassigned". 
308
309   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
310
311   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
312   TClonesArray * trackSegments = gime->TrackSegments();
313     
314   //Make arrays to mark clusters already chosen
315   Int_t * emcExist = 0;
316   if(fEmcLast > fEmcFirst)
317     emcExist = new Int_t[fEmcLast-fEmcFirst] ;
318   
319   Int_t index;
320   for(index = 0; index <fEmcLast-fEmcFirst; index ++)
321     emcExist[index] = 1 ;
322   
323   
324   Int_t * tpcExist = 0;
325   Int_t nTracks = fTPC[fModule-1]->GetSize() ;
326   if(nTracks>0)
327     tpcExist = new Int_t[nTracks] ;
328   
329   for(index = 0; index <nTracks; index ++)
330     tpcExist[index] = 1 ;
331   
332   
333   // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance 
334   TIter nextUp(fLinkUpArray) ;
335   
336   AliPHOSLink * linkUp ;
337   
338   AliPHOSCpvRecPoint * nullpointer = 0 ;
339   
340   while ( (linkUp =  static_cast<AliPHOSLink *>(nextUp()) ) ){  
341
342     if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){
343
344       if(tpcExist[linkUp->GetCpv()]){ //Track still exist
345          Float_t dx,dz ; 
346          linkUp->GetXZ(dx,dz) ;
347          new ((* trackSegments)[fNTrackSegments]) 
348            AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) , 
349                                nullpointer,
350                                linkUp->GetTrack(),dx,dz) ;
351        (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
352        fNTrackSegments++ ;
353        emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc  that Cpv was found 
354        //mark track as already used 
355        tpcExist[linkUp->GetCpv()] = kFALSE ;
356       } //if CpvUp still exist
357     } 
358   }        
359
360   //look through emc recPoints left without CPV
361   if(emcExist){ //if there is emc rec point
362     Int_t iEmcRP ;
363     for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst  ; iEmcRP++ ){
364       if(emcExist[iEmcRP] > 0 ){
365        new ((*trackSegments)[fNTrackSegments])  
366          AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)), 
367                            nullpointer) ;
368        (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
369        fNTrackSegments++;    
370       } 
371     }
372   }
373   delete [] emcExist ; 
374   delete [] tpcExist ; 
375 }
376
377 //____________________________________________________________________________
378 void  AliPHOSTrackSegmentMakerv2::Exec(Option_t *option)
379 {
380   // Steering method to perform track segment construction for events
381   // in the range from fFirstEvent to fLastEvent.
382   // This range is optionally set by SetEventRange().
383   // if fLastEvent=-1 (by default), then process events until the end.
384   
385   if(strstr(option,"tim"))
386     gBenchmark->Start("PHOSTSMaker");
387  
388   if(strstr(option,"print")) {
389     Print() ; 
390     return ; 
391   }
392   
393   AliPHOSGetter * gime = AliPHOSGetter::Instance() ;  
394  
395   const AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
396
397   if (fLastEvent == -1) 
398     fLastEvent = gime->MaxEvent() - 1 ;
399   else 
400     fLastEvent = TMath::Min(fFirstEvent,gime->MaxEvent());
401   Int_t nEvents   = fLastEvent - fFirstEvent + 1;
402
403   Int_t ievent ; 
404   for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
405     gime->Event(ievent,"DR") ;
406    //Make some initializations 
407     fNTrackSegments = 0 ;
408     fEmcFirst = 0 ;    
409     fEmcLast  = 0 ;   
410     
411     gime->TrackSegments()->Clear();
412
413     //    if(!ReadRecPoints(ievent))   continue; //reads RecPoints for event ievent
414     
415     for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ) {
416       FillOneModule() ; 
417       MakeLinks() ;
418       MakePairs() ;
419     }
420
421     WriteTrackSegments() ;
422
423     if(strstr(option,"deb"))
424       PrintTrackSegments(option);
425     
426     //increment the total number of track segments per run 
427     fTrackSegmentsInRun += gime->TrackSegments()->GetEntriesFast() ; 
428   }
429   
430   if(strstr(option,"tim")){
431     gBenchmark->Stop("PHOSTSMaker");
432     Info("Exec", "took %f seconds for making TS %f seconds per event", 
433           gBenchmark->GetCpuTime("PHOSTSMaker"), 
434           gBenchmark->GetCpuTime("PHOSTSMaker")/nEvents) ;
435    }
436   if(fWrite) //do not unload in "on flight" mode
437     Unload();
438 }
439 //____________________________________________________________________________
440 void AliPHOSTrackSegmentMakerv2::Unload() 
441 {
442   // Unloads the task from the folder
443   AliPHOSGetter * gime = AliPHOSGetter::Instance() ;  
444   gime->PhosLoader()->UnloadRecPoints() ;
445   gime->PhosLoader()->UnloadTracks() ;
446 }
447
448 //____________________________________________________________________________
449 void AliPHOSTrackSegmentMakerv2::Print(const Option_t *)const
450 {
451   //  Print TrackSegmentMaker parameters
452
453   TString message("") ;
454   if( strcmp(GetName(), "") != 0 ) {
455     message = "\n======== AliPHOSTrackSegmentMakerv2 ========\n" ; 
456     message += "Making Track segments\n" ;
457     message += "with parameters:\n" ; 
458     message += "     Maximal EMC - TPC distance (cm) %f\n" ;
459     message += "============================================\n" ;
460     Info("Print", message.Data(),fRtpc) ;
461   }
462   else
463     Info("Print", "AliPHOSTrackSegmentMakerv2 not initialized ") ;
464 }
465
466 //____________________________________________________________________________
467 void AliPHOSTrackSegmentMakerv2::WriteTrackSegments()
468 {
469   // Writes found TrackSegments to TreeR. Creates branches 
470   // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
471   // In the former branch found TrackSegments are stored, while 
472   // in the latter all parameters, with which TS were made. 
473   // ROOT does not allow overwriting existing branches, therefore
474   // first we check, if branches with the same title already exist.
475   // If yes - exits without writing.
476
477   AliPHOSGetter *gime = AliPHOSGetter::Instance() ; 
478
479   TClonesArray * trackSegments = gime->TrackSegments() ; 
480   trackSegments->Expand(trackSegments->GetEntriesFast()) ;
481
482   if(fWrite){ //We write TreeT
483     TTree * treeT = gime->TreeT();
484     
485     //First TS
486     Int_t bufferSize = 32000 ; 
487     TBranch * tsBranch = treeT->Branch("PHOSTS",&trackSegments,bufferSize);
488     tsBranch->Fill() ;  
489     
490     gime->WriteTracks("OVERWRITE");
491     gime->WriteTrackSegmentMaker("OVERWRITE");
492   }
493 }
494
495
496 //____________________________________________________________________________
497 void AliPHOSTrackSegmentMakerv2::PrintTrackSegments(Option_t * option)
498 {
499   // option deb - prints # of found TrackSegments
500   // option deb all - prints as well indexed of found RecParticles assigned to the TS
501
502   TClonesArray * trackSegments = AliPHOSGetter::Instance()->TrackSegments() ; 
503
504   Info("PrintTrackSegments", "Results from TrackSegmentMaker:") ; 
505   printf("nevent: %d\n", gAlice->GetEvNumber()) ; 
506   printf("        Found %d TrackSegments\n", trackSegments->GetEntriesFast() ); 
507   
508   if(strstr(option,"all")) {  // printing found TS
509     printf("TrackSegment #  EMC RP#  CPV RP#\n") ; 
510     Int_t index;
511     for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
512       AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ; 
513       printf("   %d           %d        %d \n", ts->GetIndexInList(), ts->GetEmcIndex(), ts->GetTrackIndex() ) ; 
514     }   
515   }
516 }