]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGHF/hfe/macros/configs/PbPb/ConfigHFEnpePbPb.C
Fix for end-of-line style
[u/mrichter/AliRoot.git] / PWGHF / hfe / macros / configs / PbPb / ConfigHFEnpePbPb.C
CommitLineData
1452614c 1TF1* GetEtaCorrection(){
2 TString list=gSystem->Getenv("LIST");
3
4 TString etaMap="$ALICE_ROOT/PWGHF/hfe/macros/configs/PbPb/EtaCorrMapsTPC.root";
5 if (gSystem->AccessPathName(gSystem->ExpandPathName(etaMap.Data()))){
6 Error("ConfigPbPb2010_Cent","Eta map not found: %s",etaMap.Data());
7 return 0;
8 }
9
10 TFile f(etaMap.Data());
11 if (!f.IsOpen()) return 0;
12 gROOT->cd();
13 TList *keys=f.GetListOfKeys();
14
15 for (Int_t i=0; i<keys->GetEntries(); ++i){
16 TString kName=keys->At(i)->GetName();
17 TPRegexp reg(kName);
18 if (reg.MatchB(list)){
19 printf("Using Eta Correction Function: %s\n",kName.Data());
20 return (TF1*)f.Get(kName.Data());
21 }
22 }
23 return 0;
24}
25
26Bool_t ReadContaminationFunctions(TString filename, TF1 **functions, double sigma){
27 TFile *in = TFile::Open(Form("$ALICE_ROOT/PWGHF/hfe/macros/configs/PbPb/%s", filename.Data()));
28 gROOT->cd();
29 int isig = static_cast<int>(sigma * 100.);
30 printf("Getting hadron background for the sigma cut: %d\n", isig);
31 bool status = kTRUE;
32 for(int icent = 0; icent < 12; icent++){
33 functions[icent] = dynamic_cast<TF1 *>(in->Get(Form("hback_%d_%d", isig, icent)));
34 if(functions[icent]) printf("Config for centrality class %d found\n", icent);
35 else{
36 printf("Config for the centrality class %d not found\n", icent);
37 status = kFALSE;
38 }
39 }
40 delete in;
41 return status;
42}
43
44AliAnalysisTaskHFE* ConfigHFEnpePbPb(Bool_t useMC, Bool_t isAOD, TString appendix,
45 UChar_t TPCcl=70, UChar_t TPCclPID = 80,
46 UChar_t ITScl=3, Double_t DCAxy=1000., Double_t DCAz=1000.,
47 Double_t* tpcdEdxcutlow=NULL, Double_t* tpcdEdxcuthigh=NULL,
48 Double_t TOFs=3., Int_t TOFmis=0,
49 Int_t itshitpixel = 0, Double_t itsChi2PerClusters, Double_t tpcClShared,
50 Bool_t etacor = kFALSE, Bool_t multicor = kFALSE,
51 Double_t etami=-0.8, Double_t etama=0.8,
52 Double_t assETAm=-0.8, Double_t assETAp=0.8,
53 Int_t assITS=2,
54 Int_t assTPCcl=100, Int_t assTPCPIDcl=80,
55 Double_t assDCAr=1.0, Double_t assDCAz=2.0,
56 Double_t *assTPCSminus=NULL, Double_t *assTPCSplus=NULL,
2b66effc 57 Bool_t useCat1Tracks = kTRUE, Bool_t useCat2Tracks = kTRUE,
58 Bool_t rejectMCFake = kFALSE)
1452614c 59{
60 Bool_t kAnalyseTaggedTracks = kFALSE;
61 Bool_t kApplyPreselection = kTRUE;
62
63 //***************************************//
64 // Setting up the HFE cuts //
65 //***************************************//
66
67 AliHFEcuts *hfecuts = new AliHFEcuts(appendix,"HFE cuts for PbPb");
68 //hfecuts->SetQAOn();
69 hfecuts->CreateStandardCuts();
70 hfecuts->SetMinNClustersTPC(TPCcl);
71 hfecuts->SetMinNClustersTPCPID(TPCclPID);
72 hfecuts->SetMinNClustersITS(ITScl);
73 hfecuts->SetMinRatioTPCclusters(0.6);
74 hfecuts->SetTPCmodes(AliHFEextraCuts::kFoundAll, AliHFEextraCuts::kFoundAllOverFindable);
75 hfecuts->SetCutITSpixel(itshitpixel);
76 hfecuts->SetCheckITSLayerStatus(kFALSE);
77 hfecuts->SetMaxChi2perClusterITS(itsChi2PerClusters);
78 hfecuts->SetEtaRange(etami,etama);
79 hfecuts->SetFractionOfSharedTPCClusters(tpcClShared);
80 hfecuts->SetAcceptKinkMothers();
81 if(isAOD) hfecuts->SetAODFilterBit(2);
82
77db0996 83 if((itshitpixel==AliHFEextraCuts::kAny) || (itshitpixel==AliHFEextraCuts::kSecond))
1452614c 84 hfecuts->SetProductionVertex(0,7,0,7);
85
86 hfecuts->SetMaxImpactParam(DCAxy,DCAz);
87 hfecuts->SetUseMixedVertex(kTRUE);
88 hfecuts->SetVertexRange(10.);
89
90 // TOF settings:
91 Int_t usetof=0;
92 Bool_t kTOFmis=kFALSE;
93 if (TOFs>0.){
94 usetof = 1;
95 printf("CONFIGURATION FILE: TOF is used \n");
96 hfecuts->SetTOFPIDStep(kTRUE);
46961971
R
97 if(useMC) hfecuts->SetMatchTOFLabel(kTRUE);
98 //if(useMC) hfecuts->SetMatchTOFLabel(kFALSE);
1452614c 99 printf("CONFIGURATION FILE: TOF PID step is requested !!!! \n");
100 if (TOFmis>0){
101 kTOFmis = kTRUE;
102 printf("CONFIGURATION FILE: TOF mismatch rejection is set ON \n");
103 }
104 }
105
106 //***************************************//
107 // Setting up the task //
108 //***************************************//
109
110 AliAnalysisTaskHFE *task = new AliAnalysisTaskHFE(Form("HFEtask%s",appendix.Data()));
111 printf("task %p\n", task);
112 task->SetPbPbAnalysis();
113 task->SetRemovePileUp(kFALSE);
114 task->SetHFECuts(hfecuts);
115 task->SetRejectKinkMother(kFALSE);
116 task->GetPIDQAManager()->SetHighResolutionHistos();
2b66effc 117 if(useMC && rejectMCFake) task->SetRejectMCFakeTracks(kTRUE); // MC label negative
1452614c 118
119 // Determine the centrality estimator
120 task->SetCentralityEstimator("V0M");
121
122 //***************************************//
123 // Prepare preselection //
124 // This mimics the ESD->AOD filter in //
125 // case of the ESD analysis and selects //
126 // only tracks which will be selected in //
127 // the AOD analysis with the given filter//
128 // bit. Not to be applied for AODS. //
129 // For pPb the cuts used are (bit 4) //
130 // esdTrackCutsHG0 from file $ALICE_ROOT///
131 // ANALYSIS/macros/AddTaskESDFilter.C //
132 //***************************************//
133 if(kApplyPreselection){
134 AliESDtrackCuts* esdTrackCutsH = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
135 esdTrackCutsH->SetClusterRequirementITS(AliESDtrackCuts::kSPD, AliESDtrackCuts::kAny);
136 task->SetHFECutsPreselect(esdTrackCutsH);
137 printf("Put a preselection cut\n");
138 task->SetFillNoCuts(kTRUE);
139 }
140
141 //***************************************//
142 // Variable manager //
143 //***************************************//
144 // Define Variables
145 Double_t ptbinning[36] = {0., 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1., 1.1, 1.2, 1.3, 1.4, 1.5, 1.75, 2., 2.25, 2.5, 2.75, 3., 3.5, 4., 4.5, 5., 5.5, 6., 7., 8., 10., 12., 14., 16., 18., 20.};
146 Double_t etabinning[17] = {-0.8, -0.7, -0.6, -0.5, -0.4, -0.3, -0.2, -0.1, 0., 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8};
147
148 Int_t sizept=(sizeof(ptbinning)/sizeof(double))-1;
149 Int_t sizeeta=(sizeof(etabinning)/sizeof(double))-1;
150
151 AliHFEvarManager *vm = task->GetVarManager();
152 vm->AddVariable("pt", sizept, ptbinning);
153 vm->AddVariable("eta", sizeeta, -0.8,0.8);
154 vm->AddVariable("phi",21, -0, 2*TMath::Pi());
155 vm->AddVariable("charge");
156 vm->AddVariable("source");
157 vm->AddVariable("centrality");
158
159 // For the moment, remove the part dedicated to the background subtraction.
160 // It will be implemented in a different way, reading it from a root file.
161
162 //***************************************//
163 // Configure the PID //
164 //***************************************//
165
166 // Define PID
167 AliHFEpid *pid = task->GetPID();
168 if(useMC) pid->SetHasMCData(kTRUE);
169
170 if (usetof){
171 pid->AddDetector("TOF", 0);
172 pid->AddDetector("TPC", 1);
173 } else {
174 pid->AddDetector("TPC", 0);
175 }
176
177 // Configure TPC PID
178 // do the identical thing in data and MC
179 Double_t paramsTPCdEdxcutlow[12] ={0.0, 0.0, 0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0};
180 if(tpcdEdxcutlow) memcpy(paramsTPCdEdxcutlow,tpcdEdxcutlow,sizeof(paramsTPCdEdxcutlow));
181
182 Double_t paramsTPCdEdxcuthigh[12] ={3.0, 3.0, 3.0,3.0,3.0,3.0,3.0,3.0,3.0,3.0,3.0,3.0};
183 if(tpcdEdxcuthigh) memcpy(paramsTPCdEdxcuthigh,tpcdEdxcuthigh,sizeof(paramsTPCdEdxcuthigh));
184
185 char *cutmodel;
186 cutmodel="pol0";
187
188 for(Int_t a=0;a<11;a++){
189 // Not necessary anymore, since the PbPb case is handled similarly to the pp case
190 // cout << a << " " << paramsTPCdEdxcut[a] << endl;
191 Double_t tpcparamlow[1]={paramsTPCdEdxcutlow[a]};
192 Float_t tpcparamhigh=paramsTPCdEdxcuthigh[a];
193 pid->ConfigureTPCcentralityCut(a,cutmodel,tpcparamlow,tpcparamhigh);
194 }
195
196 if(!useMC){
197 AliHFEpidTPC *tpcpid = pid->GetDetPID(AliHFEpid::kTOFpid);
198 if(etacor){
199 // Apply eta correction
200 TF1 *etacorrection = GetEtaCorrection();
201 if(etacorrection) tpcpid->SetEtaCorrection(etacorrection);
202 }
203 if(multicor){
204 TF1 *centralityCorrection = new TF1("centralityCorrection", "pol1", 0., 10000.);
205 centralityCorrection->SetParameter(0, 1.0);
206 centralityCorrection->SetParameter(1, -0.00002);
207 tpcpid->SetCentralityCorrection(centralityCorrection);
208 }
209 }
210
211 // Configure TOF PID
212 if (usetof){
213 pid->ConfigureTOF(TOFs);
214 AliHFEpidTOF *tofpid = pid->GetDetPID(AliHFEpid::kTOFpid);
215 if (kTOFmis){
216 tofpid->SetRejectTOFmismatch();
217 }
218 }
219
220 // To make different upper TOF cut to see contamination effect
221 // The below two lines should be removed after this check
222 //AliHFEpidTOF *tofpid = pid->GetDetPID(AliHFEpid::kTOFpid);
223 //if(TOFs<3.) tofpid->SetTOFnSigmaBand(-3,TOFs); //only to check the assymmetric tof cut
224
225 // Load hadron background
226 if(!useMC){
227 Bool_t status = kTRUE;
228 TF1 *hBackground[12];
229 status = ReadContaminationFunctions("hadronContamination_PbPbTPC.root", hBackground, tpcdEdxcutlow[0]);
230 for(Int_t a=0;a<12;a++) {
231 //printf("back %f \n",hBackground[a]);
232 if(status) task->SetBackGroundFactorsFunction(hBackground[a],a);
233 else printf("not all background functions found\n");
234 }
235 }
236
237 //***************************************//
238 // Configure NPE plugin //
239 //***************************************//
240
241 AliHFENonPhotonicElectron *backe = new AliHFENonPhotonicElectron(Form("HFEBackGroundSubtractionPID2%s",appendix.Data()),"Background subtraction"); //appendix
242 //Setting the Cuts for the Associated electron-pool
243 AliHFEcuts *hfeBackgroundCuts = new AliHFEcuts(Form("HFEBackSub%s",appendix.Data()),"Background sub Cuts");
244 // hfeBackgroundCuts->SetEtaRange(assETA);
245 hfeBackgroundCuts->SetEtaRange(assETAm,assETAp);
246 hfeBackgroundCuts->SetPtRange(0.1,1e10);
247
248 hfeBackgroundCuts->SetMaxChi2perClusterTPC(4);
249 hfeBackgroundCuts->SetMinNClustersITS(assITS);
250 hfeBackgroundCuts->SetMinNClustersTPC(assTPCcl);
251 hfeBackgroundCuts->SetMinNClustersTPCPID(assTPCPIDcl);
252 hfeBackgroundCuts->SetMaxImpactParam(assDCAr,assDCAz);
253 if(isAOD) hfeBackgroundCuts->SetAODFilterBit(0);
254 //hfeBackgroundCuts->SetQAOn(); // QA
255
256 AliHFEpid *pidbackground = backe->GetPIDBackground();
257 if(useMC) pidbackground->SetHasMCData(kTRUE);
258 pidbackground->AddDetector("TPC", 0);
259 Double_t paramsTPCdEdxcutlowAssoc[12] ={-3.0,-3.0,-3.0,-3.0,-3.0,-3.0,-3.0,-3.0,-3.0,-3.0,-3.0,-3.0};
260 if(assTPCSminus) memcpy(paramsTPCdEdxcutlowAssoc,assTPCSminus,sizeof(paramsTPCdEdxcutlowAssoc));
261
262 Double_t paramsTPCdEdxcuthighAssoc[12] ={3.0, 3.0, 3.0,3.0,3.0,3.0,3.0,3.0,3.0,3.0,3.0,3.0};
263 if(assTPCSplus) memcpy(paramsTPCdEdxcuthighAssoc,assTPCSplus,sizeof(paramsTPCdEdxcuthighAssoc));
264
265 char *cutmodelAssoc;
266 cutmodelAssoc="pol0";
267 for(Int_t a=0;a<11;a++){
268 // cout << a << " " << paramsTPCdEdxcut[a] << endl;
269 Double_t tpcparamlow[1]={paramsTPCdEdxcutlowAssoc[a]};
270 Float_t tpcparamhigh=paramsTPCdEdxcuthighAssoc[a];
271 pidbackground->ConfigureTPCcentralityCut(a,cutmodelAssoc,tpcparamlow,tpcparamhigh);
272 }
273 //backe->GetPIDBackgroundQAManager()->SetHighResolutionHistos();
274 backe->SetHFEBackgroundCuts(hfeBackgroundCuts);
275
276 // Selection of associated tracks for the pool
277 if(useCat1Tracks) backe->SelectCategory1Tracks(kTRUE);
278 if(useCat2Tracks){
279 backe->SelectCategory2Tracks(kTRUE);
280 backe->SetITSMeanShift(-0.5);
281 }
282
283 // apply opening angle cut to reduce file size
284 backe->SetMaxInvMass(0.3);
285 backe->SetPtBinning(sizept, ptbinning);
286 backe->SetEtaBinning(sizeeta, etabinning);
287
288 task->SetHFEBackgroundSubtraction(backe);
289
290 //***************************************//
291 // V0 tagged tracks //
292 //***************************************//
293
294 if(kAnalyseTaggedTracks){
295 AliHFEcuts *v0trackCuts = new AliHFEcuts("V0trackCuts", "Track Cuts for tagged track Analysis");
296 v0trackCuts->CreateStandardCuts();
297 v0trackCuts->SetMinNClustersTPC(TPCcl);
298 v0trackCuts->SetMinNClustersTPCPID(TPCclPID);
299 v0trackCuts->SetMinRatioTPCclusters(0.6);
300 v0trackCuts->SetTPCmodes(AliHFEextraCuts::kFound, AliHFEextraCuts::kFoundOverFindable);
301 v0trackCuts->SetMinNClustersITS(1);
302 v0trackCuts->SetCutITSpixel(AliHFEextraCuts::kFirst);
303 v0trackCuts->SetCheckITSLayerStatus(kFALSE);
304 v0trackCuts->UnsetVertexRequirement();
305 //hfecuts->SetSigmaToVertex(10);
306 if(usetof) v0trackCuts->SetTOFPIDStep(kTRUE);
307 v0trackCuts->SetQAOn();
308
309 task->SwitchOnPlugin(AliAnalysisTaskHFE::kTaggedTrackAnalysis);
310 task->SetTaggedTrackCuts(v0trackCuts);
311 task->SetCleanTaggedTrack(kTRUE);
312 }
313
314 // QA
315 printf("task %p\n", task);
316 task->SetQAOn(AliAnalysisTaskHFE::kPIDqa);
317 task->SetQAOn(AliAnalysisTaskHFE::kMCqa);
318 task->SwitchOnPlugin(AliAnalysisTaskHFE::kNonPhotonicElectron);
319 task->SwitchOnPlugin(AliAnalysisTaskHFE::kDEstep);
320
321 printf("*************************************\n");
322 printf("Configuring standard Task:\n");
323 task->PrintStatus();
324 pid->PrintStatus();
325 printf("*************************************\n");
326 return task;
327}