]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - PWGGA/GammaConv/AliConversionCuts.h
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGGA / GammaConv / AliConversionCuts.h
... / ...
CommitLineData
1#ifndef ALICONVERSIONCUTS_H
2#define ALICONVERSIONCUTS_H
3
4// Class handling all kinds of selection cuts for Gamma Conversion analysis
5// Authors: Svein Lindal, Daniel Lohner *
6
7#include "AliAODpidUtil.h"
8#include "AliConversionPhotonBase.h"
9#include "AliAODConversionMother.h"
10#include "AliAODTrack.h"
11#include "AliESDtrack.h"
12#include "AliVTrack.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
21class AliESDEvent;
22class AliAODEvent;
23class AliConversionPhotonBase;
24class AliKFVertex;
25class TH1F;
26class TH2F;
27class TF1;
28class AliPIDResponse;
29class AliAnalysisCuts;
30class iostream;
31class TList;
32class AliAnalysisManager;
33class AliAODMCParticle;
34
35using namespace std;
36
37class AliConversionCuts : public AliAnalysisCuts {
38
39 public:
40
41
42 enum cutIds {
43 kisHeavyIon,
44 kCentralityMin,
45 kCentralityMax,
46 kSelectSpecialTriggerAlias,
47 kSelectSubTriggerClass,
48 kremovePileUp,
49 kExtraSignals,
50 kv0FinderType,
51 ketaCut,
52 kRCut,
53 ksinglePtCut,
54 kclsTPCCut,
55 kededxSigmaCut,
56 kpidedxSigmaCut,
57 kpiMomdedxSigmaCut,
58 kpiMaxMomdedxSigmaCut,
59 kLowPRejectionSigmaCut,
60 kTOFelectronPID,
61 kQtMaxCut,
62 kchi2GammaCut,
63 kPsiPair,
64 kdoPhotonAsymmetryCut,
65 kCosPAngle,
66 kElecShare,
67 kToCloseV0s,
68 kDcaRPrimVtx,
69 kDcaZPrimVtx,
70 kInPlaneOutOfPlane,
71 kNCuts
72 };
73
74 enum photonCuts {
75 kPhotonIn=0,
76 kOnFly,
77 kNoTracks,
78 kTrackCuts,
79 kdEdxCuts,
80 kConvPointFail,
81 kPhotonCuts,
82 kEventPlane,
83 kPhotonOut
84 };
85
86
87 Bool_t SetCutIds(TString cutString);
88 Int_t fCuts[kNCuts];
89 Bool_t SetCut(cutIds cutID, Int_t cut);
90 Bool_t UpdateCutString();
91
92
93 static const char * fgkCutNames[kNCuts];
94
95 Double_t GetCosineOfPointingAngle(const AliConversionPhotonBase * photon, AliVEvent * event) const;
96
97
98 Bool_t InitializeCutsFromCutString(const TString analysisCutSelection);
99 void SelectCollisionCandidates(UInt_t offlineTriggerMask = AliVEvent::kAny) {
100 fOfflineTriggerMask = offlineTriggerMask;
101 fTriggerSelectedManually = kTRUE;
102 }
103 void SelectSpecialTrigger(UInt_t offlineTriggerMask = AliVEvent::kAny, TString TriggerClassName = "AliVEvent::kAny" ) {
104 fOfflineTriggerMask = offlineTriggerMask;
105 fSpecialTriggerName = TriggerClassName;
106 cout << fSpecialTriggerName.Data() << endl;
107
108 }
109 void FillElectonLabelArray(AliAODConversionPhoton* photon, Int_t nV0);
110 void SetAcceptedHeader(TList *HeaderList){fHeaderList = HeaderList;}
111 void SetPreSelectionCutFlag(Bool_t preSelFlag){fPreSelCut = preSelFlag;}
112 TString *GetFoundHeader(){return fGeneratorNames;}
113
114 Int_t GetEventQuality(){return fEventQuality;}
115 Bool_t GetIsFromPileup(){return fRemovePileUp;}
116
117 AliConversionCuts(const char *name="V0Cuts", const char * title="V0 Cuts");
118 AliConversionCuts(const AliConversionCuts&);
119 AliConversionCuts& operator=(const AliConversionCuts&);
120
121 virtual ~AliConversionCuts(); //virtual destructor
122
123 static AliConversionCuts * GetStandardCuts2010PbPb();
124 static AliConversionCuts * GetStandardCuts2010pp();
125
126 virtual Bool_t IsSelected(TObject* /*obj*/){return kTRUE;}
127 virtual Bool_t IsSelected(TList* /*list*/) {return kTRUE;}
128
129 TString GetCutNumber();
130
131 void GetCentralityRange(Double_t range[2]){range[0]=10*fCentralityMin;range[1]=10*fCentralityMax;};
132
133 // Cut Selection
134 Bool_t EventIsSelected(AliVEvent *fInputEvent, AliVEvent *fMCEvent);
135 Int_t IsEventAcceptedByConversionCut(AliConversionCuts *ReaderCuts, AliVEvent *InputEvent, AliMCEvent *MCEvent, Int_t isHeavyIon);
136 Bool_t PhotonIsSelected(AliConversionPhotonBase * photon, AliVEvent * event);
137 Bool_t PhotonIsSelectedMC(TParticle *particle,AliStack *fMCStack,Bool_t checkForConvertedGamma=kTRUE);
138 Bool_t PhotonIsSelectedAODMC(AliAODMCParticle *particle,TClonesArray *aodmcArray,Bool_t checkForConvertedGamma=kTRUE);
139 Bool_t ElectronIsSelectedMC(TParticle *particle,AliStack *fMCStack);
140 Bool_t TracksAreSelected(AliVTrack * negTrack, AliVTrack * posTrack);
141 Bool_t MesonIsSelected(AliAODConversionMother *pi0,Bool_t IsSignal=kTRUE);
142 Bool_t MesonIsSelectedMC(TParticle *fMCMother,AliStack *fMCStack, Bool_t bMCDaughtersInAcceptance=kFALSE);
143
144 void InitAODpidUtil(Int_t type);
145 Bool_t InitPIDResponse();
146
147 void SetPIDResponse(AliPIDResponse * pidResponse) {fPIDResponse = pidResponse;}
148 AliPIDResponse * GetPIDResponse() { return fPIDResponse;}
149
150 void PrintCuts();
151 void PrintCutsWithValues();
152
153 void InitCutHistograms(TString name="",Bool_t preCut = kTRUE);
154 void SetFillCutHistograms(TString name="",Bool_t preCut = kTRUE){if(!fHistograms){InitCutHistograms(name,preCut);};}
155 TList *GetCutHistograms(){return fHistograms;}
156 void FillPhotonCutIndex(Int_t photoncut){if(hCutIndex)hCutIndex->Fill(photoncut);}
157 void FillV0EtaBeforedEdxCuts(Float_t v0Eta){if(hEtaDistV0s)hEtaDistV0s->Fill(v0Eta);}
158 void FillV0EtaAfterdEdxCuts(Float_t v0Eta){if(hEtaDistV0sAfterdEdxCuts)hEtaDistV0sAfterdEdxCuts->Fill(v0Eta);}
159 void SetEtaShift(Double_t etaShift) {
160 fEtaShift = etaShift;
161 }
162 void SetEtaShift(TString pPbOrPbp) {
163 Double_t etaShift = 0.0;
164 if(!pPbOrPbp.CompareTo("pPb")) etaShift = -0.465;
165 else if(!pPbOrPbp.CompareTo("Pbp")) etaShift = 0.465;
166
167 fEtaShift = etaShift;
168 }
169 Double_t GetEtaShift() {return fEtaShift;}
170 Bool_t GetDoEtaShift(){return fDoEtaShift;}
171 void DoEtaShift(Bool_t doEtaShift){fDoEtaShift = doEtaShift;}
172 void GetCorrectEtaShiftFromPeriod(TString periodName);
173
174 static AliVTrack * GetTrack(AliVEvent * event, Int_t label);
175 static AliESDtrack *GetESDTrack(AliESDEvent * event, Int_t label);
176
177 ///Cut functions
178 Bool_t SpecificTrackCuts(AliAODTrack * negTrack, AliAODTrack * posTrack,Int_t &cutIndex);
179 Bool_t SpecificTrackCuts(AliESDtrack * negTrack, AliESDtrack * posTrack,Int_t &cutIndex);
180 Bool_t AcceptanceCuts(AliConversionPhotonBase *photon);
181 Bool_t AcceptanceCut(TParticle *particle, TParticle * ePos,TParticle* eNeg);
182 Bool_t dEdxCuts(AliVTrack * track);
183 Bool_t ArmenterosQtCut(AliConversionPhotonBase *photon);
184 Bool_t AsymmetryCut(AliConversionPhotonBase *photon,AliVEvent *event);
185 Bool_t PIDProbabilityCut(AliConversionPhotonBase *photon, AliVEvent * event);
186 Bool_t SelectV0Finder(Bool_t onfly){
187 if(onfly == fUseOnFlyV0Finder) return kTRUE;
188 else return kFALSE;
189 }
190 Bool_t PhotonCuts(AliConversionPhotonBase *photon,AliVEvent *event);
191 Bool_t CorrectedTPCClusterCut(AliConversionPhotonBase *photon, AliVEvent * event);
192 Bool_t PsiPairCut(const AliConversionPhotonBase * photon) const;
193 Bool_t CosinePAngleCut(const AliConversionPhotonBase * photon, AliVEvent * event) const;
194 Bool_t RejectSharedElectronV0s(AliAODConversionPhoton* photon, Int_t nV0, Int_t nV0s);
195 Bool_t RejectToCloseV0s(AliAODConversionPhoton* photon, TList *photons, Int_t nV0);
196 Int_t IsParticleFromBGEvent(Int_t index, AliStack *MCStack, AliVEvent *InputEvent = 0x0);
197 void GetNotRejectedParticles(Int_t rejection, TList *HeaderList, AliVEvent *MCEvent);
198 void SetUseReweightingWithHistogramFromFile( Bool_t pi0reweight=kTRUE, Bool_t etareweight=kFALSE, Bool_t k0sreweight=kFALSE, TString path="$ALICE_ROOT/PWGGA/GammaConv/MCSpectraInput.root",
199 TString histoNamePi0 = "", TString histoNameEta = "", TString histoNameK0s = "",
200 TString fitNamePi0 = "", TString fitNameEta = "", TString fitNameK0s ="" ) {
201 AliInfo(Form("enabled reweighting for: pi0 : %i, eta: %i, K0s: %i",pi0reweight, etareweight, k0sreweight));
202 fDoReweightHistoMCPi0 = pi0reweight;
203 fDoReweightHistoMCEta = etareweight;
204 fDoReweightHistoMCK0s = k0sreweight;
205 fPathTrFReweighting=path;
206 fNameHistoReweightingPi0 =histoNamePi0;
207 fNameHistoReweightingEta =histoNameEta;
208 fNameHistoReweightingK0s =histoNameK0s;
209 fNameFitDataPi0 =fitNamePi0;
210 fNameFitDataEta =fitNameEta;
211 fNameFitDataK0s =fitNameK0s;
212
213 }
214 void LoadReweightingHistosMCFromFile ();
215 UChar_t DeterminePhotonQualityAOD(AliAODConversionPhoton*, AliVEvent*);
216 // Event Cuts
217 Bool_t IsCentralitySelected(AliVEvent *fInputEvent, AliVEvent *fMCEvent = NULL);
218 Double_t GetCentrality(AliVEvent *event);
219 Int_t GetNumberOfContributorsVtx(AliVEvent *event);
220 Bool_t VertexZCut(AliVEvent *fInputEvent);
221 Bool_t IsTriggerSelected(AliVEvent *fInputEvent);
222 Bool_t HasV0AND(){return fHasV0AND;}
223 Bool_t IsSDDFired(){return fIsSDDFired;}
224 Int_t IsSpecialTrigger(){return fSpecialTrigger;}
225 TString GetSpecialTriggerName(){return fSpecialTriggerName;}
226 Bool_t InPlaneOutOfPlaneCut(Double_t photonPhi, Double_t eventPlaneAngle = -100, Bool_t fill = kTRUE);
227 Int_t GetInPlaneOutOfPlaneCut(){return fInPlaneOutOfPlane;}
228
229
230 // Set Individual Cuts
231 Bool_t SetRCut(Int_t RCut);
232 Bool_t SetV0Finder(Int_t v0FinderType);
233 Bool_t SetChi2GammaCut(Int_t chi2GammaCut);
234 Bool_t SetTPCdEdxCutPionLine(Int_t pidedxSigmaCut);
235 Bool_t SetTPCdEdxCutElectronLine(Int_t ededxSigmaCut);
236 Bool_t SetSinglePtCut(Int_t singlePtCut);
237 Bool_t SetTPCClusterCut(Int_t clsTPCCut);
238 Bool_t SetEtaCut(Int_t etaCut);
239 Bool_t SetMinMomPiondEdxCut(Int_t piMinMomdedxSigmaCut);
240 Bool_t SetMaxMomPiondEdxCut(Int_t piMaxMomdedxSigmaCut);
241 Bool_t SetLowPRejectionCuts(Int_t LowPRejectionSigmaCut);
242 Bool_t SetQtMaxCut(Int_t QtMaxCut);
243 Bool_t SetTOFElectronPIDCut(Int_t TOFelectronPID);
244 Bool_t SetTRDElectronCut(Int_t TRDElectronCut);
245 Bool_t SetCentralityMin(Int_t useCentrality);
246 Bool_t SetIsHeavyIon(Int_t isHeavyIon);
247 Bool_t SetCentralityMax(Int_t centralityBin);
248 Bool_t SetPhotonAsymmetryCut(Int_t doPhotonAsymmetryCut);
249 Bool_t SetRemovePileUp(Int_t removePileUp);
250 Bool_t SetMultiplicityMethod(Int_t multiplicityMethod);
251 Bool_t SetSelectSpecialTrigger(Int_t selectSpecialTrigger);
252 Bool_t SetSelectSubTriggerClass (Int_t selectSpecialSubTriggerClass);
253 Bool_t SetCosPAngleCut(Int_t cosCut);
254 Bool_t SetPsiPairCut(Int_t psiCut);
255 Bool_t SetSharedElectronCut(Int_t sharedElec);
256 Bool_t SetToCloseV0sCut(Int_t toClose);
257 Bool_t SetRejectExtraSignalsCut(Int_t extraSignal);
258 Bool_t SetDCARPhotonPrimVtxCut(Int_t DCARPhotonPrimVtx);
259 Bool_t SetDCAZPhotonPrimVtxCut(Int_t DCAZPhotonPrimVtx);
260 Bool_t SetInPlaneOutOfPlane(Int_t inOutPlane);
261 void SetAddedSignalPDGCode(Int_t addedSignalPDGcode) {fAddedSignalPDGCode = addedSignalPDGcode;}
262 // Request Flags
263
264 Int_t IsHeavyIon(){return fIsHeavyIon;}
265 Int_t GetFirstTPCRow(Double_t radius);
266 Float_t GetWeightForMeson(TString period, Int_t index, AliStack *MCStack, AliVEvent *InputEvent = 0x0);
267
268 Bool_t UseElecSharingCut(){return fDoSharedElecCut;}
269 Bool_t UseToCloseV0sCut(){return fDoToCloseV0sCut;}
270 Int_t GetMultiplicityMethod(){return fMultiplicityMethod;}
271 Double_t GetEtaCut(){return fEtaCut;}
272 Int_t GetSignalRejection(){return fRejectExtraSignals;}
273 Int_t GetNAcceptedHeaders(){return fnHeaders; }
274 TString * GetAcceptedHeaderNames(){return fGeneratorNames;}
275 Int_t * GetAcceptedHeaderStart(){return fNotRejectedStart;}
276 Int_t * GetAcceptedHeaderEnd(){return fNotRejectedEnd;}
277 TList* GetAcceptedHeader(){return fHeaderList;}
278
279
280 protected:
281 TList *fHistograms;
282 TList *fHeaderList;
283 AliPIDResponse *fPIDResponse;
284
285
286 Int_t fEventQuality; // EventQuality
287 //cuts
288 Double_t fMaxR; //r cut
289 Double_t fMinR; //r cut
290 Double_t fEtaCut; //eta cut
291 Double_t fEtaCutMin; //eta cut
292 Double_t fPtCut; // pt cut
293 Double_t fSinglePtCut; // pt cut for electron/positron
294 Double_t fMaxZ; //z cut
295 Double_t fMinClsTPC; // minimum clusters in the TPC
296 Double_t fMinClsTPCToF; // minimum clusters to findable clusters
297 Double_t fLineCutZRSlope; //linecut
298 Double_t fLineCutZValue; //linecut
299 Double_t fLineCutZRSlopeMin; //linecut
300 Double_t fLineCutZValueMin; //linecut
301 Double_t fChi2CutConversion; //chi2cut
302 Double_t fPIDProbabilityCutNegativeParticle;
303 Double_t fPIDProbabilityCutPositiveParticle;
304 Bool_t fDodEdxSigmaCut; // flag to use the dEdxCut based on sigmas
305 Bool_t fDoTOFsigmaCut; // flag to use TOF pid cut RRnewTOF
306 Double_t fPIDTRDEfficiency; // required electron efficiency for TRD PID
307 Bool_t fDoTRDPID; // flag to use TRD pid
308 Double_t fPIDnSigmaAboveElectronLine; // sigma cut
309 Double_t fPIDnSigmaBelowElectronLine; // sigma cut
310 Double_t fTofPIDnSigmaAboveElectronLine; // sigma cut RRnewTOF
311 Double_t fTofPIDnSigmaBelowElectronLine; // sigma cut RRnewTOF
312 Double_t fPIDnSigmaAbovePionLine; // sigma cut
313 Double_t fPIDnSigmaAbovePionLineHighPt; // sigma cut
314 Double_t fPIDMinPnSigmaAbovePionLine; // sigma cut
315 Double_t fPIDMaxPnSigmaAbovePionLine; // sigma cut
316 Double_t fDoKaonRejectionLowP; // Kaon rejection at low p
317 Double_t fDoProtonRejectionLowP; // Proton rejection at low p
318 Double_t fDoPionRejectionLowP; // Pion rejection at low p
319 Double_t fPIDnSigmaAtLowPAroundKaonLine; // sigma cut
320 Double_t fPIDnSigmaAtLowPAroundProtonLine; // sigma cut
321 Double_t fPIDnSigmaAtLowPAroundPionLine; // sigma cut
322 Double_t fPIDMinPKaonRejectionLowP; // Momentum limit to apply kaon rejection
323 Double_t fPIDMinPProtonRejectionLowP; // Momentum limit to apply proton rejection
324 Double_t fPIDMinPPionRejectionLowP; // Momentum limit to apply proton rejection
325 Bool_t fDoQtGammaSelection; // Select gammas using qtMax
326 Bool_t fDo2DQt; // Select gammas using ellipse cut
327 Double_t fQtMax; // Maximum Qt from Armenteros to select Gammas
328 Double_t fXVertexCut; //vertex cut
329 Double_t fYVertexCut; //vertex cut
330 Double_t fZVertexCut; // vertexcut
331 Double_t fNSigmaMass; //nsigma cut
332 Bool_t fUseEtaMinCut; //flag
333 Bool_t fUseOnFlyV0Finder; //flag
334 Bool_t fDoPhotonAsymmetryCut; // flag to use the PhotonAsymetryCut
335 Double_t fMinPPhotonAsymmetryCut; // Min Momentum for Asymmetry Cut
336 Double_t fMinPhotonAsymmetry; // Asymmetry Cut
337 Int_t fIsHeavyIon; // flag for heavy ion
338 Int_t fDetectorCentrality; // centrality detecotor V0M or CL1
339 Int_t fModCentralityClass; // allows to select smaller centrality classes
340 Double_t fMaxVertexZ; // max z offset of vertex
341 Int_t fCentralityMin; // centrality selection lower bin value
342 Int_t fCentralityMax; // centrality selection upper bin value
343 Bool_t fUseCorrectedTPCClsInfo; // flag to use corrected tpc cl info
344 Bool_t fUseTOFpid; // flag to use tof pid
345 Int_t fMultiplicityMethod; // selected multiplicity method
346 Int_t fSpecialTrigger; // flag
347 Int_t fSpecialSubTrigger; // flag
348 Bool_t fRemovePileUp; //flag
349 Float_t fOpeningAngle; // min opening angle for meson
350 Float_t fPsiPairCut;
351 Bool_t fDo2DPsiPairChi2;
352 Float_t fCosPAngleCut;
353 Bool_t fDoToCloseV0sCut; //
354 Int_t fRejectExtraSignals;//
355 Double_t fminV0Dist; //
356 Bool_t fDoSharedElecCut; //
357 Bool_t fDoPhotonQualitySelectionCut; //
358 Int_t fPhotonQualityCut; //
359 UInt_t fOfflineTriggerMask; // Task processes collision candidates only
360 Bool_t fHasV0AND; // V0AND Offline Trigger
361 Bool_t fIsSDDFired; // SDD FIRED to select with SDD events
362 TRandom3 fRandom; //
363 Int_t fElectronArraySize; // Size of electron array
364 Int_t *fElectronLabelArray; //[fElectronArraySize]
365 Double_t fDCAZPrimVtxCut; // cut value for the maximum distance in Z between the photon & the primary vertex [cm]
366 Double_t fDCARPrimVtxCut; // cut value for the maximum distance in R between the photon & the primary vertex [cm]
367 Int_t fInPlaneOutOfPlane; // In-Plane Out-Of Plane Analysis
368 Float_t fConversionPointXArray; // Array with conversion Point x
369 Float_t fConversionPointYArray; // Array with conversion Point y
370 Float_t fConversionPointZArray; // Array with conversion Point z
371 Int_t fnHeaders; // Number of Headers
372 Int_t *fNotRejectedStart; //[fnHeaders]
373 Int_t *fNotRejectedEnd; //[fnHeaders]
374 TString *fGeneratorNames; //[fnHeaders]
375 TObjString *fCutString; // cut number used for analysis
376 AliAnalysisUtils *fUtils;
377 Double_t fEtaShift;
378 Bool_t fDoEtaShift; // Flag for Etashift
379 Bool_t fDoReweightHistoMCPi0; // Flag for reweighting Pi0 input with histogram
380 Bool_t fDoReweightHistoMCEta; // Flag for reweighting Eta input with histogram
381 Bool_t fDoReweightHistoMCK0s; // Flag for reweighting K0s input with histogram
382 TString fPathTrFReweighting; // Path for file used in reweighting
383 TString fNameHistoReweightingPi0; //Histogram name for reweighting Pi0
384 TString fNameHistoReweightingEta; //Histogram name for reweighting Eta
385 TString fNameHistoReweightingK0s; //Histogram name for reweighting K0s
386 TString fNameFitDataPi0; //Fit name for fit to spectrum of pi0s in Data
387 TString fNameFitDataEta; //Fit name for fit to spectrum of etas in Data
388 TString fNameFitDataK0s; //Fit name for fit to spectrum of k0s in Data
389 // Histograms
390 TH1F* hEtaDistV0s; //eta-distribution of all V0s after Finder selection
391 TH1F* hEtaDistV0sAfterdEdxCuts; //eta-distribution of all V0s after Finder selection after dEdx cuts
392 TH1F *hdEdxCuts; // bookkeeping for dEdx cuts
393 TH2F *hTPCdEdxbefore; // TPC dEdx before cuts
394 TH2F *hTPCdEdxafter; // TPC dEdx after cuts
395 TH2F *hTPCdEdxSigbefore; // TPC Sigma dEdx before cuts
396 TH2F *hTPCdEdxSigafter; // TPC Sigm dEdx after cuts
397 TH2F *hTOFbefore; // TOF before cuts
398 TH2F *hTOFSigbefore; // TOF Sigma before cuts
399 TH2F *hTOFSigafter; // TOF Sigma after cuts
400 TH2F *hPsiPairDeltaPhiafter; // TOF Sigma after cuts
401 TH1F *hTrackCuts; // bookkeeping for track cuts
402 TH1F *hPhotonCuts; // bookkeeping for photon specific cuts
403 TH1F *hInvMassbefore; // e+e- inv mass distribution before cuts
404 TH2F *hArmenterosbefore; // armenteros podolanski plot before cuts
405 TH1F *hInvMassafter; // e+e- inv mass distribution after cuts
406 TH2F *hArmenterosafter; // armenteros podolanski plot after cuts
407 TH1F *hAcceptanceCuts; // bookkeeping for acceptance cuts
408 TH1F *hCutIndex; // bookkeeping for cuts
409 TH1F *hV0EventCuts; // bookkeeping for event selection cuts
410 TH1F *hCentrality; // centrality distribution for selected events
411 TH2F *hCentralityVsNumberOfPrimaryTracks; // centrality distribution for selected events
412 TH1F *hVertexZ; // vertex z distribution for selected events
413 TH1F *hEventPlanePhi; //EventPlaneAngle Minus Photon Angle
414 TH1F *hTriggerClass; //fired offline trigger class
415 TH1F *hTriggerClassSelected; //selected fired offline trigger class
416 TH1D *hReweightMCHistPi0; //histogram input for reweighting Pi0
417 TH1D *hReweightMCHistEta; //histogram input for reweighting Eta
418 TH1D *hReweightMCHistK0s; //histogram input for reweighting K0s
419 TF1 *fFitDataPi0; //fit to pi0 spectrum in Data
420 TF1 *fFitDataEta; //fit to eta spectrum in Data
421 TF1 *fFitDataK0s; //fit to K0s spectrum in Data
422 Int_t fAddedSignalPDGCode;
423 Bool_t fPreSelCut; // Flag for preselection cut used in V0Reader
424 Bool_t fTriggerSelectedManually; // Flag for manual trigger selection
425 TString fSpecialTriggerName; // Name of the Special Triggers
426 TString fSpecialSubTriggerName; // Name of the Special Triggers
427 Int_t fNSpecialSubTriggerOptions;
428private:
429
430 ClassDef(AliConversionCuts,9)
431};
432
433
434inline void AliConversionCuts::InitAODpidUtil(Int_t type) {
435 if (!fPIDResponse) fPIDResponse = new AliAODpidUtil();
436 Double_t alephParameters[5];
437 // simulation
438 alephParameters[0] = 2.15898e+00/50.;
439 alephParameters[1] = 1.75295e+01;
440 alephParameters[2] = 3.40030e-09;
441 alephParameters[3] = 1.96178e+00;
442 alephParameters[4] = 3.91720e+00;
443 fPIDResponse->GetTOFResponse().SetTimeResolution(80.);
444
445 // data
446 if (type==1){
447 alephParameters[0] = 0.0283086/0.97;
448 alephParameters[1] = 2.63394e+01;
449 alephParameters[2] = 5.04114e-11;
450 alephParameters[3] = 2.12543e+00;
451 alephParameters[4] = 4.88663e+00;
452 fPIDResponse->GetTOFResponse().SetTimeResolution(130.);
453 fPIDResponse->GetTPCResponse().SetMip(50.);
454 }
455
456 fPIDResponse->GetTPCResponse().SetBetheBlochParameters(
457 alephParameters[0],alephParameters[1],alephParameters[2],
458 alephParameters[3],alephParameters[4]);
459
460 fPIDResponse->GetTPCResponse().SetSigma(3.79301e-03, 2.21280e+04);
461}
462
463
464#endif