Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:29:10

0001 #ifndef SimpleVertexAnalysis_h
0002 #define SimpleVertexAnalysis_h
0003 
0004 #include <TROOT.h>
0005 #include <TChain.h>
0006 #include <TFile.h>
0007 #include <TH1F.h>
0008 #include "TString.h"
0009 #include "TCanvas.h"
0010 #include <iostream>
0011 #include <iomanip>
0012 #include <vector>
0013 
0014 #define MAXTRACK 120
0015 
0016 
0017 
0018   /**
0019    * Root analysis code for TTree produced by SimpleVertexTree
0020    */
0021 
0022 class SimpleVertexAnalysis {
0023    public :
0024 
0025    SimpleVertexAnalysis(TTree *tree=0);
0026   /**
0027    * Constructor for a single file<br>
0028    * \param filename: The file name
0029    * \param treeName: The name of the tree
0030    */
0031    SimpleVertexAnalysis(TString filename, 
0032         TString treeName = "VertexFitter");
0033   /**
0034    * Constructor for multiple files, of the type result_1.root<br>
0035    * \param filename:   The base of he file name (e.g. "result")
0036    * \param start:  The number of the first file
0037    * \param end:    The number of the last file
0038    * \param treeName:   The name of the tree
0039    */
0040    SimpleVertexAnalysis(TString base, int start, int end, 
0041         TString treeName = "VertexFitter");
0042    ~SimpleVertexAnalysis();
0043 
0044   /**
0045    * The main method to loop over all the events and fill the histos.
0046    */
0047    void vertexLoop();
0048 
0049   /**
0050    * Sets the limits for the vertex hisotgrams:
0051    *   aMaxTrans:   vertex resolutions, transverse coordinates (x, y)
0052    *   aMaxZ:   vertex resolutions, longitudinal coordinates (z)
0053    *   aMaxTracks:  number of tracks (Sim, Rec)
0054    *   aMaxWeights: sum of weights
0055    *   aMaxTime:    cpu time per fit
0056    */
0057 
0058    void vertexHistoLimits(float aMaxTrans = 250.0, float aMaxZ = 250.0, 
0059      float aMaxTracks = 60.0,  float aMaxWeights = 35.0, float aMaxTime = 0.05);
0060 
0061   /**
0062    * The main method to loop over all the events and fill the histos.
0063    */
0064    void trackLoop();
0065   /**
0066    * Method to print all the Vertx plots to a ps file.
0067    */
0068    void psVertexResult(TString name);
0069   /**
0070    * Method to print all the Vertx plots to a serie of eps files.
0071    */
0072    void epsVertexResult(TString name);
0073   /**
0074    * Output Vertx plots into 2 canvases
0075    */
0076    void plotVertexResult();
0077 
0078   /**
0079    * Fit of the residual and pull plots with a single Gaussian, and printout of the
0080    * main results
0081    * A separate output stream can be provided (parameter out - default is std::cout)
0082    */
0083    void singleGaussianVertexResult(ostream &out = std::cout);
0084   /**
0085    * Fit of the residual and pull plots with two Gaussians (one for the core, the
0086    * other for the tails, and printout of the main results
0087    * A separate output stream can be provided (parameter out - default is std::cout)
0088    */
0089    void doubleGaussianVertexResult(ostream &out = std::cout);
0090   /**
0091    * Measurement of the coverage of the residual distributions 
0092    * (50%, 90% and 95% coverage)
0093    * A separate output stream can be provided (parameter out - default is std::cout)
0094    */
0095    void vertexCoverage(ostream &out = std::cout);
0096   /**
0097    * Output of Track parameter results into 4 canvases
0098    */
0099    void plotTrackResult();
0100   /**
0101    * Method to print all the Track parameter result plots to a ps file.
0102    */
0103    void psTrackResult(TString name);
0104 
0105   /**
0106    * Method to produce the TeX line for the statistices table 
0107    * of the Tables (TeX )
0108    * A separate output stream can be provided (parameter out - default is std::cout)
0109    */
0110    void statTeXResult(ostream &out = std::cout);
0111 
0112   /**
0113    * Method to produce the TeX line for the X-coordinate table 
0114    * of the Tables (TeX )
0115    * A separate output stream can be provided (parameter out - default is std::cout)
0116    */
0117    void xTeXResult(ostream &out = std::cout);
0118 
0119   /**
0120    * Method to produce the TeX line for the Z-coordinate table 
0121    * of the Tables (TeX )
0122    * A separate output stream can be provided (parameter out - default is std::cout)
0123    */
0124    void zTeXResult(ostream &out = std::cout);
0125 
0126   /**
0127    * Method to produce the TeX line for the X- and Z-coordinate table 
0128    * of the Tables (TeX )
0129    * A separate output stream can be provided (parameter out - default is std::cout)
0130    */
0131    void resolutionTeXResult(ostream &out = std::cout);
0132 
0133 
0134    TH1F *resX, *resY, *resZ; 
0135    TH1F *pullX, *pullY, *pullZ;
0136    TH1F *chiNorm, *chiProbability, *weight, *normWeight, *downWeight;
0137    TH1F *numberUsedRecTracks, *numberRawRecTracks, *numberSimTracks, *sharedTracks, *ratioSharedTracks, *timing;
0138    TCanvas *resCanvas, *statCanvas;
0139    float x_coverage, y_coverage, z_coverage;
0140    int failure, total;
0141 
0142    TH1F *resRecPt, *resRecPhi, *resRecTheta, *resRecTimp, *resRecLimp;
0143    TH1F *pullRecPt, *pullRecPhi, *pullRecTheta, *pullRecTimp, *pullRecLimp;
0144    TH1F *resRefPt, *resRefPhi, *resRefTheta, *resRefTimp, *resRefLimp;
0145    TH1F *pullRefPt, *pullRefPhi, *pullRefTheta, *pullRefTimp, *pullRefLimp;
0146    TH1F *pTSim, *etaSim, *pTRec, *etaRec, *pTRef, *etaRef;
0147    TCanvas *resRecCanvas,*pullRecCanvas,*resRefCanvas,*pullRefCanvas, *distCanvas;
0148    bool bookedVertexC, bookedTrackC;
0149 
0150    Int_t  Cut(Int_t entry);
0151    Int_t  GetEntry(Int_t entry);
0152    Int_t  LoadTree(Int_t entry);
0153    void   Init(TTree *tree);
0154    Bool_t Notify();
0155    void   Show(Int_t entry = -1);
0156    TString theTreeName;
0157    Int_t nentries;
0158 
0159 
0160    TTree          *fChain;   //!pointer to the analyzed TTree or TChain
0161    Int_t           fCurrent; //!current Tree number in a TChain
0162    // Declaration of leave types
0163    Int_t           vertex;
0164    Float_t         simPos_X;
0165    Float_t         simPos_Y;
0166    Float_t         simPos_Z;
0167    Float_t         recPos_X;
0168    Float_t         recPos_Y;
0169    Float_t         recPos_Z;
0170    Float_t         recErr_X;
0171    Float_t         recErr_Y;
0172    Float_t         recErr_Z;
0173    Int_t           nbrTrk_Sim;
0174    Int_t           nbrTrk_Rec;
0175    Int_t           nbrTrk_Shared;
0176    Float_t         chiTot;
0177    Float_t         ndf;
0178    Float_t         chiProb;
0179    Float_t         time;
0180    Int_t           simTracks;
0181    Int_t           simTrack_recIndex[MAXTRACK];   //[simTracks]
0182    Float_t         simPar_ptinv[MAXTRACK];   //[simTracks]
0183    Float_t         simPar_theta[MAXTRACK];   //[simTracks]
0184    Float_t         simPar_phi[MAXTRACK];   //[simTracks]
0185    Float_t         simPar_timp[MAXTRACK];   //[simTracks]
0186    Float_t         simPar_limp[MAXTRACK];   //[simTracks]
0187    Int_t           recTracks;
0188    Int_t           recTrack_simIndex[MAXTRACK];   //[recTracks]
0189    Float_t         recTrack_weight[MAXTRACK];   //[recTracks]
0190    Float_t         recPar_ptinv[MAXTRACK];   //[recTracks]
0191    Float_t         recPar_theta[MAXTRACK];   //[recTracks]
0192    Float_t         recPar_phi[MAXTRACK];   //[recTracks]
0193    Float_t         recPar_timp[MAXTRACK];   //[recTracks]
0194    Float_t         recPar_limp[MAXTRACK];   //[recTracks]
0195    Float_t         refPar_ptinv[MAXTRACK];   //[recTracks]
0196    Float_t         refPar_theta[MAXTRACK];   //[recTracks]
0197    Float_t         refPar_phi[MAXTRACK];   //[recTracks]
0198    Float_t         refPar_timp[MAXTRACK];   //[recTracks]
0199    Float_t         refPar_limp[MAXTRACK];   //[recTracks]
0200    Float_t         recErr_ptinv[MAXTRACK];   //[recTracks]
0201    Float_t         recErr_theta[MAXTRACK];   //[recTracks]
0202    Float_t         recErr_phi[MAXTRACK];   //[recTracks]
0203    Float_t         recErr_timp[MAXTRACK];   //[recTracks]
0204    Float_t         recErr_limp[MAXTRACK];   //[recTracks]
0205    Float_t         refErr_ptinv[MAXTRACK];   //[recTracks]
0206    Float_t         refErr_theta[MAXTRACK];   //[recTracks]
0207    Float_t         refErr_phi[MAXTRACK];   //[recTracks]
0208    Float_t         refErr_timp[MAXTRACK];   //[recTracks]
0209    Float_t         refErr_limp[MAXTRACK];   //[recTracks]
0210 
0211    // List of branches
0212    TBranch        *b_vertex;   //!
0213    TBranch        *b_simPos;   //!
0214    TBranch        *b_recPos;   //!
0215    TBranch        *b_recErr;   //!
0216    TBranch        *b_nbrTrk;   //!
0217    TBranch        *b_chiTot;   //!
0218    TBranch        *b_ndf;   //!
0219    TBranch        *b_chiProb;   //!
0220    TBranch        *b_time;   //!
0221    TBranch        *b_simTracks;   //!
0222    TBranch        *b_simTrack_recIndex;   //!
0223    TBranch        *b_simPar_ptinv;   //!
0224    TBranch        *b_simPar_theta;   //!
0225    TBranch        *b_simPar_phi;   //!
0226    TBranch        *b_simPar_timp;   //!
0227    TBranch        *b_simPar_limp;   //!
0228    TBranch        *b_recTracks;   //!
0229    TBranch        *b_recTrack_simIndex;   //!
0230    TBranch        *b_recTrack_weight;   //!
0231    TBranch        *b_recPar_ptinv;   //!
0232    TBranch        *b_recPar_theta;   //!
0233    TBranch        *b_recPar_phi;   //!
0234    TBranch        *b_recPar_timp;   //!
0235    TBranch        *b_recPar_limp;   //!
0236    TBranch        *b_refPar_ptinv;   //!
0237    TBranch        *b_refPar_theta;   //!
0238    TBranch        *b_refPar_phi;   //!
0239    TBranch        *b_refPar_timp;   //!
0240    TBranch        *b_refPar_limp;   //!
0241    TBranch        *b_recErr_ptinv;   //!
0242    TBranch        *b_recErr_theta;   //!
0243    TBranch        *b_recErr_phi;   //!
0244    TBranch        *b_recErr_timp;   //!
0245    TBranch        *b_recErr_limp;   //!
0246    TBranch        *b_refErr_ptinv;   //!
0247    TBranch        *b_refErr_theta;   //!
0248    TBranch        *b_refErr_phi;   //!
0249    TBranch        *b_refErr_timp;   //!
0250    TBranch        *b_refErr_limp;   //!
0251 
0252    TChain* myChain;
0253 
0254 private:
0255   /**
0256    * Calculates and prints the 50, 90 and 95% coverage for the given vector.
0257    * Returns the 90% coverage.
0258    */
0259   float getCoverage(std::vector<float> &residuals, ostream &out = std::cout);
0260   void doubleGaussianFit(TH1F *plot, ostream &out = std::cout);
0261   void epsPlot(TH1F *plot, TString name);
0262 
0263    void bookTrackHisto();
0264    void deleteTrackHisto();
0265    void bookVertexHisto();
0266    void deleteVertexHisto();
0267 
0268   float maxTrans, maxZ, maxTracks, maxWeights, maxTime;
0269   bool bookedVertex, bookedTrack;
0270 };
0271 
0272 #endif
0273 
0274 #ifdef SimpleVertexAnalysis_cxx
0275 SimpleVertexAnalysis::SimpleVertexAnalysis(TTree *tree)
0276 {
0277 // if parameter tree is not specified (or zero), connect the file
0278 // used to generate this class and read the Tree.
0279    bookedVertex=false;bookedTrack=false;
0280    bookedVertexC=false;bookedTrackC=false;
0281    if (tree == 0) {
0282       TFile *f = (TFile*)gROOT->GetListOfFiles()->FindObject("data/pvr_digi_1.root");
0283       if (!f) {
0284          f = new TFile("data/pvr_digi_1.root");
0285       }
0286       tree = (TTree*)gDirectory->Get("closest");
0287 
0288    }
0289    Init(tree);
0290    vertexHistoLimits();
0291 }
0292 
0293 SimpleVertexAnalysis::~SimpleVertexAnalysis()
0294 {
0295    if (!fChain) return;
0296    delete fChain->GetCurrentFile();
0297   if (bookedVertex) deleteVertexHisto();
0298    if (bookedTrack) deleteTrackHisto();
0299 }
0300 
0301 Int_t SimpleVertexAnalysis::GetEntry(Int_t entry)
0302 {
0303 // Read contents of entry.
0304    if (!fChain) return 0;
0305    return fChain->GetEntry(entry);
0306 }
0307 Int_t SimpleVertexAnalysis::LoadTree(Int_t entry)
0308 {
0309 // Set the environment to read one entry
0310    if (!fChain) return -5;
0311    Int_t centry = fChain->LoadTree(entry);
0312    if (centry < 0) return centry;
0313    if (fChain->IsA() != TChain::Class()) return centry;
0314    TChain *chain = (TChain*)fChain;
0315    if (chain->GetTreeNumber() != fCurrent) {
0316       fCurrent = chain->GetTreeNumber();
0317       Notify();
0318    }
0319    return centry;
0320 }
0321 
0322 void SimpleVertexAnalysis::Init(TTree *tree)
0323 {
0324    // The Init() function is called when the selector needs to initialize
0325    // a new tree or chain. Typically here the branch addresses of the tree
0326    // will be set. It is normaly not necessary to make changes to the
0327    // generated code, but the routine can be extended by the user if needed.
0328    // Init() will be called many times when running with PROOF.
0329 
0330    // Set branch addresses
0331    if (tree == 0) return;
0332    fChain = tree;
0333    fCurrent = -1;
0334    fChain->SetMakeClass(1);
0335 
0336    fChain->SetBranchAddress("vertex",&vertex);
0337    fChain->SetBranchAddress("simPos",&simPos_X);
0338    fChain->SetBranchAddress("recPos",&recPos_X);
0339    fChain->SetBranchAddress("recErr",&recErr_X);
0340    fChain->SetBranchAddress("nbrTrk",&nbrTrk_Sim);
0341    fChain->SetBranchAddress("chiTot",&chiTot);
0342    fChain->SetBranchAddress("ndf",&ndf);
0343    fChain->SetBranchAddress("chiProb",&chiProb);
0344    fChain->SetBranchAddress("time",&time);
0345    fChain->SetBranchAddress("simTracks",&simTracks);
0346    fChain->SetBranchAddress("simTrack_recIndex",simTrack_recIndex);
0347    fChain->SetBranchAddress("simPar_ptinv",simPar_ptinv);
0348    fChain->SetBranchAddress("simPar_theta",simPar_theta);
0349    fChain->SetBranchAddress("simPar_phi",simPar_phi);
0350    fChain->SetBranchAddress("simPar_timp",simPar_timp);
0351    fChain->SetBranchAddress("simPar_limp",simPar_limp);
0352    fChain->SetBranchAddress("recTracks",&recTracks);
0353    fChain->SetBranchAddress("recTrack_simIndex",recTrack_simIndex);
0354    fChain->SetBranchAddress("recTrack_weight",recTrack_weight);
0355    fChain->SetBranchAddress("recPar_ptinv",recPar_ptinv);
0356    fChain->SetBranchAddress("recPar_theta",recPar_theta);
0357    fChain->SetBranchAddress("recPar_phi",recPar_phi);
0358    fChain->SetBranchAddress("recPar_timp",recPar_timp);
0359    fChain->SetBranchAddress("recPar_limp",recPar_limp);
0360    fChain->SetBranchAddress("refPar_ptinv",refPar_ptinv);
0361    fChain->SetBranchAddress("refPar_theta",refPar_theta);
0362    fChain->SetBranchAddress("refPar_phi",refPar_phi);
0363    fChain->SetBranchAddress("refPar_timp",refPar_timp);
0364    fChain->SetBranchAddress("refPar_limp",refPar_limp);
0365    fChain->SetBranchAddress("recErr_ptinv",recErr_ptinv);
0366    fChain->SetBranchAddress("recErr_theta",recErr_theta);
0367    fChain->SetBranchAddress("recErr_phi",recErr_phi);
0368    fChain->SetBranchAddress("recErr_timp",recErr_timp);
0369    fChain->SetBranchAddress("recErr_limp",recErr_limp);
0370    fChain->SetBranchAddress("refErr_ptinv",refErr_ptinv);
0371    fChain->SetBranchAddress("refErr_theta",refErr_theta);
0372    fChain->SetBranchAddress("refErr_phi",refErr_phi);
0373    fChain->SetBranchAddress("refErr_timp",refErr_timp);
0374    fChain->SetBranchAddress("refErr_limp",refErr_limp);
0375    Notify();
0376 }
0377 
0378 Bool_t SimpleVertexAnalysis::Notify()
0379 {
0380    // The Notify() function is called when a new file is opened. This
0381    // can be either for a new TTree in a TChain or when when a new TTree
0382    // is started when using PROOF. Typically here the branch pointers
0383    // will be retrieved. It is normaly not necessary to make changes
0384    // to the generated code, but the routine can be extended by the
0385    // user if needed.
0386 
0387    // Get branch pointers
0388    b_vertex = fChain->GetBranch("vertex");
0389    b_simPos = fChain->GetBranch("simPos");
0390    b_recPos = fChain->GetBranch("recPos");
0391    b_recErr = fChain->GetBranch("recErr");
0392    b_nbrTrk = fChain->GetBranch("nbrTrk");
0393    b_chiTot = fChain->GetBranch("chiTot");
0394    b_ndf = fChain->GetBranch("ndf");
0395    b_chiProb = fChain->GetBranch("chiProb");
0396    b_time = fChain->GetBranch("time");
0397    b_simTracks = fChain->GetBranch("simTracks");
0398    b_simTrack_recIndex = fChain->GetBranch("simTrack_recIndex");
0399    b_simPar_ptinv = fChain->GetBranch("simPar_ptinv");
0400    b_simPar_theta = fChain->GetBranch("simPar_theta");
0401    b_simPar_phi = fChain->GetBranch("simPar_phi");
0402    b_simPar_timp = fChain->GetBranch("simPar_timp");
0403    b_simPar_limp = fChain->GetBranch("simPar_limp");
0404    b_recTracks = fChain->GetBranch("recTracks");
0405    b_recTrack_simIndex = fChain->GetBranch("recTrack_simIndex");
0406    b_recTrack_weight = fChain->GetBranch("recTrack_weight");
0407    b_recPar_ptinv = fChain->GetBranch("recPar_ptinv");
0408    b_recPar_theta = fChain->GetBranch("recPar_theta");
0409    b_recPar_phi = fChain->GetBranch("recPar_phi");
0410    b_recPar_timp = fChain->GetBranch("recPar_timp");
0411    b_recPar_limp = fChain->GetBranch("recPar_limp");
0412    b_refPar_ptinv = fChain->GetBranch("refPar_ptinv");
0413    b_refPar_theta = fChain->GetBranch("refPar_theta");
0414    b_refPar_phi = fChain->GetBranch("refPar_phi");
0415    b_refPar_timp = fChain->GetBranch("refPar_timp");
0416    b_refPar_limp = fChain->GetBranch("refPar_limp");
0417    b_recErr_ptinv = fChain->GetBranch("recErr_ptinv");
0418    b_recErr_theta = fChain->GetBranch("recErr_theta");
0419    b_recErr_phi = fChain->GetBranch("recErr_phi");
0420    b_recErr_timp = fChain->GetBranch("recErr_timp");
0421    b_recErr_limp = fChain->GetBranch("recErr_limp");
0422    b_refErr_ptinv = fChain->GetBranch("refErr_ptinv");
0423    b_refErr_theta = fChain->GetBranch("refErr_theta");
0424    b_refErr_phi = fChain->GetBranch("refErr_phi");
0425    b_refErr_timp = fChain->GetBranch("refErr_timp");
0426    b_refErr_limp = fChain->GetBranch("refErr_limp");
0427 
0428    return kTRUE;
0429 }
0430 
0431 void SimpleVertexAnalysis::Show(Int_t entry)
0432 {
0433 // Print contents of entry.
0434 // If entry is not specified, print current entry
0435    if (!fChain) return;
0436    fChain->Show(entry);
0437 }
0438 Int_t SimpleVertexAnalysis::Cut(Int_t entry)
0439 {
0440 // This function may be called from Loop.
0441 // returns  1 if entry is accepted.
0442 // returns -1 otherwise.
0443    return entry;
0444 }
0445 #endif // #ifdef SimpleVertexAnalysis_cxx
0446