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