]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
from Redmer Bertens:
authormkrzewic <mkrzewic@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 15 Mar 2013 16:24:42 +0000 (16:24 +0000)
committermkrzewic <mkrzewic@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 15 Mar 2013 16:24:42 +0000 (16:24 +0000)
  add a stub for the jet flow analysis task with the accompanying Add... macro

PWG/CMakelibPWGflowTasks.pkg
PWG/FLOW/Tasks/AliAnalysisTaskJetFlow.cxx [new file with mode: 0644]
PWG/FLOW/Tasks/AliAnalysisTaskJetFlow.h [new file with mode: 0644]
PWG/PWGflowTasksLinkDef.h
PWGCF/FLOW/macros/AddTaskJetFlow.C [new file with mode: 0644]

index d562b29762576238a335dddc264915a0fd76bde3..969ea607f9e996a705ec60b3163cd4f9848206db 100644 (file)
@@ -57,6 +57,7 @@ set ( SRCS
     FLOW/Tasks/AliFlowVZEROResults.cxx
     FLOW/Tasks/AliFlowVZEROQA.cxx
     FLOW/Tasks/AliAnalysisTaskFlowEPCascade.cxx
+    FLOW/Tasks/AliAnalysisTaskJetFlow.cxx
     )
 
 string ( REPLACE ".cxx" ".h" HDRS "${SRCS}" )
