]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
add event plane decorrelation histogram (in eta)
authorrbertens <rbertens@cern.ch>
Wed, 1 Oct 2014 12:37:50 +0000 (14:37 +0200)
committerrbertens <rbertens@cern.ch>
Wed, 1 Oct 2014 12:38:02 +0000 (14:38 +0200)
PWGCF/FLOW/macros/AddTaskJetFlowMC.C
PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskJetV2.cxx
PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskJetV2.h

index 62a78ff3e911be5ec101995d044635609895ad17..d86bfd0b4a2e28dfdeb40cf296f83feec5632a15 100644 (file)
@@ -26,7 +26,10 @@ AliAnalysisTaskJetFlowMC* AddTaskJetFlowMC(
   const char *name              ="AliAnalysisTaskJetFlowMC",
   Bool_t doQA                   = kFALSE,
   Bool_t doDecay                = kFALSE,       // be sure to load pythia libs
-  Bool_t doEmbedding            = kFALSE        // not to be used on train
+  Bool_t doEmbedding            = kFALSE,       // not to be used on train
+  Double_t ptHardPythiaMin      = 0.,           // pt hard bin lower bound
+  Double_t ptHardPythiaMax      = 10.           // pt hard bin upper bound
+
   )
 {  
   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
@@ -49,19 +52,27 @@ AliAnalysisTaskJetFlowMC* AddTaskJetFlowMC(
   // to decay tracks using as a default pythia
   if(doDecay) task->SetDecayer(TaskJetFlowMC::GetDecayer(kTRUE));
   // to embed pythia jets to the events
-  if(doEmbedding) AliJetEmbeddingFromGenTask* eTask = TaskJetFlowMC::EmbedGeneratedJets(TaskJetFlowMC::GetPythiaGenerator(), outputTracks);
+  if(doEmbedding) AliJetEmbeddingFromGenTask* eTask = TaskJetFlowMC::EmbedGeneratedJets(TaskJetFlowMC::GetPythiaGenerator(2760., ptHardPythiaMin, ptHardPythiaMax), outputTracks);
   return task;
 }
 
 namespace TaskJetFlowMC {
 
     TF1* GetSpectrum() {
-        // thermal spectrum used for ALICE SIMULATION PLOTS ALI-SIMUL-75145 ALI-SIMUL-75171
+        // full spectrum used for ALICE SIMULATION PLOTS ALI-SIMUL-75145 ALI-SIMUL-75171
+        // combination of boltzmann spectrum and hard jet spectrum
         TF1* fspectrum = new TF1("fspectrum", "[0]*(TMath::Power([1], 2)*x*TMath::Exp(-[1]*x))+(x>1)*[2]*(1.17676e-01*TMath::Sqrt(0.1396*0.1396+x*x)*TMath::Power(1.+1./[3]/8.21795e-01*TMath::Sqrt(0.1396*0.1396+x*x),-1.*[3]))*(1/(1 + TMath::Exp(-(x - [4])/[5])))", .2, 200.);
-        fspectrum->SetParameters(2434401791.20528 ,2.98507 ,10069622.25117 ,5.50000 ,2.80000 ,0.20000 );
+        fspectrum->SetParameters(2434401791.20528, 2.98507, 10069622.25117, 5.50000, 2.80000, 0.20000);
         return fspectrum;   
     }
 
+    TF1* GetThermalSpectrum() {
+        // pure boltzmann part of thermal spectrum
+        TF1* boltzmann = new TF1("boltzmann", "[0]*(TMath::Power([1], 2)*x*TMath::Exp(-[1]*x))");
+        boltzmann->SetParameters(2434401791.20528, 2.98507);
+        return boltzmann;
+    }
+
     TVirtualDecayer* GetDecayer(Bool_t local = kTRUE) {
         // setup a decayer
         if(local) gSystem->Load("$ALICE_ROOT/lib/tgt_linuxx8664gcc/libpythia6");
index 9ef2b3480c26fe7b61e324a39caba2f17d2d812a..6c7f1c4da30ead0601acb43f15a9a4fb9cdbfa71 100644 (file)
@@ -38,6 +38,7 @@
 #include <TF2.h>
 #include <TH1F.h>
 #include <TH2F.h>
+#include <TH3F.h>
 #include <TProfile.h>
 #include <TFile.h>
 // aliroot includes
@@ -82,6 +83,7 @@ AliAnalysisTaskJetV2::AliAnalysisTaskJetV2() : AliAnalysisTaskEmcalJet("AliAnaly
         fHistPsiVZEROALeadingJet[i] = 0;  
         fHistPsiVZEROCLeadingJet[i] = 0;
         fHistPsiVZEROCombLeadingJet[i] = 0;
+        fHistPsi2Correlation[i] = 0;
         fHistRhoPackage[i] = 0;
         fHistRho[i] = 0;
         fHistRCPhiEta[i] = 0;
@@ -140,6 +142,7 @@ AliAnalysisTaskJetV2::AliAnalysisTaskJetV2(const char* name, runModeType type) :
         fHistPsiVZEROALeadingJet[i] = 0;  
         fHistPsiVZEROCLeadingJet[i] = 0;
         fHistPsiVZEROCombLeadingJet[i] = 0;
+        fHistPsi2Correlation[i] = 0;
         fHistRhoPackage[i] = 0;
         fHistRho[i] = 0;
         fHistRCPhiEta[i] = 0;
@@ -380,7 +383,7 @@ TH1F* AliAnalysisTaskJetV2::BookTH1F(const char* name, const char* x, Int_t bins
     return histogram;   
 }
 //_____________________________________________________________________________
-TH2F* AliAnalysisTaskJetV2::BookTH2F(const char* name, const char* x, const char*y, Int_t binsx, Double_t minx, Double_t maxx, Int_t binsy, Double_t miny, Double_t maxy, Int_t c, Bool_t append)
+TH2F* AliAnalysisTaskJetV2::BookTH2F(const char* name, const char* x, const char* y, Int_t binsx, Double_t minx, Double_t maxx, Int_t binsy, Double_t miny, Double_t maxy, Int_t c, Bool_t append)
 {
     // book a TH2F and connect it to the output container
     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
@@ -399,6 +402,28 @@ TH2F* AliAnalysisTaskJetV2::BookTH2F(const char* name, const char* x, const char
     return histogram;   
 }
 //_____________________________________________________________________________
+TH3F* AliAnalysisTaskJetV2::BookTH3F(const char* name, const char* x, const char* y, const char* z, Int_t binsx, Double_t minx, Double_t maxx, Int_t binsy, Double_t miny, Double_t maxy, Int_t binsz, Double_t minz, Double_t maxz, Int_t c, Bool_t append)
+{
+    // book a TH2F and connect it to the output container
+    if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
+    if(fReduceBinsXByFactor > 0 ) {
+        binsx = TMath::Nint(binsx/fReduceBinsXByFactor);
+        binsy = TMath::Nint(binsy/fReduceBinsXByFactor);
+        binsz = TMath::Nint(binsz/fReduceBinsXByFactor);
+    }
+    if(!fOutputList) return 0x0;
+    TString title(name);
+    if(c!=-1) { // format centrality dependent histograms accordingly
+        name = Form("%s_%i", name, c);
+        title += Form("_%i-%i", (int)fCentralityClasses->At(c), (int)(fCentralityClasses->At((1+c))));
+    }
+    title += Form(";%s;%s;%s", x, y, z);
+    TH3F* histogram = new TH3F(name, title.Data(), binsx, minx, maxx, binsy, miny, maxy, binsz, minz, maxz);
+    histogram->Sumw2();
+    if(append) fOutputList->Add(histogram);
+    return histogram;   
+}
+//_____________________________________________________________________________
 void AliAnalysisTaskJetV2::UserCreateOutputObjects()
 {
     // create output objects. also initializes some default values in case they aren't 
@@ -446,6 +471,7 @@ void AliAnalysisTaskJetV2::UserCreateOutputObjects()
             fHistPsiVZEROALeadingJet[i] =   BookTH2F("fHistPsiVZEROALeadingJet", "p_{t} [GeV/c]", "#Psi_{VZEROA}", 350, -100, 250, 50, -1.*TMath::Pi()/2., TMath::Pi()/2., i);
             fHistPsiVZEROCLeadingJet[i] =   BookTH2F("fHistPsiVZEROCLeadingJet", "p_{t} [GeV/c]", "#Psi_{VZEROC}", 350, -100, 250, 50, -1.*TMath::Pi()/2., TMath::Pi()/2., i);
             fHistPsiVZEROCombLeadingJet[i] = BookTH2F("fHistPsiVZEROCombLeadingJet", "p_{t} [GeV/c]", "#Psi_{VZEROComb}", 350, -100, 250, 50, -1.*TMath::Pi()/2., TMath::Pi()/2., i);
+            fHistPsi2Correlation[i] = BookTH3F("fHistPsi2Correlation", "#Psi_{TPC}", "#Psi_{VZEROA}", "#Psi_{VZEROC}",  20, -1.*TMath::Pi()/2., TMath::Pi()/2., 20, -1.*TMath::Pi()/2., TMath::Pi()/2., 20, -1.*TMath::Pi()/2., TMath::Pi()/2., i);
         }
     }
 
@@ -1758,7 +1784,7 @@ void AliAnalysisTaskJetV2::FillClusterHistograms() const
 //_____________________________________________________________________________
 void AliAnalysisTaskJetV2::FillEventPlaneHistograms(Double_t vzero[2][2], Double_t* vzeroComb, Double_t* tpc) const
 {
-    // fill event plane histograms
+    // fill event plane histograms, only called in qa mode
     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
     fHistPsiControl->Fill(0.5, vzero[0][0]);    // vzero a psi2
     fHistPsiControl->Fill(1.5, vzero[1][0]);    // vzero c psi2
@@ -1793,7 +1819,8 @@ void AliAnalysisTaskJetV2::FillEventPlaneHistograms(Double_t vzero[2][2], Double
         fHistPsiVZEROCLeadingJet[fInCentralitySelection]->Fill(pt, vzero[1][0]);
         fHistPsiVZEROCombLeadingJet[fInCentralitySelection]->Fill(pt, vzeroComb[0]);
     }
-
+    // correlation of event planes
+    fHistPsi2Correlation[fInCentralitySelection]->Fill(tpc[0], vzero[0][0], vzero[1][0]);
 }
 //_____________________________________________________________________________
 void AliAnalysisTaskJetV2::FillRhoHistograms()
index 033ed99a05a3049bd5c2116ef9c8fb6d7db48d69..2fbb4b40aaa08175ba08337729db32bbf128e00f 100644 (file)
@@ -19,8 +19,9 @@
 
 class TFile;
 class TF1;
-class THF1;
-class THF2;
+class TH1F;
+class TH2F;
+class TH3F;
 class TProfile;
 class AliLocalRhoParameter;
 class AliClusterContainer;
@@ -50,6 +51,7 @@ class AliAnalysisTaskJetV2 : public AliAnalysisTaskEmcalJet {
         virtual Bool_t          Run();
         TH1F*                   BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Int_t c = -1, Bool_t append = kTRUE);
         TH2F*                   BookTH2F(const char* name, const char* x, const char* y, Int_t binsx, Double_t minx, Double_t maxx, Int_t binsy, Double_t miny, Double_t maxy, Int_t c = -1, Bool_t append = kTRUE);
+        TH3F*                   BookTH3F(const char* name, const char* x, const char* y, const char* z, Int_t binsx, Double_t minx, Double_t maxx, Int_t binsy, Double_t miny, Double_t maxy, Int_t binsz, Double_t minz, Double_t maxz, Int_t c = -1, Bool_t append = kTRUE);
         /* inline */    static Double_t PhaseShift(Double_t x) {  
             while (x>=TMath::TwoPi())x-=TMath::TwoPi();
             while (x<0.)x+=TMath::TwoPi();
@@ -354,6 +356,7 @@ class AliAnalysisTaskJetV2 : public AliAnalysisTaskEmcalJet {
         TH2F*                   fHistPsiVZEROALeadingJet[10];   //! correlation vzeroa EP, LJ pt
         TH2F*                   fHistPsiVZEROCLeadingJet[10];   //! correlation vzeroc EP, LJ pt
         TH2F*                   fHistPsiVZEROCombLeadingJet[10];//! correlation vzerocomb EP, LJ pt
+        TH3F*                   fHistPsi2Correlation[10];       //! correlation of event planes
         // background
         TH1F*                   fHistRhoPackage[10];    //! rho as estimated by emcal jet package
         TH1F*                   fHistRho[10];           //! background