]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/EBYE/AliAnalysisTaskBF.cxx
parameter set for crude clustering and refined clustering
[u/mrichter/AliRoot.git] / PWG2 / EBYE / AliAnalysisTaskBF.cxx
1 #include "TChain.h"
2 #include "TTree.h"
3 #include "TCanvas.h"
4 #include "TLorentzVector.h"
5 #include "TGraphErrors.h"
6
7 #include "AliAnalysisTask.h"
8 #include "AliAnalysisManager.h"
9
10 #include "AliESDEvent.h"
11 #include "AliESDInputHandler.h"
12 #include "AliAODEvent.h"
13 #include "AliAODInputHandler.h"
14 #include "AliMCEventHandler.h"
15 #include "AliMCEvent.h"
16 #include "AliStack.h"
17
18 #include "AliBalance.h"
19
20 #include "AliAnalysisTaskBF.h"
21
22 // Analysis task for the BF code
23 // Authors: Panos Cristakoglou@cern.ch
24
25 ClassImp(AliAnalysisTaskBF)
26
27 //________________________________________________________________________
28 AliAnalysisTaskBF::AliAnalysisTaskBF(const char *name) 
29   : AliAnalysisTask(name, ""), fESD(0), fAOD(0), fMC(0), fBalance(0) {
30   // Constructor
31
32   // Define input and output slots here
33   // Input slot #0 works with a TChain
34   DefineInput(0, TChain::Class());
35   // Output slot #0 writes into a TH1 container
36   DefineOutput(0, AliBalance::Class());
37 }
38
39 //________________________________________________________________________
40 void AliAnalysisTaskBF::ConnectInputData(Option_t *) {
41   // Connect ESD or AOD here
42   // Called once
43   TString gAnalysisLevel = fBalance->GetAnalysisLevel();
44
45   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
46   if (!tree) {
47     Printf("ERROR: Could not read chain from input slot 0");
48   } 
49   else {
50     if(gAnalysisLevel == "ESD") {
51       AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
52       
53       if (!esdH) {
54         Printf("ERROR: Could not get ESDInputHandler");
55       } 
56       else
57         fESD = esdH->GetEvent();
58     }//ESD
59     else if(gAnalysisLevel == "AOD") {
60       AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
61       
62       if (!aodH) {
63         Printf("ERROR: Could not get AODInputHandler");
64       } else
65         fAOD = aodH->GetEvent();
66     }//AOD
67     else if(gAnalysisLevel == "MC") {
68       AliMCEventHandler* mcH = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
69       if (!mcH) {
70         Printf("ERROR: Could not retrieve MC event handler");
71       }
72       else
73         fMC = mcH->MCEvent();
74     }//MC
75     else
76       Printf("Wrong analysis type: Only ESD, AOD and MC types are allowed!");
77   }
78 }
79
80 //________________________________________________________________________
81 void AliAnalysisTaskBF::CreateOutputObjects() {
82   // Create histograms
83   // Called once
84   
85 }
86
87 //________________________________________________________________________
88 void AliAnalysisTaskBF::Exec(Option_t *) {
89   // Main loop
90   // Called for each event
91   TString gAnalysisLevel = fBalance->GetAnalysisLevel();
92
93   TObjArray *array = new TObjArray();
94   
95   //ESD analysis
96   if(gAnalysisLevel == "ESD") {
97     if (!fESD) {
98       Printf("ERROR: fESD not available");
99       return;
100     }
101     
102     Printf("There are %d tracks in this event", fESD->GetNumberOfTracks());
103     for (Int_t iTracks = 0; iTracks < fESD->GetNumberOfTracks(); iTracks++) {
104       AliESDtrack* track = fESD->GetTrack(iTracks);
105       if (!track) {
106         Printf("ERROR: Could not receive track %d", iTracks);
107         continue;
108       }
109       array->Add(track);
110     } //track loop
111   }//ESD analysis
112   //AOD analysis
113   else if(gAnalysisLevel == "AOD") {
114     if (!fAOD) {
115       Printf("ERROR: fAOD not available");
116       return;
117     }
118     
119     Printf("There are %d tracks in this event", fAOD->GetNumberOfTracks());
120     for (Int_t iTracks = 0; iTracks < fAOD->GetNumberOfTracks(); iTracks++) {
121       AliAODtrack* track = fAOD->GetTrack(iTracks);
122       if (!track) {
123         Printf("ERROR: Could not receive track %d", iTracks);
124         continue;
125       }
126       array->Add(track);
127     } //track loop
128   }//AOD analysis
129   //MC analysis
130   else if(gAnalysisLevel == "MC") {
131     if (!fMC) {
132       Printf("ERROR: fMC not available");
133       return;
134     }
135     
136     Printf("There are %d tracks in this event", fMC->GetNumberOfPrimaries());
137     for (Int_t iTracks = 0; iTracks < fMC->GetNumberOfPrimaries(); iTracks++) {
138       AliMCParticle* track = fMC->GetTrack(iTracks);
139       if (!track) {
140         Printf("ERROR: Could not receive particle %d", iTracks);
141         continue;
142       }
143       array->Add(track);
144     } //track loop
145   }//MC analysis
146
147   fBalance->CalculateBalance(array);
148   
149   delete array;
150
151   // Post output data.
152   PostData(0, fBalance);
153 }      
154
155 //________________________________________________________________________
156 void AliAnalysisTaskBF::Terminate(Option_t *) {
157   // Draw result to the screen
158   // Called once at the end of the query
159
160   fBalance = dynamic_cast<AliBalance*> (GetOutputData(0));
161   if (!fBalance) {
162     Printf("ERROR: fBalance not available");
163     return;
164   }
165   
166   TGraphErrors *gr = fBalance->DrawBalance();
167   gr->SetMarkerStyle(20);
168   gr->Draw("AP");
169
170   fBalance->PrintResults();
171 }