]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSLoader.cxx
Extracting PHOS and EMCAL trackers from the correspondig reconstructors (Yu.Belikov)
[u/mrichter/AliRoot.git] / PHOS / AliPHOSLoader.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //_________________________________________________________________________
19 //  A singleton. This class should be used in the analysis stage to get 
20 //  reconstructed objects: Digits, RecPoints, TrackSegments and RecParticles,
21 //  instead of directly reading them from galice.root file. This container 
22 //  ensures, that one reads Digits, made of these particular digits, RecPoints, 
23 //  made of these particular RecPoints, TrackSegments and RecParticles. 
24 //  This becomes non trivial if there are several identical branches, produced with
25 //  different set of parameters. 
26 //
27 //  An example of how to use (see also class AliPHOSAnalyser):
28 //  for(Int_t irecp = 0; irecp < gime->NRecParticles() ; irecp++)
29 //     AliPHOSRecParticle * part = gime->RecParticle(1) ;
30 //     ................
31 //  please->GetEvent(event) ;    // reads new event from galice.root
32 //                  
33 //*-- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (RRC KI & SUBATECH)
34 //*--         Completely redesigned by Dmitri Peressounko March 2001  
35 //
36 //*-- YS June 2001 : renamed the original AliPHOSIndexToObject and make
37 //*--         systematic usage of TFolders without changing the interface        
38 //////////////////////////////////////////////////////////////////////////////
39
40
41 // --- ROOT system ---
42
43 #include "TFile.h"
44 #include "TTree.h"
45 #include "TROOT.h"
46
47 // --- Standard library ---
48
49 // --- AliRoot header files ---
50 #include "AliLog.h"
51 #include "AliPHOSLoader.h"
52 #include "AliPHOS.h"
53 #include "AliPHOSHit.h"
54 #include "AliPHOSCalibrationDB.h"
55 #include "AliPHOSGetter.h"
56
57 ClassImp(AliPHOSLoader)
58
59
60 const TString AliPHOSLoader::fgkHitsName("HITS");//Name for TClonesArray with hits from one event
61 const TString AliPHOSLoader::fgkSDigitsName("SDIGITS");//Name for TClonesArray 
62 const TString AliPHOSLoader::fgkDigitsName("DIGITS");//Name for TClonesArray 
63 const TString AliPHOSLoader::fgkEmcRecPointsName("EMCRECPOINTS");//Name for TClonesArray 
64 const TString AliPHOSLoader::fgkCpvRecPointsName("CPVRECPOINTS");//Name for TClonesArray 
65 const TString AliPHOSLoader::fgkTracksName("TRACKS");//Name for TClonesArray 
66 const TString AliPHOSLoader::fgkRecParticlesName("RECPARTICLES");//Name for TClonesArray
67
68 const TString AliPHOSLoader::fgkEmcRecPointsBranchName("PHOSEmcRP");//Name for branch with EMC Reconstructed Points
69 const TString AliPHOSLoader::fgkCpvRecPointsBranchName("PHOSCpvRP");//Name for branch with CPV Reconstructed Points
70 const TString AliPHOSLoader::fgkTrackSegmentsBranchName("PHOSTS");//Name for branch with TrackSegments
71 const TString AliPHOSLoader::fgkRecParticlesBranchName("PHOSRP");//Name for branch with Reconstructed Particles
72 //____________________________________________________________________________ 
73 AliPHOSLoader::AliPHOSLoader()
74  {
75   fDebug = 0;
76  }
77 //____________________________________________________________________________ 
78 AliPHOSLoader::AliPHOSLoader(const Char_t *detname,const Char_t *eventfoldername):
79       AliLoader(detname,eventfoldername)
80 {
81   fDebug=0;
82 }
83 //____________________________________________________________________________ 
84
85 AliPHOSLoader::~AliPHOSLoader()
86 {
87   //remove and delete arrays
88   Clean(fgkHitsName);
89   Clean(fgkSDigitsName);
90   Clean(fgkDigitsName);
91   Clean(fgkEmcRecPointsName);
92   Clean(fgkCpvRecPointsName);
93   Clean(fgkTracksName);
94   Clean(fgkRecParticlesName);
95   CleanFolders() ;
96   // set to 0x0 the objgetter in AliGetter ... weird isn it !
97   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; // (AliLoader::GetRunLoader()->GetFileName()).Data()) ; 
98   if (gime) 
99     gime->Reset() ;
100 }
101
102 //____________________________________________________________________________ 
103 void AliPHOSLoader::CleanFolders()
104  {
105    CleanRecParticles();
106    AliLoader::CleanFolders();
107  }
108
109 //____________________________________________________________________________ 
110 Int_t AliPHOSLoader::SetEvent()
111 {
112 //Cleans loaded stuff and and sets Files and Directories
113 // do not post any data to folder/tasks
114
115
116  Int_t retval = AliLoader::SetEvent();
117   if (retval)
118    {
119      AliError("returned error");
120      return retval;
121    }
122
123
124   if (Hits()) Hits()->Clear();
125   if (SDigits()) SDigits()->Clear();
126   if (Digits()) Digits()->Clear();
127   if (EmcRecPoints()) EmcRecPoints()->Clear();
128   if (CpvRecPoints()) CpvRecPoints()->Clear();
129   if (TrackSegments()) TrackSegments()->Clear();
130   if (RecParticles()) RecParticles()->Clear();
131    
132   return 0;
133 }
134
135 //____________________________________________________________________________ 
136 Int_t AliPHOSLoader::GetEvent()
137 {
138 //Overloads GetEvent method called by AliRunLoader::GetEvent(Int_t) method
139 //to add Rec Particles specific for PHOS
140
141 //First call the original method to get whatever from std. setup is needed
142   Int_t retval;
143   
144   retval = AliLoader::GetEvent();
145   if (retval)
146    {
147      AliError("returned error");
148      return retval;
149    }
150   
151   if (GetHitsDataLoader()->GetBaseDataLoader()->IsLoaded()) ReadHits();
152   if (GetSDigitsDataLoader()->GetBaseDataLoader()->IsLoaded()) ReadSDigits();
153   if (GetDigitsDataLoader()->GetBaseDataLoader()->IsLoaded()) ReadDigits();
154   if (GetRecPointsDataLoader()->GetBaseDataLoader()->IsLoaded()) ReadRecPoints();
155   if (GetTracksDataLoader()->GetBaseDataLoader()->IsLoaded()) ReadTracks();
156   if (GetRecParticlesDataLoader()->GetBaseDataLoader()->IsLoaded()) ReadRecParticles();
157
158
159 //Now, check if RecPart were loaded  
160   return 0;
161 }
162
163 //____________________________________________________________________________ 
164 const AliPHOS * AliPHOSLoader::PHOS() 
165 {
166   // returns the PHOS object 
167   AliPHOS * phos = dynamic_cast<AliPHOS*>(GetModulesFolder()->FindObject(fDetectorName));
168   if ( phos == 0x0) 
169     if (fDebug)
170       cout << "WARNING: AliPHOSLoader::PHOS -> PHOS module not found in Folders" << endl ; 
171   return phos ; 
172 }  
173
174 //____________________________________________________________________________ 
175 const AliPHOSGeometry * AliPHOSLoader::PHOSGeometry() 
176 {
177   // Return PHOS geometry
178   AliPHOSGeometry * rv = 0 ; 
179   if (PHOS() )
180     rv =  PHOS()->GetGeometry();
181   return rv ; 
182
183
184
185 //____________________________________________________________________________ 
186 Int_t AliPHOSLoader::LoadHits(Option_t* opt)
187 {  
188 //------- Hits ----------------------
189 //Overload (extends) LoadHits implemented in AliLoader
190 //
191   Int_t res;
192   
193   //First call the AliLoader's method to send the TreeH to folder
194   res = AliLoader::LoadHits(opt);
195   
196   if (res)
197    {//oops, error
198      AliError("returned error");
199      return res;
200    }
201
202   //read the data from tree in folder and send it to folder
203   res = ReadHits();
204   return 0;
205 }
206
207
208 //____________________________________________________________________________ 
209 Int_t AliPHOSLoader::LoadSDigits(Option_t* opt)
210 {
211   //---------- SDigits -------------------------
212   Int_t res;
213   //First call the AliLoader's method to send the TreeS to folder
214   res = AliLoader::LoadSDigits(opt);
215   if (res)
216    {//oops, error
217      AliError("returned error");
218      return res;
219    }
220   return ReadSDigits();
221    
222
223 //____________________________________________________________________________ 
224 Int_t AliPHOSLoader::LoadDigits(Option_t* opt)
225
226   //---------- Digits -------------------------
227   Int_t res;
228   //First call the AliLoader's method to send the TreeS to folder
229   res = AliLoader::LoadDigits(opt);
230   if (res)
231    {//oops, error
232      AliError("returned error");
233      return res;
234    }
235   return ReadDigits();
236 }
237 //____________________________________________________________________________ 
238 Int_t AliPHOSLoader::LoadRecPoints(Option_t* opt) 
239 {
240   // -------------- RecPoints -------------------------------------------
241   Int_t res;
242   //First call the AliLoader's method to send the TreeR to folder
243   res = AliLoader::LoadRecPoints(opt);
244   if (res)
245    {//oops, error
246      AliError("returned error");
247      return res;
248    }
249
250   TFolder * phosFolder = GetDetectorDataFolder();
251   if ( phosFolder  == 0x0 ) 
252    {
253      AliError("Can not get detector data folder");
254      return 1;
255    }
256   return ReadRecPoints();
257 }
258 //____________________________________________________________________________ 
259
260 Int_t  AliPHOSLoader::LoadTracks(Option_t* opt)
261 {
262  //Loads Tracks: Open File, Reads Tree and posts, Read Data and Posts
263   AliDebug(1, Form("opt = %s",opt));
264  Int_t res;
265  res = AliLoader::LoadTracks(opt);
266  if (res)
267    {//oops, error
268       AliError("returned error");
269       return res;
270    }  
271  return ReadTracks();
272 }
273
274 //____________________________________________________________________________ 
275 Int_t AliPHOSLoader::LoadRecParticles(Option_t* opt) 
276 {
277   // -------------- RecPoints -------------------------------------------
278   Int_t res;
279   //First call the AliLoader's method to send the TreeT to folder
280   res = AliLoader::LoadRecParticles(opt);
281   if (res)
282    {//oops, error
283      AliError("returned error");
284      return res;
285    }
286   return ReadRecParticles();
287 }
288
289 //____________________________________________________________________________ 
290
291 Int_t AliPHOSLoader::PostHits()
292 {
293   // -------------- Hits -------------------------------------------
294   Int_t reval = AliLoader::PostHits();
295   if (reval)
296    {
297      AliError("returned error");
298      return reval;
299    }
300   return ReadHits();
301 }
302 //____________________________________________________________________________ 
303
304 Int_t AliPHOSLoader::PostSDigits()
305 {
306   // -------------- SDigits -------------------------------------------
307   Int_t reval = AliLoader::PostSDigits();
308   if (reval)
309    {
310      AliError("returned error");
311      return reval;
312    }
313   return ReadSDigits();
314 }
315 //____________________________________________________________________________ 
316
317 Int_t AliPHOSLoader::PostDigits()
318 {
319   // -------------- Digits -------------------------------------------
320   Int_t reval = AliLoader::PostDigits();
321   if (reval)
322    {
323      AliError("returned error");
324      return reval;
325    }
326   return ReadDigits();
327 }
328 //____________________________________________________________________________ 
329
330 Int_t AliPHOSLoader::PostRecPoints()
331 {
332   // -------------- RecPoints -------------------------------------------
333   Int_t reval = AliLoader::PostRecPoints();
334   if (reval)
335    {
336      AliError("returned error");
337      return reval;
338    }
339   return ReadRecPoints();
340 }
341
342 //____________________________________________________________________________ 
343
344 Int_t AliPHOSLoader::PostRecParticles()
345 {
346   // -------------- RecParticles -------------------------------------------
347   Int_t reval = AliLoader::PostRecParticles();
348   if (reval)
349    {
350      AliError("returned error");
351      return reval;
352    }
353   return ReadRecParticles();
354 }
355 //____________________________________________________________________________ 
356
357 Int_t AliPHOSLoader::PostTracks()
358 {
359   // -------------- Tracks -------------------------------------------
360   Int_t reval = AliLoader::PostTracks();
361   if (reval)
362    {
363      AliError("returned error");
364      return reval;
365    }
366   return ReadTracks();
367 }
368 //____________________________________________________________________________ 
369
370
371
372 //____________________________________________________________________________ 
373 Int_t AliPHOSLoader::ReadHits()
374 {
375 // If there is no Clones Array in folder creates it and sends to folder
376 // then tries to read
377 // Reads the first entry of PHOS branch in hit tree TreeH()
378 // Reads data from TreeH and stores it in TClonesArray that sits in DetectorDataFolder
379 //
380   TObject** hitref = HitsRef();
381   if(hitref == 0x0)
382    {
383      MakeHitsArray();
384      hitref = HitsRef();
385    }
386
387   TClonesArray* hits = dynamic_cast<TClonesArray*>(*hitref);
388
389   TTree* treeh = TreeH();
390   
391   if(treeh == 0)
392    {
393     AliError("Cannot read TreeH from folder");
394     return 1;
395   }
396   
397   TBranch * hitsbranch = treeh->GetBranch(fDetectorName);
398   if (hitsbranch == 0) 
399    {
400     AliError("Cannot find branch PHOS"); 
401     return 1;
402   }
403
404   AliDebug(1, "Reading Hits");
405   
406   if (hitsbranch->GetEntries() > 1)
407    {
408     TClonesArray * tempo =  new TClonesArray("AliPHOSHit",1000);
409
410     hitsbranch->SetAddress(&tempo);
411     Int_t index = 0 ; 
412     Int_t i = 0 ;
413     for (i = 0 ; i < hitsbranch->GetEntries(); i++) 
414      {
415       hitsbranch->GetEntry(i) ;
416       Int_t j = 0 ;
417       for ( j = 0 ; j < tempo->GetEntries() ; j++) 
418        {
419          AliPHOSHit* hit = (AliPHOSHit*)tempo->At(j); 
420          new((*hits)[index]) AliPHOSHit( *hit ) ;
421          index++ ; 
422        }
423      }
424     delete tempo;
425    }
426   else 
427    {
428     hitsbranch->SetAddress(hitref);
429     hitsbranch->GetEntry(0) ;
430    }
431
432   return 0;
433 }
434 //____________________________________________________________________________ 
435 Int_t AliPHOSLoader::ReadSDigits()
436 {
437 // Read the summable digits tree TreeS():
438 // Check if TClones is in folder
439 // if not create and add to folder
440 // connect to tree if available
441 // Read the data
442
443   TObject** sdref = SDigitsRef();
444   if(sdref == 0x0)
445    {
446      MakeSDigitsArray();
447      sdref = SDigitsRef();
448    }
449    
450   TTree * treeS = TreeS();
451   if(treeS==0)
452    {
453      //May happen if file is truncated or new in LoadSDigits
454      //AliError("There is no SDigit Tree");
455      return 0;
456    }
457   
458   TBranch * branch = treeS->GetBranch(fDetectorName);
459   if (branch == 0) 
460    {//easy, maybe just a new tree
461     //AliError("Cannot find branch PHOS"); 
462     return 0;
463   }
464     
465   branch->SetAddress(SDigitsRef());
466   branch->GetEntry(0);
467   return 0;
468 }
469
470 //____________________________________________________________________________ 
471 Int_t AliPHOSLoader::ReadDigits()
472 {
473 // Read the summable digits tree TreeS():
474 // Check if TClones is in folder
475 // if not create and add to folder
476 // connect to tree if available
477 // Read the data
478   
479   TObject** dref = DigitsRef();
480   if(dref == 0x0)
481    {//if there is not array in folder, create it and put it there
482      MakeDigitsArray();
483      dref = DigitsRef();
484    }
485
486   TTree * treeD = TreeD();
487   if(treeD==0)
488    {
489      //May happen if file is truncated or new in LoadSDigits
490      //AliError("There is no Digit Tree");
491      return 0;
492    }
493
494   TBranch * branch = treeD->GetBranch(fDetectorName);
495   if (branch == 0) 
496    {//easy, maybe just a new tree
497     //AliError("Cannot find branch ",fDetectorName.Data()); 
498     return 0;
499    }
500   
501   branch->SetAddress(dref);//connect branch to buffer sitting in folder
502   branch->GetEntry(0);//get first event 
503
504   return 0;  
505 }
506
507 //____________________________________________________________________________ 
508
509 void AliPHOSLoader::Track(Int_t itrack)
510 {
511 // Read the first entry of PHOS branch in hit tree gAlice->TreeH()
512   if(TreeH()== 0)
513    {
514      if (LoadHits())
515       {
516         AliError("Can not load hits.");
517         return;
518       } 
519    }
520   
521   TBranch * hitsbranch = dynamic_cast<TBranch*>(TreeH()->GetListOfBranches()->FindObject("PHOS")) ;
522   if ( !hitsbranch ) {
523     if (fDebug)
524       cout << "WARNING:  AliPHOSLoader::ReadTreeH -> Cannot find branch PHOS" << endl ; 
525     return ;
526   }  
527   if(!Hits()) PostHits();
528
529   hitsbranch->SetAddress(HitsRef());
530   hitsbranch->GetEntry(itrack);
531
532 }
533 //____________________________________________________________________________ 
534 void AliPHOSLoader::ReadTreeQA()
535 {
536   // Read the digit tree gAlice->TreeQA()
537   // so far only PHOS knows about this Tree  
538
539   if(PHOS()->TreeQA()== 0){
540     cerr <<   "ERROR: AliPHOSLoader::ReadTreeQA: can not read TreeQA " << endl ;
541     return ;
542   }
543   
544   TBranch * qabranch = PHOS()->TreeQA()->GetBranch("PHOS");
545   if (!qabranch) { 
546     if (fDebug)
547       cout << "WARNING: AliPHOSLoader::ReadTreeQA -> Cannot find QA Alarms for PHOS" << endl ;
548     return ; 
549   }   
550   
551 //  if(!Alarms()) PostQA();
552
553   qabranch->SetAddress(AlarmsRef()) ;
554
555   qabranch->GetEntry(0) ;
556   
557 }
558
559
560 //____________________________________________________________________________ 
561 Int_t AliPHOSLoader::ReadRecPoints()
562 {
563 //Creates and posts to folder an array container, 
564 //connects branch in tree (if exists), and reads data to array
565   
566   MakeRecPointsArray();
567   
568   TObjArray * cpva = 0x0 ; 
569   TObjArray * emca = 0x0 ; 
570   
571   TTree * treeR = TreeR();
572  
573   if(treeR==0)
574    {
575      //May happen if file is truncated or new in LoadSDigits
576      return 0;
577    }
578
579   Int_t retval = 0;
580   TBranch * emcbranch = treeR->GetBranch(fgkEmcRecPointsBranchName);
581
582   if (emcbranch == 0x0)
583    {
584      AliError(Form("Can not get branch with EMC Rec. Points named %s",
585                    fgkEmcRecPointsBranchName.Data()));
586      retval = 1;
587    }
588   else
589    {
590      emcbranch->SetAddress(&emca) ;
591      emcbranch->GetEntry(0) ;
592    }
593   TBranch * cpvbranch = treeR->GetBranch(fgkCpvRecPointsBranchName);
594   if (cpvbranch == 0x0)
595    {
596      AliError(Form("Can not get branch with CPV Rec. Points named %s",
597                    fgkCpvRecPointsBranchName.Data()));
598      retval = 2;
599    }
600   else
601    {
602      cpvbranch->SetAddress(&cpva);
603      cpvbranch->GetEntry(0) ;
604    }
605
606   Int_t ii ; 
607   Int_t maxemc = emca->GetEntries() ; 
608   for ( ii= 0 ; ii < maxemc ; ii++ ) 
609     EmcRecPoints()->Add(emca->At(ii)) ;
610  
611   Int_t maxcpv = cpva->GetEntries() ;
612   for ( ii= 0 ; ii < maxcpv ; ii++ )
613     CpvRecPoints()->Add(cpva->At(ii)) ; 
614
615   return retval;
616 }
617
618 //____________________________________________________________________________ 
619 Int_t AliPHOSLoader::ReadTracks()
620 {
621 //Creates and posts to folder an array container, 
622 //connects branch in tree (if exists), and reads data to arry
623
624   TObject** trkref = TracksRef();
625   if ( trkref == 0x0 )   
626    {//Create and post array
627     MakeTrackSegmentsArray();
628     trkref = TracksRef();
629    }
630
631   TTree * treeT = TreeT();
632   if(treeT==0)
633    {
634      //May happen if file is truncated or new in LoadSDigits, or the file is in update mode, 
635      //but tracking was not performed yet for a current event
636      //AliError("There is no Tree with Tracks");
637      return 0;
638    }
639   
640   TBranch * branch = treeT->GetBranch(fgkTrackSegmentsBranchName);
641   if (branch == 0) 
642    {//easy, maybe just a new tree
643     AliError(Form("Cannot find branch named %s",
644                   fgkTrackSegmentsBranchName.Data()));
645     return 0;
646   }
647
648   branch->SetAddress(trkref);//connect branch to buffer sitting in folder
649   branch->GetEntry(0);//get first event 
650   
651   return 0;
652 }
653 //____________________________________________________________________________ 
654
655 Int_t AliPHOSLoader::ReadRecParticles()
656 {
657 //Reads Reconstructed  Particles from file
658 //Creates and posts to folder an array container, 
659 //connects branch in tree (if exists), and reads data to arry
660   
661   TObject** recpartref = RecParticlesRef();
662   
663   if ( recpartref == 0x0 )   
664    {//Create and post array
665     MakeRecParticlesArray();
666     recpartref = RecParticlesRef();
667    }
668
669   TTree * treeP = TreeP();
670   if(treeP==0)
671    {
672      //May happen if file is truncated or new in LoadSDigits, 
673      //or the file is in update mode, 
674      //but tracking was not performed yet for a current event
675      //     AliError("There is no Tree with Tracks and Reconstructed Particles");
676      return 0;
677    }
678   
679   TBranch * branch = treeP->GetBranch(fgkRecParticlesBranchName);
680   if (branch == 0) 
681    {//easy, maybe just a new tree
682     AliError(Form("Cannot find branch %s",
683                   fgkRecParticlesBranchName.Data())); 
684     return 0;
685   }
686
687   branch->SetAddress(recpartref);//connect branch to buffer sitting in folder
688   branch->GetEntry(0);//get first event 
689   
690   return 0;
691 }
692
693
694 AliPHOSGeometry* AliPHOSLoader::GetPHOSGeometry()
695 {
696   //returns PHOS geometry from gAlice 
697   //static Method used by some classes where it is not convienient to pass eventfoldername
698   if (gAlice == 0x0)
699     return 0x0;
700   AliPHOS* phos=dynamic_cast<AliPHOS*>(gAlice->GetDetector("PHOS"));
701   if (phos == 0x0)
702     return 0x0;
703   return phos->GetGeometry();
704 }
705 /***************************************************************************************/
706
707 AliPHOSLoader* AliPHOSLoader::GetPHOSLoader(const  char* eventfoldername)
708 {
709   // Return PHOS loader
710   AliRunLoader* rn  = AliRunLoader::GetRunLoader(eventfoldername);
711   if (rn == 0x0) {
712     printf("Can not find Run Loader in folder %s", eventfoldername);
713     return 0x0 ; 
714   }
715   return dynamic_cast<AliPHOSLoader*>(rn->GetLoader("PHOSLoader"));
716 }
717 /***************************************************************************************/
718
719 Bool_t AliPHOSLoader::BranchExists(const TString& recName)
720 {
721   // Check if a branch named redName exists
722   if (fBranchTitle.IsNull()) return kFALSE;
723   TString dataname, zername ;
724   TTree* tree;
725   if(recName == "SDigits")
726    {
727     tree = TreeS();
728     dataname = GetDetectorName();
729     zername = "AliPHOSSDigitizer" ;
730    }
731   else
732     if(recName == "Digits"){
733       tree = TreeD();
734       dataname = GetDetectorName();
735       zername = "AliPHOSDigitizer" ;
736     }
737     else
738       if(recName == "RecPoints"){
739        tree = TreeR();
740        dataname = fgkEmcRecPointsBranchName;
741        zername = "AliPHOSClusterizer" ;
742       }
743       else
744        if(recName == "TrackSegments"){
745          tree = TreeT();
746          dataname = fgkTrackSegmentsBranchName;
747          zername = "AliPHOSTrackSegmentMaker";
748        }        
749        else
750          if(recName == "RecParticles"){
751            tree = TreeP();
752            dataname = fgkRecParticlesBranchName;
753            zername = "AliPHOSPID";
754          }
755          else
756            return kFALSE ;
757
758   
759   if(!tree ) 
760     return kFALSE ;
761
762   TObjArray * lob = static_cast<TObjArray*>(tree->GetListOfBranches()) ;
763   TIter next(lob) ; 
764   TBranch * branch = 0 ;  
765   TString titleName(fBranchTitle);
766   titleName+=":";
767
768   while ((branch = (static_cast<TBranch*>(next())))) {
769     TString branchName(branch->GetName() ) ; 
770     TString branchTitle(branch->GetTitle() ) ;  
771     if ( branchName.BeginsWith(dataname) && branchTitle.BeginsWith(fBranchTitle) ){  
772       AliWarning(Form("branch %s  with title  %s ",
773                       dataname.Data(),fBranchTitle.Data()));
774       return kTRUE ;
775     }
776     if ( branchName.BeginsWith(zername) &&  branchTitle.BeginsWith(titleName) ){
777       AliWarning(Form("branch AliPHOS... with title  %s ",
778                       branch->GetTitle()));
779       return kTRUE ; 
780     }
781   }
782   return kFALSE ;
783
784  }
785
786 void AliPHOSLoader::SetBranchTitle(const TString& btitle)
787 {
788   // Set branch title
789   if (btitle.CompareTo(fBranchTitle) == 0) return;
790   fBranchTitle = btitle;
791   ReloadAll();
792  }
793 //____________________________________________________________________________ 
794
795 void AliPHOSLoader::CleanHits()
796 {
797   // Clean Hits array
798   AliLoader::CleanHits();
799   //Clear an array 
800   TClonesArray* hits = Hits();
801   if (hits) hits->Clear();
802 }
803 //____________________________________________________________________________ 
804
805 void AliPHOSLoader::CleanSDigits()
806 {
807   // Clean SDigits array
808   AliLoader::CleanSDigits();
809   TClonesArray* sdigits = SDigits();
810   if (sdigits) sdigits->Clear();
811   
812 }
813 //____________________________________________________________________________ 
814
815 void AliPHOSLoader::CleanDigits()
816 {
817   // Clean Digits array
818   AliLoader::CleanDigits();
819   TClonesArray* digits = Digits();
820   if (digits) digits->Clear();
821 }
822 //____________________________________________________________________________ 
823
824 void AliPHOSLoader::CleanRecPoints()
825 {
826   // Clean RecPoints array
827   AliLoader::CleanRecPoints();
828   TObjArray* recpoints = EmcRecPoints();
829   if (recpoints) recpoints->Clear();
830   recpoints = CpvRecPoints();
831   if (recpoints) recpoints->Clear();
832 }
833 //____________________________________________________________________________ 
834
835 void AliPHOSLoader::CleanTracks()
836 {
837   //Cleans Tracks stuff
838   AliLoader::CleanTracks();//tree
839   
840   //and clear the array
841   TClonesArray* tracks = TrackSegments();
842   if (tracks) tracks->Clear();
843
844 }
845 //____________________________________________________________________________ 
846
847 void AliPHOSLoader::CleanRecParticles()
848  {
849    // Clean RecParticles array
850    TClonesArray *recpar = RecParticles();
851    if (recpar) recpar->Clear();
852   
853  
854  }
855 //____________________________________________________________________________ 
856
857 void AliPHOSLoader::ReadCalibrationDB(const char * database,const char * filename)
858 {
859   // Read calibration data base from file
860   if(fcdb && (strcmp(database,fcdb->GetTitle())==0))
861     return ;
862
863   TFile * file = gROOT->GetFile(filename) ;
864   if(!file)
865     file = TFile::Open(filename);
866   if(!file){
867     AliError(Form("Cannot open file %s", filename)) ;
868     return ;
869   }
870   if(fcdb)
871     fcdb->Delete() ;
872   fcdb = dynamic_cast<AliPHOSCalibrationDB *>(file->Get("AliPHOSCalibrationDB")) ;
873   if(!fcdb)
874     AliError(Form("No database %s in file %s", database, filename)) ;
875 }
876 //____________________________________________________________________________ 
877
878 // AliPHOSSDigitizer*  AliPHOSLoader::PHOSSDigitizer() 
879 // { 
880 // //return PHOS SDigitizer
881 //  return  dynamic_cast<AliPHOSSDigitizer*>(SDigitizer()) ;
882 // }
883
884 //____________________________________________________________________________ 
885 void AliPHOSLoader::MakeHitsArray()
886 {
887   // Add Hits array to the data folder
888   if (Hits()) return;
889   TClonesArray* hits = new TClonesArray("AliPHOSHit",1000);
890   hits->SetName(fgkHitsName);
891   GetDetectorDataFolder()->Add(hits);
892 }
893
894 //____________________________________________________________________________ 
895 void AliPHOSLoader::MakeSDigitsArray()
896 {
897   // Add SDigits array to the data folder
898   if ( SDigits()) return;
899   TClonesArray* sdigits = new TClonesArray("AliPHOSDigit",1);
900   sdigits->SetName(fgkSDigitsName);
901   GetDetectorDataFolder()->Add(sdigits);
902 }
903
904 //____________________________________________________________________________ 
905 void AliPHOSLoader::MakeDigitsArray()
906 {
907   // Add Digits array to the data folder
908   if ( Digits()) return;
909   TClonesArray* digits = new TClonesArray("AliPHOSDigit",1);
910   digits->SetName(fgkDigitsName);
911   GetDetectorDataFolder()->Add(digits);
912   
913 }
914
915 //____________________________________________________________________________ 
916 void AliPHOSLoader::MakeRecPointsArray()
917 {
918   // Add RecPoints array to the data folder
919   if ( EmcRecPoints() == 0x0)
920    {
921      AliDebug(9, "Making array for EMC");
922     TObjArray* emc = new TObjArray(100) ;
923     emc->SetName(fgkEmcRecPointsName) ;
924     GetDetectorDataFolder()->Add(emc);
925    }
926
927   if ( CpvRecPoints() == 0x0)
928    { 
929      AliDebug(9, "Making array for CPV");
930      TObjArray* cpv = new TObjArray(100) ;
931      cpv->SetName(fgkCpvRecPointsName);
932      GetDetectorDataFolder()->Add(cpv);
933    }
934 }
935
936 //____________________________________________________________________________ 
937 void AliPHOSLoader::MakeTrackSegmentsArray()
938 {
939   // Add TrackSegments array to the data folder
940   if ( TrackSegments()) return;
941   TClonesArray * ts = new TClonesArray("AliPHOSTrackSegment",100) ;
942   ts->SetName(fgkTracksName);
943   GetDetectorDataFolder()->Add(ts);
944
945 }
946
947 //____________________________________________________________________________ 
948 void AliPHOSLoader::MakeRecParticlesArray()
949 {
950   // Add RecParticles array to the data folder
951   if ( RecParticles()) return;
952   TClonesArray * rp = new TClonesArray("AliPHOSRecParticle",100) ;
953   rp->SetName(fgkRecParticlesName);
954   GetDetectorDataFolder()->Add(rp);
955 }