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