]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/GammaConv/AliCaloPhotonCuts.h
- added new task for pi0 reconstruction using purely calorimeter clusters
[u/mrichter/AliRoot.git] / PWGGA / GammaConv / AliCaloPhotonCuts.h
1 #ifndef ALICALOPHOTONCUTS_H
2 #define ALICALOPHOTONCUTS_H
3
4 // Class handling all kinds of selection cuts for Gamma Conversion analysis
5 // Authors: Svein Lindal, Daniel Lohner                                    *
6
7 #include "AliConversionPhotonBase.h"
8 #include "AliAODConversionMother.h"
9 #include "AliAODTrack.h"
10 #include "AliESDtrack.h"
11 #include "AliVTrack.h"
12 #include "AliVCluster.h"
13 #include "AliAODTrack.h"
14 #include "AliStack.h"
15 #include "AliAnalysisCuts.h"
16 #include "TH1F.h"
17 #include "TF1.h"
18 #include "AliAnalysisUtils.h"
19 #include "AliAnalysisManager.h"
20
21 class AliESDEvent;
22 class AliAODEvent;
23 class AliConversionPhotonBase;
24 class TH1F;
25 class TH2F;
26 class TF1;
27 class AliPIDResponse;
28 class AliAnalysisCuts;
29 class iostream;
30 class TList;
31 class AliAnalysisManager;
32 class AliAODMCParticle;
33
34 using namespace std;
35
36 class AliCaloPhotonCuts : public AliAnalysisCuts {
37                 
38         public: 
39                 enum cutIds {
40                         kClusterType,                  
41                         kEtaMin,
42                         kEtaMax,
43                         kPhiMin,
44                         kPhiMax,
45                         kDistanceToBadChannel,
46                         kTiming,
47                         kTrackMatching,
48                         kExoticCell,
49                         kMinEnery,               
50                         kNMinCells, 
51                         kMinM02,    
52                         kMaxM02,    
53                         kMinM20,
54                         kMaxM20,
55                         kDispersion,
56                         kNLM,
57                         kNCuts
58                 };
59
60                 enum photonCuts {
61                         kPhotonIn=0,
62                         kDetector,
63                         kAcceptance,
64                         kClusterQuality,
65                         kPhotonOut
66                 };
67
68                 //handeling of CutString
69                 static const char * fgkCutNames[kNCuts];
70                 Bool_t SetCutIds(TString cutString); 
71                 Int_t fCuts[kNCuts];
72                 Bool_t SetCut(cutIds cutID, Int_t cut);
73                 Bool_t UpdateCutString();
74                 void PrintCuts();
75                 void PrintCutsWithValues();
76                 
77                 Bool_t InitializeCutsFromCutString(const TString analysisCutSelection);
78                 TString GetCutNumber();
79                 
80                 //Constructors
81                 AliCaloPhotonCuts(const char *name="ClusterCuts", const char * title="Cluster Cuts");
82                 AliCaloPhotonCuts(const AliCaloPhotonCuts&);
83                 AliCaloPhotonCuts& operator=(const AliCaloPhotonCuts&);
84
85                 //virtual destructor
86                 virtual ~AliCaloPhotonCuts();                          
87
88                 virtual Bool_t IsSelected(TObject* /*obj*/){return kTRUE;}
89                 virtual Bool_t IsSelected(TList* /*list*/) {return kTRUE;}
90
91                 Bool_t ClusterIsSelected(AliVCluster* cluster, AliVEvent *event, Bool_t isMC);
92                 Bool_t ClusterIsSelectedMC(TParticle *particle,AliStack *fMCStack);
93                 Bool_t ClusterIsSelectedAODMC(AliAODMCParticle *particle,TClonesArray *aodmcArray);
94                         
95                 void InitCutHistograms(TString name="");
96                 void SetFillCutHistograms(TString name=""){if(!fHistograms){InitCutHistograms(name);};}
97                 TList *GetCutHistograms(){return fHistograms;}
98                 void FillClusterCutIndex(Int_t photoncut){if(fHistCutIndex)fHistCutIndex->Fill(photoncut);}
99                         
100                 ///Cut functions
101                 Bool_t AcceptanceCuts(AliVCluster* cluster, AliVEvent *event);
102                 Bool_t ClusterQualityCuts(AliVCluster* cluster,AliVEvent *event, Bool_t isMC);
103
104                 Bool_t MatchConvPhotonToCluster(AliAODConversionPhoton* convPhoton, AliVCluster* cluster, AliVEvent* event );
105
106                 // Set Individual Cuts
107                 Bool_t SetClusterTypeCut(Int_t);
108                 Bool_t SetMinEtaCut(Int_t);
109                 Bool_t SetMaxEtaCut(Int_t);
110                 Bool_t SetMinPhiCut(Int_t);             
111                 Bool_t SetMaxPhiCut(Int_t);             
112                 Bool_t SetDistanceToBadChannelCut(Int_t);
113                 Bool_t SetTimingCut(Int_t);
114                 Bool_t SetTrackMatchingCut(Int_t);
115                 Bool_t SetExoticCellCut(Int_t);
116                 Bool_t SetMinEnergyCut(Int_t);
117                 Bool_t SetMinNCellsCut(Int_t);
118                 Bool_t SetMaxM02(Int_t);
119                 Bool_t SetMinM02(Int_t);
120                 Bool_t SetMaxM20(Int_t);
121                 Bool_t SetMinM20(Int_t);
122                 Bool_t SetDispersion(Int_t);
123                 Bool_t SetNLM(Int_t);
124                 
125         protected:
126                 TList *fHistograms;
127                 
128                 //cuts
129                 Int_t   fClusterType;                                           // which cluster do we have
130                 Double_t fMinEtaCut;                                            // min eta cut
131                 Double_t fMaxEtaCut;                                            // max eta cut
132                 Bool_t fUseEtaCut;                                                      // flag for switching on eta cut
133                 Double_t fMinPhiCut;                                            // phi cut
134                 Double_t fMaxPhiCut;                                            // phi cut
135                 Bool_t fUsePhiCut;                                                      // flag for switching on phi cut
136                 Double_t fMinDistanceToBadChannel;                      // minimum distance to bad channel
137                 Bool_t fUseDistanceToBadChannel;                        // flag for switching on distance to bad channel cut
138                 Double_t fMaxTimeDiff;                                          // maximum time difference to triggered collision
139                 Bool_t fUseTimeDiff;                                            // flag for switching on time difference cut
140                 Double_t fMinDistTrackToCluster;                        // minimum distance between track and cluster
141                 Bool_t fUseDistTrackToCluster;                          // flag for switching on distance between track and cluster cut
142                 Double_t fExoticCell;                                           // exotic cell cut
143                 Bool_t fUseExoticCell;                                          // flag for switching on exotic cell cut
144                 Double_t fMinEnergy;                                            // minium energy per cluster
145                 Bool_t fUseMinEnergy;                                           // flag for switching on minimum energy cut
146                 Int_t fMinNCells;                                                       // minimum number of cells 
147                 Bool_t fUseNCells;                                                      // flag for switching on minimum N Cells cut
148                 Double_t fMaxM02;                                                       // maximum M02
149                 Double_t fMinM02;                                                       // minimum M02
150                 Bool_t fUseM02;                                                         // flag for switching on M02 cut
151                 Double_t fMaxM20;                                                       // maximum M20
152                 Double_t fMinM20;                                                       // minimum M20
153                 Bool_t fUseM20;                                                         // flag for switching on M20 cut
154                 Double_t fMaxDispersion;                                        // maximum dispersion
155                 Bool_t fUseDispersion;                                          // flag for switching on dispersion cut
156                 Int_t fMinNLM;                                                          // minimum number of local maxima in cluster
157                 Int_t fMaxNLM;                                                          // maximum number of local maxima in cluster
158                 Bool_t fUseNLM;                                                         // flag for switching on NLM cut
159                 
160                 // CutString
161                 TObjString *fCutString;                                         // cut number used for analysis
162                 
163                 // Histograms
164                 TH1F *fHistCutIndex;                                            // bookkeeping for cuts
165                 TH1F *fHistAcceptanceCuts;                                      // bookkeeping for acceptance cuts
166                 TH1F *fHistClusterIdentificationCuts;           // bookkeeping for cluster identification cuts
167                 
168                 TH2F* fHistClusterEtavsPhiBeforeAcc;            // eta-phi-distribution before acceptance cuts
169                 TH2F* fHistClusterEtavsPhiAfterAcc;             // eta-phi-distribution of all after acceptance cuts
170                 TH2F* fHistClusterEtavsPhiAfterQA;                      // eta-phi-distribution of all after cluster quality cuts
171                 TH1F* fHistDistanceToBadChannelBeforeAcc;   // distance to bad channel before acceptance cuts
172                 TH1F* fHistDistanceToBadChannelAfterAcc;        // distance to bad channel after acceptance cuts
173                 TH2F* fHistClusterTimevsEBeforeQA;                      // Cluster time vs E before cluster quality cuts
174                 TH2F* fHistClusterTimevsEAfterQA;                       // Cluster time vs E after cluster quality cuts
175                 TH2F* fHistExoticCellBeforeQA;                          // Exotic cell: 1-Ecross/E cell vs Ecluster before acceptance cuts
176                 TH2F* fHistExoticCellAfterQA;                           // Exotic cell: 1-Ecross/E cell vs Ecluster after cluster quality cuts
177                 TH1F* fHistNMatchedTracks;                                      // number of matched tracks
178                 TH1F* fHistDistanceTrackToClusterBeforeQA;      // distance cluster to track before acceptance cuts
179                 TH1F* fHistDistanceTrackToClusterAfterQA;       // distance cluster to track after cluster quality cuts
180                 TH1F* fHistEnergyOfClusterBeforeQA;                     // enery per cluster before acceptance cuts
181                 TH1F* fHistEnergyOfClusterAfterQA;                      // enery per cluster after cluster quality cuts
182                 TH1F* fHistNCellsBeforeQA;                                      // number of cells per cluster before acceptance cuts
183                 TH1F* fHistNCellsAfterQA;                                       // number of cells per cluster after cluster quality cuts
184                 TH1F* fHistM02BeforeQA;                                         // M02 before acceptance cuts
185                 TH1F* fHistM02AfterQA;                                          // M02 after cluster quality cuts
186                 TH1F* fHistM20BeforeQA;                                         // M20 before acceptance cuts
187                 TH1F* fHistM20AfterQA;                                          // M20 after cluster quality cuts
188                 TH1F* fHistDispersionBeforeQA;                          // dispersion before acceptance cuts
189                 TH1F* fHistDispersionAfterQA;                           // dispersion after cluster quality cuts
190                 TH1F* fHistNLMBeforeQA;                                         // number of local maxima in cluster before acceptance cuts
191                 TH1F* fHistNLMAfterQA;                                          // number of local maxima in cluster after cluster quality cuts
192                                 
193         private:
194
195                 ClassDef(AliCaloPhotonCuts,1)
196 };
197
198 #endif