]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG0/multVScentPbPb/AliTrackletTaskMulti.h
Adapted the code to run with pp data w/o bg. generation
[u/mrichter/AliRoot.git] / PWG0 / multVScentPbPb / AliTrackletTaskMulti.h
CommitLineData
a9a39f46 1#ifndef ALITRACKLETTASKMULTI_H
2#define ALITRACKLETTASKMULTI_H
3
4///////////////////////////////////////////////////////////////////////////
5// Class AliTrackletTaskMulti //
6// Analysis task to produce data and MC histos needed for tracklets //
7// dNdEta extraction in multiple bins in one go //
8// Author: ruben.shahoyan@cern.ch //
9///////////////////////////////////////////////////////////////////////////
10
11class TH1F;
12class TH2F;
6c56e693 13class TH3F;
a9a39f46 14class AliESDEvent;
15class TList;
16class TNtuple;
17
18class AliMCParticle;
19class AliITSMultRecBg;
20class AliESDTrackCuts;
21
22#include "../ITS/AliITSsegmentationSPD.h"
23#include "AliAnalysisTaskSE.h"
24#include "AliTriggerAnalysis.h"
6c56e693 25#include <TMath.h>
a9a39f46 26
27class AliTrackletTaskMulti : public AliAnalysisTaskSE {
28 public:
29 enum {kData,kBgInj,kBgRot,kBgMix,kMC};
6c56e693 30 enum {kCentV0M,kCentFMD,kCentTRK,kCentTKL,kCentCL0,kCentCL1,kCentV0MvsFMD,kCentTKLvsV0,kCentZEMvsZDC,kNCentTypes}; // what is used to define centrality
a9a39f46 31 //
32 enum { // define here id's of the standard histos in corresponding TObjArray* fHistosTr...
33 kHEtaZvCut, // histo zv vs eta for tracklets passing final selection (dist<1 or |dPhi|<narrowWindow ...)
34 kHDPhiDTheta, // measured dTheta vs dPhi
35 kHDPhiSDThetaX, // dTheta (1/sin^2 scaled if needed) vs dPhi (bending subtracted)
36 kHWDist, // Weighted distance
db1bbbe1 37 kHEtaZvSPD1, // histo zv vs eta for SPD1 single clusters
a9a39f46 38 kNStandardH // number of standard histos per centrality bin
39 };
40 enum { // define here id's of any custom histos to be added to fHistosCustom
41 kHStat, // job info (meaning of bins defined in the enum below)
42 //
43 kHStatCent, // events per centrality bin with real values on the axis
44 kHStatCentBin, // events per centrality bin
45 //
46 kHNPrimMeanMC, // <n> primaries per mult bin
47 kHNPrim2PartMC, // <n> prim per part.pair per mult bin
48 kHNPrim2BCollMC, // <n> prim per bin.coll per mult bin
49 kHNPrim2PartNpMC, // <n> prim per n part vs npart
50 kHNPrim2BCollNpMC, // <n> prim per n part vs npart
51 kHNPartMC, // n.part.pairs according to MC
52 kHNPartMeanMC, // <n> part pairs per mult bin
53 kHNBCollMC, // n.bin.colls according to MC
54 kHNBCollMeanMC,
55 //
56 kHZVtxNoSel, // Z vertex distribution before event selection
6c56e693 57 kHV0NoSel, // V0 before selection
58 kHNClSPD2NoSel, // NSPD2 before selection
59 kHZDCZEMNoSel, // ZDC ZEM before selection
a9a39f46 60 //
61 kHZVtx, // Z vertex distribution
6c56e693 62 kHV0, // V0 before selection
63 kHNClSPD2, // NSPD2 before selection
64 kHZDCZEM, // ZDC ZEM before selection
65
66
a9a39f46 67 kHZVtxMixDiff, // difference in Z vtx of mixed events
68 kHNTrMixDiff, // difference in N tracklets of mixed events
69 //
70 kHPrimPDG, // PDG code of prim tracklet
71 kHSecPDG, // PDG code of sec tracklet
72 kHPrimParPDG, // PDG code of prim tracklet parent
73 kHSecParPDG, // PDG code of sec tracklet parent
74 //
75 kHClUsedInfoL0, // used clusters of lr0
76 kHClUsedInfoL1, // used clusters of lr1
77 kHClAllInfoL0, // all clusters of lr0
78 kHClAllInfoL1, // all clusters of lr1
79 //
80 // This MUST be last one: this is just beginning of many histos (one per bin)
81 kHZVEtaPrimMC // Zv vs eta for all primary tracks (true MC multiplicity)
82 }; // custom histos
83
84 // bins for saved parameters
85 enum {kDummyBin,
6c56e693 86 kEvTot0, // events read
87 kEvTot, // events read after vertex quality selection
a9a39f46 88 kOneUnit, // just 1 to track primate merges
89 kNWorkers, // n workers
90 //
6c56e693 91 kCentVar, // cetrality var. used
a9a39f46 92 kDPhi, // dphi window
93 kDTht, // dtheta window
94 kNStd, // N.standard deviations to keep
95 kPhiShift, // bending shift
96 kThtS2, // is dtheta scaled by 1/sin^2
97 kThtCW, // on top of w.dist cut cut also on 1 sigma dThetaX
98 kPhiOvl, // overlap params
99 kZEtaOvl, // overlap params
100 kNoOvl, // flag that overlap are suppressed
101 //
102 kPhiRot, // rotation phi
103 kInjScl, // injection scaling
6c56e693 104 kEtaMin, // eta cut
105 kEtaMax, // eta cut
a9a39f46 106 kZVMin, // min ZVertex to process
107 kZVMax, // max ZVertex to process
108 //
109 kDPiSCut, // cut on dphi used to extract signal (when WDist is used in analysis, put it equal to kDPhi
110 kNStdCut, // cut on weighted distance (~1) used to extract signal
111 //
112 kMCV0Scale, // scaling value for V0 in MC
113 //
114 // here we put entries for each mult.bin
115 kBinEntries = 50,
116 kEvProcData, // events with data mult.object (ESD or reco)
117 kEvProcInj, // events Injected, total
118 kEvProcRot, // events Rotated
119 kEvProcMix, // events Mixed
120 kEntriesPerBin
121 };
122
123 //
124 AliTrackletTaskMulti(const char *name = "AliTrackletTaskMulti");
125 virtual ~AliTrackletTaskMulti();
126
127 virtual void UserCreateOutputObjects();
128 virtual void UserExec(Option_t *option);
129 virtual void Terminate(Option_t *);
130
6c56e693 131 void SetUseCentralityVar(Int_t v=kCentV0M) {fUseCentralityVar = v;}
a9a39f46 132 void SetUseMC(Bool_t mc = kFALSE) {fUseMC = mc;}
133 void SetCheckReconstructables(Bool_t c=kFALSE) {fCheckReconstructables = c;}
134 TObjArray* BookHistosSet(const char* pref, UInt_t selHistos=0xffffffff);
135 TObjArray* BookCustomHistos();
136 void AddHisto(TObjArray* histos, TObject* h, Int_t at=-1);
137 void FillHistosSet(TObjArray* histos, double eta, /*double phi,double theta,*/double dphi,double dtheta,double dthetaX,double dist);
138 // RS
139 void SetNStdDev(Float_t f=1.) {fNStdDev = f<1e-5 ? 1e-5:f;}
140 void SetScaleDThetaBySin2T(Bool_t v=kFALSE) {fScaleDTBySin2T = v;}
141 void SetCutOnDThetaX(Bool_t v=kFALSE) {fCutOnDThetaX = v;}
142 void SetPhiWindow(float w=0.08) {fDPhiWindow = w<1e-5 ? 1e-5:w;}
143 void SetThetaWindow(float w=0.025) {if (w<0) fCutOnDThetaX=kTRUE; fDThetaWindow = TMath::Abs(w)<1e-5 ? 1e-5:TMath::Abs(w);}
144 void SetPhiShift(float w=0.0045) {fDPhiShift = w;}
145 void SetPhiOverlapCut(float w=0.005) {fPhiOverlapCut = w;}
146 void SetZetaOverlapCut(float w=0.05) {fZetaOverlap = w;}
147 void SetPhiRot(float w=0) {fPhiRot = w;}
148 void SetInjScale(Float_t s=1.) {fInjScale = s>0? s:1.;}
149 void SetRemoveOverlaps(Bool_t w=kFALSE) {fRemoveOverlaps = w;}
150 //
151 void SetDPhiSCut(Float_t c=0.06) {fDPhiSCut = c;}
152 void SetNStdCut(Float_t c=1.0) {fNStdCut = c;}
153 void SetScaleMCV0(Float_t s=1.0) {fMCV0Scale = s;}
154 //
6c56e693 155 void SetEtaCut(Float_t etaCut) {fEtaMax = TMath::Abs(etaCut); fEtaMin= -fEtaMax;}
156 void SetEtaMin(Float_t etaMin) {fEtaMin = etaMin;}
157 void SetEtaMax(Float_t etaMax) {fEtaMax = etaMax;}
a9a39f46 158 void SetZVertexMin(Float_t z) {fZVertexMin = z;}
159 void SetZVertexMax(Float_t z) {fZVertexMax = z;}
160 //
161 Bool_t GetDoNormalReco() const {return fDoNormalReco;}
162 Bool_t GetDoInjection() const {return fDoInjection;}
163 Bool_t GetDoRotation() const {return fDoRotation;}
164 Bool_t GetDoMixing() const {return fDoMixing;}
165 //
166 void SetDoNormalReco(Bool_t v=kTRUE) {fDoNormalReco = v;}
167 void SetDoInjection(Bool_t v=kTRUE) {fDoInjection = v;}
168 void SetDoRotation(Bool_t v=kTRUE) {fDoRotation = v;}
169 void SetDoMixing(Bool_t v=kTRUE) {fDoMixing = v;}
170 //
a9a39f46 171 //
172 protected:
173 void InitMultReco();
174 Bool_t HaveCommonParent(const float* clLabs0,const float* clLabs1);
175 void FillHistos(Int_t type, const AliMultiplicity* mlt);
176 void FillMCPrimaries();
177 void FillSpecies(Int_t primsec, Int_t id);
178 void FillClusterInfo();
db1bbbe1 179 void FillClusterInfoFromMult(const AliMultiplicity* mlt, double zVertex);
a9a39f46 180 Int_t GetPdgBin(Int_t pdgCode);
181 void CheckReconstructables();
6c56e693 182 Int_t GetCentralityBin(Float_t percentile) const;
a9a39f46 183 //
184 protected:
185 TList* fOutput; // output list send on output slot 1
186 //
187 Bool_t fDoNormalReco; // do normal reco
188 Bool_t fDoInjection; // do injection
189 Bool_t fDoRotation; // do rotation
190 Bool_t fDoMixing; // do mixing
191 //
192 Bool_t fUseMC;
193 Bool_t fCheckReconstructables;
194 //
195 TObjArray* fHistosTrData; //! all tracklets in data
196 TObjArray* fHistosTrInj; //! injected
197 TObjArray* fHistosTrRot; //! rotated
198 TObjArray* fHistosTrMix; //! mixed
199 //
200 TObjArray* fHistosTrPrim; //! primary
201 TObjArray* fHistosTrSec; //! secondary
202 TObjArray* fHistosTrComb; //! combinatorials
203 TObjArray* fHistosTrCombU; //! combinatorials uncorrelated
204 //
205 TObjArray* fHistosTrRcblPrim; //! Primary Reconstructable
206 TObjArray* fHistosTrRcblSec; //! Secondary Reconstructable
207 TObjArray* fHistosCustom; //! custom histos
208 //
209 // Settings for the reconstruction
210 // tracklet reco settings
6c56e693 211 Float_t fEtaMin; // histos filled only for this eta range
212 Float_t fEtaMax; // histos filled only for this eta range
a9a39f46 213 Float_t fZVertexMin; // min Z vtx to process
214 Float_t fZVertexMax; // max Z vtx to process
215 //
216 Bool_t fScaleDTBySin2T; // request dTheta scaling by 1/sin^2(theta)
217 Bool_t fCutOnDThetaX; // if true, apart from NStdDev cut apply also the cut on dThetaX
218 Float_t fNStdDev; // cut on weighted distance
219 Float_t fDPhiWindow; // max dPhi
220 Float_t fDThetaWindow; // max dTheta
221 Float_t fDPhiShift; // mean bend
222 Float_t fPhiOverlapCut; // overlaps cut in phi
223 Float_t fZetaOverlap; // overlaps cut in Z
224 Float_t fPhiRot; // rotate L1 wrt L2
225 Float_t fInjScale; // scaling factor for injection
226 Bool_t fRemoveOverlaps; // request overlaps removal
227 //
228 Float_t fDPhiSCut; // cut on signal dphiS
229 Float_t fNStdCut; // cut on signal weighted distance
230 Float_t fMCV0Scale; // scaling factor for V0 in MC
231 //
232 AliITSMultRecBg *fMultReco; //! mult.reco object
233 TTree* fRPTree; //! tree of recpoints
234 TTree* fRPTreeMix; //! tree of recpoints for mixing
235 AliStack* fStack; //! MC stack
236 AliMCEvent* fMCEvent; //! MC Event
237 Float_t fESDVtx[3]; // ESD vertex
a9a39f46 238 //
a9a39f46 239 Float_t fNPart; // number of participant pairs from MC
240 Float_t fNBColl; // number of bin. collision from MC
241 Int_t fCurrCentBin; // current centrality bin
242 Int_t fNCentBins; // N of mult bins
243 Int_t fUseCentralityVar; // what is used to determine the centrality
a9a39f46 244 //
6c56e693 245 static const Float_t fgkCentPerc[]; //! centrality in percentiles
a9a39f46 246 //
6c56e693 247 static const char* fgCentSelName[]; //!centrality types
a9a39f46 248 static const char* fgkPDGNames[]; //!pdg names
249 static const Int_t fgkPDGCodes[]; //!pdg codes
250 //
251 private:
252 AliTrackletTaskMulti(const AliTrackletTaskMulti&); // not implemented
253 AliTrackletTaskMulti& operator=(const AliTrackletTaskMulti&); // not implemented
254
255 ClassDef(AliTrackletTaskMulti, 1);
256};
257
258
259#endif