verbose output of Kink tasks only when fDebug is >1
[u/mrichter/AliRoot.git] / PWG2 / KINK / AliResonanceKink.h
CommitLineData
f9afc48d 1#ifndef ALIRESONANCEKINK_H
2#define ALIRESONANCEKINK_H
10eaad41 3
4/* See cxx source for full Copyright notice */
5
6//------------------------------------------------------------------------------
f9afc48d 7// class AliResonanceKink
10eaad41 8// This task is an example of an analysis task
9// for analysing resonances having one kaon kink
10//Author: Paraskevi Ganoti, University of Athens (pganoti@phys.uoa.gr)
11//------------------------------------------------------------------------------
92adf4f6 12#include "TVector3.h"
10eaad41 13class TF1;
92adf4f6 14class TH1D;
15class TH2D;
10eaad41 16class AliESDEvent;
f27d6e67 17class AliESDtrack;
f9afc48d 18class AliESDVertex;
19class AliMCEvent;
20class TList;
21class TString;
10eaad41 22
f9afc48d 23class AliResonanceKink : public TObject {
10eaad41 24 public:
f9afc48d 25
26 enum ResonanceType {kPhi=333, kKstar0=313, kLambda1520=3124};
27 enum DaughterType {kdaughterPion=211, kdaughterKaon=321, kdaughterProton=2212};
10eaad41 28
f9afc48d 29 AliResonanceKink();
30 AliResonanceKink(Int_t nbins, Float_t nlowx, Float_t nhighx, Int_t daughter1, Int_t daughter2, Int_t resonancePDG);
31 virtual ~AliResonanceKink();
10eaad41 32
f9afc48d 33 TList* GetHistogramList();
34 void Analyse(AliESDEvent* esd, AliMCEvent* mcEvent);
10eaad41 35 Float_t GetSigmaToVertex(AliESDtrack* esdTrack) const ;
36 const AliESDVertex *GetEventVertex(const AliESDEvent* esd) const;
f9afc48d 37 void SetAnalysisType(TString type) {fAnalysisType=type;}
38 void SetPDGCodes(Int_t d1, Int_t d2, Int_t res) {fdaughter1pdg=d1; fdaughter2pdg=d2; fresonancePDGcode=res;}
39 void InitOutputHistograms(Int_t nbins, Float_t nlowx, Float_t nhighx);
92adf4f6 40 Bool_t IsAcceptedForKink(AliESDEvent *localesd, const AliESDVertex *localvertex, AliESDtrack *localtrack);
41 Bool_t IsAcceptedForTrack(AliESDEvent *localesd, const AliESDVertex *localvertex, AliESDtrack *localtrack);
42 Bool_t IsKink(AliESDEvent *localesd, Int_t kinkIndex, TVector3 trackMom);
43
44 void SetMaxNsigmaToVertex(Float_t maxNSigmaToVertex) {
45 fMaxNSigmaToVertex=maxNSigmaToVertex;
46 }
47 Float_t GetMaxNsigmaToVertex() const {return fMaxNSigmaToVertex;}
48
49 void SetPtTrackCut(Double_t minPtTrackCut) {
50 fMinPtTrackCut=minPtTrackCut;
51 }
52 Double_t GetPtTrackCut() const {return fMinPtTrackCut;}
53
54 void SetMaxDCAxy(Double_t maxDCAxy) {
55 fMaxDCAxy=maxDCAxy;
56 }
57 Double_t GetMaxDCAxy() const {return fMaxDCAxy;}
58
59 void SetMaxDCAzaxis(Double_t maxDCAzaxis) {
60 fMaxDCAzaxis=maxDCAzaxis;
61 }
62 Double_t GetMaxDCAzaxis() const {return fMaxDCAzaxis;}
63
64 void SetMinTPCclusters(Int_t minTPCclusters) {
65 fMinTPCclusters=minTPCclusters;
66 }
67 Int_t GetMinTPCclusters() const {return fMinTPCclusters;}
68
69 void SetMaxChi2PerTPCcluster(Double_t maxChi2PerTPCcluster) {
70 fMaxChi2PerTPCcluster=maxChi2PerTPCcluster;
71 }
72 Double_t GetMaxChi2PerTPCcluster() const {return fMaxChi2PerTPCcluster;}
73
74 void SetMaxCov0(Double_t maxCov0) {
75 fMaxCov0=maxCov0;
76 }
77 Double_t GetMaxCov0() const {return fMaxCov0;}
78
79 void SetMaxCov2(Double_t maxCov2) {
80 fMaxCov2=maxCov2;
81 }
82 Double_t GetMaxCov2() const {return fMaxCov2;}
83
84 void SetMaxCov5(Double_t maxCov5) {
85 fMaxCov5=maxCov5;
86 }
87 Double_t GetMaxCov5() const {return fMaxCov5;}
88
89 void SetMaxCov9(Double_t maxCov9) {
90 fMaxCov9=maxCov9;
91 }
92 Double_t GetMaxCov9() const {return fMaxCov9;}
93
94 void SetMaxCov14(Double_t maxCov14) {
95 fMaxCov14=maxCov14;
96 }
97 Double_t GetMaxCov14() const {return fMaxCov14;}
98
99 //void SetTPCrefit() {Int_t fTPCrefitFlag=kTRUE;}
100
10eaad41 101 private:
f9afc48d 102
92adf4f6 103 TList *fListOfHistos; // List
104 TH1D *fOpeningAngle; // Opening
105 TH1D *fInvariantMass; // invMass spectrum
106 TH1D *fInvMassTrue; // invMassTrue spectrum
107 TH1D *fPhiBothKinks; // bothKaonsKinks
108 TH1D *fRecPt; // pT spectrum
109 TH1D *fRecEta; // Eta spectrum
110 TH2D *fRecEtaPt; // Eta pT spectrum
111 TH1D *fSimPt; // pT Sim spectrum
112 TH1D *fSimEta; // Eta Sim spectrum
113 TH2D *fSimEtaPt; // Eta pT Sim spectrum
114 TH1D *fSimPtKink; // pT Sim one kaon kink spectrum
115 TH1D *fSimEtaKink; // Eta Sim one kaon kink spectrum spectrum
116 TH2D *fSimEtaPtKink; // Eta pT Sim one kaon kink spectrum
117 TH1D *fhdr ; // radial impact
118 TH1D *fhdz ; // z impact
10eaad41 119 TF1 *f1;
120 TF1 *f2;
121 TString fAnalysisType;//"ESD" or "MC"
92adf4f6 122 TH1D *fvtxz ; // vtx z component
123 Int_t fNbins; // bins
124 Float_t fLowX;// lowx
125 Float_t fHighX; // high x
f9afc48d 126 Int_t fdaughter1pdg;
127 Int_t fdaughter2pdg;
128 Int_t fresonancePDGcode;
92adf4f6 129 Float_t fMaxNSigmaToVertex;
130 Double_t fMinPtTrackCut;
131 Double_t fMaxDCAxy;
132 Double_t fMaxDCAzaxis;
133 Int_t fMinTPCclusters;
134 Double_t fMaxChi2PerTPCcluster;
135 Double_t fMaxCov0;
136 Double_t fMaxCov2;
137 Double_t fMaxCov5;
138 Double_t fMaxCov9;
139 Double_t fMaxCov14;
140// Bool_t fTPCrefitFlag;
10eaad41 141
f9afc48d 142 AliResonanceKink(const AliResonanceKink&); // not implemented
143 AliResonanceKink& operator=(const AliResonanceKink&); // not implemented
10eaad41 144
f9afc48d 145 ClassDef(AliResonanceKink, 1); // example of analysis
10eaad41 146};
147
148#endif