Changes in examples/ExamplePVtxFit.C [127644a:0d65492] in git
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
examples/ExamplePVtxFit.C
r127644a r0d65492 1 1 /* 2 2 Example of using vertex fitter class to fit primary vertex 3 assumed to be generated in (0,0,0) 3 assumed to be generated in (0,0,0). 4 To run compiled version: 5 root 6 root> .L examples/LoadPvtxFit.C+ 7 root> LoadPVtxFit() 8 root> ExamplePVtxFit("infile.root", Nevents) 9 10 March 4, 2021 11 F. Bedeschi, INFN-Pisa, Italy 4 12 */ 5 6 #ifdef __CLING__7 R__LOAD_LIBRARY(libDelphes)8 13 9 14 #include "classes/DelphesClasses.h" … … 13 18 #include "external/TrackCovariance/VertexFit.h" 14 19 15 #endif 16 20 #include <TClonesArray.h> 21 #include <TChain.h> 22 #include <TVectorD.h> 23 #include <TMatrixDSym.h> 24 #include <TH1.h> 25 #include <TCanvas.h> 26 #include <TStyle.h> 17 27 18 28 //------------------------------------------------------------------------------ … … 46 56 treeReader->ReadEntry(entry); 47 57 Int_t Ntr = 0; // # of tracks from primary vertex 48 Int_t NtrG = branchTrack->GetEntries(); 58 Int_t NtrG = branchTrack->GetEntries(); // Total # of tracks 49 59 TVectorD** pr = new TVectorD * [NtrG]; 50 60 TMatrixDSym** cv = new TMatrixDSym * [NtrG]; … … 58 68 Track* trk = (Track*)branchTrack->At(it); 59 69 // 60 // Get associated generated particle 70 // Get associated generated particle 61 71 GenParticle* gp = (GenParticle*)trk->Particle.GetObject(); 62 72 // 63 // Position of origin in meters64 Double_t x = 1.0e-3 *gp->X;65 Double_t y = 1.0e-3 *gp->Y;66 Double_t z = 1.0e-3 *gp->Z;73 // Position of origin (mm) 74 Double_t x = gp->X; 75 Double_t y = gp->Y; 76 Double_t z = gp->Z; 67 77 // 68 // grouptracks originating from the primary vertex78 // Select group of tracks originating from the primary vertex 69 79 if (x == 0.0 && y == 0.0) 70 80 { … … 78 88 Double_t oPar[5] = { obsD0, obsPhi, obsC, obsZ0, obsCtg }; 79 89 TVectorD obsPar(5, oPar); // Fill observed parameters 80 TVector3 xv(x, y, z);81 90 // 82 91 pr[Ntr] = new TVectorD(obsPar); … … 85 94 } 86 95 } // End loop on tracks 87 //std::cout << "Total of " << Ntr << " primary tracks out of " << NtrG << " tracks" << std::endl;88 96 } 89 97 // … … 107 115 } 108 116 109 //std::cout << "Vertex chi2/Ndof = " << Chi2 / Ndof << std::endl;110 117 // 111 118 // Cleanup … … 118 125 // Show resulting histograms 119 126 // 120 TCanvas* Cnv = new TCanvas("Cnv", "Delphes generated track plots", 50, 50, 900, 500);127 TCanvas* Cnv = new TCanvas("Cnv", "Delphes primary vertex fit pulls", 50, 50, 900, 500); 121 128 Cnv->Divide(2, 2); 122 Cnv->cd(1); gPad->SetLogy(1); gStyle->SetOptFit(1111); 129 Cnv->cd(1); gPad->SetLogy(1); gStyle->SetOptFit(1111); 123 130 hXpull->Fit("gaus"); hXpull->Draw(); 124 Cnv->cd(2); gPad->SetLogy(1); gStyle->SetOptFit(1111); 131 Cnv->cd(2); gPad->SetLogy(1); gStyle->SetOptFit(1111); 125 132 hYpull->Fit("gaus"); hYpull->Draw(); 126 133 Cnv->cd(3); gPad->SetLogy(1); gStyle->SetOptFit(1111);
Note:
See TracChangeset
for help on using the changeset viewer.