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),
162 fKineDataLoader(new AliDataLoader(fgkDefaultKineFileName,fgkKineContainerName,"Kinematics")),
163 fTrackRefsDataLoader(new AliDataLoader(fgkDefaultTrackRefsFileName,fgkTrackRefsContainerName,"Track References")),
170 TString errmsg("Parameter is NULL");
171 AliError(errmsg.Data());
176 TObject* obj = fEventFolder->FindObject(fgkRunLoaderName);
178 { //if it is, then sth. is going wrong... exits aliroot session
179 TString errmsg("In Event Folder Named ");
180 errmsg+=fEventFolder->GetName();
181 errmsg+=" object named "+fgkRunLoaderName+" already exists. I am confused ...";
183 AliError(errmsg.Data());
185 return;//never reached
188 if (!fgRunLoader) fgRunLoader = this;
190 fEventFolder->Add(this);//put myself to the folder to accessible for all
194 /**************************************************************************/
196 Int_t AliRunLoader::GetEvent(Int_t evno)
198 //Gets event number evno
199 //Reloads all data properly
200 //PH if (fCurrentEvent == evno) return 0;
204 AliError("Can not give the event with negative number");
208 if (evno >= GetNumberOfEvents())
210 AliError(Form("There is no event with number %d",evno));
214 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
215 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
216 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
217 AliDebug(1, Form(" GETTING EVENT %d",evno));
218 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
219 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
220 AliDebug(1, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>");
222 fCurrentEvent = evno;
226 //Reload header (If header was loaded)
229 retval = TreeE()->GetEvent(fCurrentEvent);
232 AliError(Form("Cannot find event: %d\n ",fCurrentEvent));
236 //Reload stack (If header was loaded)
237 if (TreeE()) fStack = GetHeader()->Stack();
238 //Set event folder in stack (it does not mean that we read kinematics from file)
239 if( GetTrigger() && TreeCT() ) {
240 retval = TreeCT()->GetEvent(fCurrentEvent);
242 AliError(Form("Error occured while GetEvent for Trigger. Event %d",evno));
250 AliError(Form("Error occured while setting event %d",evno));
254 //Post Track References
255 retval = fTrackRefsDataLoader->GetEvent();
258 AliError(Form("Error occured while GetEvent for Track References. Event %d",evno));
262 //Read Kinematics if loaded
263 retval = fKineDataLoader->GetEvent();
266 AliError(Form("Error occured while GetEvent for Kinematics. Event %d",evno));
270 if (fStack && fKineDataLoader->GetBaseLoader(0)->IsLoaded())
272 fStack->ConnectTree(TreeK());
274 if (fStack->GetEvent() == kFALSE)
276 AliError(Form("Error occured while GetEvent for Stack. Event %d",evno));
281 //Trigger data reloading in all loaders
282 TIter next(fLoaders);
284 while((loader = (AliLoader*)next()))
286 retval = loader->GetEvent();
289 AliError(Form("Error occured while getting event for %s. Event %d.",
290 loader->GetDetectorName().Data(), evno));
295 SetDetectorAddresses();
299 /**************************************************************************/
300 Int_t AliRunLoader::SetEvent()
302 //if kinematics was loaded Cleans folder data
306 retval = fKineDataLoader->SetEvent();
309 AliError("SetEvent for Kinamtics Data Loader retutned error.");
312 retval = fTrackRefsDataLoader->SetEvent();
315 AliError("SetEvent for Track References Data Loader retutned error.");
319 TIter next(fLoaders);
321 while((loader = (AliLoader*)next()))
323 retval = loader->SetEvent();
326 AliError(Form("SetEvent for %s Data Loader retutned error.",loader->GetName()));
333 /**************************************************************************/
335 Int_t AliRunLoader::SetEventNumber(Int_t evno)
337 //cleans folders and sets the root dirs in files
338 if (fCurrentEvent == evno) return 0;
339 fCurrentEvent = evno;
343 /**************************************************************************/
344 AliCDBEntry* AliRunLoader::GetCDBEntry(const char* name) const
346 //Get an AliCDBEntry from the run data storage
348 if ( !(AliCDBManager::Instance()->IsDefaultStorageSet()) ) {
349 AliError("No run data storage defined!");
352 return AliCDBManager::Instance()->GetDefaultStorage()->Get(name, GetHeader()->GetRun());
356 /**************************************************************************/
357 AliRunLoader* AliRunLoader::Open
358 (const char* filename, const char* eventfoldername, Option_t* option)
360 //Opens a desired file 'filename'
361 //gets the the run-Loader and mounts it desired folder
362 //returns the pointer to run Loader which can be further used for accessing data
363 //in case of error returns NULL
365 static const TString kwebaddress("http://alisoft.cern.ch/people/skowron/codedoc/split/index.html");
366 AliDebugClass(1,Form("\n\n\nNew I/O strcture: See more info:\n %s\n\n\n",kwebaddress.Data()));
368 AliRunLoader* result = 0x0;
370 /* ************************************************ */
371 /* Chceck if folder with given name already exists */
372 /* ************************************************ */
374 TObject* obj = AliConfig::Instance()->GetTopFolder()->FindObject(eventfoldername);
377 TFolder* fold = dynamic_cast<TFolder*>(obj);
380 AliErrorClass("Such a obejct already exists in top alice folder and it is not a folder.");
384 //check if we can get RL from that folder
385 result = AliRunLoader::GetRunLoader(eventfoldername);
388 AliErrorClass(Form("Folder %s already exists, and can not find session there. Can not mount.",eventfoldername));
392 if (result->GetFileName().CompareTo(filename) != 0)
394 AliErrorClass("Other file is mounted in demanded folder. Can not mount.");
398 //check if now is demanded (re)creation
399 if ( AliLoader::TestFileOption(option) == kFALSE)
401 AliErrorClass(Form("Session already exists in folder %s and this session option is %s. Unable to proceed.",
402 eventfoldername,option));
406 //check if demanded option is update and existing one
407 TString tmpstr(option);
408 if ( (tmpstr.CompareTo("update",TString::kIgnoreCase) == 0) &&
409 (result->fGAFile->IsWritable() == kFALSE) )
411 AliErrorClass(Form("Session already exists in folder %s and is not writable while this session option is %s. Unable to proceed.",
412 eventfoldername,option));
416 AliWarningClass("Session is already opened and mounted in demanded folder");
417 if (!fgRunLoader) fgRunLoader = result; //PH get access from any place
419 } //end of checking in case of existance of object named identically that folder session is being opened
422 TFile * gAliceFile = TFile::Open(filename,option);//open a file
424 {//null pointer returned
425 AliFatalClass(Form("Can not open file %s.",filename));
429 if (gAliceFile->IsOpen() == kFALSE)
430 {//pointer to valid object returned but file is not opened
431 AliErrorClass(Form("Can not open file %s.",filename));
435 //if file is "read" or "update" than we try to find AliRunLoader there - if not found cry and exit
436 //else create new AliRunLoader
437 if ( AliLoader::TestFileOption(option) )
439 AliDebugClass(1, "Reading RL from file");
441 result = dynamic_cast<AliRunLoader*>(gAliceFile->Get(fgkRunLoaderName));//get the run Loader from the file
444 AliErrorClass(Form("Can not find run-Loader in file %s.",filename));
445 delete gAliceFile;//close the file
448 Int_t tmp = result->SetEventFolderName(eventfoldername);//mount a event folder
449 if (tmp)//if SetEvent returned error
451 AliErrorClass(Form("Can not mount event in folder %s.",eventfoldername));
452 delete result; //delete run-Loader
453 delete gAliceFile;//close the file
459 AliDebugClass(1, Form("Creating new AliRunLoader. Folder name is %s",eventfoldername));
462 result = new AliRunLoader(eventfoldername);
464 catch (TString& errmsg)
466 AliErrorClass(Form("AliRunLoader constrcutor has thrown exception: %s\n",errmsg.Data()));
468 delete gAliceFile;//close the file
473 //procedure for extracting dir name from the file name
474 TString fname(filename);
475 Int_t nsl = fname.Last('#');//look for hash in file name
477 if (nsl < 0) {//hash not found
478 nsl = fname.Last('/');// look for slash
480 nsl = fname.Last(':');// look for colon e.g. rfio:galice.root
483 if (nsl < 0) dirname = "./"; // no directory path, use "."
484 else dirname = fname.Remove(nsl+1);// directory path
486 AliDebugClass(1, Form("Dir name is : %s",dirname.Data()));
488 result->SetDirName(dirname);
489 result->SetGAliceFile(gAliceFile);//set the pointer to gAliceFile
490 if (!fgRunLoader) fgRunLoader = result; //PH get access from any place
493 /**************************************************************************/
494 Int_t AliRunLoader::GetNumberOfEvents()
496 //returns number of events in Run
500 retval = LoadHeader();
503 AliError("Error occured while loading header");
507 return (Int_t)TreeE()->GetEntries();
509 /**************************************************************************/
511 void AliRunLoader::MakeHeader()
513 //Makes header and connects it to header tree (if it exists)
517 AliDebug(1, "Creating new Header Object");
518 fHeader= new AliHeader();
520 TTree* tree = TreeE();
523 AliDebug(1, "Got Tree from folder.");
524 TBranch* branch = tree->GetBranch(fgkHeaderBranchName);
527 AliDebug(1, "Creating new branch");
528 branch = tree->Branch(fgkHeaderBranchName, "AliHeader", &fHeader, 4000, 0);
529 branch->SetAutoDelete(kFALSE);
533 AliDebug(1, "Got Branch from Tree");
534 branch->SetAddress(&fHeader);
535 tree->GetEvent(fCurrentEvent);
536 fStack = fHeader->Stack(); //should be safe - if we created Stack, header returns pointer to the same object
540 fStack->ConnectTree(TreeK());
546 AliDebug(1, "Header does not have a stack.");
550 AliDebug(1, "Exiting MakeHeader method");
552 /**************************************************************************/
554 void AliRunLoader::MakeStack()
556 //Creates the stack object - do not connect the tree
559 fStack = new AliStack(10000);
562 /**************************************************************************/
564 void AliRunLoader::MakeTrigger()
566 // Makes trigger object and connects it to trigger tree (if it exists)
568 if( fCTrigger == 0x0 ) {
569 AliDebug( 1, "Creating new Trigger Object" );
570 fCTrigger = new AliCentralTrigger();
572 TTree* tree = TreeCT();
574 fCTrigger->MakeBranch( fgkTriggerBranchName, tree );
575 tree->GetEvent( fCurrentEvent );
578 AliDebug( 1, "Exiting MakeTrigger method" );
580 /**************************************************************************/
582 void AliRunLoader::MakeTree(Option_t *option)
585 const char *oK = strstr(option,"K"); //Kine
586 const char *oE = strstr(option,"E"); //Header
587 const char *oGG = strstr(option,"GG"); //Central TriGGer
591 if (fKineDataLoader->GetBaseLoader(0)->IsLoaded() == kFALSE)
593 AliError("Load Kinematics first");
598 fKineDataLoader->MakeTree();
601 fStack->ConnectTree(TreeK());
602 WriteKinematics("OVERWRITE");
609 TTree* tree = new TTree(fgkHeaderContainerName,"Tree with Headers");
610 GetEventFolder()->Add(tree);
612 WriteHeader("OVERWRITE");
617 // create the CTP Trigger output file and tree
618 TFile* file = gROOT->GetFile( fgkDefaultTriggerFileName );
620 file = TFile::Open( gSystem->ConcatFileName( fUnixDirName.Data(), fgkDefaultTriggerFileName.Data() ), "RECREATE" ) ;
624 TTree* tree = new TTree( fgkTriggerContainerName, "Tree with Central Trigger Mask" );
625 GetEventFolder()->Add(tree);
627 // WriteHeader("OVERWRITE");
630 TIter next(fLoaders);
632 while((loader = (AliLoader*)next()))
634 loader->MakeTree(option);
638 /**************************************************************************/
640 Int_t AliRunLoader::LoadgAlice()
642 //Loads gAlice from file
645 AliWarning("AliRun is already in folder. Unload first.");
648 AliRun* alirun = dynamic_cast<AliRun*>(fGAFile->Get(fgkGAliceName));
651 AliError(Form("Can not find gAlice in file %s",fGAFile->GetName()));
654 alirun->SetRunLoader(this);
657 AliWarning(Form("gAlice already exists. Putting retrived object in folder named %s",
658 GetEventFolder()->GetName()));
664 SetDetectorAddresses();//calls SetTreeAddress for all detectors
667 /**************************************************************************/
669 Int_t AliRunLoader::LoadHeader()
671 //loads treeE and reads header object for current event
674 AliWarning("Header is already loaded. Use ReloadHeader to force reload. Nothing done");
678 if (GetEventFolder() == 0x0)
680 AliError("Event folder not specified yet");
686 AliError("Session not opened. Use AliRunLoader::Open");
690 if (fGAFile->IsOpen() == kFALSE)
692 AliError("Session not opened. Use AliRunLoader::Open");
696 TTree* tree = dynamic_cast<TTree*>(fGAFile->Get(fgkHeaderContainerName));
699 AliError(Form("Can not find header tree named %s in file %s",
700 fgkHeaderContainerName.Data(),fGAFile->GetName()));
704 if (tree == TreeE()) return 0;
707 GetEventFolder()->Add(tree);
708 MakeHeader();//creates header object and connects to tree
712 /**************************************************************************/
714 Int_t AliRunLoader::LoadTrigger(Option_t* option)
719 AliWarning("Trigger is already loaded. Nothing done");
723 if( GetEventFolder() == 0x0 ) {
724 AliError("Event folder not specified yet");
727 // get the CTP Trigger output file and tree
728 TString trgfile = gSystem->ConcatFileName( fUnixDirName.Data(),
729 fgkDefaultTriggerFileName.Data() );
730 TFile* file = gROOT->GetFile( trgfile );
732 file = TFile::Open( trgfile, option ) ;
733 if (!file || file->IsOpen() == kFALSE ) {
734 AliError( Form( "Can not open trigger file %s", trgfile.Data() ) );
740 TTree* tree = dynamic_cast<TTree*>(file->Get( fgkTriggerContainerName ));
742 AliError( Form( "Can not find trigger tree named %s in file %s",
743 fgkTriggerContainerName.Data(), file->GetName() ) );
749 fCTrigger = dynamic_cast<AliCentralTrigger*>(file->Get( "AliCentralTrigger" ));
750 GetEventFolder()->Add( tree );
756 /**************************************************************************/
758 Int_t AliRunLoader::LoadKinematics(Option_t* option)
760 //Loads the kinematics
761 Int_t retval = fKineDataLoader->GetBaseLoader(0)->Load(option);
764 AliError("Error occured while loading kinamatics tree.");
769 fStack->ConnectTree(TreeK());
770 retval = fStack->GetEvent();
771 if ( retval == kFALSE)
773 AliError("Error occured while loading kinamatics tree.");
780 /**************************************************************************/
782 Int_t AliRunLoader::OpenDataFile(const TString& filename,TFile*& file,TDirectory*& dir,Option_t* opt,Int_t cl)
784 //Opens File with kinematics
787 if (file->IsOpen() == kFALSE)
788 {//pointer is not null but file is not opened
789 AliWarning("Pointer to file is not null, but file is not opened");//risky any way
791 file = 0x0; //proceed with opening procedure
795 AliWarning(Form("File %s already opened",filename.Data()));
799 //try to find if that file is opened somewere else
800 file = (TFile *)( gROOT->GetListOfFiles()->FindObject(filename) );
803 if(file->IsOpen() == kTRUE)
805 AliWarning(Form("File %s already opened by sombody else.",file->GetName()));
810 file = TFile::Open(filename,opt);
813 AliError(Form("Can not open file %s",filename.Data()));
816 if (file->IsOpen() == kFALSE)
817 {//file is not opened
818 AliError(Form("Can not open file %s",filename.Data()));
822 file->SetCompressionLevel(cl);
824 dir = AliLoader::ChangeDir(file,fCurrentEvent);
827 AliError(Form("Can not change to root directory in file %s",filename.Data()));
832 /**************************************************************************/
834 TTree* AliRunLoader::TreeE() const
836 //returns the tree from folder; shortcut method
837 if (AliDebugLevel() > 10) fEventFolder->ls();
838 TObject *obj = fEventFolder->FindObject(fgkHeaderContainerName);
839 return (obj)?dynamic_cast<TTree*>(obj):0x0;
841 /**************************************************************************/
843 TTree* AliRunLoader::TreeCT() const
845 //returns the tree from folder; shortcut method
846 if (AliDebugLevel() > 10) fEventFolder->ls();
847 TObject *obj = fEventFolder->FindObject(fgkTriggerContainerName);
848 return (obj)?dynamic_cast<TTree*>(obj):0x0;
850 /**************************************************************************/
852 AliHeader* AliRunLoader::GetHeader() const
854 //returns pointer header object
857 /**************************************************************************/
859 AliCentralTrigger* AliRunLoader::GetTrigger() const
861 //returns pointer trigger object
865 /**************************************************************************/
867 TTree* AliRunLoader::TreeK() const
869 //returns the tree from folder; shortcut method
870 TObject *obj = GetEventFolder()->FindObject(fgkKineContainerName);
871 return (obj)?dynamic_cast<TTree*>(obj):0x0;
873 /**************************************************************************/
875 TTree* AliRunLoader::TreeTR() const
877 //returns the tree from folder; shortcut method
878 TObject* obj = GetEventFolder()->FindObject(fgkTrackRefsContainerName);
879 return (obj)?dynamic_cast<TTree*>(obj):0x0;
881 /**************************************************************************/
883 AliRun* AliRunLoader::GetAliRun() const
885 //returns AliRun which sits in the folder
886 if (fEventFolder == 0x0) return 0x0;
887 TObject *obj = fEventFolder->FindObject(fgkGAliceName);
888 return (obj)?dynamic_cast<AliRun*>(obj):0x0;
890 /**************************************************************************/
892 Int_t AliRunLoader::WriteHeader(Option_t* opt)
895 AliDebug(1, "WRITING HEADER");
897 TTree* tree = TreeE();
900 AliWarning("Can not find Header Tree in Folder");
903 if (fGAFile->IsWritable() == kFALSE)
905 AliError(Form("File %s is not writable",fGAFile->GetName()));
909 TObject* obj = fGAFile->Get(fgkHeaderContainerName);
911 { //if they exist, see if option OVERWRITE is used
913 if(tmp.Contains("OVERWRITE",TString::kIgnoreCase) == 0)
914 {//if it is not used - give an error message and return an error code
915 AliError("Tree already exisists. Use option \"OVERWRITE\" to overwrite previous data");
920 tree->SetDirectory(fGAFile);
921 tree->Write(0,TObject::kOverwrite);
923 AliDebug(1, "WRITTEN\n\n");
928 /**************************************************************************/
930 Int_t AliRunLoader::WriteTrigger(Option_t* opt)
933 AliDebug( 1, "WRITING TRIGGER" );
935 TTree* tree = TreeCT();
937 AliWarning("Can not find Trigger Tree in Folder");
941 TFile* file = gROOT->GetFile( gSystem->ConcatFileName( fUnixDirName.Data(), fgkDefaultTriggerFileName.Data() ) ) ;
942 if( !file || !file->IsOpen() ) {
943 AliError( "can't write Trigger, file is not open" );
947 TObject* obj = file->Get( fgkTriggerContainerName );
948 if( obj ) { //if they exist, see if option OVERWRITE is used
950 if( tmp.Contains( "OVERWRITE", TString::kIgnoreCase ) == 0) {
951 //if it is not used - give an error message and return an error code
952 AliError( "Tree already exisists. Use option \"OVERWRITE\" to overwrite previous data" );
957 fCTrigger->Write( 0, TObject::kOverwrite );
958 tree->Write( 0, TObject::kOverwrite );
961 AliDebug(1, "WRITTEN\n\n");
965 /**************************************************************************/
967 Int_t AliRunLoader::WriteAliRun(Option_t* /*opt*/)
969 //writes AliRun object to the file
971 if (GetAliRun()) GetAliRun()->Write();
974 /**************************************************************************/
976 Int_t AliRunLoader::WriteKinematics(Option_t* opt)
979 return fKineDataLoader->GetBaseLoader(0)->WriteData(opt);
981 /**************************************************************************/
982 Int_t AliRunLoader::WriteTrackRefs(Option_t* opt)
984 //writes Track References tree
985 return fTrackRefsDataLoader->GetBaseLoader(0)->WriteData(opt);
987 /**************************************************************************/
989 Int_t AliRunLoader::WriteHits(Option_t* opt)
991 //Calls WriteHits for all loaders
994 TIter next(fLoaders);
996 while((loader = (AliLoader*)next()))
998 res = loader->WriteHits(opt);
1001 AliError(Form("Failed to write hits for %s (%d)",loader->GetDetectorName().Data(),res));
1007 /**************************************************************************/
1009 Int_t AliRunLoader::WriteSDigits(Option_t* opt)
1011 //Calls WriteSDigits for all loaders
1014 TIter next(fLoaders);
1016 while((loader = (AliLoader*)next()))
1018 res = loader->WriteSDigits(opt);
1021 AliError(Form("Failed to write summable digits for %s.",loader->GetDetectorName().Data()));
1027 /**************************************************************************/
1029 Int_t AliRunLoader::WriteDigits(Option_t* opt)
1031 //Calls WriteDigits for all loaders
1034 TIter next(fLoaders);
1036 while((loader = (AliLoader*)next()))
1038 res = loader->WriteDigits(opt);
1041 AliError(Form("Failed to write digits for %s.",loader->GetDetectorName().Data()));
1047 /**************************************************************************/
1049 Int_t AliRunLoader::WriteRecPoints(Option_t* opt)
1051 //Calls WriteRecPoints for all loaders
1054 TIter next(fLoaders);
1056 while((loader = (AliLoader*)next()))
1058 res = loader->WriteRecPoints(opt);
1061 AliError(Form("Failed to write Reconstructed Points for %s.",
1062 loader->GetDetectorName().Data()));
1068 /**************************************************************************/
1070 Int_t AliRunLoader::WriteTracks(Option_t* opt)
1072 //Calls WriteTracks for all loaders
1075 TIter next(fLoaders);
1077 while((loader = (AliLoader*)next()))
1079 res = loader->WriteTracks(opt);
1082 AliError(Form("Failed to write Tracks for %s.",
1083 loader->GetDetectorName().Data()));
1089 /**************************************************************************/
1091 Int_t AliRunLoader::WriteRunLoader(Option_t* /*opt*/)
1093 //Writes itself to the file
1095 this->Write(0,TObject::kOverwrite);
1098 /**************************************************************************/
1100 Int_t AliRunLoader::SetEventFolderName(const TString& name)
1102 //sets top folder name for this run; of alread
1105 AliError("Name is empty");
1109 //check if such a folder already exists - try to find it in alice top folder
1110 TObject* obj = AliConfig::Instance()->GetTopFolder()->FindObject(name);
1113 TFolder* fold = dynamic_cast<TFolder*>(obj);
1116 AliError("Such a obejct already exists in top alice folder and it is not a folder.");
1119 //folder which was found is our folder
1120 if (fEventFolder == fold)
1126 AliError("Such a folder already exists in top alice folder. Can not mount.");
1131 //event is alredy connected, just change name of the folder
1134 fEventFolder->SetName(name);
1138 if (fKineDataLoader == 0x0)
1139 fKineDataLoader = new AliDataLoader(fgkDefaultKineFileName,fgkKineContainerName,"Kinematics");
1141 if ( fTrackRefsDataLoader == 0x0)
1142 fTrackRefsDataLoader = new AliDataLoader(fgkDefaultTrackRefsFileName,fgkTrackRefsContainerName,"Track References");
1144 //build the event folder structure
1145 AliDebug(1, Form("Creating new event folder named %s",name.Data()));
1146 fEventFolder = AliConfig::Instance()->BuildEventFolder(name,"Event Folder");
1147 fEventFolder->Add(this);//put myself to the folder to accessible for all
1149 TIter next(fLoaders);
1151 while((loader = (AliLoader*)next()))
1153 loader->Register(fEventFolder);//build folder structure for this detector
1156 fKineDataLoader->SetEventFolder(GetEventFolder());
1157 fTrackRefsDataLoader->SetEventFolder(GetEventFolder());
1158 fKineDataLoader->SetFolder(GetEventFolder());
1159 fTrackRefsDataLoader->SetFolder(GetEventFolder());
1161 fEventFolder->SetOwner();
1164 /**************************************************************************/
1166 void AliRunLoader::AddLoader(AliLoader* loader)
1168 //Adds the Loader for given detector
1169 if (loader == 0x0) //if null shout and exit
1171 AliError("Parameter is NULL");
1174 loader->SetDirName(fUnixDirName);
1175 if (fEventFolder) loader->SetEventFolder(fEventFolder); //if event folder is already defined,
1176 //pass information to the Loader
1177 fLoaders->Add(loader);//add the Loader to the array
1179 /**************************************************************************/
1181 void AliRunLoader::AddLoader(AliDetector* det)
1183 //Asks module (detector) ro make a Loader and stores in the array
1184 if (det == 0x0) return;
1185 AliLoader* get = det->GetLoader();//try to get loader
1186 if (get == 0x0) get = det->MakeLoader(fEventFolder->GetName());//if did not obtain, ask to make it
1190 AliDebug(1, Form("Detector: %s Loader : %s",det->GetName(),get->GetName()));
1195 /**************************************************************************/
1197 AliLoader* AliRunLoader::GetLoader(const char* detname) const
1199 //returns loader for given detector
1200 //note that naming convention is TPCLoader not just TPC
1201 return (AliLoader*)fLoaders->FindObject(detname);
1204 /**************************************************************************/
1206 AliLoader* AliRunLoader::GetLoader(AliDetector* det) const
1208 //get loader for detector det
1209 if(det == 0x0) return 0x0;
1210 TString getname(det->GetName());
1212 AliDebug(1, Form(" Loader name is %s",getname.Data()));
1213 return GetLoader(getname);
1216 /**************************************************************************/
1218 void AliRunLoader::CleanFolders()
1220 // fEventFolder->Add(this);//put myself to the folder to accessible for all
1227 /**************************************************************************/
1229 void AliRunLoader::CleanDetectors()
1231 //Calls CleanFolders for all detectors
1232 TIter next(fLoaders);
1234 while((loader = (AliLoader*)next()))
1236 loader->CleanFolders();
1239 /**************************************************************************/
1241 void AliRunLoader::RemoveEventFolder()
1243 //remove all the tree of event
1244 //all the stuff changing EbE stays untached (PDGDB, tasks, etc.)
1246 if (fEventFolder == 0x0) return;
1247 fEventFolder->SetOwner(kFALSE);//don't we want to deleted while removing the folder that we are sitting in
1248 fEventFolder->Remove(this);//remove us drom folder
1250 AliConfig::Instance()->GetTopFolder()->SetOwner(); //brings ownership back for fEventFolder since it sits in top folder
1251 AliConfig::Instance()->GetTopFolder()->Remove(fEventFolder); //remove the event tree
1252 delete fEventFolder;
1254 /**************************************************************************/
1256 void AliRunLoader::SetGAliceFile(TFile* gafile)
1258 //sets pointer to galice.root file
1262 /**************************************************************************/
1264 Int_t AliRunLoader::LoadHits(Option_t* detectors,Option_t* opt)
1266 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1268 AliDebug(1, "Loading Hits");
1272 const char* oAll = strstr(detectors,"all");
1275 AliDebug(1, "Option is All");
1280 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1281 loaders = &arr;//get the pointer array
1284 AliDebug(1, Form("For detectors. Number of detectors chosen for loading %d",loaders->GetEntries()));
1286 TIter next(loaders);
1288 while((loader = (AliLoader*)next()))
1290 AliDebug(1, Form(" Calling LoadHits(%s) for %s",opt,loader->GetName()));
1291 loader->LoadHits(opt);
1293 AliDebug(1, "Done");
1297 /**************************************************************************/
1299 Int_t AliRunLoader::LoadSDigits(Option_t* detectors,Option_t* opt)
1301 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1306 const char* oAll = strstr(detectors,"all");
1313 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1314 loaders = &arr;//get the pointer to array
1317 TIter next(loaders);
1319 while((loader = (AliLoader*)next()))
1321 loader->LoadSDigits(opt);
1326 /**************************************************************************/
1328 Int_t AliRunLoader::LoadDigits(Option_t* detectors,Option_t* opt)
1330 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1335 const char* oAll = strstr(detectors,"all");
1342 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1343 loaders = &arr;//get the pointer array
1346 TIter next(loaders);
1348 while((loader = (AliLoader*)next()))
1350 loader->LoadDigits(opt);
1354 /**************************************************************************/
1356 Int_t AliRunLoader::LoadRecPoints(Option_t* detectors,Option_t* opt)
1358 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1363 const char* oAll = strstr(detectors,"all");
1370 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1371 loaders = &arr;//get the pointer array
1374 TIter next(loaders);
1376 while((loader = (AliLoader*)next()))
1378 loader->LoadRecPoints(opt);
1382 /**************************************************************************/
1384 Int_t AliRunLoader::LoadRecParticles(Option_t* detectors,Option_t* opt)
1386 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1391 const char* oAll = strstr(detectors,"all");
1398 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1399 loaders = &arr;//get the pointer array
1402 TIter next(loaders);
1404 while((loader = (AliLoader*)next()))
1406 loader->LoadRecParticles(opt);
1410 /**************************************************************************/
1412 Int_t AliRunLoader::LoadTracks(Option_t* detectors,Option_t* opt)
1414 //LoadHits in selected detectors i.e. detectors="ITS TPC TRD" or "all"
1419 const char* oAll = strstr(detectors,"all");
1426 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1427 loaders = &arr;//get the pointer array
1430 TIter next(loaders);
1432 while((loader = (AliLoader*)next()))
1434 loader->LoadTracks(opt);
1438 /**************************************************************************/
1440 void AliRunLoader::UnloadHits(Option_t* detectors)
1442 //unloads hits for detectors specified in parameter
1446 const char* oAll = strstr(detectors,"all");
1453 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1454 loaders = &arr;//get the pointer to array
1457 TIter next(loaders);
1459 while((loader = (AliLoader*)next()))
1461 loader->UnloadHits();
1464 /**************************************************************************/
1466 void AliRunLoader::UnloadSDigits(Option_t* detectors)
1468 //unloads SDigits for detectors specified in parameter
1472 const char* oAll = strstr(detectors,"all");
1479 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1480 loaders = &arr;//get the pointer to array
1483 TIter next(loaders);
1485 while((loader = (AliLoader*)next()))
1487 loader->UnloadSDigits();
1490 /**************************************************************************/
1492 void AliRunLoader::UnloadDigits(Option_t* detectors)
1494 //unloads Digits for detectors specified in parameter
1498 const char* oAll = strstr(detectors,"all");
1505 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1506 loaders = &arr;//get the pointer to array
1509 TIter next(loaders);
1511 while((loader = (AliLoader*)next()))
1513 loader->UnloadDigits();
1516 /**************************************************************************/
1518 void AliRunLoader::UnloadRecPoints(Option_t* detectors)
1520 //unloads RecPoints for detectors specified in parameter
1524 const char* oAll = strstr(detectors,"all");
1531 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1532 loaders = &arr;//get the pointer to array
1535 TIter next(loaders);
1537 while((loader = (AliLoader*)next()))
1539 loader->UnloadRecPoints();
1542 /**************************************************************************/
1544 void AliRunLoader::UnloadAll(Option_t* detectors)
1546 //calls UnloadAll for detectors names specified in parameter
1547 // option "all" passed can be passed
1551 const char* oAll = strstr(detectors,"all");
1558 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1559 loaders = &arr;//get the pointer to array
1562 TIter next(loaders);
1564 while((loader = (AliLoader*)next()))
1566 loader->UnloadAll();
1569 /**************************************************************************/
1571 void AliRunLoader::UnloadTracks(Option_t* detectors)
1573 //unloads Tracks for detectors specified in parameter
1577 const char* oAll = strstr(detectors,"all");
1584 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1585 loaders = &arr;//get the pointer to array
1588 TIter next(loaders);
1590 while((loader = (AliLoader*)next()))
1592 loader->UnloadTracks();
1595 /**************************************************************************/
1597 void AliRunLoader::UnloadRecParticles(Option_t* detectors)
1599 //unloads Particles for detectors specified in parameter
1603 const char* oAll = strstr(detectors,"all");
1610 GetListOfDetectors(detectors,arr);//this method looks for all Loaders corresponding to names (many) specified in detectors option
1611 loaders = &arr;//get the pointer to array
1614 TIter next(loaders);
1616 while((loader = (AliLoader*)next()))
1618 loader->UnloadRecParticles();
1621 /**************************************************************************/
1623 AliRunLoader* AliRunLoader::GetRunLoader(const char* eventfoldername)
1625 //returns RunLoader from folder named eventfoldername
1626 TFolder* evfold= dynamic_cast<TFolder*>(AliConfig::Instance()->GetTopFolder()->FindObject(eventfoldername));
1631 AliRunLoader* runget = dynamic_cast<AliRunLoader*>(evfold->FindObject(AliRunLoader::fgkRunLoaderName));
1635 /**************************************************************************/
1637 AliLoader* AliRunLoader::GetDetectorLoader(const char* detname, const char* eventfoldername)
1639 //get the loader of the detector with the given name from the global
1641 AliRunLoader* runLoader = GetRunLoader(eventfoldername);
1643 AliErrorClass("No run loader found");
1646 return runLoader->GetDetectorLoader(detname);
1648 /**************************************************************************/
1650 AliLoader* AliRunLoader::GetDetectorLoader(const char* detname)
1652 //get the loader of the detector with the given name from the global
1655 char loadername[256];
1656 sprintf(loadername, "%sLoader", detname);
1657 AliLoader* loader = GetLoader(loadername);
1659 AliError(Form("No loader for %s found", detname));
1664 /**************************************************************************/
1666 TTree* AliRunLoader::GetTreeH(const char* detname, Bool_t maketree, const char* eventfoldername)
1668 //get the tree with hits of the detector with the given name
1669 //if maketree is true and the tree does not exist, the tree is created
1670 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1671 if (!loader) return NULL;
1672 if (!loader->TreeH() && maketree) loader->MakeTree("H");
1673 return loader->TreeH();
1676 /**************************************************************************/
1678 TTree* AliRunLoader::GetTreeH(const char* detname, Bool_t maketree)
1680 //get the tree with hits of the detector with the given name
1681 //if maketree is true and the tree does not exist, the tree is created
1682 AliLoader* loader = GetDetectorLoader(detname);
1683 if (!loader) return NULL;
1684 if (!loader->TreeH() && maketree) loader->MakeTree("H");
1685 return loader->TreeH();
1687 /**************************************************************************/
1689 TTree* AliRunLoader::GetTreeS(const char* detname, Bool_t maketree,const char* eventfoldername)
1691 //get the tree with summable digits of the detector with the given name
1692 //if maketree is true and the tree does not exist, the tree is created
1693 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1694 if (!loader) return NULL;
1695 if (!loader->TreeS() && maketree) loader->MakeTree("S");
1696 return loader->TreeS();
1698 /**************************************************************************/
1700 TTree* AliRunLoader::GetTreeS(const char* detname, Bool_t maketree)
1702 //get the tree with summable digits of the detector with the given name
1703 //if maketree is true and the tree does not exist, the tree is created
1704 AliLoader* loader = GetDetectorLoader(detname);
1705 if (!loader) return NULL;
1706 if (!loader->TreeS() && maketree) loader->MakeTree("S");
1707 return loader->TreeS();
1709 /**************************************************************************/
1711 TTree* AliRunLoader::GetTreeD(const char* detname, Bool_t maketree,const char* eventfoldername)
1713 //get the tree with digits of the detector with the given name
1714 //if maketree is true and the tree does not exist, the tree is created
1715 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1716 if (!loader) return NULL;
1717 if (!loader->TreeD() && maketree) loader->MakeTree("D");
1718 return loader->TreeD();
1720 /**************************************************************************/
1722 TTree* AliRunLoader::GetTreeD(const char* detname, Bool_t maketree)
1724 //get the tree with digits of the detector with the given name
1725 //if maketree is true and the tree does not exist, the tree is created
1726 AliLoader* loader = GetDetectorLoader(detname);
1727 if (!loader) return NULL;
1728 if (!loader->TreeD() && maketree) loader->MakeTree("D");
1729 return loader->TreeD();
1731 /**************************************************************************/
1732 TTree* AliRunLoader::GetTreeR(const char* detname, Bool_t maketree,const char* eventfoldername)
1734 //get the tree with clusters of the detector with the given name
1735 //if maketree is true and the tree does not exist, the tree is created
1736 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1737 if (!loader) return NULL;
1738 if (!loader->TreeR() && maketree) loader->MakeTree("R");
1739 return loader->TreeR();
1741 /**************************************************************************/
1743 TTree* AliRunLoader::GetTreeR(const char* detname, Bool_t maketree)
1745 //get the tree with clusters of the detector with the given name
1746 //if maketree is true and the tree does not exist, the tree is created
1747 AliLoader* loader = GetDetectorLoader(detname);
1748 if (!loader) return NULL;
1749 if (!loader->TreeR() && maketree) loader->MakeTree("R");
1750 return loader->TreeR();
1752 /**************************************************************************/
1754 TTree* AliRunLoader::GetTreeT(const char* detname, Bool_t maketree,const char* eventfoldername)
1756 //get the tree with tracks of the detector with the given name
1757 //if maketree is true and the tree does not exist, the tree is created
1758 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1759 if (!loader) return NULL;
1760 if (!loader->TreeT() && maketree) loader->MakeTree("T");
1761 return loader->TreeT();
1763 /**************************************************************************/
1765 TTree* AliRunLoader::GetTreeT(const char* detname, Bool_t maketree)
1767 //get the tree with tracks of the detector with the given name
1768 //if maketree is true and the tree does not exist, the tree is created
1769 AliLoader* loader = GetDetectorLoader(detname);
1770 if (!loader) return NULL;
1771 if (!loader->TreeT() && maketree) loader->MakeTree("T");
1772 return loader->TreeT();
1774 /**************************************************************************/
1776 TTree* AliRunLoader::GetTreeP(const char* detname, Bool_t maketree,const char* eventfoldername)
1778 //get the tree with particles of the detector with the given name
1779 //if maketree is true and the tree does not exist, the tree is created
1780 AliLoader* loader = GetDetectorLoader(detname,eventfoldername);
1781 if (!loader) return NULL;
1782 if (!loader->TreeP() && maketree) loader->MakeTree("P");
1783 return loader->TreeP();
1785 /**************************************************************************/
1787 TTree* AliRunLoader::GetTreeP(const char* detname, Bool_t maketree)
1789 //get the tree with particles of the detector with the given name
1790 //if maketree is true and the tree does not exist, the tree is created
1791 AliLoader* loader = GetDetectorLoader(detname);
1792 if (!loader) return NULL;
1793 if (!loader->TreeP() && maketree) loader->MakeTree("P");
1794 return loader->TreeP();
1797 /**************************************************************************/
1799 void AliRunLoader::CdGAFile()
1801 //sets gDirectory to galice file
1803 if(fGAFile) fGAFile->cd();
1806 /**************************************************************************/
1808 void AliRunLoader::GetListOfDetectors(const char * namelist,TObjArray& pointerarray) const
1810 //this method looks for all Loaders corresponding
1811 //to names (many) specified in namelist i.e. namelist ("ITS TPC TRD")
1815 strcpy(dets,namelist);//compiler cries when char* = const Option_t*;
1816 dets[strlen(dets)+1] = '\n';//set endl at the end of string
1821 tmp = sscanf(pdet,"%s",buff);//read the string from the input string pdet into buff
1822 if ( (buff[0] == 0) || (tmp == 0) ) break; //if not read
1824 pdet = strstr(pdet,buff) + strlen(buff);//move the input pointer about number of bytes (letters) read
1825 //I am aware that is a little complicated. I don't know the number of spaces between detector names
1826 //so I read the string, than I find where it starts (strstr) and move the pointer about length of a string
1827 //If there is a better way, please write me (Piotr.Skowronski@cern.ch)
1828 //construct the Loader name
1829 TString getname(buff);
1831 AliLoader* loader = GetLoader(getname);//get the Loader
1834 pointerarray.Add(loader);
1838 AliError(Form("Can not find Loader for %s",buff));
1844 /*****************************************************************************/
1846 void AliRunLoader::Clean(const TString& name)
1848 //removes object with given name from event folder and deletes it
1849 if (GetEventFolder() == 0x0) return;
1850 TObject* obj = GetEventFolder()->FindObject(name);
1853 AliDebug(1, Form("name=%s, cleaning %s.",GetName(),name.Data()));
1854 GetEventFolder()->Remove(obj);
1860 /*****************************************************************************/
1862 TTask* AliRunLoader::GetRunDigitizer()
1864 //returns Run Digitizer from folder
1866 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1867 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetDigitizerTaskName());
1868 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1870 /*****************************************************************************/
1872 TTask* AliRunLoader::GetRunSDigitizer()
1874 //returns SDigitizer Task from folder
1876 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1877 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetSDigitizerTaskName());
1878 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1880 /*****************************************************************************/
1882 TTask* AliRunLoader::GetRunReconstructioner()
1884 //returns Reconstructioner Task from folder
1885 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1886 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetReconstructionerTaskName());
1887 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1889 /*****************************************************************************/
1891 TTask* AliRunLoader::GetRunTracker()
1893 //returns Tracker Task from folder
1894 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1895 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetTrackerTaskName());
1896 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1898 /*****************************************************************************/
1900 TTask* AliRunLoader::GetRunPIDTask()
1902 //returns Tracker Task from folder
1903 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1904 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetPIDTaskName());
1905 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1907 /*****************************************************************************/
1909 TTask* AliRunLoader::GetRunQATask()
1911 //returns Quality Assurance Task from folder
1912 TFolder* topf = AliConfig::Instance()->GetTaskFolder();
1915 AliErrorClass("Can not get task folder from AliConfig");
1918 TObject* obj = topf->FindObjectAny(AliConfig::Instance()->GetQATaskName());
1919 return (obj)?dynamic_cast<TTask*>(obj):0x0;
1922 /*****************************************************************************/
1924 void AliRunLoader::SetCompressionLevel(Int_t cl)
1926 //Sets Compression Level in all files
1927 if (fGAFile) fGAFile->SetCompressionLevel(cl);
1928 SetKineComprLevel(cl);
1929 SetTrackRefsComprLevel(cl);
1930 TIter next(fLoaders);
1932 while((loader = (AliLoader*)next()))
1934 loader->SetCompressionLevel(cl);
1937 /**************************************************************************/
1939 void AliRunLoader::SetKineComprLevel(Int_t cl)
1941 //Sets comression level in Kine File
1942 fKineDataLoader->SetCompressionLevel(cl);
1944 /**************************************************************************/
1946 void AliRunLoader::SetTrackRefsComprLevel(Int_t cl)
1948 //Sets comression level in Track Refences File
1949 fTrackRefsDataLoader->SetCompressionLevel(cl);
1951 /**************************************************************************/
1953 void AliRunLoader::UnloadHeader()
1955 //removes TreeE from folder and deletes it
1956 // as well as fHeader object
1961 /**************************************************************************/
1963 void AliRunLoader::UnloadTrigger()
1965 //removes TreeCT from folder and deletes it
1966 // as well as fHeader object
1972 /**************************************************************************/
1974 void AliRunLoader::UnloadKinematics()
1976 //Unloads Kinematics
1977 fKineDataLoader->GetBaseLoader(0)->Unload();
1979 /**************************************************************************/
1981 void AliRunLoader::UnloadTrackRefs()
1983 //Unloads Track Refernces
1984 fTrackRefsDataLoader->GetBaseLoader(0)->Unload();
1986 /**************************************************************************/
1988 void AliRunLoader::UnloadgAlice()
1991 if (gAlice == GetAliRun())
1993 AliDebug(1, "Set gAlice = 0x0");
1994 gAlice = 0x0;//if gAlice is the same that in folder (to be deleted by the way of folder)
1996 AliRun* alirun = GetAliRun();
1997 if (GetEventFolder()) GetEventFolder()->Remove(alirun);
2000 /**************************************************************************/
2002 void AliRunLoader::MakeTrackRefsContainer()
2004 // Makes a tree for Track References
2005 fTrackRefsDataLoader->MakeTree();
2007 /**************************************************************************/
2009 Int_t AliRunLoader::LoadTrackRefs(Option_t* option)
2011 //Load track references from file (opens file and posts tree to folder)
2013 return fTrackRefsDataLoader->GetBaseLoader(0)->Load(option);
2015 /**************************************************************************/
2017 void AliRunLoader::SetDirName(TString& dirname)
2019 //sets directory name
2020 if (dirname.IsNull()) return;
2021 fUnixDirName = dirname;
2022 fKineDataLoader->SetDirName(dirname);
2023 fTrackRefsDataLoader->SetDirName(dirname);
2025 TIter next(fLoaders);
2027 while((loader = (AliLoader*)next()))
2029 loader->SetDirName(dirname);
2033 /*****************************************************************************/
2035 Int_t AliRunLoader::GetFileOffset() const
2037 //returns the file number that is added to the file name for current event
2038 return Int_t(fCurrentEvent/fNEventsPerFile);
2041 /*****************************************************************************/
2042 const TString AliRunLoader::SetFileOffset(const TString& fname)
2044 //adds the the number to the file name at proper place for current event
2045 Long_t offset = (Long_t)GetFileOffset();
2046 if (offset < 1) return fname;
2048 soffset += offset;//automatic conversion to string
2049 TString dotroot(".root");
2050 const TString& offfsetdotroot = offset + dotroot;
2051 TString out = fname;
2052 out = out.ReplaceAll(dotroot,offfsetdotroot);
2053 AliDebug(1, Form(" in=%s out=%s",fname.Data(),out.Data()));
2056 /*****************************************************************************/
2058 void AliRunLoader::SetDigitsFileNameSuffix(const TString& suffix)
2060 //adds the suffix before ".root",
2061 //e.g. TPC.Digits.root -> TPC.DigitsMerged.root
2062 //made on Jiri Chudoba demand
2064 TIter next(fLoaders);
2066 while((loader = (AliLoader*)next()))
2068 loader->SetDigitsFileNameSuffix(suffix);
2071 /*****************************************************************************/
2073 TString AliRunLoader::GetFileName() const
2075 //returns name of galice file
2077 if (fGAFile == 0x0) return result;
2078 result = fGAFile->GetName();
2081 /*****************************************************************************/
2083 void AliRunLoader::SetDetectorAddresses()
2085 //calls SetTreeAddress for all detectors
2086 if (GetAliRun()==0x0) return;
2087 TIter next(GetAliRun()->Modules());
2089 while((mod = (AliModule*)next()))
2091 AliDetector* det = dynamic_cast<AliDetector*>(mod);
2092 if (det) det->SetTreeAddress();
2095 /*****************************************************************************/
2097 void AliRunLoader::Synchronize()
2099 //synchrinizes all writtable files
2100 TIter next(fLoaders);
2102 while((loader = (AliLoader*)next()))
2104 loader->Synchronize();
2107 fKineDataLoader->Synchronize();
2108 fTrackRefsDataLoader->Synchronize();
2110 TFile* file = gROOT->GetFile( gSystem->ConcatFileName( fUnixDirName.Data(), fgkDefaultTriggerFileName.Data() ) ) ;
2111 if( file ) file->Flush();
2113 if (fGAFile) fGAFile->Flush();
2115 /*****************************************************************************/
2116 /*****************************************************************************/