updates JetSample to check cluster-track matching
authormverweij <marta.verweij@cern.ch>
Mon, 17 Mar 2014 08:22:36 +0000 (09:22 +0100)
committerhristov <Peter.Hristov@cern.ch>
Thu, 27 Mar 2014 15:25:07 +0000 (16:25 +0100)
PWGJE/EMCALJetTasks/AliAnalysisTaskEmcalJetSample.cxx
PWGJE/EMCALJetTasks/AliAnalysisTaskEmcalJetSample.h

index 075bf5e..abeb7f2 100644 (file)
@@ -7,11 +7,13 @@
 #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 "AliEmcalJet.h"
 #include "AliRhoParameter.h"
@@ -19,6 +21,7 @@
 #include "AliJetContainer.h"
 #include "AliParticleContainer.h"
 #include "AliClusterContainer.h"
+#include "AliPicoTrack.h"
 
 #include "AliAnalysisTaskEmcalJetSample.h"
 
@@ -34,6 +37,8 @@ AliAnalysisTaskEmcalJetSample::AliAnalysisTaskEmcalJetSample() :
   fHistJetsPtArea(0),
   fHistJetsPtLeadHad(0),
   fHistJetsCorrPtArea(0),
+  fHistPtDEtaDPhiTrackClus(0),
+  fHistPtDEtaDPhiClusTrack(0),
   fJetsCont(0),
   fTracksCont(0),
   fCaloClustersCont(0)
@@ -72,6 +77,8 @@ AliAnalysisTaskEmcalJetSample::AliAnalysisTaskEmcalJetSample(const char *name) :
   fHistJetsPtArea(0),
   fHistJetsPtLeadHad(0),
   fHistJetsCorrPtArea(0),
+  fHistPtDEtaDPhiTrackClus(0),
+  fHistPtDEtaDPhiClusTrack(0),
   fJetsCont(0),
   fTracksCont(0),
   fCaloClustersCont(0)
@@ -184,6 +191,15 @@ void AliAnalysisTaskEmcalJetSample::UserCreateOutputObjects()
       }
     }
   }
+
+  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.
 }
 
@@ -196,27 +212,18 @@ Bool_t AliAnalysisTaskEmcalJetSample::FillHistograms()
     AliVTrack *track = static_cast<AliVTrack*>(fTracksCont->GetNextAcceptParticle(0)); 
     while(track) {
       fHistTracksPt[fCentBin]->Fill(track->Pt()); 
-      Int_t emc1 = track->GetEMCALcluster();
-      Printf("EMCAL cluster %d",emc1);
       track = static_cast<AliVTrack*>(fTracksCont->GetNextAcceptParticle());
     }
   }
   
   if (fCaloClustersCont) {
-    //    AliVCluster *cluster = fCaloClustersCont->GetNextAcceptCluster(0); 
-    AliAODCaloCluster *cluster = static_cast<AliAODCaloCluster*>(fCaloClustersCont->GetNextAcceptCluster(0));
+    AliVCluster *cluster = fCaloClustersCont->GetNextAcceptCluster(0); 
     while(cluster) {
-
       TLorentzVector nPart;
       cluster->GetMomentum(nPart, fVertex);
       fHistClustersPt[fCentBin]->Fill(nPart.Pt());
-      
-      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());
+
+      cluster = fCaloClustersCont->GetNextAcceptCluster();
     }
   }
 
@@ -240,10 +247,76 @@ Bool_t AliAnalysisTaskEmcalJetSample::FillHistograms()
     if(jet) fHistLeadingJetPt[fCentBin]->Fill(jet->Pt());
   }
 
+  CheckClusTrackMatching();
+
   return kTRUE;
 }
 
 //________________________________________________________________________
+void AliAnalysisTaskEmcalJetSample::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 AliAnalysisTaskEmcalJetSample::ExecOnce() {
 
   AliAnalysisTaskEmcalJet::ExecOnce();
index 1edd125..bd759de 100644 (file)
@@ -5,6 +5,7 @@
 
 class TH1;
 class TH2;
+class TH3;
 class AliJetContainer;
 class AliParticleContainer;
 class AliClusterContainer;
@@ -25,6 +26,7 @@ class AliAnalysisTaskEmcalJetSample : public AliAnalysisTaskEmcalJet {
   void                        ExecOnce();
   Bool_t                      FillHistograms()   ;
   Bool_t                      Run()              ;
+  void                        CheckClusTrackMatching();
 
   // General histograms
   TH1                       **fHistTracksPt;            //!Track pt spectrum
@@ -34,6 +36,8 @@ class AliAnalysisTaskEmcalJetSample : public AliAnalysisTaskEmcalJet {
   TH2                       **fHistJetsPtArea;          //!Jet pt vs. area
   TH2                       **fHistJetsPtLeadHad;       //!Jet pt vs. leading hadron
   TH2                       **fHistJetsCorrPtArea;      //!Jet pt - bkg vs. area
+  TH3                        *fHistPtDEtaDPhiTrackClus; //!track pt, delta eta, delta phi to matched cluster
+  TH3                        *fHistPtDEtaDPhiClusTrack; //!cluster pt, delta eta, delta phi to matched track
 
   AliJetContainer            *fJetsCont;                   //!Jets
   AliParticleContainer       *fTracksCont;                 //!Tracks