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