]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskQualityAssurancePA.h
54f9d88ddbe3c14e28c0cdcc813c3a6715e6b035
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / UserTasks / AliAnalysisTaskQualityAssurancePA.h
1 #ifndef ALIANALYSISTASKQUALITYASSURANCEPA_H
2 #define ALIANALYSISTASKQUALITYASSURANCEPA_H
3
4 //  #define DEBUGMODE
5
6
7 class TH1F;
8 class TH2F;
9 class TList;
10 class TClonesArray;
11 class TString;
12 class AliEmcalJet;
13 class AliRhoParameter;
14 class AliVParticle;
15 class AliLog;
16 class AliAnalysisUtils;
17
18 #ifndef ALIANALYSISTASKSE_H
19 #include <Riostream.h>
20 #include <TROOT.h>
21 #include <TString.h>
22 #include <TFile.h>
23 #include <TCint.h>
24 #include <TChain.h>
25 #include <TTree.h>
26 #include <TKey.h>
27 #include <TProfile.h>
28 #include <TH1F.h>
29 #include <TH2F.h>
30 #include <TCanvas.h>
31 #include <TList.h>
32 #include <TClonesArray.h>
33 #include <TObject.h>
34 #include <TMath.h>
35 #include <TSystem.h>
36 #include <TInterpreter.h>
37 #include <TH1.h>
38 #include "AliAnalysisTask.h"
39 #include "AliCentrality.h"
40 #include "AliStack.h"
41 #include "AliESDEvent.h"
42 #include "AliESDInputHandler.h"
43 #include "AliAODEvent.h"
44 #include "AliAODHandler.h"
45 #include "AliAnalysisManager.h"
46 #include "AliAnalysisTaskSE.h"
47 #endif
48
49 #include <TRandom3.h>
50 #include "AliGenPythiaEventHeader.h"
51 #include "AliMCEvent.h"
52 #include "AliLog.h"
53 #include <AliEmcalJet.h>
54 #include <AliRhoParameter.h>
55 #include "AliVEventHandler.h"
56 #include "AliVParticle.h"
57 #include "AliAnalysisUtils.h"
58
59
60 class AliAnalysisTaskQualityAssurancePA : public AliAnalysisTaskSE {
61  public:
62
63   AliAnalysisTaskQualityAssurancePA() : AliAnalysisTaskSE(), fOutputList(0), fAnalyzeQA(1), fAnalyzeJets(1), fAnalyzePythia(0), fHasTracks(0), fHasClusters(0), fHasJets(0), fIsMC(0), fJetArray(0), fTrackArray(0), fClusterArray(0), fJetArrayName(0), fTrackArrayName(0), fClusterArrayName(0), fRunNumbers(0), fNumPtHardBins(11), fSignalJetRadius(0.4), fNumberExcludedJets(2), fSignalJetEtaWindow(0.5), fTrackEtaWindow(0.9), fClusterEtaWindow(0.7), fVertexWindow(10.0), fVertexMaxR(1.0), fMinTrackPt(0.150), fMinClusterPt(0.300), fMinJetPt(1.0), fMinJetArea(0.4), fFirstLeadingJet(0), fSecondLeadingJet(0), fNumberSignalJets(0), fCrossSection(0.0), fTrials(0.0), fRandom(0), fHelperClass(0), fInitialized(0), fTaskInstanceCounter(0), fHistList(0), fHistCount(0) {}
64   
65   AliAnalysisTaskQualityAssurancePA(const char *name, const char* trackArrayName, const char* clusterArrayName, const char* jetArrayName);
66
67   // Standard  functions
68   virtual ~AliAnalysisTaskQualityAssurancePA();
69   virtual void     UserCreateOutputObjects();
70   virtual void     UserExec(Option_t *option);
71   virtual void     Terminate(Option_t *);
72   
73   // Setters
74   void SetAnalyzeTracks(Bool_t val) {fAnalyzeQA = val;}
75   void SetAnalyzeJets(Bool_t val) {fAnalyzeJets = val;}
76   void SetAnalyzePythia(Bool_t val) {fAnalyzePythia = val;}
77
78   void SetTrackMinPt(Double_t minPt) {fMinJetPt = minPt;}
79   void SetSignalJetMinPt(Double_t minPt) {fMinJetPt = minPt;}
80   void SetSignalJetMinArea(Double_t minArea) {fMinJetArea = minArea;}
81   void SetRunNumbers(const char* runNumbers) {*fRunNumbers = runNumbers;}
82   void SetNumberOfPtHardBins(Int_t count) {fNumPtHardBins = count;}
83
84   void SetSignalJetRadius(Double_t radius) {fSignalJetRadius = radius;}
85   void SetAcceptanceWindows(Double_t trackEta, Double_t vertexZ, Double_t vertexMaxR, Double_t signalJetRadius)
86   {
87     fVertexWindow = vertexZ;
88     fVertexMaxR = vertexMaxR;
89     fTrackEtaWindow = trackEta;
90     fSignalJetRadius = signalJetRadius;
91     fSignalJetEtaWindow = fTrackEtaWindow-fSignalJetRadius;
92   }
93
94   // Getters
95   Int_t GetInstanceCounter() {return fTaskInstanceCounter;}
96
97  private:
98   // Calculation functions
99   void      GetSignalJets();
100   Int_t     GetLeadingJets(TClonesArray* jetArray, Int_t* jetIDArray);
101   Double_t  GetPtHard();
102   Int_t     GetPtHardBin();
103
104   // Cut checks
105   Bool_t    IsTrackInAcceptance(AliVParticle* track);
106   Bool_t    IsClusterInAcceptance(AliVCluster* cluster);
107   Bool_t    IsSignalJetInAcceptance(AliEmcalJet* jet);
108   
109   // Some helpers
110   Double_t  EtaToTheta(Double_t arg){return 2.*atan(exp(-arg));} 
111   Double_t  ThetaToEta(Double_t arg)
112   {
113     if ((arg > TMath::Pi()) || (arg < 0.0))
114     {
115       AliError(Form("ThetaToEta got wrong input! (%f)", arg));
116       return 0.0;
117     }
118     return -log(tan(arg/2.));
119   }
120   Double_t  GetDeltaPhi(Double_t phi1, Double_t phi2) {return min(TMath::Abs(phi1-phi2),TMath::TwoPi()- TMath::Abs(phi1-phi2));}
121
122   // #### This functions return the ratio of a rectangle that is covered by a circle
123   Double_t MCGetOverlapCircleRectancle(Double_t cPosX, Double_t cPosY, Double_t cRadius, Double_t rPosXmin, Double_t rPosXmax, Double_t rPosYmin, Double_t rPosYmax)
124   {
125     const Int_t kTests = 1000;
126     Int_t hits = 0;
127     TRandom3 randomGen(0);
128    
129     // Loop over kTests-many tests
130     for (Int_t i=0; i<kTests; i++)
131     {
132       //Choose random position in rectangle for the tester
133       Double_t tmpTestX = randomGen.Uniform(rPosXmin, rPosXmax);
134       Double_t tmpTestY = randomGen.Uniform(rPosYmin, rPosYmax);
135
136       //Check, if tester is in circle. If yes, increment circle counter.
137       Double_t tmpDistance = TMath::Sqrt( (tmpTestX - cPosX)*(tmpTestX - cPosX) + (tmpTestY - cPosY)*(tmpTestY - cPosY) );
138       if(tmpDistance < cRadius)
139         hits++;
140     }
141
142     // return ratio
143     return (static_cast<Double_t>(hits)/static_cast<Double_t>(kTests));
144   }
145
146   void FillHistogram(const char* runNumber, const char * key, Double_t x);
147   void FillHistogram(const char* runNumber, const char * key, Double_t x, Double_t y);
148   void FillHistogram(const char* runNumber, const char * key, Double_t x, Double_t y, Double_t add);
149   const char* GetHistoName(const char* runNumber, const char* name)
150   {
151     if (fIsMC)    
152       return Form("H%d_%s_%s_MC", fTaskInstanceCounter, name, runNumber);
153
154     return Form("H%d_%s_%s", fTaskInstanceCounter, name, runNumber);
155   }
156
157   template <class T> T* AddHistogram1D(const char* runNumber, const char* name = "CustomHistogram", const char* title = "NO_TITLE", const char* options = "", Int_t xBins = 100, Double_t xMin = 0.0, Double_t xMax = 20.0, const char* xTitle = "x axis", const char* yTitle = "y axis");
158
159   template <class T> T* AddHistogram2D(const char* runNumber, const char* name = "CustomHistogram", const char* title = "NO_TITLE", const char* options = "", Int_t xBins = 100, Double_t xMin = 0.0, Double_t xMax = 20.0, Int_t yBins = 100, Double_t yMin = 0.0, Double_t yMax = 20.0, const char* xTitle = "x axis", const char* yTitle = "y axis", const char* zTitle = "z axis");
160
161
162
163   // Standard functions
164   Bool_t    Notify();
165   void      Calculate(AliVEvent* event);
166   void      ExecOnce();
167   void      Init ();
168
169   TList*              fOutputList;            //! Output list
170   // ########## TRIGGERS 
171   Bool_t              fAnalyzeQA;             // trigger if tracks should be analyzed
172   Bool_t              fAnalyzeJets;           // trigger if jets should be processed
173   Bool_t              fAnalyzePythia;         // trigger if pythia properties should be processed
174   Bool_t              fHasTracks;             // trigger if tracks are actually valid
175   Bool_t              fHasClusters;           // trigger if clusters are actually valid
176   Bool_t              fHasJets;               // trigger if jets are actually valid
177   Bool_t              fIsMC;                  // trigger if data is MC (for naming reasons)
178   // ########## SOURCE INFORMATION
179   TClonesArray*       fJetArray;              //! object containing the jets
180   TClonesArray*       fTrackArray;            //! object containing the tracks
181   TClonesArray*       fClusterArray;          //! object containing the clusters
182   TString*            fJetArrayName;          // name of object containing the jets
183   TString*            fTrackArrayName;        // name of object containing the tracks
184   TString*            fClusterArrayName;      // name of object containing the tracks
185   TString*            fRunNumbers;            // analyzed run numbers
186   Int_t               fNumPtHardBins;         // Number of used pt hard bins
187
188   // ########## JET/DIJET/RC PROPERTIES
189   Double_t            fSignalJetRadius;       // Radius for the signal jets
190   Int_t               fNumberExcludedJets;    // Number of jets to be excluded
191
192   // ########## CUTS 
193   Double_t            fSignalJetEtaWindow;    // +- window in eta for signal jets
194   Double_t            fTrackEtaWindow;        // +- window in eta for tracks
195   Double_t            fClusterEtaWindow;      // +- window in eta for clusters
196   Double_t            fVertexWindow;          // +- window in Z for the vertex
197   Double_t            fVertexMaxR;            // +- window in R for the vertex (distance in xy-plane)
198   Double_t            fMinTrackPt;            // Min track pt to be accepted
199   Double_t            fMinClusterPt;          // Min track pt to be accepted
200   Double_t            fMinJetPt;              // Min jet pt to be accepted
201   Double_t            fMinJetArea;            // Min jet area to be accepted
202
203   // ########## EVENT PROPERTIES
204   AliEmcalJet*        fFirstLeadingJet;       //! leading jet in event
205   AliEmcalJet*        fSecondLeadingJet;      //! next to leading jet in event
206   Int_t               fNumberSignalJets;      //! Number of signal jets in event
207   AliEmcalJet*        fSignalJets[1024];      //! memory for signal jet pointers
208   Double_t            fCrossSection;          //! value is filled, if pythia header is accessible
209   Double_t            fTrials;                //! value is filled, if pythia header is accessible
210
211   // ########## GENERAL VARS
212   TRandom3*           fRandom;                //! A random number
213   AliAnalysisUtils*   fHelperClass;          //! Vertex selection helper
214   Bool_t              fInitialized;           //! trigger if tracks/jets are loaded
215   Int_t               fTaskInstanceCounter;   // for naming reasons
216   TList*              fHistList;              // Histogram list
217   Int_t               fHistCount;             // Histogram count
218
219   AliAnalysisTaskQualityAssurancePA(const AliAnalysisTaskQualityAssurancePA&); // not implemented
220   AliAnalysisTaskQualityAssurancePA& operator=(const AliAnalysisTaskQualityAssurancePA&); // not implemented
221
222   ClassDef(AliAnalysisTaskQualityAssurancePA, 1); // QA helper task for pA
223
224 };
225 #endif