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