]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/FLOW/Tasks/AliAnalysisTaskJetFlow.cxx
55ef2187402cce6cb83ca7666ed4b611a687553c
[u/mrichter/AliRoot.git] / PWG / FLOW / Tasks / AliAnalysisTaskJetFlow.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* 
17  * AliAnaysisTaskJetFlow
18  * author: Redmer Alexander Bertens
19  * rbertens@cern.ch, rbertens@nikhef.nl, r.a.bertens@uu.nl 
20  *
21  * Interface task between EMCal jet framework and flow package
22  *
23  * This task expects POI's in a TClonesArray  (e.g. from running it after 
24  * $ALICE_ROOT/PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskRhoVnModulation() )
25  * and connects them into an AliFlowEvent which is filled with either VZERO tracks or 
26  * TPC trakcs. 
27  * 
28  * POI's can be supplied as AliEmcalJets or as AliVTracks, note the different behavior
29  * with respect to the Pt value: for AliEmcalJets subtracted Pt is used by default, for VParticles
30  * Pt is taken directly.
31  *
32  * */
33
34 // root includes
35 #include <TChain.h>
36 #include <TH1F.h>
37 #include <TArrayD.h>
38 #include <TProfile.h>
39 #include <TString.h>
40 #include <TList.h>
41 #include <TClonesArray.h>
42 // aliroot includes
43 #include <AliAODEvent.h>
44 #include <AliESDEvent.h>
45 #include <AliEmcalJet.h>
46 #include <AliAnalysisTask.h>
47 #include <AliAnalysisManager.h>
48 #include <AliFlowEventSimple.h>
49 #include <AliFlowEvent.h>
50 #include <AliFlowTrack.h>
51 #include <AliFlowTrackCuts.h>
52 #include <AliFlowEventCuts.h>
53 #include <AliFlowCommonConstants.h>
54 // local includes
55 #include "AliAnalysisTaskJetFlow.h"
56 // EMCAL jet framework includes
57 #include <AliRhoParameter.h>
58 #include <AliPicoTrack.h>
59 #include <AliAnalysisTaskRhoVnModulation.h>
60
61 class AliAnalysisTaskJetFlow;
62
63 using namespace std;
64
65 ClassImp(AliAnalysisTaskJetFlow)
66
67 AliAnalysisTaskJetFlow::AliAnalysisTaskJetFlow() : AliAnalysisTaskSE(), 
68     fDebug(-1), fJetsName(0), fTracksName(0), fPois(0x0), fOutputList(0), fDataType(kESD), fVParticleAnalysis(kFALSE), fMinimizeDiffBins(kTRUE), fDoVZEROFlowAnalysis(kTRUE), fDoQC2FlowAnalysis(kTRUE), fDoQC4FlowAnalysis(kFALSE), fDoQCFPAnalysis(kFALSE), fDoSPFPAnalysis(kFALSE), fDoMultWeight(kTRUE), fDoPtWeight(0), fInitialized(kFALSE), fUsePtWeight(kFALSE), fCCMinPt(1), fCCMaxPt(150), fCCBinsInPt(50), fCentralityMin(-1), fCentralityMax(-1), fPtBins(0), fCutsRP_VZERO(0), fCutsNull(0), fCutsEvent(0), fFlowEvent_TPC(0), fFlowEvent_VZERO(0), fRhoVn(0), fHistAnalysisSummary(0), fCentralitySelection(0), fVZEROAEP(0), fVZEROCEP(0), fv2VZEROA(0), fv2VZEROC(0), fRefCumulants(0), fDiffCumlantsV2(0), fDiffCumlantsV3(0), fQC2v2(0), fQC2v3(0), fTempA(0), fTempC(0)
69 { /* default constructor for ROOT IO */ }
70 //_____________________________________________________________________________
71 AliAnalysisTaskJetFlow::AliAnalysisTaskJetFlow(
72         const char* name,
73         AliAnalysisTaskRhoVnModulation* rhoTask, 
74         Bool_t VPart,
75         Bool_t VZEROEP, 
76         Bool_t QC2,
77         Bool_t QC4,
78         Bool_t FlowPackageSP,
79         Bool_t FlowPackageQC  
80         ) : AliAnalysisTaskSE(name),
81     fDebug(-1), fJetsName(0), fTracksName(0), fPois(0x0), fOutputList(0), fDataType(kESD), fVParticleAnalysis(VPart), fMinimizeDiffBins(kTRUE), fDoVZEROFlowAnalysis(VZEROEP), fDoQC2FlowAnalysis(QC2), fDoQC4FlowAnalysis(QC4), fDoQCFPAnalysis(FlowPackageQC), fDoSPFPAnalysis(FlowPackageSP), fDoMultWeight(kTRUE), fDoPtWeight(0), fInitialized(kFALSE), fUsePtWeight(kFALSE), fCCMinPt(1), fCCMaxPt(150), fCCBinsInPt(50), fCentralityMin(-1), fCentralityMax(-1), fPtBins(0), fCutsRP_VZERO(0x0), fCutsNull(0), fCutsEvent(0), fFlowEvent_TPC(0), fFlowEvent_VZERO(0), fRhoVn(rhoTask), fHistAnalysisSummary(0), fCentralitySelection(0), fVZEROAEP(0), fVZEROCEP(0), fv2VZEROA(0), fv2VZEROC(0), fRefCumulants(0), fDiffCumlantsV2(0), fDiffCumlantsV3(0), fQC2v2(0), fQC2v3(0), fTempA(0), fTempC(0)
82 {
83     // constructor
84     DefineInput(0, TChain::Class());
85     DefineOutput(1, TList::Class());
86     fJetsName = rhoTask->GetJetsName();
87     fTracksName = rhoTask->GetTracksName(); 
88     if(FlowPackageSP || FlowPackageQC)    DefineOutput(2, AliFlowEventSimple::Class());
89     if(FlowPackageSP && FlowPackageQC)    DefineOutput(3, AliFlowEventSimple::Class());
90 }
91 //_____________________________________________________________________________
92 AliAnalysisTaskJetFlow::~AliAnalysisTaskJetFlow()
93 {
94     // destructor
95     if(fOutputList)             delete fOutputList;
96     if(fFlowEvent_TPC)          delete fFlowEvent_TPC;
97     if(fFlowEvent_VZERO)        delete fFlowEvent_VZERO;
98     if(fCutsEvent)              delete fCutsEvent;
99     if(fCutsRP_VZERO)           delete fCutsRP_VZERO;
100     if(fCutsNull)               delete fCutsNull;
101     if(fPtBins)                 delete fPtBins;
102     if(fTempA)                  delete fTempA;
103     if(fTempC)                  delete fTempC;
104 }
105 //_____________________________________________________________________________
106 void AliAnalysisTaskJetFlow::UserCreateOutputObjects()
107 {
108     // create output objects
109     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
110     // create the cut objects
111     if(fDoSPFPAnalysis || fDoQCFPAnalysis) {
112         fCutsEvent = new AliFlowEventCuts();
113         fCutsEvent->SetRefMultMethod(AliESDtrackCuts::kTrackletsITSTPC); 
114         fCutsNull = new AliFlowTrackCuts("CutsNull");
115         fCutsNull->SetParamType(AliFlowTrackCuts::kGlobal);
116         fCutsNull->SetEtaRange(+1, -1);
117         fCutsNull->SetPtRange(+1, -1);
118         fCutsRP_VZERO = new AliFlowTrackCuts();
119         fCutsRP_VZERO = fCutsRP_VZERO->GetStandardVZEROOnlyTrackCuts();
120         if(fDoSPFPAnalysis) {
121             (fVParticleAnalysis) ? fFlowEvent_VZERO = new AliFlowEvent(10000) : fFlowEvent_VZERO = new AliFlowEvent(100);
122         }
123         fFlowEvent_TPC = new AliFlowEvent(10000);
124     }
125     fOutputList = new TList();
126     fOutputList->SetOwner(kTRUE);
127     // histograms
128     fHistAnalysisSummary = new TH1F("fHistAnalysisSummary", "fHistAnalysisSummary", 4, -.5, 3.5);
129     fHistAnalysisSummary->GetXaxis()->SetBinLabel(1, "fDataType");
130     fHistAnalysisSummary->SetBinContent(1, (int)fDataType);
131     fHistAnalysisSummary->GetXaxis()->SetBinLabel(2, "fCentralityMin");
132     fHistAnalysisSummary->SetBinContent(2, fCentralityMin);
133     fHistAnalysisSummary->GetXaxis()->SetBinLabel(3, "fCentralityMax");
134     fHistAnalysisSummary->SetBinContent(3, fCentralityMax);
135     fOutputList->Add(fHistAnalysisSummary);
136     if(fVParticleAnalysis && ! fPtBins) { // FIXME inherits from flow package now
137         Double_t pt[fCCBinsInPt+1];
138         Double_t width = (fCCMaxPt - fCCMinPt ) / (float)fCCBinsInPt;
139         for(Int_t i(0); i < fCCBinsInPt+1; i++) pt[i] = fCCMinPt+width*i;
140         fPtBins = new TArrayD(fCCBinsInPt+1, pt);  // assume they will be charged particles
141     } else if (!fPtBins) {
142         Double_t pt[] = {0., 10., 20., 30., 40., 50., 80., 110., 140., 170., 200.};
143         fPtBins = new TArrayD(sizeof(pt)/sizeof(pt[0]), pt);     // assuming jets
144     }
145     if(fDoVZEROFlowAnalysis) {
146         fv2VZEROA = new TProfile("Differential v_{2}^{obs} VZEROA", "Differential v_{2}^{obs} VZEROA", fPtBins->GetSize()-1, fPtBins->GetArray());
147         fv2VZEROC = new TProfile("Differential v_{2}^{obs} VZEROC", "Differential v_{2}^{obs} VZEROC", fPtBins->GetSize()-1, fPtBins->GetArray());
148         fv2VZEROA->GetXaxis()->SetTitle("Pt [GeV/c]");
149         fv2VZEROA->GetYaxis()->SetTitle("v_{2}^{obs}");
150         fv2VZEROC->GetXaxis()->SetTitle("Pt [GeV/c]");
151         fv2VZEROC->GetYaxis()->SetTitle("v_{2}^{obs}");
152         fOutputList->Add(fv2VZEROA);
153         fOutputList->Add(fv2VZEROC);
154         fVZEROAEP = new TH1F("V0A_EP", "V0A_EP", 100, -TMath::Pi()/2., TMath::Pi()/2.);
155         fVZEROCEP = new TH1F("V0C_EP", "V0C_EP", 100, -TMath::Pi()/2., TMath::Pi()/2.);
156         fOutputList->Add(fVZEROAEP);
157         fOutputList->Add(fVZEROCEP);
158         // bookkeeping histo's, will not be stored
159         fTempA = (TProfile*)fv2VZEROA->Clone("temp_a");
160         fTempC = (TProfile*)fv2VZEROC->Clone("temp_c");
161     }
162     if(fDoQC2FlowAnalysis) {
163         fRefCumulants = new TProfile("Reference cumulants", "Reference cumulants", 2, -0.5, 1.5);
164         fRefCumulants->GetXaxis()->SetBinLabel(1, "c_{2}[2]");
165         fRefCumulants->GetXaxis()->SetBinLabel(2, "c_{3}[2]");
166         fOutputList->Add(fRefCumulants);
167         fDiffCumlantsV2 = new TProfile("Differential cumulants v2", "Differential cumulants v2", fPtBins->GetSize()-1, fPtBins->GetArray());
168         fDiffCumlantsV2->GetXaxis()->SetTitle("Pt [GeV/c]");
169         fDiffCumlantsV2->GetYaxis()->SetTitle("v_{2}[2]");
170         fOutputList->Add(fDiffCumlantsV2);
171         fDiffCumlantsV3 = new TProfile("Differential cumulants v3", "Differential cumulants v3", fPtBins->GetSize()-1, fPtBins->GetArray());
172         fDiffCumlantsV3->GetXaxis()->SetTitle("Pt [GeV/c]");
173         fDiffCumlantsV3->GetYaxis()->SetTitle("v_{3}[2]");
174         fOutputList->Add(fDiffCumlantsV3);
175         fQC2v2 = new TH1F("Differential v_{2}[2]", "Differential v_{2}[2]", fPtBins->GetSize()-1, fPtBins->GetArray());
176         fQC2v2->Sumw2();
177         fQC2v3 = new TH1F("Differential v_{3}[2]", "Differential v_{3}[2]", fPtBins->GetSize()-1, fPtBins->GetArray());
178         fQC2v3->Sumw2();
179         fOutputList->Add(fQC2v2); 
180         fOutputList->Add(fQC2v3);
181     }
182     // qa
183     fCentralitySelection = new TH1F("fCentralitySelection", "fCentralitySelection", 100, 0, 100);
184     fOutputList->Add(fCentralitySelection);
185     PostData(1, fOutputList);
186     // create the flow event and configure the static cc object
187     Bool_t slotTwoFilled(kFALSE);
188     if(fFlowEvent_VZERO) {
189         PostData(2, fFlowEvent_VZERO);
190         slotTwoFilled = kTRUE;
191     }
192     if(fFlowEvent_TPC) {
193         (slotTwoFilled) ? PostData(3, fFlowEvent_TPC) : PostData(2, fFlowEvent_TPC);
194     }
195     if(fFlowEvent_VZERO || fFlowEvent_TPC) {
196         AliFlowCommonConstants* cc = AliFlowCommonConstants::GetMaster();
197         cc->SetPtMin(fCCMinPt);
198         cc->SetPtMax(fCCMaxPt);
199         cc->SetNbinsPt(fCCBinsInPt);
200         if(fMinimizeDiffBins) { // minimize differential bins to reduce the risk of numerical instability
201             cc->SetNbinsMult(1);    // only reduces output size
202             cc->SetNbinsPhi(1);     // only reduces output size
203             cc->SetNbinsEta(1);     // reduces instability
204             cc->SetNbinsQ(1);       // only reduces output size
205             cc->SetNbinsMass(1);    // reduces instability but should be one in any case ...
206         }
207     }
208 }
209 //_____________________________________________________________________________
210 void AliAnalysisTaskJetFlow::UserExec(Option_t *)
211 {
212     // user exec
213     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
214     if(!InputEvent()) return; // coverity (and sanity)
215     // see if the analysis is initialized 
216     if(!fInitialized) { 
217         if(dynamic_cast<AliAODEvent*>(InputEvent())) fDataType = kAOD; // determine the datatype
218         else if(dynamic_cast<AliESDEvent*>(InputEvent())) fDataType = kESD;
219         (fVParticleAnalysis) ? fPois = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fTracksName.Data())) : fPois = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fJetsName.Data()));
220         if(!fPois) return; // couldn't get expected input data
221         fInitialized = kTRUE;
222     }
223     if(!PassesCuts()) return; // event quality cuts and centrality determination
224     // execute the requested flow methods
225     if(fDoVZEROFlowAnalysis)            DoVZEROFlowAnalysis();
226     if(fDoQC2FlowAnalysis)              DoQC2FlowAnalysis();
227     if(fDoQC4FlowAnalysis)              DoQC4FlowAnalysis();
228     Bool_t post(0x0);   // post only when analysis succeeded
229     if(fFlowEvent_TPC || fFlowEvent_VZERO) DoFlowPackageFlowAnalysis();
230     // push the output for the different analyses to file
231     PostData(1, fOutputList);
232     if(fFlowEvent_VZERO) {
233         PostData(2, fFlowEvent_VZERO);
234         post = kTRUE;
235     }
236     if(fFlowEvent_TPC) { 
237         (post) ? PostData(3, fFlowEvent_TPC) : PostData(2, fFlowEvent_TPC);
238     }
239 }
240 //_____________________________________________________________________________
241 Bool_t AliAnalysisTaskJetFlow::PassesCuts() 
242 {
243     // passes event cuts
244     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
245     if(!fRhoVn->PassesCuts(InputEvent())) return kFALSE;          // event quality cuts
246     if(InputEvent()->GetCentrality()->GetCentralityPercentile("V0M") < fCentralityMin || InputEvent()->GetCentrality()->GetCentralityPercentile("V0M") > fCentralityMax) return kFALSE;  // cutting on centrality quality is done in event quality cuts
247     fCentralitySelection->Fill(InputEvent()->GetCentrality()->GetCentralityPercentile("V0M"));
248     return kTRUE;
249 }
250 //_____________________________________________________________________________
251 void AliAnalysisTaskJetFlow::DoVZEROFlowAnalysis()
252 {
253     // flow with the VZERO event plane method
254     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
255     Double_t vzero[2][2];
256     fRhoVn->CalculateEventPlaneVZERO(vzero);
257     Double_t Q2a(vzero[0][0]), Q2c(vzero[1][0]);        // just for readability
258     fTempA->Reset();            fTempC->Reset();        // clear the containers for a new iteration
259     fVZEROAEP->Fill(Q2a);       fVZEROCEP->Fill(Q2c);
260     Int_t iPois(fPois->GetEntriesFast());
261     if(fVParticleAnalysis) {
262         for(Int_t i(0); i < iPois; i++) {
263             AliVTrack* poi = static_cast<AliVTrack*>(fPois->At(i));
264             if(fRhoVn->PassesCuts(poi)) {
265                 if(!fDoMultWeight) {
266                     fTempA->Fill(poi->Pt(), TMath::Cos(2.*fRhoVn->PhaseShift((poi->Phi()-Q2a), 2)));
267                     fTempC->Fill(poi->Pt(), TMath::Cos(2.*fRhoVn->PhaseShift((poi->Phi()-Q2c), 2)));
268                 } else {
269                     fv2VZEROA->Fill(poi->Pt(), TMath::Cos(2.*fRhoVn->PhaseShift((poi->Phi()-Q2a), 2)));
270                     fv2VZEROC->Fill(poi->Pt(), TMath::Cos(2.*fRhoVn->PhaseShift((poi->Phi()-Q2c), 2)));
271                 }
272             }
273         }
274     } else {
275         for(Int_t i(0); i < iPois; i++) {
276             AliEmcalJet* poi = static_cast<AliEmcalJet*>(fPois->At(i));
277             if(poi && poi->PtSub()) {
278                 if(!fDoMultWeight) {
279                     fTempA->Fill(poi->PtSub(), TMath::Cos(2.*fRhoVn->PhaseShift((poi->Phi()-Q2a), 2)));
280                     fTempC->Fill(poi->PtSub(), TMath::Cos(2.*fRhoVn->PhaseShift((poi->Phi()-Q2c), 2)));
281                 } else {
282                     fv2VZEROA->Fill(poi->PtSub(), TMath::Cos(2.*fRhoVn->PhaseShift((poi->Phi()-Q2a), 2)));
283                     fv2VZEROC->Fill(poi->PtSub(), TMath::Cos(2.*fRhoVn->PhaseShift((poi->Phi()-Q2c), 2)));
284                 }
285             }
286         }
287     }
288     if(!fDoMultWeight) {
289         for(Int_t i(0); i < fv2VZEROA->GetXaxis()->GetNbins(); i++) {
290             fv2VZEROA->Fill(fPtBins->At(i)+(fPtBins->At(i)+fPtBins->At(1+i))/2., fTempA->GetBinContent(i+1));
291             fv2VZEROC->Fill(fPtBins->At(i)+(fPtBins->At(i)+fPtBins->At(1+i))/2., fTempC->GetBinContent(i+1));
292         }
293     }
294 }
295 //_____________________________________________________________________________
296 void AliAnalysisTaskJetFlow::DoQC2FlowAnalysis()
297 {
298     // flow analysis with the qc2 method
299     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
300     // reference flow is taken from the pico track selection and evaluated in 
301     // the AliAnalysisTaskiRhoVnModulation class
302     for(Int_t harm(2); harm < 4; harm++) {        // loop over harmonics
303         Double_t reQn(0), imQn(0), mQ(0);         // total q-vector
304         // get the total q-vectors and reference cumulants
305         (fDoPtWeight) ? fRhoVn->QCnQnk(harm, 1, reQn, imQn) : fRhoVn->QCnQnk(harm, 0, reQn, imQn);
306         (fDoPtWeight) ? mQ = fRhoVn->QCnM11() : mQ = fRhoVn->QCnM();
307         if(mQ < 2) continue;    // avoid division by zero
308         fRefCumulants->Fill((double)(harm-2), ((reQn*reQn+imQn*imQn)-mQ)/(mQ*(mQ-1)), mQ*(mQ-1));
309         // differential flow vectors
310         Double_t repn[fPtBins->GetSize()-1];      // real part of q-vector of all poi's
311         Double_t impn[fPtBins->GetSize()-1];      // im part of q-vector of all poi's
312         Double_t mp[fPtBins->GetSize()-1];        // poi multiplicity
313         Double_t reqn[fPtBins->GetSize()-1];      // real part of q-vectors of poi's labeled as rp
314         Double_t imqn[fPtBins->GetSize()-1];      // im part of q-vectors of poi's labeled as rp
315         Double_t mq[fPtBins->GetSize()-1];        // multiplicity of poi's labeled as rp
316         for(Int_t i(0); i < fPtBins->GetSize(); i++) {
317             repn[i] = 0;
318             impn[i] = 0;
319             mp[i] = 0;
320             reqn[i] = 0;
321             imqn[i] = 0;
322             mq[i] = 0;
323         }
324         // calculate differential q-vectors and fill the profile with cumulants
325         QCnDiffentialFlowVectors(repn, impn, mp, reqn, imqn, mq, harm);
326         // FIXME differential evnet weights
327         for(Int_t i(0); i < fPtBins->GetSize(); i++) {
328             if(mp[i]*mQ - mq[i] <= 0 ) continue;        // avoid division by zero
329             Double_t atPt(fPtBins->At(i)+0.5*(fPtBins->At(i+1)-fPtBins->At(i)));      // pt value
330             Double_t diffC(((repn[i]*reQn+impn[i]*imQn)-mq[i])/(mp[i]*mQ-mq[i]));
331             Double_t eventW(mp[i]*mQ-mq[i]);
332             (harm == 2 ) ? fDiffCumlantsV2->Fill(atPt, diffC /*weight*/) : fDiffCumlantsV3->Fill(atPt, diffC, eventW);
333         }
334     } 
335 }
336 //_____________________________________________________________________________
337 void AliAnalysisTaskJetFlow::DoQC4FlowAnalysis()
338 {
339     // flow analysis with the qc4 method - see comments at qc2 FIXME not implemented yet
340     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
341 }
342 //_____________________________________________________________________________
343 Bool_t AliAnalysisTaskJetFlow::DoFlowPackageFlowAnalysis() 
344 {
345     // name's a bit misleading: this function does anlaysis using FlowAnalysis classes in /PWG/FLOW/Base/
346     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
347     // get the jet array, which is added as an extension to the AliVEvent by the jetfinder
348     Int_t nAcceptedJets(0), iPois(fPois->GetEntriesFast());
349     if(iPois <= 0) return kFALSE;
350     if(fFlowEvent_VZERO) {
351         fCutsNull->SetEvent(InputEvent(), MCEvent());
352         fCutsRP_VZERO->SetEvent(InputEvent(), MCEvent());
353         fFlowEvent_VZERO->ClearFast();
354         // the event is filled with rp's only, poi's will be added manually
355         fFlowEvent_VZERO->Fill(fCutsRP_VZERO, fCutsNull);
356         fFlowEvent_VZERO->SetReferenceMultiplicity(fCutsEvent->RefMult(InputEvent(), MCEvent()));
357     }
358     if(fFlowEvent_TPC) {
359         fCutsNull->SetEvent(InputEvent(), MCEvent());
360         // in this case, both poi's and rp's will be added manually
361         fFlowEvent_TPC->ClearFast();
362         fFlowEvent_TPC->SetReferenceMultiplicity(fCutsEvent->RefMult(InputEvent(), MCEvent()));
363     }
364     // loop over jets and inject them as POI's
365     if(fVParticleAnalysis) {
366         for(Int_t i(0); i < iPois; i++) {
367             AliVTrack* poi = static_cast<AliVTrack*>(fPois->At(i));
368             if(fRhoVn->PassesCuts(poi)) {
369                 nAcceptedJets++;
370                 // AliFlowTracks are created on the stack 
371                 AliFlowTrack flowTrack = AliFlowTrack(poi);
372                 flowTrack.SetForPOISelection(kTRUE);
373                 if(fFlowEvent_TPC) {
374                     fFlowEvent_TPC->SetNumberOfRPs(fFlowEvent_TPC->GetNumberOfRPs()+1);
375                     fFlowEvent_TPC->SetNumberOfPOIs(fFlowEvent_TPC->GetNumberOfPOIs()+1);
376                     flowTrack.SetForRPSelection(kTRUE);
377                     fFlowEvent_TPC->InsertTrack(&flowTrack);
378                 }
379                 if(fFlowEvent_VZERO) {
380                     flowTrack.SetForRPSelection(kFALSE);
381                     fFlowEvent_VZERO->InsertTrack(&flowTrack);
382                 }
383             }
384         }
385     } else {
386         for(Int_t i(0); i < iPois; i++) {
387             AliEmcalJet* poi = static_cast<AliEmcalJet*>(fPois->At(i));
388             if(poi) {
389                 if(poi->PtSub()) {
390                     fHistAnalysisSummary->SetBinContent(4, 1);
391                     continue;
392                 }
393                 nAcceptedJets++;
394                 AliFlowTrack flowTrack = AliFlowTrack(poi);
395                 flowTrack.SetPt(poi->PtSub());
396                 flowTrack.SetForPOISelection(kTRUE);
397                 flowTrack.SetForRPSelection(kFALSE);
398                 if(fFlowEvent_TPC)      fFlowEvent_TPC->InsertTrack(&flowTrack);
399                 if(fFlowEvent_VZERO)    fFlowEvent_VZERO->InsertTrack(&flowTrack);
400             }
401         }
402     }
403     if(fFlowEvent_TPC)      fFlowEvent_TPC->TagSubeventsInEta(-10, 0, 0, 10);
404     if(fFlowEvent_VZERO)    fFlowEvent_VZERO->TagSubeventsInEta(-10, 0, 0, 10);
405     return (nAcceptedJets < 1) ? kFALSE : kTRUE;
406 }
407 //_____________________________________________________________________________
408 void AliAnalysisTaskJetFlow::QCnDiffentialFlowVectors(Double_t* repn, Double_t* impn, Double_t *mp, Double_t *reqn, Double_t *imqn, Double_t* mq, Int_t n) 
409 {
410     // get (for now) unweighted differential flow vectors
411     // FIXME move (part of) this code to AliAnalysisTaskRhoVnModulation
412     Int_t iPois(fPois->GetEntriesFast());
413     if(fVParticleAnalysis) {
414         for(Int_t i(0); i < iPois; i++) {
415             for(Int_t ptBin(0); ptBin < fPtBins->GetSize()-1; ptBin++) {
416                 AliVTrack* poi = static_cast<AliVTrack*>(fPois->At(i));
417                 if(fRhoVn->PassesCuts(poi)) {       // inherit cuts from mother task
418                     if(poi->Pt() >= fPtBins->At(ptBin) && poi->Pt() < fPtBins->At(ptBin+1)) {
419                             // fill the flow vectors assuming that all poi's are in the rp selection (true by design)
420                             repn[ptBin]+=TMath::Cos(((Double_t)n)*poi->Phi());
421                             impn[ptBin]+=TMath::Sin(((Double_t)n)*poi->Phi());
422                             mp[ptBin]++;
423                             reqn[ptBin]+=TMath::Cos(((Double_t)n)*poi->Phi());
424                             imqn[ptBin]+=TMath::Sin(((Double_t)n)*poi->Phi());
425                             mq[ptBin]++;
426                     }
427                 }
428             }
429         }
430     } else {
431         for(Int_t i(0); i < iPois; i++) {
432             for(Int_t ptBin(0); ptBin < fPtBins->GetSize()-1; ptBin++) {
433                 AliEmcalJet* poi = static_cast<AliEmcalJet*>(fPois->At(i));
434                 if(poi && poi->PtSub() > 0) {   // note here that no cuts are needed since only accepted jets have PtSub set !    
435                     if(poi->PtSub() >= fPtBins->At(ptBin) && poi->PtSub() < fPtBins->At(ptBin+1)) {    
436                             // fill the flow vectors assuming that all poi's are in the rp selection (true by design)  
437                             repn[ptBin]+=TMath::Cos(((Double_t)n)*poi->Phi());
438                             impn[ptBin]+=TMath::Sin(((Double_t)n)*poi->Phi());
439                             mp[ptBin]++;        // qn isn't filled, no overlap between poi's and rp's
440                     }
441                 }
442             }
443         }
444     }
445 }
446 //_____________________________________________________________________________
447 void AliAnalysisTaskJetFlow::Terminate(Option_t *)
448 {
449     // terminate
450     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
451 }
452 //_____________________________________________________________________________