]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/SPECTRA/AliAnalysisTaskProtons.cxx
Validated macros for the global analysis train...
[u/mrichter/AliRoot.git] / PWG2 / SPECTRA / AliAnalysisTaskProtons.cxx
1 #include "TStyle.h"
2 #include "TChain.h"
3 #include "TTree.h"
4 #include "TString.h"
5 #include "TList.h"
6 #include "TFile.h"
7 #include "TH2F.h"
8 #include "TCanvas.h"
9
10 #include "AliAnalysisTask.h"
11 #include "AliAnalysisManager.h"
12
13 #include "AliESDEvent.h"
14 #include "AliESDInputHandler.h"
15 #include "AliAODEvent.h"
16 #include "AliAODInputHandler.h"
17 #include "AliMCEventHandler.h"
18 #include "AliMCEvent.h"
19 #include "AliStack.h"
20 #include "AliESDVertex.h"
21
22 #include "AliProtonAnalysis.h"
23 #include "AliProtonAnalysisBase.h"
24 #include "AliAnalysisTaskProtons.h"
25
26 //-----------------------------------------------------------------
27 //                 AliAnalysisTakProtons class
28 //   This is the task to run the \bar{p}/p analysis
29 //   Origin: Panos Christakoglou | Panos.Christakoglou@cern.ch
30 //-----------------------------------------------------------------
31
32 ClassImp(AliAnalysisTaskProtons)
33   
34 //________________________________________________________________________ 
35 AliAnalysisTaskProtons::AliAnalysisTaskProtons()
36   : AliAnalysisTask(), fESD(0), fAOD(0), fMC(0),
37     fList(0), fProtonAnalysis(0) {
38   //Dummy constructor
39   
40 }
41
42 //________________________________________________________________________
43 AliAnalysisTaskProtons::AliAnalysisTaskProtons(const char *name) 
44   : AliAnalysisTask(name, ""), fESD(0), fAOD(0), fMC(0),
45     fList(0), fProtonAnalysis(0) {
46   // Constructor
47   
48   // Define input and output slots here
49   // Input slot #0 works with a TChain
50   DefineInput(0, TChain::Class());
51   // Output slot #0 writes into a TList container
52   DefineOutput(0, TList::Class());
53 }
54
55 //________________________________________________________________________
56 void AliAnalysisTaskProtons::ConnectInputData(Option_t *) {
57   // Connect ESD or AOD here
58   // Called once
59   TString gAnalysisLevel = dynamic_cast<AliProtonAnalysisBase*>(fProtonAnalysis->GetProtonAnalysisBaseObject())->GetAnalysisLevel(); 
60
61   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
62   if (!tree) {
63     Printf("ERROR: Could not read chain from input slot 0");
64   } else {
65     if(gAnalysisLevel == "ESD") {
66       AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
67      
68       if (!esdH) {
69         Printf("ERROR: Could not get ESDInputHandler");
70       } else
71         fESD = esdH->GetEvent();
72     }
73     else if(gAnalysisLevel == "AOD") {
74      AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
75      
76      if (!aodH) {
77        Printf("ERROR: Could not get AODInputHandler");
78      } else
79        fAOD = aodH->GetEvent();
80     }
81     else if(gAnalysisLevel == "MC") {
82      AliMCEventHandler* mcH = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
83      if (!mcH) {
84        Printf("ERROR: Could not retrieve MC event handler");
85      }
86      else
87        fMC = mcH->MCEvent();
88     }
89     else
90       Printf("Wrong analysis type: Only ESD, AOD and MC types are allowed!");
91   }
92 }
93
94 //________________________________________________________________________
95 void AliAnalysisTaskProtons::CreateOutputObjects() {
96   // Create output objects
97   // Called once
98   fList = new TList();
99   fList->Add(fProtonAnalysis->GetProtonYPtHistogram());
100   fList->Add(fProtonAnalysis->GetAntiProtonYPtHistogram());
101   fList->Add(fProtonAnalysis->GetEventHistogram());
102   fList->Add(fProtonAnalysis->GetProtonContainer());
103   fList->Add(fProtonAnalysis->GetAntiProtonContainer());
104 }
105
106 //________________________________________________________________________
107 void AliAnalysisTaskProtons::Exec(Option_t *) {
108   // Main loop
109   // Called for each event
110   TString gAnalysisLevel = dynamic_cast<AliProtonAnalysisBase*>(fProtonAnalysis->GetProtonAnalysisBaseObject())->GetAnalysisLevel(); 
111   
112   if(gAnalysisLevel == "ESD") {
113     if (!fESD) {
114       Printf("ERROR: fESD not available");
115       return;
116     }
117
118     if(dynamic_cast<AliProtonAnalysisBase*>(fProtonAnalysis->GetProtonAnalysisBaseObject())->IsEventTriggered(fESD,dynamic_cast<AliProtonAnalysisBase*>(fProtonAnalysis->GetProtonAnalysisBaseObject())->GetTriggerMode())) {
119       const AliESDVertex *vertex = dynamic_cast<AliProtonAnalysisBase*>(fProtonAnalysis->GetProtonAnalysisBaseObject())->GetVertex(fESD,dynamic_cast<AliProtonAnalysisBase*>(fProtonAnalysis->GetProtonAnalysisBaseObject())->GetAnalysisMode(),dynamic_cast<AliProtonAnalysisBase*>(fProtonAnalysis->GetProtonAnalysisBaseObject())->GetVxMax(),dynamic_cast<AliProtonAnalysisBase*>(fProtonAnalysis->GetProtonAnalysisBaseObject())->GetVyMax(),dynamic_cast<AliProtonAnalysisBase*>(fProtonAnalysis->GetProtonAnalysisBaseObject())->GetVzMax());
120       if(vertex) {
121         Printf("Proton ESD analysis task: There are %d tracks in this event", fESD->GetNumberOfTracks());
122         fProtonAnalysis->Analyze(fESD,vertex);
123       }//reconstructed vertex
124     }//triggered event
125   }//ESD analysis              
126   
127   else if(gAnalysisLevel == "AOD") {
128     if (!fAOD) {
129       Printf("ERROR: fAOD not available");
130       return;
131     }
132     
133     Printf("Proton AOD analysis task: There are %d tracks in this event", fAOD->GetNumberOfTracks());
134     fProtonAnalysis->Analyze(fAOD);
135   }//AOD analysis
136   
137   else if(gAnalysisLevel == "MC") {
138     if (!fMC) {
139       Printf("ERROR: Could not retrieve MC event");
140       return;
141     }
142
143     AliStack* stack = fMC->Stack();
144     if (!stack) {
145       Printf("ERROR: Could not retrieve the stack");
146       return;
147     }
148     Printf("Proton MC analysis task: There are %d primaries in this event", stack->GetNprimary());
149     fProtonAnalysis->Analyze(stack,kFALSE);//kTRUE in case of inclusive measurement
150   }//MC analysis                      
151
152   // Post output data.
153   PostData(0, fList);
154 }      
155
156 //________________________________________________________________________
157 void AliAnalysisTaskProtons::Terminate(Option_t *) {
158   // Draw result to the screen
159   // Called once at the end of the query
160   gStyle->SetPalette(1,0);
161
162   fList = dynamic_cast<TList*> (GetOutputData(0));
163   if (!fList) {
164     Printf("ERROR: fList not available");
165     return;
166   }
167    
168   TH2F *fHistYPtProtons = (TH2F *)fList->At(0);
169   TH2F *fHistYPtAntiProtons = (TH2F *)fList->At(1);
170     
171   TCanvas *c1 = new TCanvas("c1","p-\bar{p}",200,0,800,400);
172   c1->SetFillColor(10); c1->SetHighLightColor(10); c1->Divide(2,1);
173
174   c1->cd(1)->SetLeftMargin(0.15); c1->cd(1)->SetBottomMargin(0.15);  
175   if (fHistYPtProtons) fHistYPtProtons->DrawCopy("colz");
176   c1->cd(2)->SetLeftMargin(0.15); c1->cd(2)->SetBottomMargin(0.15);  
177   if (fHistYPtAntiProtons) fHistYPtAntiProtons->DrawCopy("colz");
178
179   TCanvas *c2 = dynamic_cast<AliProtonAnalysisBase*>(fProtonAnalysis->GetProtonAnalysisBaseObject())->GetListOfCuts();
180   TFile *flocal = TFile::Open("ListOfCuts.root","recreate");
181   c2->Write();
182   flocal->Close();
183 }
184