]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGCF/Correlations/DPhi/AliAnalysisTaskTwoPlusOne.cxx
Merge branch 'master' into TPCdev
[u/mrichter/AliRoot.git] / PWGCF / Correlations / DPhi / AliAnalysisTaskTwoPlusOne.cxx
CommitLineData
c9ae910e 1#include <TROOT.h>
2#include <TChain.h>
3#include <TFile.h>
4#include <TList.h>
5#include <TMath.h>
6#include <TTree.h>
7//#include <TParameter.h>
8
9#include "AliAnalysisTaskTwoPlusOne.h"
10#include "AliCFParticle.h"
11#include "AliAnalyseLeadingTrackUE.h"
12#include "AliTwoPlusOneContainer.h"
13#include "AliUEHist.h"
14
15#include "AliAODHandler.h"
16#include "AliAODInputHandler.h"
17#include "AliVParticle.h"
18#include "AliCFContainer.h"
19
20#include "AliEventPoolManager.h"
21
22
23ClassImp( AliAnalysisTaskTwoPlusOne )
24
25//________________________________________________________________________
26AliAnalysisTaskTwoPlusOne::AliAnalysisTaskTwoPlusOne(const char *name)
27: AliAnalysisTaskSE(name),
28 fMixingTracks(10000),
29 fAnalyseUE(0x0),
30// pointers to UE classes
31 fHistos(0x0),
32// handlers and events
33 fAOD(0x0),
34 fPoolMgr(0x0),
35// histogram settings
36 fListOfHistos(0x0),
37// event QA
38 fnTracksVertex(1), // QA tracks pointing to principal vertex (= 3 default)
39 fZVertex(7.),
40 fCentralityMethod("V0M"),
41// track cuts
42 fTrackEtaCut(0.9),
43 fTrackEtaCutMin(-1.),
44 fPtMin(0.5),
45 fDCAXYCut(0),
46 fSharedClusterCut(-1),
47 fCrossedRowsCut(-1),
48 fFoundFractionCut(-1),
49 fFilterBit(0xFF),
50 fTrackStatus(0),
48c9fd73 51 fThreeParticleMixed(0),
c9ae910e 52 fCustomBinning(),
53 fAlpha(0.2)
54{
55
56 DefineOutput(1, TList::Class());
57}
58
59//________________________________________________________________________
60AliAnalysisTaskTwoPlusOne::~AliAnalysisTaskTwoPlusOne()
61{
62 // Destructor
63
64}
65
66//________________________________________________________________________
67void AliAnalysisTaskTwoPlusOne::UserCreateOutputObjects()
68{
69 // Initialize class with main algorithms, event and track selection.
70 fAnalyseUE = new AliAnalyseLeadingTrackUE();
71 fAnalyseUE->SetParticleSelectionCriteria(fFilterBit, kFALSE, fTrackEtaCut, fTrackEtaCutMin, fPtMin);
72 fAnalyseUE->SetDCAXYCut(fDCAXYCut);
73 fAnalyseUE->SetSharedClusterCut(fSharedClusterCut);
74 fAnalyseUE->SetCrossedRowsCut(fCrossedRowsCut);
75 fAnalyseUE->SetFoundFractionCut(fFoundFractionCut);
76 fAnalyseUE->SetTrackStatus(fTrackStatus);
77 fAnalyseUE->SetDebug(fDebug);
78 fAnalyseUE->DefineESDCuts(fFilterBit);
c9ae910e 79
80 fListOfHistos = new TList();
81 fListOfHistos->SetOwner(kTRUE);
82
83 fHistos = new AliTwoPlusOneContainer("AliTwoPlusOneContainer", fCustomBinning, fAlpha);
84 fHistos->GetData()->SetTrackEtaCut(fTrackEtaCut);
85
86 fListOfHistos->Add(fHistos);
87
88 fListOfHistos->Add(new TH1F("eventStat", ";;events", 4, -0.5, 3.5));
89
90 PostData(1,fListOfHistos);
91
92 // Add task configuration to output list
93 AddSettingsTree();
94
95 // event mixing
96 Int_t poolsize = 1000; // Maximum number of events, ignored in the present implemention of AliEventPoolManager
97
98 Int_t nCentralityBins = fHistos->GetData()->GetTrackHist(AliUEHist::kToward)->GetNBins(3);
99
100
101 Double_t* centralityBins = (Double_t*) fHistos->GetData()->GetTrackHist(AliUEHist::kToward)->GetAxis(3, 0)->GetXbins()->GetArray();
102
103 Int_t nZvtxBins = fHistos->GetData()->GetTrackHist(AliUEHist::kToward)->GetNBins(5);
104
105 Double_t* zvtxbin = (Double_t*) fHistos->GetData()->GetTrackHist(AliUEHist::kToward)->GetAxis(5, 0)->GetXbins()->GetArray();
106
107 fPoolMgr = new AliEventPoolManager(poolsize, fMixingTracks, nCentralityBins, centralityBins, nZvtxBins, zvtxbin);
108 fPoolMgr->SetTargetValues(fMixingTracks, 0.1, 5);
109}
110
111//________________________________________________________________________
112void AliAnalysisTaskTwoPlusOne::UserExec(Option_t *)
113{
114 // exec (per event)
115 fAnalyseUE->NextEvent();
116
117 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(0);
118
119 // Support for AOD based analysis
120 fAOD = dynamic_cast<AliAODEvent*> (InputEvent());
121
122 Double_t centrality = 0;
123 AliCentrality *centralityObj = 0;
124
125 centralityObj = fAOD->GetHeader()->GetCentralityP();
126
127
128 if (centralityObj)
129 centrality = centralityObj->GetCentralityPercentile(fCentralityMethod);
130 else
131 centrality = -1;
132
133 if (fAOD){
134 // remove outliers
135 if (centrality == 0)
136 {
137 if (fAOD->GetVZEROData())
138 {
139 Float_t multV0 = 0;
140 for (Int_t i=0; i<64; i++)
141 multV0 += fAOD->GetVZEROData()->GetMultiplicity(i);
142 if (multV0 < 19500)
143 {
144 centrality = -1;
145 AliInfo("Rejecting event due to too small V0 multiplicity");
146 }
147 }
148 }
149 }
150
151
152 AliInfo(Form("Centrality is %f", centrality));
153
154
155 // Vertex selection *************************************************
156 if(!fAnalyseUE->VertexSelection(InputEvent(), fnTracksVertex, fZVertex)) return;
157
158 // optimization
159 if (centrality < 0)
160 return;
161
162 TObjArray* tracks = fAnalyseUE->GetAcceptedParticles(InputEvent(), 0, kTRUE, -1, kTRUE);
163 // create a list of reduced objects. This speeds up processing and reduces memory consumption for the event pool
164 TObjArray* tracksClone = CloneAndReduceTrackList(tracks);
165 delete tracks;
166
167 const AliVVertex* vertex = InputEvent()->GetPrimaryVertex();
168 Double_t zVtx = vertex->GetZ();
169
170 fHistos->FillCorrelations(centrality, zVtx, AliTwoPlusOneContainer::kSameNS, tracksClone, tracksClone, tracksClone, tracksClone, 1.0);//same event for near and away side
171
172 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(1);
173
174 // event mixing
175
176 // 1. First get an event pool corresponding in mult (cent) and
177 // zvertex to the current event. Once initialized, the pool
178 // should contain nMix (reduced) events. This routine does not
179 // pre-scan the chain. The first several events of every chain
180 // will be skipped until the needed pools are filled to the
181 // specified depth. If the pool categories are not too rare, this
182 // should not be a problem. If they are rare, you could lose
183 // statistics.
184
185 // 2. Collect the whole pool's content of tracks into one TObjArray
186 // (bgTracks), which is effectively a single background super-event.
187
188 // 3. The reduced and bgTracks arrays must both be passed into
189 // FillCorrelations(). Also nMix should be passed in, so a weight
190 // of 1./nMix can be applied.
191
192 AliEventPool* pool = fPoolMgr->GetEventPool(centrality, zVtx);
193
194 if (!pool)
195 AliFatal(Form("No pool found for centrality = %f, zVtx = %f", centrality, zVtx));
196 if (pool->IsReady()){
197 Int_t nMix = pool->GetCurrentNEvents();
198
199 ((TH1F*) fListOfHistos->FindObject("eventStat"))->Fill(2);
200
201 // Fill mixed-event histos here
202 for (Int_t jMix=0; jMix<nMix; jMix++){
203 TObjArray* bgTracks = pool->GetEvent(jMix);
204
48c9fd73 205 //standard mixed event
206 if(!fThreeParticleMixed)
207 fHistos->FillCorrelations(centrality, zVtx, AliTwoPlusOneContainer::kMixedNS, tracksClone, tracksClone, bgTracks, bgTracks, 1.0 / nMix);
208
209 //mixed combinatorics
210 fHistos->FillCorrelations(centrality, zVtx, AliTwoPlusOneContainer::kMixedCombNS, tracksClone, bgTracks, tracksClone, bgTracks, 1.0 / nMix);
c9ae910e 211 }
48c9fd73 212
213
214 // use 3 particle mixed event
215 if(fThreeParticleMixed && nMix>1){
216 TObjArray* tracks_t2 = pool->GetEvent(0);
217 for (Int_t jMix=1; jMix<nMix; jMix++){
218
219 TObjArray* bgTracks = pool->GetEvent(jMix);
220
221 fHistos->FillCorrelations(centrality, zVtx, AliTwoPlusOneContainer::kMixedNS, tracksClone, tracks_t2, bgTracks, bgTracks, 1.0 / (nMix-1));
222 }
223 }
224
225
c9ae910e 226 }
227
228 // ownership is with the pool now
229 pool->UpdatePool(tracksClone);
230}
231
232//________________________________________________________________________
233void AliAnalysisTaskTwoPlusOne::Terminate(Option_t*)
234{
235 //terminate function is called at the end
236 //can be used to draw histograms etc.
237
238
239}
240
241
242TObjArray* AliAnalysisTaskTwoPlusOne::CloneAndReduceTrackList(TObjArray* tracks)
243{
244 // clones a track list by using AliCFParticle which uses much less memory (used for event mixing)
245
246 TObjArray* tracksClone = new TObjArray;
247 tracksClone->SetOwner(kTRUE);
248
249 for (Int_t i=0; i<tracks->GetEntriesFast(); i++)
250 {
251 AliVParticle* particle = (AliVParticle*) tracks->UncheckedAt(i);
252 AliCFParticle* copy = new AliCFParticle(particle->Pt(), particle->Eta(), particle->Phi(), particle->Charge(), 0);
253 copy->SetUniqueID(particle->GetUniqueID());
254 tracksClone->Add(copy);
255 }
256
257 return tracksClone;
258}
259
260
261//____________________________________________________________________
262void AliAnalysisTaskTwoPlusOne::AddSettingsTree()
263{
264 //Write settings to output list
265 TTree *settingsTree = new TTree("UEAnalysisSettings","Analysis Settings in UE estimation");
266 settingsTree->Branch("fnTracksVertex", &fnTracksVertex,"nTracksVertex/I");
267 settingsTree->Branch("fZVertex", &fZVertex,"ZVertex/D");
268 //settingsTree->Branch("fCentralityMethod", fCentralityMethod.Data(),"CentralityMethod/C");
269 settingsTree->Branch("fTrackEtaCut", &fTrackEtaCut, "TrackEtaCut/D");
270 settingsTree->Branch("fTrackEtaCutMin", &fTrackEtaCutMin, "TrackEtaCutMin/D");
271 settingsTree->Branch("fPtMin", &fPtMin, "PtMin/D");
272 settingsTree->Branch("fFilterBit", &fFilterBit,"FilterBit/I");
273 settingsTree->Branch("fSharedClusterCut", &fSharedClusterCut,"SharedClusterCut/D");
274 settingsTree->Branch("fCrossedRowsCut", &fCrossedRowsCut,"CrossedRowsCut/I");
275 settingsTree->Branch("fFoundFractionCut", &fFoundFractionCut,"FoundFractionCut/D");
276 settingsTree->Branch("fTrackStatus", &fTrackStatus,"TrackStatus/I");
48c9fd73 277 settingsTree->Branch("fThreeParticleMixed", &fThreeParticleMixed,"fThreeParticleMixed/I");
c9ae910e 278 settingsTree->Branch("fMixingTracks", &fMixingTracks,"MixingTracks/I");
279
280 settingsTree->Fill();
281 fListOfHistos->Add(settingsTree);
282}