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