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