File indexing completed on 2024-06-22 02:24:05
0001 #include "RecoParticleFlow/PFProducer/interface/BlockElementImporterBase.h"
0002 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementTrack.h"
0003 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementGsfTrack.h"
0004 #include "DataFormats/ParticleFlowReco/interface/PFRecTrack.h"
0005 #include "DataFormats/TrackReco/interface/Track.h"
0006 #include "DataFormats/MuonReco/interface/Muon.h"
0007 #include "DataFormats/Common/interface/ValueMap.h"
0008 #include "RecoParticleFlow/PFProducer/interface/PFMuonAlgo.h"
0009 #include "RecoParticleFlow/PFTracking/interface/PFTrackAlgoTools.h"
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011
0012
0013
0014
0015 class TrackTimingImporter : public BlockElementImporterBase {
0016 public:
0017 TrackTimingImporter(const edm::ParameterSet& conf, edm::ConsumesCollector& cc)
0018 : BlockElementImporterBase(conf, cc),
0019 useTimeQuality_(conf.getParameter<bool>("useTimeQuality")),
0020 timeQualityThreshold_(useTimeQuality_ ? conf.getParameter<double>("timeQualityThreshold") : -99.),
0021 srcTime_(cc.consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("timeValueMap"))),
0022 srcTimeError_(cc.consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("timeErrorMap"))),
0023 srcTimeQuality_(useTimeQuality_
0024 ? cc.consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("timeQualityMap"))
0025 : edm::EDGetTokenT<edm::ValueMap<float>>()),
0026 srcTimeGsf_(cc.consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("timeValueMapGsf"))),
0027 srcTimeErrorGsf_(cc.consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("timeErrorMapGsf"))),
0028 srcTimeQualityGsf_(
0029 useTimeQuality_ ? cc.consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("timeQualityMapGsf"))
0030 : edm::EDGetTokenT<edm::ValueMap<float>>()),
0031 debug_(conf.getUntrackedParameter<bool>("debug", false)) {}
0032
0033 void importToBlock(const edm::Event&, ElementList&) const override;
0034
0035 private:
0036 const bool useTimeQuality_;
0037 const double timeQualityThreshold_;
0038 edm::EDGetTokenT<edm::ValueMap<float>> srcTime_, srcTimeError_, srcTimeQuality_, srcTimeGsf_, srcTimeErrorGsf_,
0039 srcTimeQualityGsf_;
0040 const bool debug_;
0041 };
0042
0043 DEFINE_EDM_PLUGIN(BlockElementImporterFactory, TrackTimingImporter, "TrackTimingImporter");
0044
0045 void TrackTimingImporter::importToBlock(const edm::Event& e, BlockElementImporterBase::ElementList& elems) const {
0046 typedef BlockElementImporterBase::ElementList::value_type ElementType;
0047
0048 auto const& time = e.get(srcTime_);
0049 auto const& timeErr = e.get(srcTimeError_);
0050 auto const& timeGsf = e.get(srcTimeGsf_);
0051 auto const& timeErrGsf = e.get(srcTimeErrorGsf_);
0052
0053 edm::Handle<edm::ValueMap<float>> timeQualH, timeQualGsfH;
0054 if (useTimeQuality_) {
0055 e.getByToken(srcTimeQuality_, timeQualH);
0056 e.getByToken(srcTimeQualityGsf_, timeQualGsfH);
0057 }
0058
0059 for (auto& elem : elems) {
0060 if (reco::PFBlockElement::TRACK == elem->type()) {
0061 const auto& ref = elem->trackRef();
0062 if (time.contains(ref.id())) {
0063 const bool assocQuality = useTimeQuality_ ? (*timeQualH)[ref] > timeQualityThreshold_ : true;
0064 if (assocQuality) {
0065 elem->setTime(time[ref], timeErr[ref]);
0066 } else {
0067 elem->setTime(0., -1.);
0068 }
0069 }
0070 if (debug_) {
0071 edm::LogInfo("TrackTimingImporter")
0072 << "Track with pT / eta " << ref->pt() << " / " << ref->eta() << " has time: " << elem->time() << " +/- "
0073 << elem->timeError() << std::endl;
0074 }
0075 } else if (reco::PFBlockElement::GSF == elem->type()) {
0076 const auto& ref = static_cast<const reco::PFBlockElementGsfTrack*>(elem.get())->GsftrackRef();
0077 if (timeGsf.contains(ref.id())) {
0078 const bool assocQuality = useTimeQuality_ ? (*timeQualGsfH)[ref] > timeQualityThreshold_ : true;
0079 if (assocQuality) {
0080 elem->setTime(timeGsf[ref], timeErrGsf[ref]);
0081 } else {
0082 elem->setTime(0., -1.);
0083 }
0084 }
0085 if (debug_) {
0086 edm::LogInfo("TrackTimingImporter")
0087 << "Track with pT / eta " << ref->pt() << " / " << ref->eta() << " has time: " << elem->time() << " +/- "
0088 << elem->timeError() << std::endl;
0089 }
0090 }
0091 }
0092 }