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