]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - ANALYSIS/AliAnalysisManager.cxx
Added separate flag to mark that the LocalInit phase was done + method AliAnalysisMan...
[u/mrichter/AliRoot.git] / ANALYSIS / AliAnalysisManager.cxx
index ed79819c22c2002c257ff7b117805cbb418a9cbd..7fbad0c67d13330f26535982386d57e153ef3832 100644 (file)
@@ -30,6 +30,7 @@
 #include <cerrno>
 #include <Riostream.h>
 #include <TError.h>
+#include <TMap.h>
 #include <TClass.h>
 #include <TFile.h>
 #include <TMath.h>
@@ -92,18 +93,22 @@ AliAnalysisManager::AliAnalysisManager(const char *name, const char *title)
                     fMaxEntries(0),
                     fStatisticsMsg(),
                     fRequestedBranches(),
-                    fStatistics(0)
+                    fStatistics(0),
+                    fGlobals(0)
 {
 // Default constructor.
    fgAnalysisManager = this;
    fgCommonFileName  = "AnalysisResults.root";
-   fTasks      = new TObjArray();
-   fTopTasks   = new TObjArray();
-   fZombies    = new TObjArray();
-   fContainers = new TObjArray();
-   fInputs     = new TObjArray();
-   fOutputs    = new TObjArray();
-   fParamCont  = new TObjArray();
+   if (TClass::IsCallingNew() != TClass::kDummyNew) {
+     fTasks      = new TObjArray();
+     fTopTasks   = new TObjArray();
+     fZombies    = new TObjArray();
+     fContainers = new TObjArray();
+     fInputs     = new TObjArray();
+     fOutputs    = new TObjArray();
+     fParamCont  = new TObjArray();
+     fGlobals    = new TMap();
+   }  
    SetEventLoop(kTRUE);
 }
 
@@ -141,7 +146,8 @@ AliAnalysisManager::AliAnalysisManager(const AliAnalysisManager& other)
                     fMaxEntries(other.fMaxEntries),
                     fStatisticsMsg(other.fStatisticsMsg),
                     fRequestedBranches(other.fRequestedBranches),
-                    fStatistics(other.fStatistics)
+                    fStatistics(other.fStatistics),
+                    fGlobals(other.fGlobals)
 {
 // Copy constructor.
    fTasks      = new TObjArray(*other.fTasks);
@@ -194,6 +200,7 @@ AliAnalysisManager& AliAnalysisManager::operator=(const AliAnalysisManager& othe
       fStatisticsMsg = other.fStatisticsMsg;
       fRequestedBranches = other.fRequestedBranches;
       fStatistics = other.fStatistics;
+      fGlobals = new TMap();
    }
    return *this;
 }
@@ -215,6 +222,7 @@ AliAnalysisManager::~AliAnalysisManager()
    if (fMCtruthEventHandler) delete fMCtruthEventHandler;
    if (fEventPool) delete fEventPool;
    if (fgAnalysisManager==this) fgAnalysisManager = NULL;
+   if (fGlobals) {fGlobals->DeleteAll(); delete fGlobals;}
 }
 
 //______________________________________________________________________________
@@ -243,7 +251,7 @@ Int_t AliAnalysisManager::GetRunFromAlienPath(const char *path)
    ind1 = s.Index("/00");
    if (ind1>0) {
       ind2 = s.Index("/",ind1+1);
-      if (ind2>0) srun = s(ind1+1, ind2-ind1-1);
+      if (ind2-ind1>8) srun = s(ind1+1, ind2-ind1-1);
    }   
    if (srun.IsNull()) {
       ind1 = s.Index("/LHC");
@@ -957,6 +965,7 @@ void AliAnalysisManager::Terminate()
    }
    gROOT->cd();
    next1.Reset();
