]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliAnalysisManager.h
Merging THbtp and HBTP in one library. Comiplation on Windows/Cygwin
[u/mrichter/AliRoot.git] / ANALYSIS / AliAnalysisManager.h
1 #ifndef ALIANALYSISMANAGER_H
2 #define ALIANALYSISMANAGER_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 /* $Id$ */
7 // Author: Andrei Gheata, 31/05/2006
8
9 //==============================================================================
10 //   AliAnalysysManager - Manager analysis class. Allows creation of several
11 // analysis tasks and data containers storing their input/output. Allows 
12 // connecting/chaining tasks via shared data containers. Serializes the current
13 // event for all tasks depending only on initial input data.
14 //==============================================================================
15
16 #ifndef ROOT_TNamed
17 #include <TNamed.h>
18 #endif
19
20 class TClass;
21 class TTree;
22 class AliAnalysisDataContainer;
23 class AliAnalysisTask;
24 class AliVEventHandler;
25
26
27 class AliAnalysisManager : public TNamed {
28
29 public:
30
31 enum EAliAnalysisContType {
32    kExchangeContainer  = 0,
33    kInputContainer   = 1,
34    kOutputContainer  = 2
35 };   
36
37 enum EAliAnalysisExecMode {
38    kLocalAnalysis    = 0,
39    kProofAnalysis    = 1,
40    kGridAnalysis     = 2
41 };
42
43 enum EAliAnalysisFlags {
44    kEventLoop        = BIT(14),
45    kDisableBranches  = BIT(15),
46    kUseDataSet       = BIT(16)
47 };   
48
49    AliAnalysisManager(const char *name = "mgr", const char *title="");
50    virtual            ~AliAnalysisManager();
51
52    AliAnalysisManager(const AliAnalysisManager& other);
53    AliAnalysisManager& operator=(const AliAnalysisManager& other);
54    
55    // Management
56    void                StartAnalysis(const char *type="local", TTree *tree=0, Long64_t nentries=1234567890, Long64_t firstentry=0);
57    void                StartAnalysis(const char *type, const char *dataset, Long64_t nentries=1234567890, Long64_t firstentry=0);
58
59    virtual void        Init(TTree *tree);   
60    virtual Bool_t      Notify();
61    virtual void        SlaveBegin(TTree *tree);
62    virtual Bool_t      ProcessCut(Long64_t entry) {return Process(entry);}
63    virtual Bool_t      Process(Long64_t entry);
64    virtual Int_t       GetEntry(Long64_t entry, Int_t getall = 0);
65    void                PackOutput(TList *target);
66    void                UnpackOutput(TList *source);
67    virtual void        Terminate();
68
69    // Getters/Setters
70    static AliAnalysisManager *GetAnalysisManager() {return fgAnalysisManager;}
71    TObjArray          *GetContainers() const {return fContainers;}
72    UInt_t              GetDebugLevel() const {return fDebug;}
73    TObjArray          *GetInputs() const     {return fInputs;}
74    TObjArray          *GetOutputs() const    {return fOutputs;}
75    TObjArray          *GetTasks() const      {return fTasks;}
76    TObjArray          *GetTopTasks() const   {return fTopTasks;}
77    TTree              *GetTree() const       {return fTree;}
78    TObjArray          *GetZombieTasks() const {return fZombies;}
79    Long64_t            GetCurrentEntry() const {return fCurrentEntry;}
80    EAliAnalysisExecMode 
81                        GetAnalysisType() const {return fMode;}
82    Bool_t              IsUsingDataSet() const  {return TObject::TestBit(kUseDataSet);}
83
84    void                SetAnalysisType(EAliAnalysisExecMode mode) {fMode = mode;}
85    void                SetCurrentEntry(Long64_t entry) {fCurrentEntry = entry;}
86    void                SetDebugLevel(UInt_t level) {fDebug = level;}
87    void                SetSpecialOutputLocation(const char *location) {fSpecialOutputLocation = location;}
88    void                SetDisableBranches(Bool_t disable=kTRUE) {TObject::SetBit(kDisableBranches,disable);}
89    void                SetCollectSysInfoEach(Int_t nevents=0) {fNSysInfo = nevents;}
90    void                SetInputEventHandler(AliVEventHandler*  handler)  {fInputEventHandler   = handler;}
91    void                SetOutputEventHandler(AliVEventHandler*  handler) {fOutputEventHandler  = handler;}
92    void                SetMCtruthEventHandler(AliVEventHandler* handler) {fMCtruthEventHandler = handler;}
93    void                SetNSysInfo(Long64_t nevents) {fNSysInfo = nevents;}
94    AliVEventHandler*   GetInputEventHandler()   {return fInputEventHandler;}
95    AliVEventHandler*   GetOutputEventHandler()  {return fOutputEventHandler;}
96    AliVEventHandler*   GetMCtruthEventHandler() {return fMCtruthEventHandler;}
97
98    // Container handling
99    AliAnalysisDataContainer *CreateContainer(const char *name, TClass *datatype, 
100                        EAliAnalysisContType type     = kExchangeContainer, 
101                        const char          *filename = NULL);
102    
103    // Including tasks and getting them
104    void                 AddTask(AliAnalysisTask *task);
105    AliAnalysisTask     *GetTask(const char *name) const;
106    
107    // Connecting data containers to task inputs/outputs
108    Bool_t               ConnectInput(AliAnalysisTask *task, Int_t islot,
109                                      AliAnalysisDataContainer *cont);
110    Bool_t               ConnectOutput(AliAnalysisTask *task, Int_t islot,
111                                      AliAnalysisDataContainer *cont);
112    // Garbage collection
113    void                 CleanContainers();
114    
115    // Analysis initialization and execution, status
116    Bool_t               InitAnalysis();
117    Bool_t               IsInitialized() const {return fInitOK;}
118    Bool_t               IsEventLoop() const {return TObject::TestBit(kEventLoop);}
119    void                 ResetAnalysis();
120    void                 ExecAnalysis(Option_t *option="");
121    void                 FinishAnalysis();
122    void                 PrintStatus(Option_t *option="all") const;
123
124 protected:
125    void                 ImportWrappers(TList *source);
126    void                 SetEventLoop(Bool_t flag=kTRUE) {TObject::SetBit(kEventLoop,flag);}
127
128 private:
129    TTree                  *fTree;                //! Input tree in case of TSelector model
130    AliVEventHandler       *fInputEventHandler;   //  Optional common input  event handler
131    AliVEventHandler       *fOutputEventHandler;  //  Optional common output event handler
132    AliVEventHandler       *fMCtruthEventHandler; //  Optional common MC Truth event handler
133    Long64_t                fCurrentEntry;        //! Current processed entry in the tree
134    Long64_t                fNSysInfo;            // Event frequency for collecting system information
135    EAliAnalysisExecMode    fMode;                // Execution mode
136    Bool_t                  fInitOK;              // Initialisation done
137    UInt_t                  fDebug;               // Debug level
138    TString                 fSpecialOutputLocation; // URL/path where the special outputs will be copied
139    TObjArray              *fTasks;               // List of analysis tasks
140    TObjArray              *fTopTasks;            // List of top tasks
141    TObjArray              *fZombies;             // List of zombie tasks
142    TObjArray              *fContainers;          // List of all containers
143    TObjArray              *fInputs;              // List of containers with input data
144    TObjArray              *fOutputs;             // List of containers with results
145
146    static AliAnalysisManager *fgAnalysisManager; //! static pointer to object instance
147    ClassDef(AliAnalysisManager,3)  // Analysis manager class
148 };   
149 #endif