]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/RESONANCES/AliRsnMiniAnalysisTask.h
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGLF / RESONANCES / AliRsnMiniAnalysisTask.h
1 #ifndef ALIRSNMINIANALYSISTASK_H
2 #define ALIRSNMINIANALYSISTASK_H
3
4 //
5 // Analysis task for 'mini' sub-package
6 // Contains all definitions needed for running an analysis:
7 // -- global event cut
8 // -- a list of track cuts (any number)
9 // -- definitions of output histograms
10 // -- values to be computed.
11 //
12
13 #include <TString.h>
14 #include <TClonesArray.h>
15
16 #include "AliAnalysisTaskSE.h"
17
18 #include "AliRsnEvent.h"
19 #include "AliRsnMiniValue.h"
20 #include "AliRsnMiniOutput.h"
21
22 class TList;
23
24 class AliTriggerAnalysis;
25 class AliRsnMiniEvent;
26 class AliRsnCutSet;
27
28 class AliRsnMiniAnalysisTask : public AliAnalysisTaskSE {
29
30 public:
31
32    AliRsnMiniAnalysisTask();
33    AliRsnMiniAnalysisTask(const char *name, Bool_t isMC = kFALSE);
34    AliRsnMiniAnalysisTask(const AliRsnMiniAnalysisTask &copy);
35    AliRsnMiniAnalysisTask &operator=(const AliRsnMiniAnalysisTask &copy);
36    virtual ~AliRsnMiniAnalysisTask();
37
38    void                UseMC(Bool_t yn = kTRUE)           {fUseMC = yn;}
39    void                UseESDTriggerMask(UInt_t trgMask)     {fTriggerMask = trgMask;}
40    void                UseCentrality(const char *type)    {fUseCentrality = kTRUE; fCentralityType = type; fCentralityType.ToUpper();}
41    void                SetUseCentralityPatch(Bool_t isAOD049) {fUseAOD049CentralityPatch = isAOD049;}
42    void                SetUseCentralityPatchPbPb2011(Int_t centralityPatchPbPb2011) {fUseCentralityPatchPbPb2011 = centralityPatchPbPb2011;}
43    void                UseMultiplicity(const char *type)  {fUseCentrality = kFALSE; fCentralityType = type; fCentralityType.ToUpper();}
44    void                UseContinuousMix()                 {fContinuousMix = kTRUE;}
45    void                UseBinnedMix()                     {fContinuousMix = kFALSE;}
46    void                SetNMix(Int_t nmix)                {fNMix = nmix;}
47    void                SetMaxDiffMult (Double_t val)      {fMaxDiffMult  = val;}
48    void                SetMaxDiffVz   (Double_t val)      {fMaxDiffVz    = val;}
49    void                SetMaxDiffAngle(Double_t val)      {fMaxDiffAngle = val;}
50    void                SetEventCuts(AliRsnCutSet *cuts)   {fEventCuts    = cuts;}
51    void                SetMixPrintRefresh(Int_t n)        {fMixPrintRefresh = n;}
52    void                SetMaxNDaughters(Short_t n)        {fMaxNDaughters = n;}
53    void                SetCheckMomentumConservation(Bool_t checkP) {fCheckP = checkP;}
54    void                SetCheckFeedDown(Bool_t checkFeedDown)      {fCheckFeedDown = checkFeedDown;}
55    void                SetDselection(UShort_t originDselection);
56    void                SetRejectCandidateIfNotFromQuark(Bool_t opt){fRejectIfNoQuark=opt;}
57    void                SetMotherAcceptanceCutMinPt(Float_t minPt)  {fMotherAcceptanceCutMinPt = minPt;}
58    void                SetMotherAcceptanceCutMaxEta(Float_t maxEta){fMotherAcceptanceCutMaxEta = maxEta;}
59    void                KeepMotherInAcceptance(Bool_t keepMotherInAcceptance) {fKeepMotherInAcceptance = keepMotherInAcceptance;}
60    Int_t               AddTrackCuts(AliRsnCutSet *cuts);
61    TClonesArray       *Outputs()                          {return &fHistograms;}
62    TClonesArray       *Values()                           {return &fValues;}
63    Short_t             GetMaxNDaughters()                 {return fMaxNDaughters;}
64    void                SetEventQAHist(TString type,TH2F *histo);
65    void                UseBigOutput(Bool_t b=kTRUE) { fBigOutput = b; }
66
67    virtual void        UserCreateOutputObjects();
68    virtual void        UserExec(Option_t *);
69    virtual void        Terminate(Option_t *);
70    virtual void        FinishTaskOutput();
71
72    Int_t               ValueID(AliRsnMiniValue::EType type, Bool_t useMC = kFALSE);
73    Int_t               CreateValue(AliRsnMiniValue::EType type, Bool_t useMC = kFALSE);
74    AliRsnMiniOutput   *CreateOutput(const char *name, AliRsnMiniOutput::EOutputType type, AliRsnMiniOutput::EComputation src);
75    AliRsnMiniOutput   *CreateOutput(const char *name, const char *outType, const char *compType);
76
77 private:
78
79    Char_t   CheckCurrentEvent();
80    void     FillMiniEvent(Char_t evType);
81    Double_t ComputeAngle();
82    Double_t ComputeCentrality(Bool_t isESD);
83    Double_t ComputeMultiplicity(Bool_t isESD,TString type);
84    Double_t ComputeTracklets();
85    Double_t ApplyCentralityPatchAOD049();
86    Double_t ApplyCentralityPatchPbPb2011();
87    void     FillTrueMotherESD(AliRsnMiniEvent *event);
88    void     FillTrueMotherAOD(AliRsnMiniEvent *event);
89    void     StoreTrueMother(AliRsnMiniPair *pair, AliRsnMiniEvent *event);
90    Bool_t   EventsMatch(AliRsnMiniEvent *event1, AliRsnMiniEvent *event2);
91
92    Bool_t               fUseMC;           //  use or not MC info
93    Int_t                fEvNum;           //! absolute event counter
94    UInt_t               fTriggerMask;   //trigger mask
95    Bool_t               fUseCentrality;   //  if true, use centrality for event, otherwise use multiplicity
96    TString              fCentralityType;  //  definition used to choose what centrality or multiplicity to use
97    Bool_t               fUseAOD049CentralityPatch; //flag to enable AOD049 centrality patch
98    Int_t                fUseCentralityPatchPbPb2011; //for PbPb 2011 centrality flattening
99
100    Bool_t               fContinuousMix;   //  mixing --> technique chosen (continuous or binned)
101    Int_t                fNMix;            //  mixing --> required number of mixes
102    Double_t             fMaxDiffMult;     //  mixing --> max difference in multiplicity
103    Double_t             fMaxDiffVz;       //  mixing --> max difference in Vz of prim vert
104    Double_t             fMaxDiffAngle;    //  mixing --> max difference in reaction plane angle
105
106    TList               *fOutput;          //  output list
107    TClonesArray         fHistograms;      //  list of histogram definitions
108    TClonesArray         fValues;          //  list of values to be computed
109    TH1F                *fHEventStat;      //  histogram of event statistics
110    TH1F                *fHAEventsVsMulti; //  histogram of event statistics
111    TH1F                *fHAEventsVsTracklets; //  histogram of event statistics
112    TH2F                *fHAEventVz;       //  histogram of vertex-z vs. multiplicity/centrality
113    TH2F                *fHAEventMultiCent;//  histogram of multiplicity vs. centrality
114    TH2F                *fHAEventPlane;    //  histogram of event plane vs. multiplicity/centrality
115
116    AliRsnCutSet        *fEventCuts;       //  cuts on events
117    TObjArray            fTrackCuts;       //  list of single track cuts
118    AliRsnEvent          fRsnEvent;        //! interface object to the event
119    TTree               *fEvBuffer;        //! mini-event buffer
120    AliTriggerAnalysis  *fTriggerAna;      //! trigger analysis
121    AliESDtrackCuts     *fESDtrackCuts;    //! quality cut for ESD tracks
122    AliRsnMiniEvent     *fMiniEvent;       //! mini-event cursor
123    Bool_t               fBigOutput;       // flag if open file for output list
124    Int_t                fMixPrintRefresh; // how often info in mixing part is printed
125    Short_t              fMaxNDaughters;   // maximum number of allowed mother's daughter
126    Bool_t               fCheckP;          // flag to set in order to check the momentum conservation for mothers
127    
128    Bool_t               fCheckFeedDown;      // flag to set in order to check the particle feed down (specific for D meson analysis)
129    UShort_t             fOriginDselection;   // flag to select D0 origins. 0 Only from charm 1 only from beauty 2 both from charm and beauty (specific for D meson analysis)
130    Bool_t               fKeepDfromB;         // flag for the feed down from b quark decay (specific for D meson analysis)                       
131    Bool_t               fKeepDfromBOnly;     // flag to keep only the charm particles that comes from beauty decays (specific for D meson analysis)
132    Bool_t               fRejectIfNoQuark;    // flag to remove events not generated with PYTHIA
133    Float_t              fMotherAcceptanceCutMinPt;              // cut value to apply when selecting the mothers inside a defined acceptance
134    Float_t              fMotherAcceptanceCutMaxEta;             // cut value to apply when selecting the mothers inside a defined acceptance
135    Bool_t               fKeepMotherInAcceptance;                // flag to keep also mothers in acceptance
136
137    ClassDef(AliRsnMiniAnalysisTask, 11);   // AliRsnMiniAnalysisTask
138 };
139
140
141 #endif