]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/AliAnalysisTaskProtons.cxx
Fix includes
[u/mrichter/AliRoot.git] / PWG2 / AliAnalysisTaskProtons.cxx
1 #include "TChain.h"
2 #include "TTree.h"
3 #include "TList.h"
4 #include "TH2F.h"
5 #include "TF1.h"
6 #include "TCanvas.h"
7 #include "TLorentzVector.h"
8
9 #include "AliAnalysisTask.h"
10 #include "AliAnalysisManager.h"
11
12 #include "AliESDEvent.h"
13 #include "AliESDInputHandler.h"
14 #include "AliAODEvent.h"
15 #include "AliAODInputHandler.h"
16
17 #include "AliProtonAnalysis.h"
18 #include "AliAnalysisTaskProtons.h"
19
20 // Analysis task creating a the 2d y-p_t spectrum of p and antip
21 // Author: Panos Cristakoglou
22
23 ClassImp(AliAnalysisTaskProtons)
24
25 //________________________________________________________________________
26 AliAnalysisTaskProtons::AliAnalysisTaskProtons(const char *name) 
27 : AliAnalysisTask(name, ""), fESD(0), fAOD(0), fAnalysisType("ESD"), 
28   fList(0), fAnalysis(0), 
29   fElectronFunction(0), fMuonFunction(0),
30   fPionFunction(0), fKaonFunction(0), fProtonFunction(0),
31   fFunctionUsed(kFALSE) { 
32   // Constructor
33
34   // Define input and output slots here
35   // Input slot #0 works with a TChain
36   DefineInput(0, TChain::Class());
37   // Output slot #0 writes into a TList container
38   DefineOutput(0, TList::Class());
39 }
40
41 //________________________________________________________________________
42 void AliAnalysisTaskProtons::ConnectInputData(Option_t *) {
43   // Connect ESD or AOD here
44   // Called once
45
46   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
47   if (!tree) {
48     Printf("ERROR: Could not read chain from input slot 0");
49   } else {
50     // Disable all branches and enable only the needed ones
51     // The next two lines are different when data produced as AliESDEvent is read
52     if(fAnalysisType == "ESD") {
53 // In train mode branches can be disabled at the level of ESD handler (M.G.)
54 //      tree->SetBranchStatus("*", kFALSE);
55       tree->SetBranchStatus("*Tracks.*", kTRUE);
56
57       AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
58       
59       if (!esdH) {
60         Printf("ERROR: Could not get ESDInputHandler");
61       } else
62         fESD = esdH->GetEvent();
63     }
64     else if(fAnalysisType == "AOD") {
65       AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
66       
67       if (!aodH) {
68         Printf("ERROR: Could not get AODInputHandler");
69       } else
70         fAOD = aodH->GetEvent();
71     }
72     else 
73       Printf("Wrong analysis type: Only ESD and AOD types are allowed!");
74   }
75 }
76
77 //________________________________________________________________________
78 void AliAnalysisTaskProtons::CreateOutputObjects() {
79   // Create histograms
80   // Called once
81   Double_t partFrac[5] = {0.01, 0.01, 0.85, 0.10, 0.05};
82
83   fAnalysis = new AliProtonAnalysis();
84   fAnalysis->InitHistograms(10,-1.0,1.0,30,0.1,3.1);
85   if(fAnalysisType == "ESD") {
86     fAnalysis->SetMinTPCClusters(50);
87     fAnalysis->SetMinITSClusters(1);
88     fAnalysis->SetMaxChi2PerTPCCluster(3.5);
89     fAnalysis->SetMaxCov11(2.0);
90     fAnalysis->SetMaxCov22(2.0);
91     fAnalysis->SetMaxCov33(0.5);
92     fAnalysis->SetMaxCov44(0.5);
93     fAnalysis->SetMaxCov55(2.0);
94     fAnalysis->SetMaxSigmaToVertex(3.);
95     fAnalysis->SetITSRefit();
96     fAnalysis->SetTPCRefit();
97   }
98   if(fFunctionUsed)
99     fAnalysis->SetPriorProbabilityFunctions(fElectronFunction,
100                                             fMuonFunction,
101                                             fPionFunction,
102                                             fKaonFunction,
103                                             fProtonFunction);
104   else
105     fAnalysis->SetPriorProbabilities(partFrac);
106
107   fList = new TList();
108   fList->Add(fAnalysis->GetProtonYPtHistogram()); 
109   fList->Add(fAnalysis->GetAntiProtonYPtHistogram()); 
110 }
111
112 //________________________________________________________________________
113 void AliAnalysisTaskProtons::Exec(Option_t *) {
114   // Main loop
115   // Called for each event
116   
117   if(fAnalysisType == "ESD") {
118     if (!fESD) {
119       Printf("ERROR: fESD not available");
120       return;
121     }
122   
123     Printf("Proton analysis task: There are %d tracks in this event", fESD->GetNumberOfTracks());
124     fAnalysis->Analyze(fESD);
125   }
126   else if(fAnalysisType == "AOD") {
127     if (!fAOD) {
128       Printf("ERROR: fAOD not available");
129       return;
130     }
131     
132     Printf("Proton analysis task: There are %d tracks in this event", fAOD->GetNumberOfTracks());
133     fAnalysis->Analyze(fAOD);
134   }
135
136   // Post output data.
137   PostData(0, fList);
138 }      
139
140 //________________________________________________________________________
141 void AliAnalysisTaskProtons::Terminate(Option_t *) {
142   // Draw result to the screen
143   // Called once at the end of the query
144   
145   fList = dynamic_cast<TList*> (GetOutputData(0));
146   if (!fList) {
147     Printf("ERROR: fList not available");
148     return;
149   }
150    
151   TH2F *fHistYPtProtons = (TH2F *)fList->At(0);
152   TH2F *fHistYPtAntiProtons = (TH2F *)fList->At(1);
153     
154   TCanvas *c2 = new TCanvas("c2","p-\bar{p}",200,0,800,400);
155   c2->SetFillColor(10); c2->SetHighLightColor(10); c2->Divide(2,1);
156
157   c2->cd(1)->SetLeftMargin(0.15); c2->cd(1)->SetBottomMargin(0.15);  
158   if (fHistYPtProtons) fHistYPtProtons->DrawCopy("colz");
159   c2->cd(2)->SetLeftMargin(0.15); c2->cd(2)->SetBottomMargin(0.15);  
160   if (fHistYPtAntiProtons) fHistYPtAntiProtons->DrawCopy("colz");
161 }
162
163