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