diff --git a/PWG/FLOW/Tasks/AliAnalysisTaskJetFlow.cxx b/PWG/FLOW/Tasks/AliAnalysisTaskJetFlow.cxx
new file mode 100644 (file)
index 0000000..fdf4d45
--- /dev/null
@@ -0,0 +1,136 @@
+/**************************************************************************
+ * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ *                                                                        *
+ * Author: The ALICE Off-line Project.                                    *
+ * Contributors are mentioned in the code where appropriate.              *
+ *                                                                        *
+ * Permission to use, copy, modify and distribute this software and its   *
+ * documentation strictly for non-commercial purposes is hereby granted   *
+ * without fee, provided that the above copyright notice appears in all   *
+ * copies and that both the copyright notice and this permission notice   *
+ * appear in the supporting documentation. The authors make no claims     *
+ * about the suitability of this software for any purpose. It is          *
+ * provided "as is" without express or implied warranty.                  *
+ **************************************************************************/
+
+/* 
+ * AliAnaysisTaskJet
+ * author: Redmer Alexander Bertens
+ * rbertens@cern.ch, rbertens@nikhef.nl, r.a.bertens@uu.nl 
+ *
+ * jet correlation task - correlates jets to the vzero reaction plane using
+ * the scalar product method 
+ *
+ * this task should be run AFTER 
+ * $ALICE_ROOT/PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskRhoVnModulation()
+ *
+ * */
+
+#include <TChain.h>
+#include <TH1F.h>
+#include <AliAnalysisTask.h>
+#include <AliAnalysisManager.h>
+#include <AliFlowEventSimple.h>
+#include <AliFlowTrack.h>
+#include <AliFlowTrackCuts.h>
+#include <AliFlowCommonConstants.h>
+#include <TString.h>
+
+#include "AliAnalysisTaskJetFlow.h"
+
+class AliAnalysisTaskJetFlow;
+
+using namespace std;
+
+ClassImp(AliAnalysisTaskJetFlow)
+
+AliAnalysisTaskJetFlow::AliAnalysisTaskJetFlow() : AliAnalysisTaskSE(), 
+    fDebug(-1), fJetsName(0), fOutputList(0), fCutsRP(0), fCutsPOI(0), fCutsNull(0), fFlowEvent(0)/* , fHistPt(0) */
+{
+    // default constructor
+}
+//_____________________________________________________________________________
+AliAnalysisTaskJetFlow::AliAnalysisTaskJetFlow(const char* name) : AliAnalysisTaskSE(name),
+    fDebug(-1), fJetsName(0), fOutputList(0), fCutsRP(0), fCutsPOI(0), fCutsNull(0), fFlowEvent(0) /*, fHistPt(0) */
+{
+    // constructor
+    DefineInput(0, TChain::Class());
+    DefineOutput(1, TList::Class());
+    DefineOutput(2, AliFlowEventSimple::Class());
+}
+//_____________________________________________________________________________
+AliAnalysisTaskJetFlow::~AliAnalysisTaskJetFlow()
+{
+    // destructor
+    if(fOutputList)     delete fOutputList;
+    if(fFlowEvent)      delete fFlowEvent;
+}
+//_____________________________________________________________________________
+void AliAnalysisTaskJetFlow::UserCreateOutputObjects()
+{
+    // create output objects
+    if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
+    fOutputList = new TList();
+    fOutputList->SetOwner(kTRUE);
+    // histograms - not that this task does not produce output itself at this moment
+    /* fHistPt = new TH1F("fHistPt", "fHistPt", 100, 0, 100); */
+    /* fOutputList->Add(fHistPt); */
+    PostData(1, fOutputList);
+    // create the flow event and configure the static cc object
+    fFlowEvent = new AliFlowEvent(1000);
+    PostData(2, fFlowEvent);
+    AliFlowCommonConstants* cc = AliFlowCommonConstants::GetMaster();
+    cc->SetPtMax(100);
+    cc->SetNbinsPt(20);
+}
+//_____________________________________________________________________________
+void AliAnalysisTaskJetFlow::UserExec(Option_t *)
+{
+    // user exec
+    if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
+    
+    if(!InputEvent() || !fCutsNull || !fCutsRP) return; // coverity (and sanity)
+    // get the jet array, which is added as an extension to the AliVEvent by the jetfinder
+    TClonesArray* jets = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fJetsName.Data()));
+    if(jets) {
+        Int_t iJets = jets->GetEntriesFast();
+        if(iJets <= 0) {
+            if(fDebug>0) printf(" > Retrieved empty AliVEvent extension, aborting ! < \n ");
+            return;
+        }
+        // prepare the track selection for RP's and the flow event
+        fCutsRP->SetEvent(InputEvent(), MCEvent());
+        fCutsNull->SetEvent(InputEvent(), MCEvent());
+        fFlowEvent->ClearFast();
+        fFlowEvent->Fill(fCutsRP, fCutsNull);
+        fFlowEvent->SetReferenceMultiplicity(64);       // hardcoded for vzero
+        fFlowEvent->DefineDeadZone(0, 0, 0, 0);
+        // loop over jets and inject them as POI's
+        for(Int_t i(0); i < iJets; i++) {
+            AliVParticle* jet = static_cast<AliVParticle*>(jets->At(i));
+            if(jet) {
+                /* fHistPt->Fill(jet->Pt()); */
+                if(jet->Pt() <= 0) continue;
+                AliFlowTrack* flowTrack = new AliFlowTrack(jet);
+                flowTrack->SetForPOISelection(kTRUE);
+                flowTrack->SetForRPSelection(kFALSE);
+                fFlowEvent->InsertTrack(flowTrack);
+            }       
+        }
+    }
+    else if(fDebug > 0 ) printf(" Failed to find TClones Jet array while using name %s \n ", fJetsName.Data());
+    PostData(1, fOutputList);
+    PostData(2, fFlowEvent);
+    if(fDebug>0) {
+        printf(" Event summary \n");
+        printf(" > number of POI's (jets) %i ", fFlowEvent->NumberOfTracks() - fFlowEvent->GetNumberOfRPs());
+        printf(" > number of RP's %i \n", fFlowEvent->GetNumberOfRPs());
+    }
+}
+//_____________________________________________________________________________
+void AliAnalysisTaskJetFlow::Terminate(Option_t *)
+{
+    // terminate
+    if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
+}
+//_____________________________________________________________________________
diff --git a/PWG/FLOW/Tasks/AliAnalysisTaskJetFlow.h b/PWG/FLOW/Tasks/AliAnalysisTaskJetFlow.h
new file mode 100644 (file)
index 0000000..07288a7
--- /dev/null
@@ -0,0 +1,51 @@
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. */
+/* See cxx source for full Copyright notice */
+/* $Id$ */
+
+#ifndef AliAnalysisTaskJetFlow_H
+#define AliAnalysisTaskJetFlow_H
+
+#include <AliAnalysisTaskSE.h>
+#include <AliFlowTrackCuts.h>
+#include <AliFlowEvent.h>
+#include <TString.h>
+
+class AliAnalysisTaskJetFlow : public AliAnalysisTaskSE
+{
+    public:
+                                AliAnalysisTaskJetFlow();
+                                AliAnalysisTaskJetFlow(const char *name);
+        virtual                 ~AliAnalysisTaskJetFlow();
+
+        virtual void            UserCreateOutputObjects();
+        virtual void            UserExec(Option_t* option);
+        virtual void            Terminate(Option_t* option);
+        // setters
+        void                    SetCutsRP(AliFlowTrackCuts* c)          {fCutsRP = c;}
+        void                    SetCutsPOI(AliFlowTrackCuts* c)         {fCutsPOI = c;}
+        void                    SetCutsNull(AliFlowTrackCuts* c)        {fCutsNull = c;}
+        void                    SetJetCollectionName(TString jets)      {fJetsName = jets;}
+        void                    SetDebugMode(Int_t d)                   {fDebug = d;}
+
+    private:
+
+        // technical stuff
+        Int_t                   fDebug;         // debug level (0 none, 1 fcn calls, 2 verbose)
+        TString                 fJetsName;      // name of jet list
+        TList*                  fOutputList;    //! output list
+        // cut objects
+        AliFlowTrackCuts*       fCutsRP;        // rp cuts
+        AliFlowTrackCuts*       fCutsPOI;       // poi cuts
+        AliFlowTrackCuts*       fCutsNull;      // empty cuts
+        // input, output
+        AliFlowEvent*           fFlowEvent;     //! flow event
+        // histograms
+        /* TH1F*                   fHistPt;        //! dummy histogram */
+
+        AliAnalysisTaskJetFlow(const AliAnalysisTaskJetFlow&); // not implemented
+        AliAnalysisTaskJetFlow& operator=(const AliAnalysisTaskJetFlow&); // not implemented
+
+        ClassDef(AliAnalysisTaskJetFlow, 1);
+};
+
+#endif
index 9ae59c453e97e8e9108426e8f20ae2b1873bcf05..0a5c5f513ee08e3aa8e2cd98f6943917090f9959 100644 (file)
@@ -37,6 +37,7 @@
 #pragma link C++ class AliFlowVZEROResults+;
 #pragma link C++ class AliFlowVZEROQA+;
 #pragma link C++ class AliAnalysisTaskFlowEPCascade+;
