]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/AnalysisMacros/Kine/AliAnalysisTaskRLPt.cxx
Adding the example on how to access the kine tree with tha TaskRL
[u/mrichter/AliRoot.git] / PWG2 / AnalysisMacros / Kine / AliAnalysisTaskRLPt.cxx
1 #define AliAnalysisTaskRLPt_cxx
2
3 #include "TChain.h"
4 #include "TH1.h"
5 #include "TCanvas.h"
6 #include "TSystem.h"
7 #include "TROOT.h"
8 #include "TParticle.h"
9
10 #include "AliESD.h"
11 #include "AliLog.h"
12 #include "AliStack.h"
13
14 #include "AliAnalysisTaskRL.h"
15 #include "AliAnalysisTaskRLPt.h"
16
17 ClassImp(AliAnalysisTaskRLPt)
18
19 //________________________________________________________________________
20 AliAnalysisTaskRLPt::AliAnalysisTaskRLPt(const char *name) :AliAnalysisTaskRL(name,""), fESD(0), fHistPt(0) {
21   // Constructor.
22   // Input slot #0 works with an Ntuple
23   DefineInput(0, TChain::Class());
24   // Output slot #0 writes into a TH1 container
25   DefineOutput(0, TH1F::Class());
26 }
27
28 //___________________________________________________________________________
29 void AliAnalysisTaskRLPt::ConnectInputData(Option_t *) {
30   // Initialize branches.
31   printf("   ConnectInputData of task %s\n", GetName());
32   if (!fESD) {
33     char ** address = (char **)GetBranchAddress(0, "ESD");
34     if (address) fESD = (AliESD*)(*address);
35     if (!fESD) {
36       fESD = new AliESD();
37       SetBranchAddress(0, "ESD", &fESD);
38     }
39   }
40 }
41
42 //___________________________________________________________________________
43 void AliAnalysisTaskRLPt::CreateOutputObjects() {
44   printf("   CreateOutputObjects of task %s\n", GetName());
45   if (!fHistPt) {
46     fHistPt = new TH1F("fHistPt","This is the Pt distribution",15,0.1,3.1);
47     fHistPt->SetStats(kTRUE);
48     fHistPt->GetXaxis()->SetTitle("P_{T} [GeV]");
49     fHistPt->GetYaxis()->SetTitle("#frac{dN}{dP_{T}}");
50     fHistPt->GetXaxis()->SetTitleColor(1);
51     fHistPt->SetMarkerStyle(kFullCircle);
52   }
53 }
54
55 //________________________________________________________________________
56 void AliAnalysisTaskRLPt::Exec(Option_t *) {
57   // Task making a pt distribution.
58   // Get input data
59   TTree *tinput = (TTree*)GetInputData(0);
60   Long64_t ientry = tinput->GetReadEntry();
61   if (AliAnalysisTaskRL::GetEntry(ientry) == kFALSE) {
62     printf("Couldn't get event from the runLoader\n");
63     return;
64   }  
65   if (!fESD) return;
66
67   AliStack* stack = GetStack();
68   if (!stack) {
69   AliDebug(AliLog::kError, "Stack not available");
70     //return kFALSE;
71   }
72   // loop over mc particles
73   Int_t nPrim = stack->GetNprimary();
74   printf("Particles: %d - Tracks: %d \n",nPrim,fESD->GetNumberOfTracks());
75   for(Int_t i = 0; i < nPrim; i++) {
76     TParticle * particle = stack->Particle(i); 
77     if(TMath::Abs(particle->Eta()) > 1.0) continue;
78     fHistPt->Fill(particle->Pt());
79   }
80
81   // Post final data. It will be written to a file with option "RECREATE"
82   PostData(0, fHistPt);
83 }      
84
85 //________________________________________________________________________
86 void AliAnalysisTaskRLPt::Terminate(Option_t *) {
87   // Draw some histogram at the end.
88   if (!gROOT->IsBatch()) {
89     TCanvas *c1 = new TCanvas("c1","Pt",10,10,310,310);
90     c1->SetFillColor(10); c1->SetHighLightColor(10);
91     c1->cd(1)->SetLeftMargin(0.15); c1->cd(1)->SetBottomMargin(0.15);  
92     c1->cd(1)->SetLogy();
93     fHistPt->DrawCopy("E");
94   }
95 }