]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/multVScentPbPb/AliTrackletTaskUni.h
Fix coding violations
[u/mrichter/AliRoot.git] / PWG0 / multVScentPbPb / AliTrackletTaskUni.h
1 #ifndef ALITRACKLETTASKUNI_H
2 #define ALITRACKLETTASKUNI_H
3
4 ///////////////////////////////////////////////////////////////////////////
5 // Class AliTrackletTask                                                 //
6 // Analysis task to study performance of tracklet reconstruction         //
7 // algorithm and combinatorial background                                //
8 // Author:  M. Nicassio (INFN Bari)                                      //
9 // Contact: Maria.Nicassio@ba.infn.it, Domenico.Elia@ba.infn.it          //
10 ///////////////////////////////////////////////////////////////////////////
11
12 class TH1F; 
13 class TH2F;
14 class AliESDEvent;
15 class TList;
16 class TNtuple;
17
18 class AliMCParticle;
19 class AliITSMultRecBg;
20
21 #include "../ITS/AliITSsegmentationSPD.h"
22 #include "AliAnalysisTaskSE.h"
23 #include "AliTriggerAnalysis.h" 
24
25 class AliTrackletTaskUni : public AliAnalysisTaskSE {
26  public:
27   enum {kData,kBgInj,kBgRot,kBgMix,kMC};
28   //
29   enum {  // define here id's of the standard histos in corresponding TObjArray* fHistosTr...
30     kHEtaZvDist,      // 3 d sparse histo with dist  (uncut) vs zv vs eta
31     kHEtaZvDPhiS,     // 3 d sparse histo with dphiS (uncut) vs zv vs eta
32     kHEtaZvCut,       // zv vs eta with strict cut on tracklets applied (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     kHEtaDPhiS,       // dPhi (bending subtracted) vs eta
36     kHEtaDThetaX,     // dTheta (1/sin^2 scaled if needed) vs eta
37     kHEtaDist,        // Weighted distance vs eta
38     kHZvDPhiS,        // dPhi (bending subtracted) vs Zv
39     kHZvDThetaX,      // dTheta (1/sin^2 scaled if needed) vs Zv
40     kHZvDist          // Weighted distance vs Zv
41   };
42   enum { // define here id's of any custom histos to be added to fHistosCustom
43     kHStat,            // job info (meaning of bins defined in the enum below)
44     kHZVEtaPrimMC,     // Zv vs eta for all primary tracks (true MC multiplicity)
45     //
46     kHZVtxNoSel,       // Z vertex distribution before event selection
47     kHNTrackletsNoSel, // N tracklets before event selection
48     kHNClSPD1NoSel,    // N clusters on SPD1 before event selection
49     kHNClSPD2NoSel,    // N clusters on SPD2 before event selection
50     kHV0NoSel,         // V0 mult before selection
51     kHV0NClSPD2NoSel,  // V0 - nspd2 correlation
52     //
53     kHZVtx,            // Z vertex distribution
54     kHNTracklets,      // N tracklets
55     kHNClSPD1,         // N clusters on SPD1
56     kHNClSPD2,         // N clusters on SPD2
57     kHV0,              // V0 mult after selection
58     //
59     kHZVtxMixDiff,     // difference in Z vtx of mixed events
60     kHNTrMixDiff,      // difference in N tracklets of mixed events
61     //
62     kHPrimPDG,         // PDG code of prim tracklet
63     kHSecPDG,          // PDG code of sec tracklet
64     kHPrimParPDG,      // PDG code of prim tracklet parent
65     kHSecParPDG        // PDG code of sec tracklet parent
66   }; // custom histos
67
68   // bins for saved parameters
69   enum {kDummyBin,
70         kEvTot,       // events read
71         kEvProcData,  // events with data mult.object (ESD or reco)
72         kEvProcInj,   // events Injected
73         kEvProcRot,   // events Rotated
74         kEvProcMix,   // events Mixed
75         //
76         kDPhi,        // dphi window
77         kDTht,        // dtheta window
78         kNStd,        // N.standard deviations to keep
79         kPhiShift,    // bending shift
80         kThtS2,       // is dtheta scaled by 1/sin^2
81         kThtCW,       // on top of w.dist cut cut also on 1 sigma dThetaX
82         kPhiOvl,      // overlap params
83         kZEtaOvl,     // overlap params
84         kNoOvl,       // flag that overlap are suppressed
85         //
86         kPhiRot,      // rotation phi
87         kInjScl,      // injection scaling
88         kEtaCut,      // eta cut
89         kZVMin,       // min ZVertex to process
90         kZVMax,       // max ZVertex to process
91         kTrcMin,      // min mult to process
92         kTrcMax,      // max mult to process
93         //
94         kOneUnit=49,  // just 1 to track mergings
95         kNWorkers=50, // n workers
96         kNStatBins
97   };
98
99   //
100   AliTrackletTaskUni(const char *name = "AliTrackletTaskUni");
101   virtual ~AliTrackletTaskUni(); 
102   
103   virtual void  UserCreateOutputObjects();
104   virtual void  UserExec(Option_t *option);
105   virtual void  Terminate(Option_t *);
106
107   void       SetUseMC(Bool_t mc = kFALSE)              {fUseMC = mc;}
108   void       SetCheckReconstructables(Bool_t c=kFALSE) {fCheckReconstructables = c;}
109   TObjArray* BookHistosSet(const char* pref, UInt_t selHistos=0xffffffff);
110   TObjArray* BookCustomHistos();
111   void       AddHisto(TObjArray* histos, TObject* h, Int_t at=-1);
112   void       FillHistosSet(TObjArray* histos, double phi,double theta,double dphi,double dtheta,double dist);
113   // RS
114   void       SetNStdDev(Float_t f=1.)           {fNStdDev = f<1e-5 ? 1e-5:f;}
115   void       SetScaleDThetaBySin2T(Bool_t v=kFALSE) {fScaleDTBySin2T = v;}
116   void       SetCutOnDThetaX(Bool_t v=kFALSE)   {fCutOnDThetaX = v;}
117   void       SetPhiWindow(float w=0.08)         {fDPhiWindow   = w<1e-5 ? 1e-5:w;}
118   void       SetThetaWindow(float w=0.025)      {if (w<0) fCutOnDThetaX=kTRUE; fDThetaWindow = TMath::Abs(w)<1e-5 ? 1e-5:TMath::Abs(w);}
119   void       SetPhiShift(float w=0.0045)        {fDPhiShift = w;}
120   void       SetPhiOverlapCut(float w=0.005)    {fPhiOverlapCut = w;}
121   void       SetZetaOverlapCut(float w=0.05)    {fZetaOverlap = w;}
122   void       SetPhiRot(float w=0)               {fPhiRot = w;}
123   void       SetInjScale(Float_t s=1.)          {fInjScale = s>0? s:1.;}
124   void       SetRemoveOverlaps(Bool_t w=kFALSE) {fRemoveOverlaps = w;}
125   //
126   void       SetEtaCut(Float_t eta)             {fEtaCut = eta;}
127   void       SetZVertexMin(Float_t z)           {fZVertexMin = z;}
128   void       SetZVertexMax(Float_t z)           {fZVertexMax = z;}
129   void       SetMultCutMin(Int_t n=0)           {fMultCutMin = n;}
130   void       SetMultCutMax(Int_t n=99999)       {fMultCutMax = n;}
131   //
132   Bool_t     GetDoNormalReco()             const {return fDoNormalReco;}
133   Bool_t     GetDoInjection()              const {return fDoInjection;}
134   Bool_t     GetDoRotation()               const {return fDoRotation;}
135   Bool_t     GetDoMixing()                 const {return fDoMixing;}
136   //
137   void       SetDoNormalReco(Bool_t v=kTRUE)    {fDoNormalReco = v;}
138   void       SetDoInjection(Bool_t v=kTRUE)     {fDoInjection = v;}
139   void       SetDoRotation(Bool_t v=kTRUE)      {fDoRotation = v;}
140   void       SetDoMixing(Bool_t v=kTRUE)        {fDoMixing = v;}
141   //
142   /*
143   void       SetTrigger(AliTriggerAnalysis::Trigger trigger)  { fTrigger = trigger; }
144   void       SetMCCentralityBin(MCCentralityBin mccentrbin)   { fMCCentralityBin = mccentrbin;}
145   void       SetCentralityLowLim(Float_t centrlowlim)         { fCentrLowLim = centrlowlim;}
146   void       SetCentralityUpLim(Float_t centruplim)           { fCentrUpLim = centruplim;}
147   void       SetCentralityEst(TString centrest)               { fCentrEst = centrest;}
148   */
149   //
150  protected:
151   void       InitMultReco();
152   Bool_t     HaveCommonParent(const float* clLabs0,const float* clLabs1);
153   void       FillHistos(Int_t type, const AliMultiplicity* mlt);
154   void       FillMCPrimaries(TH2F* hetaz);
155   void       FillSpecies(Int_t primsec, Int_t id, Double_t dist);
156   Int_t      GetPdgBin(Int_t pdgCode);
157   void       CheckReconstructables();
158   //
159  protected:
160   TList*       fOutput;                   // output list send on output slot 1 
161   //
162   Bool_t       fDoNormalReco;              // do normal reco
163   Bool_t       fDoInjection;               // do injection
164   Bool_t       fDoRotation;                // do rotation
165   Bool_t       fDoMixing;                  // do mixing
166   //
167   Bool_t       fUseMC; 
168   Bool_t       fCheckReconstructables;
169   //
170   TObjArray*   fHistosTrData;              //! all tracklets in data
171   TObjArray*   fHistosTrInj;               //! injected
172   TObjArray*   fHistosTrRot;               //! rotated
173   TObjArray*   fHistosTrMix;               //! mixed
174   //
175   TObjArray*   fHistosTrPrim;              //! primary
176   TObjArray*   fHistosTrSec;               //! secondary
177   TObjArray*   fHistosTrComb;              //! combinatorials
178   TObjArray*   fHistosTrCombU;             //! combinatorials uncorrelated
179   //
180   TObjArray*   fHistosTrRcblPrim;          //! Primary Reconstructable
181   TObjArray*   fHistosTrRcblSec;           //! Secondary Reconstructable
182   TObjArray*   fHistosCustom;              //! custom histos
183   //
184   // Settings for the reconstruction
185   // tracklet reco settings
186   Float_t      fEtaCut;                    // histos filled only for this eta range
187   Float_t      fZVertexMin;                // min Z vtx to process
188   Float_t      fZVertexMax;                // max Z vtx to process
189   Int_t        fMultCutMin;                // min mult in ESD to process?
190   Int_t        fMultCutMax;                // max mult in ESD to process?
191   //
192   Bool_t       fScaleDTBySin2T;            // request dTheta scaling by 1/sin^2(theta)
193   Bool_t       fCutOnDThetaX;              // if true, apart from NStdDev cut apply also the cut on dThetaX
194   Float_t      fNStdDev;                   // cut on weighted distance
195   Float_t      fDPhiWindow;                // max dPhi
196   Float_t      fDThetaWindow;              // max dTheta
197   Float_t      fDPhiShift;                 // mean bend
198   Float_t      fPhiOverlapCut;             // overlaps cut in phi
199   Float_t      fZetaOverlap;               // overlaps cut in Z
200   Float_t      fPhiRot;                    // rotate L1 wrt L2
201   Float_t      fInjScale;                  // scaling factor for injection
202   Bool_t       fRemoveOverlaps;            // request overlaps removal
203   //
204   AliITSMultRecBg *fMultReco;              //! mult.reco object
205   TTree*       fRPTree;                    //! tree of recpoints
206   TTree*       fRPTreeMix;                 //! tree of recpoints for mixing
207   AliStack*    fStack;                     //! MC stack
208   AliMCEvent*  fMCEvent;                   //! MC Event
209   Float_t      fESDVtx[3];                 //  ESD vertex
210   //
211   /*
212   AliTriggerAnalysis::Trigger fTrigger;    // requested trigger
213   MCCentralityBin fMCCentralityBin;        // to select MC centrality bin in which corrections are calculated
214   Float_t      fCentrLowLim;               // to select centrality bin on data
215   Float_t      fCentrUpLim;                // to select centrality bin on data
216   TString      fCentrEst;                  // to select centrality estimator
217   */
218   static const char*  fgkPDGNames[];                //!pdg names
219   static const Int_t  fgkPDGCodes[];                //!pdg codes
220   //
221  private:    
222   AliTrackletTaskUni(const AliTrackletTaskUni&); // not implemented
223   AliTrackletTaskUni& operator=(const AliTrackletTaskUni&); // not implemented 
224   
225   ClassDef(AliTrackletTaskUni, 1);  
226 };
227
228
229 #endif