]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/JetTasks/AliAnalysisTaskMinijet.h
Covarity fix
[u/mrichter/AliRoot.git] / PWG4 / JetTasks / AliAnalysisTaskMinijet.h
1 #ifndef ALIANALYSISTASKMINIJET_H\r
2 #define ALIANALYSISTASKMINIJET_H\r
3 \r
4 // Two-particle correlations using all particles over pt threshold\r
5 // Extract mini-jet yield and fragmentation properties via Delta-Phi histograms\r
6 // Can use ESD or AOD, reconstructed and Monte Carlo data as input\r
7 // Author: eva.sicking@cern.ch\r
8 \r
9 class TList;\r
10 class TH1F;\r
11 class TH2F;\r
12 class TProfile;\r
13 \r
14 class AliESDtrackCuts;\r
15 \r
16 \r
17 #include "AliAnalysisTaskSE.h"\r
18 \r
19 class AliAnalysisTaskMinijet : public AliAnalysisTaskSE {\r
20  public:\r
21   AliAnalysisTaskMinijet(const char *name="<default name>");\r
22   virtual ~AliAnalysisTaskMinijet();\r
23   \r
24   virtual void   UserCreateOutputObjects();\r
25   virtual void   UserExec(Option_t* option);\r
26   virtual void   Terminate(Option_t *);\r
27   \r
28   Int_t LoopESD  (Float_t **pt, Float_t **eta, Float_t **phi, Short_t **charge, Int_t **nTracksTracklets);\r
29   Int_t LoopESDMC(Float_t **pt, Float_t **eta, Float_t **phi, Short_t **charge,  Int_t **nTracksTracklets);\r
30   Int_t LoopAOD  (Float_t **pt, Float_t **eta, Float_t **phi, Short_t **charge,  Int_t **nTracksTracklets);\r
31   Int_t LoopAODMC(Float_t **pt, Float_t **eta, Float_t **phi, Short_t **charge,  Int_t **nTracksTracklets);\r
32   void  Analyse  (const Float_t* pt, const Float_t* eta, const Float_t* phi,  const Short_t *charge, Int_t ntacks, Int_t ntacklets=0, const Int_t nAll=0, Int_t mode=0);\r
33   const void  CleanArrays(const Float_t *pt, const Float_t *eta, const Float_t *phi, const Short_t *charge, const Int_t *nTracksTracklets=0);\r
34   const Bool_t SelectParticlePlusCharged(Short_t charge, Int_t pdg, Bool_t prim);\r
35   const Bool_t SelectParticle(Short_t charge, Int_t pdg, Bool_t prim);\r
36 \r
37 \r
38   void UseMC(Bool_t useMC=kTRUE, Bool_t mcOnly=kFALSE)    {fUseMC = useMC; fMcOnly=mcOnly;}\r
39 \r
40   virtual void   SetCuts(AliESDtrackCuts* cuts)           {fCuts = cuts;}\r
41 \r
42   void   SetRadiusCut(Float_t radiusCut)                  {fRadiusCut = radiusCut;}  \r
43   void   SetTriggerPtCut(Float_t triggerPtCut)            {fTriggerPtCut = triggerPtCut;}  \r
44   void   SetAssociatePtCut(Float_t associatePtCut)        {fAssociatePtCut = associatePtCut;} \r
45   void   SetEventAxis(Int_t leadingOrRandom)              {fLeadingOrRandom = leadingOrRandom;}  \r
46   void   SetMode(Int_t mode)                              {fMode = mode;}\r
47   void   SetMaxVertexZ(Float_t vertexZCut)                {fVertexZCut = vertexZCut;}\r
48   void   SetMaxEta(Float_t etaCut)                        {fEtaCut = etaCut;}\r
49   void   SetMaxEtaSeed(Float_t etaCutSeed)                {fEtaCutSeed = etaCutSeed;}\r
50 \r
51   void   SelectParticles(Int_t selectParticles)           {fSelectParticles = selectParticles;}\r
52   void   SelectParticlesAssoc(Int_t selectParticlesAssoc) {fSelectParticlesAssoc = selectParticlesAssoc;}\r
53 \r
54 \r
55  private:\r
56 \r
57   Bool_t       fUseMC;                      // flag for Monte Carlo usages\r
58   Bool_t       fMcOnly;                     // flag defines, if only MC data is used in analysis or also reconstructed data\r
59   AliESDtrackCuts* fCuts;                   // List of cuts for ESDs\r
60   Float_t      fRadiusCut;                  // radius cut \r
61   Float_t      fTriggerPtCut;               // cut on particle pt used as event axis\r
62   Float_t      fAssociatePtCut;             // cut on particle pt used for correlations\r
63   Int_t        fLeadingOrRandom;            // event axis:leading track or random track\r
64   Int_t        fMode;                       // ESD(=0) of AOD(=1) reading \r
65   Float_t      fVertexZCut;                 // vertex cut\r
66   Float_t      fEtaCut;                     // eta acceptance cut\r
67   Float_t      fEtaCutSeed;                 // eta acceptance cut for seed\r
68   Int_t        fSelectParticles;            // only in cas of MC: use also neutral particles or not \r
69   Int_t        fSelectParticlesAssoc;       // only in cas of MC: use also neutral particles or not \r
70 \r
71   AliESDEvent *fESDEvent;                   //! esd event\r
72   AliAODEvent *fAODEvent;                   //! aod event\r
73   Int_t        fNMcPrimAccept;              // global variable for mc multiplucity\r
74   Float_t      fVzEvent;                    // global variable for rec vertex position\r
75   \r
76   TList      *fHists;                      // output list\r
77   TH1F       *fHistPt;                     // Pt spectrum ESD\r
78   TH1F       *fHistPtMC;                   // Pt spectrum MC\r
79   TH2F       *fNmcNch;                     // N mc - N ch rec\r
80   TProfile   *fPNmcNch;                     // N mc - N ch rec\r
81   TH2F       *fChargedPi0;                 // charged versus charged+Pi0\r
82   TH1F       * fVertexZ[4];                 // z of vertex\r
83 \r
84   TH1F       * fPt[4];                      // pt\r
85   TH1F       * fEta[4];                     // et\r
86   TH1F       * fPhi[4];                     // phi\r
87   TH1F       * fDcaXY[4];                   // dca xy direction\r
88   TH1F       * fDcaZ[4];                    // dca z direction\r
89 \r
90   TH1F       * fPtSeed[4];                  // pt of seed (event axis)\r
91   TH1F       * fEtaSeed[4];                 // eta of seed \r
92   TH1F       * fPhiSeed[4];                 // phi of seed\r
93 \r
94   TH1F       * fPtOthers[4];                // pt of all other particels used in dEtadPhi\r
95   TH1F       * fEtaOthers[4];               // eta of all other particels used in dEtadPhi\r
96   TH1F       * fPhiOthers[4];               // phi of all other particels used in dEtadPhi\r
97   TH2F       * fPtEtaOthers[4];             // pt-eta of all other particels used in dEtadPhi\r
98 \r
99 \r
100   TH2F       * fPhiEta[4];                  // eta - phi\r
101   TH2F       * fDPhiDEtaEventAxis[4];       // correlation dEta-dPhi towards event axis\r
102   TH2F       * fDPhiDEtaEventAxisSeeds[4];  // correlation dEta-dPhi towards event axis of trigger particles\r
103   TH1F       * fTriggerNch[4];              // number of triggers with accepted-track number\r
104   TH2F       * fTriggerNchSeeds[4];         // number of triggers with accepted-track number\r
105   TH1F       * fTriggerTracklet[4];         // number of triggers with accepted-tracklet number\r
106   TH2F       * fNch07Nch[4];                // nCharged with pT>fTriggerPtCut vs nCharged\r
107   TProfile   * fPNch07Nch[4];                // nCharged with pT>fTriggerPtCut vs nCharged\r
108   TH2F       * fNch07Tracklet[4];           // nCharged with pT>fTriggerPtCut vs nTracklet\r
109   TH2F       * fNchTracklet[4];             // nCharged vs nTracklet\r
110   TProfile   * fPNch07Tracklet[4];           // nCharged with pT>fTriggerPtCut vs nTracklet\r
111 \r
112   TH1F       * fDPhiEventAxis[4];           // delta phi of associate tracks to event axis\r
113   TH1F       * fDPhiEventAxisNchBin[4][150];// delta phi of associate tracks to event axis per Nch bin\r
114   TH1F       * fDPhiEventAxisNchBinTrig[4][150];// "" for all possoble trigger particles\r
115 \r
116   TH1F       * fDPhiEventAxisTrackletBin[4][150]; // delta phi of associate tracks to event axis per Nch bin\r
117   TH1F       * fDPhiEventAxisTrackletBinTrig[4][150]; // "" for all possible trigger particles\r
118 \r
119   AliAnalysisTaskMinijet(const AliAnalysisTaskMinijet&); // not implemented\r
120   AliAnalysisTaskMinijet& operator=(const AliAnalysisTaskMinijet&); // not implemented\r
121   \r
122   ClassDef(AliAnalysisTaskMinijet, 1); // example of analysis\r
123 };\r
124 \r
125 #endif\r