]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/STRANGENESS/Cascades/AliAnalysisTaskCheckCascadepp276.h
Fix for end-of-line style
[u/mrichter/AliRoot.git] / PWGLF / STRANGENESS / Cascades / AliAnalysisTaskCheckCascadepp276.h
1 #ifndef ALIANALYSISTASKCHECKCASCADEPP276_H
2 #define ALIANALYSISTASKCHECKCASCADEPP276_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 //                Adapted to Cascade : A.Maire Mar2008, antonin.maire@ires.in2p3.fr
15 //                Modified :           A.Maire Mar2010, antonin.maire@ires.in2p3.fr
16 //                Modified for PbPb analysis: M. Nicassio Feb 2011, maria.nicassio@ba.infn.it
17 //                Modified for pp@2.76 analysis: D. Colella Feb2012, domenico.colella@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 AliAnalysisTaskCheckCascadepp276 : public AliAnalysisTaskSE {
37  public:
38   AliAnalysisTaskCheckCascadepp276();
39   AliAnalysisTaskCheckCascadepp276(const char *name);
40   virtual ~AliAnalysisTaskCheckCascadepp276();
41   
42   virtual void   UserCreateOutputObjects();
43   virtual void   UserExec(Option_t *option);
44   virtual Int_t  DoESDTrackWithTPCrefitMultiplicity(const AliESDEvent *lESDevent);
45          //virtual Int_t  Tracks2V0vertices(AliESDEvent *event);  
46          //virtual Int_t  V0sTracks2CascadeVertices(AliESDEvent *event); 
47          //virtual Double_t Det(Double_t a00, Double_t a01, Double_t a10, Double_t a11) const;
48          //virtual Double_t Det(Double_t a00,Double_t a01,Double_t a02,
49          //      Double_t a10,Double_t a11,Double_t a12,
50          //      Double_t a20,Double_t a21,Double_t a22) const;
51
52          //virtual Double_t PropagateToDCA(AliESDv0 *vtx,AliExternalTrackParam *trk,Double_t b);
53   virtual void   Terminate(Option_t *);
54   
55   void SetAnalysisType               (const char* analysisType          = "ESD"  ) { fAnalysisType                = analysisType;               }
56   void SetRelaunchV0CascVertexers    (Bool_t rerunV0CascVertexers       = kFALSE ) { fkRerunV0CascVertexers       = rerunV0CascVertexers;       }
57   void SetSDDSelection               (Bool_t sddOnSelection             = kTRUE  ) { fkSDDSelectionOn             = sddOnSelection;             }
58   void SetQualityCutZprimVtxPos      (Bool_t qualityCutZprimVtxPos      = kTRUE  ) { fkQualityCutZprimVtxPos      = qualityCutZprimVtxPos;      }
59   void SetQualityCutNoTPConlyPrimVtx (Bool_t qualityCutNoTPConlyPrimVtx = kTRUE  ) { fkQualityCutNoTPConlyPrimVtx = qualityCutNoTPConlyPrimVtx; }
60   void SetQualityCutTPCrefit         (Bool_t qualityCutTPCrefit         = kTRUE  ) { fkQualityCutTPCrefit         = qualityCutTPCrefit;         }
61   void SetQualityCutnTPCcls          (Bool_t qualityCutnTPCcls          = kTRUE  ) { fkQualityCutnTPCcls          = qualityCutnTPCcls;          }
62   void SetQualityCutPileup           (Bool_t qualityCutPileup           = kTRUE  ) { fkQualityCutPileup           = qualityCutPileup;           }
63   void SetWithSDDOn                  (Bool_t withsddOn                  = kTRUE  ) { fwithSDD                     = withsddOn;                  }
64   void SetQualityCutMinnTPCcls       (Int_t  minnTPCcls                 = 70     ) { fMinnTPCcls                  = minnTPCcls;                 }
65   void SetExtraSelections            (Bool_t extraSelections            = kFALSE ) { fkExtraSelections            = extraSelections;            }
66   void SetVertexRange                (Float_t vtxrange                  = 10.0   ) { fVtxRange                    = vtxrange;                   }
67   void SetVertexRangeMin             (Float_t vtxrangemin               = 0.0    ) { fVtxRangeMin                 = vtxrangemin;                }
68   void SetMinptCutOnDaughterTracks   (Float_t minptdaughtrks            = 0.0    ) { fMinPtCutOnDaughterTracks    = minptdaughtrks;             }
69   void SetEtaCutOnDaughterTracks     (Float_t etadaughtrks              = 0.8    ) { fEtaCutOnDaughterTracks      = etadaughtrks;               }
70
71  private:
72         // Note : In ROOT, "//!" means "do not stream the data from Master node to Worker node" ...
73         // your data member object is created on the worker nodes and streaming is not needed.
74         // http://root.cern.ch/download/doc/11InputOutput.pdf, page 14
75
76
77         TString         fAnalysisType;                  // "ESD" or "AOD" analysis type 
78         AliESDtrackCuts *fESDtrackCuts;                 // ESD track cuts used for primary track definition
79         AliPIDResponse  *fPIDResponse;                  //! PID response object
80
81         Bool_t          fkRerunV0CascVertexers;         // Boolean : kTRUE = relaunch both V0 + Cascade vertexers
82         Bool_t          fkSDDSelectionOn;               // Boolena : kTRUE = select events with SDD on
83         Bool_t          fkQualityCutZprimVtxPos;        // Boolean : kTRUE = cut on the prim.vtx  z-position
84         Bool_t          fkQualityCutNoTPConlyPrimVtx;   // Boolean : kTRUE = prim vtx should be SPD or Tracking vertex
85         Bool_t          fkQualityCutTPCrefit;           // Boolean : kTRUE = ask for TPCrefit for the 3 daughter tracks
86         Bool_t          fkQualityCutnTPCcls;            // Boolean : kTRUE = ask for fMinnTPCcls TPC clusters for each daughter track
87         Bool_t          fkQualityCutPileup;             // Boolean : kTRUE = ask for No Pileup events
88         Bool_t          fwithSDD;                       // Boolean : kTRUE = select events with SDD reco
89         Int_t           fMinnTPCcls;                    // Minimum number of TPC cluster for daughter tracks
90         Bool_t          fkExtraSelections;              // Boolean : kTRUE = apply tighter selections, before starting the analysis
91         Float_t         fVtxRange;                      // to select events with |zvtx|<fVtxRange cm
92         Float_t         fVtxRangeMin;                   // to select events with |zvtx|>fVtxRangeMin cm
93         Float_t         fMinPtCutOnDaughterTracks;      // minimum pt cut on daughter tracks
94         Float_t         fEtaCutOnDaughterTracks;        // pseudorapidity cut on daughter tracks
95        
96         Double_t        fV0Sels[7];                     // Array to store the 7 values for the different selections V0 related (if fkRerunV0CascVertexers)
97         Double_t        fCascSels[8];                   // Array to store the 8 values for the different selections Casc. related (if fkRerunV0CascVertexers)
98
99         TList      *fListHistCascade;                   //! List of Cascade histograms
100         
101         // Cascades multiplicity plots
102         TH1F   *fHistCascadeMultiplicityBeforeAnySel;                 //! Cascade multiplicity distribution before any evnt selection 
103         TH1F   *fHistCascadeMultiplicityAfterSDDSel;                  //! Cascade multiplicity distribution after evnt selection on the SDD
104         TH1F   *fHistCascadeMultiplicityAfterPhysicsSel;              //! Cascade multiplicity distribution after evnt Physics Selection  
105         TH1F   *fHistCascadeMultiplicityForSelEvtNoTPCOnly;           //! Cascade multiplicity distribution after evnt noTPCOnly selection
106         TH1F   *fHistCascadeMultiplicityForSelEvtNoTPCOnlyNoPileup;   //! Cascade multiplicity distribution after evnt PileUp selection
107         TH1F   *fHistCascadeMultiplicityAfterVertexCutSel;            //! Cascade multiplicity distribution after evnt selection on the Z vertex position cut
108         // Tracks multiplicity plots
109         TH1F   *fHistTrackMultiplicityBeforeAnySel;                   //! Track multiplicity distribution before any evnt selection  
110         TH1F   *fHistTrackMultiplicityAfterSDDSel;                    //! Track multiplicity distribution after evnt selection on the SDD
111         TH1F   *fHistTrackMultiplicityAfterPhysicsSel;                //! Track multiplicity distribution after evnt Physics Selection
112         TH1F   *fHistTrackMultiplicityForSelEvtNoTPCOnly;             //! Track multiplicity distribution after evnt noTPCOnly selection
113         TH1F   *fHistTrackMultiplicityForSelEvtNoTPCOnlyNoPileup;     //! Track multiplicity distributionafter evnt PileUp selection
114         TH1F   *fHistTrackMultiplicityAfterVertexCutSel;              //! Track multiplicity distribution after evnt selection on the Z vertex position cut
115         // Vertex position plots (BestVertex)
116         TH1F   *fHistPVx;                                             //! Best primary vertex X position distribution after all evnt selection
117         TH1F   *fHistPVy;                                             //! Best primary vertex Y position distribution after all evnt selection
118         TH1F   *fHistPVz;                                             //! Best primary vertex Z position distribution after all evnt selection
119         TH1F   *fHistPVxAnalysis;                                     //! Best primary vertex X position distribution after all evnt selection and |z|>10cm cut
120         TH1F   *fHistPVyAnalysis;                                     //! Best primary vertex Y position distribution after all evnt selection and |z|>10cm cut    
121         TH1F   *fHistPVzAnalysis;                                     //! Best primary vertex Z position distribution after all evnt selection and |z|>10cm cut
122         // TPC cluster distributions for daughters
123         TH1F   *fHistPosV0TPCClusters;                                //! TPC clusters distribution for Positive V0 daughter track
124         TH1F   *fHistNegV0TPCClusters;                                //! TPC clusters distribution for Negative V0 daughter track
125         TH1F   *fHistBachTPCClusters;                                 //! TPC clusters distribution for Bachelor V0 daughter track
126         // Cut's variables distributions
127         TH1F   *fHistEffMassXi;                                       //! reconstructed cascade effective mass
128         TH1F   *fHistDcaXiDaughters;                                  //! dca between Xi's daughters
129         TH1F   *fHistDcaBachToPrimVertex;                             //! dca of the bachelor track to primary vertex
130         TH1F   *fHistXiCosineOfPointingAngle;                         //! cosine of Xi pointing angle in a cascade
131         TH1F   *fHistXiRadius;                                        //! (transverse) radius of the cascade vertex
132         TH1F   *fHistMassLambdaAsCascDghter;                          //! Test Invariant Mass of Lambda coming from Cascade 
133         TH1F   *fHistDcaV0DaughtersXi;                                //! Dca between V0 daughters, for the V0 associated to a cascade
134         TH1F   *fHistDcaV0ToPrimVertexXi;                             //! Dca of V0 to primary vertex, for the V0 associated to a cascade
135         TH1F   *fHistV0CosineOfPointingAngleXi;                       //! Cosine of V0 pointing angle, for the V0 associated to a cascade
136         TH1F   *fHistV0RadiusXi;                                      //! V0 (transverse) distance distribution, for the V0 associated to a cascade
137         TH1F   *fHistDcaPosToPrimVertexXi;                            //! Dca of V0 positive daughter to primary vertex, for the V0 associated to a cascade
138         TH1F   *fHistDcaNegToPrimVertexXi;                            //! Dca of V0 negative daughter to primary vertex, for the V0 associated to a cascade
139         // Invariant mass distributions
140         TH1F   *fHistMassXiMinus;                                     //! reconstructed cascade effective mass, under Xi- hyp.
141         TH1F   *fHistMassXiPlus;                                      //! reconstructed cascade effective mass, under Xi+ hyp.
142         TH1F   *fHistMassOmegaMinus;                                  //! reconstructed cascade effective mass, under Omega- hyp.
143         TH1F   *fHistMassOmegaPlus;                                   //! reconstructed cascade effective mass, under Omega+ hyp.
144         // Transverse and total momentum distributions
145         TH1F   *fHistXiTransvMom;                                     //! Xi transverse momentum, around the mass peak of Xi-/+
146         TH1F   *fHistXiTotMom;                                        //! Xi momentum norm, around the mass peak of Xi-/+
147         TH1F   *fHistBachTransvMomXi;                                 //! bachelor transverse momentum, for cand. around the mass peak of Xi-/+
148         TH1F   *fHistBachTotMomXi;                                    //! bachelor momentum norm, for cand. around the mass peak of Xi-/+
149         // Others QA plots
150         TH1F   *fHistChargeXi;                                        //! Charge sign of the cascade candidate
151         TH1F   *fHistV0toXiCosineOfPointingAngle;                     //! Cos. of Pointing angle between the V0 mom and the Xi-V0 vtx line
152         TH1F   *fHistRapXi;                                           //! rapidity of Xi candidates, around the mass peak of Xi-/+
153         TH1F   *fHistRapOmega;                                        //! rapidity of Omega candidates, around the mass peak of Omega-/+
154         TH1F   *fHistEtaXi;                                           //! eta distrib. of all the cascade candidates, around the mass peak of Xi-/+
155         TH1F   *fHistEtaBachXi;                                       
156         TH1F   *fHistEtaPosXi;                                        
157         TH1F   *fHistEtaNegXi;                                        
158         TH1F   *fHistThetaXi;                                         //! theta distrib. of all the cascade candidates, around the mass peak of Xi-/+
159         TH1F   *fHistPhiXi;                                           //! phi distrib. of all the cascade candidates, around the mass peak of Xi-/+
160         TH2F   *f2dHistArmenteros;                                    //! alpha(casc. cand.) Vs PtArm(casc. cand.)
161         TH2F   *f2dHistEffMassLambdaVsEffMassXiMinus;                 //! Xi- Eff mass Vs V0 Eff mass, under Xi- hyp.
162         TH2F   *f2dHistEffMassXiVsEffMassOmegaMinus;                  //! Xi- Eff mass Vs Omega- Eff mass, for negative cascades
163         TH2F   *f2dHistEffMassLambdaVsEffMassXiPlus;                  //! Xi+ Eff mass Vs V0 Eff mass, under Xi+ hyp. 
164         TH2F   *f2dHistEffMassXiVsEffMassOmegaPlus;                   //! Xi+ Eff mass Vs Omega+ Eff mass, for positive cascades
165         TH2F   *f2dHistXiRadiusVsEffMassXiMinus;                      //! transv. casc. decay radius Vs Xi- Eff mass, under Xi- hyp.
166         TH2F   *f2dHistXiRadiusVsEffMassXiPlus;                       //! transv. casc. decay radius Vs Xi+ Eff mass, under Xi+ hyp.
167         TH2F   *f2dHistXiRadiusVsEffMassOmegaMinus;                   //! transv. casc. decay radius Vs Omega- Eff mass, under Omega- hyp.
168         TH2F   *f2dHistXiRadiusVsEffMassOmegaPlus;                    //! transv. casc. decay radius Vs Omega+ Eff mass, under Omega+ hyp.
169         TH2F   *f2dHistTPCdEdxOfCascDghters;                          //! TPC Bethe-Bloch curve, populated with the cascade daughters
170         TH2F   *f2dHistDcaXiDaughtersvsInvMass;                       //! cut variables vs inv. mass
171         TH2F   *f2dHistDcaBachToPrimVertexvsInvMass;                  //! cut variables vs inv. mass
172         TH2F   *f2dHistXiCosineOfPointingAnglevsInvMass;              //! cut variables vs inv. mass
173         TH2F   *f2dHistMassLambdaAsCascDghtervsInvMass;               //! cut variables vs inv. mass
174         TH2F   *f2dHistDcaV0DaughtersXivsInvMass;                     //! cut variables vs inv. mass 
175         TH2F   *f2dHistDcaV0ToPrimVertexXivsInvMass;                  //! cut variables vs inv. mass 
176         // Containers for cuts study 
177         AliCFContainer  *fCFContCascadePIDXiMinus;                       //! for Xi-   : Container to store any 3D histos with the different PID flavours
178         AliCFContainer  *fCFContCascadePIDXiPlus;                        //! for Xi+   : Container to store any 3D histos with the different PID flavours
179         AliCFContainer  *fCFContCascadePIDOmegaMinus;                    //! for Omega-: Container to store any 3D histos with the different PID flavours
180         AliCFContainer  *fCFContCascadePIDOmegaPlus;                     //! for Omega+: Container to store any 3D histos with the different PID flavours
181         AliCFContainer  *fCFContCascadeCuts;                             //! Container meant to store all the relevant distributions corresponding to the cut variables
182
183
184   AliAnalysisTaskCheckCascadepp276(const AliAnalysisTaskCheckCascadepp276&);            // not implemented
185   AliAnalysisTaskCheckCascadepp276& operator=(const AliAnalysisTaskCheckCascadepp276&); // not implemented
186   
187   ClassDef(AliAnalysisTaskCheckCascadepp276, 7);
188 };
189
190 #endif