1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 /* History of cvs commits:
21 * Revision 1.18 2006/08/28 10:01:56 kharlov
22 * Effective C++ warnings fixed (Timur Pocheptsov)
24 * Revision 1.17 2006/08/25 16:00:53 kharlov
25 * Compliance with Effective C++AliPHOSHit.cxx
27 * Revision 1.16 2006/08/01 12:15:04 cvetan
28 * Adding a constructor from TFolder. Needed by AliReconstruction plugin scheme
30 * Revision 1.15 2005/07/12 20:07:35 hristov
31 * Changes needed to run simulation and reconstrruction in the same AliRoot session
33 * Revision 1.14 2005/05/28 14:19:04 schutz
34 * Compilation warnings fixed by T.P.
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.
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) ;
51 // please->GetEvent(event) ; // reads new event from galice.root
53 //-- Author: Yves Schutz (SUBATECH) & Dmitri Peressounko (RRC KI & SUBATECH)
54 //-- Completely redesigned by Dmitri Peressounko March 2001
56 //-- YS June 2001 : renamed the original AliPHOSIndexToObject and make
57 //-- systematic usage of TFolders without changing the interface
58 //////////////////////////////////////////////////////////////////////////////
61 // --- ROOT system ---
67 // --- Standard library ---
69 // --- AliRoot header files ---
70 #include "AliObjectLoader.h"
72 #include "AliPHOSLoader.h"
74 #include "AliPHOSHit.h"
75 #include "AliPHOSCalibrationDB.h"
77 ClassImp(AliPHOSLoader)
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
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)
97 //____________________________________________________________________________
98 AliPHOSLoader::AliPHOSLoader(const Char_t *detname,const Char_t *eventfoldername) :
99 AliLoader(detname, eventfoldername),
100 fBranchTitle(), fcdb(0), fDebug(0)
104 //____________________________________________________________________________
105 AliPHOSLoader::AliPHOSLoader(const Char_t *detname,TFolder *topfolder):
106 AliLoader(detname,topfolder),
107 fBranchTitle(), fcdb(0), fDebug(0)
112 //____________________________________________________________________________
113 AliPHOSLoader::AliPHOSLoader(const AliPHOSLoader & obj):
114 AliLoader(obj),fBranchTitle(obj.GetBranchTitle()),fcdb(obj.CalibrationDB()),
115 fDebug(obj.GetDebug())
119 //____________________________________________________________________________
121 AliPHOSLoader::~AliPHOSLoader()
123 //remove and delete arrays
125 Clean(fgkSDigitsName);
126 Clean(fgkDigitsName);
127 Clean(fgkEmcRecPointsName);
128 Clean(fgkCpvRecPointsName);
129 Clean(fgkTracksName);
130 Clean(fgkRecParticlesName);
134 //____________________________________________________________________________
135 void AliPHOSLoader::CleanFolders()
138 AliLoader::CleanFolders();
141 //____________________________________________________________________________
142 Int_t AliPHOSLoader::SetEvent()
144 //Cleans loaded stuff and and sets Files and Directories
145 // do not post any data to folder/tasks
148 Int_t retval = AliLoader::SetEvent();
151 AliError("returned error");
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();
167 //____________________________________________________________________________
168 Int_t AliPHOSLoader::GetEvent()
170 //Overloads GetEvent method called by AliRunLoader::GetEvent(Int_t) method
171 //to add Rec Particles specific for PHOS
173 //First call the original method to get whatever from std. setup is needed
176 retval = AliLoader::GetEvent();
179 AliError("returned error");
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();
191 //Now, check if RecPart were loaded
195 // //____________________________________________________________________________
196 // const AliPHOS * AliPHOSLoader::PHOS()
198 // // returns the PHOS object
199 // AliPHOS * phos = dynamic_cast<AliPHOS*>(GetModulesFolder()->FindObject(fDetectorName));
202 // cout << "WARNING: AliPHOSLoader::PHOS -> PHOS module not found in Folders" << endl ;
206 //____________________________________________________________________________
207 Int_t AliPHOSLoader::LoadHits(Option_t* opt)
209 //------- Hits ----------------------
210 //Overload (extends) LoadHits implemented in AliLoader
214 //First call the AliLoader's method to send the TreeH to folder
215 res = AliLoader::LoadHits(opt);
219 AliError("returned error");
223 //read the data from tree in folder and send it to folder
229 //____________________________________________________________________________
230 Int_t AliPHOSLoader::LoadSDigits(Option_t* opt)
232 //---------- SDigits -------------------------
234 //First call the AliLoader's method to send the TreeS to folder
235 res = AliLoader::LoadSDigits(opt);
238 AliError("returned error");
241 return ReadSDigits();
244 //____________________________________________________________________________
245 Int_t AliPHOSLoader::LoadDigits(Option_t* opt)
247 //---------- Digits -------------------------
249 //First call the AliLoader's method to send the TreeS to folder
250 res = AliLoader::LoadDigits(opt);
253 AliError("returned error");
258 //____________________________________________________________________________
259 Int_t AliPHOSLoader::LoadRecPoints(Option_t* opt)
261 // -------------- RecPoints -------------------------------------------
263 //First call the AliLoader's method to send the TreeR to folder
264 res = AliLoader::LoadRecPoints(opt);
267 AliError("returned error");
271 TFolder * phosFolder = GetDetectorDataFolder();
272 if ( phosFolder == 0x0 )
274 AliError("Can not get detector data folder");
277 return ReadRecPoints();
279 //____________________________________________________________________________
281 Int_t AliPHOSLoader::LoadTracks(Option_t* opt)
283 //Loads Tracks: Open File, Reads Tree and posts, Read Data and Posts
284 AliDebug(1, Form("opt = %s",opt));
286 res = AliLoader::LoadTracks(opt);
289 AliError("returned error");
295 //____________________________________________________________________________
296 Int_t AliPHOSLoader::LoadRecParticles(Option_t* opt)
298 // -------------- RecPoints -------------------------------------------
300 //First call the AliLoader's method to send the TreeT to folder
301 res = AliLoader::LoadRecParticles(opt);
304 AliError("returned error");
307 return ReadRecParticles();
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
316 // -------------- Hits -------------------------------------------
317 Int_t reval = AliLoader::PostHits();
320 AliError("returned error");
324 return const_cast<AliPHOSLoader *>(this)->ReadHits();
326 //____________________________________________________________________________
328 Int_t AliPHOSLoader::PostSDigits()const
330 // -------------- SDigits -------------------------------------------
331 Int_t reval = AliLoader::PostSDigits();
334 AliError("returned error");
337 return const_cast<AliPHOSLoader *>(this)->ReadSDigits();
339 //____________________________________________________________________________
341 Int_t AliPHOSLoader::PostDigits()const
343 // -------------- Digits -------------------------------------------
344 Int_t reval = AliLoader::PostDigits();
347 AliError("returned error");
350 return const_cast<AliPHOSLoader *>(this)->ReadDigits();
352 //____________________________________________________________________________
354 Int_t AliPHOSLoader::PostRecPoints()const
356 // -------------- RecPoints -------------------------------------------
357 Int_t reval = AliLoader::PostRecPoints();
360 AliError("returned error");
363 return const_cast<AliPHOSLoader *>(this)->ReadRecPoints();
366 //____________________________________________________________________________
368 Int_t AliPHOSLoader::PostRecParticles()const
370 // -------------- RecParticles -------------------------------------------
371 Int_t reval = AliLoader::PostRecParticles();
374 AliError("returned error");
377 return const_cast<AliPHOSLoader *>(this)->ReadRecParticles();
379 //____________________________________________________________________________
381 Int_t AliPHOSLoader::PostTracks()const
383 // -------------- Tracks -------------------------------------------
384 Int_t reval = AliLoader::PostTracks();
387 AliError("returned error");
390 return const_cast<AliPHOSLoader *>(this)->ReadTracks();
392 //____________________________________________________________________________
396 //____________________________________________________________________________
397 Int_t AliPHOSLoader::ReadHits()
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
404 TObject** hitref = HitsRef();
411 TClonesArray* hits = dynamic_cast<TClonesArray*>(*hitref);
413 TTree* treeh = TreeH();
417 AliError("Cannot read TreeH from folder");
421 TBranch * hitsbranch = treeh->GetBranch(fDetectorName);
424 AliError("Cannot find branch PHOS");
428 AliDebug(1, "Reading Hits");
430 if (hitsbranch->GetEntries() > 1)
432 TClonesArray * tempo = new TClonesArray("AliPHOSHit",1000);
434 hitsbranch->SetAddress(&tempo);
437 for (i = 0 ; i < hitsbranch->GetEntries(); i++)
439 hitsbranch->GetEntry(i) ;
441 for ( j = 0 ; j < tempo->GetEntries() ; j++)
443 AliPHOSHit* hit = (AliPHOSHit*)tempo->At(j);
444 new((*hits)[index]) AliPHOSHit( *hit ) ;
453 hitsbranch->SetAddress(hitref);
454 hitsbranch->GetEntry(0) ;
459 //____________________________________________________________________________
460 Int_t AliPHOSLoader::ReadSDigits()
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
468 TObject** sdref = SDigitsRef();
472 sdref = SDigitsRef();
475 TTree * treeS = TreeS();
478 //May happen if file is truncated or new in LoadSDigits
479 //AliError("There is no SDigit Tree");
483 TBranch * branch = treeS->GetBranch(fDetectorName);
485 {//easy, maybe just a new tree
486 //AliError("Cannot find branch PHOS");
490 branch->SetAddress(SDigitsRef());
495 //____________________________________________________________________________
496 Int_t AliPHOSLoader::ReadDigits()
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
504 TObject** dref = DigitsRef();
506 {//if there is not array in folder, create it and put it there
511 TTree * treeD = TreeD();
514 //May happen if file is truncated or new in LoadSDigits
515 //AliError("There is no Digit Tree");
519 TBranch * branch = treeD->GetBranch(fDetectorName);
521 {//easy, maybe just a new tree
522 //AliError("Cannot find branch ",fDetectorName.Data());
526 branch->SetAddress(dref);//connect branch to buffer sitting in folder
527 branch->GetEntry(0);//get first event
532 //____________________________________________________________________________
534 void AliPHOSLoader::Track(Int_t itrack)
536 // Read the first entry of PHOS branch in hit tree gAlice->TreeH()
541 AliError("Can not load hits.");
546 TBranch * hitsbranch = dynamic_cast<TBranch*>(TreeH()->GetListOfBranches()->FindObject("PHOS")) ;
549 cout << "WARNING: AliPHOSLoader::ReadTreeH -> Cannot find branch PHOS" << endl ;
552 if(!Hits()) PostHits();
554 hitsbranch->SetAddress(HitsRef());
555 hitsbranch->GetEntry(itrack);
559 //____________________________________________________________________________
560 Int_t AliPHOSLoader::ReadRecPoints()
562 //Creates and posts to folder an array container,
563 //connects branch in tree (if exists), and reads data to array
565 MakeRecPointsArray();
567 TObjArray * cpva = 0x0 ;
568 TObjArray * emca = 0x0 ;
570 TTree * treeR = TreeR();
574 //May happen if file is truncated or new in LoadSDigits
579 TBranch * emcbranch = treeR->GetBranch(fgkEmcRecPointsBranchName);
581 if (emcbranch == 0x0)
583 AliError(Form("Can not get branch with EMC Rec. Points named %s",
584 fgkEmcRecPointsBranchName.Data()));
589 emcbranch->SetAddress(&emca) ;
590 emcbranch->GetEntry(0) ;
592 TBranch * cpvbranch = treeR->GetBranch(fgkCpvRecPointsBranchName);
593 if (cpvbranch == 0x0)
595 AliError(Form("Can not get branch with CPV Rec. Points named %s",
596 fgkCpvRecPointsBranchName.Data()));
601 cpvbranch->SetAddress(&cpva);
602 cpvbranch->GetEntry(0) ;
606 Int_t maxemc = emca->GetEntries() ;
607 for ( ii= 0 ; ii < maxemc ; ii++ )
608 EmcRecPoints()->Add(emca->At(ii)) ;
610 Int_t maxcpv = cpva->GetEntries() ;
611 for ( ii= 0 ; ii < maxcpv ; ii++ )
612 CpvRecPoints()->Add(cpva->At(ii)) ;
617 //____________________________________________________________________________
618 Int_t AliPHOSLoader::ReadTracks()
620 //Creates and posts to folder an array container,
621 //connects branch in tree (if exists), and reads data to arry
623 TObject** trkref = TracksRef();
625 {//Create and post array
626 MakeTrackSegmentsArray();
627 trkref = TracksRef();
630 TTree * treeT = TreeT();
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");
639 TBranch * branch = treeT->GetBranch(fgkTrackSegmentsBranchName);
640 // AliInfo(Form("Branch named %s is opened: 0x%z",
641 // fgkTrackSegmentsBranchName.Data(),branch));
643 {//easy, maybe just a new tree
644 AliError(Form("Cannot find branch named %s",
645 fgkTrackSegmentsBranchName.Data()));
649 branch->SetAddress(trkref);//connect branch to buffer sitting in folder
650 branch->GetEntry(0);//get first event
654 //____________________________________________________________________________
656 Int_t AliPHOSLoader::ReadRecParticles()
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
662 TObject** recpartref = RecParticlesRef();
664 if ( recpartref == 0x0 )
665 {//Create and post array
666 MakeRecParticlesArray();
667 recpartref = RecParticlesRef();
670 TTree * treeP = TreeP();
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");
680 TBranch * branch = treeP->GetBranch(fgkRecParticlesBranchName);
682 {//easy, maybe just a new tree
683 AliError(Form("Cannot find branch %s",
684 fgkRecParticlesBranchName.Data()));
688 branch->SetAddress(recpartref);//connect branch to buffer sitting in folder
689 branch->GetEntry(0);//get first event
695 /***************************************************************************************/
697 AliPHOSLoader* AliPHOSLoader::GetPHOSLoader(const char* eventfoldername)
699 // Return PHOS loader
700 AliRunLoader* rn = AliRunLoader::GetRunLoader(eventfoldername);
702 printf("Can not find Run Loader in folder %s", eventfoldername);
705 return dynamic_cast<AliPHOSLoader*>(rn->GetLoader("PHOSLoader"));
707 /***************************************************************************************/
709 Bool_t AliPHOSLoader::BranchExists(const TString& recName)
711 // Check if a branch named redName exists
712 if (fBranchTitle.IsNull()) return kFALSE;
713 TString dataname, zername ;
715 if(recName == "SDigits")
718 dataname = GetDetectorName();
719 zername = "AliPHOSSDigitizer" ;
722 if(recName == "Digits"){
724 dataname = GetDetectorName();
725 zername = "AliPHOSDigitizer" ;
728 if(recName == "RecPoints"){
730 dataname = fgkEmcRecPointsBranchName;
731 zername = "AliPHOSClusterizer" ;
734 if(recName == "TrackSegments"){
736 dataname = fgkTrackSegmentsBranchName;
737 zername = "AliPHOSTrackSegmentMaker";
740 if(recName == "RecParticles"){
742 dataname = fgkRecParticlesBranchName;
743 zername = "AliPHOSPID";
752 TObjArray * lob = static_cast<TObjArray*>(tree->GetListOfBranches()) ;
754 TBranch * branch = 0 ;
755 TString titleName(fBranchTitle);
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()));
766 if ( branchName.BeginsWith(zername) && branchTitle.BeginsWith(titleName) ){
767 AliWarning(Form("branch AliPHOS... with title %s ",
768 branch->GetTitle()));
776 void AliPHOSLoader::SetBranchTitle(const TString& btitle)
779 if (btitle.CompareTo(fBranchTitle) == 0) return;
780 fBranchTitle = btitle;
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
790 AliLoader::CleanHits();
792 TClonesArray* hits = const_cast<AliPHOSLoader *>(this)->Hits();
793 if (hits) hits->Clear();
795 //____________________________________________________________________________
797 void AliPHOSLoader::CleanSDigits()const
799 // Clean SDigits array
800 AliLoader::CleanSDigits();
801 TClonesArray* sdigits = const_cast<AliPHOSLoader *>(this)->SDigits();
802 if (sdigits) sdigits->Clear();
805 //____________________________________________________________________________
807 void AliPHOSLoader::CleanDigits()const
809 // Clean Digits array
810 AliLoader::CleanDigits();
811 TClonesArray* digits = const_cast<AliPHOSLoader *>(this)->Digits();
812 if (digits) digits->Clear();
814 //____________________________________________________________________________
816 void AliPHOSLoader::CleanRecPoints()const
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();
825 //____________________________________________________________________________
827 void AliPHOSLoader::CleanTracks()const
829 //Cleans Tracks stuff
830 AliLoader::CleanTracks();//tree
832 //and clear the array
833 TClonesArray* tracks = const_cast<AliPHOSLoader *>(this)->TrackSegments();
834 if (tracks) tracks->Clear();
837 //____________________________________________________________________________
839 void AliPHOSLoader::CleanRecParticles()
841 // Clean RecParticles array
842 TClonesArray *recpar = RecParticles();
843 if (recpar) recpar->Clear();
847 //____________________________________________________________________________
849 void AliPHOSLoader::ReadCalibrationDB(const char * database,const char * filename)
851 // Read calibration data base from file
852 if(fcdb && (strcmp(database,fcdb->GetTitle())==0))
855 TFile * file = gROOT->GetFile(filename) ;
857 file = TFile::Open(filename);
859 AliError(Form("Cannot open file %s", filename)) ;
864 fcdb = dynamic_cast<AliPHOSCalibrationDB *>(file->Get("AliPHOSCalibrationDB")) ;
866 AliError(Form("No database %s in file %s", database, filename)) ;
868 //____________________________________________________________________________
870 // AliPHOSSDigitizer* AliPHOSLoader::PHOSSDigitizer()
872 // //return PHOS SDigitizer
873 // return dynamic_cast<AliPHOSSDigitizer*>(SDigitizer()) ;
876 //____________________________________________________________________________
877 void AliPHOSLoader::MakeHitsArray()
879 // Add Hits array to the data folder
881 TClonesArray* hits = new TClonesArray("AliPHOSHit",1000);
882 hits->SetName(fgkHitsName);
883 GetDetectorDataFolder()->Add(hits);
886 //____________________________________________________________________________
887 void AliPHOSLoader::MakeSDigitsArray()
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);
896 //____________________________________________________________________________
897 void AliPHOSLoader::MakeDigitsArray()
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);
907 //____________________________________________________________________________
908 void AliPHOSLoader::MakeRecPointsArray()
910 // Add RecPoints array to the data folder
911 if ( EmcRecPoints() == 0x0)
913 AliDebug(9, "Making array for EMC");
914 TObjArray* emc = new TObjArray(100) ;
915 emc->SetName(fgkEmcRecPointsName) ;
916 GetDetectorDataFolder()->Add(emc);
919 if ( CpvRecPoints() == 0x0)
921 AliDebug(9, "Making array for CPV");
922 TObjArray* cpv = new TObjArray(100) ;
923 cpv->SetName(fgkCpvRecPointsName);
924 GetDetectorDataFolder()->Add(cpv);
928 //____________________________________________________________________________
929 void AliPHOSLoader::MakeTrackSegmentsArray()
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);
939 //____________________________________________________________________________
940 void AliPHOSLoader::MakeRecParticlesArray()
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);