]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/STRANGENESS/Cascades/AliAnalysisTaskCheckCascade.h
Cascade analysis code moved in Cascades folder
[u/mrichter/AliRoot.git] / PWGLF / STRANGENESS / Cascades / AliAnalysisTaskCheckCascade.h
1 #ifndef ALIANALYSISTASKCHECKCASCADE_H
2 #define ALIANALYSISTASKCHECKCASCADE_H
3
4 /*  See cxx source for full Copyright notice */
5
6 //-----------------------------------------------------------------
7 //                 AliAnalysisTaskCheckCascade class
8 //            (AliAnalysisTaskCheckCascade)
9 //            This task has four roles :
10 //              1. QAing the Cascades from ESD and AOD
11 //                 Origin:  AliAnalysisTaskESDCheckV0 by B.H. 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 //-----------------------------------------------------------------
18
19 class TList;
20 class TH1F;
21 class TH2F;
22 class TH3F;
23 class TVector3;
24 #include "THnSparse.h"
25  
26 class AliESDpid;
27 class AliESDEvent;
28 class AliESDtrackCuts;
29 class AliCFContainer;
30
31
32 #include "TString.h"
33
34 #include "AliAnalysisTaskSE.h"
35
36 class AliAnalysisTaskCheckCascade : public AliAnalysisTaskSE {
37  public:
38   AliAnalysisTaskCheckCascade();
39   AliAnalysisTaskCheckCascade(const char *name);
40   virtual ~AliAnalysisTaskCheckCascade();
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                                             Double_t  lRapCascade);
50   virtual Int_t  DoESDTrackWithTPCrefitMultiplicity(const AliESDEvent *lESDevent);
51   
52   virtual void   Terminate(Option_t *);
53   
54   void SetCollidingSystems           (Short_t collidingSystems          = 0    ) { fCollidingSystems = collidingSystems;}
55   void SetAnalysisType               (const char* analysisType          = "ESD") { fAnalysisType     = analysisType;}
56   void SetTriggerMaskType            (const char* triggerMaskType       = "kMB") { fTriggerMaskType  = triggerMaskType;}
57   void SetRelaunchV0CascVertexers    (Bool_t rerunV0CascVertexers       = 0    ) { fkRerunV0CascVertexers         =  rerunV0CascVertexers;      }
58   void SetQualityCutZprimVtxPos      (Bool_t qualityCutZprimVtxPos      = kTRUE) { fkQualityCutZprimVtxPos        =  qualityCutZprimVtxPos;     }
59   void SetRejectEventPileUp          (Bool_t rejectPileUp               = kTRUE) { fkRejectEventPileUp            =  rejectPileUp;              }
60   void SetQualityCutNoTPConlyPrimVtx (Bool_t qualityCutNoTPConlyPrimVtx = kTRUE) { fkQualityCutNoTPConlyPrimVtx   =  qualityCutNoTPConlyPrimVtx;}
61   void SetQualityCutTPCrefit         (Bool_t qualityCutTPCrefit         = kTRUE) { fkQualityCutTPCrefit           =  qualityCutTPCrefit;        }
62   void SetQualityCut80TPCcls         (Bool_t qualityCut80TPCcls         = kTRUE) { fkQualityCut80TPCcls           =  qualityCut80TPCcls;        }
63   void SetAlephParamFor1PadTPCCluster(Bool_t onePadTPCCluster           = kTRUE) { fkIsDataRecoWith1PadTPCCluster = onePadTPCCluster ;          }
64   void SetExtraSelections            (Bool_t extraSelections            = 0    ) { fkExtraSelections              =  extraSelections;           }
65   void SetAngularCorrelationType     (const char* correlationType       = "TrigLeadingTrck-AssoCasc") {fAngularCorrelationType = correlationType; }
66
67
68  private:
69         // Note : In ROOT, "//!" means "do not stream the data from Master node to Worker node" ...
70         // your data member object is created on the worker nodes and streaming is not needed.
71         // http://root.cern.ch/download/doc/11InputOutput.pdf, page 14
72
73
74         TString         fAnalysisType;                  // "ESD" or "AOD" analysis type
75         TString         fTriggerMaskType;               // type of trigger to use in UserExec : AliVEvent::kMB or kHighMult ?
76         Short_t         fCollidingSystems;              // 0 = pp collisions or 1 = AA collisions
77
78         AliESDpid       *fESDpid;                       // Tool data member to manage the TPC Bethe-Bloch info
79         AliESDtrackCuts *fESDtrackCuts;                 // ESD track cuts used for primary track definition
80         //TPaveText       *fPaveTextBookKeeping;          // TString to store all the relevant info necessary for book keeping (v0 cuts, cascade cuts, quality cuts, ...)
81
82         Bool_t          fkRerunV0CascVertexers;         // Boolean : kTRUE = relaunch both V0 + Cascade vertexers
83         Bool_t          fkQualityCutZprimVtxPos;        // Boolean : kTRUE = cut on the prim.vtx  z-position
84         Bool_t          fkRejectEventPileUp;            // Boolean : kTRUE = enable the rejection of events tagged as pile-up by SPD (AliESDEvent::IsPileupFromSPD)
85         Bool_t          fkQualityCutNoTPConlyPrimVtx;   // Boolean : kTRUE = prim vtx should be SPD or Tracking vertex
86         Bool_t          fkQualityCutTPCrefit;           // Boolean : kTRUE = ask for TPCrefit for the 3 daughter tracks
87         Bool_t          fkQualityCut80TPCcls;           // Boolean : kTRUE = ask for 80 TPC clusters for each daughter track
88         Bool_t          fkIsDataRecoWith1PadTPCCluster; // Boolean : kTRUE = data reconstructed with the "one pad cluster" algo in TPC (important for the ALEPH param for TPC PID)
89         Bool_t          fkExtraSelections;              // Boolean : kTRUE = apply tighter selections, before starting the analysis
90         TString         fAngularCorrelationType;        // type of correlation to compute : "TrigAnyCasc-AssoAnyPrim", "TrigCascLeading-AssoAnyPrim", "TrigLeadingTrck-AssoCasc"
91         
92         Double_t        fAlephParameters[5];            // Array to store the 5 param values for the TPC Bethe-Bloch parametrisation
93         Double_t        fV0Sels[7];                     // Array to store the 7 values for the different selections V0 related (if fkRerunV0CascVertexers)
94         Double_t        fCascSels[8];                   // Array to store the 8 values for the different selections Casc. related (if fkRerunV0CascVertexers)
95
96
97              TList      *fListHistCascade;              //! List of Cascade histograms
98         
99         // - General histos (filled before the trigger selection)
100         TH1F    *fHistCascadeMultiplicityBeforeTrigSel; //! Cascade multiplicity distribution
101          
102         // - General histos (filled for any triggered event)
103         TH1F    *fHistCascadeMultiplicityForTrigEvt;              //! Cascade multiplicity distribution
104         TH1F    *fHistTrackMultiplicityForTrigEvt;                //! Track multiplicity distribution (without any cut = include ITS stand-alone + global tracks)
105         TH1F    *fHistTPCrefitTrackMultiplicityForTrigEvt;        //! Track multiplicity distribution for tracks with TPCrefit
106         TH1F    *fHistPrimaryTrackMultiplicityForTrigEvt;         //! Track multiplicity distribution for the primary tracks (See fESDtrackCuts)
107         TH1F    *fHistEstimateITSTPCMultiplicityForTrigEvt;       //! (ITS+TPC tracks + SPD tracklets) multiplicity estimation from Ruben's work
108         
109         // - General histos (filled for any triggered event + (|z(prim. vtx)| < 10 cm ) )
110         TH1F    *fHistCascadeMultiplicityForTrigEvtAndZprimVtx;   //! Cascade multiplicity distribution
111         
112         // - General histos (filled for any triggered event + (|z(prim. vtx)| < 10 cm ) + after pile-up rejection from SPD )
113         TH1F    *fHistCascadeMultiplicityForTrigEvtNonPiledUpAndZprimVtx;   //! Cascade multiplicity distribution
114
115         // - General histos (filled for events selected in this analysis (|z(prim. vtx)| < 10 cm + prim vtx = SPD or Tracking) )
116         TH1F    *fHistCascadeMultiplicityForSelEvt;     //! Cascade multiplicity distribution
117         TH1F    *fHistPosBestPrimaryVtxXForSelEvt;      //! (best) primary vertex position distribution in x 
118         TH1F    *fHistPosBestPrimaryVtxYForSelEvt;      //! (best) primary vertex position distribution in y
119         TH1F    *fHistPosBestPrimaryVtxZForSelEvt;      //! (best) primary vertex position distribution in z
120         
121         
122         
123
124         // - Characteristics for event with >1 cascade : Track Multiplicity, TPC clusters + Prim. vertex positions
125         TH1F    *fHistTPCrefitTrackMultiplicityForCascadeEvt;   //! TPCrefit Track multiplicity distribution for event with >1 cascade candidate (NB: after quality sel. within the task)
126         TH1F    *fHistPrimaryTrackMultiplicityForCascadeEvt;    //! Track multiplicity distribution for the primary tracks for event with >1 cascade candidate (See fESDtrackCuts)
127
128         TH1F    *fHistPosV0TPCClusters;                 //! TPC clusters distribution for Positive V0 daughter track
129         TH1F    *fHistNegV0TPCClusters;                 //! TPC clusters distribution for Negative V0 daughter track
130         TH1F    *fHistBachTPCClusters;                  //! TPC clusters distribution for Bachelor             track
131
132         TH1F    *fHistVtxStatus;                        //! Is there a tracking vertex in the cascade event ?
133
134                 // Vtx coming from the full tracking, for events containing at least a cascade
135         TH1F    *fHistPosTrkgPrimaryVtxXForCascadeEvt;  //! primary vertex position distribution in x 
136         TH1F    *fHistPosTrkgPrimaryVtxYForCascadeEvt;  //! primary vertex position distribution in y
137         TH1F    *fHistPosTrkgPrimaryVtxZForCascadeEvt;  //! primary vertex position distribution in z
138         TH1F    *fHistTrkgPrimaryVtxRadius;             //! primary vertex (3D) radius distribution 
139
140                 // Best primary Vtx available, for events containing at least a cascade
141         TH1F    *fHistPosBestPrimaryVtxXForCascadeEvt;  //! (best) primary vertex position distribution in x 
142         TH1F    *fHistPosBestPrimaryVtxYForCascadeEvt;  //! (best) primary vertex position distribution in y
143         TH1F    *fHistPosBestPrimaryVtxZForCascadeEvt;  //! (best) primary vertex position distribution in z
144         TH1F    *fHistBestPrimaryVtxRadius;             //! (best) primary vertex radius distribution 
145
146                 // Correlation Best Vtx / Full Tracking Vtx
147         TH2F    *f2dHistTrkgPrimVtxVsBestPrimVtx;       //!  Radius of prim. Vtx from tracks Vs Radius of best Prim. Vtx
148
149
150 // PART 1 : Adavanced QA
151 // - Typical histos on the variables used for the selection of cascades
152         TH1F    *fHistEffMassXi;                        //! reconstructed cascade effective mass
153         TH1F    *fHistChi2Xi;                           //! chi2 value
154         TH1F    *fHistDcaXiDaughters;                   //! dca between Xi's daughters
155         TH1F    *fHistDcaBachToPrimVertex;              //! dca of the bachelor track to primary vertex
156         TH1F    *fHistXiCosineOfPointingAngle;          //! cosine of Xi pointing angle in a cascade
157         TH1F    *fHistXiRadius;                         //! (transverse) radius of the cascade vertex 
158                 
159         // - Histos about ~ the "V0 selection part" of the cascade,  coming by inheritance from AliESDv0
160         TH1F    *fHistMassLambdaAsCascDghter;           //! Test Invariant Mass of Lambda coming from Cascade
161         TH1F    *fHistV0Chi2Xi;                         //! V0 chi2 distribution, for the V0 associated to a cascade
162         TH1F    *fHistDcaV0DaughtersXi;                 //! Dca between V0 daughters, for the V0 associated to a cascade
163         TH1F    *fHistDcaV0ToPrimVertexXi;              //! Dca of V0 to primary vertex, for the V0 associated to a cascade     
164         TH1F    *fHistV0CosineOfPointingAngleXi;        //! Cosine of V0 pointing angle, for the V0 associated to a cascade
165         TH1F    *fHistV0RadiusXi;                       //! V0 (transverse) distance distribution, for the V0 associated to a cascade
166
167         TH1F    *fHistDcaPosToPrimVertexXi;             //! Dca of V0 positive daughter to primary vertex, for the V0 associated to a cascade
168         TH1F    *fHistDcaNegToPrimVertexXi;             //! Dca of V0 negative daughter to primary vertex, for the V0 associated to a cascade
169         
170
171         // - Effective mass histos for cascades.
172         TH1F    *fHistMassXiMinus;                      //! reconstructed cascade effective mass, under Xi- hyp.
173         TH1F    *fHistMassXiPlus;                       //! reconstructed cascade effective mass, under Xi+ hyp.
174         TH1F    *fHistMassOmegaMinus;                   //! reconstructed cascade effective mass, under Omega- hyp.
175         TH1F    *fHistMassOmegaPlus;                    //! reconstructed cascade effective mass, under Omega+ hyp.
176         
177         TH1F    *fHistMassWithCombPIDXiMinus;           //! reconstructed Xi- effective mass, with bach. comb PID
178         TH1F    *fHistMassWithCombPIDXiPlus;            //! reconstructed Xi+ effective mass, with bach. comb PID
179         TH1F    *fHistMassWithCombPIDOmegaMinus;        //! reconstructed Omega- effective mass, with bach. comb PID
180         TH1F    *fHistMassWithCombPIDOmegaPlus;         //! reconstructed Omega+ effective mass, with bach. comb PID
181
182         // - Complements for QA
183         TH1F    *fHistXiTransvMom;                      //! Xi transverse momentum, around the mass peak of Xi-/+ 
184         TH1F    *fHistXiTotMom;                         //! Xi momentum norm, around the mass peak of Xi-/+
185         
186         TH1F    *fHistBachTransvMomXi;                  //! bachelor transverse momentum, for cand. around the mass peak of Xi-/+
187         TH1F    *fHistBachTotMomXi;                     //! bachelor momentum norm, for cand. around the mass peak of Xi-/+
188
189         TH1F    *fHistChargeXi;                         //! Charge sign of the cascade candidate
190         TH1F    *fHistV0toXiCosineOfPointingAngle;      //! Cos. of Pointing angle between the V0 mom and the Xi-V0 vtx line
191   
192         TH1F    *fHistRapXi;                            //! rapidity of Xi candidates, around the mass peak of Xi-/+
193         TH1F    *fHistRapOmega;                         //! rapidity of Omega candidates, around the mass peak of Omega-/+
194         TH1F    *fHistEtaXi;                            //! eta distrib. of all the cascade candidates, around the mass peak of Xi-/+
195         TH1F    *fHistThetaXi;                          //! theta distrib. of all the cascade candidates, around the mass peak of Xi-/+
196         TH1F    *fHistPhiXi;                            //! phi distrib. of all the cascade candidates, around the mass peak of Xi-/+
197         
198         TH1F    *fHistcTauXiMinus;                      //! lifetime c.Tau, around the mass peak of Xi-
199         TH1F    *fHistcTauXiPlus;                       //! lifetime c.Tau, around the mass peak of Xi+
200         TH1F    *fHistcTauOmegaMinus;                   //! lifetime c.Tau, around the mass peak of Omega-
201         TH1F    *fHistcTauOmegaPlus;                    //! lifetime c.Tau, around the mass peak of Omega+
202         
203         TH2F    *f2dHistArmenteros;                     //! alpha(casc. cand.) Vs PtArm(casc. cand.)
204         
205         TH2F    *f2dHistEffMassLambdaVsEffMassXiMinus;  //! Xi- Eff mass Vs V0 Eff mass, under Xi- hyp.
206         TH2F    *f2dHistEffMassXiVsEffMassOmegaMinus;   //! Xi- Eff mass Vs Omega- Eff mass, for negative cascades
207         TH2F    *f2dHistEffMassLambdaVsEffMassXiPlus;   //! Xi+ Eff mass Vs V0 Eff mass, under Xi+ hyp.
208         TH2F    *f2dHistEffMassXiVsEffMassOmegaPlus;    //! Xi+ Eff mass Vs Omega+ Eff mass, for positive cascades
209         
210         TH2F    *f2dHistXiRadiusVsEffMassXiMinus;       //! transv. casc. decay radius Vs Xi- Eff mass, under Xi- hyp.
211         TH2F    *f2dHistXiRadiusVsEffMassXiPlus;        //! transv. casc. decay radius Vs Xi+ Eff mass, under Xi+ hyp.
212         TH2F    *f2dHistXiRadiusVsEffMassOmegaMinus;    //! transv. casc. decay radius Vs Omega- Eff mass, under Omega- hyp.
213         TH2F    *f2dHistXiRadiusVsEffMassOmegaPlus;     //! transv. casc. decay radius Vs Omega+ Eff mass, under Omega+ hyp.
214         
215         TH2F    *f2dHistTPCdEdxOfCascDghters;           //! TPC Bethe-Bloch curve, populated with the cascade daughters
216         
217         
218         // PART 2 : TH3F needed for pt spectrum and yield extraction
219         // Without any PID
220         TH3F    *f3dHistXiPtVsEffMassVsYXiMinus;        //! casc. transv. momemtum Vs Xi- Eff mass Vs Y
221         TH3F    *f3dHistXiPtVsEffMassVsYXiPlus;         //! casc. transv. momemtum Vs Xi+ Eff mass Vs Y
222         TH3F    *f3dHistXiPtVsEffMassVsYOmegaMinus;     //! casc. transv. momemtum Vs Omega- Eff mass Vs Y
223         TH3F    *f3dHistXiPtVsEffMassVsYOmegaPlus;      //! casc. transv. momemtum Vs Omega+ Eff mass Vs Y
224         
225         // Compilation of all PID plots (3D = casc. transv. momemtum Vs Casc Eff mass Vs Y), stored into an AliCFContainer
226         AliCFContainer  *fCFContCascadePIDXiMinus;      //! for Xi-   : Container to store any 3D histos with the different PID flavours
227         AliCFContainer  *fCFContCascadePIDXiPlus;       //! for Xi+   : Container to store any 3D histos with the different PID flavours
228         AliCFContainer  *fCFContCascadePIDOmegaMinus;   //! for Omega-: Container to store any 3D histos with the different PID flavours
229         AliCFContainer  *fCFContCascadePIDOmegaPlus;    //! for Omega+: Container to store any 3D histos with the different PID flavours
230         
231         
232         
233         // PART 3 : Towards the optimisation of topological selections / systematics
234         AliCFContainer  *fCFContCascadeCuts;            //! Container meant to store all the relevant distributions corresponding to the cut variables
235         
236         
237         // PART 4 :  Azimuthal correlation study
238         THnSparseF      *fHnSpAngularCorrXiMinus;       //! Delta Phi(Casc,any trck) Vs Delta Eta(Casc,any trck) Vs Casc Pt Vs Pt of the tracks Vs Eff Mass
239         THnSparseF      *fHnSpAngularCorrXiPlus;        //! Delta Phi(Casc,any trck) Vs Delta Eta(Casc,any trck) Vs Casc Pt Vs Pt of the tracks Vs Eff Mass
240         THnSparseF      *fHnSpAngularCorrOmegaMinus;    //! Delta Phi(Casc,any trck) Vs Delta Eta(Casc,any trck) Vs Casc Pt Vs Pt of the tracks Vs Eff Mass
241         THnSparseF      *fHnSpAngularCorrOmegaPlus;     //! Delta Phi(Casc,any trck) Vs Delta Eta(Casc,any trck) Vs Casc Pt Vs Pt of the tracks Vs Eff Mass
242
243
244   AliAnalysisTaskCheckCascade(const AliAnalysisTaskCheckCascade&);            // not implemented
245   AliAnalysisTaskCheckCascade& operator=(const AliAnalysisTaskCheckCascade&); // not implemented
246   
247   ClassDef(AliAnalysisTaskCheckCascade, 13);
248 };
249
250 #endif