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