Skip to content

Commit

Permalink
Add a new StartAnalysis method to work directly on a TFileCollection.
Browse files Browse the repository at this point in the history
Use case is for AFs that have locally generated data that are not on the Grid and so not accessible through named datasets using the « Find; … » syntax.
  • Loading branch information
aphecetche authored and hristov committed Nov 27, 2015
1 parent f3a55c1 commit aa5236f
Show file tree
Hide file tree
Showing 2 changed files with 58 additions and 0 deletions.
56 changes: 56 additions & 0 deletions ANALYSIS/ANALYSIS/AliAnalysisManager.cxx
Original file line number Diff line number Diff line change
Expand Up @@ -33,6 +33,7 @@
#include <TMap.h>
#include <TClass.h>
#include <TFile.h>
#include <TFileCollection.h>
#include <TTreeCache.h>
#include <TEnv.h>
#include <TMath.h>
Expand Down Expand Up @@ -2132,6 +2133,61 @@ Long64_t AliAnalysisManager::StartAnalysis(const char *type, const char *dataset
return retv;
}

//______________________________________________________________________________
Long64_t AliAnalysisManager::StartAnalysis(const char *type, TFileCollection* dataset, Long64_t nentries, Long64_t firstentry)
{
// Start analysis for this manager on a given dataset. Analysis task can be:
// LOCAL, PROOF or GRID. Process nentries starting from firstentry.

AliInfo("Using the new direct TFileCollection interface !!!!");

if (!fInitOK) {
Error("StartAnalysis","Analysis manager was not initialized !");
return -1;
}
if (!dataset) {
Error("StartAnalysis","Can not work with a NULL TFileCollection !");
return -1;
}
fIsRemote = kTRUE;
if (fDebug > 1) printf("StartAnalysis %s\n",GetName());
TString anaType = type;
anaType.ToLower();
if (!anaType.Contains("proof")) {
Error("StartAnalysis", "Cannot process datasets in %s mode. Try PROOF.", type);
return -1;
}
fMode = kProofAnalysis;
TString line;
TString proofProcessOpt;
SetEventLoop(kTRUE);
// Set the dataset flag
TObject::SetBit(kUseDataSet);
fTree = 0;

if (!gROOT->GetListOfProofs() || !gROOT->GetListOfProofs()->GetEntries()) {
Error("StartAnalysis", "No PROOF!!! Exiting.");
return -1;
}

// Initialize locally all tasks
RunLocalInit();

line.Form("gProof->AddInput((TObject*)%p);", this);
gROOT->ProcessLine(line);
Long_t retv;
line.Form("gProof->Process((TFileCollection *)%p, \"AliAnalysisSelector\", \"%s\", %lld, %lld);",
dataset, proofProcessOpt.Data(), nentries, firstentry);
char *dispDataset = new char[101];
strncpy(dispDataset, dataset->GetName(), 100);
strncpy(&dispDataset[97], "...", 3);
dispDataset[100] = '\0';
cout << "===== RUNNING PROOF ANALYSIS " << GetName() << " ON DATASET " << dispDataset << endl;
delete dispDataset;
retv = (Long_t)gROOT->ProcessLine(line);
return retv;
}

//______________________________________________________________________________
TFile *AliAnalysisManager::OpenFile(AliAnalysisDataContainer *cont, const char *option, Bool_t ignoreProof)
{
Expand Down
2 changes: 2 additions & 0 deletions ANALYSIS/ANALYSIS/AliAnalysisManager.h
Original file line number Diff line number Diff line change
Expand Up @@ -23,6 +23,7 @@
class TClass;
class TTree;
class TFile;
class TFileCollection;
class TStopwatch;
class TMap;
class AliAnalysisSelector;
Expand Down Expand Up @@ -86,6 +87,7 @@ enum EAliAnalysisFlags {
void RegisterExtraFile(const char *fname);
Long64_t StartAnalysis(const char *type, TTree * const tree, Long64_t nentries=1234567890, Long64_t firstentry=0);
Long64_t StartAnalysis(const char *type, const char *dataset, Long64_t nentries=1234567890, Long64_t firstentry=0);
Long64_t StartAnalysis(const char *type, TFileCollection* dataset, Long64_t nentries=1234567890, Long64_t firstentry=0);
Long64_t StartAnalysis(const char *type, Long64_t nentries=1234567890, Long64_t firstentry=0);
virtual void SlaveBegin(TTree *tree);
virtual void Terminate();
Expand Down

0 comments on commit aa5236f

Please sign in to comment.