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 //*-- 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.
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
42 // // One can specify the title for each branch
43 // root [2] r->SetBranchFileName("RecPoints","RecPoints1") ;
45 // // One can change parameters of reconstruction algorithms
46 // root [3] r->GetClusterizer()->SetEmcLocalMaxCut(0.02)
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
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 ///////////////////////////////////////////////////////////////////////////////////////////////////
61 // --- ROOT system ---
63 #include "TClonesArray.h"
67 // --- Standard library ---
70 // --- AliRoot header files ---
72 #include "AliPHOSReconstructioner.h"
73 #include "AliPHOSClusterizerv1.h"
74 #include "AliPHOSDigitizer.h"
75 #include "AliPHOSSDigitizer.h"
76 #include "AliPHOSTrackSegmentMakerv1.h"
77 #include "AliPHOSPIDv1.h"
78 #include "AliPHOSFastRecParticle.h"
79 #include "AliPHOSCpvRecPoint.h"
81 ClassImp(AliPHOSReconstructioner)
83 //____________________________________________________________________________
84 AliPHOSReconstructioner::AliPHOSReconstructioner():TTask("AliPHOSReconstructioner","")
92 fHeaderFileName = "galice.root" ;
94 fIsInitialized = kFALSE ;
98 //____________________________________________________________________________
99 AliPHOSReconstructioner::AliPHOSReconstructioner(const char* headerFile):TTask("AliPHOSReconstructioner","")
103 fHeaderFileName = headerFile ;
106 fSDigitizer = new AliPHOSSDigitizer(fHeaderFileName.Data(),fSDigitsBranch.Data()) ;
110 fDigitizer = new AliPHOSDigitizer(fHeaderFileName.Data(),fDigitsBranch.Data()) ;
115 fClusterizer = new AliPHOSClusterizerv1(fHeaderFileName.Data(),fRecPointBranch.Data()) ;
120 fTSMaker = new AliPHOSTrackSegmentMakerv1(fHeaderFileName.Data(),fTSBranch.Data()) ;
125 fPID = new AliPHOSPIDv1(fHeaderFileName.Data(),fRecPartBranch.Data()) ;
128 fIsInitialized = kTRUE ;
131 //____________________________________________________________________________
132 void AliPHOSReconstructioner::Exec(Option_t *option)
134 //chesk, if the names of branches, which should be made conicide with already
139 gAlice->GetEvent(0) ;
141 if(fSDigitizer->IsActive()&& gAlice->TreeS()){ //Will produce SDigits
143 TBranch * sdigitsBranch = 0;
144 TBranch * sdigitizerBranch = 0;
146 TObjArray * branches = gAlice->TreeS()->GetListOfBranches() ;
148 Bool_t phosNotFound = kTRUE ;
149 Bool_t sdigitizerNotFound = kTRUE ;
151 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
153 sdigitsBranch=(TBranch *) branches->At(ibranch) ;
154 if(( strcmp("PHOS",sdigitsBranch->GetName())==0 ) &&
155 (fSDigitsBranch.CompareTo(sdigitsBranch->GetTitle())== 0 ))
156 phosNotFound = kFALSE ;
158 if(sdigitizerNotFound){
159 sdigitizerBranch = (TBranch *) branches->At(ibranch) ;
160 if(( strcmp(sdigitizerBranch->GetName(),"AliPHOSSDigitizer") == 0) &&
161 (fSDigitsBranch.CompareTo(sdigitizerBranch->GetTitle())== 0 ) )
162 sdigitizerNotFound = kFALSE ;
166 if(!(sdigitizerNotFound && phosNotFound)){
167 cout << "AliPHOSReconstructioner error: "<< endl ;
168 cout << " Branches ''PHOS'' or ''AliPHOSSDigitizer'' with title ``" << fSDigitsBranch.Data() << "''" << endl ;
169 cout << " already exist in TreeS. ROOT does not allow updating/overwriting." << endl ;
170 cout << " Specify another title for branches or use ''StartFrom()'' method" << endl ;
172 //mark all tasks as inactive
175 while((task=(TTask*)next()))
176 task->SetActive(kFALSE) ;
182 if(fDigitizer->IsActive() && gAlice->TreeD()){ //Will produce Digits
183 TBranch * digitsBranch = 0;
184 TBranch * digitizerBranch = 0;
186 TObjArray * branches = gAlice->TreeD()->GetListOfBranches() ;
188 Bool_t phosNotFound = kTRUE ;
189 Bool_t digitizerNotFound = kTRUE ;
191 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
193 digitsBranch=(TBranch *) branches->At(ibranch) ;
194 if(( strcmp("PHOS",digitsBranch->GetName())==0 ) &&
195 (fDigitsBranch.CompareTo(digitsBranch->GetTitle())== 0 ))
196 phosNotFound = kFALSE ;
198 if(digitizerNotFound){
199 digitizerBranch = (TBranch *) branches->At(ibranch) ;
200 if(( strcmp(digitizerBranch->GetName(),"AliPHOSDigitizer") == 0) &&
201 (fDigitsBranch.CompareTo(digitizerBranch->GetTitle())== 0 ) )
202 digitizerNotFound = kFALSE ;
206 if(!(digitizerNotFound && phosNotFound)){
207 cout << "AliPHOSReconstructioner error: "<< endl ;
208 cout << " Branches ''PHOS'' or ''AliPHOSDigitizer'' with title ``" << fDigitsBranch.Data() << "''" << endl ;
209 cout << " already exist in TreeD. ROOT does not allow updating/overwriting." << endl ;
210 cout << " Specify another title for branches or use ''StartFrom()'' method" << endl ;
212 //mark all tasks as inactive
215 while((task=(TTask*)next()))
216 task->SetActive(kFALSE) ;
222 if(fClusterizer->IsActive() && gAlice->TreeR()){ //Will produce RecPoints
223 TBranch * emcBranch = 0;
224 TBranch * cpvBranch = 0;
225 TBranch * clusterizerBranch = 0;
227 TObjArray * branches = gAlice->TreeR()->GetListOfBranches() ;
229 Bool_t emcNotFound = kTRUE ;
230 Bool_t cpvNotFound = kTRUE ;
231 Bool_t clusterizerNotFound = kTRUE ;
233 for(ibranch = 0;ibranch <branches->GetEntries();ibranch++){
236 emcBranch=(TBranch *) branches->At(ibranch) ;
237 if(fRecPointBranch.CompareTo(emcBranch->GetTitle())==0 )
238 if( strcmp(emcBranch->GetName(),"PHOSEmcRP") == 0)
239 emcNotFound = kFALSE ;
242 cpvBranch=(TBranch *) branches->At(ibranch) ;
243 if(fRecPointBranch.CompareTo(cpvBranch->GetTitle())==0 )
244 if( strcmp(cpvBranch->GetName(),"PHOSCpvRP") == 0)
245 cpvNotFound = kFALSE ;
247 if(clusterizerNotFound){
248 clusterizerBranch = (TBranch *) branches->At(ibranch) ;
249 if( fRecPointBranch.CompareTo(clusterizerBranch->GetTitle()) == 0)
250 if( strcmp(clusterizerBranch->GetName(),"AliPHOSClusterizer") == 0)
251 clusterizerNotFound = kFALSE ;
255 if(!(clusterizerNotFound && emcNotFound && cpvNotFound)){
256 cout << "AliPHOSReconstructioner error: "<< endl ;
257 cout << " Branches ''PHOSEmcRP'', ''PHOSCpvRP'' or ''AliPHOSClusterizer'' with title ``"
258 << fRecPointBranch.Data() << "''" << endl ;
259 cout << " already exist in TreeR. ROOT does not allow updating/overwriting." << endl ;
260 cout << " Specify another title for branches or use ''StartFrom()'' method" << endl ;
262 //mark all tasks as inactive
265 while((task=(TTask*)next()))
266 task->SetActive(kFALSE) ;
271 if(fTSMaker->IsActive() && gAlice->TreeR()){ //Produce TrackSegments
273 TBranch * tsMakerBranch = 0;
274 TBranch * tsBranch = 0;
276 TObjArray * branches = gAlice->TreeR()->GetListOfBranches() ;
278 Bool_t tsMakerNotFound = kTRUE ;
279 Bool_t tsNotFound = kTRUE ;
281 for(ibranch = 0;(ibranch <branches->GetEntries())&&(tsMakerNotFound||tsNotFound);ibranch++){
283 tsMakerBranch=(TBranch *) branches->At(ibranch) ;
284 if( fTSBranch.CompareTo(tsMakerBranch->GetTitle())==0 )
285 if( strcmp(tsMakerBranch->GetName(),"AliPHOSTrackSegmentMaker") == 0)
286 tsMakerNotFound = kFALSE ;
289 tsBranch=(TBranch *) branches->At(ibranch) ;
290 if( fTSBranch.CompareTo(tsBranch->GetTitle())==0 )
291 if( strcmp(tsBranch->GetName(),"PHOSTS") == 0)
292 tsNotFound = kFALSE ;
296 if(!(tsMakerNotFound &&tsNotFound) ){
297 cout << "AliPHOSReconstructioner error: "<< endl ;
298 cout << " Branches ''PHOSTS'' or ''AliPHOSTrackSegmentMaker'' with title ``"
299 << fTSBranch.Data() << "''" << endl ;
300 cout << " already exist in TreeR. ROOT does not allow updating/overwriting." << endl ;
301 cout << " Specify another title for branches or use ''StartFrom()'' method" << endl ;
303 //mark all tasks as inactive
306 while((task=(TTask*)next()))
307 task->SetActive(kFALSE) ;
314 if(fPID->IsActive() && gAlice->TreeR()){ //Produce RecParticles
315 TBranch * pidBranch = 0;
316 TBranch * rpBranch = 0;
318 TObjArray * branches = gAlice->TreeR()->GetListOfBranches() ;
320 Bool_t pidNotFound = kTRUE ;
321 Bool_t rpNotFound = kTRUE ;
323 for(ibranch = 0;(ibranch <branches->GetEntries()) && pidNotFound && rpNotFound ;ibranch++){
325 pidBranch=(TBranch *) branches->At(ibranch) ;
326 if( (strcmp(fRecPartBranch,pidBranch->GetTitle())==0 ) &&
327 (strcmp(pidBranch->GetName(),"AliPHOSPID") == 0) )
328 pidNotFound = kFALSE ;
331 rpBranch=(TBranch *) branches->At(ibranch) ;
332 if( (strcmp(fRecPartBranch,rpBranch->GetTitle())==0 ) &&
333 (strcmp(rpBranch->GetName(),"PHOSRP") == 0) )
334 rpNotFound = kFALSE ;
338 if(!pidNotFound || !rpNotFound ){
339 cout << "AliPHOSReconstructioner error: "<< endl ;
340 cout << " Branches ''PHOSRP'' or ''AliPHOSPID'' with title ``"
341 << fRecPartBranch.Data() << "''" << endl ;
342 cout << " already exist in TreeR. ROOT does not allow updating/overwriting." << endl ;
343 cout << " Specify another title for branches." << endl ;
345 //mark all tasks as inactive
348 while((task=(TTask*)next()))
349 task->SetActive(kFALSE) ;
355 //____________________________________________________________________________
356 void AliPHOSReconstructioner::Init()
358 // initiliaze Reconstructioner if necessary: we can not do this in default constructor
364 fSDigitizer = new AliPHOSSDigitizer(fHeaderFileName.Data(),fSDigitsBranch.Data()) ;
368 fDigitizer = new AliPHOSDigitizer(fHeaderFileName.Data(),fDigitsBranch.Data()) ;
372 fClusterizer = new AliPHOSClusterizerv1(fHeaderFileName.Data(),fRecPointBranch.Data()) ;
376 fTSMaker = new AliPHOSTrackSegmentMakerv1(fHeaderFileName.Data(),fTSBranch.Data()) ;
381 fPID = new AliPHOSPIDv1(fHeaderFileName.Data(),fRecPartBranch.Data()) ;
384 fIsInitialized = kTRUE ;
387 //____________________________________________________________________________
388 AliPHOSReconstructioner::~AliPHOSReconstructioner()
398 delete fClusterizer ;
406 //____________________________________________________________________________
407 void AliPHOSReconstructioner::SetBranchTitle(const char* branch, const char * title){
408 //Diverge correcpoinding branch to the file "title"
410 if(strcmp(branch,"SDigits") == 0){
411 fSDigitizer->SetSDigitsBranch(title) ;
412 fDigitizer->SetSDigitsBranch(title) ;
413 fSDigitsBranch = title ;
417 if(strcmp(branch,"Digits") == 0){
418 fDigitizer->SetDigitsBranch(title) ;
419 fClusterizer->SetDigitsBranch(title) ;
420 fDigitsBranch = title ;
424 if(strcmp(branch,"RecPoints") == 0){
425 fClusterizer->SetRecPointsBranch(title) ;
426 fTSMaker->SetRecPointsBranch(title) ;
427 fRecPointBranch = title ;
431 if(strcmp(branch,"TrackSegments") == 0){
432 fTSMaker->SetTrackSegmentsBranch(title) ;
433 fPID->SetTrackSegmentsBranch(title) ;
438 if(strcmp(branch,"RecParticles") == 0){
439 fPID->SetRecParticlesBranch(title) ;
440 fRecPartBranch = title ;
444 cout << "There is no branch " << branch << "!"<< endl ;
445 cout << "Available branches `SDigits', `Digits', `RecPoints', `TrackSegments' and `RecParticles' " << endl ;
448 //____________________________________________________________________________
449 void AliPHOSReconstructioner::StartFrom(char * module,char* title)
451 // in the next pass of reconstruction (call ExecuteTask()) reconstruction will
452 // start from the module "module", and in the case of non zero title all
453 // pruduced branches will have title "title". The following "modules" are recognized
454 // "SD" - AliPHOSSDigitizer,
455 // "D" - AliPHOSDigitizer
456 // "C" - AliPHOSClusterizer
457 // "TS" - AliPHOSTrackSegmentMaker
463 char * moduleName = new char[30];
464 if(strstr(module,"SD"))
465 sprintf(moduleName,"AliPHOSSDigitizer") ;
467 if(strstr(module,"D") )
468 sprintf(moduleName,"AliPHOSDigitizer") ;
470 if(strstr(module,"C") || strstr(module,"RecPoint") )
471 sprintf(moduleName,"AliPHOSClusterizer") ;
473 if(strstr(module,"TS") || strstr(module,"Track") )
474 sprintf(moduleName,"AliPHOSTrackSegmentMaker") ;
476 if(strstr(module,"PID") || strstr(module,"Particle") || strstr(module,"RP") )
477 sprintf(moduleName,"AliPHOSPID") ;
479 cout << "Do not know such a module / Rec Object " << endl;
485 Bool_t active = kFALSE ;
486 while((task=(TTask*)next())){
487 if (strcmp(moduleName,task->GetName())==0)
489 task->SetActive(active) ;
490 if(active && title){ // set title to branches
491 switch(strlen(task->GetName()) ) {
492 case 17: // "AliPHOSSDigitizer"
493 fSDigitizer->SetSDigitsBranch(title) ;
494 fDigitizer->SetSDigitsBranch(title) ;
495 fSDigitsBranch = title ;
497 case 16: //"AliPHOSDigitizer"
498 fDigitizer->SetDigitsBranch(title) ;
499 fClusterizer->SetDigitsBranch(title) ;
500 fDigitsBranch = title ;
502 case 18: //"AliPHOSClusterizer"
503 fClusterizer->SetRecPointsBranch(title) ;
504 fTSMaker->SetRecPointsBranch(title) ;
505 fRecPointBranch = title ;
507 case 24: //"AliPHOSTrackSegmentMaker"
508 fTSMaker->SetTrackSegmentsBranch(title) ;
509 fPID->SetTrackSegmentsBranch(title) ;
512 case 10: // "AliPHOSPID"
513 fPID->SetRecParticlesBranch(title) ;
514 fRecPartBranch = title ;
523 //____________________________________________________________________________
524 void AliPHOSReconstructioner::Print(Option_t * option)const
526 // Print the parameters of the reconstructioner
527 cout << "-----------------AliPHOSReconstructioner---------------" << endl ;
528 cout << " Reconstruction of the header file " <<fHeaderFileName.Data() << endl ;
529 cout << " with the following modules: " << endl ;
531 if(fSDigitizer->IsActive()){
532 cout << " (+) " << fSDigitizer->GetName() << " to branch : " << fSDigitsBranch.Data() << endl ;
535 if(fDigitizer->IsActive()){
536 cout << " (+) " << fDigitizer->GetName() << " to branch : " << fDigitsBranch.Data() << endl ;
540 if(fClusterizer->IsActive()){
541 cout << " (+) " <<fClusterizer->GetName() << " to branch : " <<fRecPointBranch.Data() << endl ;
545 if(fTSMaker->IsActive()){
546 cout << " (+) " << fTSMaker->GetName() << " to branch : " << fTSBranch.Data() << endl ;
551 if(fPID->IsActive()){
552 cout << " (+) " << fPID->GetName() << " to branch : " <<fRecPartBranch.Data() << endl ;