]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSTrackSegmentMakerv1.cxx
cosmetics
[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 = fRcpv ; 
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       AliESDtrack *track;
160       for (Int_t iTrack=0; iTrack<nTracks; iTrack++) {
161         track = fESD->GetTrack(iTrack);
162         if (track->IsPHOS()) 
163           continue ; 
164         track->GetOuterXYZ(xyz);     // track coord on the cylinder of PHOS radius
165         if ((TMath::Abs(xyz[0])+TMath::Abs(xyz[1])+TMath::Abs(xyz[2]))<=0)
166           continue;
167         track->GetOuterPxPyPz(pxyz); // track momentum ibid.
168         vecDist = PropagateToPlane(xyz,pxyz,"CPV",cpvClu->GetPHOSMod());
169         //      Info("GetDistanceInPHOSPlane","Track %d propagation to CPV = (%f,%f,%f)",
170         //     iTrack,vecDist.X(),vecDist.Y(),vecDist.Z());
171         vecDist -= vecCpv;
172         distance2Track = TMath::Sqrt(vecDist.X()*vecDist.X() + vecDist.Z()*vecDist.Z());
173         // Find the closest track to the EMC recpoint
174         if (distance2Track < minDistance) {
175           minDistance = distance2Track;
176           iClosestTrack = iTrack;
177         }
178       }
179
180       if (iClosestTrack != -1) {
181         track = fESD->GetTrack(iClosestTrack);
182         track->GetOuterPxPyPz(pxyz); // track momentum ibid.
183         TVector3 vecCpvGlobal; // Global position of the CPV recpoint
184         AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
185         const AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
186         geom->GetGlobal((AliRecPoint*)cpvClu,vecCpvGlobal);
187         for (Int_t ixyz=0; ixyz<3; ixyz++)
188           xyz[ixyz] = vecCpvGlobal[ixyz];
189         vecDist = PropagateToPlane(xyz,pxyz,"EMC",cpvClu->GetPHOSMod());
190 //      Info("GetDistanceInPHOSPlane","Track %d propagation to EMC = (%f,%f,%f)",
191 //           iClosestTrack,vecDist.X(),vecDist.Y(),vecDist.Z());
192         vecDist -= vecEmc;
193         distance2Track = TMath::Sqrt(vecDist.X()*vecDist.X() + vecDist.Z()*vecDist.Z());
194       }
195 //     } else {
196 //       // If no ESD exists, than simply find EMC-CPV distance
197 //       distance = (vecCpv - vecEmc).Mag() ;
198     
199       if(distance2Track < fRcpv + 2*delta )
200         trackindex = iClosestTrack ; 
201       //      toofar = kFALSE ;
202     }
203     //     Info("GetDistanceInPHOSPlane","cpv-emc distance is %f cm",
204     //   distance);
205   }
206   
207   return distance2Cpv ;
208 }
209
210 //____________________________________________________________________________
211 TVector3  AliPHOSTrackSegmentMakerv1::PropagateToPlane(Double_t *x, Double_t *p,
212                                                        char *det, Int_t moduleNumber) const
213 {
214   // Propagate a straight-line track from the origin point x
215   // along the direction p to the CPV or EMC module moduleNumber
216   // Returns a local position of such a propagation
217
218   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
219   const AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
220   TVector3 moduleCenter = geom->GetModuleCenter(det,moduleNumber);
221   TVector3 vertex(x);
222   TVector3 direction(p);
223
224 //   Info("PropagateToCPV","Center of the %s module %d is (%f,%f,%f)",
225 //        det,moduleNumber,moduleCenter[0],moduleCenter[1],moduleCenter[2]);
226
227   Double_t time = (moduleCenter.Mag2() - vertex.Dot(moduleCenter)) /
228     (direction.Dot(moduleCenter));
229   TVector3 globalIntersection = vertex + direction*time;
230   return geom->Global2Local(globalIntersection,moduleNumber);
231 }
232
233 //____________________________________________________________________________
234 void  AliPHOSTrackSegmentMakerv1::Init()
235 {
236   // Make all memory allocations that are not possible in default constructor
237   
238   AliPHOSGetter* gime = AliPHOSGetter::Instance(GetTitle(), fEventFolderName.Data());
239   
240   fLinkUpArray  = new TClonesArray("AliPHOSLink", 1000); 
241   if ( !gime->TrackSegmentMaker() ) {
242     gime->PostTrackSegmentMaker(this);
243   }
244 }
245
246 //____________________________________________________________________________
247 void  AliPHOSTrackSegmentMakerv1::InitParameters()
248 {
249   //Initializes parameters
250   fRcpv      = 10. ;   
251   fEmcFirst  = 0 ;    
252   fEmcLast   = 0 ;   
253   fCpvFirst  = 0 ;   
254   fCpvLast   = 0 ;   
255   fLinkUpArray = 0 ;
256   fTrackSegmentsInRun       = 0 ; 
257   SetEventRange(0,-1) ;
258 }
259
260
261 //____________________________________________________________________________
262 void  AliPHOSTrackSegmentMakerv1::MakeLinks()const
263
264   // Finds distances (links) between all EMC and CPV clusters, 
265   // which are not further apart from each other than fRcpv 
266   // and sort them in accordance with this distance
267   
268   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
269   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
270   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
271
272   fLinkUpArray->Clear() ;    
273
274   AliPHOSCpvRecPoint * cpv ;
275   AliPHOSEmcRecPoint * emcclu ;
276
277   Int_t iLinkUp  = 0 ;
278   
279   Int_t iEmcRP;
280   for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
281     emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
282
283     //Bool_t toofar ;        
284     Int_t iCpv = 0 ;    
285     for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) { 
286       
287       cpv = dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(iCpv)) ;
288       Int_t track = -1 ; 
289       Float_t r = GetDistanceInPHOSPlane(emcclu, cpv, track) ;     
290       //      if(toofar)
291       //        continue ;       
292       if(r < fRcpv) { 
293         new ((*fLinkUpArray)[iLinkUp++])  AliPHOSLink(r, iEmcRP, iCpv, track) ;
294       }      
295     }
296   } 
297   
298   fLinkUpArray->Sort() ;  //first links with smallest distances
299 }
300
301 //____________________________________________________________________________
302 void  AliPHOSTrackSegmentMakerv1::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   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
313   TClonesArray * trackSegments = gime->TrackSegments();
314     
315   //Make arrays to mark clusters already chosen
316   Int_t * emcExist = 0;
317   if(fEmcLast > fEmcFirst)
318     emcExist = new Int_t[fEmcLast-fEmcFirst] ;
319   
320   Int_t index;
321   for(index = 0; index <fEmcLast-fEmcFirst; index ++)
322     emcExist[index] = 1 ;
323   
324   Bool_t * cpvExist = 0;
325   if(fCpvLast > fCpvFirst)
326     cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
327   for(index = 0; index <fCpvLast-fCpvFirst; index ++)
328     cpvExist[index] = kTRUE ;
329   
330   
331   // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance 
332   TIter nextUp(fLinkUpArray) ;
333   
334   AliPHOSLink * linkUp ;
335   
336   AliPHOSCpvRecPoint * nullpointer = 0 ;
337   
338   while ( (linkUp =  static_cast<AliPHOSLink *>(nextUp()) ) ){  
339
340     if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){
341
342       if(cpvExist[linkUp->GetCpv()-fCpvFirst]){ //CPV still exist
343          new ((* trackSegments)[fNTrackSegments]) 
344            AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) , 
345                                dynamic_cast<AliPHOSCpvRecPoint *>(cpvRecPoints->At(linkUp->GetCpv())) , 
346                                linkUp->GetTrack()) ;
347          
348        (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
349        fNTrackSegments++ ;
350        emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc  that Cpv was found 
351        //mark CPV recpoint as already used 
352        cpvExist[linkUp->GetCpv()-fCpvFirst] = kFALSE ;
353       } //if CpvUp still exist
354     } 
355   }        
356
357   //look through emc recPoints left without CPV
358   if(emcExist){ //if there is emc rec point
359     Int_t iEmcRP ;
360     for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst  ; iEmcRP++ ){
361       if(emcExist[iEmcRP] > 0 ){
362        new ((*trackSegments)[fNTrackSegments])  
363          AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)), 
364                            nullpointer) ;
365        (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
366        fNTrackSegments++;    
367       } 
368     }
369   }
370   delete [] emcExist ; 
371   delete [] cpvExist ; 
372 }
373
374 //____________________________________________________________________________
375 void  AliPHOSTrackSegmentMakerv1::Exec(Option_t *option)
376 {
377   // Steering method to perform track segment construction for events
378   // in the range from fFirstEvent to fLastEvent.
379   // This range is optionally set by SetEventRange().
380   // if fLastEvent=-1 (by default), then process events until the end.
381   
382   if(strstr(option,"tim"))
383     gBenchmark->Start("PHOSTSMaker");
384  
385   if(strstr(option,"print")) {
386     Print() ; 
387     return ; 
388   }
389   
390   AliPHOSGetter * gime = AliPHOSGetter::Instance(GetTitle()) ;  
391  
392   const AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
393
394   if (fLastEvent == -1) 
395     fLastEvent = gime->MaxEvent() - 1 ;
396   else 
397     fLastEvent = TMath::Min(fFirstEvent,gime->MaxEvent());
398   Int_t nEvents   = fLastEvent - fFirstEvent + 1;
399
400   Int_t ievent ; 
401   for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
402     gime->Event(ievent,"R") ;
403    //Make some initializations 
404     fNTrackSegments = 0 ;
405     fEmcFirst = 0 ;    
406     fEmcLast  = 0 ;   
407     fCpvFirst = 0 ;   
408     fCpvLast  = 0 ;   
409     
410     gime->TrackSegments()->Clear();
411
412     //    if(!ReadRecPoints(ievent))   continue; //reads RecPoints for event ievent
413     
414     for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ) {
415       FillOneModule() ; 
416       MakeLinks() ;
417       MakePairs() ;
418     }
419
420     WriteTrackSegments() ;
421
422     if(strstr(option,"deb"))
423       PrintTrackSegments(option);
424     
425     //increment the total number of track segments per run 
426     fTrackSegmentsInRun += gime->TrackSegments()->GetEntriesFast() ; 
427   }
428   
429   if(strstr(option,"tim")){
430     gBenchmark->Stop("PHOSTSMaker");
431     Info("Exec", "took %f seconds for making TS %f seconds per event", 
432           gBenchmark->GetCpuTime("PHOSTSMaker"), 
433           gBenchmark->GetCpuTime("PHOSTSMaker")/nEvents) ;
434    }
435   Unload();
436 }
437
438 //____________________________________________________________________________
439 void AliPHOSTrackSegmentMakerv1::Unload() 
440 {
441   // Unloads the task from the folder
442   AliPHOSGetter * gime = AliPHOSGetter::Instance() ;  
443   gime->PhosLoader()->UnloadRecPoints() ;
444   gime->PhosLoader()->UnloadTracks() ;
445 }
446
447 //____________________________________________________________________________
448 void AliPHOSTrackSegmentMakerv1::Print()const
449 {
450   //  Print TrackSegmentMaker parameters
451
452   TString message("") ;
453   if( strcmp(GetName(), "") != 0 ) {
454     message = "\n======== AliPHOSTrackSegmentMakerv1 ========\n" ; 
455     message += "Making Track segments\n" ;
456     message += "with parameters:\n" ; 
457     message += "     Maximal EMC - CPV distance (cm) %f\n" ;
458     message += "============================================\n" ;
459     Info("Print", message.Data(),fRcpv) ;
460   }
461   else
462     Info("Print", "AliPHOSTrackSegmentMakerv1 not initialized ") ;
463 }
464
465 //____________________________________________________________________________
466 void AliPHOSTrackSegmentMakerv1::WriteTrackSegments()
467 {
468   // Writes found TrackSegments to TreeR. Creates branches 
469   // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
470   // In the former branch found TrackSegments are stored, while 
471   // in the latter all parameters, with which TS were made. 
472   // ROOT does not allow overwriting existing branches, therefore
473   // first we check, if branches with the same title already exist.
474   // If yes - exits without writing.
475
476   AliPHOSGetter *gime = AliPHOSGetter::Instance() ; 
477
478   TClonesArray * trackSegments = gime->TrackSegments() ; 
479   trackSegments->Expand(trackSegments->GetEntriesFast()) ;
480
481   TTree * treeT = gime->TreeT();
482  
483   //First TS
484   Int_t bufferSize = 32000 ; 
485   TBranch * tsBranch = treeT->Branch("PHOSTS",&trackSegments,bufferSize);
486   tsBranch->Fill() ;  
487
488   gime->WriteTracks("OVERWRITE");
489   gime->WriteTrackSegmentMaker("OVERWRITE");
490 }
491
492
493 //____________________________________________________________________________
494 void AliPHOSTrackSegmentMakerv1::PrintTrackSegments(Option_t * option)
495 {
496   // option deb - prints # of found TrackSegments
497   // option deb all - prints as well indexed of found RecParticles assigned to the TS
498
499   TClonesArray * trackSegments = AliPHOSGetter::Instance()->TrackSegments() ; 
500
501   Info("PrintTrackSegments", "Results from TrackSegmentMaker:") ; 
502   printf("nevent: %d\n", gAlice->GetEvNumber()) ; 
503   printf("        Found %d TrackSegments\n", trackSegments->GetEntriesFast() ); 
504   
505   if(strstr(option,"all")) {  // printing found TS
506     printf("TrackSegment #  EMC RP#  CPV RP#\n") ; 
507     Int_t index;
508     for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
509       AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ; 
510       printf("   %d           %d        %d \n", ts->GetIndexInList(), ts->GetEmcIndex(), ts->GetCpvIndex() ) ; 
511     }   
512   }
513 }