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 <TObjArray.h>
59 #include "AliConfig.h"
60 #include "AliLoader.h"
61 #include "AliHeader.h"
63 #include "AliDetector.h"
64 #include "AliCDBManager.h"
65 #include "AliCDBLocal.h"
66 #include "AliCentralTrigger.h"
68 ClassImp(AliRunLoader)
70 AliRunLoader* AliRunLoader::fgRunLoader = 0x0;
72 const TString AliRunLoader::fgkRunLoaderName("RunLoader");
73 const TString AliRunLoader::fgkHeaderBranchName("Header");
74 const TString AliRunLoader::fgkTriggerBranchName("ClassMask");
75 const TString AliRunLoader::fgkHeaderContainerName("TE");
76 const TString AliRunLoader::fgkTriggerContainerName("TreeCT");
77 const TString AliRunLoader::fgkKineContainerName("TreeK");
78 const TString AliRunLoader::fgkTrackRefsContainerName("TreeTR");
79 const TString AliRunLoader::fgkKineBranchName("Particles");
80 const TString AliRunLoader::fgkDefaultKineFileName("Kinematics.root");
81 const TString AliRunLoader::fgkDefaultTrackRefsFileName("TrackRefs.root");
82 const TString AliRunLoader::fgkGAliceName("gAlice");
83 const TString AliRunLoader::fgkDefaultTriggerFileName("Trigger.root");
84 /**************************************************************************/
86 AliRunLoader::AliRunLoader():
96 fTrackRefsDataLoader(0x0),
100 AliConfig::Instance();//force to build the folder structure
101 if (!fgRunLoader) fgRunLoader = this;
103 /**************************************************************************/
105 AliRunLoader::AliRunLoader(const char* eventfoldername):
106 TNamed(fgkRunLoaderName,fgkRunLoaderName),
107 fLoaders(new TObjArray()),
115 fKineDataLoader(new AliDataLoader(fgkDefaultKineFileName,fgkKineContainerName,"Kinematics")),
116 fTrackRefsDataLoader(new AliDataLoader(fgkDefaultTrackRefsFileName,fgkTrackRefsContainerName,"Track References")),
121 SetEventFolderName(eventfoldername);
122 if (!fgRunLoader) fgRunLoader = this;
124 /**************************************************************************/
126 AliRunLoader::~AliRunLoader()
129 if (fgRunLoader == this) fgRunLoader = 0x0;
135 fLoaders->SetOwner();
139 delete fKineDataLoader;
140 delete fTrackRefsDataLoader;
145 //fEventFolder is deleted by the way of removing - TopAliceFolder owns it
146 if( fCTrigger ) delete fCTrigger;
151 /**************************************************************************/
153 AliRunLoader::AliRunLoader(TFolder* topfolder):
154 TNamed(fgkRunLoaderName,fgkRunLoaderName),
155 fLoaders(new TObjArray()),
156 fEventFolder(topfolder),
163 fKineDataLoader(new AliDataLoader(fgkDefaultKineFileName,fgkKineContainerName,"Kinematics")),
164 fTrackRefsDataLoader(new AliDataLoader(fgkDefaultTrackRefsFileName,fgkTrackRefsContainerName,"Track References")),
171 TString errmsg("Parameter is NULL");
172 AliError(errmsg.Data());
177 TObject* obj = fEventFolder->FindObject(fgkRunLoaderName);
179 { //if it is, then sth. is going wrong... exits aliroot session
180 TString errmsg("In Event Folder Named ");
181 errmsg+=fEventFolder->GetName();
182 errmsg+=" object named "+fgkRunLoaderName+" already exists. I am confused ...";
184 AliError(errmsg.Data());
186 return;//never reached
189 if (!fgRunLoader) fgRunLoader = this;
191 fEventFolder->Add(this);//put myself to the folder to accessible for all
195 /**************************************************************************/
197 Int_t AliRunLoader::GetEvent(Int_t evno)
199 //Gets event number evno
200 //Reloads all data properly
201 //PH if (fCurrentEvent == evno) return 0;
205 AliError("Can not give the event with negative number");
209 if (evno >= GetNumberOfEvents())
211 AliError(Form("There is no event with number %d",evno));
215 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
216 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
217 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
218 AliDebug(1, Form(" GETTING EVENT %d",evno));
219 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
220 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
221 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
223 fCurrentEvent = evno;
227 //Reload header (If header was loaded)
230 retval = TreeE()->GetEvent(fCurrentEvent);
233 AliError(Form("Cannot find event: %d\n ",fCurrentEvent));
237 //Reload stack (If header was loaded)
238 if (TreeE()) fStack = GetHeader()->Stack();
239 //Set event folder in stack (it does not mean that we read kinematics from file)
240 if( GetTrigger() && TreeCT() ) {
241 retval = TreeCT()->GetEvent(fCurrentEvent);
243 AliError(Form("Error occured while GetEvent for Trigger. Event %d",evno));
251 AliError(Form("Error occured while setting event %d",evno));
255 //Post Track References
256 retval = fTrackRefsDataLoader->GetEvent();
259 AliError(Form("Error occured while GetEvent for Track References. Event %d",evno));
263 //Read Kinematics if loaded
264 retval = fKineDataLoader->GetEvent();
267 AliError(Form("Error occured while GetEvent for Kinematics. Event %d",evno));
271 if (fStack && fKineDataLoader->GetBaseLoader(0)->IsLoaded())
273 fStack->ConnectTree(TreeK());
275 if (fStack->GetEvent() == kFALSE)
277 AliError(Form("Error occured while GetEvent for Stack. Event %d",evno));
282 //Trigger data reloading in all loaders
283 TIter next(fLoaders);
285 while((loader = (AliLoader*)next()))
287 retval = loader->GetEvent();
290 AliError(Form("Error occured while getting event for %s. Event %d.",
291 loader->GetDetectorName().Data(), evno));
296 SetDetectorAddresses();
300 /**************************************************************************/
301 Int_t AliRunLoader::SetEvent()
303 //if kinematics was loaded Cleans folder data
307 retval = fKineDataLoader->SetEvent();
310 AliError("SetEvent for Kinamtics Data Loader retutned error.");
313 retval = fTrackRefsDataLoader->SetEvent();
316 AliError("SetEvent for Track References Data Loader retutned error.");
320 TIter next(fLoaders);
322 while((loader = (AliLoader*)next()))
324 retval = loader->SetEvent();
327 AliError(Form("SetEvent for %s Data Loader retutned error.",loader->GetName()));
334 /**************************************************************************/
336 Int_t AliRunLoader::SetEventNumber(Int_t evno)
338 //cleans folders and sets the root dirs in files
339 if (fCurrentEvent == evno) return 0;
340 fCurrentEvent = evno;
344 /**************************************************************************/
345 AliCDBEntry* AliRunLoader::GetCDBEntry(const char* name) const
347 //Get an AliCDBEntry from the run data storage
349 if ( !(AliCDBManager::Instance()->IsDefaultStorageSet()) ) {
350 AliError("No run data storage defined!");
353 return AliCDBManager::Instance()->GetDefaultStorage()->Get(name, GetHeader()->GetRun());
357 /**************************************************************************/
358 AliRunLoader* AliRunLoader::Open
359 (const char* filename, const char* eventfoldername, Option_t* option)
361 //Opens a desired file 'filename'
362 //gets the the run-Loader and mounts it desired folder
363 //returns the pointer to run Loader which can be further used for accessing data
364 //in case of error returns NULL
366 static const TString kwebaddress("http://alisoft.cern.ch/people/skowron/codedoc/split/index.html");
367 AliDebugClass(1,Form("\n\n\nNew I/O strcture: See more info:\n %s\n\n\n",kwebaddress.Data()));
369 AliRunLoader* result = 0x0;
371 /* ************************************************ */
372 /* Chceck if folder with given name already exists */
373 /* ************************************************ */
375 TObject* obj = AliConfig::Instance()->GetTopFolder()->FindObject(eventfoldername);
378 TFolder* fold = dynamic_cast<TFolder*>(obj);
381 AliErrorClass("Such a obejct already exists in top alice folder and it is not a folder.");
385 //check if we can get RL from that folder
386 result = AliRunLoader::GetRunLoader(eventfoldername);
389 AliErrorClass(Form("Folder %s already exists, and can not find session there. Can not mount.",eventfoldername));
393 if (result->GetFileName().CompareTo(filename) != 0)
395 AliErrorClass("Other file is mounted in demanded folder. Can not mount.");
399 //check if now is demanded (re)creation
400 if ( AliLoader::TestFileOption(option) == kFALSE)
402 AliErrorClass(Form("Session already exists in folder %s and this session option is %s. Unable to proceed.",
403 eventfoldername,option));
407 //check if demanded option is update and existing one
408 TString tmpstr(option);
409 if ( (tmpstr.CompareTo("update",TString::kIgnoreCase) == 0) &&
410 (result->fGAFile->IsWritable() == kFALSE) )
412 AliErrorClass(Form("Session already exists in folder %s and is not writable while this session option is %s. Unable to proceed.",
413 eventfoldername,option));
417 AliWarningClass("Session is already opened and mounted in demanded folder");
418 if (!fgRunLoader) fgRunLoader = result; //PH get access from any place
420 } //end of checking in case of existance of object named identically that folder session is being opened
423 TFile * gAliceFile = TFile::Open(filename,option);//open a file
425 {//null pointer returned
426 AliFatalClass(Form("Can not open file %s.",filename));
430 if (gAliceFile->IsOpen() == kFALSE)
431 {//pointer to valid object returned but file is not opened
432 AliErrorClass(Form("Can not open file %s.",filename));
436 //if file is "read" or "update" than we try to find AliRunLoader there - if not found cry and exit
437 //else create new AliRunLoader
438 if ( AliLoader::TestFileOption(option) )
440 AliDebugClass(1, "Reading RL from file");
442 result = dynamic_cast<AliRunLoader*>(gAliceFile->Get(fgkRunLoaderName));//get the run Loader from the file
445 AliErrorClass(Form("Can not find run-Loader in file %s.",filename));
446 delete gAliceFile;//close the file
449 Int_t tmp = result->SetEventFolderName(eventfoldername);//mount a event folder
450 if (tmp)//if SetEvent returned error
452 AliErrorClass(Form("Can not mount event in folder %s.",eventfoldername));
453 delete result; //delete run-Loader
454 delete gAliceFile;//close the file
460 AliDebugClass(1, Form("Creating new AliRunLoader. Folder name is %s",eventfoldername));
463 result = new AliRunLoader(eventfoldername);
465 catch (TString& errmsg)
467 AliErrorClass(Form("AliRunLoader constrcutor has thrown exception: %s\n",errmsg.Data()));
469 delete gAliceFile;//close the file
474 //procedure for extracting dir name from the file name
475 TString fname(filename);
476 Int_t nsl = fname.Last('#');//look for hash in file name
478 if (nsl < 0) {//hash not found
479 nsl = fname.Last('/');// look for slash
481 nsl = fname.Last(':');// look for colon e.g. rfio:galice.root
484 if (nsl < 0) dirname = "./"; // no directory path, use "."
485 else dirname = fname.Remove(nsl+1);// directory path
487 AliDebugClass(1, Form("Dir name is : %s",dirname.Data()));
489 result->SetDirName(dirname);
490 result->SetGAliceFile(gAliceFile);//set the pointer to gAliceFile
491 if (!fgRunLoader) fgRunLoader = result; //PH get access from any place
494 /**************************************************************************/
495 Int_t AliRunLoader::GetNumberOfEvents()
497 //returns number of events in Run
501 retval = LoadHeader();
504 AliError("Error occured while loading header");
508 return (Int_t)TreeE()->GetEntries();
510 /**************************************************************************/
512 void AliRunLoader::MakeHeader()
514 //Makes header and connects it to header tree (if it exists)
518 AliDebug(1, "Creating new Header Object");
519 fHeader= new AliHeader();
521 TTree* tree = TreeE();
524 AliDebug(1, "Got Tree from folder.");
525 TBranch* branch = tree->GetBranch(fgkHeaderBranchName);
528 AliDebug(1, "Creating new branch");
529 branch = tree->Branch(fgkHeaderBranchName, "AliHeader", &fHeader, 4000, 0);
530 branch->SetAutoDelete(kFALSE);
534 AliDebug(1, "Got Branch from Tree");
535 branch->SetAddress(&fHeader);
536 tree->GetEvent(fCurrentEvent);
537 fStack = fHeader->Stack(); //should be safe - if we created Stack, header returns pointer to the same object
541 fStack->ConnectTree(TreeK());
547 AliDebug(1, "Header does not have a stack.");
551 AliDebug(1, "Exiting MakeHeader method");
553 /**************************************************************************/
555 void AliRunLoader::MakeStack()
557 //Creates the stack object - do not connect the tree
560 fStack = new AliStack(10000);
563 /**************************************************************************/
565 void AliRunLoader::MakeTrigger()
567 // Makes trigger object and connects it to trigger tree (if it exists)
569 if( fCTrigger == 0x0 ) {
570 AliDebug( 1, "Creating new Trigger Object" );
571 fCTrigger = new AliCentralTrigger();
573 TTree* tree = TreeCT();
575 fCTrigger->MakeBranch( fgkTriggerBranchName, tree );
576 tree->GetEvent( fCurrentEvent );
579 AliDebug( 1, "Exiting MakeTrigger method" );
581 /**************************************************************************/
583 void AliRunLoader::MakeTree(Option_t *option)
586 const char *oK = strstr(option,"K"); //Kine
587 const char *oE = strstr(option,"E"); //Header
588 const char *oGG = strstr(option,"GG"); //Central TriGGer
592 if (fKineDataLoader->GetBaseLoader(0)->IsLoaded() == kFALSE)
594 AliError("Load Kinematics first");
599 fKineDataLoader->MakeTree();
602 fStack->ConnectTree(TreeK());
603 WriteKinematics("OVERWRITE");
610 TTree* tree = new TTree(fgkHeaderContainerName,"Tree with Headers");
611 GetEventFolder()->Add(tree);
613 WriteHeader("OVERWRITE");
618 // create the CTP Trigger output file and tree
619 TFile* file = gROOT->GetFile( fgkDefaultTriggerFileName );
621 file = TFile::Open( gSystem->ConcatFileName( fUnixDirName.Data(), fgkDefaultTriggerFileName.Data() ), "RECREATE" ) ;
625 TTree* tree = new TTree( fgkTriggerContainerName, "Tree with Central Trigger Mask" );
626 GetEventFolder()->Add(tree);
628 // WriteHeader("OVERWRITE");
631 TIter next(fLoaders);
633 while((loader = (AliLoader*)next()))
635 loader->MakeTree(option);
639 /**************************************************************************/
641 Int_t AliRunLoader::LoadgAlice()
643 //Loads gAlice from file
646 AliWarning("AliRun is already in folder. Unload first.");
649 AliRun* alirun = dynamic_cast<AliRun*>(fGAFile->Get(fgkGAliceName));
652 AliError(Form("Can not find gAlice in file %s",fGAFile->GetName()));
655 alirun->SetRunLoader(this);
658 AliWarning(Form("gAlice already exists. Putting retrived object in folder named %s",
659 GetEventFolder()->GetName()));
665 SetDetectorAddresses();//calls SetTreeAddress for all detectors
668 /**************************************************************************/
670 Int_t AliRunLoader::LoadHeader()
672 //loads treeE and reads header object for current event
675 AliWarning("Header is already loaded. Use ReloadHeader to force reload. Nothing done");
679 if (GetEventFolder() == 0x0)
681 AliError("Event folder not specified yet");
687 AliError("Session not opened. Use AliRunLoader::Open");
691 if (fGAFile->IsOpen() == kFALSE)
693 AliError("Session not opened. Use AliRunLoader::Open");
697 TTree* tree = dynamic_cast<TTree*>(fGAFile->Get(fgkHeaderContainerName));
700 AliError(Form("Can not find header tree named %s in file %s",
701 fgkHeaderContainerName.Data(),fGAFile->GetName()));
705 if (tree == TreeE()) return 0;
708 GetEventFolder()->Add(tree);
709 MakeHeader();//creates header object and connects to tree
713 /**************************************************************************/
715 Int_t AliRunLoader::LoadTrigger(Option_t* option)
720 AliWarning("Trigger is already loaded. Nothing done");
724 if( GetEventFolder() == 0x0 ) {
725 AliError("Event folder not specified yet");
728 // get the CTP Trigger output file and tree
729 TString trgfile = gSystem->ConcatFileName( fUnixDirName.Data(),
730 fgkDefaultTriggerFileName.Data() );
731 TFile* file = gROOT->GetFile( trgfile );
733 file = TFile::Open( trgfile, option ) ;
734 if (!file || file->IsOpen() == kFALSE ) {
735 AliError( Form( "Can not open trigger file %s", trgfile.Data() ) );
741 TTree* tree = dynamic_cast<TTree*>(file->Get( fgkTriggerContainerName ));
743 AliError( Form( "Can not find trigger tree named %s in file %s",
744 fgkTriggerContainerName.Data(), file->GetName() ) );
750 fCTrigger = dynamic_cast<AliCentralTrigger*>(file->Get( "AliCentralTrigger" ));
751 GetEventFolder()->Add( tree );
757 /**************************************************************************/
759 Int_t AliRunLoader::LoadKinematics(Option_t* option)
761 //Loads the kinematics
762 Int_t retval = fKineDataLoader->GetBaseLoader(0)->Load(option);
765 AliError("Error occured while loading kinamatics tree.");
770 fStack->ConnectTree(TreeK());
771 retval = fStack->GetEvent();
772 if ( retval == kFALSE)
774 AliError("Error occured while loading kinamatics tree.");
781 /**************************************************************************/
783 Int_t AliRunLoader::OpenDataFile(const TString& filename,TFile*& file,TDirectory*& dir,Option_t* opt,Int_t cl)
785 //Opens File with kinematics
788 if (file->IsOpen() == kFALSE)
789 {//pointer is not null but file is not opened
790 AliWarning("Pointer to file is not null, but file is not opened");//risky any way
792 file = 0x0; //proceed with opening procedure
796 AliWarning(Form("File %s already opened",filename.Data()));
800 //try to find if that file is opened somewere else
801 file = (TFile *)( gROOT->GetListOfFiles()->FindObject(filename) );
804 if(file->IsOpen() == kTRUE)
806 AliWarning(Form("File %s already opened by sombody else.",file->GetName()));
811 file = TFile::Open(filename,opt);
814 AliError(Form("Can not open file %s",filename.Data()));
817 if (file->IsOpen() == kFALSE)
818 {//file is not opened
819 AliError(Form("Can not open file %s",filename.Data()));
823 file->SetCompressionLevel(cl);
825 dir = AliLoader::ChangeDir(file,fCurrentEvent);
828 AliError(Form("Can not change to root directory in file %s",filename.Data()));
833 /**************************************************************************/
835 TTree* AliRunLoader::TreeE() const
837 //returns the tree from folder; shortcut method
838 if (AliDebugLevel() > 10) fEventFolder->ls();
839 TObject *obj = fEventFolder->FindObject(fgkHeaderContainerName);
840 return (obj)?dynamic_cast<TTree*>(obj):0x0;
842 /**************************************************************************/
844 TTree* AliRunLoader::TreeCT() const
846 //returns the tree from folder; shortcut method
847 if (AliDebugLevel() > 10) fEventFolder->ls();
848 TObject *obj = fEventFolder->FindObject(fgkTriggerContainerName);
849 return (obj)?dynamic_cast<TTree*>(obj):0x0;
851 /**************************************************************************/
853 AliHeader* AliRunLoader::GetHeader() const
855 //returns pointer header object
858 /**************************************************************************/
860 AliCentralTrigger* AliRunLoader::GetTrigger() const
862 //returns pointer trigger object
866 /**************************************************************************/
868 TTree* AliRunLoader::TreeK() const
870 //returns the tree from folder; shortcut method
871 TObject *obj = GetEventFolder()->FindObject(fgkKineContainerName);
872 return (obj)?dynamic_cast<TTree*>(obj):0x0;
874 /**************************************************************************/
876 TTree* AliRunLoader::TreeTR() const
878 //returns the tree from folder; shortcut method
879 TObject* obj = GetEventFolder()->FindObject(fgkTrackRefsContainerName);
880 return (obj)?dynamic_cast<TTree*>(obj):0x0;
882 /**************************************************************************/
884 AliRun* AliRunLoader::GetAliRun() const
886 //returns AliRun which sits in the folder
887 if (fEventFolder == 0x0) return 0x0;
888 TObject *obj = fEventFolder->FindObject(fgkGAliceName);
889 return (obj)?dynamic_cast<AliRun*>(obj):0x0;
891 /**************************************************************************/
893 Int_t AliRunLoader::WriteHeader(Option_t* opt)
896 AliDebug(1, "WRITING HEADER");
898 TTree* tree = TreeE();
901 AliWarning("Can not find Header Tree in Folder");
904 if (fGAFile->IsWritable() == kFALSE)
906 AliError(Form("File %s is not writable",fGAFile->GetName()));
910 TObject* obj = fGAFile->Get(fgkHeaderContainerName);
912 { //if they exist, see if option OVERWRITE is used
914 if(tmp.Contains("OVERWRITE",TString::kIgnoreCase) == 0)
915 {//if it is not used - give an error message and return an error code
916 AliError("Tree already exisists. Use option \"OVERWRITE\" to overwrite previous data");
921 tree->SetDirectory(fGAFile);
922 tree->Write(0,TObject::kOverwrite);
924 AliDebug(1, "WRITTEN\n\n");
929 /**************************************************************************/
931 Int_t AliRunLoader::WriteTrigger(Option_t* opt)
934 AliDebug( 1, "WRITING TRIGGER" );
936 TTree* tree = TreeCT();
938 AliWarning("Can not find Trigger Tree in Folder");
942 TFile* file = gROOT->GetFile( gSystem->ConcatFileName( fUnixDirName.Data(), fgkDefaultTriggerFileName.Data() ) ) ;
943 if( !file || !file->IsOpen() ) {
944 AliError( "can't write Trigger, file is not open" );
948 TObject* obj = file->Get( fgkTriggerContainerName );
949 if( obj ) { //if they exist, see if option OVERWRITE is used
951 if( tmp.Contains( "OVERWRITE", TString::kIgnoreCase ) == 0) {
952 //if it is not used - give an error message and return an error code
953 AliError( "Tree already exisists. Use option \"OVERWRITE\" to overwrite previous data" );
958 fCTrigger->Write( 0, TObject::kOverwrite );
959 tree->Write( 0, TObject::kOverwrite );
962 AliDebug(1, "WRITTEN\n\n");
966 /**************************************************************************/
968 Int_t AliRunLoader::WriteAliRun(Option_t* /*opt*/)
970 //writes AliRun object to the file
972 if (GetAliRun()) GetAliRun()->Write();
975 /**************************************************************************/
977 Int_t AliRunLoader::WriteKinematics(Option_t* opt)
980 return fKineDataLoader->GetBaseLoader(0)->WriteData(opt);
982 /**************************************************************************/
983 Int_t AliRunLoader::WriteTrackRefs(Option_t* opt)
985 //writes Track References tree
986 return fTrackRefsDataLoader->GetBaseLoader(0)->WriteData(opt);
988 /**************************************************************************/
990 Int_t AliRunLoader::WriteHits(Option_t* opt)
992 //Calls WriteHits for all loaders
995 TIter next(fLoaders);
997 while((loader = (AliLoader*)next()))
999 res = loader->WriteHits(opt);
1002 AliError(Form("Failed to write hits for %s (%d)",loader->GetDetectorName().Data(),res));
1008 /**************************************************************************/
1010 Int_t AliRunLoader::WriteSDigits(Option_t* opt)
1012 //Calls WriteSDigits for all loaders
1015 TIter next(fLoaders);
1017 while((loader = (AliLoader*)next()))
1019 res = loader->WriteSDigits(opt);
1022 AliError(Form("Failed to write summable digits for %s.",loader->GetDetectorName().Data()));
1028 /**************************************************************************/
1030 Int_t AliRunLoader::WriteDigits(Option_t* opt)
1032 //Calls WriteDigits for all loaders
1035 TIter next(fLoaders);
1037 while((loader = (AliLoader*)next()))
1039 res = loader->WriteDigits(opt);
1042 AliError(Form("Failed to write digits for %s.",loader->GetDetectorName().Data()));
1048 /**************************************************************************/
1050 Int_t AliRunLoader::WriteRecPoints(Option_t* opt)
1052 //Calls WriteRecPoints for all loaders
1055 TIter next(fLoaders);
1057 while((loader = (AliLoader*)next()))
1059 res = loader->WriteRecPoints(opt);
1062 AliError(Form("Failed to write Reconstructed Points for %s.",
1063 loader->GetDetectorName().Data()));
1069 /**************************************************************************/
1071 Int_t AliRunLoader::WriteTracks(Option_t* opt)
1073 //Calls WriteTracks for all loaders
1076 TIter next(fLoaders);
1078 while((loader = (AliLoader*)next()))
1080 res = loader->WriteTracks(opt);
1083 AliError(Form("Failed to write Tracks for %s.",
1084 loader->GetDetectorName().Data()));
1090 /**************************************************************************/
1092 Int_t AliRunLoader::WriteRunLoader(Option_t* /*opt*/)
1094 //Writes itself to the file
1096 this->Write(0,TObject::kOverwrite);
1099 /**************************************************************************/
1101 Int_t AliRunLoader::SetEventFolderName(const TString& name)
1103 //sets top folder name for this run; of alread
1106 AliError("Name is empty");
1110 //check if such a folder already exists - try to find it in alice top folder
1111 TObject* obj = AliConfig::Instance()->GetTopFolder()->FindObject(name);
1114 TFolder* fold = dynamic_cast<TFolder*>(obj);
1117 AliError("Such a obejct already exists in top alice folder and it is not a folder.");
1120 //folder which was found is our folder
1121 if (fEventFolder == fold)
1127 AliError("Such a folder already exists in top alice folder. Can not mount.");
1132 //event is alredy connected, just change name of the folder
1135 fEventFolder->SetName(name);
1139 if (fKineDataLoader == 0x0)
1140 fKineDataLoader = new AliDataLoader(fgkDefaultKineFileName,fgkKineContainerName,"Kinematics");
1142 if ( fTrackRefsDataLoader == 0x0)
1143 fTrackRefsDataLoader = new AliDataLoader(fgkDefaultTrackRefsFileName,fgkTrackRefsContainerName,"Track References");
1145 //build the event folder structure
1146 AliDebug(1, Form("Creating new event folder named %s",name.Data()));
1147 fEventFolder = AliConfig::Instance()->BuildEventFolder(name,"Event Folder");
1148 fEventFolder->Add(this);//put myself to the folder to accessible for all
1150 TIter next(fLoaders);
1152 while((loader = (AliLoader*)next()))
1154 loader->Register(fEventFolder);//build folder structure for this detector
1157 fKineDataLoader->SetEventFolder(GetEventFolder());
1158 fTrackRefsDataLoader->SetEventFolder(GetEventFolder());
1159 fKineDataLoader->SetFolder(GetEventFolder());
1160 fTrackRefsDataLoader->SetFolder(GetEventFolder());
1162 fEventFolder->SetOwner();
1165 /**************************************************************************/
1167 void AliRunLoader::AddLoader(AliLoader* loader)
1169 //Adds the Loader for given detector
1170 if (loader == 0x0) //if null shout and exit
1172 AliError("Parameter is NULL");
1175 loader->SetDirName(fUnixDirName);
1176 if (fEventFolder) loader->SetEventFolder(fEventFolder); //if event folder is already defined,
1177 //pass information to the Loader
1178 fLoaders->Add(loader);//add the Loader to the array
1180 /**************************************************************************/
1182 void AliRunLoader::AddLoader(AliDetector* det)
1184 //Asks module (detector) ro make a Loader and stores in the array
1185 if (det == 0x0) return;
1186 AliLoader* get = det->GetLoader();//try to get loader
1187 if (get == 0x0) get = det->MakeLoader(fEventFolder->GetName());//if did not obtain, ask to make it
1191 AliDebug(1, Form("Detector: %s Loader : %s",det->GetName(),get->GetName()));
1196 /**************************************************************************/
1198 AliLoader* AliRunLoader::GetLoader(const char* detname) const
1200 //returns loader for given detector
1201 //note that naming convention is TPCLoader not just TPC
1202 return (AliLoader*)fLoaders->FindObject(detname);
1205 /**************************************************************************/
1207 AliLoader* AliRunLoader::GetLoader(AliDetector* det) const
1209 //get loader for detector det
1210 if(det == 0x0) return 0x0;
1211 TString getname(det->GetName());
1213 AliDebug(1, Form(" Loader name is %s",getname.Data()));
1214 return GetLoader(getname);
1217 /**************************************************************************/
1219 void AliRunLoader::CleanFolders()
1221 // fEventFolder->Add(this);//put myself to the folder to accessible for all
1228 /**************************************************************************/
1230 void AliRunLoader::CleanDetectors()
1232 //Calls CleanFolders for all detectors
1233 TIter next(fLoaders);
1235 while((loader = (AliLoader*)next()))
1237 loader->CleanFolders();
1240 /**************************************************************************/
1242 void AliRunLoader::RemoveEventFolder()
1244 //remove all the tree of event
1245 //all the stuff changing EbE stays untached (PDGDB, tasks, etc.)
1247 if (fEventFolder == 0x0) return;
1248 fEventFolder->SetOwner(kFALSE);//don't we want to deleted while removing the folder that we are sitting in
1249 fEventFolder->Remove(this);//remove us drom folder
1251 AliConfig::Instance()->GetTopFolder()->SetOwner(); //brings ownership back for fEventFolder since it sits in top folder
1252 AliConfig::Instance()->GetTopFolder()->Remove(fEventFolder); //remove the event tree
1253 delete fEventFolder;
1255 /**************************************************************************/
1257 void AliRunLoader::SetGAliceFile(TFile* gafile)
1259 //sets pointer to galice.root file
1263 /**************************************************************************/
1265 Int_t AliRunLoader::LoadHits(Option_t* detectors,Option_t* opt)
1267 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1269 AliDebug(1, "Loading Hits");
1273 const char* oAll = strstr(detectors,"all");
1276 AliDebug(1, "Option is All");
1281 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1282 loaders = &arr;//get the pointer array
1285 AliDebug(1, Form("For detectors. Number of detectors chosen for loading %d",loaders->GetEntries()));
1287 TIter next(loaders);
1289 while((loader = (AliLoader*)next()))
1291 AliDebug(1, Form(" Calling LoadHits(%s) for %s",opt,loader->GetName()));
1292 loader->LoadHits(opt);
1294 AliDebug(1, "Done");
1298 /**************************************************************************/
1300 Int_t AliRunLoader::LoadSDigits(Option_t* detectors,Option_t* opt)
1302 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1307 const char* oAll = strstr(detectors,"all");
1314 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1315 loaders = &arr;//get the pointer to array
1318 TIter next(loaders);
1320 while((loader = (AliLoader*)next()))
1322 loader->LoadSDigits(opt);
1327 /**************************************************************************/
1329 Int_t AliRunLoader::LoadDigits(Option_t* detectors,Option_t* opt)
1331 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1336 const char* oAll = strstr(detectors,"all");
1343 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1344 loaders = &arr;//get the pointer array
1347 TIter next(loaders);
1349 while((loader = (AliLoader*)next()))
1351 loader->LoadDigits(opt);
1355 /**************************************************************************/
1357 Int_t AliRunLoader::LoadRecPoints(Option_t* detectors,Option_t* opt)
1359 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1364 const char* oAll = strstr(detectors,"all");
1371 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1372 loaders = &arr;//get the pointer array
1375 TIter next(loaders);
1377 while((loader = (AliLoader*)next()))
1379 loader->LoadRecPoints(opt);
1383 /**************************************************************************/
1385 Int_t AliRunLoader::LoadRecParticles(Option_t* detectors,Option_t* opt)
1387 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1392 const char* oAll = strstr(detectors,"all");
1399 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1400 loaders = &arr;//get the pointer array
1403 TIter next(loaders);
1405 while((loader = (AliLoader*)next()))
1407 loader->LoadRecParticles(opt);
1411 /**************************************************************************/
1413 Int_t AliRunLoader::LoadTracks(Option_t* detectors,Option_t* opt)
1415 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1420 const char* oAll = strstr(detectors,"all");
1427 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1428 loaders = &arr;//get the pointer array
1431 TIter next(loaders);
1433 while((loader = (AliLoader*)next()))
1435 loader->LoadTracks(opt);
1439 /**************************************************************************/
1441 void AliRunLoader::UnloadHits(Option_t* detectors)
1443 //unloads hits for detectors specified in parameter
1447 const char* oAll = strstr(detectors,"all");
1454 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1455 loaders = &arr;//get the pointer to array
1458 TIter next(loaders);
1460 while((loader = (AliLoader*)next()))
1462 loader->UnloadHits();
1465 /**************************************************************************/
1467 void AliRunLoader::UnloadSDigits(Option_t* detectors)
1469 //unloads SDigits for detectors specified in parameter
1473 const char* oAll = strstr(detectors,"all");
1480 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1481 loaders = &arr;//get the pointer to array
1484 TIter next(loaders);
1486 while((loader = (AliLoader*)next()))
1488 loader->UnloadSDigits();
1491 /**************************************************************************/
1493 void AliRunLoader::UnloadDigits(Option_t* detectors)
1495 //unloads Digits for detectors specified in parameter
1499 const char* oAll = strstr(detectors,"all");
1506 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1507 loaders = &arr;//get the pointer to array
1510 TIter next(loaders);
1512 while((loader = (AliLoader*)next()))
1514 loader->UnloadDigits();
1517 /**************************************************************************/
1519 void AliRunLoader::UnloadRecPoints(Option_t* detectors)
1521 //unloads RecPoints for detectors specified in parameter
1525 const char* oAll = strstr(detectors,"all");
1532 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1533 loaders = &arr;//get the pointer to array
1536 TIter next(loaders);
1538 while((loader = (AliLoader*)next()))
1540 loader->UnloadRecPoints();
1543 /**************************************************************************/
1545 void AliRunLoader::UnloadAll(Option_t* detectors)
1547 //calls UnloadAll for detectors names specified in parameter
1548 // option "all" passed can be passed
1552 const char* oAll = strstr(detectors,"all");
1559 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1560 loaders = &arr;//get the pointer to array
1563 TIter next(loaders);
1565 while((loader = (AliLoader*)next()))
1567 loader->UnloadAll();
1570 /**************************************************************************/
1572 void AliRunLoader::UnloadTracks(Option_t* detectors)
1574 //unloads Tracks for detectors specified in parameter
1578 const char* oAll = strstr(detectors,"all");
1585 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1586 loaders = &arr;//get the pointer to array
1589 TIter next(loaders);
1591 while((loader = (AliLoader*)next()))
1593 loader->UnloadTracks();
1596 /**************************************************************************/
1598 void AliRunLoader::UnloadRecParticles(Option_t* detectors)
1600 //unloads Particles for detectors specified in parameter
1604 const char* oAll = strstr(detectors,"all");
1611 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1612 loaders = &arr;//get the pointer to array
1615 TIter next(loaders);
1617 while((loader = (AliLoader*)next()))
1619 loader->UnloadRecParticles();
1622 /**************************************************************************/
1624 AliRunLoader* AliRunLoader::GetRunLoader(const char* eventfoldername)
1626 //returns RunLoader from folder named eventfoldername
1627 TFolder* evfold= dynamic_cast<TFolder*>(AliConfig::Instance()->GetTopFolder()->FindObject(eventfoldername));
1632 AliRunLoader* runget = dynamic_cast<AliRunLoader*>(evfold->FindObject(AliRunLoader::fgkRunLoaderName));
1636 /**************************************************************************/
1638 AliLoader* AliRunLoader::GetDetectorLoader(const char* detname, const char* eventfoldername)
1640 //get the loader of the detector with the given name from the global
1642 AliRunLoader* runLoader = GetRunLoader(eventfoldername);
1644 AliErrorClass("No run loader found");
1647 return runLoader->GetDetectorLoader(detname);
1649 /**************************************************************************/
1651 AliLoader* AliRunLoader::GetDetectorLoader(const char* detname)
1653 //get the loader of the detector with the given name from the global
1656 char loadername[256];
1657 sprintf(loadername, "%sLoader", detname);
1658 AliLoader* loader = GetLoader(loadername);
1660 AliError(Form("No loader for %s found", detname));
1665 /**************************************************************************/
1667 TTree* AliRunLoader::GetTreeH(const char* detname, Bool_t maketree, const char* eventfoldername)
1669 //get the tree with hits of the detector with the given name
1670 //if maketree is true and the tree does not exist, the tree is created
1671 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1672 if (!loader) return NULL;
1673 if (!loader->TreeH() && maketree) loader->MakeTree("H");
1674 return loader->TreeH();
1677 /**************************************************************************/
1679 TTree* AliRunLoader::GetTreeH(const char* detname, Bool_t maketree)
1681 //get the tree with hits of the detector with the given name
1682 //if maketree is true and the tree does not exist, the tree is created
1683 AliLoader* loader = GetDetectorLoader(detname);
1684 if (!loader) return NULL;
1685 if (!loader->TreeH() && maketree) loader->MakeTree("H");
1686 return loader->TreeH();
1688 /**************************************************************************/
1690 TTree* AliRunLoader::GetTreeS(const char* detname, Bool_t maketree,const char* eventfoldername)
1692 //get the tree with summable digits of the detector with the given name
1693 //if maketree is true and the tree does not exist, the tree is created
1694 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1695 if (!loader) return NULL;
1696 if (!loader->TreeS() && maketree) loader->MakeTree("S");
1697 return loader->TreeS();
1699 /**************************************************************************/
1701 TTree* AliRunLoader::GetTreeS(const char* detname, Bool_t maketree)
1703 //get the tree with summable digits of the detector with the given name
1704 //if maketree is true and the tree does not exist, the tree is created
1705 AliLoader* loader = GetDetectorLoader(detname);
1706 if (!loader) return NULL;
1707 if (!loader->TreeS() && maketree) loader->MakeTree("S");
1708 return loader->TreeS();
1710 /**************************************************************************/
1712 TTree* AliRunLoader::GetTreeD(const char* detname, Bool_t maketree,const char* eventfoldername)
1714 //get the tree with digits of the detector with the given name
1715 //if maketree is true and the tree does not exist, the tree is created
1716 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1717 if (!loader) return NULL;
1718 if (!loader->TreeD() && maketree) loader->MakeTree("D");
1719 return loader->TreeD();
1721 /**************************************************************************/
1723 TTree* AliRunLoader::GetTreeD(const char* detname, Bool_t maketree)
1725 //get the tree with digits of the detector with the given name
1726 //if maketree is true and the tree does not exist, the tree is created
1727 AliLoader* loader = GetDetectorLoader(detname);
1728 if (!loader) return NULL;
1729 if (!loader->TreeD() && maketree) loader->MakeTree("D");
1730 return loader->TreeD();
1732 /**************************************************************************/
1733 TTree* AliRunLoader::GetTreeR(const char* detname, Bool_t maketree,const char* eventfoldername)
1735 //get the tree with clusters of the detector with the given name
1736 //if maketree is true and the tree does not exist, the tree is created
1737 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1738 if (!loader) return NULL;
1739 if (!loader->TreeR() && maketree) loader->MakeTree("R");
1740 return loader->TreeR();
1742 /**************************************************************************/
1744 TTree* AliRunLoader::GetTreeR(const char* detname, Bool_t maketree)
1746 //get the tree with clusters of the detector with the given name
1747 //if maketree is true and the tree does not exist, the tree is created
1748 AliLoader* loader = GetDetectorLoader(detname);
1749 if (!loader) return NULL;
1750 if (!loader->TreeR() && maketree) loader->MakeTree("R");
1751 return loader->TreeR();
1753 /**************************************************************************/
1755 TTree* AliRunLoader::GetTreeT(const char* detname, Bool_t maketree,const char* eventfoldername)
1757 //get the tree with tracks of the detector with the given name
1758 //if maketree is true and the tree does not exist, the tree is created
1759 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1760 if (!loader) return NULL;
1761 if (!loader->TreeT() && maketree) loader->MakeTree("T");
1762 return loader->TreeT();
1764 /**************************************************************************/
1766 TTree* AliRunLoader::GetTreeT(const char* detname, Bool_t maketree)
1768 //get the tree with tracks of the detector with the given name
1769 //if maketree is true and the tree does not exist, the tree is created
1770 AliLoader* loader = GetDetectorLoader(detname);
1771 if (!loader) return NULL;
1772 if (!loader->TreeT() && maketree) loader->MakeTree("T");
1773 return loader->TreeT();
1775 /**************************************************************************/
1777 TTree* AliRunLoader::GetTreeP(const char* detname, Bool_t maketree,const char* eventfoldername)
1779 //get the tree with particles of the detector with the given name
1780 //if maketree is true and the tree does not exist, the tree is created
1781 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1782 if (!loader) return NULL;
1783 if (!loader->TreeP() && maketree) loader->MakeTree("P");
1784 return loader->TreeP();
1786 /**************************************************************************/
1788 TTree* AliRunLoader::GetTreeP(const char* detname, Bool_t maketree)
1790 //get the tree with particles of the detector with the given name
1791 //if maketree is true and the tree does not exist, the tree is created
1792 AliLoader* loader = GetDetectorLoader(detname);
1793 if (!loader) return NULL;
1794 if (!loader->TreeP() && maketree) loader->MakeTree("P");
1795 return loader->TreeP();
1798 /**************************************************************************/
1800 void AliRunLoader::CdGAFile()
1802 //sets gDirectory to galice file
1804 if(fGAFile) fGAFile->cd();
1807 /**************************************************************************/
1809 void AliRunLoader::GetListOfDetectors(const char * namelist,TObjArray& pointerarray) const
1811 //this method looks for all Loaders corresponding
1812 //to names (many) specified in namelist i.e. namelist ("ITS TPC TRD")
1816 strcpy(dets,namelist);//compiler cries when char* = const Option_t*;
1817 dets[strlen(dets)+1] = '\n';//set endl at the end of string
1822 tmp = sscanf(pdet,"%s",buff);//read the string from the input string pdet into buff
1823 if ( (buff[0] == 0) || (tmp == 0) ) break; //if not read
1825 pdet = strstr(pdet,buff) + strlen(buff);//move the input pointer about number of bytes (letters) read
1826 //I am aware that is a little complicated. I don't know the number of spaces between detector names
1827 //so I read the string, than I find where it starts (strstr) and move the pointer about length of a string
1828 //If there is a better way, please write me (Piotr.Skowronski@cern.ch)
1829 //construct the Loader name
1830 TString getname(buff);
1832 AliLoader* loader = GetLoader(getname);//get the Loader
1835 pointerarray.Add(loader);
1839 AliError(Form("Can not find Loader for %s",buff));
1845 /*****************************************************************************/
1847 void AliRunLoader::Clean(const TString& name)
1849 //removes object with given name from event folder and deletes it
1850 if (GetEventFolder() == 0x0) return;
1851 TObject* obj = GetEventFolder()->FindObject(name);
1854 AliDebug(1, Form("name=%s, cleaning %s.",GetName(),name.Data()));
1855 GetEventFolder()->Remove(obj);
1861 /*****************************************************************************/
1863 TTask* AliRunLoader::GetRunDigitizer()
1865 //returns Run Digitizer from folder
1867 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1868 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetDigitizerTaskName());
1869 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1871 /*****************************************************************************/
1873 TTask* AliRunLoader::GetRunSDigitizer()
1875 //returns SDigitizer Task from folder
1877 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1878 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetSDigitizerTaskName());
1879 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1881 /*****************************************************************************/
1883 TTask* AliRunLoader::GetRunReconstructioner()
1885 //returns Reconstructioner Task from folder
1886 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1887 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetReconstructionerTaskName());
1888 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1890 /*****************************************************************************/
1892 TTask* AliRunLoader::GetRunTracker()
1894 //returns Tracker Task from folder
1895 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1896 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetTrackerTaskName());
1897 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1899 /*****************************************************************************/
1901 TTask* AliRunLoader::GetRunPIDTask()
1903 //returns Tracker Task from folder
1904 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1905 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetPIDTaskName());
1906 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1908 /*****************************************************************************/
1910 TTask* AliRunLoader::GetRunQATask()
1912 //returns Quality Assurance Task from folder
1913 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1916 AliErrorClass("Can not get task folder from AliConfig");
1919 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetQATaskName());
1920 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1923 /*****************************************************************************/
1925 void AliRunLoader::SetCompressionLevel(Int_t cl)
1927 //Sets Compression Level in all files
1928 if (fGAFile) fGAFile->SetCompressionLevel(cl);
1929 SetKineComprLevel(cl);
1930 SetTrackRefsComprLevel(cl);
1931 TIter next(fLoaders);
1933 while((loader = (AliLoader*)next()))
1935 loader->SetCompressionLevel(cl);
1938 /**************************************************************************/
1940 void AliRunLoader::SetKineComprLevel(Int_t cl)
1942 //Sets comression level in Kine File
1943 fKineDataLoader->SetCompressionLevel(cl);
1945 /**************************************************************************/
1947 void AliRunLoader::SetTrackRefsComprLevel(Int_t cl)
1949 //Sets comression level in Track Refences File
1950 fTrackRefsDataLoader->SetCompressionLevel(cl);
1952 /**************************************************************************/
1954 void AliRunLoader::UnloadHeader()
1956 //removes TreeE from folder and deletes it
1957 // as well as fHeader object
1962 /**************************************************************************/
1964 void AliRunLoader::UnloadTrigger()
1966 //removes TreeCT from folder and deletes it
1967 // as well as fHeader object
1973 /**************************************************************************/
1975 void AliRunLoader::UnloadKinematics()
1977 //Unloads Kinematics
1978 fKineDataLoader->GetBaseLoader(0)->Unload();
1980 /**************************************************************************/
1982 void AliRunLoader::UnloadTrackRefs()
1984 //Unloads Track Refernces
1985 fTrackRefsDataLoader->GetBaseLoader(0)->Unload();
1987 /**************************************************************************/
1989 void AliRunLoader::UnloadgAlice()
1992 if (gAlice == GetAliRun())
1994 AliDebug(1, "Set gAlice = 0x0");
1995 gAlice = 0x0;//if gAlice is the same that in folder (to be deleted by the way of folder)
1997 AliRun* alirun = GetAliRun();
1998 if (GetEventFolder()) GetEventFolder()->Remove(alirun);
2001 /**************************************************************************/
2003 void AliRunLoader::MakeTrackRefsContainer()
2005 // Makes a tree for Track References
2006 fTrackRefsDataLoader->MakeTree();
2008 /**************************************************************************/
2010 Int_t AliRunLoader::LoadTrackRefs(Option_t* option)
2012 //Load track references from file (opens file and posts tree to folder)
2014 return fTrackRefsDataLoader->GetBaseLoader(0)->Load(option);
2016 /**************************************************************************/
2018 void AliRunLoader::SetDirName(TString& dirname)
2020 //sets directory name
2021 if (dirname.IsNull()) return;
2022 fUnixDirName = dirname;
2023 fKineDataLoader->SetDirName(dirname);
2024 fTrackRefsDataLoader->SetDirName(dirname);
2026 TIter next(fLoaders);
2028 while((loader = (AliLoader*)next()))
2030 loader->SetDirName(dirname);
2034 /*****************************************************************************/
2036 Int_t AliRunLoader::GetFileOffset() const
2038 //returns the file number that is added to the file name for current event
2039 return Int_t(fCurrentEvent/fNEventsPerFile);
2042 /*****************************************************************************/
2043 const TString AliRunLoader::SetFileOffset(const TString& fname)
2045 //adds the the number to the file name at proper place for current event
2046 Long_t offset = (Long_t)GetFileOffset();
2047 if (offset < 1) return fname;
2049 soffset += offset;//automatic conversion to string
2050 TString dotroot(".root");
2051 const TString& offfsetdotroot = offset + dotroot;
2052 TString out = fname;
2053 out = out.ReplaceAll(dotroot,offfsetdotroot);
2054 AliDebug(1, Form(" in=%s out=%s",fname.Data(),out.Data()));
2057 /*****************************************************************************/
2059 void AliRunLoader::SetDigitsFileNameSuffix(const TString& suffix)
2061 //adds the suffix before ".root",
2062 //e.g. TPC.Digits.root -> TPC.DigitsMerged.root
2063 //made on Jiri Chudoba demand
2065 TIter next(fLoaders);
2067 while((loader = (AliLoader*)next()))
2069 loader->SetDigitsFileNameSuffix(suffix);
2072 /*****************************************************************************/
2074 TString AliRunLoader::GetFileName() const
2076 //returns name of galice file
2078 if (fGAFile == 0x0) return result;
2079 result = fGAFile->GetName();
2082 /*****************************************************************************/
2084 void AliRunLoader::SetDetectorAddresses()
2086 //calls SetTreeAddress for all detectors
2087 if (GetAliRun()==0x0) return;
2088 TIter next(GetAliRun()->Modules());
2090 while((mod = (AliModule*)next()))
2092 AliDetector* det = dynamic_cast<AliDetector*>(mod);
2093 if (det) det->SetTreeAddress();
2096 /*****************************************************************************/
2098 void AliRunLoader::Synchronize()
2100 //synchrinizes all writtable files
2101 TIter next(fLoaders);
2103 while((loader = (AliLoader*)next()))
2105 loader->Synchronize();
2108 fKineDataLoader->Synchronize();
2109 fTrackRefsDataLoader->Synchronize();
2111 TFile* file = gROOT->GetFile( gSystem->ConcatFileName( fUnixDirName.Data(), fgkDefaultTriggerFileName.Data() ) ) ;
2112 if( file ) file->Flush();
2114 if (fGAFile) fGAFile->Flush();
2116 /*****************************************************************************/
2117 /*****************************************************************************/