Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-08-23 03:25:47

0001 #include "Validation/RecoVertex/interface/PrimaryVertexAnalyzer4PUSlimmed.h"
0002 
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/EventSetup.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 
0008 // reco track and vertex
0009 #include "DataFormats/VertexReco/interface/Vertex.h"
0010 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
0011 
0012 // TrackingParticle
0013 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticleFwd.h"
0014 #include "SimDataFormats/TrackingAnalysis/interface/TrackingVertexContainer.h"
0015 
0016 // associator
0017 #include "SimTracker/VertexAssociation/interface/calculateVertexSharedTracks.h"
0018 
0019 // DQM
0020 #include "DQMServices/Core/interface/DQMStore.h"
0021 
0022 #include <numeric>
0023 
0024 namespace {
0025   template <typename T, size_t N>
0026   std::array<T, N + 1> makeLogBins(const double min, const double max) {
0027     const double minLog10 = std::log10(min);
0028     const double maxLog10 = std::log10(max);
0029     const double width = (maxLog10 - minLog10) / N;
0030     std::array<T, N + 1> ret;
0031     ret[0] = std::pow(10, minLog10);
0032     const double mult = std::pow(10, width);
0033     for (size_t i = 1; i <= N; ++i) {
0034       ret[i] = ret[i - 1] * mult;
0035     }
0036     return ret;
0037   }
0038 }  // namespace
0039 
0040 //
0041 // constructors and destructor
0042 //
0043 PrimaryVertexAnalyzer4PUSlimmed::PrimaryVertexAnalyzer4PUSlimmed(const edm::ParameterSet& iConfig)
0044     : verbose_(iConfig.getUntrackedParameter<bool>("verbose", false)),
0045       use_only_charged_tracks_(iConfig.getUntrackedParameter<bool>("use_only_charged_tracks", true)),
0046       do_generic_sim_plots_(iConfig.getUntrackedParameter<bool>("do_generic_sim_plots")),
0047       root_folder_(iConfig.getUntrackedParameter<std::string>("root_folder", "Validation/Vertices")),
0048       vecPileupSummaryInfoToken_(consumes<std::vector<PileupSummaryInfo>>(edm::InputTag(std::string("addPileupInfo")))),
0049       trackingParticleCollectionToken_(consumes<TrackingParticleCollection>(
0050           iConfig.getUntrackedParameter<edm::InputTag>("trackingParticleCollection"))),
0051       trackingVertexCollectionToken_(
0052           consumes<TrackingVertexCollection>(iConfig.getUntrackedParameter<edm::InputTag>("trackingVertexCollection"))),
0053       simToRecoAssociationToken_(
0054           consumes<reco::SimToRecoCollection>(iConfig.getUntrackedParameter<edm::InputTag>("trackAssociatorMap"))),
0055       recoToSimAssociationToken_(
0056           consumes<reco::RecoToSimCollection>(iConfig.getUntrackedParameter<edm::InputTag>("trackAssociatorMap"))),
0057       vertexAssociatorToken_(consumes<reco::VertexToTrackingVertexAssociator>(
0058           iConfig.getUntrackedParameter<edm::InputTag>("vertexAssociator"))) {
0059   reco_vertex_collections_ = iConfig.getParameter<std::vector<edm::InputTag>>("vertexRecoCollections");
0060   for (auto const& l : reco_vertex_collections_) {
0061     reco_vertex_collection_tokens_.push_back(
0062         edm::EDGetTokenT<edm::View<reco::Vertex>>(consumes<edm::View<reco::Vertex>>(l)));
0063   }
0064   errorPrintedForColl_.resize(reco_vertex_collections_.size(), false);
0065 }
0066 
0067 PrimaryVertexAnalyzer4PUSlimmed::~PrimaryVertexAnalyzer4PUSlimmed() {}
0068 
0069 //
0070 // member functions
0071 //
0072 void PrimaryVertexAnalyzer4PUSlimmed::bookHistograms(DQMStore::IBooker& i,
0073                                                      edm::Run const& iRun,
0074                                                      edm::EventSetup const& iSetup) {
0075   // TODO(rovere) make this booking method shorter and smarter,
0076   // factorizing similar histograms with different prefix in a single
0077   // method call.
0078   float log_bins[31] = {0.0,  0.0002, 0.0004, 0.0006, 0.0008, 0.001, 0.002, 0.004, 0.006, 0.008, 0.01,
0079                         0.02, 0.04,   0.06,   0.08,   0.1,    0.2,   0.4,   0.6,   0.8,   1.0,   2.0,
0080                         4.0,  6.0,    8.0,    10.0,   20.0,   40.0,  60.0,  80.0,  100.0};
0081   auto const log_mergez_bins = makeLogBins<float, 16>(1e-4, 1);  // 4 bins / 10x
0082 
0083   auto const log_pt_bins = makeLogBins<float, 100>(0.1, 1e4);
0084   float log_pt2_bins[16] = {
0085       0.0, 0.1, 0.5, 1.0, 2.0, 5.0, 10.0, 20.0, 50.0, 100.0, 200.0, 500.0, 1000.0, 2000.0, 5000.0, 10000.0};
0086   float log_ntrk_bins[25] = {0.,   2.0,  4.0,  6.0,  8.0,  10.,  12.0, 14.0, 16.0, 18.0,  22.0,  26.0, 30.0,
0087                              35.0, 40.0, 45.0, 50.0, 55.0, 60.0, 70.0, 80.0, 90.0, 100.0, 150.0, 200.0};
0088 
0089   // TODO(rovere) Possibly change or add the main DQMStore booking
0090   // interface to allow booking a TProfile with variable bin-width
0091   // using an array of floats, as done for the TH1F case, not of
0092   // doubles.
0093   double log_pt2_bins_double[16] = {
0094       0.0, 0.1, 0.5, 1.0, 2.0, 5.0, 10.0, 20.0, 50.0, 100.0, 200.0, 500.0, 1000.0, 2000.0, 5000.0, 10000.0};
0095 
0096   i.setCurrentFolder(root_folder_);
0097   if (do_generic_sim_plots_) {
0098     mes_["root_folder"]["GenVtx_vs_BX"] = i.book2D("GenVtx_vs_BX", "GenVtx_vs_BX", 16, -12.5, 3.5, 200, 0., 200.);
0099     // Generated Primary Vertex Plots
0100     mes_["root_folder"]["GenPV_X"] = i.book1D("GenPV_X", "GeneratedPV_X", 120, -0.6, 0.6);
0101     mes_["root_folder"]["GenPV_Y"] = i.book1D("GenPV_Y", "GeneratedPV_Y", 120, -0.6, 0.6);
0102     mes_["root_folder"]["GenPV_Z"] = i.book1D("GenPV_Z", "GeneratedPV_Z", 120, -60., 60.);
0103     mes_["root_folder"]["GenPV_R"] = i.book1D("GenPV_R", "GeneratedPV_R", 120, 0, 0.6);
0104     mes_["root_folder"]["GenPV_Pt2"] = i.book1D("GenPV_Pt2", "GeneratedPV_Sum-pt2", 15, &log_pt2_bins[0]);
0105     mes_["root_folder"]["GenPV_NumTracks"] =
0106         i.book1D("GenPV_NumTracks", "GeneratedPV_NumTracks", 24, &log_ntrk_bins[0]);
0107     mes_["root_folder"]["GenPV_ClosestDistanceZ"] =
0108         i.book1D("GenPV_ClosestDistanceZ", "GeneratedPV_ClosestDistanceZ", 30, &log_bins[0]);
0109 
0110     // All Generated Vertices, used for efficiency plots
0111     mes_["root_folder"]["GenAllV_NumVertices"] =
0112         i.book1D("GenAllV_NumVertices", "GeneratedAllV_NumVertices", 100, 0., 200.);
0113     mes_["root_folder"]["GenAllV_X"] = i.book1D("GenAllV_X", "GeneratedAllV_X", 120, -0.6, 0.6);
0114     mes_["root_folder"]["GenAllV_Y"] = i.book1D("GenAllV_Y", "GeneratedAllV_Y", 120, -0.6, 0.6);
0115     mes_["root_folder"]["GenAllV_Z"] = i.book1D("GenAllV_Z", "GeneratedAllV_Z", 120, -60, 60);
0116     mes_["root_folder"]["GenAllV_R"] = i.book1D("GenAllV_R", "GeneratedAllV_R", 120, 0, 0.6);
0117     mes_["root_folder"]["GenAllV_Pt2"] = i.book1D("GenAllV_Pt2", "GeneratedAllV_Sum-pt2", 15, &log_pt2_bins[0]);
0118     mes_["root_folder"]["GenAllV_NumTracks"] =
0119         i.book1D("GenAllV_NumTracks", "GeneratedAllV_NumTracks", 24, &log_ntrk_bins[0]);
0120     mes_["root_folder"]["GenAllV_ClosestDistanceZ"] =
0121         i.book1D("GenAllV_ClosestDistanceZ", "GeneratedAllV_ClosestDistanceZ", 30, &log_bins[0]);
0122     mes_["root_folder"]["GenAllV_PairDistanceZ"] =
0123         i.book1D("GenAllV_PairDistanceZ", "GeneratedAllV_PairDistanceZ", 1000, 0, 20);
0124     mes_["root_folder"]["SignalIsHighestPt2"] = i.book1D("SignalIsHighestPt2", "SignalIsHighestPt2", 2, -0.5, 1.5);
0125   }
0126 
0127   for (auto const& l : reco_vertex_collections_) {
0128     std::string label = l.label();
0129     std::string current_folder = root_folder_ + "/" + label;
0130     i.setCurrentFolder(current_folder);
0131 
0132     auto book1d = [&](const char* name, int bins, double min, double max) {
0133       mes_[label][name] = i.book1D(name, name, bins, min, max);
0134     };
0135     auto book1dlogx = [&](const char* name, int bins, const float* xbinedges) {
0136       mes_[label][name] = i.book1D(name, name, bins, xbinedges);
0137     };
0138     auto book2d = [&](const char* name, int xbins, double xmin, double xmax, int ybins, double ymin, double ymax) {
0139       mes_[label][name] = i.book2D(name, name, xbins, xmin, xmax, ybins, ymin, ymax);
0140     };
0141     auto book2dlogx = [&](const char* name, int xbins, const float* xbinedges, int ybins, double ymin, double ymax) {
0142       auto me = i.book2D(name, name, xbins, xbinedges[0], xbinedges[xbins], ybins, ymin, ymax);
0143       me->getTH2F()->GetXaxis()->Set(xbins, xbinedges);
0144       mes_[label][name] = me;
0145     };
0146 
0147     mes_[label]["RecoVtx_vs_GenVtx"] =
0148         i.bookProfile("RecoVtx_vs_GenVtx", "RecoVtx_vs_GenVtx", 125, 0., 250., 250, 0., 250.);
0149     mes_[label]["MatchedRecoVtx_vs_GenVtx"] =
0150         i.bookProfile("MatchedRecoVtx_vs_GenVtx", "MatchedRecoVtx_vs_GenVtx", 125, 0., 250., 250, 0., 250.);
0151     mes_[label]["KindOfSignalPV"] = i.book1D("KindOfSignalPV", "KindOfSignalPV", 9, -0.5, 8.5);
0152     mes_[label]["KindOfSignalPV"]->getTH1()->GetXaxis()->SetBinLabel(1, "!Highest!Assoc2Any");
0153     mes_[label]["KindOfSignalPV"]->getTH1()->GetXaxis()->SetBinLabel(2, "Highest!Assoc2Any");
0154     mes_[label]["KindOfSignalPV"]->getTH1()->GetXaxis()->SetBinLabel(3, "!HighestAssoc2First");
0155     mes_[label]["KindOfSignalPV"]->getTH1()->GetXaxis()->SetBinLabel(4, "HighestAssoc2First");
0156     mes_[label]["KindOfSignalPV"]->getTH1()->GetXaxis()->SetBinLabel(5, "!HighestAssoc2!First");
0157     mes_[label]["KindOfSignalPV"]->getTH1()->GetXaxis()->SetBinLabel(6, "HighestAssoc2!First");
0158     mes_[label]["KindOfSignalPV"]->getTH1()->GetXaxis()->SetBinLabel(7, "!HighestAssoc2First");
0159     mes_[label]["KindOfSignalPV"]->getTH1()->GetXaxis()->SetBinLabel(8, "HighestAssoc2First");
0160     mes_[label]["MisTagRate"] = i.book1D("MisTagRate", "MisTagRate", 2, -0.5, 1.5);
0161     mes_[label]["MisTagRate_vs_PU"] = i.bookProfile("MisTagRate_vs_PU", "MisTagRate_vs_PU", 125, 0., 250., 2, 0., 1.);
0162     mes_[label]["MisTagRate_vs_sum-pt2"] =
0163         i.bookProfile("MisTagRate_vs_sum-pt2", "MisTagRate_vs_sum-pt2", 15, &log_pt2_bins_double[0], 2, 0., 1.);
0164     mes_[label]["MisTagRate_vs_Z"] = i.bookProfile("MisTagRate_vs_Z", "MisTagRate_vs_Z", 120, -60., 60., 2, 0., 1.);
0165     mes_[label]["MisTagRate_vs_R"] = i.bookProfile("MisTagRate_vs_R", "MisTagRate_vs_R", 120, 0., 0.6, 2, 0., 1.);
0166     mes_[label]["MisTagRate_vs_NumTracks"] =
0167         i.bookProfile("MisTagRate_vs_NumTracks", "MisTagRate_vs_NumTracks", 100, 0., 200, 2, 0., 1.);
0168     mes_[label]["MisTagRateSignalIsHighest"] =
0169         i.book1D("MisTagRateSignalIsHighest", "MisTagRateSignalIsHighest", 2, -0.5, 1.5);
0170     mes_[label]["MisTagRateSignalIsHighest_vs_PU"] =
0171         i.bookProfile("MisTagRateSignalIsHighest_vs_PU", "MisTagRateSignalIsHighest_vs_PU", 125, 0., 250., 2, 0., 1.);
0172     mes_[label]["MisTagRateSignalIsHighest_vs_sum-pt2"] = i.bookProfile("MisTagRateSignalIsHighest_vs_sum-pt2",
0173                                                                         "MisTagRateSignalIsHighest_vs_sum-pt2",
0174                                                                         15,
0175                                                                         &log_pt2_bins_double[0],
0176                                                                         2,
0177                                                                         0.,
0178                                                                         1.);
0179     mes_[label]["MisTagRateSignalIsHighest_vs_Z"] =
0180         i.bookProfile("MisTagRateSignalIsHighest_vs_Z", "MisTagRateSignalIsHighest_vs_Z", 120, -60., 60., 2, 0., 1.);
0181     mes_[label]["MisTagRateSignalIsHighest_vs_R"] =
0182         i.bookProfile("MisTagRateSignalIsHighest_vs_R", "MisTagRateSignalIsHighest_vs_R", 120, 0., 0.6, 2, 0., 1.);
0183     mes_[label]["MisTagRateSignalIsHighest_vs_NumTracks"] = i.bookProfile(
0184         "MisTagRateSignalIsHighest_vs_NumTracks", "MisTagRateSignalIsHighest_vs_NumTracks", 100, 0., 200, 2, 0., 1.);
0185     mes_[label]["MisTagRateSignalIsNotHighest"] =
0186         i.book1D("MisTagRateSignalIsNotHighest", "MisTagRateSignalIsNotHighest", 2, -0.5, 1.5);
0187     mes_[label]["MisTagRateSignalIsNotHighest_vs_PU"] = i.bookProfile(
0188         "MisTagRateSignalIsNotHighest_vs_PU", "MisTagRateSignalIsNotHighest_vs_PU", 125, 0., 250., 2, 0., 1.);
0189     mes_[label]["MisTagRateSignalIsNotHighest_vs_sum-pt2"] = i.bookProfile("MisTagRateSignalIsNotHighest_vs_sum-pt2",
0190                                                                            "MisTagRateSignalIsNotHighest_vs_sum-pt2",
0191                                                                            15,
0192                                                                            &log_pt2_bins_double[0],
0193                                                                            2,
0194                                                                            0.,
0195                                                                            1.);
0196     mes_[label]["MisTagRateSignalIsNotHighest_vs_Z"] = i.bookProfile(
0197         "MisTagRateSignalIsNotHighest_vs_Z", "MisTagRateSignalIsNotHighest_vs_Z", 120, -60., 60., 2, 0., 1.);
0198     mes_[label]["MisTagRateSignalIsNotHighest_vs_R"] = i.bookProfile(
0199         "MisTagRateSignalIsNotHighest_vs_R", "MisTagRateSignalIsNotHighest_vs_R", 120, 0., 0.6, 2, 0., 1.);
0200     mes_[label]["MisTagRateSignalIsNotHighest_vs_NumTracks"] =
0201         i.bookProfile("MisTagRateSignalIsNotHighest_vs_NumTracks",
0202                       "MisTagRateSignalIsNotHighest_vs_NumTracks",
0203                       100,
0204                       0.,
0205                       200,
0206                       2,
0207                       0.,
0208                       1.);
0209     mes_[label]["TruePVLocationIndex"] =
0210         i.book1D("TruePVLocationIndex", "TruePVLocationIndexInRecoVertexCollection", 12, -1.5, 10.5);
0211     mes_[label]["TruePVLocationIndexCumulative"] =
0212         i.book1D("TruePVLocationIndexCumulative", "TruePVLocationIndexInRecoVertexCollectionCumulative", 3, -1.5, 1.5);
0213     mes_[label]["TruePVLocationIndexSignalIsHighest"] =
0214         i.book1D("TruePVLocationIndexSignalIsHighest",
0215                  "TruePVLocationIndexSignalIsHighestInRecoVertexCollection",
0216                  12,
0217                  -1.5,
0218                  10.5);
0219     mes_[label]["TruePVLocationIndexSignalIsNotHighest"] =
0220         i.book1D("TruePVLocationIndexSignalIsNotHighest",
0221                  "TruePVLocationIndexSignalIsNotHighestInRecoVertexCollection",
0222                  12,
0223                  -1.5,
0224                  10.5);
0225     // All Generated Vertices. Used for Efficiency plots We kind of
0226     // duplicate plots here in case we want to perform more detailed
0227     // studies on a selection of generated vertices, not on all of them.
0228     mes_[label]["GenAllAssoc2Reco_NumVertices"] =
0229         i.book1D("GenAllAssoc2Reco_NumVertices", "GeneratedAllAssoc2Reco_NumVertices", 100, 0., 200.);
0230     mes_[label]["GenAllAssoc2Reco_X"] = i.book1D("GenAllAssoc2Reco_X", "GeneratedAllAssoc2Reco_X", 120, -0.6, 0.6);
0231     mes_[label]["GenAllAssoc2Reco_Y"] = i.book1D("GenAllAssoc2Reco_Y", "GeneratedAllAssoc2Reco_Y", 120, -0.6, 0.6);
0232     mes_[label]["GenAllAssoc2Reco_Z"] = i.book1D("GenAllAssoc2Reco_Z", "GeneratedAllAssoc2Reco_Z", 120, -60, 60);
0233     mes_[label]["GenAllAssoc2Reco_R"] = i.book1D("GenAllAssoc2Reco_R", "GeneratedAllAssoc2Reco_R", 120, 0, 0.6);
0234     mes_[label]["GenAllAssoc2Reco_Pt2"] =
0235         i.book1D("GenAllAssoc2Reco_Pt2", "GeneratedAllAssoc2Reco_Sum-pt2", 15, &log_pt2_bins[0]);
0236     mes_[label]["GenAllAssoc2Reco_NumTracks"] =
0237         i.book1D("GenAllAssoc2Reco_NumTracks", "GeneratedAllAssoc2Reco_NumTracks", 24, &log_ntrk_bins[0]);
0238     mes_[label]["GenAllAssoc2Reco_ClosestDistanceZ"] =
0239         i.book1D("GenAllAssoc2Reco_ClosestDistanceZ", "GeneratedAllAssoc2Reco_ClosestDistanceZ", 30, &log_bins[0]);
0240     book1d("GenPVAssoc2RecoPV_X", 120, -0.6, 0.6);
0241     book1d("GenPVAssoc2RecoPV_Y", 120, -0.6, 0.6);
0242     book1d("GenPVAssoc2RecoPV_Z", 120, -60, 60);
0243 
0244     // All Generated Vertices Matched to a Reconstructed vertex. Used
0245     // for Efficiency plots
0246     mes_[label]["GenAllAssoc2RecoMatched_NumVertices"] =
0247         i.book1D("GenAllAssoc2RecoMatched_NumVertices", "GeneratedAllAssoc2RecoMatched_NumVertices", 100, 0., 200.);
0248     mes_[label]["GenAllAssoc2RecoMatched_X"] =
0249         i.book1D("GenAllAssoc2RecoMatched_X", "GeneratedAllAssoc2RecoMatched_X", 120, -0.6, 0.6);
0250     mes_[label]["GenAllAssoc2RecoMatched_Y"] =
0251         i.book1D("GenAllAssoc2RecoMatched_Y", "GeneratedAllAssoc2RecoMatched_Y", 120, -0.6, 0.6);
0252     mes_[label]["GenAllAssoc2RecoMatched_Z"] =
0253         i.book1D("GenAllAssoc2RecoMatched_Z", "GeneratedAllAssoc2RecoMatched_Z", 120, -60, 60);
0254     mes_[label]["GenAllAssoc2RecoMatched_R"] =
0255         i.book1D("GenAllAssoc2RecoMatched_R", "GeneratedAllAssoc2RecoMatched_R", 120, 0, 0.6);
0256     mes_[label]["GenAllAssoc2RecoMatched_Pt2"] =
0257         i.book1D("GenAllAssoc2RecoMatched_Pt2", "GeneratedAllAssoc2RecoMatched_Sum-pt2", 15, &log_pt2_bins[0]);
0258     mes_[label]["GenAllAssoc2RecoMatched_NumTracks"] =
0259         i.book1D("GenAllAssoc2RecoMatched_NumTracks", "GeneratedAllAssoc2RecoMatched_NumTracks", 24, &log_ntrk_bins[0]);
0260     mes_[label]["GenAllAssoc2RecoMatched_ClosestDistanceZ"] = i.book1D(
0261         "GenAllAssoc2RecoMatched_ClosestDistanceZ", "GeneratedAllAssoc2RecoMatched_ClosestDistanceZ", 30, &log_bins[0]);
0262     book1d("GenPVAssoc2RecoPVMatched_X", 120, -0.6, 0.6);
0263     book1d("GenPVAssoc2RecoPVMatched_Y", 120, -0.6, 0.6);
0264     book1d("GenPVAssoc2RecoPVMatched_Z", 120, -60, 60);
0265 
0266     // All Generated Vertices Multi-Matched to a Reconstructed vertex. Used
0267     // for Duplicate rate plots
0268     mes_[label]["GenAllAssoc2RecoMultiMatched_NumVertices"] = i.book1D(
0269         "GenAllAssoc2RecoMultiMatched_NumVertices", "GeneratedAllAssoc2RecoMultiMatched_NumVertices", 100, 0., 200.);
0270     mes_[label]["GenAllAssoc2RecoMultiMatched_X"] =
0271         i.book1D("GenAllAssoc2RecoMultiMatched_X", "GeneratedAllAssoc2RecoMultiMatched_X", 120, -0.6, 0.6);
0272     mes_[label]["GenAllAssoc2RecoMultiMatched_Y"] =
0273         i.book1D("GenAllAssoc2RecoMultiMatched_Y", "GeneratedAllAssoc2RecoMultiMatched_Y", 120, -0.6, 0.6);
0274     mes_[label]["GenAllAssoc2RecoMultiMatched_Z"] =
0275         i.book1D("GenAllAssoc2RecoMultiMatched_Z", "GeneratedAllAssoc2RecoMultiMatched_Z", 120, -60, 60);
0276     mes_[label]["GenAllAssoc2RecoMultiMatched_R"] =
0277         i.book1D("GenAllAssoc2RecoMultiMatched_R", "GeneratedAllAssoc2RecoMultiMatched_R", 120, 0, 0.6);
0278     mes_[label]["GenAllAssoc2RecoMultiMatched_Pt2"] = i.book1D(
0279         "GenAllAssoc2RecoMultiMatched_Pt2", "GeneratedAllAssoc2RecoMultiMatched_Sum-pt2", 15, &log_pt2_bins[0]);
0280     mes_[label]["GenAllAssoc2RecoMultiMatched_NumTracks"] = i.book1D("GenAllAssoc2RecoMultiMatched_NumTracks",
0281                                                                      "GeneratedAllAssoc2RecoMultiMatched_NumTracks",
0282                                                                      24,
0283                                                                      &log_ntrk_bins[0]);
0284     mes_[label]["GenAllAssoc2RecoMultiMatched_ClosestDistanceZ"] =
0285         i.book1D("GenAllAssoc2RecoMultiMatched_ClosestDistanceZ",
0286                  "GeneratedAllAssoc2RecoMultiMatched_ClosestDistanceZ",
0287                  30,
0288                  &log_bins[0]);
0289 
0290     // All Reco Vertices. Used for {Fake,Duplicate}-Rate plots
0291     mes_[label]["RecoAllAssoc2Gen_NumVertices"] =
0292         i.book1D("RecoAllAssoc2Gen_NumVertices", "ReconstructedAllAssoc2Gen_NumVertices", 100, 0., 200.);
0293     mes_[label]["RecoAllAssoc2Gen_X"] = i.book1D("RecoAllAssoc2Gen_X", "ReconstructedAllAssoc2Gen_X", 120, -0.6, 0.6);
0294     mes_[label]["RecoAllAssoc2Gen_Y"] = i.book1D("RecoAllAssoc2Gen_Y", "ReconstructedAllAssoc2Gen_Y", 120, -0.6, 0.6);
0295     mes_[label]["RecoAllAssoc2Gen_Z"] = i.book1D("RecoAllAssoc2Gen_Z", "ReconstructedAllAssoc2Gen_Z", 120, -60, 60);
0296     mes_[label]["RecoAllAssoc2Gen_R"] = i.book1D("RecoAllAssoc2Gen_R", "ReconstructedAllAssoc2Gen_R", 120, 0, 0.6);
0297     mes_[label]["RecoAllAssoc2Gen_Pt2"] =
0298         i.book1D("RecoAllAssoc2Gen_Pt2", "ReconstructedAllAssoc2Gen_Sum-pt2", 15, &log_pt2_bins[0]);
0299     mes_[label]["RecoAllAssoc2Gen_Ndof"] =
0300         i.book1D("RecoAllAssoc2Gen_Ndof", "ReconstructedAllAssoc2Gen_Ndof", 120, 0., 240.);
0301     mes_[label]["RecoAllAssoc2Gen_NumTracks"] =
0302         i.book1D("RecoAllAssoc2Gen_NumTracks", "ReconstructedAllAssoc2Gen_NumTracks", 24, &log_ntrk_bins[0]);
0303     mes_[label]["RecoAllAssoc2Gen_PU"] = i.book1D("RecoAllAssoc2Gen_PU", "ReconstructedAllAssoc2Gen_PU", 125, 0., 250.);
0304     mes_[label]["RecoAllAssoc2Gen_ClosestDistanceZ"] =
0305         i.book1D("RecoAllAssoc2Gen_ClosestDistanceZ", "ReconstructedAllAssoc2Gen_ClosestDistanceZ", 30, &log_bins[0]);
0306     mes_[label]["RecoAllAssoc2GenProperties"] =
0307         i.book1D("RecoAllAssoc2GenProperties", "ReconstructedAllAssoc2Gen_Properties", 8, -0.5, 7.5);
0308     mes_[label]["RecoAllAssoc2Gen_PairDistanceZ"] =
0309         i.book1D("RecoAllAssoc2Gen_PairDistanceZ", "RecoAllAssoc2Gen_PairDistanceZ", 1000, 0, 20);
0310 
0311     // All Reconstructed Vertices Matched to a Generated vertex. Used
0312     // for Fake-Rate plots
0313     mes_[label]["RecoAllAssoc2GenMatched_NumVertices"] =
0314         i.book1D("RecoAllAssoc2GenMatched_NumVertices", "ReconstructedAllAssoc2GenMatched_NumVertices", 100, 0., 200.);
0315     mes_[label]["RecoAllAssoc2GenMatched_X"] =
0316         i.book1D("RecoAllAssoc2GenMatched_X", "ReconstructedAllAssoc2GenMatched_X", 120, -0.6, 0.6);
0317     mes_[label]["RecoAllAssoc2GenMatched_Y"] =
0318         i.book1D("RecoAllAssoc2GenMatched_Y", "ReconstructedAllAssoc2GenMatched_Y", 120, -0.6, 0.6);
0319     mes_[label]["RecoAllAssoc2GenMatched_Z"] =
0320         i.book1D("RecoAllAssoc2GenMatched_Z", "ReconstructedAllAssoc2GenMatched_Z", 120, -60, 60);
0321     mes_[label]["RecoAllAssoc2GenMatched_R"] =
0322         i.book1D("RecoAllAssoc2GenMatched_R", "ReconstructedAllAssoc2GenMatched_R", 120, 0, 0.6);
0323     mes_[label]["RecoAllAssoc2GenMatched_Pt2"] =
0324         i.book1D("RecoAllAssoc2GenMatched_Pt2", "ReconstructedAllAssoc2GenMatched_Sum-pt2", 15, &log_pt2_bins[0]);
0325     mes_[label]["RecoAllAssoc2GenMatched_Ndof"] =
0326         i.book1D("RecoAllAssoc2GenMatched_Ndof", "ReconstructedAllAssoc2GenMatched_Ndof", 120, 0., 240.);
0327     mes_[label]["RecoAllAssoc2GenMatched_NumTracks"] = i.book1D(
0328         "RecoAllAssoc2GenMatched_NumTracks", "ReconstructedAllAssoc2GenMatched_NumTracks", 24, &log_ntrk_bins[0]);
0329     mes_[label]["RecoAllAssoc2GenMatched_PU"] =
0330         i.book1D("RecoAllAssoc2GenMatched_PU", "ReconstructedAllAssoc2GenMatched_PU", 125, 0., 250.);
0331     mes_[label]["RecoAllAssoc2GenMatched_ClosestDistanceZ"] =
0332         i.book1D("RecoAllAssoc2GenMatched_ClosestDistanceZ",
0333                  "ReconstructedAllAssoc2GenMatched_ClosestDistanceZ",
0334                  30,
0335                  &log_bins[0]);
0336 
0337     // All Reconstructed Vertices  Multi-Matched to a Generated vertex. Used
0338     // for Merge-Rate plots
0339     mes_[label]["RecoAllAssoc2GenMultiMatched_NumVertices"] = i.book1D(
0340         "RecoAllAssoc2GenMultiMatched_NumVertices", "ReconstructedAllAssoc2GenMultiMatched_NumVertices", 100, 0., 200.);
0341     mes_[label]["RecoAllAssoc2GenMultiMatched_X"] =
0342         i.book1D("RecoAllAssoc2GenMultiMatched_X", "ReconstructedAllAssoc2GenMultiMatched_X", 120, -0.6, 0.6);
0343     mes_[label]["RecoAllAssoc2GenMultiMatched_Y"] =
0344         i.book1D("RecoAllAssoc2GenMultiMatched_Y", "ReconstructedAllAssoc2GenMultiMatched_Y", 120, -0.6, 0.6);
0345     mes_[label]["RecoAllAssoc2GenMultiMatched_Z"] =
0346         i.book1D("RecoAllAssoc2GenMultiMatched_Z", "ReconstructedAllAssoc2GenMultiMatched_Z", 120, -60, 60);
0347     mes_[label]["RecoAllAssoc2GenMultiMatched_R"] =
0348         i.book1D("RecoAllAssoc2GenMultiMatched_R", "ReconstructedAllAssoc2GenMultiMatched_R", 120, 0, 0.6);
0349     mes_[label]["RecoAllAssoc2GenMultiMatched_Pt2"] = i.book1D(
0350         "RecoAllAssoc2GenMultiMatched_Pt2", "ReconstructedAllAssoc2GenMultiMatched_Sum-pt2", 15, &log_pt2_bins[0]);
0351     mes_[label]["RecoAllAssoc2GenMultiMatched_NumTracks"] = i.book1D("RecoAllAssoc2GenMultiMatched_NumTracks",
0352                                                                      "ReconstructedAllAssoc2GenMultiMatched_NumTracks",
0353                                                                      24,
0354                                                                      &log_ntrk_bins[0]);
0355     mes_[label]["RecoAllAssoc2GenMultiMatched_PU"] =
0356         i.book1D("RecoAllAssoc2GenMultiMatched_PU", "ReconstructedAllAssoc2GenMultiMatched_PU", 125, 0., 250.);
0357     mes_[label]["RecoAllAssoc2GenMultiMatched_ClosestDistanceZ"] =
0358         i.book1D("RecoAllAssoc2GenMultiMatched_ClosestDistanceZ",
0359                  "ReconstructedAllAssoc2GenMultiMatched_ClosestDistanceZ",
0360                  log_mergez_bins.size() - 1,
0361                  &log_mergez_bins[0]);
0362 
0363     // All Reconstructed Vertices Matched to a Multi-Matched Gen
0364     // Vertex. Used for Duplicate rate plots done w.r.t. Reco
0365     // Quantities. We basically want to ask how many times a RecoVTX
0366     // has been reconstructed and associated to a SimulatedVTX that
0367     // has been linked to at least another RecoVTX. In this sense this
0368     // RecoVTX is a duplicate of the same, real GenVTX.
0369     mes_[label]["RecoAllAssoc2MultiMatchedGen_NumVertices"] =
0370         i.book1D("RecoAllAssoc2MultiMatchedGen_NumVertices", "RecoAllAssoc2MultiMatchedGen_NumVertices", 100, 0., 200.);
0371     mes_[label]["RecoAllAssoc2MultiMatchedGen_X"] =
0372         i.book1D("RecoAllAssoc2MultiMatchedGen_X", "RecoAllAssoc2MultiMatchedGen_X", 120, -0.6, 0.6);
0373     mes_[label]["RecoAllAssoc2MultiMatchedGen_Y"] =
0374         i.book1D("RecoAllAssoc2MultiMatchedGen_Y", "RecoAllAssoc2MultiMatchedGen_Y", 120, -0.6, 0.6);
0375     mes_[label]["RecoAllAssoc2MultiMatchedGen_Z"] =
0376         i.book1D("RecoAllAssoc2MultiMatchedGen_Z", "RecoAllAssoc2MultiMatchedGen_Z", 120, -60, 60);
0377     mes_[label]["RecoAllAssoc2MultiMatchedGen_R"] =
0378         i.book1D("RecoAllAssoc2MultiMatchedGen_R", "RecoAllAssoc2MultiMatchedGen_R", 120, 0, 0.6);
0379     mes_[label]["RecoAllAssoc2MultiMatchedGen_Pt2"] =
0380         i.book1D("RecoAllAssoc2MultiMatchedGen_Pt2", "RecoAllAssoc2MultiMatchedGen_Sum-pt2", 15, &log_pt2_bins[0]);
0381     mes_[label]["RecoAllAssoc2MultiMatchedGen_NumTracks"] = i.book1D(
0382         "RecoAllAssoc2MultiMatchedGen_NumTracks", "RecoAllAssoc2MultiMatchedGen_NumTracks", 24, &log_ntrk_bins[0]);
0383     mes_[label]["RecoAllAssoc2MultiMatchedGen_PU"] =
0384         i.book1D("RecoAllAssoc2MultiMatchedGen_PU", "RecoAllAssoc2MultiMatchedGen_PU", 125, 0., 250.);
0385     mes_[label]["RecoAllAssoc2MultiMatchedGen_ClosestDistanceZ"] =
0386         i.book1D("RecoAllAssoc2MultiMatchedGen_ClosestDistanceZ",
0387                  "RecoAllAssoc2MultiMatchedGen_ClosestDistanceZ",
0388                  30,
0389                  &log_bins[0]);
0390     mes_[label]["RecoAllAssoc2GenSimForMerge_ClosestDistanceZ"] =
0391         i.book1D("RecoAllAssoc2GenSimForMerge_ClosestDistanceZ",
0392                  "RecoAllAssoc2GenSimForMerge_ClosestDistanceZ",
0393                  log_mergez_bins.size() - 1,
0394                  &log_mergez_bins[0]);
0395 
0396     // Resolution and pull histograms
0397     const double resolpt2 = 10;
0398 
0399     const double minPull = -10;
0400     const double maxPull = 10;
0401     const double nPull = 100;
0402 
0403     auto bookResolPull = [&](const std::string& prefix, const double resolx, const double resoly, const double resolz) {
0404       book1d((prefix + "_ResolX").c_str(), 100, -resolx, resolx);
0405       book1d((prefix + "_ResolY").c_str(), 100, -resoly, resoly);
0406       book1d((prefix + "_ResolZ").c_str(), 100, -resolz, resolz);
0407       book1d((prefix + "_ResolPt2").c_str(), 100, -resolpt2, resolpt2);
0408       book1d((prefix + "_PullX").c_str(), 250, -25, 25);
0409       book1d((prefix + "_PullY").c_str(), 250, -25, 25);
0410       book1d((prefix + "_PullZ").c_str(), 250, -25, 25);
0411 
0412       book2d((prefix + "_ResolX_vs_PU").c_str(), 125, 0., 250., 100, -resolx, resolx);
0413       book2d((prefix + "_ResolY_vs_PU").c_str(), 125, 0., 250., 100, -resoly, resoly);
0414       book2d((prefix + "_ResolZ_vs_PU").c_str(), 125, 0., 250., 100, -resolz, resolz);
0415       book2d((prefix + "_ResolPt2_vs_PU").c_str(), 125, 0., 250., 100, -resolpt2, resolpt2);
0416       book2d((prefix + "_PullX_vs_PU").c_str(), 125, 0., 250., nPull, minPull, maxPull);
0417       book2d((prefix + "_PullY_vs_PU").c_str(), 125, 0., 250., nPull, minPull, maxPull);
0418       book2d((prefix + "_PullZ_vs_PU").c_str(), 125, 0., 250., nPull, minPull, maxPull);
0419 
0420       book2dlogx((prefix + "_ResolX_vs_NumTracks").c_str(), 24, &log_ntrk_bins[0], 100, -resolx, resolx);
0421       book2dlogx((prefix + "_ResolY_vs_NumTracks").c_str(), 24, &log_ntrk_bins[0], 100, -resoly, resoly);
0422       book2dlogx((prefix + "_ResolZ_vs_NumTracks").c_str(), 24, &log_ntrk_bins[0], 100, -resolz, resolz);
0423       book2dlogx((prefix + "_ResolPt2_vs_NumTracks").c_str(), 24, &log_ntrk_bins[0], 100, -resolpt2, resolpt2);
0424       book2dlogx((prefix + "_PullX_vs_NumTracks").c_str(), 24, &log_ntrk_bins[0], nPull, minPull, maxPull);
0425       book2dlogx((prefix + "_PullY_vs_NumTracks").c_str(), 24, &log_ntrk_bins[0], nPull, minPull, maxPull);
0426       book2dlogx((prefix + "_PullZ_vs_NumTracks").c_str(), 24, &log_ntrk_bins[0], nPull, minPull, maxPull);
0427 
0428       book2d((prefix + "_ResolX_vs_Z").c_str(), 120, -60, 60, 100, -resolx, resolx);
0429       book2d((prefix + "_ResolY_vs_Z").c_str(), 120, -60, 60, 100, -resoly, resoly);
0430       book2d((prefix + "_ResolZ_vs_Z").c_str(), 120, -60, 60, 100, -resolz, resolz);
0431       book2d((prefix + "_ResolPt2_vs_Z").c_str(), 120, -60, 60, 100, -resolpt2, resolpt2);
0432       book2d((prefix + "_PullX_vs_Z").c_str(), 120, -60, 60, nPull, minPull, maxPull);
0433       book2d((prefix + "_PullY_vs_Z").c_str(), 120, -60, 60, nPull, minPull, maxPull);
0434       book2d((prefix + "_PullZ_vs_Z").c_str(), 120, -60, 60, nPull, minPull, maxPull);
0435 
0436       book2dlogx((prefix + "_ResolX_vs_Pt").c_str(), log_pt_bins.size() - 1, &log_pt_bins[0], 100, -resolx, resolx);
0437       book2dlogx((prefix + "_ResolY_vs_Pt").c_str(), log_pt_bins.size() - 1, &log_pt_bins[0], 100, -resoly, resoly);
0438       book2dlogx((prefix + "_ResolZ_vs_Pt").c_str(), log_pt_bins.size() - 1, &log_pt_bins[0], 100, -resolz, resolz);
0439       book2dlogx(
0440           (prefix + "_ResolPt2_vs_Pt").c_str(), log_pt_bins.size() - 1, &log_pt_bins[0], 100, -resolpt2, resolpt2);
0441       book2dlogx((prefix + "_PullX_vs_Pt").c_str(), log_pt_bins.size() - 1, &log_pt_bins[0], nPull, minPull, maxPull);
0442       book2dlogx((prefix + "_PullY_vs_Pt").c_str(), log_pt_bins.size() - 1, &log_pt_bins[0], nPull, minPull, maxPull);
0443       book2dlogx((prefix + "_PullZ_vs_Pt").c_str(), log_pt_bins.size() - 1, &log_pt_bins[0], nPull, minPull, maxPull);
0444     };
0445 
0446     bookResolPull("RecoAllAssoc2GenMatched", 0.1, 0.1, 0.1);        // Non-merged vertices
0447     bookResolPull("RecoAllAssoc2GenMatchedMerged", 0.1, 0.1, 0.1);  // Merged vertices
0448     bookResolPull(
0449         "RecoPVAssoc2GenPVMatched", 0.01, 0.01, 0.01);  // PV, when correctly matched, regardless if merged or not
0450 
0451     // Purity histograms
0452     // Reco PV (vtx0) matched to hard-scatter gen vertex
0453     book1d("RecoPVAssoc2GenPVMatched_Purity", 50, 0, 1);
0454     book1d("RecoPVAssoc2GenPVMatched_Missing", 50, 0, 1);
0455     book2d("RecoPVAssoc2GenPVMatched_Purity_vs_Index", 100, 0, 100, 50, 0, 1);
0456 
0457     // RECO PV (vtx0) not matched to hard-scatter gen vertex
0458     book1d("RecoPVAssoc2GenPVNotMatched_Purity", 50, 0, 1);
0459     book1d("RecoPVAssoc2GenPVNotMatched_Missing", 50, 0, 1);
0460     book2d("RecoPVAssoc2GenPVNotMatched_Purity_vs_Index", 100, 0, 100, 50, 0, 1);
0461 
0462     // Purity vs. fake rate
0463     book1d("RecoAllAssoc2Gen_Purity", 50, 0, 1);         // denominator
0464     book1d("RecoAllAssoc2GenMatched_Purity", 50, 0, 1);  // 1-numerator
0465 
0466     // Vertex sum(pt2)
0467     // The first two are orthogonal (i.e. their sum includes all reco vertices)
0468     book1dlogx("RecoAssoc2GenPVMatched_Pt2", 15, &log_pt2_bins[0]);
0469     book1dlogx("RecoAssoc2GenPVNotMatched_Pt2", 15, &log_pt2_bins[0]);
0470 
0471     book1dlogx("RecoAssoc2GenPVMatchedNotHighest_Pt2", 15, &log_pt2_bins[0]);
0472     book1dlogx("RecoAssoc2GenPVNotMatched_GenPVTracksRemoved_Pt2", 15, &log_pt2_bins[0]);
0473 
0474     // Shared tracks
0475     book1d("RecoAllAssoc2GenSingleMatched_SharedTrackFractionReco", 50, 0, 1);
0476     book1d("RecoAllAssoc2GenMultiMatched_SharedTrackFractionReco", 50, 0, 1);
0477     book1d("RecoAllAssoc2GenSingleMatched_SharedTrackFractionRecoMatched", 50, 0, 1);
0478     book1d("RecoAllAssoc2GenMultiMatched_SharedTrackFractionRecoMatched", 50, 0, 1);
0479 
0480     book1d("RecoAllAssoc2GenSingleMatched_SharedTrackFractionSim", 50, 0, 1);
0481     book1d("RecoAllAssoc2GenMultiMatched_SharedTrackFractionSim", 50, 0, 1);
0482     book1d("RecoAllAssoc2GenSingleMatched_SharedTrackFractionSimMatched", 50, 0, 1);
0483     book1d("RecoAllAssoc2GenMultiMatched_SharedTrackFractionSimMatched", 50, 0, 1);
0484   }
0485 }
0486 
0487 void PrimaryVertexAnalyzer4PUSlimmed::fillGenericGenVertexHistograms(const simPrimaryVertex& v) {
0488   if (v.eventId.event() == 0) {
0489     mes_["root_folder"]["GenPV_X"]->Fill(v.x);
0490     mes_["root_folder"]["GenPV_Y"]->Fill(v.y);
0491     mes_["root_folder"]["GenPV_Z"]->Fill(v.z);
0492     mes_["root_folder"]["GenPV_R"]->Fill(v.r);
0493     mes_["root_folder"]["GenPV_Pt2"]->Fill(v.ptsq);
0494     mes_["root_folder"]["GenPV_NumTracks"]->Fill(v.nGenTrk);
0495     if (v.closest_vertex_distance_z > 0.)
0496       mes_["root_folder"]["GenPV_ClosestDistanceZ"]->Fill(v.closest_vertex_distance_z);
0497   }
0498   mes_["root_folder"]["GenAllV_X"]->Fill(v.x);
0499   mes_["root_folder"]["GenAllV_Y"]->Fill(v.y);
0500   mes_["root_folder"]["GenAllV_Z"]->Fill(v.z);
0501   mes_["root_folder"]["GenAllV_R"]->Fill(v.r);
0502   mes_["root_folder"]["GenAllV_Pt2"]->Fill(v.ptsq);
0503   mes_["root_folder"]["GenAllV_NumTracks"]->Fill(v.nGenTrk);
0504   if (v.closest_vertex_distance_z > 0.)
0505     mes_["root_folder"]["GenAllV_ClosestDistanceZ"]->Fill(v.closest_vertex_distance_z);
0506 }
0507 
0508 void PrimaryVertexAnalyzer4PUSlimmed::fillRecoAssociatedGenVertexHistograms(
0509     const std::string& label, const PrimaryVertexAnalyzer4PUSlimmed::simPrimaryVertex& v) {
0510   mes_[label]["GenAllAssoc2Reco_X"]->Fill(v.x);
0511   mes_[label]["GenAllAssoc2Reco_Y"]->Fill(v.y);
0512   mes_[label]["GenAllAssoc2Reco_Z"]->Fill(v.z);
0513   mes_[label]["GenAllAssoc2Reco_R"]->Fill(v.r);
0514   mes_[label]["GenAllAssoc2Reco_Pt2"]->Fill(v.ptsq);
0515   mes_[label]["GenAllAssoc2Reco_NumTracks"]->Fill(v.nGenTrk);
0516   if (v.closest_vertex_distance_z > 0.)
0517     mes_[label]["GenAllAssoc2Reco_ClosestDistanceZ"]->Fill(v.closest_vertex_distance_z);
0518   if (!v.rec_vertices.empty()) {
0519     mes_[label]["GenAllAssoc2RecoMatched_X"]->Fill(v.x);
0520     mes_[label]["GenAllAssoc2RecoMatched_Y"]->Fill(v.y);
0521     mes_[label]["GenAllAssoc2RecoMatched_Z"]->Fill(v.z);
0522     mes_[label]["GenAllAssoc2RecoMatched_R"]->Fill(v.r);
0523     mes_[label]["GenAllAssoc2RecoMatched_Pt2"]->Fill(v.ptsq);
0524     mes_[label]["GenAllAssoc2RecoMatched_NumTracks"]->Fill(v.nGenTrk);
0525     if (v.closest_vertex_distance_z > 0.)
0526       mes_[label]["GenAllAssoc2RecoMatched_ClosestDistanceZ"]->Fill(v.closest_vertex_distance_z);
0527   }
0528   if (v.rec_vertices.size() > 1) {
0529     mes_[label]["GenAllAssoc2RecoMultiMatched_X"]->Fill(v.x);
0530     mes_[label]["GenAllAssoc2RecoMultiMatched_Y"]->Fill(v.y);
0531     mes_[label]["GenAllAssoc2RecoMultiMatched_Z"]->Fill(v.z);
0532     mes_[label]["GenAllAssoc2RecoMultiMatched_R"]->Fill(v.r);
0533     mes_[label]["GenAllAssoc2RecoMultiMatched_Pt2"]->Fill(v.ptsq);
0534     mes_[label]["GenAllAssoc2RecoMultiMatched_NumTracks"]->Fill(v.nGenTrk);
0535     if (v.closest_vertex_distance_z > 0.)
0536       mes_[label]["GenAllAssoc2RecoMultiMatched_ClosestDistanceZ"]->Fill(v.closest_vertex_distance_z);
0537   }
0538 }
0539 
0540 void PrimaryVertexAnalyzer4PUSlimmed::fillRecoAssociatedGenPVHistograms(
0541     const std::string& label, const PrimaryVertexAnalyzer4PUSlimmed::simPrimaryVertex& v, bool genPVMatchedToRecoPV) {
0542   mes_[label]["GenPVAssoc2RecoPV_X"]->Fill(v.x);
0543   mes_[label]["GenPVAssoc2RecoPV_Y"]->Fill(v.y);
0544   mes_[label]["GenPVAssoc2RecoPV_Z"]->Fill(v.z);
0545   if (genPVMatchedToRecoPV) {
0546     mes_[label]["GenPVAssoc2RecoPVMatched_X"]->Fill(v.x);
0547     mes_[label]["GenPVAssoc2RecoPVMatched_Y"]->Fill(v.y);
0548     mes_[label]["GenPVAssoc2RecoPVMatched_Z"]->Fill(v.z);
0549   }
0550 }
0551 
0552 void PrimaryVertexAnalyzer4PUSlimmed::fillGenAssociatedRecoVertexHistograms(
0553     const std::string& label, int num_pileup_vertices, PrimaryVertexAnalyzer4PUSlimmed::recoPrimaryVertex& v) {
0554   mes_[label]["RecoAllAssoc2Gen_X"]->Fill(v.x);
0555   mes_[label]["RecoAllAssoc2Gen_Y"]->Fill(v.y);
0556   mes_[label]["RecoAllAssoc2Gen_Z"]->Fill(v.z);
0557   mes_[label]["RecoAllAssoc2Gen_R"]->Fill(v.r);
0558   mes_[label]["RecoAllAssoc2Gen_Pt2"]->Fill(v.ptsq);
0559   mes_[label]["RecoAllAssoc2Gen_Ndof"]->Fill(v.recVtx->ndof());
0560   mes_[label]["RecoAllAssoc2Gen_NumTracks"]->Fill(v.nRecoTrk);
0561   mes_[label]["RecoAllAssoc2Gen_PU"]->Fill(num_pileup_vertices);
0562   mes_[label]["RecoAllAssoc2Gen_Purity"]->Fill(v.purity);
0563   if (v.closest_vertex_distance_z > 0.)
0564     mes_[label]["RecoAllAssoc2Gen_ClosestDistanceZ"]->Fill(v.closest_vertex_distance_z);
0565   if (!v.sim_vertices.empty()) {
0566     v.kind_of_vertex |= recoPrimaryVertex::MATCHED;
0567     mes_[label]["RecoAllAssoc2GenMatched_X"]->Fill(v.x);
0568     mes_[label]["RecoAllAssoc2GenMatched_Y"]->Fill(v.y);
0569     mes_[label]["RecoAllAssoc2GenMatched_Z"]->Fill(v.z);
0570     mes_[label]["RecoAllAssoc2GenMatched_R"]->Fill(v.r);
0571     mes_[label]["RecoAllAssoc2GenMatched_Pt2"]->Fill(v.ptsq);
0572     mes_[label]["RecoAllAssoc2GenMatched_Ndof"]->Fill(v.recVtx->ndof());
0573     mes_[label]["RecoAllAssoc2GenMatched_NumTracks"]->Fill(v.nRecoTrk);
0574     mes_[label]["RecoAllAssoc2GenMatched_PU"]->Fill(num_pileup_vertices);
0575     mes_[label]["RecoAllAssoc2GenMatched_Purity"]->Fill(v.purity);
0576     if (v.closest_vertex_distance_z > 0.)
0577       mes_[label]["RecoAllAssoc2GenMatched_ClosestDistanceZ"]->Fill(v.closest_vertex_distance_z);
0578 
0579     // Fill resolution and pull plots here (as in MultiTrackValidator)
0580     fillResolutionAndPullHistograms(label, num_pileup_vertices, v, false);
0581 
0582     // Now keep track of all RecoVTX associated to a SimVTX that
0583     // itself is associated to more than one RecoVTX, for
0584     // duplicate-rate plots on reco quantities.
0585     if (v.sim_vertices_internal[0]->rec_vertices.size() > 1) {
0586       v.kind_of_vertex |= recoPrimaryVertex::DUPLICATE;
0587       mes_[label]["RecoAllAssoc2MultiMatchedGen_X"]->Fill(v.x);
0588       mes_[label]["RecoAllAssoc2MultiMatchedGen_Y"]->Fill(v.y);
0589       mes_[label]["RecoAllAssoc2MultiMatchedGen_Z"]->Fill(v.z);
0590       mes_[label]["RecoAllAssoc2MultiMatchedGen_R"]->Fill(v.r);
0591       mes_[label]["RecoAllAssoc2MultiMatchedGen_Pt2"]->Fill(v.ptsq);
0592       mes_[label]["RecoAllAssoc2MultiMatchedGen_NumTracks"]->Fill(v.nRecoTrk);
0593       mes_[label]["RecoAllAssoc2MultiMatchedGen_PU"]->Fill(num_pileup_vertices);
0594       if (v.closest_vertex_distance_z > 0.)
0595         mes_[label]["RecoAllAssoc2MultiMatchedGen_ClosestDistanceZ"]->Fill(v.closest_vertex_distance_z);
0596     }
0597     // This is meant to be used as "denominator" for the merge-rate
0598     // plots produced starting from reco quantities. We   enter here
0599     // only if the reco vertex has been associated, since we need info
0600     // from the SimVTX associated to it. In this regard, the final
0601     // merge-rate plot coming from reco is not to be intended as a
0602     // pure efficiency-like plot, since the normalization is biased.
0603     if (v.sim_vertices_internal[0]->closest_vertex_distance_z > 0.)
0604       mes_[label]["RecoAllAssoc2GenSimForMerge_ClosestDistanceZ"]->Fill(
0605           v.sim_vertices_internal[0]->closest_vertex_distance_z);
0606   }
0607   // this plots are meant to be used to compute the merge rate
0608   if (v.sim_vertices.size() > 1) {
0609     v.kind_of_vertex |= recoPrimaryVertex::MERGED;
0610     mes_[label]["RecoAllAssoc2GenMultiMatched_X"]->Fill(v.x);
0611     mes_[label]["RecoAllAssoc2GenMultiMatched_Y"]->Fill(v.y);
0612     mes_[label]["RecoAllAssoc2GenMultiMatched_Z"]->Fill(v.z);
0613     mes_[label]["RecoAllAssoc2GenMultiMatched_R"]->Fill(v.r);
0614     mes_[label]["RecoAllAssoc2GenMultiMatched_Pt2"]->Fill(v.ptsq);
0615     mes_[label]["RecoAllAssoc2GenMultiMatched_NumTracks"]->Fill(v.nRecoTrk);
0616     mes_[label]["RecoAllAssoc2GenMultiMatched_PU"]->Fill(num_pileup_vertices);
0617     if (v.sim_vertices_internal[0]->closest_vertex_distance_z > 0.)
0618       mes_[label]["RecoAllAssoc2GenMultiMatched_ClosestDistanceZ"]->Fill(
0619           v.sim_vertices_internal[0]->closest_vertex_distance_z);
0620   }
0621   mes_[label]["RecoAllAssoc2GenProperties"]->Fill(v.kind_of_vertex);
0622 
0623   std::string prefix;
0624   if (v.sim_vertices.size() == 1) {
0625     prefix = "RecoAllAssoc2GenSingleMatched_SharedTrackFraction";
0626   } else if (v.sim_vertices.size() > 1) {
0627     prefix = "RecoAllAssoc2GenMultiMatched_SharedTrackFraction";
0628   }
0629 
0630   for (size_t i = 0; i < v.sim_vertices.size(); ++i) {
0631     const double sharedTracks = v.sim_vertices_num_shared_tracks[i];
0632     const simPrimaryVertex* simV = v.sim_vertices_internal[i];
0633     mes_[label][prefix + "Reco"]->Fill(sharedTracks / v.nRecoTrk);
0634     mes_[label][prefix + "RecoMatched"]->Fill(sharedTracks / v.num_matched_sim_tracks);
0635     mes_[label][prefix + "Sim"]->Fill(sharedTracks / simV->nGenTrk);
0636     mes_[label][prefix + "SimMatched"]->Fill(sharedTracks / simV->num_matched_reco_tracks);
0637   }
0638 }
0639 
0640 void PrimaryVertexAnalyzer4PUSlimmed::fillResolutionAndPullHistograms(
0641     const std::string& label,
0642     int num_pileup_vertices,
0643     PrimaryVertexAnalyzer4PUSlimmed::recoPrimaryVertex& v,
0644     bool isPV) {
0645   std::string prefix = "RecoAllAssoc2GenMatched";
0646   if (v.sim_vertices_internal.size() > 1) {
0647     prefix += "Merged";
0648   }
0649   if (isPV) {
0650     prefix = "RecoPVAssoc2GenPVMatched";
0651   }
0652 
0653   // Use the best match as defined by the vertex truth associator
0654   // reco-tracks as the best match
0655   const simPrimaryVertex& bestMatch = *(v.sim_vertices_internal[0]);
0656   const double xres = v.x - bestMatch.x;
0657   const double yres = v.y - bestMatch.y;
0658   const double zres = v.z - bestMatch.z;
0659   const double pt2res = v.ptsq - bestMatch.ptsq;
0660 
0661   const double xresol = xres;
0662   const double yresol = yres;
0663   const double zresol = zres;
0664   const double pt2resol = pt2res / v.ptsq;
0665   const double xpull = xres / v.recVtx->xError();
0666   const double ypull = yres / v.recVtx->yError();
0667   const double zpull = zres / v.recVtx->zError();
0668 
0669   mes_[label][prefix + "_ResolX"]->Fill(xresol);
0670   mes_[label][prefix + "_ResolY"]->Fill(yresol);
0671   mes_[label][prefix + "_ResolZ"]->Fill(zresol);
0672   mes_[label][prefix + "_ResolPt2"]->Fill(pt2resol);
0673   mes_[label][prefix + "_PullX"]->Fill(xpull);
0674   mes_[label][prefix + "_PullY"]->Fill(ypull);
0675   mes_[label][prefix + "_PullZ"]->Fill(zpull);
0676 
0677   mes_[label][prefix + "_ResolX_vs_PU"]->Fill(num_pileup_vertices, xresol);
0678   mes_[label][prefix + "_ResolY_vs_PU"]->Fill(num_pileup_vertices, yresol);
0679   mes_[label][prefix + "_ResolZ_vs_PU"]->Fill(num_pileup_vertices, zresol);
0680   mes_[label][prefix + "_ResolPt2_vs_PU"]->Fill(num_pileup_vertices, pt2resol);
0681   mes_[label][prefix + "_PullX_vs_PU"]->Fill(num_pileup_vertices, xpull);
0682   mes_[label][prefix + "_PullY_vs_PU"]->Fill(num_pileup_vertices, ypull);
0683   mes_[label][prefix + "_PullZ_vs_PU"]->Fill(num_pileup_vertices, zpull);
0684 
0685   mes_[label][prefix + "_ResolX_vs_NumTracks"]->Fill(v.nRecoTrk, xresol);
0686   mes_[label][prefix + "_ResolY_vs_NumTracks"]->Fill(v.nRecoTrk, yresol);
0687   mes_[label][prefix + "_ResolZ_vs_NumTracks"]->Fill(v.nRecoTrk, zresol);
0688   mes_[label][prefix + "_ResolPt2_vs_NumTracks"]->Fill(v.nRecoTrk, pt2resol);
0689   mes_[label][prefix + "_PullX_vs_NumTracks"]->Fill(v.nRecoTrk, xpull);
0690   mes_[label][prefix + "_PullY_vs_NumTracks"]->Fill(v.nRecoTrk, ypull);
0691   mes_[label][prefix + "_PullZ_vs_NumTracks"]->Fill(v.nRecoTrk, zpull);
0692 
0693   mes_[label][prefix + "_ResolX_vs_Z"]->Fill(v.z, xresol);
0694   mes_[label][prefix + "_ResolY_vs_Z"]->Fill(v.z, yresol);
0695   mes_[label][prefix + "_ResolZ_vs_Z"]->Fill(v.z, zresol);
0696   mes_[label][prefix + "_ResolPt2_vs_Z"]->Fill(v.z, pt2resol);
0697   mes_[label][prefix + "_PullX_vs_Z"]->Fill(v.z, xpull);
0698   mes_[label][prefix + "_PullY_vs_Z"]->Fill(v.z, ypull);
0699   mes_[label][prefix + "_PullZ_vs_Z"]->Fill(v.z, zpull);
0700 
0701   mes_[label][prefix + "_ResolX_vs_Pt"]->Fill(v.pt, xresol);
0702   mes_[label][prefix + "_ResolY_vs_Pt"]->Fill(v.pt, yresol);
0703   mes_[label][prefix + "_ResolZ_vs_Pt"]->Fill(v.pt, zresol);
0704   mes_[label][prefix + "_ResolPt2_vs_Pt"]->Fill(v.pt, pt2resol);
0705   mes_[label][prefix + "_PullX_vs_Pt"]->Fill(v.pt, xpull);
0706   mes_[label][prefix + "_PullY_vs_Pt"]->Fill(v.pt, ypull);
0707   mes_[label][prefix + "_PullZ_vs_Pt"]->Fill(v.pt, zpull);
0708 }
0709 
0710 bool PrimaryVertexAnalyzer4PUSlimmed::matchRecoTrack2SimSignal(const reco::TrackBaseRef& recoTrack) {
0711   auto found = r2s_->find(recoTrack);
0712 
0713   // reco track not matched to any TP
0714   if (found == r2s_->end())
0715     return false;
0716 
0717   // reco track matched to some TP from signal vertex
0718   for (const auto& tp : found->val) {
0719     if (tp.first->eventId().bunchCrossing() == 0 && tp.first->eventId().event() == 0)
0720       return true;
0721   }
0722 
0723   // reco track not matched to any TP from signal vertex
0724   return false;
0725 }
0726 
0727 void PrimaryVertexAnalyzer4PUSlimmed::calculatePurityAndFillHistograms(const std::string& label,
0728                                                                        std::vector<recoPrimaryVertex>& recopvs,
0729                                                                        int genpv_position_in_reco_collection,
0730                                                                        bool signal_is_highest_pt) {
0731   if (recopvs.empty())
0732     return;
0733 
0734   std::vector<double> vtx_sumpt_sigmatched;
0735   std::vector<double> vtx_sumpt2_sigmatched;
0736 
0737   vtx_sumpt_sigmatched.reserve(recopvs.size());
0738   vtx_sumpt2_sigmatched.reserve(recopvs.size());
0739 
0740   // Calculate purity
0741   for (auto& v : recopvs) {
0742     double sumpt_all = 0;
0743     double sumpt_sigmatched = 0;
0744     double sumpt2_sigmatched = 0;
0745     const reco::Vertex* vertex = v.recVtx;
0746     for (auto iTrack = vertex->tracks_begin(); iTrack != vertex->tracks_end(); ++iTrack) {
0747       double pt = (*iTrack)->pt();
0748       sumpt_all += pt;
0749       if (matchRecoTrack2SimSignal(*iTrack)) {
0750         sumpt_sigmatched += pt;
0751         sumpt2_sigmatched += pt * pt;
0752       }
0753     }
0754     v.purity = sumpt_sigmatched / sumpt_all;
0755 
0756     vtx_sumpt_sigmatched.push_back(sumpt_sigmatched);
0757     vtx_sumpt2_sigmatched.push_back(sumpt2_sigmatched);
0758   }
0759 
0760   const double vtxAll_sumpt_sigmatched = std::accumulate(vtx_sumpt_sigmatched.begin(), vtx_sumpt_sigmatched.end(), 0.0);
0761   const double vtxNot0_sumpt_sigmatched = vtxAll_sumpt_sigmatched - vtx_sumpt_sigmatched[0];
0762   const double missing = vtxNot0_sumpt_sigmatched / vtxAll_sumpt_sigmatched;
0763 
0764   // Fill purity
0765   std::string prefix = "RecoPVAssoc2GenPVNotMatched_";
0766   if (genpv_position_in_reco_collection == 0)
0767     prefix = "RecoPVAssoc2GenPVMatched_";
0768 
0769   mes_[label][prefix + "Purity"]->Fill(recopvs[0].purity);
0770   mes_[label][prefix + "Missing"]->Fill(missing);
0771   auto hpurity = mes_[label][prefix + "Purity_vs_Index"];
0772   for (size_t i = 0; i < recopvs.size(); ++i) {
0773     hpurity->Fill(i, recopvs[i].purity);
0774   }
0775 
0776   // Fill sumpt2
0777   for (size_t i = 0; i < recopvs.size(); ++i) {
0778     if (static_cast<int>(i) == genpv_position_in_reco_collection) {
0779       mes_[label]["RecoAssoc2GenPVMatched_Pt2"]->Fill(recopvs[i].ptsq);
0780     } else {
0781       double pt2 = recopvs[i].ptsq;
0782       mes_[label]["RecoAssoc2GenPVNotMatched_Pt2"]->Fill(pt2);
0783       // Subtract hard-scatter track pt2 from the pileup pt2
0784       double pt2_pu = pt2 - vtx_sumpt2_sigmatched[i];
0785       mes_[label]["RecoAssoc2GenPVNotMatched_GenPVTracksRemoved_Pt2"]->Fill(pt2_pu);
0786     }
0787   }
0788   if (!signal_is_highest_pt && genpv_position_in_reco_collection >= 0)
0789     mes_[label]["RecoAssoc2GenPVMatchedNotHighest_Pt2"]->Fill(recopvs[genpv_position_in_reco_collection].ptsq);
0790 }
0791 
0792 /* Extract information form TrackingParticles/TrackingVertex and fill
0793  * the helper class simPrimaryVertex with proper generation-level
0794  * information */
0795 std::vector<PrimaryVertexAnalyzer4PUSlimmed::simPrimaryVertex> PrimaryVertexAnalyzer4PUSlimmed::getSimPVs(
0796     const edm::Handle<TrackingVertexCollection>& tVC) {
0797   std::vector<PrimaryVertexAnalyzer4PUSlimmed::simPrimaryVertex> simpv;
0798   int current_event = -1;
0799 
0800   if (verbose_) {
0801     std::cout << "getSimPVs TrackingVertexCollection " << std::endl;
0802   }
0803 
0804   for (TrackingVertexCollection::const_iterator v = tVC->begin(); v != tVC->end(); ++v) {
0805     if (verbose_) {
0806       std::cout << "BunchX.EventId: " << v->eventId().bunchCrossing() << "." << (v->eventId()).event()
0807                 << " Position: " << v->position() << " G4/HepMC Vertices: " << v->g4Vertices().size() << "/"
0808                 << v->genVertices().size() << "   t = " << v->position().t() * 1.e12
0809                 << "    == 0:" << (v->position().t() > 0) << std::endl;
0810       for (TrackingVertex::g4v_iterator gv = v->g4Vertices_begin(); gv != v->g4Vertices_end(); gv++) {
0811         std::cout << *gv << std::endl;
0812       }
0813       std::cout << "----------" << std::endl;
0814 
0815     }  // end of verbose_ session
0816 
0817     // I'd rather change this and select only vertices that come from
0818     // BX=0.  We should keep only the first vertex from all the events
0819     // at BX=0.
0820     if (v->eventId().bunchCrossing() != 0)
0821       continue;
0822     if (v->eventId().event() != current_event) {
0823       current_event = v->eventId().event();
0824     } else {
0825       continue;
0826     }
0827     // TODO(rovere) is this really necessary?
0828     if (fabs(v->position().z()) > 1000)
0829       continue;  // skip funny junk vertices
0830 
0831     // could be a new vertex, check  all primaries found so far to avoid
0832     // multiple entries
0833     simPrimaryVertex sv(v->position().x(), v->position().y(), v->position().z());
0834     sv.eventId = v->eventId();
0835     sv.sim_vertex = TrackingVertexRef(tVC, std::distance(tVC->begin(), v));
0836 
0837     for (TrackingParticleRefVector::iterator iTrack = v->daughterTracks_begin(); iTrack != v->daughterTracks_end();
0838          ++iTrack) {
0839       // TODO(rovere) isn't it always the case? Is it really worth
0840       // checking this out?
0841       // sv.eventId = (**iTrack).eventId();
0842       assert((**iTrack).eventId().bunchCrossing() == 0);
0843     }
0844     // TODO(rovere) maybe get rid of this old logic completely ... ?
0845     simPrimaryVertex* vp = nullptr;  // will become non-NULL if a vertex
0846                                      // is found and then point to it
0847     for (std::vector<simPrimaryVertex>::iterator v0 = simpv.begin(); v0 != simpv.end(); v0++) {
0848       if ((sv.eventId == v0->eventId) && (fabs(sv.x - v0->x) < 1e-5) && (fabs(sv.y - v0->y) < 1e-5) &&
0849           (fabs(sv.z - v0->z) < 1e-5)) {
0850         vp = &(*v0);
0851         break;
0852       }
0853     }
0854     if (!vp) {
0855       // this is a new vertex, add it to the list of sim-vertices
0856       simpv.push_back(sv);
0857       vp = &simpv.back();
0858       if (verbose_) {
0859         std::cout << "this is a new vertex " << sv.eventId.event() << "   " << sv.x << " " << sv.y << " " << sv.z
0860                   << std::endl;
0861       }
0862     } else {
0863       if (verbose_) {
0864         std::cout << "this is not a new vertex" << sv.x << " " << sv.y << " " << sv.z << std::endl;
0865       }
0866     }
0867 
0868     // Loop over daughter track(s) as Tracking Particles
0869     for (TrackingVertex::tp_iterator iTP = v->daughterTracks_begin(); iTP != v->daughterTracks_end(); ++iTP) {
0870       auto momentum = (*(*iTP)).momentum();
0871       const reco::Track* matched_best_reco_track = nullptr;
0872       double match_quality = -1;
0873       if (use_only_charged_tracks_ && (**iTP).charge() == 0)
0874         continue;
0875       if (s2r_->find(*iTP) != s2r_->end()) {
0876         matched_best_reco_track = (*s2r_)[*iTP][0].first.get();
0877         match_quality = (*s2r_)[*iTP][0].second;
0878       }
0879       if (verbose_) {
0880         std::cout << "  Daughter momentum:      " << momentum;
0881         std::cout << "  Daughter type     " << (*(*iTP)).pdgId();
0882         std::cout << "  matched: " << (matched_best_reco_track != nullptr);
0883         std::cout << "  match-quality: " << match_quality;
0884         std::cout << std::endl;
0885       }
0886       vp->ptot.setPx(vp->ptot.x() + momentum.x());
0887       vp->ptot.setPy(vp->ptot.y() + momentum.y());
0888       vp->ptot.setPz(vp->ptot.z() + momentum.z());
0889       vp->ptot.setE(vp->ptot.e() + (**iTP).energy());
0890       vp->ptsq += ((**iTP).pt() * (**iTP).pt());
0891       // TODO(rovere) only select charged sim-particles? If so, maybe
0892       // put it as a configuration parameter?
0893       if (matched_best_reco_track) {
0894         vp->num_matched_reco_tracks++;
0895         vp->average_match_quality += match_quality;
0896       }
0897       // TODO(rovere) get rid of cuts on sim-tracks
0898       // TODO(rovere) be consistent between simulated tracks and
0899       // reconstructed tracks selection
0900       // count relevant particles
0901       if (((**iTP).pt() > 0.2) && (fabs((**iTP).eta()) < 2.5) && (**iTP).charge() != 0) {
0902         vp->nGenTrk++;
0903       }
0904     }  // End of for loop on daughters sim-particles
0905     if (vp->num_matched_reco_tracks)
0906       vp->average_match_quality /= static_cast<float>(vp->num_matched_reco_tracks);
0907     if (verbose_) {
0908       std::cout << "average number of associated tracks: "
0909                 << vp->num_matched_reco_tracks / static_cast<float>(vp->nGenTrk)
0910                 << " with average quality: " << vp->average_match_quality << std::endl;
0911     }
0912   }  // End of for loop on tracking vertices
0913 
0914   if (verbose_) {
0915     std::cout << "------- PrimaryVertexAnalyzer4PUSlimmed simPVs from "
0916                  "TrackingVertices "
0917                  "-------"
0918               << std::endl;
0919     for (std::vector<simPrimaryVertex>::iterator v0 = simpv.begin(); v0 != simpv.end(); v0++) {
0920       std::cout << "z=" << v0->z << "  event=" << v0->eventId.event() << std::endl;
0921     }
0922     std::cout << "-----------------------------------------------" << std::endl;
0923   }  // End of for summary on discovered simulated primary vertices.
0924 
0925   // In case of no simulated vertices, break here
0926   if (simpv.empty())
0927     return simpv;
0928 
0929   // Now compute the closest distance in z between all simulated vertex
0930   // first initialize
0931   auto prev_z = simpv.back().z;
0932   for (simPrimaryVertex& vsim : simpv) {
0933     vsim.closest_vertex_distance_z = std::abs(vsim.z - prev_z);
0934     prev_z = vsim.z;
0935   }
0936   // then calculate
0937   for (std::vector<simPrimaryVertex>::iterator vsim = simpv.begin(); vsim != simpv.end(); vsim++) {
0938     std::vector<simPrimaryVertex>::iterator vsim2 = vsim;
0939     vsim2++;
0940     for (; vsim2 != simpv.end(); vsim2++) {
0941       double distance = std::abs(vsim->z - vsim2->z);
0942       // need both to be complete
0943       vsim->closest_vertex_distance_z = std::min(vsim->closest_vertex_distance_z, distance);
0944       vsim2->closest_vertex_distance_z = std::min(vsim2->closest_vertex_distance_z, distance);
0945     }
0946   }
0947   return simpv;
0948 }
0949 
0950 /* Extract information form recoVertex and fill the helper class
0951  * recoPrimaryVertex with proper reco-level information */
0952 std::vector<PrimaryVertexAnalyzer4PUSlimmed::recoPrimaryVertex> PrimaryVertexAnalyzer4PUSlimmed::getRecoPVs(
0953     const edm::Handle<edm::View<reco::Vertex>>& tVC) {
0954   std::vector<PrimaryVertexAnalyzer4PUSlimmed::recoPrimaryVertex> recopv;
0955 
0956   if (verbose_) {
0957     std::cout << "getRecoPVs TrackingVertexCollection " << std::endl;
0958   }
0959 
0960   for (auto v = tVC->begin(); v != tVC->end(); ++v) {
0961     if (verbose_) {
0962       std::cout << " Position: " << v->position() << std::endl;
0963     }
0964 
0965     // Skip junk vertices
0966     if (fabs(v->z()) > 1000)
0967       continue;
0968     if (v->isFake() || !v->isValid())
0969       continue;
0970 
0971     recoPrimaryVertex sv(v->position().x(), v->position().y(), v->position().z());
0972     sv.recVtx = &(*v);
0973     sv.recVtxRef = reco::VertexBaseRef(tVC, std::distance(tVC->begin(), v));
0974     // this is a new vertex, add it to the list of reco-vertices
0975     recopv.push_back(sv);
0976     PrimaryVertexAnalyzer4PUSlimmed::recoPrimaryVertex* vp = &recopv.back();
0977 
0978     // Loop over daughter track(s)
0979     for (auto iTrack = v->tracks_begin(); iTrack != v->tracks_end(); ++iTrack) {
0980       auto momentum = (*(*iTrack)).innerMomentum();
0981       // TODO(rovere) better handle the pixelVertices, whose tracks
0982       // do not have the innerMomentum defined. This is a temporary
0983       // hack to overcome this problem.
0984       if (momentum.mag2() == 0)
0985         momentum = (*(*iTrack)).momentum();
0986       if (verbose_) {
0987         std::cout << "  Daughter momentum:      " << momentum;
0988         std::cout << std::endl;
0989       }
0990       vp->pt += std::sqrt(momentum.perp2());
0991       vp->ptsq += (momentum.perp2());
0992       vp->nRecoTrk++;
0993 
0994       auto matched = r2s_->find(*iTrack);
0995       if (matched != r2s_->end()) {
0996         vp->num_matched_sim_tracks++;
0997       }
0998 
0999     }  // End of for loop on daughters reconstructed tracks
1000   }    // End of for loop on tracking vertices
1001 
1002   if (verbose_) {
1003     std::cout << "------- PrimaryVertexAnalyzer4PUSlimmed recoPVs from "
1004                  "VertexCollection "
1005                  "-------"
1006               << std::endl;
1007     for (std::vector<recoPrimaryVertex>::iterator v0 = recopv.begin(); v0 != recopv.end(); v0++) {
1008       std::cout << "z=" << v0->z << std::endl;
1009     }
1010     std::cout << "-----------------------------------------------" << std::endl;
1011   }  // End of for summary on reconstructed primary vertices.
1012 
1013   // In case of no reco vertices, break here
1014   if (recopv.empty())
1015     return recopv;
1016 
1017   // Now compute the closest distance in z between all reconstructed vertex
1018   // first initialize
1019   auto prev_z = recopv.back().z;
1020   for (recoPrimaryVertex& vreco : recopv) {
1021     vreco.closest_vertex_distance_z = std::abs(vreco.z - prev_z);
1022     prev_z = vreco.z;
1023   }
1024   for (std::vector<recoPrimaryVertex>::iterator vreco = recopv.begin(); vreco != recopv.end(); vreco++) {
1025     std::vector<recoPrimaryVertex>::iterator vreco2 = vreco;
1026     vreco2++;
1027     for (; vreco2 != recopv.end(); vreco2++) {
1028       double distance = std::abs(vreco->z - vreco2->z);
1029       // need both to be complete
1030       vreco->closest_vertex_distance_z = std::min(vreco->closest_vertex_distance_z, distance);
1031       vreco2->closest_vertex_distance_z = std::min(vreco2->closest_vertex_distance_z, distance);
1032     }
1033   }
1034   return recopv;
1035 }
1036 
1037 void PrimaryVertexAnalyzer4PUSlimmed::resetSimPVAssociation(std::vector<simPrimaryVertex>& simpv) {
1038   for (auto& v : simpv) {
1039     v.rec_vertices.clear();
1040   }
1041 }
1042 
1043 // ------------ method called to produce the data  ------------
1044 void PrimaryVertexAnalyzer4PUSlimmed::matchSim2RecoVertices(std::vector<simPrimaryVertex>& simpv,
1045                                                             const reco::VertexSimToRecoCollection& vertex_s2r) {
1046   if (verbose_) {
1047     std::cout << "PrimaryVertexAnalyzer4PUSlimmed::matchSim2RecoVertices " << std::endl;
1048   }
1049   for (std::vector<simPrimaryVertex>::iterator vsim = simpv.begin(); vsim != simpv.end(); vsim++) {
1050     auto matched = vertex_s2r.find(vsim->sim_vertex);
1051     if (matched != vertex_s2r.end()) {
1052       for (const auto& vertexRefQuality : matched->val) {
1053         vsim->rec_vertices.push_back(&(*(vertexRefQuality.first)));
1054       }
1055     }
1056 
1057     if (verbose_) {
1058       if (!vsim->rec_vertices.empty()) {
1059         for (auto const& v : vsim->rec_vertices) {
1060           std::cout << "Found a matching vertex for genVtx " << vsim->z << " at " << v->z()
1061                     << " with sign: " << fabs(v->z() - vsim->z) / v->zError() << std::endl;
1062         }
1063       } else {
1064         std::cout << "No matching vertex for " << vsim->z << std::endl;
1065       }
1066     }
1067   }  // end for loop on simulated vertices
1068   if (verbose_) {
1069     std::cout << "Done with matching sim vertices" << std::endl;
1070   }
1071 }
1072 
1073 void PrimaryVertexAnalyzer4PUSlimmed::matchReco2SimVertices(std::vector<recoPrimaryVertex>& recopv,
1074                                                             const reco::VertexRecoToSimCollection& vertex_r2s,
1075                                                             const std::vector<simPrimaryVertex>& simpv) {
1076   for (std::vector<recoPrimaryVertex>::iterator vrec = recopv.begin(); vrec != recopv.end(); vrec++) {
1077     auto matched = vertex_r2s.find(vrec->recVtxRef);
1078     if (matched != vertex_r2s.end()) {
1079       for (const auto& vertexRefQuality : matched->val) {
1080         const auto tvPtr = &(*(vertexRefQuality.first));
1081         vrec->sim_vertices.push_back(tvPtr);
1082       }
1083 
1084       for (const TrackingVertex* tv : vrec->sim_vertices) {
1085         // Set pointers to internal simVertex objects
1086         for (const auto& vv : simpv) {
1087           if (&(*(vv.sim_vertex)) == tv) {
1088             vrec->sim_vertices_internal.push_back(&vv);
1089             continue;
1090           }
1091         }
1092 
1093         // Calculate number of shared tracks
1094         vrec->sim_vertices_num_shared_tracks.push_back(calculateVertexSharedTracks(*(vrec->recVtx), *tv, *r2s_));
1095       }
1096     }
1097 
1098     if (verbose_) {
1099       for (auto v : vrec->sim_vertices) {
1100         std::cout << "Found a matching vertex for reco: " << vrec->z << " at gen:" << v->position().z()
1101                   << " with sign: " << fabs(vrec->z - v->position().z()) / vrec->recVtx->zError() << std::endl;
1102       }
1103     }
1104   }  // end for loop on reconstructed vertices
1105 }
1106 
1107 void PrimaryVertexAnalyzer4PUSlimmed::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
1108   using edm::Handle;
1109   using edm::View;
1110   using std::cout;
1111   using std::endl;
1112   using std::vector;
1113   using namespace reco;
1114 
1115   std::vector<float> pileUpInfo_z;
1116 
1117   // get the pileup information
1118   edm::Handle<std::vector<PileupSummaryInfo>> puinfoH;
1119   if (iEvent.getByToken(vecPileupSummaryInfoToken_, puinfoH)) {
1120     for (auto const& pu_info : *puinfoH.product()) {
1121       if (do_generic_sim_plots_) {
1122         mes_["root_folder"]["GenVtx_vs_BX"]->Fill(pu_info.getBunchCrossing(), pu_info.getPU_NumInteractions());
1123       }
1124       if (pu_info.getBunchCrossing() == 0) {
1125         pileUpInfo_z = pu_info.getPU_zpositions();
1126         if (verbose_) {
1127           for (auto const& p : pileUpInfo_z) {
1128             std::cout << "PileUpInfo on Z vertex: " << p << std::endl;
1129           }
1130         }
1131         break;
1132       }
1133     }
1134   }
1135 
1136   edm::Handle<TrackingParticleCollection> TPCollectionH;
1137   iEvent.getByToken(trackingParticleCollectionToken_, TPCollectionH);
1138   if (!TPCollectionH.isValid())
1139     edm::LogWarning("PrimaryVertexAnalyzer4PUSlimmed") << "TPCollectionH is not valid";
1140 
1141   edm::Handle<TrackingVertexCollection> TVCollectionH;
1142   iEvent.getByToken(trackingVertexCollectionToken_, TVCollectionH);
1143   if (!TVCollectionH.isValid())
1144     edm::LogWarning("PrimaryVertexAnalyzer4PUSlimmed") << "TVCollectionH is not valid";
1145 
1146   // TODO(rovere) the idea is to put in case a track-selector in front
1147   // of this module and then use its label to get the selected tracks
1148   // out of the event instead of making an hard-coded selection in the
1149   // code.
1150 
1151   edm::Handle<reco::SimToRecoCollection> simToRecoH;
1152   iEvent.getByToken(simToRecoAssociationToken_, simToRecoH);
1153   if (simToRecoH.isValid())
1154     s2r_ = simToRecoH.product();
1155   else
1156     edm::LogWarning("PrimaryVertexAnalyzer4PUSlimmed") << "simToRecoH is not valid";
1157 
1158   edm::Handle<reco::RecoToSimCollection> recoToSimH;
1159   iEvent.getByToken(recoToSimAssociationToken_, recoToSimH);
1160   if (recoToSimH.isValid())
1161     r2s_ = recoToSimH.product();
1162   else
1163     edm::LogWarning("PrimaryVertexAnalyzer4PUSlimmed") << "recoToSimH is not valid";
1164 
1165   // Vertex associator
1166   edm::Handle<reco::VertexToTrackingVertexAssociator> vertexAssociatorH;
1167   iEvent.getByToken(vertexAssociatorToken_, vertexAssociatorH);
1168   if (!vertexAssociatorH.isValid()) {
1169     edm::LogWarning("PrimaryVertexAnalyzer4PUSlimmed") << "vertexAssociatorH is not valid";
1170     return;
1171   }
1172   const reco::VertexToTrackingVertexAssociator& vertexAssociator = *(vertexAssociatorH.product());
1173 
1174   std::vector<simPrimaryVertex> simpv;  // a list of simulated primary
1175                                         // MC vertices
1176   // TODO(rovere) use move semantic?
1177   simpv = getSimPVs(TVCollectionH);
1178   // TODO(rovere) 1 vertex is not, by definition, pileup, and should
1179   // probably be subtracted?
1180   int kind_of_signal_vertex = 0;
1181   int num_pileup_vertices = simpv.size();
1182   if (do_generic_sim_plots_)
1183     mes_["root_folder"]["GenAllV_NumVertices"]->Fill(simpv.size());
1184   bool signal_is_highest_pt =
1185       std::max_element(simpv.begin(), simpv.end(), [](const simPrimaryVertex& lhs, const simPrimaryVertex& rhs) {
1186         return lhs.ptsq < rhs.ptsq;
1187       }) == simpv.begin();
1188   kind_of_signal_vertex |= (signal_is_highest_pt << HIGHEST_PT);
1189   if (do_generic_sim_plots_) {
1190     mes_["root_folder"]["SignalIsHighestPt2"]->Fill(signal_is_highest_pt ? 1. : 0.);
1191     computePairDistance(simpv, mes_["root_folder"]["GenAllV_PairDistanceZ"]);
1192   }
1193 
1194   int label_index = -1;
1195   for (size_t iToken = 0, endToken = reco_vertex_collection_tokens_.size(); iToken < endToken; ++iToken) {
1196     auto const& vertex_token = reco_vertex_collection_tokens_[iToken];
1197     std::vector<recoPrimaryVertex> recopv;  // a list of reconstructed
1198                                             // primary MC vertices
1199     std::string label = reco_vertex_collections_[++label_index].label();
1200     edm::Handle<edm::View<reco::Vertex>> recVtxs;
1201     if (!iEvent.getByToken(vertex_token, recVtxs)) {
1202       if (!errorPrintedForColl_[iToken]) {
1203         edm::LogWarning("PrimaryVertexAnalyzer4PUSlimmed")
1204             << "Skipping vertex collection: " << label << " since it is missing.";
1205         errorPrintedForColl_[iToken] = true;
1206       }
1207       continue;
1208     }
1209 
1210     {
1211       // check upfront that refs to track are (likely) to be valid
1212       bool ok = true;
1213       for (const auto& v : *recVtxs) {
1214         if (v.tracksSize() > 0) {
1215           const auto& ref = v.trackRefAt(0);
1216           if (ref.isNull() || !ref.isAvailable()) {
1217             if (!errorPrintedForColl_[iToken]) {
1218               edm::LogWarning("PrimaryVertexAnalyzer4PUSlimmed")
1219                   << "Skipping vertex collection: " << label
1220                   << " since likely the track collection the vertex has refs pointing to is missing (at least the "
1221                      "first TrackBaseRef is null or not available)";
1222               errorPrintedForColl_[iToken] = true;
1223             }
1224             ok = false;
1225           }
1226         }
1227       }
1228       if (!ok)
1229         continue;
1230     }
1231 
1232     reco::VertexRecoToSimCollection vertex_r2s = vertexAssociator.associateRecoToSim(recVtxs, TVCollectionH);
1233     reco::VertexSimToRecoCollection vertex_s2r = vertexAssociator.associateSimToReco(recVtxs, TVCollectionH);
1234 
1235     resetSimPVAssociation(simpv);
1236     matchSim2RecoVertices(simpv, vertex_s2r);
1237     recopv = getRecoPVs(recVtxs);
1238     computePairDistance(recopv, mes_[label]["RecoAllAssoc2Gen_PairDistanceZ"]);
1239     matchReco2SimVertices(recopv, vertex_r2s, simpv);
1240 
1241     int num_total_gen_vertices_assoc2reco = 0;
1242     int num_total_reco_vertices_assoc2gen = 0;
1243     int num_total_gen_vertices_multiassoc2reco = 0;
1244     int num_total_reco_vertices_multiassoc2gen = 0;
1245     int num_total_reco_vertices_duplicate = 0;
1246     int genpv_position_in_reco_collection = -1;
1247     for (auto const& v : simpv) {
1248       float mistag = 1.;
1249       // TODO(rovere) put selectors here in front of fill* methods.
1250       if (v.eventId.event() == 0) {
1251         if (!recVtxs->empty() && std::find(v.rec_vertices.begin(), v.rec_vertices.end(), &((*recVtxs.product())[0])) !=
1252                                      v.rec_vertices.end()) {
1253           mistag = 0.;
1254           kind_of_signal_vertex |= (1 << IS_ASSOC2FIRST_RECO);
1255         } else {
1256           if (!v.rec_vertices.empty()) {
1257             kind_of_signal_vertex |= (1 << IS_ASSOC2ANY_RECO);
1258           }
1259         }
1260         mes_[label]["KindOfSignalPV"]->Fill(kind_of_signal_vertex);
1261         mes_[label]["MisTagRate"]->Fill(mistag);
1262         mes_[label]["MisTagRate_vs_PU"]->Fill(simpv.size(), mistag);
1263         mes_[label]["MisTagRate_vs_sum-pt2"]->Fill(v.ptsq, mistag);
1264         mes_[label]["MisTagRate_vs_Z"]->Fill(v.z, mistag);
1265         mes_[label]["MisTagRate_vs_R"]->Fill(v.r, mistag);
1266         mes_[label]["MisTagRate_vs_NumTracks"]->Fill(v.nGenTrk, mistag);
1267         if (signal_is_highest_pt) {
1268           mes_[label]["MisTagRateSignalIsHighest"]->Fill(mistag);
1269           mes_[label]["MisTagRateSignalIsHighest_vs_PU"]->Fill(simpv.size(), mistag);
1270           mes_[label]["MisTagRateSignalIsHighest_vs_sum-pt2"]->Fill(v.ptsq, mistag);
1271           mes_[label]["MisTagRateSignalIsHighest_vs_Z"]->Fill(v.z, mistag);
1272           mes_[label]["MisTagRateSignalIsHighest_vs_R"]->Fill(v.r, mistag);
1273           mes_[label]["MisTagRateSignalIsHighest_vs_NumTracks"]->Fill(v.nGenTrk, mistag);
1274         } else {
1275           mes_[label]["MisTagRateSignalIsNotHighest"]->Fill(mistag);
1276           mes_[label]["MisTagRateSignalIsNotHighest_vs_PU"]->Fill(simpv.size(), mistag);
1277           mes_[label]["MisTagRateSignalIsNotHighest_vs_sum-pt2"]->Fill(v.ptsq, mistag);
1278           mes_[label]["MisTagRateSignalIsNotHighest_vs_Z"]->Fill(v.z, mistag);
1279           mes_[label]["MisTagRateSignalIsNotHighest_vs_R"]->Fill(v.r, mistag);
1280           mes_[label]["MisTagRateSignalIsNotHighest_vs_NumTracks"]->Fill(v.nGenTrk, mistag);
1281         }
1282         // Now check at which location the Simulated PV has been
1283         // reconstructed in the primary vertex collection
1284         // at-hand. Mark it with fake index -1 if it was not
1285         // reconstructed at all.
1286 
1287         auto iv = (*recVtxs.product()).begin();
1288         for (int pv_position_in_reco_collection = 0; iv != (*recVtxs.product()).end();
1289              ++pv_position_in_reco_collection, ++iv) {
1290           if (std::find(v.rec_vertices.begin(), v.rec_vertices.end(), &(*iv)) != v.rec_vertices.end()) {
1291             mes_[label]["TruePVLocationIndex"]->Fill(pv_position_in_reco_collection);
1292             const bool genPVMatchedToRecoPV = (pv_position_in_reco_collection == 0);
1293             mes_[label]["TruePVLocationIndexCumulative"]->Fill(genPVMatchedToRecoPV ? 0 : 1);
1294 
1295             if (signal_is_highest_pt) {
1296               mes_[label]["TruePVLocationIndexSignalIsHighest"]->Fill(pv_position_in_reco_collection);
1297             } else {
1298               mes_[label]["TruePVLocationIndexSignalIsNotHighest"]->Fill(pv_position_in_reco_collection);
1299             }
1300 
1301             fillRecoAssociatedGenPVHistograms(label, v, genPVMatchedToRecoPV);
1302             if (genPVMatchedToRecoPV) {
1303               auto pv = recopv[0];
1304               assert(pv.recVtx == &(*iv));
1305               fillResolutionAndPullHistograms(label, num_pileup_vertices, pv, true);
1306             }
1307             genpv_position_in_reco_collection = pv_position_in_reco_collection;
1308             break;
1309           }
1310         }
1311 
1312         // If we reached the end, it means that the Simulated PV has not
1313         // been associated to any reconstructed vertex: mark it as
1314         // missing in the reconstructed vertex collection using the fake
1315         // index -1.
1316         if (iv == (*recVtxs.product()).end()) {
1317           mes_[label]["TruePVLocationIndex"]->Fill(-1.);
1318           mes_[label]["TruePVLocationIndexCumulative"]->Fill(-1.);
1319           if (signal_is_highest_pt)
1320             mes_[label]["TruePVLocationIndexSignalIsHighest"]->Fill(-1.);
1321           else
1322             mes_[label]["TruePVLocationIndexSignalIsNotHighest"]->Fill(-1.);
1323         }
1324       }
1325 
1326       if (!v.rec_vertices.empty())
1327         num_total_gen_vertices_assoc2reco++;
1328       if (v.rec_vertices.size() > 1)
1329         num_total_gen_vertices_multiassoc2reco++;
1330       // No need to N-tplicate the Gen-related cumulative histograms:
1331       // fill them only at the first iteration
1332       if (do_generic_sim_plots_ && label_index == 0)
1333         fillGenericGenVertexHistograms(v);
1334       fillRecoAssociatedGenVertexHistograms(label, v);
1335     }
1336     calculatePurityAndFillHistograms(label, recopv, genpv_position_in_reco_collection, signal_is_highest_pt);
1337 
1338     mes_[label]["GenAllAssoc2Reco_NumVertices"]->Fill(simpv.size(), simpv.size());
1339     mes_[label]["GenAllAssoc2RecoMatched_NumVertices"]->Fill(simpv.size(), num_total_gen_vertices_assoc2reco);
1340     mes_[label]["GenAllAssoc2RecoMultiMatched_NumVertices"]->Fill(simpv.size(), num_total_gen_vertices_multiassoc2reco);
1341     for (auto& v : recopv) {
1342       fillGenAssociatedRecoVertexHistograms(label, num_pileup_vertices, v);
1343       if (!v.sim_vertices.empty()) {
1344         num_total_reco_vertices_assoc2gen++;
1345         if (v.sim_vertices_internal[0]->rec_vertices.size() > 1) {
1346           num_total_reco_vertices_duplicate++;
1347         }
1348       }
1349       if (v.sim_vertices.size() > 1)
1350         num_total_reco_vertices_multiassoc2gen++;
1351     }
1352     mes_[label]["RecoAllAssoc2Gen_NumVertices"]->Fill(recopv.size(), recopv.size());
1353     mes_[label]["RecoAllAssoc2GenMatched_NumVertices"]->Fill(recopv.size(), num_total_reco_vertices_assoc2gen);
1354     mes_[label]["RecoAllAssoc2GenMultiMatched_NumVertices"]->Fill(recopv.size(),
1355                                                                   num_total_reco_vertices_multiassoc2gen);
1356     mes_[label]["RecoAllAssoc2MultiMatchedGen_NumVertices"]->Fill(recopv.size(), num_total_reco_vertices_duplicate);
1357     mes_[label]["RecoVtx_vs_GenVtx"]->Fill(simpv.size(), recopv.size());
1358     mes_[label]["MatchedRecoVtx_vs_GenVtx"]->Fill(simpv.size(), num_total_reco_vertices_assoc2gen);
1359   }
1360 }  // end of analyze
1361 
1362 template <class T>
1363 void PrimaryVertexAnalyzer4PUSlimmed::computePairDistance(const T& collection, MonitorElement* me) {
1364   for (unsigned int i = 0; i < collection.size(); ++i) {
1365     for (unsigned int j = i + 1; j < collection.size(); ++j) {
1366       me->Fill(std::abs(collection[i].z - collection[j].z));
1367     }
1368   }
1369 }