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