+#pragma link C++ class AliAnalysisTaskJetFlow+;
 
 #endif
 
diff --git a/PWGCF/FLOW/macros/AddTaskJetFlow.C b/PWGCF/FLOW/macros/AddTaskJetFlow.C
new file mode 100644 (file)
index 0000000..c43c404
--- /dev/null
@@ -0,0 +1,101 @@
+///////////////////////////////////////////////////////////////////
+//                                                               //            
+//                     AddTaskJetFlow                            //
+// Author: Redmer A. Bertens, Utrecht University, 2013           //
+//                                                               //
+///////////////////////////////////////////////////////////////////
+class AliAnalysisDataContainer;
+class AliFlowTrackCuts;
+
+AliAnalysisTaskJetFlow* AddTaskJetFlow( TString name = "name",
+                                        TString jets = "jets",
+                                        Bool_t debug = kTRUE  )
+
+{
+    // first check the environment (common to all tasks)
+    if(debug) printf("\n\n  >> AddTaskJetFlow <<\n");
+    TString fileName = AliAnalysisManager::GetCommonFileName();
+    fileName += ":JetFlow";
+    if(debug) printf("      - filename: %s \n",fileName.Data());
+    AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
+    if (!mgr) {
+        if(debug) cout << " Fatal error: no analysis manager found! " << endl;
+        return 0x0;
+    }
+    if (!mgr->GetInputEventHandler()) {
+        if(debug) cout << " Fatal error: no imput event handler found!" << endl;
+        return 0x0;
+    }
+   // create the task
+   AliAnalysisTaskJetFlow* task = new AliAnalysisTaskJetFlow(name.Data());   
+   if(!task) {
+        if(debug) cout << " --> Unexpected error occurred: NO TASK WAS CREATED! (could be a library problem!) " << endl;
+        return 0x0;
+    }
+   else printf(" > added task with name %s and jet collection %s < \n", name.Data(), jets);
+   task->SetJetCollectionName(jets.Data());
+
+   // pass specific objects to the task
+   AliFlowTrackCuts* CutsRP = new AliFlowTrackCuts("CutsRP");
+   CutsRP = CutsRP->GetStandardVZEROOnlyTrackCuts();
+   task->SetCutsRP(CutsRP);
+
+   AliFlowTrackCuts* CutsNull = new AliFlowTrackCuts("CutsNull");
+   CutsNull->SetParamType(AliFlowTrackCuts::kGlobal);
+   CutsNull->SetEtaRange(+1, -1);
+   CutsNull->SetPtRange(+1, -1);
+   task->SetCutsNull(CutsNull);
+   
+   mgr->AddTask(task);
+   mgr->ConnectInput(task,0,mgr->GetCommonInputContainer());
+   mgr->ConnectOutput(task,1,mgr->CreateContainer("GenericOutputContainer", TList::Class(), AliAnalysisManager::kOutputContainer, fileName.Data()));
+
+   // connect flow anaysis task
+   AliAnalysisDataContainer *flowEvent = mgr->CreateContainer("flowEvent", AliFlowEventSimple::Class(), AliAnalysisManager::kExchangeContainer);
+   mgr->ConnectOutput(task, 2, flowEvent);
+   TaskJetFlow::AddSPmethod("SP_A", -0.8, -0, 0, +0.8, "Qa", 2, flowEvent);
+   TaskJetFlow::AddSPmethod("SP_B", -0.8, -0, 0, +0.8, "Qb", 2, flowEvent);
+
+   return task;
+
+}
+//_____________________________________________________________________________
+namespace TaskJetFlow{ // use unique namespace to avoid problems in TRAIN analysis
+    void AddSPmethod(char *name, double minEtaA, double maxEtaA, double minEtaB, double maxEtaB, char *Qvector, int harmonic, AliAnalysisDataContainer *flowEvent)
+    {
+       // add sp task and filter tasks
+       AliFlowTrackSimpleCuts* a = new AliFlowTrackSimpleCuts("a");
+       a->SetEtaMin(minEtaA);
+       a->SetEtaMax(maxEtaA);
+       a->SetMassMin(-999);
+       a->SetMassMax(999);
+       AliFlowTrackSimpleCuts* b = new AliFlowTrackSimpleCuts("b");
+       b->SetEtaMin(minEtaB);
+       b->SetEtaMax(maxEtaB);
+       b->SetMassMin(-999);
+       b->SetMassMax(999);
+       TString fileName = AliAnalysisManager::GetCommonFileName();
+       fileName+=":SP";
+       AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
+       AliAnalysisDataContainer *filteredFlowEvent = mgr->CreateContainer(Form("Filter_%s", name), AliFlowEventSimple::Class(), AliAnalysisManager::kExchangeContainer);
+       AliAnalysisTaskFilterFE *tskFilter(0x0);
+       if(!strcmp("Qa", Qvector)) tskFilter = new AliAnalysisTaskFilterFE(Form("TaskFilter_%s", name), 0x0, a);
+       if(!strcmp("Qb", Qvector)) tskFilter = new AliAnalysisTaskFilterFE(Form("TaskFilter_%s", name), 0x0, b);
+       if(tskFilter) tskFilter->SetSubeventEtaRange(-10, 0, 0, 10);
+       else return;
+       mgr->AddTask(tskFilter);
+       printf( " > Added FILTER TASK < \n ");
+       mgr->ConnectInput(tskFilter, 0, flowEvent);
+       mgr->ConnectOutput(tskFilter, 1, filteredFlowEvent);
+       AliAnalysisDataContainer *outSP = mgr->CreateContainer(name, TList::Class(), AliAnalysisManager::kOutputContainer, fileName);
+       AliAnalysisTaskScalarProduct *tskSP = new AliAnalysisTaskScalarProduct(Form("TaskScalarProduct_%s", name), kFALSE);
+       tskSP->SetApplyCorrectionForNUA(kTRUE);
+       tskSP->SetHarmonic(harmonic);
+       tskSP->SetTotalQvector(Qvector);
+       tskSP->SetBookOnlyBasicCCH(kFALSE);
+       mgr->AddTask(tskSP);
+       mgr->ConnectInput(tskSP, 0, filteredFlowEvent);
+       mgr->ConnectOutput(tskSP, 1, outSP);
+    }
+}// end of namespace TaskJetFlow