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 //____________________________________________________________________
19 //////////////////////////////////////////////////////////////////////
21 // class AliRunLoader //
23 // This class aims to be the unque interface for managing data I/O. //
24 // It stores Loaders for all modules which, knows names //
25 // of the files were data are to be stored. //
27 // It aims to substitud AliRun in automatic data managing //
28 // thus there is no necessity of loading gAlice from file in order //
29 // to get access to the data. //
31 // Logical place to put the specific Loader to the given //
32 // detector is detector itself (i.e ITSLoader in ITS). //
33 // But, to load detector object one need to load gAlice, and //
34 // by the way all other detectors with their geometrieces and //
35 // so on. So, if one need to open TPC clusters there is no //
36 // principal need to read everything. //
39 // When RunLoader is read from the file it does not connect to //
40 // the folder structure automatically. It must be connected //
41 // (mounted) manualy. Default event folder is defined by //
42 // AliConfig::GetDefaultEventFolderName() //
43 // but can be mounted elsewhere. Usefull specially in merging case, //
44 // when more than pone session needs to be loaded //
46 //////////////////////////////////////////////////////////////////////
52 #include <TGeometry.h>
53 #include <TObjArray.h>
60 #include "AliConfig.h"
61 #include "AliLoader.h"
62 #include "AliHeader.h"
64 #include "AliDetector.h"
65 #include "AliRunDataStorage.h"
66 #include "AliRunDataFile.h"
68 ClassImp(AliRunLoader)
70 AliRunLoader* AliRunLoader::fgRunLoader = 0x0;
72 const TString AliRunLoader::fgkRunLoaderName("RunLoader");
74 const TString AliRunLoader::fgkHeaderBranchName("Header");
75 const TString AliRunLoader::fgkHeaderContainerName("TE");
76 const TString AliRunLoader::fgkKineContainerName("TreeK");
77 const TString AliRunLoader::fgkTrackRefsContainerName("TreeTR");
78 const TString AliRunLoader::fgkKineBranchName("Particles");
79 const TString AliRunLoader::fgkDefaultKineFileName("Kinematics.root");
80 const TString AliRunLoader::fgkDefaultTrackRefsFileName("TrackRefs.root");
81 const TString AliRunLoader::fgkGAliceName("gAlice");
82 /**************************************************************************/
84 AliRunLoader::AliRunLoader():
92 fTrackRefsDataLoader(0x0),
96 AliConfig::Instance();//force to build the folder structure
97 if (!fgRunLoader) fgRunLoader = this;
99 /**************************************************************************/
101 AliRunLoader::AliRunLoader(const char* eventfoldername):
102 TNamed(fgkRunLoaderName,fgkRunLoaderName),
103 fLoaders(new TObjArray()),
109 fKineDataLoader(new AliDataLoader(fgkDefaultKineFileName,fgkKineContainerName,"Kinematics")),
110 fTrackRefsDataLoader(new AliDataLoader(fgkDefaultTrackRefsFileName,fgkTrackRefsContainerName,"Track References")),
115 SetEventFolderName(eventfoldername);
116 if (!fgRunLoader) fgRunLoader = this;
118 /**************************************************************************/
120 AliRunLoader::AliRunLoader(const AliRunLoader &rl):
128 fKineDataLoader(0x0),
129 fTrackRefsDataLoader(0x0),
138 /**************************************************************************/
140 AliRunLoader::~AliRunLoader()
143 if (fgRunLoader == this) fgRunLoader = 0x0;
149 fLoaders->SetOwner();
153 delete fKineDataLoader;
154 delete fTrackRefsDataLoader;
159 //fEventFolder is deleted by the way of removing - TopAliceFolder owns it
164 /**************************************************************************/
166 AliRunLoader::AliRunLoader(TFolder* topfolder):
167 TNamed(fgkRunLoaderName,fgkRunLoaderName),
168 fLoaders(new TObjArray()),
169 fEventFolder(topfolder),
174 fKineDataLoader(new AliDataLoader(fgkDefaultKineFileName,fgkKineContainerName,"Kinematics")),
175 fTrackRefsDataLoader(new AliDataLoader(fgkDefaultTrackRefsFileName,fgkTrackRefsContainerName,"Track References")),
182 TString errmsg("Parameter is NULL");
183 AliError(errmsg.Data());
188 TObject* obj = fEventFolder->FindObject(fgkRunLoaderName);
190 { //if it is, then sth. is going wrong... exits aliroot session
191 TString errmsg("In Event Folder Named ");
192 errmsg+=fEventFolder->GetName();
193 errmsg+=" object named "+fgkRunLoaderName+" already exists. I am confused ...";
195 AliError(errmsg.Data());
197 return;//never reached
200 if (!fgRunLoader) fgRunLoader = this;
202 fEventFolder->Add(this);//put myself to the folder to accessible for all
205 /**************************************************************************/
207 void AliRunLoader::Copy(TObject &) const
209 AliFatal("Not implemented");
211 /**************************************************************************/
213 Int_t AliRunLoader::GetEvent(Int_t evno)
215 //Gets event number evno
216 //Reloads all data properly
217 if (fCurrentEvent == evno) return 0;
221 AliError("Can not give the event with negative number");
225 if (evno >= GetNumberOfEvents())
227 AliError(Form("There is no event with number %d",evno));
231 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
232 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
233 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
234 AliDebug(1, Form(" GETTING EVENT %d",evno));
235 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
236 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
237 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
239 fCurrentEvent = evno;
243 //Reload header (If header was loaded)
246 retval = TreeE()->GetEvent(fCurrentEvent);
249 AliError(Form("Cannot find event: %d\n ",fCurrentEvent));
253 //Reload stack (If header was loaded)
254 if (TreeE()) fStack = GetHeader()->Stack();
255 //Set event folder in stack (it does not mean that we read kinematics from file)
258 fStack->SetEventFolderName(fEventFolder->GetName());
262 AliWarning("Stack not found in header");
268 AliError(Form("Error occured while setting event %d",evno));
272 //Post Track References
273 retval = fTrackRefsDataLoader->GetEvent();
276 AliError(Form("Error occured while GetEvent for Track References. Event %d",evno));
280 //Read Kinematics if loaded
281 fKineDataLoader->GetEvent();
284 AliError(Form("Error occured while GetEvent for Kinematics. Event %d",evno));
288 if (fStack && fKineDataLoader->GetBaseLoader(0)->IsLoaded()) fStack->GetEvent();
290 //Trigger data reloading in all loaders
291 TIter next(fLoaders);
293 while((loader = (AliLoader*)next()))
295 retval = loader->GetEvent();
298 AliError(Form("Error occured while getting event for %s. Event %d.",
299 loader->GetDetectorName().Data(), evno));
304 SetDetectorAddresses();
308 /**************************************************************************/
309 Int_t AliRunLoader::SetEvent()
311 //if kinematics was loaded Cleans folder data
315 retval = fKineDataLoader->SetEvent();
318 AliError("SetEvent for Kinamtics Data Loader retutned error.");
321 retval = fTrackRefsDataLoader->SetEvent();
324 AliError("SetEvent for Track References Data Loader retutned error.");
328 TIter next(fLoaders);
330 while((loader = (AliLoader*)next()))
332 retval = loader->SetEvent();
335 AliError(Form("SetEvent for %s Data Loader retutned error.",loader->GetName()));
342 /**************************************************************************/
344 Int_t AliRunLoader::SetEventNumber(Int_t evno)
346 //cleans folders and sets the root dirs in files
347 if (fCurrentEvent == evno) return 0;
348 fCurrentEvent = evno;
352 /**************************************************************************/
353 const TObject* AliRunLoader::GetRunObject(const char* name) const
355 //Get an object from the run data storage
357 if (!AliRunDataStorage::Instance()) {
358 AliWarning("No run data storage defined. Using AliRunDataFile.");
361 return AliRunDataStorage::Instance()->Get(name, GetEventNumber());
364 /**************************************************************************/
365 AliRunLoader* AliRunLoader::Open
366 (const char* filename, const char* eventfoldername, Option_t* option)
368 //Opens a desired file 'filename'
369 //gets the the run-Loader and mounts it desired folder
370 //returns the pointer to run Loader which can be further used for accessing data
371 //in case of error returns NULL
373 static const TString kwebaddress("http://alisoft.cern.ch/people/skowron/codedoc/split/index.html");
374 AliDebugClass(1,Form("\n\n\nNew I/O strcture: See more info:\n %s\n\n\n",kwebaddress.Data()));
376 AliRunLoader* result = 0x0;
378 /* ************************************************ */
379 /* Chceck if folder with given name already exists */
380 /* ************************************************ */
382 TObject* obj = AliConfig::Instance()->GetTopFolder()->FindObject(eventfoldername);
385 TFolder* fold = dynamic_cast<TFolder*>(obj);
388 AliErrorClass("Such a obejct already exists in top alice folder and it is not a folder.");
392 //check if we can get RL from that folder
393 result = AliRunLoader::GetRunLoader(eventfoldername);
396 AliErrorClass(Form("Folder %s already exists, and can not find session there. Can not mount.",eventfoldername));
400 if (result->GetFileName().CompareTo(filename) != 0)
402 AliErrorClass("Other file is mounted in demanded folder. Can not mount.");
406 //check if now is demanded (re)creation
407 if ( AliLoader::TestFileOption(option) == kFALSE)
409 AliErrorClass(Form("Session already exists in folder %s and this session option is %s. Unable to proceed.",
410 eventfoldername,option));
414 //check if demanded option is update and existing one
415 TString tmpstr(option);
416 if ( (tmpstr.CompareTo("update",TString::kIgnoreCase) == 0) &&
417 (result->fGAFile->IsWritable() == kFALSE) )
419 AliErrorClass(Form("Session already exists in folder %s and is not writable while this session option is %s. Unable to proceed.",
420 eventfoldername,option));
424 AliWarningClass("Session is already opened and mounted in demanded folder");
426 } //end of checking in case of existance of object named identically that folder session is being opened
429 TFile * gAliceFile = TFile::Open(filename,option);//open a file
431 {//null pointer returned
432 AliFatalClass(Form("Can not open file %s.",filename));
436 if (gAliceFile->IsOpen() == kFALSE)
437 {//pointer to valid object returned but file is not opened
438 AliErrorClass(Form("Can not open file %s.",filename));
442 //if file is "read" or "update" than we try to find AliRunLoader there - if not found cry and exit
443 //else create new AliRunLoader
444 if ( AliLoader::TestFileOption(option) )
446 AliDebugClass(1, "Reading RL from file");
448 result = dynamic_cast<AliRunLoader*>(gAliceFile->Get(fgkRunLoaderName));//get the run Loader from the file
451 AliErrorClass(Form("Can not find run-Loader in file %s.",filename));
452 delete gAliceFile;//close the file
455 Int_t tmp = result->SetEventFolderName(eventfoldername);//mount a event folder
456 if (tmp)//if SetEvent returned error
458 AliErrorClass(Form("Can not mount event in folder %s.",eventfoldername));
459 delete result; //delete run-Loader
460 delete gAliceFile;//close the file
466 AliDebugClass(1, Form("Creating new AliRunLoader. Folder name is %s",eventfoldername));
469 result = new AliRunLoader(eventfoldername);
471 catch (TString& errmsg)
473 AliErrorClass(Form("AliRunLoader constrcutor has thrown exception: %s\n",errmsg.Data()));
475 delete gAliceFile;//close the file
480 //procedure for extracting dir name from the file name
481 TString fname(filename);
482 Int_t nsl = fname.Last('/');//look for slash in file name
486 Int_t nsl = fname.Last(':');//look for colon e.g. rfio:galice.root
487 if (nsl < 0) dirname = ".";//not found
488 else dirname = fname.Remove(nsl);//found
490 else dirname = fname.Remove(nsl);//slash found
492 AliDebugClass(1, Form("Dir name is : %s",dirname.Data()));
494 result->SetDirName(dirname);
495 result->SetGAliceFile(gAliceFile);//set the pointer to gAliceFile
496 if (!fgRunLoader) fgRunLoader = result; //PH get access from any place
499 /**************************************************************************/
500 Int_t AliRunLoader::GetNumberOfEvents()
502 //returns number of events in Run
506 retval = LoadHeader();
509 AliError("Error occured while loading header");
513 return (Int_t)TreeE()->GetEntries();
515 /**************************************************************************/
517 void AliRunLoader::MakeHeader()
519 //Makes header and connects it to header tree (if it exists)
523 AliDebug(1, "Creating new Header Object");
524 fHeader= new AliHeader();
526 TTree* tree = TreeE();
529 AliDebug(1, "Got Tree from folder.");
530 TBranch* branch = tree->GetBranch(fgkHeaderBranchName);
533 AliDebug(1, "Creating new branch");
534 branch = tree->Branch(fgkHeaderBranchName, "AliHeader", &fHeader, 4000, 0);
535 branch->SetAutoDelete(kFALSE);
539 AliDebug(1, "Got Branch from Tree");
540 branch->SetAddress(&fHeader);
541 tree->GetEvent(fCurrentEvent);
542 fStack = fHeader->Stack(); //should be safe - if we created Stack, header returns pointer to the same object
545 fStack->SetEventFolderName(fEventFolder->GetName());
546 if (TreeK()) fStack->GetEvent();
550 AliDebug(1, "Header does not have a stack.");
554 AliDebug(1, "Exiting MakeHeader method");
556 /**************************************************************************/
558 void AliRunLoader::MakeStack()
560 //Creates the stack object - do not connect the tree
563 fStack = new AliStack(10000);
564 fStack->SetEventFolderName(fEventFolder->GetName());
568 /**************************************************************************/
570 void AliRunLoader::MakeTree(Option_t *option)
573 const char *oK = strstr(option,"K"); //Kine
574 const char *oE = strstr(option,"E"); //Header
578 if (fKineDataLoader->GetBaseLoader(0)->IsLoaded() == kFALSE)
580 AliError("Load Kinematics first");
584 fKineDataLoader->MakeTree();
586 fStack->ConnectTree();
587 WriteKinematics("OVERWRITE");
594 TTree* tree = new TTree(fgkHeaderContainerName,"Tree with Headers");
595 GetEventFolder()->Add(tree);
597 WriteHeader("OVERWRITE");
600 TIter next(fLoaders);
602 while((loader = (AliLoader*)next()))
604 loader->MakeTree(option);
608 /**************************************************************************/
610 Int_t AliRunLoader::LoadgAlice()
612 //Loads gAlice from file
615 AliWarning("AliRun is already in folder. Unload first.");
618 AliRun* alirun = dynamic_cast<AliRun*>(fGAFile->Get(fgkGAliceName));
621 AliError(Form("Can not find gAlice in file %s",fGAFile->GetName()));
624 alirun->SetRunLoader(this);
627 AliWarning(Form("gAlice already exists. Putting retrived object in folder named %s",
628 GetEventFolder()->GetName()));
634 SetDetectorAddresses();//calls SetTreeAddress for all detectors
637 /**************************************************************************/
639 Int_t AliRunLoader::LoadHeader()
641 //loads treeE and reads header object for current event
644 AliWarning("Header is already loaded. Use ReloadHeader to force reload. Nothing done");
648 if (GetEventFolder() == 0x0)
650 AliError("Event folder not specified yet");
656 AliError("Session not opened. Use AliRunLoader::Open");
660 if (fGAFile->IsOpen() == kFALSE)
662 AliError("Session not opened. Use AliRunLoader::Open");
666 TTree* tree = dynamic_cast<TTree*>(fGAFile->Get(fgkHeaderContainerName));
669 AliError(Form("Can not find header tree named %s in file %s",
670 fgkHeaderContainerName.Data(),fGAFile->GetName()));
674 if (tree == TreeE()) return 0;
677 GetEventFolder()->Add(tree);
678 MakeHeader();//creates header object and connects to tree
682 /**************************************************************************/
684 Int_t AliRunLoader::LoadKinematics(Option_t* option)
686 //Loads the kinematics
687 Int_t retval = fKineDataLoader->GetBaseLoader(0)->Load(option);
690 AliError("Error occured while loading kinamatics tree.");
695 retval = fStack->GetEvent();
696 if ( retval == kFALSE)
698 AliError("Error occured while loading kinamatics tree.");
705 /**************************************************************************/
707 Int_t AliRunLoader::OpenDataFile(const TString& filename,TFile*& file,TDirectory*& dir,Option_t* opt,Int_t cl)
709 //Opens File with kinematics
712 if (file->IsOpen() == kFALSE)
713 {//pointer is not null but file is not opened
714 AliWarning("Pointer to file is not null, but file is not opened");//risky any way
716 file = 0x0; //proceed with opening procedure
720 AliWarning(Form("File %s already opened",filename.Data()));
724 //try to find if that file is opened somewere else
725 file = (TFile *)( gROOT->GetListOfFiles()->FindObject(filename) );
728 if(file->IsOpen() == kTRUE)
730 AliWarning(Form("File %s already opened by sombody else.",file->GetName()));
735 file = TFile::Open(filename,opt);
738 AliError(Form("Can not open file %s",filename.Data()));
741 if (file->IsOpen() == kFALSE)
742 {//file is not opened
743 AliError(Form("Can not open file %s",filename.Data()));
747 file->SetCompressionLevel(cl);
749 dir = AliLoader::ChangeDir(file,fCurrentEvent);
752 AliError(Form("Can not change to root directory in file %s",filename.Data()));
757 /**************************************************************************/
759 TTree* AliRunLoader::TreeE() const
761 //returns the tree from folder; shortcut method
762 if (AliDebugLevel() > 10) fEventFolder->ls();
763 TObject *obj = fEventFolder->FindObject(fgkHeaderContainerName);
764 return (obj)?dynamic_cast<TTree*>(obj):0x0;
766 /**************************************************************************/
768 AliHeader* AliRunLoader::GetHeader() const
770 //returns pointer header object
773 /**************************************************************************/
775 TTree* AliRunLoader::TreeK() const
777 //returns the tree from folder; shortcut method
778 TObject *obj = GetEventFolder()->FindObject(fgkKineContainerName);
779 return (obj)?dynamic_cast<TTree*>(obj):0x0;
781 /**************************************************************************/
783 TTree* AliRunLoader::TreeTR() const
785 //returns the tree from folder; shortcut method
786 TObject* obj = GetEventFolder()->FindObject(fgkTrackRefsContainerName);
787 return (obj)?dynamic_cast<TTree*>(obj):0x0;
789 /**************************************************************************/
791 AliRun* AliRunLoader::GetAliRun() const
793 //returns AliRun which sits in the folder
794 if (fEventFolder == 0x0) return 0x0;
795 TObject *obj = fEventFolder->FindObject(fgkGAliceName);
796 return (obj)?dynamic_cast<AliRun*>(obj):0x0;
798 /**************************************************************************/
800 Int_t AliRunLoader::WriteGeometry(Option_t* /*opt*/)
802 //writes geometry to the file
804 TGeometry* geo = GetAliRun()->GetGeometry();
807 AliError("Can not get geometry from gAlice");
813 /**************************************************************************/
815 Int_t AliRunLoader::WriteHeader(Option_t* opt)
818 AliDebug(1, "WRITING HEADER");
820 TTree* tree = TreeE();
823 AliWarning("Can not find Header Tree in Folder");
826 if (fGAFile->IsWritable() == kFALSE)
828 AliError(Form("File %s is not writable",fGAFile->GetName()));
832 TObject* obj = fGAFile->Get(fgkHeaderContainerName);
834 { //if they exist, see if option OVERWRITE is used
836 if(tmp.Contains("OVERWRITE",TString::kIgnoreCase) == 0)
837 {//if it is not used - give an error message and return an error code
838 AliError("Tree already exisists. Use option \"OVERWRITE\" to overwrite previous data");
843 tree->SetDirectory(fGAFile);
844 tree->Write(0,TObject::kOverwrite);
846 AliDebug(1, "WRITTEN\n\n");
850 /**************************************************************************/
852 Int_t AliRunLoader::WriteAliRun(Option_t* /*opt*/)
854 //writes AliRun object to the file
856 if (GetAliRun()) GetAliRun()->Write();
859 /**************************************************************************/
861 Int_t AliRunLoader::WriteKinematics(Option_t* opt)
864 return fKineDataLoader->GetBaseLoader(0)->WriteData(opt);
866 /**************************************************************************/
867 Int_t AliRunLoader::WriteTrackRefs(Option_t* opt)
869 //writes Track References tree
870 return fTrackRefsDataLoader->GetBaseLoader(0)->WriteData(opt);
872 /**************************************************************************/
874 Int_t AliRunLoader::WriteHits(Option_t* opt)
876 //Calls WriteHits for all loaders
879 TIter next(fLoaders);
881 while((loader = (AliLoader*)next()))
883 res = loader->WriteHits(opt);
886 AliError(Form("Failed to write hits for %s (%d)",loader->GetDetectorName().Data(),res));
892 /**************************************************************************/
894 Int_t AliRunLoader::WriteSDigits(Option_t* opt)
896 //Calls WriteSDigits for all loaders
899 TIter next(fLoaders);
901 while((loader = (AliLoader*)next()))
903 res = loader->WriteSDigits(opt);
906 AliError(Form("Failed to write summable digits for %s.",loader->GetDetectorName().Data()));
912 /**************************************************************************/
914 Int_t AliRunLoader::WriteDigits(Option_t* opt)
916 //Calls WriteDigits for all loaders
919 TIter next(fLoaders);
921 while((loader = (AliLoader*)next()))
923 res = loader->WriteDigits(opt);
926 AliError(Form("Failed to write digits for %s.",loader->GetDetectorName().Data()));
932 /**************************************************************************/
934 Int_t AliRunLoader::WriteRecPoints(Option_t* opt)
936 //Calls WriteRecPoints for all loaders
939 TIter next(fLoaders);
941 while((loader = (AliLoader*)next()))
943 res = loader->WriteRecPoints(opt);
946 AliError(Form("Failed to write Reconstructed Points for %s.",
947 loader->GetDetectorName().Data()));
953 /**************************************************************************/
955 Int_t AliRunLoader::WriteTracks(Option_t* opt)
957 //Calls WriteTracks for all loaders
960 TIter next(fLoaders);
962 while((loader = (AliLoader*)next()))
964 res = loader->WriteTracks(opt);
967 AliError(Form("Failed to write Tracks for %s.",
968 loader->GetDetectorName().Data()));
974 /**************************************************************************/
976 Int_t AliRunLoader::WriteRunLoader(Option_t* /*opt*/)
978 //Writes itself to the file
980 this->Write(0,TObject::kOverwrite);
983 /**************************************************************************/
985 Int_t AliRunLoader::SetEventFolderName(const TString& name)
987 //sets top folder name for this run; of alread
990 AliError("Name is empty");
994 //check if such a folder already exists - try to find it in alice top folder
995 TObject* obj = AliConfig::Instance()->GetTopFolder()->FindObject(name);
998 TFolder* fold = dynamic_cast<TFolder*>(obj);
1001 AliError("Such a obejct already exists in top alice folder and it is not a folder.");
1004 //folder which was found is our folder
1005 if (fEventFolder == fold)
1011 AliError("Such a folder already exists in top alice folder. Can not mount.");
1016 //event is alredy connected, just change name of the folder
1019 fEventFolder->SetName(name);
1023 if (fKineDataLoader == 0x0)
1024 fKineDataLoader = new AliDataLoader(fgkDefaultKineFileName,fgkKineContainerName,"Kinematics");
1026 if ( fTrackRefsDataLoader == 0x0)
1027 fTrackRefsDataLoader = new AliDataLoader(fgkDefaultTrackRefsFileName,fgkTrackRefsContainerName,"Track References");
1029 //build the event folder structure
1030 AliDebug(1, Form("Creating new event folder named %s",name.Data()));
1031 fEventFolder = AliConfig::Instance()->BuildEventFolder(name,"Event Folder");
1032 fEventFolder->Add(this);//put myself to the folder to accessible for all
1034 if (Stack()) Stack()->SetEventFolderName(fEventFolder->GetName());
1035 TIter next(fLoaders);
1037 while((loader = (AliLoader*)next()))
1039 loader->Register(fEventFolder);//build folder structure for this detector
1042 fKineDataLoader->SetEventFolder(GetEventFolder());
1043 fTrackRefsDataLoader->SetEventFolder(GetEventFolder());
1044 fKineDataLoader->SetFolder(GetEventFolder());
1045 fTrackRefsDataLoader->SetFolder(GetEventFolder());
1047 fEventFolder->SetOwner();
1050 /**************************************************************************/
1052 void AliRunLoader::AddLoader(AliLoader* loader)
1054 //Adds the Loader for given detector
1055 if (loader == 0x0) //if null shout and exit
1057 AliError("Parameter is NULL");
1060 loader->SetDirName(fUnixDirName);
1061 if (fEventFolder) loader->SetEventFolder(fEventFolder); //if event folder is already defined,
1062 //pass information to the Loader
1063 fLoaders->Add(loader);//add the Loader to the array
1065 /**************************************************************************/
1067 void AliRunLoader::AddLoader(AliDetector* det)
1069 //Asks module (detector) ro make a Loader and stores in the array
1070 if (det == 0x0) return;
1071 AliLoader* get = det->GetLoader();//try to get loader
1072 if (get == 0x0) get = det->MakeLoader(fEventFolder->GetName());//if did not obtain, ask to make it
1076 AliDebug(1, Form("Detector: %s Loader : %s",det->GetName(),get->GetName()));
1081 /**************************************************************************/
1083 AliLoader* AliRunLoader::GetLoader(const char* detname) const
1085 //returns loader for given detector
1086 //note that naming convention is TPCLoader not just TPC
1087 return (AliLoader*)fLoaders->FindObject(detname);
1090 /**************************************************************************/
1092 AliLoader* AliRunLoader::GetLoader(AliDetector* det) const
1094 //get loader for detector det
1095 if(det == 0x0) return 0x0;
1096 TString getname(det->GetName());
1098 AliDebug(1, Form(" Loader name is %s",getname.Data()));
1099 return GetLoader(getname);
1102 /**************************************************************************/
1104 void AliRunLoader::CleanFolders()
1106 // fEventFolder->Add(this);//put myself to the folder to accessible for all
1112 /**************************************************************************/
1114 void AliRunLoader::CleanDetectors()
1116 //Calls CleanFolders for all detectors
1117 TIter next(fLoaders);
1119 while((loader = (AliLoader*)next()))
1121 loader->CleanFolders();
1124 /**************************************************************************/
1126 void AliRunLoader::RemoveEventFolder()
1128 //remove all the tree of event
1129 //all the stuff changing EbE stays untached (PDGDB, tasks, etc.)
1131 if (fEventFolder == 0x0) return;
1132 fEventFolder->SetOwner(kFALSE);//don't we want to deleted while removing the folder that we are sitting in
1133 fEventFolder->Remove(this);//remove us drom folder
1135 AliConfig::Instance()->GetTopFolder()->SetOwner(); //brings ownership back for fEventFolder since it sits in top folder
1136 AliConfig::Instance()->GetTopFolder()->Remove(fEventFolder); //remove the event tree
1137 delete fEventFolder;
1139 /**************************************************************************/
1141 void AliRunLoader::SetGAliceFile(TFile* gafile)
1143 //sets pointer to galice.root file
1147 /**************************************************************************/
1149 Int_t AliRunLoader::LoadHits(Option_t* detectors,Option_t* opt)
1151 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1153 AliDebug(1, "Loading Hits");
1157 const char* oAll = strstr(detectors,"all");
1160 AliDebug(1, "Option is All");
1165 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1166 loaders = &arr;//get the pointer array
1169 AliDebug(1, Form("For detectors. Number of detectors chosen for loading %d",loaders->GetEntries()));
1171 TIter next(loaders);
1173 while((loader = (AliLoader*)next()))
1175 AliDebug(1, Form(" Calling LoadHits(%s) for %s",opt,loader->GetName()));
1176 loader->LoadHits(opt);
1178 AliDebug(1, "Done");
1182 /**************************************************************************/
1184 Int_t AliRunLoader::LoadSDigits(Option_t* detectors,Option_t* opt)
1186 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1191 const char* oAll = strstr(detectors,"all");
1198 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1199 loaders = &arr;//get the pointer to array
1202 TIter next(loaders);
1204 while((loader = (AliLoader*)next()))
1206 loader->LoadSDigits(opt);
1211 /**************************************************************************/
1213 Int_t AliRunLoader::LoadDigits(Option_t* detectors,Option_t* opt)
1215 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1220 const char* oAll = strstr(detectors,"all");
1227 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1228 loaders = &arr;//get the pointer array
1231 TIter next(loaders);
1233 while((loader = (AliLoader*)next()))
1235 loader->LoadDigits(opt);
1239 /**************************************************************************/
1241 Int_t AliRunLoader::LoadRecPoints(Option_t* detectors,Option_t* opt)
1243 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1248 const char* oAll = strstr(detectors,"all");
1255 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1256 loaders = &arr;//get the pointer array
1259 TIter next(loaders);
1261 while((loader = (AliLoader*)next()))
1263 loader->LoadRecPoints(opt);
1267 /**************************************************************************/
1269 Int_t AliRunLoader::LoadRecParticles(Option_t* detectors,Option_t* opt)
1271 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1276 const char* oAll = strstr(detectors,"all");
1283 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1284 loaders = &arr;//get the pointer array
1287 TIter next(loaders);
1289 while((loader = (AliLoader*)next()))
1291 loader->LoadRecParticles(opt);
1295 /**************************************************************************/
1297 Int_t AliRunLoader::LoadTracks(Option_t* detectors,Option_t* opt)
1299 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1304 const char* oAll = strstr(detectors,"all");
1311 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1312 loaders = &arr;//get the pointer array
1315 TIter next(loaders);
1317 while((loader = (AliLoader*)next()))
1319 loader->LoadTracks(opt);
1323 /**************************************************************************/
1325 void AliRunLoader::UnloadHits(Option_t* detectors)
1327 //unloads hits for detectors specified in parameter
1331 const char* oAll = strstr(detectors,"all");
1338 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1339 loaders = &arr;//get the pointer to array
1342 TIter next(loaders);
1344 while((loader = (AliLoader*)next()))
1346 loader->UnloadHits();
1349 /**************************************************************************/
1351 void AliRunLoader::UnloadSDigits(Option_t* detectors)
1353 //unloads SDigits for detectors specified in parameter
1357 const char* oAll = strstr(detectors,"all");
1364 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1365 loaders = &arr;//get the pointer to array
1368 TIter next(loaders);
1370 while((loader = (AliLoader*)next()))
1372 loader->UnloadSDigits();
1375 /**************************************************************************/
1377 void AliRunLoader::UnloadDigits(Option_t* detectors)
1379 //unloads Digits for detectors specified in parameter
1383 const char* oAll = strstr(detectors,"all");
1390 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1391 loaders = &arr;//get the pointer to array
1394 TIter next(loaders);
1396 while((loader = (AliLoader*)next()))
1398 loader->UnloadDigits();
1401 /**************************************************************************/
1403 void AliRunLoader::UnloadRecPoints(Option_t* detectors)
1405 //unloads RecPoints for detectors specified in parameter
1409 const char* oAll = strstr(detectors,"all");
1416 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1417 loaders = &arr;//get the pointer to array
1420 TIter next(loaders);
1422 while((loader = (AliLoader*)next()))
1424 loader->UnloadRecPoints();
1427 /**************************************************************************/
1429 void AliRunLoader::UnloadAll(Option_t* detectors)
1431 //calls UnloadAll for detectors names specified in parameter
1432 // option "all" passed can be passed
1436 const char* oAll = strstr(detectors,"all");
1443 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1444 loaders = &arr;//get the pointer to array
1447 TIter next(loaders);
1449 while((loader = (AliLoader*)next()))
1451 loader->UnloadAll();
1454 /**************************************************************************/
1456 void AliRunLoader::UnloadTracks(Option_t* detectors)
1458 //unloads Tracks for detectors specified in parameter
1462 const char* oAll = strstr(detectors,"all");
1469 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1470 loaders = &arr;//get the pointer to array
1473 TIter next(loaders);
1475 while((loader = (AliLoader*)next()))
1477 loader->UnloadTracks();
1480 /**************************************************************************/
1482 void AliRunLoader::UnloadRecParticles(Option_t* detectors)
1484 //unloads Particles for detectors specified in parameter
1488 const char* oAll = strstr(detectors,"all");
1495 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1496 loaders = &arr;//get the pointer to array
1499 TIter next(loaders);
1501 while((loader = (AliLoader*)next()))
1503 loader->UnloadRecParticles();
1506 /**************************************************************************/
1508 AliRunLoader* AliRunLoader::GetRunLoader(const char* eventfoldername)
1510 //returns RunLoader from folder named eventfoldername
1511 TFolder* evfold= dynamic_cast<TFolder*>(AliConfig::Instance()->GetTopFolder()->FindObject(eventfoldername));
1516 AliRunLoader* runget = dynamic_cast<AliRunLoader*>(evfold->FindObject(AliRunLoader::fgkRunLoaderName));
1520 /**************************************************************************/
1522 AliLoader* AliRunLoader::GetDetectorLoader(const char* detname, const char* eventfoldername)
1524 //get the loader of the detector with the given name from the global
1526 AliRunLoader* runLoader = GetRunLoader(eventfoldername);
1528 AliErrorClass("No run loader found");
1531 return runLoader->GetDetectorLoader(detname);
1533 /**************************************************************************/
1535 AliLoader* AliRunLoader::GetDetectorLoader(const char* detname)
1537 //get the loader of the detector with the given name from the global
1540 char loadername[256];
1541 sprintf(loadername, "%sLoader", detname);
1542 AliLoader* loader = GetLoader(loadername);
1544 AliError(Form("No loader for %s found", detname));
1549 /**************************************************************************/
1551 TTree* AliRunLoader::GetTreeH(const char* detname, Bool_t maketree, const char* eventfoldername)
1553 //get the tree with hits of the detector with the given name
1554 //if maketree is true and the tree does not exist, the tree is created
1555 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1556 if (!loader) return NULL;
1557 if (!loader->TreeH() && maketree) loader->MakeTree("H");
1558 return loader->TreeH();
1561 /**************************************************************************/
1563 TTree* AliRunLoader::GetTreeH(const char* detname, Bool_t maketree)
1565 //get the tree with hits of the detector with the given name
1566 //if maketree is true and the tree does not exist, the tree is created
1567 AliLoader* loader = GetDetectorLoader(detname);
1568 if (!loader) return NULL;
1569 if (!loader->TreeH() && maketree) loader->MakeTree("H");
1570 return loader->TreeH();
1572 /**************************************************************************/
1574 TTree* AliRunLoader::GetTreeS(const char* detname, Bool_t maketree,const char* eventfoldername)
1576 //get the tree with summable digits of the detector with the given name
1577 //if maketree is true and the tree does not exist, the tree is created
1578 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1579 if (!loader) return NULL;
1580 if (!loader->TreeS() && maketree) loader->MakeTree("S");
1581 return loader->TreeS();
1583 /**************************************************************************/
1585 TTree* AliRunLoader::GetTreeS(const char* detname, Bool_t maketree)
1587 //get the tree with summable digits of the detector with the given name
1588 //if maketree is true and the tree does not exist, the tree is created
1589 AliLoader* loader = GetDetectorLoader(detname);
1590 if (!loader) return NULL;
1591 if (!loader->TreeS() && maketree) loader->MakeTree("S");
1592 return loader->TreeS();
1594 /**************************************************************************/
1596 TTree* AliRunLoader::GetTreeD(const char* detname, Bool_t maketree,const char* eventfoldername)
1598 //get the tree with digits of the detector with the given name
1599 //if maketree is true and the tree does not exist, the tree is created
1600 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1601 if (!loader) return NULL;
1602 if (!loader->TreeD() && maketree) loader->MakeTree("D");
1603 return loader->TreeD();
1605 /**************************************************************************/
1607 TTree* AliRunLoader::GetTreeD(const char* detname, Bool_t maketree)
1609 //get the tree with digits of the detector with the given name
1610 //if maketree is true and the tree does not exist, the tree is created
1611 AliLoader* loader = GetDetectorLoader(detname);
1612 if (!loader) return NULL;
1613 if (!loader->TreeD() && maketree) loader->MakeTree("D");
1614 return loader->TreeD();
1616 /**************************************************************************/
1617 TTree* AliRunLoader::GetTreeR(const char* detname, Bool_t maketree,const char* eventfoldername)
1619 //get the tree with clusters of the detector with the given name
1620 //if maketree is true and the tree does not exist, the tree is created
1621 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1622 if (!loader) return NULL;
1623 if (!loader->TreeR() && maketree) loader->MakeTree("R");
1624 return loader->TreeR();
1626 /**************************************************************************/
1628 TTree* AliRunLoader::GetTreeR(const char* detname, Bool_t maketree)
1630 //get the tree with clusters of the detector with the given name
1631 //if maketree is true and the tree does not exist, the tree is created
1632 AliLoader* loader = GetDetectorLoader(detname);
1633 if (!loader) return NULL;
1634 if (!loader->TreeR() && maketree) loader->MakeTree("R");
1635 return loader->TreeR();
1637 /**************************************************************************/
1639 TTree* AliRunLoader::GetTreeT(const char* detname, Bool_t maketree,const char* eventfoldername)
1641 //get the tree with tracks of the detector with the given name
1642 //if maketree is true and the tree does not exist, the tree is created
1643 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1644 if (!loader) return NULL;
1645 if (!loader->TreeT() && maketree) loader->MakeTree("T");
1646 return loader->TreeT();
1648 /**************************************************************************/
1650 TTree* AliRunLoader::GetTreeT(const char* detname, Bool_t maketree)
1652 //get the tree with tracks of the detector with the given name
1653 //if maketree is true and the tree does not exist, the tree is created
1654 AliLoader* loader = GetDetectorLoader(detname);
1655 if (!loader) return NULL;
1656 if (!loader->TreeT() && maketree) loader->MakeTree("T");
1657 return loader->TreeT();
1659 /**************************************************************************/
1661 TTree* AliRunLoader::GetTreeP(const char* detname, Bool_t maketree,const char* eventfoldername)
1663 //get the tree with particles of the detector with the given name
1664 //if maketree is true and the tree does not exist, the tree is created
1665 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1666 if (!loader) return NULL;
1667 if (!loader->TreeP() && maketree) loader->MakeTree("P");
1668 return loader->TreeP();
1670 /**************************************************************************/
1672 TTree* AliRunLoader::GetTreeP(const char* detname, Bool_t maketree)
1674 //get the tree with particles of the detector with the given name
1675 //if maketree is true and the tree does not exist, the tree is created
1676 AliLoader* loader = GetDetectorLoader(detname);
1677 if (!loader) return NULL;
1678 if (!loader->TreeP() && maketree) loader->MakeTree("P");
1679 return loader->TreeP();
1682 /**************************************************************************/
1684 void AliRunLoader::CdGAFile()
1686 //sets gDirectory to galice file
1688 if(fGAFile) fGAFile->cd();
1691 /**************************************************************************/
1693 void AliRunLoader::GetListOfDetectors(const char * namelist,TObjArray& pointerarray) const
1695 //this method looks for all Loaders corresponding
1696 //to names (many) specified in namelist i.e. namelist ("ITS TPC TRD")
1700 strcpy(dets,namelist);//compiler cries when char* = const Option_t*;
1701 dets[strlen(dets)+1] = '\n';//set endl at the end of string
1706 tmp = sscanf(pdet,"%s",buff);//read the string from the input string pdet into buff
1707 if ( (buff[0] == 0) || (tmp == 0) ) break; //if not read
1709 pdet = strstr(pdet,buff) + strlen(buff);//move the input pointer about number of bytes (letters) read
1710 //I am aware that is a little complicated. I don't know the number of spaces between detector names
1711 //so I read the string, than I find where it starts (strstr) and move the pointer about length of a string
1712 //If there is a better way, please write me (Piotr.Skowronski@cern.ch)
1713 //construct the Loader name
1714 TString getname(buff);
1716 AliLoader* loader = GetLoader(getname);//get the Loader
1719 pointerarray.Add(loader);
1723 AliError(Form("Can not find Loader for %s",buff));
1729 /*****************************************************************************/
1731 void AliRunLoader::Clean(const TString& name)
1733 //removes object with given name from event folder and deletes it
1734 if (GetEventFolder() == 0x0) return;
1735 TObject* obj = GetEventFolder()->FindObject(name);
1738 AliDebug(1, Form("name=%s, cleaning %s.",GetName(),name.Data()));
1739 GetEventFolder()->Remove(obj);
1744 /*****************************************************************************/
1746 TTask* AliRunLoader::GetRunDigitizer()
1748 //returns Run Digitizer from folder
1750 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1751 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetDigitizerTaskName());
1752 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1754 /*****************************************************************************/
1756 TTask* AliRunLoader::GetRunSDigitizer()
1758 //returns SDigitizer Task from folder
1760 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1761 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetSDigitizerTaskName());
1762 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1764 /*****************************************************************************/
1766 TTask* AliRunLoader::GetRunReconstructioner()
1768 //returns Reconstructioner Task from folder
1769 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1770 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetReconstructionerTaskName());
1771 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1773 /*****************************************************************************/
1775 TTask* AliRunLoader::GetRunTracker()
1777 //returns Tracker Task from folder
1778 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1779 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetTrackerTaskName());
1780 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1782 /*****************************************************************************/
1784 TTask* AliRunLoader::GetRunPIDTask()
1786 //returns Tracker Task from folder
1787 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1788 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetPIDTaskName());
1789 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1791 /*****************************************************************************/
1793 TTask* AliRunLoader::GetRunQATask()
1795 //returns Quality Assurance Task from folder
1796 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1799 AliErrorClass("Can not get task folder from AliConfig");
1802 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetQATaskName());
1803 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1806 /*****************************************************************************/
1808 void AliRunLoader::SetCompressionLevel(Int_t cl)
1810 //Sets Compression Level in all files
1811 if (fGAFile) fGAFile->SetCompressionLevel(cl);
1812 SetKineComprLevel(cl);
1813 SetTrackRefsComprLevel(cl);
1814 TIter next(fLoaders);
1816 while((loader = (AliLoader*)next()))
1818 loader->SetCompressionLevel(cl);
1821 /**************************************************************************/
1823 void AliRunLoader::SetKineComprLevel(Int_t cl)
1825 //Sets comression level in Kine File
1826 fKineDataLoader->SetCompressionLevel(cl);
1828 /**************************************************************************/
1830 void AliRunLoader::SetTrackRefsComprLevel(Int_t cl)
1832 //Sets comression level in Track Refences File
1833 fTrackRefsDataLoader->SetCompressionLevel(cl);
1835 /**************************************************************************/
1837 void AliRunLoader::UnloadHeader()
1839 //removes TreeE from folder and deletes it
1840 // as well as fHeader object
1845 /**************************************************************************/
1847 void AliRunLoader::UnloadKinematics()
1849 //Unloads Kinematics
1850 fKineDataLoader->GetBaseLoader(0)->Unload();
1852 /**************************************************************************/
1854 void AliRunLoader::UnloadTrackRefs()
1856 //Unloads Track Refernces
1857 fTrackRefsDataLoader->GetBaseLoader(0)->Unload();
1859 /**************************************************************************/
1861 void AliRunLoader::UnloadgAlice()
1864 if (gAlice == GetAliRun())
1866 AliDebug(1, "Set gAlice = 0x0");
1867 gAlice = 0x0;//if gAlice is the same that in folder (to be deleted by the way of folder)
1869 AliRun* alirun = GetAliRun();
1870 if (GetEventFolder()) GetEventFolder()->Remove(alirun);
1873 /**************************************************************************/
1875 void AliRunLoader::MakeTrackRefsContainer()
1877 // Makes a tree for Track References
1878 fTrackRefsDataLoader->MakeTree();
1880 /**************************************************************************/
1882 Int_t AliRunLoader::LoadTrackRefs(Option_t* option)
1884 //Load track references from file (opens file and posts tree to folder)
1886 return fTrackRefsDataLoader->GetBaseLoader(0)->Load(option);
1888 /**************************************************************************/
1890 void AliRunLoader::SetDirName(TString& dirname)
1892 //sets directory name
1893 if (dirname.IsNull()) return;
1894 fUnixDirName = dirname;
1895 fKineDataLoader->SetDirName(dirname);
1896 fTrackRefsDataLoader->SetDirName(dirname);
1898 TIter next(fLoaders);
1900 while((loader = (AliLoader*)next()))
1902 loader->SetDirName(dirname);
1906 /*****************************************************************************/
1908 Int_t AliRunLoader::GetFileOffset() const
1910 //returns the file number that is added to the file name for current event
1911 return Int_t(fCurrentEvent/fNEventsPerFile);
1914 /*****************************************************************************/
1915 const TString AliRunLoader::SetFileOffset(const TString& fname)
1917 //adds the the number to the file name at proper place for current event
1918 Long_t offset = (Long_t)GetFileOffset();
1919 if (offset < 1) return fname;
1921 soffset += offset;//automatic conversion to string
1922 TString dotroot(".root");
1923 const TString& offfsetdotroot = offset + dotroot;
1924 TString out = fname;
1925 out = out.ReplaceAll(dotroot,offfsetdotroot);
1926 AliDebug(1, Form(" in=%s out=%s",fname.Data(),out.Data()));
1929 /*****************************************************************************/
1931 void AliRunLoader::SetDigitsFileNameSuffix(const TString& suffix)
1933 //adds the suffix before ".root",
1934 //e.g. TPC.Digits.root -> TPC.DigitsMerged.root
1935 //made on Jiri Chudoba demand
1937 TIter next(fLoaders);
1939 while((loader = (AliLoader*)next()))
1941 loader->SetDigitsFileNameSuffix(suffix);
1944 /*****************************************************************************/
1946 TString AliRunLoader::GetFileName() const
1948 //returns name of galice file
1950 if (fGAFile == 0x0) return result;
1951 result = fGAFile->GetName();
1954 /*****************************************************************************/
1956 void AliRunLoader::SetDetectorAddresses()
1958 //calls SetTreeAddress for all detectors
1959 if (GetAliRun()==0x0) return;
1960 TIter next(GetAliRun()->Modules());
1962 while((mod = (AliModule*)next()))
1964 AliDetector* det = dynamic_cast<AliDetector*>(mod);
1965 if (det) det->SetTreeAddress();
1968 /*****************************************************************************/
1970 void AliRunLoader::Synchronize()
1972 //synchrinizes all writtable files
1973 TIter next(fLoaders);
1975 while((loader = (AliLoader*)next()))
1977 loader->Synchronize();
1980 fKineDataLoader->Synchronize();
1981 fTrackRefsDataLoader->Synchronize();
1983 if (fGAFile) fGAFile->Flush();
1985 /*****************************************************************************/
1986 /*****************************************************************************/