]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG3/vertexingHF/AliCFTaskVertexingHF.cxx
New common interface to the CORRFW for all charm analyses (Davide, ChiaraZ)
[u/mrichter/AliRoot.git] / PWG3 / vertexingHF / AliCFTaskVertexingHF.cxx
CommitLineData
379592af 1/**************************************************************************
2 * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16//-----------------------------------------------------------------------
17// Class for HF corrections as a function of many variables
18// 6 Steps introduced: MC, MC Acc, Reco, Reco Acc, Reco Acc + ITS Cl,
19// Reco Acc + ITS Cl + PPR cuts
20// 12 variables used: pt, y, cosThetaStar, ptPi, ptK, ct,
21// dca, d0Pi, d0K, d0Pixd0K, cosPointingAngle, phi
22//
23//-----------------------------------------------------------------------
24// Author : C. Zampolli, CERN
25// D. Caffarri, Univ & INFN Padova caffarri@pd.infn.it
26//-----------------------------------------------------------------------
27//-----------------------------------------------------------------------
28// Base class for HF Unfolding (pt and eta)
29// correlation matrix filled at Acceptance and PPR level
30// Author: A.Grelli , Utrecht - agrelli@uu.nl
31//-----------------------------------------------------------------------
32#include <TCanvas.h>
33#include <TParticle.h>
34#include <TDatabasePDG.h>
35#include <TH1I.h>
36#include <TStyle.h>
37#include <TFile.h>
38
39#include "AliCFTaskVertexingHF.h"
40#include "AliStack.h"
41#include "AliMCEvent.h"
42#include "AliCFManager.h"
43#include "AliCFContainer.h"
44#include "AliLog.h"
45#include "AliAnalysisManager.h"
46#include "AliAODHandler.h"
47#include "AliAODEvent.h"
48#include "AliAODRecoDecay.h"
49#include "AliAODRecoDecayHF.h"
50#include "AliAODRecoDecayHF2Prong.h"
51#include "AliAODMCParticle.h"
52#include "AliAODMCHeader.h"
53#include "AliESDtrack.h"
54#include "TChain.h"
55#include "THnSparse.h"
56#include "TH2D.h"
57#include "AliESDtrackCuts.h"
58#include "AliRDHFCutsD0toKpi.h"
59#include "AliCFVertexingHF2Prong.h"
60#include "AliCFVertexingHF.h"
61
62//__________________________________________________________________________
63AliCFTaskVertexingHF::AliCFTaskVertexingHF() :
64AliAnalysisTaskSE(),
65fPDG(0),
66fCFManager(0x0),
67fHistEventsProcessed(0x0),
68fCorrelation(0x0),
69fCountMC(0),
70fCountAcc(0),
71fCountVertex(0),
72fCountRefit(0),
73fCountReco(0),
74fCountRecoAcc(0),
75fCountRecoITSClusters(0),
76fCountRecoPPR(0),
77fCountRecoPID(0),
78fEvents(0),
79fFillFromGenerated(kFALSE),
80fOriginDselection(0),
81fAcceptanceUnf(kTRUE),
82fCuts(0)
83{
84 //
85 //Default ctor
86 //
87}
88//___________________________________________________________________________
89AliCFTaskVertexingHF::AliCFTaskVertexingHF(const Char_t* name, AliRDHFCutsD0toKpi* cuts) :
90AliAnalysisTaskSE(name),
91fPDG(0),
92fCFManager(0x0),
93fHistEventsProcessed(0x0),
94fCorrelation(0x0),
95fCountMC(0),
96fCountAcc(0),
97fCountVertex(0),
98fCountRefit(0),
99fCountReco(0),
100fCountRecoAcc(0),
101fCountRecoITSClusters(0),
102fCountRecoPPR(0),
103fCountRecoPID(0),
104fEvents(0),
105fFillFromGenerated(kFALSE),
106fOriginDselection(0),
107fAcceptanceUnf(kTRUE),
108fCuts(cuts)
109{
110 //
111 // Constructor. Initialization of Inputs and Outputs
112 //
113 Info("AliCFTaskVertexingHF","Calling Constructor");
114 /*
115 DefineInput(0) and DefineOutput(0)
116 are taken care of by AliAnalysisTaskSE constructor
117 */
118 DefineOutput(1,TH1I::Class());
119 DefineOutput(2,AliCFContainer::Class());
120 DefineOutput(3,THnSparseD::Class());
121
122 fCuts->PrintAll();
123}
124
125//___________________________________________________________________________
126AliCFTaskVertexingHF& AliCFTaskVertexingHF::operator=(const AliCFTaskVertexingHF& c)
127{
128 //
129 // Assignment operator
130 //
131 if (this!=&c) {
132 AliAnalysisTaskSE::operator=(c) ;
133 fPDG = c.fPDG;
134 fCFManager = c.fCFManager;
135 fHistEventsProcessed = c.fHistEventsProcessed;
136 fCuts = c.fCuts;
137 }
138 return *this;
139}
140
141//___________________________________________________________________________
142AliCFTaskVertexingHF::AliCFTaskVertexingHF(const AliCFTaskVertexingHF& c) :
143AliAnalysisTaskSE(c),
144fPDG(c.fPDG),
145fCFManager(c.fCFManager),
146fHistEventsProcessed(c.fHistEventsProcessed),
147fCorrelation(c.fCorrelation),
148fCountMC(c.fCountMC),
149fCountAcc(c.fCountAcc),
150fCountVertex(c.fCountVertex),
151fCountRefit(c.fCountRefit),
152fCountReco(c.fCountReco),
153fCountRecoAcc(c.fCountRecoAcc),
154fCountRecoITSClusters(c.fCountRecoITSClusters),
155fCountRecoPPR(c.fCountRecoPPR),
156fCountRecoPID(c.fCountRecoPID),
157fEvents(c.fEvents),
158fFillFromGenerated(c.fFillFromGenerated),
159fOriginDselection(c.fOriginDselection),
160fAcceptanceUnf(c.fAcceptanceUnf),
161fCuts(c.fCuts)
162{
163 //
164 // Copy Constructor
165 //
166}
167
168//___________________________________________________________________________
169AliCFTaskVertexingHF::~AliCFTaskVertexingHF() {
170 //
171 //destructor
172 //
173 if (fCFManager) delete fCFManager ;
174 if (fHistEventsProcessed) delete fHistEventsProcessed ;
175 if (fCorrelation) delete fCorrelation ;
176 if (fCuts) delete fCuts;
177}
178
179//_________________________________________________
180void AliCFTaskVertexingHF::UserExec(Option_t *)
181{
182
183 //
184 // Main loop function
185 //
186
187 PostData(1,fHistEventsProcessed) ;
188 PostData(2,fCFManager->GetParticleContainer()) ;
189 PostData(3,fCorrelation) ;
190
191 AliESDtrackCuts* trackCuts = fCuts->GetTrackCuts();
192
193 if (fFillFromGenerated){
194 AliWarning("Flag to fill container with generated value ON ---> dca, d0pi, d0K, d0xd0, cosPointingAngle will be set as dummy!");
195 }
196
197 if (!fInputEvent) {
198 Error("UserExec","NO EVENT FOUND!");
199 return;
200 }
201
202 AliAODEvent* aodEvent = dynamic_cast<AliAODEvent*>(fInputEvent);
203
204 TClonesArray *arrayD0toKpi=0;
205 //TClonesArray *arrayBranch=0;
206
207 if(!aodEvent && AODEvent() && IsStandardAOD()) {
208 // In case there is an AOD handler writing a standard AOD, use the AOD
209 // event in memory rather than the input (ESD) event.
210 aodEvent = dynamic_cast<AliAODEvent*> (AODEvent());
211 // in this case the braches in the deltaAOD (AliAOD.VertexingHF.root)
212 // have to taken from the AOD event hold by the AliAODExtension
213 AliAODHandler* aodHandler = (AliAODHandler*)
214 ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
215 if(aodHandler->GetExtensions()) {
216 AliAODExtension *ext = (AliAODExtension*)aodHandler->GetExtensions()->FindObject("AliAOD.VertexingHF.root");
217 AliAODEvent *aodFromExt = ext->GetAOD();
218 arrayD0toKpi=(TClonesArray*)aodFromExt->GetList()->FindObject("D0toKpi");
219 }
220 } else {
221 arrayD0toKpi=(TClonesArray*)aodEvent->GetList()->FindObject("D0toKpi");
222 }
223
224 AliAODVertex *aodVtx = (AliAODVertex*)aodEvent->GetPrimaryVertex();
225 if (!aodVtx) return;
226
227 if (!arrayD0toKpi) {
228 AliError("Could not find array of HF vertices");
229 return;
230 }
231
232 fEvents++;
233 if (fEvents%10000 == 0) AliDebug(2,Form("Event %d",fEvents));
234
235 fCFManager->SetRecEventInfo(aodEvent);
236 fCFManager->SetMCEventInfo(aodEvent);
237
238
239 //******** DEFINE number of variables of the container***** for now set at 13, in the future in the config macro.
240 Int_t nVar = 13;
241
242 Double_t* containerInput = new Double_t[nVar];
243 Double_t* containerInputMC = new Double_t[nVar];
244
245 //loop on the MC event
246
247 TClonesArray* mcArray = dynamic_cast<TClonesArray*>(aodEvent->FindListObject(AliAODMCParticle::StdBranchName()));
248 if (!mcArray) {
249 AliError("Could not find Monte-Carlo in AOD");
250 return;
251 }
252 Int_t icountMC = 0;
253 Int_t icountAcc = 0;
254 Int_t icountReco = 0;
255 Int_t icountVertex = 0;
256 Int_t icountRefit = 0;
257 Int_t icountRecoAcc = 0;
258 Int_t icountRecoITSClusters = 0;
259 Int_t icountRecoPPR = 0;
260 Int_t icountRecoPID = 0;
261
262 Int_t cquarks = 0;
263
264 // These flgs will be setted i the config macro
265 UShort_t originDselection = 0;
266
267
268 AliAODMCHeader *mcHeader = dynamic_cast<AliAODMCHeader*>(aodEvent->GetList()->FindObject(AliAODMCHeader::StdBranchName()));
269 if (!mcHeader) {
270 AliError("Could not find MC Header in AOD");
271 return;
272 }
273
274
275 AliCFVertexingHF2Prong *cfVtxHF = new AliCFVertexingHF2Prong(mcArray, originDselection);
276
277 Double_t zPrimVertex = aodVtx ->GetZ();
278 Double_t zMCVertex = mcHeader->GetVtxZ();
279
280 //General settings: vertex, feed down and fill reco container with generated values.
281 cfVtxHF->SetRecoPrimVertex(zPrimVertex);
282 cfVtxHF->SetMCPrimaryVertex(zMCVertex);
283 cfVtxHF->SetFillFromGenerated(fFillFromGenerated);
284
285 // cfVtxHF->SetKeepD0fromBOnly(fKeepD0fromBOnly);
286 // cfVtxHF->SetKeepD0fromB(fKeepD0fromB);
287
288 for (Int_t iPart=0; iPart<mcArray->GetEntriesFast(); iPart++) {
289
290 AliAODMCParticle* mcPart = dynamic_cast<AliAODMCParticle*>(mcArray->At(iPart));
291
292 // check the MC-level cuts, must be a D0
293 if (!fCFManager->CheckParticleCuts(0, mcPart)) continue; // 0 stands for MC level
294
295 cfVtxHF->SetMCCandidateParam(iPart);
296 cfVtxHF->SetNVar(nVar);
297 //counting c quarks
298 cquarks += cfVtxHF->MCcquarkCounting(mcPart);
299
300 //check the candiate family at MV level
301 if (!(cfVtxHF->CheckMCPartFamily(mcPart, mcArray))) continue;
302
303 //Fill the MC container
304 Bool_t mcContainerFilled = cfVtxHF -> FillMCContainer(containerInputMC);
305
306 if (!fCuts->IsInFiducialAcceptance(containerInputMC[0],containerInputMC[1])) continue;
307 if (TMath::Abs(containerInputMC[1]) < 0.5) {
308 fCFManager->GetParticleContainer()->Fill(containerInputMC,kStepGeneratedLimAcc);
309 }
310
311 if (mcContainerFilled) {
312
313 //MC
314 fCFManager->GetParticleContainer()->Fill(containerInputMC, kStepGenerated);
315 icountMC++;
316 printf("MC cointainer filled \n");
317
318 // MC in acceptance
319 // check the MC-Acceptance level cuts
320 // since standard CF functions are not applicable, using Kine Cuts on daughters
321 Bool_t mcAccepStep = cfVtxHF-> MCAcceptanceStep();
322 if (mcAccepStep){
323 fCFManager->GetParticleContainer()->Fill(containerInputMC,kStepAcceptance);
324 printf("MC acceptance cut passed\n");
325 icountAcc++;
326
327 //MC Vertex step
328 if (fCuts->IsEventSelected(aodEvent)){
329 // filling the container if the vertex is ok
330 fCFManager->GetParticleContainer()->Fill(containerInputMC,kStepVertex) ;
331 printf("Vertex cut passed and container filled\n");
332 icountVertex++;
333
334 //mc Refit requirement
335 Bool_t mcRefitStep = cfVtxHF->MCRefitStep(aodEvent, trackCuts);
336 if (mcRefitStep){
337 fCFManager->GetParticleContainer()->Fill(containerInputMC,kStepRefit);
338 printf("MC Refit cut passed and container filled\n");
339 icountRefit++;
340 }
341 else{
342 AliDebug(3,"MC Refit cut not passed\n");
343 continue;
344 }
345
346 }
347 else{
348 AliDebug (3, "MC vertex step not passed\n");
349 continue;
350 }
351 }
352 else{
353 AliDebug (3, "MC in acceptance step not passed\n");
354 continue;
355 }
356
357 }
358 else {
359 AliDebug (3, "MC container not filled\n");
360 }
361 }
362
363 if (cquarks<2) AliDebug(2,Form("Eveith %d c-quarks", cquarks));
364 AliDebug(2,Form("Found %i MC particles that are D0!!",icountMC));
365 AliDebug(2,Form("Found %i MC particles that are D0 and satisfy Acc cuts!!",icountAcc));
366 AliDebug(2,Form("Found %i MC particles that are D0 and satisfy Vertex cuts!!",icountVertex));
367 AliDebug(2,Form("Found %i MC particles that are D0 and satisfy Refit cuts!!",icountRefit));
368
369 // Now go to rec level
370 fCountMC += icountMC;
371 fCountAcc += icountAcc;
372 fCountVertex+= icountVertex;
373 fCountRefit+= icountRefit;
374
375 AliDebug(2, Form("Found %d vertices",arrayD0toKpi->GetEntriesFast()));
376
377
378 for (Int_t iD0toKpi = 0; iD0toKpi<arrayD0toKpi->GetEntriesFast(); iD0toKpi++) {
379
380 AliAODRecoDecayHF2Prong* d0tokpi = (AliAODRecoDecayHF2Prong*)arrayD0toKpi->At(iD0toKpi);
381 Bool_t unsetvtx=kFALSE;
382 if(!d0tokpi->GetOwnPrimaryVtx()) {
383 d0tokpi->SetOwnPrimaryVtx(aodVtx); // needed to compute all variables
384 unsetvtx=kTRUE;
385 }
386
387 if (d0tokpi->GetPrimaryVtx()) printf("candidate has prim vtx \n");
388 if (aodEvent->GetPrimaryVertex()) printf ("aod event vertex\n");
389
390 Bool_t signAssociation = cfVtxHF->SetRecoCandidateParam((AliAODRecoDecayHF*)d0tokpi);
391 if (!signAssociation){
392 d0tokpi = 0x0;
393 continue;
394 }
395
396 Bool_t recoContFilled = cfVtxHF->FillRecoContainer(containerInput);
397 if (recoContFilled){
398
399 if (!fCuts->IsInFiducialAcceptance(containerInput[0],containerInput[1])) continue;
400
401 //Reco Step
402 Bool_t recoStep = cfVtxHF->RecoStep();
403 Bool_t vtxCheck = fCuts->IsEventSelected(aodEvent);
404
405 if (recoStep && recoContFilled && vtxCheck){
406 fCFManager->GetParticleContainer()->Fill(containerInput,kStepReconstructed) ;
407 icountReco++;
408 printf("Reco step passed and container filled\n");
409
410
411 //Reco in the acceptance -- take care of UNFOLDING!!!!
412 Bool_t recoAcceptanceStep = cfVtxHF->RecoAcceptStep(trackCuts);
413 if (recoAcceptanceStep) {
414 fCFManager->GetParticleContainer()->Fill(containerInput,kStepRecoAcceptance) ;
415 icountRecoAcc++;
416 printf("Reco acceptance cut passed and container filled\n");
417
418 if(fAcceptanceUnf){
419 Double_t fill[4]; //fill response matrix
420 Bool_t bUnfolding = cfVtxHF -> FillUnfoldingMatrix(fill);
421 if (bUnfolding) fCorrelation->Fill(fill);
422 }
423
424 //Number of ITS cluster requirements
425 Int_t recoITSnCluster = fCuts->IsSelected(d0tokpi, AliRDHFCuts::kTracks);
426 if (recoITSnCluster){
427 fCFManager->GetParticleContainer()->Fill(containerInput,kStepRecoITSClusters) ;
428 icountRecoITSClusters++;
429 printf("Reco n ITS cluster cut passed and container filled\n");
430
431
432 Int_t recoAnalysisCuts = fCuts->IsSelected(d0tokpi, AliRDHFCuts::kCandidate);
433 if(recoAnalysisCuts){
434 fCFManager->GetParticleContainer()->Fill(containerInput, kStepRecoPPR);
435 icountRecoPPR++;
436 printf ("Reco Analysis cuts passed and container filled \n");
437
438
439 //pid selection
440 Int_t recoPidSelection = fCuts->IsSelected(d0tokpi, AliRDHFCuts::kPID);
441 if (recoPidSelection > 0){
442 fCFManager->GetParticleContainer()->Fill(containerInput, kStepRecoPID);
443 icountRecoPID++;
444 printf ("Reco PID cuts passed and container filled \n");
445
446
447 if(!fAcceptanceUnf){
448 Double_t fill[4]; //fill response matrix
449 Bool_t bUnfolding = cfVtxHF -> FillUnfoldingMatrix(fill);
450 if (bUnfolding) fCorrelation->Fill(fill);
451 }
452 }
453 else {
454 AliDebug(3, "Analysis Cuts step not passed \n");
455 continue;
456 }
457 }
458 else {
459 AliDebug(3, "PID selection not passed \n");
460 continue;
461 }
462 }
463 else{
464 AliDebug(3, "Number of ITS cluster step not passed\n");
465 continue;
466 }
467 }
468 else{
469 AliDebug(3, "Reco acceptance step not passed\n");
470 continue;
471 }
472 }
473 else {
474 AliDebug(3, "Reco step not passed\n");
475 continue;
476 }
477 }
478
479 if(unsetvtx) d0tokpi->UnsetOwnPrimaryVtx();
480 } // end loop on D0->Kpi
481
482 AliDebug(2, Form("Found %i Reco particles that are D0!!",icountReco));
483
484 fCountReco+= icountReco;
485 fCountRecoAcc+= icountRecoAcc;
486 fCountRecoITSClusters+= icountRecoITSClusters;
487 fCountRecoPPR+= icountRecoPPR;
488 fCountRecoPID+= icountRecoPID;
489
490 fHistEventsProcessed->Fill(0);
491}
492
493//___________________________________________________________________________
494void AliCFTaskVertexingHF::Terminate(Option_t*)
495{
496 // The Terminate() function is the last function to be called during
497 // a query. It always runs on the client, it can be used to present
498 // the results graphically or save the results to file.
499
500 AliAnalysisTaskSE::Terminate();
501
502 AliInfo(Form("Found %i MC particles that are D0 in MC, in %d events",fCountMC,fEvents));
503 AliInfo(Form("Found %i MC particles that are D0 in MC and satisfy Acc cuts, in %d events",fCountAcc,fEvents));
504 AliInfo(Form("Found %i MC particles that are D0 in MC and satisfy Acc cuts, and satisfy Vertex requirement in %d events",fCountVertex,fEvents));
505 AliInfo(Form("Found %i MC particles that are D0 in MC and satisfy Acc cuts, and satisfy ITS+TPC refit requirementin %d events",fCountRefit,fEvents));
506 AliInfo(Form("Found %i reco D0 that are decaying in K+pi, in %d events",fCountReco,fEvents));
507 AliInfo(Form("Among the above, found %i reco D0 that are decaying in K+pi and are in the requested acceptance, in %d events",fCountRecoAcc,fEvents));
508 AliInfo(Form("Among the above, found %i reco D0 that are decaying in K+pi and have at least 5 clusters in ITS, in %d events",fCountRecoITSClusters,fEvents));
509 AliInfo(Form("Among the above, found %i reco D0 that are decaying in K+pi and satisfy PPR cuts, in %d events",fCountRecoPPR,fEvents));
510
511 // draw some example plots....
512
513 AliCFContainer *cont= dynamic_cast<AliCFContainer*> (GetOutputData(2));
514 if(!cont) {
515 printf("CONTAINER NOT FOUND\n");
516 return;
517 }
518 // projecting the containers to obtain histograms
519 // first argument = variable, second argument = step
520
521 // MC-level
522 TH1D* h00 = cont->ShowProjection(0,0) ; // pt
523 TH1D* h10 = cont->ShowProjection(1,0) ; // rapidity
524 TH1D* h20 = cont->ShowProjection(2,0) ; // cosThetaStar
525 TH1D* h30 = cont->ShowProjection(3,0) ; // pTpi
526 TH1D* h40 = cont->ShowProjection(4,0) ; // pTK
527 TH1D* h50 = cont->ShowProjection(5,0) ; // cT
528 TH1D* h60 = cont->ShowProjection(6,0) ; // dca
529 TH1D* h70 = cont->ShowProjection(7,0) ; // d0pi
530 TH1D* h80 = cont->ShowProjection(8,0) ; // d0K
531 TH1D* h90 = cont->ShowProjection(9,0) ; // d0xd0
532 TH1D* h100 = cont->ShowProjection(10,0) ; // cosPointingAngle
533 TH1D* h110 = cont->ShowProjection(11,0) ; // phi
534
535 // MC-Acceptance level
536 TH1D* h01 = cont->ShowProjection(0,1) ; // pt
537 TH1D* h11 = cont->ShowProjection(1,1) ; // rapidity
538 TH1D* h21 = cont->ShowProjection(2,1) ; // cosThetaStar
539 TH1D* h31 = cont->ShowProjection(3,1) ; // pTpi
540 TH1D* h41 = cont->ShowProjection(4,1) ; // pTK
541 TH1D* h51 = cont->ShowProjection(5,1) ; // cT
542 TH1D* h61 = cont->ShowProjection(6,1) ; // dca
543 TH1D* h71 = cont->ShowProjection(7,1) ; // d0pi
544 TH1D* h81 = cont->ShowProjection(8,1) ; // d0K
545 TH1D* h91 = cont->ShowProjection(9,1) ; // d0xd0
546 TH1D* h101 = cont->ShowProjection(10,1) ; // cosPointingAngle
547 TH1D* h111 = cont->ShowProjection(11,1) ; // phi
548
549 // Reco-level
550 TH1D* h04 = cont->ShowProjection(0,4) ; // pt
551 TH1D* h14 = cont->ShowProjection(1,4) ; // rapidity
552 TH1D* h24 = cont->ShowProjection(2,4) ; // cosThetaStar
553 TH1D* h34 = cont->ShowProjection(3,4) ; // pTpi
554 TH1D* h44 = cont->ShowProjection(4,4) ; // pTK
555 TH1D* h54 = cont->ShowProjection(5,4) ; // cT
556 TH1D* h64 = cont->ShowProjection(6,4) ; // dca
557 TH1D* h74 = cont->ShowProjection(7,4) ; // d0pi
558 TH1D* h84 = cont->ShowProjection(8,4) ; // d0K
559 TH1D* h94 = cont->ShowProjection(9,4) ; // d0xd0
560 TH1D* h104 = cont->ShowProjection(10,4) ; // cosPointingAngle
561 TH1D* h114 = cont->ShowProjection(11,4) ; // phi
562
563 h00->SetTitle("pT_D0 (GeV/c)");
564 h10->SetTitle("rapidity");
565 h20->SetTitle("cosThetaStar");
566 h30->SetTitle("pT_pi (GeV/c)");
567 h40->SetTitle("pT_K (Gev/c)");
568 h50->SetTitle("cT (#mum)");
569 h60->SetTitle("dca (#mum)");
570 h70->SetTitle("d0_pi (#mum)");
571 h80->SetTitle("d0_K (#mum)");
572 h90->SetTitle("d0xd0 (#mum^2)");
573 h100->SetTitle("cosPointingAngle");
574 h100->SetTitle("phi (rad)");
575
576 h00->GetXaxis()->SetTitle("pT_D0 (GeV/c)");
577 h10->GetXaxis()->SetTitle("rapidity");
578 h20->GetXaxis()->SetTitle("cosThetaStar");
579 h30->GetXaxis()->SetTitle("pT_pi (GeV/c)");
580 h40->GetXaxis()->SetTitle("pT_K (Gev/c)");
581 h50->GetXaxis()->SetTitle("cT (#mum)");
582 h60->GetXaxis()->SetTitle("dca (#mum)");
583 h70->GetXaxis()->SetTitle("d0_pi (#mum)");
584 h80->GetXaxis()->SetTitle("d0_K (#mum)");
585 h90->GetXaxis()->SetTitle("d0xd0 (#mum^2)");
586 h100->GetXaxis()->SetTitle("cosPointingAngle");
587 h110->GetXaxis()->SetTitle("phi (rad)");
588
589 h01->SetTitle("pT_D0 (GeV/c)");
590 h11->SetTitle("rapidity");
591 h21->SetTitle("cosThetaStar");
592 h31->SetTitle("pT_pi (GeV/c)");
593 h41->SetTitle("pT_K (Gev/c)");
594 h51->SetTitle("cT (#mum)");
595 h61->SetTitle("dca (#mum)");
596 h71->SetTitle("d0_pi (#mum)");
597 h81->SetTitle("d0_K (#mum)");
598 h91->SetTitle("d0xd0 (#mum^2)");
599 h101->SetTitle("cosPointingAngle");
600 h111->GetXaxis()->SetTitle("phi (rad)");
601
602 h01->GetXaxis()->SetTitle("pT_D0 (GeV/c)");
603 h11->GetXaxis()->SetTitle("rapidity");
604 h21->GetXaxis()->SetTitle("cosThetaStar");
605 h31->GetXaxis()->SetTitle("pT_pi (GeV/c)");
606 h41->GetXaxis()->SetTitle("pT_K (Gev/c)");
607 h51->GetXaxis()->SetTitle("cT (#mum)");
608 h61->GetXaxis()->SetTitle("dca (#mum)");
609 h71->GetXaxis()->SetTitle("d0_pi (#mum)");
610 h81->GetXaxis()->SetTitle("d0_K (#mum)");
611 h91->GetXaxis()->SetTitle("d0xd0 (#mum^2)");
612 h101->GetXaxis()->SetTitle("cosPointingAngle");
613 h111->GetXaxis()->SetTitle("phi (rad)");
614
615 h04->SetTitle("pT_D0 (GeV/c)");
616 h14->SetTitle("rapidity");
617 h24->SetTitle("cosThetaStar");
618 h34->SetTitle("pT_pi (GeV/c)");
619 h44->SetTitle("pT_K (Gev/c)");
620 h54->SetTitle("cT (#mum)");
621 h64->SetTitle("dca (#mum)");
622 h74->SetTitle("d0_pi (#mum)");
623 h84->SetTitle("d0_K (#mum)");
624 h94->SetTitle("d0xd0 (#mum^2)");
625 h104->SetTitle("cosPointingAngle");
626 h114->GetXaxis()->SetTitle("phi (rad)");
627
628 h04->GetXaxis()->SetTitle("pT_D0 (GeV/c)");
629 h14->GetXaxis()->SetTitle("rapidity");
630 h24->GetXaxis()->SetTitle("cosThetaStar");
631 h34->GetXaxis()->SetTitle("pT_pi (GeV/c)");
632 h44->GetXaxis()->SetTitle("pT_K (Gev/c)");
633 h54->GetXaxis()->SetTitle("cT (#mum)");
634 h64->GetXaxis()->SetTitle("dca (#mum)");
635 h74->GetXaxis()->SetTitle("d0_pi (#mum)");
636 h84->GetXaxis()->SetTitle("d0_K (#mum)");
637 h94->GetXaxis()->SetTitle("d0xd0 (#mum^2)");
638 h104->GetXaxis()->SetTitle("cosPointingAngle");
639 h114->GetXaxis()->SetTitle("phi (rad)");
640
641 Double_t max0 = h00->GetMaximum();
642 Double_t max1 = h10->GetMaximum();
643 Double_t max2 = h20->GetMaximum();
644 Double_t max3 = h30->GetMaximum();
645 Double_t max4 = h40->GetMaximum();
646 Double_t max5 = h50->GetMaximum();
647 Double_t max6 = h60->GetMaximum();
648 Double_t max7 = h70->GetMaximum();
649 Double_t max8 = h80->GetMaximum();
650 Double_t max9 = h90->GetMaximum();
651 Double_t max10 = h100->GetMaximum();
652 Double_t max11 = h110->GetMaximum();
653
654 h00->GetYaxis()->SetRangeUser(0,max0*1.2);
655 h10->GetYaxis()->SetRangeUser(0,max1*1.2);
656 h20->GetYaxis()->SetRangeUser(0,max2*1.2);
657 h30->GetYaxis()->SetRangeUser(0,max3*1.2);
658 h40->GetYaxis()->SetRangeUser(0,max4*1.2);
659 h50->GetYaxis()->SetRangeUser(0,max5*1.2);
660 h60->GetYaxis()->SetRangeUser(0,max6*1.2);
661 h70->GetYaxis()->SetRangeUser(0,max7*1.2);
662 h80->GetYaxis()->SetRangeUser(0,max8*1.2);
663 h90->GetYaxis()->SetRangeUser(0,max9*1.2);
664 h100->GetYaxis()->SetRangeUser(0,max10*1.2);
665 h110->GetYaxis()->SetRangeUser(0,max11*1.2);
666
667 h01->GetYaxis()->SetRangeUser(0,max0*1.2);
668 h11->GetYaxis()->SetRangeUser(0,max1*1.2);
669 h21->GetYaxis()->SetRangeUser(0,max2*1.2);
670 h31->GetYaxis()->SetRangeUser(0,max3*1.2);
671 h41->GetYaxis()->SetRangeUser(0,max4*1.2);
672 h51->GetYaxis()->SetRangeUser(0,max5*1.2);
673 h61->GetYaxis()->SetRangeUser(0,max6*1.2);
674 h71->GetYaxis()->SetRangeUser(0,max7*1.2);
675 h81->GetYaxis()->SetRangeUser(0,max8*1.2);
676 h91->GetYaxis()->SetRangeUser(0,max9*1.2);
677 h101->GetYaxis()->SetRangeUser(0,max10*1.2);
678 h111->GetYaxis()->SetRangeUser(0,max11*1.2);
679
680 h00->SetMarkerStyle(20);
681 h10->SetMarkerStyle(24);
682 h20->SetMarkerStyle(21);
683 h30->SetMarkerStyle(25);
684 h40->SetMarkerStyle(27);
685 h50->SetMarkerStyle(28);
686 h60->SetMarkerStyle(20);
687 h70->SetMarkerStyle(24);
688 h80->SetMarkerStyle(21);
689 h90->SetMarkerStyle(25);
690 h100->SetMarkerStyle(27);
691 h110->SetMarkerStyle(28);
692
693 h00->SetMarkerColor(2);
694 h10->SetMarkerColor(2);
695 h20->SetMarkerColor(2);
696 h30->SetMarkerColor(2);
697 h40->SetMarkerColor(2);
698 h50->SetMarkerColor(2);
699 h60->SetMarkerColor(2);
700 h70->SetMarkerColor(2);
701 h80->SetMarkerColor(2);
702 h90->SetMarkerColor(2);
703 h100->SetMarkerColor(2);
704 h110->SetMarkerColor(2);
705
706 h01->SetMarkerStyle(20) ;
707 h11->SetMarkerStyle(24) ;
708 h21->SetMarkerStyle(21) ;
709 h31->SetMarkerStyle(25) ;
710 h41->SetMarkerStyle(27) ;
711 h51->SetMarkerStyle(28) ;
712 h61->SetMarkerStyle(20);
713 h71->SetMarkerStyle(24);
714 h81->SetMarkerStyle(21);
715 h91->SetMarkerStyle(25);
716 h101->SetMarkerStyle(27);
717 h111->SetMarkerStyle(28);
718
719 h01->SetMarkerColor(8);
720 h11->SetMarkerColor(8);
721 h21->SetMarkerColor(8);
722 h31->SetMarkerColor(8);
723 h41->SetMarkerColor(8);
724 h51->SetMarkerColor(8);
725 h61->SetMarkerColor(8);
726 h71->SetMarkerColor(8);
727 h81->SetMarkerColor(8);
728 h91->SetMarkerColor(8);
729 h101->SetMarkerColor(8);
730 h111->SetMarkerColor(8);
731
732 h04->SetMarkerStyle(20) ;
733 h14->SetMarkerStyle(24) ;
734 h24->SetMarkerStyle(21) ;
735 h34->SetMarkerStyle(25) ;
736 h44->SetMarkerStyle(27) ;
737 h54->SetMarkerStyle(28) ;
738 h64->SetMarkerStyle(20);
739 h74->SetMarkerStyle(24);
740 h84->SetMarkerStyle(21);
741 h94->SetMarkerStyle(25);
742 h104->SetMarkerStyle(27);
743 h114->SetMarkerStyle(28);
744
745 h04->SetMarkerColor(4);
746 h14->SetMarkerColor(4);
747 h24->SetMarkerColor(4);
748 h34->SetMarkerColor(4);
749 h44->SetMarkerColor(4);
750 h54->SetMarkerColor(4);
751 h64->SetMarkerColor(4);
752 h74->SetMarkerColor(4);
753 h84->SetMarkerColor(4);
754 h94->SetMarkerColor(4);
755 h104->SetMarkerColor(4);
756 h114->SetMarkerColor(4);
757
758 gStyle->SetCanvasColor(0);
759 gStyle->SetFrameFillColor(0);
760 gStyle->SetTitleFillColor(0);
761 gStyle->SetStatColor(0);
762
763 // drawing in 2 separate canvas for a matter of clearity
764 TCanvas * c1 =new TCanvas("c1New","pT, rapidiy, cosThetaStar",1100,1600);
765 c1->Divide(3,3);
766
767 c1->cd(1);
768 h00->Draw("p");
769 c1->cd(1);
770 c1->cd(2);
771 h01->Draw("p");
772 c1->cd(2);
773 c1->cd(3);
774 h04->Draw("p");
775 c1->cd(3);
776 c1->cd(4);
777 h10->Draw("p");
778 c1->cd(4);
779 c1->cd(5);
780 h11->Draw("p");
781 c1->cd(5);
782 c1->cd(6);
783 h14->Draw("p");
784 c1->cd(6);
785 c1->cd(7);
786 h20->Draw("p");
787 c1->cd(7);
788 c1->cd(8);
789 h21->Draw("p");
790 c1->cd(8);
791 c1->cd(9);
792 h24->Draw("p");
793 c1->cd(9);
794 c1->cd();
795
796 TCanvas * c2 =new TCanvas("c2New","pTpi, pTK, cT",1100,1600);
797 c2->Divide(3,3);
798 c2->cd(1);
799 h30->Draw("p");
800 c2->cd(1);
801 c2->cd(2);
802 h31->Draw("p");
803 c2->cd(2);
804 c2->cd(3);
805 h34->Draw("p");
806 c2->cd(3);
807 c2->cd(4);
808 h40->Draw("p");
809 c2->cd(4);
810 c2->cd(5);
811 h41->Draw("p");
812 c2->cd(5);
813 c2->cd(6);
814 h44->Draw("p");
815 c2->cd(6);
816 c2->cd(7);
817 h50->Draw("p");
818 c2->cd(7);
819 c2->cd(8);
820 h51->Draw("p");
821 c2->cd(8);
822 c2->cd(9);
823 h54->Draw("p");
824 c2->cd(9);
825 c2->cd();
826
827 TCanvas * c3 =new TCanvas("c3New","dca, d0pi, d0K",1100,1600);
828 c3->Divide(3,3);
829 c3->cd(1);
830 h60->Draw("p");
831 c3->cd(1);
832 c3->cd(2);
833 h61->Draw("p");
834 c3->cd(2);
835 c3->cd(3);
836 h64->Draw("p");
837 c3->cd(3);
838 c3->cd(4);
839 h70->Draw("p");
840 c3->cd(4);
841 c3->cd(5);
842 h71->Draw("p");
843 c3->cd(5);
844 c3->cd(6);
845 h74->Draw("p");
846 c3->cd(6);
847 c3->cd(7);
848 h80->Draw("p");
849 c3->cd(7);
850 c3->cd(8);
851 h81->Draw("p");
852 c3->cd(8);
853 c3->cd(9);
854 h84->Draw("p");
855 c3->cd(9);
856 c3->cd();
857
858 TCanvas * c4 =new TCanvas("c4New","d0xd0, cosPointingAngle, phi",1100,1600);
859 c4->Divide(3,3);
860 c4->cd(1);
861 h90->Draw("p");
862 c4->cd(1);
863 c4->cd(2);
864 h91->Draw("p");
865 c4->cd(2);
866 c4->cd(3);
867 h94->Draw("p");
868 c4->cd(3);
869 c4->cd(4);
870 h100->Draw("p");
871 c4->cd(4);
872 c4->cd(5);
873 h101->Draw("p");
874 c4->cd(5);
875 c4->cd(6);
876 h104->Draw("p");
877 c4->cd(6);
878 c4->cd(7);
879 h110->Draw("p");
880 c4->cd(7);
881 c4->cd(8);
882 h111->Draw("p");
883 c4->cd(8);
884 c4->cd(9);
885 h114->Draw("p");
886 c4->cd(9);
887 c4->cd();
888
889 THnSparseD* hcorr = dynamic_cast<THnSparseD*> (GetOutputData(3));
890
891 TH2D* corr1 =hcorr->Projection(0,2);
892 TH2D* corr2 = hcorr->Projection(1,3);
893
894 TCanvas * c7 =new TCanvas("c7New","",800,400);
895 c7->Divide(2,1);
896 c7->cd(1);
897 corr1->Draw("text");
898 c7->cd(2);
899 corr2->Draw("text");
900
901
902 TFile* file_projection = new TFile("CFtaskHFprojectionNew.root","RECREATE");
903
904 corr1->Write();
905 corr2->Write();
906 h00->Write("pT_D0_step0");
907 h10->Write("rapidity_step0");
908 h20->Write("cosThetaStar_step0");
909 h30->Write("pT_pi_step0");
910 h40->Write("pT_K_step0");
911 h50->Write("cT_step0");
912 h60->Write("dca_step0");
913 h70->Write("d0_pi_step0");
914 h80->Write("d0_K_step0");
915 h90->Write("d0xd0_step0");
916 h100->Write("cosPointingAngle_step0");
917 h110->Write("phi_step0");
918
919 h01->Write("pT_D0_step1");
920 h11->Write("rapidity_step1");
921 h21->Write("cosThetaStar_step1");
922 h31->Write("pT_pi_step1");
923 h41->Write("pT_K_step1");
924 h51->Write("cT_step1");
925 h61->Write("dca_step1");
926 h71->Write("d0_pi_step1");
927 h81->Write("d0_K_step1");
928 h91->Write("d0xd0_step1");
929 h101->Write("cosPointingAngle_step1");
930 h111->Write("phi_step1");
931
932 h04->Write("pT_D0_step2");
933 h14->Write("rapidity_step2");
934 h24->Write("cosThetaStar_step2");
935 h34->Write("pT_pi_step2");
936 h44->Write("pT_K_step2");
937 h54->Write("cT_step2");
938 h64->Write("dca_step2");
939 h74->Write("d0_pi_step2");
940 h80->Write("d0_K_step2");
941 h94->Write("d0xd0_step2");
942 h104->Write("cosPointingAngle_step2");
943 h114->Write("phi_step2");
944
945 file_projection->Close();
946
947
948
949 /*
950 c1->SaveAs("Plots/pT_rapidity_cosThetaStar.eps");
951 c2->SaveAs("Plots/pTpi_pTK_cT.eps");
952 c3->SaveAs("Plots/dca_d0pi_d0TK.eps");
953 c4->SaveAs("Plots/d0xd0_cosPointingAngle.eps");
954
955 c1->SaveAs("Plots/pT_rapidity_cosThetaStar.gif");
956 c2->SaveAs("Plots/pTpi_pTK_cT.gif");
957 c3->SaveAs("Plots/dca_d0pi_d0TK.gif");
958 c4->SaveAs("Plots/d0xd0_cosPointingAngle.gif");
959 */
960}
961
962//___________________________________________________________________________
963void AliCFTaskVertexingHF::UserCreateOutputObjects() {
964 //HERE ONE CAN CREATE OUTPUT OBJECTS, IN PARTICULAR IF THE OBJECT PARAMETERS DON'T NEED
965 //TO BE SET BEFORE THE EXECUTION OF THE TASK
966 //
967 Info("UserCreateOutputObjects","CreateOutputObjects of task %s\n", GetName());
968
969 //slot #1
970 OpenFile(1);
971 fHistEventsProcessed = new TH1I("CFHFchist0","",1,0,1) ;
972}
973