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