]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisTaskPhiCorr.cxx
Disable Evchar if no TR read.
[u/mrichter/AliRoot.git] / ANALYSIS / AliAnalysisTaskPhiCorr.cxx
1
2 #include "TH1F.h"
3 #include "TList.h"
4 #include "TCanvas.h"
5 #include "AliAODEvent.h"
6 #include "AliAODTrack.h"
7 #include "AliAnalysisTaskPhiCorr.h"
8 #include "AliMultiEventInputHandler.h"
9
10 ClassImp(AliAnalysisTaskPhiCorr)
11
12 //________________________________________________________________________
13 AliAnalysisTaskPhiCorr::AliAnalysisTaskPhiCorr(const char *name) 
14     : AliAnalysisTaskME(name), fHists(0), fHistDphiCO(0),  fHistDphiUC(0)
15 {
16   // Constructor
17   // Define input and output slots here
18   DefineOutput(1, TList::Class());
19 }
20
21
22 //________________________________________________________________________
23 void AliAnalysisTaskPhiCorr::UserCreateOutputObjects()
24 {
25   // Create histograms
26   // Called once
27
28   fHistDphiCO = new TH1F("fHistDPhiCO", "#Delta #Phi distribution", 64, 0., 3.2);
29   fHistDphiCO->GetXaxis()->SetTitle("#Delta#Phi [rad]");
30   fHistDphiCO->GetYaxis()->SetTitle("dN/d#Phi");
31   fHistDphiCO->SetMarkerStyle(kFullCircle);
32
33   fHistDphiUC = new TH1F("fHistDPhiUC", "#Delta #Phi distribution", 64, 0., 3.2);
34   fHistDphiUC->GetXaxis()->SetTitle("#Delta#Phi [rad]");
35   fHistDphiUC->GetYaxis()->SetTitle("dN/d#Phi");
36   fHistDphiUC->SetMarkerStyle(kOpenCircle);
37   
38   fHists = new TList();
39   fHists->Add(fHistDphiCO);
40   fHists->Add(fHistDphiUC);  
41 }
42
43 //________________________________________________________________________
44 void AliAnalysisTaskPhiCorr::UserExec(Option_t *) 
45 {
46   // Main loop
47   // Called for each event
48     // Uncorrelated tracks
49     Int_t nev = fInputHandler->GetBufferSize();
50     Float_t wgt = 1./(nev*(nev-1));
51     fMixedEvent.Reset();
52     
53     for (Int_t iev = 0; iev < nev; iev++) {
54         AliAODEvent* aod = (AliAODEvent*) GetEvent(iev);
55         fMixedEvent.AddEvent(aod);
56     }
57     fMixedEvent.Init();
58     Int_t ntrack = fMixedEvent.GetNumberOfTracks();
59     if (ntrack > 1) {
60       for (Int_t itr = 0; itr < ntrack -1; itr++) {
61         for (Int_t jtr = itr+1; jtr < ntrack; jtr++) {
62           AliVParticle* track1 = fMixedEvent.GetTrack(itr);
63           AliVParticle* track2 = fMixedEvent.GetTrack(jtr);
64           Int_t iev1 =  fMixedEvent.EventIndex(itr);
65           Int_t iev2 =  fMixedEvent.EventIndex(jtr);
66
67           Float_t phi1 = track1->Phi();
68           Float_t phi2 = track2->Phi();
69           Float_t dphi = TMath::Abs(phi1 - phi2);
70           if (dphi > TMath::Pi()) dphi = 2. * TMath::Pi() - dphi;
71           if (iev1 != iev2) {
72             fHistDphiUC->Fill(dphi, wgt);
73           } else {
74             fHistDphiCO->Fill(dphi, 0.5);           
75           }
76         } // tarcks
77       } // tracks
78     } // more than 1 
79   PostData(1, fHists);
80 }      
81
82 //________________________________________________________________________
83 void AliAnalysisTaskPhiCorr::Terminate(Option_t *) 
84 {
85   // Draw result to the screen
86   // Called once at the end of the query
87
88   TCanvas *c1 = new TCanvas("AliAnalysisTaskPt","Pt",10,10,510,510);
89   c1->cd(1)->SetLogy();
90   fHistDphiCO->DrawCopy("E");
91   fHistDphiUC->DrawCopy("Esame");
92 }