]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
updates for JetSample example task
authormverweij <marta.verweij@cern.ch>
Sun, 16 Mar 2014 11:29:15 +0000 (12:29 +0100)
committerhristov <Peter.Hristov@cern.ch>
Thu, 27 Mar 2014 15:25:06 +0000 (16:25 +0100)
PWGJE/EMCALJetTasks/AliAnalysisTaskEmcalJetSample.cxx
PWGJE/EMCALJetTasks/AliAnalysisTaskEmcalJetSample.h
PWGJE/EMCALJetTasks/macros/AddTaskEmcalJetSample.C

index 2364691694595cf801f4f673cc9958227c1f4376..075bf5e8ff4f20344aa7c426e3cfae823287b856 100644 (file)
@@ -11,6 +11,7 @@
 #include <TLorentzVector.h>
 
 #include "AliVCluster.h"
+#include "AliAODCaloCluster.h"
 #include "AliVTrack.h"
 #include "AliEmcalJet.h"
 #include "AliRhoParameter.h"
@@ -26,6 +27,13 @@ ClassImp(AliAnalysisTaskEmcalJetSample)
 //________________________________________________________________________
 AliAnalysisTaskEmcalJetSample::AliAnalysisTaskEmcalJetSample() : 
   AliAnalysisTaskEmcalJet("AliAnalysisTaskEmcalJetSample", kTRUE),
+  fHistTracksPt(0),
+  fHistClustersPt(0),
+  fHistLeadingJetPt(0),
+  fHistJetsPhiEta(0),
+  fHistJetsPtArea(0),
+  fHistJetsPtLeadHad(0),
+  fHistJetsCorrPtArea(0),
   fJetsCont(0),
   fTracksCont(0),
   fCaloClustersCont(0)
