]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/STRANGENESS/Cascades/AliAnalysisTaskCheckCascadePbPb.h
Min number of TPC clusters as parameter
[u/mrichter/AliRoot.git] / PWGLF / STRANGENESS / Cascades / AliAnalysisTaskCheckCascadePbPb.h
1 #ifndef ALIANALYSISTASKCHECKCASCADEPBPB_H
2 #define ALIANALYSISTASKCHECKCASCADEPBPB_H
3
4 /*  See cxx source for full Copyright notice */
5
6 //-----------------------------------------------------------------
7 //            AliAnalysisTaskCheckCascadePbPb class
8 //              Origin AliAnalysisTaskCheckCascade
9 //              This task has four roles :
10 //                1. QAing the Cascades from ESD and AOD
11 //                   Origin:  AliAnalysisTaskESDCheckV0 by Boris Hippolyte Nov2007, hippolyt@in2p3.fr
12 //                2. Prepare the plots which stand as raw material for yield extraction (wi/wo PID)
13 //                3. Supply an AliCFContainer meant to define the optimised topological selections
14 //                4. Rough azimuthal correlation study (Eta, Phi)
15 //                Adapted to Cascade : A.Maire Mar2008, antonin.maire@ires.in2p3.fr
16 //                Modified :           A.Maire Mar2010, antonin.maire@ires.in2p3.fr
17 //                Modified for PbPb analysis: M. Nicassio Feb 2011, maria.nicassio@ba.infn.it
18 //-----------------------------------------------------------------
19
20 class TList;
21 class TH1F;
22 class TH2F;
23 class TH3F;
24 class TVector3;
25 class THnSparse;
26  
27 class AliESDEvent;
28 class AliPhysicsSelection;
29 class AliCFContainer;
30 class AliPIDResponse;
31
32 #include "TString.h"
33
34 #include "AliAnalysisTaskSE.h"
35
36 class AliAnalysisTaskCheckCascadePbPb : public AliAnalysisTaskSE {
37  public:
38   AliAnalysisTaskCheckCascadePbPb();
39   AliAnalysisTaskCheckCascadePbPb(const char *name);
40   virtual ~AliAnalysisTaskCheckCascadePbPb();
41   
42   virtual void   UserCreateOutputObjects();
43   virtual void   UserExec(Option_t *option);
44 /*        void   DoAngularCorrelation(const Char_t   *lCascType, 
45                                             Double_t  lInvMassCascade, 
46                                       const Int_t    *lArrTrackID,
47                                             TVector3 &lTVect3MomXi, 
48                                             Double_t  lEtaXi);*/
49   virtual Int_t  DoESDTrackWithTPCrefitMultiplicity(const AliESDEvent *lESDevent);
50   virtual void   Terminate(Option_t *);
51   
52   void SetAnalysisType               (const char* analysisType          = "ESD") { fAnalysisType                = analysisType;               }
53   void SetRelaunchV0CascVertexers    (Bool_t rerunV0CascVertexers       = 0    ) { fkRerunV0CascVertexers       = rerunV0CascVertexers;       }
54   void SetQualityCutZprimVtxPos      (Bool_t qualityCutZprimVtxPos      = kTRUE) { fkQualityCutZprimVtxPos      = qualityCutZprimVtxPos;      }
55   void SetQualityCutNoTPConlyPrimVtx (Bool_t qualityCutNoTPConlyPrimVtx = kTRUE) { fkQualityCutNoTPConlyPrimVtx = qualityCutNoTPConlyPrimVtx; }
56   void SetQualityCutTPCrefit         (Bool_t qualityCutTPCrefit         = kTRUE) { fkQualityCutTPCrefit         = qualityCutTPCrefit;         }
57   void SetQualityCutnTPCcls          (Bool_t qualityCutnTPCcls          = kTRUE) { fkQualityCutnTPCcls          = qualityCutnTPCcls;          }
58   void SetQualityCutMinnTPCcls       (Int_t minnTPCcls                  = 70   ) { fMinnTPCcls                  = minnTPCcls;                 }
59   void SetExtraSelections            (Bool_t extraSelections            = 0    ) { fkExtraSelections            = extraSelections;            }
60   void SetCentralityLowLim           (Float_t centrlowlim               = 0.   ) { fCentrLowLim                 = centrlowlim;                }  
61   void SetCentralityUpLim            (Float_t centruplim                = 100. ) { fCentrUpLim                  = centruplim;                 }
62   void SetCentralityEst              (TString   centrest                = "V0M") { fCentrEstimator              = centrest;                   }
63   void SetVertexRange                (Float_t vtxrange                  = 0.   ) { fVtxRange                    = vtxrange;                   }
64
65  private:
66         // Note : In ROOT, "//!" means "do not stream the data from Master node to Worker node" ...
67         // your data member object is created on the worker nodes and streaming is not needed.
68         // http://root.cern.ch/download/doc/11InputOutput.pdf, page 14
69
70
71         TString         fAnalysisType;                  // "ESD" or "AOD" analysis type 
72         AliESDtrackCuts *fESDtrackCuts;                 // ESD track cuts used for primary track definition
73         //TPaveText       *fPaveTextBookKeeping;          // TString to store all the relevant info necessary for book keeping (v0 cuts, cascade cuts, quality cuts, ...)
74         AliPIDResponse *fPIDResponse;                   //! PID response object
75
76         Bool_t          fkRerunV0CascVertexers;         // Boolean : kTRUE = relaunch both V0 + Cascade vertexers
77         Bool_t          fkQualityCutZprimVtxPos;        // Boolean : kTRUE = cut on the prim.vtx  z-position
78         Bool_t          fkQualityCutNoTPConlyPrimVtx;   // Boolean : kTRUE = prim vtx should be SPD or Tracking vertex
79         Bool_t          fkQualityCutTPCrefit;           // Boolean : kTRUE = ask for TPCrefit for the 3 daughter tracks
80         Bool_t          fkQualityCutnTPCcls;            // Boolean : kTRUE = ask for at least n TPC clusters for each daughter track
81         Int_t           fMinnTPCcls;                    // minimum number of TPC cluster for daughter tracks
82         Bool_t          fkExtraSelections;              // Boolean : kTRUE = apply tighter selections, before starting the analysis
83         Float_t         fCentrLowLim;                   // Lower limit for centrality percentile selection
84         Float_t         fCentrUpLim;                    // Upper limit for centrality percentile selection
85         TString         fCentrEstimator;                // string for the centrality estimator
86         Float_t         fVtxRange;                      // to select events with |zvtx|<fVtxRange cm
87
88        
89         Double_t        fV0Sels[7];                     // Array to store the 7 values for the different selections V0 related (if fkRerunV0CascVertexers)
90         Double_t        fCascSels[8];                   // Array to store the 8 values for the different selections Casc. related (if fkRerunV0CascVertexers)
91
92         TList      *fListHistCascade;                   //! List of Cascade histograms
93         
94         // - General histos (filled before the trigger selection)
95         TH2F    *fHistEvtsInCentralityBinsvsNtracks;  //! Events in centrality bins vs N ESDtracks
96         TH1F    *fHistCascadeMultiplicityBeforeEvSel; //! Cascade multiplicity distribution
97          
98         // - General histos (filled for any triggered event)
99         TH1F    *fHistCascadeMultiplicityForCentrEvt;              //! Cascade multiplicity distribution
100         TH1F    *fHistTrackMultiplicityForCentrEvt;                //! Track multiplicity distribution (without any cut = include ITS stand-alone + global tracks)
101         TH1F    *fHistTPCrefitTrackMultiplicityForCentrEvt;        //! Track multiplicity distribution for tracks with TPCrefit
102         
103         // - General histos (filled for events selected in this analysis (|z(prim. vtx)| < 10 cm + prim vtx = SPD or Tracking) )
104         TH1F    *fHistCascadeMultiplicityForSelEvt;     //! Cascade multiplicity distribution
105         TH1F    *fHistPosBestPrimaryVtxXForSelEvt;      //! (best) primary vertex position distribution in x 
106         TH1F    *fHistPosBestPrimaryVtxYForSelEvt;      //! (best) primary vertex position distribution in y
107         TH1F    *fHistPosBestPrimaryVtxZForSelEvt;      //! (best) primary vertex position distribution in z
108         
109         
110         
111
112         // - Characteristics for event with >1 cascade : Track Multiplicity, TPC clusters + Prim. vertex positions
113         TH1F    *fHistTPCrefitTrackMultiplicityForCascadeEvt;   //! TPCrefit Track multiplicity distribution for event with >1 cascade candidate (NB: after quality sel. within the task)
114
115         TH1F    *fHistPosV0TPCClusters;                 //! TPC clusters distribution for Positive V0 daughter track
116         TH1F    *fHistNegV0TPCClusters;                 //! TPC clusters distribution for Negative V0 daughter track
117         TH1F    *fHistBachTPCClusters;                  //! TPC clusters distribution for Bachelor             track
118
119         TH1F    *fHistVtxStatus;                        //! Is there a tracking vertex in the cascade event ?
120
121                 // Vtx coming from the full tracking, for events containing at least a cascade
122         TH1F    *fHistPosTrkgPrimaryVtxXForCascadeEvt;  //! primary vertex position distribution in x 
123         TH1F    *fHistPosTrkgPrimaryVtxYForCascadeEvt;  //! primary vertex position distribution in y
124         TH1F    *fHistPosTrkgPrimaryVtxZForCascadeEvt;  //! primary vertex position distribution in z
125         TH1F    *fHistTrkgPrimaryVtxRadius;             //! primary vertex (3D) radius distribution 
126
127                 // Best primary Vtx available, for events containing at least a cascade
128         TH1F    *fHistPosBestPrimaryVtxXForCascadeEvt;  //! (best) primary vertex position distribution in x 
129         TH1F    *fHistPosBestPrimaryVtxYForCascadeEvt;  //! (best) primary vertex position distribution in y
130         TH1F    *fHistPosBestPrimaryVtxZForCascadeEvt;  //! (best) primary vertex position distribution in z
131         TH1F    *fHistBestPrimaryVtxRadius;             //! (best) primary vertex radius distribution 
132
133                 // Correlation Best Vtx / Full Tracking Vtx
134         TH2F    *f2dHistTrkgPrimVtxVsBestPrimVtx;       //!  Radius of prim. Vtx from tracks Vs Radius of best Prim. Vtx
135
136
137 // PART 1 : Adavanced QA
138 // - Typical histos on the variables used for the selection of cascades
139         TH1F    *fHistEffMassXi;                        //! reconstructed cascade effective mass
140         TH1F    *fHistChi2Xi;                           //! chi2 value
141         TH1F    *fHistDcaXiDaughters;                   //! dca between Xi's daughters
142         TH1F    *fHistDcaBachToPrimVertex;              //! dca of the bachelor track to primary vertex
143         TH1F    *fHistXiCosineOfPointingAngle;          //! cosine of Xi pointing angle in a cascade
144         TH1F    *fHistXiRadius;                         //! (transverse) radius of the cascade vertex 
145                 
146         // - Histos about ~ the "V0 selection part" of the cascade,  coming by inheritance from AliESDv0
147         TH1F    *fHistMassLambdaAsCascDghter;           //! Test Invariant Mass of Lambda coming from Cascade
148         TH1F    *fHistV0Chi2Xi;                         //! V0 chi2 distribution, for the V0 associated to a cascade
149         TH1F    *fHistDcaV0DaughtersXi;                 //! Dca between V0 daughters, for the V0 associated to a cascade
150         TH1F    *fHistDcaV0ToPrimVertexXi;              //! Dca of V0 to primary vertex, for the V0 associated to a cascade     
151         TH1F    *fHistV0CosineOfPointingAngle;          //! Cosine of V0 pointing angle, for the V0 associated to a cascade
152         TH1F    *fHistV0RadiusXi;                       //! V0 (transverse) distance distribution, for the V0 associated to a cascade
153
154         TH1F    *fHistDcaPosToPrimVertexXi;             //! Dca of V0 positive daughter to primary vertex, for the V0 associated to a cascade
155         TH1F    *fHistDcaNegToPrimVertexXi;             //! Dca of V0 negative daughter to primary vertex, for the V0 associated to a cascade
156         
157
158         // - Effective mass histos for cascades.
159         TH1F    *fHistMassXiMinus;                      //! reconstructed cascade effective mass, under Xi- hyp.
160         TH1F    *fHistMassXiPlus;                       //! reconstructed cascade effective mass, under Xi+ hyp.
161         TH1F    *fHistMassOmegaMinus;                   //! reconstructed cascade effective mass, under Omega- hyp.
162         TH1F    *fHistMassOmegaPlus;                    //! reconstructed cascade effective mass, under Omega+ hyp.
163         
164         TH1F    *fHistMassWithCombPIDXiMinus;           //! reconstructed Xi- effective mass, with bach. comb PID
165         TH1F    *fHistMassWithCombPIDXiPlus;            //! reconstructed Xi+ effective mass, with bach. comb PID
166         TH1F    *fHistMassWithCombPIDOmegaMinus;        //! reconstructed Omega- effective mass, with bach. comb PID
167         TH1F    *fHistMassWithCombPIDOmegaPlus;         //! reconstructed Omega+ effective mass, with bach. comb PID
168
169         // - Complements for QA
170         TH1F    *fHistXiTransvMom;                      //! Xi transverse momentum, around the mass peak of Xi-/+ 
171         TH1F    *fHistXiTotMom;                         //! Xi momentum norm, around the mass peak of Xi-/+
172         
173         TH1F    *fHistBachTransvMomXi;                  //! bachelor transverse momentum, for cand. around the mass peak of Xi-/+
174         TH1F    *fHistBachTotMomXi;                     //! bachelor momentum norm, for cand. around the mass peak of Xi-/+
175
176         TH1F    *fHistChargeXi;                         //! Charge sign of the cascade candidate
177         TH1F    *fHistV0toXiCosineOfPointingAngle;      //! Cos. of Pointing angle between the V0 mom and the Xi-V0 vtx line
178   
179         TH1F    *fHistRapXi;                            //! rapidity of Xi candidates, around the mass peak of Xi-/+
180         TH1F    *fHistRapOmega;                         //! rapidity of Omega candidates, around the mass peak of Omega-/+
181         TH1F    *fHistEtaXi;                            //! eta distrib. of all the cascade candidates, around the mass peak of Xi-/+
182         TH1F    *fHistThetaXi;                          //! theta distrib. of all the cascade candidates, around the mass peak of Xi-/+
183         TH1F    *fHistPhiXi;                            //! phi distrib. of all the cascade candidates, around the mass peak of Xi-/+
184         
185         TH2F    *f2dHistArmenteros;                     //! alpha(casc. cand.) Vs PtArm(casc. cand.)
186         
187         TH2F    *f2dHistEffMassLambdaVsEffMassXiMinus;  //! Xi- Eff mass Vs V0 Eff mass, under Xi- hyp.
188         TH2F    *f2dHistEffMassXiVsEffMassOmegaMinus;   //! Xi- Eff mass Vs Omega- Eff mass, for negative cascades
189         TH2F    *f2dHistEffMassLambdaVsEffMassXiPlus;   //! Xi+ Eff mass Vs V0 Eff mass, under Xi+ hyp.
190         TH2F    *f2dHistEffMassXiVsEffMassOmegaPlus;    //! Xi+ Eff mass Vs Omega+ Eff mass, for positive cascades
191         
192         TH2F    *f2dHistXiRadiusVsEffMassXiMinus;       //! transv. casc. decay radius Vs Xi- Eff mass, under Xi- hyp.
193         TH2F    *f2dHistXiRadiusVsEffMassXiPlus;        //! transv. casc. decay radius Vs Xi+ Eff mass, under Xi+ hyp.
194         TH2F    *f2dHistXiRadiusVsEffMassOmegaMinus;    //! transv. casc. decay radius Vs Omega- Eff mass, under Omega- hyp.
195         TH2F    *f2dHistXiRadiusVsEffMassOmegaPlus;     //! transv. casc. decay radius Vs Omega+ Eff mass, under Omega+ hyp.
196         
197         TH2F    *f2dHistTPCdEdxOfCascDghters;           //! TPC Bethe-Bloch curve, populated with the cascade daughters
198         
199         
200         // PART 2 : TH3F needed for pt spectrum and yield extraction
201         // Without any PID
202 /*      TH3F    *f3dHistXiPtVsEffMassVsYXiMinus;        //! casc. transv. momemtum Vs Xi- Eff mass Vs Y
203         TH3F    *f3dHistXiPtVsEffMassVsYXiPlus;         //! casc. transv. momemtum Vs Xi+ Eff mass Vs Y
204         TH3F    *f3dHistXiPtVsEffMassVsYOmegaMinus;     //! casc. transv. momemtum Vs Omega- Eff mass Vs Y
205         TH3F    *f3dHistXiPtVsEffMassVsYOmegaPlus;      //! casc. transv. momemtum Vs Omega+ Eff mass Vs Y
206 */      
207         // Compilation of all PID plots (3D = casc. transv. momemtum Vs Casc Eff mass Vs Y), stored into an AliCFContainer
208         AliCFContainer  *fCFContCascadePIDXiMinus;      //! for Xi-   : Container to store any 3D histos with the different PID flavours
209         AliCFContainer  *fCFContCascadePIDXiPlus;       //! for Xi+   : Container to store any 3D histos with the different PID flavours
210         AliCFContainer  *fCFContCascadePIDOmegaMinus;   //! for Omega-: Container to store any 3D histos with the different PID flavours
211         AliCFContainer  *fCFContCascadePIDOmegaPlus;    //! for Omega+: Container to store any 3D histos with the different PID flavours
212         
213         
214         
215         // PART 3 : Towards the optimisation of topological selections / systematics
216         AliCFContainer  *fCFContCascadeCuts;            //! Container meant to store all the relevant distributions corresponding to the cut variables
217         
218         
219         // PART 4 :  Azimuthal correlation study
220 /*      THnSparseF      *fHnSpAngularCorrXiMinus;       //! Delta Phi(Casc,any trck) Vs Delta Eta(Casc,any trck) Vs Casc Pt Vs Pt of the tracks Vs Eff Mass
221         THnSparseF      *fHnSpAngularCorrXiPlus;        //! Delta Phi(Casc,any trck) Vs Delta Eta(Casc,any trck) Vs Casc Pt Vs Pt of the tracks Vs Eff Mass
222         THnSparseF      *fHnSpAngularCorrOmegaMinus;    //! Delta Phi(Casc,any trck) Vs Delta Eta(Casc,any trck) Vs Casc Pt Vs Pt of the tracks Vs Eff Mass
223         THnSparseF      *fHnSpAngularCorrOmegaPlus;     //! Delta Phi(Casc,any trck) Vs Delta Eta(Casc,any trck) Vs Casc Pt Vs Pt of the tracks Vs Eff Mass
224 */
225         TH1F *fV0Ampl;                                  //! histo to check the V0 amplitude distribution  
226
227         TH2F    *fHistDcaXiDaughtersvsInvMass;          //! cut variables vs inv. mass
228         TH2F    *fHistDcaBachToPrimVertexvsInvMass;     //! cut variables vs inv. mass
229         TH2F    *fHistXiCosineOfPointingAnglevsInvMass; //! cut variables vs inv. mass
230         TH2F    *fHistMassLambdaAsCascDghtervsInvMass;  //! cut variables vs inv. mass
231         TH2F    *fHistDcaV0DaughtersXivsInvMass;        //! cut variables vs inv. mass
232         TH2F    *fHistDcaV0ToPrimVertexXivsInvMass;     //! cut variables vs inv. mass
233
234
235
236   AliAnalysisTaskCheckCascadePbPb(const AliAnalysisTaskCheckCascadePbPb&);            // not implemented
237   AliAnalysisTaskCheckCascadePbPb& operator=(const AliAnalysisTaskCheckCascadePbPb&); // not implemented
238   
239   ClassDef(AliAnalysisTaskCheckCascadePbPb, 4);
240 };
241
242 #endif