]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/muon/AliCFMuonSingleTask1.C
Added LHC10h run list for flow analysis (Giacomo)
[u/mrichter/AliRoot.git] / PWG3 / muon / AliCFMuonSingleTask1.C
1 Bool_t AliCFMuonSingleTask1(Int_t runmin = 17, Int_t runmax = 17)
2 {
3     Bool_t isMC=kTRUE;
4     
5     TBenchmark benchmark;
6     benchmark.Start("AliMuonSingleTask1");
7     AliLog::SetGlobalDebugLevel(0);
8     Load() ; // load the required libraries
9     TChain * analysisChain ;
10     analysisChain = new TChain("esdTree");
11
12     Char_t RunFile[256];
13     for(Int_t i=runmin; i<=runmax; i++){
14         sprintf(RunFile,"/home/lopez/alice/data/TEST/run%d-10/AliESDs.root",i);
15         analysisChain->Add(RunFile);
16     }
17
18     enum             { kEta,  kY, kPhi, kPt, kP3, kHit, kChi2Fit, kTrM, kChi2TrM,  kContrN,  kVt,  kVz, kTrig, kDCA, kZcoor, kRabs, kCharge, kTheta,kNVars };
19     Int_t  nBins[] = {   5 , 5  ,  45 , 60 ,150 ,  20 ,      20 ,  4  ,      20 ,     202 , 100 ,  100 ,  10 , 500 ,   1000,    7 ,     3 , 100      };
20     Double_t min[] = {  -4.,-4. ,-180.,  0.,  0.,   0.,       0., -0.5,       0.,     -2.5,   0., -100.,   0.,   0., -3000.,  171.,  -1.5 , 2.95      };
21     Double_t max[] = {-2.5.,-2.5, 180., 30.,150.,  20.,      20.,  3.5,      10.,    199.5, 200.,  100.,  10., 500.,  1000.,  178.,   1.5 , 3.15     };
22     
23     Double_t *binLimits = 0;
24     Int_t nSteps=1; if (isMC) nSteps=2;
25     AliCFContainer* contCF = new AliCFContainer("container", "", nSteps, kNVars, nBins);
26     for (Int_t var=kNVars; var--;) {
27         binLimits = new Double_t[nBins[var]+1];
28         for (Int_t i=0; i<=nBins[var]; i++) binLimits[i]=min[var]+i*(max[var]-min[var])/nBins[var];
29         contCF->SetBinLimits(var, binLimits);
30         delete [] binLimits; binLimits=0;
31     }
32     
33     TList *qaList = new TList();
34     TObjArray *genList = new TObjArray(0);
35     AliCFTrackKineCuts *genKineCuts = new AliCFTrackKineCuts("genKineCuts", "GenKineCuts");
36     genKineCuts->SetPtRange(min[kPt], max[kPt]);
37     genKineCuts->SetRapidityRange(min[kY], max[kY]);
38     genKineCuts->SetQAOn(qaList);
39     genList->AddLast(genKineCuts);
40
41     TObjArray *recList = new TObjArray(0);
42     AliCFTrackKineCuts *recKineCuts = new AliCFTrackKineCuts("recKineCuts", "RecKineCuts");
43     recKineCuts->SetPtRange(min[kPt], max[kPt]);
44     recKineCuts->SetRapidityRange(min[kY], max[kY]);
45     recKineCuts->SetQAOn(qaList);
46     recList->AddLast(recKineCuts);
47 //__
48
49     AliCFManager* managerCF = new AliCFManager() ;
50     managerCF->SetParticleContainer(contCF);
51     managerCF->SetParticleCutsList(AliCFManager::kPartGenCuts, genList);
52     managerCF->SetParticleCutsList(AliCFManager::kPartAccCuts, recList);
53
54     AliCFMuonSingleTask1 *taskMuonCF = new AliCFMuonSingleTask1("AliMuonSingleTask1");
55     taskMuonCF->SetCFManager(managerCF);
56     taskMuonCF->SetQAList(qaList);
57     taskMuonCF->SetUseMC(isMC);
58
59     AliAnalysisManager *mgr = new AliAnalysisManager("TestManager");
60     mgr->SetAnalysisType(AliAnalysisManager::kLocalAnalysis);
61     
62     AliMCEventHandler*  mcHandler = new AliMCEventHandler();
63     mgr->SetMCtruthEventHandler(mcHandler);
64
65 //  AliInputEventHandler* dataHandler ;
66 //  dataHandler = new AliESDInputHandler();
67 //  mgr->SetInputEventHandler(dataHandler);
68
69     AliESDInputHandler *esdHandler = new AliESDInputHandler();
70     esdHandler->SetReadFriends(kFALSE);
71     mgr->SetInputEventHandler(esdHandler);
72
73
74     AliAnalysisDataContainer *cinput0  = mgr->CreateContainer("cchain0",TChain::Class(),AliAnalysisManager::kInputContainer);
75     cinput0->SetData(analysisChain);
76
77     mgr->AddTask(taskMuonCF);
78     mgr->ConnectInput(taskMuonCF, 0, mgr->GetCommonInputContainer());  
79
80     Char_t fileName[256];
81     sprintf(fileName,"muonCF_%d_%d.root",runmin,runmax);
82     printf("Analysis output in %s \n",fileName);
83
84     mgr->ConnectOutput(taskMuonCF,1,mgr->CreateContainer("chist",TH1I::Class(),AliAnalysisManager::kOutputContainer,fileName));
85     mgr->ConnectOutput(taskMuonCF,2,mgr->CreateContainer("ccont",AliCFContainer::Class(),AliAnalysisManager::kOutputContainer,fileName));
86  
87     printf("READY TO RUN\n");
88     //RUN !!!
89     if (mgr->InitAnalysis()) {
90         mgr->PrintStatus();
91         mgr->StartAnalysis("local",analysisChain);
92   }
93     benchmark.Stop("AliMuonSingleTask1");
94     benchmark.Show("AliMuonSingleTask1");
95
96     return kTRUE ;
97 }
98
99 void Load() {
100
101   //load the required aliroot libraries
102     gSystem->Load("libANALYSIS") ;
103     gSystem->Load("libANALYSISalice") ;
104     gSystem->Load("$ALICE_ROOT/lib/tgt_linux/libCORRFW.so") ;
105
106   //compile online the task class
107     gSystem->SetIncludePath("-I. -I$ALICE_ROOT/include -I$ALICE_ROOT/MUON -I$ALICE_ROOT/STEER -I$ROOTSYS/include");
108     gROOT->LoadMacro("./AliCFMuonSingleTask1.cxx+");
109 }