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