]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/testAna.C
JYMMY includes.
[u/mrichter/AliRoot.git] / ANALYSIS / testAna.C
1 #include "TNtuple.h"
2 #include "TH1.h"
3 #include "TFile.h"
4 #include "TCanvas.h"
5 #include "TSystem.h"
6 #include "AliAnalysisTask.h"
7 #include "AliAnalysisManager.h"
8 #include "AliAnalysisDataContainer.h"
9
10 //______________________________________________________________________________
11 class TaskPt : public AliAnalysisTask {
12
13 public:
14    TaskPt(const char *name);
15    virtual ~TaskPt() {;}
16    
17    virtual void              Exec(Option_t *option);
18    
19    ClassDef(TaskPt, 0); // example of analysis
20 };
21
22 //______________________________________________________________________________
23 class TaskCombine : public AliAnalysisTask {
24
25 public:
26    TaskCombine(const char *name);
27    virtual ~TaskCombine() {;}
28    
29    virtual void              Exec(Option_t *option);
30    
31    ClassDef(TaskCombine, 0); // example of combination
32 };
33
34 ClassImp(TaskPt)
35
36 //______________________________________________________________________________
37 TaskPt::TaskPt(const char *name)
38        :AliAnalysisTask(name,"")
39 {
40 // Constructor.
41    // Input slot #0 works with an Ntuple
42    DefineInput(0, TNtuple::Class());
43    // Output slot #0 writes into a TH1 container
44    DefineOutput(0, TH1F::Class());
45 }
46
47 //______________________________________________________________________________
48 void TaskPt::Exec(Option_t *)
49 {
50 // Task making a pt distribution.
51    static Int_t icalls=0;
52    icalls++;
53    TCanvas *c1 = new TCanvas(Form("c%i",icalls),"Dynamic Filling Example",200,10,700,500);
54    TH1F *h1 = new TH1F(Form("hpt%i",icalls),"This is the pt distribution",100,0,5);
55    // Get input data
56    TNtuple *ntuple = (TNtuple*)GetInputData(0);
57    if (!ntuple) {
58       printf("WOOPS ! Where is input 0 for %s ?\n", GetName());
59       return;
60    }
61    Float_t px, py, pt;
62    const Int_t kUPDATE = 1000;
63    ntuple->SetBranchAddress("px", &px);
64    ntuple->SetBranchAddress("py", &py);
65    Long64_t nentries = ntuple->GetEntries();
66    for (Long64_t i=0; i<nentries; i++) {
67       ntuple->GetEntry(i);
68       pt = TMath::Sqrt(px*px+py*py);
69       h1->Fill(pt);
70       if (i && (i%kUPDATE) == 0) {
71          if (i == kUPDATE) h1->Draw();
72          c1->Modified();
73          c1->Update();
74          if (gSystem->ProcessEvents())
75             break;
76       }
77    }
78    // Post final data. It will be written to a file with option "RECREATE"
79 //   h1->Draw();
80    PostData(0, h1, "RECREATE");
81 }      
82
83 ClassImp(TaskCombine)
84
85 //______________________________________________________________________________
86 TaskCombine::TaskCombine(const char *name)
87        :AliAnalysisTask(name,"")
88 {
89 // Constructor.
90    // Input slot #0 works with a TH1
91    DefineInput(0, TH1F::Class());
92    // Input slot #1 works with a TH1
93    DefineInput(1, TH1F::Class());
94    // Output slot #0 writes into a TH1 container
95    DefineOutput(0, TH1F::Class());
96 }
97    
98 //______________________________________________________________________________
99 void TaskCombine::Exec(Option_t *)
100 {
101 // Task combining 2 histograms.
102    new TCanvas("c2","h1+h2",200,10,700,500);
103    TH1F *h1 = (TH1F*)GetInputData(0);
104    TH1F *h2 = (TH1F*)GetInputData(1);
105    TH1F *hsum = new TH1F("hsum","h1+h2",100,0,5);
106    hsum->Add(h1);
107    hsum->Add(h2);
108    hsum->Draw();
109    PostData(0,hsum);
110 }   
111
112
113 //______________________________________________________________________________
114 void testAna()
115 {
116    // Make the analysis manager
117    TFile *f = 0;
118    AliAnalysisManager *mgr = new AliAnalysisManager();
119    // Make a task
120    AliAnalysisTask *task1 = new TaskPt("PtTask1");
121    mgr->AddTask(task1);
122    AliAnalysisTask *task2 = new TaskPt("PtTask2");
123    mgr->AddTask(task2);
124    AliAnalysisTask *task = new TaskCombine("TaskCombine");
125    mgr->AddTask(task);
126    // Create containers for input/output
127    AliAnalysisDataContainer *cinput = mgr->CreateContainer("cntuple", 
128                       TNtuple::Class(),AliAnalysisManager::kInputContainer);
129    AliAnalysisDataContainer *coutput1 = mgr->CreateContainer("chist1", TH1::Class());
130    AliAnalysisDataContainer *coutput2 = mgr->CreateContainer("chist2", TH1::Class());
131    AliAnalysisDataContainer *coutput = mgr->CreateContainer("chistsum", 
132                       TH1::Class(),AliAnalysisManager::kOutputContainer);
133    coutput1->SetFileName("output1.root");                   
134    coutput2->SetFileName("output2.root");                   
135    mgr->ConnectInput(task1,0,cinput);
136    mgr->ConnectInput(task2,0,cinput);
137    mgr->ConnectOutput(task1,0,coutput1);
138    mgr->ConnectOutput(task2,0,coutput2);
139    mgr->ConnectInput(task,0,coutput1);
140    mgr->ConnectInput(task,1,coutput2);
141    mgr->ConnectOutput(task,0,coutput);
142    // Open data
143    if (!gSystem->AccessPathName("hsimple.root")) {
144       f = new TFile("hsimple.root");
145       TNtuple *ntpl = (TNtuple*)f->Get("ntuple");
146       cinput->SetData(ntpl);
147    } else {
148       printf("FIRST run $ROOTSYS/tutorials/hsimple.C\n");
149       return;
150    }   
151    
152    if (mgr->InitAnalysis()) {
153       mgr->PrintStatus();
154       mgr->ExecAnalysis();
155    }   
156 }                         
157