]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisTaskPhiCorr.cxx
- Info message corrected for common output container
[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 "AliMultiAODInputHandler.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     
52     for (Int_t iev = 0; iev < nev; iev++) {
53         for (Int_t jev = (iev + 1); jev < nev; jev++) {
54             AliAODEvent* aod1 = GetEvent(iev);
55             AliAODEvent* aod2 = GetEvent(jev);
56             
57             Int_t ntracks1 =  aod1->GetNumberOfTracks();
58             Int_t ntracks2 =  aod2->GetNumberOfTracks();
59             
60             printf("Number of tracks %5d:%5d %5d:%5d\n", iev, ntracks1, jev, ntracks2);
61             
62             for (Int_t iTracks = 0; iTracks < ntracks1; iTracks++) {
63                 for (Int_t jTracks = 0; jTracks < ntracks2; jTracks++) {
64                     
65                     
66                     AliAODTrack* track1 = aod1->GetTrack(iTracks);
67                     AliAODTrack* track2 = aod2->GetTrack(jTracks);
68                     
69                     Float_t phi1 = track1->Phi();
70                     Float_t phi2 = track2->Phi();
71                     Float_t dphi = TMath::Abs(phi1 - phi2);
72                     if (dphi > TMath::Pi()) dphi = 2. * TMath::Pi() - dphi;
73                     fHistDphiUC->Fill(dphi, wgt);
74                 } // tracks
75             } // tracks 
76         } // event loop
77     } // event loop
78     
79 //    Correlated 
80     AliAODEvent* aod = fInputHandler->GetLatestEvent();
81     
82     Int_t ntracks = aod->GetNumberOfTracks();
83     for (Int_t iTracks = 0; iTracks < ntracks; iTracks++) {
84         for (Int_t jTracks = (iTracks+1); jTracks < ntracks; jTracks++) {
85               AliAODTrack* track1 = aod->GetTrack(iTracks);
86               AliAODTrack* track2 = aod->GetTrack(jTracks);
87               Float_t phi1 = track1->Phi();
88               Float_t phi2 = track2->Phi();
89               Float_t dphi = TMath::Abs(phi1 - phi2);
90               if (dphi > TMath::Pi()) dphi = 2. * TMath::Pi() - dphi;
91               fHistDphiCO->Fill(dphi);
92         } // tracks
93     }
94     
95   // Post output data.
96   PostData(1, fHists);
97 }      
98
99 //________________________________________________________________________
100 void AliAnalysisTaskPhiCorr::Terminate(Option_t *) 
101 {
102   // Draw result to the screen
103   // Called once at the end of the query
104
105   TCanvas *c1 = new TCanvas("AliAnalysisTaskPt","Pt",10,10,510,510);
106   c1->cd(1)->SetLogy();
107   fHistDphiUC->DrawCopy("E");
108   fHistDphiCO->DrawCopy("Esame");
109 }