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