Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-05-10 02:20:38

0001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0002 #include "DQM/TrackingMonitorSource/interface/TrackToTrackComparisonHists.h"
0003 
0004 #include "DQMServices/Core/interface/DQMStore.h"
0005 
0006 #include "DataFormats/TrackReco/interface/Track.h"
0007 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0008 #include "DataFormats/VertexReco/interface/Vertex.h"
0009 #include "DataFormats/Luminosity/interface/LumiSummary.h"
0010 
0011 #include "DataFormats/Math/interface/deltaR.h"
0012 #include "CommonTools/TriggerUtils/interface/GenericTriggerEventFlag.h"
0013 
0014 //
0015 // constructors and destructor
0016 //
0017 TrackToTrackComparisonHists::TrackToTrackComparisonHists(const edm::ParameterSet& iConfig)
0018     : monitoredTrackInputTag_(iConfig.getParameter<edm::InputTag>("monitoredTrack")),
0019       referenceTrackInputTag_(iConfig.getParameter<edm::InputTag>("referenceTrack")),
0020       topDirName_(iConfig.getParameter<std::string>("topDirName")),
0021       dRmin_(iConfig.getParameter<double>("dRmin")),
0022       pTCutForPlateau_(iConfig.getParameter<double>("pTCutForPlateau")),
0023       dxyCutForPlateau_(iConfig.getParameter<double>("dxyCutForPlateau")),
0024       dzWRTPvCut_(iConfig.getParameter<double>("dzWRTPvCut")),
0025       requireValidHLTPaths_(iConfig.getParameter<bool>("requireValidHLTPaths")),
0026       genTriggerEventFlag_(new GenericTriggerEventFlag(
0027           iConfig.getParameter<edm::ParameterSet>("genericTriggerEventPSet"), consumesCollector(), *this))
0028 
0029 {
0030   initialize_parameter(iConfig);
0031 
0032   //now do what ever initialization is needed
0033   monitoredTrackToken_ = consumes<edm::View<reco::Track>>(monitoredTrackInputTag_);
0034   referenceTrackToken_ = consumes<edm::View<reco::Track>>(referenceTrackInputTag_);
0035   monitoredBSToken_ = consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("monitoredBeamSpot"));
0036   referenceBSToken_ = consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("referenceBeamSpot"));
0037   monitoredPVToken_ = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("monitoredPrimaryVertices"));
0038   referencePVToken_ = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("referencePrimaryVertices"));
0039   lumiScalersToken_ = consumes<LumiScalersCollection>(iConfig.getParameter<edm::InputTag>("scalers"));
0040   onlineMetaDataDigisToken_ =
0041       consumes<OnlineLuminosityRecord>(iConfig.getParameter<edm::InputTag>("onlineMetaDataDigis"));
0042 
0043   referenceTracksMEs_.label = referenceTrackInputTag_.label();
0044   matchedReferenceTracksMEs_.label = referenceTrackInputTag_.label() + "_matched";
0045 
0046   monitoredTracksMEs_.label = monitoredTrackInputTag_.label();
0047   unMatchedMonitoredTracksMEs_.label = monitoredTrackInputTag_.label() + "_unMatched";
0048 
0049   matchTracksMEs_.label = "matches";
0050 }
0051 
0052 TrackToTrackComparisonHists::~TrackToTrackComparisonHists() {
0053   if (genTriggerEventFlag_)
0054     genTriggerEventFlag_.reset();
0055 }
0056 
0057 void TrackToTrackComparisonHists::beginJob(const edm::EventSetup& iSetup) {}
0058 
0059 void TrackToTrackComparisonHists::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0060   LogDebug("TrackToTrackComparisonHists")
0061       << " requireValidHLTPaths_ " << requireValidHLTPaths_ << " hltPathsAreValid_  " << hltPathsAreValid_ << "\n";
0062   // if valid HLT paths are required,
0063   // analyze event only if paths are valid
0064   if (requireValidHLTPaths_ and (not hltPathsAreValid_)) {
0065     return;
0066   }
0067 
0068   LogDebug("TrackToTrackComparisonHists") << " genTriggerEventFlag_->on() " << genTriggerEventFlag_->on()
0069                                           << "  accept:  " << genTriggerEventFlag_->accept(iEvent, iSetup) << "\n";
0070   // Filter out events if Trigger Filtering is requested
0071   if (genTriggerEventFlag_->on() && !genTriggerEventFlag_->accept(iEvent, iSetup)) {
0072     return;
0073   }
0074 
0075   //
0076   //  Get Lumi/LS Info
0077   //
0078 
0079   unsigned int ls = iEvent.id().luminosityBlock();
0080 
0081   double onlinelumi = -1.f;
0082   double PU = -1.f;
0083 
0084   auto const lumiScalersHandle = iEvent.getHandle(lumiScalersToken_);
0085   auto const onlineMetaDataDigisHandle = iEvent.getHandle(onlineMetaDataDigisToken_);
0086 
0087   if (onlineMetaDataDigisHandle.isValid()) {
0088     onlinelumi = onlineMetaDataDigisHandle->instLumi();
0089     PU = onlineMetaDataDigisHandle->avgPileUp();
0090   } else if (lumiScalersHandle.isValid() and not lumiScalersHandle->empty()) {
0091     edm::LogError("TrackToTrackComparisonHists") << "onlineMetaDataDigisHandle not found, trying SCAL";
0092     auto const scalit = lumiScalersHandle->begin();
0093     onlinelumi = scalit->instantLumi();
0094     PU = scalit->pileup();
0095   } else {
0096     edm::LogError("TrackToTrackComparisonHists") << "lumiScalersHandle not found or empty, skipping event";
0097     return;
0098   }
0099 
0100   //
0101   //  Get Reference Track Info
0102   //
0103   edm::Handle<edm::View<reco::Track>> referenceTracksHandle;
0104   iEvent.getByToken(referenceTrackToken_, referenceTracksHandle);
0105   if (!referenceTracksHandle.isValid()) {
0106     edm::LogError("TrackToTrackComparisonHists")
0107         << "referenceTracksHandle with input tag " << referenceTrackInputTag_.encode() << " not found, skipping event";
0108     return;
0109   }
0110   const edm::View<reco::Track>& referenceTracks = *referenceTracksHandle;
0111 
0112   edm::Handle<reco::BeamSpot> referenceBSHandle;
0113   iEvent.getByToken(referenceBSToken_, referenceBSHandle);
0114   if (!referenceBSHandle.isValid()) {
0115     edm::LogError("TrackToTrackComparisonHists") << "referenceBSHandle not found, skipping event";
0116     return;
0117   }
0118   reco::BeamSpot referenceBS = *referenceBSHandle;
0119 
0120   edm::Handle<reco::VertexCollection> referencePVHandle;
0121   iEvent.getByToken(referencePVToken_, referencePVHandle);
0122   if (!referencePVHandle.isValid()) {
0123     edm::LogError("TrackToTrackComparisonHists") << "referencePVHandle not found, skipping event";
0124     return;
0125   }
0126   if (referencePVHandle->empty()) {
0127     edm::LogInfo("TrackToTrackComparisonHists") << "referencePVHandle->size is 0 ";
0128     return;
0129   }
0130   reco::Vertex referencePV = referencePVHandle->at(0);
0131 
0132   //
0133   //  Get Monitored Track Info
0134   //
0135   edm::Handle<edm::View<reco::Track>> monitoredTracksHandle;
0136   iEvent.getByToken(monitoredTrackToken_, monitoredTracksHandle);
0137   if (!monitoredTracksHandle.isValid()) {
0138     edm::LogError("TrackToTrackComparisonHists")
0139         << "monitoredTracksHandle with input tag " << monitoredTrackInputTag_.encode() << " not found, skipping event";
0140     return;
0141   }
0142   const edm::View<reco::Track>& monitoredTracks = *monitoredTracksHandle;
0143 
0144   edm::Handle<reco::BeamSpot> monitoredBSHandle;
0145   iEvent.getByToken(monitoredBSToken_, monitoredBSHandle);
0146   if (!monitoredTracksHandle.isValid()) {
0147     edm::LogError("TrackToTrackComparisonHists") << "monitoredBSHandle not found, skipping event";
0148     return;
0149   }
0150   reco::BeamSpot monitoredBS = *monitoredBSHandle;
0151 
0152   edm::Handle<reco::VertexCollection> monitoredPVHandle;
0153   iEvent.getByToken(monitoredPVToken_, monitoredPVHandle);
0154   if (!monitoredPVHandle.isValid()) {
0155     edm::LogError("TrackToTrackComparisonHists") << "monitoredPVHandle not found, skipping event";
0156     return;
0157   }
0158   if (monitoredPVHandle->empty()) {
0159     edm::LogInfo("TrackToTrackComparisonHists") << "monitoredPVHandle->size is 0 ";
0160     return;
0161   }
0162   reco::Vertex monitoredPV = monitoredPVHandle->at(0);
0163 
0164   edm::LogInfo("TrackToTrackComparisonHists")
0165       << "analyzing " << monitoredTrackInputTag_.process() << ":" << monitoredTrackInputTag_.label() << ":"
0166       << monitoredTrackInputTag_.instance() << " w.r.t. " << referenceTrackInputTag_.process() << ":"
0167       << referenceTrackInputTag_.label() << ":" << referenceTrackInputTag_.instance() << " \n";
0168 
0169   //
0170   // Build the dR maps
0171   //
0172   idx2idxByDoubleColl monitored2referenceColl;
0173   fillMap(monitoredTracks, referenceTracks, monitored2referenceColl, dRmin_);
0174 
0175   idx2idxByDoubleColl reference2monitoredColl;
0176   fillMap(referenceTracks, monitoredTracks, reference2monitoredColl, dRmin_);
0177 
0178   unsigned int nReferenceTracks(0);           // Counts the number of refernce tracks
0179   unsigned int nMatchedReferenceTracks(0);    // Counts the number of matched refernce tracks
0180   unsigned int nMonitoredTracks(0);           // Counts the number of monitored tracks
0181   unsigned int nUnmatchedMonitoredTracks(0);  // Counts the number of unmatched monitored tracks
0182 
0183   //
0184   // loop over reference tracks
0185   //
0186   LogDebug("TrackToTrackComparisonHists") << "\n# of tracks (reference): " << referenceTracks.size() << "\n";
0187   for (idx2idxByDoubleColl::const_iterator pItr = reference2monitoredColl.begin(), eItr = reference2monitoredColl.end();
0188        pItr != eItr;
0189        ++pItr) {
0190     nReferenceTracks++;
0191     int trackIdx = pItr->first;
0192     reco::Track track = referenceTracks.at(trackIdx);
0193 
0194     float dzWRTpv = track.dz(referencePV.position());
0195     if (fabs(dzWRTpv) > dzWRTPvCut_)
0196       continue;
0197 
0198     fill_generic_tracks_histos(*&referenceTracksMEs_, &track, &referenceBS, &referencePV, ls, onlinelumi, PU);
0199 
0200     std::map<double, int> trackDRmap = pItr->second;
0201     if (trackDRmap.empty()) {
0202       (matchedReferenceTracksMEs_.h_dRmin)->Fill(-1.);
0203       (matchedReferenceTracksMEs_.h_dRmin_l)->Fill(-1.);
0204       continue;
0205     }
0206 
0207     double dRmin = trackDRmap.begin()->first;
0208     (referenceTracksMEs_.h_dRmin)->Fill(dRmin);
0209     (referenceTracksMEs_.h_dRmin_l)->Fill(dRmin);
0210 
0211     bool matched = false;
0212     if (dRmin < dRmin_)
0213       matched = true;
0214 
0215     if (matched) {
0216       nMatchedReferenceTracks++;
0217       fill_generic_tracks_histos(*&matchedReferenceTracksMEs_, &track, &referenceBS, &referencePV, ls, onlinelumi, PU);
0218       (matchedReferenceTracksMEs_.h_dRmin)->Fill(dRmin);
0219       (matchedReferenceTracksMEs_.h_dRmin_l)->Fill(dRmin);
0220 
0221       int matchedTrackIndex = trackDRmap[dRmin];
0222       reco::Track matchedTrack = monitoredTracks.at(matchedTrackIndex);
0223       fill_matching_tracks_histos(*&matchTracksMEs_, &track, &matchedTrack, &referenceBS, &referencePV);
0224     }
0225 
0226   }  // Over reference tracks
0227 
0228   //
0229   // loop over monitoed tracks
0230   //
0231   LogDebug("TrackToTrackComparisonHists") << "\n# of tracks (monitored): " << monitoredTracks.size() << "\n";
0232   for (idx2idxByDoubleColl::const_iterator pItr = monitored2referenceColl.begin(), eItr = monitored2referenceColl.end();
0233        pItr != eItr;
0234        ++pItr) {
0235     nMonitoredTracks++;
0236     int trackIdx = pItr->first;
0237     reco::Track track = monitoredTracks.at(trackIdx);
0238 
0239     float dzWRTpv = track.dz(monitoredPV.position());
0240     if (fabs(dzWRTpv) > dzWRTPvCut_)
0241       continue;
0242 
0243     fill_generic_tracks_histos(*&monitoredTracksMEs_, &track, &monitoredBS, &monitoredPV, ls, onlinelumi, PU);
0244 
0245     std::map<double, int> trackDRmap = pItr->second;
0246     if (trackDRmap.empty()) {
0247       (unMatchedMonitoredTracksMEs_.h_dRmin)->Fill(-1.);
0248       (unMatchedMonitoredTracksMEs_.h_dRmin_l)->Fill(-1.);
0249       continue;
0250     }
0251 
0252     double dRmin = trackDRmap.begin()->first;
0253     (monitoredTracksMEs_.h_dRmin)->Fill(dRmin);
0254     (monitoredTracksMEs_.h_dRmin_l)->Fill(dRmin);
0255 
0256     bool matched = false;
0257     if (dRmin < dRmin_)
0258       matched = true;
0259 
0260     if (!matched) {
0261       nUnmatchedMonitoredTracks++;
0262       fill_generic_tracks_histos(
0263           *&unMatchedMonitoredTracksMEs_, &track, &monitoredBS, &monitoredPV, ls, onlinelumi, PU);
0264       (unMatchedMonitoredTracksMEs_.h_dRmin)->Fill(dRmin);
0265       (unMatchedMonitoredTracksMEs_.h_dRmin_l)->Fill(dRmin);
0266     }
0267 
0268   }  // over monitoed tracks
0269 
0270   edm::LogInfo("TrackToTrackComparisonHists")
0271       << "Total reference tracks: " << nReferenceTracks << "\n"
0272       << "Total matched reference tracks: " << nMatchedReferenceTracks << "\n"
0273       << "Total monitored tracks: " << nMonitoredTracks << "\n"
0274       << "Total unMatched monitored tracks: " << nUnmatchedMonitoredTracks << "\n";
0275 }
0276 
0277 void TrackToTrackComparisonHists::bookHistograms(DQMStore::IBooker& ibooker,
0278                                                  edm::Run const& iRun,
0279                                                  edm::EventSetup const& iSetup) {
0280   if (genTriggerEventFlag_ && genTriggerEventFlag_->on())
0281     genTriggerEventFlag_->initRun(iRun, iSetup);
0282 
0283   // check if every HLT path specified has a valid match in the HLT Menu
0284   hltPathsAreValid_ =
0285       (genTriggerEventFlag_ && genTriggerEventFlag_->on() && genTriggerEventFlag_->allHLTPathsAreValid());
0286 
0287   // if valid HLT paths are required,
0288   // create DQM outputs only if all paths are valid
0289   if (requireValidHLTPaths_ and (not hltPathsAreValid_)) {
0290     return;
0291   }
0292 
0293   std::string dir = topDirName_;
0294 
0295   bookHistos(ibooker, referenceTracksMEs_, "ref", dir);
0296   bookHistos(ibooker, matchedReferenceTracksMEs_, "ref_matched", dir);
0297 
0298   bookHistos(ibooker, monitoredTracksMEs_, "mon", dir);
0299   bookHistos(ibooker, unMatchedMonitoredTracksMEs_, "mon_unMatched", dir);
0300 
0301   book_matching_tracks_histos(ibooker, matchTracksMEs_, "matches", dir);
0302 }
0303 
0304 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0305 void TrackToTrackComparisonHists::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0306   edm::ParameterSetDescription desc;
0307 
0308   desc.add<bool>("requireValidHLTPaths", true);
0309 
0310   desc.add<edm::InputTag>("monitoredTrack", edm::InputTag("hltMergedTracks"));
0311   desc.add<edm::InputTag>("monitoredBeamSpot", edm::InputTag("hltOnlineBeamSpot"));
0312   desc.add<edm::InputTag>("monitoredPrimaryVertices", edm::InputTag("hltVerticesPFSelector"));
0313 
0314   desc.add<edm::InputTag>("referenceTrack", edm::InputTag("generalTracks"));
0315   desc.add<edm::InputTag>("referenceBeamSpot", edm::InputTag("offlineBeamSpot"));
0316   desc.add<edm::InputTag>("referencePrimaryVertices", edm::InputTag("offlinePrimaryVertices"));
0317   desc.add<edm::InputTag>("scalers", edm::InputTag("scalersRawToDigi"));
0318   desc.add<edm::InputTag>("onlineMetaDataDigis", edm::InputTag("onlineMetaDataDigis"));
0319 
0320   desc.add<std::string>("topDirName", "HLT/Tracking/ValidationWRTOffline");
0321   desc.add<double>("dRmin", 0.002);
0322 
0323   desc.add<double>("pTCutForPlateau", 0.9);
0324   desc.add<double>("dxyCutForPlateau", 2.5);
0325   desc.add<double>("dzWRTPvCut", 1e6);
0326 
0327   edm::ParameterSetDescription genericTriggerEventPSet;
0328   GenericTriggerEventFlag::fillPSetDescription(genericTriggerEventPSet);
0329   desc.add<edm::ParameterSetDescription>("genericTriggerEventPSet", genericTriggerEventPSet);
0330 
0331   edm::ParameterSetDescription histoPSet;
0332   fillHistoPSetDescription(histoPSet);
0333   desc.add<edm::ParameterSetDescription>("histoPSet", histoPSet);
0334 
0335   descriptions.add("trackToTrackComparisonHists", desc);
0336 }
0337 
0338 void TrackToTrackComparisonHists::fillMap(const edm::View<reco::Track>& tracks1,
0339                                           const edm::View<reco::Track>& tracks2,
0340                                           idx2idxByDoubleColl& map,
0341                                           float dRMin) {
0342   //
0343   // loop on tracks1
0344   //
0345   int i = 0;
0346   for (const auto& track1 : tracks1) {
0347     std::map<double, int> tmp;
0348     int j = 0;
0349     float smallest_dR = 1e9;
0350     int smallest_dR_j = -1;
0351 
0352     //
0353     // loop on tracks2
0354     //
0355     for (const auto& track2 : tracks2) {
0356       double dR = reco::deltaR(track1.eta(), track1.phi(), track2.eta(), track2.phi());
0357 
0358       if (dR < smallest_dR) {
0359         smallest_dR = dR;
0360         smallest_dR_j = j;
0361       }
0362 
0363       if (dR < dRMin) {
0364         tmp[dR] = j;
0365       }
0366 
0367       j++;
0368     }
0369 
0370     //
0371     // If there are no tracks that pass the dR store the smallest (for debugging/validating matching)
0372     //
0373     if (tmp.empty())
0374       tmp[smallest_dR] = smallest_dR_j;
0375 
0376     map.push_back(std::make_pair(i, tmp));
0377     i++;
0378   }
0379 }
0380 
0381 void TrackToTrackComparisonHists::bookHistos(DQMStore::IBooker& ibooker,
0382                                              generalME& mes,
0383                                              TString label,
0384                                              std::string& dir) {
0385   book_generic_tracks_histos(ibooker, mes, label, dir);
0386 }
0387 
0388 void TrackToTrackComparisonHists::book_generic_tracks_histos(DQMStore::IBooker& ibooker,
0389                                                              generalME& mes,
0390                                                              TString label,
0391                                                              std::string& dir) {
0392   ibooker.cd();
0393   ibooker.setCurrentFolder(dir);
0394   (mes.h_pt) = ibooker.book1D(label + "_pt", "track p_{T}", Pt_nbin, Pt_rangeMin, Pt_rangeMax);
0395   (mes.h_eta) = ibooker.book1D(label + "_eta", "track pseudorapidity", Eta_nbin, Eta_rangeMin, Eta_rangeMax);
0396   (mes.h_phi) = ibooker.book1D(label + "_phi", "track #phi", Phi_nbin, Phi_rangeMin, Phi_rangeMax);
0397   (mes.h_dxy) =
0398       ibooker.book1D(label + "_dxy", "track transverse dca to beam spot", Dxy_nbin, Dxy_rangeMin, Dxy_rangeMax);
0399   (mes.h_dz) = ibooker.book1D(label + "_dz", "track longitudinal dca to beam spot", Dz_nbin, Dz_rangeMin, Dz_rangeMax);
0400   (mes.h_dxyWRTpv) = ibooker.book1D(
0401       label + "_dxyWRTpv", "track transverse dca to primary vertex", Dxy_nbin, Dxy_rangeMin, Dxy_rangeMax);
0402   (mes.h_dzWRTpv) = ibooker.book1D(
0403       label + "_dzWRTpv", "track longitudinal dca to primary vertex", Dz_nbin, 0.1 * Dz_rangeMin, 0.1 * Dz_rangeMax);
0404   (mes.h_charge) = ibooker.book1D(label + "_charge", "track charge", 5, -2, 2);
0405   (mes.h_hits) = ibooker.book1D(label + "_hits", "track number of hits", 35, -0.5, 34.5);
0406   (mes.h_dRmin) = ibooker.book1D(label + "_dRmin", "track min dR", 100, 0., 0.01);
0407   (mes.h_dRmin_l) = ibooker.book1D(label + "_dRmin_l", "track min dR", 100, 0., 0.4);
0408 
0409   (mes.h_pt_vs_eta) = ibooker.book2D(label + "_ptVSeta",
0410                                      "track p_{T} vs #eta",
0411                                      Eta_nbin,
0412                                      Eta_rangeMin,
0413                                      Eta_rangeMax,
0414                                      Pt_nbin,
0415                                      Pt_rangeMin,
0416                                      Pt_rangeMax);
0417 
0418   // counts of tracks vs lumi
0419   // for this moment, xmin,xmax and binning are hardcoded, maybe in future in a config file!
0420   // have to add (declare) this in the .h file as well
0421   (mes.h_onlinelumi) = ibooker.book1D(label + "_onlinelumi",
0422                                       "number of tracks vs onlinelumi",
0423                                       onlinelumi_nbin,
0424                                       onlinelumi_rangeMin,
0425                                       onlinelumi_rangeMax);
0426   (mes.h_ls) = ibooker.book1D(label + "_ls", "number of tracks vs ls", ls_nbin, ls_rangeMin, ls_rangeMax);
0427   (mes.h_PU) = ibooker.book1D(label + "_PU", "number of tracks vs PU", PU_nbin, PU_rangeMin, PU_rangeMax);
0428 }
0429 
0430 void TrackToTrackComparisonHists::book_matching_tracks_histos(DQMStore::IBooker& ibooker,
0431                                                               matchingME& mes,
0432                                                               TString label,
0433                                                               std::string& dir) {
0434   ibooker.cd();
0435   ibooker.setCurrentFolder(dir);
0436 
0437   (mes.h_hits_vs_hits) = ibooker.book2D(
0438       label + "_hits_vs_hits", "monitored track # hits vs reference track # hits", 35, -0.5, 34.5, 35, -0.5, 34.5);
0439   (mes.h_pt_vs_pt) = ibooker.book2D(label + "_pt_vs_pt",
0440                                     "monitored track p_{T} vs reference track p_{T}",
0441                                     Pt_nbin,
0442                                     Pt_rangeMin,
0443                                     Pt_rangeMax,
0444                                     Pt_nbin,
0445                                     Pt_rangeMin,
0446                                     Pt_rangeMax);
0447   (mes.h_eta_vs_eta) = ibooker.book2D(label + "_eta_vs_eta",
0448                                       "monitored track #eta vs reference track #eta",
0449                                       Eta_nbin,
0450                                       Eta_rangeMin,
0451                                       Eta_rangeMax,
0452                                       Eta_nbin,
0453                                       Eta_rangeMin,
0454                                       Eta_rangeMax);
0455   (mes.h_phi_vs_phi) = ibooker.book2D(label + "_phi_vs_phi",
0456                                       "monitored track #phi vs reference track #phi",
0457                                       Phi_nbin,
0458                                       Phi_rangeMin,
0459                                       Phi_rangeMax,
0460                                       Phi_nbin,
0461                                       Phi_rangeMin,
0462                                       Phi_rangeMax);
0463 
0464   (mes.h_dPt) = ibooker.book1D(label + "_dPt", "#Delta track #P_T", ptRes_nbin, ptRes_rangeMin, ptRes_rangeMax);
0465   (mes.h_dEta) = ibooker.book1D(label + "_dEta", "#Delta track #eta", etaRes_nbin, etaRes_rangeMin, etaRes_rangeMax);
0466   (mes.h_dPhi) = ibooker.book1D(label + "_dPhi", "#Delta track #phi", phiRes_nbin, phiRes_rangeMin, phiRes_rangeMax);
0467   (mes.h_dDxy) = ibooker.book1D(
0468       label + "_dDxy", "#Delta track transverse dca to beam spot", dxyRes_nbin, dxyRes_rangeMin, dxyRes_rangeMax);
0469   (mes.h_dDz) = ibooker.book1D(
0470       label + "_dDz", "#Delta track longitudinal dca to beam spot", dzRes_nbin, dzRes_rangeMin, dzRes_rangeMax);
0471   (mes.h_dDxyWRTpv) = ibooker.book1D(label + "_dDxyWRTpv",
0472                                      "#Delta track transverse dca to primary vertex ",
0473                                      dxyRes_nbin,
0474                                      dxyRes_rangeMin,
0475                                      dxyRes_rangeMax);
0476   (mes.h_dDzWRTpv) = ibooker.book1D(label + "_dDzWRTpv",
0477                                     "#Delta track longitudinal dca to primary vertex",
0478                                     dzRes_nbin,
0479                                     dzRes_rangeMin,
0480                                     dzRes_rangeMax);
0481   (mes.h_dCharge) = ibooker.book1D(label + "_dCharge", "#Delta track charge", 5, -2.5, 2.5);
0482   (mes.h_dHits) = ibooker.book1D(label + "_dHits", "#Delta track number of hits", 39, -19.5, 19.5);
0483 }
0484 
0485 void TrackToTrackComparisonHists::fill_generic_tracks_histos(generalME& mes,
0486                                                              reco::Track* trk,
0487                                                              reco::BeamSpot* bs,
0488                                                              reco::Vertex* pv,
0489                                                              unsigned int ls,
0490                                                              double onlinelumi,
0491                                                              double PU,
0492                                                              bool requirePlateau) {
0493   float pt = trk->pt();
0494   float eta = trk->eta();
0495   float phi = trk->phi();
0496   float dxy = trk->dxy(bs->position());
0497   float dz = trk->dz(bs->position());
0498   float dxyWRTpv = trk->dxy(pv->position());
0499   float dzWRTpv = trk->dz(pv->position());
0500   float charge = trk->charge();
0501   float nhits = trk->hitPattern().numberOfValidHits();
0502 
0503   bool dxyOnPlateau = (fabs(dxyWRTpv) < dxyCutForPlateau_);
0504   bool pTOnPlateau = (pt > pTCutForPlateau_);
0505 
0506   if (dxyOnPlateau || !requirePlateau) {
0507     (mes.h_pt)->Fill(pt);
0508   }
0509 
0510   if ((pTOnPlateau && dxyOnPlateau) || !requirePlateau) {
0511     (mes.h_eta)->Fill(eta);
0512     (mes.h_phi)->Fill(phi);
0513     (mes.h_dz)->Fill(dz);
0514     (mes.h_dzWRTpv)->Fill(dzWRTpv);
0515     (mes.h_charge)->Fill(charge);
0516     (mes.h_hits)->Fill(nhits);
0517     (mes.h_onlinelumi)->Fill(onlinelumi);
0518     (mes.h_ls)->Fill(ls);
0519     (mes.h_PU)->Fill(PU);
0520   }
0521 
0522   if (pTOnPlateau || !requirePlateau) {
0523     (mes.h_dxy)->Fill(dxy);
0524     (mes.h_dxyWRTpv)->Fill(dxyWRTpv);
0525   }
0526 
0527   (mes.h_pt_vs_eta)->Fill(eta, pt);
0528 }
0529 
0530 void TrackToTrackComparisonHists::fill_matching_tracks_histos(
0531     matchingME& mes, reco::Track* mon, reco::Track* ref, reco::BeamSpot* bs, reco::Vertex* pv) {
0532   float mon_pt = mon->pt();
0533   float mon_eta = mon->eta();
0534   float mon_phi = mon->phi();
0535   float mon_dxy = mon->dxy(bs->position());
0536   float mon_dz = mon->dz(bs->position());
0537   float mon_dxyWRTpv = mon->dxy(pv->position());
0538   float mon_dzWRTpv = mon->dz(pv->position());
0539   float mon_charge = mon->charge();
0540   float mon_nhits = mon->hitPattern().numberOfValidHits();
0541 
0542   float ref_pt = ref->pt();
0543   float ref_eta = ref->eta();
0544   float ref_phi = ref->phi();
0545   float ref_dxy = ref->dxy(bs->position());
0546   float ref_dz = ref->dz(bs->position());
0547   float ref_dxyWRTpv = ref->dxy(pv->position());
0548   float ref_dzWRTpv = ref->dz(pv->position());
0549   float ref_charge = ref->charge();
0550   float ref_nhits = ref->hitPattern().numberOfValidHits();
0551 
0552   (mes.h_hits_vs_hits)->Fill(ref_nhits, mon_nhits);
0553   (mes.h_pt_vs_pt)->Fill(ref_pt, mon_pt);
0554   (mes.h_eta_vs_eta)->Fill(ref_eta, mon_eta);
0555   (mes.h_phi_vs_phi)->Fill(ref_phi, mon_phi);
0556 
0557   (mes.h_dPt)->Fill(ref_pt - mon_pt);
0558   (mes.h_dEta)->Fill(ref_eta - mon_eta);
0559   (mes.h_dPhi)->Fill(ref_phi - mon_phi);
0560   (mes.h_dDxy)->Fill(ref_dxy - mon_dxy);
0561   (mes.h_dDz)->Fill(ref_dz - mon_dz);
0562   (mes.h_dDxyWRTpv)->Fill(ref_dxyWRTpv - mon_dxyWRTpv);
0563   (mes.h_dDzWRTpv)->Fill(ref_dzWRTpv - mon_dzWRTpv);
0564   (mes.h_dCharge)->Fill(ref_charge - mon_charge);
0565   (mes.h_dHits)->Fill(ref_nhits - mon_nhits);
0566 }
0567 
0568 void TrackToTrackComparisonHists::initialize_parameter(const edm::ParameterSet& iConfig) {
0569   const edm::ParameterSet& pset = iConfig.getParameter<edm::ParameterSet>("histoPSet");
0570 
0571   Eta_rangeMin = pset.getParameter<double>("Eta_rangeMin");
0572   Eta_rangeMax = pset.getParameter<double>("Eta_rangeMax");
0573   Eta_nbin = pset.getParameter<unsigned int>("Eta_nbin");
0574 
0575   Pt_rangeMin = pset.getParameter<double>("Pt_rangeMin");
0576   Pt_rangeMax = pset.getParameter<double>("Pt_rangeMax");
0577   Pt_nbin = pset.getParameter<unsigned int>("Pt_nbin");
0578 
0579   Phi_rangeMin = pset.getParameter<double>("Phi_rangeMin");
0580   Phi_rangeMax = pset.getParameter<double>("Phi_rangeMax");
0581   Phi_nbin = pset.getParameter<unsigned int>("Phi_nbin");
0582 
0583   Dxy_rangeMin = pset.getParameter<double>("Dxy_rangeMin");
0584   Dxy_rangeMax = pset.getParameter<double>("Dxy_rangeMax");
0585   Dxy_nbin = pset.getParameter<unsigned int>("Dxy_nbin");
0586 
0587   Dz_rangeMin = pset.getParameter<double>("Dz_rangeMin");
0588   Dz_rangeMax = pset.getParameter<double>("Dz_rangeMax");
0589   Dz_nbin = pset.getParameter<unsigned int>("Dz_nbin");
0590 
0591   ptRes_rangeMin = pset.getParameter<double>("ptRes_rangeMin");
0592   ptRes_rangeMax = pset.getParameter<double>("ptRes_rangeMax");
0593   ptRes_nbin = pset.getParameter<unsigned int>("ptRes_nbin");
0594 
0595   phiRes_rangeMin = pset.getParameter<double>("phiRes_rangeMin");
0596   phiRes_rangeMax = pset.getParameter<double>("phiRes_rangeMax");
0597   phiRes_nbin = pset.getParameter<unsigned int>("phiRes_nbin");
0598 
0599   etaRes_rangeMin = pset.getParameter<double>("etaRes_rangeMin");
0600   etaRes_rangeMax = pset.getParameter<double>("etaRes_rangeMax");
0601   etaRes_nbin = pset.getParameter<unsigned int>("etaRes_nbin");
0602 
0603   dxyRes_rangeMin = pset.getParameter<double>("dxyRes_rangeMin");
0604   dxyRes_rangeMax = pset.getParameter<double>("dxyRes_rangeMax");
0605   dxyRes_nbin = pset.getParameter<unsigned int>("dxyRes_nbin");
0606 
0607   dzRes_rangeMin = pset.getParameter<double>("dzRes_rangeMin");
0608   dzRes_rangeMax = pset.getParameter<double>("dzRes_rangeMax");
0609   dzRes_nbin = pset.getParameter<unsigned int>("dzRes_nbin");
0610 
0611   ls_rangeMin = pset.getParameter<unsigned int>("ls_rangeMin");
0612   ls_rangeMax = pset.getParameter<unsigned int>("ls_rangeMax");
0613   ls_nbin = pset.getParameter<unsigned int>("ls_nbin");
0614 
0615   onlinelumi_rangeMin = pset.getParameter<double>("onlinelumi_rangeMin");
0616   onlinelumi_rangeMax = pset.getParameter<double>("onlinelumi_rangeMax");
0617   onlinelumi_nbin = pset.getParameter<unsigned int>("onlinelumi_nbin");
0618 
0619   PU_rangeMin = pset.getParameter<double>("PU_rangeMin");
0620   PU_rangeMax = pset.getParameter<double>("PU_rangeMax");
0621   PU_nbin = pset.getParameter<unsigned int>("PU_nbin");
0622 }
0623 
0624 void TrackToTrackComparisonHists::fillHistoPSetDescription(edm::ParameterSetDescription& pset) {
0625   pset.add<double>("Eta_rangeMin", -2.5);
0626   pset.add<double>("Eta_rangeMax", 2.5);
0627   pset.add<unsigned int>("Eta_nbin", 50);
0628 
0629   pset.add<double>("Pt_rangeMin", 0.1);
0630   pset.add<double>("Pt_rangeMax", 100.0);
0631   pset.add<unsigned int>("Pt_nbin", 1000);
0632 
0633   pset.add<double>("Phi_rangeMin", -3.1416);
0634   pset.add<double>("Phi_rangeMax", 3.1416);
0635   pset.add<unsigned int>("Phi_nbin", 36);
0636 
0637   pset.add<double>("Dxy_rangeMin", -1.0);
0638   pset.add<double>("Dxy_rangeMax", 1.0);
0639   pset.add<unsigned int>("Dxy_nbin", 300);
0640 
0641   pset.add<double>("Dz_rangeMin", -30.0);
0642   pset.add<double>("Dz_rangeMax", 30.0);
0643   pset.add<unsigned int>("Dz_nbin", 60);
0644 
0645   pset.add<double>("ptRes_rangeMin", -0.1);
0646   pset.add<double>("ptRes_rangeMax", 0.1);
0647   pset.add<unsigned int>("ptRes_nbin", 100);
0648 
0649   pset.add<double>("phiRes_rangeMin", -0.01);
0650   pset.add<double>("phiRes_rangeMax", 0.01);
0651   pset.add<unsigned int>("phiRes_nbin", 300);
0652 
0653   pset.add<double>("etaRes_rangeMin", -0.01);
0654   pset.add<double>("etaRes_rangeMax", 0.01);
0655   pset.add<unsigned int>("etaRes_nbin", 300);
0656 
0657   pset.add<double>("dxyRes_rangeMin", -0.05);
0658   pset.add<double>("dxyRes_rangeMax", 0.05);
0659   pset.add<unsigned int>("dxyRes_nbin", 500);
0660 
0661   pset.add<double>("dzRes_rangeMin", -0.05);
0662   pset.add<double>("dzRes_rangeMax", 0.05);
0663   pset.add<unsigned int>("dzRes_nbin", 150);
0664 
0665   pset.add<unsigned int>("ls_rangeMin", 0);
0666   pset.add<unsigned int>("ls_rangeMax", 3000);
0667   pset.add<unsigned int>("ls_nbin", 300);
0668 
0669   pset.add<double>("onlinelumi_rangeMin", 0.0);
0670   pset.add<double>("onlinelumi_rangeMax", 20000.0);
0671   pset.add<unsigned int>("onlinelumi_nbin", 200);
0672 
0673   pset.add<double>("PU_rangeMin", 0.0);
0674   pset.add<double>("PU_rangeMax", 120.0);
0675   pset.add<unsigned int>("PU_nbin", 120);
0676 }