AliPhysics  a60a912 (a60a912)
AddTaskForwardFlowQC.C
Go to the documentation of this file.
1 
37 void AddTaskForwardFlowQC(Int_t maxMom = 5,
38  TString fwdDet = "FMD",
39  Bool_t useEtaGap = kFALSE,
40  Bool_t use3cor = kFALSE,
41  Bool_t mc = kFALSE,
42  Double_t outlierCutFMD = 4.0,
43  Double_t outlierCutSPD = 4.0,
44  Double_t etaGap = 2.0,
45  UInt_t useTracksForRef = 0,
46  Bool_t useCent = kFALSE,
47  Bool_t ispA = kFALSE,
48  Bool_t useMCVtx = kFALSE,
49  Bool_t satVtx = kFALSE,
50  Bool_t addFlow = kFALSE)
51 {
52  // --- Get analysis manager ----------------------------------------
53  AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
54  if (!mgr)
55  Fatal("","No analysis manager to connect to.");
56 
57  // --- Check that we have an AOD input handler ---------------------
58  UShort_t aodInput = 0;
59  if (!(aodInput = AliForwardUtil::CheckForAOD()))
60  Fatal("","Cannot proceed without and AOD handler");
61  if (aodInput == 2 &&
62  (!AliForwardUtil::CheckForTask("AliForwardMultiplicityBase") ||
63  !AliForwardUtil::CheckForTask("AliCentralMultiplicityTask")))
64  Fatal("","The relevant tasks weren't added to the train");
65 
66  // --- For the selected flow tasks the input and output is set -----
67  fwdDet.ToUpper();
68 
69  // --- Set flow flags --------------------------------------------
70  if (useEtaGap && use3cor)
71  Fatal("", "You're doing it wrong! Cannot do both eta-gap and 3-sub");
73  if (useEtaGap) flags = AliForwardFlowTaskQC::kEtaGap;
74  if (use3cor) flags = AliForwardFlowTaskQC::k3Cor;
75  if (satVtx) flags |= AliForwardFlowTaskQC::kSatVtx;
76  if (fwdDet.Contains("FMD")) flags |= AliForwardFlowTaskQC::kNUAcorr;
77  if (fwdDet.Contains("FMD")) flags |= AliForwardFlowTaskQC::kFMD;
78  else if (fwdDet.Contains("VZERO")) flags |= AliForwardFlowTaskQC::kVZERO;
79  if (useTracksForRef == 1) flags |= AliForwardFlowTaskQC::kTPC;
80  if (useTracksForRef == 2) flags |= AliForwardFlowTaskQC::kHybrid;
81 
82  const char* name = Form("ForwardFlowQC%s%s", fwdDet.Data(), AliForwardFlowTaskQC::GetQCType(flags, false));
83  AliForwardFlowTaskQC* task = 0;
84  // --- Set up adding flow to MC input ----------------------------
85  if (mc) {
87  mcTask->SetUseImpactParameter(!useCent);
88  mcTask->SetUseMCHeaderVertex(useMCVtx);
89  if (addFlow) {
90  mcTask->SetUseFlowWeights(true);
91  }
92  task = mcTask;
93  }
94  // --- Or use normal task ----------------------------------------
95  else
96  task = new AliForwardFlowTaskQC(name);
97 
98  mgr->AddTask(task);
99 // mgr->SetSkipTerminate(true);
100 // task->SelectCollisionCandidates(AliVEvent::kCentral);
101  task->SetFlowFlags(flags);
102 
103  // --- Set eta gap value -----------------------------------------
104  if (useEtaGap) {
105  if (useTracksForRef) task->SetEtaGapValue(0.4);
106  else task->SetEtaGapValue(etaGap);
107  }
108  else if (useTracksForRef && fwdDet.Contains("FMD")) task->SetEtaGapValue(0.0);
109 
110  // --- Check which harmonics to calculate --------------------------
111  task->SetMaxFlowMoment(maxMom);
112 
113  // --- Set non-default axis for vertices ---------------------------
114  TAxis* a = 0;
115  if (satVtx) {
116  a = new TAxis(6, 93.75, 318.75);
117  }
118  else
119  a = new TAxis(20, -10, 10);
120 // a = new TAxis(10, -5, 5);
121  task->SetVertexAxis(a);
122 
123  // --- Set non-default axis for centrality -------------------------
124  TAxis* centAxis = 0;
125  if (ispA) {
126  Double_t cent[] = {0, 5, 10, 15, 20, 30, 40, 50, 60, 70, 80, 90, 100});
127  Int_t nBins = sizeof(cent)/sizeof(Double_t) -1;
128  centAxis = new TAxis(nBins, cent);
129  } else {
130  Double_t cent[] = {0, 2.5, 5, 7.5, 10, 12.5, 15, 20, 25, 30, 35, 40, 45, 50, 60, 70, 80, 90, 95, 100 };
131  Int_t nBins = sizeof(cent)/sizeof(Double_t) -1;
132  centAxis = new TAxis(nBins, cent);
133  }
134  task->SetCentralityAxis(centAxis);
135 
136  // --- Set sigma cuts for outliers ---------------------------------
137  task->SetDetectorCuts(outlierCutFMD, outlierCutSPD);
138 
139  // --- Setup track cuts --------------------------------------------
140  if (useTracksForRef > 0) {
141  AliAnalysisFilter* trackCuts = new AliAnalysisFilter("trackFilter");
142  if (!ispA) {
143  if (useTracksForRef == 1) { // tpc tracks
144  AliESDtrackCuts* tpcTrackCuts = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
145  tpcTrackCuts->SetPtRange(0.2, 5.0);
146  tpcTrackCuts->SetEtaRange(-0.8, 0.8);
147  tpcTrackCuts->SetMinNClustersTPC(70);
148  trackCuts->AddCuts(tpcTrackCuts);
149  task->SetTrackCuts(trackCuts);
150  } else if (useTracksForRef == 2) { // hybrid tracks
151  // Basic cuts for global tracks - working for 10h!
152  AliESDtrackCuts* baseCuts = new AliESDtrackCuts("base");
153  TFormula *f1NClustersTPCLinearPtDep = new TFormula("f1NClustersTPCLinearPtDep","70.+30./20.*x");
154  baseCuts->SetMinNClustersTPCPtDep(f1NClustersTPCLinearPtDep,20.);
155  baseCuts->SetMinNClustersTPC(70);
156  baseCuts->SetMaxChi2PerClusterTPC(4);
157  baseCuts->SetRequireTPCStandAlone(kTRUE); //cut on NClustersTPC and chi2TPC Iter1
158  baseCuts->SetAcceptKinkDaughters(kFALSE);
159  baseCuts->SetRequireTPCRefit(kTRUE);
160  baseCuts->SetMaxFractionSharedTPCClusters(0.4);
161  // ITS
162  baseCuts->SetRequireITSRefit(kTRUE);
163  //accept secondaries
164  baseCuts->SetMaxDCAToVertexXY(2.4);
165  baseCuts->SetMaxDCAToVertexZ(3.2);
166  baseCuts->SetDCAToVertex2D(kTRUE);
167  //reject fakes
168  baseCuts->SetMaxChi2PerClusterITS(36);
169  baseCuts->SetMaxChi2TPCConstrainedGlobal(36);
170  baseCuts->SetRequireSigmaToVertex(kFALSE);
171  // baseCuts->SetEtaRange(-0.9,0.9);
172  // baseCuts->SetPtRange(0.15, 1E+15.);
173  baseCuts->SetEtaRange(-0.8, 0.8);
174  baseCuts->SetPtRange(0.2, 5.);
175 
176  // Good global tracks
177  AliESDtrackCuts* globalCuts = baseCuts->Clone("global");
178  globalCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD, AliESDtrackCuts::kAny);
179 
180  // tracks with vertex constraint
181  AliESDtrackCuts* constrainedCuts = baseCuts->Clone("vertexConstrained");
182  constrainedCuts->SetRequireITSRefit(kFALSE);
183 
184  // Add
185  trackCuts->AddCuts(globalCuts);
186  trackCuts->AddCuts(constrainedCuts);
187  task->SetTrackCuts(trackCuts);
188  } // end of hybrid
189  } // end of Pb-Pb
190  if (ispA) {
191  if (useTracksForRef == 1) { // tpc tracks
192  AliESDtrackCuts* tpcTrackCuts = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
193  tpcTrackCuts->SetPtRange(0.2, 5.0);
194  tpcTrackCuts->SetEtaRange(-0.8, 0.8);
195  tpcTrackCuts->SetMinNClustersTPC(70);
196  trackCuts->AddCuts(tpcTrackCuts);
197  task->SetTrackCuts(trackCuts);
198  } else if (useTracksForRef == 2) { // hybrid tracks
199  AliESDtrackCuts* baseCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2011(kFALSE);
200  baseCuts->SetMaxDCAToVertexXY(2.4);
201  baseCuts->SetMaxDCAToVertexZ(3.2);
202  baseCuts->SetDCAToVertex2D(kTRUE);
203  baseCuts->SetMaxChi2TPCConstrainedGlobal(36);
204  baseCuts->SetMaxFractionSharedTPCClusters(0.4);
205  // extra
206  baseCuts->SetEtaRange(-0.8, 0.8);
207  baseCuts->SetPtRange(0.2, 5.);
208 
209  // Global good tracks
210  AliESDtrackCuts* globalCuts = baseCuts->Clone("global");
211 
212  // tracks with vertex constraint
213  AliESDtrackCuts* constrainedCuts = baseCuts->Clone("vertexConstrained");
214  constrainedCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD, AliESDtrackCuts::kOff);
215  constrainedCuts->SetRequireITSRefit(kTRUE);
216 
217  // Add
218  trackCuts->AddCuts(globalCuts);
219  trackCuts->AddCuts(constrainedCuts);
220  task->SetTrackCuts(trackCuts);
221  } // end of hybrid
222  } // end of pA
223  } // end of tracks
224 
225  // --- Create containers for output --------------------------------
226  const char* sumName = Form("FlowQCSums%s%s", fwdDet.Data(), AliForwardFlowTaskQC::GetQCType(flags, false));
227  const char* resName = Form("FlowQCResults%s%s", fwdDet.Data(), AliForwardFlowTaskQC::GetQCType(flags, false));
228  AliAnalysisDataContainer* sums =
229  mgr->CreateContainer(sumName, TList::Class(),
230  AliAnalysisManager::kOutputContainer,
231  AliAnalysisManager::GetCommonFileName());
232  AliAnalysisDataContainer* output =
233  mgr->CreateContainer(resName, TList::Class(),
234  AliAnalysisManager::kParamContainer,
235  AliAnalysisManager::GetCommonFileName());
236  mgr->ConnectInput(task, 0, mgr->GetCommonInputContainer());
237  mgr->ConnectOutput(task, 1, sums);
238  mgr->ConnectOutput(task, 2, output);
239 }
240 /*
241  * EOF
242  */
double Double_t
Definition: External.C:58
int Int_t
Definition: External.C:63
void AddTaskForwardFlowQC(Int_t maxMom=5, TString fwdDet="FMD", Bool_t useEtaGap=kFALSE, Bool_t use3cor=kFALSE, Bool_t mc=kFALSE, Double_t outlierCutFMD=4.0, Double_t outlierCutSPD=4.0, Double_t etaGap=2.0, UInt_t useTracksForRef=0, Bool_t useCent=kFALSE, Bool_t ispA=kFALSE, Bool_t useMCVtx=kFALSE, Bool_t satVtx=kFALSE, Bool_t addFlow=kFALSE)
unsigned int UInt_t
Definition: External.C:33
void SetUseFlowWeights(Bool_t use=kTRUE)
static UShort_t CheckForAOD()
unsigned short UShort_t
Definition: External.C:28
void SetUseMCHeaderVertex(Bool_t use=kTRUE)
static const Char_t * GetQCType(UShort_t flags, Bool_t prependUS=kTRUE)
bool Bool_t
Definition: External.C:53
static Bool_t CheckForTask(const char *clsOrName, Bool_t cls=true)
void SetUseImpactParameter(Bool_t use=kTRUE)