+   TString copiedFiles;
    while ((output=(AliAnalysisDataContainer*)next1())) {
       // Close all files at output
       TDirectory *opwd = gDirectory;
@@ -964,14 +973,20 @@ void AliAnalysisManager::Terminate()
          // Clear file list to release object ownership to user.
 //         output->GetFile()->Clear();
          output->GetFile()->Close();
-         output->SetFile(NULL);
          // Copy merged outputs in alien if requested
-         if (fSpecialOutputLocation.Length() && 
-             fSpecialOutputLocation.BeginsWith("alien://")) {
+         if (fSpecialOutputLocation.BeginsWith("alien://")) {
+            if (copiedFiles.Contains(output->GetFile()->GetName())) {
+               if (opwd) opwd->cd();
+               output->SetFile(NULL);
+               continue;
+            } 
             Info("Terminate", "Copy file %s to %s", output->GetFile()->GetName(),fSpecialOutputLocation.Data()); 
+            gROOT->ProcessLine("if (!gGrid) TGrid::Connect(\"alien:\");");
             TFile::Cp(output->GetFile()->GetName(), 
                       Form("%s/%s", fSpecialOutputLocation.Data(), output->GetFile()->GetName()));
+            copiedFiles += output->GetFile()->GetName();
          }             
+         output->SetFile(NULL);
       }   
       if (opwd) opwd->cd();
    }   
@@ -1497,6 +1512,22 @@ void AliAnalysisManager::ResetAnalysis()
    CleanContainers();
 }
 
