]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Update to include local rho
authorrreed <rosi.jan.reed@cern.ch>
Thu, 13 Feb 2014 05:36:58 +0000 (00:36 -0500)
committermverweij <marta.verweij@cern.ch>
Thu, 13 Feb 2014 08:53:59 +0000 (09:53 +0100)
PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskEmcalJetSpectra.cxx
PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskEmcalJetSpectra.h
PWGJE/EMCALJetTasks/macros/AddTaskEmcalJetSpectra.C

index c851908a777a67e62dab9576935b05ac755f811b..ca18194d5391ddd2e775bade4dfa562fe650176b 100644 (file)
@@ -29,6 +29,8 @@
 #include "AliVCluster.h"
 #include "AliRhoParameter.h"
 #include "AliEmcalParticle.h"
+#include "AliLocalRhoParameter.h"
+#include "AliAnalysisTaskLocalRho.h"
 
 ClassImp(AliAnalysisTaskEmcalJetSpectra)
 
@@ -52,9 +54,12 @@ AliAnalysisTaskEmcalJetSpectra::AliAnalysisTaskEmcalJetSpectra() :
     fHistJetPtvsEP[i]           = 0;
     fHistJetPtvsEPBias[i]       = 0;
     fHistRhovsEP[i]             = 0;
-
+    fHistCorJetPtfromLocalRho[i]= 0;
+    fHistCorJetPtfromGlobalRho[i] = 0;
   }
+  fLocalRhoVal = 0;
   SetMakeGeneralHistograms(kTRUE);
+    
 }
 
 //________________________________________________________________________
@@ -76,7 +81,10 @@ AliAnalysisTaskEmcalJetSpectra::AliAnalysisTaskEmcalJetSpectra(const char *name)
     fHistJetPtvsEP[i]           = 0;
     fHistJetPtvsEPBias[i]       = 0;
     fHistRhovsEP[i]             = 0;
+    fHistCorJetPtfromLocalRho[i]= 0;
+    fHistCorJetPtfromGlobalRho[i] = 0;
    }
+   fLocalRhoVal = 0;
    SetMakeGeneralHistograms(kTRUE);
  }
 
@@ -142,6 +150,16 @@ void AliAnalysisTaskEmcalJetSpectra::UserCreateOutputObjects()
     title = TString(Form("Rho vs EP cent bin %i",i));
     fHistRhovsEP[i] = new TH2F(name,title,500,0,500,400,-2*TMath::Pi(),2*TMath::Pi());
     fOutput->Add(fHistRhovsEP[i]);
+      
+      name = TString(Form("NjetvsCorrJetPtfromLocalRho_%i",i));
+      title = TString(Form("Njets vs Corrected jet pT from Local Rho cent bin %i",i));
+      fHistCorJetPtfromLocalRho[i] = new TH1F(name,title, 500, -250,250);
+      fOutput->Add(fHistCorJetPtfromLocalRho[i]);
+  
+      name = TString(Form("NjetvsCorrJetPtfromGlobalRho_%i",i));
+      title = TString(Form("Njets vs Corrected jet pT from Global Rho cent bin %i",i));
+      fHistCorJetPtfromGlobalRho[i] = new TH1F(name,title, 500, -250,250);
+      fOutput->Add(fHistCorJetPtfromGlobalRho[i]);
   }
 
   
@@ -207,6 +225,13 @@ Bool_t AliAnalysisTaskEmcalJetSpectra::Run()
       continue;
     fHistTrackPt[centbin]->Fill(track->Pt());
   }
+    
+    if(!fLocalRho) {
+        cout<<"name: "<<fLocalRhoName.Data()<<endl;
+        cout<<"found no LocalRho, try to get it from Event based on name"<<endl;
+        fLocalRho = GetLocalRhoFromEvent(fLocalRhoName);
+    }
+
 
   fHistEP0[centbin]->Fill(fEPV0);
   fHistEP0A[centbin]->Fill(fEPV0A);
