]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGHF/vertexingHF/macros/AddTaskDplus.C
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGHF / vertexingHF / macros / AddTaskDplus.C
CommitLineData
2fbb7748 1AliAnalysisTaskSEDplus *AddTaskDplus(Int_t system=0/*0=pp,1=PbPb*/,
2 Float_t minC=0, Float_t maxC=100,
3 Bool_t storeNtuple=kFALSE,
aa72cb06 4 Bool_t doSparse=kFALSE,
504df8ba 5 Bool_t readMC=kFALSE,
2fbb7748 6 TString finDirname="Loose",
0c22e2ac 7 TString filename="",
cee75703 8 TString finAnObjname="AnalysisCuts",
d94bafb6 9 Int_t etaRange=0,
10 Bool_t cutsDistr=kFALSE)
95e5b6b5 11{
12 //
13 // Test macro for the AliAnalysisTaskSE for D+ candidates
14
15 //Invariant mass histogram and
16 // association with MC truth (using MC info in AOD)
17 // R. Bala, bala@to.infn.it
18 // Get the pointer to the existing analysis manager via the static access method.
19 //==============================================================================
20 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
21 if (!mgr) {
22 ::Error("AddTaskDplus", "No analysis manager to connect to.");
95e5b6b5 23 }
24
ae39ad32 25 Bool_t stdcuts=kFALSE;
0c22e2ac 26 TFile* filecuts;
27 if( filename.EqualTo("") ) {
28 stdcuts=kTRUE;
29 } else {
30 filecuts=TFile::Open(filename.Data());
31 if(!filecuts ||(filecuts&& !filecuts->IsOpen())){
32 AliFatal("Input file not found : check your cut object");
33 }
95e5b6b5 34 }
fad29dfd 35
36
37 //Analysis Task
95e5b6b5 38
fad29dfd 39
ae39ad32 40 AliRDHFCutsDplustoKpipi* analysiscuts=new AliRDHFCutsDplustoKpipi();
2fbb7748 41 if(stdcuts) {
42 if(system==0) analysiscuts->SetStandardCutsPP2010();
43 else if(system==1){
0c22e2ac 44 analysiscuts->SetStandardCutsPbPb2011();
2fbb7748 45 analysiscuts->SetMinCentrality(minC);
46 analysiscuts->SetMaxCentrality(maxC);
47 // analysiscuts->SetUseAOD049(kTRUE);
48 analysiscuts->SetUseCentrality(AliRDHFCuts::kCentV0M);
49 }
50 }
51 else analysiscuts = (AliRDHFCutsDplustoKpipi*)filecuts->Get(finAnObjname);
fad29dfd 52
9e08baa8 53 AliAnalysisTaskSEDplus *dplusTask = new AliAnalysisTaskSEDplus("DplusAnalysis",analysiscuts,storeNtuple);
95e5b6b5 54 dplusTask->SetReadMC(readMC);
5fc4893f 55 dplusTask->SetDoLikeSign(kFALSE);
73173a6a 56 // dplusTask->SetUseTPCpid(kTRUE);
57 //dplusTask->SetUseTOFpid(kTRUE);
95e5b6b5 58 dplusTask->SetDebugLevel(0);
fad29dfd 59 dplusTask->SetMassLimits(0.2);
5fc4893f 60 dplusTask->SetUseBit(kTRUE);
d94bafb6 61 dplusTask->SetCutsDistr(cutsDistr);
aa7302d4 62 dplusTask->SetSystem(system);
aa72cb06 63 if (doSparse) dplusTask->SetDoImpactParameterHistos(kTRUE);
cee75703 64 if(etaRange==1) dplusTask->SetUseOnlyPositiveEta();
65 if(etaRange==-1) dplusTask->SetUseOnlyNegativeEta();
2fbb7748 66
95e5b6b5 67 mgr->AddTask(dplusTask);
fad29dfd 68
69 // Create containers for input/output
2fbb7748 70
71 TString inname = "cinputDplus";
72 TString outname = "coutputDplus";
73 TString cutsname = "coutputDplusCuts";
74 TString normname = "coutputDplusNorm";
75 TString ntuplename = "coutputDplus2";
76 inname += finDirname.Data();
77 outname += finDirname.Data();
78 cutsname += finDirname.Data();
79 normname += finDirname.Data();
80 ntuplename += finDirname.Data();
81 TString centr=Form("%.0f%.0f",analysiscuts->GetMinCentrality(),analysiscuts->GetMaxCentrality());
82 inname += centr;
83 outname += centr;
84 cutsname += centr;
85 normname += centr;
86 ntuplename += centr;
87
88
89 AliAnalysisDataContainer *cinputDplus = mgr->CreateContainer(inname,TChain::Class(),
fad29dfd 90 AliAnalysisManager::kInputContainer);
95e5b6b5 91 TString outputfile = AliAnalysisManager::GetCommonFileName();
92 outputfile += ":PWG3_D2H_InvMassDplus";
fad29dfd 93
2fbb7748 94 AliAnalysisDataContainer *coutputDplusCuts = mgr->CreateContainer(cutsname,TList::Class(),
fad29dfd 95 AliAnalysisManager::kOutputContainer,
96 outputfile.Data());
97
2fbb7748 98 AliAnalysisDataContainer *coutputDplus = mgr->CreateContainer(outname,TList::Class(),
fad29dfd 99 AliAnalysisManager::kOutputContainer,
95e5b6b5 100 outputfile.Data());
2fbb7748 101 AliAnalysisDataContainer *coutputDplusNorm = mgr->CreateContainer(normname,AliNormalizationCounter::Class(),
a96083b9 102 AliAnalysisManager::kOutputContainer,
103 outputfile.Data());
fad29dfd 104
95e5b6b5 105 if(storeNtuple){
2fbb7748 106 AliAnalysisDataContainer *coutputDplus2 = mgr->CreateContainer(ntuplename,TNtuple::Class(),
fad29dfd 107 AliAnalysisManager::kOutputContainer,
45bbb192 108 outputfile.Data());
fad29dfd 109
95e5b6b5 110 coutputDplus2->SetSpecialOutput();
111 }
112 mgr->ConnectInput(dplusTask,0,mgr->GetCommonInputContainer());
fad29dfd 113
95e5b6b5 114 mgr->ConnectOutput(dplusTask,1,coutputDplus);
fad29dfd 115
95e5b6b5 116 mgr->ConnectOutput(dplusTask,2,coutputDplusCuts);
a96083b9 117
118 mgr->ConnectOutput(dplusTask,3,coutputDplusNorm);
95e5b6b5 119 if(storeNtuple){
a96083b9 120 mgr->ConnectOutput(dplusTask,4,coutputDplus2);
95e5b6b5 121 }
122 return dplusTask;
123}