]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/SPECTRA/AliAnalysisTaskCheckV0.h
Update for complying with train (first round)
[u/mrichter/AliRoot.git] / PWG2 / SPECTRA / AliAnalysisTaskCheckV0.h
1 #ifndef AliAnalysisTaskCheckV0_cxx
2 #define AliAnalysisTaskCheckV0_cxx
3
4 /*  See cxx source for full Copyright notice */
5
6 //-----------------------------------------------------------------
7 //                 AliAnalysisTaskCheckV0 class
8 //            This task is for QAing the V0s from ESD/AOD
9 //              Origin: B.H. Nov2007, hippolyt@in2p3.fr
10 //-----------------------------------------------------------------
11
12 class TString;
13 class TList;
14 class TH1F;
15 class AliESDEvent;
16 class AliAODEvent;
17
18 #include "AliAnalysisTaskSE.h"
19
20 class AliAnalysisTaskCheckV0 : public AliAnalysisTaskSE {
21  public:
22   AliAnalysisTaskCheckV0();
23   AliAnalysisTaskCheckV0(const char *name);
24   virtual ~AliAnalysisTaskCheckV0() {}
25   
26   virtual void   ConnectInputData(Option_t *);
27   virtual void   CreateOutputObjects();
28   virtual void   Exec(Option_t *option);
29   virtual void   Terminate(Option_t *);
30
31   void   SetCollidingSystems(Int_t collidingSystems = 0) {fCollidingSystems = collidingSystems;}
32   void   SetAnalysisType(const char* analysisType) {fAnalysisType = analysisType;}
33   
34  private:
35   AliESDEvent *fESD;                            //! ESD object
36   AliAODEvent *fAOD;                            //! AOD object
37   TString      fAnalysisType;                   //  ESD or AOD
38   Int_t        fCollidingSystems;               //  Colliding systems 0/1 for pp/PbPb  
39   TList       *fListHist;                       //! List of histograms
40   TH1F        *fHistPrimaryVertexPosX;          //! Primary vertex position in X
41   TH1F        *fHistPrimaryVertexPosY;          //! Primary vertex position in Y
42   TH1F        *fHistPrimaryVertexPosZ;          //! Primary vertex position in Z
43   TH1F        *fHistTrackMultiplicity;          //! Track multiplicity distribution
44   TH1F        *fHistV0Multiplicity;             //! V0 multiplicity distribution
45   TH1F        *fHistV0OnFlyStatus;              //! V0 on fly status distribution
46
47               // V0 offline distributions
48   TH1F        *fHistV0MultiplicityOff;          //! V0 multiplicity distribution offline
49   TH1F        *fHistV0Chi2Off;                  //! V0 chi2 distribution
50   TH1F        *fHistDcaV0DaughtersOff;          //! Dca between V0 daughters
51   TH1F        *fHistV0CosineOfPointingAngleOff; //! Cosine of V0 pointing angle
52   TH1F        *fHistV0RadiusOff;                //! V0 radial distance distribution
53   TH1F        *fHistDcaV0ToPrimVertexOff;       //! Dca of V0 to primary vertex
54   TH1F        *fHistDcaPosToPrimVertexOff;      //! Dca of V0 positive daughter to primary vertex
55   TH1F        *fHistDcaNegToPrimVertexOff;      //! Dca of V0 negative daughter to primary vertex
56
57   TH1F        *fHistMassK0Off;                  //! Invariant Mass of K0s
58   TH1F        *fHistMassLambdaOff;              //! Invariant Mass of Lambda
59   TH1F        *fHistMassAntiLambdaOff;          //! Invariant Mass of Anti-Lambda
60
61               // V0 on-the-fly distributions
62   TH1F        *fHistV0MultiplicityOn;           //! V0 multiplicity distribution offline
63   TH1F        *fHistV0Chi2On;                   //! V0 chi2 distribution
64   TH1F        *fHistDcaV0DaughtersOn;           //! Dca between V0 daughters
65   TH1F        *fHistV0CosineOfPointingAngleOn;  //! Cosine of V0 pointing angle
66   TH1F        *fHistV0RadiusOn;                 //! V0 radial distance distribution
67   TH1F        *fHistDcaV0ToPrimVertexOn;        //! Dca of V0 to primary vertex
68   TH1F        *fHistDcaPosToPrimVertexOn;       //! Dca of V0 positive daughter to primary vertex
69   TH1F        *fHistDcaNegToPrimVertexOn;       //! Dca of V0 negative daughter to primary vertex
70
71   TH1F        *fHistMassK0On;                   //! Invariant Mass of K0s
72   TH1F        *fHistMassLambdaOn;               //! Invariant Mass of Lambda
73   TH1F        *fHistMassAntiLambdaOn;           //! Invariant Mass of Anti-Lambda
74    
75   AliAnalysisTaskCheckV0(const AliAnalysisTaskCheckV0&);            // not implemented
76   AliAnalysisTaskCheckV0& operator=(const AliAnalysisTaskCheckV0&); // not implemented
77   
78   ClassDef(AliAnalysisTaskCheckV0, 0);
79 };
80
81 #endif