]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSGetter.cxx
Trigger pattern getter added
[u/mrichter/AliRoot.git] / PHOS / AliPHOSGetter.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 /* $Log:
19    29.05.2001 Yuri Kharlov:
20               Everywhere reading the treese TTree->GetEvent(i)
21               is replaced by reading the branches TBranch->GetEntry(0)
22 */
23 /* $Log:
24    08.2002 Dmitri Peressounko:
25
26 */
27
28 //_________________________________________________________________________
29 //  A singleton. This class should be used in the analysis stage to get 
30 //  reconstructed objects: Digits, RecPoints, TrackSegments and RecParticles,
31 //  instead of directly reading them from galice.root file. This container 
32 //  ensures, that one reads Digits, made of these particular digits, RecPoints, 
33 //  made of these particular RecPoints, TrackSegments and RecParticles. 
34 //  This becomes non trivial if there are several identical branches, produced with
35 //  different set of parameters. 
36 //
37 //  An example of how to use (see also class AliPHOSAnalyser):
38 //  AliPHOSGetter * gime = AliPHOSGetter::GetInstance("galice.root","test") ;
39 //  for(Int_t irecp = 0; irecp < gime->NRecParticles() ; irecp++)
40 //     AliPHOSRecParticle * part = gime->RecParticle(1) ;
41 //     ................
42 //  gime->Event(event) ;    // reads new event from galice.root
43 //                  
44 //*-- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (RRC KI & SUBATECH)
45 //*--         Completely redesigned by Dmitri Peressounko March 2001  
46 //
47 //*-- YS June 2001 : renamed the original AliPHOSIndexToObject and make
48 //*--         systematic usage of TFolders without changing the interface        
49 //////////////////////////////////////////////////////////////////////////////
50
51 // --- ROOT system ---
52 #include "TFile.h"
53 #include "TTree.h"
54 #include "TROOT.h"
55 #include "TObjString.h"
56 #include "TFolder.h"
57 #include "TParticle.h"
58
59 // --- Standard library ---
60
61 // --- AliRoot header files ---
62 #include "AliRun.h"
63 #include "AliConfig.h"
64 #include "AliPHOSGetter.h"
65 #include "AliPHOS.h"
66 #include "AliPHOSDigitizer.h"
67 #include "AliPHOSSDigitizer.h"
68 #include "AliPHOSClusterizerv1.h"
69 #include "AliPHOSTrackSegmentMakerv1.h"
70 #include "AliPHOSTrackSegment.h"
71 #include "AliPHOSPIDv1.h" 
72 #include "AliPHOSGeometry.h"
73 #include "AliPHOSRaw2Digits.h"
74 //#include "AliPHOSCalibrationDB.h"
75 #include "AliPHOSBeamTestEvent.h"
76 ClassImp(AliPHOSGetter)
77   
78   AliPHOSGetter * AliPHOSGetter::fgObjGetter = 0 ; 
79   TFile * AliPHOSGetter::fFile = 0 ; 
80
81 //____________________________________________________________________________ 
82 AliPHOSGetter::AliPHOSGetter(const char* headerFile, const char* branchTitle, const Bool_t toSplit )
83 {
84   // This is the ctor called by GetInstance and the only one that can be used 
85
86   if( fHeaderFile.Contains("_") ) {
87     Fatal("AliPHOSGetter", "Invalid file name (_ not allowed) %s", fHeaderFile.Data() ) ;
88   }
89
90   //Initialize  all data
91
92   fFailed = kFALSE ;   
93   fDebug  = 0 ; 
94   fAlice  = 0 ; 
95   fBTE    = 0 ;
96
97   fToSplit    = toSplit ;
98   fHeaderFile = headerFile ; 
99
100   fPrimaries = new TObjArray(1) ;
101
102   fModuleFolder    = dynamic_cast<TFolder*>(gROOT->FindObjectAny("Folders/Run/Configuration/Modules")); 
103   fPrimariesFolder = dynamic_cast<TFolder*>(gROOT->FindObjectAny("Folders/RunMC/Event/Data")); 
104   fHitsFolder      = dynamic_cast<TFolder*>(gROOT->FindObjectAny("Folders/RunMC/Event/Data/Hits")); 
105   fSDigitsFolder   = dynamic_cast<TFolder*>(gROOT->FindObjectAny("Folders/RunMC/Event/Data/SDigits")); 
106   fDigitsFolder    = dynamic_cast<TFolder*>(gROOT->FindObjectAny("Folders/Run/Event/Data")); 
107   fRecoFolder      = dynamic_cast<TFolder*>(gROOT->FindObjectAny("Folders/Run/Event/RecData")); 
108   fQAFolder        = dynamic_cast<TFolder*>(gROOT->FindObjectAny("Folders/Run/Conditions/QA")); 
109   fTasksFolder     = dynamic_cast<TFolder*>(gROOT->FindObjectAny("Folders/Tasks")) ; 
110
111   //Set titles to branches and create PHOS specific folders
112   SetTitle(branchTitle) ;
113
114   if ( fHeaderFile != "aliroot"  ) { // to call the getter without a file
115     //open headers file
116     fFile = static_cast<TFile*>(gROOT->GetFile(fHeaderFile.Data() ) ) ;
117
118     if(!fFile) {    //if file was not opened yet, read gAlice
119       fFile = TFile::Open(fHeaderFile.Data(),"update") ;
120       if (!fFile->IsOpen()) {
121         Error("AliPHOSGetter", "Cannot open %s", fHeaderFile.Data() ) ; 
122         fFailed = kTRUE ;
123         return ;  
124       }
125     }
126     gAlice = dynamic_cast<AliRun *>(fFile->Get("gAlice")) ;
127   }
128   
129   if (!gAlice) {
130     Error("AliPHOSGetter", "Cannot find gAlice in %s", fHeaderFile.Data() ) ; 
131     fFailed = kTRUE ;
132     return ; 
133   }
134   if (!PHOS()) {
135     if (fDebug)
136       Info("AliPHOSGetter", "-> Posting PHOS to Folders") ; 
137     if (gAlice->GetDetector("PHOS")) {
138       AliConfig * conf = AliConfig::Instance() ; 
139       conf->Add(static_cast<AliDetector*>(gAlice->GetDetector("PHOS"))) ; 
140       conf->Add(static_cast<AliModule*>(gAlice->GetDetector("PHOS"))) ; 
141     }
142     else 
143       Error("AliPHOSGetter", "detector PHOS not found") ;    
144   }
145
146   fDebug=0;
147 }
148
149 //____________________________________________________________________________ 
150 AliPHOSGetter::~AliPHOSGetter(){
151
152   if (fPrimaries) {
153     fPrimaries->Delete() ; 
154     delete fPrimaries ; 
155   }
156
157   TFolder * phosF = dynamic_cast<TFolder *>(fSDigitsFolder->FindObject("PHOS")) ;
158   TCollection * folderslist = phosF->GetListOfFolders() ; 
159   TIter next(folderslist) ; 
160   TFolder * folder = 0 ; 
161   while ( (folder = static_cast<TFolder*>(next())) ) 
162     phosF->Remove(folder) ; 
163   
164   if (fFile) { 
165     fFile->Close() ;  
166     delete fFile ;
167     fFile = 0 ;
168   }
169   fgObjGetter = 0 ; 
170   
171 }
172
173 //____________________________________________________________________________ 
174 void AliPHOSGetter::CloseFile()
175 {
176   if(gAlice)
177     delete gAlice ;  
178   gAlice = 0 ; 
179   if(fAlice)
180     delete fAlice ; 
181   fAlice = 0 ; 
182 }
183
184 //____________________________________________________________________________ 
185 const TFolder * AliPHOSGetter::Folder(const TString what) const {
186
187   // returns the PHOS folder required by what
188   // what = hits, sdigits, digits
189
190   if ( what == "hits" ) 
191     return dynamic_cast<const TFolder *>(fHitsFolder->FindObject("PHOS")) ; 
192   else if ( what == "sdigits" ) 
193     return  dynamic_cast<const TFolder *>(fSDigitsFolder->FindObject("PHOS")) ; 
194   else if ( what == "digits" ) 
195     return  dynamic_cast<const TFolder *>(fDigitsFolder->FindObject("PHOS")) ; 
196   else {
197     Error("GetFolder", " %s illegal option (hits, sdigits, digits)", what.Data() ) ; 
198     return 0 ; 
199   }
200 }
201
202 //____________________________________________________________________________ 
203 AliPHOSGetter * AliPHOSGetter::GetInstance()
204 {
205   // Returns the pointer of the unique instance already defined
206   
207   if ( fgObjGetter ) {
208     return fgObjGetter ;
209   }
210   else {
211     //Warning("GetInstance", "not yet initialized") ;
212     return 0 ; 
213   }
214 }
215
216 //____________________________________________________________________________ 
217 AliPHOSGetter * AliPHOSGetter::GetInstance(const char* headerFile,
218                                            const char* branchTitle,
219                                            const Bool_t toSplit)
220 {
221   // Creates and returns the pointer of the unique instance
222   // Must be called only when the environment has changed
223
224   if(!fgObjGetter){
225     fgObjGetter = new AliPHOSGetter(headerFile,branchTitle,toSplit) ;
226     if(fgObjGetter->fFailed)
227       return 0;
228     else
229       return fgObjGetter ;
230   }
231
232   //First checks, if header file already opened
233   if(!fgObjGetter->fFile){
234      fgObjGetter = new AliPHOSGetter(headerFile,branchTitle,toSplit) ;
235     if(fgObjGetter->fFailed)
236       return 0;
237     else
238       return fgObjGetter ;
239   }
240
241   if(fgObjGetter->fHeaderFile.CompareTo(headerFile)==0){ //Opened the same header file
242     if((fgObjGetter->fBranchTitle.CompareTo(branchTitle) == 0)&&   //Open the same branch title
243        (toSplit==fgObjGetter->fToSplit)){                          //Nothing should be cleaned
244     }
245     else{ //Clean all data and AliPHOS...zers
246       if(fgObjGetter->fToSplit)
247         fgObjGetter->CloseSplitFiles() ;
248       fgObjGetter->CleanWhiteBoard() ;
249       fgObjGetter->fToSplit = toSplit ;
250       fgObjGetter->SetTitle(branchTitle) ;
251     }
252   }
253   else{  //Close already opened files, clean memory and open new header file
254     if(gAlice){ //should first delete gAlice, then close file
255       //Should be in dtor of PHOS, but if one changes path ...
256       fgObjGetter->fModuleFolder->Remove(fgObjGetter->fModuleFolder->FindObject("PHOS")) ; 
257       delete gAlice ;      
258     }
259     if(fgObjGetter->fFile){
260       fgObjGetter->fFile->Close() ;
261       fgObjGetter->fFile=0;
262     }
263     if(fgObjGetter->fToSplit)
264       fgObjGetter->CloseSplitFiles() ;
265     fgObjGetter->CleanWhiteBoard() ;    
266     fgObjGetter = new AliPHOSGetter(headerFile,branchTitle,toSplit) ;
267
268   }
269   return fgObjGetter ; 
270   
271 }
272
273 //____________________________________________________________________________ 
274 const Bool_t AliPHOSGetter::BranchExists(const TString recName) const
275 {
276   //Looks in the tree Tree"name" if branch with current name olready exists
277
278   TString filename("") ;
279   TString name, dataname, zername ;
280   if(recName == "SDigits"){
281     filename=fSDigitsFileName ;
282     name = "TreeS0" ;
283     dataname = "PHOS" ;
284     zername = "AliPHOSSDigitizer" ;
285   }
286   else
287     if(recName == "Digits"){
288       filename=fDigitsFileName ;
289       name = "TreeD0" ;
290       dataname = "PHOS" ;
291       zername = "AliPHOSDigitizer" ;
292     }
293     else
294       if(recName == "RecPoints"){
295         filename=fRecPointsFileName ;
296         name = "TreeR0" ;
297         dataname = "PHOSEmcRP" ;
298         zername = "AliPHOSClusterizer" ;
299       }
300       else
301         if(recName == "TrackSegments"){
302           filename=fTrackSegmentsFileName ;
303           name = "TreeR0" ;
304           dataname = "PHOSTS" ;
305           zername = "AliPHOSTrackSegmentMaker" ;
306         }        
307         else
308           if(recName == "RecParticles"){
309             filename= fRecParticlesFileName ;
310             name = "TreeR0" ;
311             dataname = "PHOSRP" ;
312             zername = "AliPHOSPID" ;
313           }
314           else
315             return kFALSE ;
316
317   TFile * file ;
318   TTree * tree ;
319   if(fToSplit){
320     file = static_cast<TFile*>(gROOT->GetFile(filename.Data() ) ) ;
321     if(!file)
322       file = TFile::Open(fSDigitsFileName.Data(),"update");
323   }
324   else
325     file = fFile ;
326
327   tree = (TTree *)file->Get(name.Data()) ;
328   
329   if(!tree ) 
330     return kFALSE ;
331
332   TObjArray * lob = static_cast<TObjArray*>(tree->GetListOfBranches()) ;
333   TIter next(lob) ; 
334   TBranch * branch = 0 ;  
335   TString titleName(fBranchTitle);
336   titleName+=":";
337   while ((branch = (static_cast<TBranch*>(next())))) {
338     TString branchName(branch->GetName() ) ; 
339     TString branchTitle(branch->GetTitle() ) ;  
340     if ( branchName.BeginsWith(dataname) && branchTitle.BeginsWith(fBranchTitle) ){  
341       Warning("BranchExists", "branch %s with title %s already exists in %s", dataname.Data(), fBranchTitle.Data(), name.Data() ) ;
342       return kTRUE ;
343     }
344     if ( branchName.BeginsWith(zername) &&  branchTitle.BeginsWith(titleName) ){
345       Warning("BranchExists", "branch AliPHOS... with title %s already exists in %s", branch->GetTitle(), name.Data() ) ;     
346       return kTRUE ; 
347     }
348   }
349   //We can't delete three if gAlice points to it... To be redisigned somehow???!!!
350   if(!fToSplit){
351     if(name.Contains("TreeS"))
352       if(tree!=gAlice->TreeS())
353         tree->Delete();
354     if(name.Contains("TreeD"))
355       if(tree!=gAlice->TreeD())
356         tree->Delete();
357     if(name.Contains("TreeR"))
358       if(tree!=gAlice->TreeR())
359         tree->Delete();    
360   }
361   return kFALSE ;
362   
363 }
364
365 //____________________________________________________________________________ 
366 void AliPHOSGetter::ListBranches(Int_t event) const  
367 {
368   
369   TBranch * branch = 0 ; 
370   if (gAlice->GetEvent(event) == -1)
371     return ; 
372   
373   TTree * t =  gAlice->TreeH() ; 
374   if(t){
375     Info("ListBranches", "-> ****** Hits    : ") ; 
376     TObjArray * lob = t->GetListOfBranches() ;
377     TIter next(lob) ; 
378     while ( (branch = static_cast<TBranch*>(next())) )
379       Info("ListBranches", "             %s", branch->GetName()) ; 
380   } else
381     Warning("ListBranches", "TreeH not found for event %d", event ) ;  
382   
383   t = gAlice->TreeS() ;
384   if(t){
385     Info("ListBranches", "-> ****** SDigits : ") ; 
386     TObjArray * lob = t->GetListOfBranches() ;
387     TIter next(lob) ; 
388     while ( (branch = static_cast<TBranch*>(next())) )
389       Info("ListBranches", "             %s %s", branch->GetName(), branch->GetTitle()) ; 
390   } else 
391     Warning("ListBranches", "TreeS not found for event %d", event)  ;  
392   
393   
394   t = gAlice->TreeD() ;
395   if(t){
396     Info("ListBranches", "-> ****** Digits  : ") ; 
397     TObjArray * lob = t->GetListOfBranches() ;
398     TIter next(lob) ; 
399     while ( (branch = static_cast<TBranch*>(next())) )
400       Info("ListBranches", "             %s %s", branch->GetName(), branch->GetTitle()) ; 
401   } else 
402     Warning("ListBranches", "TreeD not found for event %d", event) ;  
403   
404
405   t = gAlice->TreeR() ;
406   if(t){
407     Info("ListBranches", "-> ****** Recon   : ") ; 
408     TObjArray * lob = t->GetListOfBranches() ;
409     TIter next(lob) ; 
410     while ( (branch = static_cast<TBranch*>(next())) )
411      Info("ListBranches", "             %s %s", branch->GetName(), branch->GetTitle()) ; 
412   } else 
413     Warning("ListBranches", "TreeR not found for event %d", event) ;  
414   
415 }
416
417 //____________________________________________________________________________ 
418 void AliPHOSGetter::NewBranch(TString name, Int_t event)  
419 {
420   fBranchTitle = fSDigitsTitle = fDigitsTitle = fRecPointsTitle = fTrackSegmentsTitle = fRecParticlesTitle =  name ; 
421   Event(event) ; 
422 }
423
424 //____________________________________________________________________________ 
425 Bool_t AliPHOSGetter::NewFile(TString name)  
426 {
427   fHeaderFile = name ; 
428   fFile->Close() ; 
429   fFailed = kFALSE; 
430
431   fFile = static_cast<TFile*>(gROOT->GetFile(fHeaderFile.Data() ) ) ;
432   if(!fFile) {    //if file was not opened yet, read gAlice
433     fFile = TFile::Open(fHeaderFile.Data(),"update") ;
434     if (!fFile->IsOpen()) {
435       Error("NewFile", "Cannot open %s", fHeaderFile.Data() ) ; 
436       fFailed = kTRUE ;
437       return fFailed ;  
438     }
439     gAlice = static_cast<AliRun *>(fFile->Get("gAlice")) ;
440   } 
441   
442   if (!gAlice) {
443     Error("AliPHOSGetter", "Cannot find gAlice in %s", fHeaderFile.Data() ) ; 
444     fFailed = kTRUE ;
445     return fFailed ; 
446   }
447   return fFailed ; 
448 }
449
450 //____________________________________________________________________________ 
451 const AliPHOS * AliPHOSGetter::PHOS() 
452 {
453   // returns the PHOS object 
454   AliPHOS * phos = dynamic_cast<AliPHOS*>(fModuleFolder->FindObject("PHOS")) ;  
455   if (!phos) 
456     if (fDebug)
457       Warning("PHOS", "-> PHOS module not found in Folders") ; 
458   return phos ; 
459 }  
460
461 //____________________________________________________________________________ 
462 const AliPHOSGeometry * AliPHOSGetter::PHOSGeometry() 
463 {
464   AliPHOSGeometry * rv = 0 ; 
465   if (PHOS() )
466     rv =  PHOS()->GetGeometry() ;
467   return rv ; 
468
469
470 //____________________________________________________________________________ 
471 const Bool_t AliPHOSGetter::PostPrimaries(void) const 
472 {  //------- Primaries ----------------------
473
474   // the hierarchy is //Folders/RunMC/Event/Data/Primaries
475   
476   TFolder * primariesFolder = dynamic_cast<TFolder*>(fPrimariesFolder->FindObject("Primaries")) ; 
477   if ( !primariesFolder ) {
478     if (fDebug) {
479       Warning("PostPrimaries", "-> Folder //%s/Primaries/ not found!", fPrimariesFolder->GetName()) ;
480       Info("PostPrimaries", "-> Adding Folder //%s/Primaries", fPrimariesFolder->GetName()) ;
481     }
482     primariesFolder = fPrimariesFolder->AddFolder("Primaries", "Primaries particles from TreeK") ; 
483   }    
484   TClonesArray *primaries=  new TClonesArray("TParticle",1000) ;
485   primaries->SetName("Primaries") ;
486   primariesFolder->Add(primaries) ; 
487   
488   return kTRUE;
489
490
491 //____________________________________________________________________________ 
492 TObject** AliPHOSGetter::PrimariesRef(void) const 
493 {  //------- Primaries ----------------------
494
495   
496   // the hierarchy is //Folders/RunMC/Event/Data/Primaries
497   if ( !fPrimariesFolder ) {
498     Fatal("PrimariesRef", "Folder //%s not found", fPrimariesFolder) ;
499   }    
500  
501   TFolder * primariesFolder = dynamic_cast<TFolder *>(fPrimariesFolder->FindObject("Primaries")) ;
502   if ( !primariesFolder ) {
503     Fatal("PrimariesRef", "Folder //%s/Primaries/ not found", fPrimariesFolder) ;  
504   }
505  
506   TObject * p = primariesFolder->FindObject("Primaries") ;
507   if(!p) {
508     Fatal("PrimariesRef","%s /Primaries not found !", primariesFolder->GetName() ) ; 
509   }
510
511   return primariesFolder->GetListOfFolders()->GetObjectRef(p) ;
512 }
513
514 //____________________________________________________________________________ 
515 const Bool_t AliPHOSGetter::PostHits(void) const 
516 {  //------- Hits ----------------------
517
518   // the hierarchy is //Folders/RunMC/Event/Data/PHOS/Hits
519   
520   TFolder * phosFolder = dynamic_cast<TFolder*>(fHitsFolder->FindObject("PHOS")) ; 
521   if ( !phosFolder ) {
522     if (fDebug) {
523       Warning("PostHits", "-> Folder //%s/PHOS/ not found!", fHitsFolder) ;
524       Info("PostHits", "-> Adding Folder //%s/PHOS/", fHitsFolder) ;
525     }
526     phosFolder = fHitsFolder->AddFolder("PHOS", "Hits from PHOS") ; 
527   }    
528   TClonesArray *hits=  new TClonesArray("AliPHOSHit",1000) ;
529   hits->SetName("Hits") ;
530   phosFolder->Add(hits) ; 
531   
532   return kTRUE;
533
534
535 //____________________________________________________________________________ 
536 TObject** AliPHOSGetter::HitsRef(void) const 
537 {  //------- Hits ----------------------
538
539   
540   // the hierarchy is //Folders/RunMC/Event/Data/PHOS/Hits
541   if ( !fHitsFolder ) {
542     Fatal("HitsRef", "Folder //%s not found !", fHitsFolder) ;
543   }    
544  
545   TFolder * phosFolder = dynamic_cast<TFolder *>(fHitsFolder->FindObject("PHOS")) ;
546   if ( !phosFolder ) {
547     Fatal("HitsRef", "Folder //%s/PHOS/ not found !", fHitsFolder) ;  
548   }
549  
550   TObject * h = phosFolder->FindObject("Hits") ;
551   if(!h) {
552     Fatal("HitsRef", "%s/Hits not fount !", phosFolder->GetName() ) ; 
553   }
554
555   return phosFolder->GetListOfFolders()->GetObjectRef(h) ;
556 }
557
558 //____________________________________________________________________________ 
559 const Bool_t AliPHOSGetter::PostSDigits(const char * name, const char * headerFile) const 
560 {  //---------- SDigits -------------------------
561
562   
563   // the hierarchy is //Folders/RunMC/Event/Data/PHOS/SDigits/headerFile/sdigitsname
564   // because you can have sdigits from several hit files for mixing
565   
566   TFolder * phosFolder = dynamic_cast<TFolder*>(fSDigitsFolder->FindObject("PHOS")) ;
567   if ( !phosFolder ) {
568     if (fDebug) {
569       Warning("PostSDigits", "-> Folder //%s/PHOS/ not found!", fSDigitsFolder) ;
570       Info("PostSDigits", "-> Adding Folder //%s/PHOS/", fHitsFolder) ;
571     }
572     phosFolder = fSDigitsFolder->AddFolder("PHOS", "SDigits from PHOS") ; 
573   }    
574
575   
576   TString subdir(headerFile) ;
577   if(fToSplit){
578     subdir.Remove(subdir.Last('/')+1,subdir.Length()) ;
579     subdir.ReplaceAll("/","_") ; 
580     subdir+="PHOS.SDigits." ;
581     if(name && (strcmp(name,"Default")!=0)){
582       subdir+=name ;
583       subdir+="." ;
584     }
585     subdir+="root" ;
586   }
587     
588   TFolder * phosSubFolder = dynamic_cast<TFolder*>(phosFolder->FindObject(subdir)) ; 
589   if ( !phosSubFolder ) 
590     phosSubFolder = phosFolder->AddFolder(subdir, ""); 
591   
592
593   TObject * sd  = phosSubFolder->FindObject(name); 
594   if ( !sd ) {
595     TClonesArray * sdigits = new TClonesArray("AliPHOSDigit",1) ;
596     sdigits->SetName(name) ;
597     phosSubFolder->Add(sdigits) ;
598   }
599   
600   return kTRUE;
601
602 //____________________________________________________________________________ 
603 TObject** AliPHOSGetter::SDigitsRef(const char * name, const char * foldername) const 
604 {  //------- SDigits ----------------------
605   
606   // the hierarchy is //Folders/RunMC/Event/Data/PHOS/SDigits/filename/SDigits
607
608   if ( !fSDigitsFolder ) {
609     Fatal("SDigitsRef", "Folder //%s not found !", fSDigitsFolder) ;
610   }    
611  
612   TFolder * phosFolder = static_cast<TFolder *>(fSDigitsFolder->FindObject("PHOS")) ;
613   if ( !phosFolder ) {
614     Fatal("SDigitsRef", "Folder //%s/PHOS not found !", fSDigitsFolder) ;
615   }
616
617   TFolder * phosSubFolder = 0 ;
618
619   phosSubFolder = dynamic_cast<TFolder *>(phosFolder->FindObject(foldername)) ;
620   
621   if(!phosSubFolder) {
622     Fatal("SDigitsRef", "Folder //Folders/RunMC/Event/Data/PHOS/%s not found !", foldername) ;
623   }
624
625   TObject * dis = phosSubFolder->FindObject(name) ;
626   if(!dis){
627     Fatal("SDigitsRef", "object %s not found !", name) ; 
628   }
629
630   return phosSubFolder->GetListOfFolders()->GetObjectRef(dis) ;
631
632 }
633
634 //____________________________________________________________________________ 
635 const Bool_t AliPHOSGetter::PostSDigitizer(AliPHOSSDigitizer * sdigitizer) const 
636 {  //---------- SDigitizer -------------------------
637     
638   // the hierarchy is //Folders/Tasks/SDigitizer/PHOS/sdigitsname
639
640
641   TTask * sd  = dynamic_cast<TTask*>(fTasksFolder->FindObject("SDigitizer")) ; 
642
643   if ( !sd ) {
644     Error("PostDigitizer", "Task //%s/SDigitizer not found !", fTasksFolder) ;
645     return kFALSE ;
646   }        
647   TTask * phos = dynamic_cast<TTask*>(sd->GetListOfTasks()->FindObject("PHOS")) ; 
648   if ( !phos )  {
649     if (fDebug) {
650       Warning("PostSDigitizer", "->//%s/SDigitizer/PHOS/ not found!", fTasksFolder) ;  
651       Info("PostSDigitizer", "-> Adding //%s/SDigitizer/PHOS/", fTasksFolder) ;
652     }
653     phos = new TTask("PHOS", "") ; 
654     sd->Add(phos) ; 
655   } 
656   AliPHOSSDigitizer * phossd  = dynamic_cast<AliPHOSSDigitizer *>(phos->GetListOfTasks()->FindObject( sdigitizer->GetName() )); 
657   if (phossd) { 
658     if (fDebug)
659       Info("PostSDigitizer", "-> Task %s already exists", sdigitizer->GetName()) ; 
660     phos->GetListOfTasks()->Remove(phossd) ;
661   }
662   phos->Add(sdigitizer) ;       
663   return kTRUE; 
664   
665 }
666
667 //____________________________________________________________________________ 
668 TObject** AliPHOSGetter::SDigitizerRef(const char * name) const 
669 {  
670
671   TTask * sd  = dynamic_cast<TTask*>(fTasksFolder->FindObject("SDigitizer")) ; 
672   if ( !sd ) {
673     Fatal("SDigitizerRef", "Task //%s/SDigitizer not found !", fTasksFolder) ;
674   }        
675
676   TTask * phos = dynamic_cast<TTask*>(sd->GetListOfTasks()->FindObject("PHOS")) ; 
677   if ( !phos )  {
678     Fatal("SDigitizerRef", "//%s/SDigitizer/PHOS not found !", fTasksFolder) ;
679   }        
680
681   TTask * task = dynamic_cast<TTask*>(phos->GetListOfTasks()->FindObject(name)) ; 
682
683   return phos->GetListOfTasks()->GetObjectRef(task) ;
684
685 }
686
687 //____________________________________________________________________________ 
688 const Bool_t AliPHOSGetter::PostSDigitizer(const char * name, const char * file) const 
689 {  //---------- SDigitizer -------------------------
690   
691  // the hierarchy is //Folders/Tasks/SDigitizer/PHOS/sdigitsname
692
693   TTask * sd  = dynamic_cast<TTask*>(fTasksFolder->FindObject("SDigitizer")) ; 
694   if ( !sd ) {
695     Error("PostSDigitizer", "Task //%s/SDigitizer not found !", fTasksFolder) ;
696     return kFALSE ;
697   }        
698
699   TTask * phos = dynamic_cast<TTask*>(sd->GetListOfTasks()->FindObject("PHOS")) ; 
700   if ( !phos )  {
701     if (fDebug) {
702       Error("PostSDigitizer", "->  //%s/SDigitizer/PHOS/ not found!", fTasksFolder) ;
703       Info("PostSDigitizer", "-> Adding  //%s/SDigitizer/PHOS", fTasksFolder) ;
704     }
705     phos = new TTask("PHOS", "") ; 
706     sd->Add(phos) ; 
707   } 
708
709   TString sdname(name) ;
710   sdname.Append(":") ;
711   sdname.Append(file);
712   sdname.ReplaceAll("/","_") ; 
713   AliPHOSSDigitizer * phossd  = dynamic_cast<AliPHOSSDigitizer *>(phos->GetListOfTasks()->FindObject( sdname )); 
714   if (!phossd) {
715     phossd = new AliPHOSSDigitizer() ;  
716     //Note, we can not call constructor with parameters: it will call Getter and scew up everething
717     phossd->SetName(sdname) ;
718     phossd->SetTitle(file) ;
719     phos->Add(phossd) ; 
720   }
721   return kTRUE; 
722   
723 }
724
725 //____________________________________________________________________________ 
726 const Bool_t AliPHOSGetter::PostDigits(const char * name) const 
727 {  //---------- Digits -------------------------
728
729   // the hierarchy is //Folders/Run/Event/Data/PHOS/SDigits/name
730
731   TFolder * phosFolder  = dynamic_cast<TFolder*>(fDigitsFolder->FindObject("PHOS")) ;
732
733   if ( !phosFolder ) {
734     if (fDebug) {
735       Warning("PostDigitizer", "-> Folder //%s/PHOS/ not found!", fDigitsFolder) ;
736       Info("PostDigitizer", "-> Adding Folder //%s/PHOS/", fDigitsFolder) ;
737     }
738     phosFolder = fDigitsFolder->AddFolder("PHOS", "Digits from PHOS") ;  
739   }    
740  
741   TObject*  dig = phosFolder->FindObject( name ) ;
742   if ( !dig ) {
743     TClonesArray * digits = new TClonesArray("AliPHOSDigit",1000) ;
744     digits->SetName(name) ;
745     phosFolder->Add(digits) ;  
746   }
747   return kTRUE; 
748 }
749
750 //____________________________________________________________________________ 
751 TObject** AliPHOSGetter::DigitsRef(const char * name) const 
752 { //------- Digits ----------------------
753   
754   // the hierarchy is //Folders/Run/Event/Data/PHOS/Digits/name
755
756   if ( !fDigitsFolder ) {
757     Fatal("DigitsRef", "Folder //%s not found !", fDigitsFolder) ;
758   }    
759   
760   TFolder * phosFolder  = dynamic_cast<TFolder*>(fDigitsFolder->FindObject("PHOS")) ; 
761   if ( !phosFolder ) {
762     Fatal("DigitsRef", "Folder //%s/PHOS/ not found !", fDigitsFolder) ;
763   }    
764
765   TObject * d = phosFolder->FindObject(name) ;
766   if(!d) {
767     Fatal("DigitsRef", "object %s not found !", name) ; 
768   }
769
770   return phosFolder->GetListOfFolders()->GetObjectRef(d) ;
771
772 }
773
774 //____________________________________________________________________________ 
775 const Bool_t AliPHOSGetter::PostDigitizer(AliPHOSDigitizer * digitizer) const 
776 {  //---------- Digitizer -------------------------
777   
778   TTask * sd  = dynamic_cast<TTask*>(fTasksFolder->FindObject("Digitizer")) ; 
779
780   if ( !sd ) {
781     Error("PostDigitizer", "Task //%s/Digitizer not found !", fTasksFolder) ;
782     return kFALSE ;
783   }        
784   TTask * phos = dynamic_cast<TTask*>(sd->GetListOfTasks()->FindObject("PHOS")) ; 
785   if ( !phos )  {
786     if (fDebug) {
787       Error("PostDigitizer", "//%s/Digitizer/PHOS not found!", fTasksFolder) ;
788       Info("PostDigitizer", "Adding //%s/Digitizer/PHOS", fTasksFolder) ; 
789     }
790     phos = new TTask("PHOS", "") ; 
791     sd->Add(phos) ; 
792   } 
793
794     AliPHOSDigitizer * phosd = dynamic_cast<AliPHOSDigitizer*>(phos->GetListOfTasks()->FindObject(digitizer->GetName())) ; 
795     if (phosd) { 
796       phosd->Delete() ;
797       phos->GetListOfTasks()->Remove(phosd) ;
798     }
799     phos->Add(digitizer) ; 
800     return kTRUE; 
801 }  
802
803 //____________________________________________________________________________ 
804 const Bool_t AliPHOSGetter::PostDigitizer(const char * name) const 
805 {  //---------- Digitizer -------------------------
806   
807  // the hierarchy is //Folders/Tasks/SDigitizer/PHOS/sdigitsname
808
809   TTask * d  = dynamic_cast<TTask*>(fTasksFolder->FindObject("Digitizer")) ; 
810   if ( !d ) {
811     Error("PostDigitizer", "Task //%s/Digitizer not found !", fTasksFolder) ;
812     return kFALSE ;
813   }        
814
815   TTask * phos = dynamic_cast<TTask*>(d->GetListOfTasks()->FindObject("PHOS")) ; 
816   if ( !phos )  {
817     if (fDebug) {
818       Warning("PostDigitizer", "//%s/Digitizer/PHOS not found!", fTasksFolder) ; 
819       Info("PostDigitizer", "Adding //%s/Digitizer/PHOS", fTasksFolder) ;
820     }
821     phos = new TTask("PHOS", "") ; 
822     d->Add(phos) ; 
823
824
825   TTask * phosd = dynamic_cast<TTask*>(phos->GetListOfTasks()->FindObject(fDigitsTitle)) ; 
826   if (!phosd) { 
827     if(strcmp(name, "Digitizer")==0){
828       phosd = new AliPHOSDigitizer() ;
829       phosd->SetName(fDigitsTitle) ;
830       phosd->SetTitle(fHeaderFile) ;
831       phos->Add(phosd) ;
832     } 
833     else{
834       phosd = new AliPHOSRaw2Digits() ;
835       phosd->SetName(fDigitsTitle) ;
836       phosd->SetTitle(fHeaderFile) ;
837       phos->Add(phosd) ;
838     }      
839   }
840   return kTRUE;  
841 }
842
843 //____________________________________________________________________________ 
844 TObject** AliPHOSGetter::DigitizerRef(const char * name) const 
845 {  
846   TTask * sd  = dynamic_cast<TTask*>(fTasksFolder->FindObject("Digitizer")) ; 
847   if ( !sd ) {
848     Fatal("DigitizerRef", "Task //%s/Digitizer not found !", fTasksFolder) ;
849   }        
850
851   TTask * phos = dynamic_cast<TTask*>(sd->GetListOfTasks()->FindObject("PHOS")) ; 
852   if ( !phos )  {
853     Fatal("DigitizerRef", "//%s/Digitizer/PHOS", fTasksFolder) ;
854   }        
855
856   TTask * task = dynamic_cast<TTask*>(phos->GetListOfTasks()->FindObject(name)) ; 
857
858   return phos->GetListOfTasks()->GetObjectRef(task) ;
859
860 }
861  
862 //____________________________________________________________________________ 
863 const Bool_t AliPHOSGetter::PostRecPoints(const char * name) const 
864 { // -------------- RecPoints -------------------------------------------
865   
866   // the hierarchy is //Folders/Run/Event/RecData/PHOS/EMCARecPoints/name
867   // the hierarchy is //Folders/Run/Event/RecData/PHOS/CPVRecPoints/name
868
869   TFolder * phosFolder  = dynamic_cast<TFolder*>(fRecoFolder->FindObject("PHOS")) ; 
870   
871   if ( !phosFolder ) {
872     if (fDebug) {
873       Warning("PostRecPoints", "-> Folder //%s/PHOS/ not found!", fRecoFolder->GetName()) ;
874       Info("PostRecPoints", "-> Adding Folder //%s/PHOS/", fRecoFolder->GetName()) ;
875     }
876     phosFolder = fRecoFolder->AddFolder("PHOS", "Reconstructed data from PHOS") ;  
877   }    
878   
879   // EMCA RecPoints 
880   TFolder * phosRPoEMCAFolder  = dynamic_cast<TFolder*>(phosFolder->FindObject("EMCARecPoints")) ;
881   if ( !phosRPoEMCAFolder ) {
882     if (fDebug) {
883       Warning("PostRecPoints", "-> Folder //%s/PHOS/EMCARecPoints/ not found!", fRecoFolder->GetName()) ;
884       Info("PostRecPoints", "-> Adding Folder //%s/PHOS/EMCARecPoints", fRecoFolder->GetName()) ;
885     }
886     phosRPoEMCAFolder = phosFolder->AddFolder("EMCARecPoints", "EMCA RecPoints from PHOS") ;  
887   }    
888   
889   TObject * erp = phosFolder->FindObject( name ) ;
890   if ( !erp )   {
891     TObjArray * emcrp = new TObjArray(100) ;
892     emcrp->SetName(name) ;
893     phosRPoEMCAFolder->Add(emcrp) ;  
894   }
895
896   // CPV RecPoints 
897   TFolder * phosRPoCPVFolder  = dynamic_cast<TFolder*>(phosFolder->FindObject("CPVRecPoints")) ;
898   if ( !phosRPoCPVFolder ) {
899     if (fDebug) {
900       Warning("PostRecPoints", "-> Folder //%s/PHOS/CPVRecPoints/ not found!", fRecoFolder->GetName()) ;
901       Info("PostRecPoints", "-> Adding Folder //%s/PHOS/CPVRecPoints/", fRecoFolder->GetName()) ;
902     }
903     phosRPoCPVFolder = phosFolder->AddFolder("CPVRecPoints", "CPV RecPoints from PHOS") ;  
904   }    
905   
906   TObject * crp =  phosRPoCPVFolder->FindObject( name ) ;
907   if ( !crp )   {
908     TObjArray * cpvrp = new TObjArray(100) ;
909     cpvrp->SetName(name) ;
910     phosRPoCPVFolder->Add(cpvrp) ;  
911   }
912   return kTRUE; 
913 }
914
915 //____________________________________________________________________________ 
916 TObject** AliPHOSGetter::EmcRecPointsRef(const char * name) const 
917 { // -------------- RecPoints -------------------------------------------
918   
919   // the hierarchy is //Folders/Run/Event/RecData/PHOS/EMCARecPoints/name
920    
921   if ( !fRecoFolder ) {
922     Fatal("EmcRecPointsRef", "Folder //%s not found !", fRecoFolder->GetName() ) ;
923   }    
924
925   TFolder * phosFolder  = dynamic_cast<TFolder*>(fRecoFolder->FindObject("PHOS/EMCARecPoints")) ; 
926   if ( !phosFolder ) {
927      Fatal("EmcRecPointsRef", "Folder //%s/PHOS/EMCARecPoints/ not found !", fRecoFolder->GetName() ) ;
928   }    
929
930
931   TObject * erp = phosFolder->FindObject(name ) ;
932   if ( !erp )   {
933     Fatal("EmcRecPointsRef", "object %s not found !", name) ; 
934   }
935   return phosFolder->GetListOfFolders()->GetObjectRef(erp) ;
936   
937
938
939 //____________________________________________________________________________ 
940 TObject** AliPHOSGetter::CpvRecPointsRef(const char * name) const 
941 { // -------------- RecPoints -------------------------------------------
942   
943   // the hierarchy is //Folders/Run/Event/RecData/PHOS/CPVRecPoints/name
944    
945   if ( !fRecoFolder ) {
946     Fatal("CpvRecPointsRef", "Folder //%s not found !", fRecoFolder->GetName() ) ;
947   }    
948
949   TFolder * phosFolder  = dynamic_cast<TFolder*>(fRecoFolder->FindObject("PHOS/CPVRecPoints")) ; 
950   if ( !phosFolder ) {
951     Fatal("CpvRecPointsRef", "Folder //%s/PHOS/CPVRecPoints/ not found !", fRecoFolder->GetName() ) ;
952   }    
953
954   TObject * crp = phosFolder->FindObject(name ) ;
955   if ( !crp )   {
956     Fatal("CpvRecPointsRef", "object %s nott found", name ) ; 
957   }
958   return phosFolder->GetListOfFolders()->GetObjectRef(crp) ;
959   
960
961
962 //____________________________________________________________________________ 
963 const Bool_t AliPHOSGetter::PostClusterizer(AliPHOSClusterizer * clu) const 
964 { // ------------------ AliPHOSClusterizer ------------------------
965   
966   // the hierarchy is //Folders/Tasks/Reconstructioner/PHOS/sdigitsname
967
968   TTask * tasks  = dynamic_cast<TTask*>(fTasksFolder->FindObject("Reconstructioner")) ; 
969
970   if ( !tasks ) {
971     Error("PostClusterizer", "Task //%s/Reconstructioner not found !", fTasksFolder) ;
972     return kFALSE ;
973   }        
974         
975   TTask * phos = dynamic_cast<TTask*>(tasks->GetListOfTasks()->FindObject("PHOS")) ; 
976   if ( !phos )  {
977     if (fDebug) {
978       Warning("PostClusterizer", "//%s/Reconstructioner/PHOS not found!", fTasksFolder) ; 
979       Info("PostClusterizer", "Adding //%s/Reconstructioner/PHOS", fTasksFolder) ; 
980     }
981     phos = new TTask("PHOS", "") ; 
982     tasks->Add(phos) ; 
983   } 
984
985   AliPHOSClusterizer * phoscl = dynamic_cast<AliPHOSClusterizer*>(phos->GetListOfTasks()->FindObject(clu->GetName())) ; 
986   if (phoscl) { 
987     if (fDebug)
988       Info("PostClusterizer", "Task %s already exists", clu->GetName()) ; 
989     phoscl->Delete() ; 
990     phos->GetListOfTasks()->Remove(phoscl) ;
991   }
992   phos->Add(clu) ;      
993   return kTRUE; 
994
995
996 //____________________________________________________________________________ 
997 TObject** AliPHOSGetter::ClusterizerRef(const char * name) const 
998 { // ------------------ AliPHOSClusterizer ------------------------
999   
1000   TTask * tasks  = dynamic_cast<TTask*>(fTasksFolder->FindObject("Reconstructioner")) ; 
1001
1002   if ( !tasks ) {
1003     Fatal("ClusterizerRef", "Task //%s/Reconstructioner not found !", fTasksFolder->GetName() ) ;
1004   }        
1005         
1006   TTask * phos = dynamic_cast<TTask*>(tasks->GetListOfTasks()->FindObject("PHOS")) ; 
1007   if ( !phos )  {
1008     Fatal("ClusterizerRef", " //%s/Reconstructioner/PHOS not founf !", fTasksFolder->GetName() ) ; 
1009   }   
1010
1011   TList * l = phos->GetListOfTasks() ; 
1012   TIter it(l) ;
1013   TTask * task ;
1014   TTask * clu = 0 ;
1015   TString cluname(name) ;
1016   cluname+=":clu" ;
1017   while((task = static_cast<TTask *>(it.Next()) )){
1018     TString taskname(task->GetName()) ;
1019     if(taskname.BeginsWith(cluname)){
1020       clu = task ;
1021       break ;
1022     }
1023   }
1024
1025   if(!clu) {
1026     Fatal("ClusterizerRef", "Task //%s/Reconstructioner/clusterizer/%s not found", fTasksFolder->GetName(), name) ;
1027   }
1028
1029   return l->GetObjectRef(clu) ;
1030
1031 }
1032
1033 //____________________________________________________________________________ 
1034 const Bool_t AliPHOSGetter::PostClusterizer(const char * name) const 
1035 { // ------------------ AliPHOSClusterizer ------------------------
1036
1037   // the hierarchy is //Folders/Tasks/Reconstructioner/PHOS/sdigitsname
1038   
1039   TTask * tasks  = dynamic_cast<TTask*>(fTasksFolder->FindObject("Reconstructioner")) ; 
1040
1041   if ( !tasks ) {
1042     Error("PostClusterizer", "Task//%s/Reconstructioner not found !", fTasksFolder) ; 
1043     return kFALSE ;
1044   }        
1045   
1046   TTask * phos = dynamic_cast<TTask*>(tasks->GetListOfTasks()->FindObject("PHOS")) ; 
1047   if ( !phos )  {
1048     if (fDebug) {
1049       Warning("PostClusterizer", "//%s/Reconstructioner/PHOS not found!", fTasksFolder) ;
1050       Info("PostClusterizer", "Adding //%s/Reconstructioner/PHOS", fTasksFolder) ;
1051     }
1052     phos = new TTask("PHOS", "") ; 
1053     tasks->Add(phos) ; 
1054   } 
1055
1056   TList * l = phos->GetListOfTasks() ;   
1057   TIter it(l) ;
1058   TString clun(name) ;
1059   clun+=":clu" ; 
1060   TTask * task ;
1061   while((task = static_cast<TTask *>(it.Next()) )){
1062     TString taskname(task->GetName()) ;
1063     if(taskname.BeginsWith(clun))
1064       return kTRUE ;
1065   }
1066
1067   AliPHOSClusterizerv1 * phoscl = new AliPHOSClusterizerv1() ;
1068   clun+="-v1" ; 
1069   phoscl->SetName(clun) ;
1070   phoscl->SetTitle(fHeaderFile) ;
1071   phos->Add(phoscl) ;
1072   return kTRUE; 
1073   
1074 }
1075
1076 //____________________________________________________________________________ 
1077 const Bool_t AliPHOSGetter::PostTrackSegments(const char * name) const 
1078 { // ---------------TrackSegments -----------------------------------
1079   
1080   // the hierarchy is //Folders/Run/Event/RecData/PHOS/TrackSegments/name
1081
1082   TFolder * phosFolder  = dynamic_cast<TFolder*>(fRecoFolder->FindObject("PHOS")) ; 
1083   
1084   if ( !phosFolder ) {
1085     if (fDebug) {
1086       Warning("PostTrackSegments", "-> Folder //%s/PHOS/ not found", fRecoFolder->GetName()) ;
1087       Info("PostTrackSegments", "-> Adding Folder //%s/PHOS", fRecoFolder->GetName()) ;
1088     }
1089     phosFolder = fRecoFolder->AddFolder("PHOS", "Reconstructed data from PHOS") ;  
1090   }    
1091
1092   TFolder * phosTSFolder  = dynamic_cast<TFolder*>(phosFolder->FindObject("TrackSegments")) ;
1093   if ( !phosTSFolder ) {
1094     if (fDebug) {
1095       Warning("PostTrackSegments", "-> Folder //%s/PHOS/TrackSegments/ not found!", fRecoFolder->GetName() ) ; 
1096       Info("PostTrackSegments", "-> Adding Folder //%s/PHOS/TrackSegments/", fRecoFolder->GetName()) ; 
1097     }
1098     phosTSFolder = phosFolder->AddFolder("TrackSegments", "TrackSegments from PHOS") ;  
1099   }    
1100   
1101   TObject * tss =  phosTSFolder->FindObject( name ) ;
1102   if (!tss) {
1103     TClonesArray * ts = new TClonesArray("AliPHOSTrackSegment",100) ;
1104     ts->SetName(name) ;
1105     phosTSFolder->Add(ts) ;  
1106   }
1107   return kTRUE; 
1108
1109
1110 //____________________________________________________________________________ 
1111 TObject** AliPHOSGetter::TrackSegmentsRef(const char * name) const 
1112 { // ---------------TrackSegments -----------------------------------
1113   
1114   // the hierarchy is //Folders/Run/Event/RecData/PHOS/TrackSegments/name
1115
1116  if ( !fRecoFolder ) {
1117     Fatal("TrackSegmentsRef", "Folder //%s not found !", fRecoFolder->GetName() ) ;
1118   }    
1119
1120   TFolder * phosFolder  = dynamic_cast<TFolder*>(fRecoFolder->FindObject("PHOS/TrackSegments")) ; 
1121   if ( !phosFolder ) {
1122     Fatal("TrackSegmentsRef", "Folder //%s/PHOS/TrackSegments/ not found !", fRecoFolder->GetName() ) ;
1123   }    
1124   
1125   TObject * tss =  phosFolder->FindObject(name) ;
1126   if (!tss) {
1127     Fatal("TrackSegmentsRef", "object %s not found !", name) ;  
1128   }
1129   return phosFolder->GetListOfFolders()->GetObjectRef(tss) ;
1130
1131
1132 //____________________________________________________________________________ 
1133 const Bool_t AliPHOSGetter::PostTrackSegmentMaker(AliPHOSTrackSegmentMaker * tsmaker) const 
1134 { //------------Track Segment Maker ------------------------------
1135   
1136   // the hierarchy is //Folders/Tasks/Reconstructioner/PHOS/sdigitsname
1137
1138   TTask * tasks  = dynamic_cast<TTask*>(fTasksFolder->FindObject("Reconstructioner")) ; 
1139
1140   if ( !tasks ) {
1141     Error("PostTrackSegmentMaker", "Task //%s/Reconstructioner not found !", fTasksFolder) ;
1142     return kFALSE ;
1143   }        
1144         
1145   TTask * phos = dynamic_cast<TTask*>(tasks->GetListOfTasks()->FindObject("PHOS")) ; 
1146   if ( !phos )  {
1147     if (fDebug) {
1148       Warning("PostTrackSegmentMaker", "//%s/Reconstructioner/PHOS not found!", fTasksFolder) ; 
1149       Info("PostTrackSegmentMaker", "Adding //%s/Reconstructioner/PHOS", fTasksFolder) ;
1150     }
1151     phos = new TTask("PHOS", "") ; 
1152     tasks->Add(phos) ; 
1153   } 
1154
1155   AliPHOSTrackSegmentMaker * phosts = 
1156     dynamic_cast<AliPHOSTrackSegmentMaker*>(phos->GetListOfTasks()->FindObject(tsmaker->GetName())) ; 
1157   if (phosts) { 
1158     phosts->Delete() ;
1159     phos->GetListOfTasks()->Remove(phosts) ;
1160   }
1161   phos->Add(tsmaker) ;      
1162   return kTRUE; 
1163   
1164
1165 //____________________________________________________________________________ 
1166 const Bool_t AliPHOSGetter::PostTrackSegmentMaker(const char * name) const 
1167 { //------------Track Segment Maker ------------------------------
1168   
1169   // the hierarchy is //Folders/Tasks/Reconstructioner/PHOS/sdigitsname
1170   
1171   TTask * tasks  = dynamic_cast<TTask*>(fTasksFolder->FindObject("Reconstructioner")) ; 
1172   
1173   if ( !tasks ) {
1174     Error("PostTrackSegmentMaker", "Task //%s/Reconstructioner not found !", fTasksFolder->GetName() ) ;
1175     return kFALSE ;
1176   }        
1177   
1178   TTask * phos = dynamic_cast<TTask*>(tasks->GetListOfTasks()->FindObject("PHOS")) ; 
1179   if ( !phos )  {
1180     if (fDebug) {
1181       Warning("PostTrackSegmentMaker", "//%s/Reconstructioner/PHOS not found!", fTasksFolder->GetName() ) ; 
1182       Info("PostTrackSegmentMaker", "Adding //%s/Reconstructioner/PHOS", fTasksFolder->GetName()) ;
1183     }
1184     phos = new TTask("PHOS", "") ; 
1185     tasks->Add(phos) ; 
1186   } 
1187
1188   TList * l = phos->GetListOfTasks() ;   
1189   TIter it(l) ;
1190   TString tsn(name);
1191   tsn+=":tsm" ; 
1192   TTask * task ;
1193   while((task = static_cast<TTask *>(it.Next()) )){
1194     TString taskname(task->GetName()) ;
1195     if(taskname.BeginsWith(tsn))
1196       return kTRUE ;
1197   }
1198   
1199   AliPHOSTrackSegmentMakerv1 * phosts = new AliPHOSTrackSegmentMakerv1() ;
1200   tsn+="-v1" ;
1201   phosts->SetName(tsn) ;
1202   phosts->SetTitle(fHeaderFile) ;
1203   phos->Add(phosts) ;      
1204   return kTRUE; 
1205   
1206
1207
1208 //____________________________________________________________________________ 
1209 TObject** AliPHOSGetter::TSMakerRef(const char * name) const 
1210 { //------------Track Segment Maker ------------------------------
1211   
1212   TTask * tasks  = dynamic_cast<TTask*>(fTasksFolder->FindObject("Reconstructioner")) ; 
1213
1214   if ( !tasks ) {
1215     Fatal("TSMakerRef", "Task //%s/Reconstructioner not found !", fTasksFolder->GetName() ) ;
1216   }        
1217         
1218   TTask * phos = dynamic_cast<TTask*>(tasks->GetListOfTasks()->FindObject("PHOS")) ; 
1219   if ( !phos )  {
1220     Fatal("TSMakerRef", "//%s/Reconstructioner/PHOS not found !", fTasksFolder->GetName() ) ; 
1221   }   
1222
1223   TList * l = phos->GetListOfTasks() ; 
1224   TIter it(l) ;
1225   TTask * task ;
1226   TTask * tsm = 0 ;
1227   TString tsmname(name) ;
1228   tsmname+=":tsm" ;
1229   while((task = static_cast<TTask *>(it.Next()) )){
1230     TString taskname(task->GetName()) ;
1231     if(taskname.BeginsWith(tsmname)){
1232       tsm = task ;
1233       break ;
1234     }
1235   }
1236   
1237   if(!tsm) {
1238    Fatal("TSMakerRef", "Task //%s/Reconstructioner/PHOS/TrackSegmentMarker/%s not found !", fTasksFolder->GetName(),  name) ;
1239   }
1240  
1241   return l->GetObjectRef(tsm) ;
1242
1243
1244
1245 //____________________________________________________________________________ 
1246 const Bool_t AliPHOSGetter::PostRecParticles(const char * name) const 
1247 {  // -------------------- RecParticles ------------------------
1248   
1249   // the hierarchy is //Folders/Run/Event/RecData/PHOS/RecParticles/name
1250
1251   TFolder * phosFolder  = dynamic_cast<TFolder*>(fRecoFolder->FindObject("PHOS")) ; 
1252   
1253   if ( !phosFolder ) {
1254     if (fDebug) {
1255       Warning("PostRecParticles", "-> Folder //%s/PHOS/ not found!", fRecoFolder->GetName()) ;
1256       Info("PostRecParticles", "-> Adding Folder //%s/PHOS/", fRecoFolder->GetName()) ;
1257     }
1258     phosFolder = fRecoFolder->AddFolder("PHOS", "Reconstructed data from PHOS") ;  
1259   }    
1260
1261  TFolder * phosRPaFolder  = dynamic_cast<TFolder*>(phosFolder->FindObject("RecParticles")) ;
1262   if ( !phosRPaFolder ) {
1263     if (fDebug) {
1264       Warning("PostRecParticles", "-> Folder //%s/PHOS/RecParticles/ not found!", fRecoFolder->GetName()) ;
1265       Info("PostRecParticles", "-> Adding Folder //%s/PHOS/RecParticles/", fRecoFolder->GetName()) ;
1266     }
1267     phosRPaFolder = phosFolder->AddFolder("RecParticles", "RecParticles from PHOS") ;  
1268   } 
1269
1270   TObject * rps = phosRPaFolder->FindObject( name )  ;
1271   if ( !rps ) {
1272     TClonesArray * rp = new TClonesArray("AliPHOSRecParticle",100) ;
1273     rp->SetName(name) ;    
1274     phosRPaFolder->Add(rp) ;  
1275   }
1276   return kTRUE; 
1277
1278
1279 //____________________________________________________________________________ 
1280 TObject** AliPHOSGetter::RecParticlesRef(const char * name) const 
1281 { // ---------------RecParticles -----------------------------------
1282   
1283   // the hierarchy is //Folders/Run/Event/RecData/PHOS/TrackSegments/name
1284
1285  if ( !fRecoFolder ) {
1286     Fatal("RecParticlesRef", "Folder//%s not found !", fRecoFolder->GetName() ) ; 
1287   }    
1288
1289   TFolder * phosFolder  = dynamic_cast<TFolder*>(fRecoFolder->FindObject("PHOS/RecParticles")) ; 
1290   if ( !phosFolder ) {
1291     Fatal("RecParticlesRef", "Folder //%s/PHOS/RecParticles/ not found !", fRecoFolder->GetName() ) ;
1292   }    
1293
1294   TObject * tss =  phosFolder->FindObject(name  ) ;
1295   if (!tss) {
1296     Fatal("RecParticlesRef", "object %s not found !", name) ; 
1297   }
1298   return phosFolder->GetListOfFolders()->GetObjectRef(tss) ;
1299 }
1300 //____________________________________________________________________________ 
1301 const UShort_t AliPHOSGetter::EventPattern(void){
1302   if(fBTE)
1303     return fBTE->GetPattern() ;
1304   else
1305     return 0 ;
1306 }
1307 //____________________________________________________________________________ 
1308 const Float_t AliPHOSGetter::BeamEnergy(void){
1309   if(fBTE)
1310     return fBTE->GetBeamEnergy() ;
1311   else
1312     return 0 ;
1313 }
1314 //____________________________________________________________________________ 
1315 const Bool_t AliPHOSGetter::PostPID(AliPHOSPID * pid) const 
1316 {      // ------------AliPHOS PID -----------------------------
1317
1318   TTask * tasks  = dynamic_cast<TTask*>(fTasksFolder->FindObject("Reconstructioner")) ; 
1319
1320   if ( !tasks ) {
1321     Error("PostPID", "Task //%s/Reconstructioner not found !", fTasksFolder) ;
1322     return kFALSE ;
1323   }        
1324   
1325   TTask * phos = dynamic_cast<TTask*>(tasks->GetListOfTasks()->FindObject("PHOS")) ; 
1326   if ( !phos )  {
1327     if (fDebug) {
1328       Warning("PostPID", "//%s/Reconstructioner/PHOS not found!", fTasksFolder) ; 
1329       Info("PostPID", "Adding //%s/Reconstructioner/PHOS", fTasksFolder) ;
1330     }
1331     phos = new TTask("PHOS", "") ; 
1332     tasks->Add(phos) ; 
1333   } 
1334
1335   AliPHOSPID * phospid = dynamic_cast<AliPHOSPID*>(phos->GetListOfTasks()->FindObject(pid->GetName())) ; 
1336   if (phospid) { 
1337     if (fDebug)
1338       Info("PostPID", "-> Task %s qlready exists", pid->GetName()) ; 
1339     phos->GetListOfTasks()->Remove(phospid) ;
1340   }
1341   
1342   phos->Add(pid) ;      
1343   return kTRUE; 
1344
1345
1346 //____________________________________________________________________________ 
1347 const Bool_t AliPHOSGetter::PostPID(const char * name) const 
1348 {     
1349   // the hierarchy is //Folders/Tasks/Reconstructioner/PHOS/sdigitsname
1350   
1351   TTask * tasks  = dynamic_cast<TTask*>(fTasksFolder->FindObject("Reconstructioner")) ; 
1352
1353   if ( !tasks ) {
1354     Error("PostPID", "Task //%s/Reconstructioner not found !", fTasksFolder->GetName() ) ;
1355     return kFALSE ;
1356   }        
1357   
1358   TTask * phos = dynamic_cast<TTask*>(tasks->GetListOfTasks()->FindObject("PHOS")) ; 
1359   if ( !phos )  {
1360     if (fDebug) {
1361       Warning("PostPID", "//%s/Reconstructioner/PHOS not found!", fTasksFolder->GetName()) ; 
1362       Info("PostPID", "Adding //%s/Reconstructioner/PHOS", fTasksFolder->GetName()) ;
1363     }
1364     phos = new TTask("PHOS", "") ; 
1365     tasks->Add(phos) ; 
1366   } 
1367
1368   TList * l = phos->GetListOfTasks() ;   
1369   TIter it(l) ;
1370   TString pidname(name) ;
1371   pidname+=":pid" ;
1372   TTask * task ;
1373   while((task = static_cast<TTask *>(it.Next()) )){
1374     TString taskname(task->GetName()) ;
1375     if(taskname.BeginsWith(pidname))
1376       return kTRUE ;
1377   }
1378  
1379   AliPHOSPIDv1 * phospid = new AliPHOSPIDv1() ;
1380   pidname+="-v1" ;
1381   phospid->SetName(pidname) ; 
1382   phospid->SetTitle(fHeaderFile) ;
1383   phos->Add(phospid) ;      
1384   
1385   return kTRUE; 
1386
1387
1388 //____________________________________________________________________________ 
1389 TObject** AliPHOSGetter::PIDRef(const char * name) const 
1390 { //------------PID ------------------------------
1391
1392   TTask * tasks  = dynamic_cast<TTask*>(fTasksFolder->FindObject("Reconstructioner")) ; 
1393
1394   if ( !tasks ) {
1395     Fatal("PIDRef", "Task //%s/Reconstructioner not found !", fTasksFolder->GetName() ) ;
1396   }        
1397         
1398   TTask * phos = dynamic_cast<TTask*>(tasks->GetListOfTasks()->FindObject("PHOS")) ; 
1399   if ( !phos )  {
1400     Fatal("PIDRef", "//%s/Reconstructioner/PHOS not found !", fTasksFolder->GetName() ) ; 
1401   }   
1402   
1403   TList * l = phos->GetListOfTasks() ; 
1404   TIter it(l) ;
1405   TTask * task ;
1406   TTask * pid = 0 ;
1407   TString pidname(name) ;
1408   pidname+=":pid" ;
1409   while((task = static_cast<TTask *>(it.Next()) )){
1410     TString taskname(task->GetName()) ;
1411     if(taskname.BeginsWith(pidname)){
1412       pid = task ;
1413       break ;
1414     }
1415   }
1416   
1417   if(!pid) {
1418     Fatal("PIDRef", "Task //%s/Reconstructioner/PHOS/PID/%s not found !", fTasksFolder->GetName(), name) ;
1419   }
1420   
1421     return l->GetObjectRef(pid) ;
1422
1423
1424 //____________________________________________________________________________ 
1425 const Bool_t AliPHOSGetter::PostQA(void) const 
1426 { // ------------------ QA ---------------------------------
1427
1428   // the hierarchy is //Folders/Run/Conditions/QA/PHOS/alarmsName
1429
1430   TFolder * phosFolder = dynamic_cast<TFolder*>(fQAFolder->FindObject("PHOS")) ; 
1431   if ( !phosFolder ) {
1432     if (fDebug) {
1433       Warning("PostQA", "-> Folder //%s/PHOS/ not found!", fQAFolder) ;
1434       Info("PostQA", "-> Adding Folder //%s/PHOS", fQAFolder) ;
1435     }
1436     phosFolder = fQAFolder->AddFolder("PHOS", "QA from PHOS") ; 
1437   }      
1438
1439   return kTRUE;
1440 }
1441
1442 //____________________________________________________________________________ 
1443 TObject** AliPHOSGetter::AlarmsRef(void) const 
1444 {  //------- Alarms ----------------------
1445
1446   
1447   // the hierarchy is //Folders/Run/Conditions/QA/PHOS
1448   if ( !fQAFolder ) {
1449     Fatal("AlarmsRef", "Folder //%s not found !", fQAFolder) ;
1450   }    
1451  
1452   TFolder * phosFolder = dynamic_cast<TFolder *>(fQAFolder->FindObject("PHOS")) ;
1453   if ( !phosFolder ) {
1454     Fatal("AlarmsRef", "Folder //%s/PHOS/ not found !", fQAFolder) ;
1455   }
1456    
1457   return fQAFolder->GetListOfFolders()->GetObjectRef(phosFolder) ;
1458 }
1459
1460
1461 //____________________________________________________________________________ 
1462 TTree * AliPHOSGetter::TreeK(TString filename)  
1463 {
1464
1465   // returns TreeK from file filename
1466   // usefull in case of split file
1467
1468   if ( filename.IsNull() ) 
1469     filename = fHeaderFile ; 
1470
1471   TFile * file = 0 ; 
1472   file = static_cast<TFile*>(gROOT->GetFile(filename.Data() ) ) ;
1473 //   if (file && (filename != fHeaderFile) ) {  // file already open 
1474 //     file->Close() ; 
1475 //     //delete fAlice ; 
1476 //   }
1477   if(!file || !file->IsOpen())    
1478     file = TFile::Open(filename.Data(), "read") ;
1479   if(filename != fHeaderFile ){
1480     fAlice = dynamic_cast<AliRun *>(file->Get("gAlice")) ;
1481   } 
1482   TString treeName("TreeK") ; 
1483   treeName += EventNumber()  ; 
1484   TTree * tree = dynamic_cast<TTree *>(file->Get(treeName.Data())) ;
1485   if (!tree && fDebug)  
1486     Warning("TreeK", "-> %s not found in %s", treeName.Data(), filename.Data()) ; 
1487   
1488   return tree ;                       
1489 }
1490
1491 //____________________________________________________________________________ 
1492 TTree * AliPHOSGetter::TreeH(TString filename)  
1493 {
1494
1495   // returns TreeH from file filename
1496   // usefull in case of split file
1497
1498   if ( filename.IsNull() ) 
1499     filename = fHeaderFile ; 
1500
1501   TFile * file = 0 ; 
1502   file = static_cast<TFile*>(gROOT->GetFile(filename.Data() ) ) ;
1503   if (!file) { // file not open yet
1504     file = TFile::Open(filename.Data(), "read") ; 
1505   }
1506   TString treeName("TreeH") ; 
1507   treeName += EventNumber()  ; 
1508   TTree * tree = static_cast<TTree *>(file->Get(treeName.Data())) ;
1509   if (!tree && fDebug)  
1510     Warning("TreeH", "-> %s not found in %s", treeName.Data(), filename.Data()) ; 
1511   
1512   return tree ;                       
1513 }
1514
1515 //____________________________________________________________________________ 
1516 TTree * AliPHOSGetter::TreeS(TString filename)  
1517 {
1518
1519   // returns TreeS from file filename
1520   // usefull in case of split file
1521
1522   if ( filename.IsNull() ) 
1523     filename = fHeaderFile ; 
1524
1525   TFile * file = 0 ; 
1526   file = static_cast<TFile*>(gROOT->GetFile(filename.Data() ) ) ;
1527   if (!file) { // file not open yet
1528     file = TFile::Open(filename.Data(), "read") ; 
1529   }
1530   TString treeName("TreeS") ; 
1531   treeName += EventNumber()  ; 
1532   TTree * tree = static_cast<TTree *>(file->Get(treeName.Data())) ;
1533   if (!tree && fDebug)  
1534     Warning("TreeS", "-> %s not found in %s", treeName.Data(), filename.Data() ); 
1535   
1536   return tree ;                       
1537 }
1538
1539 //____________________________________________________________________________ 
1540 TTree * AliPHOSGetter::TreeD(TString filename)  
1541 {
1542
1543   // returns TreeD from file filename
1544   // usefull in case of split file
1545
1546   if ( filename.IsNull() ) 
1547     filename = fHeaderFile ; 
1548
1549   TFile * file = 0 ; 
1550   file = static_cast<TFile*>(gROOT->GetFile(filename.Data() ) ) ;
1551   if (!file) { // file not open yet
1552     file = TFile::Open(filename.Data(), "read") ; 
1553   }
1554   TString treeName("TreeD") ; 
1555   treeName += EventNumber()  ; 
1556   TTree * tree = static_cast<TTree *>(file->Get(treeName.Data())) ;
1557   if (!tree && fDebug)  
1558     Warning("TreeD", "-> %s not found in %s", treeName.Data(), filename.Data()) ; 
1559   
1560   return tree ;                       
1561 }
1562
1563 //____________________________________________________________________________ 
1564 const TParticle * AliPHOSGetter::Primary(Int_t index) const 
1565 {
1566   // Return primary particle numbered by <index>
1567
1568   if(index < 0) 
1569     return 0 ;
1570   TParticle *  p = 0 ;
1571   if (fAlice) 
1572     p = fAlice->Particle(index) ; 
1573   else 
1574     p = gAlice->Particle(index) ; 
1575   
1576   return p ; 
1577     
1578 }
1579
1580 //____________________________________________________________________________ 
1581 const TParticle * AliPHOSGetter::Secondary(TParticle* p, Int_t index) const
1582 {
1583   // Return first (index=1) or second (index=2) secondary particle of primary particle p 
1584
1585   if(index <= 0) 
1586     return 0 ;
1587   if(index > 2)
1588     return 0 ;
1589
1590   if(p) {
1591   Int_t daughterIndex = p->GetDaughter(index-1) ; 
1592   return  gAlice->Particle(daughterIndex) ; 
1593   }
1594   else
1595     return 0 ;
1596 }
1597
1598 //____________________________________________________________________________ 
1599 Int_t AliPHOSGetter::ReadTreeD(const Int_t event)
1600 {
1601   // Read the digit tree gAlice->TreeD()  
1602   
1603   TTree * treeD ;
1604   if(fToSplit){
1605     TFile * file = static_cast<TFile*>(gROOT->GetFile(fDigitsFileName)); 
1606     if(!file) 
1607       file = TFile::Open(fDigitsFileName) ;      
1608     // Get Digits Tree header from file
1609     TString treeName("TreeD") ;
1610     treeName += event ; 
1611     treeD = dynamic_cast<TTree*>(file->Get(treeName.Data()));
1612     if(!treeD){ // TreeD not found in header file
1613       if (fDebug)
1614         Warning("ReadTreeD", "-> Cannot find TreeD in %s", fDigitsFileName.Data()) ;
1615       return 1;
1616     }
1617   }
1618   else
1619     treeD = gAlice->TreeD() ;
1620   
1621   TObjArray * lob = static_cast<TObjArray*>(treeD->GetListOfBranches()) ;
1622   TIter next(lob) ; 
1623   TBranch * branch = 0 ; 
1624   TBranch * digitsbranch = 0 ; 
1625   TBranch * digitizerbranch = 0 ; 
1626   Bool_t phosfound = kFALSE, digitizerfound = kFALSE ; 
1627   
1628   while ( (branch = static_cast<TBranch*>(next())) && (!phosfound || !digitizerfound) ) {
1629     if ( (strcmp(branch->GetName(), "PHOS")==0) && (strcmp(branch->GetTitle(), fDigitsTitle)==0) ) {
1630       digitsbranch = branch ; 
1631       phosfound = kTRUE ;
1632     }
1633     else if ( ((strcmp(branch->GetName(), "AliPHOSDigitizer")==0)||
1634                (strcmp(branch->GetName(), "AliPHOSRaw2Digits")==0)) &&
1635               (strcmp(branch->GetTitle(), fDigitsTitle)==0) ) {
1636       digitizerbranch = branch ; 
1637       digitizerfound = kTRUE ; 
1638     }
1639   }
1640   
1641   if ( !phosfound || !digitizerfound ) {
1642     if (fDebug)
1643       Warning("ReadTreeD", "-> Cannot find Digits and/or Digitizer with name %s", fDigitsTitle.Data()) ;
1644     return 2; 
1645   }   
1646   
1647   //read digits
1648   if(!Digits(fDigitsTitle) ) 
1649     PostDigits(fDigitsTitle);
1650   digitsbranch->SetAddress(DigitsRef(fDigitsTitle)) ;
1651   digitsbranch->GetEntry(0) ;
1652   
1653   // read  the Digitizer
1654   if(Digitizer()){
1655     if(strcmp(Digitizer()->IsA()->GetName(),digitizerbranch->GetName())!=0){
1656       RemoveTask("D", fDigitsTitle) ;
1657       if(strcmp(digitizerbranch->GetName(), "AliPHOSDigitizer")==0)
1658         PostDigitizer("Digitizer") ;
1659       else
1660         PostDigitizer("Raw2Digits") ;
1661     }
1662   }
1663   else{
1664     if(strcmp(digitizerbranch->GetName(), "AliPHOSDigitizer")==0)
1665       PostDigitizer("Digitizer") ;
1666     else
1667       PostDigitizer("Raw2Digits") ;
1668   }
1669     
1670
1671   digitizerbranch->SetAddress(DigitizerRef(fDigitsTitle)) ;
1672   digitizerbranch->GetEntry(0) ;
1673   
1674
1675 //   if((!fcdb)&&(strcmp(digitizerbranch->GetName(), "AliPHOSRaw2Digits")==0))
1676 //     ReadCalibrationDB("Primordial","beamtest.root") ;
1677
1678
1679   if(gAlice->TreeD()!=treeD)
1680     treeD->Delete();
1681
1682   return 0 ; 
1683 }
1684 //____________________________________________________________________________ 
1685 //void AliPHOSGetter::ReadCalibrationDB(const char * database,const char * filename){
1686 //
1687 //  if(fcdb && (strcmp(database,fcdb->GetTitle())==0))
1688 //    return ;
1689 //
1690 //  TFile * file = gROOT->GetFile(filename) ;
1691 //  if(!file)
1692 //    file = TFile::Open(filename);
1693 //  if(!file){
1694 //    Error ("ReadCalibrationDB", "Cannot open file %s", filename) ;
1695 //    return ;
1696 //  }
1697 //  if(fcdb)
1698 //    fcdb->Delete() ;
1699 //  fcdb = dynamic_cast<AliPHOSCalibrationDB *>(file->Get("AliPHOSCalibrationDB")) ;
1700 //  if(!fcdb)
1701 //    Error ("ReadCalibrationDB", "No database %s in file %s", database, filename) ;
1702 //}
1703
1704 //____________________________________________________________________________ 
1705 Int_t AliPHOSGetter::ReadTreeH()
1706 {
1707   // Read the first entry of PHOS branch in hit tree gAlice->TreeH()
1708   
1709   TTree * treeH = gAlice->TreeH() ;
1710
1711   if(!treeH) {// TreeH not found in header file
1712  
1713     if (fDebug) 
1714       Warning("ReadTreeH", "-> Cannot find TreeH in %s", fHeaderFile.Data() ) ;
1715     
1716     TString searchFileName("PHOS.Hits") ; 
1717     if((strcmp(fBranchTitle.Data(),"Default")!=0)&&(strcmp(fBranchTitle.Data(),"")!=0)){
1718       searchFileName+="." ;
1719       searchFileName += fBranchTitle ;
1720     }
1721     searchFileName+=".root" ;
1722     
1723     if ( (treeH = TreeH(searchFileName)) ) { //found TreeH in the file which contains the hits
1724       if (fDebug) 
1725         Info("ReadTreeH", "-> TreeH found in %s", searchFileName.Data()) ; 
1726       
1727     } else {
1728       Error("ReadTreeH", "TreeH not found") ; 
1729       return 1;
1730     }  
1731   }
1732   
1733   TBranch * hitsbranch = static_cast<TBranch*>(treeH->GetBranch("PHOS")) ;
1734   if ( !hitsbranch ) {
1735     if (fDebug)
1736       Warning("ReadTreeH", "-> Cannot find branch PHOS") ; 
1737     return 2;
1738   }
1739   if(!Hits())
1740     PostHits() ;
1741
1742   if (hitsbranch->GetEntries() > 1 ) {
1743     (dynamic_cast<TClonesArray*> (*HitsRef()))->Clear() ;
1744     TClonesArray * tempo =  new TClonesArray("AliPHOSHit",1000) ;
1745     TClonesArray * hits = dynamic_cast<TClonesArray*>(*HitsRef()) ; 
1746     hitsbranch->SetAddress(&tempo) ;
1747     Int_t index = 0 ; 
1748     Int_t i = 0 ;
1749     for (i = 0 ; i < hitsbranch->GetEntries() ; i++) {
1750       hitsbranch->GetEntry(i) ;
1751       Int_t j = 0 ; 
1752       for ( j = 0 ; j < tempo->GetEntries() ; j++) { 
1753         const AliPHOSHit * hit = static_cast<const AliPHOSHit *>(tempo->At(j)) ; 
1754         new((*hits)[index]) AliPHOSHit( *hit ) ;
1755         index++ ; 
1756       }
1757     }
1758     delete tempo ; 
1759   }
1760   else {
1761     (dynamic_cast<TClonesArray*> (*HitsRef()))->Clear() ;
1762     hitsbranch->SetAddress(HitsRef()) ;
1763     hitsbranch->GetEntry(0) ;
1764   }
1765   return 0 ; 
1766 }
1767
1768 //____________________________________________________________________________ 
1769 void AliPHOSGetter::Track(const Int_t itrack) 
1770 {
1771   // Read the first entry of PHOS branch in hit tree gAlice->TreeH()
1772
1773   if(gAlice->TreeH()== 0){
1774     Error("Track", "Cannot read TreeH") ;
1775     return ;
1776   }
1777   
1778   TBranch * hitsbranch = dynamic_cast<TBranch*>(gAlice->TreeH()->GetListOfBranches()->FindObject("PHOS")) ;
1779   if ( !hitsbranch ) {
1780     if (fDebug)
1781       Warning("Track", "Cannot find branch PHOS") ; 
1782     return ;
1783   }  
1784   if(!Hits())
1785     PostHits() ;
1786
1787   (dynamic_cast<TClonesArray*> (*HitsRef()))->Clear() ;
1788   hitsbranch->SetAddress(HitsRef()) ;
1789   hitsbranch->GetEntry(itrack) ;
1790
1791 }
1792
1793 //____________________________________________________________________________ 
1794 void AliPHOSGetter::ReadTreeQA()
1795 {
1796   // Read the digit tree gAlice->TreeQA()
1797   // so far only PHOS knows about this Tree  
1798
1799   if(PHOS()->TreeQA()== 0){
1800     Error("ReadTreeQA", "Cannot read TreeQA") ;
1801     return ;
1802   }
1803   
1804   TBranch * qabranch = PHOS()->TreeQA()->GetBranch("PHOS") ; 
1805   if (!qabranch) { 
1806     if (fDebug)
1807       Warning("ReadTreeQA", "Cannot find QA Alarms for PHOS");
1808     return ; 
1809   }   
1810   
1811   if(!Alarms())
1812     PostQA() ; 
1813
1814   qabranch->SetAddress(AlarmsRef()) ;
1815
1816   qabranch->GetEntry(0) ;
1817  
1818 //   PostQA("PHOS") ; 
1819 //   TFolder * alarmsF = Alarms() ; 
1820 //   alarmsF->Clear() ; 
1821 //   qabranch->SetAddress(&alarmsF) ;
1822 //   qabranch->GetEntry(0) ;
1823   
1824 }
1825
1826 //____________________________________________________________________________ 
1827 Int_t AliPHOSGetter::ReadTreeR(const Int_t event)
1828 {
1829   // Read the reconstrunction tree gAlice->TreeR()
1830   // A particularity has been introduced here :
1831   //  if gime->Event(ievent,"R") is called branches with the current title are read, the current title
1832   //   being for example give in AliPHOSPID(fileName, title)
1833   //  if gime(Event(ievent, "RA") is called the title of the branches is not checked anymore, "A" stands for any
1834   // This is a feature needed by PID to be able to reconstruct several times particles (each time a ther title is given)
1835   // from a given set of TrackSegments (with a given name)
1836   // This is why any is NOT used to read the branch of RecParticles
1837   // any migh have become obsolete : to be checked
1838   // See AliPHOSPIDv1    
1839
1840   TTree * treeR ;
1841   if(fToSplit){
1842     TFile * file = static_cast<TFile*>(gROOT->GetFile(fRecPointsFileName)); 
1843     if(!file) 
1844       file = TFile::Open(fRecPointsFileName) ;      
1845     // Get Digits Tree header from file
1846     TString treeName("TreeR") ;
1847     treeName += event ; 
1848     treeR = dynamic_cast<TTree*>(file->Get(treeName.Data()));
1849     if(!treeR){ // TreeR not found in header file
1850       if (fDebug)
1851         Warning("ReadTreeD", "-> Cannot find TreeR in %s", fRecPointsFileName.Data()) ;
1852       return 1;
1853     }
1854   }
1855   else
1856     treeR = gAlice->TreeR() ;
1857   
1858   // RecPoints 
1859   TObjArray * lob = static_cast<TObjArray*>(treeR->GetListOfBranches()) ;
1860   TIter next(lob) ; 
1861   TBranch * branch = 0 ; 
1862   TBranch * emcbranch = 0 ; 
1863   TBranch * cpvbranch = 0 ; 
1864   TBranch * clusterizerbranch = 0 ; 
1865   Bool_t phosemcrpfound = kFALSE, phoscpvrpfound = kFALSE, clusterizerfound = kFALSE ; 
1866
1867   
1868   while ( (branch = static_cast<TBranch*>(next())) && (!phosemcrpfound || !phoscpvrpfound || !clusterizerfound) ) {
1869     if(strcmp(branch->GetTitle(), fRecPointsTitle)==0 ) {
1870       if ( strcmp(branch->GetName(), "PHOSEmcRP")==0) {
1871         emcbranch = branch ; 
1872         phosemcrpfound = kTRUE ;
1873       }
1874       else if ( strcmp(branch->GetName(), "PHOSCpvRP")==0) {
1875         cpvbranch = branch ; 
1876         phoscpvrpfound = kTRUE ;
1877       }
1878       else if(strcmp(branch->GetName(), "AliPHOSClusterizer")==0){
1879         clusterizerbranch = branch ; 
1880         clusterizerfound = kTRUE ; 
1881       }
1882     }
1883   }
1884
1885   if ( !phosemcrpfound || !phoscpvrpfound || !clusterizerfound) {
1886     if (fDebug)
1887       Warning("ReadTreeR", "-> Cannot find RecPoints and/or Clusterizer with name %s", fRecPointsTitle.Data() ) ;
1888  
1889   } else { 
1890     if(!EmcRecPoints(fRecPointsTitle) ) 
1891       PostRecPoints(fRecPointsTitle) ;
1892     
1893     emcbranch->SetAddress(EmcRecPointsRef(fRecPointsTitle)) ;
1894     emcbranch->GetEntry(0) ;
1895
1896     cpvbranch->SetAddress(CpvRecPointsRef(fRecPointsTitle)) ; 
1897     cpvbranch->GetEntry(0) ;  
1898     
1899     if(!Clusterizer(fRecPointsTitle) )
1900       PostClusterizer(fRecPointsTitle) ;
1901     
1902     clusterizerbranch->SetAddress(ClusterizerRef(fRecPointsTitle)) ;
1903     clusterizerbranch->GetEntry(0) ;
1904   }
1905   
1906   //------------------- TrackSegments ---------------------
1907   next.Reset() ; 
1908   TBranch * tsbranch = 0 ; 
1909   TBranch * tsmakerbranch = 0 ; 
1910   Bool_t phostsfound = kFALSE, tsmakerfound = kFALSE ; 
1911   while ( (branch = static_cast<TBranch*>(next())) && (!phostsfound || !tsmakerfound) ) {
1912     if(strcmp(branch->GetTitle(), fTrackSegmentsTitle)==0 )  {
1913       if ( strcmp(branch->GetName(), "PHOSTS")==0){
1914         tsbranch = branch ; 
1915         phostsfound = kTRUE ;
1916       }
1917       else if(strcmp(branch->GetName(), "AliPHOSTrackSegmentMaker")==0) {
1918         tsmakerbranch = branch ; 
1919         tsmakerfound  = kTRUE ; 
1920       }
1921     }
1922   }
1923
1924   if ( !phostsfound || !tsmakerfound ) {
1925     if (fDebug)
1926       Warning("ReadTreeR", "-> Cannot find TrackSegments and/or TrackSegmentMaker with name %s", fTrackSegmentsTitle.Data() ) ;
1927   } else { 
1928     // Read and Post the TrackSegments
1929     if(!TrackSegments(fTrackSegmentsTitle))
1930       PostTrackSegments(fTrackSegmentsTitle) ;
1931     tsbranch->SetAddress(TrackSegmentsRef(fTrackSegmentsTitle)) ;
1932     tsbranch->GetEntry(0) ;
1933
1934     // Read and Post the TrackSegment Maker
1935     if(!TrackSegmentMaker(fTrackSegmentsTitle))
1936       PostTrackSegmentMaker(fTrackSegmentsTitle) ;
1937     tsmakerbranch->SetAddress(TSMakerRef(fTrackSegmentsTitle)) ;
1938     tsmakerbranch->GetEntry(0) ;
1939  }
1940   
1941   
1942   //------------ RecParticles ----------------------------
1943   next.Reset() ; 
1944   TBranch * rpabranch = 0 ; 
1945   TBranch * pidbranch = 0 ; 
1946   Bool_t phosrpafound = kFALSE, pidfound = kFALSE ; 
1947   
1948   while ( (branch = static_cast<TBranch*>(next())) && (!phosrpafound || !pidfound) ) 
1949     if(strcmp(branch->GetTitle(), fRecParticlesTitle)==0) {   
1950       if ( strcmp(branch->GetName(), "PHOSRP")==0) {   
1951         rpabranch = branch ; 
1952         phosrpafound = kTRUE ;
1953       }
1954       else if (strcmp(branch->GetName(), "AliPHOSPID")==0) {
1955         pidbranch = branch ; 
1956         pidfound  = kTRUE ; 
1957       }
1958     }
1959   
1960   if ( !phosrpafound || !pidfound ) {
1961     if (fDebug)
1962       Warning("ReadTreeR", "-> Cannot find RecParticles and/or PID with name %s", fRecParticlesTitle.Data() ) ; 
1963   } else { 
1964     // Read and Post the RecParticles
1965     if(!RecParticles(fRecParticlesTitle)) 
1966       PostRecParticles(fRecParticlesTitle) ;
1967     rpabranch->SetAddress(RecParticlesRef(fRecParticlesTitle)) ;
1968     rpabranch->GetEntry(0) ;
1969     // Read and Post the PID
1970     if(!PID(fRecParticlesTitle))
1971       PostPID(fRecParticlesTitle) ;
1972     pidbranch->SetAddress(PIDRef(fRecParticlesTitle)) ;
1973     pidbranch->GetEntry(0) ;
1974   }
1975
1976   if(gAlice->TreeR()!=treeR)
1977     treeR->Delete();
1978   return 0 ; 
1979 }
1980
1981 //____________________________________________________________________________ 
1982 Int_t AliPHOSGetter::ReadTreeS(const Int_t event)
1983 {
1984   // Reads the SDigits treeS from all files  
1985   // Files, which should be opened are listed in phosF
1986   // So, first get list of files
1987   TFolder * phosF = dynamic_cast<TFolder *>(fSDigitsFolder->FindObject("PHOS")) ;
1988   if (!phosF) 
1989     phosF = fSDigitsFolder->AddFolder("PHOS", "SDigits from PHOS") ; 
1990   TCollection * folderslist = phosF->GetListOfFolders() ; 
1991   
1992   // Now iterate over the list of files and read TreeS into Whiteboard
1993   TIter next(folderslist) ; 
1994   TFolder * folder = 0 ; 
1995   TFile * file; 
1996   TTree * treeS = 0;
1997   while ( (folder = static_cast<TFolder*>(next())) ) {
1998     TString fileName("") ;
1999     fileName = folder->GetName() ; 
2000     fileName.ReplaceAll("_","/") ; 
2001     file = static_cast<TFile*>(gROOT->GetFile(fileName)); 
2002     if(!file) 
2003       file = TFile::Open(fileName) ;      
2004     // Get SDigits Tree header from file
2005     TString treeName("TreeS") ;
2006     treeName += event ; 
2007     treeS = dynamic_cast<TTree*>(file->Get(treeName.Data()));
2008
2009     if(!treeS){ // TreeS not found in header file
2010       if (fDebug)
2011         Warning("ReadTreeS", "-> Cannot find TreeS in %s", fileName.Data()) ;
2012       return 1;
2013     }
2014     
2015     //set address of the SDigits and SDigitizer
2016     TBranch   * sdigitsBranch    = 0;
2017     TBranch   * sdigitizerBranch = 0;
2018     TBranch   * branch           = 0 ;  
2019     TObjArray * lob = static_cast<TObjArray*>(treeS->GetListOfBranches()) ;
2020     TIter next(lob) ; 
2021     Bool_t phosfound = kFALSE, sdigitizerfound = kFALSE ; 
2022
2023     while ( (branch = static_cast<TBranch*>(next())) && (!phosfound || !sdigitizerfound) ) {
2024       if ( (strcmp(branch->GetName(), "PHOS")==0) && (strcmp(branch->GetTitle(), fSDigitsTitle)==0) ) {
2025         phosfound = kTRUE ;
2026         sdigitsBranch = branch ; 
2027       }
2028       
2029       else if ( (strcmp(branch->GetName(), "AliPHOSSDigitizer")==0) && 
2030                 (strcmp(branch->GetTitle(), fSDigitsTitle)==0) ) {
2031         sdigitizerfound = kTRUE ; 
2032         sdigitizerBranch = branch ;
2033       }
2034     }
2035     if ( !phosfound || !sdigitizerfound ) {
2036       if (fDebug)
2037         Warning("ReadSDigits", "-> Digits and/or Digitizer branch with name %s not found", GetName()) ;
2038       return 2; 
2039     }   
2040     
2041     if ( !folder->FindObject(fSDigitsTitle) ){  
2042       TClonesArray * sdigits = new TClonesArray("AliPHOSDigit",1) ;
2043       sdigits->SetName(fSDigitsTitle) ;
2044       folder->Add(sdigits) ;
2045     }
2046
2047     ((TClonesArray*) (*SDigitsRef(fSDigitsTitle,folder->GetName())))->Clear() ;
2048     sdigitsBranch->SetAddress(SDigitsRef(fSDigitsTitle,folder->GetName())) ;
2049     sdigitsBranch->GetEntry(0) ;
2050     
2051     TString sdname(fSDigitsTitle) ;
2052     sdname+=":" ;
2053     sdname+=folder->GetName() ;
2054     if(!SDigitizer(sdname) ) 
2055       PostSDigitizer(fSDigitsTitle,folder->GetName()) ;
2056     sdigitizerBranch->SetAddress(SDigitizerRef(sdname)) ;
2057     sdigitizerBranch->GetEntry(0) ; 
2058     if(gAlice->TreeS()!=treeS)
2059       treeS->Delete();
2060   }    
2061   return 0 ; 
2062 }
2063
2064 //____________________________________________________________________________ 
2065 void AliPHOSGetter::ReadTreeS(TTree * treeS, Int_t input)
2066
2067 {  // Read the summable digits fron treeS()  
2068
2069
2070   TString filename("mergefile") ;
2071   filename+= input ;
2072
2073   TFolder * phosFolder = dynamic_cast<TFolder*>(fSDigitsFolder->FindObject("PHOS")) ; 
2074   if ( !phosFolder ) { 
2075    phosFolder = fSDigitsFolder->AddFolder("PHOS", "SDigits from PHOS") ; 
2076   } 
2077   TFolder * folder=(TFolder*)phosFolder->FindObject(filename) ;
2078   //set address of the SDigits and SDigitizer
2079   TBranch   * sdigitsBranch    = 0;
2080   TBranch   * sdigitizerBranch = 0;
2081   TBranch   * branch           = 0 ;  
2082   TObjArray * lob = (TObjArray*)treeS->GetListOfBranches() ;
2083   TIter next(lob) ; 
2084   Bool_t phosfound = kFALSE, sdigitizerfound = kFALSE ; 
2085   
2086   while ( (branch = (TBranch*)next()) && (!phosfound || !sdigitizerfound) ) {
2087     if ( strcmp(branch->GetName(), "PHOS")==0) {
2088       phosfound = kTRUE ;
2089       sdigitsBranch = branch ; 
2090     }
2091     
2092     else if ( strcmp(branch->GetName(), "AliPHOSSDigitizer")==0) {
2093       sdigitizerfound = kTRUE ; 
2094       sdigitizerBranch = branch ;
2095     }
2096   }
2097   if ( !phosfound || !sdigitizerfound ) {
2098     if (fDebug)
2099       Warning("ReadTreeS", "-> Digits and/or Digitizer branch not found") ;
2100     return ; 
2101   }   
2102   
2103   if (!folder || !(folder->FindObject(sdigitsBranch->GetTitle()) ) )
2104     PostSDigits(sdigitsBranch->GetTitle(),filename) ;
2105
2106   sdigitsBranch->SetAddress(SDigitsRef(sdigitsBranch->GetTitle(),folder->GetName())) ;
2107   sdigitsBranch->GetEntry(0) ;
2108   
2109   TString sdname(sdigitsBranch->GetTitle()) ;
2110   sdname+=":" ;
2111   sdname+=filename ;
2112   
2113   if(!SDigitizer(sdigitsBranch->GetTitle()) )
2114     PostSDigitizer(sdigitsBranch->GetTitle(),filename) ;
2115   sdigitizerBranch->SetAddress(SDigitizerRef(sdname)) ;
2116   sdigitizerBranch->GetEntry(0) ;
2117   if(gAlice->TreeS()!=treeS)
2118     treeS->Delete();
2119 }    
2120
2121
2122 //____________________________________________________________________________ 
2123 void AliPHOSGetter::ReadPrimaries()
2124 {
2125   // a lot simplified.... if 2 files are opened then we have a problem
2126
2127   TClonesArray * ar = 0  ; 
2128   if(! (ar = Primaries()) ) { 
2129     PostPrimaries() ;
2130     ar = Primaries() ; 
2131   }
2132   ar->Delete() ; 
2133   
2134   if (TreeK(fHeaderFile)) { // treeK found in header file
2135     if (fDebug) 
2136       Info("ReadPrimaries", "-> TreeK found in %s", fHeaderFile.Data() ); 
2137     fNPrimaries = gAlice->GetNtrack() ; 
2138     fAlice = 0 ; 
2139   
2140   } else { // treeK not found in header file
2141     Error("ReadPrimaries", "TreeK not found") ; 
2142     return ;  
2143   }
2144
2145   Int_t index = 0 ; 
2146   for (index = 0 ; index < fNPrimaries; index++) { 
2147     new ((*ar)[index]) TParticle(*(Primary(index)));
2148   }
2149 }
2150
2151 //____________________________________________________________________________ 
2152 void AliPHOSGetter::Event(const Int_t event, const char* opt)  
2153 {
2154   // Reads the content of all Tree's S, D and R
2155
2156   if (event >= gAlice->TreeE()->GetEntries() ) {
2157     Error("Event", "%d not found in TreeE !", event) ; 
2158     return ; 
2159   }
2160
2161   TBranch * btb = gAlice->TreeE()->GetBranch("AliPHOSBeamTestEvent") ;
2162   if(btb){
2163     if(!fBTE)
2164       fBTE = new AliPHOSBeamTestEvent() ;
2165     btb->SetAddress(&fBTE) ;
2166     btb->GetEntry(event) ;
2167   }
2168   else{
2169     if(fBTE){
2170       delete fBTE ;
2171       fBTE = 0 ;
2172     }
2173   }
2174
2175   Bool_t any = kFALSE ; 
2176   if (strstr(opt,"A") ) // do not check the title of the branches
2177     any = kTRUE; 
2178
2179   gAlice->GetEvent(event) ; 
2180
2181   if( strstr(opt,"R") )
2182     ReadTreeR(event) ;
2183
2184   if( strstr(opt,"D") )
2185     ReadTreeD(event) ;
2186
2187   if(strstr(opt,"S") )
2188     ReadTreeS(event) ;
2189
2190   if(strstr(opt,"H") )
2191     ReadTreeH() ;
2192    
2193   if( strstr(opt,"Q") )
2194     ReadTreeQA() ;
2195
2196   if( strstr(opt,"P") || (strcmp(opt,"")==0) )
2197     ReadPrimaries() ;
2198
2199  
2200 }
2201
2202 //____________________________________________________________________________ 
2203 TObject * AliPHOSGetter::ReturnO(TString what, TString name, TString file) const 
2204 {
2205   // get the object named "what" from the folder
2206   // folders are named like //Folders
2207
2208   if ( file.IsNull() ) 
2209     file = fHeaderFile ; 
2210   if( name.IsNull() )
2211     name = fBranchTitle ;
2212
2213   TFolder * folder = 0 ;
2214   TObject * phosO  = 0 ; 
2215
2216   if ( what.CompareTo("Primaries") == 0 ) {
2217     folder = dynamic_cast<TFolder *>(fPrimariesFolder->FindObject("Primaries")) ; 
2218     if (folder) 
2219       phosO  = dynamic_cast<TObject *>(folder->FindObject("Primaries")) ;  
2220     else 
2221       return 0 ; 
2222   }
2223   else if ( what.CompareTo("Hits") == 0 ) {
2224     folder = dynamic_cast<TFolder *>(fHitsFolder->FindObject("PHOS")) ; 
2225     if (folder) 
2226       phosO  = dynamic_cast<TObject *>(folder->FindObject("Hits")) ;  
2227   }
2228   else if ( what.CompareTo("SDigits") == 0 ) {
2229     if(fToSplit){
2230       file.Remove(file.Last('/')+1,file.Length()-file.Last('/')-1) ;
2231       file.ReplaceAll("/","_") ; 
2232       file+="PHOS.SDigits." ;
2233       if(name && (strcmp(name,"Default")!=0)){
2234         file+=name ;
2235         file+="." ;
2236       }
2237       file+="root" ;
2238     }
2239     TString path = "PHOS/" + file  ;
2240     folder = dynamic_cast<TFolder *>(fSDigitsFolder->FindObject(path.Data())) ; 
2241     if (folder) { 
2242       if (name.IsNull())
2243         name = fSDigitsTitle ; 
2244       phosO  = dynamic_cast<TObject *>(folder->FindObject(name)) ; 
2245     }
2246   }
2247   else if ( what.CompareTo("Digits") == 0 ){
2248     folder = dynamic_cast<TFolder *>(fDigitsFolder->FindObject("PHOS")) ; 
2249     if (folder) { 
2250       if (name.IsNull())
2251         name = fDigitsTitle ; 
2252       phosO  = dynamic_cast<TObject *>(folder->FindObject(name)) ; 
2253     } 
2254   }
2255   else if ( what.CompareTo("EmcRecPoints") == 0 ) {
2256     folder = dynamic_cast<TFolder *>(fRecoFolder->FindObject("PHOS/EMCARecPoints")) ; 
2257     if (folder) { 
2258       if (name.IsNull())
2259         name = fRecPointsTitle ; 
2260       phosO  = dynamic_cast<TObject *>(folder->FindObject(name)) ;
2261     } 
2262   }
2263   else if ( what.CompareTo("CpvRecPoints") == 0 ) {
2264     folder = dynamic_cast<TFolder *>(fRecoFolder->FindObject("PHOS/CPVRecPoints")) ; 
2265     if (folder) { 
2266       if (name.IsNull())
2267         name = fRecPointsTitle ; 
2268       phosO  = dynamic_cast<TObject *>(folder->FindObject(name)) ; 
2269     }   
2270   }
2271   else if ( what.CompareTo("TrackSegments") == 0 ) {
2272     folder = dynamic_cast<TFolder *>(fRecoFolder->FindObject("PHOS/TrackSegments")) ; 
2273     if (folder) { 
2274       if (name.IsNull())
2275         name = fTrackSegmentsTitle ; 
2276       phosO  = dynamic_cast<TObject *>(folder->FindObject(name)) ; 
2277     }   
2278   }
2279   else if ( what.CompareTo("RecParticles") == 0 ) {
2280     folder = dynamic_cast<TFolder *>(fRecoFolder->FindObject("PHOS/RecParticles")) ; 
2281    if (folder) { 
2282       if (name.IsNull())
2283         name = fRecParticlesTitle ; 
2284       phosO  = dynamic_cast<TObject *>(folder->FindObject(name)) ;
2285     }   
2286  }
2287   else if ( what.CompareTo("Alarms") == 0 ){ 
2288     if (name.IsNull() ) 
2289       phosO = dynamic_cast<TObject *>(fQAFolder->FindObject("PHOS")) ;  
2290     else {
2291       folder = dynamic_cast<TFolder *>(fQAFolder->FindObject("PHOS")) ; 
2292       if (!folder) 
2293         phosO = 0 ; 
2294       else 
2295         phosO = dynamic_cast<TObject *>(folder->FindObject(name)) ;  
2296     }
2297   }
2298   if (!phosO) {
2299     if(fDebug)
2300       Warning("ReturnO", "Object %s not found in PHOS", what.Data() ) ; 
2301     return 0 ;
2302   }
2303
2304   return phosO ;
2305 }
2306   
2307 //____________________________________________________________________________ 
2308 const TTask * AliPHOSGetter::ReturnT(TString what, TString name) const 
2309 {
2310   // get the TTask named "what" from the folder
2311   // folders are named like //Folders/Tasks/what/PHOS/name
2312
2313   TString search(what) ; 
2314   if ( what.CompareTo("Clusterizer") == 0 ) 
2315     search = "Reconstructioner" ; 
2316   else if ( what.CompareTo("TrackSegmentMaker") == 0 ) 
2317     search = "Reconstructioner" ; 
2318   else if ( what.CompareTo("PID") == 0 ) 
2319     search = "Reconstructioner" ; 
2320   else if ( what.CompareTo("QATasks") == 0 ) 
2321     search = "QA" ; 
2322
2323   TTask * tasks = dynamic_cast<TTask*>(fTasksFolder->FindObject(search)) ; 
2324
2325   if (!tasks) {
2326     Error("ReturnT", "Task %s not found !", what.Data() ) ;  
2327     return 0 ; 
2328   }
2329
2330   TTask * phosT = dynamic_cast<TTask*>(tasks->GetListOfTasks()->FindObject("PHOS")) ; 
2331   if (!phosT) { 
2332      Error("ReturnT", "Task %s/PHOS not found !", what.Data() ) ;  
2333     return 0 ; 
2334   }
2335   
2336   TList * list = phosT->GetListOfTasks() ; 
2337  
2338   if (what.CompareTo("SDigitizer") == 0) {  
2339     if ( name.IsNull() )
2340       name =  fSDigitsTitle ; 
2341   } else  if (what.CompareTo("Digitizer") == 0){ 
2342     if ( name.IsNull() )
2343       name =  fDigitsTitle ;
2344   } else  if (what.CompareTo("Clusterizer") == 0){ 
2345     if ( name.IsNull() )
2346       name =  fRecPointsTitle ;
2347     name.Append(":clu") ;
2348   }
2349   else  if (what.CompareTo("TrackSegmentMaker") == 0){ 
2350     if ( name.IsNull() )
2351       name =  fTrackSegmentsTitle ;
2352     name.Append(":tsm") ;
2353   }
2354   else  if (what.CompareTo("PID") == 0){ 
2355     if ( name.IsNull() )
2356       name =  fRecParticlesTitle ;
2357     name.Append(":pid") ;
2358   }
2359   else  if (what.CompareTo("QATasks") == 0){ 
2360     if ( name.IsNull() )
2361       return phosT ;
2362   }
2363   
2364   TIter it(list) ;
2365   TTask * task = 0 ; 
2366   while((task = static_cast<TTask *>(it.Next()) )){
2367     TString taskname(task->GetName()) ;
2368     if(taskname.BeginsWith(name))
2369       return task ;
2370   }
2371   
2372   if(fDebug)
2373     Warning("ReturnT", "-> Task %s/PHOS/%s not found", search.Data(), name.Data() ) ; 
2374   return 0 ;
2375 }
2376
2377 //____________________________________________________________________________ 
2378 void AliPHOSGetter::RemoveTask(TString opt, TString name) const 
2379 {
2380   // remove a task from the folder
2381   // path is fTasksFolder/SDigitizer/PHOS/name
2382   
2383   TTask * task = 0 ; 
2384   TTask * phos = 0 ; 
2385   TList * lofTasks = 0 ; 
2386
2387   if (opt == "S") { // SDigitizer
2388     task = dynamic_cast<TTask*>(fTasksFolder->FindObject("SDigitizer")) ;
2389     if (!task) 
2390       return ; 
2391   }
2392   else if (opt == "D") { // Digitizer
2393     task = dynamic_cast<TTask*>(fTasksFolder->FindObject("Digitizer")) ;
2394     if (!task) 
2395       return ; 
2396   }
2397   else if (opt == "C" || opt == "T" || opt == "P"  ) { // Clusterizer, TrackSegmentMaker, PID
2398     task = dynamic_cast<TTask*>(fTasksFolder->FindObject("Reconstructioner")) ;
2399     if (!task) 
2400       return ; 
2401   }
2402   else {
2403     Warning("RemoveTask", "Unknown option %s", opt.Data() ); 
2404     return ; 
2405   }
2406   phos =  dynamic_cast<TTask*>(task->GetListOfTasks()->FindObject("PHOS")) ;
2407   if (!phos)
2408     return ; 
2409   lofTasks = phos->GetListOfTasks() ;
2410   if (!lofTasks) 
2411     return ; 
2412   TObject * obj = lofTasks->FindObject(name) ; 
2413   if (obj) 
2414     lofTasks->Remove(obj) ;
2415    
2416 }
2417
2418 //____________________________________________________________________________ 
2419 void AliPHOSGetter::RemoveObjects(TString opt, TString name) const 
2420 {
2421   // remove SDigits from the folder
2422   // path is fSDigitsFolder/fHeaderFileName/name
2423
2424   TFolder * phos     = 0 ; 
2425   TFolder * phosmain = 0 ; 
2426
2427   if (opt == "H") { // Hits
2428     phos = dynamic_cast<TFolder*>(fHitsFolder->FindObject("PHOS")) ;
2429     if (!phos) 
2430       return ;
2431     name = "Hits" ; 
2432   }
2433
2434   else if ( opt == "S") { // SDigits
2435     phosmain = dynamic_cast<TFolder*>(fSDigitsFolder->FindObject("PHOS")) ;
2436     if (!phosmain) 
2437       return ;
2438     phos = dynamic_cast<TFolder*>(phosmain->FindObject(fHeaderFile)) ;
2439     if (!phos) 
2440       return ;
2441   }
2442   
2443   else if (opt == "D") { // Digits
2444     phos = dynamic_cast<TFolder*>(fDigitsFolder->FindObject("PHOS")) ;
2445     if (!phos) 
2446       return ;
2447   }
2448
2449   else if (opt == "RE") { // EMCARecPoints
2450     phos = dynamic_cast<TFolder*>(fRecoFolder->FindObject("PHOS/EMCARecPoints")) ;
2451     if (!phos) 
2452       return ;
2453   }
2454
2455   else if (opt == "RC") { // CPVRecPoints
2456     phos = dynamic_cast<TFolder*>(fRecoFolder->FindObject("PHOS/CPVRecPoints")) ;
2457     if (!phos) 
2458       return ;
2459   }  
2460
2461   else if (opt == "T") { // TrackSegments
2462     phos = dynamic_cast<TFolder*>(fRecoFolder->FindObject("PHOS/TrackSegments")) ;
2463     if (!phos) 
2464       return ;
2465   }
2466
2467   else if (opt == "P") { // RecParticles
2468     phos = dynamic_cast<TFolder*>(fRecoFolder->FindObject("PHOS/RecParticles")) ;
2469     if (!phos) 
2470       return ;
2471   }
2472   
2473   else {
2474     Warning("RemoveObjects", "Unknown option %s", opt.Data() ) ; 
2475     return ; 
2476   }
2477   
2478   TObjArray * ar  = dynamic_cast<TObjArray*>(phos->FindObject(name)) ; 
2479   if (ar) { 
2480     phos->Remove(ar) ;
2481     ar->Delete() ; 
2482     delete ar ; 
2483   }
2484
2485   if (opt == "S") 
2486     phosmain->Remove(phos) ; 
2487   
2488 }
2489
2490 //____________________________________________________________________________ 
2491 void AliPHOSGetter::RemoveSDigits() const 
2492 {
2493   TFolder * phos= dynamic_cast<TFolder*>(fSDigitsFolder->FindObject("PHOS")) ;
2494   if (!phos) 
2495     return ;
2496   
2497   phos->SetOwner() ; 
2498   phos->Clear() ; 
2499 }
2500
2501 //____________________________________________________________________________ 
2502 void AliPHOSGetter::CleanWhiteBoard(void){
2503
2504   TFolder * phosmain = 0 ; 
2505   TFolder * phos ;
2506   TObjArray * ar ;
2507   TList * lofTasks = 0 ; 
2508   TTask * task = 0 ; 
2509   TTask * phost = 0 ; 
2510   
2511   // Hits  
2512   phos = dynamic_cast<TFolder*>(fHitsFolder->FindObject("PHOS")) ;
2513   if (phos){  
2514     TObjArray * ar  = dynamic_cast<TObjArray*>(phos->FindObject("Hits")) ; 
2515     if (ar) { 
2516       phos->Remove(ar) ;
2517       ar->Delete() ; 
2518       delete ar ; 
2519     }
2520   }
2521   
2522   // SDigits
2523   phosmain = dynamic_cast<TFolder*>(fSDigitsFolder->FindObject("PHOS")) ;
2524   if (phosmain){ 
2525     phos = dynamic_cast<TFolder*>(phosmain->FindObject(fHeaderFile)) ;
2526     if (phos) {
2527       ar  = dynamic_cast<TObjArray*>(phos->FindObject(fSDigitsTitle)) ; 
2528       if (ar) { 
2529         phos->Remove(ar) ;
2530         ar->Delete() ; 
2531         delete ar ; 
2532       }
2533     }
2534     phosmain->Remove(phos) ; 
2535   }
2536
2537   
2538   // Digits
2539   phos = dynamic_cast<TFolder*>(fDigitsFolder->FindObject("PHOS")) ;
2540   if (phos){ 
2541     ar  = dynamic_cast<TObjArray*>(phos->FindObject(fDigitsTitle)) ; 
2542     if (ar) { 
2543       phos->Remove(ar) ;
2544       ar->Delete() ; 
2545       delete ar ; 
2546     }
2547   }
2548
2549
2550   // EMCARecPoints
2551   phos = dynamic_cast<TFolder*>(fRecoFolder->FindObject("PHOS/EMCARecPoints")) ;
2552   if (phos){ 
2553     ar  = dynamic_cast<TObjArray*>(phos->FindObject(fRecPointsTitle)) ; 
2554     if (ar) { 
2555       phos->Remove(ar) ;
2556       ar->Delete() ; 
2557       delete ar ; 
2558     }
2559   }
2560
2561   
2562   // CPVRecPoints
2563   phos = dynamic_cast<TFolder*>(fRecoFolder->FindObject("PHOS/CPVRecPoints")) ;
2564   if (phos){ 
2565     ar  = dynamic_cast<TObjArray*>(phos->FindObject(fRecPointsTitle)) ; 
2566     if (ar) { 
2567       phos->Remove(ar) ;
2568       ar->Delete() ; 
2569       delete ar ; 
2570     }
2571   }  
2572
2573   
2574   // TrackSegments
2575   phos = dynamic_cast<TFolder*>(fRecoFolder->FindObject("PHOS/TrackSegments")) ;
2576   if (phos) { 
2577     ar  = dynamic_cast<TObjArray*>(phos->FindObject(fTrackSegmentsTitle)) ; 
2578     if (ar) { 
2579       phos->Remove(ar) ;
2580       ar->Delete() ; 
2581       delete ar ; 
2582     }
2583   }
2584   
2585
2586
2587   // RecParticles
2588   phos = dynamic_cast<TFolder*>(fRecoFolder->FindObject("PHOS/RecParticles")) ;
2589   if (phos){ 
2590     ar  = dynamic_cast<TObjArray*>(phos->FindObject(fRecParticlesTitle)) ; 
2591     if (ar) { 
2592       phos->Remove(ar) ;
2593       ar->Delete() ; 
2594       delete ar ; 
2595     }
2596   }
2597
2598
2599   //---- Now Tasks ----------- 
2600
2601   TObject * obj ;
2602   TString sdname(fSDigitsTitle);
2603   
2604   // Digitizer
2605   task = dynamic_cast<TTask*>(fTasksFolder->FindObject("Digitizer")) ;
2606   if (task){ 
2607     phost =  dynamic_cast<TTask*>(task->GetListOfTasks()->FindObject("PHOS")) ;
2608     if (phost){
2609       lofTasks = phost->GetListOfTasks() ;
2610       if (lofTasks){ 
2611         obj = lofTasks->FindObject(sdname.Data()) ; 
2612         if (obj) 
2613           lofTasks->Remove(obj) ;
2614       }
2615     }      
2616   }
2617   
2618
2619   sdname.Append(":") ;
2620   // Clusterizer, TrackSegmentMaker, PID
2621   task = dynamic_cast<TTask*>(fTasksFolder->FindObject("Reconstructioner")) ;
2622   if (task){ 
2623     phost =  dynamic_cast<TTask*>(task->GetListOfTasks()->FindObject("PHOS")) ;
2624     if (phost){
2625       lofTasks = phost->GetListOfTasks() ;
2626       TIter next(lofTasks);
2627       while((obj=next())){ 
2628         TString oname(obj->GetName()) ;
2629         if (oname.BeginsWith(sdname)){ 
2630           lofTasks->Remove(obj) ;
2631         }
2632       }
2633     }  
2634   }
2635
2636
2637   // SDigitizer
2638   sdname.Append(fHeaderFile) ;
2639   task = dynamic_cast<TTask*>(fTasksFolder->FindObject("SDigitizer")) ;
2640   if (task) {
2641     phost =  dynamic_cast<TTask*>(task->GetListOfTasks()->FindObject("PHOS")) ;
2642     if (phost){
2643       lofTasks = phost->GetListOfTasks() ;
2644       if (lofTasks){ 
2645         obj = lofTasks->FindObject(sdname.Data()) ; 
2646         if (obj) 
2647           lofTasks->Remove(obj) ;
2648       }
2649     }
2650   }  
2651
2652 }
2653 //____________________________________________________________________________ 
2654 void AliPHOSGetter::SetTitle(const char * branchTitle ) 
2655 {
2656   fBranchTitle        = branchTitle ;
2657   fSDigitsTitle       = branchTitle ; 
2658   fDigitsTitle        = branchTitle ; 
2659   fRecPointsTitle     = branchTitle ; 
2660   fRecParticlesTitle  = branchTitle ; 
2661   fTrackSegmentsTitle = branchTitle ; 
2662   if(fToSplit){
2663     //First - extract full path if necessary
2664     TString sFileName(fHeaderFile) ;
2665     Ssiz_t islash = sFileName.Last('/') ;
2666     if(islash<sFileName.Length())
2667       sFileName.Remove(islash+1,sFileName.Length()) ;
2668     else
2669       sFileName="" ;
2670     //Now construct file names
2671     fSDigitsFileName       = sFileName ;
2672     fDigitsFileName        = sFileName ; 
2673     fRecPointsFileName     = sFileName ; 
2674     fRecParticlesFileName  = sFileName ; 
2675     fTrackSegmentsFileName = sFileName ; 
2676     fSDigitsFileName      += "PHOS.SDigits." ;
2677     fDigitsFileName       += "PHOS.Digits." ; 
2678     fRecPointsFileName    += "PHOS.RecData." ; 
2679     fTrackSegmentsFileName+= "PHOS.RecData." ; 
2680     fRecParticlesFileName += "PHOS.RecData." ; 
2681     if((strcmp(fBranchTitle.Data(),"Default")!=0)&&(strcmp(fBranchTitle.Data(),"")!=0)){
2682       fSDigitsFileName      += fBranchTitle ;
2683       fSDigitsFileName      += "." ;
2684       fDigitsFileName       += fBranchTitle ; 
2685       fDigitsFileName       += "." ; 
2686       fRecPointsFileName    += fBranchTitle ; 
2687       fRecPointsFileName    += "." ; 
2688       fRecParticlesFileName += fBranchTitle ; 
2689       fRecParticlesFileName += "." ; 
2690       fTrackSegmentsFileName+= fBranchTitle ; 
2691       fTrackSegmentsFileName+= "." ; 
2692     }
2693     fSDigitsFileName      += "root" ;
2694     fDigitsFileName       += "root" ; 
2695     fRecPointsFileName    += "root" ; 
2696     fRecParticlesFileName += "root" ; 
2697     fTrackSegmentsFileName+= "root" ; 
2698   }else{
2699     fSDigitsFileName       = fHeaderFile ;
2700
2701     fDigitsFileName        = "" ; 
2702     fRecPointsFileName     = "" ; 
2703     fRecParticlesFileName  = "" ; 
2704     fTrackSegmentsFileName = "" ; 
2705   }
2706   TFolder * phosFolder ; 
2707   phosFolder = dynamic_cast<TFolder*>(fHitsFolder->FindObject("PHOS")) ; 
2708   if ( !phosFolder ) 
2709     phosFolder = fHitsFolder->AddFolder("PHOS", "Hits from PHOS") ; 
2710
2711   phosFolder = dynamic_cast<TFolder*>(fSDigitsFolder->FindObject("PHOS")) ;
2712   if ( !phosFolder ) 
2713     phosFolder = fSDigitsFolder->AddFolder("PHOS", "SDigits from PHOS") ; 
2714
2715   //Make folder for SDigits
2716   fSDigitsFileName.ReplaceAll("/","_") ;
2717   phosFolder->AddFolder(fSDigitsFileName.Data(),"");
2718
2719   phosFolder  = dynamic_cast<TFolder*>(fDigitsFolder->FindObject("PHOS")) ;
2720   if ( !phosFolder ) 
2721     phosFolder = fDigitsFolder->AddFolder("PHOS", "Digits from PHOS") ;  
2722
2723   phosFolder  = dynamic_cast<TFolder*>(fRecoFolder->FindObject("PHOS")) ; 
2724   if ( !phosFolder )
2725     phosFolder = fRecoFolder->AddFolder("PHOS", "Reconstructed data from PHOS") ;  
2726   
2727 }
2728 //____________________________________________________________________________ 
2729 void AliPHOSGetter::CloseSplitFiles(void){
2730   TFile * file ;
2731   file = static_cast<TFile*>(gROOT->GetFile(fSDigitsFileName.Data() ) ) ;
2732   if(file)
2733     file->Close() ;
2734   file = static_cast<TFile*>(gROOT->GetFile(fDigitsFileName.Data() ) ) ;
2735   if(file)
2736     file->Close() ;
2737   file = static_cast<TFile*>(gROOT->GetFile(fRecPointsFileName.Data() ) ) ;
2738   if(file)
2739     file->Close() ;
2740   file = static_cast<TFile*>(gROOT->GetFile(fTrackSegmentsFileName.Data() ) ) ;
2741   if(file)
2742     file->Close() ;
2743   file = static_cast<TFile*>(gROOT->GetFile(fRecParticlesFileName.Data() ) ) ;
2744   if(file)
2745     file->Close() ;
2746
2747 }