]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/correlationHF/macros/AddTaskD0Correlations.C
39e68c91228df65aef36f0dc4905e9763dafb5a2
[u/mrichter/AliRoot.git] / PWGHF / correlationHF / macros / AddTaskD0Correlations.C
1 AliAnalysisTaskSED0Correlations *AddTaskD0Correlations(Bool_t readMC=kFALSE, Int_t system=0/*0=pp,1=PbPb*/, Int_t flagD0D0bar=0,
2                                        Float_t minC=0, Float_t maxC=0, TString finDirname="Output", TString cutsfilename="D0toKpiCuts.root", 
3                                         TString cutsfilename2="AssocPartCuts.root", TString cutsD0name="D0toKpiCuts", TString cutsTrkname="AssociatedTrkCuts",                                          Bool_t flagAOD049=kFALSE, Int_t standardbins=1, Bool_t stdcuts=kFALSE)
4 {
5   //
6   // AddTask for the AliAnalysisTaskSE for D0 candidates
7   // invariant mass histogram and association with MC truth 
8   // (using MC info in AOD) and cut variables distributions
9   // C.Bianchin  chiara.bianchin@pd.infn.it
10   // F.Colamaria fabio.colamaria@ba.infn.it
11
12
13   // Get the pointer to the existing analysis manager via the static access method.
14   //==============================================================================
15   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
16   if (!mgr) {
17     ::Error("AddTaskD0Distr", "No analysis manager to connect to.");
18     return NULL;
19   }   
20
21   TString filename="",out1name="",out2name="",out3name="",out4name="",out5name="",out6name="",out7name="",inname="";
22   filename = AliAnalysisManager::GetCommonFileName();
23   filename += ":PWG3_D2H_";
24   filename+="D0InvMass";
25   if(flagD0D0bar==1)filename+="D0";
26   if(flagD0D0bar==2)filename+="D0bar";
27   //list mass
28   out1name="coutputmassD0Mass";
29   if(flagD0D0bar==1)out1name+="D0";
30   if(flagD0D0bar==2)out1name+="D0bar";
31   //hist entries
32   out2name="nEntriesD0";
33   if(flagD0D0bar==1)out2name+="D0";
34   if(flagD0D0bar==2)out2name+="D0bar";
35   //cuts object
36   out3name="cutsD0"; 
37   if(flagD0D0bar==1)out3name+="D0";
38   if(flagD0D0bar==2)out3name+="D0bar";
39   //AliNormalizationCounter
40   out4name="normalizationCounter";
41   if(flagD0D0bar==1)out4name+="D0";
42   if(flagD0D0bar==2)out4name+="D0bar";
43   //correlation outputs
44   out5name ="correlations";
45   if(flagD0D0bar==1)out5name+="D0";
46   if(flagD0D0bar==2)out5name+="D0bar";
47   //correlation further studies
48   out6name ="MCStudyPlots";
49   if(flagD0D0bar==1)out6name+="D0";
50   if(flagD0D0bar==2)out6name+="D0bar";
51   //correlated trk cuts
52   out7name ="cutsTracks";
53   if(flagD0D0bar==1)out7name+="D0";
54   if(flagD0D0bar==2)out7name+="D0bar";
55   inname="cinputmassD0_0";
56   if(flagD0D0bar==1)inname+="D0";
57   if(flagD0D0bar==2)inname+="D0bar";
58
59   filename += finDirname.Data();
60   out1name += finDirname.Data();
61   out2name += finDirname.Data(); 
62   out3name += finDirname.Data(); 
63   out4name += finDirname.Data(); 
64   out5name += finDirname.Data();
65   out6name += finDirname.Data();
66   out7name += finDirname.Data();
67   inname += finDirname.Data();
68
69     //setting my cut values
70
71     //cuts order
72     //     printf("    |M-MD0| [GeV]    < %f\n",fD0CorrCuts[0]);
73     //     printf("    dca    [cm]  < %f\n",fD0CorrCuts[1]);
74     //     printf("    cosThetaStar     < %f\n",fD0CorrCuts[2]);
75     //     printf("    pTK     [GeV/c]    > %f\n",fD0CorrCuts[3]);
76     //     printf("    pTpi    [GeV/c]    > %f\n",fD0CorrCuts[4]);
77     //     printf("    |d0K|  [cm]  < %f\n",fD0CorrCuts[5]);
78     //     printf("    |d0pi| [cm]  < %f\n",fD0CorrCuts[6]);
79     //     printf("    d0d0  [cm^2] < %f\n",fD0CorrCuts[7]);
80     //     printf("    cosThetaPoint    > %f\n",fD0CorrCuts[8]);
81
82   TFile* filecuts=new TFile(cutsfilename.Data());
83   if(!filecuts->IsOpen()){
84     cout<<"Input file not found for D0 cuts: using std cut object"<<endl;
85     stdcuts=kTRUE;
86   }
87   TFile* filecuts2=new TFile(cutsfilename2.Data());
88   if(!filecuts2->IsOpen()){
89     cout<<"Input file not found for tracks cuts!"<<endl;
90     return;
91   }
92
93   //Cuts for D0
94   AliRDHFCutsD0toKpi* RDHFD0Corrs=new AliRDHFCutsD0toKpi();
95 //RDHFD0Corrs->SetUsePhysicsSelection(kFALSE);
96   if(stdcuts) {
97     RDHFD0Corrs->SetLowPt(kFALSE); //low-pt special PID disabled
98     if(system==0) RDHFD0Corrs->SetStandardCutsPP2010();
99     else {
100       RDHFD0Corrs->SetStandardCutsPbPb2010();
101       if(minC!=0 && maxC!=0) { //if centrality 0 and 0 leave the values in the cut object
102         RDHFD0Corrs->SetMinCentrality(minC);
103         RDHFD0Corrs->SetMaxCentrality(maxC);
104       }
105       if(flagAOD049)RDHFD0Corrs->SetUseAOD049(kTRUE);
106       RDHFD0Corrs->SetUseCentrality(AliRDHFCuts::kCentV0M);
107     }
108   }
109   else {
110     RDHFD0Corrs = (AliRDHFCutsD0toKpi*)filecuts->Get(cutsD0name.Data());
111     if(!RDHFD0Corrs){
112       cout<<"Specific AliRDHFCuts not found"<<endl;
113       return;
114     }
115     if(flagAOD049)RDHFD0Corrs->SetUseAOD049(kTRUE);
116     if(minC!=0 && maxC!=0) { //if centrality 0 and 0 leave the values in the cut object
117       RDHFD0Corrs->SetMinCentrality(minC);
118       RDHFD0Corrs->SetMaxCentrality(maxC);
119     } 
120   }
121
122   //Cuts for correlated tracks/K0
123   AliHFAssociatedTrackCuts* corrCuts=new AliHFAssociatedTrackCuts();
124   corrCuts = (AliHFAssociatedTrackCuts*)filecuts2->Get(cutsTrkname.Data());
125   if(!corrCuts){
126       cout<<"Specific AliHFAssociatedTrackCuts not found"<<endl;
127       return;
128   }
129   corrCuts->PrintAll();
130
131   TString centr="";
132   if(minC!=0 && maxC!=0) centr = Form("%.0f%.0f",minC,maxC);
133   else centr = Form("%.0f%.0f",RDHFD0Corrs->GetMinCentrality(),RDHFD0Corrs->GetMaxCentrality());
134   out1name+=centr;
135   out2name+=centr;
136   out3name+=centr;
137   out4name+=centr;
138   out5name+=centr;
139   out6name+=centr;
140   out7name+=centr;
141   inname+=centr;
142
143   // Aanalysis task    
144   TString taskname="MassAndDistrAnalysis"; taskname.Prepend("D0");
145   AliAnalysisTaskSED0Correlations *massD0Task = new AliAnalysisTaskSED0Correlations(taskname.Data(),RDHFD0Corrs,corrCuts);
146   massD0Task->SetDebugLevel(2);
147   massD0Task->SetReadMC(readMC);
148   massD0Task->SetFillOnlyD0D0bar(flagD0D0bar);
149   massD0Task->SetSystem(system); //0=pp, 1=PbPb
150
151 //*********************
152 //correlation settings
153 //*********************
154
155   if(standardbins==1) {
156     printf("Standard bins (from D0Mass cuts object)\n");
157     massD0Task->SetNPtBinsCorr(RDHFD0Corrs->GetNPtBins()); 
158     massD0Task->SetPtBinsLimsCorr(RDHFD0Corrs->GetPtBinLimits());
159   } else {
160     Double_t ptlimits[15] = {0,0.5,1,2,3,4,5,6,7,8,12,16,20,24,9999};
161     massD0Task->SetNPtBinsCorr(15);
162     massD0Task->SetPtBinsLimsCorr(ptlimits);
163   }
164   Double_t pttreshlow[15] = {0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
165   Double_t pttreshup[15] = {999.,999.,999.,999.,999.,999.,999.,999.,999.,999.,999.,999.,999.,999.,999.};
166   massD0Task->SetPtTreshLow(pttreshlow);
167   massD0Task->SetPtTreshUp(pttreshup);
168   massD0Task->PrintBinsAndLimits();
169
170   //  massD0Task->SetRejectSDDClusters(kTRUE);
171   //  massD0Task->SetWriteVariableTree(kTRUE);
172
173   mgr->AddTask(massD0Task);
174   
175   //
176   // Create containers for input/output
177   AliAnalysisDataContainer *cinputmassD0 = mgr->CreateContainer(inname,TChain::Class(),AliAnalysisManager::kInputContainer);
178
179   AliAnalysisDataContainer *coutputmassD01 = mgr->CreateContainer(out1name,TList::Class(),AliAnalysisManager::kOutputContainer, filename.Data()); //mass
180   AliAnalysisDataContainer *coutputmassD02 = mgr->CreateContainer(out2name,TH1F::Class(),AliAnalysisManager::kOutputContainer, filename.Data()); //nev
181   AliAnalysisDataContainer *coutputmassD03 = mgr->CreateContainer(out3name,AliRDHFCutsD0toKpi::Class(),AliAnalysisManager::kOutputContainer, filename.Data()); //cuts
182   AliAnalysisDataContainer *coutputmassD04 = mgr->CreateContainer(out4name,AliNormalizationCounter::Class(),AliAnalysisManager::kOutputContainer, filename.Data()); //counter
183   AliAnalysisDataContainer *coutputmassD05 = mgr->CreateContainer(out5name,TList::Class(),AliAnalysisManager::kOutputContainer, filename.Data()); //correlations
184   AliAnalysisDataContainer *coutputmassD06 = mgr->CreateContainer(out6name,TList::Class(),AliAnalysisManager::kOutputContainer, filename.Data()); //MC study plots (for corrs)
185   AliAnalysisDataContainer *coutputmassD07 = mgr->CreateContainer(out7name,AliHFAssociatedTrackCuts::Class(),AliAnalysisManager::kOutputContainer, filename.Data()); //cuts for tracks/K0
186
187   mgr->ConnectInput(massD0Task,0,mgr->GetCommonInputContainer());
188
189   mgr->ConnectOutput(massD0Task,1,coutputmassD01);
190   mgr->ConnectOutput(massD0Task,2,coutputmassD02);
191   mgr->ConnectOutput(massD0Task,3,coutputmassD03);
192   mgr->ConnectOutput(massD0Task,4,coutputmassD04);
193   mgr->ConnectOutput(massD0Task,5,coutputmassD05);
194   mgr->ConnectOutput(massD0Task,6,coutputmassD06);
195   mgr->ConnectOutput(massD0Task,7,coutputmassD07);
196
197   return massD0Task;
198 }