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