]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGLF/STRANGENESS/LambdaK0PbPb/AliV0CutVariations.C
A messy macro to generate preliminary plots with K0s and Lambda spectra.
[u/mrichter/AliRoot.git] / PWGLF / STRANGENESS / LambdaK0PbPb / AliV0CutVariations.C
CommitLineData
dc9b43d9 1#define AliV0CutVariations_cxx
2// The class definition in AliV0CutVariations.h has been generated automatically
3// by the ROOT utility TTree::MakeSelector(). This class is derived
4// from the ROOT class TSelector. For more information on the TSelector
5// framework see $ROOTSYS/README/README.SELECTOR or the ROOT User Manual.
6
7// The following methods are defined in this file:
8// Begin(): called every time a loop on the tree starts,
9// a convenient place to create your histograms.
10// SlaveBegin(): called after Begin(), when on PROOF called only on the
11// slave servers.
12// Process(): called for each event, in this function you decide what
13// to read and fill your histograms.
14// SlaveTerminate: called at the end of the loop on the tree, when on PROOF
15// called only on the slave servers.
16// Terminate(): called at the end of the loop on the tree,
17// a convenient place to draw/fit your histograms.
18//
19// To use this file, try the following session on your Tree T:
20//
21// Root > T->Process("AliV0CutVariations.C")
22// Root > T->Process("AliV0CutVariations.C","some options")
23// Root > T->Process("AliV0CutVariations.C+")
24//
25
26#include "AliV0CutVariations.h"
27
28#include "Riostream.h"
29#include "TCanvas.h"
30#include "TPDGCode.h"
31#include <TDatabasePDG.h>
32#include <TH2.h>
33#include <TStyle.h>
34#include <TH1F.h>
35#include <TH2F.h>
36#include <TH3F.h>
37
38AliV0CutVariations::AliV0CutVariations(TTree * /*tree*/ ) :
39fChain(0),
40fIsMC(kFALSE),
41fCMin(0.),
42fCMax(90.),
43fCPA(0.9975),
44fDCA(1.0),
45fTPCcr(70.),
46fTPCcrfd(0.8),
47fDCApv(0.1),
48fMult(0),
49fdEdx(0),
50fdEdxPid(0),
51fCosPA(0),
52fDtrDCA(0),
53fTPCrows(0),
54fTPCratio(0),
55fPrimDCA(0),
56
57fK0sM(0),
58fK0sSi(0),
59fK0sMC(0),
60fK0sAs(0),
61
62fLambdaM(0),
63fLambdaSi(0),
64fLambdaMC(0),
65fLambdaAs(0),
66
67fLambdaEff(0),
68fLambdaPt(0),
69
70fLambdaBarM(0),
71fLambdaBarSi(0),
72fLambdaBarMC(0),
73fLambdaBarAs(0),
74
75fLambdaBarEff(0),
76fLambdaBarPt(0),
77
78fLambdaFromXi(0),
79fXiM(0),
80fXiSiP(0),
81
82fLambdaBarFromXiBar(0),
83fXiBarM(0),
84fXiBarSiP(0)
85{ }
86
87void AliV0CutVariations::Begin(TTree * /*tree*/)
88{
89 // The Begin() function is called at the start of the query.
90 // When running with PROOF Begin() is only called on the client.
91 // The tree argument is deprecated (on PROOF 0 is passed).
92
93 TString option = GetOption();
94
95}
96
97void AliV0CutVariations::SlaveBegin(TTree * tree)
98{
99 // The SlaveBegin() function is called after the Begin() function.
100 // When running with PROOF SlaveBegin() is called on each slave server.
101 // The tree argument is deprecated (on PROOF 0 is passed).
102
103 Init(tree);
104
105 cout<<tree->GetEntries()<<endl;
106
107 TString option = GetOption();
108
109 Int_t lbins=100; // number of bins in lt
110 Int_t kbins=33; // number of bins in pt of Xi
111 Double_t ltMax=100.;
112 Double_t ptMax=12.;
113
114 const Double_t xBins[]={
115 0.0,0.1,0.2,0.3,0.4,0.5,0.6,0.7,0.8,0.9,
116 1.0,1.1,1.2,1.3,1.4,1.5,1.6,1.7,1.8,1.9,2.0,
117 2.2,2.4,2.6,2.8,3.0,3.2,3.4,3.6,3.8,4.0,
118 4.5,5.0,5.5,6.5,8.0,10.0,12.0
119 };
120 const Int_t nBins=sizeof(xBins)/sizeof(Double_t) - 1;
121
122 Double_t yBins[lbins+1];
123 for (Int_t i=0; i<=(lbins+1); i++) yBins[i]=i*ltMax/lbins;
124 Double_t zBins[kbins+1];
125 for (Int_t i=0; i<=(kbins+1); i++) zBins[i]=i*(ptMax+2)/kbins;
126
127
128 //fMult=new TH1F("fMult","Multiplicity",1100,0.,3300);
129 //fMult->GetXaxis()->SetTitle("N tracks");
130 fMult=new TH1F("fMult","Centrality",100,0.,100);
131 fMult->GetXaxis()->SetTitle("Centrality (%)");
132 fOutput->Add(fMult);
133
134
135 fCosPA=new TH1F("fCosPA","Cos(PA) distribution",50,0.9975,1.0005);
136 fCosPA->GetXaxis()->SetTitle("Cos(PA)");
137 fOutput->Add(fCosPA);
138
139 fDtrDCA=new TH1F("fDtrDCA","DCA between V0 daughters",50,0.0,1.5);
140 fDtrDCA->GetXaxis()->SetTitle("DCA (rel. u.)");
141 fOutput->Add(fDtrDCA);
142
143 fTPCrows=new TH1F("fTPCrows","TPC crossed pad rows",180,0.,180.);
144 fTPCrows->GetXaxis()->SetTitle("TPC crossed pad rows");
145 fOutput->Add(fTPCrows);
146
147 fTPCratio=new TH1F("fTPCratio","TPC crossed/findable pad rows",50,0.0,1.5);
148 fTPCratio->GetXaxis()->SetTitle("TPC crossed/findable pad rows");
149 fOutput->Add(fTPCratio);
150
151 fPrimDCA=new TH1F("fPrimDCA","DCA wrt the primary vertex",50,0.0,1.5);
152 fPrimDCA->GetXaxis()->SetTitle("DCA wrt the PV (cm)");
153 fOutput->Add(fPrimDCA);
154
155
156 fdEdx=new TH2F("fdEdx","dE/dx",50,0.2,3,50,0.,6.);
157 fOutput->Add(fdEdx);
158
159 fdEdxPid=new TH2F("fdEdxPid","dE/dx with PID",50,0.2,3,50,0.,6.);
160 fOutput->Add(fdEdxPid);
161
162 fK0sM =
163 new TH2F("fK0sM", "Mass for K^{0}_{s}", 50,0.448,0.548,nBins,xBins);
164 fK0sM->GetXaxis()->SetTitle("Mass (GeV/c)");
165 fOutput->Add(fK0sM);
166
167 fK0sSi =
168 new TH2F("fK0sSi","L_{T} vs p_{T} for K^{0}_{s}, side-band subtracted",
169 nBins,xBins,lbins,0.,ltMax);
170 fK0sSi->GetXaxis()->SetTitle("p_{T} (GeV/c)");
171 fK0sSi->GetYaxis()->SetTitle("L_{T} (cm)");
172 fOutput->Add(fK0sSi);
173
174 fK0sMC =
175 new TH2F("fK0sMC","L_{T} vs p_{T} for K^{0}_{s}, from MC stack",
176 nBins,xBins,lbins,0.,ltMax);
177 fK0sMC->GetXaxis()->SetTitle("p_{T} (GeV/c)");
178 fK0sMC->GetYaxis()->SetTitle("L_{T} (cm)");
179 fOutput->Add(fK0sMC);
180
181 fK0sAs =
182 new TH2F("fK0sAs", "L_{T} vs p_{T} for K^{0}_{s}, associated",
183 nBins,xBins,lbins,0.,ltMax);
184 fK0sAs->GetXaxis()->SetTitle("p_{T} (GeV/c)");
185 fK0sAs->GetYaxis()->SetTitle("L_{T} (cm)");
186 fOutput->Add(fK0sAs);
187
188 //----------------------
189
190 fLambdaM =
191 new TH2F("fLambdaM","Mass for \\Lambda", 100, 1.065, 1.165,nBins,xBins);
192 fLambdaM->GetXaxis()->SetTitle("Mass (GeV/c)");
193 fOutput->Add(fLambdaM);
194
195 fLambdaSi =
196 new TH2F("fLambdaSi","L_{T} vs p_{T} for \\Lambda, side-band subtructed",
197 nBins,xBins,lbins,0.,ltMax);
198 fLambdaSi->GetXaxis()->SetTitle("p_{T} (GeV/c)");
199 fLambdaSi->GetYaxis()->SetTitle("L_{T} (cm)");
200 fOutput->Add(fLambdaSi);
201
202 fLambdaMC =
203 new TH2F("fLambdaMC","L_{T} vs p_{T} for \\Lambda, from MC stack",
204 nBins,xBins,lbins,0.,ltMax);
205 fLambdaMC->GetXaxis()->SetTitle("p_{T} (GeV/c)");
206 fLambdaMC->GetYaxis()->SetTitle("L_{T} (cm)");
207 fOutput->Add(fLambdaMC);
208
209 fLambdaAs =
210 new TH2F("fLambdaAs","L_{T} vs p_{T} for \\Lambda, associated",
211 nBins,xBins,lbins,0.,ltMax);
212 fLambdaAs->GetXaxis()->SetTitle("p_{T} (GeV/c)");
213 fLambdaAs->GetYaxis()->SetTitle("L_{T} (cm)");
214 fOutput->Add(fLambdaAs);
215
216 //----------------------
217
218 fLambdaEff=fLambdaAs->ProjectionX();
219 fLambdaEff->SetName("fLambdaEff");
220 fLambdaEff->SetTitle("Efficiency for #Lambda");
221 fOutput->Add(fLambdaEff);
222
223 fLambdaPt=fLambdaAs->ProjectionX();
224 fLambdaPt->SetName("fLambdaPt");
225 fLambdaPt->SetTitle("Raw #Lambda pT spectrum");
226 fOutput->Add(fLambdaPt);
227 //----------------------
228
229 fLambdaBarM =
230 new TH2F("fLambdaBarM","Mass for anti-\\Lambda",100,1.065,1.165,nBins,xBins);
231 fLambdaBarM->GetXaxis()->SetTitle("Mass (GeV/c)");
232 fOutput->Add(fLambdaBarM);
233
234 fLambdaBarSi =
235 new TH2F("fLambdaBarSi","L_{T} vs p_{T} for anti-\\Lambda, side-band subtructed",
236 nBins,xBins,lbins,0.,ltMax);
237 fLambdaBarSi->GetXaxis()->SetTitle("p_{T} (GeV/c)");
238 fLambdaBarSi->GetYaxis()->SetTitle("L_{T} (cm)");
239 fOutput->Add(fLambdaBarSi);
240
241 fLambdaBarMC =
242 new TH2F("fLambdaBarMC","L_{T} vs p_{T} for anti-\\Lambda, from MC stack",
243 nBins,xBins,lbins,0.,ltMax);
244 fLambdaBarMC->GetXaxis()->SetTitle("p_{T} (GeV/c)");
245 fLambdaBarMC->GetYaxis()->SetTitle("L_{T} (cm)");
246 fOutput->Add(fLambdaBarMC);
247
248 fLambdaBarAs =
249 new TH2F("fLambdaBarAs","L_{T} vs p_{T} for anti-\\Lambda, associated",
250 nBins,xBins,lbins,0.,ltMax);
251 fLambdaBarAs->GetXaxis()->SetTitle("p_{T} (GeV/c)");
252 fLambdaBarAs->GetYaxis()->SetTitle("L_{T} (cm)");
253 fOutput->Add(fLambdaBarAs);
254
255
256 //----------------------
257
258 fLambdaBarEff=fLambdaBarAs->ProjectionX();
259 fLambdaBarEff->SetName("fLambdaBarEff");
260 fLambdaBarEff->SetTitle("Efficiency for anti-#Lambda");
261 fOutput->Add(fLambdaBarEff);
262
263 fLambdaBarPt=fLambdaBarAs->ProjectionX();
264 fLambdaBarPt->SetName("fLambdaBarPt");
265 fLambdaBarPt->SetTitle("Raw anti-#Lambda pT spectrum");
266 fOutput->Add(fLambdaBarPt);
267
268 //----------------------
269 fLambdaFromXi=new TH3F("fLambdaFromXi","L_{T} vs p_{T} vs p_{T} of \\Xi for \\Lambda from \\Xi",
270 nBins,xBins,lbins,yBins,kbins,zBins);
271 fOutput->Add(fLambdaFromXi);
272
273 fXiM =
274 new TH2F("fXiM", "\\Xi mass distribution", 50, 1.271, 1.371,33,0.,ptMax+2);
275 fOutput->Add(fXiM);
276
277 fXiSiP = new TH1F("fXiSiP", "Pt for \\Xi, side-band subracted",
278 33,0.,ptMax+2);
279 fOutput->Add(fXiSiP);
280
281
282 fLambdaBarFromXiBar=new TH3F("fLambdaBarFromXiBar","L_{T} vs p_{T} vs p_{T} of anti-\\Xi for anti-\\Lambda from anti-\\Xi",
283 nBins,xBins,lbins,yBins,kbins,zBins);
284 fOutput->Add(fLambdaBarFromXiBar);
285
286 fXiBarM =
287 new TH2F("fXiBarM", "anti-\\Xi mass distribution", 50, 1.271, 1.371,33,0.,ptMax+2);
288 fOutput->Add(fXiBarM);
289
290 fXiBarSiP = new TH1F("fXiBarSiP", "Pt for anti-\\Xi, side-band subracted",
291 33,0.,ptMax+2);
292 fOutput->Add(fXiBarSiP);
293
294}
295
296Bool_t AliV0CutVariations::AcceptTracks() {
297 //if (!t->IsOn(AliAODTrack::kTPCrefit)) return kFALSE;
298 //if (t->GetKinkIndex(0)>0) return kFALSE;
299
300 if (fTreeVariableLeastNbrCrossedRows < fTPCcr) return kFALSE;
301 if (fTreeVariableLeastRatioCrossedRowsOverFindable < fTPCcrfd) return kFALSE;
302
303 if (TMath::Abs(fTreeVariableNegEta) > 0.8) return kFALSE;
304 if (TMath::Abs(fTreeVariablePosEta) > 0.8) return kFALSE;
305
306 fTPCrows->Fill(fTreeVariableLeastNbrCrossedRows);
307 fTPCratio->Fill(fTreeVariableLeastRatioCrossedRowsOverFindable);
308
309 return kTRUE;
310}
311
312Bool_t AliV0CutVariations::AcceptV0()
313{
314 //if (v0->GetOnFlyStatus()) return kFALSE;
315
316 if (fTreeVariableV0CosineOfPointingAngle < fCPA) return kFALSE;
317 if (fTreeVariableDcaV0Daughters > fDCA) return kFALSE;
318
319 if (!AcceptTracks()) return kFALSE;
320
321 if (TMath::Abs(fTreeVariableDcaNegToPrimVertex) < fDCApv) return kFALSE;
322 if (TMath::Abs(fTreeVariableDcaPosToPrimVertex) < fDCApv) return kFALSE;
323
324 if (fTreeVariableV0Radius < 0.9) return kFALSE;
325 if (fTreeVariableV0Radius > 100) return kFALSE;
326
327 fCosPA->Fill(fTreeVariableV0CosineOfPointingAngle);
328 fDtrDCA->Fill(fTreeVariableDcaV0Daughters);
329 fPrimDCA->Fill(fTreeVariableDcaNegToPrimVertex);
330 fPrimDCA->Fill(fTreeVariableDcaPosToPrimVertex);
331
332 return kTRUE;
333}
334
335Bool_t AliV0CutVariations::AcceptPID(Int_t code)
336{
337 //const AliAODPid *pid=ptrack->GetDetPid();
338 //if (!pid) return kTRUE;
339 if (code > 0) {
340 if (fTreeVariablePosTransvMomentum > 1.) return kTRUE;
341 } else {
342 if (fTreeVariableNegTransvMomentum > 1.) return kTRUE;
343 }
344
345
346 if (fIsMC) {
347 // MC PID
348 if (code > 0) {
349 if (fTreeVariablePIDPositive == code) return kTRUE;
350 } else {
351 if (fTreeVariablePIDNegative == code) return kTRUE;
352 }
353 } else {
354 // Real PID
355 if (code > 0) {
356 if (TMath::Abs(fTreeVariableNSigmasPosProton) < 3.) return kTRUE;
357 } else {
358 if (TMath::Abs(fTreeVariableNSigmasNegProton) < 3.) return kTRUE;
359 }
360 }
361
362 return kFALSE;
363}
364
365
366Bool_t AliV0CutVariations::Process(Long64_t entry)
367{
368 // The Process() function is called for each entry in the tree (or possibly
369 // keyed object in the case of PROOF) to be processed. The entry argument
370 // specifies which entry in the currently loaded tree is to be processed.
371 // It can be passed to either AliV0CutVariations::GetEntry() or TBranch::GetEntry()
372 // to read either all or the required parts of the data. When processing
373 // keyed objects with PROOF, the object is already loaded and is available
374 // via the fObject pointer.
375 //
376 // This function should contain the "body" of the analysis. It can contain
377 // simple or elaborate selection criteria, run algorithms on the data
378 // of the event and typically fill histograms.
379 //
380 // The processing can be stopped by calling Abort().
381 //
382 // Use fStatus to set the return value of TTree::Process().
383 //
384 // The return value is currently not used.
385
386 const Double_t yMax=0.5;
387 const Double_t pMin=0.0;
388 //const Double_t lMax=0.001;
389
390 fChain->GetTree()->GetEntry(entry);
391
392 cout<<entry<<'\r';
393
394 if (fTreeVariableMultiplicity<fCMin) return kFALSE;
395 if (fTreeVariableMultiplicity>fCMax) return kFALSE;
396
397 fMult->Fill(fTreeVariableMultiplicity);
398
399 Double_t pt=0;
400 Double_t lt=0;
401 //+++++++ MC
402 if (fIsMC) {
403 Int_t code=fTreeVariablePID;
404
405 if (code != kK0Short)
406 if (code != kLambda0)
407 if (code != kLambda0Bar) goto real;
408
409 pt=fTreeVariablePtMC;
410 if (pt < pMin) goto real;
411
412 if (TMath::Abs(fTreeVariableRapMC) > yMax) goto real;
413
414 //if (TMath::Abs(fTreeVariableV0CreationRadius) > lMax) goto real;
415 if (fTreeVariablePrimaryStatus != 1) goto real;
416
417 lt=fTreeVariableV0Radius;
418 switch (code) {
419 case kK0Short:
420 fK0sMC->Fill(pt,lt);
421 break;
422 case kLambda0:
423 fLambdaMC->Fill(pt,lt);
424 break;
425 case kLambda0Bar:
426 fLambdaBarMC->Fill(pt,lt);
427 break;
428 default: break;
429 }
430
431 }
432 //+++++++
433
434 real:
435
436 pt=fTreeVariablePt;
437 if (pt < pMin) return kFALSE;
438
439 if (!AcceptV0()) return kFALSE;
440
441 lt=fTreeVariableV0Radius;
442 if (lt/pt > 3*7.89/1.1157) return kFALSE;
443
444 //--- V0 switches
445 Bool_t isK0s=kTRUE;
446 Bool_t isLambda=kTRUE;
447 Bool_t isLambdaBar=kTRUE;
448
449 if (0.4977*lt/pt > 3*2.68) isK0s=kFALSE;
450 if (1.1157*lt/pt > 3*7.89) isLambdaBar=isLambda=kFALSE;
451
452 if (!AcceptPID(kProton)) isLambda=kFALSE;
453 if (!AcceptPID(kProtonBar)) isLambdaBar=kFALSE;
454
455 Double_t yK0s=fTreeVariableRapK0Short;
456 Double_t yLam=fTreeVariableRapLambda;
457 if (yK0s > yMax) isK0s=kFALSE;
458 if (yLam > yMax) isLambda=isLambdaBar=kFALSE;
459 //---
460
461 Double_t mass=0., m=0., s=0.;
462 if (isK0s) {
463 mass=fTreeVariableInvMassK0s;
464 fK0sM->Fill(mass,pt);
465
466 m=TDatabasePDG::Instance()->GetParticle(kK0Short)->Mass();
467 s=0.0044 + (0.008-0.0044)/(10-1)*(pt - 1.);
468 if (TMath::Abs(m-mass) < 3*s) {
469 fK0sSi->Fill(pt,lt);
470 } else {
471 isK0s=kFALSE;
472 }
473 if (TMath::Abs(m-mass + 4.5*s) < 1.5*s) {
474 fK0sSi->Fill(pt,lt,-1);
475 }
476 if (TMath::Abs(m-mass - 4.5*s) < 1.5*s) {
477 fK0sSi->Fill(pt,lt,-1);
478 }
479 }
480
481 if (isLambda) {
482 mass=fTreeVariableInvMassLambda;
483 fLambdaM->Fill(mass,pt);
484
485 m=TDatabasePDG::Instance()->GetParticle(kLambda0)->Mass();
486 //s=0.0027 + (0.004-0.0027)/(10-1)*(pt-1);
487 //s=0.0015 + (0.002-0.0015)/(2.6-1)*(pt-1);
488 s=0.0023 + (0.004-0.0023)/(6-1)*(pt-1);
489 if (TMath::Abs(m-mass) < 3*s) {
490 fLambdaSi->Fill(pt,lt);
491 } else {
492 isLambda=kFALSE;
493 }
494 if (TMath::Abs(m-mass + 4.5*s) < 1.5*s) {
495 fLambdaSi->Fill(pt,lt,-1);
496 }
497 if (TMath::Abs(m-mass - 4.5*s) < 1.5*s) {
498 fLambdaSi->Fill(pt,lt,-1);
499 }
500 }
501
502 if (isLambdaBar) {
503 mass=fTreeVariableInvMassAntiLambda;
504 fLambdaBarM->Fill(mass,pt);
505
506 m=TDatabasePDG::Instance()->GetParticle(kLambda0Bar)->Mass();
507 //s=0.0027 + (0.004-0.0027)/(10-1)*(pt-1);
508 //s=0.0015 + (0.002-0.0015)/(2.6-1)*(pt-1);
509 s=0.0023 + (0.004-0.0023)/(6-1)*(pt-1);
510 if (TMath::Abs(m-mass) < 3*s) {
511 fLambdaBarSi->Fill(pt,lt);
512 } else {
513 isLambdaBar=kFALSE;
514 }
515 if (TMath::Abs(m-mass + 4.5*s) < 1.5*s) {
516 fLambdaBarSi->Fill(pt,lt,-1);
517 }
518 if (TMath::Abs(m-mass - 4.5*s) < 1.5*s) {
519 fLambdaBarSi->Fill(pt,lt,-1);
520 }
521 }
522
523 if (!fIsMC) return kFALSE;
524
525 //++++++ MC
526 if (!isK0s)
527 if (!isLambda)
528 if (!isLambdaBar) return kFALSE;//check MC only for the accepted V0s
529
530 Int_t code=fTreeVariablePID;
531 if (code != kK0Short)
532 if (code != kLambda0)
533 if (code != kLambda0Bar) return kFALSE;
534
535 Double_t ptAs=fTreeVariablePtMC;
536 if (ptAs < pMin) return kFALSE;
537
538 Double_t yAs=fTreeVariableRapMC;
539 if (TMath::Abs(yAs) > yMax) return kFALSE;
540
541 Double_t ltAs=fTreeVariableV0Radius;
542 if (fTreeVariablePrimaryStatus == 1) {
543 switch (code) {
544 case kK0Short:
545 if (isK0s) fK0sAs->Fill(ptAs,ltAs);
546 break;
547 case kLambda0:
548 if (isLambda) fLambdaAs->Fill(ptAs,ltAs);
549 break;
550 case kLambda0Bar:
551 if (isLambdaBar) fLambdaBarAs->Fill(ptAs,ltAs);
552 break;
553 default: break;
554 }
555 } else {
556 if (code==kK0Short) return kTRUE;
557
558 Int_t xcode=fTreeVariablePIDMother;
559 Double_t xpt=fTreeVariablePtXiMother;
560
561 switch (code) {
562 case kLambda0:
563 if (isLambda)
564 if ((xcode==kXiMinus) || (xcode==3322))
565 fLambdaFromXi->Fill(ptAs,ltAs,xpt);
566 break;
567 case kLambda0Bar:
568 if (isLambdaBar)
569 if ((xcode==kXiPlusBar)||(xcode==-3322))
570 fLambdaBarFromXiBar->Fill(ptAs,ltAs,xpt);
571 break;
572 default: break;
573 }
574 }
575
576 return kTRUE;
577}
578
579void AliV0CutVariations::SlaveTerminate()
580{
581 // The SlaveTerminate() function is called after all entries or objects
582 // have been processed. When running with PROOF SlaveTerminate() is called
583 // on each slave server.
584
585}
586
587void AliV0CutVariations::Terminate()
588{
589 // The Terminate() function is the last function to be called during
590 // a query. It always runs on the client, it can be used to present
591 // the results graphically or save the results to file.
592
593 cout<<endl;
594
595 if (!fOutput) {
596 Printf("ERROR: fOutput not available");
597 return;
598 }
599
600 fMult = dynamic_cast<TH1F*>(fOutput->FindObject("fMult")) ;
601 if (!fMult) {
602 Printf("ERROR: fMult not available");
603 return;
604 }
605
606 fdEdx = dynamic_cast<TH2F*>(fOutput->FindObject("fdEdx")) ;
607 if (!fdEdx) {
608 Printf("ERROR: fdEdx not available");
609 return;
610 }
611
612 fdEdxPid = dynamic_cast<TH2F*>(fOutput->FindObject("fdEdxPid")) ;
613 if (!fdEdxPid) {
614 Printf("ERROR: fdEdxPid not available");
615 return;
616 }
617
618
619 fK0sMC = dynamic_cast<TH2F*>(fOutput->FindObject("fK0sMC")) ;
620 if (!fK0sMC) {
621 Printf("ERROR: fK0sMC not available");
622 return;
623 }
624 TH1D *k0sMcPx=fK0sMC->ProjectionX(); k0sMcPx->Sumw2();
625 fK0sAs = dynamic_cast<TH2F*>(fOutput->FindObject("fK0sAs")) ;
626 if (!fK0sAs) {
627 Printf("ERROR: fK0sAs not available");
628 return;
629 }
630 TH1D *k0sAsPx=fK0sAs->ProjectionX();
631 k0sAsPx->Sumw2(); //k0sAsPx->Scale(0.69);
632
633
634
635 // Lambda histograms
636 fLambdaFromXi = dynamic_cast<TH3F*>(fOutput->FindObject("fLambdaFromXi")) ;
637 if (!fLambdaFromXi) {
638 Printf("ERROR: fLambdaFromXi not available");
639 return;
640 }
641 TH1D *lambdaFromXiPx=fLambdaFromXi->ProjectionX(); lambdaFromXiPx->Sumw2();
642
643 fLambdaMC = dynamic_cast<TH2F*>(fOutput->FindObject("fLambdaMC")) ;
644 if (!fLambdaMC) {
645 Printf("ERROR: fLambdaMC not available");
646 return;
647 }
648 TH1D *lambdaMcPx=fLambdaMC->ProjectionX(); lambdaMcPx->Sumw2();
649
650 fLambdaAs = dynamic_cast<TH2F*>(fOutput->FindObject("fLambdaAs")) ;
651 if (!fLambdaAs) {
652 Printf("ERROR: fLambdaAs not available");
653 return;
654 }
655 TH1D *lambdaAsPx=fLambdaAs->ProjectionX();
656 lambdaAsPx->Sumw2(); //lambdaAsPx->Scale(0.64);
657
658 fLambdaSi = dynamic_cast<TH2F*>(fOutput->FindObject("fLambdaSi")) ;
659 if (!fLambdaSi) {
660 Printf("ERROR: fLambdaSi not available");
661 return;
662 }
663 TH1D *lambdaSiPx=fLambdaSi->ProjectionX();
664 lambdaSiPx->SetName("fLambdaPt");
665 lambdaSiPx->Sumw2();
666
667 fLambdaEff = dynamic_cast<TH1D*>(fOutput->FindObject("fLambdaEff")) ;
668 if (!fLambdaEff) {
669 Printf("ERROR: fLambdaEff not available");
670 return;
671 }
672 fLambdaPt = dynamic_cast<TH1D*>(fOutput->FindObject("fLambdaPt")) ;
673 if (!fLambdaPt) {
674 Printf("ERROR: fLambdaPt not available");
675 return;
676 }
677
678
679 // anti-Lambda histograms
680 fLambdaBarFromXiBar =
681 dynamic_cast<TH3F*>(fOutput->FindObject("fLambdaBarFromXiBar")) ;
682 if (!fLambdaBarFromXiBar) {
683 Printf("ERROR: fLambdaBarFromXiBar not available");
684 return;
685 }
686 TH1D *lambdaBarFromXiBarPx=
687 fLambdaBarFromXiBar->ProjectionX("Bar"); lambdaBarFromXiBarPx->Sumw2();
688
689 fLambdaBarMC = dynamic_cast<TH2F*>(fOutput->FindObject("fLambdaBarMC")) ;
690 if (!fLambdaBarMC) {
691 Printf("ERROR: fLambdaBarMC not available");
692 return;
693 }
694 TH1D *lambdaBarMcPx=fLambdaBarMC->ProjectionX(); lambdaBarMcPx->Sumw2();
695
696 fLambdaBarAs = dynamic_cast<TH2F*>(fOutput->FindObject("fLambdaBarAs")) ;
697 if (!fLambdaBarAs) {
698 Printf("ERROR: fLambdaBarAs not available");
699 return;
700 }
701 TH1D *lambdaBarAsPx=fLambdaBarAs->ProjectionX();
702 lambdaBarAsPx->Sumw2(); //lambdaBarAsPx->Scale(0.64);
703
704 fLambdaBarSi = dynamic_cast<TH2F*>(fOutput->FindObject("fLambdaBarSi")) ;
705 if (!fLambdaBarSi) {
706 Printf("ERROR: fLambdaBarSi not available");
707 return;
708 }
709 TH1D *lambdaBarSiPx=fLambdaBarSi->ProjectionX();
710 lambdaBarSiPx->SetName("fLambdaBarPt");
711 lambdaBarSiPx->Sumw2();
712
713 fLambdaBarEff = dynamic_cast<TH1D*>(fOutput->FindObject("fLambdaBarEff")) ;
714 if (!fLambdaBarEff) {
715 Printf("ERROR: fLambdaBarEff not available");
716 return;
717 }
718 fLambdaBarPt = dynamic_cast<TH1D*>(fOutput->FindObject("fLambdaBarPt")) ;
719 if (!fLambdaBarPt) {
720 Printf("ERROR: fLambdaBarPt not available");
721 return;
722 }
723
724
725 if (!gROOT->IsBatch()) {
726
727 TCanvas *c1 = new TCanvas("c1","Mulitplicity");
728 c1->SetLogy();
729 fMult->DrawCopy() ;
730
731 new TCanvas("c2","dE/dx");
732 fdEdx->DrawCopy() ;
733
734 new TCanvas("c3","dE/dx with PID");
735 fdEdxPid->DrawCopy() ;
736
737 if (fIsMC) {
738 /*
739 TH1D effK(*k0sAsPx); effK.SetTitle("Efficiency for K0s");
740 effK.Divide(k0sAsPx,k0sMcPx,1,1,"b");
741 new TCanvas("c4","Efficiency for K0s");
742 effK.DrawCopy("E") ;
743 */
744
745 //+++ Lambdas
746 fLambdaEff->Divide(lambdaAsPx,lambdaMcPx,1,1,"b");
747 new TCanvas("c5","Efficiency for #Lambda");
748 fLambdaEff->DrawCopy("E") ;
749
750 lambdaSiPx->Add(lambdaFromXiPx,-1);
751 lambdaSiPx->Divide(fLambdaEff);
752
753 new TCanvas("c6","Corrected #Lambda pt");
754 lambdaSiPx->SetTitle("Corrected #Lambda pt");
755 *fLambdaPt = *lambdaSiPx;
756 fLambdaPt->SetLineColor(2);
757 fLambdaPt->DrawCopy("E");
758
759 lambdaMcPx->DrawCopy("same");
760
761
762 //+++ anti-Lambdas
763 fLambdaBarEff->Divide(lambdaBarAsPx,lambdaBarMcPx,1,1,"b");
764 new TCanvas("c7","Efficiency for anti-#Lambda");
765 fLambdaBarEff->DrawCopy("E") ;
766
767 lambdaBarSiPx->Add(lambdaBarFromXiBarPx,-1);
768 lambdaBarSiPx->Divide(fLambdaBarEff);
769
770 new TCanvas("c8","Corrected anti-#Lambda pt");
771 lambdaBarSiPx->SetTitle("Corrected anti-#Lambda pt");
772 *fLambdaBarPt = *lambdaBarSiPx;
773 fLambdaBarPt->SetLineColor(2);
774 fLambdaBarPt->DrawCopy("E");
775
776 lambdaBarMcPx->DrawCopy("same");
777 } else {
778 new TCanvas("c6","Raw #Lambda pt");
779 lambdaSiPx->SetTitle("Raw #Lambda pt");
780 *fLambdaPt = *lambdaSiPx;
781 fLambdaPt->SetLineColor(2);
782 fLambdaPt->DrawCopy("E");
783
784
785 new TCanvas("c7","Raw anti-#Lambda pt");
786 lambdaBarSiPx->SetTitle("Raw anti-#Lambda pt");
787 *fLambdaBarPt = *lambdaBarSiPx;
788 fLambdaBarPt->SetLineColor(2);
789 fLambdaBarPt->DrawCopy("E");
790 }
791 }
792
793 TFile *file = 0;
794 if (fIsMC)
795 file=TFile::Open("AliV0CutVariationsMC.root", "RECREATE");
796 else
797 file=TFile::Open("AliV0CutVariations.root", "RECREATE");
798
799 fOutput->Write() ;
800 file->Close() ;
801 delete file ;
802
803}