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