]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG3/ReadAODVertexingHF.C
Example macro to read AOD with HF candidates (Andrea)
[u/mrichter/AliRoot.git] / PWG3 / ReadAODVertexingHF.C
CommitLineData
fcca9870 1void ReadAODVertexingHF(const char *aodFileName="AliAOD.root",
2 const char *aodHFFileName="AliAOD.VertexingHF.root")
3{
4 //
5 // Example macro to read D0->Kpi candidates from AOD and apply cuts
6 // Origin: A.Dainese
7 //
8 gSystem->Load("libANALYSIS.so");
9 gSystem->Load("libANALYSISalice.so");
10 gSystem->Load("libAOD.so");
11 gSystem->Load("libPWG3base.so");
12
13 // create a test histogram
14 TH2F *hCPtaVSd0d0 = new TH2F("hCPtaVSd0d0","D^{0} correlation plot",1000,-50000,50000,1000,-1,1);
15 hCPtaVSd0d0->SetXTitle("Product of impact parameters [#mu m^{2}]");
16 hCPtaVSd0d0->SetYTitle("Cosine of pointing angle");
17 TH1F *hMass = new TH1F("hMass","D^{0} mass plot",100,1.7,2);
18 hMass->SetXTitle("Invariant mass [GeV]");
19 hMass->SetYTitle("Entries");
20 TH1F *hSecVtxZ = new TH1F("hSecVtxZ","D^{0} decay vertex z",1000,-10,10);
21 hSecVtxZ->SetXTitle("z of decay vertex [cm]");
22 hSecVtxZ->SetYTitle("Entries");
23
24 // open input file and get the TTree
25 TFile inFile(aodFileName,"READ");
26 if (!inFile.IsOpen()) return;
27
28 TTree *aodTree = (TTree*)inFile.Get("aodTree");
29 aodTree->AddFriend("aodTree",aodHFFileName);
30
31 AliAODEvent *aod = new AliAODEvent();
32
33 aod->ReadFromTree(aodTree);
34
35 // load heavy flavour vertices
36 TClonesArray *arrayVerticesHF =
37 (TClonesArray*)aod->GetList()->FindObject("VerticesHF");
38
39 // load D0->Kpi candidates
40 TClonesArray *arrayD0toKpi =
41 (TClonesArray*)aod->GetList()->FindObject("D0toKpi");
42
43 // load 3prong candidates
44 TClonesArray *array3Prong =
45 (TClonesArray*)aod->GetList()->FindObject("Charm3Prong");
46
47
48 Double_t cutsD0[9]=
49 // cutsD0[0] = inv. mass half width [GeV]
50 // cutsD0[1] = dca [cm]
51 // cutsD0[2] = cosThetaStar
52 // cutsD0[3] = pTK [GeV/c]
53 // cutsD0[4] = pTPi [GeV/c]
54 // cutsD0[5] = d0K [cm] upper limit!
55 // cutsD0[6] = d0Pi [cm] upper limit!
56 // cutsD0[7] = d0d0 [cm^2]
57 // cutsD0[8] = cosThetaPoint
58 {1000.,
59 100000.,
60 1.1,
61 0.,
62 0.,
63 100000.,
64 100000.,
65 100000000.,
66 -1.1};
67
68 Int_t nTotHF=0,nTotD0toKpi=0,nTot3Prong=0;
69 AliAODVertex *vtx1=0;
70
71 // loop over events
72 Int_t nEvents = aodTree->GetEntries();
73 for (Int_t nEv = 0; nEv < nEvents; nEv++) {
74 cout<<"\n------------ Event: "<<nEv<<" ------------------"<<endl;
75
76 // read event
77 aodTree->GetEvent(nEv);
78 //aodTree->BranchRef();
79
80 //print event info
81 aod->GetHeader()->Print();
82
83 // primary vertex
84 vtx1 = (AliAODVertex*)aod->GetPrimaryVertex();
85 vtx1->Print();
86
87 // make trkIDtoEntry register (temporary)
88 Int_t trkIDtoEntry[100000];
89 for(Int_t it=0;it<aod->GetNumberOfTracks();it++) {
90 AliAODTrack *track = aod->GetTrack(it);
91 trkIDtoEntry[track->GetID()]=it;
92 }
93
94 // loop over D0->Kpi candidates
95 Int_t nD0toKpi = arrayD0toKpi->GetEntriesFast();
96 nTotD0toKpi += nD0toKpi;
97 cout<<"Number of D0->Kpi: "<<nD0toKpi<<endl;
98
99 for (Int_t iD0toKpi = 0; iD0toKpi < nD0toKpi; iD0toKpi++) {
100 AliAODRecoDecayHF2Prong *d = (AliAODRecoDecayHF2Prong*)arrayD0toKpi->UncheckedAt(iD0toKpi);
101 Bool_t unsetvtx=kFALSE;
102 if(!d->GetOwnPrimaryVtx()) {
103 d->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
104 unsetvtx=kTRUE;
105 }
106 Int_t okD0=0,okD0bar=0;
107 if(d->SelectD0(cutsD0,okD0,okD0bar)) {
108 //cout<<1e8*d->Prodd0d0()<<endl;
109 hMass->Fill(d->InvMassD0(),0.5);
110 hMass->Fill(d->InvMassD0bar(),0.5);
111 hCPtaVSd0d0->Fill(1e8*d->Prodd0d0(),d->CosPointingAngle());
112 hSecVtxZ->Fill(d->GetSecVtxZ());
113 //cout<<d->GetSecVtxZ() <<endl;
114
115 // get daughter AOD tracks
116 AliAODTrack *trk0=aod->GetTrack(trkIDtoEntry[d->GetProngID(0)]);
117 AliAODTrack *trk1=aod->GetTrack(trkIDtoEntry[d->GetProngID(1)]);
118 cout<<"pt of positive track: "<<trk0->Pt()<<endl;
119 // this will be replaced by
120 //AliAODTrack *trk0 = (AliAODTrack*)(d->GetSecondaryVtx()->GetDaughter(0));
121 }
122 if(unsetvtx) d->UnsetOwnPrimaryVtx();
123 }
124 // Instanciates vHF and loads its parameters
125
126 // count 3prong candidates
127 Int_t n3Prong = array3Prong->GetEntriesFast();
128 nTot3Prong += n3Prong;
129 cout<<"Number of Charm->3Prong: "<<n3Prong<<endl;
130
131 // loop over HF vertices
132 Int_t nVtxsHF = arrayVerticesHF->GetEntriesFast();
133 nTotHF += nVtxsHF;
134 cout<<"Number of heavy-flavour vertices: "<<nVtxsHF<<endl;
135 for (Int_t iVtx = 0; iVtx < nVtxsHF; iVtx++) {
136 AliAODVertex *vtxHF = (AliAODVertex*)arrayVerticesHF->UncheckedAt(iVtx);
137 // print info
138 //cout << iVtx << ": vertex z position: " << vtxHF->GetZ() << endl;
139 }
140
141 }
142
143 printf("\n Total HF vertices: %d\n",nTotHF);
144 printf("\n Total D0->Kpi: %d\n",nTotD0toKpi);
145 printf("\n Total Charm->3Prong: %d\n",nTot3Prong);
146
147 TCanvas *c = new TCanvas("c","c",0,0,1000,1000);
148 c->Divide(2,2);
149 c->cd(1);
150 hCPtaVSd0d0->Draw("colz");
151 c->cd(2);
152 hMass->SetFillColor(4);
153 hMass->Draw();
154 c->cd(3);
155 hSecVtxZ->SetFillColor(2);
156 hSecVtxZ->Draw();
157
158 return;
159}