Fork me on GitHub

source: git/examples/ExamplePVtxFit.C@ c41f48e

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

Vertex fit example

  • Property mode set to 100644
File size: 4.2 KB
Line 
1/*
2Example of using vertex fitter class to fit primary vertex
3assumed to be generated in (0,0,0)
4*/
5
6#include "classes/DelphesClasses.h"
7#include "external/ExRootAnalysis/ExRootTreeReader.h"
8#include "modules/TrackCovariance.h"
9#include "external/TrackCovariance/TrkUtil.h"
10#include "external/TrackCovariance/VertexFit.h"
11
12#include <TClonesArray.h>
13#include <TChain.h>
14#include <TVector3.h>
15#include <TVectorD.h>
16#include <TMatrixDSym.h>
17#include <TH1.h>
18#include <TCanvas.h>
19#include <TStyle.h>
20
21//------------------------------------------------------------------------------
22
23void ExamplePVtxFit(const char* inputFile, Int_t Nevent = 5)
24{
25 // Create chain of root trees
26 TChain chain("Delphes");
27 chain.Add(inputFile);
28
29 // Create object of class ExRootTreeReader
30 ExRootTreeReader* treeReader = new ExRootTreeReader(&chain);
31 Long64_t numberOfEntries = treeReader->GetEntries();
32
33 // Get pointers to branches used in this analysis
34 TClonesArray* branchGenPart = treeReader->UseBranch("Particle");
35 TClonesArray* branchTrack = treeReader->UseBranch("Track");
36
37 // Book histograms
38 Int_t Nbin = 100;
39 TH1D* hXpull = new TH1D("hXpull", "Pull X vertex component", Nbin, -10., 10.);
40 TH1D* hYpull = new TH1D("hYpull", "Pull Y vertex component", Nbin, -10., 10.);
41 TH1D* hZpull = new TH1D("hZpull", "Pull Z vertex component", Nbin, -10., 10.);
42 TH1D* hChi2 = new TH1D("hChi2", "Vertex #chi^{2}/N_{dof}", Nbin, 0., 10.);
43 //
44 // Loop over all events
45 Int_t Nev = TMath::Min(Nevent, (Int_t) numberOfEntries);
46 for (Int_t entry = 0; entry < Nev; ++entry)
47 {
48 // Load selected branches with data from specified event
49 treeReader->ReadEntry(entry);
50 Int_t Ntr = 0; // # of tracks from primary vertex
51 Int_t NtrG = branchTrack->GetEntries();
52 TVectorD** pr = new TVectorD * [NtrG];
53 TMatrixDSym** cv = new TMatrixDSym * [NtrG];
54 // If event contains at least 1 track
55 //
56 if (branchTrack->GetEntries() > 0)
57 {
58 // Loop on tracks
59 for (Int_t it = 0; it < branchTrack->GetEntries(); it++)
60 {
61 Track* trk = (Track*)branchTrack->At(it);
62 //
63 // Get associated generated particle
64 GenParticle* gp = (GenParticle*)trk->Particle.GetObject();
65 //
66 // Position of origin in meters
67 Double_t x = 1.0e-3 * gp->X;
68 Double_t y = 1.0e-3 * gp->Y;
69 Double_t z = 1.0e-3 * gp->Z;
70 //
71 // group tracks originating from the primary vertex
72 if (x == 0.0 && y == 0.0)
73 {
74 //
75 // Reconstructed track parameters
76 Double_t obsD0 = trk->D0;
77 Double_t obsPhi = trk->Phi;
78 Double_t obsC = trk->C;
79 Double_t obsZ0 = trk->DZ;
80 Double_t obsCtg = trk->CtgTheta;
81 Double_t oPar[5] = { obsD0, obsPhi, obsC, obsZ0, obsCtg };
82 TVectorD obsPar(5, oPar); // Fill observed parameters
83 TVector3 xv(x, y, z);
84 //
85 pr[Ntr] = new TVectorD(obsPar);
86 cv[Ntr] = new TMatrixDSym(TrkUtil::CovToMm(trk->CovarianceMatrix()));
87 Ntr++;
88 }
89 } // End loop on tracks
90 //std::cout << "Total of " << Ntr << " primary tracks out of " << NtrG << " tracks" << std::endl;
91 }
92 //
93 // Fit primary vertex
94 Int_t MinTrk = 2; // Minumum # tracks for vertex fit
95 if (Ntr >= MinTrk) {
96 VertexFit* Vtx = new VertexFit(Ntr, pr, cv);
97 TVectorD xvtx = Vtx->GetVtx();
98 TMatrixDSym covX = Vtx->GetVtxCov();
99 Double_t Chi2 = Vtx->GetVtxChi2();
100 Double_t Ndof = 2 * (Double_t)Ntr - 3;
101 Double_t PullX = xvtx(0) / TMath::Sqrt(covX(0, 0));
102 Double_t PullY = xvtx(1) / TMath::Sqrt(covX(1, 1));
103 Double_t PullZ = xvtx(2) / TMath::Sqrt(covX(2, 2));
104 //
105 // Fill histograms
106 hXpull->Fill(PullX);
107 hYpull->Fill(PullY);
108 hZpull->Fill(PullZ);
109 hChi2->Fill(Chi2 / Ndof);
110 }
111
112 //std::cout << "Vertex chi2/Ndof = " << Chi2 / Ndof << std::endl;
113 //
114 // Cleanup
115 for (Int_t i = 0; i < Ntr; i++) delete pr[i];
116 for (Int_t i = 0; i < Ntr; i++) delete cv[i];
117 delete[] pr;
118 delete[] cv;
119 }
120 //
121 // Show resulting histograms
122 //
123 TCanvas* Cnv = new TCanvas("Cnv", "Delphes generated track plots", 50, 50, 900, 500);
124 Cnv->Divide(2, 2);
125 Cnv->cd(1); gPad->SetLogy(1); gStyle->SetOptFit(1111);
126 hXpull->Fit("gaus"); hXpull->Draw();
127 Cnv->cd(2); gPad->SetLogy(1); gStyle->SetOptFit(1111);
128 hYpull->Fit("gaus"); hYpull->Draw();
129 Cnv->cd(3); gPad->SetLogy(1); gStyle->SetOptFit(1111);
130 hZpull->Fit("gaus"); hZpull->Draw();
131 Cnv->cd(4); hChi2->Draw();
132}
Note: See TracBrowser for help on using the repository browser.