]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG3/AliAnalysisVertexingHFTest.C
Move AliD0toKpi* and AliBtoJPSI* from libPWG3base to libPWG3 (Andrea)
[u/mrichter/AliRoot.git] / PWG3 / AliAnalysisVertexingHFTest.C
CommitLineData
6a213b59 1//--------------------------------------------------------------------------
2// Test macro for reconstruction of heavy-flavour vertexing candidates
3//
7456a9e7 4// mode: 0 -> input ESD, output AOD
5// 1 -> input ESD, output TTree
6//
6a213b59 7// Andrea Dainese, andrea.dainese@lnl.infn.it
8//--------------------------------------------------------------------------
9
7456a9e7 10void AliAnalysisVertexingHFTest(Int_t mode=0)
11{
6a213b59 12
7456a9e7 13 gSystem->Load("libANALYSIS.so");
6a213b59 14 gSystem->Load("libAOD.so");
15 gSystem->Load("libPWG3base.so");
16
17 AliAnalysisVertexingHF *vHF = new AliAnalysisVertexingHF();
7456a9e7 18 //vHF->SetUseTRef();
6a213b59 19 //--- switch-off candidates finding (default: all on)
20 //vHF->SetD0toKpiOff();
7456a9e7 21 vHF->SetJPSItoEleOff();
b82f6d67 22 vHF->Set3ProngOff();
23 vHF->Set4ProngOff();
24 //--- secondary vertex with KF?
7456a9e7 25 vHF->SetSecVtxWithKF();
6a213b59 26 //--- set cuts for single-track selection
27 vHF->SetITSrefitRequired();
28 vHF->SetBothSPDNotRequired();
29 vHF->SetMinITSCls(5);
7456a9e7 30 vHF->SetMinPtCut(0.);
6a213b59 31 vHF->SetMind0Cut(0.);
32 //--- set cuts for candidates selection
33 //vHF->SetD0toKpiCuts();
34 //vHF->SetBtoJPSICuts();
35 //vHF->SetDplusCuts();
36 //--- set this if you want to reconstruct primary vertex candidate by
37 // candidate using other tracks in the event (for pp, broad
38 // interaction region)
39 //vHF->SetRecoPrimVtxSkippingTrks();
40 //--- OR set this if you want to remove the candidate daughters from
41 // the primary vertex, without recostructing it from scratch
42 //vHF->SetRmTrksFromPrimVtx();
43
44 //--- check the settings
45 vHF->PrintStatus();
46 //--- verbose
7456a9e7 47 //vHF->SetDebug(1);
6a213b59 48
7456a9e7 49 switch(mode) {
50 case 0:
51 vHF->FindCandidatesESDtoAOD();
52 break;
53 case 1:
54 TTree *trees = new TTree[1];
55 AliAODRecoDecayHF2Prong *rd2=0;
56 // AliAODRecoDecayHF3Prong *rd3=0;
57 //AliAODRecoDecayHF4Prong *rd4=0;
58 trees[0].Branch("D0toKpi","AliAODRecoDecayHF2Prong",&rd2);
59 //trees[1].Branch("JPSItoEle","AliAODRecoDecayHF2Prong",&rd2);
60 //trees[2].Branch("Charmto3Prong","AliAODRecoDecayHF3Prong",&rd3);
61 //trees[3].Branch("D0to4Prong","AliAODRecoDecayHF4Prong",&rd4);
62
63 TFile *inesd = new TFile("AliESDs.root");
64 AliESDEvent *esd = new AliESDEvent();
65 TTree *esdTree = (TTree*)inesd->Get("esdTree");
66 esd->ReadFromTree(esdTree);
67
68 for(Int_t i=0; i<esdTree->GetEntries(); i++) {
69 esdTree->GetEvent(i);
70 vHF->FindCandidates(esd,trees);
71 //if(i==evLast) i=evFirst;
72 }
73
74 delete esdTree;
75 inesd->Close();
76 delete inesd;
77
78 // Write trees with candidates
79 TFile *fout = new TFile("AliAnalysisVertexingHF.root","recreate");
80 trees[0].Write("TreeD0toKpi");
81 //trees[1].Write("TreeJPSItoEle");
82 //trees[2].Write("TreeCharm3Prong");
83 //trees[3].Write("TreeD0to4Prong");
84 fout->Close();
85 delete fout;
86 break;
6a213b59 87 }
88
6a213b59 89 return;
90}