]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSLoader.cxx
Forgotten: in reconstruct to get branchname from RunLoader. Clean of the task is...
[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   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      Error("SetEvent","AliLoader::SetEvent 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      Error("GetEvent","AliLoader::GetEvent 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      Error("LoadHits","AliLoader::LoadHits 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      Error("PostSDigits","AliLoader::LoadSDigits 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      Error("LoadDigits","AliLoader::LoadDigits 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      Error("LoadRecPoints","AliLoader::LoadRecPoints returned error");
247      return res;
248    }
249
250   TFolder * phosFolder = GetDetectorDataFolder();
251   if ( phosFolder  == 0x0 ) 
252    {
253      Error("PostDigits","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  if (GetDebug()) Info("LoadTracks","opt = %s",opt);
264  Int_t res;
265  res = AliLoader::LoadTracks(opt);
266  if (res)
267    {//oops, error
268       Error("LoadTracks","AliLoader::LoadTracks 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      Error("LoadRecParticles","AliLoader::LoadRecParticles 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      Error("PostHits","AliLoader::  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      Error("PostSDigits","AliLoader::PostSDigits  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      Error("PostDigits","AliLoader::PostDigits  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      Error("PostRecPoints","AliLoader::PostRecPoints  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      Error("PostRecParticles","AliLoader::PostRecParticles  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      Error("PostTracks","AliLoader::PostTracks  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     Error("ReadHits"," Cannot read TreeH from folder");
394     return 1;
395   }
396   
397   TBranch * hitsbranch = treeh->GetBranch(fDetectorName);
398   if (hitsbranch == 0) 
399    {
400     Error("ReadHits"," Cannot find branch PHOS"); 
401     return 1;
402   }
403   
404   if (GetDebug()) Info("ReadHits","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      //Error("ReadSDigits","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     //Error("ReadSDigits"," 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      //Error("ReadDigits","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     //Error("ReadDigits"," 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         Error("Track","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      Error("ReadRecPoints","Can not get branch with EMC Rec. Points named %s",fgkEmcRecPointsBranchName.Data());
585      retval = 1;
586    }
587   else
588    {
589      emcbranch->SetAddress(&emca) ;
590      emcbranch->GetEntry(0) ;
591    }
592   TBranch * cpvbranch = treeR->GetBranch(fgkCpvRecPointsBranchName);
593   if (cpvbranch == 0x0)
594    {
595      Error("ReadRecPoints","Can not get branch with CPV Rec. Points named %s",fgkCpvRecPointsBranchName.Data());
596      retval = 2;
597    }
598   else
599    {
600      cpvbranch->SetAddress(&cpva);
601      cpvbranch->GetEntry(0) ;
602    }
603
604   Int_t ii ; 
605   Int_t maxemc = emca->GetEntries() ; 
606   for ( ii= 0 ; ii < maxemc ; ii++ ) 
607     EmcRecPoints()->Add(emca->At(ii)) ;
608  
609   Int_t maxcpv = cpva->GetEntries() ;
610   for ( ii= 0 ; ii < maxcpv ; ii++ )
611     CpvRecPoints()->Add(cpva->At(ii)) ; 
612
613   return retval;
614 }
615
616 //____________________________________________________________________________ 
617 Int_t AliPHOSLoader::ReadTracks()
618 {
619 //Creates and posts to folder an array container, 
620 //connects branch in tree (if exists), and reads data to arry
621
622   TObject** trkref = TracksRef();
623   if ( trkref == 0x0 )   
624    {//Create and post array
625     MakeTrackSegmentsArray();
626     trkref = TracksRef();
627    }
628
629   TTree * treeT = TreeT();
630   if(treeT==0)
631    {
632      //May happen if file is truncated or new in LoadSDigits, or the file is in update mode, 
633      //but tracking was not performed yet for a current event
634      //Error("ReadTracks","There is no Tree with Tracks");
635      return 0;
636    }
637   
638   TBranch * branch = treeT->GetBranch(fgkTrackSegmentsBranchName);
639   if (branch == 0) 
640    {//easy, maybe just a new tree
641     Error("ReadTracks"," Cannot find branch named %s",fgkTrackSegmentsBranchName.Data());
642     return 0;
643   }
644
645   branch->SetAddress(trkref);//connect branch to buffer sitting in folder
646   branch->GetEntry(0);//get first event 
647   
648   return 0;
649 }
650 //____________________________________________________________________________ 
651
652 Int_t AliPHOSLoader::ReadRecParticles()
653 {
654 //Reads Reconstructed  Particles from file
655 //Creates and posts to folder an array container, 
656 //connects branch in tree (if exists), and reads data to arry
657   
658   TObject** recpartref = RecParticlesRef();
659   
660   if ( recpartref == 0x0 )   
661    {//Create and post array
662     MakeRecParticlesArray();
663     recpartref = RecParticlesRef();
664    }
665
666   TTree * treeP = TreeP();
667   if(treeP==0)
668    {
669      //May happen if file is truncated or new in LoadSDigits, 
670      //or the file is in update mode, 
671      //but tracking was not performed yet for a current event
672      //     Error("ReadRecParticles","There is no Tree with Tracks and Reconstructed Particles");
673      return 0;
674    }
675   
676   TBranch * branch = treeP->GetBranch(fgkRecParticlesBranchName);
677   if (branch == 0) 
678    {//easy, maybe just a new tree
679     Error("ReadRecParticles"," Cannot find branch %s",fgkRecParticlesBranchName.Data()); 
680     return 0;
681   }
682
683   branch->SetAddress(recpartref);//connect branch to buffer sitting in folder
684   branch->GetEntry(0);//get first event 
685   
686   return 0;
687 }
688
689
690 AliPHOSGeometry* AliPHOSLoader::GetPHOSGeometry()
691 {
692   //returns PHOS geometry from gAlice 
693   //static Method used by some classes where it is not convienient to pass eventfoldername
694   if (gAlice == 0x0)
695     return 0x0;
696   AliPHOS* phos=dynamic_cast<AliPHOS*>(gAlice->GetDetector("PHOS"));
697   if (phos == 0x0)
698     return 0x0;
699   return phos->GetGeometry();
700 }
701 /***************************************************************************************/
702
703 AliPHOSLoader* AliPHOSLoader::GetPHOSLoader(const  char* eventfoldername)
704 {
705   // Return PHOS loader
706   AliRunLoader* rn  = AliRunLoader::GetRunLoader(eventfoldername);
707   if (rn == 0x0)
708    {
709      cerr<<"Error: <AliPHOSLoader::GetPHOSLoader>: "
710          << "Can not find Run Loader in folder "<<eventfoldername<<endl;
711      return 0x0;
712    }
713   return dynamic_cast<AliPHOSLoader*>(rn->GetLoader("PHOSLoader"));
714 }
715 /***************************************************************************************/
716
717 Bool_t AliPHOSLoader::BranchExists(const TString& recName)
718 {
719   // Check if a branch named redName exists
720   if (fBranchTitle.IsNull()) return kFALSE;
721   TString dataname, zername ;
722   TTree* tree;
723   if(recName == "SDigits")
724    {
725     tree = TreeS();
726     dataname = GetDetectorName();
727     zername = "AliPHOSSDigitizer" ;
728    }
729   else
730     if(recName == "Digits"){
731       tree = TreeD();
732       dataname = GetDetectorName();
733       zername = "AliPHOSDigitizer" ;
734     }
735     else
736       if(recName == "RecPoints"){
737        tree = TreeR();
738        dataname = fgkEmcRecPointsBranchName;
739        zername = "AliPHOSClusterizer" ;
740       }
741       else
742        if(recName == "TrackSegments"){
743          tree = TreeT();
744          dataname = fgkTrackSegmentsBranchName;
745          zername = "AliPHOSTrackSegmentMaker";
746        }        
747        else
748          if(recName == "RecParticles"){
749            tree = TreeP();
750            dataname = fgkRecParticlesBranchName;
751            zername = "AliPHOSPID";
752          }
753          else
754            return kFALSE ;
755
756   
757   if(!tree ) 
758     return kFALSE ;
759
760   TObjArray * lob = static_cast<TObjArray*>(tree->GetListOfBranches()) ;
761   TIter next(lob) ; 
762   TBranch * branch = 0 ;  
763   TString titleName(fBranchTitle);
764   titleName+=":";
765
766   while ((branch = (static_cast<TBranch*>(next())))) {
767     TString branchName(branch->GetName() ) ; 
768     TString branchTitle(branch->GetTitle() ) ;  
769     if ( branchName.BeginsWith(dataname) && branchTitle.BeginsWith(fBranchTitle) ){  
770       Warning("BranchExists","branch %s  with title  %s ",dataname.Data(),fBranchTitle.Data());
771       return kTRUE ;
772     }
773     if ( branchName.BeginsWith(zername) &&  branchTitle.BeginsWith(titleName) ){
774       Warning("BranchExists","branch AliPHOS... with title  %s ",branch->GetTitle());
775       return kTRUE ; 
776     }
777   }
778   return kFALSE ;
779
780  }
781
782 void AliPHOSLoader::SetBranchTitle(const TString& btitle)
783 {
784   // Set branch title
785   if (btitle.CompareTo(fBranchTitle) == 0) return;
786   fBranchTitle = btitle;
787   ReloadAll();
788  }
789 //____________________________________________________________________________ 
790
791 void AliPHOSLoader::CleanHits()
792 {
793   // Clean Hits array
794   AliLoader::CleanHits();
795   //Clear an array 
796   TClonesArray* hits = Hits();
797   if (hits) hits->Clear();
798 }
799 //____________________________________________________________________________ 
800
801 void AliPHOSLoader::CleanSDigits()
802 {
803   // Clean SDigits array
804   AliLoader::CleanSDigits();
805   TClonesArray* sdigits = SDigits();
806   if (sdigits) sdigits->Clear();
807   
808 }
809 //____________________________________________________________________________ 
810
811 void AliPHOSLoader::CleanDigits()
812 {
813   // Clean Digits array
814   AliLoader::CleanDigits();
815   TClonesArray* digits = Digits();
816   if (digits) digits->Clear();
817 }
818 //____________________________________________________________________________ 
819
820 void AliPHOSLoader::CleanRecPoints()
821 {
822   // Clean RecPoints array
823   AliLoader::CleanRecPoints();
824   TObjArray* recpoints = EmcRecPoints();
825   if (recpoints) recpoints->Clear();
826   recpoints = CpvRecPoints();
827   if (recpoints) recpoints->Clear();
828 }
829 //____________________________________________________________________________ 
830
831 void AliPHOSLoader::CleanTracks()
832 {
833   //Cleans Tracks stuff
834   AliLoader::CleanTracks();//tree
835   
836   //and clear the array
837   TClonesArray* tracks = TrackSegments();
838   if (tracks) tracks->Clear();
839
840 }
841 //____________________________________________________________________________ 
842
843 void AliPHOSLoader::CleanRecParticles()
844  {
845    // Clean RecParticles array
846    TClonesArray *recpar = RecParticles();
847    if (recpar) recpar->Clear();
848   
849  
850  }
851 //____________________________________________________________________________ 
852
853 void AliPHOSLoader::ReadCalibrationDB(const char * database,const char * filename)
854 {
855   // Read calibration data base from file
856   if(fcdb && (strcmp(database,fcdb->GetTitle())==0))
857     return ;
858
859   TFile * file = gROOT->GetFile(filename) ;
860   if(!file)
861     file = TFile::Open(filename);
862   if(!file){
863     Error ("ReadCalibrationDB", "Cannot open file %s", filename) ;
864     return ;
865   }
866   if(fcdb)
867     fcdb->Delete() ;
868   fcdb = dynamic_cast<AliPHOSCalibrationDB *>(file->Get("AliPHOSCalibrationDB")) ;
869   if(!fcdb)
870     Error ("ReadCalibrationDB", "No database %s in file %s", database, filename) ;
871 }
872 //____________________________________________________________________________ 
873
874 // AliPHOSSDigitizer*  AliPHOSLoader::PHOSSDigitizer() 
875 // { 
876 // //return PHOS SDigitizer
877 //  return  dynamic_cast<AliPHOSSDigitizer*>(SDigitizer()) ;
878 // }
879
880 //____________________________________________________________________________ 
881 void AliPHOSLoader::MakeHitsArray()
882 {
883   // Add Hits array to the data folder
884   if (Hits()) return;
885   TClonesArray* hits = new TClonesArray("AliPHOSHit",1000);
886   hits->SetName(fgkHitsName);
887   GetDetectorDataFolder()->Add(hits);
888 }
889
890 //____________________________________________________________________________ 
891 void AliPHOSLoader::MakeSDigitsArray()
892 {
893   // Add SDigits array to the data folder
894   if ( SDigits()) return;
895   TClonesArray* sdigits = new TClonesArray("AliPHOSDigit",1);
896   sdigits->SetName(fgkSDigitsName);
897   GetDetectorDataFolder()->Add(sdigits);
898 }
899
900 //____________________________________________________________________________ 
901 void AliPHOSLoader::MakeDigitsArray()
902 {
903   // Add Digits array to the data folder
904   if ( Digits()) return;
905   TClonesArray* digits = new TClonesArray("AliPHOSDigit",1);
906   digits->SetName(fgkDigitsName);
907   GetDetectorDataFolder()->Add(digits);
908   
909 }
910
911 //____________________________________________________________________________ 
912 void AliPHOSLoader::MakeRecPointsArray()
913 {
914   // Add RecPoints array to the data folder
915   if ( EmcRecPoints() == 0x0)
916    {
917     if (GetDebug()>9) Info("MakeRecPointsArray","Making array for EMC");
918     TObjArray* emc = new TObjArray(100) ;
919     emc->SetName(fgkEmcRecPointsName) ;
920     GetDetectorDataFolder()->Add(emc);
921    }
922
923   if ( CpvRecPoints() == 0x0)
924    {
925     if (GetDebug()>9) Info("MakeRecPointsArray","Making array for CPV");
926     TObjArray* cpv = new TObjArray(100) ;
927     cpv->SetName(fgkCpvRecPointsName);
928     GetDetectorDataFolder()->Add(cpv);
929    }
930 }
931
932 //____________________________________________________________________________ 
933 void AliPHOSLoader::MakeTrackSegmentsArray()
934 {
935   // Add TrackSegments array to the data folder
936   if ( TrackSegments()) return;
937   TClonesArray * ts = new TClonesArray("AliPHOSTrackSegment",100) ;
938   ts->SetName(fgkTracksName);
939   GetDetectorDataFolder()->Add(ts);
940
941 }
942
943 //____________________________________________________________________________ 
944 void AliPHOSLoader::MakeRecParticlesArray()
945 {
946   // Add RecParticles array to the data folder
947   if ( RecParticles()) return;
948   TClonesArray * rp = new TClonesArray("AliPHOSRecParticle",100) ;
949   rp->SetName(fgkRecParticlesName);
950   GetDetectorDataFolder()->Add(rp);
951 }