]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/AliAnalysisTaskProtons.cxx
Using AliMagWrapCheb
[u/mrichter/AliRoot.git] / PWG2 / AliAnalysisTaskProtons.cxx
1 #include "TChain.h"
2 #include "TTree.h"
3 #include "TString.h"
4 #include "TList.h"
5 #include "TH2F.h"
6 #include "TH1I.h"
7 #include "TF1.h"
8 #include "TCanvas.h"
9 #include "TLorentzVector.h"
10
11 #include "AliAnalysisTask.h"
12 #include "AliAnalysisManager.h"
13
14 #include "AliESDEvent.h"
15 #include "AliESDInputHandler.h"
16 #include "AliAODEvent.h"
17 #include "AliAODInputHandler.h"
18 #include "AliMCEventHandler.h"
19 #include "AliMCEvent.h"
20 #include "AliStack.h"
21 #include "AliCFContainer.h"
22
23 #include "AliProtonAnalysis.h"
24 #include "AliAnalysisTaskProtons.h"
25
26 // Analysis task creating a the 2d y-p_t spectrum of p and antip
27 // Author: Panos Cristakoglou
28
29 ClassImp(AliAnalysisTaskProtons)
30
31 //________________________________________________________________________ 
32 AliAnalysisTaskProtons::AliAnalysisTaskProtons()
33   : AliAnalysisTask(), fESD(0), fAOD(0), fMC(0), fAnalysisType("ESD"),
34     fList(0), fAnalysis(0),
35     fElectronFunction(0), fMuonFunction(0),
36     fPionFunction(0), fKaonFunction(0), fProtonFunction(0),
37     fFunctionUsed(kFALSE) {
38   //Dummy constructor
39                                                                                                    
40 }
41
42 //________________________________________________________________________
43 AliAnalysisTaskProtons::AliAnalysisTaskProtons(const char *name) 
44 : AliAnalysisTask(name, ""), fESD(0), fAOD(0), fMC(0), fAnalysisType("ESD"),
45   fList(0), fAnalysis(0), 
46   fElectronFunction(0), fMuonFunction(0),
47   fPionFunction(0), fKaonFunction(0), fProtonFunction(0),
48   fFunctionUsed(kFALSE) { 
49   // Constructor
50
51   // Define input and output slots here
52   // Input slot #0 works with a TChain
53   DefineInput(0, TChain::Class());
54   // Output slot #0 writes into a TList container
55   DefineOutput(0, TList::Class());
56 }
57
58 //________________________________________________________________________
59 void AliAnalysisTaskProtons::ConnectInputData(Option_t *) {
60   // Connect ESD or AOD here
61   // Called once
62
63   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
64   if (!tree) {
65     Printf("ERROR: Could not read chain from input slot 0");
66   } else {
67     if(fAnalysisType == "ESD") {
68       tree->SetBranchStatus("*", kFALSE);
69       tree->SetBranchStatus("Tracks.*", kTRUE);
70       
71       AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
72      
73       if (!esdH) {
74         Printf("ERROR: Could not get ESDInputHandler");
75       } else
76         fESD = esdH->GetEvent();
77     }
78     else if(fAnalysisType == "AOD") {
79      AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
80      
81      if (!aodH) {
82        Printf("ERROR: Could not get AODInputHandler");
83      } else
84        fAOD = aodH->GetEvent();
85     }
86     else if(fAnalysisType == "MC") {
87      AliMCEventHandler* mcH = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
88      if (!mcH) {
89        Printf("ERROR: Could not retrieve MC event handler");
90      }
91      else
92        fMC = mcH->MCEvent();
93     }
94     else
95       Printf("Wrong analysis type: Only ESD, AOD and MC types are allowed!");
96   }
97 }
98
99 //________________________________________________________________________
100 void AliAnalysisTaskProtons::CreateOutputObjects() {
101   // Create histograms
102   // Called once
103   //Prior probabilities
104   Double_t partFrac[5] = {0.01, 0.01, 0.85, 0.10, 0.05};
105   
106   //proton analysis object
107   fAnalysis = new AliProtonAnalysis();
108   fAnalysis->InitAnalysisHistograms(10,-1.0,1.0,26,0.5,3.1);
109
110   if(fAnalysisType == "ESD") {
111     //Use of TPConly tracks
112     fAnalysis->UseTPCOnly();
113
114     //TPC related cuts       
115     fAnalysis->SetMinTPCClusters(50);
116     fAnalysis->SetMaxChi2PerTPCCluster(3.5);
117     fAnalysis->SetMaxCov11(2.0);
118     fAnalysis->SetMaxCov22(2.0);
119     fAnalysis->SetMaxCov33(0.5);
120     fAnalysis->SetMaxCov44(0.5);
121     fAnalysis->SetMaxCov55(2.0);
122     fAnalysis->SetMaxSigmaToVertex(2.5);
123     fAnalysis->SetTPCRefit();
124
125     //ITS related cuts - to be used for the analysis of global tracking 
126     //fAnalysis->SetMinITSClusters(5);
127     //fAnalysis->SetITSRefit();
128   }
129
130   if(fFunctionUsed)
131     fAnalysis->SetPriorProbabilityFunctions(fElectronFunction,
132                                             fMuonFunction,
133                                             fPionFunction,
134                                             fKaonFunction,
135                                             fProtonFunction);
136   else
137     fAnalysis->SetPriorProbabilities(partFrac);
138
139   fList = new TList();
140   fList->Add(fAnalysis->GetProtonYPtHistogram());
141   fList->Add(fAnalysis->GetAntiProtonYPtHistogram());
142   fList->Add(fAnalysis->GetEventHistogram());
143   fList->Add(fAnalysis->GetProtonContainer());
144   fList->Add(fAnalysis->GetAntiProtonContainer());
145 }
146
147 //________________________________________________________________________
148 void AliAnalysisTaskProtons::Exec(Option_t *) {
149   // Main loop
150   // Called for each event
151   
152   if(fAnalysisType == "ESD") {
153     if (!fESD) {
154       Printf("ERROR: fESD not available");
155       return;
156     }
157
158     Printf("Proton ESD analysis task: There are %d tracks in this event", fESD->GetNumberOfTracks());
159     fAnalysis->Analyze(fESD);
160   }//ESD analysis              
161   
162   else if(fAnalysisType == "AOD") {
163     if (!fAOD) {
164       Printf("ERROR: fAOD not available");
165       return;
166     }
167     
168     Printf("Proton AOD analysis task: There are %d tracks in this event", fAOD->GetNumberOfTracks());
169     fAnalysis->Analyze(fAOD);
170   }//AOD analysis                                                                                                                                                                
171   
172   else if(fAnalysisType == "MC") {
173     if (!fMC) {
174       Printf("ERROR: Could not retrieve MC event");
175       return;
176     }
177
178     AliStack* stack = fMC->Stack();
179     if (!stack) {
180       Printf("ERROR: Could not retrieve the stack");
181       return;
182     }
183     Printf("Proton MC analysis task: There are %d primaries in this event", stack->GetNprimary());
184     fAnalysis->Analyze(stack);
185   }//MC analysis                      
186
187   // Post output data.
188   PostData(0, fList);
189 }      
190
191 //________________________________________________________________________
192 void AliAnalysisTaskProtons::Terminate(Option_t *) {
193   // Draw result to the screen
194   // Called once at the end of the query
195   
196   fList = dynamic_cast<TList*> (GetOutputData(0));
197   if (!fList) {
198     Printf("ERROR: fList not available");
199     return;
200   }
201    
202   TH2F *fHistYPtProtons = (TH2F *)fList->At(0);
203   TH2F *fHistYPtAntiProtons = (TH2F *)fList->At(1);
204     
205   TCanvas *c1 = new TCanvas("c1","p-\bar{p}",200,0,800,400);
206   c1->SetFillColor(10); c1->SetHighLightColor(10); c1->Divide(2,1);
207
208   c1->cd(1)->SetLeftMargin(0.15); c1->cd(1)->SetBottomMargin(0.15);  
209   if (fHistYPtProtons) fHistYPtProtons->DrawCopy("colz");
210   c1->cd(2)->SetLeftMargin(0.15); c1->cd(2)->SetBottomMargin(0.15);  
211   if (fHistYPtAntiProtons) fHistYPtAntiProtons->DrawCopy("colz");
212 }
213
214