1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //_________________________________________________________________________
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.
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
41 // // One can specify the title for each branch
42 // root [2] r->SetBranchFileName("RecPoints","RecPoints1") ;
44 // // One can change parameters of reconstruction algorithms
45 // root [3] r->GetClusterizer()->SetEmcLocalMaxCut(0.02)
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
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 ///////////////////////////////////////////////////////////////////////////////////////////////////
60 // --- ROOT system ---
62 #include "TClonesArray.h"
66 // --- Standard library ---
69 // --- AliRoot header files ---
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"
80 ClassImp(AliPHOSReconstructioner)
82 //____________________________________________________________________________
83 AliPHOSReconstructioner::AliPHOSReconstructioner():TTask("AliPHOSReconstructioner","")
91 fHeaderFileName = "galice.root" ;
93 fIsInitialized = kFALSE ;
97 //____________________________________________________________________________
98 AliPHOSReconstructioner::AliPHOSReconstructioner(const char* headerFile):TTask("AliPHOSReconstructioner","")
102 fHeaderFileName = headerFile ;
105 fSDigitizer = new AliPHOSSDigitizer(fHeaderFileName.Data(),fSDigitsBranch.Data()) ;
109 fDigitizer = new AliPHOSDigitizer(fHeaderFileName.Data(),fDigitsBranch.Data()) ;
114 fClusterizer = new AliPHOSClusterizerv1(fHeaderFileName.Data(),fRecPointBranch.Data()) ;
119 fTSMaker = new AliPHOSTrackSegmentMakerv1(fHeaderFileName.Data(),fTSBranch.Data()) ;
124 fPID = new AliPHOSPIDv1(fHeaderFileName.Data(),fRecPartBranch.Data()) ;
127 fIsInitialized = kTRUE ;
130 //____________________________________________________________________________
131 void AliPHOSReconstructioner::Exec(Option_t *option){
132 //chesk, if the names of branches, which should be made conicide with already
137 gAlice->GetEvent(0) ;
139 if(fSDigitizer->IsActive()&& gAlice->TreeS()){ //Will produce SDigits
141 TBranch * sdigitsBranch = 0;
142 TBranch * sdigitizerBranch = 0;
144 TObjArray * branches = gAlice->TreeS()->GetListOfBranches() ;
146 Bool_t phosNotFound = kTRUE ;
147 Bool_t sdigitizerNotFound = kTRUE ;
149 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
151 sdigitsBranch=(TBranch *) branches->At(ibranch) ;
152 if(( strcmp("PHOS",sdigitsBranch->GetName())==0 ) &&
153 (fSDigitsBranch.CompareTo(sdigitsBranch->GetTitle())== 0 ))
154 phosNotFound = kFALSE ;
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 ;
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 ;
170 //mark all tasks as inactive
173 while((task=(TTask*)next()))
174 task->SetActive(kFALSE) ;
180 if(fDigitizer->IsActive() && gAlice->TreeD()){ //Will produce Digits
181 TBranch * digitsBranch = 0;
182 TBranch * digitizerBranch = 0;
184 TObjArray * branches = gAlice->TreeD()->GetListOfBranches() ;
186 Bool_t phosNotFound = kTRUE ;
187 Bool_t digitizerNotFound = kTRUE ;
189 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
191 digitsBranch=(TBranch *) branches->At(ibranch) ;
192 if(( strcmp("PHOS",digitsBranch->GetName())==0 ) &&
193 (fDigitsBranch.CompareTo(digitsBranch->GetTitle())== 0 ))
194 phosNotFound = kFALSE ;
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 ;
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 ;
210 //mark all tasks as inactive
213 while((task=(TTask*)next()))
214 task->SetActive(kFALSE) ;
220 if(fClusterizer->IsActive() && gAlice->TreeR()){ //Will produce RecPoints
221 TBranch * emcBranch = 0;
222 TBranch * cpvBranch = 0;
223 TBranch * clusterizerBranch = 0;
225 TObjArray * branches = gAlice->TreeR()->GetListOfBranches() ;
227 Bool_t emcNotFound = kTRUE ;
228 Bool_t cpvNotFound = kTRUE ;
229 Bool_t clusterizerNotFound = kTRUE ;
231 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
234 emcBranch=(TBranch *) branches->At(ibranch) ;
235 if(fRecPointBranch.CompareTo(emcBranch->GetTitle())==0 )
236 if( strcmp(emcBranch->GetName(),"PHOSEmcRP") == 0)
237 emcNotFound = kFALSE ;
240 cpvBranch=(TBranch *) branches->At(ibranch) ;
241 if(fRecPointBranch.CompareTo(cpvBranch->GetTitle())==0 )
242 if( strcmp(cpvBranch->GetName(),"PHOSCpvRP") == 0)
243 cpvNotFound = kFALSE ;
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 ;
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 ;
260 //mark all tasks as inactive
263 while((task=(TTask*)next()))
264 task->SetActive(kFALSE) ;
269 if(fTSMaker->IsActive() && gAlice->TreeR()){ //Produce TrackSegments
271 TBranch * tsMakerBranch = 0;
272 TBranch * tsBranch = 0;
274 TObjArray * branches = gAlice->TreeR()->GetListOfBranches() ;
276 Bool_t tsMakerNotFound = kTRUE ;
277 Bool_t tsNotFound = kTRUE ;
279 for(ibranch = 0;(ibranch <branches->GetEntries())&&(tsMakerNotFound||tsNotFound);ibranch++){
281 tsMakerBranch=(TBranch *) branches->At(ibranch) ;
282 if( fTSBranch.CompareTo(tsMakerBranch->GetTitle())==0 )
283 if( strcmp(tsMakerBranch->GetName(),"AliPHOSTrackSegmentMaker") == 0)
284 tsMakerNotFound = kFALSE ;
287 tsBranch=(TBranch *) branches->At(ibranch) ;
288 if( fTSBranch.CompareTo(tsBranch->GetTitle())==0 )
289 if( strcmp(tsBranch->GetName(),"PHOSTS") == 0)
290 tsNotFound = kFALSE ;
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 ;
301 //mark all tasks as inactive
304 while((task=(TTask*)next()))
305 task->SetActive(kFALSE) ;
312 if(fPID->IsActive() && gAlice->TreeR()){ //Produce RecParticles
313 TBranch * pidBranch = 0;
314 TBranch * rpBranch = 0;
316 TObjArray * branches = gAlice->TreeR()->GetListOfBranches() ;
318 Bool_t pidNotFound = kTRUE ;
319 Bool_t rpNotFound = kTRUE ;
321 for(ibranch = 0;(ibranch <branches->GetEntries()) && pidNotFound && rpNotFound ;ibranch++){
323 pidBranch=(TBranch *) branches->At(ibranch) ;
324 if( (strcmp(fRecPartBranch,pidBranch->GetTitle())==0 ) &&
325 (strcmp(pidBranch->GetName(),"AliPHOSPID") == 0) )
326 pidNotFound = kFALSE ;
329 rpBranch=(TBranch *) branches->At(ibranch) ;
330 if( (strcmp(fRecPartBranch,rpBranch->GetTitle())==0 ) &&
331 (strcmp(rpBranch->GetName(),"PHOSRP") == 0) )
332 rpNotFound = kFALSE ;
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 ;
343 //mark all tasks as inactive
346 while((task=(TTask*)next()))
347 task->SetActive(kFALSE) ;
353 //____________________________________________________________________________
354 void AliPHOSReconstructioner::Init()
356 //initiase Reconstructioner if necessary: we can not do this in default constructor
362 fSDigitizer = new AliPHOSSDigitizer(fHeaderFileName.Data(),fSDigitsBranch.Data()) ;
366 fDigitizer = new AliPHOSDigitizer(fHeaderFileName.Data(),fDigitsBranch.Data()) ;
370 fClusterizer = new AliPHOSClusterizerv1(fHeaderFileName.Data(),fRecPointBranch.Data()) ;
374 fTSMaker = new AliPHOSTrackSegmentMakerv1(fHeaderFileName.Data(),fTSBranch.Data()) ;
379 fPID = new AliPHOSPIDv1(fHeaderFileName.Data(),fRecPartBranch.Data()) ;
382 fIsInitialized = kTRUE ;
385 //____________________________________________________________________________
386 AliPHOSReconstructioner::~AliPHOSReconstructioner()
396 delete fClusterizer ;
404 //____________________________________________________________________________
405 void AliPHOSReconstructioner::SetBranchTitle(const char* branch, const char * title){
406 //Diverge correcpoinding branch to the file "title"
408 if(strcmp(branch,"SDigits") == 0){
409 fSDigitizer->SetSDigitsBranch(title) ;
410 fDigitizer->SetSDigitsBranch(title) ;
411 fSDigitsBranch = title ;
415 if(strcmp(branch,"Digits") == 0){
416 fDigitizer->SetDigitsBranch(title) ;
417 fClusterizer->SetDigitsBranch(title) ;
418 fDigitsBranch = title ;
422 if(strcmp(branch,"RecPoints") == 0){
423 fClusterizer->SetRecPointsBranch(title) ;
424 fTSMaker->SetRecPointsBranch(title) ;
425 fRecPointBranch = title ;
429 if(strcmp(branch,"TrackSegments") == 0){
430 fTSMaker->SetTrackSegmentsBranch(title) ;
431 fPID->SetTrackSegmentsBranch(title) ;
436 if(strcmp(branch,"RecParticles") == 0){
437 fPID->SetRecParticlesBranch(title) ;
438 fRecPartBranch = title ;
442 cout << "There is no branch " << branch << "!"<< endl ;
443 cout << "Available branches `SDigits', `Digits', `RecPoints', `TrackSegments' and `RecParticles' " << endl ;
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
460 char * moduleName = new char[30];
461 if(strstr(module,"SD"))
462 sprintf(moduleName,"AliPHOSSDigitizer") ;
464 if(strstr(module,"D") )
465 sprintf(moduleName,"AliPHOSDigitizer") ;
467 if(strstr(module,"C") || strstr(module,"RecPoint") )
468 sprintf(moduleName,"AliPHOSClusterizer") ;
470 if(strstr(module,"TS") || strstr(module,"Track") )
471 sprintf(moduleName,"AliPHOSTrackSegmentMaker") ;
473 if(strstr(module,"PID") || strstr(module,"Particle") || strstr(module,"RP") )
474 sprintf(moduleName,"AliPHOSPID") ;
476 cout << "Do not know such a module / Rec Object " << endl;
482 Bool_t active = kFALSE ;
483 while((task=(TTask*)next())){
484 if (strcmp(moduleName,task->GetName())==0)
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 ;
494 case 16: //"AliPHOSDigitizer"
495 fDigitizer->SetDigitsBranch(title) ;
496 fClusterizer->SetDigitsBranch(title) ;
497 fDigitsBranch = title ;
499 case 18: //"AliPHOSClusterizer"
500 fClusterizer->SetRecPointsBranch(title) ;
501 fTSMaker->SetRecPointsBranch(title) ;
502 fRecPointBranch = title ;
504 case 24: //"AliPHOSTrackSegmentMaker"
505 fTSMaker->SetTrackSegmentsBranch(title) ;
506 fPID->SetTrackSegmentsBranch(title) ;
509 case 10: // "AliPHOSPID"
510 fPID->SetRecParticlesBranch(title) ;
511 fRecPartBranch = title ;
520 //____________________________________________________________________________
521 void AliPHOSReconstructioner::Print(Option_t * option)const {
523 cout << "-----------------AliPHOSReconstructioner---------------" << endl ;
524 cout << " Reconstruction of the header file " <<fHeaderFileName.Data() << endl ;
525 cout << " with the following modules: " << endl ;
527 if(fSDigitizer->IsActive()){
528 cout << " (+) " << fSDigitizer->GetName() << " to branch : " << fSDigitsBranch.Data() << endl ;
531 if(fDigitizer->IsActive()){
532 cout << " (+) " << fDigitizer->GetName() << " to branch : " << fDigitsBranch.Data() << endl ;
536 if(fClusterizer->IsActive()){
537 cout << " (+) " <<fClusterizer->GetName() << " to branch : " <<fRecPointBranch.Data() << endl ;
541 if(fTSMaker->IsActive()){
542 cout << " (+) " << fTSMaker->GetName() << " to branch : " << fTSBranch.Data() << endl ;
547 if(fPID->IsActive()){
548 cout << " (+) " << fPID->GetName() << " to branch : " <<fRecPartBranch.Data() << endl ;