]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/SPECTRA/ParticleEfficiency/SteerAnalysisTaskParticleEfficiency.C
Split: fix refs to AddTaskCentrality.C
[u/mrichter/AliRoot.git] / PWGLF / SPECTRA / ParticleEfficiency / SteerAnalysisTaskParticleEfficiency.C
1 SteerAnalysisTaskParticleEfficiency(const Char_t *inputfilename, Int_t maxFiles = kMaxInt, Int_t maxEv = kMaxInt)
2 {
3
4   /* include path for ACLic */
5   gSystem->AddIncludePath("-I$ALICE_ROOT/include");
6   gSystem->AddIncludePath("-I$ALICE_ROOT/TOF");
7   /* load libraries */
8   gSystem->Load("libANALYSIS");
9   gSystem->Load("libANALYSISalice");
10   /* build analysis task class */
11   gROOT->LoadMacro("AliAnalysisTaskParticleEfficiency.cxx+g");
12
13   /* setup input chain */
14   TString str = inputfilename;
15   const Char_t *filename;
16   TChain *chain = new TChain("esdTree");
17   if (str.EndsWith(".xml")) {
18     TGrid::Connect("alien://");
19     Info("SteerTaskParticleEfficiency", "reading data list from collection:");
20     TAlienCollection coll(inputfilename, maxFiles);
21     coll.Reset();
22     while (coll.Next()) {
23       filename = coll.GetTURL();
24       Info("SteerTaskParticleEfficiency", Form("%s", filename));
25       chain->Add(filename);
26     }
27   }
28   else if (str.EndsWith(".txt")) {
29     Info("SteerTaskParticleEfficiency", "reading data list from text file:");
30     ifstream is(inputfilename);
31     Char_t buf[4096];
32     while(!is.eof()) {
33       is.getline(buf, 4096);
34       if (is.eof()) break;
35       chain->Add(buf);
36       Info("SteerTaskParticleEfficiency", Form("%s", buf));
37     }
38     is.close();
39   }
40   else {
41     Info("SteerTaskParticleEfficiency", "single file:");
42     filename = inputfilename;
43     Info("SteerTaskParticleEfficiency", Form("%s", filename));
44     chain->Add(filename);
45   }
46   Info("SteerTaskParticleEfficiency", Form("chain is ready: %d events", chain->GetEntries()));
47
48   /* create analysis manager */
49   AliAnalysisManager *mgr = new AliAnalysisManager("ParticleEfficiency");
50
51   /* define input event handler */
52   AliESDInputHandler *esdh = new AliESDInputHandler();
53   esdh->SetReadFriends(kFALSE);
54   mgr->SetInputEventHandler(esdh);
55
56   /* define MC truth event handler */
57   AliMCEventHandler *mch = new AliMCEventHandler();
58   mgr->SetMCtruthEventHandler(mch);
59   
60   /* define output handler */
61   AliAODHandler *outputh = new AliAODHandler();
62   mgr->SetOutputEventHandler(outputh);
63   
64   /* add tasks */
65   gROOT->LoadMacro("$ALICE_ROOT/ANALYSIS/macros/AddTaskPhysicsSelection.C");
66   AliPhysicsSelectionTask* physSelTask = AddTaskPhysicsSelection(kTRUE);
67   gROOT->LoadMacro("$ALICE_ROOT/OADB/macros/AddTaskCentrality.C");
68   AliCentralitySelectionTask *centralityTask = AddTaskCentrality(); 
69   centralityTask->SetMCInput();
70   gROOT->LoadMacro("AddAnalysisTaskParticleEfficiency.C");
71   
72   AddAnalysisTaskParticleEfficiency("pi+");
73   AddAnalysisTaskParticleEfficiency("pi-");
74   AddAnalysisTaskParticleEfficiency("K+");
75   AddAnalysisTaskParticleEfficiency("K-");
76   AddAnalysisTaskParticleEfficiency("proton");
77   AddAnalysisTaskParticleEfficiency("antiproton");
78   AddAnalysisTaskParticleEfficiency("K*0");
79   AddAnalysisTaskParticleEfficiency("K*0_bar");
80   AddAnalysisTaskParticleEfficiency("phi");
81
82   /* start analysis */
83   mgr->SetDebugLevel(0);
84   if (!mgr->InitAnalysis()) return;
85   mgr->PrintStatus();
86   mgr->StartAnalysis("local", chain, maxEv);
87
88   /* create dummy file to tell we are done */
89   gSystem->Exec("touch done");
90
91 }