]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSReconstructioner.cxx
Changed drwaing routines.
[u/mrichter/AliRoot.git] / PHOS / AliPHOSReconstructioner.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
16 /* $Id$ */
17
18 //_________________________________________________________________________
19 //*--
20 //*-- Author: Gines Martinez & Yves Schutz (SUBATECH) 
21 //*-- Complitely redisigned by Dmitri Peressounko (SUBATECH & RRC KI) March 2001
22 /////////////////////////////////////////////////////////////////////////////////////
23 //  Wrapping class for reconstruction. Allows to produce reconstruction from 
24 //  different steps: from previously produced hits,sdigits, etc. Each new reconstruction
25 //  flow (e.g. gigits, made from them RecPoints, made from them TrackSegments, made from 
26 //  them RecParticles) are distinguished by the title of created branches. One can 
27 //  use this title as a comment, see user case below. Thanks to getters, one can set 
28 //  parameters to reconstruction bricks. The full set of parameters is saved in the 
29 //  correspoinding branch: e.g. parameters of clusterizer are stored in branch 
30 //  TreeR::AliPHOSClusterizer with the same title as the branch contaning RecPoints 
31 //  themself. TTree does not support overwriting, therefore one can not produce several 
32 //  branches with the same names and titles - use different titles.
33 //
34 //  User case: 
35 //
36 //  root [0] AliPHOSReconstructioner * r = new AliPHOSReconstructioner("galice.root")
37 //              //  Set the header file
38 //  root [1] r->ExecuteTask() 
39 //              //  Make full cheine of reconstruction
40 //
41 //              // One can specify the title for each branch 
42 //  root [2] r->SetBranchFileName("RecPoints","RecPoints1") ;
43 //      
44 //             // One can change parameters of reconstruction algorithms
45 //  root [3] r->GetClusterizer()->SetEmcLocalMaxCut(0.02)
46 //
47 //             // One can specify the starting point of the reconstruction and title of all 
48 //             // branches produced in this pass
49 //  root [4] r->StartFrom("AliPHOSClusterizer","Local max cut 0.02") 
50 //             // means that will use already generated Digits and produce only RecPoints, 
51 //             // TS and RecParticles 
52 //
53 //             // And finally one can call ExecuteTask() with the following options
54 //  root [5] r->ExecuteTask("debug all timing")
55 //             // deb     - prints the numbers of produced SDigits, Digits etc.
56 //             // deb all - prints in addition list of made SDigits, digits etc.
57 //             // timing  - prints benchmarking results
58 ///////////////////////////////////////////////////////////////////////////////////////////////////
59
60 // --- ROOT system ---
61
62 #include "TClonesArray.h"
63 #include "TROOT.h"
64 #include "TTree.h"
65
66 // --- Standard library ---
67 #include <iostream.h>   
68
69 // --- AliRoot header files ---
70 #include "AliRun.h"
71 #include "AliPHOSReconstructioner.h"
72 #include "AliPHOSClusterizerv1.h"
73 #include "AliPHOSDigitizer.h"
74 #include "AliPHOSSDigitizer.h"
75 #include "AliPHOSTrackSegmentMakerv1.h"
76 #include "AliPHOSPIDv1.h"
77 #include "AliPHOSFastRecParticle.h"
78 #include "AliPHOSCpvRecPoint.h"
79
80 ClassImp(AliPHOSReconstructioner)
81
82 //____________________________________________________________________________
83   AliPHOSReconstructioner::AliPHOSReconstructioner():TTask("AliPHOSReconstructioner","")
84 {
85   // ctor
86   fDigitizer   = 0 ;
87   fClusterizer = 0 ;
88   fTSMaker     = 0 ;
89   fPID         = 0 ; 
90   fSDigitizer  = 0 ;
91   fHeaderFileName = "galice.root" ;
92
93   fIsInitialized = kFALSE ;
94
95
96
97 //____________________________________________________________________________
98 AliPHOSReconstructioner::AliPHOSReconstructioner(const char* headerFile):TTask("AliPHOSReconstructioner","")
99 {
100   // ctor
101   
102   fHeaderFileName = headerFile ;
103
104   fSDigitsBranch="" ; 
105   fSDigitizer  = new AliPHOSSDigitizer(fHeaderFileName.Data(),fSDigitsBranch.Data()) ; 
106   Add(fSDigitizer) ;
107
108   fDigitsBranch="" ; 
109   fDigitizer   = new AliPHOSDigitizer(fHeaderFileName.Data(),fDigitsBranch.Data()) ; 
110   Add(fDigitizer) ;
111
112
113   fRecPointBranch="" ; 
114   fClusterizer = new AliPHOSClusterizerv1(fHeaderFileName.Data(),fRecPointBranch.Data()) ; 
115   Add(fClusterizer) ;
116   
117
118   fTSBranch="" ; 
119   fTSMaker     = new AliPHOSTrackSegmentMakerv1(fHeaderFileName.Data(),fTSBranch.Data()) ;
120   Add(fTSMaker) ;
121   
122   
123   fRecPartBranch="" ; 
124   fPID         = new AliPHOSPIDv1(fHeaderFileName.Data(),fRecPartBranch.Data()) ;
125   Add(fPID) ;
126   
127   fIsInitialized = kTRUE ;
128   
129
130 //____________________________________________________________________________
131 void AliPHOSReconstructioner::Exec(Option_t *option){
132   //chesk, if the names of branches, which should be made conicide with already
133   //existing
134   if(!fIsInitialized)
135     Init() ;
136
137   gAlice->GetEvent(0) ;
138
139   if(fSDigitizer->IsActive()&& gAlice->TreeS()){ //Will produce SDigits
140
141     TBranch * sdigitsBranch = 0;
142     TBranch * sdigitizerBranch = 0;
143
144     TObjArray * branches = gAlice->TreeS()->GetListOfBranches() ;
145     Int_t ibranch;
146     Bool_t phosNotFound = kTRUE ;
147     Bool_t sdigitizerNotFound = kTRUE ;
148
149     for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){            
150       if(phosNotFound){
151         sdigitsBranch=(TBranch *) branches->At(ibranch) ;
152         if(( strcmp("PHOS",sdigitsBranch->GetName())==0 ) &&
153            (fSDigitsBranch.CompareTo(sdigitsBranch->GetTitle())== 0 ))
154           phosNotFound = kFALSE ;
155       }
156       if(sdigitizerNotFound){
157         sdigitizerBranch = (TBranch *) branches->At(ibranch) ;
158         if(( strcmp(sdigitizerBranch->GetName(),"AliPHOSSDigitizer") == 0) &&
159            (fSDigitsBranch.CompareTo(sdigitizerBranch->GetTitle())== 0 ) )
160           sdigitizerNotFound = kFALSE ;
161       }
162     }
163     
164     if(!(sdigitizerNotFound && phosNotFound)){
165       cout << "AliPHOSReconstructioner error: "<< endl ;
166       cout << "       Branches ''PHOS'' or ''AliPHOSSDigitizer'' with title ``" << fSDigitsBranch.Data() << "''" << endl ;
167       cout << "       already exist in TreeS. ROOT does not allow updating/overwriting." << endl ;
168       cout << "       Specify another title for branches or use ''StartFrom()'' method" << endl ;
169       
170       //mark all tasks as inactive
171       TIter next(fTasks);
172       TTask *task;
173       while((task=(TTask*)next()))
174         task->SetActive(kFALSE) ;
175       
176       return ;
177     }
178   }
179
180   if(fDigitizer->IsActive() && gAlice->TreeD()){ //Will produce Digits
181     TBranch * digitsBranch = 0;
182     TBranch * digitizerBranch = 0;
183     
184     TObjArray * branches = gAlice->TreeD()->GetListOfBranches() ;
185     Int_t ibranch;
186     Bool_t phosNotFound = kTRUE ;
187     Bool_t digitizerNotFound = kTRUE ;
188     
189     for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){            
190       if(phosNotFound){
191         digitsBranch=(TBranch *) branches->At(ibranch) ;
192         if(( strcmp("PHOS",digitsBranch->GetName())==0 ) &&
193            (fDigitsBranch.CompareTo(digitsBranch->GetTitle())== 0 ))
194           phosNotFound = kFALSE ;
195       }
196       if(digitizerNotFound){
197         digitizerBranch = (TBranch *) branches->At(ibranch) ;
198         if(( strcmp(digitizerBranch->GetName(),"AliPHOSDigitizer") == 0) &&
199            (fDigitsBranch.CompareTo(digitizerBranch->GetTitle())== 0 ) )
200           digitizerNotFound = kFALSE ;
201       }
202     }
203     
204     if(!(digitizerNotFound && phosNotFound)){
205       cout << "AliPHOSReconstructioner error: "<< endl ;
206       cout << "       Branches ''PHOS'' or ''AliPHOSDigitizer'' with title ``" << fDigitsBranch.Data() << "''" << endl ;
207       cout << "       already exist in TreeD. ROOT does not allow updating/overwriting." << endl ;
208       cout << "       Specify another title for branches or use ''StartFrom()'' method" << endl ;
209       
210       //mark all tasks as inactive
211       TIter next(fTasks);
212       TTask *task;
213       while((task=(TTask*)next()))
214         task->SetActive(kFALSE) ;
215       
216       return ;
217     }
218   }
219
220   if(fClusterizer->IsActive() && gAlice->TreeR()){ //Will produce RecPoints
221     TBranch * emcBranch = 0;
222     TBranch * cpvBranch = 0;
223     TBranch * clusterizerBranch = 0;
224     
225     TObjArray * branches = gAlice->TreeR()->GetListOfBranches() ;
226     Int_t ibranch;
227     Bool_t emcNotFound = kTRUE ;
228     Bool_t cpvNotFound = kTRUE ;  
229     Bool_t clusterizerNotFound = kTRUE ;
230     
231     for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
232       
233       if(emcNotFound){
234         emcBranch=(TBranch *) branches->At(ibranch) ;
235         if(fRecPointBranch.CompareTo(emcBranch->GetTitle())==0 )
236           if( strcmp(emcBranch->GetName(),"PHOSEmcRP") == 0) 
237             emcNotFound = kFALSE ;
238       }
239       if(cpvNotFound){
240         cpvBranch=(TBranch *) branches->At(ibranch) ;
241         if(fRecPointBranch.CompareTo(cpvBranch->GetTitle())==0 )
242           if( strcmp(cpvBranch->GetName(),"PHOSCpvRP") == 0) 
243             cpvNotFound = kFALSE ;
244       }
245       if(clusterizerNotFound){
246         clusterizerBranch = (TBranch *) branches->At(ibranch) ;
247         if( fRecPointBranch.CompareTo(clusterizerBranch->GetTitle()) == 0)
248           if( strcmp(clusterizerBranch->GetName(),"AliPHOSClusterizer") == 0) 
249             clusterizerNotFound = kFALSE ;
250       }
251     }
252
253     if(!(clusterizerNotFound && emcNotFound && cpvNotFound)){
254       cout << "AliPHOSReconstructioner error: "<< endl ;
255       cout << "       Branches ''PHOSEmcRP'', ''PHOSCpvRP'' or ''AliPHOSClusterizer'' with title ``" 
256            << fRecPointBranch.Data() << "''" << endl ;
257       cout << "       already exist in TreeR. ROOT does not allow updating/overwriting." << endl ;
258       cout << "       Specify another title for branches or use ''StartFrom()'' method" << endl ;
259       
260       //mark all tasks as inactive
261       TIter next(fTasks);
262       TTask *task;
263       while((task=(TTask*)next()))
264         task->SetActive(kFALSE) ;
265       return ;
266     }
267   }
268   
269   if(fTSMaker->IsActive() && gAlice->TreeR()){ //Produce TrackSegments
270
271     TBranch * tsMakerBranch = 0;
272     TBranch * tsBranch = 0;
273     
274     TObjArray * branches = gAlice->TreeR()->GetListOfBranches() ;
275     Int_t ibranch;
276     Bool_t tsMakerNotFound = kTRUE ;
277     Bool_t tsNotFound = kTRUE ;
278     
279     for(ibranch = 0;(ibranch <branches->GetEntries())&&(tsMakerNotFound||tsNotFound);ibranch++){
280       if(tsMakerNotFound){
281         tsMakerBranch=(TBranch *) branches->At(ibranch) ;
282         if( fTSBranch.CompareTo(tsMakerBranch->GetTitle())==0 )
283           if( strcmp(tsMakerBranch->GetName(),"AliPHOSTrackSegmentMaker") == 0) 
284             tsMakerNotFound = kFALSE ;
285       }
286       if(tsNotFound){
287         tsBranch=(TBranch *) branches->At(ibranch) ;
288         if( fTSBranch.CompareTo(tsBranch->GetTitle())==0 )
289           if( strcmp(tsBranch->GetName(),"PHOSTS") == 0) 
290             tsNotFound = kFALSE ;
291       }
292     }
293     
294     if(!(tsMakerNotFound &&tsNotFound) ){
295       cout << "AliPHOSReconstructioner error: "<< endl ;
296       cout << "       Branches ''PHOSTS'' or ''AliPHOSTrackSegmentMaker'' with title ``" 
297            << fTSBranch.Data() << "''" << endl ;
298       cout << "       already exist in TreeR. ROOT does not allow updating/overwriting." << endl ;
299       cout << "       Specify another title for branches or use ''StartFrom()'' method" << endl ;
300       
301       //mark all tasks as inactive
302       TIter next(fTasks);
303       TTask *task;
304       while((task=(TTask*)next()))
305         task->SetActive(kFALSE) ;
306       return ;
307       
308     }
309     
310   }
311
312   if(fPID->IsActive() && gAlice->TreeR()){ //Produce RecParticles
313     TBranch * pidBranch = 0;
314     TBranch * rpBranch = 0;
315     
316     TObjArray * branches = gAlice->TreeR()->GetListOfBranches() ;
317     Int_t ibranch;
318     Bool_t pidNotFound = kTRUE ;
319     Bool_t rpNotFound = kTRUE ;
320     
321     for(ibranch = 0;(ibranch <branches->GetEntries()) && pidNotFound && rpNotFound ;ibranch++){
322       if(pidNotFound){
323         pidBranch=(TBranch *) branches->At(ibranch) ;
324         if( (strcmp(fRecPartBranch,pidBranch->GetTitle())==0 ) &&
325             (strcmp(pidBranch->GetName(),"AliPHOSPID") == 0) )
326           pidNotFound = kFALSE ;
327       }
328       if(rpNotFound){
329         rpBranch=(TBranch *) branches->At(ibranch) ;
330         if( (strcmp(fRecPartBranch,rpBranch->GetTitle())==0 ) &&
331             (strcmp(rpBranch->GetName(),"PHOSRP") == 0) )
332           rpNotFound = kFALSE ;
333       }
334     }
335     
336     if(!pidNotFound  || !rpNotFound ){
337       cout << "AliPHOSReconstructioner error: "<< endl ;
338       cout << "       Branches ''PHOSRP'' or ''AliPHOSPID'' with title ``" 
339            << fRecPartBranch.Data() << "''" << endl ;
340       cout << "       already exist in TreeR. ROOT does not allow updating/overwriting." << endl ;
341       cout << "       Specify another title for branches." << endl ;
342       
343       //mark all tasks as inactive
344       TIter next(fTasks);
345       TTask *task;
346       while((task=(TTask*)next()))
347         task->SetActive(kFALSE) ;
348       return ;
349     }
350     
351   }
352 }
353 //____________________________________________________________________________
354  void AliPHOSReconstructioner::Init()
355 {
356   //initiase Reconstructioner if necessary: we can not do this in default constructor
357
358   if(!fIsInitialized){
359     // Initialisation
360
361     fSDigitsBranch="" ; 
362     fSDigitizer  = new AliPHOSSDigitizer(fHeaderFileName.Data(),fSDigitsBranch.Data()) ; 
363     Add(fSDigitizer) ;
364
365     fDigitsBranch="" ; 
366     fDigitizer   = new AliPHOSDigitizer(fHeaderFileName.Data(),fDigitsBranch.Data()) ; 
367     Add(fDigitizer) ;
368
369     fRecPointBranch="" ; 
370     fClusterizer = new AliPHOSClusterizerv1(fHeaderFileName.Data(),fRecPointBranch.Data()) ; 
371     Add(fClusterizer) ;
372
373     fTSBranch="" ; 
374     fTSMaker     = new AliPHOSTrackSegmentMakerv1(fHeaderFileName.Data(),fTSBranch.Data()) ;
375     Add(fTSMaker) ;
376
377
378     fRecPartBranch="" ; 
379     fPID         = new AliPHOSPIDv1(fHeaderFileName.Data(),fRecPartBranch.Data()) ;
380     Add(fPID) ;
381     
382     fIsInitialized = kTRUE ;
383   }
384
385 //____________________________________________________________________________
386 AliPHOSReconstructioner::~AliPHOSReconstructioner()
387 {
388   
389   if(fSDigitizer)
390     delete fSDigitizer ;
391   
392   if(fDigitizer)
393     delete fDigitizer ;
394   
395   if(fClusterizer)
396     delete fClusterizer ;
397   
398   if(fTSMaker)
399     delete fTSMaker ;
400   
401   if(fPID)
402     delete fPID ;
403
404 //____________________________________________________________________________
405 void AliPHOSReconstructioner::SetBranchTitle(const char* branch, const char * title){
406   //Diverge correcpoinding branch to the file "title"
407
408   if(strcmp(branch,"SDigits") == 0){ 
409     fSDigitizer->SetSDigitsBranch(title) ;
410     fDigitizer->SetSDigitsBranch(title) ;
411     fSDigitsBranch = title ;
412     return ;
413   }
414   
415   if(strcmp(branch,"Digits") == 0){ 
416     fDigitizer->SetDigitsBranch(title) ;
417     fClusterizer->SetDigitsBranch(title) ;
418     fDigitsBranch = title ;
419     return ;
420   }
421
422   if(strcmp(branch,"RecPoints") == 0){ 
423     fClusterizer->SetRecPointsBranch(title) ;
424     fTSMaker->SetRecPointsBranch(title) ;
425     fRecPointBranch = title ;
426     return ;
427   }
428
429   if(strcmp(branch,"TrackSegments") == 0){
430     fTSMaker->SetTrackSegmentsBranch(title) ;
431     fPID->SetTrackSegmentsBranch(title) ;
432     fTSBranch = title ;
433     return ;
434   }
435
436   if(strcmp(branch,"RecParticles") == 0){ 
437     fPID->SetRecParticlesBranch(title) ;
438     fRecPartBranch = title ;
439     return ;
440   }
441
442   cout << "There is no branch " << branch << "!"<< endl ;
443   cout << "Available branches `SDigits', `Digits', `RecPoints', `TrackSegments' and `RecParticles' " << endl ;
444   
445 }
446 //____________________________________________________________________________
447 void AliPHOSReconstructioner::StartFrom(char * module,char* title){
448   // in the next pass of reconstruction (call ExecuteTask()) reconstruction will 
449   // start from the module "module", and in the case of nonsero title all 
450   // pruduced branches will have title "title". Recognize the following "modules"
451   // "SD" - AliPHOSSDigitizer,
452   // "D"  - AliPHOSDigitizer
453   // "C"  - AliPHOSClusterizer
454   // "TS" - AliPHOSTrackSegmentMaker
455   // "RP" - AliPHOSPID
456
457   if(!fIsInitialized)
458     Init() ;
459
460   char * moduleName = new char[30];
461   if(strstr(module,"SD"))
462     sprintf(moduleName,"AliPHOSSDigitizer") ;
463   else
464     if(strstr(module,"D") )
465       sprintf(moduleName,"AliPHOSDigitizer") ;
466     else
467       if(strstr(module,"C") || strstr(module,"RecPoint") )
468         sprintf(moduleName,"AliPHOSClusterizer") ;
469       else
470         if(strstr(module,"TS") || strstr(module,"Track") )
471           sprintf(moduleName,"AliPHOSTrackSegmentMaker") ;
472         else
473           if(strstr(module,"PID") || strstr(module,"Particle") || strstr(module,"RP") )
474             sprintf(moduleName,"AliPHOSPID") ;
475           else{
476             cout << "Do not know such a module / Rec Object " << endl;
477             return ;
478           }
479   
480   TIter next(fTasks);
481   TTask *task;
482   Bool_t active = kFALSE ;
483   while((task=(TTask*)next())){ 
484     if (strcmp(moduleName,task->GetName())==0)  
485       active = kTRUE;
486     task->SetActive(active) ;
487     if(active && title){ // set title to branches
488       switch(strlen(task->GetName()) ) {
489       case 17:   // "AliPHOSSDigitizer"
490         fSDigitizer->SetSDigitsBranch(title) ;
491         fDigitizer->SetSDigitsBranch(title) ;
492         fSDigitsBranch = title ;
493         break ;
494       case 16:   //"AliPHOSDigitizer"
495         fDigitizer->SetDigitsBranch(title) ;
496         fClusterizer->SetDigitsBranch(title) ;
497         fDigitsBranch = title ;
498         break ;
499       case 18:   //"AliPHOSClusterizer"
500         fClusterizer->SetRecPointsBranch(title) ;
501         fTSMaker->SetRecPointsBranch(title) ;
502         fRecPointBranch = title ;
503         break ;
504       case 24:   //"AliPHOSTrackSegmentMaker"
505         fTSMaker->SetTrackSegmentsBranch(title) ;
506         fPID->SetTrackSegmentsBranch(title) ;
507         fTSBranch = title ;
508         break ;
509       case 10:   // "AliPHOSPID"
510         fPID->SetRecParticlesBranch(title) ;
511         fRecPartBranch = title ;
512         break ;
513       }
514       
515     }
516   }
517   
518   delete moduleName;
519 }
520 //____________________________________________________________________________
521 void AliPHOSReconstructioner::Print(Option_t * option)const {
522   
523   cout << "-----------------AliPHOSReconstructioner---------------" << endl ;
524   cout << " Reconstruction of the header file " <<fHeaderFileName.Data() << endl ;
525   cout << " with the following modules: " << endl ;
526
527   if(fSDigitizer->IsActive()){
528     cout << "   (+)   " << fSDigitizer->GetName() << " to branch : " << fSDigitsBranch.Data() << endl ; 
529     cout << endl ;
530   }
531   if(fDigitizer->IsActive()){
532     cout << "   (+)   " << fDigitizer->GetName() << " to branch : " << fDigitsBranch.Data() << endl ;  
533     cout <<  endl ;
534   }
535   
536   if(fClusterizer->IsActive()){
537     cout << "   (+)   " <<fClusterizer->GetName() << " to branch : " <<fRecPointBranch.Data()  << endl ;  
538     cout <<  endl ;
539   }
540
541   if(fTSMaker->IsActive()){
542     cout << "   (+)   " << fTSMaker->GetName() << " to branch : " << fTSBranch.Data() << endl ;  
543     cout <<  endl ;
544   }
545
546
547   if(fPID->IsActive()){
548     cout << "   (+)   " << fPID->GetName() << " to branch : " <<fRecPartBranch.Data()  << endl ;  
549     cout <<  endl ;
550   }
551
552
553 }