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