CMS 3D CMS Logo

TrackTimingImporter.cc
Go to the documentation of this file.
11 
12 // this doesn't actually import anything,
13 // but rather applies time stamps to tracks after they are all inserted
14 
16 public:
18  : BlockElementImporterBase(conf, cc),
19  useTimeQuality_(conf.existsAs<edm::InputTag>("timeQualityMap")),
20  timeQualityThreshold_(useTimeQuality_ ? conf.getParameter<double>("timeQualityThreshold") : -99.),
21  srcTime_(cc.consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("timeValueMap"))),
22  srcTimeError_(cc.consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("timeErrorMap"))),
24  ? cc.consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("timeQualityMap"))
25  : edm::EDGetTokenT<edm::ValueMap<float>>()),
26  srcTimeGsf_(cc.consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("timeValueMapGsf"))),
27  srcTimeErrorGsf_(cc.consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("timeErrorMapGsf"))),
29  useTimeQuality_ ? cc.consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("timeQualityMapGsf"))
30  : edm::EDGetTokenT<edm::ValueMap<float>>()),
31  debug_(conf.getUntrackedParameter<bool>("debug", false)) {}
32 
33  void importToBlock(const edm::Event&, ElementList&) const override;
34 
35 private:
36  const bool useTimeQuality_;
37  const double timeQualityThreshold_;
40  const bool debug_;
41 };
42 
44 
47 
48  auto const& time = e.get(srcTime_);
49  auto const& timeErr = e.get(srcTimeError_);
50  auto const& timeGsf = e.get(srcTimeGsf_);
51  auto const& timeErrGsf = e.get(srcTimeErrorGsf_);
52 
53  edm::Handle<edm::ValueMap<float>> timeQualH, timeQualGsfH;
54  if (useTimeQuality_) {
55  e.getByToken(srcTimeQuality_, timeQualH);
56  e.getByToken(srcTimeQualityGsf_, timeQualGsfH);
57  }
58 
59  for (auto& elem : elems) {
60  if (reco::PFBlockElement::TRACK == elem->type()) {
61  const auto& ref = elem->trackRef();
62  if (time.contains(ref.id())) {
63  const bool assocQuality = useTimeQuality_ ? (*timeQualH)[ref] > timeQualityThreshold_ : true;
64  if (assocQuality) {
65  elem->setTime(time[ref], timeErr[ref]);
66  } else {
67  elem->setTime(0., -1.);
68  }
69  }
70  if (debug_) {
71  edm::LogInfo("TrackTimingImporter")
72  << "Track with pT / eta " << ref->pt() << " / " << ref->eta() << " has time: " << elem->time() << " +/- "
73  << elem->timeError() << std::endl;
74  }
75  } else if (reco::PFBlockElement::GSF == elem->type()) {
76  const auto& ref = static_cast<const reco::PFBlockElementGsfTrack*>(elem.get())->GsftrackRef();
77  if (timeGsf.contains(ref.id())) {
78  const bool assocQuality = useTimeQuality_ ? (*timeQualGsfH)[ref] > timeQualityThreshold_ : true;
79  if (assocQuality) {
80  elem->setTime(timeGsf[ref], timeErrGsf[ref]);
81  } else {
82  elem->setTime(0., -1.);
83  }
84  }
85  if (debug_) {
86  edm::LogInfo("TrackTimingImporter")
87  << "Track with pT / eta " << ref->pt() << " / " << ref->eta() << " has time: " << elem->time() << " +/- "
88  << elem->timeError() << std::endl;
89  }
90  }
91  }
92 }
edm::EDGetTokenT< edm::ValueMap< float > > srcTimeQuality_
void importToBlock(const edm::Event &, ElementList &) const override
edm::EDGetTokenT< edm::ValueMap< float > > srcTimeErrorGsf_
TrackTimingImporter(const edm::ParameterSet &conf, edm::ConsumesCollector &cc)
edm::EDGetTokenT< edm::ValueMap< float > > srcTimeError_
edm::EDGetTokenT< edm::ValueMap< float > > srcTime_
const double timeQualityThreshold_
Log< level::Info, false > LogInfo
HLT enums.
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< std::unique_ptr< reco::PFBlockElement > > ElementList
edm::EDGetTokenT< edm::ValueMap< float > > srcTimeGsf_
edm::EDGetTokenT< edm::ValueMap< float > > srcTimeQualityGsf_