]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSTrackSegmentMakerv1.cxx
907881e130a4bb4695987be03f6f479529d6f6aa
[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 + (possibly) PPSD RecPoint
20 // To find TrackSegments we do the following: 
21 //  for each EMC RecPoints we look at
22 //   CPV/PPSD RecPoints in the radious fR0. 
23 //  If there is such a CPV RecPoint, 
24 //   we make "Link" it is just indexes of EMC and CPV/PPSD 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/PPSD RecPoint we make TrackSegment 
30 // consisting from EMC alone. There is no TrackSegments without EMC RecPoint.
31 //
32 // In principle this class should be called from AliPHOSReconstructioner, but 
33 // one can use it as well in standalone mode.
34 // Use  case:
35 //  root [0] AliPHOSTrackSegmentMakerv1 * t = new AliPHOSTrackSegmentMaker("galice.root", "tracksegmentsname", "recpointsname")
36 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
37 //               // reads gAlice from header file "galice.root", uses recpoints stored in the branch names "recpointsname" (default = "Default")
38 //               // and saves recpoints in branch named "tracksegmentsname" (default = "recpointsname")                       
39 //  root [1] t->ExecuteTask()
40 //  root [2] t->SetMaxEmcPpsdDistance(5)
41 //  root [3] t->SetTrackSegmentsBranch("max distance 5 cm")
42 //  root [4] t->ExecuteTask("deb all time") 
43 //                 
44 //*-- Author: Dmitri Peressounko (RRC Ki & SUBATECH) & Yves Schutz (SUBATECH) 
45 //
46
47 // --- ROOT system ---
48 #include "TROOT.h"
49 #include "TFile.h"
50 #include "TFolder.h"
51 #include "TTree.h"
52 #include "TSystem.h"
53 #include "TBenchmark.h"
54 // --- Standard library ---
55
56 #include <iostream.h>
57 #include <iomanip.h>
58
59 // --- AliRoot header files ---
60
61 #include "AliPHOSTrackSegmentMakerv1.h"
62 #include "AliPHOSClusterizerv1.h"
63 #include "AliPHOSTrackSegment.h"
64 #include "AliPHOSCpvRecPoint.h"
65 #include "AliPHOSLink.h"
66 #include "AliPHOSGetter.h"
67 #include "AliPHOS.h"
68 #include "AliRun.h"
69
70 ClassImp( AliPHOSTrackSegmentMakerv1) 
71
72
73 //____________________________________________________________________________
74   AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1() : AliPHOSTrackSegmentMaker()
75 {
76   // default ctor (to be used mainly by Streamer)
77
78   InitParameters() ; 
79   fHeaderFileName           = "" ;
80   fRecPointsBranchTitle     = "" ;
81   fTrackSegmentsBranchTitle = "" ; 
82   fFrom                     = "" ; 
83
84   fTrackSegmentsInRun       = 0 ; 
85
86   fDefaultInit = kTRUE ; 
87 }
88
89 //____________________________________________________________________________
90  AliPHOSTrackSegmentMakerv1::AliPHOSTrackSegmentMakerv1(const char * headerFile, const char * name, const char * from) : AliPHOSTrackSegmentMaker(headerFile, name)
91 {
92   // ctor
93
94   InitParameters() ; 
95   fHeaderFileName           = GetTitle() ;
96   fRecPointsBranchTitle     = GetName() ;
97   fTrackSegmentsBranchTitle = GetName() ; 
98   fTrackSegmentsInRun       = 0 ; 
99
100   if ( from == 0 ) 
101     fFrom = name ; 
102   else
103     fFrom = from ; 
104   Init() ;
105
106   fDefaultInit = kFALSE ; 
107
108 }
109
110 //____________________________________________________________________________
111  AliPHOSTrackSegmentMakerv1::~AliPHOSTrackSegmentMakerv1()
112
113   // dtor
114   // fDefaultInit = kTRUE if TrackSegmentMaker created by default ctor (to get just the parameters)
115   
116   if (!fDefaultInit) {
117     delete fLinkUpArray  ;
118     
119     AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
120     
121     // remove the task from the folder list
122     gime->RemoveTask("T",GetName()) ;
123     TString name(GetName()) ; 
124     name.ReplaceAll("tsm", "clu") ; 
125     gime->RemoveTask("C",name) ;
126     
127     // remove the data from the folder list
128     name = GetName() ; 
129     name.Remove(name.Index(":")) ; 
130     gime->RemoveObjects("RE", name) ; // EMCARecPoints
131     gime->RemoveObjects("RC", name) ; // CPVRecPoints
132     gime->RemoveObjects("T", name) ;  // TrackSegments
133     
134     // Delete gAlice
135     gime->CloseFile() ; 
136     
137     fSplitFile = 0 ; 
138   }
139 }
140
141
142 //____________________________________________________________________________
143 const TString AliPHOSTrackSegmentMakerv1::BranchName() const 
144 {  
145   TString branchName(GetName() ) ;
146   branchName.Remove(branchName.Index(Version())-1) ;
147   return branchName ;
148 }
149
150 //____________________________________________________________________________
151 void  AliPHOSTrackSegmentMakerv1::FillOneModule()
152 {
153   // Finds first and last indexes between which 
154   // clusters from one PHOS module are
155   
156   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
157   TObjArray * emcRecPoints = gime->EmcRecPoints(fFrom) ; 
158   TObjArray * cpvRecPoints = gime->CpvRecPoints(fFrom) ; 
159  
160   //First EMC clusters
161   Int_t totalEmc = emcRecPoints->GetEntriesFast() ;
162   for(fEmcFirst = fEmcLast; (fEmcLast < totalEmc) &&  
163         ((dynamic_cast<AliPHOSRecPoint *>(emcRecPoints->At(fEmcLast)))->GetPHOSMod() == fModule ); 
164       fEmcLast ++)  ;
165   
166   //Now CPV clusters
167   Int_t totalCpv = cpvRecPoints->GetEntriesFast() ;
168
169     for(fCpvFirst = fCpvLast; (fCpvLast < totalCpv) && 
170           ((dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(fCpvLast)))->GetPHOSMod() == fModule ); 
171         fCpvLast ++) ;
172       
173 }
174
175 //____________________________________________________________________________
176 Float_t  AliPHOSTrackSegmentMakerv1::GetDistanceInPHOSPlane(AliPHOSEmcRecPoint * emcClu,AliPHOSRecPoint * cpvClu, Bool_t &toofar)const
177 {
178   // Calculates the distance between the EMC RecPoint and the PPSD RecPoint
179   // Clusters are sorted in "rows" and "columns" of width 1 cm
180
181   Float_t delta = 1 ;  // Width of the rows in sorting of RecPoints (in cm)
182                        // if you change this value, change it as well in xxxRecPoint::Compare()
183   Float_t r = fR0 ;
184  
185   TVector3 vecEmc ;
186   TVector3 vecCpv ;
187   
188   emcClu->GetLocalPosition(vecEmc) ;
189   cpvClu->GetLocalPosition(vecCpv)  ; 
190
191   if(emcClu->GetPHOSMod() == cpvClu->GetPHOSMod()){ 
192     if(vecCpv.X() <= vecEmc.X() + fR0 + 2*delta ){ 
193
194       vecCpv = vecCpv  - vecEmc ; 
195       r = vecCpv.Mag() ;
196       toofar = kFALSE ;
197
198     } // if  xPpsd >= xEmc + ...
199     else 
200       toofar = kTRUE ;
201   } 
202   else 
203     toofar = kTRUE ;
204
205   //toofar = kFALSE ;
206  
207   
208   return r ;
209 }
210
211 //____________________________________________________________________________
212 void  AliPHOSTrackSegmentMakerv1::Init()
213 {
214   // Make all memory allocations that are not possible in default constructor
215   
216   if ( strcmp(GetTitle(), "") == 0 )
217     SetTitle("galice.root") ;
218     
219   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), fFrom.Data()) ; 
220   if ( gime == 0 ) {
221     cerr << "ERROR: AliPHOSTrackSegmentMakerv1::Init -> Could not obtain the Getter object !" << endl ; 
222     return ;
223   } 
224   
225   fLinkUpArray  = new TClonesArray("AliPHOSLink", 1000); 
226   
227   //add Task to //YSAlice/tasks/Reconstructioner/PHOS
228   gime->PostTrackSegmentMaker(this) ;
229
230   // create a folder on the white board //YSAlice/WhiteBoard/RecPoints/PHOS/trackSegmentsName
231   gime->PostTrackSegments(BranchName()) ; 
232
233 }
234
235 //____________________________________________________________________________
236 void  AliPHOSTrackSegmentMakerv1::InitParameters()
237 {
238   fR0        = 10. ;   
239   fEmcFirst  = 0 ;    
240   fEmcLast   = 0 ;   
241   fCpvFirst  = 0 ;   
242   fCpvLast   = 0 ;   
243   fLinkUpArray = 0 ;
244   TString tsmName( GetName()) ; 
245   if (tsmName.IsNull() ) 
246     tsmName = "Default" ; 
247   tsmName.Append(":") ; 
248   tsmName.Append(Version()) ; 
249   SetName(tsmName) ;
250 }
251
252
253 //____________________________________________________________________________
254 void  AliPHOSTrackSegmentMakerv1::MakeLinks()const
255
256   // Finds distances (links) between all EMC and PPSD clusters, 
257   // which are not further apart from each other than fR0 
258   // and sort them in accordance with this distance
259   
260   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
261   TObjArray * emcRecPoints = gime->EmcRecPoints(fFrom) ; 
262   TObjArray * cpvRecPoints = gime->CpvRecPoints(fFrom) ; 
263
264   fLinkUpArray->Clear() ;    
265
266   AliPHOSRecPoint * cpv ;
267   AliPHOSEmcRecPoint * emcclu ;
268
269   Int_t iLinkUp  = 0 ;
270   
271   Int_t iEmcRP;
272   for(iEmcRP = fEmcFirst; iEmcRP < fEmcLast; iEmcRP++ ) {
273     emcclu = dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP)) ;
274
275     Bool_t toofar ;        
276     Int_t iCpv = 0 ;    
277     for(iCpv = fCpvFirst; iCpv < fCpvLast;iCpv++ ) { 
278       
279       cpv = dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(iCpv)) ;
280       Float_t r = GetDistanceInPHOSPlane(emcclu, cpv, toofar) ;
281       
282       if(toofar)
283         break ;  
284       if(r < fR0) { 
285         new ((*fLinkUpArray)[iLinkUp++])  AliPHOSLink(r, iEmcRP, iCpv) ;
286       }      
287     }
288   } 
289   
290   fLinkUpArray->Sort() ;  //first links with smallest distances
291 }
292
293 //____________________________________________________________________________
294 void  AliPHOSTrackSegmentMakerv1::MakePairs()
295
296   // Using the previously made list of "links", we found the smallest link - i.e. 
297   // link with the least distance between EMC and CPV and pointing to still 
298   // unassigned RecParticles. We assign these RecPoints to TrackSegment and 
299   // remove them from the list of "unassigned". 
300   
301   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
302   TObjArray * emcRecPoints = gime->EmcRecPoints(fFrom) ; 
303   TObjArray * cpvRecPoints = gime->CpvRecPoints(fFrom) ; 
304   TClonesArray * trackSegments = gime->TrackSegments(BranchName()) ;   
305     
306   //Make arrays to mark clusters already chosen
307   Int_t * emcExist = 0;
308   if(fEmcLast > fEmcFirst)
309     emcExist = new Int_t[fEmcLast-fEmcFirst] ;
310   
311   Int_t index;
312   for(index = 0; index <fEmcLast-fEmcFirst; index ++)
313     emcExist[index] = 1 ;
314   
315   Bool_t * cpvExist = 0;
316   if(fCpvLast > fCpvFirst)
317     cpvExist = new Bool_t[fCpvLast-fCpvFirst] ;
318   for(index = 0; index <fCpvLast-fCpvFirst; index ++)
319     cpvExist[index] = kTRUE ;
320   
321   
322   // Finds the smallest links and makes pairs of CPV and EMC clusters with smallest distance 
323   TIter nextUp(fLinkUpArray) ;
324   
325   AliPHOSLink * linkUp ;
326   
327   AliPHOSRecPoint * nullpointer = 0 ;
328   
329   while ( (linkUp =  static_cast<AliPHOSLink *>(nextUp()) ) ){  
330
331     if(emcExist[linkUp->GetEmc()-fEmcFirst] != -1){ //without ppsd Up yet 
332
333       if(cpvExist[linkUp->GetPpsd()-fCpvFirst]){ //CPV still exist
334         
335         new ((* trackSegments)[fNTrackSegments]) 
336           AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(linkUp->GetEmc())) , 
337                               dynamic_cast<AliPHOSRecPoint *>(cpvRecPoints->At(linkUp->GetPpsd()))) ;
338         (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
339         fNTrackSegments++ ;
340         
341         emcExist[linkUp->GetEmc()-fEmcFirst] = -1 ; //Mark emc  that Cpv was found 
342         //mark CPV recpoint as already used 
343         cpvExist[linkUp->GetPpsd()-fCpvFirst] = kFALSE ;
344       } //if ppsdUp still exist
345     } 
346   }      
347
348   //look through emc recPoints left without CPV/PPSD
349   if(emcExist){ //if there is emc rec point
350     Int_t iEmcRP ;
351     for(iEmcRP = 0; iEmcRP < fEmcLast-fEmcFirst  ; iEmcRP++ ){
352       if(emcExist[iEmcRP] > 0 ){
353         new ((*trackSegments)[fNTrackSegments])  
354           AliPHOSTrackSegment(dynamic_cast<AliPHOSEmcRecPoint *>(emcRecPoints->At(iEmcRP+fEmcFirst)), 
355                               nullpointer) ;
356         (dynamic_cast<AliPHOSTrackSegment *>(trackSegments->At(fNTrackSegments)))->SetIndexInList(fNTrackSegments);
357         fNTrackSegments++;    
358       } 
359     }
360   }
361   delete [] emcExist ; 
362   delete [] cpvExist ; 
363 }
364
365 //____________________________________________________________________________
366 void  AliPHOSTrackSegmentMakerv1::Exec(Option_t * option)
367 {
368   // STEERing method
369
370   if( strcmp(GetName(), "")== 0 ) 
371     Init() ;
372
373   if(strstr(option,"tim"))
374     gBenchmark->Start("PHOSTSMaker");
375  
376   if(strstr(option,"print")) {
377     Print("") ; 
378     return ; 
379   }
380
381   gAlice->GetEvent(0) ;
382   //check, if the branch with name of this" already exits?
383   if (gAlice->TreeR()) { 
384     TObjArray * lob = static_cast<TObjArray*>(gAlice->TreeR()->GetListOfBranches()) ;
385     TIter next(lob) ; 
386     TBranch * branch = 0 ;  
387     Bool_t phostsfound = kFALSE, tracksegmentmakerfound = kFALSE ; 
388     
389     TString branchname = GetName() ;
390     branchname.Remove(branchname.Index(Version())-1) ;
391     
392     while ( (branch = static_cast<TBranch*>(next())) && (!phostsfound || !tracksegmentmakerfound) ) {
393       if ( (strcmp(branch->GetName(), "PHOSTS")==0) && (strcmp(branch->GetTitle(), branchname.Data())==0) ) 
394         phostsfound = kTRUE ;
395       
396       else if ( (strcmp(branch->GetName(), "AliPHOSTrackSegmentMaker")==0) && (strcmp(branch->GetTitle(), GetName())==0) ) 
397         tracksegmentmakerfound = kTRUE ; 
398     }
399     
400     if ( phostsfound || tracksegmentmakerfound ) {
401       cerr << "WARNING: AliPHOSTrackSegmentMakerv1::Exec -> TrackSegments and/or TrackSegmentMaker branch with name " 
402            << branchname.Data() << " already exits" << endl ;
403       return ; 
404     }       
405   }
406
407   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
408   const AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
409   Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
410   Int_t ievent ;
411
412   for(ievent = 0; ievent < nevents; ievent++){
413
414     gime->Event(ievent,"R") ;
415     //Make some initializations 
416     fNTrackSegments = 0 ;
417     fEmcFirst = 0 ;    
418     fEmcLast  = 0 ;   
419     fCpvFirst = 0 ;   
420     fCpvLast  = 0 ;   
421     gime->TrackSegments(BranchName())->Clear() ; 
422
423     //    if(!ReadRecPoints(ievent))   continue; //reads RecPoints for event ievent
424     
425     for(fModule = 1; fModule <= geom->GetNModules() ; fModule++ ){
426       
427       FillOneModule() ; 
428       
429       MakeLinks() ;
430       
431       MakePairs() ;
432       
433     }
434
435     WriteTrackSegments(ievent) ;
436
437     if(strstr(option,"deb"))
438       PrintTrackSegments(option) ;
439     
440     //increment the total number of track segments per run 
441     fTrackSegmentsInRun += gime->TrackSegments(BranchName())->GetEntriesFast() ; 
442
443   }
444
445   if(strstr(option,"tim")){
446     gBenchmark->Stop("PHOSTSMaker");
447     cout << "AliPHOSTSMaker:" << endl ;
448     cout << "  took " << gBenchmark->GetCpuTime("PHOSTSMaker") << " seconds for making TS " 
449          <<  gBenchmark->GetCpuTime("PHOSTSMaker")/nevents << " seconds per event " << endl ;
450     cout << endl ;
451   }
452     
453 }
454
455 //____________________________________________________________________________
456 void AliPHOSTrackSegmentMakerv1::Print(Option_t * option)const
457 {
458   //  Print TrackSegmentMaker parameters
459
460   if( strcmp(GetName(), "") != 0 ) {
461     cout <<  "======== AliPHOSTrackSegmentMakerv1 ========" << endl ;
462     cout <<  "Making Track segments "<< endl ;
463     cout <<  "    Headers file:                   " << fHeaderFileName.Data() << endl ;
464     cout <<  "    RecPoints branch file name:     " << fRecPointsBranchTitle.Data() << endl ;
465     cout <<  "    TrackSegments Branch file name: " << fTrackSegmentsBranchTitle.Data() << endl ;
466     cout <<  "with parameters: " << endl ;
467     cout <<  "    Maximal EMC - CPV (PPSD) distance (cm)" << fR0 << endl ;
468     cout <<  "============================================" << endl ;
469   }
470   else
471     cout << "AliPHOSTrackSegmentMakerv1 not initialized " << endl ;
472 }
473
474 //____________________________________________________________________________
475 void AliPHOSTrackSegmentMakerv1::WriteTrackSegments(Int_t event)
476 {
477   // Writes found TrackSegments to TreeR. Creates branches 
478   // "PHOSTS" and "AliPHOSTrackSegmentMaker" with the same title.
479   // In the former branch found TrackSegments are stored, while 
480   // in the latter all parameters, with which TS were made. 
481   // ROOT does not allow overwriting existing branches, therefore
482   // first we check, if branches with the same title already exist.
483   // If yes - exits without writing.
484   
485   AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ; 
486
487   TClonesArray * trackSegments = gime->TrackSegments(BranchName()) ; 
488   trackSegments->Expand(trackSegments->GetEntriesFast()) ;
489   TTree * treeR = gAlice->TreeR();
490
491   if (!treeR) 
492     gAlice->MakeTree("R", fSplitFile);
493   treeR = gAlice->TreeR(); 
494
495   //First TS
496   Int_t bufferSize = 32000 ;    
497   TBranch * tsBranch = treeR->Branch("PHOSTS",&trackSegments,bufferSize);
498   tsBranch->SetTitle(BranchName());
499
500   //Second -TSMaker
501   Int_t splitlevel = 0 ;
502   AliPHOSTrackSegmentMakerv1 * ts = this ;
503   TBranch * tsMakerBranch = treeR->Branch("AliPHOSTrackSegmentMaker","AliPHOSTrackSegmentMakerv1",
504                                           &ts,bufferSize,splitlevel);
505   tsMakerBranch->SetTitle(BranchName());
506
507   tsBranch->Fill() ;  
508   tsMakerBranch->Fill() ;
509
510   treeR->AutoSave() ; //Write(0,kOverwrite) ;  
511   
512 }
513
514
515 //____________________________________________________________________________
516 void AliPHOSTrackSegmentMakerv1::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   TString taskName(GetName()) ; 
521   taskName.Remove(taskName.Index(Version())-1) ;
522   
523   TClonesArray * trackSegments = AliPHOSGetter::GetInstance()->TrackSegments(taskName) ; 
524
525   
526   cout << "AliPHOSTrackSegmentMakerv1: event "<<gAlice->GetEvNumber()  << endl ;
527   cout << "       Found " << trackSegments->GetEntriesFast() << "  trackSegments " << endl ;
528   
529   if(strstr(option,"all")) {  // printing found TS
530     cout << "TrackSegment # " << "    EMC RP#    " << "    CPV RP#    " << endl ; 
531     
532     Int_t index;
533     for (index = 0 ; index <trackSegments->GetEntriesFast() ; index++) {
534       AliPHOSTrackSegment * ts = (AliPHOSTrackSegment * )trackSegments->At(index) ; 
535       cout<<"   "<< setw(4) << ts->GetIndexInList() << "            " 
536           <<setw(4) << ts->GetEmcIndex()<< "            " 
537           <<setw(4) << ts->GetCpvIndex()<< "            " << endl ;
538     }   
539     
540     cout << "-------------------------------------------------------"<< endl ;
541   }
542 }