Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:33:08

0001 #include "Validation/RecoMuon/plugins/MuonTrackValidator.h"
0002 #include "FWCore/Framework/interface/MakerMacros.h"
0003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0004 
0005 #include "DataFormats/TrackReco/interface/Track.h"
0006 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0007 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
0008 #include "SimDataFormats/Vertex/interface/SimVertexContainer.h"
0009 #include "SimDataFormats/PileupSummaryInfo/interface/PileupSummaryInfo.h"
0010 #include "SimTracker/TrackerHitAssociation/interface/TrackerHitAssociator.h"
0011 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
0012 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
0013 #include "SimTracker/TrackAssociation/interface/TrackingParticleIP.h"
0014 
0015 #include "TMath.h"
0016 #include <TF1.h>
0017 
0018 using namespace std;
0019 using namespace edm;
0020 
0021 void MuonTrackValidator::bookHistograms(DQMEDAnalyzer::DQMStore::IBooker& ibooker,
0022                                         edm::Run const&,
0023                                         edm::EventSetup const& setup) {
0024   for (unsigned int ww = 0; ww < associators.size(); ww++) {
0025     for (unsigned int www = 0; www < label.size(); www++) {
0026       ibooker.cd();
0027       InputTag algo = label[www];
0028       string dirName = dirName_;
0029 
0030       auto setBinLogX = [this](TH1* th1) {
0031         if (this->useLogPt) {
0032           BinLogX(th1);
0033         }
0034       };
0035 
0036       if (!algo.process().empty())
0037         dirName += algo.process() + "_";
0038       if (!algo.label().empty())
0039         dirName += algo.label();
0040       if (!algo.instance().empty())
0041         dirName += ("_" + algo.instance());
0042       if (dirName.find("Tracks") < dirName.length()) {
0043         dirName.replace(dirName.find("Tracks"), 6, "Trks");
0044       }
0045       if (dirName.find("UpdatedAtVtx") < dirName.length()) {
0046         dirName.replace(dirName.find("UpdatedAtVtx"), 12, "UpdAtVtx");
0047       }
0048       string assoc = associators[ww];
0049       if (assoc.find("tpToTkmuTrackAssociation") < assoc.length()) {
0050         dirName += "_TkAsso";
0051       }
0052       std::replace(dirName.begin(), dirName.end(), ':', '_');
0053       ibooker.setCurrentFolder(dirName);
0054 
0055       h_tracks.push_back(
0056           ibooker.book1D("Ntracks", "Number of reconstructed tracks", nintNTracks, minNTracks, maxNTracks));
0057       h_fakes.push_back(ibooker.book1D("Nfakes", "Number of fake reco tracks", nintFTracks, minFTracks, maxFTracks));
0058       h_charge.push_back(ibooker.book1D("Ncharge", "track charge", 3, -1.5, 1.5));
0059 
0060       h_recoeta.push_back(ibooker.book1D("num_reco_eta", "N of reco track vs eta", nintEta, minEta, maxEta));
0061       h_assoceta.push_back(ibooker.book1D(
0062           "num_assoSimToReco_eta", "N of associated tracks (simToReco) vs eta", nintEta, minEta, maxEta));
0063       h_assoc2eta.push_back(ibooker.book1D(
0064           "num_assoRecoToSim_eta", "N of associated (recoToSim) tracks vs eta", nintEta, minEta, maxEta));
0065       h_simuleta.push_back(ibooker.book1D("num_simul_eta", "N of simulated tracks vs eta", nintEta, minEta, maxEta));
0066       h_misideta.push_back(ibooker.book1D("num_chargemisid_eta",
0067                                           "N of associated (simToReco) tracks with charge misID vs eta",
0068                                           nintEta,
0069                                           minEta,
0070                                           maxEta));
0071 
0072       h_recopT.push_back(ibooker.book1D("num_reco_pT", "N of reco track vs pT", nintPt, minPt, maxPt, setBinLogX));
0073       h_assocpT.push_back(ibooker.book1D(
0074           "num_assoSimToReco_pT", "N of associated tracks (simToReco) vs pT", nintPt, minPt, maxPt, setBinLogX));
0075       h_assoc2pT.push_back(ibooker.book1D(
0076           "num_assoRecoToSim_pT", "N of associated (recoToSim) tracks vs pT", nintPt, minPt, maxPt, setBinLogX));
0077       h_simulpT.push_back(
0078           ibooker.book1D("num_simul_pT", "N of simulated tracks vs pT", nintPt, minPt, maxPt, setBinLogX));
0079       h_misidpT.push_back(ibooker.book1D("num_chargemisid_pT",
0080                                          "N of associated (simToReco) tracks with charge misID vs pT",
0081                                          nintPt,
0082                                          minPt,
0083                                          maxPt,
0084                                          setBinLogX));
0085 
0086       h_recophi.push_back(ibooker.book1D("num_reco_phi", "N of reco track vs phi", nintPhi, minPhi, maxPhi));
0087       h_assocphi.push_back(ibooker.book1D(
0088           "num_assoSimToReco_phi", "N of associated tracks (simToReco) vs phi", nintPhi, minPhi, maxPhi));
0089       h_assoc2phi.push_back(ibooker.book1D(
0090           "num_assoRecoToSim_phi", "N of associated (recoToSim) tracks vs phi", nintPhi, minPhi, maxPhi));
0091       h_simulphi.push_back(ibooker.book1D("num_simul_phi", "N of simulated tracks vs phi", nintPhi, minPhi, maxPhi));
0092       h_misidphi.push_back(ibooker.book1D("num_chargemisid_phi",
0093                                           "N of associated (simToReco) tracks with charge misID vs phi",
0094                                           nintPhi,
0095                                           minPhi,
0096                                           maxPhi));
0097 
0098       h_recohit.push_back(ibooker.book1D("num_reco_hit", "N of reco tracks vs N SimHits", nintNHit, minNHit, maxNHit));
0099       h_assochit.push_back(ibooker.book1D(
0100           "num_assoSimToReco_hit", "N of associated tracks (simToReco) vs N SimHits", nintNHit, minNHit, maxNHit));
0101       h_assoc2hit.push_back(ibooker.book1D(
0102           "num_assoRecoToSim_hit", "N of associated (recoToSim) tracks vs N Rechits", nintNHit, minNHit, maxNHit));
0103       h_simulhit.push_back(
0104           ibooker.book1D("num_simul_hit", "N of simulated tracks vs N SimHits", nintNHit, minNHit, maxNHit));
0105       h_misidhit.push_back(ibooker.book1D("num_chargemisid_hit",
0106                                           "N of associated (recoToSim) tracks with charge misID vs N RecHits",
0107                                           nintNHit,
0108                                           minNHit,
0109                                           maxNHit));
0110 
0111       h_recodxy.push_back(ibooker.book1D("num_reco_dxy", "N of reco track vs dxy", nintDxy, minDxy, maxDxy));
0112       h_assocdxy.push_back(ibooker.book1D(
0113           "num_assoSimToReco_dxy", "N of associated tracks (simToReco) vs dxy", nintDxy, minDxy, maxDxy));
0114       h_assoc2dxy.push_back(ibooker.book1D(
0115           "num_assoRecoToSim_dxy", "N of associated (recoToSim) tracks vs dxy", nintDxy, minDxy, maxDxy));
0116       h_simuldxy.push_back(ibooker.book1D("num_simul_dxy", "N of simulated tracks vs dxy", nintDxy, minDxy, maxDxy));
0117       h_misiddxy.push_back(ibooker.book1D("num_chargemisid_dxy",
0118                                           "N of associated (simToReco) tracks with charge misID vs dxy",
0119                                           nintDxy,
0120                                           minDxy,
0121                                           maxDxy));
0122       h_recodz.push_back(ibooker.book1D("num_reco_dz", "N of reco track vs dz", nintDz, minDz, maxDz));
0123       h_assocdz.push_back(
0124           ibooker.book1D("num_assoSimToReco_dz", "N of associated tracks (simToReco) vs dz", nintDz, minDz, maxDz));
0125       h_assoc2dz.push_back(
0126           ibooker.book1D("num_assoRecoToSim_dz", "N of associated (recoToSim) tracks vs dz", nintDz, minDz, maxDz));
0127       h_simuldz.push_back(ibooker.book1D("num_simul_dz", "N of simulated tracks vs dz", nintDz, minDz, maxDz));
0128       h_misiddz.push_back(ibooker.book1D(
0129           "num_chargemisid_dz", "N of associated (simToReco) tracks with charge misID vs dz", nintDz, minDz, maxDz));
0130 
0131       h_assocRpos.push_back(ibooker.book1D(
0132           "num_assoSimToReco_Rpos", "N of associated tracks (simToReco) vs Radius", nintRpos, minRpos, maxRpos));
0133       h_simulRpos.push_back(
0134           ibooker.book1D("num_simul_Rpos", "N of simulated tracks vs Radius", nintRpos, minRpos, maxRpos));
0135 
0136       h_assocZpos.push_back(ibooker.book1D(
0137           "num_assoSimToReco_Zpos", "N of associated tracks (simToReco) vs Z", nintZpos, minZpos, maxZpos));
0138       h_simulZpos.push_back(ibooker.book1D("num_simul_Zpos", "N of simulated tracks vs Z", nintZpos, minZpos, maxZpos));
0139 
0140       h_recopu.push_back(ibooker.book1D("num_reco_pu", "N of reco track vs pu", nintPU, minPU, maxPU));
0141       h_assocpu.push_back(
0142           ibooker.book1D("num_assoSimToReco_pu", "N of associated tracks (simToReco) vs pu", nintPU, minPU, maxPU));
0143       h_assoc2pu.push_back(
0144           ibooker.book1D("num_assoRecoToSim_pu", "N of associated (recoToSim) tracks vs pu", nintPU, minPU, maxPU));
0145       h_simulpu.push_back(ibooker.book1D("num_simul_pu", "N of simulated tracks vs pu", nintPU, minPU, maxPU));
0146       h_misidpu.push_back(ibooker.book1D(
0147           "num_chargemisid_pu", "N of associated (simToReco) charge misIDed tracks vs pu", nintPU, minPU, maxPU));
0148 
0149       h_nchi2.push_back(ibooker.book1D("chi2", "Track normalized #chi^{2}", 80, 0., 20.));
0150       h_nchi2_prob.push_back(ibooker.book1D("chi2prob", "Probability of track normalized #chi^{2}", 100, 0., 1.));
0151 
0152       chi2_vs_nhits.push_back(
0153           ibooker.book2D("chi2_vs_nhits", "#chi^{2} vs nhits", nintNHit, minNHit, maxNHit, 20, 0., 10.));
0154       chi2_vs_eta.push_back(ibooker.book2D("chi2_vs_eta", "chi2_vs_eta", nintEta, minEta, maxEta, 40, 0., 20.));
0155       chi2_vs_phi.push_back(ibooker.book2D("chi2_vs_phi", "#chi^{2} vs #phi", nintPhi, minPhi, maxPhi, 40, 0., 20.));
0156 
0157       h_nhits.push_back(ibooker.book1D("nhits", "Number of hits per track", nintNHit, minNHit, maxNHit));
0158       nhits_vs_eta.push_back(
0159           ibooker.book2D("nhits_vs_eta", "Number of Hits vs eta", nintEta, minEta, maxEta, nintNHit, minNHit, maxNHit));
0160       nhits_vs_phi.push_back(
0161           ibooker.book2D("nhits_vs_phi", "#hits vs #phi", nintPhi, minPhi, maxPhi, nintNHit, minNHit, maxNHit));
0162 
0163       if (do_MUOhitsPlots) {
0164         nDThits_vs_eta.push_back(ibooker.book2D(
0165             "nDThits_vs_eta", "Number of DT hits vs eta", nintEta, minEta, maxEta, nintDTHit, minDTHit, maxDTHit));
0166         nCSChits_vs_eta.push_back(ibooker.book2D(
0167             "nCSChits_vs_eta", "Number of CSC hits vs eta", nintEta, minEta, maxEta, nintCSCHit, minCSCHit, maxCSCHit));
0168         nRPChits_vs_eta.push_back(ibooker.book2D(
0169             "nRPChits_vs_eta", "Number of RPC hits vs eta", nintEta, minEta, maxEta, nintRPCHit, minRPCHit, maxRPCHit));
0170         if (useGEMs_)
0171           nGEMhits_vs_eta.push_back(ibooker.book2D(
0172               "nGEMhits_vs_eta", "Number of GEM hits vs eta", nintEta, minEta, maxEta, nintNHit, minNHit, maxNHit));
0173         if (useME0_)
0174           nME0hits_vs_eta.push_back(ibooker.book2D(
0175               "nME0hits_vs_eta", "Number of ME0 hits vs eta", nintEta, minEta, maxEta, nintNHit, minNHit, maxNHit));
0176       }
0177 
0178       if (do_TRKhitsPlots) {
0179         nTRK_LayersWithMeas_vs_eta.push_back(ibooker.book2D("nTRK_LayersWithMeas_vs_eta",
0180                                                             "# TRK Layers with measurement vs eta",
0181                                                             nintEta,
0182                                                             minEta,
0183                                                             maxEta,
0184                                                             nintLayers,
0185                                                             minLayers,
0186                                                             maxLayers));
0187         nPixel_LayersWithMeas_vs_eta.push_back(ibooker.book2D("nPixel_LayersWithMeas_vs_eta",
0188                                                               "Number of Pixel Layers with measurement vs eta",
0189                                                               nintEta,
0190                                                               minEta,
0191                                                               maxEta,
0192                                                               nintPixels,
0193                                                               minPixels,
0194                                                               maxPixels));
0195         h_nmisslayers_inner.push_back(ibooker.book1D(
0196             "nTRK_misslayers_inner", "Number of missing inner TRK layers", nintLayers, minLayers, maxLayers));
0197         h_nmisslayers_outer.push_back(ibooker.book1D(
0198             "nTRK_misslayers_outer", "Number of missing outer TRK layers", nintLayers, minLayers, maxLayers));
0199         h_nlosthits.push_back(ibooker.book1D("nlosthits", "Number of lost hits per track", 6, -0.5, 5.5));
0200         nlosthits_vs_eta.push_back(ibooker.book2D(
0201             "nlosthits_vs_eta", "Number of lost hits per track vs eta", nintEta, minEta, maxEta, 6, -0.5, 5.5));
0202       }
0203 
0204       ptres_vs_eta.push_back(ibooker.book2D("ptres_vs_eta",
0205                                             "p_{T} Relative Residual vs #eta",
0206                                             nintEta,
0207                                             minEta,
0208                                             maxEta,
0209                                             ptRes_nbin,
0210                                             ptRes_rangeMin,
0211                                             ptRes_rangeMax));
0212       ptres_vs_phi.push_back(ibooker.book2D("ptres_vs_phi",
0213                                             "p_{T} Relative Residual vs #phi",
0214                                             nintPhi,
0215                                             minPhi,
0216                                             maxPhi,
0217                                             ptRes_nbin,
0218                                             ptRes_rangeMin,
0219                                             ptRes_rangeMax));
0220       ptres_vs_pt.push_back(ibooker.book2D("ptres_vs_pt",
0221                                            "p_{T} Relative Residual vs p_{T}",
0222                                            nintPt,
0223                                            minPt,
0224                                            maxPt,
0225                                            ptRes_nbin,
0226                                            ptRes_rangeMin,
0227                                            ptRes_rangeMax,
0228                                            setBinLogX));
0229       h_ptpull.push_back(ibooker.book1D("ptpull", "p_{T} Pull", 100, -10., 10.));
0230       ptpull_vs_eta.push_back(
0231           ibooker.book2D("ptpull_vs_eta", "p_{T} Pull vs #eta", nintEta, minEta, maxEta, 100, -10., 10.));
0232       ptpull_vs_phi.push_back(
0233           ibooker.book2D("ptpull_vs_phi", "p_{T} Pull vs #phi", nintPhi, minPhi, maxPhi, 100, -10., 10.));
0234       h_qoverppull.push_back(ibooker.book1D("qoverppull", "q/p Pull", 100, -10., 10.));
0235 
0236       h_etaRes.push_back(ibooker.book1D("etaRes", "#eta residual", etaRes_nbin, etaRes_rangeMin, etaRes_rangeMax));
0237       etares_vs_eta.push_back(ibooker.book2D("etares_vs_eta",
0238                                              "#eta Residual vs #eta",
0239                                              nintEta,
0240                                              minEta,
0241                                              maxEta,
0242                                              etaRes_nbin,
0243                                              etaRes_rangeMin,
0244                                              etaRes_rangeMax));
0245 
0246       thetaCotres_vs_eta.push_back(ibooker.book2D("thetaCotres_vs_eta",
0247                                                   "cot(#theta) Residual vs #eta",
0248                                                   nintEta,
0249                                                   minEta,
0250                                                   maxEta,
0251                                                   cotThetaRes_nbin,
0252                                                   cotThetaRes_rangeMin,
0253                                                   cotThetaRes_rangeMax));
0254       thetaCotres_vs_pt.push_back(ibooker.book2D("thetaCotres_vs_pt",
0255                                                  "cot(#theta) Residual vs p_{T}",
0256                                                  nintPt,
0257                                                  minPt,
0258                                                  maxPt,
0259                                                  cotThetaRes_nbin,
0260                                                  cotThetaRes_rangeMin,
0261                                                  cotThetaRes_rangeMax,
0262                                                  setBinLogX));
0263       h_thetapull.push_back(ibooker.book1D("thetapull", "#theta Pull", 100, -10., 10.));
0264       thetapull_vs_eta.push_back(
0265           ibooker.book2D("thetapull_vs_eta", "#theta Pull vs #eta", nintEta, minEta, maxEta, 100, -10, 10));
0266       thetapull_vs_phi.push_back(
0267           ibooker.book2D("thetapull_vs_phi", "#theta Pull vs #phi", nintPhi, minPhi, maxPhi, 100, -10, 10));
0268 
0269       phires_vs_eta.push_back(ibooker.book2D("phires_vs_eta",
0270                                              "#phi Residual vs #eta",
0271                                              nintEta,
0272                                              minEta,
0273                                              maxEta,
0274                                              phiRes_nbin,
0275                                              phiRes_rangeMin,
0276                                              phiRes_rangeMax));
0277       phires_vs_pt.push_back(ibooker.book2D("phires_vs_pt",
0278                                             "#phi Residual vs p_{T}",
0279                                             nintPt,
0280                                             minPt,
0281                                             maxPt,
0282                                             phiRes_nbin,
0283                                             phiRes_rangeMin,
0284                                             phiRes_rangeMax,
0285                                             setBinLogX));
0286       phires_vs_phi.push_back(ibooker.book2D("phires_vs_phi",
0287                                              "#phi Residual vs #phi",
0288                                              nintPhi,
0289                                              minPhi,
0290                                              maxPhi,
0291                                              phiRes_nbin,
0292                                              phiRes_rangeMin,
0293                                              phiRes_rangeMax));
0294       h_phipull.push_back(ibooker.book1D("phipull", "#phi Pull", 100, -10., 10.));
0295       phipull_vs_eta.push_back(
0296           ibooker.book2D("phipull_vs_eta", "#phi Pull vs #eta", nintEta, minEta, maxEta, 100, -10, 10));
0297       phipull_vs_phi.push_back(
0298           ibooker.book2D("phipull_vs_phi", "#phi Pull vs #phi", nintPhi, minPhi, maxPhi, 100, -10, 10));
0299 
0300       dxyres_vs_eta.push_back(ibooker.book2D("dxyres_vs_eta",
0301                                              "dxy Residual vs #eta",
0302                                              nintEta,
0303                                              minEta,
0304                                              maxEta,
0305                                              dxyRes_nbin,
0306                                              dxyRes_rangeMin,
0307                                              dxyRes_rangeMax));
0308       dxyres_vs_pt.push_back(ibooker.book2D("dxyres_vs_pt",
0309                                             "dxy Residual vs p_{T}",
0310                                             nintPt,
0311                                             minPt,
0312                                             maxPt,
0313                                             dxyRes_nbin,
0314                                             dxyRes_rangeMin,
0315                                             dxyRes_rangeMax,
0316                                             setBinLogX));
0317       h_dxypull.push_back(ibooker.book1D("dxypull", "dxy Pull", 100, -10., 10.));
0318       dxypull_vs_eta.push_back(
0319           ibooker.book2D("dxypull_vs_eta", "dxy Pull vs #eta", nintEta, minEta, maxEta, 100, -10, 10));
0320 
0321       dzres_vs_eta.push_back(ibooker.book2D(
0322           "dzres_vs_eta", "dz Residual vs #eta", nintEta, minEta, maxEta, dzRes_nbin, dzRes_rangeMin, dzRes_rangeMax));
0323       dzres_vs_pt.push_back(ibooker.book2D("dzres_vs_pt",
0324                                            "dz Residual vs p_{T}",
0325                                            nintPt,
0326                                            minPt,
0327                                            maxPt,
0328                                            dzRes_nbin,
0329                                            dzRes_rangeMin,
0330                                            dzRes_rangeMax,
0331                                            setBinLogX));
0332       h_dzpull.push_back(ibooker.book1D("dzpull", "dz Pull", 100, -10., 10.));
0333       dzpull_vs_eta.push_back(
0334           ibooker.book2D("dzpull_vs_eta", "dz Pull vs #eta", nintEta, minEta, maxEta, 100, -10, 10));
0335 
0336       nRecHits_vs_nSimHits.push_back(ibooker.book2D(
0337           "nRecHits_vs_nSimHits", "nRecHits vs nSimHits", nintNHit, minNHit, maxNHit, nintNHit, minNHit, maxNHit));
0338 
0339       if (MABH) {
0340         h_PurityVsQuality.push_back(
0341             ibooker.book2D("PurityVsQuality", "Purity vs Quality (MABH)", 20, 0.01, 1.01, 20, 0.01, 1.01));
0342       }
0343 
0344       if (associators[ww] == "trackAssociatorByChi2") {
0345         h_assochi2.push_back(ibooker.book1D("assocChi2", "track association #chi^{2}", 1000, 0., 100.));
0346         h_assochi2_prob.push_back(ibooker.book1D("assocChi2_prob", "probability of association #chi^{2}", 100, 0., 1.));
0347       } else if (associators[ww] == "trackAssociatorByHits") {
0348         h_assocFraction.push_back(ibooker.book1D("assocFraction", "fraction of shared hits", 22, 0., 1.1));
0349         h_assocSharedHit.push_back(ibooker.book1D("assocSharedHit", "number of shared hits", 41, -0.5, 40.5));
0350       }
0351 
0352     }  //for (unsigned int www=0;www<label.size();www++)
0353   }    //for (unsigned int ww=0;ww<associators.size();ww++)
0354 }
0355 
0356 void MuonTrackValidator::analyze(const edm::Event& event, const edm::EventSetup& setup) {
0357   using namespace reco;
0358 
0359   edm::LogInfo("MuonTrackValidator") << "\n===================================================="
0360                                      << "\n"
0361                                      << "Analyzing new event"
0362                                      << "\n"
0363                                      << "====================================================\n"
0364                                      << "\n";
0365 
0366   edm::Handle<std::vector<PileupSummaryInfo> > puinfoH;
0367   int PU_NumInteractions(-1);
0368 
0369   if (parametersDefiner == "LhcParametersDefinerForTP") {
0370     // PileupSummaryInfo is contained only in collision events
0371     event.getByToken(pileupinfo_Token, puinfoH);
0372     for (std::vector<PileupSummaryInfo>::const_iterator puInfoIt = puinfoH->begin(); puInfoIt != puinfoH->end();
0373          ++puInfoIt) {
0374       if (puInfoIt->getBunchCrossing() == 0) {
0375         PU_NumInteractions = puInfoIt->getPU_NumInteractions();
0376         break;
0377       }
0378     }
0379 
0380   } else if (parametersDefiner == "CosmicParametersDefinerForTP") {
0381     edm::Handle<SimHitTPAssociationProducer::SimHitTPAssociationList> simHitsTPAssoc;
0382     //warning: make sure the TP collection used in the map is the same used here
0383     event.getByToken(_simHitTpMapTag, simHitsTPAssoc);
0384     cosmicParametersDefinerTP_->initEvent(simHitsTPAssoc);
0385     cosmictpSelector.initEvent(simHitsTPAssoc);
0386   }
0387 
0388   TrackingParticleRefVector TPrefV;
0389   const TrackingParticleRefVector* ptr_TPrefV = nullptr;
0390   edm::Handle<TrackingParticleCollection> TPCollection_H;
0391   edm::Handle<TrackingParticleRefVector> TPCollectionRefVector_H;
0392 
0393   if (label_tp_refvector) {
0394     event.getByToken(tp_refvector_Token, TPCollectionRefVector_H);
0395     ptr_TPrefV = TPCollectionRefVector_H.product();
0396   } else {
0397     event.getByToken(tp_Token, TPCollection_H);
0398     size_t nTP = TPCollection_H->size();
0399     for (size_t i = 0; i < nTP; ++i) {
0400       TPrefV.push_back(TrackingParticleRef(TPCollection_H, i));
0401     }
0402     ptr_TPrefV = &TPrefV;
0403   }
0404   TrackingParticleRefVector const& tPC = *ptr_TPrefV;
0405 
0406   edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
0407   bool bs_Available = event.getByToken(bsSrc_Token, recoBeamSpotHandle);
0408   reco::BeamSpot bs;
0409   if (bs_Available)
0410     bs = *recoBeamSpotHandle;
0411   edm::LogVerbatim("MuonTrackValidator") << bs;
0412 
0413   std::vector<const reco::TrackToTrackingParticleAssociator*> associator;
0414   if (UseAssociators) {
0415     edm::Handle<reco::TrackToTrackingParticleAssociator> theAssociator;
0416     for (auto const& associatorName : associators) {
0417       event.getByLabel(associatorName, theAssociator);
0418       associator.push_back(theAssociator.product());
0419     }
0420   }
0421 
0422   int w = 0;
0423   for (unsigned int ww = 0; ww < associators.size(); ww++) {
0424     for (unsigned int www = 0; www < label.size(); www++) {
0425       //
0426       //get collections from the event
0427       //
0428       edm::Handle<edm::View<Track> > trackCollection;
0429       unsigned int trackCollectionSize = 0;
0430 
0431       reco::RecoToSimCollection recSimColl;
0432       reco::SimToRecoCollection simRecColl;
0433 
0434       // account for missing track collections (HLT case)
0435       if (!event.getByToken(track_Collection_Token[www], trackCollection) && ignoremissingtkcollection_) {
0436         recSimColl.post_insert();
0437         simRecColl.post_insert();
0438       }
0439 
0440       //associate tracks to TrackingParticles
0441       else {
0442         trackCollectionSize = trackCollection->size();
0443 
0444         if (UseAssociators) {
0445           edm::LogVerbatim("MuonTrackValidator")
0446               << "Analyzing " << label[www].process() << ":" << label[www].label() << ":" << label[www].instance()
0447               << " with " << associators[ww].c_str() << "\n";
0448 
0449           LogTrace("MuonTrackValidator") << "Calling associateRecoToSim method"
0450                                          << "\n";
0451           recSimColl = associator[ww]->associateRecoToSim(trackCollection, TPCollection_H);
0452           LogTrace("MuonTrackValidator") << "Calling associateSimToReco method"
0453                                          << "\n";
0454           simRecColl = associator[ww]->associateSimToReco(trackCollection, TPCollection_H);
0455         } else {
0456           edm::LogVerbatim("MuonTrackValidator")
0457               << "Analyzing " << label[www].process() << ":" << label[www].label() << ":" << label[www].instance()
0458               << " with " << associatormap.process() << ":" << associatormap.label() << ":" << associatormap.instance()
0459               << "\n";
0460 
0461           Handle<reco::SimToRecoCollection> simtorecoCollectionH;
0462           event.getByToken(simToRecoCollection_Token, simtorecoCollectionH);
0463           simRecColl = *simtorecoCollectionH.product();
0464 
0465           Handle<reco::RecoToSimCollection> recotosimCollectionH;
0466           event.getByToken(recoToSimCollection_Token, recotosimCollectionH);
0467           recSimColl = *recotosimCollectionH.product();
0468         }
0469       }
0470 
0471       //
0472       //fill simulation histograms
0473       //
0474       edm::LogVerbatim("MuonTrackValidator") << "\n# of TrackingParticles: " << tPC.size() << "\n";
0475       int ats = 0;
0476       int st = 0;
0477 
0478       for (size_t i = 0; i < tPC.size(); i++) {
0479         bool TP_is_matched = false;
0480         bool isChargeOK = true;
0481         double quality = 0.;
0482 
0483         const TrackingParticleRef& tpr = tPC[i];
0484         const TrackingParticle& tp = *tpr;
0485 
0486         TrackingParticle::Vector momentumTP;
0487         TrackingParticle::Point vertexTP;
0488         double dxySim = 0;
0489         double dzSim = 0;
0490 
0491         //If the TrackingParticle is collision-like, get the momentum and vertex at production state
0492         //and the impact parameters w.r.t. PCA
0493         if (parametersDefiner == "LhcParametersDefinerForTP") {
0494           LogTrace("MuonTrackValidator") << "TrackingParticle " << i;
0495           if (!tpSelector(tp))
0496             continue;
0497           momentumTP = tp.momentum();
0498           vertexTP = tp.vertex();
0499           TrackingParticle::Vector momentum = lhcParametersDefinerTP_->momentum(event, setup, tpr);
0500           TrackingParticle::Point vertex = lhcParametersDefinerTP_->vertex(event, setup, tpr);
0501           dxySim = TrackingParticleIP::dxy(vertex, momentum, bs.position());
0502           dzSim = TrackingParticleIP::dz(vertex, momentum, bs.position());
0503         }
0504         //for cosmics get the momentum and vertex at PCA
0505         else if (parametersDefiner == "CosmicParametersDefinerForTP") {
0506           edm::LogVerbatim("MuonTrackValidator") << "TrackingParticle " << i;
0507           if (!cosmictpSelector(tpr, &bs, event, setup))
0508             continue;
0509           momentumTP = cosmicParametersDefinerTP_->momentum(event, setup, tpr);
0510           vertexTP = cosmicParametersDefinerTP_->vertex(event, setup, tpr);
0511           dxySim = TrackingParticleIP::dxy(vertexTP, momentumTP, bs.position());
0512           dzSim = TrackingParticleIP::dz(vertexTP, momentumTP, bs.position());
0513         }
0514 
0515         // Number of counted SimHits depend on the selection of tracker and muon detectors (via cfg parameters)
0516         int nSimHits = 0;
0517         if (usetracker && usemuon) {
0518           nSimHits = tpr.get()->numberOfHits();
0519         } else if (!usetracker && usemuon) {
0520           nSimHits = tpr.get()->numberOfHits() - tpr.get()->numberOfTrackerHits();
0521         } else if (usetracker && !usemuon) {
0522           nSimHits = tpr.get()->numberOfTrackerHits();
0523         }
0524 
0525         edm::LogVerbatim("MuonTrackValidator") << "--------------------Selected TrackingParticle #" << tpr.key()
0526                                                << "  (N counted simhits = " << nSimHits << ")";
0527         edm::LogVerbatim("MuonTrackValidator")
0528             << "momentumTP: pt = " << sqrt(momentumTP.perp2()) << ", pz = " << momentumTP.z()
0529             << ", \t vertexTP: radius = " << sqrt(vertexTP.perp2()) << ",  z = " << vertexTP.z();
0530         st++;
0531 
0532         double TPeta = momentumTP.eta();
0533         double xTPeta = getEta(TPeta);  // may be |eta| in histos according to useFabsEta
0534         double TPpt = sqrt(momentumTP.perp2());
0535         double xTPpt = getPt(TPpt);  // may be 1/pt in histos according to useInvPt
0536         double TPphi = momentumTP.phi();
0537         double TPrpos = sqrt(vertexTP.perp2());
0538         double TPzpos = vertexTP.z();
0539 
0540         int assoc_recoTrack_NValidHits = 0;
0541         if (simRecColl.find(tpr) != simRecColl.end()) {
0542           auto const& rt = simRecColl[tpr];
0543           if (!rt.empty()) {
0544             RefToBase<Track> assoc_recoTrack = rt.begin()->first;
0545             TP_is_matched = true;
0546             ats++;
0547             if (assoc_recoTrack->charge() != tpr->charge())
0548               isChargeOK = false;
0549             quality = rt.begin()->second;
0550             assoc_recoTrack_NValidHits = assoc_recoTrack->numberOfValidHits();
0551             edm::LogVerbatim("MuonTrackValidator") << "-----------------------------associated to Track #"
0552                                                    << assoc_recoTrack.key() << " with quality:" << quality << "\n";
0553           }
0554         } else {
0555           edm::LogVerbatim("MuonTrackValidator")
0556               << "TrackingParticle #" << tpr.key() << " with pt,eta,phi: " << sqrt(momentumTP.perp2()) << " , "
0557               << momentumTP.eta() << " , " << momentumTP.phi() << " , "
0558               << " NOT associated to any reco::Track"
0559               << "\n";
0560         }
0561 
0562         // histos for efficiency vs eta
0563         fillPlotNoFlow(h_simuleta[w], xTPeta);
0564         if (TP_is_matched) {
0565           fillPlotNoFlow(h_assoceta[w], xTPeta);
0566           if (!isChargeOK)
0567             fillPlotNoFlow(h_misideta[w], xTPeta);
0568         }
0569 
0570         // histos for efficiency vs phi
0571         fillPlotNoFlow(h_simulphi[w], TPphi);
0572         if (TP_is_matched) {
0573           fillPlotNoFlow(h_assocphi[w], TPphi);
0574           if (!isChargeOK)
0575             fillPlotNoFlow(h_misidphi[w], TPphi);
0576         }
0577 
0578         // histos for efficiency vs pT
0579         fillPlotNoFlow(h_simulpT[w], xTPpt);
0580         if (TP_is_matched) {
0581           fillPlotNoFlow(h_assocpT[w], xTPpt);
0582           if (!isChargeOK)
0583             fillPlotNoFlow(h_misidpT[w], xTPpt);
0584         }
0585 
0586         // histos for efficiency vs dxy
0587         fillPlotNoFlow(h_simuldxy[w], dxySim);
0588         if (TP_is_matched) {
0589           fillPlotNoFlow(h_assocdxy[w], dxySim);
0590           if (!isChargeOK)
0591             fillPlotNoFlow(h_misiddxy[w], dxySim);
0592         }
0593 
0594         // histos for efficiency vs dz
0595         fillPlotNoFlow(h_simuldz[w], dzSim);
0596         if (TP_is_matched) {
0597           fillPlotNoFlow(h_assocdz[w], dzSim);
0598           if (!isChargeOK)
0599             fillPlotNoFlow(h_misiddz[w], dzSim);
0600         }
0601 
0602         // histos for efficiency vs Radius
0603         fillPlotNoFlow(h_simulRpos[w], TPrpos);
0604         if (TP_is_matched)
0605           fillPlotNoFlow(h_assocRpos[w], TPrpos);
0606 
0607         // histos for efficiency vs z position
0608         fillPlotNoFlow(h_simulZpos[w], TPzpos);
0609         if (TP_is_matched)
0610           fillPlotNoFlow(h_assocZpos[w], TPzpos);
0611 
0612         // histos for efficiency vs Number of Hits
0613         fillPlotNoFlow(h_simulhit[w], nSimHits);
0614         if (TP_is_matched) {
0615           fillPlotNoFlow(h_assochit[w], nSimHits);
0616           nRecHits_vs_nSimHits[w]->Fill(nSimHits, assoc_recoTrack_NValidHits);
0617 
0618           // charge misid is more useful w.r.t. nRecHits (filled after)
0619           //if (!isChargeOK) fillPlotNoFlow(h_misidhit[w], nSimHits);
0620         }
0621 
0622         // histos for efficiency vs PileUp
0623         fillPlotNoFlow(h_simulpu[w], PU_NumInteractions);
0624         if (TP_is_matched) {
0625           fillPlotNoFlow(h_assocpu[w], PU_NumInteractions);
0626           if (!isChargeOK)
0627             fillPlotNoFlow(h_misidpu[w], PU_NumInteractions);
0628         }
0629 
0630       }  // End for (size_t i = 0; i < tPCeff.size(); i++) {
0631 
0632       //
0633       //fill reconstructed track histograms
0634       //
0635       edm::LogVerbatim("MuonTrackValidator")
0636           << "\n# of reco::Tracks with " << label[www].process() << ":" << label[www].label() << ":"
0637           << label[www].instance() << ": " << trackCollectionSize << "\n";
0638 
0639       int at = 0;
0640       int rT = 0;
0641       for (edm::View<Track>::size_type i = 0; i < trackCollectionSize; ++i) {
0642         bool Track_is_matched = false;
0643         bool isChargeOK = true;
0644         RefToBase<Track> track(trackCollection, i);
0645         int nRecHits = track->numberOfValidHits();
0646         rT++;
0647 
0648         std::vector<std::pair<TrackingParticleRef, double> > tp;
0649         TrackingParticleRef tpr;
0650 
0651         // new logic (bidirectional)
0652         if (BiDirectional_RecoToSim_association) {
0653           edm::LogVerbatim("MuonTrackValidator") << "----------------------------------------Track #" << track.key()
0654                                                  << " (N valid rechits = " << nRecHits << ")";
0655 
0656           if (recSimColl.find(track) != recSimColl.end()) {
0657             tp = recSimColl[track];
0658             if (!tp.empty()) {
0659               tpr = tp.begin()->first;
0660               // RtS and StR must associate the same pair !
0661               if (simRecColl.find(tpr) != simRecColl.end()) {
0662                 auto const& assoc_track_checkback = simRecColl[tpr].begin()->first;
0663 
0664                 if (assoc_track_checkback.key() == track.key()) {
0665                   Track_is_matched = true;
0666                   at++;
0667                   if (track->charge() != tpr->charge())
0668                     isChargeOK = false;
0669                   double Purity = tp.begin()->second;
0670                   double Quality = simRecColl[tpr].begin()->second;
0671                   edm::LogVerbatim("MuonTrackValidator")
0672                       << "with pt=" << track->pt() << " associated with purity:" << Purity << " to TrackingParticle #"
0673                       << tpr.key() << "\n";
0674                   if (MABH)
0675                     h_PurityVsQuality[w]->Fill(Quality, Purity);
0676                 }
0677               }
0678             }
0679           }
0680 
0681           if (!Track_is_matched)
0682             edm::LogVerbatim("MuonTrackValidator")
0683                 << "with pt=" << track->pt() << " NOT associated to any TrackingParticle"
0684                 << "\n";
0685         }
0686         // old logic, valid for cosmics 2-legs reco (bugged for collision scenario)
0687         else {
0688           if (recSimColl.find(track) != recSimColl.end()) {
0689             tp = recSimColl[track];
0690             if (!tp.empty()) {
0691               tpr = tp.begin()->first;
0692               Track_is_matched = true;
0693               at++;
0694               if (track->charge() != tpr->charge())
0695                 isChargeOK = false;
0696               edm::LogVerbatim("MuonTrackValidator") << "reco::Track #" << track.key() << " with pt=" << track->pt()
0697                                                      << " associated with quality:" << tp.begin()->second << "\n";
0698             }
0699           } else {
0700             edm::LogVerbatim("MuonTrackValidator") << "reco::Track #" << track.key() << " with pt=" << track->pt()
0701                                                    << " NOT associated to any TrackingParticle"
0702                                                    << "\n";
0703           }
0704         }
0705 
0706         double etaRec = track->eta();
0707         double xetaRec = getEta(etaRec);
0708 
0709         double ptRec = track->pt();
0710         double xptRec = getPt(ptRec);
0711 
0712         double qoverpRec = track->qoverp();
0713         double phiRec = track->phi();
0714         double thetaRec = track->theta();
0715         double dxyRec = track->dxy(bs.position());
0716         double dzRec = track->dz(bs.position());
0717 
0718         double qoverpError = track->qoverpError();
0719         double ptError = track->ptError();
0720         double thetaError = track->thetaError();
0721         double phiError = track->phiError();
0722         double dxyError = track->dxyError();
0723         double dzError = track->dzError();
0724 
0725         // histos for fake rate vs eta
0726         fillPlotNoFlow(h_recoeta[w], xetaRec);
0727         if (Track_is_matched) {
0728           fillPlotNoFlow(h_assoc2eta[w], xetaRec);
0729         }
0730 
0731         // histos for fake rate vs phi
0732         fillPlotNoFlow(h_recophi[w], phiRec);
0733         if (Track_is_matched) {
0734           fillPlotNoFlow(h_assoc2phi[w], phiRec);
0735         }
0736 
0737         // histos for fake rate vs pT
0738         fillPlotNoFlow(h_recopT[w], xptRec);
0739         if (Track_is_matched) {
0740           fillPlotNoFlow(h_assoc2pT[w], xptRec);
0741         }
0742 
0743         // histos for fake rate vs dxy
0744         fillPlotNoFlow(h_recodxy[w], dxyRec);
0745         if (Track_is_matched) {
0746           fillPlotNoFlow(h_assoc2dxy[w], dxyRec);
0747         }
0748 
0749         // histos for fake rate vs dz
0750         fillPlotNoFlow(h_recodz[w], dzRec);
0751         if (Track_is_matched) {
0752           fillPlotNoFlow(h_assoc2dz[w], dzRec);
0753         }
0754 
0755         // histos for fake rate vs Number of RecHits in track
0756         fillPlotNoFlow(h_recohit[w], nRecHits);
0757         if (Track_is_matched) {
0758           fillPlotNoFlow(h_assoc2hit[w], nRecHits);
0759           // charge misid w.r.t. nRecHits
0760           if (!isChargeOK)
0761             fillPlotNoFlow(h_misidhit[w], nRecHits);
0762         }
0763 
0764         // histos for fake rate vs Number of PU interactions
0765         fillPlotNoFlow(h_recopu[w], PU_NumInteractions);
0766         if (Track_is_matched) {
0767           fillPlotNoFlow(h_assoc2pu[w], PU_NumInteractions);
0768         }
0769 
0770         // Fill other histos
0771         TrackingParticle* tpp = const_cast<TrackingParticle*>(tpr.get());
0772         // TrackingParticle parameters at point of closest approach to the beamline
0773         TrackingParticle::Vector momentumTP;
0774         TrackingParticle::Point vertexTP;
0775 
0776         if (parametersDefiner == "LhcParametersDefinerForTP") {
0777           // following reco plots are made only from tracks associated to selected signal TPs
0778           if (!(Track_is_matched && tpSelector(*tpp)))
0779             continue;
0780           else {
0781             momentumTP = lhcParametersDefinerTP_->momentum(event, setup, tpr);
0782             vertexTP = lhcParametersDefinerTP_->vertex(event, setup, tpr);
0783           }
0784         } else if (parametersDefiner == "CosmicParametersDefinerForTP") {
0785           // following reco plots are made only from tracks associated to selected signal TPs
0786           if (!(Track_is_matched && cosmictpSelector(tpr, &bs, event, setup)))
0787             continue;
0788           else {
0789             momentumTP = cosmicParametersDefinerTP_->momentum(event, setup, tpr);
0790             vertexTP = cosmicParametersDefinerTP_->vertex(event, setup, tpr);
0791           }
0792         }
0793 
0794         if (associators[ww] == "trackAssociatorByChi2") {
0795           //association chi2
0796           double assocChi2 = -tp.begin()->second;  //in association map is stored -chi2
0797           h_assochi2[www]->Fill(assocChi2);
0798           h_assochi2_prob[www]->Fill(TMath::Prob((assocChi2)*5, 5));
0799         } else if (associators[ww] == "trackAssociatorByHits") {
0800           double fraction = tp.begin()->second;
0801           h_assocFraction[www]->Fill(fraction);
0802           h_assocSharedHit[www]->Fill(fraction * nRecHits);
0803         }
0804 
0805         h_charge[w]->Fill(track->charge());
0806 
0807         // Hits
0808         h_nhits[w]->Fill(nRecHits);
0809         nhits_vs_eta[w]->Fill(xetaRec, nRecHits);
0810         nhits_vs_phi[w]->Fill(phiRec, nRecHits);
0811 
0812         if (do_MUOhitsPlots) {
0813           nDThits_vs_eta[w]->Fill(xetaRec, track->hitPattern().numberOfValidMuonDTHits());
0814           nCSChits_vs_eta[w]->Fill(xetaRec, track->hitPattern().numberOfValidMuonCSCHits());
0815           nRPChits_vs_eta[w]->Fill(xetaRec, track->hitPattern().numberOfValidMuonRPCHits());
0816           if (useGEMs_)
0817             nGEMhits_vs_eta[w]->Fill(xetaRec, track->hitPattern().numberOfValidMuonGEMHits());
0818           if (useME0_)
0819             nME0hits_vs_eta[w]->Fill(xetaRec, track->hitPattern().numberOfValidMuonME0Hits());
0820         }
0821 
0822         if (do_TRKhitsPlots) {
0823           nTRK_LayersWithMeas_vs_eta[w]->Fill(xetaRec, track->hitPattern().trackerLayersWithMeasurement());
0824           nPixel_LayersWithMeas_vs_eta[w]->Fill(xetaRec, track->hitPattern().pixelLayersWithMeasurement());
0825           h_nlosthits[w]->Fill(track->numberOfLostHits());
0826           h_nmisslayers_inner[w]->Fill(track->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS));
0827           h_nmisslayers_outer[w]->Fill(track->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_OUTER_HITS));
0828           nlosthits_vs_eta[w]->Fill(xetaRec, track->numberOfLostHits());
0829         }
0830 
0831         // normalized chi2
0832         h_nchi2[w]->Fill(track->normalizedChi2());
0833         h_nchi2_prob[w]->Fill(TMath::Prob(track->chi2(), (int)track->ndof()));
0834         chi2_vs_nhits[w]->Fill(nRecHits, track->normalizedChi2());
0835         chi2_vs_eta[w]->Fill(xetaRec, track->normalizedChi2());
0836         chi2_vs_phi[w]->Fill(phiRec, track->normalizedChi2());
0837 
0838         double ptSim = sqrt(momentumTP.perp2());
0839         double xptSim = getPt(ptSim);
0840         double qoverpSim = tpr->charge() / sqrt(momentumTP.x() * momentumTP.x() + momentumTP.y() * momentumTP.y() +
0841                                                 momentumTP.z() * momentumTP.z());
0842         double etaSim = momentumTP.eta();
0843         double thetaSim = momentumTP.theta();
0844         double phiSim = momentumTP.phi();
0845         double dxySim = TrackingParticleIP::dxy(vertexTP, momentumTP, bs.position());
0846         double dzSim = TrackingParticleIP::dz(vertexTP, momentumTP, bs.position());
0847 
0848         double etares = etaRec - etaSim;
0849         double ptRelRes = (ptRec - ptSim) / ptSim;  // relative residual -> resolution
0850         double ptPull = (ptRec - ptSim) / ptError;
0851         double qoverpPull = (qoverpRec - qoverpSim) / qoverpError;
0852         double thetaPull = (thetaRec - thetaSim) / thetaError;
0853         double phiDiff = phiRec - phiSim;
0854         if (abs(phiDiff) > M_PI) {
0855           if (phiDiff > 0.)
0856             phiDiff = phiDiff - 2. * M_PI;
0857           else
0858             phiDiff = phiDiff + 2. * M_PI;
0859         }
0860         double phiPull = phiDiff / phiError;
0861         double dxyPull = (dxyRec - dxySim) / dxyError;
0862         double dzPull = (dzRec - dzSim) / dzError;
0863 
0864         h_etaRes[w]->Fill(etares);
0865         etares_vs_eta[w]->Fill(xetaRec, etares);
0866 
0867         ptres_vs_eta[w]->Fill(xetaRec, ptRelRes);
0868         ptres_vs_pt[w]->Fill(xptSim, ptRelRes);
0869         ptres_vs_phi[w]->Fill(phiRec, ptRelRes);
0870         h_ptpull[w]->Fill(ptPull);
0871         ptpull_vs_eta[w]->Fill(xetaRec, ptPull);
0872         ptpull_vs_phi[w]->Fill(phiRec, ptPull);
0873         h_qoverppull[w]->Fill(qoverpPull);
0874 
0875         thetaCotres_vs_eta[w]->Fill(xetaRec, cos(thetaRec) / sin(thetaRec) - cos(thetaSim) / sin(thetaSim));
0876         thetaCotres_vs_pt[w]->Fill(xptSim, cos(thetaRec) / sin(thetaRec) - cos(thetaSim) / sin(thetaSim));
0877         h_thetapull[w]->Fill(thetaPull);
0878         thetapull_vs_eta[w]->Fill(xetaRec, thetaPull);
0879         thetapull_vs_phi[w]->Fill(phiRec, thetaPull);
0880 
0881         phires_vs_eta[w]->Fill(xetaRec, phiDiff);
0882         phires_vs_pt[w]->Fill(xptSim, phiDiff);
0883         phires_vs_phi[w]->Fill(phiRec, phiDiff);
0884         h_phipull[w]->Fill(phiPull);
0885         phipull_vs_eta[w]->Fill(xetaRec, phiPull);
0886         phipull_vs_phi[w]->Fill(phiRec, phiPull);
0887 
0888         dxyres_vs_eta[w]->Fill(xetaRec, dxyRec - dxySim);
0889         dxyres_vs_pt[w]->Fill(xptSim, dxyRec - dxySim);
0890         h_dxypull[w]->Fill(dxyPull);
0891         dxypull_vs_eta[w]->Fill(xetaRec, dxyPull);
0892 
0893         dzres_vs_eta[w]->Fill(xetaRec, dzRec - dzSim);
0894         dzres_vs_pt[w]->Fill(xptSim, dzRec - dzSim);
0895         h_dzpull[w]->Fill(dzPull);
0896         dzpull_vs_eta[w]->Fill(xetaRec, dzPull);
0897 
0898         double contrib_Qoverp = qoverpPull * qoverpPull / 5;
0899         double contrib_dxy = dxyPull * dxyPull / 5;
0900         double contrib_dz = dzPull * dzPull / 5;
0901         double contrib_theta = thetaPull * thetaPull / 5;
0902         double contrib_phi = phiPull * phiPull / 5;
0903         double assoChi2 = contrib_Qoverp + contrib_dxy + contrib_dz + contrib_theta + contrib_phi;
0904 
0905         LogTrace("MuonTrackValidator") << "normalized Chi2 (track 5-dofs matching) = " << assoChi2 << "\n"
0906                                        << "\t contrib_Qoverp = " << contrib_Qoverp << "\n"
0907                                        << "\t contrib_theta = " << contrib_theta << "\n"
0908                                        << "\t contrib_phi = " << contrib_phi << "\n"
0909                                        << "\t contrib_dxy = " << contrib_dxy << "\n"
0910                                        << "\t contrib_dz = " << contrib_dz << "\n";
0911 
0912         LogTrace("MuonTrackValidator") << "ptRec = " << ptRec << "\n"
0913                                        << "etaRec = " << etaRec << "\n"
0914                                        << "qoverpRec = " << qoverpRec << "\n"
0915                                        << "thetaRec = " << thetaRec << "\n"
0916                                        << "phiRec = " << phiRec << "\n"
0917                                        << "dxyRec = " << dxyRec << "\n"
0918                                        << "dzRec = " << dzRec << "\n"
0919                                        << ""
0920                                        << "\n"
0921                                        << "qoverpError = " << qoverpError << "\n"
0922                                        << "thetaError = " << thetaError << "\n"
0923                                        << "phiError = " << phiError << "\n"
0924                                        << "dxyError = " << dxyError << "\n"
0925                                        << "dzError = " << dzError << "\n"
0926                                        << ""
0927                                        << "\n"
0928                                        << "ptSim = " << ptSim << "\n"
0929                                        << "etaSim = " << etaSim << "\n"
0930                                        << "qoverpSim = " << qoverpSim << "\n"
0931                                        << "thetaSim = " << thetaSim << "\n"
0932                                        << "phiSim = " << phiSim << "\n"
0933                                        << "dxySim = " << dxySim << "\n"
0934                                        << "dzSim = " << dzSim << "\n";
0935       }  // End of for(edm::View<Track>::size_type i=0; i<trackCollectionSize; ++i) {
0936 
0937       h_tracks[w]->Fill(at);
0938       h_fakes[w]->Fill(rT - at);
0939       edm::LogVerbatim("MuonTrackValidator") << "Total Simulated: " << st << "\n"
0940                                              << "Total Associated (simToReco): " << ats << "\n"
0941                                              << "Total Reconstructed: " << rT << "\n"
0942                                              << "Total Associated (recoToSim): " << at << "\n"
0943                                              << "Total Fakes: " << rT - at << "\n";
0944       w++;
0945     }  // End of for (unsigned int www=0;www<label.size();www++){
0946   }    //END of for (unsigned int ww=0;ww<associators.size();ww++){
0947 }