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():
95 fTrackRefsDataLoader(0x0),
99 AliConfig::Instance();//force to build the folder structure
100 if (!fgRunLoader) fgRunLoader = this;
102 /**************************************************************************/
104 AliRunLoader::AliRunLoader(const char* eventfoldername):
105 TNamed(fgkRunLoaderName,fgkRunLoaderName),
106 fLoaders(new TObjArray()),
113 fKineDataLoader(new AliDataLoader(fgkDefaultKineFileName,fgkKineContainerName,"Kinematics")),
114 fTrackRefsDataLoader(new AliDataLoader(fgkDefaultTrackRefsFileName,fgkTrackRefsContainerName,"Track References")),
119 SetEventFolderName(eventfoldername);
120 if (!fgRunLoader) fgRunLoader = this;
122 /**************************************************************************/
124 AliRunLoader::~AliRunLoader()
127 if (fgRunLoader == this) fgRunLoader = 0x0;
133 fLoaders->SetOwner();
137 delete fKineDataLoader;
138 delete fTrackRefsDataLoader;
143 //fEventFolder is deleted by the way of removing - TopAliceFolder owns it
144 if( fCTrigger ) delete fCTrigger;
149 /**************************************************************************/
151 AliRunLoader::AliRunLoader(TFolder* topfolder):
152 TNamed(fgkRunLoaderName,fgkRunLoaderName),
153 fLoaders(new TObjArray()),
154 fEventFolder(topfolder),
160 fKineDataLoader(new AliDataLoader(fgkDefaultKineFileName,fgkKineContainerName,"Kinematics")),
161 fTrackRefsDataLoader(new AliDataLoader(fgkDefaultTrackRefsFileName,fgkTrackRefsContainerName,"Track References")),
168 TString errmsg("Parameter is NULL");
169 AliError(errmsg.Data());
174 TObject* obj = fEventFolder->FindObject(fgkRunLoaderName);
176 { //if it is, then sth. is going wrong... exits aliroot session
177 TString errmsg("In Event Folder Named ");
178 errmsg+=fEventFolder->GetName();
179 errmsg+=" object named "+fgkRunLoaderName+" already exists. I am confused ...";
181 AliError(errmsg.Data());
183 return;//never reached
186 if (!fgRunLoader) fgRunLoader = this;
188 fEventFolder->Add(this);//put myself to the folder to accessible for all
192 /**************************************************************************/
194 Int_t AliRunLoader::GetEvent(Int_t evno)
196 //Gets event number evno
197 //Reloads all data properly
198 //PH if (fCurrentEvent == evno) return 0;
202 AliError("Can not give the event with negative number");
206 if (evno >= GetNumberOfEvents())
208 AliError(Form("There is no event with number %d",evno));
212 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
213 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
214 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
215 AliDebug(1, Form(" GETTING EVENT %d",evno));
216 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
217 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
218 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
220 fCurrentEvent = evno;
224 //Reload header (If header was loaded)
227 retval = TreeE()->GetEvent(fCurrentEvent);
230 AliError(Form("Cannot find event: %d\n ",fCurrentEvent));
234 //Reload stack (If header was loaded)
235 if (TreeE()) fStack = GetHeader()->Stack();
236 //Set event folder in stack (it does not mean that we read kinematics from file)
237 if( GetTrigger() && TreeCT() ) {
238 retval = TreeCT()->GetEvent(fCurrentEvent);
240 AliError(Form("Error occured while GetEvent for Trigger. Event %d",evno));
248 AliError(Form("Error occured while setting event %d",evno));
252 //Post Track References
253 retval = fTrackRefsDataLoader->GetEvent();
256 AliError(Form("Error occured while GetEvent for Track References. Event %d",evno));
260 //Read Kinematics if loaded
261 retval = fKineDataLoader->GetEvent();
264 AliError(Form("Error occured while GetEvent for Kinematics. Event %d",evno));
268 if (fStack && fKineDataLoader->GetBaseLoader(0)->IsLoaded())
270 fStack->ConnectTree(TreeK());
272 if (fStack->GetEvent() == kFALSE)
274 AliError(Form("Error occured while GetEvent for Stack. Event %d",evno));
279 //Trigger data reloading in all loaders
280 TIter next(fLoaders);
282 while((loader = (AliLoader*)next()))
284 retval = loader->GetEvent();
287 AliError(Form("Error occured while getting event for %s. Event %d.",
288 loader->GetDetectorName().Data(), evno));
293 SetDetectorAddresses();
297 /**************************************************************************/
298 Int_t AliRunLoader::SetEvent()
300 //if kinematics was loaded Cleans folder data
304 retval = fKineDataLoader->SetEvent();
307 AliError("SetEvent for Kinamtics Data Loader retutned error.");
310 retval = fTrackRefsDataLoader->SetEvent();
313 AliError("SetEvent for Track References Data Loader retutned error.");
317 TIter next(fLoaders);
319 while((loader = (AliLoader*)next()))
321 retval = loader->SetEvent();
324 AliError(Form("SetEvent for %s Data Loader retutned error.",loader->GetName()));
331 /**************************************************************************/
333 Int_t AliRunLoader::SetEventNumber(Int_t evno)
335 //cleans folders and sets the root dirs in files
336 if (fCurrentEvent == evno) return 0;
337 fCurrentEvent = evno;
341 /**************************************************************************/
342 AliCDBEntry* AliRunLoader::GetCDBEntry(const char* name) const
344 //Get an AliCDBEntry from the run data storage
346 if ( !(AliCDBManager::Instance()->IsDefaultStorageSet()) ) {
347 AliError("No run data storage defined!");
350 return AliCDBManager::Instance()->GetDefaultStorage()->Get(name, GetHeader()->GetRun());
354 /**************************************************************************/
355 AliRunLoader* AliRunLoader::Open
356 (const char* filename, const char* eventfoldername, Option_t* option)
358 //Opens a desired file 'filename'
359 //gets the the run-Loader and mounts it desired folder
360 //returns the pointer to run Loader which can be further used for accessing data
361 //in case of error returns NULL
363 static const TString kwebaddress("http://alisoft.cern.ch/people/skowron/codedoc/split/index.html");
364 AliDebugClass(1,Form("\n\n\nNew I/O strcture: See more info:\n %s\n\n\n",kwebaddress.Data()));
366 AliRunLoader* result = 0x0;
368 /* ************************************************ */
369 /* Chceck if folder with given name already exists */
370 /* ************************************************ */
372 TObject* obj = AliConfig::Instance()->GetTopFolder()->FindObject(eventfoldername);
375 TFolder* fold = dynamic_cast<TFolder*>(obj);
378 AliErrorClass("Such a obejct already exists in top alice folder and it is not a folder.");
382 //check if we can get RL from that folder
383 result = AliRunLoader::GetRunLoader(eventfoldername);
386 AliErrorClass(Form("Folder %s already exists, and can not find session there. Can not mount.",eventfoldername));
390 if (result->GetFileName().CompareTo(filename) != 0)
392 AliErrorClass("Other file is mounted in demanded folder. Can not mount.");
396 //check if now is demanded (re)creation
397 if ( AliLoader::TestFileOption(option) == kFALSE)
399 AliErrorClass(Form("Session already exists in folder %s and this session option is %s. Unable to proceed.",
400 eventfoldername,option));
404 //check if demanded option is update and existing one
405 TString tmpstr(option);
406 if ( (tmpstr.CompareTo("update",TString::kIgnoreCase) == 0) &&
407 (result->fGAFile->IsWritable() == kFALSE) )
409 AliErrorClass(Form("Session already exists in folder %s and is not writable while this session option is %s. Unable to proceed.",
410 eventfoldername,option));
414 AliWarningClass("Session is already opened and mounted in demanded folder");
415 if (!fgRunLoader) fgRunLoader = result; //PH get access from any place
417 } //end of checking in case of existance of object named identically that folder session is being opened
420 TFile * gAliceFile = TFile::Open(filename,option);//open a file
422 {//null pointer returned
423 AliFatalClass(Form("Can not open file %s.",filename));
427 if (gAliceFile->IsOpen() == kFALSE)
428 {//pointer to valid object returned but file is not opened
429 AliErrorClass(Form("Can not open file %s.",filename));
433 //if file is "read" or "update" than we try to find AliRunLoader there - if not found cry and exit
434 //else create new AliRunLoader
435 if ( AliLoader::TestFileOption(option) )
437 AliDebugClass(1, "Reading RL from file");
439 result = dynamic_cast<AliRunLoader*>(gAliceFile->Get(fgkRunLoaderName));//get the run Loader from the file
442 AliErrorClass(Form("Can not find run-Loader in file %s.",filename));
443 delete gAliceFile;//close the file
446 Int_t tmp = result->SetEventFolderName(eventfoldername);//mount a event folder
447 if (tmp)//if SetEvent returned error
449 AliErrorClass(Form("Can not mount event in folder %s.",eventfoldername));
450 delete result; //delete run-Loader
451 delete gAliceFile;//close the file
457 AliDebugClass(1, Form("Creating new AliRunLoader. Folder name is %s",eventfoldername));
460 result = new AliRunLoader(eventfoldername);
462 catch (TString& errmsg)
464 AliErrorClass(Form("AliRunLoader constrcutor has thrown exception: %s\n",errmsg.Data()));
466 delete gAliceFile;//close the file
471 //procedure for extracting dir name from the file name
472 TString fname(filename);
473 Int_t nsl = fname.Last('#');//look for hash in file name
475 if (nsl < 0) {//hash not found
476 nsl = fname.Last('/');// look for slash
478 nsl = fname.Last(':');// look for colon e.g. rfio:galice.root
481 if (nsl < 0) dirname = "./"; // no directory path, use "."
482 else dirname = fname.Remove(nsl+1);// directory path
484 AliDebugClass(1, Form("Dir name is : %s",dirname.Data()));
486 result->SetDirName(dirname);
487 result->SetGAliceFile(gAliceFile);//set the pointer to gAliceFile
488 if (!fgRunLoader) fgRunLoader = result; //PH get access from any place
491 /**************************************************************************/
492 Int_t AliRunLoader::GetNumberOfEvents()
494 //returns number of events in Run
498 retval = LoadHeader();
501 AliError("Error occured while loading header");
505 return (Int_t)TreeE()->GetEntries();
507 /**************************************************************************/
509 void AliRunLoader::MakeHeader()
511 //Makes header and connects it to header tree (if it exists)
515 AliDebug(1, "Creating new Header Object");
516 fHeader= new AliHeader();
518 TTree* tree = TreeE();
521 AliDebug(1, "Got Tree from folder.");
522 TBranch* branch = tree->GetBranch(fgkHeaderBranchName);
525 AliDebug(1, "Creating new branch");
526 branch = tree->Branch(fgkHeaderBranchName, "AliHeader", &fHeader, 4000, 0);
527 branch->SetAutoDelete(kFALSE);
531 AliDebug(1, "Got Branch from Tree");
532 branch->SetAddress(&fHeader);
533 tree->GetEvent(fCurrentEvent);
534 fStack = fHeader->Stack(); //should be safe - if we created Stack, header returns pointer to the same object
538 fStack->ConnectTree(TreeK());
544 AliDebug(1, "Header does not have a stack.");
548 AliDebug(1, "Exiting MakeHeader method");
550 /**************************************************************************/
552 void AliRunLoader::MakeStack()
554 //Creates the stack object - do not connect the tree
557 fStack = new AliStack(10000);
560 /**************************************************************************/
562 void AliRunLoader::MakeTrigger()
564 // Makes trigger object and connects it to trigger tree (if it exists)
566 if( fCTrigger == 0x0 ) {
567 AliDebug( 1, "Creating new Trigger Object" );
568 fCTrigger = new AliCentralTrigger();
570 TTree* tree = TreeCT();
572 fCTrigger->MakeBranch( fgkTriggerBranchName, tree );
573 tree->GetEvent( fCurrentEvent );
576 AliDebug( 1, "Exiting MakeTrigger method" );
578 /**************************************************************************/
580 void AliRunLoader::MakeTree(Option_t *option)
583 const char *oK = strstr(option,"K"); //Kine
584 const char *oE = strstr(option,"E"); //Header
585 const char *oGG = strstr(option,"GG"); //Central TriGGer
589 if (fKineDataLoader->GetBaseLoader(0)->IsLoaded() == kFALSE)
591 AliError("Load Kinematics first");
596 fKineDataLoader->MakeTree();
599 fStack->ConnectTree(TreeK());
600 WriteKinematics("OVERWRITE");
607 TTree* tree = new TTree(fgkHeaderContainerName,"Tree with Headers");
608 GetEventFolder()->Add(tree);
610 WriteHeader("OVERWRITE");
615 // create the CTP Trigger output file and tree
616 TFile* file = gROOT->GetFile( fgkDefaultTriggerFileName );
618 file = TFile::Open( gSystem->ConcatFileName( fUnixDirName.Data(), fgkDefaultTriggerFileName.Data() ), "RECREATE" ) ;
622 TTree* tree = new TTree( fgkTriggerContainerName, "Tree with Central Trigger Mask" );
623 GetEventFolder()->Add(tree);
625 // WriteHeader("OVERWRITE");
628 TIter next(fLoaders);
630 while((loader = (AliLoader*)next()))
632 loader->MakeTree(option);
636 /**************************************************************************/
638 Int_t AliRunLoader::LoadgAlice()
640 //Loads gAlice from file
643 AliWarning("AliRun is already in folder. Unload first.");
646 AliRun* alirun = dynamic_cast<AliRun*>(fGAFile->Get(fgkGAliceName));
649 AliError(Form("Can not find gAlice in file %s",fGAFile->GetName()));
652 alirun->SetRunLoader(this);
655 AliWarning(Form("gAlice already exists. Putting retrived object in folder named %s",
656 GetEventFolder()->GetName()));
662 SetDetectorAddresses();//calls SetTreeAddress for all detectors
665 /**************************************************************************/
667 Int_t AliRunLoader::LoadHeader()
669 //loads treeE and reads header object for current event
672 AliWarning("Header is already loaded. Use ReloadHeader to force reload. Nothing done");
676 if (GetEventFolder() == 0x0)
678 AliError("Event folder not specified yet");
684 AliError("Session not opened. Use AliRunLoader::Open");
688 if (fGAFile->IsOpen() == kFALSE)
690 AliError("Session not opened. Use AliRunLoader::Open");
694 TTree* tree = dynamic_cast<TTree*>(fGAFile->Get(fgkHeaderContainerName));
697 AliError(Form("Can not find header tree named %s in file %s",
698 fgkHeaderContainerName.Data(),fGAFile->GetName()));
702 if (tree == TreeE()) return 0;
705 GetEventFolder()->Add(tree);
706 MakeHeader();//creates header object and connects to tree
710 /**************************************************************************/
712 Int_t AliRunLoader::LoadTrigger(Option_t* option)
717 AliWarning("Trigger is already loaded. Nothing done");
721 if( GetEventFolder() == 0x0 ) {
722 AliError("Event folder not specified yet");
725 // get the CTP Trigger output file and tree
726 TString trgfile = gSystem->ConcatFileName( fUnixDirName.Data(),
727 fgkDefaultTriggerFileName.Data() );
728 TFile* file = gROOT->GetFile( trgfile );
730 file = TFile::Open( trgfile, option ) ;
731 if (!file || file->IsOpen() == kFALSE ) {
732 AliError( Form( "Can not open trigger file %s", trgfile.Data() ) );
738 TTree* tree = dynamic_cast<TTree*>(file->Get( fgkTriggerContainerName ));
740 AliError( Form( "Can not find trigger tree named %s in file %s",
741 fgkTriggerContainerName.Data(), file->GetName() ) );
747 fCTrigger = dynamic_cast<AliCentralTrigger*>(file->Get( "AliCentralTrigger" ));
748 GetEventFolder()->Add( tree );
754 /**************************************************************************/
756 Int_t AliRunLoader::LoadKinematics(Option_t* option)
758 //Loads the kinematics
759 Int_t retval = fKineDataLoader->GetBaseLoader(0)->Load(option);
762 AliError("Error occured while loading kinamatics tree.");
767 fStack->ConnectTree(TreeK());
768 retval = fStack->GetEvent();
769 if ( retval == kFALSE)
771 AliError("Error occured while loading kinamatics tree.");
778 /**************************************************************************/
780 Int_t AliRunLoader::OpenDataFile(const TString& filename,TFile*& file,TDirectory*& dir,Option_t* opt,Int_t cl)
782 //Opens File with kinematics
785 if (file->IsOpen() == kFALSE)
786 {//pointer is not null but file is not opened
787 AliWarning("Pointer to file is not null, but file is not opened");//risky any way
789 file = 0x0; //proceed with opening procedure
793 AliWarning(Form("File %s already opened",filename.Data()));
797 //try to find if that file is opened somewere else
798 file = (TFile *)( gROOT->GetListOfFiles()->FindObject(filename) );
801 if(file->IsOpen() == kTRUE)
803 AliWarning(Form("File %s already opened by sombody else.",file->GetName()));
808 file = TFile::Open(filename,opt);
811 AliError(Form("Can not open file %s",filename.Data()));
814 if (file->IsOpen() == kFALSE)
815 {//file is not opened
816 AliError(Form("Can not open file %s",filename.Data()));
820 file->SetCompressionLevel(cl);
822 dir = AliLoader::ChangeDir(file,fCurrentEvent);
825 AliError(Form("Can not change to root directory in file %s",filename.Data()));
830 /**************************************************************************/
832 TTree* AliRunLoader::TreeE() const
834 //returns the tree from folder; shortcut method
835 if (AliDebugLevel() > 10) fEventFolder->ls();
836 TObject *obj = fEventFolder->FindObject(fgkHeaderContainerName);
837 return (obj)?dynamic_cast<TTree*>(obj):0x0;
839 /**************************************************************************/
841 TTree* AliRunLoader::TreeCT() const
843 //returns the tree from folder; shortcut method
844 if (AliDebugLevel() > 10) fEventFolder->ls();
845 TObject *obj = fEventFolder->FindObject(fgkTriggerContainerName);
846 return (obj)?dynamic_cast<TTree*>(obj):0x0;
848 /**************************************************************************/
850 AliHeader* AliRunLoader::GetHeader() const
852 //returns pointer header object
855 /**************************************************************************/
857 AliCentralTrigger* AliRunLoader::GetTrigger() const
859 //returns pointer trigger object
863 /**************************************************************************/
865 TTree* AliRunLoader::TreeK() const
867 //returns the tree from folder; shortcut method
868 TObject *obj = GetEventFolder()->FindObject(fgkKineContainerName);
869 return (obj)?dynamic_cast<TTree*>(obj):0x0;
871 /**************************************************************************/
873 TTree* AliRunLoader::TreeTR() const
875 //returns the tree from folder; shortcut method
876 TObject* obj = GetEventFolder()->FindObject(fgkTrackRefsContainerName);
877 return (obj)?dynamic_cast<TTree*>(obj):0x0;
879 /**************************************************************************/
881 AliRun* AliRunLoader::GetAliRun() const
883 //returns AliRun which sits in the folder
884 if (fEventFolder == 0x0) return 0x0;
885 TObject *obj = fEventFolder->FindObject(fgkGAliceName);
886 return (obj)?dynamic_cast<AliRun*>(obj):0x0;
888 /**************************************************************************/
890 Int_t AliRunLoader::WriteHeader(Option_t* opt)
893 AliDebug(1, "WRITING HEADER");
895 TTree* tree = TreeE();
898 AliWarning("Can not find Header Tree in Folder");
901 if (fGAFile->IsWritable() == kFALSE)
903 AliError(Form("File %s is not writable",fGAFile->GetName()));
907 TObject* obj = fGAFile->Get(fgkHeaderContainerName);
909 { //if they exist, see if option OVERWRITE is used
911 if(tmp.Contains("OVERWRITE",TString::kIgnoreCase) == 0)
912 {//if it is not used - give an error message and return an error code
913 AliError("Tree already exisists. Use option \"OVERWRITE\" to overwrite previous data");
918 tree->SetDirectory(fGAFile);
919 tree->Write(0,TObject::kOverwrite);
921 AliDebug(1, "WRITTEN\n\n");
926 /**************************************************************************/
928 Int_t AliRunLoader::WriteTrigger(Option_t* opt)
931 AliDebug( 1, "WRITING TRIGGER" );
933 TTree* tree = TreeCT();
935 AliWarning("Can not find Trigger Tree in Folder");
939 TFile* file = gROOT->GetFile( gSystem->ConcatFileName( fUnixDirName.Data(), fgkDefaultTriggerFileName.Data() ) ) ;
940 if( !file || !file->IsOpen() ) {
941 AliError( "can't write Trigger, file is not open" );
945 TObject* obj = file->Get( fgkTriggerContainerName );
946 if( obj ) { //if they exist, see if option OVERWRITE is used
948 if( tmp.Contains( "OVERWRITE", TString::kIgnoreCase ) == 0) {
949 //if it is not used - give an error message and return an error code
950 AliError( "Tree already exisists. Use option \"OVERWRITE\" to overwrite previous data" );
955 fCTrigger->Write( 0, TObject::kOverwrite );
956 tree->Write( 0, TObject::kOverwrite );
959 AliDebug(1, "WRITTEN\n\n");
963 /**************************************************************************/
965 Int_t AliRunLoader::WriteAliRun(Option_t* /*opt*/)
967 //writes AliRun object to the file
969 if (GetAliRun()) GetAliRun()->Write();
972 /**************************************************************************/
974 Int_t AliRunLoader::WriteKinematics(Option_t* opt)
977 return fKineDataLoader->GetBaseLoader(0)->WriteData(opt);
979 /**************************************************************************/
980 Int_t AliRunLoader::WriteTrackRefs(Option_t* opt)
982 //writes Track References tree
983 return fTrackRefsDataLoader->GetBaseLoader(0)->WriteData(opt);
985 /**************************************************************************/
987 Int_t AliRunLoader::WriteHits(Option_t* opt)
989 //Calls WriteHits for all loaders
992 TIter next(fLoaders);
994 while((loader = (AliLoader*)next()))
996 res = loader->WriteHits(opt);
999 AliError(Form("Failed to write hits for %s (%d)",loader->GetDetectorName().Data(),res));
1005 /**************************************************************************/
1007 Int_t AliRunLoader::WriteSDigits(Option_t* opt)
1009 //Calls WriteSDigits for all loaders
1012 TIter next(fLoaders);
1014 while((loader = (AliLoader*)next()))
1016 res = loader->WriteSDigits(opt);
1019 AliError(Form("Failed to write summable digits for %s.",loader->GetDetectorName().Data()));
1025 /**************************************************************************/
1027 Int_t AliRunLoader::WriteDigits(Option_t* opt)
1029 //Calls WriteDigits for all loaders
1032 TIter next(fLoaders);
1034 while((loader = (AliLoader*)next()))
1036 res = loader->WriteDigits(opt);
1039 AliError(Form("Failed to write digits for %s.",loader->GetDetectorName().Data()));
1045 /**************************************************************************/
1047 Int_t AliRunLoader::WriteRecPoints(Option_t* opt)
1049 //Calls WriteRecPoints for all loaders
1052 TIter next(fLoaders);
1054 while((loader = (AliLoader*)next()))
1056 res = loader->WriteRecPoints(opt);
1059 AliError(Form("Failed to write Reconstructed Points for %s.",
1060 loader->GetDetectorName().Data()));
1066 /**************************************************************************/
1068 Int_t AliRunLoader::WriteTracks(Option_t* opt)
1070 //Calls WriteTracks for all loaders
1073 TIter next(fLoaders);
1075 while((loader = (AliLoader*)next()))
1077 res = loader->WriteTracks(opt);
1080 AliError(Form("Failed to write Tracks for %s.",
1081 loader->GetDetectorName().Data()));
1087 /**************************************************************************/
1089 Int_t AliRunLoader::WriteRunLoader(Option_t* /*opt*/)
1091 //Writes itself to the file
1093 this->Write(0,TObject::kOverwrite);
1096 /**************************************************************************/
1098 Int_t AliRunLoader::SetEventFolderName(const TString& name)
1100 //sets top folder name for this run; of alread
1103 AliError("Name is empty");
1107 //check if such a folder already exists - try to find it in alice top folder
1108 TObject* obj = AliConfig::Instance()->GetTopFolder()->FindObject(name);
1111 TFolder* fold = dynamic_cast<TFolder*>(obj);
1114 AliError("Such a obejct already exists in top alice folder and it is not a folder.");
1117 //folder which was found is our folder
1118 if (fEventFolder == fold)
1124 AliError("Such a folder already exists in top alice folder. Can not mount.");
1129 //event is alredy connected, just change name of the folder
1132 fEventFolder->SetName(name);
1136 if (fKineDataLoader == 0x0)
1137 fKineDataLoader = new AliDataLoader(fgkDefaultKineFileName,fgkKineContainerName,"Kinematics");
1139 if ( fTrackRefsDataLoader == 0x0)
1140 fTrackRefsDataLoader = new AliDataLoader(fgkDefaultTrackRefsFileName,fgkTrackRefsContainerName,"Track References");
1142 //build the event folder structure
1143 AliDebug(1, Form("Creating new event folder named %s",name.Data()));
1144 fEventFolder = AliConfig::Instance()->BuildEventFolder(name,"Event Folder");
1145 fEventFolder->Add(this);//put myself to the folder to accessible for all
1147 TIter next(fLoaders);
1149 while((loader = (AliLoader*)next()))
1151 loader->Register(fEventFolder);//build folder structure for this detector
1154 fKineDataLoader->SetEventFolder(GetEventFolder());
1155 fTrackRefsDataLoader->SetEventFolder(GetEventFolder());
1156 fKineDataLoader->SetFolder(GetEventFolder());
1157 fTrackRefsDataLoader->SetFolder(GetEventFolder());
1159 fEventFolder->SetOwner();
1162 /**************************************************************************/
1164 void AliRunLoader::AddLoader(AliLoader* loader)
1166 //Adds the Loader for given detector
1167 if (loader == 0x0) //if null shout and exit
1169 AliError("Parameter is NULL");
1172 loader->SetDirName(fUnixDirName);
1173 if (fEventFolder) loader->SetEventFolder(fEventFolder); //if event folder is already defined,
1174 //pass information to the Loader
1175 fLoaders->Add(loader);//add the Loader to the array
1177 /**************************************************************************/
1179 void AliRunLoader::AddLoader(AliDetector* det)
1181 //Asks module (detector) ro make a Loader and stores in the array
1182 if (det == 0x0) return;
1183 AliLoader* get = det->GetLoader();//try to get loader
1184 if (get == 0x0) get = det->MakeLoader(fEventFolder->GetName());//if did not obtain, ask to make it
1188 AliDebug(1, Form("Detector: %s Loader : %s",det->GetName(),get->GetName()));
1193 /**************************************************************************/
1195 AliLoader* AliRunLoader::GetLoader(const char* detname) const
1197 //returns loader for given detector
1198 //note that naming convention is TPCLoader not just TPC
1199 return (AliLoader*)fLoaders->FindObject(detname);
1202 /**************************************************************************/
1204 AliLoader* AliRunLoader::GetLoader(AliDetector* det) const
1206 //get loader for detector det
1207 if(det == 0x0) return 0x0;
1208 TString getname(det->GetName());
1210 AliDebug(1, Form(" Loader name is %s",getname.Data()));
1211 return GetLoader(getname);
1214 /**************************************************************************/
1216 void AliRunLoader::CleanFolders()
1218 // fEventFolder->Add(this);//put myself to the folder to accessible for all
1225 /**************************************************************************/
1227 void AliRunLoader::CleanDetectors()
1229 //Calls CleanFolders for all detectors
1230 TIter next(fLoaders);
1232 while((loader = (AliLoader*)next()))
1234 loader->CleanFolders();
1237 /**************************************************************************/
1239 void AliRunLoader::RemoveEventFolder()
1241 //remove all the tree of event
1242 //all the stuff changing EbE stays untached (PDGDB, tasks, etc.)
1244 if (fEventFolder == 0x0) return;
1245 fEventFolder->SetOwner(kFALSE);//don't we want to deleted while removing the folder that we are sitting in
1246 fEventFolder->Remove(this);//remove us drom folder
1248 AliConfig::Instance()->GetTopFolder()->SetOwner(); //brings ownership back for fEventFolder since it sits in top folder
1249 AliConfig::Instance()->GetTopFolder()->Remove(fEventFolder); //remove the event tree
1250 delete fEventFolder;
1252 /**************************************************************************/
1254 void AliRunLoader::SetGAliceFile(TFile* gafile)
1256 //sets pointer to galice.root file
1260 /**************************************************************************/
1262 Int_t AliRunLoader::LoadHits(Option_t* detectors,Option_t* opt)
1264 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1266 AliDebug(1, "Loading Hits");
1270 const char* oAll = strstr(detectors,"all");
1273 AliDebug(1, "Option is All");
1278 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1279 loaders = &arr;//get the pointer array
1282 AliDebug(1, Form("For detectors. Number of detectors chosen for loading %d",loaders->GetEntries()));
1284 TIter next(loaders);
1286 while((loader = (AliLoader*)next()))
1288 AliDebug(1, Form(" Calling LoadHits(%s) for %s",opt,loader->GetName()));
1289 loader->LoadHits(opt);
1291 AliDebug(1, "Done");
1295 /**************************************************************************/
1297 Int_t AliRunLoader::LoadSDigits(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 to array
1315 TIter next(loaders);
1317 while((loader = (AliLoader*)next()))
1319 loader->LoadSDigits(opt);
1324 /**************************************************************************/
1326 Int_t AliRunLoader::LoadDigits(Option_t* detectors,Option_t* opt)
1328 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1333 const char* oAll = strstr(detectors,"all");
1340 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1341 loaders = &arr;//get the pointer array
1344 TIter next(loaders);
1346 while((loader = (AliLoader*)next()))
1348 loader->LoadDigits(opt);
1352 /**************************************************************************/
1354 Int_t AliRunLoader::LoadRecPoints(Option_t* detectors,Option_t* opt)
1356 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1361 const char* oAll = strstr(detectors,"all");
1368 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1369 loaders = &arr;//get the pointer array
1372 TIter next(loaders);
1374 while((loader = (AliLoader*)next()))
1376 loader->LoadRecPoints(opt);
1380 /**************************************************************************/
1382 Int_t AliRunLoader::LoadRecParticles(Option_t* detectors,Option_t* opt)
1384 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1389 const char* oAll = strstr(detectors,"all");
1396 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1397 loaders = &arr;//get the pointer array
1400 TIter next(loaders);
1402 while((loader = (AliLoader*)next()))
1404 loader->LoadRecParticles(opt);
1408 /**************************************************************************/
1410 Int_t AliRunLoader::LoadTracks(Option_t* detectors,Option_t* opt)
1412 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1417 const char* oAll = strstr(detectors,"all");
1424 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1425 loaders = &arr;//get the pointer array
1428 TIter next(loaders);
1430 while((loader = (AliLoader*)next()))
1432 loader->LoadTracks(opt);
1436 /**************************************************************************/
1438 void AliRunLoader::UnloadHits(Option_t* detectors)
1440 //unloads hits for detectors specified in parameter
1444 const char* oAll = strstr(detectors,"all");
1451 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1452 loaders = &arr;//get the pointer to array
1455 TIter next(loaders);
1457 while((loader = (AliLoader*)next()))
1459 loader->UnloadHits();
1462 /**************************************************************************/
1464 void AliRunLoader::UnloadSDigits(Option_t* detectors)
1466 //unloads SDigits for detectors specified in parameter
1470 const char* oAll = strstr(detectors,"all");
1477 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1478 loaders = &arr;//get the pointer to array
1481 TIter next(loaders);
1483 while((loader = (AliLoader*)next()))
1485 loader->UnloadSDigits();
1488 /**************************************************************************/
1490 void AliRunLoader::UnloadDigits(Option_t* detectors)
1492 //unloads Digits for detectors specified in parameter
1496 const char* oAll = strstr(detectors,"all");
1503 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1504 loaders = &arr;//get the pointer to array
1507 TIter next(loaders);
1509 while((loader = (AliLoader*)next()))
1511 loader->UnloadDigits();
1514 /**************************************************************************/
1516 void AliRunLoader::UnloadRecPoints(Option_t* detectors)
1518 //unloads RecPoints for detectors specified in parameter
1522 const char* oAll = strstr(detectors,"all");
1529 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1530 loaders = &arr;//get the pointer to array
1533 TIter next(loaders);
1535 while((loader = (AliLoader*)next()))
1537 loader->UnloadRecPoints();
1540 /**************************************************************************/
1542 void AliRunLoader::UnloadAll(Option_t* detectors)
1544 //calls UnloadAll for detectors names specified in parameter
1545 // option "all" passed can be passed
1549 const char* oAll = strstr(detectors,"all");
1556 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1557 loaders = &arr;//get the pointer to array
1560 TIter next(loaders);
1562 while((loader = (AliLoader*)next()))
1564 loader->UnloadAll();
1567 /**************************************************************************/
1569 void AliRunLoader::UnloadTracks(Option_t* detectors)
1571 //unloads Tracks for detectors specified in parameter
1575 const char* oAll = strstr(detectors,"all");
1582 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1583 loaders = &arr;//get the pointer to array
1586 TIter next(loaders);
1588 while((loader = (AliLoader*)next()))
1590 loader->UnloadTracks();
1593 /**************************************************************************/
1595 void AliRunLoader::UnloadRecParticles(Option_t* detectors)
1597 //unloads Particles for detectors specified in parameter
1601 const char* oAll = strstr(detectors,"all");
1608 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1609 loaders = &arr;//get the pointer to array
1612 TIter next(loaders);
1614 while((loader = (AliLoader*)next()))
1616 loader->UnloadRecParticles();
1619 /**************************************************************************/
1621 AliRunLoader* AliRunLoader::GetRunLoader(const char* eventfoldername)
1623 //returns RunLoader from folder named eventfoldername
1624 TFolder* evfold= dynamic_cast<TFolder*>(AliConfig::Instance()->GetTopFolder()->FindObject(eventfoldername));
1629 AliRunLoader* runget = dynamic_cast<AliRunLoader*>(evfold->FindObject(AliRunLoader::fgkRunLoaderName));
1633 /**************************************************************************/
1635 AliLoader* AliRunLoader::GetDetectorLoader(const char* detname, const char* eventfoldername)
1637 //get the loader of the detector with the given name from the global
1639 AliRunLoader* runLoader = GetRunLoader(eventfoldername);
1641 AliErrorClass("No run loader found");
1644 return runLoader->GetDetectorLoader(detname);
1646 /**************************************************************************/
1648 AliLoader* AliRunLoader::GetDetectorLoader(const char* detname)
1650 //get the loader of the detector with the given name from the global
1653 char loadername[256];
1654 sprintf(loadername, "%sLoader", detname);
1655 AliLoader* loader = GetLoader(loadername);
1657 AliError(Form("No loader for %s found", detname));
1662 /**************************************************************************/
1664 TTree* AliRunLoader::GetTreeH(const char* detname, Bool_t maketree, const char* eventfoldername)
1666 //get the tree with hits of the detector with the given name
1667 //if maketree is true and the tree does not exist, the tree is created
1668 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1669 if (!loader) return NULL;
1670 if (!loader->TreeH() && maketree) loader->MakeTree("H");
1671 return loader->TreeH();
1674 /**************************************************************************/
1676 TTree* AliRunLoader::GetTreeH(const char* detname, Bool_t maketree)
1678 //get the tree with hits of the detector with the given name
1679 //if maketree is true and the tree does not exist, the tree is created
1680 AliLoader* loader = GetDetectorLoader(detname);
1681 if (!loader) return NULL;
1682 if (!loader->TreeH() && maketree) loader->MakeTree("H");
1683 return loader->TreeH();
1685 /**************************************************************************/
1687 TTree* AliRunLoader::GetTreeS(const char* detname, Bool_t maketree,const char* eventfoldername)
1689 //get the tree with summable digits of the detector with the given name
1690 //if maketree is true and the tree does not exist, the tree is created
1691 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1692 if (!loader) return NULL;
1693 if (!loader->TreeS() && maketree) loader->MakeTree("S");
1694 return loader->TreeS();
1696 /**************************************************************************/
1698 TTree* AliRunLoader::GetTreeS(const char* detname, Bool_t maketree)
1700 //get the tree with summable digits of the detector with the given name
1701 //if maketree is true and the tree does not exist, the tree is created
1702 AliLoader* loader = GetDetectorLoader(detname);
1703 if (!loader) return NULL;
1704 if (!loader->TreeS() && maketree) loader->MakeTree("S");
1705 return loader->TreeS();
1707 /**************************************************************************/
1709 TTree* AliRunLoader::GetTreeD(const char* detname, Bool_t maketree,const char* eventfoldername)
1711 //get the tree with digits of the detector with the given name
1712 //if maketree is true and the tree does not exist, the tree is created
1713 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1714 if (!loader) return NULL;
1715 if (!loader->TreeD() && maketree) loader->MakeTree("D");
1716 return loader->TreeD();
1718 /**************************************************************************/
1720 TTree* AliRunLoader::GetTreeD(const char* detname, Bool_t maketree)
1722 //get the tree with digits of the detector with the given name
1723 //if maketree is true and the tree does not exist, the tree is created
1724 AliLoader* loader = GetDetectorLoader(detname);
1725 if (!loader) return NULL;
1726 if (!loader->TreeD() && maketree) loader->MakeTree("D");
1727 return loader->TreeD();
1729 /**************************************************************************/
1730 TTree* AliRunLoader::GetTreeR(const char* detname, Bool_t maketree,const char* eventfoldername)
1732 //get the tree with clusters of the detector with the given name
1733 //if maketree is true and the tree does not exist, the tree is created
1734 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1735 if (!loader) return NULL;
1736 if (!loader->TreeR() && maketree) loader->MakeTree("R");
1737 return loader->TreeR();
1739 /**************************************************************************/
1741 TTree* AliRunLoader::GetTreeR(const char* detname, Bool_t maketree)
1743 //get the tree with clusters of the detector with the given name
1744 //if maketree is true and the tree does not exist, the tree is created
1745 AliLoader* loader = GetDetectorLoader(detname);
1746 if (!loader) return NULL;
1747 if (!loader->TreeR() && maketree) loader->MakeTree("R");
1748 return loader->TreeR();
1750 /**************************************************************************/
1752 TTree* AliRunLoader::GetTreeT(const char* detname, Bool_t maketree,const char* eventfoldername)
1754 //get the tree with tracks of the detector with the given name
1755 //if maketree is true and the tree does not exist, the tree is created
1756 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1757 if (!loader) return NULL;
1758 if (!loader->TreeT() && maketree) loader->MakeTree("T");
1759 return loader->TreeT();
1761 /**************************************************************************/
1763 TTree* AliRunLoader::GetTreeT(const char* detname, Bool_t maketree)
1765 //get the tree with tracks of the detector with the given name
1766 //if maketree is true and the tree does not exist, the tree is created
1767 AliLoader* loader = GetDetectorLoader(detname);
1768 if (!loader) return NULL;
1769 if (!loader->TreeT() && maketree) loader->MakeTree("T");
1770 return loader->TreeT();
1772 /**************************************************************************/
1774 TTree* AliRunLoader::GetTreeP(const char* detname, Bool_t maketree,const char* eventfoldername)
1776 //get the tree with particles of the detector with the given name
1777 //if maketree is true and the tree does not exist, the tree is created
1778 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1779 if (!loader) return NULL;
1780 if (!loader->TreeP() && maketree) loader->MakeTree("P");
1781 return loader->TreeP();
1783 /**************************************************************************/
1785 TTree* AliRunLoader::GetTreeP(const char* detname, Bool_t maketree)
1787 //get the tree with particles of the detector with the given name
1788 //if maketree is true and the tree does not exist, the tree is created
1789 AliLoader* loader = GetDetectorLoader(detname);
1790 if (!loader) return NULL;
1791 if (!loader->TreeP() && maketree) loader->MakeTree("P");
1792 return loader->TreeP();
1795 /**************************************************************************/
1797 void AliRunLoader::CdGAFile()
1799 //sets gDirectory to galice file
1801 if(fGAFile) fGAFile->cd();
1804 /**************************************************************************/
1806 void AliRunLoader::GetListOfDetectors(const char * namelist,TObjArray& pointerarray) const
1808 //this method looks for all Loaders corresponding
1809 //to names (many) specified in namelist i.e. namelist ("ITS TPC TRD")
1813 strcpy(dets,namelist);//compiler cries when char* = const Option_t*;
1814 dets[strlen(dets)+1] = '\n';//set endl at the end of string
1819 tmp = sscanf(pdet,"%s",buff);//read the string from the input string pdet into buff
1820 if ( (buff[0] == 0) || (tmp == 0) ) break; //if not read
1822 pdet = strstr(pdet,buff) + strlen(buff);//move the input pointer about number of bytes (letters) read
1823 //I am aware that is a little complicated. I don't know the number of spaces between detector names
1824 //so I read the string, than I find where it starts (strstr) and move the pointer about length of a string
1825 //If there is a better way, please write me (Piotr.Skowronski@cern.ch)
1826 //construct the Loader name
1827 TString getname(buff);
1829 AliLoader* loader = GetLoader(getname);//get the Loader
1832 pointerarray.Add(loader);
1836 AliError(Form("Can not find Loader for %s",buff));
1842 /*****************************************************************************/
1844 void AliRunLoader::Clean(const TString& name)
1846 //removes object with given name from event folder and deletes it
1847 if (GetEventFolder() == 0x0) return;
1848 TObject* obj = GetEventFolder()->FindObject(name);
1851 AliDebug(1, Form("name=%s, cleaning %s.",GetName(),name.Data()));
1852 GetEventFolder()->Remove(obj);
1858 /*****************************************************************************/
1860 TTask* AliRunLoader::GetRunDigitizer()
1862 //returns Run Digitizer from folder
1864 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1865 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetDigitizerTaskName());
1866 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1868 /*****************************************************************************/
1870 TTask* AliRunLoader::GetRunSDigitizer()
1872 //returns SDigitizer Task from folder
1874 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1875 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetSDigitizerTaskName());
1876 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1878 /*****************************************************************************/
1880 TTask* AliRunLoader::GetRunReconstructioner()
1882 //returns Reconstructioner Task from folder
1883 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1884 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetReconstructionerTaskName());
1885 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1887 /*****************************************************************************/
1889 TTask* AliRunLoader::GetRunTracker()
1891 //returns Tracker Task from folder
1892 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1893 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetTrackerTaskName());
1894 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1896 /*****************************************************************************/
1898 TTask* AliRunLoader::GetRunPIDTask()
1900 //returns Tracker Task from folder
1901 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1902 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetPIDTaskName());
1903 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1905 /*****************************************************************************/
1907 TTask* AliRunLoader::GetRunQATask()
1909 //returns Quality Assurance Task from folder
1910 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1913 AliErrorClass("Can not get task folder from AliConfig");
1916 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetQATaskName());
1917 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1920 /*****************************************************************************/
1922 void AliRunLoader::SetCompressionLevel(Int_t cl)
1924 //Sets Compression Level in all files
1925 if (fGAFile) fGAFile->SetCompressionLevel(cl);
1926 SetKineComprLevel(cl);
1927 SetTrackRefsComprLevel(cl);
1928 TIter next(fLoaders);
1930 while((loader = (AliLoader*)next()))
1932 loader->SetCompressionLevel(cl);
1935 /**************************************************************************/
1937 void AliRunLoader::SetKineComprLevel(Int_t cl)
1939 //Sets comression level in Kine File
1940 fKineDataLoader->SetCompressionLevel(cl);
1942 /**************************************************************************/
1944 void AliRunLoader::SetTrackRefsComprLevel(Int_t cl)
1946 //Sets comression level in Track Refences File
1947 fTrackRefsDataLoader->SetCompressionLevel(cl);
1949 /**************************************************************************/
1951 void AliRunLoader::UnloadHeader()
1953 //removes TreeE from folder and deletes it
1954 // as well as fHeader object
1959 /**************************************************************************/
1961 void AliRunLoader::UnloadTrigger()
1963 //removes TreeCT from folder and deletes it
1964 // as well as fHeader object
1970 /**************************************************************************/
1972 void AliRunLoader::UnloadKinematics()
1974 //Unloads Kinematics
1975 fKineDataLoader->GetBaseLoader(0)->Unload();
1977 /**************************************************************************/
1979 void AliRunLoader::UnloadTrackRefs()
1981 //Unloads Track Refernces
1982 fTrackRefsDataLoader->GetBaseLoader(0)->Unload();
1984 /**************************************************************************/
1986 void AliRunLoader::UnloadgAlice()
1989 if (gAlice == GetAliRun())
1991 AliDebug(1, "Set gAlice = 0x0");
1992 gAlice = 0x0;//if gAlice is the same that in folder (to be deleted by the way of folder)
1994 AliRun* alirun = GetAliRun();
1995 if (GetEventFolder()) GetEventFolder()->Remove(alirun);
1998 /**************************************************************************/
2000 void AliRunLoader::MakeTrackRefsContainer()
2002 // Makes a tree for Track References
2003 fTrackRefsDataLoader->MakeTree();
2005 /**************************************************************************/
2007 Int_t AliRunLoader::LoadTrackRefs(Option_t* option)
2009 //Load track references from file (opens file and posts tree to folder)
2011 return fTrackRefsDataLoader->GetBaseLoader(0)->Load(option);
2013 /**************************************************************************/
2015 void AliRunLoader::SetDirName(TString& dirname)
2017 //sets directory name
2018 if (dirname.IsNull()) return;
2019 fUnixDirName = dirname;
2020 fKineDataLoader->SetDirName(dirname);
2021 fTrackRefsDataLoader->SetDirName(dirname);
2023 TIter next(fLoaders);
2025 while((loader = (AliLoader*)next()))
2027 loader->SetDirName(dirname);
2031 /*****************************************************************************/
2033 Int_t AliRunLoader::GetFileOffset() const
2035 //returns the file number that is added to the file name for current event
2036 return Int_t(fCurrentEvent/fNEventsPerFile);
2039 /*****************************************************************************/
2040 const TString AliRunLoader::SetFileOffset(const TString& fname)
2042 //adds the the number to the file name at proper place for current event
2043 Long_t offset = (Long_t)GetFileOffset();
2044 if (offset < 1) return fname;
2046 soffset += offset;//automatic conversion to string
2047 TString dotroot(".root");
2048 const TString& offfsetdotroot = offset + dotroot;
2049 TString out = fname;
2050 out = out.ReplaceAll(dotroot,offfsetdotroot);
2051 AliDebug(1, Form(" in=%s out=%s",fname.Data(),out.Data()));
2054 /*****************************************************************************/
2056 void AliRunLoader::SetDigitsFileNameSuffix(const TString& suffix)
2058 //adds the suffix before ".root",
2059 //e.g. TPC.Digits.root -> TPC.DigitsMerged.root
2060 //made on Jiri Chudoba demand
2062 TIter next(fLoaders);
2064 while((loader = (AliLoader*)next()))
2066 loader->SetDigitsFileNameSuffix(suffix);
2069 /*****************************************************************************/
2071 TString AliRunLoader::GetFileName() const
2073 //returns name of galice file
2075 if (fGAFile == 0x0) return result;
2076 result = fGAFile->GetName();
2079 /*****************************************************************************/
2081 void AliRunLoader::SetDetectorAddresses()
2083 //calls SetTreeAddress for all detectors
2084 if (GetAliRun()==0x0) return;
2085 TIter next(GetAliRun()->Modules());
2087 while((mod = (AliModule*)next()))
2089 AliDetector* det = dynamic_cast<AliDetector*>(mod);
2090 if (det) det->SetTreeAddress();
2093 /*****************************************************************************/
2095 void AliRunLoader::Synchronize()
2097 //synchrinizes all writtable files
2098 TIter next(fLoaders);
2100 while((loader = (AliLoader*)next()))
2102 loader->Synchronize();
2105 fKineDataLoader->Synchronize();
2106 fTrackRefsDataLoader->Synchronize();
2108 TFile* file = gROOT->GetFile( gSystem->ConcatFileName( fUnixDirName.Data(), fgkDefaultTriggerFileName.Data() ) ) ;
2109 if( file ) file->Flush();
2111 if (fGAFile) fGAFile->Flush();
2113 /*****************************************************************************/
2114 /*****************************************************************************/