]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/vertexingHF/charmFlow/AddTaskFlowD2H.C
* Fix mass band for dstar analysis
[u/mrichter/AliRoot.git] / PWGHF / vertexingHF / charmFlow / AddTaskFlowD2H.C
1 class AliAnalysisDataContainer;
2 class AliFlowEventCuts;
3 class AliFlowEventTrackCuts;
4 class AliFlowEventTrackSimpleCuts;
5 class AliRDHFCutsD0toKpi;
6 class AliRDHFCutsDStartoKpipi; 
7
8 void AddTaskFlowD2H(TString fileNameCuts, TString folderName, Int_t nDmeson, Int_t myHarmonic, 
9                     Bool_t bDoQC, Bool_t bDoSPTPC, Bool_t bDoSPVZERO, Bool_t bDoEPTPC, Bool_t bDoEPVZERO, 
10                     Int_t ptBinWidth, Double_t gapTPC, Double_t etaVZERO1, Double_t etaVZERO2, Double_t etaVZERO3, Double_t etaVZERO4,
11                     Bool_t bOldApproach=kFALSE, Bool_t shrinkSP=kFALSE ) {
12   TFile *filecuts = TFile::Open( fileNameCuts.Data() );
13   if( (!filecuts) || ( filecuts && !filecuts->IsOpen()) ){
14     AliFatal("Could not open cuts file.");
15   }
16
17   TString fileName = AliAnalysisManager::GetCommonFileName();
18   fileName.ReplaceAll(".root","");
19
20   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
21   AliAnalysisDataContainer *cinput1 = mgr->GetCommonInputContainer();
22
23   TString thecuts = folderName;
24   TString sgn = DMesonName( nDmeson ) + Form("w%d",ptBinWidth);
25
26   //**********************************************************************
27   // FLOW TRACK CUTS
28   AliFlowTrackCuts* cutsRFPTPC = new AliFlowTrackCuts( "GlobalRFPTPC" );
29   cutsRFPTPC->SetParamType(AliFlowTrackCuts::kGlobal);
30   cutsRFPTPC->SetPtRange(0.2,5.);
31   cutsRFPTPC->SetEtaRange(-0.8,0.8);
32   cutsRFPTPC->SetMinNClustersTPC(70);
33   cutsRFPTPC->SetMinChi2PerClusterTPC(0.2);
34   cutsRFPTPC->SetMaxChi2PerClusterTPC(4.0);
35   cutsRFPTPC->SetAcceptKinkDaughters(kFALSE);
36   cutsRFPTPC->SetMinimalTPCdedx(10.);
37   cutsRFPTPC->SetAODfilterBit(1);
38   //  cutsRFPTPC->SetQA(kTRUE);
39   AliFlowTrackCuts* cutsRFPVZE = new AliFlowTrackCuts( "GlobalRFPVZE" );
40   cutsRFPVZE->SetParamType(AliFlowTrackCuts::kV0);
41   cutsRFPVZE->SetEtaRange(-10,+10);
42   cutsRFPVZE->SetPhiMin(0);
43   cutsRFPVZE->SetPhiMax(TMath::TwoPi());
44   //  cutsRFPVZE->SetQA(kTRUE);
45   AliFlowTrackCuts* cutsPOI = new AliFlowTrackCuts( "GlobalPOI" );
46   cutsPOI->SetParamType(AliFlowTrackCuts::kGlobal);
47   cutsPOI->SetPtRange(0.2,-5.);
48   cutsPOI->SetEtaRange(0.8,-0.8);
49
50   //**********************************************************************
51   // POI FILTER CUTS
52   AliFlowTrackSimpleCuts *filterPOIQC[50]; // MASS BANDS
53   AliFlowTrackSimpleCuts *filterPOISP[50][2]; // MASS BANDS || ETA || SUBEVENT GAP
54   int myMassBands = MassBands( nDmeson, bOldApproach );
55   for(int mb=0; mb!=myMassBands; ++mb) {
56     filterPOISP[mb][0] = new AliFlowTrackSimpleCuts( Form("FilterPOISP_MB%d_ETANEG",mb) );
57     filterPOISP[mb][0]->SetEtaMin( -0.8 );
58     filterPOISP[mb][0]->SetEtaMax( 0. );
59     filterPOISP[mb][0]->SetMassMin( MassBandLowEdge(nDmeson,mb,bOldApproach) );
60     filterPOISP[mb][0]->SetMassMax( MassBandLowEdge(nDmeson,mb+1,bOldApproach) );
61     filterPOISP[mb][1] = new AliFlowTrackSimpleCuts( Form("FilterPOISP_MB%d_ETAPOS",mb) );
62     filterPOISP[mb][1]->SetEtaMin( 0. );
63     filterPOISP[mb][1]->SetEtaMax( +0.8 );
64     filterPOISP[mb][1]->SetMassMin( MassBandLowEdge(nDmeson,mb,bOldApproach) );
65     filterPOISP[mb][1]->SetMassMax( MassBandLowEdge(nDmeson,mb+1,bOldApproach) );
66     filterPOIQC[mb] = new AliFlowTrackSimpleCuts( Form("FilterPOIQC_MB%d",mb) );
67     filterPOIQC[mb]->SetEtaMin( -0.8 );
68     filterPOIQC[mb]->SetEtaMax( +0.8 );
69     filterPOIQC[mb]->SetMassMin( MassBandLowEdge(nDmeson,mb,bOldApproach) );
70     filterPOIQC[mb]->SetMassMax( MassBandLowEdge(nDmeson,mb+1,bOldApproach) );
71   }
72
73   // * DMESON SELECTOR *************************************************
74   AliAnalysisTaskFlowD2H *taskSel;
75   switch (nDmeson) {
76   case ( AliRDHFCuts::kD0toKpiCuts ):
77     AliRDHFCutsD0toKpi *myCutsD0 = (AliRDHFCutsD0toKpi*)filecuts->Get("D0toKpiCuts");
78     if(!myCutsD0) {
79       AliFatal("Problems reaching D0toKpiCuts");
80     }
81     taskSel = new AliAnalysisTaskFlowD2H( Form("TaskD0Selector_%s",thecuts.Data()), 
82                                           cutsRFPTPC, cutsRFPVZE, myCutsD0, nDmeson );
83     break;
84   case ( AliRDHFCuts::kDstarCuts ):
85     AliRDHFCutsDStartoKpipi *myCutsDStar = (AliRDHFCutsDStartoKpipi*)filecuts->Get("DStartoKpipiCuts");
86     if(!myCutsDStar) {
87       AliFatal("Problems reaching DStarToKpipiCuts");
88     }
89     taskSel = new AliAnalysisTaskFlowD2H( Form("TaskDStarSelector_%s",thecuts.Data()),
90                                           cutsRFPTPC, cutsRFPVZE, myCutsDStar, nDmeson);
91     break;
92   case (AliRDHFCuts::kDplusCuts):
93     AliRDHFCutsDplustoKpipi *myCutsDplus = (AliRDHFCutsDplustoKpipi*)filecuts->Get("AnalysisCuts");
94     if(!myCutsDplus) {
95       AliFatal("Problems reaching AnalysisCuts");
96     }
97     taskSel = new AliAnalysisTaskFlowD2H( Form("TaskDplusSelector_%s",thecuts.Data()), 
98                                           cutsRFPTPC, cutsRFPVZE, myCutsDplus, nDmeson );
99     break;
100   }
101   taskSel->SetCommonConstants( MassBins(nDmeson), MinMass(nDmeson), MaxMass(nDmeson), ptBinWidth );
102
103   //taskSel->SelectCollisionCandidates(trigger);
104   //taskSel->SetDebug();
105   AliAnalysisDataContainer *coutHist = mgr->CreateContainer( Form("%sSelector_%s",sgn.Data(),thecuts.Data()),
106                                                              TList::Class(),AliAnalysisManager::kOutputContainer,
107                                                              Form("%s.root:FlowD2H_%s",fileName.Data(),thecuts.Data()) );
108   AliAnalysisDataContainer *exc_TPC = mgr->CreateContainer( Form("TPCEventWithCandidates_%s_%s",sgn.Data(),thecuts.Data()),
109                                                             AliFlowEventSimple::Class(),
110                                                             AliAnalysisManager::kExchangeContainer );
111   AliAnalysisDataContainer *exc_VZE = mgr->CreateContainer( Form("VZEEventWithCandidates_%s_%s",sgn.Data(),thecuts.Data()),
112                                                             AliFlowEventSimple::Class(),
113                                                             AliAnalysisManager::kExchangeContainer );
114   mgr->AddTask(taskSel);
115   mgr->ConnectOutput(taskSel,1,coutHist);
116   mgr->ConnectOutput(taskSel,2,exc_TPC);
117   mgr->ConnectOutput(taskSel,3,exc_VZE);
118   mgr->ConnectInput (taskSel,0,cinput1);
119
120   // * HANGING ANALYSIS TASKS ********************************************
121   int harm = myHarmonic;
122   for(int mb=0; mb!=myMassBands; ++mb) {
123     if(bDoQC) {
124       AddQCmethod( Form("%sQCTPCMB%d",sgn.Data(),mb), fileName.Data(), thecuts.Data(), harm, exc_TPC, filterPOIQC[mb]);
125     }
126     if(bDoSPVZERO) {
127       AddSPmethod( Form("%sSPVZEMB%d",sgn.Data(),mb), fileName.Data(), thecuts.Data(), etaVZERO1, etaVZERO2, etaVZERO3, etaVZERO4,
128                    "Qa", harm, exc_VZE, 0, filterPOIQC[mb], NULL, false, shrinkSP );
129       AddSPmethod( Form("%sSPVZEMB%d",sgn.Data(),mb), fileName.Data(), thecuts.Data(), etaVZERO1, etaVZERO2, etaVZERO3, etaVZERO4,
130                    "Qb", harm, exc_VZE, 0, filterPOIQC[mb], NULL, false, shrinkSP );
131     }
132     if(bDoEPVZERO) {
133       AddSPmethod( Form("%sEPVZEMB%d",sgn.Data(),mb), fileName.Data(), thecuts.Data(), etaVZERO1, etaVZERO2, etaVZERO3, etaVZERO4,
134                    "QaQb", harm, exc_VZE, 0, filterPOIQC[mb], NULL, true, shrinkSP );
135     }
136     if(bDoSPTPC) {
137       for(int eg=0; eg!=2; ++eg) {
138         AddSPmethod( Form("%sSPTPCMB%d",sgn.Data(),mb), fileName.Data(), thecuts.Data(), -0.8, -gapTPC*eg, +gapTPC*eg, +0.8,
139                      "Qa", harm, exc_TPC, eg, filterPOISP[mb][1], NULL, false, shrinkSP );
140         AddSPmethod( Form("%sSPTPCMB%d",sgn.Data(),mb), fileName.Data(), thecuts.Data(), -0.8, -gapTPC*eg, +gapTPC*eg, +0.8,
141                      "Qb", harm, exc_TPC, eg, filterPOISP[mb][0], NULL, false, shrinkSP );
142       }
143     }
144     if(bDoEPTPC) {
145       AddSPmethod( Form("%sEPTPCMB%d",sgn.Data(),mb), fileName.Data(), thecuts.Data(), -0.8, -0.0, +0.0, +0.8,
146                    "QaQb", harm, exc_TPC, 0, filterPOIQC[mb], NULL, true, shrinkSP );
147     }
148   }
149 }
150
151 void AddSPmethod(char *name, char *fileName, char *thecuts,
152                  double minEtaA, double maxEtaA, double minEtaB, double maxEtaB,
153                  char *Qvector, int harmonic, AliAnalysisDataContainer *flowEvent, int eg,
154                  AliFlowTrackSimpleCuts *cutsPOI=NULL, AliFlowTrackSimpleCuts *cutsRFP=NULL,
155                  bool bEP, bool shrink=false ) {
156   TString myFolder = Form("v%d_%s",harmonic,thecuts);
157   TString myNameSP = Form("%sSPv%d%sGAP%d_%s",name,harmonic,Qvector,eg,thecuts);
158   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
159   AliAnalysisDataContainer *flowEvent2 = mgr->CreateContainer( Form("Filter_%s", myNameSP.Data()),
160                                                                AliFlowEventSimple::Class(),
161                                                                AliAnalysisManager::kExchangeContainer );
162   AliAnalysisTaskFilterFE *tskFilter = new AliAnalysisTaskFilterFE( Form("TaskFilter_%s",myNameSP.Data()),
163                                                                     cutsRFP, cutsPOI);
164   tskFilter->SetSubeventEtaRange( minEtaA, maxEtaA, minEtaB, maxEtaB );
165   mgr->AddTask(tskFilter);
166   mgr->ConnectInput( tskFilter,0,flowEvent);
167   mgr->ConnectOutput(tskFilter,1,flowEvent2);
168   AliAnalysisDataContainer *outSP = mgr->CreateContainer( myNameSP.Data(),
169                                                           TList::Class(),AliAnalysisManager::kOutputContainer,
170                                                           Form("%s.root:FlowD2H_SP_%s",fileName,myFolder.Data()) );
171   AliAnalysisTaskScalarProduct *tskSP = new AliAnalysisTaskScalarProduct( Form("TaskScalarProduct_%s",
172                                                                                myNameSP.Data()),kFALSE);
173   tskSP->SetApplyCorrectionForNUA(kTRUE);
174   tskSP->SetHarmonic(harmonic);
175   tskSP->SetTotalQvector(Qvector);
176   if(bEP) tskSP->SetBehaveAsEP();
177   if(shrink) tskSP->SetBookOnlyBasicCCH(kTRUE);
178   mgr->AddTask(tskSP);
179   mgr->ConnectInput( tskSP,0,flowEvent2);
180   mgr->ConnectOutput(tskSP,1,outSP);
181 }
182
183 void AddQCmethod(char *name, char *fileName, char *thecuts,
184                  int harmonic, AliAnalysisDataContainer *flowEvent,
185                  AliFlowTrackSimpleCuts *cutsPOI=NULL, AliFlowTrackSimpleCuts *cutsRFP=NULL) {
186   TString myFolder = Form("v%d_%s",harmonic,thecuts);
187   TString myName = Form("%sv%d_%s",name,harmonic,thecuts);
188   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
189   AliAnalysisDataContainer *flowEvent2 = mgr->CreateContainer( Form("Filter_%s", myName.Data()),
190                                                                AliFlowEventSimple::Class(),
191                                                                AliAnalysisManager::kExchangeContainer );
192   AliAnalysisTaskFilterFE *tskFilter = new AliAnalysisTaskFilterFE( Form("TaskFilter_%s",myName.Data()),
193                                                                     cutsRFP, cutsPOI);
194   mgr->AddTask(tskFilter);
195   mgr->ConnectInput( tskFilter,0,flowEvent);
196   mgr->ConnectOutput(tskFilter,1,flowEvent2);
197   AliAnalysisDataContainer *outQC = mgr->CreateContainer( myName.Data(),
198                                                           TList::Class(),AliAnalysisManager::kOutputContainer,
199                                                           Form("%s.root:FlowD2H_QC_%s",fileName,myFolder.Data()) );
200   AliAnalysisTaskQCumulants *tskQC = new AliAnalysisTaskQCumulants( Form("TaskQCumulants_%s",
201                                                                          myName.Data()),kFALSE);
202   tskQC->SetApplyCorrectionForNUA(kTRUE);
203   tskQC->SetHarmonic(harmonic);
204   tskQC->SetBookOnlyBasicCCH(kTRUE);
205   mgr->AddTask(tskQC);
206   mgr->ConnectInput( tskQC,0,flowEvent2);
207   mgr->ConnectOutput(tskQC,1,outQC);
208 }
209
210 int MassBands( int nDmeson, bool bOldApproach=false ) {
211   switch (nDmeson) {
212   case ( AliRDHFCuts::kDplusCuts ):
213   case ( AliRDHFCuts::kD0toKpiCuts ):
214     if(bOldApproach) return 5;
215     else return 26;
216   case ( AliRDHFCuts::kDstarCuts ):
217     return 20;
218   }
219 }
220
221 double MassBandLowEdge( int nDmeson, int mb, bool bOldApproach=false ) {
222   switch (nDmeson) {
223   case ( AliRDHFCuts::kDplusCuts ):
224   case ( AliRDHFCuts::kD0toKpiCuts ): // 2 + 20 + 4
225     double lowEdgeMinimal[5+1] = {1.75,1.80,1.83,1.90,1.93,2.03};
226     double lowEdge[26+1] = { 1.66, 1.71, 1.76, 1.77, 1.78, 1.79, 1.80, 1.81, 1.82, 1.83,
227                              1.84, 1.85, 1.86, 1.87, 1.88, 1.89, 1.90, 1.91, 1.92, 1.93,
228                              1.94, 1.95, 1.96, 2.01, 2.06, 2.11, 2.16 };
229     if(bOldApproach) return lowEdgeMinimal[mb];
230     else return lowEdge[mb];
231   case ( AliRDHFCuts::kDstarCuts ): // 2 + 10 + 3
232     double lowEdge[20+1] = {0.1380, 0.1396, 0.1412, 0.1420, 0.1428, 0.1436, 0.1444, 0.1452, 0.1460, 0.1468, 0.1476, 0.1484, 0.1492, 0.1500,
233                             0.1508, 0.1516, 0.1524, 0.1532, 0.1548, 0.1564, 0.1580};
234     return lowEdge[mb];
235   }
236 }
237
238 double MinMass( int nDmeson ) {
239   switch (nDmeson) {
240   case ( AliRDHFCuts::kDplusCuts ):
241   case ( AliRDHFCuts::kD0toKpiCuts ):
242     return 1.66;
243   case ( AliRDHFCuts::kDstarCuts ):
244     return 0.138;
245   }
246 }
247
248 double MaxMass( int nDmeson ) {
249   switch (nDmeson) {
250   case ( AliRDHFCuts::kDplusCuts ):
251   case ( AliRDHFCuts::kD0toKpiCuts ):
252     return 2.16;
253   case ( AliRDHFCuts::kDstarCuts ):
254     return 0.158;
255   }
256 }
257
258 int MassBins( int nDmeson ) {
259   switch (nDmeson) {
260   case ( AliRDHFCuts::kDplusCuts ):
261   case ( AliRDHFCuts::kD0toKpiCuts ):
262     return 50;
263   case ( AliRDHFCuts::kDstarCuts ):
264     return 25;
265   }
266 }
267
268 TString DMesonName( int nDmeson ) {
269   TString toReturn;
270   switch (nDmeson) {
271   case ( AliRDHFCuts::kDplusCuts ):
272     toReturn = "DPlus"; break;
273   case ( AliRDHFCuts::kD0toKpiCuts ):
274     toReturn = "D0"; break;
275   case ( AliRDHFCuts::kDstarCuts ):
276     toReturn = "DStar"; break;
277   }
278   return toReturn;
279 }