+//______________________________________________________________________________
+void AliAnalysisManager::RunLocalInit()
+{
+// Run LocalInit method for all tasks.
+   TDirectory *cdir = gDirectory;
+   if (IsTrainInitialized()) return;
+   TIter nextTask(fTasks);
+   AliAnalysisTask *task;
+   while ((task=(AliAnalysisTask*)nextTask())) {
+      gROOT->cd();
+      task->LocalInit();
+   }
+   cdir->cd();
+   TObject::SetBit(kTasksInitialized, kTRUE);
+}   
+
 //______________________________________________________________________________
 Long64_t AliAnalysisManager::StartAnalysis(const char *type, Long64_t nentries, Long64_t firstentry)
 {
@@ -1534,11 +1565,7 @@ Long64_t AliAnalysisManager::StartAnalysis(const char *type, TTree * const tree,
    TString anaType = type;
    anaType.ToLower();
    fMode = kLocalAnalysis;
-   Bool_t runlocalinit = kTRUE;
-   if (anaType.Contains("file")) {
-      runlocalinit = kFALSE;
-      fIsRemote = kTRUE;
-   }   
+   if (anaType.Contains("file"))      fIsRemote = kTRUE;
    if (anaType.Contains("proof"))     fMode = kProofAnalysis;
    else if (anaType.Contains("grid")) fMode = kGridAnalysis;
    else if (anaType.Contains("mix"))  fMode = kMixingAnalysis;
@@ -1555,12 +1582,7 @@ Long64_t AliAnalysisManager::StartAnalysis(const char *type, TTree * const tree,
          // Write analysis manager in the analysis file
          cout << "===== RUNNING GRID ANALYSIS: " << GetName() << endl;
          // run local task configuration
-         TIter nextTask(fTasks);
-         AliAnalysisTask *task;
-         while ((task=(AliAnalysisTask*)nextTask())) {
-            task->LocalInit();
-            gROOT->cd();
-         }
+         RunLocalInit();
          if (!fGridHandler->StartAnalysis(nentries, firstentry)) {
             Info("StartAnalysis", "Grid analysis was stopped and cannot be terminated");
             cdir->cd();
@@ -1606,13 +1628,7 @@ Long64_t AliAnalysisManager::StartAnalysis(const char *type, TTree * const tree,
    // Initialize locally all tasks (happens for all modes)
    TIter next(fTasks);
    AliAnalysisTask *task;
-   if (runlocalinit) {
-      while ((task=(AliAnalysisTask*)next())) {
-         task->LocalInit();
-         gROOT->cd();
-      }
-      if (getsysInfo) AliSysInfo::AddStamp("LocalInit_all", 0);
-   }   
+   RunLocalInit();
    
    switch (fMode) {
       case kLocalAnalysis:
@@ -1794,12 +1810,8 @@ Long64_t AliAnalysisManager::StartAnalysis(const char *type, const char *dataset
    }   
 
    // Initialize locally all tasks
-   TIter next(fTasks);
-   AliAnalysisTask *task;
-   while ((task=(AliAnalysisTask*)next())) {
-      task->LocalInit();
-   }
-   
+   RunLocalInit();
+      
    line = Form("gProof->AddInput((TObject*)%p);", this);
    gROOT->ProcessLine(line);
    Long_t retv;
@@ -2439,3 +2451,101 @@ const char* AliAnalysisManager::GetOADBPath()
       
    return oadbPath;
 }
+
+//______________________________________________________________________________
+void AliAnalysisManager::SetGlobalStr(const char *key, const char *value)
+{
+// Define a custom string variable mapped to a global unique name. The variable
+// can be then retrieved by a given analysis macro via GetGlobalStr(key).
+   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
+   if (!mgr) {
+      ::Error("AliAnalysisManager::SetGlobalStr", "No analysis manager defined");
+      return;
+   }   
+   Bool_t valid = kFALSE;
+   TString existing = AliAnalysisManager::GetGlobalStr(key, valid);
+   if (valid) {
+      ::Error("AliAnalysisManager::SetGlobalStr", "Global %s = %s already defined.", key, existing.Data());
+      return;
+   }
+   mgr->GetGlobals()->Add(new TObjString(key), new TObjString(value));
+}
+
+//______________________________________________________________________________
+const char *AliAnalysisManager::GetGlobalStr(const char *key, Bool_t &valid)
+{
+// Static method to retrieve a global variable defined via SetGlobalStr.
+   valid = kFALSE;
+   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
+   if (!mgr) return 0;
+   TObject *value = mgr->GetGlobals()->GetValue(key);
+   if (!value) return 0;
+   valid = kTRUE;
+   return value->GetName();
+}
+
+//______________________________________________________________________________
+void AliAnalysisManager::SetGlobalInt(const char *key, Int_t value)
+{
+// Define a custom integer variable mapped to a global unique name. The variable
+// can be then retrieved by a given analysis macro via GetGlobalInt(key).
+   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
+   if (!mgr) {
+      ::Error("AliAnalysisManager::SetGlobalStr", "No analysis manager defined");
+      return;
+   }   
+   Bool_t valid = kFALSE;
+   Int_t existing = AliAnalysisManager::GetGlobalInt(key, valid);
+   if (valid) {
+      ::Error("AliAnalysisManager::SetGlobalInt", "Global %s = %i already defined.", key, existing);
+      return;
+   }
+   mgr->GetGlobals()->Add(new TObjString(key), new TObjString(TString::Format("%i",value)));
+}
+
+//______________________________________________________________________________
+Int_t AliAnalysisManager::GetGlobalInt(const char *key, Bool_t &valid)
+{
+// Static method to retrieve a global variable defined via SetGlobalInt.
+   valid = kFALSE;
+   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
+   if (!mgr) return 0;
+   TObject *value = mgr->GetGlobals()->GetValue(key);
+   if (!value) return 0;
+   valid = kTRUE;
+   TString s = value->GetName();
+   return s.Atoi();
+}
+
+//______________________________________________________________________________
+void AliAnalysisManager::SetGlobalDbl(const char *key, Double_t value)
+{
+// Define a custom double precision variable mapped to a global unique name. The variable
+// can be then retrieved by a given analysis macro via GetGlobalInt(key).
+   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
+   if (!mgr) {
+      ::Error("AliAnalysisManager::SetGlobalStr", "No analysis manager defined");
+      return;
+   }   
+   Bool_t valid = kFALSE;
+   Double_t existing = AliAnalysisManager::GetGlobalDbl(key, valid);
+   if (valid) {
+      ::Error("AliAnalysisManager::SetGlobalInt", "Global %s = %g already defined.", key, existing);
+      return;
+   }
+   mgr->GetGlobals()->Add(new TObjString(key), new TObjString(TString::Format("%f.16",value)));
+}
+
+//______________________________________________________________________________
+Double_t AliAnalysisManager::GetGlobalDbl(const char *key, Bool_t &valid)
+{
+// Static method to retrieve a global variable defined via SetGlobalDbl.
+   valid = kFALSE;
+   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
+   if (!mgr) return 0;
+   TObject *value = mgr->GetGlobals()->GetValue(key);
+   if (!value) return 0;
+   valid = kTRUE;
+   TString s = value->GetName();
+   return s.Atof();
+}