Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-06-12 23:29:47

0001 // user includes
0002 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0003 #include "DQMServices/Core/interface/DQMStore.h"
0004 #include "DQMServices/Core/interface/MonitorElement.h"
0005 #include "DataFormats/TrackReco/interface/Track.h"
0006 #include "DataFormats/VertexReco/interface/Vertex.h"
0007 #include "FWCore/Framework/interface/Event.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0010 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0011 #include "FWCore/Utilities/interface/transform.h"  // for edm::vector_transform
0012 
0013 // ROOT includes
0014 #include "TLorentzVector.h"
0015 
0016 // standard includes
0017 #include <fmt/printf.h>
0018 
0019 class ShortenedTrackResolution : public DQMEDAnalyzer {
0020 public:
0021   ShortenedTrackResolution(const edm::ParameterSet &);
0022   static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0023 
0024 protected:
0025   void analyze(edm::Event const &iEvent, edm::EventSetup const &iSetup) override;
0026   void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
0027 
0028 private:
0029   const std::string folderName_;
0030   const std::vector<std::string> hitsRemain_;
0031   const double minTracksEta_;
0032   const double maxTracksEta_;
0033   const double minTracksPt_;
0034   const double maxTracksPt_;
0035 
0036   const double maxDr_;
0037   const edm::InputTag tracksTag_;
0038   const std::vector<edm::InputTag> tracksRerecoTag_;
0039   const edm::EDGetTokenT<std::vector<reco::Track>> tracksToken_;
0040   const std::vector<edm::EDGetTokenT<std::vector<reco::Track>>> tracksRerecoToken_;
0041 
0042   std::vector<MonitorElement *> histsPtRatioAll_;
0043   std::vector<MonitorElement *> histsPtDiffAll_;
0044   std::vector<MonitorElement *> histsEtaDiffAll_;
0045   std::vector<MonitorElement *> histsPhiDiffAll_;
0046   std::vector<MonitorElement *> histsPtRatioVsDeltaRAll_;
0047   std::vector<MonitorElement *> histsDeltaPtOverPtAll_;
0048   std::vector<MonitorElement *> histsPtAll_;
0049   std::vector<MonitorElement *> histsNhitsAll_;
0050   std::vector<MonitorElement *> histsDeltaRAll_;
0051 };
0052 
0053 // -----------------------------
0054 // constructors and destructor
0055 // -----------------------------
0056 ShortenedTrackResolution::ShortenedTrackResolution(const edm::ParameterSet &ps)
0057     : folderName_(ps.getUntrackedParameter<std::string>("folderName", "TrackRefitting")),
0058       hitsRemain_(ps.getUntrackedParameter<std::vector<std::string>>("hitsRemainInput")),
0059       minTracksEta_(ps.getUntrackedParameter<double>("minTracksEtaInput", 0.0)),
0060       maxTracksEta_(ps.getUntrackedParameter<double>("maxTracksEtaInput", 2.2)),
0061       minTracksPt_(ps.getUntrackedParameter<double>("minTracksPtInput", 15.0)),
0062       maxTracksPt_(ps.getUntrackedParameter<double>("maxTracksPtInput", 99999.9)),
0063       maxDr_(ps.getUntrackedParameter<double>("maxDrInput", 0.01)),
0064       tracksTag_(ps.getUntrackedParameter<edm::InputTag>("tracksInputTag", edm::InputTag("generalTracks", "", "DQM"))),
0065       tracksRerecoTag_(ps.getUntrackedParameter<std::vector<edm::InputTag>>("tracksRerecoInputTag")),
0066       tracksToken_(consumes<std::vector<reco::Track>>(tracksTag_)),
0067       tracksRerecoToken_(edm::vector_transform(
0068           tracksRerecoTag_, [this](edm::InputTag const &tag) { return consumes<std::vector<reco::Track>>(tag); })) {
0069   histsPtRatioAll_.clear();
0070   histsPtDiffAll_.clear();
0071   histsEtaDiffAll_.clear();
0072   histsPhiDiffAll_.clear();
0073   histsPtRatioVsDeltaRAll_.clear();
0074   histsDeltaPtOverPtAll_.clear();
0075   histsPtAll_.clear();
0076   histsNhitsAll_.clear();
0077   histsDeltaRAll_.clear();
0078 
0079   const size_t n = hitsRemain_.size();
0080   histsPtRatioAll_.reserve(n);
0081   histsPtDiffAll_.reserve(n);
0082   histsEtaDiffAll_.reserve(n);
0083   histsPhiDiffAll_.reserve(n);
0084   histsPtRatioVsDeltaRAll_.reserve(n);
0085   histsDeltaPtOverPtAll_.reserve(n);
0086   histsPtAll_.reserve(n);
0087   histsNhitsAll_.reserve(n);
0088   histsDeltaRAll_.reserve(n);
0089 }
0090 
0091 //__________________________________________________________________________________
0092 void ShortenedTrackResolution::bookHistograms(DQMStore::IBooker &iBook,
0093                                               edm::Run const &iRun,
0094                                               edm::EventSetup const &iSetup) {
0095   auto book1D = [&](const std::string &name, const std::string &title, int bins, double min, double max) {
0096     return iBook.book1D(name.c_str(), title.c_str(), bins, min, max);
0097   };
0098 
0099   auto book2D = [&](const std::string &name,
0100                     const std::string &title,
0101                     int binsX,
0102                     double minX,
0103                     double maxX,
0104                     int binsY,
0105                     double minY,
0106                     double maxY) {
0107     return iBook.book2D(name.c_str(), title.c_str(), binsX, minX, maxX, binsY, minY, maxY);
0108   };
0109 
0110   std::string currentFolder = folderName_ + "/Resolutions";
0111   iBook.setCurrentFolder(currentFolder);
0112 
0113   for (const auto &label : hitsRemain_) {
0114     std::string name, title;
0115 
0116     name = fmt::sprintf("trackPtRatio_%s", label);
0117     title =
0118         fmt::sprintf("Short Track p_{T} / Full Track p_{T} - %s layers;p_{T}^{short}/p_{T}^{full};n. tracks", label);
0119     histsPtRatioAll_.push_back(book1D(name, title, 100, 0.5, 1.5));
0120 
0121     name = fmt::sprintf("trackPtDiff_%s", label);
0122     title = fmt::sprintf(
0123         "Short Track p_{T} - Full Track p_{T} - %s layers;p_{T}^{short} - p_{T}^{full} [GeV];n. tracks", label);
0124     histsPtDiffAll_.push_back(book1D(name, title, 100, -10., 10.));
0125 
0126     name = fmt::sprintf("trackEtaDiff_%s", label);
0127     title = fmt::sprintf("Short Track #eta - Full Track #eta - %s layers;#eta^{short} - #eta^{full};n. tracks", label);
0128     histsEtaDiffAll_.push_back(book1D(name, title, 100, -0.001, 0.001));
0129 
0130     name = fmt::sprintf("trackPhiDiff_%s", label);
0131     title = fmt::sprintf("Short Track #phi - Full Track #phi - %s layers;#phi^{short} - #phi^{full};n. tracks", label);
0132     histsPhiDiffAll_.push_back(book1D(name, title, 100, -0.001, 0.001));
0133 
0134     name = fmt::sprintf("trackPtRatioVsDeltaR_%s", label);
0135     title = fmt::sprintf(
0136         "Short Track p_{T} / Full Track p_{T} - %s layers vs "
0137         "#DeltaR;#DeltaR(short,full);p_{T}^{short}/p_{T}^{full} [GeV];n. tracks",
0138         label);
0139     histsPtRatioVsDeltaRAll_.push_back(book2D(name, title, 100, 0., 0.01, 101, -0.05, 2.05));
0140 
0141     name = fmt::sprintf("trackDeltaPtOverPt_%s", label);
0142     title = fmt::sprintf(
0143         "Short Track p_{T} - Full Track p_{T} / Full Track p_{T} - %s layers;"
0144         "p_{T}^{short} - p_{T}^{full} / p^{full}_{T};n. tracks",
0145         label);
0146     histsDeltaPtOverPtAll_.push_back(book1D(name, title, 101, -5., 5.));
0147   }
0148 
0149   currentFolder = folderName_ + "/TrackProperties";
0150   iBook.setCurrentFolder(currentFolder);
0151 
0152   for (const auto &label : hitsRemain_) {
0153     std::string name, title;
0154 
0155     name = fmt::sprintf("trackPt_%s", label);
0156     title = fmt::sprintf("Short Track p_{T} - %s layers;p_{T}^{short} [GeV];n. tracks", label);
0157     histsPtAll_.push_back(book1D(name, title, 100, 0., 100.));
0158 
0159     name = fmt::sprintf("trackNhits_%s", label);
0160     title = fmt::sprintf("Short Track n. hits - %s layers;n. hits per track;n. tracks", label);
0161     histsNhitsAll_.push_back(book1D(name, title, 20, -0.5, 19.5));
0162 
0163     name = fmt::sprintf("trackDeltaR_%s", label);
0164     title = fmt::sprintf("Short Track / Full Track #DeltaR - %s layers;#DeltaR(short,full);n. tracks", label);
0165     histsDeltaRAll_.push_back(book1D(name, title, 100, 0., 0.005));
0166   }
0167 }
0168 
0169 //__________________________________________________________________________________
0170 void ShortenedTrackResolution::analyze(edm::Event const &iEvent, edm::EventSetup const &iSetup) {
0171   const auto &tracks = iEvent.getHandle(tracksToken_);
0172 
0173   if (!tracks.isValid()) {
0174     edm::LogError("ShortenedTrackResolution") << "Missing input track collection " << tracksTag_.encode() << std::endl;
0175     return;
0176   }
0177 
0178   for (const auto &track : *tracks) {
0179     const reco::HitPattern &hp = track.hitPattern();
0180     if (int(int(hp.numberOfValidHits()) - int(hp.numberOfAllHits(reco::HitPattern::TRACK_HITS))) != 0) {
0181       break;
0182     }
0183 
0184     TLorentzVector tvec;
0185     tvec.SetPtEtaPhiM(track.pt(), track.eta(), track.phi(), 0.0);
0186 
0187     int i = 0;  // token index
0188     for (const auto &token : tracksRerecoToken_) {
0189       const auto &tracks_rereco = iEvent.getHandle(token);
0190 
0191       for (const auto &track_rereco : *tracks_rereco) {
0192         TLorentzVector trerecovec;
0193         trerecovec.SetPtEtaPhiM(track_rereco.pt(), track_rereco.eta(), track_rereco.phi(), 0.0);
0194         double deltaR = tvec.DeltaR(trerecovec);
0195 
0196         if (deltaR < maxDr_) {
0197           if (track_rereco.pt() >= minTracksPt_ && track_rereco.pt() <= maxTracksPt_ &&
0198               std::abs(track_rereco.eta()) >= minTracksEta_ && std::abs(track_rereco.eta()) <= maxTracksEta_) {
0199             histsPtRatioAll_[i]->Fill(1.0 * track_rereco.pt() / track.pt());
0200             histsPtDiffAll_[i]->Fill(track_rereco.pt() - track.pt());
0201             histsDeltaPtOverPtAll_[i]->Fill((track_rereco.pt() - track.pt()) / track.pt());
0202             histsEtaDiffAll_[i]->Fill(track_rereco.eta() - track.eta());
0203             histsPhiDiffAll_[i]->Fill(track_rereco.phi() - track.phi());
0204             histsPtRatioVsDeltaRAll_[i]->Fill(deltaR, track_rereco.pt() / track.pt());
0205             histsPtAll_[i]->Fill(track_rereco.pt());
0206             histsNhitsAll_[i]->Fill(track_rereco.numberOfValidHits());
0207             histsDeltaRAll_[i]->Fill(deltaR);
0208           }
0209         }
0210       }
0211       ++i;
0212     }
0213   }
0214 }
0215 
0216 //__________________________________________________________________________________
0217 void ShortenedTrackResolution::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0218   edm::ParameterSetDescription desc;
0219   desc.addUntracked<std::string>("folderName", "TrackRefitting");
0220   desc.addUntracked<std::vector<std::string>>("hitsRemainInput", {});
0221   desc.addUntracked<double>("minTracksEtaInput", 0.0);
0222   desc.addUntracked<double>("maxTracksEtaInput", 2.2);
0223   desc.addUntracked<double>("minTracksPtInput", 15.0);
0224   desc.addUntracked<double>("maxTracksPtInput", 99999.9);
0225   desc.addUntracked<double>("maxDrInput", 0.01);
0226   desc.addUntracked<edm::InputTag>("tracksInputTag", edm::InputTag("generalTracks", "", "DQM"));
0227   desc.addUntracked<std::vector<edm::InputTag>>("tracksRerecoInputTag", {});
0228   descriptions.addWithDefaultLabel(desc);
0229 }
0230 
0231 // Define this as a plug-in
0232 #include "FWCore/Framework/interface/MakerMacros.h"
0233 DEFINE_FWK_MODULE(ShortenedTrackResolution);