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