]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliPriorsTask.C
Merge branch 'feature-doxygen'
[u/mrichter/AliRoot.git] / ANALYSIS / AliPriorsTask.C
1
2 Bool_t AliPriorsTask()
3 {
4   
5   TBenchmark benchmark;
6   benchmark.Start("AliPriorsTask");
7
8   AliLog::SetGlobalDebugLevel(0);
9
10   Load() ; //load the required libraries
11
12   TChain * analysisChain = new TChain("esdTree");
13   analysisChain->Add("AliESDs.root");
14
15   Int_t nIter = 5;
16
17   // create the task
18   AliPriorsTask *task = new AliPriorsTask("Task for Priors");
19   Double_t priors[5]={0.2,0.2,0.2,0.2,0.2};
20   task->SetPriors(priors);
21   task->SetNiterations(nIter);
22
23   // Make the analysis manager
24   AliAnalysisManager *mgr = new AliAnalysisManager("TestIteration");
25
26   AliMCEventHandler*  mcHandler = new AliMCEventHandler();
27   mgr->SetMCtruthEventHandler(mcHandler);
28
29
30    AliESDInputHandler* esdHandler = new AliESDInputHandler();
31    mgr->SetInputEventHandler(esdHandler);
32
33
34   // Create and connect containers for input/output
35   AliAnalysisDataContainer *cinput0  = mgr->GetCommonInputContainer();
36   //slot 0 : default output tree (by default handled by AliAnalysisTaskSE)
37   AliAnalysisDataContainer *coutput0 = mgr->CreateContainer("ctree0", TTree::Class(),AliAnalysisManager::kOutputContainer,"output.root");
38   // output TH1I for event counting
39   AliAnalysisDataContainer *coutput1 = mgr->CreateContainer("chist0", TH1I::Class(),AliAnalysisManager::kOutputContainer,"output.root");
40   AliAnalysisDataContainer *coutput2 = mgr->CreateContainer("cgraph0", TH2D::Class(),AliAnalysisManager::kOutputContainer,"output.root");
41
42
43   cinput0->SetData(analysisChain);
44
45   mgr->AddTask(task);
46   mgr->ConnectInput(task,0,cinput0);
47   mgr->ConnectOutput(task,0,coutput0);
48   mgr->ConnectOutput(task,1,coutput1);
49   mgr->ConnectOutput(task,2,coutput2);
50
51
52
53   AliEventPoolLoop* loop = new AliEventPoolLoop(nIter);
54   loop->SetChain(analysisChain);
55   mgr->SetEventPool(loop);
56
57   //mgr->SetDebugLevel(1); 
58   printf("READY TO RUN\n");
59   //RUN !!!
60   if (mgr->InitAnalysis()) {
61     mgr->PrintStatus();
62     mgr->StartAnalysis("mixing",analysisChain);
63     //mgr->StartAnalysis("local",analysisChain);
64     
65   }
66
67   benchmark.Stop("AliPriorsTask");
68   benchmark.Show("AliPriorsTask");
69
70   return kTRUE ;
71 }
72
73 void Load() {
74
75   //load the required aliroot libraries
76   gSystem->Load("libANALYSIS") ;
77   gSystem->Load("libANALYSISalice") ;
78
79   //compile online the task class
80   gSystem->SetIncludePath("-I. -I$ALICE_ROOT/include -I$ROOTSYS/include");
81   gROOT->LoadMacro("./AliPriorsTask.cxx++g");
82 }