]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/AnalysisMacros/Local/AliAnalysisTaskPt.cxx
2191befe83d8a255dd57efb9591cfccf0961e471
[u/mrichter/AliRoot.git] / PWG2 / AnalysisMacros / Local / AliAnalysisTaskPt.cxx
1 #define AliAnalysisTaskPt_cxx
2 #include "TROOT.h"
3 #include "TChain.h"
4 #include "TH1.h"
5 #include "TCanvas.h"
6 #include "TSystem.h"
7
8 #include "AliAnalysisTask.h"
9
10 #include "AliESD.h"
11
12 #include "AliAnalysisTaskPt.h"
13
14 ClassImp(AliAnalysisTaskPt)
15
16 //________________________________________________________________________
17 AliAnalysisTaskPt::AliAnalysisTaskPt(const char *name) :AliAnalysisTask(name,""), fESD(0), fHistPt(0) {
18   // Constructor.
19   // Input slot #0 works with an Ntuple
20   DefineInput(0, TChain::Class());
21   // Output slot #0 writes into a TH1 container
22   DefineOutput(0, TH1F::Class());
23 }
24
25 //________________________________________________________________________
26 void AliAnalysisTaskPt::ConnectInputData(Option_t *) {
27   printf("   ConnectInputData %s\n", GetName());
28
29   char ** address = (char **)GetBranchAddress(0, "ESD");
30   if (address) {
31     fESD = (AliESD*)(*address);
32   }
33   else  {
34     fESD = new AliESD();
35     SetBranchAddress(0, "ESD", &fESD);
36   }
37 }
38
39 //________________________________________________________________________
40 void AliAnalysisTaskPt::CreateOutputObjects() {
41   if (!fHistPt) {
42     fHistPt = new TH1F("fHistPt","This is the Pt distribution",15,0.1,3.1);
43     fHistPt->SetStats(kTRUE);
44     fHistPt->GetXaxis()->SetTitle("P_{T} [GeV]");
45     fHistPt->GetYaxis()->SetTitle("#frac{dN}{dP_{T}}");
46     fHistPt->GetXaxis()->SetTitleColor(1);
47     fHistPt->SetMarkerStyle(kFullCircle);
48   }
49 }
50
51 //________________________________________________________________________
52 void AliAnalysisTaskPt::Exec(Option_t *) {
53   // Task making a pt distribution.
54   // Get input data
55   TChain *chain = (TChain*)GetInputData(0);
56   Long64_t ientry = chain->GetReadEntry();
57   if (!fESD) return;
58
59   printf("Tracks: %d \n",fESD->GetNumberOfTracks());
60
61   for(Int_t iTracks = 0; iTracks < fESD->GetNumberOfTracks(); iTracks++) {
62     AliESDtrack * ESDTrack = fESD->GetTrack(iTracks);
63     //UInt_t status = ESDTrack->GetStatus();
64     Double_t momentum[3];
65     ESDTrack->GetPxPyPz(momentum);
66     Double_t Pt = sqrt(pow(momentum[0],2) + pow(momentum[1],2));
67     fHistPt->Fill(Pt);
68   }//track loop 
69   // Post final data. It will be written to a file with option "RECREATE"
70   PostData(0, fHistPt);
71 }      
72
73 //________________________________________________________________________
74 void AliAnalysisTaskPt::Terminate(Option_t *) {
75   // Draw some histogram at the end.
76   if (!gROOT->IsBatch()) {
77     TCanvas *c1 = new TCanvas("c1","Pt",10,10,310,310);
78     c1->SetFillColor(10);
79     c1->SetHighLightColor(10);
80     
81     c1->cd(1)->SetLeftMargin(0.15);
82     c1->cd(1)->SetBottomMargin(0.15);  
83     c1->cd(1)->SetLogy();
84     //fHistPt = (TH1F*)GetOutputData(0);
85     if (fHistPt) fHistPt->DrawCopy("E");
86   }
87 }