@@ -239,6 +264,10 @@ Bool_t AliAnalysisTaskEmcalJetSpectra::Run()
      fHistRawJetPtvsTrackPt[centbin]->Fill(jet->Pt(),jet->MaxTrackPt());
      fHistJetPtvsdEP[centbin]->Fill(jetPt,RelativePhi((fEPV0+TMath::Pi()),jet->Phi()));
      fHistJetPtvsEP[centbin]->Fill(jetPt,fEPV0);
+    fLocalRhoVal = fLocalRho->GetLocalVal(jet->Phi(), 0.2);
+      Double_t jetPtLocal = jet->Pt() - jet->Area()*fLocalRhoVal;
+       fHistCorJetPtfromLocalRho[centbin]->Fill(jetPtLocal);
+        fHistCorJetPtfromGlobalRho[centbin]->Fill(jetPt);
      if (jet->MaxTrackPt()>5.0){
        fHistJetPtvsdEPBias[centbin]->Fill(jetPt,RelativePhi((fEPV0+TMath::Pi()),jet->Phi()));
        fHistJetPtvsEPBias[centbin]->Fill(jetPt,fEPV0);
index 8c29edd13e7c3a5a08328f2885404b77aa94f2ab..c7849f4bb292163eb457b4f8564121c647c21f5e 100644 (file)
@@ -7,6 +7,7 @@
 class TH1F;
 class TH2F;
 class THnSparse;
+//class AliLocalRhoParameter;
 
 #include "AliAnalysisTaskEmcalJet.h"
 
@@ -15,6 +16,7 @@ class AliAnalysisTaskEmcalJetSpectra : public AliAnalysisTaskEmcalJet {
   AliAnalysisTaskEmcalJetSpectra();
   AliAnalysisTaskEmcalJetSpectra(const char *name);
   virtual ~AliAnalysisTaskEmcalJetSpectra() {}
   
   
   virtual void           UserCreateOutputObjects();
@@ -23,6 +25,8 @@ class AliAnalysisTaskEmcalJetSpectra : public AliAnalysisTaskEmcalJet {
   Bool_t                 Run();
   virtual Int_t          GetCentBin(Double_t cent) const;
   Float_t                RelativePhi(Double_t mphi,Double_t vphi) const;
+  Float_t                RelativeEPJET(Double_t jetAng, Double_t EPAng) const;
+  Double_t             fLocalRhoVal;
 
  private:
   TH2F                  *fHistRhovsCent; //!
@@ -39,6 +43,8 @@ class AliAnalysisTaskEmcalJetSpectra : public AliAnalysisTaskEmcalJet {
   TH2F                  *fHistJetPtvsEP[6];//!
   TH2F                  *fHistJetPtvsEPBias[6];//!
   TH2F                  *fHistRhovsEP[6]; //!
+  TH1F                  *fHistCorJetPtfromLocalRho[6]; //!
+  TH1F                  *fHistCorJetPtfromGlobalRho[6]; //!
 
 
 
index e389be50bdedfcf06ba588fe6bc97dd9d56c4ec0..2d4381411932a283efb1c8f8ddb7d5a8e2a88652 100644 (file)
@@ -3,8 +3,9 @@
 AliAnalysisTaskEmcalJetSpectra* AddTaskEmcalJetSpectra(
    const char *outfilename    = "AnalysisOutput.root",
    const char *nJets          = "Jets",
-   UInt_t type                = AliAnalysisTaskEmcal::kTPC,
+   UInt_t type                = 0,//AliAnalysisTaskEmcal::kTPC,
    const char *nRhosChEm      = "rhoChEm",
+   const char *lrho           = "lrho",
    const Double_t minPhi      = 1.8,
    const Double_t maxPhi      = 2.74,
    const Double_t minEta      = -0.3,
@@ -39,6 +40,7 @@ AliAnalysisTaskEmcalJetSpectra* AddTaskEmcalJetSpectra(
   spectratask->SetJetsName(nJets);
   spectratask->SetAnaType(type);
   spectratask->SetRhoName(nRhosChEm);
+  spectratask->SetLocalRhoName(lrho);
   spectratask->SetJetPhiLimits(minPhi,maxPhi);
   spectratask->SetJetEtaLimits(minEta,maxEta);
   spectratask->SetJetAreaCut(minArea);