Fork me on GitHub

Ignore:
File:
1 edited

Legend:

Unmodified
Added
Removed
  • examples/ExamplePVtxFit.C

    r127644a r0d65492  
    11/*
    22Example of using vertex fitter class to fit primary vertex
    3 assumed to be generated in (0,0,0)
     3assumed to be generated in (0,0,0).
     4To run compiled version:
     5root
     6root> .L examples/LoadPvtxFit.C+
     7root> LoadPVtxFit()
     8root> ExamplePVtxFit("infile.root", Nevents)
     9
     10March 4, 2021
     11F. Bedeschi, INFN-Pisa, Italy
    412*/
    5 
    6 #ifdef __CLING__
    7 R__LOAD_LIBRARY(libDelphes)
    813
    914#include "classes/DelphesClasses.h"
     
    1318#include "external/TrackCovariance/VertexFit.h"
    1419
    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>
    1727
    1828//------------------------------------------------------------------------------
     
    4656                treeReader->ReadEntry(entry);
    4757                Int_t Ntr = 0;  // # of tracks from primary vertex
    48                 Int_t NtrG = branchTrack->GetEntries();
     58                Int_t NtrG = branchTrack->GetEntries(); // Total # of tracks
    4959                TVectorD** pr = new TVectorD * [NtrG];
    5060                TMatrixDSym** cv = new TMatrixDSym * [NtrG];
     
    5868                                Track* trk = (Track*)branchTrack->At(it);
    5969                                //
    60                                 // Get associated generated particle
     70                                // Get associated generated particle 
    6171                                GenParticle* gp = (GenParticle*)trk->Particle.GetObject();
    6272                                //
    63                                 // Position of origin in meters
    64                                 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;
    6777                                //
    68                                 // group tracks originating from the primary vertex
     78                                // Select group of tracks originating from the primary vertex
    6979                                if (x == 0.0 && y == 0.0)
    7080                                {
     
    7888                                        Double_t oPar[5] = { obsD0, obsPhi, obsC, obsZ0, obsCtg };
    7989                                        TVectorD obsPar(5, oPar);       // Fill observed parameters
    80                                         TVector3 xv(x, y, z);
    8190                                        //
    8291                                        pr[Ntr] = new TVectorD(obsPar);
     
    8594                                }
    8695                        }               // End loop on tracks
    87                         //std::cout << "Total of " << Ntr << " primary tracks out of " << NtrG << " tracks" << std::endl;
    8896                }
    8997                //
     
    107115                }
    108116
    109                 //std::cout << "Vertex chi2/Ndof = " << Chi2 / Ndof << std::endl;
    110117                //
    111118                // Cleanup
     
    118125        // Show resulting histograms
    119126        //
    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);
    121128        Cnv->Divide(2, 2);
    122         Cnv->cd(1); gPad->SetLogy(1); gStyle->SetOptFit(1111);
     129        Cnv->cd(1); gPad->SetLogy(1); gStyle->SetOptFit(1111); 
    123130        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); 
    125132        hYpull->Fit("gaus"); hYpull->Draw();
    126133        Cnv->cd(3); gPad->SetLogy(1); gStyle->SetOptFit(1111);
Note: See TracChangeset for help on using the changeset viewer.