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