AliPhysics  695988a (695988a)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliFlowEventCuts.h
Go to the documentation of this file.
1 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. */
2 /* See cxx source for full Copyright notice */
3 /* $Id$ */
4 
5 // AliFlowEventCuts:
6 // An event cut class
7 // origin: Mikolaj Krzewicki (mikolaj.krzewicki@cern.ch)
8 
9 #ifndef ALIFLOWEVENTCUTS_H
10 #define ALIFLOWEVENTCUTS_H
11 
12 #include <float.h>
13 #include <limits.h>
14 #include "TNamed.h"
15 
16 class AliVEvent;
17 class AliMCEvent;
18 class TBrowser;
19 class AliAnalysisUtils;
20 #include "TList.h"
21 #include "TH1.h"
22 #include "AliTriggerAnalysis.h"
23 #include "AliFlowTrackCuts.h"
24 #include "AliFlowEventSimpleCuts.h"
25 
27 
28  public:
30 
32  AliFlowEventCuts(const char* name, const char* title = "AliFlowEventCuts");
33  AliFlowEventCuts(const AliFlowEventCuts& someCuts);
35  virtual ~AliFlowEventCuts();
36 
37  virtual Bool_t IsSelected(TObject* obj, TObject *objmc);
38 
39  Bool_t PassesCuts(AliVEvent* event, AliMCEvent *mcevent);
40 
42 
46  void SetRefMultMax(Int_t value) {fRefMultMax=value;fCutRefMult=kTRUE;}
47  void SetRefMultMin(Int_t value) {fRefMultMin=value;fCutRefMult=kTRUE;}
48  void SetRefMultRange(Int_t min, Int_t max) {fRefMultMin=min;fRefMultMax=max;fCutRefMult=kTRUE;}
51  void SetImpactParameterRange(Double_t min, Double_t max) {fImpactParameterMin=min;fImpactParameterMax=max;fCutImpactParameter=kTRUE;}
52  void SetPrimaryVertexXrange(Double_t min, Double_t max)
54  void SetPrimaryVertexYrange(Double_t min, Double_t max)
56  void SetPrimaryVertexZrange(Double_t min, Double_t max)
58  void SetNContributorsRange(Int_t min, Int_t max=INT_MAX)
60  void SetMeanPtRange(Double_t min, Double_t max) {fCutMeanPt=kTRUE; fMeanPtMax=max; fMeanPtMin=min;}
62  void SetCutZDCtiming(Bool_t c=kTRUE) {fCutZDCtiming=c;}
63  void SetCutSPDTRKVtxZ(Bool_t b=kTRUE) {fCutSPDTRKVtxZ=b;}
66 
67  Int_t GetNumberOfTracksMax() const {return fNumberOfTracksMax;}
68  Int_t GetNumberOfTracksMin() const {return fNumberOfTracksMin;}
69  Int_t GetRefMultMax() const {return fRefMultMax;}
70  Int_t GetRefMultMin() const {return fRefMultMin;}
72  void SetRefMultMethod(AliESDtrackCuts::MultEstTrackType m) { fRefMultMethodAliESDtrackCuts=m;
75  void SetRefMultCuts( AliFlowTrackCuts* cuts ) {fRefMultCuts=static_cast<AliFlowTrackCuts*>(cuts->Clone());}
76  void SetMeanPtCuts( AliFlowTrackCuts* cuts ) {fMeanPtCuts=static_cast<AliFlowTrackCuts*>(cuts->Clone());}
78  void DefineHistograms();
79  void SetQA(Bool_t b=kTRUE) {if (b) DefineHistograms();}
80  TList* GetQA() const {return fQA;}
81  TH1* QAbefore(Int_t i) {return static_cast<TH1*>(static_cast<TList*>(fQA->At(0))->At(i));}
82  TH1* QAafter(Int_t i) {return static_cast<TH1*>(static_cast<TList*>(fQA->At(1))->At(i));}
83 
84  Int_t RefMult(AliVEvent* event, AliMCEvent *mcEvent = 0x0);
85  //Int_t GetRefMult() {return fRefMult;}
86  Int_t GetReferenceMultiplicity(AliVEvent* event, AliMCEvent *mcEvent) {return RefMult(event,mcEvent);}
87  const char* CentrMethName(refMultMethod method) const;
90 
91  Float_t GetCentrality(AliVEvent* event, AliMCEvent* mcEvent);
92  void SetUsedDataset(Bool_t b=kTRUE) {fData2011=b;} // confusing name, better use different interface
93  void SetLHC10h(Bool_t b=kTRUE) {fData2011=(!b);} // TODO let cut object determine runnumber and period
94  void SetLHC11h(Bool_t b=kTRUE) {fData2011=b;} // use this only as 'manual override'
95  void SetCheckPileup(Bool_t b=kFALSE) {fCheckPileUp = b;} // In case a pile-up rejection is required
96 
97  void Browse(TBrowser* b);
98  Long64_t Merge(TCollection* list);
100 
101 
102  private:
103  TList* fQA; //QA
104  Bool_t fCutNumberOfTracks;//cut on # of tracks
105  Int_t fNumberOfTracksMax; //limits
106  Int_t fNumberOfTracksMin; //limits
107  Bool_t fCutRefMult; //cut on refmult
108  refMultMethod fRefMultMethod; //how do we calculate refmult?
109  Bool_t fUseAliESDtrackCutsRefMult; //use AliESDtrackCuts for refmult calculation
110  AliESDtrackCuts::MultEstTrackType fRefMultMethodAliESDtrackCuts;
111  Int_t fRefMultMax; //max refmult
112  Int_t fRefMultMin; //min refmult
114  AliFlowTrackCuts* fMeanPtCuts; //mean pt cuts
115  AliFlowTrackCuts* fStandardTPCcuts; //Standard TPC cuts
116  AliFlowTrackCuts* fStandardGlobalCuts; //StandardGlobalCuts
117  AliAnalysisUtils* fUtils;
118  Bool_t fCutPrimaryVertexX; //cut on x of prim vtx
119  Double_t fPrimaryVertexXmax; //max x prim vtx
120  Double_t fPrimaryVertexXmin; //min x prim vtx
121  Bool_t fCutPrimaryVertexY; //cut on y of prim vtx
122  Double_t fPrimaryVertexYmax; //max y prim vtx
123  Double_t fPrimaryVertexYmin; //min y prim vtx
124  Bool_t fCutPrimaryVertexZ; //cut on z of prim vtx
125  Double_t fPrimaryVertexZmax; //max z prim vtx
126  Double_t fPrimaryVertexZmin; //min z prim vtx
127  Bool_t fCutNContributors; //cut on number of contributors
128  Int_t fNContributorsMax; //maximal number of contrib
129  Int_t fNContributorsMin; //minimal number of contrib
130  Bool_t fCutMeanPt; //cut on mean pt
131  Double_t fMeanPtMax; //max mean pt
132  Double_t fMeanPtMin; //min mean pt
133  Bool_t fCutSPDvertexerAnomaly; //cut on the spd vertexer anomaly
134  Bool_t fCutSPDTRKVtxZ; //require compatibility between SPDvertexz TRKvertexz
135  Bool_t fCutTPCmultiplicityOutliers; //cut TPC multiplicity outliers
136  Bool_t fCutTPCmultiplicityOutliersAOD; // cut TPC outliers in 10h or 11h aod
137  Bool_t fUseCentralityUnchecked; //use the unchecked method
138  refMultMethod fCentralityPercentileMethod; //where to get the percentile from
139  Bool_t fCutZDCtiming; //cut on ZDC timing
140  AliTriggerAnalysis fTrigAna; //trigger analysis object
141  Bool_t fCutImpactParameter; //cut on impact parameter (MC header)
142  Double_t fImpactParameterMin; // min impact parameter
143  Double_t fImpactParameterMax; // max impact parameter
145  Bool_t fData2011; //2011 data is used
146  Bool_t fCheckPileUp; //pile-up
147  ClassDef(AliFlowEventCuts,6)
148 };
149 
150 #endif
151 
152 
Double_t fPrimaryVertexZmax
AliESDtrackCuts::MultEstTrackType fRefMultMethodAliESDtrackCuts
Float_t GetCentrality(AliVEvent *event, AliMCEvent *mcEvent)
void SetMeanPtCuts(AliFlowTrackCuts *cuts)
Double_t fPrimaryVertexXmax
void SetRefMultRange(Int_t min, Int_t max)
Bool_t fData2011
correlation between TPCMult and GlobalMult
void SetPrimaryVertexYrange(Double_t min, Double_t max)
void SetRefMultMax(Int_t value)
const char * title
Definition: MakeQAPdf.C:26
Long64_t Merge(TCollection *list)
TH1 * QAbefore(Int_t i)
Bool_t fCutTPCmultiplicityOutliersAOD
Int_t GetReferenceMultiplicity(AliVEvent *event, AliMCEvent *mcEvent)
AliFlowTrackCuts * fStandardTPCcuts
AliFlowTrackCuts * fMeanPtCuts
Double_t fPrimaryVertexYmin
refMultMethod GetRefMultMethod() const
Double_t fPrimaryVertexYmax
TList * list
Int_t GetNumberOfTracksMin() const
TH1 * QAafter(Int_t i)
AliFlowEventCuts & operator=(const AliFlowEventCuts &someCuts)
void SetRefMultMethod(refMultMethod m)
AliAnalysisUtils * fUtils
const char * CentrMethName(refMultMethod method) const
AliFlowTrackCuts * GetRefMultCuts() const
refMultMethod fRefMultMethod
Bool_t PassesCuts(AliVEvent *event, AliMCEvent *mcevent)
void SetImpactParameterMax(Double_t value)
void SetImpactParameterRange(Double_t min, Double_t max)
void SetCutSPDvertexerAnomaly(Bool_t b=kTRUE)
Int_t RefMult(AliVEvent *event, AliMCEvent *mcEvent=0x0)
void SetNContributorsRange(Int_t min, Int_t max=INT_MAX)
void SetUsedDataset(Bool_t b=kTRUE)
void SetRefMultCuts(AliFlowTrackCuts *cuts)
refMultMethod fCentralityPercentileMethod
void SetPrimaryVertexZrange(Double_t min, Double_t max)
void SetRefMultMethod(AliESDtrackCuts::MultEstTrackType m)
Int_t GetRefMultMin() const
void Browse(TBrowser *b)
Bool_t fCutPrimaryVertexX
analysis utils object
void SetNumberOfTracksRange(Int_t min, Int_t max)
AliTriggerAnalysis fTrigAna
Double_t fPrimaryVertexZmin
void SetLHC11h(Bool_t b=kTRUE)
Bool_t fUseAliESDtrackCutsRefMult
void SetNumberOfTracksMax(Int_t value)
Double_t fImpactParameterMax
void SetUseCentralityUnchecked(Bool_t b=kTRUE)
Int_t GetRefMultMax() const
AliFlowTrackCuts * fStandardGlobalCuts
void SetRefMultMin(Int_t value)
Bool_t fUseCentralityUnchecked
AliFlowTrackCuts * fRefMultCuts
TList * GetQA() const
void SetCutTPCmultiplicityOutliersAOD(Bool_t b=kTRUE)
static AliFlowEventCuts * StandardCuts()
Double_t fImpactParameterMin
void SetCentralityPercentileMethod(refMultMethod m)
TH2F * GetCorrelationTPCvsGlobalMultiplicity()
void SetCutTPCmultiplicityOutliers(Bool_t b=kTRUE)
void SetImpactParameterMin(Double_t value)
Double_t fPrimaryVertexXmin
virtual Bool_t IsSelected(TObject *obj, TObject *objmc)
void SetMeanPtRange(Double_t min, Double_t max)
void SetCutZDCtiming(Bool_t c=kTRUE)
Int_t GetNumberOfTracksMax() const
void SetCutSPDTRKVtxZ(Bool_t b=kTRUE)
void SetLHC10h(Bool_t b=kTRUE)
void SetPrimaryVertexXrange(Double_t min, Double_t max)
Bool_t fCutTPCmultiplicityOutliers
void SetCheckPileup(Bool_t b=kFALSE)
void SetNumberOfTracksMin(Int_t value)
void SetQA(Bool_t b=kTRUE)