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