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