Fork me on GitHub

source: git/examples/ExamplePVtxFind.C@ 193c267c

Last change on this file since 193c267c was 82db145, checked in by Franco BEDESCHI <bed@…>, 4 years ago

Vertex fit add/remove tracks - TrkUtil has cluster counting info

  • Property mode set to 100644
File size: 6.8 KB
Line 
1/*
2Example of using vertex fitter class to fit primary vertex
3assumed to be generated in (0,0,0)
4*/
5
6#ifdef __CLING__
7R__LOAD_LIBRARY(libDelphes)
8
9#include "classes/DelphesClasses.h"
10#include "external/ExRootAnalysis/ExRootTreeReader.h"
11#include "modules/TrackCovariance.h"
12#include "external/TrackCovariance/TrkUtil.h"
13#include "external/TrackCovariance/VertexFit.h"
14
15#endif
16
17
18//------------------------------------------------------------------------------
19
20void ExamplePVtxFind(const char* inputFile, Int_t Nevent = 5)
21{
22 // Create chain of root trees
23 TChain chain("Delphes");
24 chain.Add(inputFile);
25
26 // Create object of class ExRootTreeReader
27 ExRootTreeReader* treeReader = new ExRootTreeReader(&chain);
28 Long64_t numberOfEntries = treeReader->GetEntries();
29
30 // Get pointers to branches used in this analysis
31 TClonesArray* branchGenPart = treeReader->UseBranch("Particle");
32 TClonesArray* branchTrack = treeReader->UseBranch("Track");
33
34 // Book histograms
35 // Vertex fit pulls
36 Int_t Nbin = 100;
37 TH1D* hXpull = new TH1D("hXpull", "Pull X vertex component", Nbin, -10., 10.);
38 TH1D* hYpull = new TH1D("hYpull", "Pull Y vertex component", Nbin, -10., 10.);
39 TH1D* hZpull = new TH1D("hZpull", "Pull Z vertex component", Nbin, -10., 10.);
40 TH1D* hChi2 = new TH1D("hChi2", "Vertex #chi^{2}/N_{dof}", Nbin, 0., 10.);
41 // Track vertex associations
42 TH1D* hTrPrim = new TH1D("hTrPrim", "Available primary tracks", 41, -0.5, 40.5);
43 TH1D* hTrFound = new TH1D("hTrFound", "Found primary tracks", 41, -0.5, 40.5);
44 TH1D* hTrFoundG = new TH1D("hTrFoundG", "Found GOOD primary tracks", 41, -0.5, 40.5);
45 TH1D* hTrFoundB = new TH1D("hTrFoundB", "Found BAD primary tracks", 41, -0.5, 40.5);
46 //
47 // Loop over all events
48 Int_t Nev = TMath::Min(Nevent, (Int_t)numberOfEntries);
49 for (Int_t entry = 0; entry < Nev; ++entry)
50 {
51 // Load selected branches with data from specified event
52 treeReader->ReadEntry(entry);
53 Int_t Ntr = 0; // # of starting tracks used in primary vertex
54 Int_t NtrG = branchTrack->GetEntries();
55 //std::cout << "Event opened containing " << NtrG << " tracks" << std::endl;
56 TVectorD** pr = new TVectorD * [NtrG];
57 TMatrixDSym** cv = new TMatrixDSym * [NtrG];
58 //
59 // test Particle branch
60 Int_t Ngen = branchGenPart->GetEntries();
61 //std::cout << "Nr. of generated particles: " << Ngen << std::endl;
62 // If event contains at least 1 track
63 //
64 Double_t Nprim = 0.0;
65 if (branchTrack->GetEntries() > 0)
66 {
67 // Loop on tracks
68 for (Int_t it = 0; it < branchTrack->GetEntries(); it++)
69 {
70 Track* trk = (Track*)branchTrack->At(it);
71 //
72 // Start fitting all available tracks
73 //
74 // Reconstructed track parameters
75 Double_t obsD0 = trk->D0;
76 Double_t obsPhi = trk->Phi;
77 Double_t obsC = trk->C;
78 Double_t obsZ0 = trk->DZ;
79 Double_t obsCtg = trk->CtgTheta;
80 //std::cout << "Got track parameters for track " << it << std::endl;
81 //
82 // Load tracks for vertex fit if impact parameters is not ridiculous
83 Double_t Dmax = 1.0; // max is 1 mm
84 if (TMath::Abs(obsD0) < Dmax) {
85 Double_t oPar[5] = { obsD0, obsPhi, obsC, obsZ0, obsCtg };
86 TVectorD obsPar(5, oPar); // Fill observed parameters
87 pr[Ntr] = new TVectorD(obsPar);
88 cv[Ntr] = new TMatrixDSym(trk->CovarianceMatrix());
89 Ntr++;
90 //std::cout << "Track loaded Ntr= " << Ntr << std::endl;
91 }
92 //
93 // Get associated generated particle
94 GenParticle* gp = (GenParticle*)trk->Particle.GetObject();
95 //std::cout << "GenParticle pointer "<<gp << std::endl;
96 //
97 // Position of origin in mm
98 Double_t x = gp->X;
99 Double_t y = gp->Y;
100 Double_t z = gp->Z;
101 //std::cout << "Got position of origin " << std::endl;
102 //
103 // Count tracks originating from the primary vertex
104 if (x == 0.0 && y == 0.0) Nprim++;
105
106 } // End loop on tracks
107 }
108 //
109 // Fit primary vertex
110 //
111 Int_t Nfound = Ntr;
112 //std::cout << "Found tracks "<<Nfound << std::endl;
113 Int_t MinTrk = 2; // Minumum # tracks for vertex fit
114 Double_t MaxChi2 = 6.;
115 if (Ntr >= MinTrk) {
116 VertexFit* Vtx = new VertexFit(Ntr, pr, cv);
117 //std::cout << "Vertex fit created " << std::endl;
118 //
119 // Remove tracks with large chi2
120 Bool_t Done = kFALSE;
121 while (!Done){
122 //std::cout << "After while " << std::endl;
123 // Find largest Chi2 contribution
124 TVectorD Chi2List = Vtx->GetVtxChi2List(); // Get contributions to Chi2
125 //std::cout << "After Chi2List. " << std::endl; Chi2List.Print();
126 Double_t* Chi2L = new Double_t[Nfound];
127 Chi2L = Chi2List.GetMatrixArray();
128 Int_t iMax = TMath::LocMax(Nfound, Chi2L);
129 //std::cout << "iMax = "<<iMax << std::endl;
130 Double_t Chi2Mx = Chi2L[iMax];
131 //std::cout << "Chi2Mx "<<Chi2Mx << std::endl;
132 if (Chi2Mx > MaxChi2 && Nfound > 2) {
133 Vtx->RemoveTrk(iMax);
134 Nfound--;
135 }
136 else {
137 Done = kTRUE;
138 }
139 //std::cout << "Done = " << Done << std::endl;
140 //delete[] Chi2L;
141 //std::cout << "Array Chi2L removed " << std::endl;
142 }
143 //
144 //std::cout << "Before getting vertex " << std::endl;
145 //
146 // Require minimum number of tracks in vertex
147 Int_t Nmin = 4;
148 if (Nfound >= Nmin) {
149 TVectorD xvtx = Vtx->GetVtx();
150 //std::cout << "Found vertex " << xvtx(0)<<", "<<xvtx(1)<<", "<<xvtx(2) << std::endl;
151 TMatrixDSym covX = Vtx->GetVtxCov();
152 Double_t Chi2 = Vtx->GetVtxChi2();
153 Double_t Ndof = 2 * (Double_t)Nfound - 3;
154 Double_t PullX = xvtx(0) / TMath::Sqrt(covX(0, 0));
155 Double_t PullY = xvtx(1) / TMath::Sqrt(covX(1, 1));
156 Double_t PullZ = xvtx(2) / TMath::Sqrt(covX(2, 2));
157 //
158 // Fill histograms
159 hXpull->Fill(PullX);
160 hYpull->Fill(PullY);
161 hZpull->Fill(PullZ);
162 hChi2->Fill(Chi2 / Ndof);
163 //
164 hTrPrim->Fill(Nprim);
165 hTrFound->Fill((Double_t)Nfound);
166 //std::cout << "Histograms filled " << std::endl;
167 }
168 //
169 delete Vtx;
170 }
171
172 //std::cout << "Vertex chi2/Ndof = " << Chi2 / Ndof << std::endl;
173 //
174 // Cleanup
175 for (Int_t i = 0; i < Ntr; i++) delete pr[i];
176 for (Int_t i = 0; i < Ntr; i++) delete cv[i];
177 delete[] pr;
178 delete[] cv;
179 }
180 //
181 // Show resulting histograms
182 //
183 TCanvas* Cnv = new TCanvas("Cnv", "Delphes primary vertex pulls", 50, 50, 900, 500);
184 Cnv->Divide(2, 2);
185 Cnv->cd(1); gPad->SetLogy(1); gStyle->SetOptFit(1111);
186 hXpull->Fit("gaus"); hXpull->Draw();
187 Cnv->cd(2); gPad->SetLogy(1); gStyle->SetOptFit(1111);
188 hYpull->Fit("gaus"); hYpull->Draw();
189 Cnv->cd(3); gPad->SetLogy(1); gStyle->SetOptFit(1111);
190 hZpull->Fit("gaus"); hZpull->Draw();
191 Cnv->cd(4); hChi2->Draw();
192 //
193 TCanvas* CnvN = new TCanvas("CnvN", "Primary tracks found", 100, 100, 900, 500);
194 CnvN->Divide(2, 1);
195 CnvN->cd(1);
196 hTrPrim->Draw();
197 CnvN-> cd(2);
198 hTrFound->SetLineColor(kRed);
199 hTrFound->Draw();
200}
Note: See TracBrowser for help on using the repository browser.