]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGUD/multVScentPbPb/AliTrackletTaskMulti.h
Update PR task: drathee
[u/mrichter/AliRoot.git] / PWGUD / multVScentPbPb / AliTrackletTaskMulti.h
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
11 class TH1F; 
12 class TH2F;
13 class TH3F;
14 class AliESDEvent;
15 class TList;
16 class TNtuple;
17
18 class AliMCParticle;
19 class AliITSMultRecBg;
20 class AliESDTrackCuts;
21
22 #include "../ITS/AliITSsegmentationSPD.h"
23 #include "AliAnalysisTaskSE.h"
24 #include "AliTriggerAnalysis.h" 
25 #include <TMath.h>
26
27 class AliTrackletTaskMulti : public AliAnalysisTaskSE {
28  public:
29   enum {kData,kBgInj,kBgRot,kBgMix,kMC};
30   enum {kCentV0M,kCentFMD,kCentTRK,kCentTKL,kCentCL0,kCentCL1,kCentV0MvsFMD,kCentTKLvsV0,kCentZEMvsZDC,kNCentTypes}; // what is used to define centrality
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 
37     kHEtaZvSPD1,      // histo zv vs eta for SPD1 single clusters
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
57     kHV0NoSel,         // V0 before selection
58     kHNClSPD2NoSel,    // NSPD2 before selection
59     kHZDCZEMNoSel,     // ZDC ZEM before selection
60     //
61     kHZVtx,            // Z vertex distribution
62     kHV0,              // V0 before selection
63     kHNClSPD2,         // NSPD2 before selection
64     kHZDCZEM,          // ZDC ZEM before selection
65     
66
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,
86         kEvTot0,      // events read
87         kEvTot,       // events read after vertex quality selection
88         kOneUnit,     // just 1 to track primate merges
89         kNWorkers,    // n workers
90         //
91         kCentVar,     // cetrality var. used
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
104         kEtaMin,      // eta cut
105         kEtaMax,      // eta cut
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   void       RegisterStat();
131
132   void       SetUseCentralityVar(Int_t v=kCentV0M)     {fUseCentralityVar = v;}
133   void       SetUseMC(Bool_t mc = kFALSE)              {fUseMC = mc;}
134   void       SetCheckReconstructables(Bool_t c=kFALSE) {fCheckReconstructables = c;}
135   TObjArray* BookHistosSet(const char* pref, UInt_t selHistos=0xffffffff);
136   TObjArray* BookCustomHistos();
137   void       AddHisto(TObjArray* histos, TObject* h, Int_t at=-1);
138   void       FillHistosSet(TObjArray* histos, double eta, /*double phi,double theta,*/double dphi,double dtheta,double dthetaX,double dist);
139   // RS
140   void       SetNStdDev(Float_t f=1.)           {fNStdDev = f<1e-5 ? 1e-5:f;}
141   void       SetScaleDThetaBySin2T(Bool_t v=kFALSE) {fScaleDTBySin2T = v;}
142   void       SetCutOnDThetaX(Bool_t v=kFALSE)   {fCutOnDThetaX = v;}
143   void       SetPhiWindow(float w=0.08)         {fDPhiWindow   = w<1e-5 ? 1e-5:w;}
144   void       SetThetaWindow(float w=0.025)      {if (w<0) fCutOnDThetaX=kTRUE; fDThetaWindow = TMath::Abs(w)<1e-5 ? 1e-5:TMath::Abs(w);}
145   void       SetPhiShift(float w=0.0045)        {fDPhiShift = w;}
146   void       SetPhiOverlapCut(float w=0.005)    {fPhiOverlapCut = w;}
147   void       SetZetaOverlapCut(float w=0.05)    {fZetaOverlap = w;}
148   void       SetPhiRot(float w=0)               {fPhiRot = w;}
149   void       SetInjScale(Float_t s=1.)          {fInjScale = s>0? s:1.;}
150   void       SetRemoveOverlaps(Bool_t w=kFALSE) {fRemoveOverlaps = w;}
151   //
152   void       SetDPhiSCut(Float_t c=0.06)        {fDPhiSCut = c;}
153   void       SetNStdCut(Float_t c=1.0)          {fNStdCut = c;}
154   void       SetScaleMCV0(Float_t s=1.0)        {fMCV0Scale = s;}  
155   //
156   void       SetEtaCut(Float_t etaCut)          {fEtaMax = TMath::Abs(etaCut); fEtaMin= -fEtaMax;}
157   void       SetEtaMin(Float_t etaMin)          {fEtaMin = etaMin;}
158   void       SetEtaMax(Float_t etaMax)          {fEtaMax = etaMax;}
159   void       SetZVertexMin(Float_t z)           {fZVertexMin = z;}
160   void       SetZVertexMax(Float_t z)           {fZVertexMax = z;}
161   //
162   Bool_t     GetDoNormalReco()             const {return fDoNormalReco;}
163   Bool_t     GetDoInjection()              const {return fDoInjection;}
164   Bool_t     GetDoRotation()               const {return fDoRotation;}
165   Bool_t     GetDoMixing()                 const {return fDoMixing;}
166   //
167   void       SetDoNormalReco(Bool_t v=kTRUE)    {fDoNormalReco = v;}
168   void       SetDoInjection(Bool_t v=kTRUE)     {fDoInjection = v;}
169   void       SetDoRotation(Bool_t v=kTRUE)      {fDoRotation = v;}
170   void       SetDoMixing(Bool_t v=kTRUE)        {fDoMixing = v;}
171   //
172   //
173  protected:
174   void       InitMultReco();
175   Bool_t     HaveCommonParent(const float* clLabs0,const float* clLabs1);
176   void       FillHistos(Int_t type, const AliMultiplicity* mlt);
177   void       FillMCPrimaries();
178   void       FillSpecies(Int_t primsec, Int_t id);
179   void       FillClusterInfo();
180   void       FillClusterInfoFromMult(const AliMultiplicity* mlt, double zVertex);
181   Int_t      GetPdgBin(Int_t pdgCode);
182   void       CheckReconstructables();
183   Int_t      GetCentralityBin(Float_t percentile) const;
184   //
185  protected:
186   TList*       fOutput;                   // output list send on output slot 1 
187   //
188   Bool_t       fDoNormalReco;              // do normal reco
189   Bool_t       fDoInjection;               // do injection
190   Bool_t       fDoRotation;                // do rotation
191   Bool_t       fDoMixing;                  // do mixing
192   //
193   Bool_t       fUseMC; 
194   Bool_t       fCheckReconstructables;
195   //
196   TObjArray*   fHistosTrData;              //! all tracklets in data
197   TObjArray*   fHistosTrInj;               //! injected
198   TObjArray*   fHistosTrRot;               //! rotated
199   TObjArray*   fHistosTrMix;               //! mixed
200   //
201   TObjArray*   fHistosTrPrim;              //! primary
202   TObjArray*   fHistosTrSec;               //! secondary
203   TObjArray*   fHistosTrComb;              //! combinatorials
204   TObjArray*   fHistosTrCombU;             //! combinatorials uncorrelated
205   //
206   TObjArray*   fHistosTrRcblPrim;          //! Primary Reconstructable
207   TObjArray*   fHistosTrRcblSec;           //! Secondary Reconstructable
208   TObjArray*   fHistosCustom;              //! custom histos
209   //
210   // Settings for the reconstruction
211   // tracklet reco settings
212   Float_t      fEtaMin;                    // histos filled only for this eta range
213   Float_t      fEtaMax;                    // histos filled only for this eta range
214   Float_t      fZVertexMin;                // min Z vtx to process
215   Float_t      fZVertexMax;                // max Z vtx to process
216   //
217   Bool_t       fScaleDTBySin2T;            // request dTheta scaling by 1/sin^2(theta)
218   Bool_t       fCutOnDThetaX;              // if true, apart from NStdDev cut apply also the cut on dThetaX
219   Float_t      fNStdDev;                   // cut on weighted distance
220   Float_t      fDPhiWindow;                // max dPhi
221   Float_t      fDThetaWindow;              // max dTheta
222   Float_t      fDPhiShift;                 // mean bend
223   Float_t      fPhiOverlapCut;             // overlaps cut in phi
224   Float_t      fZetaOverlap;               // overlaps cut in Z
225   Float_t      fPhiRot;                    // rotate L1 wrt L2
226   Float_t      fInjScale;                  // scaling factor for injection
227   Bool_t       fRemoveOverlaps;            // request overlaps removal
228   //
229   Float_t      fDPhiSCut;                  // cut on signal dphiS
230   Float_t      fNStdCut;                   // cut on signal weighted distance
231   Float_t      fMCV0Scale;                 // scaling factor for V0 in MC
232   //
233   AliITSMultRecBg *fMultReco;              //! mult.reco object
234   TTree*       fRPTree;                    //! tree of recpoints
235   TTree*       fRPTreeMix;                 //! tree of recpoints for mixing
236   AliStack*    fStack;                     //! MC stack
237   AliMCEvent*  fMCEvent;                   //! MC Event
238   Float_t      fESDVtx[3];                 //  ESD vertex
239   //
240   Float_t fNPart;                          // number of participant pairs from MC
241   Float_t fNBColl;                         // number of bin. collision from MC
242   Int_t  fCurrCentBin;                     // current centrality bin
243   Int_t  fNCentBins;                       // N of mult bins
244   Int_t  fUseCentralityVar;                // what is used to determine the centrality
245   //
246   static const Float_t fgkCentPerc[];               //! centrality in percentiles
247   //
248   static const char*  fgCentSelName[];              //!centrality types
249   static const char*  fgkPDGNames[];                //!pdg names
250   static const Int_t  fgkPDGCodes[];                //!pdg codes
251   //
252  private:    
253   AliTrackletTaskMulti(const AliTrackletTaskMulti&); // not implemented
254   AliTrackletTaskMulti& operator=(const AliTrackletTaskMulti&); // not implemented 
255   
256   ClassDef(AliTrackletTaskMulti, 1);  
257 };
258
259
260 #endif