]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/RESONANCES/AliRsnMiniAnalysisTask.h
Recovered charge and phi values in AliRsnValueDaughter (overwritten by wrong commit)
[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                UseCentrality(const char *type)    {fUseCentrality = kTRUE; fCentralityType = type; fCentralityType.ToUpper();}
40    void                SetUseCentralityPatch(Bool_t isAOD049) {fUseAOD049CentralityPatch = isAOD049;}
41    void                UseMultiplicity(const char *type)  {fUseCentrality = kFALSE; fCentralityType = type; fCentralityType.ToUpper();}
42    void                UseContinuousMix()                 {fContinuousMix = kTRUE;}
43    void                UseBinnedMix()                     {fContinuousMix = kFALSE;}
44    void                SetNMix(Int_t nmix)                {fNMix = nmix;}
45    void                SetMaxDiffMult (Double_t val)      {fMaxDiffMult  = val;}
46    void                SetMaxDiffVz   (Double_t val)      {fMaxDiffVz    = val;}
47    void                SetMaxDiffAngle(Double_t val)      {fMaxDiffAngle = val;}
48    void                SetEventCuts(AliRsnCutSet *cuts)   {fEventCuts    = cuts;}
49    void                SetMixPrintRefresh(Int_t n)        {fMixPrintRefresh = n;}
50    Int_t               AddTrackCuts(AliRsnCutSet *cuts);
51    TClonesArray       *Outputs()                          {return &fHistograms;}
52    TClonesArray       *Values()                           {return &fValues;}
53    void                SetEventQAHist(TString type,TH2F *histo);
54
55    virtual void        UserCreateOutputObjects();
56    virtual void        UserExec(Option_t *);
57    virtual void        Terminate(Option_t *);
58    virtual void        FinishTaskOutput();
59
60    Int_t               ValueID(AliRsnMiniValue::EType type, Bool_t useMC = kFALSE);
61    Int_t               CreateValue(AliRsnMiniValue::EType type, Bool_t useMC = kFALSE);
62    AliRsnMiniOutput   *CreateOutput(const char *name, AliRsnMiniOutput::EOutputType type, AliRsnMiniOutput::EComputation src);
63    AliRsnMiniOutput   *CreateOutput(const char *name, const char *outType, const char *compType);
64
65 private:
66
67    Char_t   CheckCurrentEvent();
68    void     FillMiniEvent(Char_t evType);
69    Double_t ComputeAngle();
70    Double_t ComputeCentrality(Bool_t isESD);
71    Double_t ComputeMultiplicity(Bool_t isESD,TString type);
72    Double_t ApplyCentralityPatchAOD049();
73    void     FillTrueMotherESD(AliRsnMiniEvent *event);
74    void     FillTrueMotherAOD(AliRsnMiniEvent *event);
75    void     StoreTrueMother(AliRsnMiniPair *pair, AliRsnMiniEvent *event);
76    Bool_t   EventsMatch(AliRsnMiniEvent *event1, AliRsnMiniEvent *event2);
77
78    Bool_t               fUseMC;           //  use or not MC info
79    Int_t                fEvNum;           //! absolute event counter
80    Bool_t               fUseCentrality;   //  if true, use centrality for event, otherwise use multiplicity
81    TString              fCentralityType;  //  definition used to choose what centrality or multiplicity to use
82    Bool_t               fUseAOD049CentralityPatch; //flag to enable AOD049 centrality patch
83
84    Bool_t               fContinuousMix;   //  mixing --> technique chosen (continuous or binned)
85    Int_t                fNMix;            //  mixing --> required number of mixes
86    Double_t             fMaxDiffMult;     //  mixing --> max difference in multiplicity
87    Double_t             fMaxDiffVz;       //  mixing --> max difference in Vz of prim vert
88    Double_t             fMaxDiffAngle;    //  mixing --> max difference in reaction plane angle
89
90    TList               *fOutput;          //  output list
91    TClonesArray         fHistograms;      //  list of histogram definitions
92    TClonesArray         fValues;          //  list of values to be computed
93    TH1F                *fHEventStat;      //  histogram of event statistics
94    TH1F                *fHAEventsVsMulti; //  histogram of event statistics
95    TH2F                *fHAEventVz;       //  histogram of vertex-z vs. multiplicity/centrality
96    TH2F                *fHAEventMultiCent;//  histogram of multiplicity vs. centrality
97    TH2F                *fHAEventPlane;    //  histogram of event plane vs. multiplicity/centrality
98
99    AliRsnCutSet        *fEventCuts;       //  cuts on events
100    TObjArray            fTrackCuts;       //  list of single track cuts
101    AliRsnEvent          fRsnEvent;        //! interface object to the event
102    TTree               *fEvBuffer;        //! mini-event buffer
103    AliTriggerAnalysis  *fTriggerAna;      //! trigger analysis
104    AliESDtrackCuts     *fESDtrackCuts;    //! quality cut for ESD tracks
105    AliRsnMiniEvent     *fMiniEvent;       //! mini-event cursor
106    Bool_t               fBigOutput;       // flag if open file for output list
107    Int_t                fMixPrintRefresh; // how often info in mixing part is printed
108
109    ClassDef(AliRsnMiniAnalysisTask, 5);   // AliRsnMiniAnalysisTask
110 };
111
112 inline Int_t AliRsnMiniAnalysisTask::CreateValue(AliRsnMiniValue::EType type, Bool_t useMC)
113 {
114 //
115 // Create a new value in the task,
116 // and returns its ID, which is needed for setting up histograms.
117 // If that value was already initialized, returns its ID and does not recreate it.
118 //
119
120    Int_t valID = ValueID(type, useMC);
121    if (valID >= 0 && valID < fValues.GetEntries()) {
122       AliInfo(Form("Value '%s' is already created in slot #%d", AliRsnMiniValue::ValueName(type, useMC), valID));
123    } else {
124       valID = fValues.GetEntries();
125       AliInfo(Form("Creating value '%s' in slot #%d", AliRsnMiniValue::ValueName(type, useMC), valID));
126       new (fValues[valID]) AliRsnMiniValue(type, useMC);
127    }
128
129    return valID;
130 }
131
132 inline Int_t AliRsnMiniAnalysisTask::ValueID(AliRsnMiniValue::EType type, Bool_t useMC)
133 {
134 //
135 // Searches if a value computation is initialized
136 //
137
138    const char *name = AliRsnMiniValue::ValueName(type, useMC);
139    TObject *obj = fValues.FindObject(name);
140    if (obj)
141       return fValues.IndexOf(obj);
142    else
143       return -1;
144 }
145
146 inline AliRsnMiniOutput *AliRsnMiniAnalysisTask::CreateOutput
147 (const char *name, AliRsnMiniOutput::EOutputType type, AliRsnMiniOutput::EComputation src)
148 {
149 //
150 // Create a new histogram definition in the task,
151 // which is then returned to the user for its configuration
152 //
153
154    Int_t n = fHistograms.GetEntries();
155    AliRsnMiniOutput *newDef = new (fHistograms[n]) AliRsnMiniOutput(name, type, src);
156
157    return newDef;
158 }
159
160 inline AliRsnMiniOutput *AliRsnMiniAnalysisTask::CreateOutput
161 (const char *name, const char *outType, const char *compType)
162 {
163 //
164 // Create a new histogram definition in the task,
165 // which is then returned to the user for its configuration
166 //
167
168    Int_t n = fHistograms.GetEntries();
169    AliRsnMiniOutput *newDef = new (fHistograms[n]) AliRsnMiniOutput(name, outType, compType);
170
171    return newDef;
172 }
173 #endif