]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/SPECTRA/IdentifiedHighPt/grid/AliAnalysisTaskHighPtDeDx.h
TPC rdEdx analysis code, macros and more (P.Christiansen)
[u/mrichter/AliRoot.git] / PWGLF / SPECTRA / IdentifiedHighPt / grid / AliAnalysisTaskHighPtDeDx.h
1 #ifndef ALIANALYSISTASKHIGHPTDEDX_H
2 #define ALIANALYSISTASKHIGHPTDEDX_H
3
4 // ROOT includes
5 #include <TList.h>
6 #include <TH1.h>
7 #include <TTreeStream.h>
8 #include <TRandom.h>
9 #include <TObject.h>
10
11 // AliRoot includes
12 #include <AliAnalysisTaskSE.h>
13 #include <AliESDEvent.h>
14 #include <AliAODEvent.h>
15 #include <AliMCEvent.h>
16 #include <AliAnalysisFilter.h>
17 #include <AliStack.h>
18 #include <AliGenEventHeader.h>
19 #include <AliVHeader.h>
20 #include <AliAODMCParticle.h> 
21
22 #include "DebugClasses.C"
23
24
25
26 class AliAnalysisTaskHighPtDeDx : public AliAnalysisTaskSE {
27  public:
28   enum AnalysisMode { kInvalid = -1, kGlobalTrk = 0x1, kTPCTrk = 0x2 }; 
29   AliAnalysisTaskHighPtDeDx();
30   AliAnalysisTaskHighPtDeDx(const char *name);
31
32   virtual ~AliAnalysisTaskHighPtDeDx();
33
34   virtual void   UserCreateOutputObjects();
35   virtual void   UserExec(Option_t *option);
36
37   Bool_t   GetAnalysisMC() { return fAnalysisMC; }   
38   Double_t GetVtxCut() { return fVtxCut; }   
39   Double_t GetEtaCut() { return fEtaCut; }     
40   Double_t GetMinPt() { return fMinPt; }   
41   Int_t    GetTreeOption() { return fTreeOption; }  
42
43   virtual void  SetTrigger1(UInt_t ktriggerInt1) {ftrigBit1 = ktriggerInt1;}
44   virtual void  SetTrigger2(UInt_t ktriggerInt2) {ftrigBit2 = ktriggerInt2;}
45   virtual void  SetTrackFilter(AliAnalysisFilter* trackF) {fTrackFilter = trackF;}
46   virtual void  SetTrackFilterGolden(AliAnalysisFilter* trackF) {fTrackFilterGolden = trackF;}
47   virtual void  SetTrackFilterTPC(AliAnalysisFilter* trackF) {fTrackFilterTPC = trackF;}
48   virtual void  SetProduceTPCBranch(Bool_t prodtpcb) {fTPCBranch = prodtpcb;}
49   virtual void  SetAnalysisType(const char* analysisType) {fAnalysisType = analysisType;}
50   virtual void  SetAnalysisMC(Bool_t isMC) {fAnalysisMC = isMC;}
51   virtual void  SetVtxCut(Double_t vtxCut){fVtxCut = vtxCut;}
52   virtual void  SetEtaCut(Double_t etaCut){fEtaCut = etaCut;}
53   virtual void  SetMinPt(Double_t value) {fMinPt = value;}   
54   virtual void  SetMinCent(Float_t minvalc) {fMinCent = minvalc;}
55   virtual void  SetMaxCent(Float_t maxvalc) {fMaxCent = maxvalc;}
56   virtual void  SetLowPtFraction(Double_t value) {fLowPtFraction = value;}   
57   virtual void  SetTreeOption(Int_t value) {fTreeOption = value;}    
58   virtual void  SetAnalysisPbPb(Bool_t isanaPbPb) {fAnalysisPbPb = isanaPbPb;}
59
60  private:
61   virtual Float_t GetVertex(const AliVEvent* event) const;
62   virtual void AnalyzeESD(AliESDEvent* esd); 
63   virtual void AnalyzeAOD(AliAODEvent* aod); 
64   virtual void ProduceArrayTrksESD(AliESDEvent* event, AnalysisMode anamode );
65   virtual void ProduceArrayTrksAOD(AliAODEvent* event, AnalysisMode anamode );
66   Short_t   GetPidCode(Int_t pdgCode) const;
67   void      ProcessMCTruthESD();
68   void      ProcessMCTruthAOD(); 
69   void      Sort(TClonesArray* array, Bool_t isMC);
70   Short_t   GetPythiaEventProcessType(Int_t pythiaType);
71   Short_t   GetDPMjetEventProcessType(Int_t dpmJetType);
72   ULong64_t GetEventIdAsLong(AliVHeader* header) const;
73
74   TParticle* FindPrimaryMother(AliStack* stack, Int_t label);
75   Int_t      FindPrimaryMotherLabel(AliStack* stack, Int_t label);
76
77   AliAODMCParticle* FindPrimaryMotherAOD(AliAODMCParticle* startParticle);
78
79   static const Double_t fgkClight;   // Speed of light (cm/ps)
80
81   AliESDEvent* fESD;                  //! ESD object
82   AliAODEvent* fAOD;                  //! AOD object
83   AliMCEvent*  fMC;                   //! MC object
84   AliStack*    fMCStack;              //! MC ESD stack
85   TClonesArray* fMCArray;             //! MC array for AOD
86   AliAnalysisFilter* fTrackFilter;    //  Track Filter, old cuts 2010
87   AliAnalysisFilter* fTrackFilterGolden;    //  Track Filter, set 2010 with golden cuts
88   AliAnalysisFilter* fTrackFilterTPC; // track filter for TPC only tracks
89   TString       fAnalysisType;        //  "ESD" or "AOD"
90   Bool_t        fAnalysisMC;          //  Real(kFALSE) or MC(kTRUE) flag
91   Bool_t        fAnalysisPbPb;        //  true you want to analyze PbPb data, false for pp
92   Bool_t        fTPCBranch;           //tru if you want to produce the TPC branch
93   TRandom*      fRandom;              //! random number generator
94   DeDxEvent*    fEvent;               //! event pointer
95   TClonesArray* fTrackArrayGlobalPar;          //! track array pointer, global tracks
96   TClonesArray* fTrackArrayTPCPar;          //! track array pointer, tpc track parameters
97   TClonesArray* fTrackArrayMC;        //! MC track array pointer
98   TClonesArray* fVZEROArray;          //! array of the v0 cells.
99
100
101   //
102   // Cuts and options
103   //
104   UInt_t       ftrigBit1;
105   UInt_t       ftrigBit2;
106   Double_t     fVtxCut;             // Vtx cut on z position in cm
107   Double_t     fEtaCut;             // Eta cut used to select particles
108   Double_t     fMinPt;              // Min pt - for histogram limits
109   Double_t     fLowPtFraction;      // Fraction of tracks below min pt to keep
110   Int_t        fTreeOption;         // 0: no tree, >0: enable debug tree
111   Float_t      fMinCent; //minimum centrality
112   Float_t      fMaxCent; //maximum centrality
113   //
114   // Help variables
115   //
116   Short_t      fMcProcessType;      // -1=invalid, 0=data, 1=ND, 2=SD, 3=DD
117   Short_t      fTriggeredEventMB;   // 1 = triggered, 0 = not trigged (MC only)
118   Short_t      fVtxStatus;          // -1 = no vtx, 0 = outside cut, 1 = inside cut
119   Float_t      fZvtx;               // z vertex
120   Float_t      fZvtxMC;             // z vertex MC (truth)
121   Int_t        fRun;                // run no
122   ULong64_t    fEventId;            // unique event id
123               
124   //
125   // Output objects
126   //
127   TList*        fListOfObjects;     //! Output list of objects
128   TH1I*         fEvents;            //! No of accepted events
129   TH1I*         fVtx;               //! Event vertex info
130   TH1I*         fVtxMC;             //! Event vertex info for ALL MC events
131   TH1F*         fVtxBeforeCuts;     //! Vertex z dist before cuts
132   TH1F*         fVtxAfterCuts;      //! Vertex z dist after cuts
133   TH1F* fn1;
134   TH1F* fn2;
135   TH1F* fcent;
136   TTree*        fTree;              //! Debug tree 
137
138   ClassDef(AliAnalysisTaskHighPtDeDx, 1);    //Analysis task for high pt analysis 
139 };
140
141 #endif