@@ -33,7 +41,15 @@ AliAnalysisTaskEmcalJetSample::AliAnalysisTaskEmcalJetSample() :
 {
   // Default constructor.
 
-  for (Int_t i = 0; i < 4; i++) {
+  fHistTracksPt       = new TH1*[fNcentBins];
+  fHistClustersPt     = new TH1*[fNcentBins];
+  fHistLeadingJetPt   = new TH1*[fNcentBins];
+  fHistJetsPhiEta     = new TH2*[fNcentBins];
+  fHistJetsPtArea     = new TH2*[fNcentBins];
+  fHistJetsPtLeadHad  = new TH2*[fNcentBins];
+  fHistJetsCorrPtArea = new TH2*[fNcentBins];
+
+  for (Int_t i = 0; i < fNcentBins; i++) {
     fHistTracksPt[i] = 0;
     fHistClustersPt[i] = 0;
     fHistLeadingJetPt[i] = 0;
@@ -49,14 +65,28 @@ AliAnalysisTaskEmcalJetSample::AliAnalysisTaskEmcalJetSample() :
 //________________________________________________________________________
 AliAnalysisTaskEmcalJetSample::AliAnalysisTaskEmcalJetSample(const char *name) : 
   AliAnalysisTaskEmcalJet(name, kTRUE),
+  fHistTracksPt(0),
+  fHistClustersPt(0),
+  fHistLeadingJetPt(0),
+  fHistJetsPhiEta(0),
+  fHistJetsPtArea(0),
+  fHistJetsPtLeadHad(0),
+  fHistJetsCorrPtArea(0),
   fJetsCont(0),
   fTracksCont(0),
   fCaloClustersCont(0)
 {
   // Standard constructor.
 
+  fHistTracksPt       = new TH1*[fNcentBins];
+  fHistClustersPt     = new TH1*[fNcentBins];
+  fHistLeadingJetPt   = new TH1*[fNcentBins];
+  fHistJetsPhiEta     = new TH2*[fNcentBins];
+  fHistJetsPtArea     = new TH2*[fNcentBins];
+  fHistJetsPtLeadHad  = new TH2*[fNcentBins];
+  fHistJetsCorrPtArea = new TH2*[fNcentBins];
 
-  for (Int_t i = 0; i < 4; i++) {
+  for (Int_t i = 0; i < fNcentBins; i++) {
     fHistTracksPt[i] = 0;
     fHistClustersPt[i] = 0;
     fHistLeadingJetPt[i] = 0;
@@ -82,13 +112,20 @@ void AliAnalysisTaskEmcalJetSample::UserCreateOutputObjects()
 
   AliAnalysisTaskEmcalJet::UserCreateOutputObjects();
 
-  fJetsCont         = GetJetContainer(0);
-  fTracksCont       = fJetsCont->GetParticleContainer();
-  fCaloClustersCont = fJetsCont->GetClusterContainer();
+  fJetsCont           = GetJetContainer(0);
+  if(fJetsCont) { //get particles and clusters connected to jets
+    fTracksCont       = fJetsCont->GetParticleContainer();
+    fCaloClustersCont = fJetsCont->GetClusterContainer();
+  } else {        //no jets, just analysis tracks and clusters
+    fTracksCont       = GetParticleContainer(0);
+    fCaloClustersCont = GetClusterContainer(0);
+  }
+  fTracksCont->SetClassName("AliVTrack");
+  fCaloClustersCont->SetClassName("AliAODCaloCluster");
 
   TString histname;
 
-  for (Int_t i = 0; i < 4; i++) {
+  for (Int_t i = 0; i < fNcentBins; i++) {
     if (fParticleCollArray.GetEntriesFast()>0) {
       histname = "fHistTracksPt_";
       histname += i;
@@ -156,23 +193,30 @@ Bool_t AliAnalysisTaskEmcalJetSample::FillHistograms()
   // Fill histograms.
 
   if (fTracksCont) {
-    AliVParticle *track = fTracksCont->GetNextAcceptParticle(0); 
+    AliVTrack *track = static_cast<AliVTrack*>(fTracksCont->GetNextAcceptParticle(0)); 
     while(track) {
       fHistTracksPt[fCentBin]->Fill(track->Pt()); 
-      
-      track = fTracksCont->GetNextAcceptParticle(); 
+      Int_t emc1 = track->GetEMCALcluster();
+      Printf("EMCAL cluster %d",emc1);
+      track = static_cast<AliVTrack*>(fTracksCont->GetNextAcceptParticle());
     }
   }
   
   if (fCaloClustersCont) {
-    AliVCluster *cluster = fCaloClustersCont->GetNextAcceptCluster(0); 
+    //    AliVCluster *cluster = fCaloClustersCont->GetNextAcceptCluster(0); 
+    AliAODCaloCluster *cluster = static_cast<AliAODCaloCluster*>(fCaloClustersCont->GetNextAcceptCluster(0));
     while(cluster) {
 
       TLorentzVector nPart;
       cluster->GetMomentum(nPart, fVertex);
       fHistClustersPt[fCentBin]->Fill(nPart.Pt());
       
-      cluster = fCaloClustersCont->GetNextAcceptCluster(); 
+      AliVTrack *mt = NULL;
+      Printf("N matched tracks: %d",cluster->GetNTracksMatched());
+      if(cluster->GetNTracksMatched()>1)
+       mt = static_cast<AliVTrack*>(cluster->GetTrackMatched(0));
+      if(mt) Printf("matched track pt: %f eta: %f phi: %f",mt->Pt(),mt->Eta(),mt->Phi());
+      cluster = static_cast<AliAODCaloCluster*>(fCaloClustersCont->GetNextAcceptCluster());
     }
   }
 
index cf7279771665163d5e29f8f856e818282f8b124f..1edd125f30b0e207773680c598da179089e79758 100644 (file)
@@ -27,13 +27,13 @@ class AliAnalysisTaskEmcalJetSample : public AliAnalysisTaskEmcalJet {
   Bool_t                      Run()              ;
 
   // General histograms
-  TH1                        *fHistTracksPt[4];            //!Track pt spectrum
-  TH1                        *fHistClustersPt[4];          //!Cluster pt spectrum
-  TH1                        *fHistLeadingJetPt[4];        //!Leading jet pt spectrum
-  TH2                        *fHistJetsPhiEta[4];          //!Phi-Eta distribution of jets
-  TH2                        *fHistJetsPtArea[4];          //!Jet pt vs. area
-  TH2                        *fHistJetsPtLeadHad[4];       //!Jet pt vs. leading hadron
-  TH2                        *fHistJetsCorrPtArea[4];      //!Jet pt - bkg vs. area
+  TH1                       **fHistTracksPt;            //!Track pt spectrum
+  TH1                       **fHistClustersPt;          //!Cluster pt spectrum
+  TH1                       **fHistLeadingJetPt;        //!Leading jet pt spectrum
+  TH2                       **fHistJetsPhiEta;          //!Phi-Eta distribution of jets
+  TH2                       **fHistJetsPtArea;          //!Jet pt vs. area
+  TH2                       **fHistJetsPtLeadHad;       //!Jet pt vs. leading hadron
+  TH2                       **fHistJetsCorrPtArea;      //!Jet pt - bkg vs. area
 
   AliJetContainer            *fJetsCont;                   //!Jets
   AliParticleContainer       *fTracksCont;                 //!Tracks
@@ -43,6 +43,6 @@ class AliAnalysisTaskEmcalJetSample : public AliAnalysisTaskEmcalJet {
   AliAnalysisTaskEmcalJetSample(const AliAnalysisTaskEmcalJetSample&);            // not implemented
   AliAnalysisTaskEmcalJetSample &operator=(const AliAnalysisTaskEmcalJetSample&); // not implemented
 
-  ClassDef(AliAnalysisTaskEmcalJetSample, 2) // jet sample analysis task
+  ClassDef(AliAnalysisTaskEmcalJetSample, 3) // jet sample analysis task
 };
 #endif
index fc5df3ca84e809ed56f737cafb6c37e56f7d3594..f9838555a5cfab3cc3a6d904609fe5b09aa1699a 100644 (file)
@@ -5,6 +5,7 @@ AliAnalysisTaskEmcalJetSample* AddTaskEmcalJetSample(
   const char *nclusters          = "CaloClusters",
   const char *njets              = "Jets",
   const char *nrho               = "Rho",
+  Int_t       nCentBins          = 1,
   Double_t    jetradius          = 0.2,
   Double_t    jetptcut           = 1,
   Double_t    jetareacut         = 0.6,
@@ -51,8 +52,11 @@ AliAnalysisTaskEmcalJetSample* AddTaskEmcalJetSample(
   Printf("name: %s",name.Data());
 
   AliAnalysisTaskEmcalJetSample* jetTask = new AliAnalysisTaskEmcalJetSample(name);
+  jetTask->SetCentRange(0.,100.);
+  jetTask->SetNCentBins(nCentBins);
 
   AliParticleContainer *trackCont  = jetTask->AddParticleContainer(ntracks);
+  trackCont->SetClassName("AliVTrack");
   AliClusterContainer *clusterCont = jetTask->AddClusterContainer(nclusters);
 
   TString strType(type);