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