1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 * author: Redmer Alexander Bertens
19 * rbertens@cern.ch, rbertens@nikhef.nl, r.a.bertens@uu.nl
21 * jet correlation task - correlates jets to the vzero reaction plane using
22 * the scalar product method
24 * this task should be run AFTER
25 * $ALICE_ROOT/PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskRhoVnModulation()
31 #include <AliAnalysisTask.h>
32 #include <AliAnalysisManager.h>
33 #include <AliFlowEventSimple.h>
34 #include <AliFlowTrack.h>
35 #include <AliFlowTrackCuts.h>
36 #include <AliFlowCommonConstants.h>
38 #include <AliAODEvent.h>
39 #include <AliESDEvent.h>
41 #include "AliAnalysisTaskJetFlow.h"
43 class AliAnalysisTaskJetFlow;
47 ClassImp(AliAnalysisTaskJetFlow)
49 AliAnalysisTaskJetFlow::AliAnalysisTaskJetFlow() : AliAnalysisTaskSE(),
50 fDebug(-1), fJetsName(0), fOutputList(0), fDataType(kESD), fPtBump(0), fCentralityMin(-1), fCentralityMax(-1), fCutsRP(0), fCutsPOI(0), fCutsNull(0), fFlowEvent(0), fHistAnalysisSummary(0)
51 { /* default constructor */ }
52 //_____________________________________________________________________________
53 AliAnalysisTaskJetFlow::AliAnalysisTaskJetFlow(const char* name) : AliAnalysisTaskSE(name),
54 fDebug(-1), fJetsName(0), fOutputList(0), fDataType(kESD), fPtBump(0), fCentralityMin(-1), fCentralityMax(-1), fCutsRP(0), fCutsPOI(0), fCutsNull(0), fFlowEvent(0), fHistAnalysisSummary(0)
57 DefineInput(0, TChain::Class());
58 DefineOutput(1, TList::Class());
59 DefineOutput(2, AliFlowEventSimple::Class());
61 //_____________________________________________________________________________
62 AliAnalysisTaskJetFlow::~AliAnalysisTaskJetFlow()
65 if(fOutputList) delete fOutputList;
66 if(fFlowEvent) delete fFlowEvent;
68 //_____________________________________________________________________________
69 void AliAnalysisTaskJetFlow::LocalInit()
72 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
73 if(dynamic_cast<AliAODEvent*>(InputEvent())) fDataType = kAOD; // determine the datatype
74 else if(dynamic_cast<AliESDEvent*>(InputEvent())) fDataType = kESD;
76 //_____________________________________________________________________________
77 void AliAnalysisTaskJetFlow::UserCreateOutputObjects()
79 // create output objects
80 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
81 fOutputList = new TList();
82 fOutputList->SetOwner(kTRUE);
84 fHistAnalysisSummary = new TH1F("fHistAnalysisSummary", "fHistAnalysisSummary", 4, -.5, 3.5);
85 fHistAnalysisSummary->GetXaxis()->SetBinLabel(1, "fDataType");
86 fHistAnalysisSummary->SetBinContent(1, (int)fDataType);
87 fHistAnalysisSummary->GetXaxis()->SetBinLabel(2, "fCentralityMin");
88 fHistAnalysisSummary->SetBinContent(2, fCentralityMin);
89 fHistAnalysisSummary->GetXaxis()->SetBinLabel(3, "fCentralityMax");
90 fHistAnalysisSummary->SetBinContent(3, fCentralityMax);
91 fHistAnalysisSummary->GetXaxis()->SetBinLabel(4, "pt bias");
92 fOutputList->Add(fHistAnalysisSummary);
93 PostData(1, fOutputList);
94 // create the flow event and configure the static cc object
95 fFlowEvent = new AliFlowEvent(1000);
96 PostData(2, fFlowEvent);
97 AliFlowCommonConstants* cc = AliFlowCommonConstants::GetMaster();
98 cc->SetPtMax(100+fPtBump);
101 //_____________________________________________________________________________
102 void AliAnalysisTaskJetFlow::UserExec(Option_t *)
105 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
107 if(!InputEvent() || !fCutsNull || !fCutsRP) return; // coverity (and sanity)
108 // get the jet array, which is added as an extension to the AliVEvent by the jetfinder
109 TClonesArray* jets = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fJetsName.Data()));
111 Int_t iJets = jets->GetEntriesFast();
113 if(fDebug>0) printf(" > Retrieved empty AliVEvent extension, aborting ! < \n ");
116 // prepare the track selection for RP's and the flow event
117 fCutsRP->SetEvent(InputEvent(), MCEvent());
118 fCutsNull->SetEvent(InputEvent(), MCEvent());
119 fFlowEvent->ClearFast();
120 fFlowEvent->Fill(fCutsRP, fCutsNull);
121 fFlowEvent->SetReferenceMultiplicity(64); // hardcoded for vzero
122 fFlowEvent->DefineDeadZone(0, 0, 0, 0);
123 // loop over jets and inject them as POI's
124 for(Int_t i(0); i < iJets; i++) {
125 AliVParticle* jet = static_cast<AliVParticle*>(jets->At(i));
127 if(jet->Pt() + fPtBump <= 0) {
128 fHistAnalysisSummary->SetBinContent(4, 1);
131 AliFlowTrack* flowTrack = new AliFlowTrack(jet);
132 flowTrack->SetPt(flowTrack->Pt() + fPtBump);
133 flowTrack->SetForPOISelection(kTRUE);
134 flowTrack->SetForRPSelection(kFALSE);
135 fFlowEvent->InsertTrack(flowTrack);
139 else if(fDebug > 0 ) printf(" Failed to find TClones Jet array while using name %s \n ", fJetsName.Data());
140 PostData(1, fOutputList);
141 PostData(2, fFlowEvent);
143 printf(" Event summary \n");
144 printf(" > number of POI's (jets) %i ", fFlowEvent->NumberOfTracks() - fFlowEvent->GetNumberOfRPs());
145 printf(" > number of RP's %i \n", fFlowEvent->GetNumberOfRPs());
148 //_____________________________________________________________________________
149 Bool_t AliAnalysisTaskJetFlow::PassesCuts(AliVEvent* event)
152 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
153 if(!event) return kFALSE;
154 if(TMath::Abs(InputEvent()->GetPrimaryVertex()->GetZ()) > 10.) return kFALSE;
155 // aod and esd specific checks
158 AliESDEvent* esdEvent = static_cast<AliESDEvent*>(InputEvent());
159 if( (!esdEvent) || (TMath::Abs(esdEvent->GetPrimaryVertexSPD()->GetZ() - esdEvent->GetPrimaryVertex()->GetZ()) > .5) ) return kFALSE;
162 AliAODEvent* aodEvent = static_cast<AliAODEvent*>(InputEvent());
163 if( (!aodEvent) || (TMath::Abs(aodEvent->GetPrimaryVertexSPD()->GetZ() - aodEvent->GetPrimaryVertex()->GetZ()) > .5) ) return kFALSE;
167 Float_t cent(InputEvent()->GetCentrality()->GetCentralityPercentile("V0M"));
168 return (cent <= fCentralityMin || cent > fCentralityMax || TMath::Abs(cent-InputEvent()->GetCentrality()->GetCentralityPercentile("TRK")) > 5.) ? kFALSE : kTRUE;
170 //_____________________________________________________________________________
171 void AliAnalysisTaskJetFlow::Terminate(Option_t *)
174 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
176 //_____________________________________________________________________________