]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Add EMCal only (no jets) example task
authormverweij <marta.verweij@cern.ch>
Tue, 24 Jun 2014 12:55:39 +0000 (14:55 +0200)
committermverweij <marta.verweij@cern.ch>
Tue, 24 Jun 2014 13:04:46 +0000 (15:04 +0200)
activate clus-track matching checker

PWG/CMakelibPWGEMCAL.pkg
PWG/EMCAL/AliAnalysisTaskEmcalSample.cxx [new file with mode: 0644]
PWG/EMCAL/AliAnalysisTaskEmcalSample.h [new file with mode: 0644]
PWG/EMCAL/macros/AddTaskEmcalSample.C [new file with mode: 0644]
PWG/EMCAL/macros/runEMCalAnalysis.C [new file with mode: 0644]
PWG/PWGEMCALLinkDef.h

index 2131d799aeab6b269aa8e2486583aae4219f99ba..5cadd245258de49d8532debd877bedc93e0df5fb 100644 (file)
@@ -29,6 +29,7 @@
 set ( SRCS 
  EMCAL/AliAnalysisTaskEMCALClusterizeFast.cxx
  EMCAL/AliAnalysisTaskEmcal.cxx
+ EMCAL/AliAnalysisTaskEmcalSample.cxx
  EMCAL/AliClusterContainer.cxx
  EMCAL/AliEMCALClusterParams.cxx
  EMCAL/AliEmcalAodTrackFilterTask.cxx
diff --git a/PWG/EMCAL/AliAnalysisTaskEmcalSample.cxx b/PWG/EMCAL/AliAnalysisTaskEmcalSample.cxx
new file mode 100644 (file)
index 0000000..e28188d
--- /dev/null
@@ -0,0 +1,238 @@
+// $Id$
+//
+// Emcal sample analysis task.
+//
+// Author: S.Aiola, M. Verweij
+
+#include <TClonesArray.h>
+#include <TH1F.h>
+#include <TH2F.h>
+#include <TH3F.h>
+#include <TList.h>
+#include <TLorentzVector.h>
+
+#include "AliVCluster.h"
+#include "AliAODCaloCluster.h"
+#include "AliESDCaloCluster.h"
+#include "AliVTrack.h"
+#include "AliLog.h"
+#include "AliParticleContainer.h"
+#include "AliClusterContainer.h"
+#include "AliPicoTrack.h"
+
+#include "AliAnalysisTaskEmcalSample.h"
+
+ClassImp(AliAnalysisTaskEmcalSample)
+
+//________________________________________________________________________
+AliAnalysisTaskEmcalSample::AliAnalysisTaskEmcalSample() : 
+  AliAnalysisTaskEmcal("AliAnalysisTaskEmcalSample", kTRUE),
+  fHistTracksPt(0),
+  fHistClustersPt(0),
+  fHistPtDEtaDPhiTrackClus(0),
+  fHistPtDEtaDPhiClusTrack(0),
+  fTracksCont(0),
+  fCaloClustersCont(0)
+{
+  // Default constructor.
+
+  fHistTracksPt       = new TH1*[fNcentBins];
+  fHistClustersPt     = new TH1*[fNcentBins];
+
+  for (Int_t i = 0; i < fNcentBins; i++) {
+    fHistTracksPt[i] = 0;
+    fHistClustersPt[i] = 0;
+  }
+
+  SetMakeGeneralHistograms(kTRUE);
+}
+
+//________________________________________________________________________
+AliAnalysisTaskEmcalSample::AliAnalysisTaskEmcalSample(const char *name) : 
+  AliAnalysisTaskEmcal(name, kTRUE),
+  fHistTracksPt(0),
+  fHistClustersPt(0),
+  fHistPtDEtaDPhiTrackClus(0),
+  fHistPtDEtaDPhiClusTrack(0),
+  fTracksCont(0),
+  fCaloClustersCont(0)
+{
+  // Standard constructor.
+
+  fHistTracksPt       = new TH1*[fNcentBins];
+  fHistClustersPt     = new TH1*[fNcentBins];
+
+  for (Int_t i = 0; i < fNcentBins; i++) {
+    fHistTracksPt[i] = 0;
+    fHistClustersPt[i] = 0;
+  }
+
+  SetMakeGeneralHistograms(kTRUE);
+}
+
+//________________________________________________________________________
+AliAnalysisTaskEmcalSample::~AliAnalysisTaskEmcalSample()
+{
+  // Destructor.
+}
+
+//________________________________________________________________________
+void AliAnalysisTaskEmcalSample::UserCreateOutputObjects()
+{
+  // Create user output.
+
+  AliAnalysisTaskEmcal::UserCreateOutputObjects();
+
+  fTracksCont       = GetParticleContainer(0);
+  fCaloClustersCont = GetClusterContainer(0);
+  fTracksCont->SetClassName("AliVTrack");
+  fCaloClustersCont->SetClassName("AliAODCaloCluster");
+
+  TString histname;
+
+  for (Int_t i = 0; i < fNcentBins; i++) {
+    if (fParticleCollArray.GetEntriesFast()>0) {
+      histname = "fHistTracksPt_";
+      histname += i;
+      fHistTracksPt[i] = new TH1F(histname.Data(), histname.Data(), fNbins / 2, fMinBinPt, fMaxBinPt / 2);
+      fHistTracksPt[i]->GetXaxis()->SetTitle("p_{T,track} (GeV/c)");
+      fHistTracksPt[i]->GetYaxis()->SetTitle("counts");
+      fOutput->Add(fHistTracksPt[i]);
+    }
+
+    if (fClusterCollArray.GetEntriesFast()>0) {
+      histname = "fHistClustersPt_";
+      histname += i;
+      fHistClustersPt[i] = new TH1F(histname.Data(), histname.Data(), fNbins / 2, fMinBinPt, fMaxBinPt / 2);
+      fHistClustersPt[i]->GetXaxis()->SetTitle("p_{T,clus} (GeV/c)");
+      fHistClustersPt[i]->GetYaxis()->SetTitle("counts");
+      fOutput->Add(fHistClustersPt[i]);
+    }
+  }
+
+  histname = "fHistPtDEtaDPhiTrackClus";
+  fHistPtDEtaDPhiTrackClus = new TH3F(histname.Data(),Form("%s;#it{p}_{T}^{track};#Delta#eta;#Delta#varphi",histname.Data()),100,0.,100.,100,-0.1,0.1,100,-0.1,0.1);
+  fOutput->Add(fHistPtDEtaDPhiTrackClus);
+
+  histname = "fHistPtDEtaDPhiClusTrack";
+  fHistPtDEtaDPhiClusTrack = new TH3F(histname.Data(),Form("%s;#it{p}_{T}^{clus};#Delta#eta;#Delta#varphi",histname.Data()),100,0.,100.,100,-0.1,0.1,100,-0.1,0.1);
+  fOutput->Add(fHistPtDEtaDPhiClusTrack);
+
+  PostData(1, fOutput); // Post data for ALL output slots > 0 here.
+}
+
+//________________________________________________________________________
+Bool_t AliAnalysisTaskEmcalSample::FillHistograms()
+{
+  // Fill histograms.
+
+  if (fTracksCont) {
+    AliVTrack *track = static_cast<AliVTrack*>(fTracksCont->GetNextAcceptParticle(0)); 
+    while(track) {
+      fHistTracksPt[fCentBin]->Fill(track->Pt()); 
+      track = static_cast<AliVTrack*>(fTracksCont->GetNextAcceptParticle());
+    }
+  }
+  
+  if (fCaloClustersCont) {
+    AliVCluster *cluster = fCaloClustersCont->GetNextAcceptCluster(0); 
+    while(cluster) {
+      TLorentzVector nPart;
+      cluster->GetMomentum(nPart, fVertex);
+      fHistClustersPt[fCentBin]->Fill(nPart.Pt());
+      cluster = fCaloClustersCont->GetNextAcceptCluster();
+    }
+  }
+
+  CheckClusTrackMatching();
+
+  return kTRUE;
+}
+
+//________________________________________________________________________
+void AliAnalysisTaskEmcalSample::CheckClusTrackMatching()
+{
+  
+  if(!fTracksCont || !fCaloClustersCont)
+    return;
+
+  Double_t deta = 999;
+  Double_t dphi = 999;
+
+  //Get closest cluster to track
+  AliVTrack *track = static_cast<AliVTrack*>(fTracksCont->GetNextAcceptParticle(0)); 
+  while(track) {
+    //Get matched cluster
+    Int_t emc1 = track->GetEMCALcluster();
+    if(fCaloClustersCont && emc1>=0) {
+      AliVCluster *clusMatch = fCaloClustersCont->GetCluster(emc1);
+      if(clusMatch) {
+       AliPicoTrack::GetEtaPhiDiff(track, clusMatch, dphi, deta);
+       fHistPtDEtaDPhiTrackClus->Fill(track->Pt(),deta,dphi);
+      }
+    }
+    track = static_cast<AliVTrack*>(fTracksCont->GetNextAcceptParticle());
+  }
+  
+  //Get closest track to cluster
+  AliVCluster *cluster = fCaloClustersCont->GetNextAcceptCluster(0); 
+  while(cluster) {
+    TLorentzVector nPart;
+    cluster->GetMomentum(nPart, fVertex);
+    fHistClustersPt[fCentBin]->Fill(nPart.Pt());
+    
+    //Get matched track
+    AliVTrack *mt = NULL;      
+    AliAODCaloCluster *acl = dynamic_cast<AliAODCaloCluster*>(cluster);
+    if(acl) {
+      if(acl->GetNTracksMatched()>1)
+       mt = static_cast<AliVTrack*>(acl->GetTrackMatched(0));
+    }
+    else {
+      AliESDCaloCluster *ecl = dynamic_cast<AliESDCaloCluster*>(cluster);
+      Int_t im = ecl->GetTrackMatchedIndex();
+      if(fTracksCont && im>=0) {
+       mt = static_cast<AliVTrack*>(fTracksCont->GetParticle(im));
+      }
+    }
+    if(mt) {
+      AliPicoTrack::GetEtaPhiDiff(mt, cluster, dphi, deta);
+      fHistPtDEtaDPhiClusTrack->Fill(nPart.Pt(),deta,dphi);
+      
+      /* //debugging
+        if(mt->IsEMCAL()) {
+        Int_t emc1 = mt->GetEMCALcluster();
+        Printf("current id: %d  emc1: %d",fCaloClustersCont->GetCurrentID(),emc1);
+        AliVCluster *clm = fCaloClustersCont->GetCluster(emc1);
+        AliPicoTrack::GetEtaPhiDiff(mt, clm, dphi, deta);
+        Printf("deta: %f dphi: %f",deta,dphi);
+        }
+      */
+    }
+    cluster = fCaloClustersCont->GetNextAcceptCluster();
+  }
+}
+
+//________________________________________________________________________
+void AliAnalysisTaskEmcalSample::ExecOnce() {
+
+  AliAnalysisTaskEmcal::ExecOnce();
+
+  if (fTracksCont && fTracksCont->GetArray() == 0) fTracksCont = 0;
+  if (fCaloClustersCont && fCaloClustersCont->GetArray() == 0) fCaloClustersCont = 0;
+
+}
+
+//________________________________________________________________________
+Bool_t AliAnalysisTaskEmcalSample::Run()
+{
+  // Run analysis code here, if needed. It will be executed before FillHistograms().
+
+  return kTRUE;  // If return kFALSE FillHistogram() will NOT be executed.
+}
+
+//________________________________________________________________________
+void AliAnalysisTaskEmcalSample::Terminate(Option_t *) 
+{
+  // Called once at the end of the analysis.
+}
diff --git a/PWG/EMCAL/AliAnalysisTaskEmcalSample.h b/PWG/EMCAL/AliAnalysisTaskEmcalSample.h
new file mode 100644 (file)
index 0000000..b3bf22a
--- /dev/null
@@ -0,0 +1,45 @@
+#ifndef ALIANALYSISTASKEMCALSAMPLE_H
+#define ALIANALYSISTASKEMCALSAMPLE_H
+
+// $Id$
+
+class TH1;
+class TH2;
+class TH3;
+class AliParticleContainer;
+class AliClusterContainer;
+
+#include "AliAnalysisTaskEmcal.h"
+
+class AliAnalysisTaskEmcalSample : public AliAnalysisTaskEmcal {
+ public:
+
+  AliAnalysisTaskEmcalSample();
+  AliAnalysisTaskEmcalSample(const char *name);
+  virtual ~AliAnalysisTaskEmcalSample();
+
+  void                        UserCreateOutputObjects();
+  void                        Terminate(Option_t *option);
+
+ protected:
+  void                        ExecOnce();
+  Bool_t                      FillHistograms()   ;
+  Bool_t                      Run()              ;
+  void                        CheckClusTrackMatching();
+
+  // General histograms
+  TH1                       **fHistTracksPt;            //!Track pt spectrum
+  TH1                       **fHistClustersPt;          //!Cluster pt spectrum
+  TH3                        *fHistPtDEtaDPhiTrackClus; //!track pt, delta eta, delta phi to matched cluster
+  TH3                        *fHistPtDEtaDPhiClusTrack; //!cluster pt, delta eta, delta phi to matched track
+
+  AliParticleContainer       *fTracksCont;                 //!Tracks
+  AliClusterContainer        *fCaloClustersCont;           //!Clusters  
+
+ private:
+  AliAnalysisTaskEmcalSample(const AliAnalysisTaskEmcalSample&);            // not implemented
+  AliAnalysisTaskEmcalSample &operator=(const AliAnalysisTaskEmcalSample&); // not implemented
+
+  ClassDef(AliAnalysisTaskEmcalSample, 1) // emcal sample analysis task
+};
+#endif
diff --git a/PWG/EMCAL/macros/AddTaskEmcalSample.C b/PWG/EMCAL/macros/AddTaskEmcalSample.C
new file mode 100644 (file)
index 0000000..2bde65a
--- /dev/null
@@ -0,0 +1,59 @@
+// $Id$
+
+AliAnalysisTaskEmcalSample* AddTaskEmcalSample(
+  const char *ntracks            = "Tracks",
+  const char *nclusters          = "CaloClusters",
+  Int_t       nCentBins          = 1,
+  const char *taskname           = "AliAnalysisTaskEmcalSample"
+)
+{  
+  // Get the pointer to the existing analysis manager via the static access method.
+  //==============================================================================
+  AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
+  if (!mgr)
+  {
+    ::Error("AddTaskEmcalSample", "No analysis manager to connect to.");
+    return NULL;
+  }  
+  
+  // Check the analysis type using the event handlers connected to the analysis manager.
+  //==============================================================================
+  if (!mgr->GetInputEventHandler())
+  {
+    ::Error("AddTaskEmcalSample", "This task requires an input event handler");
+    return NULL;
+  }
+  
+  //-------------------------------------------------------
+  // Init the task and do settings
+  //-------------------------------------------------------
+
+  TString name(taskname);
+  Printf("name: %s",name.Data());
+
+  AliAnalysisTaskEmcalSample* emcTask = new AliAnalysisTaskEmcalSample(name);
+  emcTask->SetCentRange(0.,100.);
+  emcTask->SetNCentBins(nCentBins);
+
+  AliParticleContainer *trackCont  = emcTask->AddParticleContainer(ntracks);
+  trackCont->SetClassName("AliVTrack");
+  AliClusterContainer *clusterCont = emcTask->AddClusterContainer(nclusters);
+  
+  //-------------------------------------------------------
+  // Final settings, pass to manager and set the containers
+  //-------------------------------------------------------
+  
+  mgr->AddTask(emcTask);
+  
+  // Create containers for input/output
+  AliAnalysisDataContainer *cinput1  = mgr->GetCommonInputContainer()  ;
+  TString contname(name);
+  contname += "_histos";
+  AliAnalysisDataContainer *coutput1 = mgr->CreateContainer(contname.Data(), 
+                                                           TList::Class(),AliAnalysisManager::kOutputContainer,
+                                                           Form("%s", AliAnalysisManager::GetCommonFileName()));
+  mgr->ConnectInput  (emcTask, 0,  cinput1 );
+  mgr->ConnectOutput (emcTask, 1, coutput1 );
+  
+  return emcTask;
+}
diff --git a/PWG/EMCAL/macros/runEMCalAnalysis.C b/PWG/EMCAL/macros/runEMCalAnalysis.C
new file mode 100644 (file)
index 0000000..3f3e8c8
--- /dev/null
@@ -0,0 +1,350 @@
+// runEMCalAnalysis.C
+// =====================
+// This macro can be used to run EMCal analysis within the EMCal part of the Jet Framework.
+//
+// Examples:
+// -> Analyze ESDs from the pA pilot run on the AliEn grid with your task in AnaClass.cxx/.h
+//     dataType = "ESD", useGrid = kTRUE, pattern = "*ESDs/pass2/*ESDs.root", addCXXs = "AnaClass.cxx", 
+//     addHs = "AnaClass.h", gridDir = "/alice/data/2012/LHC12g", gridMode = "full", runNumbers = "188359 188362"
+//     
+// -> Analyze AODs (up to 96 files) locally given in files_aod.txt
+//     dataType = "AOD", useGrid = kFALSE, numLocalFiles = 96
+//
+// MERGING ON ALIEN
+// ++++++++++++++++
+// If you run on the grid, you can monitor the jobs with alimonitor.cern.ch. When enough of them are in DONE state,
+// you have to merge the output. This can be done automatically, if you just change the gridMode to "terminate" and
+// give the EXACT name of the task whose output should be merged in uniqueName.
+// 
+//
+// Authors: R. Haake, S. Aiola, M. Verweij
+
+#include <ctime>
+#include "TGrid.h"
+
+AliAnalysisGrid* CreateAlienHandler(const char* uniqueName, const char* gridDir, const char* gridMode, const char* runNumbers, 
+                                     const char* pattern, TString additionalCode, TString additionalHeaders, Int_t maxFilesPerWorker, 
+                                     Int_t workerTTL, Bool_t isMC);
+                                    
+//______________________________________________________________________________
+void runEMCalAnalysis(
+         Bool_t         useGrid             = kTRUE,                      // local or grid
+         const char*    gridMode            = "test",                      // set the grid run mode (can be "full", "test", "offline", "submit" or "terminate")
+         const char*    dataType            = "AOD",                       // set the analysis type, AOD, ESD or sESD
+         const char*    pattern             = "*ESDs/pass2/AOD145/*AOD.root",    // file pattern (here one can specify subdirs like passX etc.) (used on grid)
+         const char*    gridDir             = "/alice/data/2011/LHC11h_2",   // dir on alien, where the files live (used on grid)
+         const char*    runNumbers          = "167903 167915",             // considered run numbers (used on grid)
+         UInt_t         numLocalFiles       = 1,                          // number of files analyzed locally  
+         const char*    runPeriod           = "LHC11h",                    // set the run period (used on grid)
+         const char*    uniqueName          = "EMCal_LEGOTrainTest",     // sets base string for the name of the task on the grid
+         UInt_t         pSel                = AliVEvent::kAny,             // used event selection for every task except for the analysis tasks
+         Bool_t         useTender           = kTRUE,                       // trigger, if tender, track and cluster selection should be used (always)
+         Bool_t         isMC                = kFALSE,                      // trigger, if MC handler should be used
+         // Here you have to specify additional code files you want to use but that are not in aliroot
+         const char*    addCXXs             = "",
+         const char*    addHs               = "",
+
+         // These two settings depend on the dataset and your quotas on the AliEN services
+         Int_t          maxFilesPerWorker   = 4,
+         Int_t          workerTTL           = 7200
+
+         )
+{
+
+  // Some pre-settings and constants
+  enum AlgoType {kKT, kANTIKT};
+  enum JetType  {kFULLJETS, kCHARGEDJETS, kNEUTRALJETS};
+  gSystem->SetFPEMask();
+  gSystem->Setenv("ETRAIN_ROOT", ".");
+  gSystem->Setenv("ETRAIN_PERIOD", runPeriod);
+  // change this objects to strings
+  TString usedData(dataType);
+  TString additionalCXXs(addCXXs);
+  TString additionalHs(addHs);
+  cout << dataType << " analysis chosen" << endl;
+  if (useGrid)  
+  {
+    cout << "-- using AliEn grid.\n";
+    if (usedData == "sESD") 
+    {
+      cout << "Skimmed ESD analysis not available on the grid!" << endl;
+      return;
+    }
+  }
+  else
+    cout << "-- using local analysis.\n";
+  
+
+  // Load necessary libraries
+  LoadLibs();
+
+  // Create analysis manager
+  AliAnalysisManager* mgr = new AliAnalysisManager(uniqueName);
+
+  // Check type of input and create handler for it
+  TString localFiles("-1");
+  if(usedData == "AOD")
+  {
+    localFiles = "files_aod.txt";
+    gROOT->LoadMacro("$ALICE_ROOT/ANALYSIS/macros/train/AddAODHandler.C");
+    AliAODInputHandler* aodH = AddAODHandler();
+  }
+  else if((usedData == "ESD") || (usedData == "sESD"))
+  {
+    if (usedData == "ESD")
+      localFiles = "files_esd.txt";
+    else
+      localFiles = "files_sesd.txt";
+    
+    gROOT->LoadMacro("$ALICE_ROOT/ANALYSIS/macros/train/AddESDHandler.C");
+    AliESDInputHandler* esdH = AddESDHandler();
+  }
+  else
+  {
+    cout << "Data type not recognized! You have to specify ESD, AOD, or sESD!\n";
+  }
+
+  if(!useGrid)
+    cout << "Using " << localFiles.Data() << " as input file list.\n";
+
+  // Create MC handler, if MC is demanded
+  if (isMC && (usedData != "AOD"))
+  {
+    AliMCEventHandler* mcH = new AliMCEventHandler();
+    mcH->SetPreReadMode(AliMCEventHandler::kLmPreRead);
+    mcH->SetReadTR(kTRUE);
+    mgr->SetMCtruthEventHandler(mcH); 
+  }
+  
+  // ################# Now: Add some basic tasks
+
+  // Physics selection task
+  gROOT->LoadMacro("$ALICE_ROOT/PWG/EMCAL/macros/AddTaskEmcalPhysicsSelection.C");
+  AliPhysicsSelectionTask *physSelTask = AddTaskEmcalPhysicsSelection(kTRUE, kTRUE, pSel, 5, 5, 10, kTRUE, -1, -1, -1, -1);
+
+  if (!physSelTask) 
+  {
+    cout << "no physSelTask"; 
+    return; 
+  }
+
+  // Centrality task
+  if (usedData == "ESD") {
+    gROOT->LoadMacro("$ALICE_ROOT/ANALYSIS/macros/AddTaskCentrality.C");
+    AliCentralitySelectionTask *centralityTask = AddTaskCentrality(kTRUE);
+  }
+
+  // Compatibility task, only needed for skimmed ESD
+  if (usedData == "sESD") {
+    gROOT->LoadMacro("$ALICE_ROOT/PWG/EMCAL/macros/AddTaskEmcalCompat.C");
+    AliEmcalCompatTask *comptask = AddTaskEmcalCompat();
+  }
+
+  // Setup task
+  gROOT->LoadMacro("$ALICE_ROOT/PWG/EMCAL/macros/AddTaskEmcalSetup.C");
+  AliEmcalSetupTask *setupTask = AddTaskEmcalSetup();
+  setupTask->SetGeoPath("$ALICE_ROOT/OADB/EMCAL");
+  
+  // Tender Supplies
+  if (useTender) {
+    gROOT->LoadMacro("$ALICE_ROOT/PWG/EMCAL/macros/AddTaskEmcalPreparation.C");
+    //adjust pass when running locally. On grid give empty string, will be picked up automatically from path to ESD/AOD file
+    AliAnalysisTaskSE *clusm = AddTaskEmcalPreparation(runPeriod,"pass1"); 
+  }
+
+  // Names of the different objects passed around; these are the default names; added here mostly for documentation purposes
+  // rhoName is only set if the background subtraction is switched on (doBkg)
+  TString tracksName = "AODFilterTracks";
+  TString clustersName = "EmcCaloClusters";
+
+  // ################# Now: Call jet preparation macro (picotracks, hadronic corrected caloclusters, ...) 
+  gROOT->LoadMacro("$ALICE_ROOT/PWG/EMCAL/macros/AddTaskMatchingChain.C");
+  AddTaskMatchingChain(runPeriod,pSel,
+                      clustersName,1.,kTRUE,
+                      0.1,kTRUE,kTRUE);
+  
+  // ################# Now: Add analysis task
+  // Here you can put in your AddTaskMacro for your task
+  gROOT->LoadMacro("$ALICE_ROOT/PWG/EMCAL/macros/AddTaskEmcalSample.C");
+  AliAnalysisTaskEmcalSample* anaTask = AddTaskEmcalSample(tracksName, clustersName, 4);
+
+  // Set the physics selection for all given tasks
+  TObjArray *toptasks = mgr->GetTasks();
+  for (Int_t i=0; i<toptasks->GetEntries(); ++i) 
+  {
+    AliAnalysisTaskSE *task = dynamic_cast<AliAnalysisTaskSE*>(toptasks->At(i));
+    if (!task)
+      continue;
+    if (task->InheritsFrom("AliPhysicsSelectionTask"))
+      continue;
+    ::Info("setPSel", "Set physics selection for %s (%s)", task->GetName(), task->ClassName());
+    task->SelectCollisionCandidates(pSel);
+  }
+
+  mgr->SetUseProgressBar(1, 25);
+        
+  if (!mgr->InitAnalysis()) 
+    return;
+  mgr->PrintStatus();
+
+  if (useGrid) 
+  {  // GRID CALCULATION
+
+    AliAnalysisGrid *plugin = CreateAlienHandler(uniqueName, gridDir, gridMode, runNumbers, pattern, additionalCXXs, additionalHs, maxFilesPerWorker, workerTTL, isMC);
+    mgr->SetGridHandler(plugin);
+
+    // start analysis
+    cout << "Starting GRID Analysis...";
+    mgr->SetDebugLevel(0);
+    mgr->StartAnalysis("grid",10);
+  }
+  else
+  {  // LOCAL CALCULATION
+
+    TChain* chain = 0;
+    if (usedData == "AOD") 
+    {
+      gROOT->LoadMacro("$ALICE_ROOT/PWG/EMCAL/macros/CreateAODChain.C");
+      chain = CreateAODChain(localFiles.Data(), numLocalFiles);
+    }
+    else
+    {  // ESD or skimmed ESD
+      gROOT->LoadMacro("$ALICE_ROOT/PWG/EMCAL/macros/CreateESDChain.C");
+      chain = CreateESDChain(localFiles.Data(), numLocalFiles);
+    }
+    
+    // start analysis
+    cout << "Starting LOCAL Analysis...";
+    mgr->SetDebugLevel(2);
+    mgr->StartAnalysis("local", chain);
+  }
+}
+
+//______________________________________________________________________________
+void LoadLibs()
+{
+  // Load common libraries (better too many than too few)
+  gSystem->Load("libTree");
+  gSystem->Load("libVMC");
+  gSystem->Load("libGeom");
+  gSystem->Load("libGui");
+  gSystem->Load("libXMLParser");
+  gSystem->Load("libMinuit");
+  gSystem->Load("libMinuit2");
+  gSystem->Load("libProof");
+  gSystem->Load("libPhysics");
+  gSystem->Load("libSTEERBase");
+  gSystem->Load("libESD");
+  gSystem->Load("libAOD");
+  gSystem->Load("libOADB");
+  gSystem->Load("libANALYSIS");
+  gSystem->Load("libCDB");
+  gSystem->Load("libRAWDatabase");
+  gSystem->Load("libSTEER");
+  gSystem->Load("libEVGEN");
+  gSystem->Load("libANALYSISalice");
+  gSystem->Load("libCORRFW");
+  gSystem->Load("libTOFbase");
+  //gSystem->Load("libTOFrec");
+  gSystem->Load("libRAWDatabase.so");
+  gSystem->Load("libRAWDatarec.so");
+  gSystem->Load("libTPCbase.so");
+  gSystem->Load("libTPCrec.so");
+  gSystem->Load("libITSbase.so");
+  gSystem->Load("libITSrec.so");
+  gSystem->Load("libTRDbase.so");
+  gSystem->Load("libTENDER.so");
+  gSystem->Load("libSTAT.so");
+  gSystem->Load("libTRDrec.so");
+  gSystem->Load("libHMPIDbase.so");
+  gSystem->Load("libPWGPP.so");
+  gSystem->Load("libPWGHFbase");
+  gSystem->Load("libPWGDQdielectron");
+  gSystem->Load("libPWGHFhfe");
+  gSystem->Load("libEMCALUtils");
+  gSystem->Load("libPHOSUtils");
+  gSystem->Load("libPWGCaloTrackCorrBase");
+  gSystem->Load("libEMCALraw");
+  gSystem->Load("libEMCALbase");
+  gSystem->Load("libEMCALrec");
+  gSystem->Load("libTRDbase");
+  gSystem->Load("libVZERObase");
+  gSystem->Load("libVZEROrec");
+  gSystem->Load("libTENDER");
+  gSystem->Load("libTENDERSupplies");
+  gSystem->Load("libPWGTools");
+  gSystem->Load("libPWGEMCAL");
+  gSystem->Load("libESDfilter");
+  gSystem->Load("libPWGGAEMCALTasks");
+  gSystem->Load("libPWGCFCorrelationsBase");
+  gSystem->Load("libPWGCFCorrelationsDPhi");
+
+  // include paths
+  gSystem->AddIncludePath("-Wno-deprecated");
+  gSystem->AddIncludePath("-I$ALICE_ROOT -I$ALICE_ROOT/include -I$ALICE_ROOT/EMCAL");
+  gSystem->AddIncludePath("-I$ALICE_ROOT/PWGDQ/dielectron -I$ALICE_ROOT/PWGHF/hfe");
+}
+
+AliAnalysisGrid* CreateAlienHandler(const char* uniqueName, const char* gridDir, const char* gridMode, const char* runNumbers, 
+                                     const char* pattern, TString additionalCode, TString additionalHeaders, Int_t maxFilesPerWorker, 
+                                     Int_t workerTTL, Bool_t isMC)
+{
+  TDatime currentTime;
+  TString tmpName(uniqueName);
+
+  // Only add current date and time when not in terminate mode! In this case the exact name has to be supplied by the user
+  if(strcmp(gridMode, "terminate"))
+  {
+    tmpName += "_";
+    tmpName += currentTime.GetDate();
+    tmpName += "_";
+    tmpName += currentTime.GetTime();
+  }
+
+  TString tmpAdditionalLibs("");
+  tmpAdditionalLibs = Form("libTree.so libVMC.so libGeom.so libGui.so libXMLParser.so libMinuit.so libMinuit2.so libProof.so libPhysics.so libSTEERBase.so libESD.so libAOD.so libOADB.so libANALYSIS.so libCDB.so libRAWDatabase.so libSTEER.so libANALYSISalice.so libCORRFW.so libTOFbase.so libRAWDatabase.so libRAWDatarec.so libTPCbase.so libTPCrec.so libITSbase.so libITSrec.so libTRDbase.so libTENDER.so libSTAT.so libTRDrec.so libHMPIDbase.so libPWGPP.so libPWGHFbase.so libPWGDQdielectron.so libPWGHFhfe.so libEMCALUtils.so libPHOSUtils.so libPWGCaloTrackCorrBase.so libEMCALraw.so libEMCALbase.so libEMCALrec.so libTRDbase.so libVZERObase.so libVZEROrec.so libTENDER.so libTENDERSupplies.so libESDfilter.so libPWGTools.so libPWGEMCAL.so libPWGGAEMCALTasks.so libPWGCFCorrelationsBase.so libPWGCFCorrelationsDPhi.so  %s %s",additionalCode.Data(),additionalHeaders.Data());
+
+
+  TString macroName("");
+  TString execName("");
+  TString jdlName("");
+  macroName = Form("%s.C", tmpName.Data());
+  execName = Form("%s.sh", tmpName.Data());
+  jdlName = Form("%s.jdl", tmpName.Data());
+
+  AliAnalysisAlien *plugin = new AliAnalysisAlien();
+  plugin->SetOverwriteMode();
+  plugin->SetRunMode(gridMode);
+     
+  // Here you can set the (Ali)ROOT version you want to use
+  plugin->SetAPIVersion("V1.1x");
+  plugin->SetROOTVersion("v5-34-08-6");
+  plugin->SetAliROOTVersion("vAN-20140525");
+
+  plugin->SetGridDataDir(gridDir); // e.g. "/alice/sim/LHC10a6"
+  plugin->SetDataPattern(pattern); //dir structure in run directory
+  if (!isMC)
+   plugin->SetRunPrefix("000");
+
+  plugin->AddRunList(runNumbers);
+
+  plugin->SetGridWorkingDir(Form("work/%s",tmpName.Data()));
+  plugin->SetGridOutputDir("output"); // In this case will be $HOME/work/output
+
+  plugin->SetAnalysisSource(additionalCode.Data());
+  plugin->SetAdditionalLibs(tmpAdditionalLibs.Data());
+
+  plugin->SetDefaultOutputs(kTRUE);
+  //plugin->SetMergeExcludes("");
+  plugin->SetAnalysisMacro(macroName.Data());
+  plugin->SetSplitMaxInputFileNumber(maxFilesPerWorker);
+  plugin->SetExecutable(execName.Data());
+  plugin->SetTTL(workerTTL);
+  plugin->SetInputFormat("xml-single");
+  plugin->SetJDLName(jdlName.Data());
+  plugin->SetPrice(1);      
+  plugin->SetSplitMode("se");
+  plugin->SetNtestFiles(1);
+
+  return plugin;
+}
index 965837f27e38b9142a44c0fa165af345895f9bfb..fe58f6d58befc6517cd416c8ece127e68e8bbcbe 100644 (file)
@@ -6,6 +6,7 @@
 
 #pragma link C++ class AliAnalysisTaskEMCALClusterizeFast+;
 #pragma link C++ class AliAnalysisTaskEmcal+;
+#pragma link C++ class AliAnalysisTaskEmcalSample+;
 #pragma link C++ class AliClusterContainer+;
 #pragma link C++ class AliEMCALClusterParams+;
 #pragma link C++ class AliEmcalAodTrackFilterTask+;