]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/Macros/RunAliComparisonTask.C
0937503a12af193968dc64cfb906eaf91a2bf7fb
[u/mrichter/AliRoot.git] / PWG1 / Macros / RunAliComparisonTask.C
1 void RunAliComparisonTask(TChain  *chain = 0, Bool_t aProof = kTRUE, Bool_t aDebug = kFALSE)
2 {
3   //
4   // Set mag field map (needed to propagate track to the DCA)
5   //
6   Int_t magField = 2;  // 0 - 0.2 T, 1 = 0.4 T, 2  - 0.5 T
7   magFMap = new AliMagFMaps("Maps","Maps", 2, 1., 10., magField);
8   AliTracker::SetFieldMap(magFMap,kFALSE);
9
10   //
11   // Create global cuts objects 
12   //
13
14   // Create ESD track reconstruction cuts
15   AliRecInfoCuts *pRecInfoCuts = new AliRecInfoCuts(); 
16   if(pRecInfoCuts) {
17     pRecInfoCuts->SetPtRange(0.15,200.0);
18     pRecInfoCuts->SetMaxAbsTanTheta(1.0);
19     pRecInfoCuts->SetMinNClustersTPC(10);
20     pRecInfoCuts->SetMinTPCsignalN(50);
21
22         pRecInfoCuts->SetHistogramsOn(kFALSE); 
23   } else {
24     AliDebug(AliLog::kError, "ERROR: Cannot create AliRecInfoCuts object");
25   }
26
27   // Create MC track reconstruction cuts
28   AliMCInfoCuts  *pMCInfoCuts = new AliMCInfoCuts();
29   if(pMCInfoCuts) {
30     pMCInfoCuts->SetMinRowsWithDigits(50);
31     pMCInfoCuts->SetMaxR(0.001);  
32     pMCInfoCuts->SetMaxVz(0.001); 
33     pMCInfoCuts->SetRangeTPCSignal(0.5,1.4); 
34   } else {
35     AliDebug(AliLog::kError, "ERROR: Cannot AliMCInfoCuts object");
36   }
37
38   //
39   // Create comparison objects and set cuts 
40   //
41
42   // Resolution
43   AliComparisonRes *pCompRes = new AliComparisonRes(); 
44   if(!pCompRes) {
45     AliDebug(AliLog::kError, "ERROR: Cannot create AliComparisonRes object");
46   }
47   pCompRes->SetAliRecInfoCuts(pRecInfoCuts);
48   pCompRes->SetAliMCInfoCuts(pMCInfoCuts);
49
50   // Efficiency
51   AliComparisonEff *pCompEff =  new AliComparisonEff();
52   if(!pCompEff) {
53     AliDebug(AliLog::kError, "ERROR: Cannot create AliComparisonEff object");
54   }
55   pCompEff->SetAliRecInfoCuts(pRecInfoCuts);
56   pCompEff->SetAliMCInfoCuts(pMCInfoCuts);
57
58   // dE/dx
59   AliComparisonDEdx *pCompDEdx = new AliComparisonDEdx();
60   if(!pCompDEdx) {
61     AliDebug(AliLog::kError, "ERROR: Cannot create AliComparisonDEdx object");
62   }
63   pCompDEdx->SetAliRecInfoCuts(pRecInfoCuts);
64   pCompDEdx->SetAliMCInfoCuts(pMCInfoCuts);
65   pCompDEdx->SetMCPtMin(0.5);
66   pCompDEdx->SetMCAbsTanThetaMax(0.5);
67   pCompDEdx->SetMCPdgCode(pMCInfoCuts->GetPiP()); // only pi+ particles
68
69   // DCA
70   AliComparisonDCA *pCompDCA = new AliComparisonDCA();
71   if(!pCompDCA) {
72     AliDebug(AliLog::kError, "ERROR: Cannot create AliComparisonDCA object");
73   }
74   pCompDCA->SetAliRecInfoCuts(pRecInfoCuts);
75   pCompDCA->SetAliMCInfoCuts(pMCInfoCuts);
76
77   // Create the analysis manager
78   mgr = new AliAnalysisManager("testAnalysis");
79
80   // Create, add task
81   task = new AliComparisonTask;
82
83   //task->SetMagField(magField);
84
85   task->AddComparisonObject( pCompRes );
86   task->AddComparisonObject( pCompEff );
87   task->AddComparisonObject( pCompDEdx );
88   task->AddComparisonObject( pCompDCA );
89   
90   mgr->AddTask(task);
91
92   // Attach input
93   cInput  = mgr->CreateContainer("cInput", TChain::Class(), AliAnalysisManager::kInputContainer);
94   mgr->ConnectInput(task, 0, cInput);
95
96   // Attach output
97   cOutput = mgr->CreateContainer("cOutput", TList::Class(), AliAnalysisManager::kOutputContainer,"Output.root");
98   mgr->ConnectOutput(task, 0, cOutput);
99
100   // Enable debug printouts
101   if (aDebug)
102     mgr->SetDebugLevel(2);
103
104   // Run analysis
105   mgr->InitAnalysis();
106   mgr->PrintStatus();
107
108   if(chain) {
109     mgr->StartAnalysis((aProof) ? "proof" : "local", chain);
110   } else {
111     AliDebug(AliLog::kError, "ERROR: No chain available");
112   }
113 }