Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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 *> histsPtAll_;
0043 };
0044 
0045 // -----------------------------
0046 // constructors and destructor
0047 // -----------------------------
0048 ShortenedTrackResolution::ShortenedTrackResolution(const edm::ParameterSet &ps)
0049     : folderName_(ps.getUntrackedParameter<std::string>("folderName", "TrackRefitting")),
0050       hitsRemain_(ps.getUntrackedParameter<std::vector<std::string>>("hitsRemainInput")),
0051       minTracksEta_(ps.getUntrackedParameter<double>("minTracksEtaInput", 0.0)),
0052       maxTracksEta_(ps.getUntrackedParameter<double>("maxTracksEtaInput", 2.2)),
0053       minTracksPt_(ps.getUntrackedParameter<double>("minTracksPtInput", 15.0)),
0054       maxTracksPt_(ps.getUntrackedParameter<double>("maxTracksPtInput", 99999.9)),
0055       maxDr_(ps.getUntrackedParameter<double>("maxDrInput", 0.01)),
0056       tracksTag_(ps.getUntrackedParameter<edm::InputTag>("tracksInputTag", edm::InputTag("generalTracks", "", "DQM"))),
0057       tracksRerecoTag_(ps.getUntrackedParameter<std::vector<edm::InputTag>>("tracksRerecoInputTag")),
0058       tracksToken_(consumes<std::vector<reco::Track>>(tracksTag_)),
0059       tracksRerecoToken_(edm::vector_transform(
0060           tracksRerecoTag_, [this](edm::InputTag const &tag) { return consumes<std::vector<reco::Track>>(tag); })) {
0061   histsPtAll_.clear();
0062 }
0063 
0064 //__________________________________________________________________________________
0065 void ShortenedTrackResolution::bookHistograms(DQMStore::IBooker &iBook,
0066                                               edm::Run const &iRun,
0067                                               edm::EventSetup const &iSetup) {
0068   std::string currentFolder = folderName_ + "/";
0069   iBook.setCurrentFolder(currentFolder);
0070 
0071   for (int i = 0; i < int(hitsRemain_.size()); ++i) {
0072     histsPtAll_.push_back(iBook.book1D(
0073         fmt::sprintf("trackPtRatio_%s", hitsRemain_[i]).c_str(),
0074         fmt::sprintf("Short Track p_{T} / Full Track p_{T} - %s layers;p_{T}^{short}/p_{T}^{full};n. tracks",
0075                      hitsRemain_[i])
0076             .c_str(),
0077         101,
0078         -0.05,
0079         2.05));
0080   }
0081 }
0082 
0083 //__________________________________________________________________________________
0084 void ShortenedTrackResolution::analyze(edm::Event const &iEvent, edm::EventSetup const &iSetup) {
0085   const auto &tracks = iEvent.getHandle(tracksToken_);
0086 
0087   if (!tracks.isValid()) {
0088     edm::LogError("ShortenedTrackResolution") << "Missing input track collection " << tracksTag_.encode() << std::endl;
0089     return;
0090   }
0091 
0092   for (const auto &track : *tracks) {
0093     const reco::HitPattern &hp = track.hitPattern();
0094     if (int(int(hp.numberOfValidHits()) - int(hp.numberOfAllHits(reco::HitPattern::TRACK_HITS))) != 0) {
0095       break;
0096     }
0097 
0098     TLorentzVector tvec;
0099     tvec.SetPtEtaPhiM(track.pt(), track.eta(), track.phi(), 0.0);
0100 
0101     int i = 0;  // token index
0102     for (const auto &token : tracksRerecoToken_) {
0103       const auto &tracks_rereco = iEvent.getHandle(token);
0104 
0105       for (const auto &track_rereco : *tracks_rereco) {
0106         TLorentzVector trerecovec;
0107         trerecovec.SetPtEtaPhiM(track_rereco.pt(), track_rereco.eta(), track_rereco.phi(), 0.0);
0108         double deltaR = tvec.DeltaR(trerecovec);
0109 
0110         if (deltaR < maxDr_) {
0111           if (track_rereco.pt() >= minTracksPt_ && track_rereco.pt() <= maxTracksPt_ &&
0112               std::abs(track_rereco.eta()) >= minTracksEta_ && std::abs(track_rereco.eta()) <= maxTracksEta_) {
0113             histsPtAll_[i]->Fill(1.0 * track_rereco.pt() / track.pt());
0114           }
0115         }
0116       }
0117       ++i;
0118     }
0119   }
0120 }
0121 
0122 //__________________________________________________________________________________
0123 void ShortenedTrackResolution::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0124   edm::ParameterSetDescription desc;
0125   desc.addUntracked<std::string>("folderName", "TrackRefitting");
0126   desc.addUntracked<std::vector<std::string>>("hitsRemainInput", {});
0127   desc.addUntracked<double>("minTracksEtaInput", 0.0);
0128   desc.addUntracked<double>("maxTracksEtaInput", 2.2);
0129   desc.addUntracked<double>("minTracksPtInput", 15.0);
0130   desc.addUntracked<double>("maxTracksPtInput", 99999.9);
0131   desc.addUntracked<double>("maxDrInput", 0.01);
0132   desc.addUntracked<edm::InputTag>("tracksInputTag", edm::InputTag("generalTracks", "", "DQM"));
0133   desc.addUntracked<std::vector<edm::InputTag>>("tracksRerecoInputTag", {});
0134   descriptions.addWithDefaultLabel(desc);
0135 }
0136 
0137 // Define this as a plug-in
0138 #include "FWCore/Framework/interface/MakerMacros.h"
0139 DEFINE_FWK_MODULE(ShortenedTrackResolution);