Fork me on GitHub

source: git/examples/ExamplePVtxFind.C@ 46b3e01

Last change on this file since 46b3e01 was b750b0a, checked in by Franco BEDESCHI <bed@…>, 3 years ago

Fix CovarianceMatrix scaling error. Vertexing improvements.

  • Property mode set to 100644
File size: 9.5 KB
Line 
1/*
2Example of using vertex fitter class to fit primary vertex
3assumed to be generated with Pythia8/ee_zh_smr-shf.cmn
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 TH1D* hTrDiff = new TH1D("hTrDiff", "Found - available tracks", 21, -10.5, 10.5);
47 //
48 // Loop over all events
49 Int_t Nev = TMath::Min(Nevent, (Int_t)numberOfEntries);
50 for (Int_t entry = 0; entry < Nev; ++entry)
51 {
52 // Load selected branches with data from specified event
53 treeReader->ReadEntry(entry);
54 Int_t NtrG = branchTrack->GetEntries();
55 if(entry%500 ==0)std::cout << "Event "<<entry<<" 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 Double_t xpv = 0.0; // Init true primary vertex position
66 Double_t ypv = 0.0;
67 Double_t zpv = 0.0;
68 if (branchTrack->GetEntries() > 0)
69 {
70 // Loop on tracks
71 for (Int_t it = 0; it < branchTrack->GetEntries(); it++)
72 {
73 Track* trk = (Track*)branchTrack->At(it);
74 //
75 // Start fitting all available tracks
76 //
77 // Reconstructed track parameters
78 Double_t obsD0 = trk->D0;
79 Double_t obsPhi = trk->Phi;
80 Double_t obsC = trk->C;
81 Double_t obsZ0 = trk->DZ;
82 Double_t obsCtg = trk->CtgTheta;
83 //std::cout << "Got track parameters for track " << it << std::endl;
84 //
85 // Load all tracks for vertex fit
86 Double_t oPar[5] = { obsD0, obsPhi, obsC, obsZ0, obsCtg };
87 TVectorD obsPar(5, oPar); // Fill observed parameters
88 pr[it] = new TVectorD(obsPar);
89 cv[it] = new TMatrixDSym(trk->CovarianceMatrix());
90 // Check covariance
91 /*
92 TMatrixDSymEigen Eign(*cv[it]);
93 TVectorD lambda = Eign.GetEigenValues();
94 Bool_t good = kTRUE;
95 for(Int_t i=0; i<5; i++) if(lambda(i)<0.0) good = kFALSE;
96 if(!good){
97 std::cout<<"Cov["<<it<<"]"<<" eigenvalues: "
98 <<lambda(0)<<", "<<lambda(1)<<", "<<lambda(2)<<", "
99 <<lambda(3)<<", "<<lambda(4)<<std::endl;
100 TMatrixDSym Ctmp = *cv[it];
101 std::cout<<"Error D (dir - cov) "<<trk->ErrorD0
102 <<" - "<<TMath::Sqrt(Ctmp(0,0))<<std::endl;
103 }
104 */
105 //std::cout << "Track loaded it= " << it << std::endl;
106 //
107 // Find true primary vertex
108 GenParticle* gp = (GenParticle*)trk->Particle.GetObject();
109 //std::cout << "GenParticle pointer "<<gp << std::endl;
110 //
111 // Position of origin in mm
112 Double_t x = gp->X;
113 Double_t y = gp->Y;
114 Double_t z = gp->Z;
115 Bool_t prim = kTRUE; // Is primary?
116 Int_t mp = gp->M1; // Mother
117 while (mp > 0) {
118 GenParticle* gm =
119 (GenParticle*)branchGenPart->At(mp);
120 Double_t xm = gm->X;
121 Double_t ym = gm->Y;
122 Double_t zm = gm->Z;
123 if (x != xm || y != ym || z != zm) {
124 prim = kFALSE;
125 break;
126 }
127 else mp = gm->M1;
128 }
129 if (prim) { // It's a primary track
130 Nprim++;
131 xpv = x; // Store true primary
132 ypv = y;
133 zpv = z;
134 }
135
136 } // End loop on tracks
137 }
138 std::cout<<"PVtxFind true vertex: Nprim= "<<Nprim<<", x,y,z= "<<xpv<<", "<<ypv<<", "<<zpv<<std::endl;
139 //
140 // Find primary vertex
141 //
142 //Beam constraint
143 TVectorD xpvc(3);
144 xpvc(0) = 1.0;
145 xpvc(1) = -2.0;
146 xpvc(2) = 10.0;
147 TMatrixDSym covpvc(3); covpvc.Zero();
148 covpvc(0, 0) = 0.0097 * 0.0097;
149 covpvc(1, 1) = 2.55e-05 * 2.55e-05;
150 covpvc(2, 2) = 0.64 * 0.64;
151 if(Nprim == 0){
152 xpv = xpvc(0);
153 ypv = xpvc(1);
154 zpv = xpvc(2);
155 }
156 //
157 //
158 // Skim tracks
159 Int_t nSkim = 0;
160 Int_t* nSkimmed = new Int_t[NtrG];
161 TVectorD** PrSk = new TVectorD * [1];
162 TMatrixDSym** CvSk = new TMatrixDSym * [1];
163 Double_t MaxChi2 = 2.0+3*2.;
164 for (Int_t n = 0; n < NtrG; n++) {
165 PrSk[0] = new TVectorD(*pr[n]);
166 CvSk[0] = new TMatrixDSym(*cv[n]);
167 VertexFit* Vskim = new VertexFit(1,PrSk, CvSk);
168 Vskim->AddVtxConstraint(xpvc, covpvc);
169 Double_t Chi2One = Vskim->GetVtxChi2();
170 //std::cout<<"Track "<<n<<", Chi2 = "<<Chi2One<<std::endl;
171 if (Chi2One < MaxChi2) {
172 nSkimmed[nSkim] = n;
173 //std::cout << "nSkimmed[" << nSkim << "] = " << n << std::endl;
174 nSkim++;
175 }
176 // Cleanup
177 delete Vskim;
178 }
179 delete PrSk[0];
180 delete CvSk[0];
181 delete[] PrSk;
182 delete[] CvSk;
183 //
184 // Load tracks for primary fit
185 Int_t MinTrk = 1; // Minumum # tracks for vertex fit
186 if (nSkim >= MinTrk) {
187 TVectorD** PrFit = new TVectorD * [nSkim];
188 TMatrixDSym** CvFit = new TMatrixDSym * [nSkim];
189 for (Int_t n = 0; n < nSkim; n++) {
190 PrFit[n] = new TVectorD(*pr[nSkimmed[n]]);
191 CvFit[n] = new TMatrixDSym(*cv[nSkimmed[n]]);
192 }
193 delete[] nSkimmed;
194 Int_t Nfound = nSkim;
195 if(entry%500 ==0)std::cout << "Found tracks "<<Nfound << std::endl;
196 const Int_t MaxFound = 100; Double_t Chi2LL[MaxFound]; Double_t *Chi2L = Chi2LL;
197 VertexFit* Vtx = new VertexFit(nSkim, PrFit, CvFit);
198 //std::cout << "Vertex fit created " << std::endl;
199 Vtx->AddVtxConstraint(xpvc, covpvc);
200 //
201 // Remove tracks with large chi2
202 Double_t MaxChi2Fit = 4.5;
203 Bool_t Done = kFALSE;
204 while (!Done) {
205 //std::cout << "After while " << std::endl;
206 // Find largest Chi2 contribution
207 TVectorD Chi2List = Vtx->GetVtxChi2List(); // Get contributions to Chi2
208 //std::cout << "After Chi2List. " << std::endl; Chi2List.Print();
209 //Double_t* Chi2L = new Double_t[Nfound];
210 Chi2L = Chi2List.GetMatrixArray();
211 Int_t iMax = TMath::LocMax(Nfound, Chi2L);
212 //std::cout << "iMax = "<<iMax << std::endl;
213 Double_t Chi2Mx = Chi2L[iMax];
214 //std::cout << "Chi2Mx "<<Chi2Mx << std::endl;
215 if (Chi2Mx > MaxChi2Fit && Nfound > 1) {
216 //std::cout << "Before remove. Nfound = "<<Nfound << std::endl;
217 Vtx->RemoveTrk(iMax);
218 //std::cout << "After remove." << std::endl;
219 Nfound--;
220 }
221 else {
222 Done = kTRUE;
223 }
224 }
225 //
226 //std::cout << "Before getting vertex " << std::endl;
227 //
228 // Require minimum number of tracks in vertex
229 Int_t Nmin = 1;
230 if (Nfound >= Nmin) {
231 TVectorD xvtx = Vtx->GetVtx();
232 std::cout << "Found vertex " << xvtx(0)<<", "<<xvtx(1)<<", "<<xvtx(2) << std::endl;
233 TMatrixDSym covX = Vtx->GetVtxCov();
234 Double_t Chi2 = Vtx->GetVtxChi2();
235 Double_t Ndof = 2 * (Double_t)Nfound;
236 Double_t PullX = (xvtx(0) - xpv) / TMath::Sqrt(covX(0, 0));
237 Double_t PullY = (xvtx(1) - ypv) / TMath::Sqrt(covX(1, 1));
238 Double_t PullZ = (xvtx(2) - zpv) / TMath::Sqrt(covX(2, 2));
239 //
240 // Fill histograms
241 hXpull->Fill(PullX);
242 hYpull->Fill(PullY);
243 hZpull->Fill(PullZ);
244 hChi2->Fill(Chi2 / Ndof);
245 //
246 hTrPrim->Fill(Nprim);
247 hTrFound->Fill((Double_t)Nfound);
248 hTrDiff->Fill((Double_t)Nfound-Nprim);
249 //std::cout << "Histograms filled " << std::endl;
250 }
251 //
252 // Clean
253 delete Vtx;
254 for (Int_t i = 0; i < nSkim; i++) delete PrFit[i];
255 for (Int_t i = 0; i < nSkim; i++) delete CvFit[i];
256 delete[] PrFit;
257 delete[] CvFit;
258 }
259
260 //std::cout << "Vertex chi2/Ndof = " << Chi2 / Ndof << std::endl;
261 //
262 // Cleanup
263 for (Int_t i = 0; i < NtrG; i++) delete pr[i];
264 for (Int_t i = 0; i < NtrG; i++) delete cv[i];
265 delete[] pr;
266 delete[] cv;
267 }
268 //
269 // Show resulting histograms
270 //
271 TCanvas* Cnv = new TCanvas("Cnv", "Delphes primary vertex pulls", 50, 50, 900, 500);
272 Cnv->Divide(2, 2);
273 Cnv->cd(1); gPad->SetLogy(1); gStyle->SetOptFit(1111);
274 hXpull->Fit("gaus"); hXpull->Draw();
275 Cnv->cd(2); gPad->SetLogy(1); gStyle->SetOptFit(1111);
276 hYpull->Fit("gaus"); hYpull->Draw();
277 Cnv->cd(3); gPad->SetLogy(1); gStyle->SetOptFit(1111);
278 hZpull->Fit("gaus"); hZpull->Draw();
279 Cnv->cd(4); hChi2->Draw();
280 //
281 TCanvas* CnvN = new TCanvas("CnvN", "Primary tracks found", 100, 100, 900, 500);
282 CnvN->Divide(2, 2);
283 CnvN->cd(1);
284 hTrPrim->Draw();
285 CnvN->cd(2);
286 hTrFound->SetLineColor(kRed);
287 hTrFound->Draw();
288 CnvN->cd(3);
289 hTrDiff->Draw();
290}
Note: See TracBrowser for help on using the repository browser.