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