]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/SPECTRA/AliProtonQAAnalysis.h
Adding more QA histograms: eta vs phi vs Nclusters
[u/mrichter/AliRoot.git] / PWG2 / SPECTRA / AliProtonQAAnalysis.h
1 #ifndef ALIPROTONQAANALYSIS_H
2 #define ALIPROTONQAANALYSIS_H
3
4 /*  See cxx source for full Copyright notice */
5
6
7 /* $Id: AliProtonQAAnalysis.h 29114 2008-10-03 16:49:02Z pchrist $ */
8
9 //-------------------------------------------------------------------------
10 //                       Class AliProtonQAAnalysis
11 //   This is the class for the baryon (proton) analysis
12 //
13 //    Origin: Panos Christakoglou | Panos.Christakoglou@cern.ch
14 //-------------------------------------------------------------------------
15
16 #include "TObject.h"
17 #include "TH1I.h"
18 #include "TList.h"
19 #include "TArrayI.h"
20
21 #include "AliPID.h"
22
23 class TF1;
24 class TH1F;
25 class TH3F;
26
27 class AliESDEvent;
28 class AliESDtrack;
29 class AliStack;
30 class AliGenEventHeader;
31 class AliESDVertex;
32
33 class AliProtonQAAnalysis : public TObject {
34  public:
35   AliProtonQAAnalysis();
36   virtual ~AliProtonQAAnalysis();
37
38   void SetEtaMode() {fAnalysisEtaMode = kTRUE;}
39
40   void UseTPCOnly() {fUseTPCOnly = kTRUE;}
41   void UseHybridTPC() {fUseTPCOnly = kTRUE; fUseHybridTPC = kTRUE;}
42
43   //Cut functions
44   void    SetPointOnITSLayer1() {fPointOnITSLayer1Flag = kTRUE;}
45   void    SetPointOnITSLayer2() {fPointOnITSLayer2Flag = kTRUE;}
46   void    SetPointOnITSLayer3() {fPointOnITSLayer3Flag = kTRUE;}
47   void    SetPointOnITSLayer4() {fPointOnITSLayer4Flag = kTRUE;}
48   void    SetPointOnITSLayer5() {fPointOnITSLayer5Flag = kTRUE;}
49   void    SetPointOnITSLayer6() {fPointOnITSLayer6Flag = kTRUE;}
50   void    SetMinITSClusters(Int_t minITSClusters) {
51     fMinITSClusters = minITSClusters;
52     fMinITSClustersFlag = kTRUE;
53   }
54   void    SetMaxChi2PerITSCluster(Double_t maxChi2PerITSCluster) {
55     fMaxChi2PerITSCluster = maxChi2PerITSCluster;
56     fMaxChi2PerITSClusterFlag = kTRUE;
57   }
58   void    SetMinTPCClusters(Int_t minTPCClusters) {
59     fMinTPCClusters = minTPCClusters;
60     fMinTPCClustersFlag = kTRUE;
61   }
62   void    SetMaxChi2PerTPCCluster(Double_t maxChi2PerTPCCluster) {
63     fMaxChi2PerTPCCluster = maxChi2PerTPCCluster;
64     fMaxChi2PerTPCClusterFlag = kTRUE;
65   }
66   void    SetMaxCov11(Double_t maxCov11) {
67     fMaxCov11 = maxCov11; fMaxCov11Flag = kTRUE;}
68   void    SetMaxCov22(Double_t maxCov22) {
69     fMaxCov22 = maxCov22; fMaxCov22Flag = kTRUE;}
70   void    SetMaxCov33(Double_t maxCov33) {
71     fMaxCov33 = maxCov33; fMaxCov33Flag = kTRUE;}
72   void    SetMaxCov44(Double_t maxCov44) {
73     fMaxCov44 = maxCov44; fMaxCov44Flag = kTRUE;}
74   void    SetMaxCov55(Double_t maxCov55) {
75     fMaxCov55 = maxCov55; fMaxCov55Flag = kTRUE;}
76   void    SetMaxSigmaToVertex(Double_t maxSigmaToVertex) {
77     fMaxSigmaToVertex = maxSigmaToVertex;
78     fMaxSigmaToVertexFlag = kTRUE;
79   }
80   void    SetMaxSigmaToVertexTPC(Double_t maxSigmaToVertex) {
81     fMaxSigmaToVertexTPC = maxSigmaToVertex;
82     fMaxSigmaToVertexTPCFlag = kTRUE;
83   }
84   void    SetMaxDCAXY(Double_t maxDCAXY) {
85     fMaxDCAXY = maxDCAXY;
86     fMaxDCAXYFlag = kTRUE;
87   }
88   void    SetMaxDCAXYTPC(Double_t maxDCAXY) {
89     fMaxDCAXYTPC = maxDCAXY;
90     fMaxDCAXYTPCFlag = kTRUE;
91   }
92   void    SetMaxDCAZ(Double_t maxDCAZ) {
93     fMaxDCAZ = maxDCAZ;
94     fMaxDCAZFlag = kTRUE;
95   }
96   void    SetMaxDCAZTPC(Double_t maxDCAZ) {
97     fMaxDCAZTPC = maxDCAZ;
98     fMaxDCAZTPCFlag = kTRUE;
99   }
100   void    SetMaxConstrainChi2(Double_t maxConstrainChi2) {
101     fMaxConstrainChi2 = maxConstrainChi2;
102     fMaxConstrainChi2Flag = kTRUE;
103   }
104   void    SetITSRefit() {fITSRefitFlag = kTRUE;}
105   void    SetTPCRefit() {fTPCRefitFlag = kTRUE;}
106   void    SetESDpid() {fESDpidFlag = kTRUE;}
107   void    SetTPCpid() {fTPCpidFlag = kTRUE;}
108
109   //Prior probabilities
110   void SetPriorProbabilities(Double_t *partFrac) {
111     for(Int_t i = 0; i < AliPID::kSPECIESN; i++) fPartFrac[i] = partFrac[i];} 
112   void SetPriorProbabilityFunctions(TF1 *felectron, TF1 *fmuon, TF1 *fpion, TF1 *fkaon, TF1 *fproton) {
113     fFunctionProbabilityFlag = kTRUE;
114     fElectronFunction = felectron; 
115     fMuonFunction = fmuon; 
116     fPionFunction = fpion;
117     fKaonFunction = fkaon;
118     fProtonFunction = fproton;
119   } 
120   Double_t GetParticleFraction(Int_t i, Double_t p);
121
122   //Vertex QA
123   void RunVertexQA(AliGenEventHeader *header,
124                    AliStack *stack,
125                    AliESDEvent *esd);
126   TList *GetVertexQAList() {return fQAVertexList;}
127
128   //QA histograms
129   void SetQAYPtBins(Int_t nbinsY, Double_t minY, Double_t maxY,
130                     Int_t nbinsPt, Double_t minPt, Double_t maxPt);
131   void RunQAAnalysis(AliStack *stack, 
132                      AliESDEvent *esd,
133                      const AliESDVertex *vertex);
134   void SetRunQAAnalysis();
135   TList *GetGlobalQAList() {return fGlobalQAList;}
136
137   //Efficiency plots (reconstruction & PID)
138   void RunEfficiencyAnalysis(AliStack *stack, 
139                              AliESDEvent *esd,
140                              const AliESDVertex *vertex);
141   void SetRunEfficiencyAnalysis(Bool_t gUseCuts) {
142     fRunEfficiencyAnalysis = kTRUE;
143     fUseCutsInEfficiency = gUseCuts;
144   }
145   TList *GetEfficiencyQAList() {return fEfficiencyList;}
146
147   //MC analysis
148   void RunMCAnalysis(AliStack* stack);
149   void SetRunMCAnalysis() {fRunMCAnalysis = kTRUE;}
150   void SetMCProcessId(Int_t id) {
151     fMCProcessIdFlag = kTRUE;
152     fMCProcessId = id;
153   }
154   void SetMotherParticlePDGCode(Int_t pdgCode) {
155     fMotherParticlePDGCodeFlag = kTRUE;
156     fMotherParticlePDGCode = pdgCode;
157   }
158   TList *GetPDGList() {return fPDGList;}
159   TList *GetMCProcessesList() {return fMCProcessesList;}
160
161   TList *GetAcceptedCutList() {return fAcceptedCutList;}
162   TList *GetRejectedCutList() {return fRejectedCutList;}
163   TList *GetAcceptedDCAList() {return fAcceptedDCAList;}
164   TList *GetRejectedDCAList() {return fRejectedDCAList;}
165
166  private:
167   AliProtonQAAnalysis(const AliProtonQAAnalysis&); // Not implemented
168   AliProtonQAAnalysis& operator=(const AliProtonQAAnalysis&);// Not implemented
169
170   Bool_t   IsAccepted(AliESDEvent *esd,
171                       const AliESDVertex *vertex, 
172                       AliESDtrack *track);
173   Bool_t   IsInPhaseSpace(AliESDtrack *track);
174
175   void     FillQA(AliStack *stack,
176                   AliESDEvent *esd,
177                   const AliESDVertex *vertex,
178                   AliESDtrack *track);
179
180   void     InitVertexQA();
181   void     InitQA();
182   void     InitMCAnalysis();
183   void     InitCutLists();
184   void     InitEfficiencyAnalysis();
185
186   Bool_t   IsLabelUsed(TArrayI array, Int_t label);
187   Int_t    ConvertPDGToInt(Int_t pdgCode);
188   Float_t  GetSigmaToVertex(AliESDtrack* esdTrack) const; 
189   Double_t Rapidity(Double_t Px, Double_t Py, Double_t Pz) const;
190   
191   Bool_t fAnalysisEtaMode; //run the QA in eta or y
192
193   Int_t fNBinsY; //number of bins in eta or y
194   Float_t fMinY, fMaxY; //min & max value of eta or y
195   Int_t fNBinsPt;  //number of bins in pT
196   Float_t fMinPt, fMaxPt; //min & max value of pT
197   
198   //cuts
199   Int_t fMinTPCClusters, fMinITSClusters; //min TPC & ITS clusters
200   Double_t fMaxChi2PerTPCCluster, fMaxChi2PerITSCluster; //max chi2 per TPC & ITS cluster
201   Double_t fMaxCov11, fMaxCov22, fMaxCov33, fMaxCov44, fMaxCov55; //max values of cov. matrix
202   Double_t fMaxSigmaToVertex; //max sigma to vertex cut
203   Double_t fMaxSigmaToVertexTPC; //max sigma to vertex cut
204   Double_t fMaxDCAXY, fMaxDCAXYTPC; //max DCA xy
205   Double_t fMaxDCAZ, fMaxDCAZTPC; //max DCA z
206   Double_t fMaxConstrainChi2; //max constrain chi2 - vertex
207   Bool_t fMinTPCClustersFlag, fMinITSClustersFlag; //shows if this cut is used or not
208   Bool_t fMaxChi2PerTPCClusterFlag, fMaxChi2PerITSClusterFlag; //shows if this cut is used or not
209   Bool_t fMaxCov11Flag, fMaxCov22Flag, fMaxCov33Flag, fMaxCov44Flag, fMaxCov55Flag; //shows if this cut is used or not
210   Bool_t fMaxSigmaToVertexFlag; //shows if this cut is used or not
211   Bool_t fMaxSigmaToVertexTPCFlag; //shows if this cut is used or not
212   Bool_t fMaxDCAXYFlag, fMaxDCAXYTPCFlag; //shows if this cut is used or not
213   Bool_t fMaxDCAZFlag, fMaxDCAZTPCFlag; //shows if this cut is used or not
214   Bool_t fMaxConstrainChi2Flag; //shows if this cut is used or not
215   Bool_t fITSRefitFlag, fTPCRefitFlag; //shows if this cut is used or not
216   Bool_t fESDpidFlag, fTPCpidFlag; //shows if this cut is used or not
217   Bool_t fPointOnITSLayer1Flag, fPointOnITSLayer2Flag; //shows if this cut is used or not
218   Bool_t fPointOnITSLayer3Flag, fPointOnITSLayer4Flag; //shows if this cut is used or not
219   Bool_t fPointOnITSLayer5Flag, fPointOnITSLayer6Flag; //shows if this cut is used or not
220
221   //QA histograms
222   //Bool_t fQAHistograms; //Boolean to activate the QA histograms
223   TList *fGlobalQAList; //TList storing the directories for the QA histograms
224   TList *fQAVertexList; //TList storing the vertex QA plots
225   TList *fQA2DList; //TList storing the accepted primary/secondary (anti)protons
226   TList *fQAPrimaryProtonsAcceptedList; //list of the QA histos for accepted primary protons
227   TList *fQAPrimaryProtonsRejectedList; //list of the QA histos for rejected primary protons
228   TList *fQASecondaryProtonsAcceptedList; //list of the QA histos for accepted secondary protons
229   TList *fQASecondaryProtonsRejectedList; //list of the QA histos for rejected secondary protons
230   TList *fQAPrimaryAntiProtonsAcceptedList; //list of the QA histos for accepted primary antiprotons
231   TList *fQAPrimaryAntiProtonsRejectedList; //list of the QA histos for rejected primary antiprotons
232   TList *fQASecondaryAntiProtonsAcceptedList; //list of the QA histos for accepted secondary antiprotons
233   TList *fQASecondaryAntiProtonsRejectedList; //list of the QA histos for rejected secondary antiprotons
234
235   //pid
236   Bool_t fFunctionProbabilityFlag; //flag: kTRUE if functions used
237   Double_t fPartFrac[10]; //prior probabilities
238   TF1  *fElectronFunction; //momentum dependence of the prior probs
239   TF1  *fMuonFunction; //momentum dependence of the prior probs
240   TF1  *fPionFunction; //momentum dependence of the prior probs
241   TF1  *fKaonFunction; //momentum dependence of the prior probs
242   TF1  *fProtonFunction; //momentum dependence of the prior probs
243
244   //Detectors
245   Bool_t fUseTPCOnly; //kTRUE if TPC only information is used
246   Bool_t fUseHybridTPC; //kTRUE if TPC info is used for momentum - PID and ITS for vertex & points
247
248   //MC analysis
249   TList *fPDGList; //list with the 3D histograms: y-pt-pdg (anti)protons
250   TList *fMCProcessesList; //list with the MC processes for every secondary (anti)proton
251   Bool_t fRunMCAnalysis; //run this part or not
252   Bool_t fMCProcessIdFlag; //flag to see if we should check the process id
253   UInt_t fMCProcessId; //process id based on the TMCProcess
254   Bool_t fMotherParticlePDGCodeFlag; //flag to see if we should check the pdg code of the mother particle
255   Int_t  fMotherParticlePDGCode; //pdg code of the mother particle
256
257   TList *fAcceptedCutList;// list of the cut parameters' histograms
258   TList *fRejectedCutList;// list of the cut parameters' histograms
259   TList *fAcceptedDCAList;// list of the DCA histograms
260   TList *fRejectedDCAList;// list of the DCA histograms
261
262   //Efficiency (reconstruction & PID)
263   Bool_t fRunEfficiencyAnalysis; //run this part or not
264   Bool_t fUseCutsInEfficiency;//use the cuts in the reco and pid efficiency
265
266   TList *fEfficiencyList;// list of the efficiency histograms
267
268   ClassDef(AliProtonQAAnalysis,0);
269 };
270
271 #endif